Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Member Functions | Friends | List of all members
Nektar::CFLtester Class Reference

#include <CFLtester.h>

Inheritance diagram for Nektar::CFLtester:
[legend]

Public Member Functions

virtual ~CFLtester ()
 
- Public Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
SOLVER_UTILS_EXPORT AdvectionSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
virtual SOLVER_UTILS_EXPORT ~AdvectionSystem ()
 
SOLVER_UTILS_EXPORT AdvectionSharedPtr GetAdvObject ()
 Returns the advection object held by this instance. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleGetElmtCFLVals (void)
 
SOLVER_UTILS_EXPORT NekDouble GetCFLEstimate (int &elmtid)
 
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep (const Array< OneD, const Array< OneD, NekDouble >> &inarray)
 Calculate the larger time-step mantaining the problem stable. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void SetUpTraceNormals (void)
 
SOLVER_UTILS_EXPORT void InitObject ()
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise ()
 Perform any initialisation necessary before solving the problem. More...
 
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem. More...
 
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space. More...
 
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space. More...
 
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve. More...
 
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. More...
 
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name. More...
 
template<class T >
std::shared_ptr< T > as ()
 
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name. More...
 
SOLVER_UTILS_EXPORT LibUtilities::SessionReaderSharedPtr GetSession ()
 Get Session name. More...
 
SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr GetPressure ()
 Get pressure field if available. More...
 
SOLVER_UTILS_EXPORT void ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics. More...
 
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda. More...
 
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction (std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
 Get a SessionFunction by name. More...
 
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields. More...
 
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution. More...
 
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. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf]. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write checkpoint file of custom data fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write input fields to the given filename. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file. More...
 
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. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, std::vector< std::string > &fieldStr, Array< OneD, Array< OneD, NekDouble > > &coeffs)
 Output a field. Input field data into array from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, MultiRegions::ExpListSharedPtr &pField, std::string &pFieldName)
 Output a field. Input field data into ExpList from the given file. More...
 
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary. More...
 
SOLVER_UTILS_EXPORT Array< OneD, MultiRegions::ExpListSharedPtr > & UpdateFields ()
 
SOLVER_UTILS_EXPORT LibUtilities::FieldMetaDataMapUpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated. More...
 
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time. More...
 
SOLVER_UTILS_EXPORT int GetNcoeffs ()
 
SOLVER_UTILS_EXPORT int GetNcoeffs (const int eid)
 
SOLVER_UTILS_EXPORT int GetNumExpModes ()
 
SOLVER_UTILS_EXPORT const Array< OneD, int > GetNumExpModesPerExp ()
 
SOLVER_UTILS_EXPORT int GetNvariables ()
 
SOLVER_UTILS_EXPORT const std::string GetVariable (unsigned int i)
 
SOLVER_UTILS_EXPORT int GetTraceTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTraceNpoints ()
 
SOLVER_UTILS_EXPORT int GetExpSize ()
 
SOLVER_UTILS_EXPORT int GetPhys_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetCoeff_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTotPoints (int n)
 
SOLVER_UTILS_EXPORT int GetNpoints ()
 
SOLVER_UTILS_EXPORT int GetSteps ()
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep ()
 
SOLVER_UTILS_EXPORT void CopyFromPhysField (const int i, Array< OneD, NekDouble > &output)
 
SOLVER_UTILS_EXPORT void CopyToPhysField (const int i, Array< OneD, NekDouble > &output)
 
SOLVER_UTILS_EXPORT void SetSteps (const int steps)
 
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
 
SOLVER_UTILS_EXPORT void FwdTransFields ()
 
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
 
SOLVER_UTILS_EXPORT int GetCheckpointNumber ()
 
SOLVER_UTILS_EXPORT void SetCheckpointNumber (int num)
 
SOLVER_UTILS_EXPORT int GetCheckpointSteps ()
 
SOLVER_UTILS_EXPORT void SetCheckpointSteps (int num)
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to identify if operator is negated in DoSolve. More...
 

Static Public Member Functions

static EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 

Static Public Attributes

static std::string className = GetEquationSystemFactory().RegisterCreatorFunction("CFLtester", CFLtester::create, "Testing CFL restriction")
 

