Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EquationSystem.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File EquationSystem.h
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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Base class for individual solvers.
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #ifndef NEKTAR_SOLVERUTILS_EQUATIONSYSTEM_H
37 #define NEKTAR_SOLVERUTILS_EQUATIONSYSTEM_H
38 
48 #include <MultiRegions/ExpList.h>
51 #include <SolverUtils/Core/Misc.h>
52 
53 namespace Nektar
54 {
55  namespace SolverUtils
56  {
58 
59  /// A shared pointer to an EquationSystem object
60  typedef boost::shared_ptr<EquationSystem> EquationSystemSharedPtr;
61  /// Datatype of the NekFactory used to instantiate classes derived from
62  /// the EquationSystem class.
64  std::string, EquationSystem,
68 
69  struct loadedFldField {
70  std::vector<LibUtilities::FieldDefinitionsSharedPtr> fieldDef;
71  std::vector<std::vector<NekDouble> > fieldData;
72  } ;
73 
74  /// A base class for describing how to solve specific equations.
75  class EquationSystem : public boost::enable_shared_from_this<EquationSystem>
76  {
77  public:
78  /// Destructor
80 
81  // Set up trace normals if required
83 
84  /// Initialises the members of this object.
85  SOLVER_UTILS_EXPORT inline void InitObject();
86 
87  /// Perform any initialisation necessary before solving the problem.
88  SOLVER_UTILS_EXPORT inline void DoInitialise();
89 
90  /// Solve the problem.
91  SOLVER_UTILS_EXPORT inline void DoSolve();
92 
93  /// Transform from coefficient to physical space.
95 
96  /// Transform from physical to coefficient space.
98 
99  /// Perform output operations after solve.
100  SOLVER_UTILS_EXPORT inline void Output();
101 
102  /// Linf error computation
103  SOLVER_UTILS_EXPORT inline NekDouble LinfError(unsigned int field, const Array<OneD,NekDouble> &exactsoln = NullNekDouble1DArray);
104 
105  /// Get Session name
107  {
108  return m_sessionName;
109  }
110 
111  template<class T>
112  boost::shared_ptr<T> as()
113  {
114 #if defined __INTEL_COMPILER && BOOST_VERSION > 105200
115  typedef typename boost::shared_ptr<T>::element_type E;
116  E * p = dynamic_cast< E* >( shared_from_this().get() );
117  ASSERTL1(p, "Cannot perform cast");
118  return boost::shared_ptr<T>( shared_from_this(), p );
119 #else
120  return boost::dynamic_pointer_cast<T>( shared_from_this() );
121 #endif
122  }
123 
124  /// Reset Session name
125  SOLVER_UTILS_EXPORT void ResetSessionName(std::string newname)
126  {
127  m_sessionName = newname;
128  }
129 
130  /// Get Session name
132  {
133  return m_session;
134  }
135 
136  /// Get pressure field if available
138 
139  /// Print a summary of parameters and solver characteristics.
140  SOLVER_UTILS_EXPORT inline void PrintSummary(std::ostream &out);
141 
142  /// Set parameter m_lambda
143  SOLVER_UTILS_EXPORT inline void SetLambda(NekDouble lambda);
144 
145  /// Evaluates a function as specified in the session file.
147  Array<OneD, Array<OneD, NekDouble> >& pArray,
148  std::string pFunctionName,
149  const NekDouble pTime = 0.0,
150  const int domain = 0);
151 
152  /// Populate given fields with the function from session.
154  std::vector<std::string> pFieldNames,
155  Array<OneD, Array<OneD, NekDouble> > &pFields,
156  const std::string& pName,
157  const NekDouble& pTime = 0.0,
158  const int domain = 0);
159 
160  /// Populate given fields with the function from session.
162  std::vector<std::string> pFieldNames,
164  const std::string& pName,
165  const NekDouble& pTime = 0.0,
166  const int domain = 0);
167 
168  // Populate an array with a function variable from session.
170  std::string pFieldName,
171  Array<OneD, NekDouble>& pArray,
172  const std::string& pFunctionName,
173  const NekDouble& pTime = 0.0,
174  const int domain = 0);
175 
176  // Describe a function.
178  std::string pFieldName,
179  const std::string &pFunctionName,
180  const int domain);
181 
182  /// Perform initialisation of the base flow.
185 
186  /// Initialise the data in the dependent fields.
188  NekDouble initialtime = 0.0,
189  bool dumpInitialConditions = true,
190  const int domain = 0);
191 
192  /// Evaluates an exact solution
194  int field,
195  Array<OneD, NekDouble> &outfield,
196  const NekDouble time);
197 
198  /// Compute the L2 error between fields and a given exact
199  /// solution.
201  unsigned int field,
202  const Array<OneD,NekDouble> &exactsoln,
203  bool Normalised = false);
204 
205  /// Compute the L2 error of the fields
207  unsigned int field,
208  bool Normalised = false)
209  {
210  return L2Error(field,NullNekDouble1DArray,Normalised);
211  }
212 
213  /// Compute error (L2 and L_inf) over an larger set of quadrature
214  /// points return [L2 Linf]
216  unsigned int field);
217 
218  /// Compute the inner product \f$ (\nabla \phi \cdot F) \f$.
220  const Array<OneD, Array<OneD, NekDouble> > &F,
221  Array<OneD, NekDouble> &outarray);
222 
223  /// Compute the inner product \f$ (\phi, \nabla \cdot F) \f$.
225  const Array<OneD, Array<OneD, NekDouble> > &F,
226  Array<OneD, NekDouble> &outarray);
227 
228  /// Compute the inner product \f$ (\phi, V\cdot \nabla u) \f$.
230  const Array<OneD, Array<OneD, NekDouble> > &V,
232  Array<OneD, NekDouble> &outarray,
233  bool UseContCoeffs = false);
234 
235  /// Compute the non-conservative advection \f$ (V \cdot \nabla u)
236  /// \f$.
238  const Array<OneD, Array<OneD, NekDouble> > &V,
240  Array<OneD, NekDouble> &outarray,
242 
243  /// Calculate the weak discontinuous Galerkin advection.
245  const Array<OneD, Array<OneD, NekDouble> >& InField,
246  Array<OneD, Array<OneD, NekDouble> >& OutField,
247  bool NumericalFluxIncludesNormal = true,
248  bool InFieldIsInPhysSpace = false,
249  int nvariables = 0);
250 
251  /// Calculate weak DG Diffusion in the LDG form.
253  const Array<OneD, Array<OneD, NekDouble> >& InField,
254  Array<OneD, Array<OneD, NekDouble> >& OutField,
255  bool NumericalFluxIncludesNormal = true,
256  bool InFieldIsInPhysSpace = false);
257 
258  /// Write checkpoint file of #m_fields.
259  SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n);
260 
261  /// Write checkpoint file of custom data fields.
263  const int n,
265  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
266  std::vector<std::string> &variables);
267 
268  /// Write base flow file of #m_fields.
269  SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow(const int n);
270 
271  /// Write field data to the given filename.
272  SOLVER_UTILS_EXPORT void WriteFld(const std::string &outname);
273 
274  /// Write input fields to the given filename.
276  const std::string &outname,
278  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
279  std::vector<std::string> &variables);
280 
281  /// Input field data from the given file.
283  const std::string &infile,
285 
286  /// Input field data from the given file to multiple domains
288  const std::string &infile,
290  const int ndomains);
291 
292  /// Output a field.
293  /// Input field data into array from the given file.
295  const std::string &infile,
296  std::vector<std::string> &fieldStr,
297  Array<OneD, Array<OneD, NekDouble> > &coeffs);
298 
299  /// Output a field.
300  /// Input field data into ExpList from the given file.
302  const std::string &infile,
304  std::string &pFieldName);
305 
306  /// Builds map of which element holds each history point.
308 
309  /// Probe each history point and write to file.
310  SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out);
311 
312  /// Write out a session summary.
314 
315  SOLVER_UTILS_EXPORT inline Array<
317 
318 
319  /// Get hold of FieldInfoMap so it can be updated
322 
323  /// Return final time
325 
326  SOLVER_UTILS_EXPORT inline int GetNcoeffs();
327 
328  SOLVER_UTILS_EXPORT inline int GetNcoeffs(const int eid);
329 
330  SOLVER_UTILS_EXPORT inline int GetNumExpModes();
331 
334 
335  SOLVER_UTILS_EXPORT inline int GetNvariables();
336 
337  SOLVER_UTILS_EXPORT inline const std::string
338  GetVariable(unsigned int i);
339 
341 
343 
344  SOLVER_UTILS_EXPORT inline int GetExpSize();
345 
346  SOLVER_UTILS_EXPORT inline int GetPhys_Offset(int n);
347 
348  SOLVER_UTILS_EXPORT inline int GetCoeff_Offset(int n);
349 
350  SOLVER_UTILS_EXPORT inline int GetTotPoints();
351 
352  SOLVER_UTILS_EXPORT inline int GetTotPoints(int n);
353 
354  SOLVER_UTILS_EXPORT inline int GetNpoints();
355 
357 
358  SOLVER_UTILS_EXPORT inline int GetSteps();
359 
361 
362  SOLVER_UTILS_EXPORT inline void CopyFromPhysField(const int i,
363  Array<OneD, NekDouble> &output);
364 
365  SOLVER_UTILS_EXPORT inline void CopyToPhysField(const int i,
366  Array<OneD, NekDouble> &output);
367 
368  SOLVER_UTILS_EXPORT inline void SetSteps(const int steps);
369 
371 
373 
375  const int i,
376  Array<OneD, Array<OneD, NekDouble> >&physfield,
377  Array<OneD, Array<OneD, NekDouble> >&flux);
378 
380  const int i,
381  Array<OneD, Array<OneD, NekDouble> >&physfield,
382  Array<OneD, Array<OneD, NekDouble> >&fluxX,
383  Array<OneD, Array<OneD, NekDouble> > &fluxY);
384 
386  const int i,
387  const int j,
388  Array<OneD, Array<OneD, NekDouble> > &physfield,
389  Array<OneD, Array<OneD, NekDouble> > &flux);
390 
392  Array<OneD, Array<OneD, NekDouble> > &physfield,
393  Array<OneD, Array<OneD, NekDouble> > &numflux);
394 
396  Array<OneD, Array<OneD, NekDouble> > &physfield,
397  Array<OneD, Array<OneD, NekDouble> > &numfluxX,
398  Array<OneD, Array<OneD, NekDouble> > &numfluxY);
399 
401  const Array<OneD, Array<OneD, NekDouble> > &ufield,
402  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &uflux);
403 
405  const Array<OneD, Array<OneD, NekDouble> > &ufield,
406  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &qfield,
407  Array<OneD, Array<OneD, NekDouble> > &qflux);
408 
410  const bool modbasis);
411 
412  /// Perform a case-insensitive string comparison.
414  const std::string & s1, const std::string& s2) ;
415 
417  {
418  return m_nchk;
419  }
420 
422  {
423  m_nchk = num;
424  }
425 
427  {
428  return m_checksteps;
429  }
430 
432  {
433  m_checksteps = num;
434  }
435 
437  const NekDouble time)
438  {
439  m_time = time;
440  }
441 
443  const int step)
444  {
445  m_initialStep = step;
446  }
447 
448  /// Evaluates the boundary conditions at the given time.
450 
451  /// Virtual function to identify if operator is negated in DoSolve
452  SOLVER_UTILS_EXPORT virtual bool v_NegatedOp();
453 
454  protected:
455  /// Communicator
457  /// The session reader
459  /// Field input/output
461  /// Map of interpolator objects
462  std::map<std::string, FieldUtils::Interpolator > m_interpolators;
463  /// pts fields we already read from disk: {funcFilename: (filename, ptsfield)}
464  std::map<std::string, std::pair<std::string, LibUtilities::PtsFieldSharedPtr> > m_loadedPtsFields;
465  // fld fiels already loaded from disk: {funcFilename: (filename, loadedFldField)}
466  std::map<std::string, std::pair<std::string, loadedFldField> > m_loadedFldFields;
467  /// Array holding all dependent variables.
469  /// Base fields.
471  /// Array holding all dependent variables.
473  /// Pointer to boundary conditions object.
475  /// Pointer to graph defining mesh.
477  /// Name of the session.
478  std::string m_sessionName;
479  /// Current time of simulation.
481  /// Number of the step where the simulation should begin
483  /// Finish time of the simulation.
485  /// Time step size
487  /// Lambda constant in real system if one required.
489  /// Time between checkpoints.
491  /// Number of checkpoints written so far
492  int m_nchk;
493  /// Number of steps to take.
494  int m_steps;
495  /// Number of steps between checkpoints.
497  /// Spatial dimension (>= expansion dim).
499  /// Expansion dimension.
500  int m_expdim;
501  /// Flag to determine if single homogeneous mode is used.
503  /// Flag to determine if half homogeneous mode is used.
505  /// Flag to determine if use multiple homogenenous modes are used.
507  /// Flag to determine if FFT is used for homogeneous transform.
508  bool m_useFFT;
509  /**
510  * \brief Flag to determine if dealiasing is used for
511  * homogeneous simulations.
512  */
514  /**
515  * \brief Flag to determine if dealisising is usde for the
516  * Spectral/hp element discretisation.
517  */
519  /// Type of projection; e.g continuous or discontinuous.
521  /// Array holding trace normals for DG simulations in the forwards direction.
523  /// 1 x nvariable x nq
525  /// 2 x m_spacedim x nq
527  /// Flag to indicate if the fields should be checked for singularity.
529  /// Map to identify relevant solver info to dump in output fields
531 
532  /// Number of Quadrature points used to work out the error
534 
535  /// Parameter for homogeneous expansions
537  {
542  };
543 
544 
545 
547 
548  NekDouble m_LhomX; ///< physical length in X direction (if homogeneous)
549  NekDouble m_LhomY; ///< physical length in Y direction (if homogeneous)
550  NekDouble m_LhomZ; ///< physical length in Z direction (if homogeneous)
551 
552  int m_npointsX; ///< number of points in X direction (if homogeneous)
553  int m_npointsY; ///< number of points in Y direction (if homogeneous)
554  int m_npointsZ; ///< number of points in Z direction (if homogeneous)
555 
556  int m_HomoDirec; ///< number of homogenous directions
557 
558 
559  /// Initialises EquationSystem class members.
561 
562  // Here for consistency purposes with old version
563  int nocase_cmp(const std::string & s1, const std::string& s2)
564  {
565  return NoCaseStringCompare(s1,s2);
566  }
567 
568  SOLVER_UTILS_EXPORT virtual void v_InitObject();
569 
570  /// Virtual function for initialisation implementation.
571  SOLVER_UTILS_EXPORT virtual void v_DoInitialise();
572 
573  /// Virtual function for solve implementation.
574  SOLVER_UTILS_EXPORT virtual void v_DoSolve();
575 
576 
577  /// Virtual function for the L_inf error computation between fields and a given exact solution.
579  unsigned int field,
580  const Array<OneD, NekDouble> &exactsoln = NullNekDouble1DArray);
581 
582  /// Virtual function for the L_2 error computation between fields and a given exact solution.
584  unsigned int field,
585  const Array<OneD, NekDouble> &exactsoln = NullNekDouble1DArray,
586  bool Normalised = false);
587 
588  /// Virtual function for transformation to physical space.
590 
591  /// Virtual function for transformation to coefficient space.
593 
594  /// Virtual function for generating summary information.
596 
598  NekDouble initialtime = 0.0,
599  bool dumpInitialConditions = true,
600  const int domain = 0);
601 
603  unsigned int field,
604  Array<OneD, NekDouble> &outfield,
605  const NekDouble time);
606 
607  // Populate an array with a function variable from session.
609  std::string pFieldName,
610  Array<OneD, NekDouble>& pArray,
611  const std::string& pFunctionName,
612  const NekDouble& pTime = 0.0,
613  const int domain = 0);
614 
615  // Populate an array with a function variable from session.
617  std::string pFieldName,
618  Array<OneD, NekDouble>& pArray,
619  const std::string& pFunctionName,
620  const NekDouble& pTime = 0.0,
621  const int domain = 0);
622 
624  std::string pFieldName,
625  Array<OneD, NekDouble>& pArray,
626  const std::string& pFunctionName,
627  const NekDouble& pTime = 0.0,
628  const int domain = 0);
629 
631  std::string funcFilename,
632  std::string filename,
634 
635  //Initialise m_base in order to store the base flow from a file
637 
638  // Fill m_base with the values stored in a fld file
640  std::string pInfile,
642 
643  // Ouptut field information
644  SOLVER_UTILS_EXPORT virtual void v_Output(void);
645 
646  // Get pressure field if available
648 
650  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
651  std::vector<std::string> &variables);
652 
653  private:
654 
657  const int i, Array<OneD,
658  Array<OneD, NekDouble> >&physfield,
659  Array<OneD, Array<OneD, NekDouble> >&flux);
660 
662  const int i, const int j,
663  Array<OneD, Array<OneD, NekDouble> >&physfield,
664  Array<OneD, Array<OneD, NekDouble> >&flux);
665 
667  const int i, Array<OneD,
668  Array<OneD, NekDouble> >&physfield,
669  Array<OneD, Array<OneD, NekDouble> >&fluxX,
670  Array<OneD, Array<OneD, NekDouble> > &fluxY);
671 
673  Array<OneD, Array<OneD, NekDouble> > &physfield,
674  Array<OneD, Array<OneD, NekDouble> > &numflux);
675 
677  Array<OneD, Array<OneD, NekDouble> > &physfield,
678  Array<OneD, Array<OneD, NekDouble> > &numfluxX,
679  Array<OneD, Array<OneD, NekDouble> > &numfluxY);
680 
682  const Array<OneD, Array<OneD, NekDouble> > &ufield,
683  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &uflux);
684 
686  const Array<OneD, Array<OneD, NekDouble> > &ufield,
687  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &qfield,
688  Array<OneD, Array<OneD, NekDouble > > &qflux);
689 
690  SOLVER_UTILS_EXPORT void PrintProgressbar(const int position,
691  const int goal) const
692  {
693  LibUtilities::PrintProgressbar(position, goal, "Interpolating");
694  }
695  };
696 
697 
698  /**
699  * This is the second part of the two-phase initialisation process.
700  * Calls to virtual functions will correctly resolve to the derived class
701  * during this phase of the construction.
702  */
704  {
705  v_InitObject();
706  }
707 
708 
709  /**
710  * This allows initialisation of the solver which cannot be completed
711  * during object construction (such as setting of initial conditions).
712  *
713  * Public interface routine to virtual function implementation.
714  */
716  {
717  v_DoInitialise();
718  }
719 
720 
721  /**
722  * Performs the transformation from coefficient to physical space.
723  *
724  * Public interface routine to virtual function implementation.
725  */
727  {
729  }
730 
731  /**
732  * Performs the transformation from physical to coefficient space.
733  *
734  * Public interface routine to virtual function implementation.
735  */
737  {
739  }
740 
741 
742  /**
743  * Performs the actual solve.
744  *
745  * Public interface routine to virtual function implementation.
746  */
747  inline void EquationSystem::DoSolve(void)
748  {
749  v_DoSolve();
750  }
751 
752 
753  /**
754  * Perform output operations after solve.
755  */
756  inline void EquationSystem::Output(void)
757  {
758  v_Output();
759  }
760 
761  /**
762  * L_inf Error computation
763  * Public interface routine to virtual function implementation.
764  */
765  inline NekDouble EquationSystem::LinfError(unsigned int field, const Array<OneD,NekDouble> &exactsoln)
766  {
767  return v_LinfError(field, exactsoln);
768  }
769 
770  /**
771  * L_2 Error computation
772  * Public interface routine to virtual function implementation.
773  */
774  inline NekDouble EquationSystem::L2Error(unsigned int field, const Array<OneD,NekDouble> &exactsoln, bool Normalised)
775  {
776  return v_L2Error(field, exactsoln, Normalised);
777  }
778 
779  /**
780  * Get Pressure field if available
781  */
783  {
784  return v_GetPressure();
785  }
786 
787  /**
788  * Prints a summary of variables and problem parameters.
789  *
790  * Public interface routine to virtual function implementation.
791  *
792  * @param out The ostream object to write to.
793  */
794  inline void EquationSystem::PrintSummary(std::ostream &out)
795  {
796  if (m_session->GetComm()->GetRank() == 0)
797  {
798  std::vector<std::pair<std::string, std::string> > vSummary;
799  v_GenerateSummary(vSummary);
800 
801  out << "=======================================================================" << std::endl;
802  SummaryList::const_iterator x;
803  for (x = vSummary.begin(); x != vSummary.end(); ++x)
804  {
805  out << "\t";
806  out.width(20);
807  out << x->first << ": " << x->second << std::endl;
808  }
809  out << "=======================================================================" << std::endl;
810  }
811  }
812 
814  {
815  m_lambda = lambda;
816  }
817 
819  bool dumpInitialConditions,
820  const int domain)
821  {
822  v_SetInitialConditions(initialtime,dumpInitialConditions,domain);
823  }
824 
825  /// Evaluates an exact solution
827  Array<OneD, NekDouble> &outfield,
828  const NekDouble time)
829  {
830  v_EvaluateExactSolution(field, outfield, time);
831  }
832 
834  {
835  return m_fields;
836  }
837 
838  /// Return final time
840  {
841  return m_time;
842  }
843 
844  inline int EquationSystem::GetNcoeffs(void)
845  {
846  return m_fields[0]->GetNcoeffs();
847  }
848 
849  inline int EquationSystem::GetNcoeffs(const int eid)
850  {
851  return m_fields[0]->GetNcoeffs(eid);
852  }
853 
855  {
856  return m_graph->GetExpansions().begin()->second->m_basisKeyVector[0]
857  .GetNumModes();
858  }
859 
861  {
862  return m_fields[0]->EvalBasisNumModesMaxPerExp();
863  }
864 
866  {
867  return m_session->GetVariables().size();
868  }
869 
870  inline const std::string EquationSystem::GetVariable(unsigned int i)
871  {
872  return m_session->GetVariable(i);
873  }
874 
876  {
877  return GetTraceNpoints();
878  }
879 
881  {
882  return m_fields[0]->GetTrace()->GetNpoints();
883  }
884 
885  inline int EquationSystem::GetExpSize(void)
886  {
887  return m_fields[0]->GetExpSize();
888  }
889 
891  {
892  return m_fields[0]->GetPhys_Offset(n);
893  }
894 
896  {
897  return m_fields[0]->GetCoeff_Offset(n);
898  }
899 
901  {
902  return m_fields[0]->GetNpoints();
903  }
904 
906  {
907  return m_fields[0]->GetTotPoints(n);
908  }
909 
910  inline int EquationSystem::GetNpoints(void)
911  {
912  return m_fields[0]->GetNpoints();
913  }
914 
916  {
917  return (m_fields.num_elements() - 1);
918  }
919 
920  inline int EquationSystem::GetSteps(void)
921  {
922  return m_steps;
923  }
924 
926  {
927  return m_timestep;
928  }
929 
930  inline void EquationSystem::SetSteps(const int steps)
931  {
932  m_steps = steps;
933  }
934 
935  inline void EquationSystem::CopyFromPhysField(const int i,
936  Array<OneD, NekDouble> &output)
937  {
938  Vmath::Vcopy(output.num_elements(), m_fields[i]->GetPhys(), 1, output, 1 );
939  }
940 
941  inline void EquationSystem::CopyToPhysField(const int i,
942  Array<OneD, NekDouble> &output)
943  {
944  Vmath::Vcopy(output.num_elements(), output, 1, m_fields[i]->UpdatePhys(), 1 );
945  }
946 
947  inline void EquationSystem::GetFluxVector(const int i,
948  Array<OneD, Array<OneD, NekDouble> >&physfield,
950  {
951  v_GetFluxVector(i,physfield, flux);
952  }
953 
954  inline void EquationSystem::GetFluxVector(const int i,
955  Array<OneD, Array<OneD, NekDouble> >&physfield,
958  {
959  v_GetFluxVector(i,physfield, fluxX, fluxY);
960  }
961 
962  inline void EquationSystem::GetFluxVector(const int i, const int j,
963  Array<OneD, Array<OneD, NekDouble> > &physfield,
965  {
966  v_GetFluxVector(i,j,physfield,flux);
967  }
968 
970  Array<OneD, Array<OneD, NekDouble> > &numflux)
971  {
972  v_NumericalFlux(physfield, numflux);
973  }
974 
976  Array<OneD, Array<OneD, NekDouble> > &numfluxX,
977  Array<OneD, Array<OneD, NekDouble> > &numfluxY)
978  {
979  v_NumericalFlux(physfield, numfluxX, numfluxY);
980  }
981 
983  const Array<OneD, Array<OneD, NekDouble> > &ufield,
985  {
986  v_NumFluxforScalar(ufield, uflux);
987  }
988 
990  const Array<OneD, Array<OneD, NekDouble> > &ufield,
993  {
994  v_NumFluxforVector(ufield, qfield, qflux);
995  }
996  }
997 }
998 
999 #endif
SOLVER_UTILS_EXPORT void EvaluateFunctionPts(std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
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.
SOLVER_UTILS_EXPORT void NumFluxforScalar(const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
f SOLVER_UTILS_EXPORT void AdvectionNonConservativeForm(const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray)
Compute the non-conservative advection.
SOLVER_UTILS_EXPORT void InitialiseBaseFlow(Array< OneD, Array< OneD, NekDouble > > &base)
Perform initialisation of the base flow.
LibUtilities::NekFactory< std::string, EquationSystem, const LibUtilities::SessionReaderSharedPtr & > EquationSystemFactory
Datatype of the NekFactory used to instantiate classes derived from the EquationSystem class...
SOLVER_UTILS_EXPORT std::string DescribeFunction(std::string pFieldName, const std::string &pFunctionName, const int domain)
Provide a description of a function for a given field name.
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_gradtan
1 x nvariable x nq
A base class for describing how to solve specific equations.
static Array< OneD, NekDouble > NullNekDouble1DArray
int PrintProgressbar(const int position, const int goal, const string message, int lastprogress=-1)
Prints a progressbar.
Definition: Progressbar.hpp:69
SOLVER_UTILS_EXPORT NekDouble L2Error(unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
Compute the L2 error between fields and a given exact solution.
SOLVER_UTILS_EXPORT void LoadPts(std::string funcFilename, std::string filename, Nektar::LibUtilities::PtsFieldSharedPtr &outPts)
NekDouble m_time
Current time of simulation.
SOLVER_UTILS_EXPORT void WeakAdvectionNonConservativeForm(const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, bool UseContCoeffs=false)
Compute the inner product .
std::map< std::string, std::pair< std::string, loadedFldField > > m_loadedFldFields
NekDouble m_timestep
Time step size.
SOLVER_UTILS_EXPORT void SetUpTraceNormals(void)
SOLVER_UTILS_EXPORT int NoCaseStringCompare(const std::string &s1, const std::string &s2)
Perform a case-insensitive string comparison.
SOLVER_UTILS_EXPORT void SetSteps(const int steps)
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:47
NekDouble m_lambda
Lambda constant in real system if one required.
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.
SOLVER_UTILS_EXPORT const Array< OneD, int > GetNumExpModesPerExp()
Array< OneD, bool > m_checkIfSystemSingular
Flag to indicate if the fields should be checked for singularity.
int m_expdim
Expansion dimension.
SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n)
Write checkpoint file of m_fields.
std::vector< std::vector< NekDouble > > fieldData
NekDouble m_LhomX
physical length in X direction (if homogeneous)
SOLVER_UTILS_EXPORT Array< OneD, MultiRegions::ExpListSharedPtr > & UpdateFields()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
std::map< std::string, std::pair< std::string, LibUtilities::PtsFieldSharedPtr > > m_loadedPtsFields
pts fields we already read from disk: {funcFilename: (filename, ptsfield)}
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].
ProjectionType
Type of Galerkin projection.
SOLVER_UTILS_EXPORT void WeakAdvectionGreensDivergenceForm(const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
Compute the inner product .
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()
SOLVER_UTILS_EXPORT void TransCoeffToPhys()
Transform from coefficient to physical space.
SOLVER_UTILS_EXPORT void InitObject()
Initialises the members of this object.
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput(std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
SOLVER_UTILS_EXPORT void PrintProgressbar(const int position, const int goal) const
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)
virtual SOLVER_UTILS_EXPORT ~EquationSystem()
Destructor.
std::string m_sessionName
Name of the session.
int m_nchk
Number of checkpoints written so far.
SOLVER_UTILS_EXPORT void WeakAdvectionDivergenceForm(const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
Compute the inner product .
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
std::map< std::string, std::string > FieldMetaDataMap
Definition: FieldIO.h:54
int m_checksteps
Number of steps between checkpoints.
SOLVER_UTILS_EXPORT int GetCoeff_Offset(int n)
SOLVER_UTILS_EXPORT void EvaluateFunctionExp(std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
boost::shared_ptr< PtsField > PtsFieldSharedPtr
Definition: PtsField.h:178
LibUtilities::CommSharedPtr m_comm
Communicator.
SOLVER_UTILS_EXPORT void SetModifiedBasis(const bool modbasis)
SOLVER_UTILS_EXPORT LibUtilities::FieldMetaDataMap & UpdateFieldMetaDataMap()
Get hold of FieldInfoMap so it can be updated.
SOLVER_UTILS_EXPORT void ResetSessionName(std::string newname)
Reset Session name.
SOLVER_UTILS_EXPORT int GetTotPoints()
SOLVER_UTILS_EXPORT const std::string GetVariable(unsigned int i)
SOLVER_UTILS_EXPORT void Output()
Perform output operations after solve.
NekDouble m_fintime
Finish time of the simulation.
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Initialisation object for EquationSystem.
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_tanbasis
2 x m_spacedim x nq
SOLVER_UTILS_EXPORT void CopyToPhysField(const int i, Array< OneD, NekDouble > &output)
boost::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.
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
SOLVER_UTILS_EXPORT void FwdTransFields()
Array< OneD, MultiRegions::ExpListSharedPtr > m_derivedfields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT void WeakDGDiffusion(const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false)
Calculate weak DG Diffusion in the LDG form.
SOLVER_UTILS_EXPORT void EvaluateFunctionFld(std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
boost::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
SOLVER_UTILS_EXPORT int GetCheckpointSteps()
SOLVER_UTILS_EXPORT int GetSteps()
bool m_multipleModes
Flag to determine if use multiple homogenenous modes are used.
virtual SOLVER_UTILS_EXPORT Array< OneD, bool > v_GetSystemSingularChecks()
SOLVER_UTILS_EXPORT void TransPhysToCoeff()
Transform from physical to coefficient space.
SOLVER_UTILS_EXPORT void NumericalFlux(Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
SOLVER_UTILS_EXPORT void SetUpBaseFields(SpatialDomains::MeshGraphSharedPtr &mesh)
SOLVER_UTILS_EXPORT void CopyFromPhysField(const int i, Array< OneD, NekDouble > &output)
SOLVER_UTILS_EXPORT void SetCheckpointSteps(int num)
int nocase_cmp(const std::string &s1, const std::string &s2)
int m_npointsY
number of points in Y direction (if homogeneous)
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff()
Virtual function for transformation to coefficient space.
virtual SOLVER_UTILS_EXPORT void v_GetFluxVector(const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
virtual SOLVER_UTILS_EXPORT void v_Output(void)
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
Pointer to boundary conditions object.
boost::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
int m_spacedim
Spatial dimension (>= expansion dim).
SOLVER_UTILS_EXPORT NekDouble L2Error(unsigned int field, bool Normalised=false)
Compute the L2 error of the fields.
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
boost::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:309
double NekDouble
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
SOLVER_UTILS_EXPORT void DoInitialise()
Perform any initialisation necessary before solving the problem.
SOLVER_UTILS_EXPORT LibUtilities::SessionReaderSharedPtr GetSession()
Get Session name.
SOLVER_UTILS_EXPORT NekDouble GetTimeStep()
SOLVER_UTILS_EXPORT void GetFluxVector(const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
SOLVER_UTILS_EXPORT void EvaluateFunction(Array< OneD, Array< OneD, NekDouble > > &pArray, std::string pFunctionName, const NekDouble pTime=0.0, const int domain=0)
Evaluates a function as specified in the session file.
SOLVER_UTILS_EXPORT void ScanForHistoryPoints()
Builds map of which element holds each history point.
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.
SOLVER_UTILS_EXPORT std::string GetSessionName()
Get Session name.
EquationSystemFactory & GetEquationSystemFactory()
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.
int m_npointsX
number of points in X direction (if homogeneous)
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
SOLVER_UTILS_EXPORT int GetPhys_Offset(int n)
SOLVER_UTILS_EXPORT NekDouble GetFinalTime()
Return final time.
SOLVER_UTILS_EXPORT void PrintSummary(std::ostream &out)
Print a summary of parameters and solver characteristics.
SOLVER_UTILS_EXPORT void WriteFld(const std::string &outname)
Write field data to the given filename.
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetExpSize()
#define SOLVER_UTILS_EXPORT
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
Array< OneD, MultiRegions::ExpListSharedPtr > m_base
Base fields.
SOLVER_UTILS_EXPORT NekDouble LinfError(unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
Linf error computation.
boost::shared_ptr< BoundaryConditions > BoundaryConditionsSharedPtr
Definition: Conditions.h:284
SOLVER_UTILS_EXPORT int GetTraceNpoints()
SOLVER_UTILS_EXPORT int GetNcoeffs()
LibUtilities::FieldIOSharedPtr m_fld
Field input/output.
HomogeneousType
Parameter for homogeneous expansions.
virtual SOLVER_UTILS_EXPORT void v_NumFluxforVector(const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
int m_initialStep
Number of the step where the simulation should begin.
SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr GetPressure()
Get pressure field if available.
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 NumFluxforVector(const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
SOLVER_UTILS_EXPORT void ZeroPhysFields()
SOLVER_UTILS_EXPORT void EvaluateExactSolution(int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
Evaluates an exact solution.
SOLVER_UTILS_EXPORT void SessionSummary(SummaryList &vSummary)
Write out a session summary.
SOLVER_UTILS_EXPORT void SetTime(const NekDouble time)
SpatialDomains::MeshGraphSharedPtr m_graph
Pointer to graph defining mesh.
virtual SOLVER_UTILS_EXPORT void v_NumFluxforScalar(const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
SOLVER_UTILS_EXPORT void WriteHistoryData(std::ostream &out)
Probe each history point and write to file.
SOLVER_UTILS_EXPORT void SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
Initialise the data in the dependent fields.
SOLVER_UTILS_EXPORT void WeakDGAdvection(const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false, int nvariables=0)
Calculate the weak discontinuous Galerkin advection.
SOLVER_UTILS_EXPORT int GetNumElmVelocity()
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow(const int n)
Write base flow file of m_fields.
SOLVER_UTILS_EXPORT void ImportFldBase(std::string pInfile, SpatialDomains::MeshGraphSharedPtr pGraph)
virtual SOLVER_UTILS_EXPORT void v_DoSolve()
Virtual function for solve implementation.
std::map< std::string, FieldUtils::Interpolator > m_interpolators
Map of interpolator objects.
int m_NumQuadPointsError
Number of Quadrature points used to work out the error.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:228
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:442
SOLVER_UTILS_EXPORT void DoSolve()
Solve the problem.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
std::vector< LibUtilities::FieldDefinitionsSharedPtr > fieldDef
SOLVER_UTILS_EXPORT void SetInitialStep(const int step)
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.
SOLVER_UTILS_EXPORT void SetLambda(NekDouble lambda)
Set parameter m_lambda.
SOLVER_UTILS_EXPORT EquationSystem(const LibUtilities::SessionReaderSharedPtr &pSession)
Initialises EquationSystem class members.
virtual SOLVER_UTILS_EXPORT void v_NumericalFlux(Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
enum HomogeneousType m_HomogeneousType
SOLVER_UTILS_EXPORT void SetCheckpointNumber(int num)
SOLVER_UTILS_EXPORT int GetCheckpointNumber()
Provides a generic Factory class.
Definition: NekFactory.hpp:116
SOLVER_UTILS_EXPORT void ImportFld(const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
Input field data from the given file.