Nektar++
Public Member Functions | Static 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:
[legend]

Public Member Functions

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

Static Public Attributes

static std::string cmdSetStartTime
 
static std::string cmdSetStartChkNum
 

Protected Member Functions

SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareField=true) override
 Init object for UnsteadySystem class. More...
 
SOLVER_UTILS_EXPORT void v_DoSolve () override
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintStatusInformation (const int step, const NekDouble cpuTime)
 Print Status Information. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintSummaryStatistics (const NekDouble intTime)
 Print Summary Statistics. More...
 
SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true) override
 Sets up initial conditions. More...
 
SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &s) override
 Print a summary of time stepping parameters. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Return the timestep to be used for the next step in the time-marching loop. More...
 
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_RequireFwdTrans ()
 
virtual SOLVER_UTILS_EXPORT void v_SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
virtual SOLVER_UTILS_EXPORT bool v_UpdateTimeStepCheck ()
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time, int &nchk)
 
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff (const Array< OneD, Array< OneD, NekDouble > > vel, StdRegions::VarCoeffMap &varCoeffMap)
 Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h time velocity. More...
 
SOLVER_UTILS_EXPORT void DoDummyProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Perform dummy projection. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises EquationSystem class members. More...
 
virtual SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareFeld=true)
 Initialisation object for EquationSystem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true)
 Virtual function for initialisation implementation. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Virtual function for solve implementation. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space. More...
 
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &l)
 Virtual function for generating summary information. More...
 
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
virtual SOLVER_UTILS_EXPORT void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure (void)
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp (void)
 Virtual function to identify if operator is negated in DoSolve. More...
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 

Protected Attributes

LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::vector< int > m_intVariables
 
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 Maximun cfl in cfl growth. More...
 
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
std::ofstream m_errFile
 
NekDouble m_epsilon
 Diffusion coefficient. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
bool m_verbose
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh. More...
 
std::string m_sessionName
 Name of the session. More...
 
NekDouble m_time
 Current time of simulation. More...
 
int m_initialStep
 Number of the step where the simulation should begin. More...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
NekDouble m_lastCheckTime
 
NekDouble m_TimeIncrementFactor
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_infosteps
 Number of time steps between outputting status information. More...
 
int m_iterPIT = 0
 Number of parallel-in-time time iteration. More...
 
int m_windowPIT = 0
 Index of windows for parallel-in-time time iteration. More...
 
int m_spacedim
 Spatial dimension (>= expansion dim). More...
 
int m_expdim
 Expansion dimension. More...
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used. More...
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used. More...
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used. More...
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform. More...
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations. More...
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation. More...
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction. More...
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity. More...
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields. More...
 
Array< OneD, NekDoublem_movingFrameVelsxyz
 Moving frame of reference velocities (u, v, w, omega_x, omega_y, omega_z, a_x, a_y, a_z, domega_x, domega_y, domega_z) More...
 
Array< OneD, NekDoublem_movingFrameData
 Moving frame of reference angles with respect to the. More...
 
boost::numeric::ublas::matrix< NekDoublem_movingFrameProjMat
 Projection matrix for transformation between inertial and moving. More...
 
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error. More...
 
enum HomogeneousType m_HomogeneousType
 
NekDouble m_LhomX
 physical length in X direction (if homogeneous) More...
 
NekDouble m_LhomY
 physical length in Y direction (if homogeneous) More...
 
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous) More...
 
int m_npointsX
 number of points in X direction (if homogeneous) More...
 
int m_npointsY
 number of points in Y direction (if homogeneous) More...
 
int m_npointsZ
 number of points in Z direction (if homogeneous) More...
 
int m_HomoDirec
 number of homogenous directions More...
 

Private Member Functions

SOLVER_UTILS_EXPORT void AppendOutput1D (void)
 Print the solution at each solution point in a txt file. More...
 
void InitializeSteadyState ()
 
bool CheckSteadyState (int step, const NekDouble &totCPUTime=0.0)
 Calculate whether the system has reached a steady state by observing residuals to a user-defined tolerance. More...
 

Additional Inherited Members

- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D , eHomogeneous2D , eHomogeneous3D , eNotHomogeneous }
 Parameter for homogeneous expansions. More...
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 
static std::string projectionTypeLookupIds []
 

Detailed Description

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 45 of file UnsteadySystem.h.

Constructor & Destructor Documentation

◆ ~UnsteadySystem()

Nektar::SolverUtils::UnsteadySystem::~UnsteadySystem ( )
override

Destructor.

Destructor for the class UnsteadyAdvection.

Definition at line 176 of file UnsteadySystem.cpp.

177{
178}

◆ UnsteadySystem()

Nektar::SolverUtils::UnsteadySystem::UnsteadySystem ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
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 74 of file UnsteadySystem.cpp.

77 : EquationSystem(pSession, pGraph)
78
79{
80}
SOLVER_UTILS_EXPORT EquationSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises EquationSystem class members.

Member Function Documentation

◆ AppendOutput1D()

void Nektar::SolverUtils::UnsteadySystem::AppendOutput1D ( void  )
private

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 635 of file UnsteadySystem.cpp.

636{
637 // Coordinates of the quadrature points in the real physical space.
638 Array<OneD, NekDouble> x(GetNpoints());
639 Array<OneD, NekDouble> y(GetNpoints());
640 Array<OneD, NekDouble> z(GetNpoints());
641 m_fields[0]->GetCoords(x, y, z);
642
643 // Print out the solution in a txt file.
644 ofstream outfile;
645 outfile.open("solution1D.txt");
646 for (int i = 0; i < GetNpoints(); i++)
647 {
648 outfile << scientific << setw(17) << setprecision(16) << x[i] << " "
649 << m_fields[m_intVariables[0]]->GetPhys()[i] << endl;
650 }
651 outfile << endl << endl;
652 outfile.close();
653}
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetNpoints()
std::vector< double > z(NPUPPER)

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

Referenced by v_DoSolve().

◆ CheckForRestartTime()

void Nektar::SolverUtils::UnsteadySystem::CheckForRestartTime ( NekDouble time,
int &  nchk 
)
protected

Definition at line 658 of file UnsteadySystem.cpp.

