Nektar++
EquationSystem.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: EquationSystem.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Main wrapper class for Advection Diffusion Reaction Solver
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <boost/core/ignore_unused.hpp>
36 
39 
41 #include <LocalRegions/MatrixKey.h>
42 #include <MultiRegions/ContField.h>
45 
46 #include <MultiRegions/ExpList.h>
49 
52 
53 #include <GlobalMapping/Mapping.h>
54 
55 #include <boost/format.hpp>
56 
57 #include <iostream>
58 #include <string>
59 
60 using namespace std;
61 
62 namespace Nektar
63 {
64 namespace SolverUtils
65 {
66 
67 std::string EquationSystem::equationSystemTypeLookupIds[2] = {
68  LibUtilities::SessionReader::RegisterEnumValue("DEALIASING", "True", 0),
69  LibUtilities::SessionReader::RegisterEnumValue("DEALIASING", "False", 1)};
70 
71 /**
72  * @class EquationSystem
73  *
74  * This class is a base class for all solver implementations. It
75  * provides the underlying generic functionality and interface for
76  * solving equations.
77  *
78  * To solve a steady-state equation, create a derived class from this
79  * class and reimplement the virtual functions to provide custom
80  * implementation for the problem.
81  *
82  * To solve unsteady problems, derive from the UnsteadySystem class
83  * instead which provides general time integration.
84  */
86 {
87  static EquationSystemFactory instance;
88  return instance;
89 }
90 
91 /**
92  * This constructor is protected as the objects of this class are never
93  * instantiated directly.
94  * @param pSession The session reader holding problem parameters.
95  */
96 EquationSystem::EquationSystem(
99  : m_comm(pSession->GetComm()), m_session(pSession), m_graph(pGraph),
100  m_lambda(0), m_infosteps(10),
101  m_fieldMetaDataMap(LibUtilities::NullFieldMetaDataMap)
102 {
103  // set up session names in fieldMetaDataMap
104  const vector<std::string> filenames = m_session->GetFilenames();
105 
106  for (int i = 0; i < filenames.size(); ++i)
107  {
108  string sessionname = "SessionName";
109  sessionname += boost::lexical_cast<std::string>(i);
110  m_fieldMetaDataMap[sessionname] = filenames[i];
111  m_fieldMetaDataMap["ChkFileNum"] = boost::lexical_cast<std::string>(0);
112  }
113 }
114 
115 /**
116  * @brief Initialisation object for EquationSystem.
117  */
118 void EquationSystem::v_InitObject(bool DeclareFields)
119 {
120  // Save the basename of input file name for output details
121  m_sessionName = m_session->GetSessionName();
122 
123  // Instantiate a field reader/writer
125 
126  // Also read and store the boundary conditions
129  m_session, m_graph);
130 
131  // Set space dimension for use in class
132  m_spacedim = m_graph->GetSpaceDimension();
133 
134  // Setting parameteres for homogenous problems
135  m_HomoDirec = 0;
136  m_useFFT = false;
137  m_homogen_dealiasing = false;
138  m_singleMode = false;
139  m_halfMode = false;
140  m_multipleModes = false;
142 
143  m_verbose = m_session->DefinesCmdLineArgument("verbose");
144 
145  if (m_session->DefinesSolverInfo("HOMOGENEOUS"))
146  {
147  std::string HomoStr = m_session->GetSolverInfo("HOMOGENEOUS");
148  m_spacedim = 3;
149 
150  if ((HomoStr == "HOMOGENEOUS1D") || (HomoStr == "Homogeneous1D") ||
151  (HomoStr == "1D") || (HomoStr == "Homo1D"))
152  {
154  m_session->LoadParameter("LZ", m_LhomZ);
155  m_HomoDirec = 1;
156 
157  if (m_session->DefinesSolverInfo("ModeType"))
158  {
159  m_session->MatchSolverInfo("ModeType", "SingleMode",
160  m_singleMode, false);
161  m_session->MatchSolverInfo("ModeType", "HalfMode", m_halfMode,
162  false);
163  m_session->MatchSolverInfo("ModeType", "MultipleModes",
164  m_multipleModes, false);
165  }
166 
167  // Stability Analysis flags
168  if (m_session->DefinesSolverInfo("ModeType"))
169  {
170  if (m_singleMode)
171  {
172  m_npointsZ = 2;
173  }
174  else if (m_halfMode)
175  {
176  m_npointsZ = 1;
177  }
178  else if (m_multipleModes)
179  {
180  m_npointsZ = m_session->GetParameter("HomModesZ");
181  }
182  else
183  {
184  ASSERTL0(false, "SolverInfo ModeType not valid");
185  }
186  }
187  else
188  {
189  m_npointsZ = m_session->GetParameter("HomModesZ");
190  }
191  }
192 
193  if ((HomoStr == "HOMOGENEOUS2D") || (HomoStr == "Homogeneous2D") ||
194  (HomoStr == "2D") || (HomoStr == "Homo2D"))
195  {
197  m_session->LoadParameter("HomModesY", m_npointsY);
198  m_session->LoadParameter("LY", m_LhomY);
199  m_session->LoadParameter("HomModesZ", m_npointsZ);
200  m_session->LoadParameter("LZ", m_LhomZ);
201  m_HomoDirec = 2;
202  }
203 
204  if ((HomoStr == "HOMOGENEOUS3D") || (HomoStr == "Homogeneous3D") ||
205  (HomoStr == "3D") || (HomoStr == "Homo3D"))
206  {
208  m_session->LoadParameter("HomModesY", m_npointsY);
209  m_session->LoadParameter("LY", m_LhomY);
210  m_session->LoadParameter("HomModesZ", m_npointsZ);
211  m_session->LoadParameter("LZ", m_LhomZ);
212  m_HomoDirec = 2;
213  }
214 
215  m_session->MatchSolverInfo("USEFFT", "FFTW", m_useFFT, false);
216 
217  m_session->MatchSolverInfo("DEALIASING", "True", m_homogen_dealiasing,
218  false);
219  }
220  else
221  {
222  // set to default value so can use to identify 2d or 3D
223  // (homogeneous) expansions
224  m_npointsZ = 1;
225  }
226 
227  m_session->MatchSolverInfo("SPECTRALHPDEALIASING", "True",
228  m_specHP_dealiasing, false);
229  if (m_specHP_dealiasing == false)
230  {
231  m_session->MatchSolverInfo("SPECTRALHPDEALIASING", "On",
232  m_specHP_dealiasing, false);
233  }
234 
235  // Options to determine type of projection from file or directly
236  // from constructor
237  if (m_session->DefinesSolverInfo("PROJECTION"))
238  {
239  std::string ProjectStr = m_session->GetSolverInfo("PROJECTION");
240 
241  if ((ProjectStr == "Continuous") || (ProjectStr == "Galerkin") ||
242  (ProjectStr == "CONTINUOUS") || (ProjectStr == "GALERKIN"))
243  {
245  }
246  else if ((ProjectStr == "MixedCGDG") ||
247  (ProjectStr == "Mixed_CG_Discontinuous"))
248  {
250  }
251  else if (ProjectStr == "DisContinuous")
252  {
254  }
255  else
256  {
257  ASSERTL0(false, "PROJECTION value not recognised");
258  }
259  }
260  else
261  {
262  cerr << "Projection type not specified in SOLVERINFO,"
263  "defaulting to continuous Galerkin"
264  << endl;
266  }
267 
268  // Enforce singularity check for some problems
270 
271  int i;
272  int nvariables = m_session->GetVariables().size();
273  bool DeclareCoeffPhysArrays = true;
274 
276  m_spacedim = m_graph->GetSpaceDimension() + m_HomoDirec;
277  m_expdim = m_graph->GetMeshDimension();
278 
279  if (DeclareFields) // declare field if required
280  {
281  /// Continuous field
284  {
285  switch (m_expdim)
286  {
287  case 1:
288  {
291  {
292  const LibUtilities::PointsKey PkeyY(
294  const LibUtilities::BasisKey BkeyY(
296  const LibUtilities::PointsKey PkeyZ(
298  const LibUtilities::BasisKey BkeyZ(
300 
301  for (i = 0; i < m_fields.size(); i++)
302  {
303  m_fields[i] = MemoryManager<
305  AllocateSharedPtr(m_session, BkeyY, BkeyZ,
308  m_session->GetVariable(i));
309  }
310  }
311  else
312  {
313  for (i = 0; i < m_fields.size(); i++)
314  {
315  m_fields[i] =
319  m_session->GetVariable(i));
320  }
321  }
322  break;
323  }
324  case 2:
325  {
327  {
328  // Fourier single mode stability analysis
329  if (m_singleMode)
330  {
331  const LibUtilities::PointsKey PkeyZ(
332  m_npointsZ,
334 
335  const LibUtilities::BasisKey BkeyZ(
337  PkeyZ);
338 
339  for (i = 0; i < m_fields.size(); i++)
340  {
341  m_fields[i] = MemoryManager<
343  AllocateSharedPtr(
344  m_session, BkeyZ, m_LhomZ, m_useFFT,
346  m_session->GetVariable(i),
348  }
349  }
350  // Half mode stability analysis
351  else if (m_halfMode)
352  {
353  const LibUtilities::PointsKey PkeyZ(
354  m_npointsZ,
356 
357  const LibUtilities::BasisKey BkeyZR(
359  PkeyZ);
360 
361  const LibUtilities::BasisKey BkeyZI(
363  PkeyZ);
364 
365  for (i = 0; i < m_fields.size(); i++)
366  {
367  if (m_session->GetVariable(i).compare("w") == 0)
368  {
369  m_fields[i] = MemoryManager<
370  MultiRegions::
371  ContField3DHomogeneous1D>::
372  AllocateSharedPtr(
373  m_session, BkeyZI, m_LhomZ,
375  m_graph, m_session->GetVariable(i),
377  }
378  else
379  {
380  m_fields[i] = MemoryManager<
381  MultiRegions::
382  ContField3DHomogeneous1D>::
383  AllocateSharedPtr(
384  m_session, BkeyZR, m_LhomZ,
386  m_graph, m_session->GetVariable(i),
388  }
389  }
390  }
391  // Normal homogeneous 1D
392  else
393  {
394  const LibUtilities::PointsKey PkeyZ(
396  const LibUtilities::BasisKey BkeyZ(
398 
399  for (i = 0; i < m_fields.size(); i++)
400  {
401  m_fields[i] = MemoryManager<
403  AllocateSharedPtr(
404  m_session, BkeyZ, m_LhomZ, m_useFFT,
406  m_session->GetVariable(i),
408  }
409  }
410  }
411  else
412  {
413  i = 0;
417  m_session->GetVariable(i),
418  DeclareCoeffPhysArrays,
420  m_fields[0] = firstfield;
421  for (i = 1; i < m_fields.size(); i++)
422  {
423  if (m_graph->SameExpansionInfo(
424  m_session->GetVariable(0),
425  m_session->GetVariable(i)))
426  {
427  m_fields[i] =
430  *firstfield, m_graph,
431  m_session->GetVariable(i),
432  DeclareCoeffPhysArrays,
434  }
435  else
436  {
437  m_fields[i] =
441  m_session->GetVariable(i),
442  DeclareCoeffPhysArrays,
444  }
445  }
446 
447  if (m_projectionType ==
449  {
450  /// Setting up the normals
453 
454  for (i = 0; i < m_spacedim; ++i)
455  {
456  m_traceNormals[i] =
458  }
459 
460  m_fields[0]->GetTrace()->GetNormals(m_traceNormals);
461  }
462  }
463 
464  break;
465  }
466  case 3:
467  {
468  i = 0;
472  m_session->GetVariable(i),
473  DeclareCoeffPhysArrays,
475 
476  m_fields[0] = firstfield;
477  for (i = 1; i < m_fields.size(); i++)
478  {
479  if (m_graph->SameExpansionInfo(
480  m_session->GetVariable(0),
481  m_session->GetVariable(i)))
482  {
483  m_fields[i] =
486  *firstfield, m_graph,
487  m_session->GetVariable(i),
488  DeclareCoeffPhysArrays,
490  }
491  else
492  {
493  m_fields[i] =
497  m_session->GetVariable(i),
498  DeclareCoeffPhysArrays,
500  }
501  }
502 
503  if (m_projectionType ==
505  {
506  /// Setting up the normals
509  for (i = 0; i < m_spacedim; ++i)
510  {
511  m_traceNormals[i] =
513  }
514 
515  m_fields[0]->GetTrace()->GetNormals(m_traceNormals);
516  // Call the trace on all fields to ensure DG setup.
517  for (i = 1; i < m_fields.size(); ++i)
518  {
519  m_fields[i]->GetTrace();
520  }
521  }
522  break;
523  }
524  default:
525  ASSERTL0(false, "Expansion dimension not recognised");
526  break;
527  }
528  }
529  // Discontinuous field
530  else
531  {
532  switch (m_expdim)
533  {
534  case 1:
535  {
538  {
539  const LibUtilities::PointsKey PkeyY(
541  const LibUtilities::BasisKey BkeyY(
543  const LibUtilities::PointsKey PkeyZ(
545  const LibUtilities::BasisKey BkeyZ(
547 
548  for (i = 0; i < m_fields.size(); i++)
549  {
550  m_fields[i] = MemoryManager<
552  AllocateSharedPtr(m_session, BkeyY, BkeyZ,
555  m_session->GetVariable(i));
556  }
557  }
558  else
559  {
560  for (i = 0; i < m_fields.size(); i++)
561  {
562  m_fields[i] =
566  m_session->GetVariable(i));
567  }
568  }
569 
570  break;
571  }
572  case 2:
573  {
575  {
576  const LibUtilities::PointsKey PkeyZ(
578  const LibUtilities::BasisKey BkeyZ(
580 
581  for (i = 0; i < m_fields.size(); i++)
582  {
583  m_fields[i] = MemoryManager<
585  AllocateSharedPtr(m_session, BkeyZ, m_LhomZ,
586  m_useFFT,
588  m_session->GetVariable(i));
589  }
590  }
591  else
592  {
593  i = 0;
597  m_session->GetVariable(i));
598  m_fields[0] = firstfield;
599  for (i = 1; i < m_fields.size(); i++)
600  {
601  if (m_graph->SameExpansionInfo(
602  m_session->GetVariable(0),
603  m_session->GetVariable(i)))
604  {
605  m_fields[i] =
608  *firstfield, m_graph,
609  m_session->GetVariable(i));
610  }
611  else
612  {
613  m_fields[i] =
617  m_session->GetVariable(i));
618  }
619  }
620  }
621 
622  break;
623  }
624  case 3:
625  {
627  {
628  ASSERTL0(
629  false,
630  "3D fully periodic problems not implemented yet");
631  }
632  else
633  {
634  i = 0;
638  m_session->GetVariable(i));
639  m_fields[0] = firstfield;
640  for (i = 1; i < m_fields.size(); i++)
641  {
642  if (m_graph->SameExpansionInfo(
643  m_session->GetVariable(0),
644  m_session->GetVariable(i)))
645  {
646  m_fields[i] =
649  *firstfield, m_graph,
650  m_session->GetVariable(i));
651  }
652  else
653  {
654  m_fields[i] =
658  m_session->GetVariable(i));
659  }
660  }
661  }
662  break;
663  }
664  default:
665  ASSERTL0(false, "Expansion dimension not recognised");
666  break;
667  }
668 
669  // Setting up the normals
671 
672  for (i = 0; i < m_spacedim; ++i)
673  {
674  m_traceNormals[i] =
676  }
677 
678  m_fields[0]->GetTrace()->GetNormals(m_traceNormals);
679  }
680  // Zero all physical fields initially
681  ZeroPhysFields();
682  }
683 
684  // Set Default Parameter
685  m_session->LoadParameter("Time", m_time, 0.0);
686  m_session->LoadParameter("TimeStep", m_timestep, 0.0);
687  m_session->LoadParameter("NumSteps", m_steps, 0);
688  m_session->LoadParameter("IO_CheckSteps", m_checksteps, 0);
689  m_session->LoadParameter("IO_CheckTime", m_checktime, 0.0);
690  m_session->LoadParameter("FinTime", m_fintime, 0);
691  m_session->LoadParameter("NumQuadPointsError", m_NumQuadPointsError, 0);
692 
693  // Check uniqueness of checkpoint output
694  ASSERTL0((m_checktime == 0.0 && m_checksteps == 0) ||
695  (m_checktime > 0.0 && m_checksteps == 0) ||
696  (m_checktime == 0.0 && m_checksteps > 0),
697  "Only one of IO_CheckTime and IO_CheckSteps "
698  "should be set!");
699  m_session->LoadParameter("TimeIncrementFactor", m_TimeIncrementFactor, 1.0);
700 
701  m_nchk = 0;
702  m_pararealIter = 0;
703  m_useInitialCondition = true;
704 }
705 
706 /**
707  * @brief Destructor for class EquationSystem.
708  */
710 {
712  LocalRegions::MatrixKey::opLess>::ClearManager();
714  LocalRegions::MatrixKey::opLess>::ClearManager();
715 }
716 
718  std::string name, const MultiRegions::ExpListSharedPtr &field, bool cache)
719 {
720  MultiRegions::ExpListSharedPtr vField = field;
721  if (!field)
722  {
723  vField = m_fields[0];
724  }
725 
726  if (cache)
727  {
728  if ((m_sessionFunctions.find(name) == m_sessionFunctions.end()) ||
729  (m_sessionFunctions[name]->GetSession() != m_session) ||
730  (m_sessionFunctions[name]->GetExpansion() != vField))
731  {
734  m_session, vField, name, cache);
735  }
736 
737  return m_sessionFunctions[name];
738  }
739  else
740  {
742  new SessionFunction(m_session, vField, name, cache));
743  }
744 }
745 
746 /**
747  * If boundary conditions are time-dependent, they will be evaluated at
748  * the time specified.
749  * @param time The time at which to evaluate the BCs
750  */
752 {
753  std::string varName;
754  int nvariables = m_fields.size();
755  for (int i = 0; i < nvariables; ++i)
756  {
757  varName = m_session->GetVariable(i);
758  m_fields[i]->EvaluateBoundaryConditions(time, varName);
759  }
760 }
761 
762 /**
763  * Compute the error in the L2-norm.
764  * @param field The field to compare.
765  * @param exactsoln The exact solution to compare with.
766  * @param Normalised Normalise L2-error.
767  * @returns Error in the L2-norm.
768  */
770  const Array<OneD, NekDouble> &exactsoln,
771  bool Normalised)
772 {
773  NekDouble L2error = -1.0;
774 
775  if (m_NumQuadPointsError == 0)
776  {
777  if (m_fields[field]->GetPhysState() == false)
778  {
779  m_fields[field]->BwdTrans(m_fields[field]->GetCoeffs(),
780  m_fields[field]->UpdatePhys());
781  }
782 
783  if (exactsoln.size())
784  {
785  L2error =
786  m_fields[field]->L2(m_fields[field]->GetPhys(), exactsoln);
787  }
788  else if (m_session->DefinesFunction("ExactSolution"))
789  {
790  Array<OneD, NekDouble> exactsoln(m_fields[field]->GetNpoints());
791 
792  GetFunction("ExactSolution")
793  ->Evaluate(m_session->GetVariable(field), exactsoln, m_time);
794 
795  L2error =
796  m_fields[field]->L2(m_fields[field]->GetPhys(), exactsoln);
797  }
798  else
799  {
800  L2error = m_fields[field]->L2(m_fields[field]->GetPhys());
801  }
802 
803  if (Normalised == true)
804  {
805  Array<OneD, NekDouble> one(m_fields[field]->GetNpoints(), 1.0);
806 
807  NekDouble Vol = m_fields[field]->Integral(one);
808  L2error = sqrt(L2error * L2error / Vol);
809  }
810  }
811  else
812  {
813  Array<OneD, NekDouble> L2INF(2);
814  L2INF = ErrorExtraPoints(field);
815  L2error = L2INF[0];
816  }
817  return L2error;
818 }
819 
820 /**
821  * Compute the error in the L_inf-norm
822  * @param field The field to compare.
823  * @param exactsoln The exact solution to compare with.
824  * @returns Error in the L_inft-norm.
825  */
827  const Array<OneD, NekDouble> &exactsoln)
828 {
829  NekDouble Linferror = -1.0;
830 
831  if (m_NumQuadPointsError == 0)
832  {
833  if (m_fields[field]->GetPhysState() == false)
834  {
835  m_fields[field]->BwdTrans(m_fields[field]->GetCoeffs(),
836  m_fields[field]->UpdatePhys());
837  }
838 
839  if (exactsoln.size())
840  {
841  Linferror =
842  m_fields[field]->Linf(m_fields[field]->GetPhys(), exactsoln);
843  }
844  else if (m_session->DefinesFunction("ExactSolution"))
845  {
846  Array<OneD, NekDouble> exactsoln(m_fields[field]->GetNpoints());
847 
848  GetFunction("ExactSolution")
849  ->Evaluate(m_session->GetVariable(field), exactsoln, m_time);
850 
851  Linferror =
852  m_fields[field]->Linf(m_fields[field]->GetPhys(), exactsoln);
853  }
854  else
855  {
856  Linferror = m_fields[field]->Linf(m_fields[field]->GetPhys());
857  }
858  }
859  else
860  {
861  Array<OneD, NekDouble> L2INF(2);
862  L2INF = ErrorExtraPoints(field);
863  Linferror = L2INF[1];
864  }
865 
866  return Linferror;
867 }
868 
869 /**
870  * Compute the error in the L2-norm, L-inf for a larger number of
871  * quadrature points.
872  * @param field The field to compare.
873  * @returns Error in the L2-norm and L-inf norm.
874  */
876 {
877  int NumModes = GetNumExpModes();
878  Array<OneD, NekDouble> L2INF(2);
879 
883  LibUtilities::eGaussRadauMAlpha1Beta0);
888  const LibUtilities::BasisKey BkeyT1(LibUtilities::eModified_A, NumModes,
889  PkeyT1);
890  const LibUtilities::BasisKey BkeyT2(LibUtilities::eModified_B, NumModes,
891  PkeyT2);
892  const LibUtilities::BasisKey BkeyQ1(LibUtilities::eModified_A, NumModes,
893  PkeyQ1);
894  const LibUtilities::BasisKey BkeyQ2(LibUtilities::eModified_A, NumModes,
895  PkeyQ2);
896 
897  LibUtilities::BasisKeyVector Tkeys, Qkeys;
898 
899  // make a copy of the ExpansionInfoMap
900  SpatialDomains::ExpansionInfoMap NewExpInfo = m_graph->GetExpansionInfo();
903  NewExpInfo);
904 
905  // reset new graph with new keys
906  Tkeys.push_back(BkeyT1);
907  Tkeys.push_back(BkeyT2);
908  m_graph->ResetExpansionInfoToBasisKey(ExpInfo, LibUtilities::eTriangle,
909  Tkeys);
910  Qkeys.push_back(BkeyQ1);
911  Qkeys.push_back(BkeyQ2);
912  m_graph->ResetExpansionInfoToBasisKey(ExpInfo, LibUtilities::eQuadrilateral,
913  Qkeys);
914 
917  NewExpInfo);
918 
919  int ErrorCoordim = ErrorExp->GetCoordim(0);
920  int ErrorNq = ErrorExp->GetTotPoints();
921 
922  Array<OneD, NekDouble> ErrorXc0(ErrorNq, 0.0);
923  Array<OneD, NekDouble> ErrorXc1(ErrorNq, 0.0);
924  Array<OneD, NekDouble> ErrorXc2(ErrorNq, 0.0);
925 
926  switch (ErrorCoordim)
927  {
928  case 1:
929  ErrorExp->GetCoords(ErrorXc0);
930  break;
931  case 2:
932  ErrorExp->GetCoords(ErrorXc0, ErrorXc1);
933  break;
934  case 3:
935  ErrorExp->GetCoords(ErrorXc0, ErrorXc1, ErrorXc2);
936  break;
937  }
939  m_session->GetFunction("ExactSolution", field);
940 
941  // Evaluate the exact solution
942  Array<OneD, NekDouble> ErrorSol(ErrorNq);
943 
944  exSol->Evaluate(ErrorXc0, ErrorXc1, ErrorXc2, m_time, ErrorSol);
945 
946  // Calcualte spectral/hp approximation on the quadrature points
947  // of this new expansion basis
948  ErrorExp->BwdTrans(m_fields[field]->GetCoeffs(), ErrorExp->UpdatePhys());
949 
950  L2INF[0] = ErrorExp->L2(ErrorExp->GetPhys(), ErrorSol);
951  L2INF[1] = ErrorExp->Linf(ErrorExp->GetPhys(), ErrorSol);
952 
953  return L2INF;
954 }
955 
956 /**
957  * Set the physical fields based on a restart file, or a function
958  * describing the initial condition given in the session.
959  * @param initialtime Time at which to evaluate the function.
960  * @param dumpInitialConditions Write the initial condition to file?
961  */
963  bool dumpInitialConditions,
964  const int domain)
965 {
966  boost::ignore_unused(initialtime);
967 
969  {
970  ++m_nchk;
971  return;
972  }
973 
974  if (m_session->GetComm()->GetRank() == 0)
975  {
976  cout << "Initial Conditions:" << endl;
977  }
978 
979  if (m_session->DefinesFunction("InitialConditions"))
980  {
981  GetFunction("InitialConditions")
982  ->Evaluate(m_session->GetVariables(), m_fields, m_time, domain);
983  // Enforce C0 Continutiy of initial condiiton
986  {
987  for (int i = 0; i < m_fields.size(); ++i)
988  {
989  m_fields[i]->LocalToGlobal();
990  m_fields[i]->GlobalToLocal();
991  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
992  m_fields[i]->UpdatePhys());
993  }
994  }
995 
996  if (m_session->GetComm()->GetRank() == 0)
997  {
998  for (int i = 0; i < m_fields.size(); ++i)
999  {
1000  std::string varName = m_session->GetVariable(i);
1001  cout << " - Field " << varName << ": "
1002  << GetFunction("InitialConditions")
1003  ->Describe(varName, domain)
1004  << endl;
1005  }
1006  }
1007  }
1008  else
1009  {
1010  int nq = m_fields[0]->GetNpoints();
1011  for (int i = 0; i < m_fields.size(); i++)
1012  {
1013  Vmath::Zero(nq, m_fields[i]->UpdatePhys(), 1);
1014  m_fields[i]->SetPhysState(true);
1015  Vmath::Zero(m_fields[i]->GetNcoeffs(), m_fields[i]->UpdateCoeffs(),
1016  1);
1017  if (m_session->GetComm()->GetRank() == 0)
1018  {
1019  cout << " - Field " << m_session->GetVariable(i)
1020  << ": 0 (default)" << endl;
1021  }
1022  }
1023  }
1024 
1025  if (dumpInitialConditions && m_checksteps && m_nchk == 0 &&
1026  !ParallelInTime())
1027  {
1029  }
1030  else if (dumpInitialConditions && m_nchk == 0 && ParallelInTime())
1031  {
1032  std::string newdir = m_sessionName + ".pit";
1033  if (!fs::is_directory(newdir))
1034  {
1035  fs::create_directory(newdir);
1036  }
1037  if (m_comm->GetTimeComm()->GetRank() == 0)
1038  {
1039  WriteFld(newdir + "/" + m_sessionName + "_0.fld");
1040  }
1041  }
1042  ++m_nchk;
1043 }
1044 
1046  Array<OneD, NekDouble> &outfield,
1047  const NekDouble time)
1048 {
1049  ASSERTL0(outfield.size() == m_fields[field]->GetNpoints(),
1050  "ExactSolution array size mismatch.");
1051  Vmath::Zero(outfield.size(), outfield, 1);
1052  if (m_session->DefinesFunction("ExactSolution"))
1053  {
1054  GetFunction("ExactSolution")
1055  ->Evaluate(m_session->GetVariable(field), outfield, time);
1056  }
1057 }
1058 
1059 /**
1060  * By default, nothing needs initialising at the EquationSystem level.
1061  */
1063 {
1064 }
1065 
1066 /**
1067  *
1068  */
1070 {
1071 }
1072 
1073 /**
1074  * Virtual function to define if operator in DoSolve is
1075  * negated with regard to the strong form. This is currently
1076  * only used in Arnoldi solves. Default is false.
1077  */
1079 {
1080  return false;
1081 }
1082 
1083 /**
1084  *
1085  */
1087 {
1088 }
1089 
1090 /**
1091  *
1092  */
1094 {
1095 }
1096 
1097 /// Virtual function for generating summary information.
1099 {
1100  SessionSummary(l);
1101 }
1102 
1103 /**
1104  * Write the field data to file. The file is named according to the
1105  * session name with the extension .fld appended.
1106  */
1108 {
1109  if (!ParallelInTime())
1110  {
1111  // Serial-in-time
1112  WriteFld(m_sessionName + ".fld");
1113  }
1114  else
1115  {
1116  // Parallel-in-time
1117  std::string newdir = m_sessionName + ".pit";
1118  if (!fs::is_directory(newdir))
1119  {
1120  fs::create_directory(newdir);
1121  }
1122  WriteFld(newdir + "/" + m_sessionName + "_" +
1123  boost::lexical_cast<std::string>(
1124  m_comm->GetTimeComm()->GetRank() + 1) +
1125  ".fld");
1126  }
1127 }
1128 
1129 /**
1130  * Zero the physical fields.
1131  */
1133 {
1134  for (int i = 0; i < m_fields.size(); i++)
1135  {
1136  Vmath::Zero(m_fields[i]->GetNpoints(), m_fields[i]->UpdatePhys(), 1);
1137  }
1138 }
1139 
1140 /**
1141  * FwdTrans the m_fields members
1142  */
1144 {
1145  for (int i = 0; i < m_fields.size(); i++)
1146  {
1147  m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
1148  m_fields[i]->UpdateCoeffs());
1149  m_fields[i]->SetPhysState(false);
1150  }
1151 }
1152 
1153 /**
1154  * Write the n-th checkpoint file.
1155  * @param n The index of the checkpoint file.
1156  */
1158 {
1159  if (!ParallelInTime())
1160  {
1161  // Serial-in-time
1162  std::string outname =
1163  m_sessionName + "_" + boost::lexical_cast<std::string>(n);
1164  WriteFld(outname + ".chk");
1165  }
1166  else
1167  {
1168  // Parallel-in-time
1169  std::string paradir = m_sessionName + "_" +
1170  boost::lexical_cast<std::string>(m_pararealIter) +
1171  ".pit";
1172  if (!fs::is_directory(paradir))
1173  {
1174  fs::create_directory(paradir);
1175  }
1176  std::string outname = paradir + "/" + m_sessionName + "_" +
1177  boost::lexical_cast<std::string>(n);
1178  WriteFld(outname + ".chk");
1179  }
1180 }
1181 
1182 /**
1183  * Write the n-th checkpoint file.
1184  * @param n The index of the checkpoint file.
1185  */
1187  const int n, MultiRegions::ExpListSharedPtr &field,
1188  std::vector<Array<OneD, NekDouble>> &fieldcoeffs,
1189  std::vector<std::string> &variables)
1190 {
1191  if (!ParallelInTime())
1192  {
1193  // Serial-in-time
1194  std::string outname =
1195  m_sessionName + "_" + boost::lexical_cast<std::string>(n);
1196  WriteFld(outname, field, fieldcoeffs, variables);
1197  }
1198  else
1199  {
1200  // Parallel-in-time
1201  std::string paradir = m_sessionName + "_" +
1202  boost::lexical_cast<std::string>(m_pararealIter) +
1203  ".pit";
1204  if (!fs::is_directory(paradir))
1205  {
1206  fs::create_directory(paradir);
1207  }
1208  std::string outname = paradir + "/" + m_sessionName + "_" +
1209  boost::lexical_cast<std::string>(n);
1210  WriteFld(outname, field, fieldcoeffs, variables);
1211  }
1212 }
1213 
1214 /**
1215  * Write the n-th base flow into a .chk file
1216  * @param n The index of the base flow file.
1217  */
1219 {
1220  std::string outname =
1221  m_sessionName + "_BaseFlow_" + boost::lexical_cast<std::string>(n);
1222 
1223  WriteFld(outname + ".chk");
1224 }
1225 
1226 /**
1227  * Writes the field data to a file with the given filename.
1228  * @param outname Filename to write to.
1229  */
1230 void EquationSystem::WriteFld(const std::string &outname)
1231 {
1232  std::vector<Array<OneD, NekDouble>> fieldcoeffs(m_fields.size());
1233  std::vector<std::string> variables(m_fields.size());
1234 
1235  for (int i = 0; i < m_fields.size(); ++i)
1236  {
1237  if (m_fields[i]->GetNcoeffs() == m_fields[0]->GetNcoeffs())
1238  {
1239  fieldcoeffs[i] = m_fields[i]->UpdateCoeffs();
1240  }
1241  else
1242  {
1243  fieldcoeffs[i] = Array<OneD, NekDouble>(m_fields[0]->GetNcoeffs());
1244  m_fields[0]->ExtractCoeffsToCoeffs(
1245  m_fields[i], m_fields[i]->GetCoeffs(), fieldcoeffs[i]);
1246  }
1247  variables[i] = m_boundaryConditions->GetVariable(i);
1248  }
1249 
1250  ExtraFldOutput(fieldcoeffs, variables);
1251 
1252  WriteFld(outname, m_fields[0], fieldcoeffs, variables);
1253 }
1254 
1255 /**
1256  * Writes the field data to a file with the given filename.
1257  * @param outname Filename to write to.
1258  * @param field ExpList on which data is based.
1259  * @param fieldcoeffs An array of array of expansion coefficients.
1260  * @param variables An array of variable names.
1261  */
1262 void EquationSystem::WriteFld(const std::string &outname,
1264  std::vector<Array<OneD, NekDouble>> &fieldcoeffs,
1265  std::vector<std::string> &variables)
1266 {
1267  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef =
1268  field->GetFieldDefinitions();
1269  std::vector<std::vector<NekDouble>> FieldData(FieldDef.size());
1270 
1271  // Copy Data into FieldData and set variable
1272  for (int j = 0; j < fieldcoeffs.size(); ++j)
1273  {
1274  for (int i = 0; i < FieldDef.size(); ++i)
1275  {
1276  // Could do a search here to find correct variable
1277  FieldDef[i]->m_fields.push_back(variables[j]);
1278  field->AppendFieldData(FieldDef[i], FieldData[i], fieldcoeffs[j]);
1279  }
1280  }
1281 
1282  // Update time in field info if required
1283  if (m_fieldMetaDataMap.find("Time") != m_fieldMetaDataMap.end())
1284  {
1285  m_fieldMetaDataMap["Time"] = boost::lexical_cast<std::string>(m_time);
1286  }
1287 
1288  // Update step in field info if required
1289  if (m_fieldMetaDataMap.find("ChkFileNum") != m_fieldMetaDataMap.end())
1290  {
1291  m_fieldMetaDataMap["ChkFileNum"] =
1292  boost::lexical_cast<std::string>(m_nchk);
1293  }
1294 
1295  // If necessary, add mapping information to metadata
1296  // and output mapping coordinates
1298  fields[0] = field;
1302  mapping->Output(fieldMetaDataMap, outname);
1303 
1304  // If necessary, add informaton for moving frame reference to metadata
1305  if (m_fieldMetaDataMap.find("Theta_x") != m_fieldMetaDataMap.end())
1306  {
1307  // if one theta exists, add all three thetas
1308  std::vector<std::string> vSuffix = {"_x", "_y", "_z"};
1309  for (int i = 0; i < 3; ++i)
1310  {
1311  std::string sTheta = "Theta" + vSuffix[i];
1312  m_fieldMetaDataMap[sTheta] =
1313  boost::lexical_cast<std::string>(m_movingFrameTheta[i]);
1314  }
1315  }
1316 
1317  m_fld->Write(outname, FieldDef, FieldData, fieldMetaDataMap,
1318  m_session->GetBackups());
1319 }
1320 
1321 /**
1322  * Import field from infile and load into \a m_fields. This routine will
1323  * also perform a \a BwdTrans to ensure data is in both the physical and
1324  * coefficient storage.
1325  * @param infile Filename to read.
1326  */
1328  const std::string &infile,
1330 {
1331  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
1332  std::vector<std::vector<NekDouble>> FieldData;
1333  LibUtilities::FieldIOSharedPtr field_fld =
1335  field_fld->Import(infile, FieldDef, FieldData);
1336 
1337  // Copy FieldData into m_fields
1338  for (int j = 0; j < pFields.size(); ++j)
1339  {
1340  Vmath::Zero(pFields[j]->GetNcoeffs(), pFields[j]->UpdateCoeffs(), 1);
1341 
1342  for (int i = 0; i < FieldDef.size(); ++i)
1343  {
1344  ASSERTL1(FieldDef[i]->m_fields[j] == m_session->GetVariable(j),
1345  std::string("Order of ") + infile +
1346  std::string(" data and that defined in "
1347  "m_boundaryconditions differs"));
1348 
1349  pFields[j]->ExtractDataToCoeffs(FieldDef[i], FieldData[i],
1350  FieldDef[i]->m_fields[j],
1351  pFields[j]->UpdateCoeffs());
1352  }
1353  pFields[j]->BwdTrans(pFields[j]->GetCoeffs(), pFields[j]->UpdatePhys());
1354  }
1355 }
1356 
1357 /**
1358  * Import field from infile and load into \a m_fields. This routine will
1359  * also perform a \a BwdTrans to ensure data is in both the physical and
1360  * coefficient storage.
1361  * @param infile Filename to read.
1362  * If optionan \a ndomains is specified it assumes we loop over nodmains
1363  * for each nvariables.
1364  */
1366  const std::string &infile,
1367  Array<OneD, MultiRegions::ExpListSharedPtr> &pFields, const int ndomains)
1368 {
1369  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
1370  std::vector<std::vector<NekDouble>> FieldData;
1371 
1372  LibUtilities::Import(infile, FieldDef, FieldData);
1373 
1374  int nvariables = GetNvariables();
1375 
1376  ASSERTL0(
1377  ndomains * nvariables == pFields.size(),
1378  "Number of fields does not match the number of variables and domains");
1379 
1380  // Copy FieldData into m_fields
1381  for (int j = 0; j < ndomains; ++j)
1382  {
1383  for (int i = 0; i < nvariables; ++i)
1384  {
1385  Vmath::Zero(pFields[j * nvariables + i]->GetNcoeffs(),
1386  pFields[j * nvariables + i]->UpdateCoeffs(), 1);
1387 
1388  for (int n = 0; n < FieldDef.size(); ++n)
1389  {
1390  ASSERTL1(FieldDef[n]->m_fields[i] == m_session->GetVariable(i),
1391  std::string("Order of ") + infile +
1392  std::string(" data and that defined in "
1393  "m_boundaryconditions differs"));
1394 
1395  pFields[j * nvariables + i]->ExtractDataToCoeffs(
1396  FieldDef[n], FieldData[n], FieldDef[n]->m_fields[i],
1397  pFields[j * nvariables + i]->UpdateCoeffs());
1398  }
1399  pFields[j * nvariables + i]->BwdTrans(
1400  pFields[j * nvariables + i]->GetCoeffs(),
1401  pFields[j * nvariables + i]->UpdatePhys());
1402  }
1403  }
1404 }
1405 
1406 /**
1407  * Import field from infile and load into \a pField. This routine will
1408  * also perform a \a BwdTrans to ensure data is in both the physical and
1409  * coefficient storage.
1410  */
1411 void EquationSystem::ImportFld(const std::string &infile,
1413  std::string &pFieldName)
1414 {
1415  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
1416  std::vector<std::vector<NekDouble>> FieldData;
1417 
1418  LibUtilities::FieldIOSharedPtr field_fld =
1420  field_fld->Import(infile, FieldDef, FieldData);
1421  int idx = -1;
1422 
1423  Vmath::Zero(pField->GetNcoeffs(), pField->UpdateCoeffs(), 1);
1424 
1425  for (int i = 0; i < FieldDef.size(); ++i)
1426  {
1427  // find the index of the required field in the file.
1428  for (int j = 0; j < FieldData.size(); ++j)
1429  {
1430  if (FieldDef[i]->m_fields[j] == pFieldName)
1431  {
1432  idx = j;
1433  }
1434  }
1435  ASSERTL1(idx >= 0, "Field " + pFieldName + " not found.");
1436 
1437  pField->ExtractDataToCoeffs(FieldDef[i], FieldData[i],
1438  FieldDef[i]->m_fields[idx],
1439  pField->UpdateCoeffs());
1440  }
1441  pField->BwdTrans(pField->GetCoeffs(), pField->UpdatePhys());
1442 }
1443 
1444 /**
1445  * Import field from infile and load into the array \a coeffs.
1446  *
1447  * @param infile Filename to read.
1448  * @param fieldStr an array of string identifying fields to be imported
1449  * @param coeffs array of array of coefficients to store imported data
1450  */
1451 void EquationSystem::ImportFld(const std::string &infile,
1452  std::vector<std::string> &fieldStr,
1453  Array<OneD, Array<OneD, NekDouble>> &coeffs)
1454 {
1455 
1456  ASSERTL0(fieldStr.size() <= coeffs.size(),
1457  "length of fieldstr should be the same as pFields");
1458 
1459  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
1460  std::vector<std::vector<NekDouble>> FieldData;
1461 
1462  LibUtilities::FieldIOSharedPtr field_fld =
1464  field_fld->Import(infile, FieldDef, FieldData);
1465 
1466  // Copy FieldData into m_fields
1467  for (int j = 0; j < fieldStr.size(); ++j)
1468  {
1469  Vmath::Zero(coeffs[j].size(), coeffs[j], 1);
1470  for (int i = 0; i < FieldDef.size(); ++i)
1471  {
1472  m_fields[0]->ExtractDataToCoeffs(FieldDef[i], FieldData[i],
1473  fieldStr[j], coeffs[j]);
1474  }
1475  }
1476 }
1477 
1478 /**
1479  * Write out a summary of the session data.
1480  * @param out Output stream to write data to.
1481  */
1483 {
1484  AddSummaryItem(s, "EquationType", m_session->GetSolverInfo("EQTYPE"));
1485  AddSummaryItem(s, "Session Name", m_sessionName);
1486  AddSummaryItem(s, "Spatial Dim.", m_spacedim);
1487  AddSummaryItem(s, "Max SEM Exp. Order",
1488  m_fields[0]->EvalBasisNumModesMax());
1489 
1490  if (m_session->GetComm()->GetSize() > 1)
1491  {
1492  AddSummaryItem(s, "Num. Processes", m_session->GetComm()->GetSize());
1493  }
1494 
1496  {
1497  AddSummaryItem(s, "Quasi-3D", "Homogeneous in z-direction");
1498  AddSummaryItem(s, "Expansion Dim.", m_expdim + 1);
1499  AddSummaryItem(s, "Num. Hom. Modes (z)", m_npointsZ);
1500  AddSummaryItem(s, "Hom. length (LZ)", m_LhomZ);
1501  AddSummaryItem(s, "FFT Type", m_useFFT ? "FFTW" : "MVM");
1502  if (m_halfMode)
1503  {
1504  AddSummaryItem(s, "ModeType", "Half Mode");
1505  }
1506  else if (m_singleMode)
1507  {
1508  AddSummaryItem(s, "ModeType", "Single Mode");
1509  }
1510  else if (m_multipleModes)
1511  {
1512  AddSummaryItem(s, "ModeType", "Multiple Modes");
1513  }
1514  }
1515  else if (m_HomogeneousType == eHomogeneous2D)
1516  {
1517  AddSummaryItem(s, "Quasi-3D", "Homogeneous in yz-plane");
1518  AddSummaryItem(s, "Expansion Dim.", m_expdim + 2);
1519  AddSummaryItem(s, "Num. Hom. Modes (y)", m_npointsY);
1520  AddSummaryItem(s, "Num. Hom. Modes (z)", m_npointsZ);
1521  AddSummaryItem(s, "Hom. length (LY)", m_LhomY);
1522  AddSummaryItem(s, "Hom. length (LZ)", m_LhomZ);
1523  AddSummaryItem(s, "FFT Type", m_useFFT ? "FFTW" : "MVM");
1524  }
1525  else
1526  {
1527  AddSummaryItem(s, "Expansion Dim.", m_expdim);
1528  }
1529 
1530  if (m_session->DefinesSolverInfo("UpwindType"))
1531  {
1532  AddSummaryItem(s, "Riemann Solver",
1533  m_session->GetSolverInfo("UpwindType"));
1534  }
1535 
1536  if (m_session->DefinesSolverInfo("AdvectionType"))
1537  {
1538  std::string AdvectionType;
1539  AdvectionType = m_session->GetSolverInfo("AdvectionType");
1541  s, "Advection Type",
1542  GetAdvectionFactory().GetClassDescription(AdvectionType));
1543  }
1544 
1546  {
1547  AddSummaryItem(s, "Projection Type", "Continuous Galerkin");
1548  }
1550  {
1551  AddSummaryItem(s, "Projection Type", "Discontinuous Galerkin");
1552  }
1554  {
1555  AddSummaryItem(s, "Projection Type",
1556  "Mixed Continuous Galerkin and Discontinuous");
1557  }
1558 
1559  if (m_session->DefinesSolverInfo("DiffusionType"))
1560  {
1561  std::string DiffusionType;
1562  DiffusionType = m_session->GetSolverInfo("DiffusionType");
1564  s, "Diffusion Type",
1565  GetDiffusionFactory().GetClassDescription(DiffusionType));
1566  }
1567 }
1568 
1570 {
1571  return Array<OneD, bool>(m_session->GetVariables().size(), false);
1572 }
1573 
1575 {
1576  ASSERTL0(false, "This function is not valid for the Base class");
1578  return null;
1579 }
1580 
1582  std::vector<Array<OneD, NekDouble>> &fieldcoeffs,
1583  std::vector<std::string> &variables)
1584 {
1585  boost::ignore_unused(fieldcoeffs, variables);
1586 }
1587 
1588 } // namespace SolverUtils
1589 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
static GLOBAL_MAPPING_EXPORT MappingSharedPtr Load(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
Return a pointer to the mapping, creating it on first call.
Definition: Mapping.cpp:270
Describes the specification for a Basis.
Definition: Basis.h:50
static std::shared_ptr< FieldIO > CreateForFile(const LibUtilities::SessionReaderSharedPtr session, const std::string &filename)
Construct a FieldIO object for a given input filename.
Definition: FieldIO.cpp:226
static std::shared_ptr< FieldIO > CreateDefault(const LibUtilities::SessionReaderSharedPtr session)
Returns an object for the default FieldIO method.
Definition: FieldIO.cpp:197
Provides a generic Factory class.
Definition: NekFactory.hpp:105
Defines a specification for a set of points.
Definition: Points.h:59
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
virtual SOLVER_UTILS_EXPORT void v_DoInitialise()
Virtual function for initialisation implementation.
int m_spacedim
Spatial dimension (>= expansion dim).
SOLVER_UTILS_EXPORT int GetTraceNpoints()
bool m_useFFT
Flag to determine if FFT is used for homogeneous transform.
SOLVER_UTILS_EXPORT int GetNvariables()
int m_expdim
Expansion dimension.
LibUtilities::FieldIOSharedPtr m_fld
Field input/output.
SOLVER_UTILS_EXPORT void SessionSummary(SummaryList &vSummary)
Write out a session summary.
SOLVER_UTILS_EXPORT bool ParallelInTime()
Check if solver use Parallel-in-Time.
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
virtual SOLVER_UTILS_EXPORT ~EquationSystem()
Destructor.
SpatialDomains::MeshGraphSharedPtr m_graph
Pointer to graph defining mesh.
LibUtilities::CommSharedPtr m_comm
Communicator.
NekDouble m_timestep
Time step size.
SOLVER_UTILS_EXPORT void ImportFld(const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
Input field data from the given file.
NekDouble m_time
Current time of simulation.
bool m_multipleModes
Flag to determine if use multiple homogenenous modes are used.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff()
Virtual function for transformation to coefficient space.
NekDouble m_fintime
Finish time of the simulation.
virtual SOLVER_UTILS_EXPORT NekDouble v_L2Error(unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
Virtual function for the L_2 error computation between fields and a given exact solution.
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure(void)
int m_npointsZ
number of points in Z direction (if homogeneous)
SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n)
Write checkpoint file of m_fields.
std::map< std::string, SolverUtils::SessionFunctionSharedPtr > m_sessionFunctions
Map of known SessionFunctions.
NekDouble m_checktime
Time between checkpoints.
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareFeld=true)
Initialisation object for EquationSystem.
virtual SOLVER_UTILS_EXPORT NekDouble v_LinfError(unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
Virtual function for the L_inf error computation between fields and a given exact solution.
bool m_useInitialCondition
Flag to determine if IC are used.
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution(unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
SOLVER_UTILS_EXPORT void WriteFld(const std::string &outname)
Write field data to the given filename.
SOLVER_UTILS_EXPORT void ImportFldToMultiDomains(const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int ndomains)
Input field data from the given file to multiple domains.
NekDouble m_LhomY
physical length in Y direction (if homogeneous)
int m_npointsY
number of points in Y direction (if homogeneous)
std::string m_sessionName
Name of the session.
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow(const int n)
Write base flow file of m_fields.
bool m_specHP_dealiasing
Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
bool m_singleMode
Flag to determine if single homogeneous mode is used.
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp()
Virtual function to identify if operator is negated in DoSolve.
SOLVER_UTILS_EXPORT void ZeroPhysFields()
int m_pararealIter
Number of parareal time iteration.
int m_HomoDirec
number of homogenous directions
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput(std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables)
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
enum HomogeneousType m_HomogeneousType
SOLVER_UTILS_EXPORT Array< OneD, NekDouble > ErrorExtraPoints(unsigned int field)
Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf].
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, bool > m_checkIfSystemSingular
Flag to indicate if the fields should be checked for singularity.
bool m_homogen_dealiasing
Flag to determine if dealiasing is used for homogeneous simulations.
SOLVER_UTILS_EXPORT int GetNcoeffs()
int m_nchk
Number of checkpoints written so far.
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
virtual SOLVER_UTILS_EXPORT void v_DoSolve()
Virtual function for solve implementation.
virtual SOLVER_UTILS_EXPORT Array< OneD, bool > v_GetSystemSingularChecks()
int m_steps
Number of steps to take.
int m_NumQuadPointsError
Number of Quadrature points used to work out the error.
Array< OneD, NekDouble > m_movingFrameTheta
Moving frame of reference angles with respect to the.
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
SOLVER_UTILS_EXPORT int GetNumExpModes()
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
Pointer to boundary conditions object.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
NekDouble m_LhomZ
physical length in Z direction (if homogeneous)
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &l)
Virtual function for generating summary information.
SOLVER_UTILS_EXPORT void ExtraFldOutput(std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables)
bool m_halfMode
Flag to determine if half homogeneous mode is used.
SOLVER_UTILS_EXPORT void FwdTransFields()
int m_checksteps
Number of steps between checkpoints.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
virtual SOLVER_UTILS_EXPORT void v_Output(void)
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys()
Virtual function for transformation to physical space.
GLOBAL_MAPPING_EXPORT typedef std::shared_ptr< Mapping > MappingSharedPtr
A shared pointer to a Mapping object.
Definition: Mapping.h:50
std::vector< BasisKey > BasisKeyVector
Name for a vector of BasisKeys.
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:327
void Import(const std::string &infilename, std::vector< FieldDefinitionsSharedPtr > &fielddefs, std::vector< std::vector< NekDouble >> &fielddata, FieldMetaDataMap &fieldinfomap, const Array< OneD, int > &ElementIDs)
This function allows for data to be imported from an FLD file when a session and/or communicator is n...
Definition: FieldIO.cpp:290
std::map< std::string, std::string > FieldMetaDataMap
Definition: FieldIO.h:52
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:129
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:53
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:76
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:53
@ eFourierSingleModeSpaced
1D Non Evenly-spaced points for Single Mode analysis
Definition: PointsType.h:77
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:51
@ eFourierSingleMode
Fourier ModifiedExpansion with just the first mode .
Definition: BasisType.h:66
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50
@ eFourierHalfModeIm
Fourier Modified expansions with just the imaginary part of the first mode .
Definition: BasisType.h:70
@ eFourierHalfModeRe
Fourier Modified expansions with just the real part of the first mode .
Definition: BasisType.h:68
@ eFourier
Fourier Expansion .
Definition: BasisType.h:57
std::shared_ptr< DisContField > DisContFieldSharedPtr
Definition: DisContField.h:341
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< ContField > ContFieldSharedPtr
Definition: ContField.h:289
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:48
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41
EquationSystemFactory & GetEquationSystemFactory()
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:49
std::shared_ptr< SessionFunction > SessionFunctionSharedPtr
std::shared_ptr< ExpansionInfoMap > ExpansionInfoMapShPtr
Definition: MeshGraph.h:145
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
std::map< int, ExpansionInfoShPtr > ExpansionInfoMap
Definition: MeshGraph.h:143
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
NekMatrix< NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag >, BlockMatrixTag > DNekScalBlkMat
Definition: NekTypeDefs.hpp:68
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
double NekDouble
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294
Used to lookup the create function in NekManager.
Definition: MatrixKey.h:71