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

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

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

Protected Attributes

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

Private Member Functions

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

Constructor & Destructor Documentation

◆ ~UnsteadySystem()

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

Destructor.

Destructor for the class UnsteadyAdvection.

Definition at line 178 of file UnsteadySystem.cpp.

179{
180}

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

79 : EquationSystem(pSession, pGraph)
80
81{
82}
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 636 of file UnsteadySystem.cpp.

637{
638 // Coordinates of the quadrature points in the real physical space.
639 Array<OneD, NekDouble> x(GetNpoints());
640 Array<OneD, NekDouble> y(GetNpoints());
641 Array<OneD, NekDouble> z(GetNpoints());
642 m_fields[0]->GetCoords(x, y, z);
643
644 // Print out the solution in a txt file.
645 ofstream outfile;
646 outfile.open("solution1D.txt");
647 for (int i = 0; i < GetNpoints(); i++)
648 {
649 outfile << scientific << setw(17) << setprecision(16) << x[i] << " "
650 << m_fields[m_intVariables[0]]->GetPhys()[i] << endl;
651 }
652 outfile << endl << endl;
653 outfile.close();
654}
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 659 of file UnsteadySystem.cpp.

660{
661 if (m_session->DefinesFunction("InitialConditions"))
662 {
663 for (int i = 0; i < m_fields.size(); ++i)
664 {
666
667 vType = m_session->GetFunctionType("InitialConditions",
668 m_session->GetVariable(i));
669
671 {
672 std::string filename = m_session->GetFunctionFilename(
673 "InitialConditions", m_session->GetVariable(i));
674
675 fs::path pfilename(filename);
676
677 // Redefine path for parallel file which is in directory.
678 if (fs::is_directory(pfilename))
679 {
680 fs::path metafile("Info.xml");
681 fs::path fullpath = pfilename / metafile;
682 filename = LibUtilities::PortablePath(fullpath);
683 }
686 fld->ImportFieldMetaData(filename, m_fieldMetaDataMap);
687
688 // Check to see if time defined.
690 {
691 auto iter = m_fieldMetaDataMap.find("Time");
692 if (iter != m_fieldMetaDataMap.end())
693 {
694 time = boost::lexical_cast<NekDouble>(iter->second);
695 }
696
697 iter = m_fieldMetaDataMap.find("ChkFileNum");
698 if (iter != m_fieldMetaDataMap.end())
699 {
700 nchk = boost::lexical_cast<NekDouble>(iter->second);
701 }
702 }
703
704 break;
705 }
706 }
707 }
708 if (m_session->DefinesCmdLineArgument("set-start-time"))
709 {
710 time = boost::lexical_cast<NekDouble>(
711 m_session->GetCmdLineArgument<std::string>("set-start-time")
712 .c_str());
713 }
714 if (m_session->DefinesCmdLineArgument("set-start-chknumber"))
715 {
716 nchk = boost::lexical_cast<int>(
717 m_session->GetCmdLineArgument<std::string>("set-start-chknumber"));
718 }
719 ASSERTL0(time >= 0 && nchk >= 0,
720 "Starting time and checkpoint number should be >= 0");
721}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
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:226
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:327
std::string PortablePath(const boost::filesystem::path &path)
create portable path on different platforms for boost::filesystem path
Definition: FileSystem.cpp:45
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:53

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

858{
859 const int nPoints = GetTotPoints();
860 const int nFields = m_fields.size();
861
862 // Holds L2 errors.
863 Array<OneD, NekDouble> L2(nFields);
864
865 SteadyStateResidual(step, L2);
866
867 if (m_infosteps && m_comm->GetRank() == 0 &&
868 (((step + 1) % m_infosteps == 0) || ((step == m_initialStep))))
869 {
870 // Output time.
871 m_errFile << boost::format("%25.19e") % m_time;
872
873 m_errFile << " " << boost::format("%25.19e") % totCPUTime;
874
875 int stepp = step + 1;
876
877 m_errFile << " " << boost::format("%25.19e") % stepp;
878
879 // Output residuals.
880 for (int i = 0; i < nFields; ++i)
881 {
882 m_errFile << " " << boost::format("%25.19e") % L2[i];
883 }
884
885 m_errFile << endl;
886 }
887
888 // Calculate maximum L2 error.
889 NekDouble maxL2 = Vmath::Vmax(nFields, L2, 1);
890
891 if (m_infosteps && m_session->DefinesCmdLineArgument("verbose") &&
892 m_comm->GetRank() == 0 && ((step + 1) % m_infosteps == 0))
893 {
894 cout << "-- Maximum L^2 residual: " << maxL2 << endl;
895 }
896
897 if (maxL2 <= m_steadyStateTol)
898 {
899 return true;
900 }
901
902 for (int i = 0; i < m_fields.size(); ++i)
903 {
904 Vmath::Vcopy(nPoints, m_fields[i]->GetPhys(), 1, m_previousSolution[i],
905 1);
906 }
907
908 return false;
909}
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.cpp:940
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191

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

956{
957 boost::ignore_unused(time);
958
959 if (&inarray != &outarray)
960 {
961 for (int i = 0; i < inarray.size(); ++i)
962 {
963 Vmath::Vcopy(GetNpoints(), inarray[i], 1, outarray[i], 1);
964 }
965 }
966}

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

195{
196 return m_intScheme;
197}
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 202 of file UnsteadySystem.cpp.

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

References m_ode.

◆ GetTimeStep()

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

56 {
57 return v_GetTimeStep(inarray);
58 }
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 817 of file UnsteadySystem.cpp.

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

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