659{
660 if (m_session->DefinesFunction("InitialConditions"))
661 {
662 for (int i = 0; i < m_fields.size(); ++i)
663 {
665
666 vType = m_session->GetFunctionType("InitialConditions",
667 m_session->GetVariable(i));
668
670 {
671 std::string filename = m_session->GetFunctionFilename(
672 "InitialConditions", m_session->GetVariable(i));
673
674 fs::path pfilename(filename);
675
676 // Redefine path for parallel file which is in directory.
677 if (fs::is_directory(pfilename))
678 {
679 fs::path metafile("Info.xml");
680 fs::path fullpath = pfilename / metafile;
681 filename = LibUtilities::PortablePath(fullpath);
682 }
685 fld->ImportFieldMetaData(filename, m_fieldMetaDataMap);
686
687 // Check to see if time defined.
689 {
690 auto iter = m_fieldMetaDataMap.find("Time");
691 if (iter != m_fieldMetaDataMap.end())
692 {
693 time = boost::lexical_cast<NekDouble>(iter->second);
694 }
695
696 iter = m_fieldMetaDataMap.find("ChkFileNum");
697 if (iter != m_fieldMetaDataMap.end())
698 {
699 nchk = boost::lexical_cast<NekDouble>(iter->second);
700 }
701 }
702
703 break;
704 }
705 }
706 }
707 if (m_session->DefinesCmdLineArgument("set-start-time"))
708 {
709 time = boost::lexical_cast<NekDouble>(
710 m_session->GetCmdLineArgument<std::string>("set-start-time")
711 .c_str());
712 }
713 if (m_session->DefinesCmdLineArgument("set-start-chknumber"))
714 {
715 nchk = boost::lexical_cast<int>(
716 m_session->GetCmdLineArgument<std::string>("set-start-chknumber"));
717 }
718 ASSERTL0(time >= 0 && nchk >= 0,
719 "Starting time and checkpoint number should be >= 0");
720}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
static std::shared_ptr< FieldIO > CreateForFile(const LibUtilities::SessionReaderSharedPtr session, const std::string &filename)
Construct a FieldIO object for a given input filename.
Definition: FieldIO.cpp:224
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:322
static std::string PortablePath(const fs::path &path)
create portable path on different platforms for std::filesystem path.
Definition: Filesystem.hpp:56
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:51

References ASSERTL0, Nektar::LibUtilities::FieldIO::CreateForFile(), Nektar::LibUtilities::eFunctionTypeFile, Nektar::SolverUtils::EquationSystem::m_fieldMetaDataMap, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_session, Nektar::LibUtilities::NullFieldMetaDataMap, and Nektar::LibUtilities::PortablePath().

Referenced by v_DoInitialise().

◆ CheckSteadyState()

bool Nektar::SolverUtils::UnsteadySystem::CheckSteadyState ( int  step,
const NekDouble totCPUTime = 0.0 
)
private

Calculate whether the system has reached a steady state by observing residuals to a user-defined tolerance.

Definition at line 853 of file UnsteadySystem.cpp.

854{
855 const int nPoints = GetTotPoints();
856 const int nFields = m_fields.size();
857
858 // Holds L2 errors.
859 Array<OneD, NekDouble> L2(nFields);
860
861 SteadyStateResidual(step, L2);
862
863 if (m_infosteps && m_comm->GetRank() == 0 &&
864 (((step + 1) % m_infosteps == 0) || ((step == m_initialStep))))
865 {
866 // Output time.
867 m_errFile << boost::format("%25.19e") % m_time;
868
869 m_errFile << " " << boost::format("%25.19e") % totCPUTime;
870
871 int stepp = step + 1;
872
873 m_errFile << " " << boost::format("%25.19e") % stepp;
874
875 // Output residuals.
876 for (int i = 0; i < nFields; ++i)
877 {
878 m_errFile << " " << boost::format("%25.19e") % L2[i];
879 }
880
881 m_errFile << endl;
882 }
883
884 // Calculate maximum L2 error.
885 NekDouble maxL2 = Vmath::Vmax(nFields, L2, 1);
886
887 if (m_infosteps && m_session->DefinesCmdLineArgument("verbose") &&
888 m_comm->GetRank() == 0 && ((step + 1) % m_infosteps == 0))
889 {
890 cout << "-- Maximum L^2 residual: " << maxL2 << endl;
891 }
892
893 if (maxL2 <= m_steadyStateTol)
894 {
895 return true;
896 }
897
898 for (int i = 0; i < m_fields.size(); ++i)
899 {
900 Vmath::Vcopy(nPoints, m_fields[i]->GetPhys(), 1, m_previousSolution[i],
901 1);
902 }
903
904 return false;
905}
LibUtilities::CommSharedPtr m_comm
Communicator.
int m_infosteps
Number of time steps between outputting status information.
NekDouble m_time
Current time of simulation.
int m_initialStep
Number of the step where the simulation should begin.
SOLVER_UTILS_EXPORT int GetTotPoints()
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
Storage for previous solution for steady-state check.
NekDouble m_steadyStateTol
Tolerance to which steady state should be evaluated at.
SOLVER_UTILS_EXPORT void SteadyStateResidual(int step, Array< OneD, NekDouble > &L2)
double NekDouble
T Vmax(int n, const T *x, const int incx)
Return the maximum element in x – called vmax to avoid conflict with max.
Definition: Vmath.hpp:644
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References CellMLToNektar.pycml::format, Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::SolverUtils::EquationSystem::m_comm, m_errFile, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_infosteps, Nektar::SolverUtils::EquationSystem::m_initialStep, m_previousSolution, Nektar::SolverUtils::EquationSystem::m_session, m_steadyStateTol, Nektar::SolverUtils::EquationSystem::m_time, SteadyStateResidual(), Vmath::Vcopy(), and Vmath::Vmax().

Referenced by v_DoSolve().

◆ DoDummyProjection()

void Nektar::SolverUtils::UnsteadySystem::DoDummyProjection ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Perform dummy projection.

Definition at line 949 of file UnsteadySystem.cpp.

953{
954
955 if (&inarray != &outarray)
956 {
957 for (int i = 0; i < inarray.size(); ++i)
958 {
959 Vmath::Vcopy(GetNpoints(), inarray[i], 1, outarray[i], 1);
960 }
961 }
962}

References Nektar::SolverUtils::EquationSystem::GetNpoints(), and Vmath::Vcopy().

Referenced by Nektar::Bidomain::v_InitObject(), Nektar::BidomainRoth::v_InitObject(), and Nektar::Monodomain::v_InitObject().

◆ GetTimeIntegrationScheme()

LibUtilities::TimeIntegrationSchemeSharedPtr & Nektar::SolverUtils::UnsteadySystem::GetTimeIntegrationScheme ( )

Returns the time integration scheme.

Definition at line 191 of file UnsteadySystem.cpp.

193{
194 return m_intScheme;
195}
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
Wrapper to the time integration scheme.

References m_intScheme.

◆ GetTimeIntegrationSchemeOperators()

LibUtilities::TimeIntegrationSchemeOperators & Nektar::SolverUtils::UnsteadySystem::GetTimeIntegrationSchemeOperators ( )

Returns the time integration scheme operators.

Definition at line 200 of file UnsteadySystem.cpp.

202{
203 return m_ode;
204}
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.

References m_ode.

◆ GetTimeStep() [1/2]

SOLVER_UTILS_EXPORT NekDouble Nektar::SolverUtils::UnsteadySystem::GetTimeStep ( void  )
inline

Definition at line 58 of file UnsteadySystem.h.

59 {
61 }
SOLVER_UTILS_EXPORT NekDouble GetTimeStep()

