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

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

662  {
663  if (m_session->DefinesFunction("InitialConditions"))
664  {
665  for (int i = 0; i < m_fields.num_elements(); ++i)
666  {
668 
669  vType = m_session->GetFunctionType(
670  "InitialConditions", m_session->GetVariable(i));
671 
672  if (vType == LibUtilities::eFunctionTypeFile)
673  {
674  std::string filename
675  = m_session->GetFunctionFilename(
676  "InitialConditions", m_session->GetVariable(i));
677 
678  fs::path pfilename(filename);
679 
680  // redefine path for parallel file which is in directory
681  if(fs::is_directory(pfilename))
682  {
683  fs::path metafile("Info.xml");
684  fs::path fullpath = pfilename / metafile;
685  filename = LibUtilities::PortablePath(fullpath);
686  }
687  m_fld->ImportFieldMetaData(filename, m_fieldMetaDataMap);
688 
689  // check to see if time defined
690  if (m_fieldMetaDataMap !=
692  {
694 
695  iter = m_fieldMetaDataMap.find("Time");
696  if (iter != m_fieldMetaDataMap.end())
697  {
698  time = boost::lexical_cast<NekDouble>(
699  iter->second);
700  }
701  }
702 
703  break;
704  }
705  }
706  }
707  }
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 868 of file UnsteadySystem.cpp.

References v_GetTimeStep().

870  {
871  return v_GetTimeStep(inarray);
872  }
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::eRungeKutta2_ImprovedEuler, Nektar::LibUtilities::eRungeKutta2_ModifiedEuler, and m_intScheme.

130  {
131  NekDouble TimeStability = 0.0;
132  switch(m_intScheme->GetIntegrationMethod())
133  {
136  {
137  TimeStability = 2.784;
138  break;
139  }
143  {
144  TimeStability = 2.0;
145  break;
146  }
148  {
149  TimeStability = 1.0;
150  break;
151  }
152  default:
153  {
154  ASSERTL0(
155  false,
156  "No CFL control implementation for this time"
157  "integration scheme");
158  }
159  }
160  return TimeStability;
161  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
Adams-Bashforth Forward multi-step scheme of order 2.
Runge-Kutta multi-stage scheme 4th order explicit.
Adams-Bashforth Forward multi-step scheme of order 1.
Runge-Kutta multi-stage scheme 2nd order explicit.
double NekDouble
Runge-Kutta multi-stage scheme 2nd order explicit.
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
Wrapper to the time integration scheme.
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 465 of file UnsteadySystem.cpp.

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

Referenced by v_DoSolve().

467  {
468  // Coordinates of the quadrature points in the real physical space
472  m_fields[0]->GetCoords(x, y, z);
473 
474  // Print out the solution in a txt file
475  ofstream outfile;
476  outfile.open("solution1D.txt");
477  for(int i = 0; i < GetNpoints(); i++)
478  {
479  outfile << scientific << setw (17) << setprecision(16) << x[i]
480  << " " << solution1D[0][i] << endl;
481  }
482  outfile << endl << endl;
483  outfile.close();
484  }
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 427 of file UnsteadySystem.cpp.

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

428  {
432  }
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 167 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, Nektar::Timer::Start(), Nektar::Timer::Stop(), Nektar::Timer::TimePerTest(), v_AppendOutput1D(), v_PostIntegrate(), and v_PreIntegrate().

168  {
169  ASSERTL0(m_intScheme != 0, "No time integration scheme.");
170 
171  int i, nchk = 1;
172  int nvariables = 0;
173  int nfields = m_fields.num_elements();
174 
175  if (m_intVariables.empty())
176  {
177  for (i = 0; i < nfields; ++i)
178  {
179  m_intVariables.push_back(i);
180  }
181  nvariables = nfields;
182  }
183  else
184  {
185  nvariables = m_intVariables.size();
186  }
187 
188  // Integrate in wave-space if using homogeneous1D
190  {
191  for(i = 0; i < nfields; ++i)
192  {
193  m_fields[i]->HomogeneousFwdTrans(m_fields[i]->GetPhys(),
194  m_fields[i]->UpdatePhys());
195  m_fields[i]->SetWaveSpace(true);
196  m_fields[i]->SetPhysState(false);
197  }
198  }
199 
200  // Set up wrapper to fields data storage.
201  Array<OneD, Array<OneD, NekDouble> > fields(nvariables);
202  Array<OneD, Array<OneD, NekDouble> > tmp (nvariables);
203 
204  // Order storage to list time-integrated fields first.
205  for(i = 0; i < nvariables; ++i)
206  {
207  fields[i] = m_fields[m_intVariables[i]]->GetPhys();
208  m_fields[m_intVariables[i]]->SetPhysState(false);
209  }
210 
211  // Initialise time integration scheme
212  m_intSoln = m_intScheme->InitializeScheme(
213  m_timestep, fields, m_time, m_ode);
214 
215  // Initialise filters
217  for (x = m_filters.begin(); x != m_filters.end(); ++x)
218  {
219  (*x)->Initialise(m_fields, m_time);
220  }
221 
222  // Ensure that there is no conflict of parameters
223  if(m_cflSafetyFactor > 0.0)
224  {
225  // Check final condition
226  ASSERTL0(m_fintime == 0.0 || m_steps == 0,
227  "Final condition not unique: "
228  "fintime > 0.0 and Nsteps > 0");
229 
230  // Check timestep condition
231  ASSERTL0(m_timestep == 0.0,
232  "Timestep not unique: timestep > 0.0 & CFL > 0.0");
233  }
234 
235  // Check uniqueness of checkpoint output
236  ASSERTL0((m_checktime == 0.0 && m_checksteps == 0) ||
237  (m_checktime > 0.0 && m_checksteps == 0) ||
238  (m_checktime == 0.0 && m_checksteps > 0),
239  "Only one of IO_CheckTime and IO_CheckSteps "
240  "should be set!");
241 
242  Timer timer;
243  bool doCheckTime = false;
244  int step = 0;
245  NekDouble intTime = 0.0;
246  NekDouble lastCheckTime = 0.0;
247  NekDouble cpuTime = 0.0;
248  NekDouble elapsed = 0.0;
249 
250  while (step < m_steps ||
252  {
253  if (m_cflSafetyFactor)
254  {
255  m_timestep = GetTimeStep(fields);
256 
257  // Ensure that the final timestep finishes at the final
258  // time, or at a prescribed IO_CheckTime.
259  if (m_time + m_timestep > m_fintime && m_fintime > 0.0)
260  {
262  }
263  else if (m_checktime &&
264  m_time + m_timestep - lastCheckTime >= m_checktime)
265  {
266  lastCheckTime += m_checktime;
267  m_timestep = lastCheckTime - m_time;
268  doCheckTime = true;
269  }
270  }
271 
272  // Perform any solver-specific pre-integration steps
273  if (v_PreIntegrate(step))
274  {
275  break;
276  }
277 
278  timer.Start();
279  fields = m_intScheme->TimeIntegrate(
280  step, m_timestep, m_intSoln, m_ode);
281  timer.Stop();
282 
283  m_time += m_timestep;
284  elapsed = timer.TimePerTest(1);
285  intTime += elapsed;
286  cpuTime += elapsed;
287 
288  // Write out status information
289  if (m_session->GetComm()->GetRank() == 0 &&
290  !((step+1) % m_infosteps))
291  {
292  cout << "Steps: " << setw(8) << left << step+1 << " "
293  << "Time: " << setw(12) << left << m_time;
294 
295  if (m_cflSafetyFactor)
296  {
297  cout << " Time-step: " << setw(12)
298  << left << m_timestep;
299  }
300 
301  stringstream ss;
302  ss << cpuTime << "s";
303  cout << " CPU Time: " << setw(8) << left
304  << ss.str() << endl;
305  cpuTime = 0.0;
306  }
307 
308  // Transform data into coefficient space
309  for (i = 0; i < nvariables; ++i)
310  {
311  m_fields[m_intVariables[i]]->SetPhys(fields[i]);
312  m_fields[m_intVariables[i]]->FwdTrans_IterPerExp(
313  fields[i],
314  m_fields[m_intVariables[i]]->UpdateCoeffs());
315  m_fields[m_intVariables[i]]->SetPhysState(false);
316  }
317 
318  // Perform any solver-specific post-integration steps
319  if (v_PostIntegrate(step))
320  {
321  break;
322  }
323 
324  // Update filters
326  for (x = m_filters.begin(); x != m_filters.end(); ++x)
327  {
328  (*x)->Update(m_fields, m_time);
329  }
330 
331  // Write out checkpoint files
332  if ((m_checksteps && step && !((step + 1) % m_checksteps)) ||
333  doCheckTime)
334  {
336  {
337  vector<bool> transformed(nfields, false);
338  for(i = 0; i < nfields; i++)
339  {
340  if (m_fields[i]->GetWaveSpace())
341  {
342  m_fields[i]->SetWaveSpace(false);
343  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
344  m_fields[i]->UpdatePhys());
345  m_fields[i]->SetPhysState(true);
346  transformed[i] = true;
347  }
348  }
349  Checkpoint_Output(nchk++);
350  for(i = 0; i < nfields; i++)
351  {
352  if (transformed[i])
353  {
354  m_fields[i]->SetWaveSpace(true);
355  m_fields[i]->HomogeneousFwdTrans(
356  m_fields[i]->GetPhys(),
357  m_fields[i]->UpdatePhys());
358  m_fields[i]->SetPhysState(false);
359  }
360  }
361  }
362  else
363  {
364  Checkpoint_Output(nchk++);
365  }
366  doCheckTime = false;
367  }
368 
369  // Step advance
370  ++step;
371  }
372 
373  // Print out summary statistics
374  if (m_session->GetComm()->GetRank() == 0)
375  {
376  if (m_cflSafetyFactor > 0.0)
377  {
378  cout << "CFL safety factor : " << m_cflSafetyFactor << endl
379  << "CFL time-step : " << m_timestep << endl;
380  }
381 
382  if (m_session->GetSolverInfo("Driver") != "SteadyState")
383  {
384  cout << "Time-integration : " << intTime << "s" << endl;
385  }
386  }
387 
388  // If homogeneous, transform back into physical space if necessary.
390  {
391  for(i = 0; i < nfields; i++)
392  {
393  if (m_fields[i]->GetWaveSpace())
394  {
395  m_fields[i]->SetWaveSpace(false);
396  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
397  m_fields[i]->UpdatePhys());
398  m_fields[i]->SetPhysState(true);
399  }
400  }
401  }
402  else
403  {
404  for(i = 0; i < nvariables; ++i)
405  {
406  m_fields[m_intVariables[i]]->SetPhys(fields[i]);
407  m_fields[m_intVariables[i]]->SetPhysState(true);
408  }
409  }
410 
411  // Finalise filters
412  for (x = m_filters.begin(); x != m_filters.end(); ++x)
413  {
414  (*x)->Finalise(m_fields, m_time);
415  }
416 
417  // Print for 1D problems
418  if(m_spacedim == 1)
419  {
420  v_AppendOutput1D(fields);
421  }
422  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
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 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 438 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().

439  {
441  AddSummaryItem(s, "Advection",
442  m_explicitAdvection ? "explicit" : "implicit");
443  AddSummaryItem(s, "Diffusion",
444  m_explicitDiffusion ? "explicit" : "implicit");
445 
446  if (m_session->GetSolverInfo("EQTYPE")
447  == "SteadyAdvectionDiffusionReaction")
448  {
449  AddSummaryItem(s, "Reaction",
450  m_explicitReaction ? "explicit" : "implicit");
451  }
452 
453  AddSummaryItem(s, "Time Step", m_timestep);
454  AddSummaryItem(s, "No. of Steps", m_steps);
455  AddSummaryItem(s, "Checkpoints (steps)", m_checksteps);
456  AddSummaryItem(s, "Integration Type",
458  m_intScheme->GetIntegrationMethod()]);
459  }
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 880 of file UnsteadySystem.cpp.

References ASSERTL0.

Referenced by GetTimeStep().

882  {
883  ASSERTL0(false, "Not defined for this class");
884  return 0.0;
885  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
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 486 of file UnsteadySystem.cpp.

References ASSERTL0.

489  {
490  ASSERTL0(false,
491  "This function is not implemented for this equation.");
492  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
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 494 of file UnsteadySystem.cpp.

References ASSERTL0.

498  {
499  ASSERTL0(false,
500  "This function is not implemented for this equation.");
501  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
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 503 of file UnsteadySystem.cpp.

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

506  {
507  int i, j;
508  int nTraceNumPoints = GetTraceNpoints();
509  int nvariables = m_fields.num_elements();
510  int nqvar = uflux.num_elements();
511 
512  Array<OneD, NekDouble > Fwd (nTraceNumPoints);
513  Array<OneD, NekDouble > Bwd (nTraceNumPoints);
514  Array<OneD, NekDouble > Vn (nTraceNumPoints, 0.0);
515  Array<OneD, NekDouble > fluxtemp(nTraceNumPoints, 0.0);
516 
517  // Get the sign of (v \cdot n), v = an arbitrary vector
518 
519  // Evaulate upwind flux:
520  // uflux = \hat{u} \phi \cdot u = u^{(+,-)} n
521  for (j = 0; j < nqvar; ++j)
522  {
523  for (i = 0; i < nvariables ; ++i)
524  {
525  // Compute Fwd and Bwd value of ufield of i direction
526  m_fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
527 
528  // if Vn >= 0, flux = uFwd, i.e.,
529  // edge::eForward, if V*n>=0 <=> V*n_F>=0, pick uflux = uFwd
530  // edge::eBackward, if V*n>=0 <=> V*n_B<0, pick uflux = uFwd
531 
532  // else if Vn < 0, flux = uBwd, i.e.,
533  // edge::eForward, if V*n<0 <=> V*n_F<0, pick uflux = uBwd
534  // edge::eBackward, if V*n<0 <=> V*n_B>=0, pick uflux = uBwd
535 
536  m_fields[i]->GetTrace()->Upwind(m_traceNormals[j],
537  Fwd, Bwd, fluxtemp);
538 
539  // Imposing weak boundary condition with flux
540  // if Vn >= 0, uflux = uBwd at Neumann, i.e.,
541  // edge::eForward, if V*n>=0 <=> V*n_F>=0, pick uflux = uBwd
542  // edge::eBackward, if V*n>=0 <=> V*n_B<0, pick uflux = uBwd
543 
544  // if Vn >= 0, uflux = uFwd at Neumann, i.e.,
545  // edge::eForward, if V*n<0 <=> V*n_F<0, pick uflux = uFwd
546  // edge::eBackward, if V*n<0 <=> V*n_B>=0, pick uflux = uFwd
547 
548  if(m_fields[0]->GetBndCondExpansions().num_elements())
549  {
550  WeakPenaltyforScalar(i, ufield[i], fluxtemp);
551  }
552 
553  // if Vn >= 0, flux = uFwd*(tan_{\xi}^- \cdot \vec{n}),
554  // i.e,
555  // edge::eForward, uFwd \(\tan_{\xi}^Fwd \cdot \vec{n})
556  // edge::eBackward, uFwd \(\tan_{\xi}^Bwd \cdot \vec{n})
557 
558  // else if Vn < 0, flux = uBwd*(tan_{\xi}^- \cdot \vec{n}),
559  // i.e,
560  // edge::eForward, uBwd \(\tan_{\xi}^Fwd \cdot \vec{n})
561  // edge::eBackward, uBwd \(\tan_{\xi}^Bwd \cdot \vec{n})
562 
563  Vmath::Vmul(nTraceNumPoints,
564  m_traceNormals[j], 1,
565  fluxtemp, 1,
566  uflux[j][i], 1);
567  }
568  }
569  }
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 573 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().

577  {
578  int nTraceNumPoints = GetTraceNpoints();
579  int nvariables = m_fields.num_elements();
580  int nqvar = qfield.num_elements();
581 
582  NekDouble C11 = 1.0;
583  Array<OneD, NekDouble > Fwd(nTraceNumPoints);
584  Array<OneD, NekDouble > Bwd(nTraceNumPoints);
585  Array<OneD, NekDouble > Vn (nTraceNumPoints, 0.0);
586 
587  Array<OneD, NekDouble > qFwd (nTraceNumPoints);
588  Array<OneD, NekDouble > qBwd (nTraceNumPoints);
589  Array<OneD, NekDouble > qfluxtemp(nTraceNumPoints, 0.0);
590 
591  Array<OneD, NekDouble > uterm(nTraceNumPoints);
592 
593  // Evaulate upwind flux:
594  // qflux = \hat{q} \cdot u = q \cdot n - C_(11)*(u^+ - u^-)
595  for (int i = 0; i < nvariables; ++i)
596  {
597  qflux[i] = Array<OneD, NekDouble> (nTraceNumPoints, 0.0);
598  for (int j = 0; j < nqvar; ++j)
599  {
600  // Compute Fwd and Bwd value of ufield of jth direction
601  m_fields[i]->GetFwdBwdTracePhys(qfield[j][i],qFwd,qBwd);
602 
603  // if Vn >= 0, flux = uFwd, i.e.,
604  // edge::eForward, if V*n>=0 <=> V*n_F>=0, pick
605  // qflux = qBwd = q+
606  // edge::eBackward, if V*n>=0 <=> V*n_B<0, pick
607  // qflux = qBwd = q-
608 
609  // else if Vn < 0, flux = uBwd, i.e.,
610  // edge::eForward, if V*n<0 <=> V*n_F<0, pick
611  // qflux = qFwd = q-
612  // edge::eBackward, if V*n<0 <=> V*n_B>=0, pick
613  // qflux = qFwd = q+
614 
615  m_fields[i]->GetTrace()->Upwind(m_traceNormals[j],
616  qBwd, qFwd,
617  qfluxtemp);
618 
619  Vmath::Vmul(nTraceNumPoints,
620  m_traceNormals[j], 1,
621  qfluxtemp, 1,
622  qfluxtemp, 1);
623 
624  // Generate Stability term = - C11 ( u- - u+ )
625  m_fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
626 
627  Vmath::Vsub(nTraceNumPoints,
628  Fwd, 1, Bwd, 1,
629  uterm, 1);
630 
631  Vmath::Smul(nTraceNumPoints,
632  -1.0 * C11, uterm, 1,
633  uterm, 1);
634 
635  // Flux = {Fwd, Bwd} * (nx, ny, nz) + uterm * (nx, ny)
636  Vmath::Vadd(nTraceNumPoints,
637  uterm, 1,
638  qfluxtemp, 1,
639  qfluxtemp, 1);
640 
641  // Imposing weak boundary condition with flux
642  if (m_fields[0]->GetBndCondExpansions().num_elements())
643  {
644  WeakPenaltyforVector(i, j,
645  qfield[j][i],
646  qfluxtemp,
647  C11);
648  }
649 
650  // q_hat \cdot n = (q_xi \cdot n_xi) or (q_eta \cdot n_eta)
651  // n_xi = n_x * tan_xi_x + n_y * tan_xi_y + n_z * tan_xi_z
652  // n_xi = n_x * tan_eta_x + n_y * tan_eta_y + n_z*tan_eta_z
653  Vmath::Vadd(nTraceNumPoints,
654  qfluxtemp, 1,
655  qflux[i], 1,
656  qflux[i], 1);
657  }
658  }
659  }
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 892 of file UnsteadySystem.cpp.

Referenced by v_DoSolve().

893  {
894  return false;
895  }
bool Nektar::SolverUtils::UnsteadySystem::v_PreIntegrate ( int  step)
protectedvirtual

Reimplemented in Nektar::IncNavierStokes.

Definition at line 887 of file UnsteadySystem.cpp.

Referenced by v_DoSolve().

888  {
889  return false;
890  }
bool Nektar::SolverUtils::UnsteadySystem::v_SteadyStateCheck ( int  step)
protectedvirtual

Definition at line 897 of file UnsteadySystem.cpp.

898  {
899  return false;
900  }
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 709 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().

714  {
715  int i, e, npoints, id1, id2;
716 
717  // Number of boundary regions
718  int nbnd = m_fields[var]->GetBndCondExpansions().num_elements();
719  int Nfps, numBDEdge;
720  int nTraceNumPoints = GetTraceNpoints();
721  int cnt = 0;
722 
723  Array<OneD, NekDouble > uplus(nTraceNumPoints);
724 
725  m_fields[var]->ExtractTracePhys(physfield, uplus);
726  for (i = 0; i < nbnd; ++i)
727  {
728  // Number of boundary expansion related to that region
729  numBDEdge = m_fields[var]->
730  GetBndCondExpansions()[i]->GetExpSize();
731 
732  // Evaluate boundary values g_D or g_N from input files
734  m_session->GetFunction("InitialConditions", 0);
735 
736  npoints = m_fields[var]->
737  GetBndCondExpansions()[i]->GetNpoints();
738 
739  Array<OneD,NekDouble> BDphysics(npoints);
740  Array<OneD,NekDouble> x0(npoints,0.0);
741  Array<OneD,NekDouble> x1(npoints,0.0);
742  Array<OneD,NekDouble> x2(npoints,0.0);
743 
744  m_fields[var]->GetBndCondExpansions()[i]->GetCoords(x0,x1,x2);
745  ifunc->Evaluate(x0,x1,x2,time,BDphysics);
746 
747  // Weakly impose boundary conditions by modifying flux values
748  for (e = 0; e < numBDEdge ; ++e)
749  {
750  // Number of points on the expansion
751  Nfps = m_fields[var]->
752  GetBndCondExpansions()[i]->GetExp(e)->GetNumPoints(0);
753 
754  id1 = m_fields[var]->
755  GetBndCondExpansions()[i]->GetPhys_Offset(e);
756 
757  id2 = m_fields[0]->GetTrace()->
758  GetPhys_Offset(m_fields[0]->GetTraceMap()->
759  GetBndCondTraceToGlobalTraceMap(cnt++));
760 
761  // For Dirichlet boundary condition: uflux = g_D
762  if (m_fields[var]->GetBndConditions()[i]->
763  GetBoundaryConditionType() == SpatialDomains::eDirichlet)
764  {
765  Vmath::Vcopy(Nfps,
766  &BDphysics[id1], 1,
767  &penaltyflux[id2], 1);
768  }
769  // For Neumann boundary condition: uflux = u+
770  else if ((m_fields[var]->GetBndConditions()[i])->
771  GetBoundaryConditionType() == SpatialDomains::eNeumann)
772  {
773  Vmath::Vcopy(Nfps,
774  &uplus[id2], 1,
775  &penaltyflux[id2], 1);
776  }
777  }
778  }
779  }
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:1016
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 786 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().

793  {
794  int i, e, npoints, id1, id2;
795  int nbnd = m_fields[var]->GetBndCondExpansions().num_elements();
796  int numBDEdge, Nfps;
797  int nTraceNumPoints = GetTraceNpoints();
798  Array<OneD, NekDouble > uterm(nTraceNumPoints);
799  Array<OneD, NekDouble > qtemp(nTraceNumPoints);
800  int cnt = 0;
801 
802  m_fields[var]->ExtractTracePhys(physfield,qtemp);
803 
804  for (i = 0; i < nbnd; ++i)
805  {
806  numBDEdge = m_fields[var]->
807  GetBndCondExpansions()[i]->GetExpSize();
808 
809  // Evaluate boundary values g_D or g_N from input files
811  m_session->GetFunction("InitialConditions", 0);
812 
813  npoints = m_fields[var]->
814  GetBndCondExpansions()[i]->GetNpoints();
815 
816  Array<OneD,NekDouble> BDphysics(npoints);
817  Array<OneD,NekDouble> x0(npoints,0.0);
818  Array<OneD,NekDouble> x1(npoints,0.0);
819  Array<OneD,NekDouble> x2(npoints,0.0);
820 
821  m_fields[var]->GetBndCondExpansions()[i]->GetCoords(x0,x1,x2);
822  ifunc->Evaluate(x0,x1,x2,time,BDphysics);
823 
824  // Weakly impose boundary conditions by modifying flux values
825  for (e = 0; e < numBDEdge ; ++e)
826  {
827  Nfps = m_fields[var]->
828  GetBndCondExpansions()[i]->GetExp(e)->GetNumPoints(0);
829 
830  id1 = m_fields[var]->
831  GetBndCondExpansions()[i]->GetPhys_Offset(e);
832 
833  id2 = m_fields[0]->GetTrace()->
834  GetPhys_Offset(m_fields[0]->GetTraceMap()->
835  GetBndCondTraceToGlobalTraceMap(cnt++));
836 
837  // For Dirichlet boundary condition:
838  //qflux = q+ - C_11 (u+ - g_D) (nx, ny)
839  if(m_fields[var]->GetBndConditions()[i]->
840  GetBoundaryConditionType() == SpatialDomains::eDirichlet)
841  {
842  Vmath::Vmul(Nfps,
843  &m_traceNormals[dir][id2], 1,
844  &qtemp[id2], 1,
845  &penaltyflux[id2], 1);
846  }
847  // For Neumann boundary condition: qflux = g_N
848  else if((m_fields[var]->GetBndConditions()[i])->
849  GetBoundaryConditionType() == SpatialDomains::eNeumann)
850  {
851  Vmath::Vmul(Nfps,
852  &m_traceNormals[dir][id2], 1,
853  &BDphysics[id1], 1,
854  &penaltyflux[id2], 1);
855  }
856  }
857  }
858  }
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