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>
50 #include <SolverUtils/Core/Misc.h>
51 
52 namespace Nektar
53 {
54 namespace FieldUtils {
55 class Interpolator;
56 }
57  namespace SolverUtils
58  {
60 
61  /// A shared pointer to an EquationSystem object
62  typedef boost::shared_ptr<EquationSystem> EquationSystemSharedPtr;
63  /// Datatype of the NekFactory used to instantiate classes derived from
64  /// the EquationSystem class.
66  std::string, EquationSystem,
70 
71  struct loadedFldField {
72  std::vector<LibUtilities::FieldDefinitionsSharedPtr> fieldDef;
73  std::vector<std::vector<NekDouble> > fieldData;
74  } ;
75 
76  /// A base class for describing how to solve specific equations.
77  class EquationSystem : public boost::enable_shared_from_this<EquationSystem>
78  {
79  public:
80  /// Destructor
82 
83  // Set up trace normals if required
85 
86  /// Initialises the members of this object.
87  SOLVER_UTILS_EXPORT inline void InitObject();
88 
89  /// Perform any initialisation necessary before solving the problem.
90  SOLVER_UTILS_EXPORT inline void DoInitialise();
91 
92  /// Solve the problem.
93  SOLVER_UTILS_EXPORT inline void DoSolve();
94 
95  /// Transform from coefficient to physical space.
97 
98  /// Transform from physical to coefficient space.
100 
101  /// Perform output operations after solve.
102  SOLVER_UTILS_EXPORT inline void Output();
103 
104  /// Linf error computation
105  SOLVER_UTILS_EXPORT inline NekDouble LinfError(unsigned int field, const Array<OneD,NekDouble> &exactsoln = NullNekDouble1DArray);
106 
107  /// Get Session name
109  {
110  return m_sessionName;
111  }
112 
113  template<class T>
114  boost::shared_ptr<T> as()
115  {
116 #if defined __INTEL_COMPILER && BOOST_VERSION > 105200
117  typedef typename boost::shared_ptr<T>::element_type E;
118  E * p = dynamic_cast< E* >( shared_from_this().get() );
119  ASSERTL1(p, "Cannot perform cast");
120  return boost::shared_ptr<T>( shared_from_this(), p );
121 #else
122  return boost::dynamic_pointer_cast<T>( shared_from_this() );
123 #endif
124  }
125 
126  /// Reset Session name
127  SOLVER_UTILS_EXPORT void ResetSessionName(std::string newname)
128  {
129  m_sessionName = newname;
130  }
131 
132  /// Get Session name
134  {
135  return m_session;
136  }
137 
138  /// Get pressure field if available
140 
141  /// Print a summary of parameters and solver characteristics.
142  SOLVER_UTILS_EXPORT inline void PrintSummary(std::ostream &out);
143 
144  /// Set parameter m_lambda
145  SOLVER_UTILS_EXPORT inline void SetLambda(NekDouble lambda);
146 
147  /// Evaluates a function as specified in the session file.
149  Array<OneD, Array<OneD, NekDouble> >& pArray,
150  std::string pFunctionName,
151  const NekDouble pTime = 0.0,
152  const int domain = 0);
153 
154  /// Populate given fields with the function from session.
156  std::vector<std::string> pFieldNames,
157  Array<OneD, Array<OneD, NekDouble> > &pFields,
158  const std::string& pName,
159  const NekDouble& pTime = 0.0,
160  const int domain = 0);
161 
162  /// Populate given fields with the function from session.
164  std::vector<std::string> pFieldNames,
166  const std::string& pName,
167  const NekDouble& pTime = 0.0,
168  const int domain = 0);
169 
170  // Populate an array with a function variable from session.
172  std::string pFieldName,
173  Array<OneD, NekDouble>& pArray,
174  const std::string& pFunctionName,
175  const NekDouble& pTime = 0.0,
176  const int domain = 0);
177 
178  // Describe a function.
180  std::string pFieldName,
181  const std::string &pFunctionName,
182  const int domain);
183 
184  /// Perform initialisation of the base flow.
187 
188  /// Initialise the data in the dependent fields.
190  NekDouble initialtime = 0.0,
191  bool dumpInitialConditions = true,
192  const int domain = 0);
193 
194  /// Evaluates an exact solution
196  int field,
197  Array<OneD, NekDouble> &outfield,
198  const NekDouble time);
199 
200  /// Compute the L2 error between fields and a given exact
201  /// solution.
203  unsigned int field,
204  const Array<OneD,NekDouble> &exactsoln,
205  bool Normalised = false);
206 
207  /// Compute the L2 error of the fields
209  unsigned int field,
210  bool Normalised = false)
211  {
212  return L2Error(field,NullNekDouble1DArray,Normalised);
213  }
214 
215  /// Compute error (L2 and L_inf) over an larger set of quadrature
216  /// points return [L2 Linf]
218  unsigned int field);
219 
220  /// Compute the inner product \f$ (\nabla \phi \cdot F) \f$.
222  const Array<OneD, Array<OneD, NekDouble> > &F,
223  Array<OneD, NekDouble> &outarray);
224 
225  /// Compute the inner product \f$ (\phi, \nabla \cdot F) \f$.
227  const Array<OneD, Array<OneD, NekDouble> > &F,
228  Array<OneD, NekDouble> &outarray);
229 
230  /// Compute the inner product \f$ (\phi, V\cdot \nabla u) \f$.
232  const Array<OneD, Array<OneD, NekDouble> > &V,
234  Array<OneD, NekDouble> &outarray,
235  bool UseContCoeffs = false);
236 
237  /// Compute the non-conservative advection \f$ (V \cdot \nabla u)
238  /// \f$.
240  const Array<OneD, Array<OneD, NekDouble> > &V,
242  Array<OneD, NekDouble> &outarray,
244 
245  /// Calculate the weak discontinuous Galerkin advection.
247  const Array<OneD, Array<OneD, NekDouble> >& InField,
248  Array<OneD, Array<OneD, NekDouble> >& OutField,
249  bool NumericalFluxIncludesNormal = true,
250  bool InFieldIsInPhysSpace = false,
251  int nvariables = 0);
252 
253  /// Calculate weak DG Diffusion in the LDG form.
255  const Array<OneD, Array<OneD, NekDouble> >& InField,
256  Array<OneD, Array<OneD, NekDouble> >& OutField,
257  bool NumericalFluxIncludesNormal = true,
258  bool InFieldIsInPhysSpace = false);
259 
260  /// Write checkpoint file of #m_fields.
261  SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n);
262 
263  /// Write checkpoint file of custom data fields.
265  const int n,
267  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
268  std::vector<std::string> &variables);
269 
270  /// Write base flow file of #m_fields.
271  SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow(const int n);
272 
273  /// Write field data to the given filename.
274  SOLVER_UTILS_EXPORT void WriteFld(const std::string &outname);
275 
276  /// Write input fields to the given filename.
278  const std::string &outname,
280  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
281  std::vector<std::string> &variables);
282 
283  /// Input field data from the given file.
285  const std::string &infile,
287 
288  /// Input field data from the given file to multiple domains
290  const std::string &infile,
292  const int ndomains);
293 
294  /// Output a field.
295  /// Input field data into array from the given file.
297  const std::string &infile,
298  std::vector<std::string> &fieldStr,
299  Array<OneD, Array<OneD, NekDouble> > &coeffs);
300 
301  /// Output a field.
302  /// Input field data into ExpList from the given file.
304  const std::string &infile,
306  std::string &pFieldName);
307 
308  /// Builds map of which element holds each history point.
310 
311  /// Probe each history point and write to file.
312  SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out);
313 
314  /// Write out a session summary.
316 
317  SOLVER_UTILS_EXPORT inline Array<
319 
320 
321  /// Get hold of FieldInfoMap so it can be updated
324 
325  /// Return final time
327 
328  SOLVER_UTILS_EXPORT inline int GetNcoeffs();
329 
330  SOLVER_UTILS_EXPORT inline int GetNcoeffs(const int eid);
331 
332  SOLVER_UTILS_EXPORT inline int GetNumExpModes();
333 
336 
337  SOLVER_UTILS_EXPORT inline int GetNvariables();
338 
339  SOLVER_UTILS_EXPORT inline const std::string
340  GetVariable(unsigned int i);
341 
343 
345 
346  SOLVER_UTILS_EXPORT inline int GetExpSize();
347 
348  SOLVER_UTILS_EXPORT inline int GetPhys_Offset(int n);
349 
350  SOLVER_UTILS_EXPORT inline int GetCoeff_Offset(int n);
351 
352  SOLVER_UTILS_EXPORT inline int GetTotPoints();
353 
354  SOLVER_UTILS_EXPORT inline int GetTotPoints(int n);
355 
356  SOLVER_UTILS_EXPORT inline int GetNpoints();
357 
359 
360  SOLVER_UTILS_EXPORT inline int GetSteps();
361 
363 
364  SOLVER_UTILS_EXPORT inline void CopyFromPhysField(const int i,
365  Array<OneD, NekDouble> &output);
366 
367  SOLVER_UTILS_EXPORT inline void CopyToPhysField(const int i,
368  Array<OneD, NekDouble> &output);
369 
370  SOLVER_UTILS_EXPORT inline void SetSteps(const int steps);
371 
373 
375 
377  const int i,
378  Array<OneD, Array<OneD, NekDouble> >&physfield,
379  Array<OneD, Array<OneD, NekDouble> >&flux);
380 
382  const int i,
383  Array<OneD, Array<OneD, NekDouble> >&physfield,
384  Array<OneD, Array<OneD, NekDouble> >&fluxX,
385  Array<OneD, Array<OneD, NekDouble> > &fluxY);
386 
388  const int i,
389  const int j,
390  Array<OneD, Array<OneD, NekDouble> > &physfield,
391  Array<OneD, Array<OneD, NekDouble> > &flux);
392 
394  Array<OneD, Array<OneD, NekDouble> > &physfield,
395  Array<OneD, Array<OneD, NekDouble> > &numflux);
396 
398  Array<OneD, Array<OneD, NekDouble> > &physfield,
399  Array<OneD, Array<OneD, NekDouble> > &numfluxX,
400  Array<OneD, Array<OneD, NekDouble> > &numfluxY);
401 
403  const Array<OneD, Array<OneD, NekDouble> > &ufield,
404  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &uflux);
405 
407  const Array<OneD, Array<OneD, NekDouble> > &ufield,
408  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &qfield,
409  Array<OneD, Array<OneD, NekDouble> > &qflux);
410 
412  const bool modbasis);
413 
414  /// Perform a case-insensitive string comparison.
416  const std::string & s1, const std::string& s2) ;
417 
419  {
420  return m_nchk;
421  }
422 
424  {
425  m_nchk = num;
426  }
427 
429  {
430  return m_checksteps;
431  }
432 
434  {
435  m_checksteps = num;
436  }
437 
439  const NekDouble time)
440  {
441  m_time = time;
442  }
443 
445  const int step)
446  {
447  m_initialStep = step;
448  }
449 
450  /// Evaluates the boundary conditions at the given time.
452 
453  /// Virtual function to identify if operator is negated in DoSolve
454  SOLVER_UTILS_EXPORT virtual bool v_NegatedOp();
455 
456  protected:
457  /// Communicator
459  /// The session reader
461  /// Field input/output
463  /// Map of interpolator objects
464  std::map<std::string, FieldUtils::Interpolator > m_interpolators;
465  /// pts fields we already read from disk: {funcFilename: (filename, ptsfield)}
466  std::map<std::string, std::pair<std::string, LibUtilities::PtsFieldSharedPtr> > m_loadedPtsFields;
467  // fld fiels already loaded from disk: {funcFilename: (filename, loadedFldField)}
468  std::map<std::string, std::pair<std::string, loadedFldField> > m_loadedFldFields;
469  /// Array holding all dependent variables.
471  /// Base fields.
473  /// Array holding all dependent variables.
475  /// Pointer to boundary conditions object.
477  /// Pointer to graph defining mesh.
479  /// Name of the session.
480  std::string m_sessionName;
481  /// Current time of simulation.
483  /// Number of the step where the simulation should begin
485  /// Finish time of the simulation.
487  /// Time step size
489  /// Lambda constant in real system if one required.
491  /// Time between checkpoints.
493  /// Number of checkpoints written so far
494  int m_nchk;
495  /// Number of steps to take.
496  int m_steps;
497  /// Number of steps between checkpoints.
499  /// Spatial dimension (>= expansion dim).
501  /// Expansion dimension.
502  int m_expdim;
503  /// Flag to determine if single homogeneous mode is used.
505  /// Flag to determine if half homogeneous mode is used.
507  /// Flag to determine if use multiple homogenenous modes are used.
509  /// Flag to determine if FFT is used for homogeneous transform.
510  bool m_useFFT;
511  /**
512  * \brief Flag to determine if dealiasing is used for
513  * homogeneous simulations.
514  */
516  /**
517  * \brief Flag to determine if dealisising is usde for the
518  * Spectral/hp element discretisation.
519  */
521  /// Type of projection; e.g continuous or discontinuous.
523  /// Array holding trace normals for DG simulations in the forwards direction.
525  /// 1 x nvariable x nq
527  /// 2 x m_spacedim x nq
529  /// Flag to indicate if the fields should be checked for singularity.
531  /// Map to identify relevant solver info to dump in output fields
533 
534  /// Number of Quadrature points used to work out the error
536 
537  /// Parameter for homogeneous expansions
539  {
544  };
545 
546 
547 
549 
550  NekDouble m_LhomX; ///< physical length in X direction (if homogeneous)
551  NekDouble m_LhomY; ///< physical length in Y direction (if homogeneous)
552  NekDouble m_LhomZ; ///< physical length in Z direction (if homogeneous)
553 
554  int m_npointsX; ///< number of points in X direction (if homogeneous)
555  int m_npointsY; ///< number of points in Y direction (if homogeneous)
556  int m_npointsZ; ///< number of points in Z direction (if homogeneous)
557 
558  int m_HomoDirec; ///< number of homogenous directions
559 
560 
561  /// Initialises EquationSystem class members.
563 
564  // Here for consistency purposes with old version
565  int nocase_cmp(const std::string & s1, const std::string& s2)
566  {
567  return NoCaseStringCompare(s1,s2);
568  }
569 
570  SOLVER_UTILS_EXPORT virtual void v_InitObject();
571 
572  /// Virtual function for initialisation implementation.
573  SOLVER_UTILS_EXPORT virtual void v_DoInitialise();
574 
575  /// Virtual function for solve implementation.
576  SOLVER_UTILS_EXPORT virtual void v_DoSolve();
577 
578 
579  /// Virtual function for the L_inf error computation between fields and a given exact solution.
581  unsigned int field,
582  const Array<OneD, NekDouble> &exactsoln = NullNekDouble1DArray);
583 
584  /// Virtual function for the L_2 error computation between fields and a given exact solution.
586  unsigned int field,
587  const Array<OneD, NekDouble> &exactsoln = NullNekDouble1DArray,
588  bool Normalised = false);
589 
590  /// Virtual function for transformation to physical space.
592 
593  /// Virtual function for transformation to coefficient space.
595 
596  /// Virtual function for generating summary information.
598 
600  NekDouble initialtime = 0.0,
601  bool dumpInitialConditions = true,
602  const int domain = 0);
603 
605  unsigned int field,
606  Array<OneD, NekDouble> &outfield,
607  const NekDouble time);
608 
609  // Populate an array with a function variable from session.
611  std::string pFieldName,
612  Array<OneD, NekDouble>& pArray,
613  const std::string& pFunctionName,
614  const NekDouble& pTime = 0.0,
615  const int domain = 0);
616 
617  // Populate an array with a function variable from session.
619  std::string pFieldName,
620  Array<OneD, NekDouble>& pArray,
621  const std::string& pFunctionName,
622  const NekDouble& pTime = 0.0,
623  const int domain = 0);
624 
626  std::string pFieldName,
627  Array<OneD, NekDouble>& pArray,
628  const std::string& pFunctionName,
629  const NekDouble& pTime = 0.0,
630  const int domain = 0);
631 
633  std::string funcFilename,
634  std::string filename,
636 
637  //Initialise m_base in order to store the base flow from a file
639 
640  // Fill m_base with the values stored in a fld file
642  std::string pInfile,
644 
645  // Ouptut field information
646  SOLVER_UTILS_EXPORT virtual void v_Output(void);
647 
648  // Get pressure field if available
650 
652  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
653  std::vector<std::string> &variables);
654 
655  private:
656 
659  const int i, Array<OneD,
660  Array<OneD, NekDouble> >&physfield,
661  Array<OneD, Array<OneD, NekDouble> >&flux);
662 
664  const int i, const int j,
665  Array<OneD, Array<OneD, NekDouble> >&physfield,
666  Array<OneD, Array<OneD, NekDouble> >&flux);
667 
669  const int i, Array<OneD,
670  Array<OneD, NekDouble> >&physfield,
671  Array<OneD, Array<OneD, NekDouble> >&fluxX,
672  Array<OneD, Array<OneD, NekDouble> > &fluxY);
673 
675  Array<OneD, Array<OneD, NekDouble> > &physfield,
676  Array<OneD, Array<OneD, NekDouble> > &numflux);
677 
679  Array<OneD, Array<OneD, NekDouble> > &physfield,
680  Array<OneD, Array<OneD, NekDouble> > &numfluxX,
681  Array<OneD, Array<OneD, NekDouble> > &numfluxY);
682 
684  const Array<OneD, Array<OneD, NekDouble> > &ufield,
685  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &uflux);
686 
688  const Array<OneD, Array<OneD, NekDouble> > &ufield,
689  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &qfield,
690  Array<OneD, Array<OneD, NekDouble > > &qflux);
691 
692  SOLVER_UTILS_EXPORT void PrintProgressbar(const int position,
693  const int goal) const
694  {
695  LibUtilities::PrintProgressbar(position, goal, "Interpolating");
696  }
697  };
698 
699 
700  /**
701  * This is the second part of the two-phase initialisation process.
702  * Calls to virtual functions will correctly resolve to the derived class
703  * during this phase of the construction.
704  */
706  {
707  v_InitObject();
708  }
709 
710 
711  /**
712  * This allows initialisation of the solver which cannot be completed
713  * during object construction (such as setting of initial conditions).
714  *
715  * Public interface routine to virtual function implementation.
716  */
718  {
719  v_DoInitialise();
720  }
721 
722 
723  /**
724  * Performs the transformation from coefficient to physical space.
725  *
726  * Public interface routine to virtual function implementation.
727  */
729  {
731  }
732 
733  /**
734  * Performs the transformation from physical to coefficient space.
735  *
736  * Public interface routine to virtual function implementation.
737  */
739  {
741  }
742 
743 
744  /**
745  * Performs the actual solve.
746  *
747  * Public interface routine to virtual function implementation.
748  */
749  inline void EquationSystem::DoSolve(void)
750  {
751  v_DoSolve();
752  }
753 
754 
755  /**
756  * Perform output operations after solve.
757  */
758  inline void EquationSystem::Output(void)
759  {
760  v_Output();
761  }
762 
763  /**
764  * L_inf Error computation
765  * Public interface routine to virtual function implementation.
766  */
767  inline NekDouble EquationSystem::LinfError(unsigned int field, const Array<OneD,NekDouble> &exactsoln)
768  {
769  return v_LinfError(field, exactsoln);
770  }
771 
772  /**
773  * L_2 Error computation
774  * Public interface routine to virtual function implementation.
775  */
776  inline NekDouble EquationSystem::L2Error(unsigned int field, const Array<OneD,NekDouble> &exactsoln, bool Normalised)
777  {
778  return v_L2Error(field, exactsoln, Normalised);
779  }
780 
781  /**
782  * Get Pressure field if available
783  */
785  {
786  return v_GetPressure();
787  }
788 
789  /**
790  * Prints a summary of variables and problem parameters.
791  *
792  * Public interface routine to virtual function implementation.
793  *
794  * @param out The ostream object to write to.
795  */
796  inline void EquationSystem::PrintSummary(std::ostream &out)
797  {
798  if (m_session->GetComm()->GetRank() == 0)
799  {
800  std::vector<std::pair<std::string, std::string> > vSummary;
801  v_GenerateSummary(vSummary);
802 
803  out << "=======================================================================" << std::endl;
804  SummaryList::const_iterator x;
805  for (x = vSummary.begin(); x != vSummary.end(); ++x)
806  {
807  out << "\t";
808  out.width(20);
809  out << x->first << ": " << x->second << std::endl;
810  }
811  out << "=======================================================================" << std::endl;
812  }
813  }
814 
816  {
817  m_lambda = lambda;
818  }
819 
821  bool dumpInitialConditions,
822  const int domain)
823  {
824  v_SetInitialConditions(initialtime,dumpInitialConditions,domain);
825  }
826 
827  /// Evaluates an exact solution
829  Array<OneD, NekDouble> &outfield,
830  const NekDouble time)
831  {
832  v_EvaluateExactSolution(field, outfield, time);
833  }
834 
836  {
837  return m_fields;
838  }
839 
840  /// Return final time
842  {
843  return m_time;
844  }
845 
846  inline int EquationSystem::GetNcoeffs(void)
847  {
848  return m_fields[0]->GetNcoeffs();
849  }
850 
851  inline int EquationSystem::GetNcoeffs(const int eid)
852  {
853  return m_fields[0]->GetNcoeffs(eid);
854  }
855 
857  {
858  return m_graph->GetExpansions().begin()->second->m_basisKeyVector[0]
859  .GetNumModes();
860  }
861 
863  {
864  return m_fields[0]->EvalBasisNumModesMaxPerExp();
865  }
866 
868  {
869  return m_session->GetVariables().size();
870  }
871 
872  inline const std::string EquationSystem::GetVariable(unsigned int i)
873  {
874  return m_session->GetVariable(i);
875  }
876 
878  {
879  return GetTraceNpoints();
880  }
881 
883  {
884  return m_fields[0]->GetTrace()->GetNpoints();
885  }
886 
887  inline int EquationSystem::GetExpSize(void)
888  {
889  return m_fields[0]->GetExpSize();
890  }
891 
893  {
894  return m_fields[0]->GetPhys_Offset(n);
895  }
896 
898  {
899  return m_fields[0]->GetCoeff_Offset(n);
900  }
901 
903  {
904  return m_fields[0]->GetNpoints();
905  }
906 
908  {
909  return m_fields[0]->GetTotPoints(n);
910  }
911 
912  inline int EquationSystem::GetNpoints(void)
913  {
914  return m_fields[0]->GetNpoints();
915  }
916 
918  {
919  return (m_fields.num_elements() - 1);
920  }
921 
922  inline int EquationSystem::GetSteps(void)
923  {
924  return m_steps;
925  }
926 
928  {
929  return m_timestep;
930  }
931 
932  inline void EquationSystem::SetSteps(const int steps)
933  {
934  m_steps = steps;
935  }
936 
937  inline void EquationSystem::CopyFromPhysField(const int i,
938  Array<OneD, NekDouble> &output)
939  {
940  Vmath::Vcopy(output.num_elements(), m_fields[i]->GetPhys(), 1, output, 1 );
941  }
942 
943  inline void EquationSystem::CopyToPhysField(const int i,
944  Array<OneD, NekDouble> &output)
945  {
946  Vmath::Vcopy(output.num_elements(), output, 1, m_fields[i]->UpdatePhys(), 1 );
947  }
948 
949  inline void EquationSystem::GetFluxVector(const int i,
950  Array<OneD, Array<OneD, NekDouble> >&physfield,
952  {
953  v_GetFluxVector(i,physfield, flux);
954  }
955 
956  inline void EquationSystem::GetFluxVector(const int i,
957  Array<OneD, Array<OneD, NekDouble> >&physfield,
960  {
961  v_GetFluxVector(i,physfield, fluxX, fluxY);
962  }
963 
964  inline void EquationSystem::GetFluxVector(const int i, const int j,
965  Array<OneD, Array<OneD, NekDouble> > &physfield,
967  {
968  v_GetFluxVector(i,j,physfield,flux);
969  }
970 
972  Array<OneD, Array<OneD, NekDouble> > &numflux)
973  {
974  v_NumericalFlux(physfield, numflux);
975  }
976 
978  Array<OneD, Array<OneD, NekDouble> > &numfluxX,
979  Array<OneD, Array<OneD, NekDouble> > &numfluxY)
980  {
981  v_NumericalFlux(physfield, numfluxX, numfluxY);
982  }
983 
985  const Array<OneD, Array<OneD, NekDouble> > &ufield,
987  {
988  v_NumFluxforScalar(ufield, uflux);
989  }
990 
992  const Array<OneD, Array<OneD, NekDouble> > &ufield,
995  {
996  v_NumFluxforVector(ufield, qfield, qflux);
997  }
998  }
999 }
1000 
1001 #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.