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 ()
 
virtual SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareField=true) override
 Init object for UnsteadySystem class. More...
 
SOLVER_UTILS_EXPORT AdvectionSharedPtr GetAdvObject ()
 Returns the advection object held by this instance. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleGetElmtCFLVals (const bool FlagAcousticCFL=true)
 
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...
 
SOLVER_UTILS_EXPORT void SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeSharedPtrGetTimeIntegrationScheme ()
 Returns the time integration scheme. More...
 
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeOperatorsGetTimeIntegrationSchemeOperators ()
 Returns the time integration scheme operators. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void InitObject (bool DeclareField=true)
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise (bool dumpInitialConditions=true)
 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 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 NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. 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, const Array< OneD, const NekDouble > &input)
 
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 int GetInfoSteps ()
 
SOLVER_UTILS_EXPORT void SetInfoSteps (int num)
 
SOLVER_UTILS_EXPORT void SetIterationNumberPIT (int num)
 
SOLVER_UTILS_EXPORT void SetWindowNumberPIT (int num)
 
SOLVER_UTILS_EXPORT Array< OneD, const Array< OneD, NekDouble > > GetTraceNormals ()
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetTimeStep (const NekDouble timestep)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
SOLVER_UTILS_EXPORT bool NegatedOp ()
 Identify if operator is negated in DoSolve. More...
 
SOLVER_UTILS_EXPORT bool ParallelInTime ()
 Check if solver use Parallel-in-Time. More...
 

Static Public Member Functions

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

Static Public Attributes

static std::string className
 
- Static Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
static std::string cmdSetStartTime
 
static std::string cmdSetStartChkNum
 

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 (bool DeclareFields=true) override
 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) override
 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) override
 
virtual SOLVER_UTILS_EXPORT Array< OneD, NekDoublev_GetMaxStdVelocity (const NekDouble SpeedSoundFactor=1.0)
 
- 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...
 
virtual SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareField=true) override
 Init object for UnsteadySystem class. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve () override
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintStatusInformation (const int step, const NekDouble cpuTime)
 Print Status Information. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintSummaryStatistics (const NekDouble intTime)
 Print Summary Statistics. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true) override
 Sets up initial conditions. More...
 
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &s) override
 Print a summary of time stepping parameters. 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_PostIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_RequireFwdTrans ()
 
virtual SOLVER_UTILS_EXPORT void v_SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
virtual SOLVER_UTILS_EXPORT bool v_UpdateTimeStepCheck ()
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
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...
 
SOLVER_UTILS_EXPORT void DoDummyProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Perform dummy projection. 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 void v_InitObject (bool DeclareFeld=true)
 Initialisation object for EquationSystem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true)
 Virtual function for initialisation implementation. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Virtual function for solve implementation. 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_GenerateSummary (SummaryList &l)
 Virtual function for generating summary information. 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 bool v_NegatedOp (void)
 Virtual function to identify if operator is negated in DoSolve. More...
 
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
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::vector< int > m_intVariables
 
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 Maximun cfl in cfl growth. More...
 
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
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...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
std::ofstream m_errFile
 
NekDouble m_epsilon
 Diffusion coefficient. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
bool m_verbose
 
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...
 
NekDouble m_lastCheckTime
 
NekDouble m_TimeIncrementFactor
 
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_infosteps
 Number of time steps between outputting status information. More...
 
int m_iterPIT = 0
 Number of parallel-in-time time iteration. More...
 
int m_windowPIT = 0
 Index of windows for parallel-in-time time iteration. 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...
 
Array< OneD, NekDoublem_movingFrameVelsxyz
 Moving frame of reference velocities. More...
 
Array< OneD, NekDoublem_movingFrameTheta
 Moving frame of reference angles with respect to the. More...
 
boost::numeric::ublas::matrix< NekDoublem_movingFrameProjMat
 Projection matrix for transformation between inertial and moving. 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) override
 
