Nektar++
Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | Private Member Functions | List of all members
Nektar::SolverUtils::UnsteadySystem Class Reference

Base class for unsteady solvers. More...

#include <UnsteadySystem.h>

Inheritance diagram for Nektar::SolverUtils::UnsteadySystem:
Inheritance graph
[legend]
Collaboration diagram for Nektar::SolverUtils::UnsteadySystem:
Collaboration graph
[legend]

Public Member Functions

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 >
boost::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 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 void EvaluateFunction (Array< OneD, Array< OneD, NekDouble > > &pArray, std::string pFunctionName, const NekDouble pTime=0.0, const int domain=0)
 Evaluates a function as specified in the session file. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, Array< OneD, NekDouble > > &pFields, const std::string &pName, const int domain=0)
 Populate given fields with the function from session. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const std::string &pName, const int domain=0)
 Populate given fields with the function from session. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
 
SOLVER_UTILS_EXPORT std::string DescribeFunction (std::string pFieldName, const std::string &pFunctionName, const int domain)
 Provide a description of a function for a given field name. More...
 
SOLVER_UTILS_EXPORT void InitialiseBaseFlow (Array< OneD, Array< OneD, NekDouble > > &base)
 Perform initialisation of the base flow. 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 WeakAdvectionGreensDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\nabla \phi \cdot F) $. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\phi, \nabla \cdot F) $. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, bool UseContCoeffs=false)
 Compute the inner product $ (\phi, V\cdot \nabla u) $. More...
 
f SOLVER_UTILS_EXPORT void AdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray)
 Compute the non-conservative advection. More...
 
SOLVER_UTILS_EXPORT void WeakDGAdvection (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false, int nvariables=0)
 Calculate the weak discontinuous Galerkin advection. More...
 
SOLVER_UTILS_EXPORT void WeakDGDiffusion (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false)
 Calculate weak DG Diffusion in the LDG form. 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 ScanForHistoryPoints ()
 Builds map of which element holds each history point. More...
 
SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out)
 Probe each history point and write to 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 GetNumElmVelocity ()
 
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 SetStepsToOne ()
 
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
 
SOLVER_UTILS_EXPORT void FwdTransFields ()
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &fluxX, Array< OneD, Array< OneD, NekDouble > > &fluxY)
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, const int j, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
 
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
 
SOLVER_UTILS_EXPORT void NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
 
SOLVER_UTILS_EXPORT void NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
 
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
 
SOLVER_UTILS_EXPORT int NoCaseStringCompare (const string &s1, const string &s2)
 Perform a case-insensitive string comparison. More...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to identify if operator is negated in DoSolve. More...
 

Public Attributes

NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 

Protected Member Functions

SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises UnsteadySystem class members. More...
 
virtual SOLVER_UTILS_EXPORT void v_InitObject ()
 Init object for UnsteadySystem class. 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_GenerateSummary (SummaryList &s)
 Print a summary of time stepping parameters. 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 void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
 
virtual SOLVER_UTILS_EXPORT void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
 
virtual SOLVER_UTILS_EXPORT void v_NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
 
virtual SOLVER_UTILS_EXPORT void v_NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
 
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_SteadyStateCheck (int step)
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time)
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises EquationSystem class members. More...
 
int nocase_cmp (const string &s1, const string &s2)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. 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)
 
SOLVER_UTILS_EXPORT void SetUpBaseFields (SpatialDomains::MeshGraphSharedPtr &mesh)
 
SOLVER_UTILS_EXPORT void ImportFldBase (std::string pInfile, SpatialDomains::MeshGraphSharedPtr pGraph)
 
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

int m_infosteps
 Number of time steps between outputting status information. More...
 
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...
 
std::vector< int > m_intVariables
 
std::vector< FilterSharedPtrm_filters
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
map< std::string, Array< OneD, Array< OneD, float > > > m_interpWeights
 Map of the interpolation weights for a specific filename. More...
 
map< std::string, Array< OneD, Array< OneD, unsigned int > > > m_interpInds
 Map of the interpolation indices for a specific filename. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_base
 Base fields. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_derivedfields
 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...
 
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_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, Array< OneD, Array< OneD, NekDouble > > > m_gradtan
 1 x nvariable x nq More...
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_tanbasis
 2 x m_spacedim x nq 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...
 
int m_NumMode
 Mode to use in case of single mode analysis. More...
 

Private Member Functions

void WeakPenaltyforScalar (const int var, const Array< OneD, const NekDouble > &physfield, Array< OneD, NekDouble > &penaltyflux, NekDouble time=0.0)
 
void WeakPenaltyforVector (const int var, const int dir, const Array< OneD, const NekDouble > &physfield, Array< OneD, NekDouble > &penaltyflux, NekDouble C11, NekDouble time=0.0)
 

Additional Inherited Members

- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D, eHomogeneous2D, eHomogeneous3D, eNotHomogeneous }
 Parameter for homogeneous expansions. More...
 

Detailed Description

Base class for unsteady solvers.

Provides the underlying timestepping framework for unsteady solvers including the general timestepping routines. This class is not intended to be directly instantiated, but rather is a base class on which to define unsteady solvers.

For details on implementing unsteady solvers see sectionADRSolverModuleImplementation here

Definition at line 48 of file UnsteadySystem.h.

Constructor & Destructor Documentation

Nektar::SolverUtils::UnsteadySystem::~UnsteadySystem ( )
virtual

Destructor.

Destructor for the class UnsteadyAdvection.

Definition at line 122 of file UnsteadySystem.cpp.

123  {
124  }
Nektar::SolverUtils::UnsteadySystem::UnsteadySystem ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Initialises UnsteadySystem class members.

Processes SolverInfo parameters from the session file and sets up timestepping-specific code.

Parameters
pSessionSession object to read parameters from.

Definition at line 65 of file UnsteadySystem.cpp.

67  : EquationSystem(pSession),
68  m_infosteps(10)
69 
70  {
71  }
int m_infosteps
Number of time steps between outputting status information.
SOLVER_UTILS_EXPORT EquationSystem(const LibUtilities::SessionReaderSharedPtr &pSession)
Initialises EquationSystem class members.

Member Function Documentation