References Nektar::SolverUtils::EquationSystem::GetTimeStep().

Referenced by v_DoSolve().

◆ GetTimeStep() [2/2]

SOLVER_UTILS_EXPORT NekDouble Nektar::SolverUtils::UnsteadySystem::GetTimeStep ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray)
inline

Calculate the larger time-step mantaining the problem stable.

Definition at line 53 of file UnsteadySystem.h.

54 {
55 return v_GetTimeStep(inarray);
56 }
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.

References v_GetTimeStep().

◆ InitializeSteadyState()

void Nektar::SolverUtils::UnsteadySystem::InitializeSteadyState ( )
private

Definition at line 813 of file UnsteadySystem.cpp.

814{
815 if (m_steadyStateTol > 0.0)
816 {
817 const int nPoints = m_fields[0]->GetTotPoints();
819 Array<OneD, Array<OneD, NekDouble>>(m_fields.size());
820
821 for (int i = 0; i < m_fields.size(); ++i)
822 {
823 m_previousSolution[i] = Array<OneD, NekDouble>(nPoints);
824 Vmath::Vcopy(nPoints, m_fields[i]->GetPhys(), 1,
825 m_previousSolution[i], 1);
826 }
827
828 if (m_comm->GetRank() == 0)
829 {
830 std::string fName =
831 m_session->GetSessionName() + std::string(".resd");
832 m_errFile.open(fName.c_str());
833 m_errFile << setw(26) << left << "# Time";
834
835 m_errFile << setw(26) << left << "CPU_Time";
836
837 m_errFile << setw(26) << left << "Step";
838
839 for (int i = 0; i < m_fields.size(); ++i)
840 {
841 m_errFile << setw(26) << m_session->GetVariables()[i];
842 }
843
844 m_errFile << endl;
845 }
846 }
847}

References Nektar::SolverUtils::EquationSystem::m_comm, m_errFile, Nektar::SolverUtils::EquationSystem::m_fields, m_previousSolution, Nektar::SolverUtils::EquationSystem::m_session, m_steadyStateTol, and Vmath::Vcopy().

Referenced by v_DoInitialise().

◆ MaxTimeStepEstimator()

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 183 of file UnsteadySystem.cpp.

184{
185 return m_intScheme->GetTimeStability();
186}

References m_intScheme.

Referenced by Nektar::CompressibleFlowSystem::GetElmtTimeStep().

◆ SetTimeStep()

SOLVER_UTILS_EXPORT void Nektar::SolverUtils::UnsteadySystem::SetTimeStep ( const NekDouble  timestep)
inline

Definition at line 63 of file UnsteadySystem.h.

64 {
66 }
SOLVER_UTILS_EXPORT void SetTimeStep(const NekDouble timestep)

References Nektar::SolverUtils::EquationSystem::SetTimeStep().

◆ SteadyStateResidual()

SOLVER_UTILS_EXPORT void Nektar::SolverUtils::UnsteadySystem::SteadyStateResidual ( int  step,
Array< OneD, NekDouble > &  L2 
)
inline

Definition at line 68 of file UnsteadySystem.h.

70 {
71 v_SteadyStateResidual(step, L2);
72 }
virtual SOLVER_UTILS_EXPORT void v_SteadyStateResidual(int step, Array< OneD, NekDouble > &L2)

References v_SteadyStateResidual().

Referenced by CheckSteadyState().

◆ SVVVarDiffCoeff()

void Nektar::SolverUtils::UnsteadySystem::SVVVarDiffCoeff ( const Array< OneD, Array< OneD, NekDouble > >  vel,
StdRegions::VarCoeffMap varCoeffMap 
)
protected

Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h time velocity.

Definition at line 770 of file UnsteadySystem.cpp.

773{
774 int phystot = m_fields[0]->GetTotPoints();
775 int nvel = vel.size();
776
777 Array<OneD, NekDouble> varcoeff(phystot), tmp;
778
779 // Calculate magnitude of v.
780 Vmath::Vmul(phystot, vel[0], 1, vel[0], 1, varcoeff, 1);
781 for (int n = 1; n < nvel; ++n)
782 {
783 Vmath::Vvtvp(phystot, vel[n], 1, vel[n], 1, varcoeff, 1, varcoeff, 1);
784 }
785 Vmath::Vsqrt(phystot, varcoeff, 1, varcoeff, 1);
786
787 for (int i = 0; i < m_fields[0]->GetNumElmts(); ++i)
788 {
789 int offset = m_fields[0]->GetPhys_Offset(i);
790 int nq = m_fields[0]->GetExp(i)->GetTotPoints();
791 Array<OneD, NekDouble> unit(nq, 1.0);
792
793 int nmodes = 0;
794
795 for (int n = 0; n < m_fields[0]->GetExp(i)->GetNumBases(); ++n)
796 {
797 nmodes = max(nmodes, m_fields[0]->GetExp(i)->GetBasisNumModes(n));
798 }
799
800 NekDouble h = m_fields[0]->GetExp(i)->Integral(unit);
801 h = pow(h, (NekDouble)(1.0 / nvel)) / ((NekDouble)nmodes);
802
803 Vmath::Smul(nq, h, varcoeff + offset, 1, tmp = varcoeff + offset, 1);
804 }
805
806 // Set up map with eVarCoffLaplacian key.
807 varCoeffMap[StdRegions::eVarCoeffLaplacian] = varcoeff;
808}
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.hpp:340
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.hpp:72
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100

References Nektar::StdRegions::eVarCoeffLaplacian, Nektar::SolverUtils::EquationSystem::m_fields, Vmath::Smul(), Vmath::Vmul(), Vmath::Vsqrt(), and Vmath::Vvtvp().

Referenced by Nektar::UnsteadyViscousBurgers::DoImplicitSolve().

◆ v_DoInitialise()

void Nektar::SolverUtils::UnsteadySystem::v_DoInitialise ( bool  dumpInitialConditions = true)
overrideprotectedvirtual

Sets up initial conditions.

Sets the initial conditions.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Reimplemented in Nektar::VCSImplicit, Nektar::SolverUtils::FileSolution, Nektar::PulseWaveSystem, Nektar::MMFSWE, Nektar::CoupledLinearNS, Nektar::VCSMapping, and Nektar::VelocityCorrectionScheme.

Definition at line 594 of file UnsteadySystem.cpp.

595{
598 SetInitialConditions(m_time, dumpInitialConditions);
600}
int m_nchk
Number of checkpoints written so far.
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.
SOLVER_UTILS_EXPORT void CheckForRestartTime(NekDouble &time, int &nchk)

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

Referenced by Nektar::VCSMapping::v_DoInitialise(), and Nektar::VelocityCorrectionScheme::v_DoInitialise().

◆ v_DoSolve()

void Nektar::SolverUtils::UnsteadySystem::v_DoSolve ( void  )
overrideprotectedvirtual

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::MMFAdvection, Nektar::CFSImplicit, Nektar::MMFMaxwell, Nektar::PulseWaveSystem, Nektar::MMFSWE, and Nektar::CoupledLinearNS.