186{
187 return m_intScheme->GetTimeStability();
188}

References m_intScheme.

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

◆ SteadyStateResidual()

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

Definition at line 60 of file UnsteadySystem.h.

62 {
63 v_SteadyStateResidual(step, L2);
64 }
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 774 of file UnsteadySystem.cpp.

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

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

Referenced by Nektar::UnsteadyViscousBurgers::DoImplicitSolve(), and Nektar::UnsteadyViscousBurgers::v_InitObject().

◆ 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::PulseWaveSystem, Nektar::MMFSWE, Nektar::CoupledLinearNS, Nektar::VCSMapping, and Nektar::VelocityCorrectionScheme.

Definition at line 595 of file UnsteadySystem.cpp.

596{
599 SetInitialConditions(m_time, dumpInitialConditions);
601}
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 212 of file UnsteadySystem.cpp.

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

References AppendOutput1D(), ASSERTL0, Nektar::SolverUtils::EquationSystem::Checkpoint_Output(), CheckSteadyState(), Nektar::SolverUtils::EquationSystem::eNotHomogeneous, Nektar::SolverUtils::EquationSystem::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::SolverUtils::EquationSystem::ParallelInTime(), 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::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::SolverUtils::MMFSystem, Nektar::CFLtester, Nektar::UnsteadyDiffusion, Nektar::Bidomain, Nektar::BidomainRoth, and Nektar::Monodomain.

Definition at line 607 of file UnsteadySystem.cpp.

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

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::UnsteadyAdvectionDiffusion::v_GenerateSummary(), Nektar::UnsteadyViscousBurgers::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::CFLtester::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 729 of file UnsteadySystem.cpp.

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

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::SolverUtils::AdvectionSystem, Nektar::Bidomain, Nektar::BidomainRoth, Nektar::Monodomain, Nektar::NavierStokesCFE, Nektar::MMFDiffusion, Nektar::ImageWarpingSystem, Nektar::CoupledLinearNS, Nektar::IncNavierStokes, Nektar::SmoothedProfileMethod, Nektar::VCSMapping, Nektar::VelocityCorrectionScheme, Nektar::PulseWaveSystemOutput, Nektar::AcousticSystem, Nektar::APE, Nektar::LEE, Nektar::CFLtester, Nektar::MMFAdvection, Nektar::UnsteadyAdvection, Nektar::UnsteadyAdvectionDiffusion, Nektar::UnsteadyDiffusion, Nektar::UnsteadyInviscidBurger, 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 87 of file UnsteadySystem.cpp.

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

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::UnsteadyReactionDiffusion::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

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

Definition at line 749 of file UnsteadySystem.cpp.

750{
751 boost::ignore_unused(step);
752 return false;
753}

Referenced by v_DoSolve(), Nektar::SolverUtils::AdvectionSystem::v_PostIntegrate(), and Nektar::Dummy::v_PostIntegrate().

◆ v_PreIntegrate()

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

Reimplemented in Nektar::AcousticSystem, Nektar::UnsteadyAdvectionDiffusion, Nektar::Dummy, and Nektar::IncNavierStokes.

Definition at line 740 of file UnsteadySystem.cpp.

741{
742 boost::ignore_unused(step);
743 return false;
744}

Referenced by v_DoSolve(), Nektar::AcousticSystem::v_PreIntegrate(), and Nektar::Dummy::v_PreIntegrate().

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

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

References m_cflSafetyFactor, Nektar::SolverUtils::EquationSystem::m_infosteps, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_time, Nektar::SolverUtils::EquationSystem::m_timestep, and Nektar::SolverUtils::EquationSystem::ParallelInTime().

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

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

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, and Nektar::VelocityCorrectionScheme.

Definition at line 758 of file UnsteadySystem.cpp.

759{
760 return true;
761}

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

915{
916 boost::ignore_unused(step);
917 const int nPoints = GetTotPoints();
918 const int nFields = m_fields.size();
919
920 // Holds L2 errors.
921 Array<OneD, NekDouble> RHSL2(nFields);
922 Array<OneD, NekDouble> residual(nFields);
923 Array<OneD, NekDouble> reference(nFields);
924
925 for (int i = 0; i < nFields; ++i)
926 {
927 Array<OneD, NekDouble> tmp(nPoints);
928
929 Vmath::Vsub(nPoints, m_fields[i]->GetPhys(), 1, m_previousSolution[i],
930 1, tmp, 1);
931 Vmath::Vmul(nPoints, tmp, 1, tmp, 1, tmp, 1);
932 residual[i] = Vmath::Vsum(nPoints, tmp, 1);
933
935 tmp, 1);
936 reference[i] = Vmath::Vsum(nPoints, tmp, 1);
937 }
938
939 m_comm->GetSpaceComm()->AllReduce(residual, LibUtilities::ReduceSum);
940 m_comm->GetSpaceComm()->AllReduce(reference, LibUtilities::ReduceSum);
941
942 // L2 error.
943 for (int i = 0; i < nFields; ++i)
944 {
945 reference[i] = (reference[i] == 0) ? 1 : reference[i];
946 L2[i] = sqrt(residual[i] / reference[i]);
947 }
948}
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:890
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:414
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 766 of file UnsteadySystem.cpp.

767{
768 return true;
769}

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

◆ m_abortSteps

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

Number of steps between checks for abort conditions.

Definition at line 94 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 91 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 89 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 120 of file UnsteadySystem.h.

◆ m_errFile

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

Definition at line 117 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 101 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 109 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 82 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 104 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 106 of file UnsteadySystem.h.

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