void Nektar::SolverUtils::UnsteadySystem::CheckForRestartTime ( NekDouble time)
protected

Definition at line 683 of file UnsteadySystem.cpp.

References Nektar::LibUtilities::eFunctionTypeFile, Nektar::iterator, Nektar::SolverUtils::EquationSystem::m_fieldMetaDataMap, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_fld, Nektar::SolverUtils::EquationSystem::m_session, Nektar::LibUtilities::NullFieldMetaDataMap, and Nektar::LibUtilities::PortablePath().

Referenced by v_DoInitialise().

684  {
685  if (m_session->DefinesFunction("InitialConditions"))
686  {
687  for (int i = 0; i < m_fields.num_elements(); ++i)
688  {
690 
691  vType = m_session->GetFunctionType(
692  "InitialConditions", m_session->GetVariable(i));
693 
694  if (vType == LibUtilities::eFunctionTypeFile)
695  {
696  std::string filename
697  = m_session->GetFunctionFilename(
698  "InitialConditions", m_session->GetVariable(i));
699 
700  fs::path pfilename(filename);
701 
702  // redefine path for parallel file which is in directory
703  if(fs::is_directory(pfilename))
704  {
705  fs::path metafile("Info.xml");
706  fs::path fullpath = pfilename / metafile;
707  filename = LibUtilities::PortablePath(fullpath);
708  }
709  m_fld->ImportFieldMetaData(filename, m_fieldMetaDataMap);
710 
711  // check to see if time defined
712  if (m_fieldMetaDataMap !=
714  {
716 
717  iter = m_fieldMetaDataMap.find("Time");
718  if (iter != m_fieldMetaDataMap.end())
719  {
720  time = boost::lexical_cast<NekDouble>(
721  iter->second);
722  }
723  }
724 
725  break;
726  }
727  }
728  }
729  }
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
std::string PortablePath(const boost::filesystem::path &path)
create portable path on different platforms for boost::filesystem path
Definition: FileSystem.cpp:41
double NekDouble
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
LibUtilities::FieldIOSharedPtr m_fld
Field input/output.
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:66
NekDouble Nektar::SolverUtils::UnsteadySystem::GetTimeStep ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray)

Calculate the larger time-step mantaining the problem stable.

Return the timestep to be used for the next step in the time-marching loop.

This function can be overloaded to facilitate solver which utilise a CFL (or other) parameter to determine a maximum timestep under which the problem remains stable.

Definition at line 890 of file UnsteadySystem.cpp.

References v_GetTimeStep().

892  {
893  return v_GetTimeStep(inarray);
894  }
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.
NekDouble Nektar::SolverUtils::UnsteadySystem::MaxTimeStepEstimator ( )
protected

Get the maximum timestep estimator for cfl control.

Returns the maximum time estimator for CFL control.

Definition at line 129 of file UnsteadySystem.cpp.

References ASSERTL0, Nektar::LibUtilities::eAdamsBashforthOrder1, Nektar::LibUtilities::eAdamsBashforthOrder2, Nektar::LibUtilities::eClassicalRungeKutta4, Nektar::LibUtilities::eForwardEuler, Nektar::LibUtilities::eMidpoint, Nektar::LibUtilities::eRungeKutta2, Nektar::LibUtilities::eRungeKutta2_ImprovedEuler, Nektar::LibUtilities::eRungeKutta2_SSP, Nektar::LibUtilities::eRungeKutta3_SSP, Nektar::LibUtilities::eRungeKutta4, and m_intScheme.