Definition at line 210 of file UnsteadySystem.cpp.

211{
212 ASSERTL0(m_intScheme != nullptr, "No time integration scheme.");
213
214 int i = 1;
215 int nvariables = 0;
216 int nfields = m_fields.size();
217 if (m_intVariables.empty())
218 {
219 for (i = 0; i < nfields; ++i)
220 {
221 m_intVariables.push_back(i);
222 }
223 nvariables = nfields;
224 }
225 else
226 {
227 nvariables = m_intVariables.size();
228 }
229
230 // Integrate in wave-space if using homogeneous1D.
232 {
233 for (i = 0; i < nfields; ++i)
234 {
235 m_fields[i]->HomogeneousFwdTrans(m_fields[i]->GetTotPoints(),
236 m_fields[i]->GetPhys(),
237 m_fields[i]->UpdatePhys());
238 m_fields[i]->SetWaveSpace(true);
239 m_fields[i]->SetPhysState(false);
240 }
241 }
242
243 // Set up wrapper to fields data storage.
244 Array<OneD, Array<OneD, NekDouble>> fields(nvariables);
245
246 // Order storage to list time-integrated fields first.
247 for (i = 0; i < nvariables; ++i)
248 {
249 fields[i] = m_fields[m_intVariables[i]]->UpdatePhys();
250 m_fields[m_intVariables[i]]->SetPhysState(false);
251 }
252
253 // Initialise time integration scheme.
254 m_intScheme->InitializeScheme(m_timestep, fields, m_time, m_ode);
255
256 // Initialise filters.
257 for (auto &x : m_filters)
258 {
259 x.second->Initialise(m_fields, m_time);
260 }
261
262 LibUtilities::Timer timer;
263 bool doCheckTime = false;
264 int step = m_initialStep;
265 int stepCounter = 0;
266 NekDouble intTime = 0.0;
267 NekDouble cpuTime = 0.0;
268 NekDouble cpuPrevious = 0.0;
269 NekDouble elapsed = 0.0;
270 NekDouble totFilterTime = 0.0;
271
272 m_lastCheckTime = 0.0;
273
274 Array<OneD, int> abortFlags(2, 0);
275 string abortFile = "abort";
276 if (m_session->DefinesSolverInfo("CheckAbortFile"))
277 {
278 abortFile = m_session->GetSolverInfo("CheckAbortFile");
279 }
280
281 NekDouble tmp_cflSafetyFactor = m_cflSafetyFactor;
282 while ((step < m_steps || m_time < m_fintime - NekConstants::kNekZeroTol) &&
283 abortFlags[1] == 0)
284 {
286 {
287 tmp_cflSafetyFactor =
288 min(m_CFLEnd, m_CFLGrowth * tmp_cflSafetyFactor);
289 }
290
291 // Frozen preconditioner checks.
292 if (!m_comm->IsParallelInTime())
293 {
295 {
296 m_cflSafetyFactor = tmp_cflSafetyFactor;
297
299 {
300 m_timestep = GetTimeStep(fields);
301 }
302
303 // Ensure that the final timestep finishes at the final
304 // time, or at a prescribed IO_CheckTime.
305 if (m_time + m_timestep > m_fintime && m_fintime > 0.0)
306 {
308 }
309 else if (m_checktime &&
311 {
314 doCheckTime = true;
315 }
316 }
317 }
318
319 if (m_TimeIncrementFactor > 1.0)
320 {
321 NekDouble timeincrementFactor = m_TimeIncrementFactor;
322 m_timestep *= timeincrementFactor;
323
324 if (m_time + m_timestep > m_fintime && m_fintime > 0.0)
325 {
327 }
328 }
329
330 // Perform any solver-specific pre-integration steps.
331 timer.Start();
332 if (v_PreIntegrate(step))
333 {
334 break;
335 }
336
337 ASSERTL0(m_timestep > 0, "m_timestep < 0");
338
339 fields = m_intScheme->TimeIntegrate(stepCounter, m_timestep);
340 timer.Stop();
341
343 elapsed = timer.TimePerTest(1);
344 intTime += elapsed;
345 cpuTime += elapsed;
346
347 // Write out status information.
348 v_PrintStatusInformation(step, cpuTime);
349 if (m_infosteps &&
350 m_session->GetComm()->GetSpaceComm()->GetRank() == 0 &&
351 !((step + 1) % m_infosteps))
352 {
353 cpuPrevious = cpuTime;
354 cpuTime = 0.0;
355 }
356
357 // Transform data into coefficient space.
358 for (i = 0; i < nvariables; ++i)
359 {
360 // Copy fields into ExpList::m_phys.
361 m_fields[m_intVariables[i]]->SetPhys(fields[i]);
362 if (v_RequireFwdTrans())
363 {
364 m_fields[m_intVariables[i]]->FwdTransLocalElmt(
365 m_fields[m_intVariables[i]]->GetPhys(),
366 m_fields[m_intVariables[i]]->UpdateCoeffs());
367 }
368 m_fields[m_intVariables[i]]->SetPhysState(false);
369 }
370
371 // Perform any solver-specific post-integration steps.
372 if (v_PostIntegrate(step))
373 {
374 break;
375 }
376
377 // Check for steady-state.
379 (!((step + 1) % m_steadyStateSteps)))
380 {
381 if (CheckSteadyState(step, intTime))
382 {
383 if (m_comm->GetRank() == 0)
384 {
385 cout << "Reached Steady State to tolerance "
386 << m_steadyStateTol << endl;
387 }
388 break;
389 }
390 }
391
392 // Test for abort conditions (nan, or abort file).
393 if (m_abortSteps && !((step + 1) % m_abortSteps))
394 {
395 abortFlags[0] = 0;
396 for (i = 0; i < nvariables; ++i)
397 {
398 if (Vmath::Nnan(m_fields[m_intVariables[i]]->GetPhys().size(),
399 m_fields[m_intVariables[i]]->GetPhys(), 1) > 0)
400 {
401 abortFlags[0] = 1;
402 }
403 }
404
405 // Rank zero looks for abort file and deletes it
406 // if it exists. Communicates the abort.
407 if (m_session->GetComm()->GetSpaceComm()->GetRank() == 0)
408 {
409 if (fs::exists(abortFile))
410 {
411 fs::remove(abortFile);
412 abortFlags[1] = 1;
413 }
414 }
415
416 m_session->GetComm()->GetSpaceComm()->AllReduce(
417 abortFlags, LibUtilities::ReduceMax);
418
419 ASSERTL0(!abortFlags[0], "NaN found during time integration.");
420 }
421
422 // Update filters.
423 for (auto &x : m_filters)
424 {
425 timer.Start();
426 x.second->Update(m_fields, m_time);
427 timer.Stop();
428 elapsed = timer.TimePerTest(1);
429 totFilterTime += elapsed;
430
431 // Write out individual filter status information.
432 if (m_filtersInfosteps && m_session->GetComm()->GetRank() == 0 &&
433 !((step + 1) % m_filtersInfosteps) && !m_filters.empty() &&
434 m_session->DefinesCmdLineArgument("verbose"))
435 {
436 stringstream s0;
437 s0 << x.first << ":";
438 stringstream s1;
439 s1 << elapsed << "s";
440 stringstream s2;
441 s2 << elapsed / cpuPrevious * 100 << "%";
442 cout << "CPU time for filter " << setw(25) << left << s0.str()
443 << setw(12) << left << s1.str() << endl
444 << "\t Percentage of time integration: " << setw(10)
445 << left << s2.str() << endl;
446 }
447 }
448
449 // Write out overall filter status information.
450 if (m_filtersInfosteps && m_session->GetComm()->GetRank() == 0 &&
451 !((step + 1) % m_filtersInfosteps) && !m_filters.empty())
452 {
453 stringstream ss;
454 ss << totFilterTime << "s";
455 cout << "Total filters CPU Time:\t\t\t " << setw(10) << left
456 << ss.str() << endl;
457 }
458 totFilterTime = 0.0;
459
460 // Write out checkpoint files.
461 if ((m_checksteps && !((step + 1) % m_checksteps)) || doCheckTime)
462 {
464 {
465 // Transform to physical space for output.
466 vector<bool> transformed(nfields, false);
467 for (i = 0; i < nfields; i++)
468 {
469 if (m_fields[i]->GetWaveSpace())
470 {
471 m_fields[i]->SetWaveSpace(false);
472 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
473 m_fields[i]->UpdatePhys());
474 transformed[i] = true;
475 }
476 }
478 m_nchk++;
479
480 // Transform back to wave-space after output.
481 for (i = 0; i < nfields; i++)
482 {
483 if (transformed[i])
484 {
485 m_fields[i]->SetWaveSpace(true);
486 m_fields[i]->HomogeneousFwdTrans(
487 m_fields[i]->GetTotPoints(), m_fields[i]->GetPhys(),
488 m_fields[i]->UpdatePhys());
489 m_fields[i]->SetPhysState(false);
490 }
491 }
492 }
493 else
494 {
496 m_nchk++;
497 }
498 doCheckTime = false;
499 }
500
501 // Step advance.
502 ++step;
503 ++stepCounter;
504 }
505
506 // Print out summary statistics.
508
509 // If homogeneous, transform back into physical space if necessary.
511 {
512 for (i = 0; i < nfields; i++)
513 {
514 if (m_fields[i]->GetWaveSpace())
515 {
516 m_fields[i]->SetWaveSpace(false);
517 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
518 m_fields[i]->UpdatePhys());
519 }
520 }
521 }
522 else
523 {
524 for (i = 0; i < nvariables; ++i)
525 {
526 m_fields[m_intVariables[i]]->SetPhysState(true);
527 }
528 }
529
530 // Finalise filters.
531 for (auto &x : m_filters)
532 {
533 x.second->Finalise(m_fields, m_time);
534 }
535
536 // Print for 1D problems.
537 if (m_spacedim == 1)
538 {
540 }
541}
int m_spacedim
Spatial dimension (>= expansion dim).
NekDouble m_timestep
Time step size.
NekDouble m_fintime
Finish time of the simulation.
SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n)
Write checkpoint file of m_fields.
NekDouble m_checktime
Time between checkpoints.
enum HomogeneousType m_HomogeneousType
int m_steps
Number of steps to take.
int m_checksteps
Number of steps between checkpoints.
NekDouble m_CFLGrowth
CFL growth rate.
virtual SOLVER_UTILS_EXPORT bool v_UpdateTimeStepCheck()
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
NekDouble m_cflSafetyFactor
CFL safety factor (comprise between 0 to 1).
int m_abortSteps
Number of steps between checks for abort conditions.
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate(int step)
NekDouble m_CFLEnd
Maximun cfl in cfl growth.
virtual SOLVER_UTILS_EXPORT void v_PrintSummaryStatistics(const NekDouble intTime)
Print Summary Statistics.
virtual SOLVER_UTILS_EXPORT void v_PrintStatusInformation(const int step, const NekDouble cpuTime)
Print Status Information.
bool CheckSteadyState(int step, const NekDouble &totCPUTime=0.0)
Calculate whether the system has reached a steady state by observing residuals to a user-defined tole...
SOLVER_UTILS_EXPORT void AppendOutput1D(void)
Print the solution at each solution point in a txt file.
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate(int step)
virtual SOLVER_UTILS_EXPORT bool v_RequireFwdTrans()
int m_steadyStateSteps
Check for steady state at step interval.
int m_filtersInfosteps
Number of time steps between outputting filters information.
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
SOLVER_UTILS_EXPORT NekDouble GetTimeStep()
static const NekDouble kNekZeroTol
int Nnan(int n, const T *x, const int incx)
Return number of NaN elements of x.
Definition: Vmath.hpp:743