Protected Member Functions

 CFLtester (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
void DoOdeRhs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 
void DoOdeProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 
Array< OneD, NekDouble > & GetNormalVelocity ()
 Get the normal velocity. More...
 
virtual void v_InitObject ()
 Init object for UnsteadySystem class. More...
 
void GetFluxVector (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
 
virtual void v_GenerateSummary (SummaryList &s)
 Print a summary of time stepping parameters. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT Array< OneD, NekDoublev_GetMaxStdVelocity ()
 
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise ()
 Sets up initial conditions. More...
 
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D (Array< OneD, Array< OneD, NekDouble >> &solution1D)
 Print the solution at each solution point in a txt file. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble >> &inarray)
 Return the timestep to be used for the next step in the time-marching loop. More...
 
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_RequireFwdTrans ()
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time, int &nchk)
 
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff (const Array< OneD, Array< OneD, NekDouble >> vel, StdRegions::VarCoeffMap &varCoeffMap)
 Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h time velocity. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises EquationSystem class members. More...
 
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. More...
 
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. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space. More...
 
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
virtual SOLVER_UTILS_EXPORT void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure (void)
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 

Protected Attributes

SolverUtils::RiemannSolverSharedPtr m_riemannSolver
 
Array< OneD, Array< OneD, NekDouble > > m_velocity
 
Array< OneD, NekDoublem_traceVn
 
- Protected Attributes inherited from Nektar::SolverUtils::AdvectionSystem
SolverUtils::AdvectionSharedPtr m_advObject
 Advection term. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
int m_infosteps
 Number of time steps between outputting status information. More...
 
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
int m_nanSteps
 
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
 
NekDouble m_epsilon
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::ofstream m_errFile
 
std::vector< int > m_intVariables
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
NekDouble m_filterTimeWarning
 Number of time steps between outputting status information. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh. More...
 
std::string m_sessionName
 Name of the session. More...
 
NekDouble m_time
 Current time of simulation. More...
 
int m_initialStep
 Number of the step where the simulation should begin. More...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_spacedim
 Spatial dimension (>= expansion dim). More...
 
int m_expdim
 Expansion dimension. More...
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used. More...
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used. More...
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used. More...
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform. More...
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations. More...
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation. More...
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction. More...
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity. More...
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields. More...
 
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error. More...
 
enum HomogeneousType m_HomogeneousType
 
NekDouble m_LhomX
 physical length in X direction (if homogeneous) More...
 
NekDouble m_LhomY
 physical length in Y direction (if homogeneous) More...
 
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous) More...
 
int m_npointsX
 number of points in X direction (if homogeneous) More...
 
int m_npointsY
 number of points in Y direction (if homogeneous) More...
 
int m_npointsZ
 number of points in Z direction (if homogeneous) More...
 
int m_HomoDirec
 number of homogenous directions More...
 

Private Member Functions

virtual NekDouble v_GetTimeStep (const Array< OneD, int > ExpOrder, const Array< OneD, NekDouble > CFL, NekDouble timeCFL)
 
virtual NekDouble v_GetTimeStep (int ExpOrder, NekDouble CFL, NekDouble TimeStability)
 
Array< OneD, NekDoubleGetStdVelocity (const Array< OneD, Array< OneD, NekDouble > > inarray)
 

Friends

class MemoryManager< CFLtester >
 

Additional Inherited Members

- Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D, eHomogeneous2D, eHomogeneous3D, eNotHomogeneous }
 Parameter for homogeneous expansions. More...
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 

Detailed Description

Definition at line 81 of file CFLtester.h.

Constructor & Destructor Documentation

◆ ~CFLtester()

Nektar::CFLtester::~CFLtester ( )
virtual

Definition at line 134 of file CFLtester.cpp.

135  {
136  }

◆ CFLtester()

Nektar::CFLtester::CFLtester ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
protected

Definition at line 47 of file CFLtester.cpp.

50  : UnsteadySystem(pSession, pGraph),
51  AdvectionSystem(pSession, pGraph)
52  {
53  }
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.
SOLVER_UTILS_EXPORT AdvectionSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)

Member Function Documentation

◆ create()

static EquationSystemSharedPtr Nektar::CFLtester::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
inlinestatic

