Nektar++
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>
50 #include <SolverUtils/Core/Misc.h>
51 
52 namespace Nektar
53 {
54  namespace SolverUtils
55  {
57 
58  /// A shared pointer to an EquationSystem object
59  typedef boost::shared_ptr<EquationSystem> EquationSystemSharedPtr;
60  /// Datatype of the NekFactory used to instantiate classes derived from
61  /// the EquationSystem class.
63  std::string, EquationSystem,
67 
68  /// A base class for describing how to solve specific equations.
69  class EquationSystem : public boost::enable_shared_from_this<EquationSystem>
70  {
71  public:
72  /// Destructor
74 
75  // Set up trace normals if required
77 
78  /// Initialises the members of this object.
79  SOLVER_UTILS_EXPORT inline void InitObject();
80 
81  /// Perform any initialisation necessary before solving the problem.
82  SOLVER_UTILS_EXPORT inline void DoInitialise();
83 
84  /// Solve the problem.
85  SOLVER_UTILS_EXPORT inline void DoSolve();
86 
87  /// Transform from coefficient to physical space.
89 
90  /// Transform from physical to coefficient space.
92 
93  /// Perform output operations after solve.
94  SOLVER_UTILS_EXPORT inline void Output();
95 
96  /// Linf error computation
97  SOLVER_UTILS_EXPORT inline NekDouble LinfError(unsigned int field, const Array<OneD,NekDouble> &exactsoln = NullNekDouble1DArray);
98 
99  /// Get Session name
101  {
102  return m_sessionName;
103  }
104 
105  template<class T>
106  boost::shared_ptr<T> as()
107  {
108 #if defined __INTEL_COMPILER && BOOST_VERSION > 105200
109  typedef typename boost::shared_ptr<T>::element_type E;
110  E * p = dynamic_cast< E* >( shared_from_this().get() );
111  ASSERTL1(p, "Cannot perform cast");
112  return boost::shared_ptr<T>( shared_from_this(), p );
113 #else
114  return boost::dynamic_pointer_cast<T>( shared_from_this() );
115 #endif
116  }
117 
118  /// Reset Session name
119  SOLVER_UTILS_EXPORT void ResetSessionName(std::string newname)
120  {
121  m_sessionName = newname;
122  }
123 
124  /// Get Session name
126  {
127  return m_session;
128  }
129 
130  /// Get pressure field if available
132 
133  /// Print a summary of parameters and solver characteristics.
134  SOLVER_UTILS_EXPORT inline void PrintSummary(std::ostream &out);
135 
136  /// Set parameter m_lambda
137  SOLVER_UTILS_EXPORT inline void SetLambda(NekDouble lambda);
138 
139  /// Evaluates a function as specified in the session file.
141  Array<OneD, Array<OneD, NekDouble> >& pArray,
142  std::string pFunctionName,
143  const NekDouble pTime = 0.0,
144  const int domain = 0);
145 
146  /// Populate given fields with the function from session.
148  std::vector<std::string> pFieldNames,
149  Array<OneD, Array<OneD, NekDouble> > &pFields,
150  const std::string& pName,
151  const int domain = 0);
152 
153  /// Populate given fields with the function from session.
155  std::vector<std::string> pFieldNames,
157  const std::string& pName,
158  const int domain = 0);
159 
160  // Populate an array with a function variable from session.
162  std::string pFieldName,
163  Array<OneD, NekDouble>& pArray,
164  const std::string& pFunctionName,
165  const NekDouble& pTime = 0.0,
166  const int domain = 0);
167 
168  // Describe a function.
170  std::string pFieldName,
171  const std::string &pFunctionName,
172  const int domain);
173 
174  /// Perform initialisation of the base flow.
177 
178  /// Initialise the data in the dependent fields.
180  NekDouble initialtime = 0.0,
181  bool dumpInitialConditions = true,
182  const int domain = 0);
183 
184  /// Evaluates an exact solution
186  int field,
187  Array<OneD, NekDouble> &outfield,
188  const NekDouble time);
189 
190  /// Compute the L2 error between fields and a given exact
191  /// solution.
193  unsigned int field,
194  const Array<OneD,NekDouble> &exactsoln,
195  bool Normalised = false);
196 
197  /// Compute the L2 error of the fields
199  unsigned int field,
200  bool Normalised = false)
201  {
202  return L2Error(field,NullNekDouble1DArray,Normalised);
203  }
204 
205  /// Compute error (L2 and L_inf) over an larger set of quadrature
206  /// points return [L2 Linf]
208  unsigned int field);
209 
210  /// Compute the inner product \f$ (\nabla \phi \cdot F) \f$.
212  const Array<OneD, Array<OneD, NekDouble> > &F,
213  Array<OneD, NekDouble> &outarray);
214 
215  /// Compute the inner product \f$ (\phi, \nabla \cdot F) \f$.
217  const Array<OneD, Array<OneD, NekDouble> > &F,
218  Array<OneD, NekDouble> &outarray);
219 
220  /// Compute the inner product \f$ (\phi, V\cdot \nabla u) \f$.
222  const Array<OneD, Array<OneD, NekDouble> > &V,
224  Array<OneD, NekDouble> &outarray,
225  bool UseContCoeffs = false);
226 
227  /// Compute the non-conservative advection \f$ (V \cdot \nabla u)
228  /// \f$.
230  const Array<OneD, Array<OneD, NekDouble> > &V,
232  Array<OneD, NekDouble> &outarray,
234 
235  /// Calculate the weak discontinuous Galerkin advection.
237  const Array<OneD, Array<OneD, NekDouble> >& InField,
238  Array<OneD, Array<OneD, NekDouble> >& OutField,
239  bool NumericalFluxIncludesNormal = true,
240  bool InFieldIsInPhysSpace = false,
241  int nvariables = 0);
242 
243  /// Calculate weak DG Diffusion in the LDG form.
245  const Array<OneD, Array<OneD, NekDouble> >& InField,
246  Array<OneD, Array<OneD, NekDouble> >& OutField,
247  bool NumericalFluxIncludesNormal = true,
248  bool InFieldIsInPhysSpace = false);
249 
250  /// Write checkpoint file of #m_fields.
251  SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n);
252 
253  /// Write checkpoint file of custom data fields.
255  const int n,
257  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
258  std::vector<std::string> &variables);
259 
260  /// Write base flow file of #m_fields.
261  SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow(const int n);
262 
263  /// Write field data to the given filename.
264  SOLVER_UTILS_EXPORT void WriteFld(const std::string &outname);
265 
266  /// Write input fields to the given filename.
268  const std::string &outname,
270  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
271  std::vector<std::string> &variables);
272 
273  /// Input field data from the given file.
275  const std::string &infile,
277 
278  /// Input field data from the given file to multiple domains
280  const std::string &infile,
282  const int ndomains);
283 
284  /// Output a field.
285  /// Input field data into array from the given file.
287  const std::string &infile,
288  std::vector<std::string> &fieldStr,
289  Array<OneD, Array<OneD, NekDouble> > &coeffs);
290 
291  /// Output a field.
292  /// Input field data into ExpList from the given file.
294  const std::string &infile,
296  std::string &pFieldName);
297 
298  /// Builds map of which element holds each history point.
300 
301  /// Probe each history point and write to file.
302  SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out);
303 
304  /// Write out a session summary.
306 
307  SOLVER_UTILS_EXPORT inline Array<
309 
310 
311  /// Get hold of FieldInfoMap so it can be updated
314 
315  /// Return final time
317 
318  SOLVER_UTILS_EXPORT inline int GetNcoeffs();
319 
320  SOLVER_UTILS_EXPORT inline int GetNcoeffs(const int eid);
321 
322  SOLVER_UTILS_EXPORT inline int GetNumExpModes();
323 
326 
327  SOLVER_UTILS_EXPORT inline int GetNvariables();
328 
329  SOLVER_UTILS_EXPORT inline const std::string
330  GetVariable(unsigned int i);
331 
333 
335 
336  SOLVER_UTILS_EXPORT inline int GetExpSize();
337 
338  SOLVER_UTILS_EXPORT inline int GetPhys_Offset(int n);
339 
340  SOLVER_UTILS_EXPORT inline int GetCoeff_Offset(int n);
341 
342  SOLVER_UTILS_EXPORT inline int GetTotPoints();
343 
344  SOLVER_UTILS_EXPORT inline int GetTotPoints(int n);
345 
346  SOLVER_UTILS_EXPORT inline int GetNpoints();
347 
349 
350  SOLVER_UTILS_EXPORT inline int GetSteps();
351 
353 
354  SOLVER_UTILS_EXPORT inline void CopyFromPhysField(const int i,
355  Array<OneD, NekDouble> &output);
356 
357  SOLVER_UTILS_EXPORT inline void CopyToPhysField(const int i,
358  Array<OneD, NekDouble> &output);
359 
360  SOLVER_UTILS_EXPORT inline void SetStepsToOne();
361 
363 
365 
367  const int i,
368  Array<OneD, Array<OneD, NekDouble> >&physfield,
369  Array<OneD, Array<OneD, NekDouble> >&flux);
370 
372  const int i,
373  Array<OneD, Array<OneD, NekDouble> >&physfield,
374  Array<OneD, Array<OneD, NekDouble> >&fluxX,
375  Array<OneD, Array<OneD, NekDouble> > &fluxY);
376 
378  const int i,
379  const int j,
380  Array<OneD, Array<OneD, NekDouble> > &physfield,
381  Array<OneD, Array<OneD, NekDouble> > &flux);
382 
384  Array<OneD, Array<OneD, NekDouble> > &physfield,
385  Array<OneD, Array<OneD, NekDouble> > &numflux);
386 
388  Array<OneD, Array<OneD, NekDouble> > &physfield,
389  Array<OneD, Array<OneD, NekDouble> > &numfluxX,
390  Array<OneD, Array<OneD, NekDouble> > &numfluxY);
391 
393  const Array<OneD, Array<OneD, NekDouble> > &ufield,
394  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &uflux);
395 
397  const Array<OneD, Array<OneD, NekDouble> > &ufield,
398  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &qfield,
399  Array<OneD, Array<OneD, NekDouble> > &qflux);
400 
402  const bool modbasis);
403 
404  /// Perform a case-insensitive string comparison.
406  const string & s1, const string& s2) ;
407 
408  protected:
409  /// Communicator
411  /// The session reader
413  /// Field input/output
415  /// Map of the interpolation weights for a specific filename.
416  map<std::string, Array<OneD, Array<OneD, float> > > m_interpWeights;
417  /// Map of the interpolation indices for a specific filename.
418  map<std::string, Array<OneD, Array<OneD, unsigned int> > > m_interpInds;
419  /// Array holding all dependent variables.
421  /// Base fields.
423  /// Array holding all dependent variables.
425  /// Pointer to boundary conditions object.
427  /// Pointer to graph defining mesh.
429  /// Name of the session.
430  std::string m_sessionName;
431  /// Current time of simulation.
433  /// Finish time of the simulation.
435  /// Time step size
437  /// Lambda constant in real system if one required.
439  /// Time between checkpoints.
441  /// Number of steps to take.
442  int m_steps;
443  /// Number of steps between checkpoints.
445  /// Spatial dimension (>= expansion dim).
447  /// Expansion dimension.
448  int m_expdim;
449  /// Flag to determine if single homogeneous mode is used.
451  /// Flag to determine if half homogeneous mode is used.
453  /// Flag to determine if use multiple homogenenous modes are used.
455  /// Flag to determine if FFT is used for homogeneous transform.
456  bool m_useFFT;
457  /**
458  * \brief Flag to determine if dealiasing is used for
459  * homogeneous simulations.
460  */
462  /**
463  * \brief Flag to determine if dealisising is usde for the
464  * Spectral/hp element discretisation.
465  */
467  /// Type of projection; e.g continuous or discontinuous.
469  /// Array holding trace normals for DG simulations in the forwards direction.
471  /// 1 x nvariable x nq
473  /// 2 x m_spacedim x nq
475  /// Flag to indicate if the fields should be checked for singularity.
477  /// Map to identify relevant solver info to dump in output fields
479 
480  /// Number of Quadrature points used to work out the error
482 
483  /// Parameter for homogeneous expansions
485  {
490  };
491 
492 
493 
495 
496  NekDouble m_LhomX; ///< physical length in X direction (if homogeneous)
497  NekDouble m_LhomY; ///< physical length in Y direction (if homogeneous)
498  NekDouble m_LhomZ; ///< physical length in Z direction (if homogeneous)
499 
500  int m_npointsX; ///< number of points in X direction (if homogeneous)
501  int m_npointsY; ///< number of points in Y direction (if homogeneous)
502  int m_npointsZ; ///< number of points in Z direction (if homogeneous)
503 
504  int m_HomoDirec; ///< number of homogenous directions
505 
506  int m_NumMode; ///< Mode to use in case of single mode analysis
507 
508 
509  /// Initialises EquationSystem class members.
511 
512  // Here for consistency purposes with old version
513  int nocase_cmp(const string & s1, const string& s2)
514  {
515  return NoCaseStringCompare(s1,s2);
516  }
517 
518  SOLVER_UTILS_EXPORT virtual void v_InitObject();
519 
520  /// Evaluates the boundary conditions at the given time.
522 
523  /// Virtual function for initialisation implementation.
524  SOLVER_UTILS_EXPORT virtual void v_DoInitialise();
525 
526  /// Virtual function for solve implementation.
527  SOLVER_UTILS_EXPORT virtual void v_DoSolve();
528 
529  /// Virtual function for the L_inf error computation between fields and a given exact solution.
531  unsigned int field,
532  const Array<OneD, NekDouble> &exactsoln = NullNekDouble1DArray);
533 
534  /// Virtual function for the L_2 error computation between fields and a given exact solution.
536  unsigned int field,
537  const Array<OneD, NekDouble> &exactsoln = NullNekDouble1DArray,
538  bool Normalised = false);
539 
540  /// Virtual function for transformation to physical space.
542 
543  /// Virtual function for transformation to coefficient space.
545 
546  /// Virtual function for generating summary information.
548 
550  NekDouble initialtime = 0.0,
551  bool dumpInitialConditions = true,
552  const int domain = 0);
553 
555  unsigned int field,
556  Array<OneD, NekDouble> &outfield,
557  const NekDouble time);
558 
559  //Initialise m_base in order to store the base flow from a file
561 
562  // Fill m_base with the values stored in a fld file
564  std::string pInfile,
566 
567  // Ouptut field information
568  SOLVER_UTILS_EXPORT virtual void v_Output(void);
569 
570  // Get pressure field if available
572 
574  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
575  std::vector<std::string> &variables);
576 
577  private:
578 
581  const int i, Array<OneD,
582  Array<OneD, NekDouble> >&physfield,
583  Array<OneD, Array<OneD, NekDouble> >&flux);
584 
586  const int i, const int j,
587  Array<OneD, Array<OneD, NekDouble> >&physfield,
588  Array<OneD, Array<OneD, NekDouble> >&flux);
589 
591  const int i, Array<OneD,
592  Array<OneD, NekDouble> >&physfield,
593  Array<OneD, Array<OneD, NekDouble> >&fluxX,
594  Array<OneD, Array<OneD, NekDouble> > &fluxY);
595 
597  Array<OneD, Array<OneD, NekDouble> > &physfield,
598  Array<OneD, Array<OneD, NekDouble> > &numflux);
599 
601  Array<OneD, Array<OneD, NekDouble> > &physfield,
602  Array<OneD, Array<OneD, NekDouble> > &numfluxX,
603  Array<OneD, Array<OneD, NekDouble> > &numfluxY);
604 
606  const Array<OneD, Array<OneD, NekDouble> > &ufield,
607  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &uflux);
608 
610  const Array<OneD, Array<OneD, NekDouble> > &ufield,
611  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &qfield,
612  Array<OneD, Array<OneD, NekDouble > > &qflux);
613 
614  SOLVER_UTILS_EXPORT void PrintProgressbar(const int position,
615  const int goal) const
616  {
617  LibUtilities::PrintProgressbar(position, goal, "Interpolating");
618  }
619  };
620 
621 
622  /**
623  * This is the second part of the two-phase initialisation process.
624  * Calls to virtual functions will correctly resolve to the derived class
625  * during this phase of the construction.
626  */
628  {
629  v_InitObject();
630  }
631 
632 
633  /**
634  * This allows initialisation of the solver which cannot be completed
635  * during object construction (such as setting of initial conditions).
636  *
637  * Public interface routine to virtual function implementation.
638  */
640  {
641  v_DoInitialise();
642  }
643 
644 
645  /**
646  * Performs the transformation from coefficient to physical space.
647  *
648  * Public interface routine to virtual function implementation.
649  */
651  {
653  }
654 
655  /**
656  * Performs the transformation from physical to coefficient space.
657  *
658  * Public interface routine to virtual function implementation.
659  */
661  {
663  }
664 
665 
666  /**
667  * Performs the actual solve.
668  *
669  * Public interface routine to virtual function implementation.
670  */
671  inline void EquationSystem::DoSolve(void)
672  {
673  v_DoSolve();
674  }
675 
676 
677  /**
678  * Perform output operations after solve.
679  */
680  inline void EquationSystem::Output(void)
681  {
682  v_Output();
683  }
684 
685  /**
686  * L_inf Error computation
687  * Public interface routine to virtual function implementation.
688  */
689  inline NekDouble EquationSystem::LinfError(unsigned int field, const Array<OneD,NekDouble> &exactsoln)
690  {
691  return v_LinfError(field, exactsoln);
692  }
693 
694  /**
695  * L_2 Error computation
696  * Public interface routine to virtual function implementation.
697  */
698  inline NekDouble EquationSystem::L2Error(unsigned int field, const Array<OneD,NekDouble> &exactsoln, bool Normalised)
699  {
700  return v_L2Error(field, exactsoln, Normalised);
701  }
702 
703  /**
704  * Get Pressure field if available
705  */
707  {
708  return v_GetPressure();
709  }
710 
711  /**
712  * Prints a summary of variables and problem parameters.
713  *
714  * Public interface routine to virtual function implementation.
715  *
716  * @param out The ostream object to write to.
717  */
718  inline void EquationSystem::PrintSummary(std::ostream &out)
719  {
720  if (m_session->GetComm()->GetRank() == 0)
721  {
722  std::vector<std::pair<std::string, std::string> > vSummary;
723  v_GenerateSummary(vSummary);
724 
725  out << "=======================================================================" << endl;
726  SummaryList::const_iterator x;
727  for (x = vSummary.begin(); x != vSummary.end(); ++x)
728  {
729  out << "\t";
730  out.width(20);
731  out << x->first << ": " << x->second << endl;
732  }
733  out << "=======================================================================" << endl;
734  }
735  }
736 
738  {
739  m_lambda = lambda;
740  }
741 
743  bool dumpInitialConditions,
744  const int domain)
745  {
746  v_SetInitialConditions(initialtime,dumpInitialConditions,domain);
747  }
748 
749  /// Evaluates an exact solution
751  Array<OneD, NekDouble> &outfield,
752  const NekDouble time)
753  {
754  v_EvaluateExactSolution(field, outfield, time);
755  }
756 
758  {
759  return m_fields;
760  }
761 
762  /// Return final time
764  {
765  return m_time;
766  }
767 
768  inline int EquationSystem::GetNcoeffs(void)
769  {
770  return m_fields[0]->GetNcoeffs();
771  }
772 
773  inline int EquationSystem::GetNcoeffs(const int eid)
774  {
775  return m_fields[0]->GetNcoeffs(eid);
776  }
777 
779  {
780  return m_graph->GetExpansions().begin()->second->m_basisKeyVector[0]
781  .GetNumModes();
782  }
783 
785  {
786  return m_fields[0]->EvalBasisNumModesMaxPerExp();
787  }
788 
790  {
791  return m_session->GetVariables().size();
792  }
793 
794  inline const std::string EquationSystem::GetVariable(unsigned int i)
795  {
796  return m_session->GetVariable(i);
797  }
798 
800  {
801  return GetTraceNpoints();
802  }
803 
805  {
806  return m_fields[0]->GetTrace()->GetNpoints();
807  }
808 
809  inline int EquationSystem::GetExpSize(void)
810  {
811  return m_fields[0]->GetExpSize();
812  }
813 
815  {
816  return m_fields[0]->GetPhys_Offset(n);
817  }
818 
820  {
821  return m_fields[0]->GetCoeff_Offset(n);
822  }
823 
825  {
826  return m_fields[0]->GetNpoints();
827  }
828 
830  {
831  return m_fields[0]->GetTotPoints(n);
832  }
833 
834  inline int EquationSystem::GetNpoints(void)
835  {
836  return m_fields[0]->GetNpoints();
837  }
838 
840  {
841  return (m_fields.num_elements() - 1);
842  }
843 
844  inline int EquationSystem::GetSteps(void)
845  {
846  return m_steps;
847  }
848 
850  {
851  return m_timestep;
852  }
853 
855  {
856  m_steps=1;
857  }
858 
859  inline void EquationSystem::CopyFromPhysField(const int i,
860  Array<OneD, NekDouble> &output)
861  {
862  Vmath::Vcopy(output.num_elements(), m_fields[i]->GetPhys(), 1, output, 1 );
863  }
864 
865  inline void EquationSystem::CopyToPhysField(const int i,
866  Array<OneD, NekDouble> &output)
867  {
868  Vmath::Vcopy(output.num_elements(), output, 1, m_fields[i]->UpdatePhys(), 1 );
869  }
870 
871  inline void EquationSystem::GetFluxVector(const int i,
872  Array<OneD, Array<OneD, NekDouble> >&physfield,
874  {
875  v_GetFluxVector(i,physfield, flux);
876  }
877 
878  inline void EquationSystem::GetFluxVector(const int i,
879  Array<OneD, Array<OneD, NekDouble> >&physfield,
882  {
883  v_GetFluxVector(i,physfield, fluxX, fluxY);
884  }
885 
886  inline void EquationSystem::GetFluxVector(const int i, const int j,
887  Array<OneD, Array<OneD, NekDouble> > &physfield,
889  {
890  v_GetFluxVector(i,j,physfield,flux);
891  }
892 
894  Array<OneD, Array<OneD, NekDouble> > &numflux)
895  {
896  v_NumericalFlux(physfield, numflux);
897  }
898 
900  Array<OneD, Array<OneD, NekDouble> > &numfluxX,
901  Array<OneD, Array<OneD, NekDouble> > &numfluxY)
902  {
903  v_NumericalFlux(physfield, numfluxX, numfluxY);
904  }
905 
907  const Array<OneD, Array<OneD, NekDouble> > &ufield,
909  {
910  v_NumFluxforScalar(ufield, uflux);
911  }
912 
914  const Array<OneD, Array<OneD, NekDouble> > &ufield,
917  {
918  v_NumFluxforVector(ufield, qfield, qflux);
919  }
920  }
921 }
922 
923 #endif
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &l)
Virtual function for generating summary information.
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.
bool m_SingleMode
Flag to determine if single homogeneous mode is used.
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.
int nocase_cmp(const string &s1, const string &s2)
static Array< OneD, NekDouble > NullNekDouble1DArray
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.
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 .
NekDouble m_timestep
Time step size.
SOLVER_UTILS_EXPORT void SetUpTraceNormals(void)
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 void SetStepsToOne()
SOLVER_UTILS_EXPORT int GetNvariables()
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.
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.
map< std::string, Array< OneD, Array< OneD, unsigned int > > > m_interpInds
Map of the interpolation indices for a specific filename.
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_NumMode
Mode to use in case of single mode analysis.
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:50
std::map< std::string, std::string > FieldMetaDataMap
Definition: FieldIO.h:65
int m_checksteps
Number of steps between checkpoints.
SOLVER_UTILS_EXPORT int GetCoeff_Offset(int n)
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.
SOLVER_UTILS_EXPORT int NoCaseStringCompare(const string &s1, const string &s2)
Perform a case-insensitive string comparison.
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.
boost::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:53
SOLVER_UTILS_EXPORT int GetSteps()
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)
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)
boost::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
Definition: ExpList.h:1340
virtual SOLVER_UTILS_EXPORT void v_Output(void)
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
Pointer to boundary conditions object.
bool m_HalfMode
Flag to determine if half homogeneous mode is used.
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:236
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:269
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)
SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr GetPressure()
Get pressure field if available.
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.
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.
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:165
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:432
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:1016
void PrintProgressbar(const int position, const int goal, const string message)
Prints a progressbar.
Definition: Progressbar.hpp:69
bool m_MultipleModes
Flag to determine if use multiple homogenenous modes are used.
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.
map< std::string, Array< OneD, Array< OneD, float > > > m_interpWeights
Map of the interpolation weights for a specific filename.
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
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.