References AppendOutput1D(), ASSERTL0, Nektar::SolverUtils::EquationSystem::Checkpoint_Output(), CheckSteadyState(), Nektar::SolverUtils::EquationSystem::eNotHomogeneous, GetTimeStep(), Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::NekConstants::kNekZeroTol, m_abortSteps, m_CFLEnd, m_CFLGrowth, m_cflSafetyFactor, Nektar::SolverUtils::EquationSystem::m_checksteps, Nektar::SolverUtils::EquationSystem::m_checktime, Nektar::SolverUtils::EquationSystem::m_comm, Nektar::SolverUtils::EquationSystem::m_fields, m_filters, m_filtersInfosteps, Nektar::SolverUtils::EquationSystem::m_fintime, Nektar::SolverUtils::EquationSystem::m_HomogeneousType, m_homoInitialFwd, Nektar::SolverUtils::EquationSystem::m_infosteps, Nektar::SolverUtils::EquationSystem::m_initialStep, m_intScheme, m_intVariables, Nektar::SolverUtils::EquationSystem::m_lastCheckTime, Nektar::SolverUtils::EquationSystem::m_nchk, m_ode, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_spacedim, m_steadyStateSteps, m_steadyStateTol, Nektar::SolverUtils::EquationSystem::m_steps, Nektar::SolverUtils::EquationSystem::m_time, Nektar::SolverUtils::EquationSystem::m_TimeIncrementFactor, Nektar::SolverUtils::EquationSystem::m_timestep, Vmath::Nnan(), Nektar::LibUtilities::ReduceMax, Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), Nektar::LibUtilities::Timer::TimePerTest(), v_PostIntegrate(), v_PreIntegrate(), v_PrintStatusInformation(), v_PrintSummaryStatistics(), v_RequireFwdTrans(), and v_UpdateTimeStepCheck().

Referenced by Nektar::CFSImplicit::v_DoSolve(), and Nektar::CoupledLinearNS::v_DoSolve().