virtual NekDouble v_GetTimeStep (int ExpOrder, NekDouble CFL, NekDouble TimeStability) override
 
Array< OneD, NekDoubleGetStdVelocity (const Array< OneD, Array< OneD, NekDouble > > inarray)
 

Friends

class MemoryManager< CFLtester >
 

Additional Inherited Members

- 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 []
 
static std::string projectionTypeLookupIds []
 

Detailed Description

Definition at line 137 of file CFLtester.h.

Constructor & Destructor Documentation

◆ ~CFLtester()

Nektar::CFLtester::~CFLtester ( )
virtual

Definition at line 130 of file CFLtester.cpp.

131{
132}

◆ CFLtester()

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

Definition at line 49 of file CFLtester.cpp.

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

Member Function Documentation

◆ create()

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

Definition at line 142 of file CFLtester.h.

145 {
148 p->InitObject();
149 return p;
150 }
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.

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

◆ 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 149 of file CFLtester.cpp.

152{
153 int nvariables = inarray.size();
155
156 switch (m_projectionType)
157 {
159 {
160 // Just copy over array
161 if (inarray != outarray)
162 {
163 int npoints = GetNpoints();
164
165 for (int j = 0; j < nvariables; ++j)
166 {
167 Vmath::Vcopy(npoints, inarray[j], 1, outarray[j], 1);
168 }
169 }
170 }
171 break;
174 {
175 Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs());
176
177 for (int j = 0; j < nvariables; ++j)
178 {
179 m_fields[j]->FwdTrans(inarray[j], coeffs);
180 m_fields[j]->BwdTrans_IterPerExp(coeffs, outarray[j]);
181 }
182 break;
183 }
184 default:
185 ASSERTL0(false, "Unknown projection scheme");
186 break;
187 }
188}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetNpoints()
SOLVER_UTILS_EXPORT int GetNcoeffs()
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.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191

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().

◆ 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 134 of file CFLtester.cpp.

137{
138 int nvariables = inarray.size();
139 int npoints = GetNpoints();
140
141 m_advObject->Advect(nvariables, m_fields, m_velocity, inarray, outarray,
142 time);
143 for (int j = 0; j < nvariables; ++j)
144 {
145 Vmath::Neg(npoints, outarray[j], 1);
146 }
147}
Array< OneD, Array< OneD, NekDouble > > m_velocity
Definition: CFLtester.h:158
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:513

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

Referenced by v_InitObject().

◆ GetFluxVector()

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

Definition at line 215 of file CFLtester.cpp.

218{
219 ASSERTL1(flux[0].size() == m_velocity.size(),
220 "Dimension of flux array and velocity array do not match");
221
222 int nq = physfield[0].size();
223
224 for (int i = 0; i < flux.size(); ++i)
225 {
226 for (int j = 0; j < flux[0].size(); ++j)
227 {
228 Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1, flux[i][j], 1);
229 }
230 }
231}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
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:207

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

Referenced by v_InitObject().

◆ GetNormalVelocity()

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

Get the normal velocity.

Definition at line 193 of file CFLtester.cpp.

194{
195 // Number of trace (interface) points
196 int nTracePts = GetTraceNpoints();
197
198 // Auxiliary variable to compute the normal velocity
199 Array<OneD, NekDouble> tmp(nTracePts);
200
201 // Reset the normal velocity
202 Vmath::Zero(nTracePts, m_traceVn, 1);
203
204 for (int i = 0; i < m_velocity.size(); ++i)
205 {
206 m_fields[0]->ExtractTracePhys(m_velocity[i], tmp);
207
208 Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, tmp, 1, m_traceVn, 1,
209 m_traceVn, 1);
210 }
211
212 return m_traceVn;
213}
Array< OneD, NekDouble > m_traceVn
Definition: CFLtester.h:159
SOLVER_UTILS_EXPORT int GetTraceNpoints()
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
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:569
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487

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().

◆ GetStdVelocity()

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

