Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator 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 
45 #include <MultiRegions/ExpList.h>
47 #include <SolverUtils/Core/Misc.h>
48 
49 namespace Nektar
50 {
51  namespace SolverUtils
52  {
54 
55  /// A shared pointer to an EquationSystem object
56  typedef boost::shared_ptr<EquationSystem> EquationSystemSharedPtr;
57  /// Datatype of the NekFactory used to instantiate classes derived from
58  /// the EquationSystem class.
60  std::string, EquationSystem,
64 
65  /// A base class for describing how to solve specific equations.
66  class EquationSystem : public boost::enable_shared_from_this<EquationSystem>
67  {
68  public:
69  /// Destructor
71 
72  // Set up trace normals if required
74 
75  /// Initialises the members of this object.
76  SOLVER_UTILS_EXPORT inline void InitObject();
77 
78  /// Perform any initialisation necessary before solving the problem.
79  SOLVER_UTILS_EXPORT inline void DoInitialise();
80 
81  /// Solve the problem.
82  SOLVER_UTILS_EXPORT inline void DoSolve();
83 
84  /// Transform from coefficient to physical space.
86 
87  /// Transform from physical to coefficient space.
89 
90  /// Perform output operations after solve.
91  SOLVER_UTILS_EXPORT inline void Output();
92 
93  /// Linf error computation
94  SOLVER_UTILS_EXPORT inline NekDouble LinfError(unsigned int field, const Array<OneD,NekDouble> &exactsoln = NullNekDouble1DArray);
95 
96  /// Get Session name
98  {
99  return m_sessionName;
100  }
101 
102  template<class T>
103  boost::shared_ptr<T> as()
104  {
105 #if defined __INTEL_COMPILER && BOOST_VERSION > 105200
106  typedef typename boost::shared_ptr<T>::element_type E;
107  E * p = dynamic_cast< E* >( shared_from_this().get() );
108  ASSERTL1(p, "Cannot perform cast");
109  return boost::shared_ptr<T>( shared_from_this(), p );
110 #else
111  return boost::dynamic_pointer_cast<T>( shared_from_this() );
112 #endif
113  }
114 
115  /// Reset Session name
116  SOLVER_UTILS_EXPORT void ResetSessionName(std::string newname)
117  {
118  m_sessionName = newname;
119  }
120 
121  /// Get Session name
123  {
124  return m_session;
125  }
126 
127  /// Get pressure field if available
129 
130  /// Print a summary of parameters and solver characteristics.
131  SOLVER_UTILS_EXPORT inline void PrintSummary(std::ostream &out);
132 
133  /// Set parameter m_lambda
134  SOLVER_UTILS_EXPORT inline void SetLambda(NekDouble lambda);
135 
136  /// Evaluates a function as specified in the session file.
138  Array<OneD, Array<OneD, NekDouble> >& pArray,
139  std::string pFunctionName,
140  const NekDouble pTime = 0.0,
141  const int domain = 0);
142 
143  /// Populate given fields with the function from session.
145  std::vector<std::string> pFieldNames,
146  Array<OneD, Array<OneD, NekDouble> > &pFields,
147  const std::string& pName,
148  const int domain = 0);
149 
150  /// Populate given fields with the function from session.
152  std::vector<std::string> pFieldNames,
153  Array<OneD, MultiRegions::ExpListSharedPtr> &pFields,
154  const std::string& pName,
155  const int domain = 0);
156 
157  // Populate an array with a function variable from session.
159  std::string pFieldName,
160  Array<OneD, NekDouble>& pArray,
161  const std::string& pFunctionName,
162  const NekDouble& pTime = 0.0,
163  const int domain = 0);
164 
165  // Describe a function.
167  std::string pFieldName,
168  const std::string &pFunctionName,
169  const int domain);
170 
171  /// Perform initialisation of the base flow.
173  Array<OneD, Array<OneD, NekDouble> > &base);
174 
175  /// Initialise the data in the dependent fields.
177  NekDouble initialtime = 0.0,
178  bool dumpInitialConditions = true,
179  const int domain = 0);
180 
181  /// Evaluates an exact solution
183  int field,
184  Array<OneD, NekDouble> &outfield,
185  const NekDouble time);
186 
187  /// Compute the L2 error between fields and a given exact
188  /// solution.
190  unsigned int field,
191  const Array<OneD,NekDouble> &exactsoln,
192  bool Normalised = false);
193 
194  /// Compute the L2 error of the fields
196  unsigned int field,
197  bool Normalised = false)
198  {
199  return L2Error(field,NullNekDouble1DArray,Normalised);
200  }
201 
202  /// Compute error (L2 and L_inf) over an larger set of quadrature
203  /// points return [L2 Linf]
204  SOLVER_UTILS_EXPORT Array<OneD,NekDouble> ErrorExtraPoints(
205  unsigned int field);
206 
207  /// Compute the inner product \f$ (\nabla \phi \cdot F) \f$.
209  const Array<OneD, Array<OneD, NekDouble> > &F,
210  Array<OneD, NekDouble> &outarray);
211 
212  /// Compute the inner product \f$ (\phi, \nabla \cdot F) \f$.
214  const Array<OneD, Array<OneD, NekDouble> > &F,
215  Array<OneD, NekDouble> &outarray);
216 
217  /// Compute the inner product \f$ (\phi, V\cdot \nabla u) \f$.
219  const Array<OneD, Array<OneD, NekDouble> > &V,
220  const Array<OneD, const NekDouble> &u,
221  Array<OneD, NekDouble> &outarray,
222  bool UseContCoeffs = false);
223 
224  /// Compute the non-conservative advection \f$ (V \cdot \nabla u)
225  /// \f$.
227  const Array<OneD, Array<OneD, NekDouble> > &V,
228  const Array<OneD, const NekDouble> &u,
229  Array<OneD, NekDouble> &outarray,
230  Array<OneD, NekDouble> &wk = NullNekDouble1DArray);
231 
232  /// Calculate the weak discontinuous Galerkin advection.
234  const Array<OneD, Array<OneD, NekDouble> >& InField,
235  Array<OneD, Array<OneD, NekDouble> >& OutField,
236  bool NumericalFluxIncludesNormal = true,
237  bool InFieldIsInPhysSpace = false,
238  int nvariables = 0);
239 
240  /// Calculate weak DG Diffusion in the LDG form.
242  const Array<OneD, Array<OneD, NekDouble> >& InField,
243  Array<OneD, Array<OneD, NekDouble> >& OutField,
244  bool NumericalFluxIncludesNormal = true,
245  bool InFieldIsInPhysSpace = false);
246 
247  /// Write checkpoint file of #m_fields.
248  SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n);
249 
250  /// Write checkpoint file of custom data fields.
252  const int n,
254  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
255  std::vector<std::string> &variables);
256 
257  /// Write base flow file of #m_fields.
258  SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow(const int n);
259 
260  /// Write field data to the given filename.
261  SOLVER_UTILS_EXPORT void WriteFld(const std::string &outname);
262 
263  /// Write input fields to the given filename.
265  const std::string &outname,
267  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
268  std::vector<std::string> &variables);
269 
270  /// Input field data from the given file.
272  const std::string &infile,
273  Array<OneD, MultiRegions::ExpListSharedPtr> &pFields);
274 
275  /// Input field data from the given file to multiple domains
277  const std::string &infile,
278  Array<OneD, MultiRegions::ExpListSharedPtr> &pFields,
279  const int ndomains);
280 
281  /// Output a field.
282  /// Input field data into array from the given file.
284  const std::string &infile,
285  std::vector<std::string> &fieldStr,
286  Array<OneD, Array<OneD, NekDouble> > &coeffs);
287 
288  /// Output a field.
289  /// Input field data into ExpList from the given file.
291  const std::string &infile,
293  std::string &pFieldName);
294 
295  /// Builds map of which element holds each history point.
297 
298  /// Probe each history point and write to file.
299  SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out);
300 
301  /// Write out a session summary.
303 
304  SOLVER_UTILS_EXPORT inline Array<
306 
307 
308  /// Get hold of FieldInfoMap so it can be updated
311 
312  /// Return final time
314 
315  SOLVER_UTILS_EXPORT inline int GetNcoeffs();
316 
317  SOLVER_UTILS_EXPORT inline int GetNcoeffs(const int eid);
318 
319  SOLVER_UTILS_EXPORT inline int GetNumExpModes();
320 
321  SOLVER_UTILS_EXPORT inline const Array<OneD,int>
323 
324  SOLVER_UTILS_EXPORT inline int GetNvariables();
325 
326  SOLVER_UTILS_EXPORT inline const std::string
327  GetVariable(unsigned int i);
328 
330 
332 
333  SOLVER_UTILS_EXPORT inline int GetExpSize();
334 
335  SOLVER_UTILS_EXPORT inline int GetPhys_Offset(int n);
336 
337  SOLVER_UTILS_EXPORT inline int GetCoeff_Offset(int n);
338 
339  SOLVER_UTILS_EXPORT inline int GetTotPoints();
340 
341  SOLVER_UTILS_EXPORT inline int GetTotPoints(int n);
342 
343  SOLVER_UTILS_EXPORT inline int GetNpoints();
344 
346 
347  SOLVER_UTILS_EXPORT inline int GetSteps();
348 
350 
351  SOLVER_UTILS_EXPORT inline void CopyFromPhysField(const int i,
352  Array<OneD, NekDouble> &output);
353 
354  SOLVER_UTILS_EXPORT inline void CopyToPhysField(const int i,
355  Array<OneD, NekDouble> &output);
356 
357  SOLVER_UTILS_EXPORT inline void SetStepsToOne();
358 
360 
362 
364  const int i,
365  Array<OneD, Array<OneD, NekDouble> >&physfield,
366  Array<OneD, Array<OneD, NekDouble> >&flux);
367 
369  const int i,
370  Array<OneD, Array<OneD, NekDouble> >&physfield,
371  Array<OneD, Array<OneD, NekDouble> >&fluxX,
372  Array<OneD, Array<OneD, NekDouble> > &fluxY);
373 
375  const int i,
376  const int j,
377  Array<OneD, Array<OneD, NekDouble> > &physfield,
378  Array<OneD, Array<OneD, NekDouble> > &flux);
379 
381  Array<OneD, Array<OneD, NekDouble> > &physfield,
382  Array<OneD, Array<OneD, NekDouble> > &numflux);
383 
385  Array<OneD, Array<OneD, NekDouble> > &physfield,
386  Array<OneD, Array<OneD, NekDouble> > &numfluxX,
387  Array<OneD, Array<OneD, NekDouble> > &numfluxY);
388 
390  const Array<OneD, Array<OneD, NekDouble> > &ufield,
391  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &uflux);
392 
394  const Array<OneD, Array<OneD, NekDouble> > &ufield,
395  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &qfield,
396  Array<OneD, Array<OneD, NekDouble> > &qflux);
397 
399  const bool modbasis);
400 
401  /// Perform a case-insensitive string comparison.
403  const string & s1, const string& s2) ;
404 
405  protected:
406  /// Communicator
408  /// The session reader
410  /// Field input/output
412  /// Array holding all dependent variables.
413  Array<OneD, MultiRegions::ExpListSharedPtr> m_fields;
414  /// Base fields.
415  Array<OneD, MultiRegions::ExpListSharedPtr> m_base;
416  /// Array holding all dependent variables.
417  Array<OneD, MultiRegions::ExpListSharedPtr> m_derivedfields;
418  /// Pointer to boundary conditions object.
420  /// Pointer to graph defining mesh.
422  /// Name of the session.
423  std::string m_sessionName;
424  /// Current time of simulation.
426  /// Finish time of the simulation.
428  /// Time step size
430  /// Lambda constant in real system if one required.
432  /// Time between checkpoints.
434  /// Number of steps to take.
435  int m_steps;
436  /// Number of steps between checkpoints.
438  /// Spatial dimension (>= expansion dim).
440  /// Expansion dimension.
441  int m_expdim;
442  /// Flag to determine if single homogeneous mode is used.
444  /// Flag to determine if half homogeneous mode is used.
446  /// Flag to determine if use multiple homogenenous modes are used.
448  /// Flag to determine if FFT is used for homogeneous transform.
449  bool m_useFFT;
450  /**
451  * \brief Flag to determine if dealiasing is used for
452  * homogeneous simulations.
453  */
455  /**
456  * \brief Flag to determine if dealisising is usde for the
457  * Spectral/hp element discretisation.
458  */
460  /// Type of projection; e.g continuous or discontinuous.
462  /// Array holding trace normals for DG simulations in the forwards direction.
463  Array<OneD, Array<OneD, NekDouble> > m_traceNormals;
464  /// 1 x nvariable x nq
465  Array<OneD, Array<OneD, Array<OneD,NekDouble> > > m_gradtan;
466  /// 2 x m_spacedim x nq
467  Array<OneD, Array<OneD, Array<OneD,NekDouble> > > m_tanbasis;
468  /// Flag to indicate if the fields should be checked for singularity.
469  Array<OneD, bool> m_checkIfSystemSingular;
470  /// Map to identify relevant solver info to dump in output fields
472 
473  /// Number of Quadrature points used to work out the error
475 
476  /// Parameter for homogeneous expansions
478  {
483  };
484 
485 
486 
488 
489  NekDouble m_LhomX; ///< physical length in X direction (if homogeneous)
490  NekDouble m_LhomY; ///< physical length in Y direction (if homogeneous)
491  NekDouble m_LhomZ; ///< physical length in Z direction (if homogeneous)
492 
493  int m_npointsX; ///< number of points in X direction (if homogeneous)
494  int m_npointsY; ///< number of points in Y direction (if homogeneous)
495  int m_npointsZ; ///< number of points in Z direction (if homogeneous)
496 
497  int m_HomoDirec; ///< number of homogenous directions
498 
499  int m_NumMode; ///< Mode to use in case of single mode analysis
500 
501 
502  /// Initialises EquationSystem class members.
504 
505  // Here for consistency purposes with old version
506  int nocase_cmp(const string & s1, const string& s2)
507  {
508  return NoCaseStringCompare(s1,s2);
509  }
510 
511  SOLVER_UTILS_EXPORT virtual void v_InitObject();
512 
513  /// Evaluates the boundary conditions at the given time.
515 
516  /// Virtual function for initialisation implementation.
517  SOLVER_UTILS_EXPORT virtual void v_DoInitialise();
518 
519  /// Virtual function for solve implementation.
520  SOLVER_UTILS_EXPORT virtual void v_DoSolve();
521 
522  /// Virtual function for the L_inf error computation between fields and a given exact solution.
524  unsigned int field,
525  const Array<OneD, NekDouble> &exactsoln = NullNekDouble1DArray);
526 
527  /// Virtual function for the L_2 error computation between fields and a given exact solution.
529  unsigned int field,
530  const Array<OneD, NekDouble> &exactsoln = NullNekDouble1DArray,
531  bool Normalised = false);
532 
533  /// Virtual function for transformation to physical space.
535 
536  /// Virtual function for transformation to coefficient space.
538 
539  /// Virtual function for generating summary information.
541 
543  NekDouble initialtime = 0.0,
544  bool dumpInitialConditions = true,
545  const int domain = 0);
546 
548  unsigned int field,
549  Array<OneD, NekDouble> &outfield,
550  const NekDouble time);
551 
552  //Initialise m_base in order to store the base flow from a file
554 
555  // Fill m_base with the values stored in a fld file
557  std::string pInfile,
559 
560  // Ouptut field information
561  SOLVER_UTILS_EXPORT virtual void v_Output(void);
562 
563  // Get pressure field if available
565 
567  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
568  std::vector<std::string> &variables);
569 
570  private:
571 
572  SOLVER_UTILS_EXPORT virtual Array<OneD, bool> v_GetSystemSingularChecks();
574  const int i, Array<OneD,
575  Array<OneD, NekDouble> >&physfield,
576  Array<OneD, Array<OneD, NekDouble> >&flux);
577 
579  const int i, const int j,
580  Array<OneD, Array<OneD, NekDouble> >&physfield,
581  Array<OneD, Array<OneD, NekDouble> >&flux);
582 
584  const int i, Array<OneD,
585  Array<OneD, NekDouble> >&physfield,
586  Array<OneD, Array<OneD, NekDouble> >&fluxX,
587  Array<OneD, Array<OneD, NekDouble> > &fluxY);
588 
590  Array<OneD, Array<OneD, NekDouble> > &physfield,
591  Array<OneD, Array<OneD, NekDouble> > &numflux);
592 
594  Array<OneD, Array<OneD, NekDouble> > &physfield,
595  Array<OneD, Array<OneD, NekDouble> > &numfluxX,
596  Array<OneD, Array<OneD, NekDouble> > &numfluxY);
597 
599  const Array<OneD, Array<OneD, NekDouble> > &ufield,
600  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &uflux);
601 
603  const Array<OneD, Array<OneD, NekDouble> > &ufield,
604  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &qfield,
605  Array<OneD, Array<OneD, NekDouble > > &qflux);
606  };
607 
608 
609  /**
610  * This is the second part of the two-phase initialisation process.
611  * Calls to virtual functions will correctly resolve to the derived class
612  * during this phase of the construction.
613  */
615  {
616  v_InitObject();
617  }
618 
619 
620  /**
621  * This allows initialisation of the solver which cannot be completed
622  * during object construction (such as setting of initial conditions).
623  *
624  * Public interface routine to virtual function implementation.
625  */
627  {
628  v_DoInitialise();
629  }
630 
631 
632  /**
633  * Performs the transformation from coefficient to physical space.
634  *
635  * Public interface routine to virtual function implementation.
636  */
638  {
640  }
641 
642  /**
643  * Performs the transformation from physical to coefficient space.
644  *
645  * Public interface routine to virtual function implementation.
646  */
648  {
650  }
651 
652 
653  /**
654  * Performs the actual solve.
655  *
656  * Public interface routine to virtual function implementation.
657  */
658  inline void EquationSystem::DoSolve(void)
659  {
660  v_DoSolve();
661  }
662 
663 
664  /**
665  * Perform output operations after solve.
666  */
667  inline void EquationSystem::Output(void)
668  {
669  v_Output();
670  }
671 
672  /**
673  * L_inf Error computation
674  * Public interface routine to virtual function implementation.
675  */
676  inline NekDouble EquationSystem::LinfError(unsigned int field, const Array<OneD,NekDouble> &exactsoln)
677  {
678  return v_LinfError(field, exactsoln);
679  }
680 
681  /**
682  * L_2 Error computation
683  * Public interface routine to virtual function implementation.
684  */
685  inline NekDouble EquationSystem::L2Error(unsigned int field, const Array<OneD,NekDouble> &exactsoln, bool Normalised)
686  {
687  return v_L2Error(field, exactsoln, Normalised);
688  }
689 
690  /**
691  * Get Pressure field if available
692  */
694  {
695  return v_GetPressure();
696  }
697 
698  /**
699  * Prints a summary of variables and problem parameters.
700  *
701  * Public interface routine to virtual function implementation.
702  *
703  * @param out The ostream object to write to.
704  */
705  inline void EquationSystem::PrintSummary(std::ostream &out)
706  {
707  if (m_session->GetComm()->GetRank() == 0)
708  {
709  std::vector<std::pair<std::string, std::string> > vSummary;
710  v_GenerateSummary(vSummary);
711 
712  out << "=======================================================================" << endl;
713  SummaryList::const_iterator x;
714  for (x = vSummary.begin(); x != vSummary.end(); ++x)
715  {
716  out << "\t";
717  out.width(20);
718  out << x->first << ": " << x->second << endl;
719  }
720  out << "=======================================================================" << endl;
721  }
722  }
723 
725  {
726  m_lambda = lambda;
727  }
728 
730  bool dumpInitialConditions,
731  const int domain)
732  {
733  v_SetInitialConditions(initialtime,dumpInitialConditions,domain);
734  }
735 
736  /// Evaluates an exact solution
738  Array<OneD, NekDouble> &outfield,
739  const NekDouble time)
740  {
741  v_EvaluateExactSolution(field, outfield, time);
742  }
743 
744  inline Array<OneD, MultiRegions::ExpListSharedPtr> &EquationSystem::UpdateFields(void)
745  {
746  return m_fields;
747  }
748 
749  /// Return final time
751  {
752  return m_time;
753  }
754 
755  inline int EquationSystem::GetNcoeffs(void)
756  {
757  return m_fields[0]->GetNcoeffs();
758  }
759 
760  inline int EquationSystem::GetNcoeffs(const int eid)
761  {
762  return m_fields[0]->GetNcoeffs(eid);
763  }
764 
766  {
767  return m_graph->GetExpansions().begin()->second->m_basisKeyVector[0]
768  .GetNumModes();
769  }
770 
771  inline const Array<OneD,int> EquationSystem::GetNumExpModesPerExp(void)
772  {
773  return m_fields[0]->EvalBasisNumModesMaxPerExp();
774  }
775 
777  {
778  return m_session->GetVariables().size();
779  }
780 
781  inline const std::string EquationSystem::GetVariable(unsigned int i)
782  {
783  return m_session->GetVariable(i);
784  }
785 
787  {
788  return GetTraceNpoints();
789  }
790 
792  {
793  return m_fields[0]->GetTrace()->GetNpoints();
794  }
795 
796  inline int EquationSystem::GetExpSize(void)
797  {
798  return m_fields[0]->GetExpSize();
799  }
800 
802  {
803  return m_fields[0]->GetPhys_Offset(n);
804  }
805 
807  {
808  return m_fields[0]->GetCoeff_Offset(n);
809  }
810 
812  {
813  return m_fields[0]->GetNpoints();
814  }
815 
817  {
818  return m_fields[0]->GetTotPoints(n);
819  }
820 
821  inline int EquationSystem::GetNpoints(void)
822  {
823  return m_fields[0]->GetNpoints();
824  }
825 
827  {
828  return (m_fields.num_elements() - 1);
829  }
830 
831  inline int EquationSystem::GetSteps(void)
832  {
833  return m_steps;
834  }
835 
837  {
838  return m_timestep;
839  }
840 
842  {
843  m_steps=1;
844  }
845 
846  inline void EquationSystem::CopyFromPhysField(const int i,
847  Array<OneD, NekDouble> &output)
848  {
849  Vmath::Vcopy(output.num_elements(), m_fields[i]->GetPhys(), 1, output, 1 );
850  }
851 
852  inline void EquationSystem::CopyToPhysField(const int i,
853  Array<OneD, NekDouble> &output)
854  {
855  Vmath::Vcopy(output.num_elements(), output, 1, m_fields[i]->UpdatePhys(), 1 );
856  }
857 
858  inline void EquationSystem::GetFluxVector(const int i,
859  Array<OneD, Array<OneD, NekDouble> >&physfield,
860  Array<OneD, Array<OneD, NekDouble> >&flux)
861  {
862  v_GetFluxVector(i,physfield, flux);
863  }
864 
865  inline void EquationSystem::GetFluxVector(const int i,
866  Array<OneD, Array<OneD, NekDouble> >&physfield,
867  Array<OneD, Array<OneD, NekDouble> >&fluxX,
868  Array<OneD, Array<OneD, NekDouble> > &fluxY)
869  {
870  v_GetFluxVector(i,physfield, fluxX, fluxY);
871  }
872 
873  inline void EquationSystem::GetFluxVector(const int i, const int j,
874  Array<OneD, Array<OneD, NekDouble> > &physfield,
875  Array<OneD, Array<OneD, NekDouble> > &flux)
876  {
877  v_GetFluxVector(i,j,physfield,flux);
878  }
879 
880  inline void EquationSystem::NumericalFlux(Array<OneD, Array<OneD, NekDouble> > &physfield,
881  Array<OneD, Array<OneD, NekDouble> > &numflux)
882  {
883  v_NumericalFlux(physfield, numflux);
884  }
885 
886  inline void EquationSystem::NumericalFlux(Array<OneD, Array<OneD, NekDouble> > &physfield,
887  Array<OneD, Array<OneD, NekDouble> > &numfluxX,
888  Array<OneD, Array<OneD, NekDouble> > &numfluxY)
889  {
890  v_NumericalFlux(physfield, numfluxX, numfluxY);
891  }
892 
894  const Array<OneD, Array<OneD, NekDouble> > &ufield,
895  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &uflux)
896  {
897  v_NumFluxforScalar(ufield, uflux);
898  }
899 
901  const Array<OneD, Array<OneD, NekDouble> > &ufield,
902  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > &qfield,
903  Array<OneD, Array<OneD, NekDouble> > &qflux)
904  {
905  v_NumFluxforVector(ufield, qfield, qflux);
906  }
907  }
908 }
909 
910 #endif