◆ v_GenerateSummary()

void Nektar::SolverUtils::UnsteadySystem::v_GenerateSummary ( SummaryList s)
overrideprotectedvirtual

Print a summary of time stepping parameters.

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

Reimplemented from Nektar::SolverUtils::EquationSystem.

Reimplemented in Nektar::MMFAdvection, Nektar::UnsteadyAdvection, Nektar::UnsteadyAdvectionDiffusion, Nektar::UnsteadyInviscidBurgers, Nektar::UnsteadyViscousBurgers, Nektar::CompressibleFlowSystem, Nektar::MMFDiffusion, Nektar::ImageWarpingSystem, Nektar::CoupledLinearNS, Nektar::SmoothedProfileMethod, Nektar::VelocityCorrectionScheme, Nektar::VCSImplicit, Nektar::VCSWeakPressure, Nektar::MMFMaxwell, Nektar::PulseWavePropagation, Nektar::LinearSWE, Nektar::MMFSWE, Nektar::NonlinearPeregrine, Nektar::NonlinearSWE, Nektar::ShallowWaterSystem, Nektar::CFLtester, Nektar::UnsteadyDiffusion, Nektar::Bidomain, Nektar::BidomainRoth, and Nektar::Monodomain.

Definition at line 606 of file UnsteadySystem.cpp.

607{
609 AddSummaryItem(s, "Advect. advancement",
610 m_explicitAdvection ? "explicit" : "implicit");
611
612 AddSummaryItem(s, "Diffuse. advancement",
613 m_explicitDiffusion ? "explicit" : "implicit");
614
615 if (m_session->GetSolverInfo("EQTYPE") ==
616 "SteadyAdvectionDiffusionReaction")
617 {
618 AddSummaryItem(s, "React. advancement",
619 m_explicitReaction ? "explicit" : "implicit");
620 }
621
622 AddSummaryItem(s, "Time Step", m_timestep);
623 AddSummaryItem(s, "No. of Steps", m_steps);
624 AddSummaryItem(s, "Checkpoints (steps)", m_checksteps);
625 if (m_intScheme)
626 {
627 AddSummaryItem(s, "Integration Type", m_intScheme->GetName());
628 }
629}
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_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion 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:47

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, and Nektar::SolverUtils::EquationSystem::v_GenerateSummary().

Referenced by Nektar::UnsteadyAdvection::v_GenerateSummary(), Nektar::UnsteadyInviscidBurgers::v_GenerateSummary(), Nektar::CompressibleFlowSystem::v_GenerateSummary(), Nektar::ImageWarpingSystem::v_GenerateSummary(), Nektar::VelocityCorrectionScheme::v_GenerateSummary(), Nektar::VCSImplicit::v_GenerateSummary(), Nektar::VCSWeakPressure::v_GenerateSummary(), Nektar::PulseWavePropagation::v_GenerateSummary(), Nektar::ShallowWaterSystem::v_GenerateSummary(), Nektar::SolverUtils::MMFSystem::v_GenerateSummary(), Nektar::UnsteadyDiffusion::v_GenerateSummary(), Nektar::Bidomain::v_GenerateSummary(), Nektar::BidomainRoth::v_GenerateSummary(), and Nektar::Monodomain::v_GenerateSummary().

◆ v_GetTimeStep()

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 728 of file UnsteadySystem.cpp.

730{
731 NEKERROR(ErrorUtil::efatal, "Not defined for this class");
732 return 0.0;
733}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by GetTimeStep().

◆ v_InitObject()

void Nektar::SolverUtils::UnsteadySystem::v_InitObject ( bool  DeclareField = true)
overrideprotectedvirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Reimplemented in Nektar::PulseWavePropagation, Nektar::PulseWaveSystem, Nektar::Bidomain, Nektar::BidomainRoth, Nektar::Monodomain, Nektar::NavierStokesCFE, Nektar::MMFDiffusion, Nektar::ImageWarpingSystem, Nektar::CoupledLinearNS, Nektar::IncNavierStokes, Nektar::SmoothedProfileMethod, Nektar::VCSMapping, Nektar::VelocityCorrectionScheme, Nektar::SolverUtils::FileSolution, Nektar::AcousticSystem, Nektar::APE, Nektar::LEE, Nektar::CFLtester, Nektar::MMFAdvection, Nektar::UnsteadyAdvection, Nektar::UnsteadyAdvectionDiffusion, Nektar::UnsteadyDiffusion, Nektar::UnsteadyInviscidBurgers, Nektar::UnsteadyReactionDiffusion, Nektar::UnsteadyViscousBurgers, Nektar::CompressibleFlowSystem, Nektar::CFSImplicit, Nektar::EulerCFE, Nektar::EulerImplicitCFE, Nektar::NavierStokesCFEAxisym, Nektar::NavierStokesImplicitCFE, Nektar::Dummy, Nektar::MMFMaxwell, Nektar::LinearSWE, Nektar::MMFSWE, Nektar::NonlinearPeregrine, Nektar::NonlinearSWE, and Nektar::ShallowWaterSystem.

Definition at line 85 of file UnsteadySystem.cpp.

