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 
40 #include <LocalRegions/MatrixKey.h>
47 
48 #include <MultiRegions/ExpList2D.h> // for ExpList2D, etc
49 #include <MultiRegions/ExpList3D.h> // for ExpList3D
52 
55 
56 #include <GlobalMapping/Mapping.h>
57 
58 #include <boost/format.hpp>
59 
60 #include <iostream>
61 #include <string>
62 
63 using namespace std;
64 
65 namespace Nektar
66 {
67  namespace SolverUtils
68  {
69 
70  std::string EquationSystem::equationSystemTypeLookupIds[2] = {
71  LibUtilities::SessionReader::RegisterEnumValue("DEALIASING",
72  "True", 0),
73  LibUtilities::SessionReader::RegisterEnumValue("DEALIASING",
74  "False", 1)};
75 
76  /**
77  * @class EquationSystem
78  *
79  * This class is a base class for all solver implementations. It
80  * provides the underlying generic functionality and interface for
81  * solving equations.
82  *
83  * To solve a steady-state equation, create a derived class from this
84  * class and reimplement the virtual functions to provide custom
85  * implementation for the problem.
86  *
87  * To solve unsteady problems, derive from the UnsteadySystem class
88  * instead which provides general time integration.
89  */
91  {
92  static EquationSystemFactory instance;
93  return instance;
94  }
95 
96  /**
97  * This constructor is protected as the objects of this class are never
98  * instantiated directly.
99  * @param pSession The session reader holding problem parameters.
100  */
101  EquationSystem::EquationSystem(
102  const LibUtilities::SessionReaderSharedPtr& pSession,
104  : m_comm (pSession->GetComm()),
105  m_session (pSession),
106  m_graph (pGraph),
107  m_lambda (0),
108  m_fieldMetaDataMap(LibUtilities::NullFieldMetaDataMap)
109  {
110  // set up session names in fieldMetaDataMap
111  const vector<std::string> filenames = m_session->GetFilenames();
112 
113  for(int i = 0; i < filenames.size(); ++i)
114  {
115  string sessionname = "SessionName";
116  sessionname += boost::lexical_cast<std::string>(i);
117  m_fieldMetaDataMap[sessionname] = filenames[i];
118  m_fieldMetaDataMap["ChkFileNum"] =
119  boost::lexical_cast<std::string>(0);
120  }
121 
122  }
123 
124  /**
125  * @brief Initialisation object for EquationSystem.
126  */
128  {
129  // Save the basename of input file name for output details
130  m_sessionName = m_session->GetSessionName();
131 
132  // Instantiate a field reader/writer
134 
135  // Also read and store the boundary conditions
139 
140  // Set space dimension for use in class
141  m_spacedim = m_graph->GetSpaceDimension();
142 
143  // Setting parameteres for homogenous problems
144  m_HomoDirec = 0;
145  m_useFFT = false;
146  m_homogen_dealiasing = false;
147  m_singleMode = false;
148  m_halfMode = false;
149  m_multipleModes = false;
151 
152  if (m_session->DefinesSolverInfo("HOMOGENEOUS"))
153  {
154  std::string HomoStr = m_session->GetSolverInfo("HOMOGENEOUS");
155  m_spacedim = 3;
156 
157  if ((HomoStr == "HOMOGENEOUS1D") || (HomoStr == "Homogeneous1D")
158  || (HomoStr == "1D") || (HomoStr == "Homo1D"))
159  {
161  m_session->LoadParameter("LZ", m_LhomZ);
162  m_HomoDirec = 1;
163 
164  if(m_session->DefinesSolverInfo("ModeType"))
165  {
166  m_session->MatchSolverInfo("ModeType", "SingleMode",
167  m_singleMode, false);
168  m_session->MatchSolverInfo("ModeType", "HalfMode",
169  m_halfMode, false);
170  m_session->MatchSolverInfo("ModeType", "MultipleModes",
171  m_multipleModes, false);
172  }
173 
174  // Stability Analysis flags
175  if (m_session->DefinesSolverInfo("ModeType"))
176  {
177  if(m_singleMode)
178  {
179  m_npointsZ = 2;
180  }
181  else if(m_halfMode)
182  {
183  m_npointsZ = 1;
184  }
185  else if(m_multipleModes)
186  {
187  m_npointsZ = m_session->GetParameter("HomModesZ");
188  }
189  else
190  {
191  ASSERTL0(false, "SolverInfo ModeType not valid");
192  }
193  }
194  else
195  {
196  m_npointsZ = m_session->GetParameter("HomModesZ");
197  }
198  }
199 
200  if ((HomoStr == "HOMOGENEOUS2D") || (HomoStr == "Homogeneous2D")
201  || (HomoStr == "2D") || (HomoStr == "Homo2D"))
202  {
204  m_session->LoadParameter("HomModesY", m_npointsY);
205  m_session->LoadParameter("LY", m_LhomY);
206  m_session->LoadParameter("HomModesZ", m_npointsZ);
207  m_session->LoadParameter("LZ", m_LhomZ);
208  m_HomoDirec = 2;
209  }
210 
211  if ((HomoStr == "HOMOGENEOUS3D") || (HomoStr == "Homogeneous3D")
212  || (HomoStr == "3D") || (HomoStr == "Homo3D"))
213  {
215  m_session->LoadParameter("HomModesY", m_npointsY);
216  m_session->LoadParameter("LY", m_LhomY);
217  m_session->LoadParameter("HomModesZ", m_npointsZ);
218  m_session->LoadParameter("LZ", m_LhomZ);
219  m_HomoDirec = 2;
220  }
221 
222  m_session->MatchSolverInfo("USEFFT", "FFTW", m_useFFT, false);
223 
224  m_session->MatchSolverInfo("DEALIASING", "True",
225  m_homogen_dealiasing, false);
226  }
227  else
228  {
229  // set to default value so can use to identify 2d or 3D
230  // (homogeneous) expansions
231  m_npointsZ = 1;
232  }
233 
234  m_session->MatchSolverInfo("SPECTRALHPDEALIASING", "True",
235  m_specHP_dealiasing, false);
236  if (m_specHP_dealiasing == false)
237  {
238  m_session->MatchSolverInfo("SPECTRALHPDEALIASING", "On",
239  m_specHP_dealiasing, false);
240  }
241 
242  // Options to determine type of projection from file or directly
243  // from constructor
244  if (m_session->DefinesSolverInfo("PROJECTION"))
245  {
246  std::string ProjectStr = m_session->GetSolverInfo("PROJECTION");
247 
248  if ((ProjectStr == "Continuous") || (ProjectStr == "Galerkin") ||
249  (ProjectStr == "CONTINUOUS") || (ProjectStr == "GALERKIN"))
250  {
252  }
253  else if ((ProjectStr == "MixedCGDG") ||
254  (ProjectStr == "Mixed_CG_Discontinuous"))
255  {
257  }
258  else if(ProjectStr == "DisContinuous")
259  {
261  }
262  else
263  {
264  ASSERTL0(false,"PROJECTION value not recognised");
265  }
266  }
267  else
268  {
269  cerr << "Projection type not specified in SOLVERINFO,"
270  "defaulting to continuous Galerkin" << endl;
272  }
273 
274  // Enforce singularity check for some problems
276 
277  int i;
278  int nvariables = m_session->GetVariables().size();
279  bool DeclareCoeffPhysArrays = true;
280 
281 
283  m_spacedim = m_graph->GetSpaceDimension()+m_HomoDirec;
284  m_expdim = m_graph->GetMeshDimension();
285 
286  /// Continuous field
289  {
290  switch(m_expdim)
291  {
292  case 1:
293  {
296  {
297  const LibUtilities::PointsKey PkeyY
299  const LibUtilities::BasisKey BkeyY
301  const LibUtilities::PointsKey PkeyZ
304  BkeyZ(LibUtilities::eFourier, m_npointsZ, PkeyZ);
305 
306  for (i = 0; i < m_fields.num_elements(); i++)
307  {
308  m_fields[i] = MemoryManager<MultiRegions
309  ::ContField3DHomogeneous2D>
310  ::AllocateSharedPtr(
311  m_session, BkeyY, BkeyZ, m_LhomY,
312  m_LhomZ, m_useFFT,
314  m_session->GetVariable(i));
315  }
316  }
317  else
318  {
319  for (i = 0; i < m_fields.num_elements(); i++)
320  {
323  AllocateSharedPtr(
325  m_session->GetVariable(i));
326  }
327  }
328  break;
329  }
330  case 2:
331  {
333  {
334  // Fourier single mode stability analysis
335  if (m_singleMode)
336  {
337  const LibUtilities::PointsKey PkeyZ(
338  m_npointsZ,
340 
341  const LibUtilities::BasisKey BkeyZ(
343  m_npointsZ,
344  PkeyZ);
345 
346  for(i = 0; i < m_fields.num_elements(); i++)
347  {
348  m_fields[i] = MemoryManager<MultiRegions
349  ::ContField3DHomogeneous1D>
350  ::AllocateSharedPtr(
351  m_session, BkeyZ, m_LhomZ,
353  m_graph,
354  m_session->GetVariable(i),
356  }
357  }
358  // Half mode stability analysis
359  else if(m_halfMode)
360  {
361  const LibUtilities::PointsKey PkeyZ(
362  m_npointsZ,
364 
365  const LibUtilities::BasisKey BkeyZR(
367  m_npointsZ, PkeyZ);
368 
369  const LibUtilities::BasisKey BkeyZI(
371  m_npointsZ, PkeyZ);
372 
373 
374  for (i = 0; i < m_fields.num_elements(); i++)
375  {
376  if(m_session->GetVariable(i).compare("w")
377  == 0)
378  {
379  m_fields[i] = MemoryManager<MultiRegions
380  ::ContField3DHomogeneous1D>
381  ::AllocateSharedPtr(
382  m_session, BkeyZI, m_LhomZ,
383  m_useFFT,
385  m_graph,
386  m_session->GetVariable(i),
388  }
389  else
390  {
391  m_fields[i] = MemoryManager<MultiRegions
392  ::ContField3DHomogeneous1D>
393  ::AllocateSharedPtr(
394  m_session, BkeyZR, m_LhomZ,
396  m_graph,
397  m_session->GetVariable(i),
399  }
400 
401 
402  }
403  }
404  // Normal homogeneous 1D
405  else
406  {
407  const LibUtilities::PointsKey PkeyZ(
408  m_npointsZ,
410  const LibUtilities::BasisKey BkeyZ(
412 
413  for (i = 0; i < m_fields.num_elements(); i++)
414  {
415  m_fields[i] = MemoryManager<MultiRegions
416  ::ContField3DHomogeneous1D>
417  ::AllocateSharedPtr(
418  m_session, BkeyZ, m_LhomZ,
420  m_graph,
421  m_session->GetVariable(i),
423  }
424  }
425  }
426  else
427  {
428  i = 0;
430  firstfield = MemoryManager<MultiRegions::
431  ContField2D>::AllocateSharedPtr(
433  m_session->GetVariable(i),
434  DeclareCoeffPhysArrays,
436  m_fields[0] = firstfield;
437  for (i = 1; i < m_fields.num_elements(); i++)
438  {
439  if (m_graph->
440  SameExpansions(m_session->GetVariable(0),
441  m_session->GetVariable(i)))
442  {
443  m_fields[i] = MemoryManager<MultiRegions::
444  ContField2D>::AllocateSharedPtr(
445  *firstfield, m_graph,
446  m_session->GetVariable(i),
447  DeclareCoeffPhysArrays,
449  }
450  else
451  {
452  m_fields[i] = MemoryManager<MultiRegions
453  ::ContField2D>::AllocateSharedPtr(
454  m_session, m_graph,
455  m_session->GetVariable(i),
456  DeclareCoeffPhysArrays,
458  }
459  }
460 
461  if (m_projectionType ==
463  {
464  /// Setting up the normals
467  (m_spacedim);
468 
469  for (i = 0; i < m_spacedim; ++i)
470  {
472  (GetTraceNpoints());
473  }
474 
475  m_fields[0]->GetTrace()->
476  GetNormals(m_traceNormals);
477  }
478 
479  }
480 
481  break;
482  }
483  case 3:
484  {
485  i = 0;
489  m_session->GetVariable(i),
491 
492  m_fields[0] = firstfield;
493  for (i = 1; i < m_fields.num_elements(); i++)
494  {
495  if(m_graph->SameExpansions(
496  m_session->GetVariable(0),
497  m_session->GetVariable(i)))
498  {
499  m_fields[i] = MemoryManager<MultiRegions
500  ::ContField3D>::AllocateSharedPtr(
501  *firstfield, m_graph,
502  m_session->GetVariable(i),
504  }
505  else
506  {
507  m_fields[i] = MemoryManager<MultiRegions
508  ::ContField3D>::AllocateSharedPtr(
509  m_session, m_graph,
510  m_session->GetVariable(i),
512  }
513  }
514 
515  if (m_projectionType ==
517  {
518  /// Setting up the normals
521  (m_spacedim);
522  for(i = 0; i < m_spacedim; ++i)
523  {
524  m_traceNormals[i] =
526  }
527 
528  m_fields[0]->GetTrace()->GetNormals(m_traceNormals);
529  // Call the trace on all fields to ensure DG setup.
530  for(i = 1; i < m_fields.num_elements(); ++i)
531  {
532  m_fields[i]->GetTrace();
533  }
534  }
535  break;
536  }
537  default:
538  ASSERTL0(false,"Expansion dimension not recognised");
539  break;
540  }
541  }
542  // Discontinuous field
543  else
544  {
545  switch(m_expdim)
546  {
547  case 1:
548  {
551  {
552  const LibUtilities::PointsKey PkeyY(
554  const LibUtilities::BasisKey BkeyY(
556  const LibUtilities::PointsKey PkeyZ(
558  const LibUtilities::BasisKey BkeyZ(
560 
561  for (i = 0; i < m_fields.num_elements(); i++)
562  {
563  m_fields[i] = MemoryManager<MultiRegions
564  ::DisContField3DHomogeneous2D>
565  ::AllocateSharedPtr(
566  m_session, BkeyY, BkeyZ, m_LhomY,
567  m_LhomZ, m_useFFT,
569  m_session->GetVariable(i));
570  }
571  }
572  else
573  {
574  for (i = 0; i < m_fields.num_elements(); i++)
575  {
576  m_fields[i] = MemoryManager<MultiRegions::
577  DisContField1D>::AllocateSharedPtr(
579  m_session->GetVariable(i));
580  }
581  }
582 
583  break;
584  }
585  case 2:
586  {
588  {
589  const LibUtilities::PointsKey PkeyZ(
591  const LibUtilities::BasisKey BkeyZ(
593 
594  for (i = 0; i < m_fields.num_elements(); i++)
595  {
596  m_fields[i] = MemoryManager<MultiRegions
597  ::DisContField3DHomogeneous1D>
598  ::AllocateSharedPtr(
599  m_session, BkeyZ, m_LhomZ, m_useFFT,
601  m_session->GetVariable(i));
602  }
603  }
604  else
605  {
606  for (i = 0; i < m_fields.num_elements(); i++)
607  {
608  m_fields[i] = MemoryManager<MultiRegions::
609  DisContField2D>::AllocateSharedPtr(
611  m_session->GetVariable(i));
612  }
613  }
614 
615  break;
616  }
617  case 3:
618  {
620  {
621  ASSERTL0(false,
622  "3D fully periodic problems not implemented yet");
623  }
624  else
625  {
626  for (i = 0; i < m_fields.num_elements(); i++)
627  {
628  m_fields[i] = MemoryManager<MultiRegions::
629  DisContField3D>::AllocateSharedPtr(
631  m_session->GetVariable(i));
632  }
633  }
634  break;
635  }
636  default:
637  ASSERTL0(false, "Expansion dimension not recognised");
638  break;
639  }
640 
641  // Setting up the normals
644 
645  for (i = 0; i < m_spacedim; ++i)
646  {
647  m_traceNormals[i] =
649  }
650 
651  m_fields[0]->GetTrace()->GetNormals(m_traceNormals);
652  }
653 
654  // Set Default Parameter
655  m_session->LoadParameter("Time", m_time, 0.0);
656  m_session->LoadParameter("TimeStep", m_timestep, 0.0);
657  m_session->LoadParameter("NumSteps", m_steps, 0);
658  m_session->LoadParameter("IO_CheckSteps", m_checksteps, 0);
659  m_session->LoadParameter("IO_CheckTime", m_checktime, 0.0);
660  m_session->LoadParameter("FinTime", m_fintime, 0);
661  m_session->LoadParameter("NumQuadPointsError",
663 
664  // Check uniqueness of checkpoint output
665  ASSERTL0((m_checktime == 0.0 && m_checksteps == 0) ||
666  (m_checktime > 0.0 && m_checksteps == 0) ||
667  (m_checktime == 0.0 && m_checksteps > 0),
668  "Only one of IO_CheckTime and IO_CheckSteps "
669  "should be set!");
670 
671  m_nchk = 0;
672 
673  // Zero all physical fields initially
674  ZeroPhysFields();
675  }
676 
677  /**
678  * @brief Destructor for class EquationSystem.
679  */
681  {
684  LibUtilities::NekManager<LocalRegions::MatrixKey,
686  }
687 
689  std::string name,
690  const MultiRegions::ExpListSharedPtr &field,
691  bool cache)
692  {
693  MultiRegions::ExpListSharedPtr vField = field;
694  if (!field)
695  {
696  vField = m_fields[0];
697  }
698 
699  if (cache)
700  {
701  if ((m_sessionFunctions.find(name) == m_sessionFunctions.end())
702  || (m_sessionFunctions[name]->GetSession() != m_session)
703  || (m_sessionFunctions[name]->GetExpansion() != vField)
704  )
705  {
708  m_session, vField, name, cache);
709  }
710 
711  return m_sessionFunctions[name];
712  }
713  else
714  {
716  new SessionFunction(m_session,vField, name, cache));
717  }
718  }
719 
720  /**
721  * If boundary conditions are time-dependent, they will be evaluated at
722  * the time specified.
723  * @param time The time at which to evaluate the BCs
724  */
726  {
727  std::string varName;
728  int nvariables = m_fields.num_elements();
729  for (int i = 0; i < nvariables; ++i)
730  {
731  varName = m_session->GetVariable(i);
732  m_fields[i]->EvaluateBoundaryConditions(time, varName);
733  }
734  }
735 
736  /**
737  * Compute the error in the L2-norm.
738  * @param field The field to compare.
739  * @param exactsoln The exact solution to compare with.
740  * @param Normalised Normalise L2-error.
741  * @returns Error in the L2-norm.
742  */
744  unsigned int field,
745  const Array<OneD, NekDouble> &exactsoln,
746  bool Normalised)
747  {
748  NekDouble L2error = -1.0;
749 
750  if (m_NumQuadPointsError == 0)
751  {
752  if (m_fields[field]->GetPhysState() == false)
753  {
754  m_fields[field]->BwdTrans(m_fields[field]->GetCoeffs(),
755  m_fields[field]->UpdatePhys());
756  }
757 
758  if (exactsoln.num_elements())
759  {
760  L2error = m_fields[field]->L2(
761  m_fields[field]->GetPhys(), exactsoln);
762  }
763  else if (m_session->DefinesFunction("ExactSolution"))
764  {
766  exactsoln(m_fields[field]->GetNpoints());
767 
768  GetFunction("ExactSolution")->Evaluate(
769  m_session->GetVariable(field), exactsoln, m_time);
770 
771  L2error = m_fields[field]->L2(
772  m_fields[field]->GetPhys(), exactsoln);
773  }
774  else
775  {
776  L2error = m_fields[field]->L2(m_fields[field]->GetPhys());
777  }
778 
779  if (Normalised == true)
780  {
782  1.0);
783 
784  NekDouble Vol = m_fields[field]->PhysIntegral(one);
785  m_comm->AllReduce(Vol, LibUtilities::ReduceSum);
786 
787  L2error = sqrt(L2error*L2error/Vol);
788  }
789  }
790  else
791  {
792  Array<OneD,NekDouble> L2INF(2);
793  L2INF = ErrorExtraPoints(field);
794  L2error = L2INF[0];
795  }
796  return L2error;
797  }
798 
799  /**
800  * Compute the error in the L_inf-norm
801  * @param field The field to compare.
802  * @param exactsoln The exact solution to compare with.
803  * @returns Error in the L_inft-norm.
804  */
806  unsigned int field,
807  const Array<OneD, NekDouble> &exactsoln)
808  {
809  NekDouble Linferror = -1.0;
810 
811  if (m_NumQuadPointsError == 0)
812  {
813  if (m_fields[field]->GetPhysState() == false)
814  {
815  m_fields[field]->BwdTrans(m_fields[field]->GetCoeffs(),
816  m_fields[field]->UpdatePhys());
817  }
818 
819  if (exactsoln.num_elements())
820  {
821  Linferror = m_fields[field]->Linf(
822  m_fields[field]->GetPhys(), exactsoln);
823  }
824  else if (m_session->DefinesFunction("ExactSolution"))
825  {
827  exactsoln(m_fields[field]->GetNpoints());
828 
829  GetFunction("ExactSolution")->Evaluate(
830  m_session->GetVariable(field), exactsoln, m_time);
831 
832  Linferror = m_fields[field]->Linf(
833  m_fields[field]->GetPhys(), exactsoln);
834  }
835  else
836  {
837  Linferror = m_fields[field]->Linf(
838  m_fields[field]->GetPhys());
839  }
840  }
841  else
842  {
843  Array<OneD,NekDouble> L2INF(2);
844  L2INF = ErrorExtraPoints(field);
845  Linferror = L2INF[1];
846  }
847 
848  return Linferror;
849  }
850 
851  /**
852  * Compute the error in the L2-norm, L-inf for a larger number of
853  * quadrature points.
854  * @param field The field to compare.
855  * @returns Error in the L2-norm and L-inf norm.
856  */
858  unsigned int field)
859  {
860  int NumModes = GetNumExpModes();
861  Array<OneD,NekDouble> L2INF(2);
862 
863  const LibUtilities::PointsKey PkeyT1(
865  const LibUtilities::PointsKey PkeyT2(
867  const LibUtilities::PointsKey PkeyQ1(
869  const LibUtilities::PointsKey PkeyQ2(
871  const LibUtilities::BasisKey BkeyT1(
872  LibUtilities::eModified_A,NumModes, PkeyT1);
873  const LibUtilities::BasisKey BkeyT2(
874  LibUtilities::eModified_B, NumModes, PkeyT2);
875  const LibUtilities::BasisKey BkeyQ1(
876  LibUtilities::eModified_A, NumModes, PkeyQ1);
877  const LibUtilities::BasisKey BkeyQ2(
878  LibUtilities::eModified_A, NumModes, PkeyQ2);
879 
882  m_session, BkeyT1, BkeyT2, BkeyQ1, BkeyQ2, m_graph);
883 
884  int ErrorCoordim = ErrorExp->GetCoordim(0);
885  int ErrorNq = ErrorExp->GetTotPoints();
886 
887  Array<OneD,NekDouble> ErrorXc0(ErrorNq, 0.0);
888  Array<OneD,NekDouble> ErrorXc1(ErrorNq, 0.0);
889  Array<OneD,NekDouble> ErrorXc2(ErrorNq, 0.0);
890 
891  switch(ErrorCoordim)
892  {
893  case 1:
894  ErrorExp->GetCoords(ErrorXc0);
895  break;
896  case 2:
897  ErrorExp->GetCoords(ErrorXc0, ErrorXc1);
898  break;
899  case 3:
900  ErrorExp->GetCoords(ErrorXc0, ErrorXc1, ErrorXc2);
901  break;
902  }
904  m_session->GetFunction("ExactSolution", field);
905 
906  // Evaluate the exact solution
907  Array<OneD,NekDouble> ErrorSol(ErrorNq);
908 
909  exSol->Evaluate(ErrorXc0,ErrorXc1,ErrorXc2,m_time,ErrorSol);
910 
911  // Calcualte spectral/hp approximation on the quadrature points
912  // of this new expansion basis
913  ErrorExp->BwdTrans_IterPerExp(m_fields[field]->GetCoeffs(),
914  ErrorExp->UpdatePhys());
915 
916  L2INF[0] = ErrorExp->L2 (ErrorExp->GetPhys(), ErrorSol);
917  L2INF[1] = ErrorExp->Linf(ErrorExp->GetPhys(), ErrorSol);
918 
919  return L2INF;
920  }
921 
922 
923  /**
924  * Set the physical fields based on a restart file, or a function
925  * describing the initial condition given in the session.
926  * @param initialtime Time at which to evaluate the function.
927  * @param dumpInitialConditions Write the initial condition to file?
928  */
930  bool dumpInitialConditions,
931  const int domain)
932  {
933  boost::ignore_unused(initialtime);
934 
935  if (m_session->GetComm()->GetRank() == 0)
936  {
937  cout << "Initial Conditions:" << endl;
938  }
939 
940  if (m_session->DefinesFunction("InitialConditions"))
941  {
942  GetFunction("InitialConditions")->Evaluate(
943  m_session->GetVariables(), m_fields, m_time, domain);
944 
945  if (m_session->GetComm()->GetRank() == 0)
946  {
947 
948  for (int i = 0; i < m_fields.num_elements(); ++i)
949  {
950  std::string varName = m_session->GetVariable(i);
951  cout << " - Field " << varName << ": "
952  << GetFunction("InitialConditions")->Describe(varName, domain)
953  << endl;
954  }
955  }
956  }
957  else
958  {
959  int nq = m_fields[0]->GetNpoints();
960  for (int i = 0; i < m_fields.num_elements(); i++)
961  {
962  Vmath::Zero(nq, m_fields[i]->UpdatePhys(), 1);
963  m_fields[i]->SetPhysState(true);
965  m_fields[i]->UpdateCoeffs(), 1);
966  if (m_session->GetComm()->GetRank() == 0)
967  {
968  cout << " - Field " << m_session->GetVariable(i)
969  << ": 0 (default)" << endl;
970  }
971  }
972 
973  }
974 
975  if (dumpInitialConditions && m_checksteps)
976  {
978  m_nchk++;
979  }
980  }
981 
983  unsigned int field,
984  Array<OneD, NekDouble> &outfield,
985  const NekDouble time)
986  {
987  ASSERTL0 (outfield.num_elements() == m_fields[field]->GetNpoints(),
988  "ExactSolution array size mismatch.");
989  Vmath::Zero(outfield.num_elements(), outfield, 1);
990  if (m_session->DefinesFunction("ExactSolution"))
991  {
992  GetFunction("ExactSolution")->Evaluate(
993  m_session->GetVariable(field), outfield, time);
994  }
995  }
996 
997 
998  /**
999  * By default, nothing needs initialising at the EquationSystem level.
1000  */
1002  {
1003 
1004  }
1005 
1006  /**
1007  *
1008  */
1010  {
1011 
1012  }
1013 
1014  /**
1015  * Virtual function to define if operator in DoSolve is
1016  * negated with regard to the strong form. This is currently
1017  * only used in Arnoldi solves. Default is false.
1018  */
1020  {
1021  return false;
1022  }
1023 
1024  /**
1025  *
1026  */
1028  {
1029 
1030  }
1031 
1032  /**
1033  *
1034  */
1036  {
1037 
1038  }
1039 
1040 
1041  /// Virtual function for generating summary information.
1043  {
1044  SessionSummary(l);
1045  }
1046 
1047 
1048  /**
1049  * Write the field data to file. The file is named according to the
1050  * session name with the extension .fld appended.
1051  */
1053  {
1054  WriteFld(m_sessionName + ".fld");
1055  }
1056 
1057  /**
1058  * Zero the physical fields.
1059  */
1061  {
1062  for (int i = 0; i < m_fields.num_elements(); i++)
1063  {
1065  m_fields[i]->UpdatePhys(),1);
1066  }
1067  }
1068 
1069  /**
1070  * FwdTrans the m_fields members
1071  */
1073  {
1074  for (int i = 0; i < m_fields.num_elements(); i++)
1075  {
1076  m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
1077  m_fields[i]->UpdateCoeffs());
1078  m_fields[i]->SetPhysState(false);
1079  }
1080  }
1081 
1082  /**
1083  * Write the n-th checkpoint file.
1084  * @param n The index of the checkpoint file.
1085  */
1087  {
1088  std::string outname = m_sessionName + "_" +
1089  boost::lexical_cast<std::string>(n);
1090  WriteFld(outname + ".chk");
1091  }
1092 
1093  /**
1094  * Write the n-th checkpoint file.
1095  * @param n The index of the checkpoint file.
1096  */
1098  const int n,
1100  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
1101  std::vector<std::string> &variables)
1102  {
1103  std::string outname = m_sessionName + "_" +
1104  boost::lexical_cast<std::string>(n);
1105  WriteFld(outname, field, fieldcoeffs, variables);
1106  }
1107 
1108  /**
1109  * Write the n-th base flow into a .chk file
1110  * @param n The index of the base flow file.
1111  */
1113  {
1114  std::string outname = m_sessionName + "_BaseFlow_" +
1115  boost::lexical_cast<std::string>(n);
1116 
1117  WriteFld(outname + ".chk");
1118  }
1119 
1120  /**
1121  * Writes the field data to a file with the given filename.
1122  * @param outname Filename to write to.
1123  */
1124  void EquationSystem::WriteFld(const std::string &outname)
1125  {
1126  std::vector<Array<OneD, NekDouble> > fieldcoeffs(
1127  m_fields.num_elements());
1128  std::vector<std::string> variables(m_fields.num_elements());
1129 
1130  for (int i = 0; i < m_fields.num_elements(); ++i)
1131  {
1132  if (m_fields[i]->GetNcoeffs() == m_fields[0]->GetNcoeffs())
1133  {
1134  fieldcoeffs[i] = m_fields[i]->UpdateCoeffs();
1135  }
1136  else
1137  {
1138  fieldcoeffs[i] = Array<OneD,NekDouble>(m_fields[0]->
1139  GetNcoeffs());
1140  m_fields[0]->ExtractCoeffsToCoeffs(m_fields[i],
1141  m_fields[i]->GetCoeffs(),
1142  fieldcoeffs[i]);
1143  }
1144  variables[i] = m_boundaryConditions->GetVariable(i);
1145  }
1146 
1147  ExtraFldOutput(fieldcoeffs, variables);
1148 
1149  WriteFld(outname, m_fields[0], fieldcoeffs, variables);
1150  }
1151 
1152 
1153 
1154  /**
1155  * Writes the field data to a file with the given filename.
1156  * @param outname Filename to write to.
1157  * @param field ExpList on which data is based.
1158  * @param fieldcoeffs An array of array of expansion coefficients.
1159  * @param variables An array of variable names.
1160  */
1162  const std::string &outname,
1164  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
1165  std::vector<std::string> &variables)
1166  {
1167  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef
1168  = field->GetFieldDefinitions();
1169  std::vector<std::vector<NekDouble> > FieldData(FieldDef.size());
1170 
1171  // Copy Data into FieldData and set variable
1172  for(int j = 0; j < fieldcoeffs.size(); ++j)
1173  {
1174  for(int i = 0; i < FieldDef.size(); ++i)
1175  {
1176  // Could do a search here to find correct variable
1177  FieldDef[i]->m_fields.push_back(variables[j]);
1178  field->AppendFieldData(FieldDef[i], FieldData[i],
1179  fieldcoeffs[j]);
1180  }
1181  }
1182 
1183  // Update time in field info if required
1184  if(m_fieldMetaDataMap.find("Time") != m_fieldMetaDataMap.end())
1185  {
1186  m_fieldMetaDataMap["Time"] =
1187  boost::lexical_cast<std::string>(m_time);
1188  }
1189 
1190  // Update step in field info if required
1191  if(m_fieldMetaDataMap.find("ChkFileNum") != m_fieldMetaDataMap.end())
1192  {
1193  m_fieldMetaDataMap["ChkFileNum"] =
1194  boost::lexical_cast<std::string>(m_nchk);
1195  }
1196 
1197  // If necessary, add mapping information to metadata
1198  // and output mapping coordinates
1200  fields[0] = field;
1204  mapping->Output( fieldMetaDataMap, outname);
1205 
1206 #ifdef NEKTAR_DISABLE_BACKUPS
1207  bool backup = false;
1208 #else
1209  bool backup = true;
1210 #endif
1211 
1212  m_fld->Write(outname, FieldDef, FieldData, fieldMetaDataMap, backup);
1213  }
1214 
1215 
1216  /**
1217  * Import field from infile and load into \a m_fields. This routine will
1218  * also perform a \a BwdTrans to ensure data is in both the physical and
1219  * coefficient storage.
1220  * @param infile Filename to read.
1221  */
1223  const std::string &infile,
1225  {
1226  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
1227  std::vector<std::vector<NekDouble> > FieldData;
1228  LibUtilities::FieldIOSharedPtr field_fld =
1230  field_fld->Import(infile,FieldDef,FieldData);
1231 
1232  // Copy FieldData into m_fields
1233  for(int j = 0; j < pFields.num_elements(); ++j)
1234  {
1235  Vmath::Zero(pFields[j]->GetNcoeffs(),
1236  pFields[j]->UpdateCoeffs(),1);
1237 
1238  for(int i = 0; i < FieldDef.size(); ++i)
1239  {
1240  ASSERTL1(FieldDef[i]->m_fields[j] ==
1241  m_session->GetVariable(j),
1242  std::string("Order of ") + infile
1243  + std::string(" data and that defined in "
1244  "m_boundaryconditions differs"));
1245 
1246  pFields[j]->ExtractDataToCoeffs(FieldDef[i], FieldData[i],
1247  FieldDef[i]->m_fields[j],
1248  pFields[j]->UpdateCoeffs());
1249  }
1250  pFields[j]->BwdTrans(pFields[j]->GetCoeffs(),
1251  pFields[j]->UpdatePhys());
1252  }
1253  }
1254 
1255 
1256 
1257  /**
1258  * Import field from infile and load into \a m_fields. This routine will
1259  * also perform a \a BwdTrans to ensure data is in both the physical and
1260  * coefficient storage.
1261  * @param infile Filename to read.
1262  * If optionan \a ndomains is specified it assumes we loop over nodmains
1263  * for each nvariables.
1264  */
1266  const std::string &infile,
1268  const int ndomains)
1269  {
1270  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
1271  std::vector<std::vector<NekDouble> > FieldData;
1272 
1273  LibUtilities::Import(infile,FieldDef,FieldData);
1274 
1275  int nvariables = GetNvariables();
1276 
1277  ASSERTL0(ndomains*nvariables == pFields.num_elements(),
1278  "Number of fields does not match the number of variables and domains");
1279 
1280  // Copy FieldData into m_fields
1281  for(int j = 0; j < ndomains; ++j)
1282  {
1283  for(int i = 0; i < nvariables; ++i)
1284  {
1285  Vmath::Zero(pFields[j*nvariables+i]->GetNcoeffs(),
1286  pFields[j*nvariables+i]->UpdateCoeffs(),1);
1287 
1288  for(int n = 0; n < FieldDef.size(); ++n)
1289  {
1290  ASSERTL1(FieldDef[n]->m_fields[i] == m_session->GetVariable(i),
1291  std::string("Order of ") + infile
1292  + std::string(" data and that defined in "
1293  "m_boundaryconditions differs"));
1294 
1295  pFields[j*nvariables+i]->ExtractDataToCoeffs(
1296  FieldDef[n], FieldData[n],
1297  FieldDef[n]->m_fields[i],
1298  pFields[j*nvariables+i]->UpdateCoeffs());
1299  }
1300  pFields[j*nvariables+i]->BwdTrans(
1301  pFields[j*nvariables+i]->GetCoeffs(),
1302  pFields[j*nvariables+i]->UpdatePhys());
1303  }
1304  }
1305  }
1306 
1307  /**
1308  * Import field from infile and load into \a pField. This routine will
1309  * also perform a \a BwdTrans to ensure data is in both the physical and
1310  * coefficient storage.
1311  */
1313  const std::string &infile,
1315  std::string &pFieldName)
1316  {
1317  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
1318  std::vector<std::vector<NekDouble> > FieldData;
1319 
1320  LibUtilities::FieldIOSharedPtr field_fld =
1322  field_fld->Import(infile,FieldDef,FieldData);
1323  int idx = -1;
1324 
1325  Vmath::Zero(pField->GetNcoeffs(),pField->UpdateCoeffs(),1);
1326 
1327  for(int i = 0; i < FieldDef.size(); ++i)
1328  {
1329  // find the index of the required field in the file.
1330  for(int j = 0; j < FieldData.size(); ++j)
1331  {
1332  if (FieldDef[i]->m_fields[j] == pFieldName)
1333  {
1334  idx = j;
1335  }
1336  }
1337  ASSERTL1(idx >= 0, "Field " + pFieldName + " not found.");
1338 
1339  pField->ExtractDataToCoeffs(FieldDef[i], FieldData[i],
1340  FieldDef[i]->m_fields[idx],
1341  pField->UpdateCoeffs());
1342  }
1343  pField->BwdTrans(pField->GetCoeffs(), pField->UpdatePhys());
1344  }
1345 
1346  /**
1347  * Import field from infile and load into the array \a coeffs.
1348  *
1349  * @param infile Filename to read.
1350  * @param fieldStr an array of string identifying fields to be imported
1351  * @param coeffs array of array of coefficients to store imported data
1352  */
1354  const std::string &infile,
1355  std::vector< std::string> &fieldStr,
1356  Array<OneD, Array<OneD, NekDouble> > &coeffs)
1357  {
1358 
1359  ASSERTL0(fieldStr.size() <= coeffs.num_elements(),
1360  "length of fieldstr should be the same as pFields");
1361 
1362  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
1363  std::vector<std::vector<NekDouble> > FieldData;
1364 
1365  LibUtilities::FieldIOSharedPtr field_fld =
1367  field_fld->Import(infile,FieldDef,FieldData);
1368 
1369  // Copy FieldData into m_fields
1370  for(int j = 0; j < fieldStr.size(); ++j)
1371  {
1372  Vmath::Zero(coeffs[j].num_elements(),coeffs[j],1);
1373  for(int i = 0; i < FieldDef.size(); ++i)
1374  {
1375  m_fields[0]->ExtractDataToCoeffs(FieldDef[i], FieldData[i],
1376  fieldStr[j], coeffs[j]);
1377  }
1378  }
1379  }
1380 
1381  /**
1382  * Write out a summary of the session data.
1383  * @param out Output stream to write data to.
1384  */
1386  {
1387  AddSummaryItem(s, "EquationType",
1388  m_session->GetSolverInfo("EQTYPE"));
1389  AddSummaryItem(s, "Session Name", m_sessionName);
1390  AddSummaryItem(s, "Spatial Dim.", m_spacedim);
1391  AddSummaryItem(s, "Max SEM Exp. Order",
1392  m_fields[0]->EvalBasisNumModesMax());
1393 
1394  if (m_session->GetComm()->GetSize() > 1)
1395  {
1396  AddSummaryItem(s, "Num. Processes",
1397  m_session->GetComm()->GetSize());
1398  }
1399 
1401  {
1402  AddSummaryItem(s, "Quasi-3D", "Homogeneous in z-direction");
1403  AddSummaryItem(s, "Expansion Dim.", m_expdim + 1);
1404  AddSummaryItem(s, "Num. Hom. Modes (z)", m_npointsZ);
1405  AddSummaryItem(s, "Hom. length (LZ)", m_LhomZ);
1406  AddSummaryItem(s, "FFT Type", m_useFFT ? "FFTW" : "MVM");
1407  if (m_halfMode)
1408  {
1409  AddSummaryItem(s, "ModeType", "Half Mode");
1410  }
1411  else if (m_singleMode)
1412  {
1413  AddSummaryItem(s, "ModeType", "Single Mode");
1414  }
1415  else if (m_multipleModes)
1416  {
1417  AddSummaryItem(s, "ModeType", "Multiple Modes");
1418  }
1419  }
1420  else if(m_HomogeneousType == eHomogeneous2D)
1421  {
1422  AddSummaryItem(s, "Quasi-3D", "Homogeneous in yz-plane");
1423  AddSummaryItem(s, "Expansion Dim.", m_expdim + 2);
1424  AddSummaryItem(s, "Num. Hom. Modes (y)", m_npointsY);
1425  AddSummaryItem(s, "Num. Hom. Modes (z)", m_npointsZ);
1426  AddSummaryItem(s, "Hom. length (LY)", m_LhomY);
1427  AddSummaryItem(s, "Hom. length (LZ)", m_LhomZ);
1428  AddSummaryItem(s, "FFT Type", m_useFFT ? "FFTW" : "MVM");
1429  }
1430  else
1431  {
1432  AddSummaryItem(s, "Expansion Dim.", m_expdim);
1433  }
1434 
1435  if (m_session->DefinesSolverInfo("UpwindType"))
1436  {
1437  AddSummaryItem(s, "Riemann Solver",
1438  m_session->GetSolverInfo("UpwindType"));
1439  }
1440 
1441  if (m_session->DefinesSolverInfo("AdvectionType"))
1442  {
1443  std::string AdvectionType;
1444  AdvectionType = m_session->GetSolverInfo("AdvectionType");
1445  AddSummaryItem(s, "Advection Type", GetAdvectionFactory().
1446  GetClassDescription(AdvectionType));
1447  }
1448 
1450  {
1451  AddSummaryItem(s, "Projection Type", "Continuous Galerkin");
1452  }
1454  {
1455  AddSummaryItem(s, "Projection Type", "Discontinuous Galerkin");
1456  }
1458  {
1459  AddSummaryItem(s, "Projection Type",
1460  "Mixed Continuous Galerkin and Discontinuous");
1461  }
1462 
1463  if (m_session->DefinesSolverInfo("DiffusionType"))
1464  {
1465  std::string DiffusionType;
1466  DiffusionType = m_session->GetSolverInfo("DiffusionType");
1467  AddSummaryItem(s, "Diffusion Type", GetDiffusionFactory().
1468  GetClassDescription(DiffusionType));
1469  }
1470  }
1471 
1473  {
1474  return Array<OneD, bool>(m_session->GetVariables().size(), false);
1475  }
1476 
1478  {
1479  ASSERTL0(false, "This function is not valid for the Base class");
1481  return null;
1482  }
1483 
1485  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
1486  std::vector<std::string> &variables)
1487  {
1488  boost::ignore_unused(fieldcoeffs, variables);
1489  }
1490 
1491  }
1492 }
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &l)
Virtual function for generating summary information.
bool m_singleMode
Flag to determine if single homogeneous mode is used.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:163
NekDouble m_time
Current time of simulation.
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
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:293
NekDouble m_timestep
Time step size.
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:46
SOLVER_UTILS_EXPORT int GetNvariables()
bool m_halfMode
Flag to determine if half homogeneous mode is used.
NekDouble m_LhomZ
physical length in Z direction (if homogeneous)
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.
Array< OneD, bool > m_checkIfSystemSingular
Flag to indicate if the fields should be checked for singularity.
int m_expdim
Expansion dimension.
std::shared_ptr< ContField2D > ContField2DSharedPtr
Definition: ContField2D.h:289
SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n)
Write checkpoint file of m_fields.
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
Principle Modified Functions .
Definition: BasisType.h:48
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].
STL namespace.
NekDouble m_checktime
Time between checkpoints.
bool m_useFFT
Flag to determine if FFT is used for homogeneous transform.
NekDouble m_LhomY
physical length in Y direction (if homogeneous)
bool m_specHP_dealiasing
Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
SOLVER_UTILS_EXPORT int GetNumExpModes()
Fourier Expansion .
Definition: BasisType.h:53
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput(std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
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...
int m_npointsZ
number of points in Z direction (if homogeneous)
std::shared_ptr< ExpList2D > ExpList2DSharedPtr
Shared pointer to an ExpList2D object.
Definition: ExpList2D.h:48
virtual SOLVER_UTILS_EXPORT ~EquationSystem()
Destructor.
std::string m_sessionName
Name of the session.
int m_nchk
Number of checkpoints written so far.
std::map< std::string, std::string > FieldMetaDataMap
Definition: FieldIO.h:52
int m_checksteps
Number of steps between checkpoints.
LibUtilities::CommSharedPtr m_comm
Communicator.
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:58
NekDouble m_fintime
Finish time of the simulation.
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Initialisation object for EquationSystem.
int m_steps
Number of steps to take.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
int m_HomoDirec
number of homogenous directions
std::shared_ptr< ContField3D > ContField3DSharedPtr
Definition: ContField3D.h:208
SOLVER_UTILS_EXPORT void FwdTransFields()
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:65
bool m_multipleModes
Flag to determine if use multiple homogenenous modes are used.
Fourier Modified expansions with just the real part of the first mode .
Definition: BasisType.h:60
virtual SOLVER_UTILS_EXPORT Array< OneD, bool > v_GetSystemSingularChecks()
int m_npointsY
number of points in Y direction (if homogeneous)
Abstraction of a global continuous one-dimensional spectral/hp element expansion which approximates t...
Definition: ContField1D.h:55
std::map< std::string, SolverUtils::SessionFunctionSharedPtr > m_sessionFunctions
Map of known SessionFunctions.
NekMatrix< NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag >, BlockMatrixTag > DNekScalBlkMat
Definition: NekTypeDefs.hpp:65
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff()
Virtual function for transformation to coefficient space.
Principle Modified Functions .
Definition: BasisType.h:49
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
virtual SOLVER_UTILS_EXPORT void v_Output(void)
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
Pointer to boundary conditions object.
int m_spacedim
Spatial dimension (>= expansion dim).
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
Defines a specification for a set of points.
Definition: Points.h:59
double NekDouble
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
SOLVER_UTILS_EXPORT LibUtilities::SessionReaderSharedPtr GetSession()
Get Session name.
virtual SOLVER_UTILS_EXPORT void v_DoInitialise()
Virtual function for initialisation implementation.
bool m_homogen_dealiasing
Flag to determine if dealiasing is used for homogeneous simulations.
EquationSystemFactory & GetEquationSystemFactory()
Fourier Modified expansions with just the imaginary part of the first mode .
Definition: BasisType.h:61
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.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:131
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
static std::shared_ptr< FieldIO > CreateDefault(const LibUtilities::SessionReaderSharedPtr session)
Returns an object for the default FieldIO method.
Definition: FieldIO.cpp:195
SOLVER_UTILS_EXPORT void WriteFld(const std::string &outname)
Write field data to the given filename.
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
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
Fourier ModifiedExpansion with just the first mode .
Definition: BasisType.h:59
SOLVER_UTILS_EXPORT int GetTraceNpoints()
SOLVER_UTILS_EXPORT int GetNcoeffs()
1D Non Evenly-spaced points for Single Mode analysis
Definition: PointsType.h:66
LibUtilities::FieldIOSharedPtr m_fld
Field input/output.
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp()
Virtual function to identify if operator is negated in DoSolve.
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure(void)
SOLVER_UTILS_EXPORT void ZeroPhysFields()
SOLVER_UTILS_EXPORT void SessionSummary(SummaryList &vSummary)
Write out a session summary.
SpatialDomains::MeshGraphSharedPtr m_graph
Pointer to graph defining mesh.
SOLVER_UTILS_EXPORT void ExtraFldOutput(std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
Used to lookup the create function in NekManager.
Definition: MatrixKey.h:67
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow(const int n)
Write base flow file of m_fields.
virtual SOLVER_UTILS_EXPORT void v_DoSolve()
Virtual function for solve implementation.
std::shared_ptr< SessionFunction > SessionFunctionSharedPtr
GLOBAL_MAPPING_EXPORT typedef std::shared_ptr< Mapping > MappingSharedPtr
A shared pointer to a Mapping object.
Definition: Mapping.h:50
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376
int m_NumQuadPointsError
Number of Quadrature points used to work out the error.
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:306
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:53
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
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:268
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution(unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys()
Virtual function for transformation to physical space.
Describes the specification for a Basis.
Definition: Basis.h:49
std::shared_ptr< SessionReader > SessionReaderSharedPtr
enum HomogeneousType m_HomogeneousType
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
Provides a generic Factory class.
Definition: NekFactory.hpp:103
SOLVER_UTILS_EXPORT void ImportFld(const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
Input field data from the given file.