Definition at line 86 of file CFLtester.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

89  {
91  ::AllocateSharedPtr(pSession, pGraph);
92  p->InitObject();
93  return p;
94  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.

◆ DoOdeProjection()

void Nektar::CFLtester::DoOdeProjection ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Definition at line 154 of file CFLtester.cpp.

References ASSERTL0, Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::SolverUtils::EquationSystem::GetNcoeffs(), Nektar::SolverUtils::EquationSystem::GetNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::SolverUtils::EquationSystem::SetBoundaryConditions(), and Vmath::Vcopy().

Referenced by v_InitObject().

158  {
159  int nvariables = inarray.num_elements();
160  SetBoundaryConditions(time);
161 
162  switch(m_projectionType)
163  {
165  {
166  // Just copy over array
167  int npoints = GetNpoints();
168 
169  for(int j = 0; j < nvariables; ++j)
170  {
171  Vmath::Vcopy(npoints,inarray[j],1,outarray[j],1);
172  }
173  }
174  break;
177  {
178  Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs());
179 
180  for(int j = 0; j < nvariables; ++j)
181  {
182  m_fields[j]->FwdTrans(inarray[j],coeffs);
183  m_fields[j]->BwdTrans_IterPerExp(coeffs,outarray[j]);
184  }
185  break;
186  }
187  default:
188  ASSERTL0(false,"Unknown projection scheme");
189  break;
190  }
191  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetNcoeffs()
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

◆ DoOdeRhs()

void Nektar::CFLtester::DoOdeRhs ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Definition at line 138 of file CFLtester.cpp.

References Nektar::SolverUtils::EquationSystem::GetNpoints(), Nektar::SolverUtils::AdvectionSystem::m_advObject, Nektar::SolverUtils::EquationSystem::m_fields, m_velocity, and Vmath::Neg().

Referenced by v_InitObject().

142  {
143  int nvariables = inarray.num_elements();
144  int npoints = GetNpoints();
145 
146  m_advObject->Advect(nvariables, m_fields, m_velocity, inarray,
147  outarray, time);
148  for(int j = 0; j < nvariables; ++j)
149  {
150  Vmath::Neg(npoints,outarray[j],1);
151  }
152  }
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:399
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
Array< OneD, Array< OneD, NekDouble > > m_velocity
Definition: CFLtester.h:102

◆ GetFluxVector()

void Nektar::CFLtester::GetFluxVector ( const Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  flux 
)
protected

Definition at line 221 of file CFLtester.cpp.

References ASSERTL1, m_velocity, and Vmath::Vmul().

Referenced by v_InitObject().

224  {
225  ASSERTL1(flux[0].num_elements() == m_velocity.num_elements(),
226  "Dimension of flux array and velocity array do not match");
227 
228  int nq = physfield[0].num_elements();
229 
230  for (int i = 0; i < flux.num_elements(); ++i)
231  {
232  for (int j = 0; j < flux[0].num_elements(); ++j)
233  {
234  Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1,
235  flux[i][j], 1);
236  }
237  }
238  }
Array< OneD, Array< OneD, NekDouble > > m_velocity
Definition: CFLtester.h:102
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:186

◆ GetNormalVelocity()

Array< OneD, NekDouble > & Nektar::CFLtester::GetNormalVelocity ( )
protected

Get the normal velocity.

Definition at line 196 of file CFLtester.cpp.

References Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_traceNormals, m_traceVn, m_velocity, Vmath::Vvtvp(), and Vmath::Zero().

Referenced by v_InitObject().

197  {
198  // Number of trace (interface) points
199  int nTracePts = GetTraceNpoints();
200 
201  // Auxiliary variable to compute the normal velocity
202  Array<OneD, NekDouble> tmp(nTracePts);
203 
204  // Reset the normal velocity
205  Vmath::Zero(nTracePts, m_traceVn, 1);
206 
207  for (int i = 0; i < m_velocity.num_elements(); ++i)
208  {
209  m_fields[0]->ExtractTracePhys(m_velocity[i], tmp);
210 
211  Vmath::Vvtvp(nTracePts,
212  m_traceNormals[i], 1,
213  tmp, 1,
214  m_traceVn, 1,
215  m_traceVn, 1);
216  }
217 
218  return m_traceVn;
219  }
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:445
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetTraceNpoints()
Array< OneD, Array< OneD, NekDouble > > m_velocity
Definition: CFLtester.h:102
Array< OneD, NekDouble > m_traceVn
Definition: CFLtester.h:103
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376