86{
87 EquationSystem::v_InitObject(DeclareField);
88
89 m_initialStep = 0;
90
91 // Load SolverInfo parameters.
92 m_session->MatchSolverInfo("DIFFUSIONADVANCEMENT", "Explicit",
94 m_session->MatchSolverInfo("ADVECTIONADVANCEMENT", "Explicit",
96 m_session->MatchSolverInfo("REACTIONADVANCEMENT", "Explicit",
97 m_explicitReaction, true);
98 m_session->LoadParameter("CheckAbortSteps", m_abortSteps, 1);
99
100 // Steady state tolerance.
101 m_session->LoadParameter("SteadyStateTol", m_steadyStateTol, 0.0);
102
103 // Frequency for checking steady state.
104 m_session->LoadParameter("SteadyStateSteps", m_steadyStateSteps, 1);
105
106 // For steady problems, we do not initialise the time integration.
107 if (m_session->DefinesSolverInfo("TimeIntegrationMethod") ||
108 m_session->DefinesTimeIntScheme())
109 {
110 LibUtilities::TimeIntScheme timeInt;
111 if (m_session->DefinesTimeIntScheme())
112 {
113 timeInt = m_session->GetTimeIntScheme();
114 }
115 else
116 {
117 timeInt.method = m_session->GetSolverInfo("TimeIntegrationMethod");
118 }
119
122 timeInt.method, timeInt.variant, timeInt.order,
123 timeInt.freeParams);
124
125 // Load generic input parameters.
126 m_session->LoadParameter("IO_InfoSteps", m_infosteps, 0);
127 m_session->LoadParameter("IO_FiltersInfoSteps", m_filtersInfosteps,
128 10 * m_infosteps);
129 m_session->LoadParameter("CFL", m_cflSafetyFactor, 0.0);
130 m_session->LoadParameter("CFLEnd", m_CFLEnd, 0.0);
131 m_session->LoadParameter("CFLGrowth", m_CFLGrowth, 1.0);
132
133 // Ensure that there is no conflict of parameters.
134 if (m_cflSafetyFactor > 0.0)
135 {
136 // Check final condition.
137 ASSERTL0(m_fintime == 0.0 || m_steps == 0,
138 "Final condition not unique: "
139 "fintime > 0.0 and Nsteps > 0");
140 // Check timestep condition.
141 ASSERTL0(m_timestep == 0.0,
142 "Timestep not unique: timestep > 0.0 & CFL > 0.0");
143 }
144 else
145 {
146 ASSERTL0(m_timestep != 0.0, "Need to set either TimeStep or CFL");
147 }
148
149 // Ensure that there is no conflict of parameters.
150 if (m_CFLGrowth > 1.0)
151 {
152 // Check final condition.
154 "m_CFLEnd >= m_cflSafetyFactor required");
155 }
156
157 // Set up time to be dumped in field information.
158 m_fieldMetaDataMap["Time"] = boost::lexical_cast<std::string>(m_time);
159 }
160
161 // By default attempt to forward transform initial condition.
162 m_homoInitialFwd = true;
163
164 // Set up filters.
165 for (auto &x : m_session->GetFilters())
166 {
167 m_filters.push_back(make_pair(
168 x.first, GetFilterFactory().CreateInstance(
169 x.first, m_session, shared_from_this(), x.second)));
170 }
171}
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareFeld=true)
Initialisation object for EquationSystem.
TimeIntegrationSchemeFactory & GetTimeIntegrationSchemeFactory()
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:39

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::TimeIntScheme::freeParams, Nektar::SolverUtils::GetFilterFactory(), Nektar::LibUtilities::GetTimeIntegrationSchemeFactory(), m_abortSteps, m_CFLEnd, m_CFLGrowth, m_cflSafetyFactor, m_explicitAdvection, m_explicitDiffusion, m_explicitReaction, Nektar::SolverUtils::EquationSystem::m_fieldMetaDataMap, m_filters, m_filtersInfosteps, Nektar::SolverUtils::EquationSystem::m_fintime, m_homoInitialFwd, Nektar::SolverUtils::EquationSystem::m_infosteps, Nektar::SolverUtils::EquationSystem::m_initialStep, m_intScheme, Nektar::SolverUtils::EquationSystem::m_session, m_steadyStateSteps, m_steadyStateTol, Nektar::SolverUtils::EquationSystem::m_steps, Nektar::SolverUtils::EquationSystem::m_time, Nektar::SolverUtils::EquationSystem::m_timestep, Nektar::LibUtilities::TimeIntScheme::method, Nektar::LibUtilities::TimeIntScheme::order, Nektar::SolverUtils::EquationSystem::v_InitObject(), and Nektar::LibUtilities::TimeIntScheme::variant.

Referenced by Nektar::PulseWaveSystem::v_InitObject(), Nektar::SolverUtils::AdvectionSystem::v_InitObject(), Nektar::Bidomain::v_InitObject(), Nektar::BidomainRoth::v_InitObject(), Nektar::Monodomain::v_InitObject(), Nektar::MMFDiffusion::v_InitObject(), Nektar::MMFAdvection::v_InitObject(), Nektar::UnsteadyDiffusion::v_InitObject(), Nektar::Dummy::v_InitObject(), Nektar::MMFMaxwell::v_InitObject(), Nektar::MMFSWE::v_InitObject(), and Nektar::ShallowWaterSystem::v_InitObject().

◆ v_PostIntegrate()

bool Nektar::SolverUtils::UnsteadySystem::v_PostIntegrate ( int  step)
protectedvirtual

◆ v_PreIntegrate()

bool Nektar::SolverUtils::UnsteadySystem::v_PreIntegrate ( int  step)
protectedvirtual

◆ v_PrintStatusInformation()

void Nektar::SolverUtils::UnsteadySystem::v_PrintStatusInformation ( const int  step,
const NekDouble  cpuTime 
)
protectedvirtual

Print Status Information.

Reimplemented in Nektar::CFSImplicit.

Definition at line 543 of file UnsteadySystem.cpp.

545{
546 if (m_infosteps && m_session->GetComm()->GetSpaceComm()->GetRank() == 0 &&
547 !((step + 1) % m_infosteps))
548 {
549 if (m_comm->IsParallelInTime())
550 {
551 cout << "RANK " << m_session->GetComm()->GetTimeComm()->GetRank()
552 << " Steps: " << setw(8) << left << step + 1 << " "
553 << "Time: " << setw(12) << left << m_time;
554 }
555 else
556 {
557 cout << "Steps: " << setw(8) << left << step + 1 << " "
558 << "Time: " << setw(12) << left << m_time;
559 }
560
562 {
563 cout << " Time-step: " << setw(12) << left << m_timestep;
564 }
565
566 stringstream ss;
567 ss << cpuTime << "s";
568 cout << " CPU Time: " << setw(8) << left << ss.str() << endl;
569 }
570}

References m_cflSafetyFactor, Nektar::SolverUtils::EquationSystem::m_comm, Nektar::SolverUtils::EquationSystem::m_infosteps, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_time, and Nektar::SolverUtils::EquationSystem::m_timestep.

Referenced by v_DoSolve(), and Nektar::CFSImplicit::v_PrintStatusInformation().

◆ v_PrintSummaryStatistics()

void Nektar::SolverUtils::UnsteadySystem::v_PrintSummaryStatistics ( const NekDouble  intTime)
protectedvirtual

Print Summary Statistics.

Reimplemented in Nektar::CFSImplicit.

Definition at line 572 of file UnsteadySystem.cpp.

573{
574 if (m_session->GetComm()->GetRank() == 0)
575 {
576 if (m_cflSafetyFactor > 0.0)
577 {
578 cout << "CFL safety factor : " << m_cflSafetyFactor << endl
579 << "CFL time-step : " << m_timestep << endl;
580 }
581
582 if (m_session->GetSolverInfo("Driver") != "SteadyState" &&
583 m_session->GetSolverInfo("Driver") != "Parareal" &&
584 m_session->GetSolverInfo("Driver") != "PFASST")
585 {
586 cout << "Time-integration : " << intTime << "s" << endl;
587 }
588 }
589}

References m_cflSafetyFactor, Nektar::SolverUtils::EquationSystem::m_session, and Nektar::SolverUtils::EquationSystem::m_timestep.

Referenced by v_DoSolve(), and Nektar::CFSImplicit::v_PrintSummaryStatistics().

◆ v_RequireFwdTrans()

bool Nektar::SolverUtils::UnsteadySystem::v_RequireFwdTrans ( void  )
protectedvirtual

Reimplemented in Nektar::Dummy, Nektar::VelocityCorrectionScheme, and Nektar::SolverUtils::FileSolution.