130  {
131  NekDouble TimeStability = 0.0;
132  switch(m_intScheme->GetIntegrationMethod())
133  {
137  {
138  TimeStability = 2.784;
139  break;
140  }
147  {
148  TimeStability = 2.0;
149  break;
150  }
152  {
153  TimeStability = 1.0;
154  break;
155  }
156  default:
157  {
158  ASSERTL0(
159  false,
160  "No CFL control implementation for this time"
161  "integration scheme");
162  }
163  }
164  return TimeStability;
165  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
Adams-Bashforth Forward multi-step scheme of order 2.
Runge-Kutta multi-stage scheme 4th order explicit (old name)
Classical RungeKutta4 method (new name for eClassicalRungeKutta4)
Nonlinear SSP RungeKutta3 explicit.
Classical RungeKutta2 method (new name for eMidpoint)
Adams-Bashforth Forward multi-step scheme of order 1.
Nonlinear SSP RungeKutta2 explicit (surrogate for eRungeKutta2_ImprovedEuler)
Improved RungeKutta2 explicit (old name meaning Heun's method)
double NekDouble
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
Wrapper to the time integration scheme.
midpoint method (old name)
void Nektar::SolverUtils::UnsteadySystem::v_AppendOutput1D ( Array< OneD, Array< OneD, NekDouble > > &  solution1D)
protectedvirtual

Print the solution at each solution point in a txt file.

Stores the solution in a file for 1D problems only. This method has been implemented to facilitate the post-processing for 1D problems.

Definition at line 487 of file UnsteadySystem.cpp.

References Nektar::SolverUtils::EquationSystem::GetNpoints(), and Nektar::SolverUtils::EquationSystem::m_fields.

Referenced by v_DoSolve().

489  {
490  // Coordinates of the quadrature points in the real physical space
494  m_fields[0]->GetCoords(x, y, z);
495 
496  // Print out the solution in a txt file
497  ofstream outfile;
498  outfile.open("solution1D.txt");
499  for(int i = 0; i < GetNpoints(); i++)
500  {
501  outfile << scientific << setw (17) << setprecision(16) << x[i]
502  << " " << solution1D[0][i] << endl;
503  }
504  outfile << endl << endl;
505  outfile.close();
506  }
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Nektar::SolverUtils::UnsteadySystem::v_DoInitialise ( void  )
protectedvirtual

Sets up initial conditions.

Sets the initial conditions.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Reimplemented in Nektar::CoupledLinearNS, Nektar::PulseWaveSystem, Nektar::VelocityCorrectionScheme, and Nektar::APE.

Definition at line 449 of file UnsteadySystem.cpp.

References CheckForRestartTime(), Nektar::SolverUtils::EquationSystem::m_time, Nektar::SolverUtils::EquationSystem::SetBoundaryConditions(), and Nektar::SolverUtils::EquationSystem::SetInitialConditions().

450  {
454  }
SOLVER_UTILS_EXPORT void CheckForRestartTime(NekDouble &time)
NekDouble m_time
Current time of simulation.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
SOLVER_UTILS_EXPORT void SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
Initialise the data in the dependent fields.
void Nektar::SolverUtils::UnsteadySystem::v_DoSolve ( void  )
protectedvirtual

Solves an unsteady problem.

Initialises the time integration scheme (as specified in the session file), and perform the time integration.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Reimplemented in Nektar::CoupledLinearNS, and Nektar::PulseWaveSystem.

Definition at line 171 of file UnsteadySystem.cpp.

References ASSERTL0, Nektar::SolverUtils::EquationSystem::Checkpoint_Output(), Nektar::SolverUtils::EquationSystem::eHomogeneous1D, Nektar::SolverUtils::EquationSystem::GetTimeStep(), Nektar::iterator, Nektar::NekConstants::kNekZeroTol, m_cflSafetyFactor, Nektar::SolverUtils::EquationSystem::m_checksteps, Nektar::SolverUtils::EquationSystem::m_checktime, Nektar::SolverUtils::EquationSystem::m_fields, m_filters, Nektar::SolverUtils::EquationSystem::m_fintime, Nektar::SolverUtils::EquationSystem::m_HomogeneousType, m_homoInitialFwd, m_infosteps, m_intScheme, m_intSoln, m_intVariables, m_ode, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_spacedim, Nektar::SolverUtils::EquationSystem::m_steps, Nektar::SolverUtils::EquationSystem::m_time, Nektar::SolverUtils::EquationSystem::m_timestep, Vmath::Nnan(), Nektar::Timer::Start(), Nektar::Timer::Stop(), Nektar::Timer::TimePerTest(), v_AppendOutput1D(), v_PostIntegrate(), and v_PreIntegrate().

172  {
173  ASSERTL0(m_intScheme != 0, "No time integration scheme.");
174 
175  int i, nchk = 1;
176  int nvariables = 0;
177  int nfields = m_fields.num_elements();
178 
179  if (m_intVariables.empty())
180  {
181  for (i = 0; i < nfields; ++i)
182  {
183  m_intVariables.push_back(i);
184  }
185  nvariables = nfields;
186  }
187  else
188  {
189  nvariables = m_intVariables.size();
190  }
191 
192  // Integrate in wave-space if using homogeneous1D
194  {
195  for(i = 0; i < nfields; ++i)
196  {
197  m_fields[i]->HomogeneousFwdTrans(m_fields[i]->GetPhys(),
198  m_fields[i]->UpdatePhys());
199  m_fields[i]->SetWaveSpace(true);
200  m_fields[i]->SetPhysState(false);
201  }
202  }
203 
204  // Set up wrapper to fields data storage.
205  Array<OneD, Array<OneD, NekDouble> > fields(nvariables);
206  Array<OneD, Array<OneD, NekDouble> > tmp (nvariables);
207 
208  // Order storage to list time-integrated fields first.
209  for(i = 0; i < nvariables; ++i)
210  {
211  fields[i] = m_fields[m_intVariables[i]]->GetPhys();
212  m_fields[m_intVariables[i]]->SetPhysState(false);
213  }
214 
215  // Initialise time integration scheme
216  m_intSoln = m_intScheme->InitializeScheme(
217  m_timestep, fields, m_time, m_ode);
218 
219  // Initialise filters
221  for (x = m_filters.begin(); x != m_filters.end(); ++x)
222  {
223  (*x)->Initialise(m_fields, m_time);
224  }
225 
226  // Ensure that there is no conflict of parameters
227  if(m_cflSafetyFactor > 0.0)
228  {
229  // Check final condition
230  ASSERTL0(m_fintime == 0.0 || m_steps == 0,
231  "Final condition not unique: "
232  "fintime > 0.0 and Nsteps > 0");
233 
234  // Check timestep condition
235  ASSERTL0(m_timestep == 0.0,
236  "Timestep not unique: timestep > 0.0 & CFL > 0.0");
237  }
238 
239  // Check uniqueness of checkpoint output
240  ASSERTL0((m_checktime == 0.0 && m_checksteps == 0) ||
241  (m_checktime > 0.0 && m_checksteps == 0) ||
242  (m_checktime == 0.0 && m_checksteps > 0),
243  "Only one of IO_CheckTime and IO_CheckSteps "
244  "should be set!");
245 
246  Timer timer;
247  bool doCheckTime = false;
248  int step = 0;
249  NekDouble intTime = 0.0;
250  NekDouble lastCheckTime = 0.0;
251  NekDouble cpuTime = 0.0;
252  NekDouble elapsed = 0.0;
253 
254  while (step < m_steps ||
256  {
257  if (m_cflSafetyFactor)
258  {
259  m_timestep = GetTimeStep(fields);
260 
261  // Ensure that the final timestep finishes at the final
262  // time, or at a prescribed IO_CheckTime.
263  if (m_time + m_timestep > m_fintime && m_fintime > 0.0)
264  {
266  }
267  else if (m_checktime &&
268  m_time + m_timestep - lastCheckTime >= m_checktime)
269  {
270  lastCheckTime += m_checktime;
271  m_timestep = lastCheckTime - m_time;
272  doCheckTime = true;
273  }
274  }
275 
276  // Perform any solver-specific pre-integration steps
277  if (v_PreIntegrate(step))
278  {
279  break;
280  }
281 
282  timer.Start();
283  fields = m_intScheme->TimeIntegrate(
284  step, m_timestep, m_intSoln, m_ode);
285  timer.Stop();
286 
287  m_time += m_timestep;
288  elapsed = timer.TimePerTest(1);
289  intTime += elapsed;
290  cpuTime += elapsed;
291 
292  // Write out status information
293  if (m_session->GetComm()->GetRank() == 0 &&
294  !((step+1) % m_infosteps))
295  {
296  cout << "Steps: " << setw(8) << left << step+1 << " "
297  << "Time: " << setw(12) << left << m_time;
298 
299  if (m_cflSafetyFactor)
300  {
301  cout << " Time-step: " << setw(12)
302  << left << m_timestep;
303  }
304 
305  stringstream ss;
306  ss << cpuTime << "s";
307  cout << " CPU Time: " << setw(8) << left
308  << ss.str() << endl;
309  cpuTime = 0.0;
310  }
311 
312  // Transform data into coefficient space
313  for (i = 0; i < nvariables; ++i)
314  {
315  m_fields[m_intVariables[i]]->SetPhys(fields[i]);
316  m_fields[m_intVariables[i]]->FwdTrans_IterPerExp(
317  fields[i],
318  m_fields[m_intVariables[i]]->UpdateCoeffs());
319  m_fields[m_intVariables[i]]->SetPhysState(false);
320  }
321 
322  // Perform any solver-specific post-integration steps
323  if (v_PostIntegrate(step))
324  {
325  break;
326  }
327 
328  // search for NaN and quit if found
329  bool nanFound = false;
330  for (i = 0; i < nvariables; ++i)
331  {
332  if (Vmath::Nnan(fields[i].num_elements(), fields[i], 1) > 0)
333  {
334  cout << "NaN found in variable \""
335  << m_session->GetVariable(i)
336  << "\", terminating" << endl;
337  nanFound = true;
338  }
339  }
340 
341  if (nanFound)
342  {
343  break;
344  }
345 
346  // Update filters
348  for (x = m_filters.begin(); x != m_filters.end(); ++x)
349  {
350  (*x)->Update(m_fields, m_time);
351  }
352 
353  // Write out checkpoint files
354  if ((m_checksteps && step && !((step + 1) % m_checksteps)) ||
355  doCheckTime)
356  {
358  {
359  vector<bool> transformed(nfields, false);
360  for(i = 0; i < nfields; i++)
361  {
362  if (m_fields[i]->GetWaveSpace())
363  {
364  m_fields[i]->SetWaveSpace(false);
365  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
366  m_fields[i]->UpdatePhys());
367  m_fields[i]->SetPhysState(true);
368  transformed[i] = true;
369  }
370  }
371  Checkpoint_Output(nchk++);
372  for(i = 0; i < nfields; i++)
373  {
374  if (transformed[i])
375  {
376  m_fields[i]->SetWaveSpace(true);
377  m_fields[i]->HomogeneousFwdTrans(
378  m_fields[i]->GetPhys(),
379  m_fields[i]->UpdatePhys());
380  m_fields[i]->SetPhysState(false);
381  }
382  }
383  }
384  else
385  {
386  Checkpoint_Output(nchk++);
387  }
388  doCheckTime = false;
389  }
390 
391  // Step advance
392  ++step;
393  }
394 
395  // Print out summary statistics
396  if (m_session->GetComm()->GetRank() == 0)
397  {
398  if (m_cflSafetyFactor > 0.0)
399  {
400  cout << "CFL safety factor : " << m_cflSafetyFactor << endl
401  << "CFL time-step : " << m_timestep << endl;
402  }
403 
404  if (m_session->GetSolverInfo("Driver") != "SteadyState")
405  {
406  cout << "Time-integration : " << intTime << "s" << endl;
407  }
408  }
409 
410  // If homogeneous, transform back into physical space if necessary.
412  {
413  for(i = 0; i < nfields; i++)
414  {
415  if (m_fields[i]->GetWaveSpace())
416  {
417  m_fields[i]->SetWaveSpace(false);
418  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
419  m_fields[i]->UpdatePhys());
420  m_fields[i]->SetPhysState(true);
421  }
422  }
423  }
424  else
425  {
426  for(i = 0; i < nvariables; ++i)
427  {
428  m_fields[m_intVariables[i]]->SetPhys(fields[i]);
429  m_fields[m_intVariables[i]]->SetPhysState(true);
430  }
431  }
432 
433  // Finalise filters
434  for (x = m_filters.begin(); x != m_filters.end(); ++x)
435  {
436  (*x)->Finalise(m_fields, m_time);
437  }
438 
439  // Print for 1D problems
440  if(m_spacedim == 1)
441  {
442  v_AppendOutput1D(fields);
443  }
444  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate(int step)
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
NekDouble m_time
Current time of simulation.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
NekDouble m_timestep
Time step size.
SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n)
Write checkpoint file of m_fields.
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D(Array< OneD, Array< OneD, NekDouble > > &solution1D)
Print the solution at each solution point in a txt file.
NekDouble m_checktime
Time between checkpoints.
int m_checksteps
Number of steps between checkpoints.
NekDouble m_fintime
Finish time of the simulation.
int m_steps
Number of steps to take.
void Stop()
Definition: Timer.cpp:62
static const NekDouble kNekZeroTol
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate(int step)
int Nnan(int n, const T *x, const int incx)
Return number of NaN elements of x.
Definition: Vmath.cpp:869
int m_spacedim
Spatial dimension (>= expansion dim).
double NekDouble
SOLVER_UTILS_EXPORT NekDouble GetTimeStep()
NekDouble m_cflSafetyFactor
CFL safety factor (comprise between 0 to 1).
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
Wrapper to the time integration scheme.
std::vector< FilterSharedPtr > m_filters
int m_infosteps
Number of time steps between outputting status information.
void Start()
Definition: Timer.cpp:51
NekDouble TimePerTest(unsigned int n)
Returns amount of seconds per iteration in a test with n iterations.
Definition: Timer.cpp:108
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
enum HomogeneousType m_HomogeneousType
void Nektar::SolverUtils::UnsteadySystem::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::EquationSystem.

Reimplemented in Nektar::CoupledLinearNS, Nektar::CompressibleFlowSystem, Nektar::UnsteadyAdvectionDiffusion, Nektar::CFLtester, Nektar::UnsteadyAdvection, Nektar::NonlinearPeregrine, Nektar::VelocityCorrectionScheme, Nektar::PulseWavePropagation, Nektar::Monodomain, Nektar::Bidomain, Nektar::BidomainRoth, Nektar::ShallowWaterSystem, Nektar::LinearSWE, Nektar::EulerADCFE, Nektar::EulerCFE, Nektar::NonlinearSWE, Nektar::NavierStokesCFE, Nektar::ImageWarpingSystem, and Nektar::UnsteadyDiffusion.

Definition at line 460 of file UnsteadySystem.cpp.

References Nektar::SolverUtils::AddSummaryItem(), Nektar::SolverUtils::EquationSystem::m_checksteps, m_explicitAdvection, m_explicitDiffusion, m_explicitReaction, m_intScheme, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_steps, Nektar::SolverUtils::EquationSystem::m_timestep, Nektar::LibUtilities::TimeIntegrationMethodMap, and Nektar::SolverUtils::EquationSystem::v_GenerateSummary().

Referenced by Nektar::PulseWavePropagation::v_GenerateSummary().

461  {
463  AddSummaryItem(s, "Advection",
464  m_explicitAdvection ? "explicit" : "implicit");
465  AddSummaryItem(s, "Diffusion",
466  m_explicitDiffusion ? "explicit" : "implicit");
467 
468  if (m_session->GetSolverInfo("EQTYPE")
469  == "SteadyAdvectionDiffusionReaction")
470  {
471  AddSummaryItem(s, "Reaction",
472  m_explicitReaction ? "explicit" : "implicit");
473  }
474 
475  AddSummaryItem(s, "Time Step", m_timestep);
476  AddSummaryItem(s, "No. of Steps", m_steps);
477  AddSummaryItem(s, "Checkpoints (steps)", m_checksteps);
478  AddSummaryItem(s, "Integration Type",
480  m_intScheme->GetIntegrationMethod()]);
481  }
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &l)
Virtual function for generating summary information.
bool m_explicitReaction
Indicates if explicit or implicit treatment of reaction is used.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
NekDouble m_timestep
Time step size.
const char *const TimeIntegrationMethodMap[]
int m_checksteps
Number of steps between checkpoints.
int m_steps
Number of steps to take.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:50
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
Wrapper to the time integration scheme.
NekDouble Nektar::SolverUtils::UnsteadySystem::v_GetTimeStep ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray)
protectedvirtual