◆ GetStdVelocity()

Array< OneD, NekDouble > Nektar::CFLtester::GetStdVelocity ( const Array< OneD, Array< OneD, NekDouble > >  inarray)
private

Definition at line 347 of file CFLtester.cpp.

References ASSERTL0, Nektar::SpatialDomains::eDeformed, Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::SolverUtils::EquationSystem::m_expdim, and Nektar::SolverUtils::EquationSystem::m_fields.

Referenced by v_GetTimeStep().

349  {
350  // Checking if the problem is 2D
351  ASSERTL0(m_expdim >= 2, "Method not implemented for 1D");
352 
353  int nTotQuadPoints = GetTotPoints();
354  int n_element = m_fields[0]->GetExpSize();
355  int nvel = inarray.num_elements();
356 
357  NekDouble pntVelocity;
358 
359  // Getting the standard velocity vector on the 2D normal space
360  Array<OneD, Array<OneD, NekDouble> > stdVelocity(nvel);
361  Array<OneD, NekDouble> stdV(n_element, 0.0);
362 
363  for (int i = 0; i < nvel; ++i)
364  {
365  stdVelocity[i] = Array<OneD, NekDouble>(nTotQuadPoints);
366  }
367 
368  if (nvel == 2)
369  {
370  for (int el = 0; el < n_element; ++el)
371  {
373  m_fields[0]->GetExp(el)->as<LocalRegions::Expansion2D>();
374  int n_points = el2D->GetTotPoints();
375 
376  Array<OneD, const NekDouble> jac =
377  el2D->GetGeom2D()->GetMetricInfo()->GetJac();
378  Array<TwoD, const NekDouble> gmat =
379  el2D->GetGeom2D()->GetMetricInfo()->GetDerivFactors();
380 
381  if (el2D->GetGeom2D()->GetMetricInfo()->GetGtype()
383  {
384  for (int i = 0; i < n_points; i++)
385  {
386  stdVelocity[0][i] = gmat[0][i]*inarray[0][i]
387  + gmat[2][i]*inarray[1][i];
388 
389  stdVelocity[1][i] = gmat[1][i]*inarray[0][i]
390  + gmat[3][i]*inarray[1][i];
391  }
392  }
393  else
394  {
395  for (int i = 0; i < n_points; i++)
396  {
397  stdVelocity[0][i] = gmat[0][0]*inarray[0][i]
398  + gmat[2][0]*inarray[1][i];
399 
400  stdVelocity[1][i] = gmat[1][0]*inarray[0][i]
401  + gmat[3][0]*inarray[1][i];
402  }
403  }
404 
405 
406  for (int i = 0; i < n_points; i++)
407  {
408  pntVelocity = sqrt(stdVelocity[0][i]*stdVelocity[0][i]
409  + stdVelocity[1][i]*stdVelocity[1][i]);
410 
411  if (pntVelocity>stdV[el])
412  {
413  stdV[el] = pntVelocity;
414  }
415  }
416  }
417  }
418  else
419  {
420  for (int el = 0; el < n_element; ++el)
421  {
423  m_fields[0]->GetExp(el)->as<LocalRegions::Expansion3D>();
424 
425  int n_points = el3D->GetTotPoints();
426 
427  Array<OneD, const NekDouble> jac =
428  el3D->GetGeom3D()->GetMetricInfo()->GetJac();
429  Array<TwoD, const NekDouble> gmat =
430  el3D->GetGeom3D()->GetMetricInfo()->GetDerivFactors();
431 
432  if (el3D->GetGeom3D()->GetMetricInfo()->GetGtype()
434  {
435  for (int i = 0; i < n_points; i++)
436  {
437  stdVelocity[0][i] = gmat[0][i]*inarray[0][i]
438  + gmat[3][i]*inarray[1][i]
439  + gmat[6][i]*inarray[2][i];
440 
441  stdVelocity[1][i] = gmat[1][i]*inarray[0][i]
442  + gmat[4][i]*inarray[1][i]
443  + gmat[7][i]*inarray[2][i];
444 
445  stdVelocity[2][i] = gmat[2][i]*inarray[0][i]
446  + gmat[5][i]*inarray[1][i]
447  + gmat[8][i]*inarray[2][i];
448  }
449  }
450  else
451  {
452  Array<OneD, const NekDouble> jac =
453  el3D->GetGeom3D()->GetMetricInfo()->GetJac();
454  Array<TwoD, const NekDouble> gmat =
455  el3D->GetGeom3D()->GetMetricInfo()->GetDerivFactors();
456 
457  for (int i = 0; i < n_points; i++)
458  {
459  stdVelocity[0][i] = gmat[0][0]*inarray[0][i]
460  + gmat[3][0]*inarray[1][i]
461  + gmat[6][0]*inarray[2][i];
462 
463  stdVelocity[1][i] = gmat[1][0]*inarray[0][i]
464  + gmat[4][0]*inarray[1][i]
465  + gmat[7][0]*inarray[2][i];
466 
467  stdVelocity[2][i] = gmat[2][0]*inarray[0][i]
468  + gmat[5][0]*inarray[1][i]
469  + gmat[8][0]*inarray[2][i];
470  }
471  }
472 
473  for (int i = 0; i < n_points; i++)
474  {
475  pntVelocity = sqrt(stdVelocity[0][i]*stdVelocity[0][i]
476  + stdVelocity[1][i]*stdVelocity[1][i]
477  + stdVelocity[2][i]*stdVelocity[2][i]);
478 
479  if (pntVelocity > stdV[el])
480  {
481  stdV[el] = pntVelocity;
482  }
483  }
484  }
485  }
486 
487  return stdV;
488  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
int m_expdim
Expansion dimension.
SOLVER_UTILS_EXPORT int GetTotPoints()
std::shared_ptr< Expansion3D > Expansion3DSharedPtr
Definition: Expansion2D.h:49
double NekDouble
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
std::shared_ptr< Expansion2D > Expansion2DSharedPtr
Definition: Expansion1D.h:47
Geometry is curved or has non-constant factors.

◆ v_GenerateSummary()

void Nektar::CFLtester::v_GenerateSummary ( SummaryList s)
protectedvirtual

Print a summary of time stepping parameters.

Prints a summary with some information regards the time-stepping.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 240 of file CFLtester.cpp.

References Nektar::SolverUtils::UnsteadySystem::v_GenerateSummary().

241  {
243  }
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.

◆ v_GetTimeStep() [1/2]

NekDouble Nektar::CFLtester::v_GetTimeStep ( const Array< OneD, int >  ExpOrder,
const Array< OneD, NekDouble CFL,
NekDouble  timeCFL 
)
privatevirtual

Definition at line 247 of file CFLtester.cpp.

References GetStdVelocity(), Nektar::SolverUtils::EquationSystem::m_fields, m_velocity, and Vmath::Vmin().

251  {
252 
253  int n_element = m_fields[0]->GetExpSize();
254 
255  //const NekDouble minLengthStdTri = 1.414213;
256  //const NekDouble minLengthStdQuad = 2.0;
257  //const NekDouble cLambda = 0.2; // Spencer book pag. 317
258 
259  Array<OneD, NekDouble> tstep (n_element, 0.0);
260  Array<OneD, NekDouble> stdVelocity(n_element, 0.0);
261  stdVelocity = GetStdVelocity(m_velocity);
262 
263  for (int el = 0; el < n_element; ++el)
264  {
265  int npoints = m_fields[0]->GetExp(el)->GetTotPoints();
266  Array<OneD, NekDouble> one2D(npoints, 1.0);
267  //NekDouble Area = m_fields[0]->GetExp(el)->Integral(one2D);
268  if(std::dynamic_pointer_cast<LocalRegions::TriExp>(m_fields[0]->GetExp(el)))
269  {
270  //tstep[el] = timeCFL/(stdVelocity[el]*cLambda*(ExpOrder[el]-1)*(ExpOrder[el]-1));
271  //tstep[el] = timeCFL*minLengthStdTri/(stdVelocity[el]*cLambda*(ExpOrder[el]-1)*(ExpOrder[el]-1));
272  tstep[el] = CFL[el]/(stdVelocity[el]);
273  }
274  else if(std::dynamic_pointer_cast<LocalRegions::QuadExp>(m_fields[0]->GetExp(el)))
275  {
276  //tstep[el] = timeCFL/(stdVelocity[el]*cLambda*(ExpOrder[el]-1)*(ExpOrder[el]-1));
277  //tstep[el] = timeCFL*minLengthStdQuad/(stdVelocity[el]*cLambda*(ExpOrder[el]-1)*(ExpOrder[el]-1));
278  tstep[el] = CFL[el]/(stdVelocity[el]);
279  }
280  }
281 
282  NekDouble TimeStep = Vmath::Vmin(n_element, tstep, 1);
283 
284  return TimeStep;
285  }
T Vmin(int n, const T *x, const int incx)
Return the minimum element in x - called vmin to avoid conflict with min.
Definition: Vmath.cpp:874
Array< OneD, NekDouble > GetStdVelocity(const Array< OneD, Array< OneD, NekDouble > > inarray)
Definition: CFLtester.cpp:347
double NekDouble
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
Array< OneD, Array< OneD, NekDouble > > m_velocity
Definition: CFLtester.h:102

◆ v_GetTimeStep() [2/2]

NekDouble Nektar::CFLtester::v_GetTimeStep ( int  ExpOrder,
NekDouble  CFL,
NekDouble  TimeStability 
)
privatevirtual

Definition at line 289 of file CFLtester.cpp.

References ASSERTL0, Nektar::EigenvaluesAnaMeshesAB2, Nektar::EigenvaluesAnaMeshesRK2, Nektar::EigenvaluesAnaMeshesRK4, Nektar::SolverUtils::EquationSystem::m_fields, and class_topology::P.

293  {
294  //================================================================
295  // This function has been created just to test specific problems, hence is not general
296  // and it has been implemented in a rude fashion, as the full CFLtester class.
297  // For real CFL calculations refer to the general implementation above. (A.Bolis)
298  //================================================================
299 
300  NekDouble TimeStep;
301  NekDouble SpatialStability;
302  int n_elements = m_fields[0]->GetExpSize();
303 
304  //solve ambiguity in windows
305  NekDouble n_elem = n_elements;
306  NekDouble DH = sqrt(n_elem);
307 
308  int H = (int)DH;
309  int P = ExpOrder-1;
310 
311  //================================================================
312  // Regular meshes
313 
314  //SpatialStability = EigenvaluesRegMeshes[H-1][P-1];
315 
316  //================================================================
317  // Anisotropic meshes
318 
319  if (TimeStability == 1.0)
320  {
321  SpatialStability = EigenvaluesAnaMeshesAB2[H/2][P-1];
322  }
323  else if (TimeStability == 2.0)
324  {
325  SpatialStability = EigenvaluesAnaMeshesRK2[H/2][P-1];
326  }
327  else if (TimeStability == 2.784)
328  {
329  SpatialStability = EigenvaluesAnaMeshesRK4[H/2][P-1];
330  }
331  else
332  {
333  ASSERTL0(false,"Dominant eigenvalues database not present for this time-stepping scheme")
334  }
335 
336  //================================================================
337 
338  TimeStep = (TimeStability/SpatialStability)*CFL;
339 
340  //================================================================
341 
342  return TimeStep;
343  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
static NekDouble EigenvaluesAnaMeshesRK2[6][14]
Definition: CFLtester.h:65
static NekDouble EigenvaluesAnaMeshesRK4[6][14]
Definition: CFLtester.h:73
double NekDouble
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
static NekDouble EigenvaluesAnaMeshesAB2[6][14]
Definition: CFLtester.h:57

◆ v_InitObject()

void Nektar::CFLtester::v_InitObject ( )
protectedvirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::AdvectionSystem.

Definition at line 55 of file CFLtester.cpp.

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), DoOdeProjection(), DoOdeRhs(), Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::SolverUtils::GetAdvectionFactory(), GetFluxVector(), Nektar::SolverUtils::EquationSystem::GetFunction(), GetNormalVelocity(), Nektar::SolverUtils::GetRiemannSolverFactory(), Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::SolverUtils::AdvectionSystem::m_advObject, Nektar::SolverUtils::UnsteadySystem::m_explicitAdvection, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::UnsteadySystem::m_homoInitialFwd, Nektar::SolverUtils::UnsteadySystem::m_ode, Nektar::SolverUtils::EquationSystem::m_projectionType, m_riemannSolver, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_spacedim, m_traceVn, m_velocity, and Nektar::SolverUtils::AdvectionSystem::v_InitObject().

56  {
58 
59  m_velocity = Array<OneD, Array<OneD, NekDouble> >(m_spacedim);
60  std::vector<std::string> vel;
61  vel.push_back("Vx");
62  vel.push_back("Vy");
63  vel.push_back("Vz");
64  vel.resize(m_spacedim);
65 
66  GetFunction( "AdvectionVelocity")->Evaluate(vel, m_velocity);
67 
68  // Type of advection class to be used
69  switch(m_projectionType)
70  {
71  // Discontinuous field
73  {
74  // Do not forwards transform initial condition
75  m_homoInitialFwd = false;
76 
77  // Define the normal velocity fields
78  if (m_fields[0]->GetTrace())
79  {
80  m_traceVn = Array<OneD, NekDouble>(GetTraceNpoints());
81  }
82 
83  string advName;
84  string riemName;
85  m_session->LoadSolverInfo(
86  "AdvectionType", advName, "WeakDG");
88  GetAdvectionFactory().CreateInstance(advName, advName);
89  m_advObject->SetFluxVector(
91  m_session->LoadSolverInfo(
92  "UpwindType", riemName, "Upwind");
95  riemName, m_session);
96  m_riemannSolver->SetScalar(
97  "Vn", &CFLtester::GetNormalVelocity, this);
98 
99  m_advObject->SetRiemannSolver(m_riemannSolver);
100  m_advObject->InitObject(m_session, m_fields);
101  break;
102  }
103  // Continuous field
106  {
107  string advName;
108  m_session->LoadSolverInfo(
109  "AdvectionType", advName, "NonConservative");
111  GetAdvectionFactory().CreateInstance(advName, advName);
112  m_advObject->SetFluxVector(
113  &CFLtester::GetFluxVector, this);
114  break;
115  }
116  default:
117  {
118  ASSERTL0(false, "Unsupported projection type.");
119  break;
120  }
121  }
122 
124  {
127  }
128  else
129  {
130  ASSERTL0(false, "Implicit unsteady Advection not set up.");
131  }
132  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Definition: CFLtester.cpp:221
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Definition: CFLtester.cpp:138
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
RiemannSolverFactory & GetRiemannSolverFactory()
int m_spacedim
Spatial dimension (>= expansion dim).
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Definition: CFLtester.cpp:154
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT int GetTraceNpoints()
Array< OneD, Array< OneD, NekDouble > > m_velocity
Definition: CFLtester.h:102
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
Definition: CFLtester.h:101
Array< OneD, NekDouble > m_traceVn
Definition: CFLtester.h:103
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Init object for UnsteadySystem class.
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity.
Definition: CFLtester.cpp:196

Friends And Related Function Documentation

◆ MemoryManager< CFLtester >

friend class MemoryManager< CFLtester >
friend

Definition at line 84 of file CFLtester.h.

Member Data Documentation

◆ className

string Nektar::CFLtester::className = GetEquationSystemFactory().RegisterCreatorFunction("CFLtester", CFLtester::create, "Testing CFL restriction")
static

Definition at line 96 of file CFLtester.h.

◆ m_riemannSolver

SolverUtils::RiemannSolverSharedPtr Nektar::CFLtester::m_riemannSolver
protected

Definition at line 101 of file CFLtester.h.

Referenced by v_InitObject().

◆ m_traceVn

Array<OneD, NekDouble> Nektar::CFLtester::m_traceVn
protected

Definition at line 103 of file CFLtester.h.

Referenced by GetNormalVelocity(), and v_InitObject().

◆ m_velocity

Array<OneD, Array<OneD, NekDouble> > Nektar::CFLtester::m_velocity
protected

Definition at line 102 of file CFLtester.h.

Referenced by DoOdeRhs(), GetFluxVector(), GetNormalVelocity(), v_GetTimeStep(), and v_InitObject().