Definition at line 754 of file UnsteadySystem.cpp.

755{
756 return true;
757}

Referenced by v_DoSolve().

◆ v_SteadyStateResidual()

void Nektar::SolverUtils::UnsteadySystem::v_SteadyStateResidual ( int  step,
Array< OneD, NekDouble > &  L2 
)
protectedvirtual

Reimplemented in Nektar::CompressibleFlowSystem.

Definition at line 910 of file UnsteadySystem.cpp.

912{
913 const int nPoints = GetTotPoints();
914 const int nFields = m_fields.size();
915
916 // Holds L2 errors.
917 Array<OneD, NekDouble> RHSL2(nFields);
918 Array<OneD, NekDouble> residual(nFields);
919 Array<OneD, NekDouble> reference(nFields);
920
921 for (int i = 0; i < nFields; ++i)
922 {
923 Array<OneD, NekDouble> tmp(nPoints);
924
925 Vmath::Vsub(nPoints, m_fields[i]->GetPhys(), 1, m_previousSolution[i],
926 1, tmp, 1);
927 Vmath::Vmul(nPoints, tmp, 1, tmp, 1, tmp, 1);
928 residual[i] = Vmath::Vsum(nPoints, tmp, 1);
929
931 tmp, 1);
932 reference[i] = Vmath::Vsum(nPoints, tmp, 1);
933 }
934
935 m_comm->GetSpaceComm()->AllReduce(residual, LibUtilities::ReduceSum);
936 m_comm->GetSpaceComm()->AllReduce(reference, LibUtilities::ReduceSum);
937
938 // L2 error.
939 for (int i = 0; i < nFields; ++i)
940 {
941 reference[i] = (reference[i] == 0) ? 1 : reference[i];
942 L2[i] = sqrt(residual[i] / reference[i]);
943 }
944}
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.hpp:608
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.hpp:220
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::SolverUtils::EquationSystem::m_comm, Nektar::SolverUtils::EquationSystem::m_fields, m_previousSolution, Nektar::LibUtilities::ReduceSum, tinysimd::sqrt(), Vmath::Vmul(), Vmath::Vsub(), and Vmath::Vsum().

Referenced by SteadyStateResidual().

◆ v_UpdateTimeStepCheck()

bool Nektar::SolverUtils::UnsteadySystem::v_UpdateTimeStepCheck ( void  )
protectedvirtual

Reimplemented in Nektar::CFSImplicit.

Definition at line 762 of file UnsteadySystem.cpp.

763{
764 return true;
765}

Referenced by v_DoSolve().

Member Data Documentation

◆ cmdSetStartChkNum

std::string Nektar::SolverUtils::UnsteadySystem::cmdSetStartChkNum
static
Initial value:
=
"set-start-chknumber", "",
"Set the starting number of the checkpoint file.")
static std::string RegisterCmdLineArgument(const std::string &pName, const std::string &pShortName, const std::string &pDescription)
Registers a command-line argument with the session reader.

Definition at line 81 of file UnsteadySystem.h.

◆ cmdSetStartTime

std::string Nektar::SolverUtils::UnsteadySystem::cmdSetStartTime
static
Initial value:
=
"set-start-time", "", "Set the starting time of the simulation.")

Definition at line 80 of file UnsteadySystem.h.

◆ m_abortSteps

int Nektar::SolverUtils::UnsteadySystem::m_abortSteps
protected

Number of steps between checks for abort conditions.

Definition at line 103 of file UnsteadySystem.h.

Referenced by v_DoSolve(), and v_InitObject().

◆ m_CFLEnd

NekDouble Nektar::SolverUtils::UnsteadySystem::m_CFLEnd
protected

Maximun cfl in cfl growth.

Definition at line 100 of file UnsteadySystem.h.

Referenced by v_DoSolve(), and v_InitObject().

◆ m_CFLGrowth

NekDouble Nektar::SolverUtils::UnsteadySystem::m_CFLGrowth
protected

CFL growth rate.

Definition at line 98 of file UnsteadySystem.h.

Referenced by v_DoSolve(), and v_InitObject().

◆ m_cflSafetyFactor

NekDouble Nektar::SolverUtils::UnsteadySystem::m_cflSafetyFactor
protected

◆ m_epsilon

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

Diffusion coefficient.

Definition at line 129 of file UnsteadySystem.h.

◆ m_errFile

std::ofstream Nektar::SolverUtils::UnsteadySystem::m_errFile
protected

Definition at line 126 of file UnsteadySystem.h.

Referenced by CheckSteadyState(), and InitializeSteadyState().

◆ m_explicitAdvection

bool Nektar::SolverUtils::UnsteadySystem::m_explicitAdvection
protected

◆ m_explicitDiffusion

bool Nektar::SolverUtils::UnsteadySystem::m_explicitDiffusion
protected

◆ m_explicitReaction

bool Nektar::SolverUtils::UnsteadySystem::m_explicitReaction
protected

Indicates if explicit or implicit treatment of reaction is used.

Definition at line 110 of file UnsteadySystem.h.

Referenced by v_GenerateSummary(), and v_InitObject().

◆ m_filters

std::vector<std::pair<std::string, FilterSharedPtr> > Nektar::SolverUtils::UnsteadySystem::m_filters
protected

◆ m_filtersInfosteps

int Nektar::SolverUtils::UnsteadySystem::m_filtersInfosteps
protected

Number of time steps between outputting filters information.

Definition at line 118 of file UnsteadySystem.h.

Referenced by v_DoSolve(), and v_InitObject().

◆ m_homoInitialFwd

bool Nektar::SolverUtils::UnsteadySystem::m_homoInitialFwd
protected

◆ m_intScheme

LibUtilities::TimeIntegrationSchemeSharedPtr Nektar::SolverUtils::UnsteadySystem::m_intScheme
protected

◆ m_intVariables

std::vector<int> Nektar::SolverUtils::UnsteadySystem::m_intVariables
protected

◆ m_ode

LibUtilities::TimeIntegrationSchemeOperators Nektar::SolverUtils::UnsteadySystem::m_ode
protected

◆ m_previousSolution

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::UnsteadySystem::m_previousSolution
protected

Storage for previous solution for steady-state check.

Definition at line 91 of file UnsteadySystem.h.

Referenced by CheckSteadyState(), InitializeSteadyState(), and v_SteadyStateResidual().

◆ m_steadyStateSteps

int Nektar::SolverUtils::UnsteadySystem::m_steadyStateSteps
protected

Check for steady state at step interval.

Definition at line 113 of file UnsteadySystem.h.

Referenced by v_DoSolve(), and v_InitObject().

◆ m_steadyStateTol

NekDouble Nektar::SolverUtils::UnsteadySystem::m_steadyStateTol
protected

Tolerance to which steady state should be evaluated at.

Definition at line 115 of file UnsteadySystem.h.

Referenced by CheckSteadyState(), InitializeSteadyState(), v_DoSolve(), and v_InitObject().