Definition at line 326 of file CFLtester.cpp.

328{
329 // Checking if the problem is 2D
330 ASSERTL0(m_expdim >= 2, "Method not implemented for 1D");
331
332 int nTotQuadPoints = GetTotPoints();
333 int n_element = m_fields[0]->GetExpSize();
334 int nvel = inarray.size();
335
336 NekDouble pntVelocity;
337
338 // Getting the standard velocity vector on the 2D normal space
339 Array<OneD, Array<OneD, NekDouble>> stdVelocity(nvel);
340 Array<OneD, NekDouble> stdV(n_element, 0.0);
341
342 for (int i = 0; i < nvel; ++i)
343 {
344 stdVelocity[i] = Array<OneD, NekDouble>(nTotQuadPoints);
345 }
346
347 if (nvel == 2)
348 {
349 for (int el = 0; el < n_element; ++el)
350 {
352 m_fields[0]->GetExp(el)->as<LocalRegions::Expansion2D>();
353 int n_points = el2D->GetTotPoints();
354
355 Array<OneD, const NekDouble> jac =
356 el2D->GetGeom2D()->GetMetricInfo()->GetJac();
357 Array<TwoD, const NekDouble> gmat =
358 el2D->GetGeom2D()->GetMetricInfo()->GetDerivFactors();
359
360 if (el2D->GetGeom2D()->GetMetricInfo()->GetGtype() ==
362 {
363 for (int i = 0; i < n_points; i++)
364 {
365 stdVelocity[0][i] =
366 gmat[0][i] * inarray[0][i] + gmat[2][i] * inarray[1][i];
367
368 stdVelocity[1][i] =
369 gmat[1][i] * inarray[0][i] + gmat[3][i] * inarray[1][i];
370 }
371 }
372 else
373 {
374 for (int i = 0; i < n_points; i++)
375 {
376 stdVelocity[0][i] =
377 gmat[0][0] * inarray[0][i] + gmat[2][0] * inarray[1][i];
378
379 stdVelocity[1][i] =
380 gmat[1][0] * inarray[0][i] + gmat[3][0] * inarray[1][i];
381 }
382 }
383
384 for (int i = 0; i < n_points; i++)
385 {
386 pntVelocity = sqrt(stdVelocity[0][i] * stdVelocity[0][i] +
387 stdVelocity[1][i] * stdVelocity[1][i]);
388
389 if (pntVelocity > stdV[el])
390 {
391 stdV[el] = pntVelocity;
392 }
393 }
394 }
395 }
396 else
397 {
398 for (int el = 0; el < n_element; ++el)
399 {
401 m_fields[0]->GetExp(el)->as<LocalRegions::Expansion3D>();
402
403 int n_points = el3D->GetTotPoints();
404
405 Array<OneD, const NekDouble> jac =
406 el3D->GetGeom3D()->GetMetricInfo()->GetJac();
407 Array<TwoD, const NekDouble> gmat =
408 el3D->GetGeom3D()->GetMetricInfo()->GetDerivFactors();
409
410 if (el3D->GetGeom3D()->GetMetricInfo()->GetGtype() ==
412 {
413 for (int i = 0; i < n_points; i++)
414 {
415 stdVelocity[0][i] = gmat[0][i] * inarray[0][i] +
416 gmat[3][i] * inarray[1][i] +
417 gmat[6][i] * inarray[2][i];
418
419 stdVelocity[1][i] = gmat[1][i] * inarray[0][i] +
420 gmat[4][i] * inarray[1][i] +
421 gmat[7][i] * inarray[2][i];
422
423 stdVelocity[2][i] = gmat[2][i] * inarray[0][i] +
424 gmat[5][i] * inarray[1][i] +
425 gmat[8][i] * inarray[2][i];
426 }
427 }
428 else
429 {
430 Array<OneD, const NekDouble> jac =
431 el3D->GetGeom3D()->GetMetricInfo()->GetJac();
432 Array<TwoD, const NekDouble> gmat =
433 el3D->GetGeom3D()->GetMetricInfo()->GetDerivFactors();
434
435 for (int i = 0; i < n_points; i++)
436 {
437 stdVelocity[0][i] = gmat[0][0] * inarray[0][i] +
438 gmat[3][0] * inarray[1][i] +
439 gmat[6][0] * inarray[2][i];
440
441 stdVelocity[1][i] = gmat[1][0] * inarray[0][i] +
442 gmat[4][0] * inarray[1][i] +
443 gmat[7][0] * inarray[2][i];
444
445 stdVelocity[2][i] = gmat[2][0] * inarray[0][i] +
446 gmat[5][0] * inarray[1][i] +
447 gmat[8][0] * inarray[2][i];
448 }
449 }
450
451 for (int i = 0; i < n_points; i++)
452 {
453 pntVelocity = sqrt(stdVelocity[0][i] * stdVelocity[0][i] +
454 stdVelocity[1][i] * stdVelocity[1][i] +
455 stdVelocity[2][i] * stdVelocity[2][i]);
456
457 if (pntVelocity > stdV[el])
458 {
459 stdV[el] = pntVelocity;
460 }
461 }
462 }
463 }
464
465 return stdV;
466}
int m_expdim
Expansion dimension.
SOLVER_UTILS_EXPORT int GetTotPoints()
std::shared_ptr< Expansion2D > Expansion2DSharedPtr
Definition: Expansion1D.h:48
std::shared_ptr< Expansion3D > Expansion3DSharedPtr
Definition: Expansion2D.h:51
@ eDeformed
Geometry is curved or has non-constant factors.
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

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

Referenced by v_GetTimeStep().

◆ v_GenerateSummary()

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

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 233 of file CFLtester.cpp.

234{
236}
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.

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

◆ v_GetTimeStep() [1/2]

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

Definition at line 238 of file CFLtester.cpp.

241{
242
243 int n_element = m_fields[0]->GetExpSize();
244
245 Array<OneD, NekDouble> tstep(n_element, 0.0);
246 Array<OneD, NekDouble> stdVelocity(n_element, 0.0);
247 stdVelocity = GetStdVelocity(m_velocity);
248
249 for (int el = 0; el < n_element; ++el)
250 {
251 int npoints = m_fields[0]->GetExp(el)->GetTotPoints();
252 Array<OneD, NekDouble> one2D(npoints, 1.0);
253 if (std::dynamic_pointer_cast<LocalRegions::TriExp>(
254 m_fields[0]->GetExp(el)))
255 {
256 tstep[el] = CFL[el] / (stdVelocity[el]);
257 }
258 else if (std::dynamic_pointer_cast<LocalRegions::QuadExp>(
259 m_fields[0]->GetExp(el)))
260 {
261 tstep[el] = CFL[el] / (stdVelocity[el]);
262 }
263 }
264
265 NekDouble TimeStep = Vmath::Vmin(n_element, tstep, 1);
266
267 return TimeStep;
268}
Array< OneD, NekDouble > GetStdVelocity(const Array< OneD, Array< OneD, NekDouble > > inarray)
Definition: CFLtester.cpp:326
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:1045

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

◆ v_GetTimeStep() [2/2]

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

Definition at line 270 of file CFLtester.cpp.

272{
273 //================================================================
274 // This function has been created just to test specific problems, hence is
275 // not general and it has been implemented in a rude fashion, as the full
276 // CFLtester class. For real CFL calculations refer to the general
277 // implementation above. (A.Bolis)
278 //================================================================
279
280 NekDouble TimeStep;
281 NekDouble SpatialStability;
282 int n_elements = m_fields[0]->GetExpSize();
283
284 // solve ambiguity in windows
285 NekDouble n_elem = n_elements;
286 NekDouble DH = sqrt(n_elem);
287
288 int H = (int)DH;
289 int P = ExpOrder - 1;
290
291 //================================================================
292 // Regular meshes
293
294 // SpatialStability = EigenvaluesRegMeshes[H-1][P-1];
295
296 //================================================================
297 // Anisotropic meshes
298
299 if (TimeStability == 1.0)
300 {
301 SpatialStability = EigenvaluesAnaMeshesAB2[H / 2][P - 1];
302 }
303 else if (TimeStability == 2.0)
304 {
305 SpatialStability = EigenvaluesAnaMeshesRK2[H / 2][P - 1];
306 }
307 else if (TimeStability == 2.784)
308 {
309 SpatialStability = EigenvaluesAnaMeshesRK4[H / 2][P - 1];
310 }
311 else
312 {
313 ASSERTL0(false, "Dominant eigenvalues database not present for this "
314 "time-stepping scheme")
315 }
316
317 //================================================================
318
319 TimeStep = (TimeStability / SpatialStability) * CFL;
320
321 //================================================================
322
323 return TimeStep;
324}
@ P
Monomial polynomials .
Definition: BasisType.h:64
static NekDouble EigenvaluesAnaMeshesRK2[6][14]
Definition: CFLtester.h:97
static NekDouble EigenvaluesAnaMeshesRK4[6][14]
Definition: CFLtester.h:117
static NekDouble EigenvaluesAnaMeshesAB2[6][14]
Definition: CFLtester.h:77

References ASSERTL0, Nektar::EigenvaluesAnaMeshesAB2, Nektar::EigenvaluesAnaMeshesRK2, Nektar::EigenvaluesAnaMeshesRK4, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::LibUtilities::P, and tinysimd::sqrt().

◆ v_InitObject()

void Nektar::CFLtester::v_InitObject ( bool  DeclareField = true)
overrideprotectedvirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::AdvectionSystem.

Definition at line 55 of file CFLtester.cpp.

56{
57 AdvectionSystem::v_InitObject(DeclareFields);
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("AdvectionType", advName, "WeakDG");
87 advName, advName);
88 m_advObject->SetFluxVector(&CFLtester::GetFluxVector, this);
89 m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
92 riemName, m_session);
94 this);
95
96 m_advObject->SetRiemannSolver(m_riemannSolver);
97 m_advObject->InitObject(m_session, m_fields);
98 break;
99 }
100 // Continuous field
103 {
104 string advName;
105 m_session->LoadSolverInfo("AdvectionType", advName,
106 "NonConservative");
108 advName, advName);
109 m_advObject->SetFluxVector(&CFLtester::GetFluxVector, this);
110 break;
111 }
112 default:
113 {
114 ASSERTL0(false, "Unsupported projection type.");
115 break;
116 }
117 }
118
120 {
123 }
124 else
125 {
126 ASSERTL0(false, "Implicit unsteady Advection not set up.");
127 }
128}
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
Definition: CFLtester.h:157
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Definition: CFLtester.cpp:149
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity.
Definition: CFLtester.cpp:193
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Definition: CFLtester.cpp:215
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Definition: CFLtester.cpp:134
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
int m_spacedim
Spatial dimension (>= expansion dim).
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
RiemannSolverFactory & GetRiemannSolverFactory()

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().

Friends And Related Function Documentation

◆ MemoryManager< CFLtester >

friend class MemoryManager< CFLtester >
friend

Definition at line 117 of file CFLtester.h.

Member Data Documentation

◆ className

string Nektar::CFLtester::className
static
Initial value:
=
"CFLtester", CFLtester::create, "Testing CFL restriction")
static EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Definition: CFLtester.h:142
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
EquationSystemFactory & GetEquationSystemFactory()

Definition at line 152 of file CFLtester.h.

◆ m_riemannSolver

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

Definition at line 157 of file CFLtester.h.

Referenced by v_InitObject().

◆ m_traceVn

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

Definition at line 159 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 158 of file CFLtester.h.

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