Return the timestep to be used for the next step in the time-marching loop.

See also
UnsteadySystem::GetTimeStep

Reimplemented in Nektar::CompressibleFlowSystem.

Definition at line 902 of file UnsteadySystem.cpp.

References ASSERTL0.

Referenced by GetTimeStep().

904  {
905  ASSERTL0(false, "Not defined for this class");
906  return 0.0;
907  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void Nektar::SolverUtils::UnsteadySystem::v_InitObject ( )
protectedvirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Reimplemented in Nektar::CoupledLinearNS, Nektar::CompressibleFlowSystem, Nektar::PulseWaveSystem, Nektar::UnsteadyAdvectionDiffusion, Nektar::IncNavierStokes, Nektar::CFLtester, Nektar::UnsteadyAdvection, Nektar::UnsteadyInviscidBurger, Nektar::ShallowWaterSystem, Nektar::NonlinearPeregrine, Nektar::APE, Nektar::EulerADCFE, Nektar::NavierStokesCFE, Nektar::EulerCFE, Nektar::PulseWavePropagation, Nektar::ImageWarpingSystem, Nektar::UnsteadyDiffusion, Nektar::Monodomain, Nektar::Bidomain, Nektar::BidomainRoth, Nektar::LinearSWE, Nektar::NonlinearSWE, Nektar::PulseWaveSystemOutput, Nektar::VelocityCorrectionScheme, and Nektar::SolverUtils::AdvectionSystem.

Definition at line 76 of file UnsteadySystem.cpp.

References Nektar::SolverUtils::GetFilterFactory(), Nektar::LibUtilities::GetTimeIntegrationWrapperFactory(), m_cflSafetyFactor, m_explicitAdvection, m_explicitDiffusion, m_explicitReaction, Nektar::SolverUtils::EquationSystem::m_fieldMetaDataMap, m_filters, m_homoInitialFwd, m_infosteps, m_intScheme, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_time, and Nektar::SolverUtils::EquationSystem::v_InitObject().

Referenced by Nektar::SolverUtils::AdvectionSystem::v_InitObject().

77  {
79 
80  // Load SolverInfo parameters
81  m_session->MatchSolverInfo("DIFFUSIONADVANCEMENT","Explicit",
82  m_explicitDiffusion,true);
83  m_session->MatchSolverInfo("ADVECTIONADVANCEMENT","Explicit",
84  m_explicitAdvection,true);
85  m_session->MatchSolverInfo("REACTIONADVANCEMENT", "Explicit",
86  m_explicitReaction, true);
87 
88  // For steady problems, we do not initialise the time integration
89  if (m_session->DefinesSolverInfo("TIMEINTEGRATIONMETHOD"))
90  {
92  CreateInstance(m_session->GetSolverInfo(
93  "TIMEINTEGRATIONMETHOD"));
94 
95  // Load generic input parameters
96  m_session->LoadParameter("IO_InfoSteps", m_infosteps, 0);
97  m_session->LoadParameter("CFL", m_cflSafetyFactor, 0.0);
98 
99  // Set up time to be dumped in field information
100  m_fieldMetaDataMap["Time"] =
101  boost::lexical_cast<std::string>(m_time);
102  }
103 
104  // By default attempt to forward transform initial condition.
105  m_homoInitialFwd = true;
106 
107  // Set up filters
108  LibUtilities::FilterMap::const_iterator x;
109  LibUtilities::FilterMap f = m_session->GetFilters();
110  for (x = f.begin(); x != f.end(); ++x)
111  {
112  m_filters.push_back(GetFilterFactory().CreateInstance(
113  x->first,
114  m_session,
115  x->second));
116  }
117  }
bool m_explicitReaction
Indicates if explicit or implicit treatment of reaction is used.
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
NekDouble m_time
Current time of simulation.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Initialisation object for EquationSystem.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
std::vector< std::pair< std::string, FilterParams > > FilterMap
Definition: SessionReader.h:66
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
TimeIntegrationWrapperFactory & GetTimeIntegrationWrapperFactory()
NekDouble m_cflSafetyFactor
CFL safety factor (comprise between 0 to 1).
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
Wrapper to the time integration scheme.
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:42
std::vector< FilterSharedPtr > m_filters
int m_infosteps
Number of time steps between outputting status information.
void Nektar::SolverUtils::UnsteadySystem::v_NumericalFlux ( Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  numflux 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Reimplemented in Nektar::IncNavierStokes, Nektar::CFLtester, Nektar::PulseWavePropagation, and Nektar::ImageWarpingSystem.

Definition at line 508 of file UnsteadySystem.cpp.

References ASSERTL0.

511  {
512  ASSERTL0(false,
513  "This function is not implemented for this equation.");
514  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void Nektar::SolverUtils::UnsteadySystem::v_NumericalFlux ( Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  numfluxX,
Array< OneD, Array< OneD, NekDouble > > &  numfluxY 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 516 of file UnsteadySystem.cpp.

References ASSERTL0.

520  {
521  ASSERTL0(false,
522  "This function is not implemented for this equation.");
523  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void Nektar::SolverUtils::UnsteadySystem::v_NumFluxforScalar ( const Array< OneD, Array< OneD, NekDouble > > &  ufield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  uflux 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 525 of file UnsteadySystem.cpp.

References Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_traceNormals, Vmath::Vmul(), and WeakPenaltyforScalar().

528  {
529  int i, j;
530  int nTraceNumPoints = GetTraceNpoints();
531  int nvariables = m_fields.num_elements();
532  int nqvar = uflux.num_elements();
533 
534  Array<OneD, NekDouble > Fwd (nTraceNumPoints);
535  Array<OneD, NekDouble > Bwd (nTraceNumPoints);
536  Array<OneD, NekDouble > Vn (nTraceNumPoints, 0.0);
537  Array<OneD, NekDouble > fluxtemp(nTraceNumPoints, 0.0);
538 
539  // Get the sign of (v \cdot n), v = an arbitrary vector
540 
541  // Evaulate upwind flux:
542  // uflux = \hat{u} \phi \cdot u = u^{(+,-)} n
543  for (j = 0; j < nqvar; ++j)
544  {
545  for (i = 0; i < nvariables ; ++i)
546  {
547  // Compute Fwd and Bwd value of ufield of i direction
548  m_fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
549 
550  // if Vn >= 0, flux = uFwd, i.e.,
551  // edge::eForward, if V*n>=0 <=> V*n_F>=0, pick uflux = uFwd
552  // edge::eBackward, if V*n>=0 <=> V*n_B<0, pick uflux = uFwd
553 
554  // else if Vn < 0, flux = uBwd, i.e.,
555  // edge::eForward, if V*n<0 <=> V*n_F<0, pick uflux = uBwd
556  // edge::eBackward, if V*n<0 <=> V*n_B>=0, pick uflux = uBwd
557 
558  m_fields[i]->GetTrace()->Upwind(m_traceNormals[j],
559  Fwd, Bwd, fluxtemp);
560 
561  // Imposing weak boundary condition with flux
562  // if Vn >= 0, uflux = uBwd at Neumann, i.e.,
563  // edge::eForward, if V*n>=0 <=> V*n_F>=0, pick uflux = uBwd
564  // edge::eBackward, if V*n>=0 <=> V*n_B<0, pick uflux = uBwd
565 
566  // if Vn >= 0, uflux = uFwd at Neumann, i.e.,
567  // edge::eForward, if V*n<0 <=> V*n_F<0, pick uflux = uFwd
568  // edge::eBackward, if V*n<0 <=> V*n_B>=0, pick uflux = uFwd
569 
570  if(m_fields[0]->GetBndCondExpansions().num_elements())
571  {
572  WeakPenaltyforScalar(i, ufield[i], fluxtemp);
573  }
574 
575  // if Vn >= 0, flux = uFwd*(tan_{\xi}^- \cdot \vec{n}),
576  // i.e,
577  // edge::eForward, uFwd \(\tan_{\xi}^Fwd \cdot \vec{n})
578  // edge::eBackward, uFwd \(\tan_{\xi}^Bwd \cdot \vec{n})
579 
580  // else if Vn < 0, flux = uBwd*(tan_{\xi}^- \cdot \vec{n}),
581  // i.e,
582  // edge::eForward, uBwd \(\tan_{\xi}^Fwd \cdot \vec{n})
583  // edge::eBackward, uBwd \(\tan_{\xi}^Bwd \cdot \vec{n})
584 
585  Vmath::Vmul(nTraceNumPoints,
586  m_traceNormals[j], 1,
587  fluxtemp, 1,
588  uflux[j][i], 1);
589  }
590  }
591  }
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()
void WeakPenaltyforScalar(const int var, const Array< OneD, const NekDouble > &physfield, Array< OneD, NekDouble > &penaltyflux, NekDouble time=0.0)
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:169
void Nektar::SolverUtils::UnsteadySystem::v_NumFluxforVector ( const Array< OneD, Array< OneD, NekDouble > > &  ufield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  qfield,
Array< OneD, Array< OneD, NekDouble > > &  qflux 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 595 of file UnsteadySystem.cpp.

References Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_traceNormals, Vmath::Smul(), Vmath::Vadd(), Vmath::Vmul(), Vmath::Vsub(), and WeakPenaltyforVector().

599  {
600  int nTraceNumPoints = GetTraceNpoints();
601  int nvariables = m_fields.num_elements();
602  int nqvar = qfield.num_elements();
603 
604  NekDouble C11 = 1.0;
605  Array<OneD, NekDouble > Fwd(nTraceNumPoints);
606  Array<OneD, NekDouble > Bwd(nTraceNumPoints);
607  Array<OneD, NekDouble > Vn (nTraceNumPoints, 0.0);
608 
609  Array<OneD, NekDouble > qFwd (nTraceNumPoints);
610  Array<OneD, NekDouble > qBwd (nTraceNumPoints);
611  Array<OneD, NekDouble > qfluxtemp(nTraceNumPoints, 0.0);
612 
613  Array<OneD, NekDouble > uterm(nTraceNumPoints);
614 
615  // Evaulate upwind flux:
616  // qflux = \hat{q} \cdot u = q \cdot n - C_(11)*(u^+ - u^-)
617  for (int i = 0; i < nvariables; ++i)
618  {
619  qflux[i] = Array<OneD, NekDouble> (nTraceNumPoints, 0.0);
620  for (int j = 0; j < nqvar; ++j)
621  {
622  // Compute Fwd and Bwd value of ufield of jth direction
623  m_fields[i]->GetFwdBwdTracePhys(qfield[j][i],qFwd,qBwd);
624 
625  // if Vn >= 0, flux = uFwd, i.e.,
626  // edge::eForward, if V*n>=0 <=> V*n_F>=0, pick
627  // qflux = qBwd = q+
628  // edge::eBackward, if V*n>=0 <=> V*n_B<0, pick
629  // qflux = qBwd = q-
630 
631  // else if Vn < 0, flux = uBwd, i.e.,
632  // edge::eForward, if V*n<0 <=> V*n_F<0, pick
633  // qflux = qFwd = q-
634  // edge::eBackward, if V*n<0 <=> V*n_B>=0, pick
635  // qflux = qFwd = q+
636 
637  m_fields[i]->GetTrace()->Upwind(m_traceNormals[j],
638  qBwd, qFwd,
639  qfluxtemp);
640 
641  Vmath::Vmul(nTraceNumPoints,
642  m_traceNormals[j], 1,
643  qfluxtemp, 1,
644  qfluxtemp, 1);
645 
646  // Generate Stability term = - C11 ( u- - u+ )
647  m_fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
648 
649  Vmath::Vsub(nTraceNumPoints,
650  Fwd, 1, Bwd, 1,
651  uterm, 1);
652 
653  Vmath::Smul(nTraceNumPoints,
654  -1.0 * C11, uterm, 1,
655  uterm, 1);
656 
657  // Flux = {Fwd, Bwd} * (nx, ny, nz) + uterm * (nx, ny)
658  Vmath::Vadd(nTraceNumPoints,
659  uterm, 1,
660  qfluxtemp, 1,
661  qfluxtemp, 1);
662 
663  // Imposing weak boundary condition with flux
664  if (m_fields[0]->GetBndCondExpansions().num_elements())
665  {
666  WeakPenaltyforVector(i, j,
667  qfield[j][i],
668  qfluxtemp,
669  C11);
670  }
671 
672  // q_hat \cdot n = (q_xi \cdot n_xi) or (q_eta \cdot n_eta)
673  // n_xi = n_x * tan_xi_x + n_y * tan_xi_y + n_z * tan_xi_z
674  // n_xi = n_x * tan_eta_x + n_y * tan_eta_y + n_z*tan_eta_z
675  Vmath::Vadd(nTraceNumPoints,
676  qfluxtemp, 1,
677  qflux[i], 1,
678  qflux[i], 1);
679  }
680  }
681  }
void WeakPenaltyforVector(const int var, const int dir, const Array< OneD, const NekDouble > &physfield, Array< OneD, NekDouble > &penaltyflux, NekDouble C11, NekDouble time=0.0)
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
double NekDouble
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:329
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetTraceNpoints()
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:285
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:169
bool Nektar::SolverUtils::UnsteadySystem::v_PostIntegrate ( int  step)
protectedvirtual

Reimplemented in Nektar::CompressibleFlowSystem, and Nektar::IncNavierStokes.

Definition at line 914 of file UnsteadySystem.cpp.

Referenced by v_DoSolve().

915  {
916  return false;
917  }
bool Nektar::SolverUtils::UnsteadySystem::v_PreIntegrate ( int  step)
protectedvirtual

Reimplemented in Nektar::IncNavierStokes.

Definition at line 909 of file UnsteadySystem.cpp.

Referenced by v_DoSolve().

910  {
911  return false;
912  }
bool Nektar::SolverUtils::UnsteadySystem::v_SteadyStateCheck ( int  step)
protectedvirtual

Definition at line 919 of file UnsteadySystem.cpp.

920  {
921  return false;
922  }
void Nektar::SolverUtils::UnsteadySystem::WeakPenaltyforScalar ( const int  var,
const Array< OneD, const NekDouble > &  physfield,
Array< OneD, NekDouble > &  penaltyflux,
NekDouble  time = 0.0 
)
private

Definition at line 731 of file UnsteadySystem.cpp.

References Nektar::SpatialDomains::eDirichlet, Nektar::SpatialDomains::eNeumann, Nektar::SolverUtils::EquationSystem::GetPhys_Offset(), Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_session, and Vmath::Vcopy().

Referenced by v_NumFluxforScalar().

736  {
737  int i, e, npoints, id1, id2;
738 
739  // Number of boundary regions
740  int nbnd = m_fields[var]->GetBndCondExpansions().num_elements();
741  int Nfps, numBDEdge;
742  int nTraceNumPoints = GetTraceNpoints();
743  int cnt = 0;
744 
745  Array<OneD, NekDouble > uplus(nTraceNumPoints);
746 
747  m_fields[var]->ExtractTracePhys(physfield, uplus);
748  for (i = 0; i < nbnd; ++i)
749  {
750  // Number of boundary expansion related to that region
751  numBDEdge = m_fields[var]->
752  GetBndCondExpansions()[i]->GetExpSize();
753 
754  // Evaluate boundary values g_D or g_N from input files
756  m_session->GetFunction("InitialConditions", 0);
757 
758  npoints = m_fields[var]->
759  GetBndCondExpansions()[i]->GetNpoints();
760 
761  Array<OneD,NekDouble> BDphysics(npoints);
762  Array<OneD,NekDouble> x0(npoints,0.0);
763  Array<OneD,NekDouble> x1(npoints,0.0);
764  Array<OneD,NekDouble> x2(npoints,0.0);
765 
766  m_fields[var]->GetBndCondExpansions()[i]->GetCoords(x0,x1,x2);
767  ifunc->Evaluate(x0,x1,x2,time,BDphysics);
768 
769  // Weakly impose boundary conditions by modifying flux values
770  for (e = 0; e < numBDEdge ; ++e)
771  {
772  // Number of points on the expansion
773  Nfps = m_fields[var]->
774  GetBndCondExpansions()[i]->GetExp(e)->GetNumPoints(0);
775 
776  id1 = m_fields[var]->
777  GetBndCondExpansions()[i]->GetPhys_Offset(e);
778 
779  id2 = m_fields[0]->GetTrace()->
780  GetPhys_Offset(m_fields[0]->GetTraceMap()->
781  GetBndCondTraceToGlobalTraceMap(cnt++));
782 
783  // For Dirichlet boundary condition: uflux = g_D
784  if (m_fields[var]->GetBndConditions()[i]->
785  GetBoundaryConditionType() == SpatialDomains::eDirichlet)
786  {
787  Vmath::Vcopy(Nfps,
788  &BDphysics[id1], 1,
789  &penaltyflux[id2], 1);
790  }
791  // For Neumann boundary condition: uflux = u+
792  else if ((m_fields[var]->GetBndConditions()[i])->
793  GetBoundaryConditionType() == SpatialDomains::eNeumann)
794  {
795  Vmath::Vcopy(Nfps,
796  &uplus[id2], 1,
797  &penaltyflux[id2], 1);
798  }
799  }
800  }
801  }
boost::shared_ptr< Equation > EquationSharedPtr
SOLVER_UTILS_EXPORT int GetPhys_Offset(int n)
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT int GetTraceNpoints()
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1038
void Nektar::SolverUtils::UnsteadySystem::WeakPenaltyforVector ( const int  var,
const int  dir,
const Array< OneD, const NekDouble > &  physfield,
Array< OneD, NekDouble > &  penaltyflux,
NekDouble  C11,
NekDouble  time = 0.0 
)
private

Diffusion: Imposing weak boundary condition for q with flux uflux = g_D on Dirichlet boundary condition uflux = u_Fwd on Neumann boundary condition

Definition at line 808 of file UnsteadySystem.cpp.

References Nektar::SpatialDomains::eDirichlet, Nektar::SpatialDomains::eNeumann, Nektar::SolverUtils::EquationSystem::GetPhys_Offset(), Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_traceNormals, and Vmath::Vmul().

Referenced by v_NumFluxforVector().

815  {
816  int i, e, npoints, id1, id2;
817  int nbnd = m_fields[var]->GetBndCondExpansions().num_elements();
818  int numBDEdge, Nfps;
819  int nTraceNumPoints = GetTraceNpoints();
820  Array<OneD, NekDouble > uterm(nTraceNumPoints);
821  Array<OneD, NekDouble > qtemp(nTraceNumPoints);
822  int cnt = 0;
823 
824  m_fields[var]->ExtractTracePhys(physfield,qtemp);
825 
826  for (i = 0; i < nbnd; ++i)
827  {
828  numBDEdge = m_fields[var]->
829  GetBndCondExpansions()[i]->GetExpSize();
830 
831  // Evaluate boundary values g_D or g_N from input files
833  m_session->GetFunction("InitialConditions", 0);
834 
835  npoints = m_fields[var]->
836  GetBndCondExpansions()[i]->GetNpoints();
837 
838  Array<OneD,NekDouble> BDphysics(npoints);
839  Array<OneD,NekDouble> x0(npoints,0.0);
840  Array<OneD,NekDouble> x1(npoints,0.0);
841  Array<OneD,NekDouble> x2(npoints,0.0);
842 
843  m_fields[var]->GetBndCondExpansions()[i]->GetCoords(x0,x1,x2);
844  ifunc->Evaluate(x0,x1,x2,time,BDphysics);
845 
846  // Weakly impose boundary conditions by modifying flux values
847  for (e = 0; e < numBDEdge ; ++e)
848  {
849  Nfps = m_fields[var]->
850  GetBndCondExpansions()[i]->GetExp(e)->GetNumPoints(0);
851 
852  id1 = m_fields[var]->
853  GetBndCondExpansions()[i]->GetPhys_Offset(e);
854 
855  id2 = m_fields[0]->GetTrace()->
856  GetPhys_Offset(m_fields[0]->GetTraceMap()->
857  GetBndCondTraceToGlobalTraceMap(cnt++));
858 
859  // For Dirichlet boundary condition:
860  //qflux = q+ - C_11 (u+ - g_D) (nx, ny)
861  if(m_fields[var]->GetBndConditions()[i]->
862  GetBoundaryConditionType() == SpatialDomains::eDirichlet)
863  {
864  Vmath::Vmul(Nfps,
865  &m_traceNormals[dir][id2], 1,
866  &qtemp[id2], 1,
867  &penaltyflux[id2], 1);
868  }
869  // For Neumann boundary condition: qflux = g_N
870  else if((m_fields[var]->GetBndConditions()[i])->
871  GetBoundaryConditionType() == SpatialDomains::eNeumann)
872  {
873  Vmath::Vmul(Nfps,
874  &m_traceNormals[dir][id2], 1,
875  &BDphysics[id1], 1,
876  &penaltyflux[id2], 1);
877  }
878  }
879  }
880  }
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
boost::shared_ptr< Equation > EquationSharedPtr
SOLVER_UTILS_EXPORT int GetPhys_Offset(int n)
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT int GetTraceNpoints()
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:169

Member Data Documentation

NekDouble Nektar::SolverUtils::UnsteadySystem::m_cflSafetyFactor

CFL safety factor (comprise between 0 to 1).

Definition at line 59 of file UnsteadySystem.h.

Referenced by v_DoSolve(), and v_InitObject().

NekDouble Nektar::SolverUtils::UnsteadySystem::m_epsilon
protected

Definition at line 71 of file UnsteadySystem.h.

bool Nektar::SolverUtils::UnsteadySystem::m_explicitAdvection
protected
bool Nektar::SolverUtils::UnsteadySystem::m_explicitDiffusion
protected
bool Nektar::SolverUtils::UnsteadySystem::m_explicitReaction
protected

Indicates if explicit or implicit treatment of reaction is used.

Definition at line 77 of file UnsteadySystem.h.

Referenced by v_GenerateSummary(), and v_InitObject().

std::vector<FilterSharedPtr> Nektar::SolverUtils::UnsteadySystem::m_filters
protected
bool Nektar::SolverUtils::UnsteadySystem::m_homoInitialFwd
protected
int Nektar::SolverUtils::UnsteadySystem::m_infosteps
protected
LibUtilities::TimeIntegrationWrapperSharedPtr Nektar::SolverUtils::UnsteadySystem::m_intScheme
protected
LibUtilities::TimeIntegrationSolutionSharedPtr Nektar::SolverUtils::UnsteadySystem::m_intSoln
protected
std::vector<int> Nektar::SolverUtils::UnsteadySystem::m_intVariables
protected
LibUtilities::TimeIntegrationSchemeOperators Nektar::SolverUtils::UnsteadySystem::m_ode
protected