Nektar++
Public Member Functions | Public Attributes | 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)
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void SetUpTraceNormals (void)
 
SOLVER_UTILS_EXPORT void InitObject (bool DeclareField=true)
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise ()
 Perform any initialisation necessary before solving the problem. More...
 
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem. More...
 
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space. More...
 
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space. More...
 
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve. More...
 
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. More...
 
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name. More...
 
template<class T >
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 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, Array< OneD, NekDouble > &output)
 
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 Array< OneD, const Array< OneD, NekDouble > > GetTraceNormals ()
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to identify if operator is negated in DoSolve. More...
 

Public Attributes

NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
NekDouble m_cflNonAcoustic
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 maximun cfl in cfl growth 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)
 Init object for UnsteadySystem class. More...
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise ()
 Sets up initial conditions. More...
 
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &s)
 Print a summary of time stepping parameters. More...
 
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D (Array< OneD, Array< OneD, NekDouble >> &solution1D)
 Print the solution at each solution point in a txt file. More...
 
virtual SOLVER_UTILS_EXPORT 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)
 
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...
 
virtual SOLVER_UTILS_EXPORT bool UpdateTimeStepCheck ()
 
- 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 NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space. More...
 
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
virtual SOLVER_UTILS_EXPORT void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure (void)
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables)
 

Protected Attributes

int m_infosteps
 Number of time steps between outputting status information. More...
 
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
int m_nanSteps
 
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
NekDouble m_epsilon
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateRes = 1.0
 
NekDouble m_steadyStateRes0 = 1.0
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::ofstream m_errFile
 
std::vector< int > m_intVariables
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
NekDouble m_filterTimeWarning
 Number of time steps between outputting status information. More...
 
NekDouble m_TimeIntegLambda = 0.0
 coefff of spacial derivatives(rhs or m_F in GLM) in calculating the residual of the whole equation(used in unsteady time integrations) More...
 
bool m_flagImplicitItsStatistics
 
bool m_flagImplicitSolver = false
 
Array< OneD, NekDoublem_magnitdEstimat
 estimate the magnitude of each conserved varibles More...
 
Array< OneD, NekDoublem_locTimeStep
 local time step(notice only for jfnk other see m_cflSafetyFactor) More...
 
NekDouble m_inArrayNorm = -1.0
 
int m_TotLinItePerStep = 0
 
int m_StagesPerStep = 1
 
bool m_flagUpdatePreconMat
 
int m_maxLinItePerNewton
 
int m_TotNewtonIts = 0
 
int m_TotLinIts = 0
 
int m_TotImpStages = 0
 
bool m_CalcPhysicalAV = true
 flag to update artificial viscosity More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
bool m_verbose
 
bool m_root
 
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_timestepMax = -1.0
 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_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

void InitializeSteadyState ()
 
bool CheckSteadyState (int step)
 Calculate whether the system has reached a steady state by observing residuals to a user-defined tolerance. More...
 
bool CheckSteadyState (int step, NekDouble totCPUTime)
 

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 []
 

Detailed Description

Base class for unsteady solvers.

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

For details on implementing unsteady solvers see sectionADRSolverModuleImplementation here

Definition at line 48 of file UnsteadySystem.h.

Constructor & Destructor Documentation

◆ ~UnsteadySystem()

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

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

71  : EquationSystem(pSession, pGraph), m_infosteps(10)
72 
73 {
74 }
SOLVER_UTILS_EXPORT EquationSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises EquationSystem class members.
int m_infosteps
Number of time steps between outputting status information.

Member Function Documentation

◆ CheckForRestartTime()

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

Definition at line 635 of file UnsteadySystem.cpp.

636 {
637  if (m_session->DefinesFunction("InitialConditions"))
638  {
639  for (int i = 0; i < m_fields.size(); ++i)
640  {
642 
643  vType = m_session->GetFunctionType("InitialConditions",
644  m_session->GetVariable(i));
645 
646  if (vType == LibUtilities::eFunctionTypeFile)
647  {
648  std::string filename = m_session->GetFunctionFilename(
649  "InitialConditions", m_session->GetVariable(i));
650 
651  fs::path pfilename(filename);
652 
653  // redefine path for parallel file which is in directory
654  if (fs::is_directory(pfilename))
655  {
656  fs::path metafile("Info.xml");
657  fs::path fullpath = pfilename / metafile;
658  filename = LibUtilities::PortablePath(fullpath);
659  }
662  fld->ImportFieldMetaData(filename, m_fieldMetaDataMap);
663 
664  // check to see if time defined
666  {
667  auto iter = m_fieldMetaDataMap.find("Time");
668  if (iter != m_fieldMetaDataMap.end())
669  {
670  time = boost::lexical_cast<NekDouble>(iter->second);
671  }
672 
673  iter = m_fieldMetaDataMap.find("ChkFileNum");
674  if (iter != m_fieldMetaDataMap.end())
675  {
676  nchk = boost::lexical_cast<NekDouble>(iter->second);
677  }
678  }
679 
680  break;
681  }
682  }
683  }
684  if (m_session->DefinesCmdLineArgument("set-start-time"))
685  {
686  time = boost::lexical_cast<NekDouble>(
687  m_session->GetCmdLineArgument<std::string>("set-start-time")
688  .c_str());
689  }
690  if (m_session->DefinesCmdLineArgument("set-start-chknumber"))
691  {
692  nchk = boost::lexical_cast<int>(
693  m_session->GetCmdLineArgument<std::string>("set-start-chknumber"));
694  }
695  ASSERTL0(time >= 0 && nchk >= 0,
696  "Starting time and checkpoint number should be >= 0");
697 }
#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:227
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
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:301
std::string PortablePath(const boost::filesystem::path &path)
create portable path on different platforms for boost::filesystem path
Definition: FileSystem.cpp:41
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() [1/2]

bool Nektar::SolverUtils::UnsteadySystem::CheckSteadyState ( int  step)
private

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

Definition at line 819 of file UnsteadySystem.cpp.

820 {
821  return CheckSteadyState(step, 0.0);
822 }
bool CheckSteadyState(int step)
Calculate whether the system has reached a steady state by observing residuals to a user-defined tole...

Referenced by v_DoSolve().

◆ CheckSteadyState() [2/2]

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

Definition at line 824 of file UnsteadySystem.cpp.

825 {
826  const int nPoints = GetTotPoints();
827  const int nFields = m_fields.size();
828 
829  // Holds L2 errors.
830  Array<OneD, NekDouble> L2(nFields);
831 
832  SteadyStateResidual(step, L2);
833 
834  if (m_comm->GetRank() == 0 &&
835  (((step + 1) % m_infosteps == 0) || ((step == m_initialStep))))
836  {
837  // Output time
838  m_errFile << boost::format("%25.19e") % m_time;
839 
840  m_errFile << " " << boost::format("%25.19e") % totCPUTime;
841 
842  int stepp = step + 1;
843 
844  m_errFile << " " << boost::format("%25.19e") % stepp;
845 
846  // Output residuals
847  for (int i = 0; i < nFields; ++i)
848  {
849  m_errFile << " " << boost::format("%25.19e") % L2[i];
850  }
851 
852  m_errFile << endl;
853  }
854 
855  // Calculate maximum L2 error
856  NekDouble maxL2 = Vmath::Vmax(nFields, L2, 1);
857 
858  if (m_session->DefinesCmdLineArgument("verbose") &&
859  m_comm->GetRank() == 0 && ((step + 1) % m_infosteps == 0))
860  {
861  cout << "-- Maximum L^2 residual: " << maxL2 << endl;
862  }
863 
864  if (maxL2 <= m_steadyStateTol)
865  {
866  return true;
867  }
868 
869  for (int i = 0; i < m_fields.size(); ++i)
870  {
871  Vmath::Vcopy(nPoints, m_fields[i]->GetPhys(), 1, m_previousSolution[i],
872  1);
873  }
874 
876  m_steadyStateRes = maxL2;
877 
878  return false;
879 }
LibUtilities::CommSharedPtr m_comm
Communicator.
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:945
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

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

◆ GetTimeStep()

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

Calculate the larger time-step mantaining the problem stable.

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

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

Definition at line 707 of file UnsteadySystem.cpp.

709 {
710  return v_GetTimeStep(inarray);
711 }
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 779 of file UnsteadySystem.cpp.

780 {
781  if (m_steadyStateTol > 0.0)
782  {
783  const int nPoints = m_fields[0]->GetTotPoints();
785  Array<OneD, Array<OneD, NekDouble>>(m_fields.size());
786 
787  for (int i = 0; i < m_fields.size(); ++i)
788  {
789  m_previousSolution[i] = Array<OneD, NekDouble>(nPoints);
790  Vmath::Vcopy(nPoints, m_fields[i]->GetPhys(), 1,
791  m_previousSolution[i], 1);
792  }
793 
794  if (m_comm->GetRank() == 0)
795  {
796  std::string fName =
797  m_session->GetSessionName() + std::string(".resd");
798  m_errFile.open(fName.c_str());
799  m_errFile << setw(26) << left << "# Time";
800 
801  m_errFile << setw(26) << left << "CPU_Time";
802 
803  m_errFile << setw(26) << left << "Step";
804 
805  for (int i = 0; i < m_fields.size(); ++i)
806  {
807  m_errFile << setw(26) << m_session->GetVariables()[i];
808  }
809 
810  m_errFile << endl;
811  }
812  }
813 }

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 }
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
Wrapper to the time integration scheme.

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

60  {
61  v_SteadyStateResidual(step, L2);
62  }
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 739 of file UnsteadySystem.cpp.

742 {
743  int phystot = m_fields[0]->GetTotPoints();
744  int nvel = vel.size();
745 
746  Array<OneD, NekDouble> varcoeff(phystot), tmp;
747 
748  // calculate magnitude of v
749  Vmath::Vmul(phystot, vel[0], 1, vel[0], 1, varcoeff, 1);
750  for (int n = 1; n < nvel; ++n)
751  {
752  Vmath::Vvtvp(phystot, vel[n], 1, vel[n], 1, varcoeff, 1, varcoeff, 1);
753  }
754  Vmath::Vsqrt(phystot, varcoeff, 1, varcoeff, 1);
755 
756  for (int i = 0; i < m_fields[0]->GetNumElmts(); ++i)
757  {
758  int offset = m_fields[0]->GetPhys_Offset(i);
759  int nq = m_fields[0]->GetExp(i)->GetTotPoints();
760  Array<OneD, NekDouble> unit(nq, 1.0);
761 
762  int nmodes = 0;
763 
764  for (int n = 0; n < m_fields[0]->GetExp(i)->GetNumBases(); ++n)
765  {
766  nmodes = max(nmodes, m_fields[0]->GetExp(i)->GetBasisNumModes(n));
767  }
768 
769  NekDouble h = m_fields[0]->GetExp(i)->Integral(unit);
770  h = pow(h, (NekDouble)(1.0 / nvel)) / ((NekDouble)nmodes);
771 
772  Vmath::Smul(nq, h, varcoeff + offset, 1, tmp = varcoeff + offset, 1);
773  }
774 
775  // set up map with eVarCoffLaplacian key
776  varCoeffMap[StdRegions::eVarCoeffLaplacian] = varcoeff;
777 }
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.cpp:534
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:209
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:574
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:248

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

◆ UpdateTimeStepCheck()

bool Nektar::SolverUtils::UnsteadySystem::UpdateTimeStepCheck ( )
inlineprotectedvirtual

Reimplemented in Nektar::CFSImplicit.

Definition at line 204 of file UnsteadySystem.h.

205 {
206  return true;
207 }

Referenced by v_DoSolve().

◆ v_AppendOutput1D()

void Nektar::SolverUtils::UnsteadySystem::v_AppendOutput1D ( Array< OneD, Array< OneD, NekDouble >> &  solution1D)
protectedvirtual

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

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

Definition at line 614 of file UnsteadySystem.cpp.

616 {
617  // Coordinates of the quadrature points in the real physical space
618  Array<OneD, NekDouble> x(GetNpoints());
619  Array<OneD, NekDouble> y(GetNpoints());
620  Array<OneD, NekDouble> z(GetNpoints());
621  m_fields[0]->GetCoords(x, y, z);
622 
623  // Print out the solution in a txt file
624  ofstream outfile;
625  outfile.open("solution1D.txt");
626  for (int i = 0; i < GetNpoints(); i++)
627  {
628  outfile << scientific << setw(17) << setprecision(16) << x[i] << " "
629  << solution1D[0][i] << endl;
630  }
631  outfile << endl << endl;
632  outfile.close();
633 }
SOLVER_UTILS_EXPORT int GetNpoints()

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

Referenced by v_DoSolve().

◆ v_DoInitialise()

void Nektar::SolverUtils::UnsteadySystem::v_DoInitialise ( void  )
protectedvirtual

Sets up initial conditions.

Sets the initial conditions.

Reimplemented from Nektar::SolverUtils::EquationSystem.

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

Definition at line 576 of file UnsteadySystem.cpp.

577 {
582 }
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  )
protectedvirtual

Solves an unsteady problem.

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

Reimplemented from Nektar::SolverUtils::EquationSystem.

Reimplemented in Nektar::CoupledLinearNS, Nektar::MMFSWE, Nektar::PulseWaveSystem, Nektar::MMFMaxwell, and Nektar::MMFAdvection.

Definition at line 192 of file UnsteadySystem.cpp.

193 {
194  ASSERTL0(m_intScheme != 0, "No time integration scheme.");
195 
196  int i = 1;
197  int nvariables = 0;
198  int nfields = m_fields.size();
199 
200  if (m_intVariables.empty())
201  {
202  for (i = 0; i < nfields; ++i)
203  {
204  m_intVariables.push_back(i);
205  }
206  nvariables = nfields;
207  }
208  else
209  {
210  nvariables = m_intVariables.size();
211  }
212 
213  // Integrate in wave-space if using homogeneous1D
215  {
216  for (i = 0; i < nfields; ++i)
217  {
218  m_fields[i]->HomogeneousFwdTrans(m_fields[i]->GetPhys(),
219  m_fields[i]->UpdatePhys());
220  m_fields[i]->SetWaveSpace(true);
221  m_fields[i]->SetPhysState(false);
222  }
223  }
224 
225  // Set up wrapper to fields data storage.
226  Array<OneD, Array<OneD, NekDouble>> fields(nvariables);
227 
228  // Order storage to list time-integrated fields first.
229  for (i = 0; i < nvariables; ++i)
230  {
231  fields[i] = m_fields[m_intVariables[i]]->GetPhys();
232  m_fields[m_intVariables[i]]->SetPhysState(false);
233  }
234 
235  // Initialise time integration scheme
236  m_intScheme->InitializeScheme(m_timestep, fields, m_time, m_ode);
237 
238  // Initialise filters
239  for (auto &x : m_filters)
240  {
241  x.second->Initialise(m_fields, m_time);
242  }
243 
244  LibUtilities::Timer timer;
245  bool doCheckTime = false;
246  int step = m_initialStep;
247  int stepCounter = 0;
248  int restartStep = -1;
249  NekDouble intTime = 0.0;
250  NekDouble cpuTime = 0.0;
251  NekDouble cpuPrevious = 0.0;
252  NekDouble elapsed = 0.0;
253  NekDouble totFilterTime = 0.0;
254 
255  m_lastCheckTime = 0.0;
256 
257  m_TotNewtonIts = 0;
258  m_TotLinIts = 0;
259  m_TotImpStages = 0;
260 
261  Array<OneD, int> abortFlags(2, 0);
262  string abortFile = "abort";
263  if (m_session->DefinesSolverInfo("CheckAbortFile"))
264  {
265  abortFile = m_session->GetSolverInfo("CheckAbortFile");
266  }
267 
268  NekDouble tmp_cflSafetyFactor = m_cflSafetyFactor;
269 
271  while ((step < m_steps || m_time < m_fintime - NekConstants::kNekZeroTol) &&
272  abortFlags[1] == 0)
273  {
274  restartStep++;
275 
276  if (m_CFLGrowth > 1.0 && m_cflSafetyFactor < m_CFLEnd)
277  {
278  tmp_cflSafetyFactor =
279  min(m_CFLEnd, m_CFLGrowth * tmp_cflSafetyFactor);
280  }
281 
282  m_flagUpdatePreconMat = true;
283 
284  // Flag to update AV
285  m_CalcPhysicalAV = true;
286  // Frozen preconditioner checks
287  if (UpdateTimeStepCheck())
288  {
289  m_cflSafetyFactor = tmp_cflSafetyFactor;
290 
291  if (m_cflSafetyFactor)
292  {
293  m_timestep = GetTimeStep(fields);
294  }
295 
296  // Ensure that the final timestep finishes at the final
297  // time, or at a prescribed IO_CheckTime.
298  if (m_time + m_timestep > m_fintime && m_fintime > 0.0)
299  {
301  }
302  else if (m_checktime &&
304  {
307  doCheckTime = true;
308  }
309  }
310 
311  if (m_TimeIncrementFactor > 1.0)
312  {
313  NekDouble timeincrementFactor = m_TimeIncrementFactor;
314  m_timestep *= timeincrementFactor;
315 
316  if (m_time + m_timestep > m_fintime && m_fintime > 0.0)
317  {
319  }
320  }
321 
322  // Perform any solver-specific pre-integration steps
323  timer.Start();
324  if (v_PreIntegrate(step))
325  {
326  break;
327  }
328 
329  m_StagesPerStep = 0;
330  m_TotLinItePerStep = 0;
331 
332  ASSERTL0(m_timestep > 0, "m_timestep < 0");
333 
334  fields = m_intScheme->TimeIntegrate(stepCounter, m_timestep, m_ode);
335  timer.Stop();
336 
337  m_time += m_timestep;
338  elapsed = timer.TimePerTest(1);
339  intTime += elapsed;
340  cpuTime += elapsed;
341 
342  // Write out status information
343  if (m_session->GetComm()->GetRank() == 0 && !((step + 1) % m_infosteps))
344  {
345  cout << "Steps: " << setw(8) << left << step + 1 << " "
346  << "Time: " << setw(12) << left << m_time;
347 
348  if (m_cflSafetyFactor)
349  {
350  cout << " Time-step: " << setw(12) << left << m_timestep;
351  }
352 
353  stringstream ss;
354  ss << cpuTime << "s";
355  cout << " CPU Time: " << setw(8) << left << ss.str() << endl;
356  cpuPrevious = cpuTime;
357  cpuTime = 0.0;
358 
360  {
361  cout << " &&"
362  << " TotImpStages= " << m_TotImpStages
363  << " TotNewtonIts= " << m_TotNewtonIts
364  << " TotLinearIts = " << m_TotLinIts << endl;
365  }
366  }
367 
368  // Transform data into coefficient space
369  for (i = 0; i < nvariables; ++i)
370  {
371  // copy fields into ExpList::m_phys and assign the new
372  // array to fields
373  m_fields[m_intVariables[i]]->SetPhys(fields[i]);
374  fields[i] = m_fields[m_intVariables[i]]->UpdatePhys();
375  if (v_RequireFwdTrans())
376  {
377  m_fields[m_intVariables[i]]->FwdTransLocalElmt(
378  fields[i], m_fields[m_intVariables[i]]->UpdateCoeffs());
379  }
380  m_fields[m_intVariables[i]]->SetPhysState(false);
381  }
382 
383  // Perform any solver-specific post-integration steps
384  if (v_PostIntegrate(step))
385  {
386  break;
387  }
388 
389  // Check for steady-state
390  if (m_steadyStateTol > 0.0 && (!((step + 1) % m_steadyStateSteps)))
391  {
392  if (CheckSteadyState(step, intTime))
393  {
394  if (m_comm->GetRank() == 0)
395  {
396  cout << "Reached Steady State to tolerance "
397  << m_steadyStateTol << endl;
398  }
399  break;
400  }
401  }
402 
403  // test for abort conditions (nan, or abort file)
404  if (m_abortSteps && !((step + 1) % m_abortSteps))
405  {
406  abortFlags[0] = 0;
407  for (i = 0; i < nvariables; ++i)
408  {
409  if (Vmath::Nnan(fields[i].size(), fields[i], 1) > 0)
410  {
411  abortFlags[0] = 1;
412  }
413  }
414 
415  // rank zero looks for abort file and deltes it
416  // if it exists. The communicates the abort
417  if (m_session->GetComm()->GetRank() == 0)
418  {
419  if (boost::filesystem::exists(abortFile))
420  {
421  boost::filesystem::remove(abortFile);
422  abortFlags[1] = 1;
423  }
424  }
425 
426  m_session->GetComm()->AllReduce(abortFlags,
428 
429  ASSERTL0(!abortFlags[0], "NaN found during time integration.");
430  }
431 
432  // Update filters
433  for (auto &x : m_filters)
434  {
435  timer.Start();
436  x.second->Update(m_fields, m_time);
437  timer.Stop();
438  elapsed = timer.TimePerTest(1);
439  totFilterTime += elapsed;
440 
441  // Write out individual filter status information
442  if (m_session->GetComm()->GetRank() == 0 &&
443  !((step + 1) % m_filtersInfosteps) && !m_filters.empty() &&
444  m_session->DefinesCmdLineArgument("verbose"))
445  {
446  stringstream s0;
447  s0 << x.first << ":";
448  stringstream s1;
449  s1 << elapsed << "s";
450  stringstream s2;
451  s2 << elapsed / cpuPrevious * 100 << "%";
452  cout << "CPU time for filter " << setw(25) << left << s0.str()
453  << setw(12) << left << s1.str() << endl
454  << "\t Percentage of time integration: " << setw(10)
455  << left << s2.str() << endl;
456  }
457  }
458 
459  // Write out overall filter status information
460  if (m_session->GetComm()->GetRank() == 0 &&
461  !((step + 1) % m_filtersInfosteps) && !m_filters.empty())
462  {
463  stringstream ss;
464  ss << totFilterTime << "s";
465  cout << "Total filters CPU Time:\t\t\t " << setw(10) << left
466  << ss.str() << endl;
467  }
468  totFilterTime = 0.0;
469 
470  // Write out checkpoint files
471  if ((m_checksteps && !((step + 1) % m_checksteps)) || doCheckTime)
472  {
474  {
475  vector<bool> transformed(nfields, false);
476  for (i = 0; i < nfields; i++)
477  {
478  if (m_fields[i]->GetWaveSpace())
479  {
480  m_fields[i]->SetWaveSpace(false);
481  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
482  m_fields[i]->UpdatePhys());
483  m_fields[i]->SetPhysState(true);
484  transformed[i] = true;
485  }
486  }
488  m_nchk++;
489  for (i = 0; i < nfields; i++)
490  {
491  if (transformed[i])
492  {
493  m_fields[i]->SetWaveSpace(true);
494  m_fields[i]->HomogeneousFwdTrans(
495  m_fields[i]->GetPhys(), m_fields[i]->UpdatePhys());
496  m_fields[i]->SetPhysState(false);
497  }
498  }
499  }
500  else
501  {
503  m_nchk++;
504  }
505  doCheckTime = false;
506  }
507 
508  // Step advance
509  ++step;
510  ++stepCounter;
511  }
512 
513  // Print out summary statistics
514  if (m_session->GetComm()->GetRank() == 0)
515  {
516  if (m_cflSafetyFactor > 0.0)
517  {
518  cout << "CFL safety factor : " << m_cflSafetyFactor << endl
519  << "CFL time-step : " << m_timestep << endl;
520  }
521 
522  if (m_session->GetSolverInfo("Driver") != "SteadyState")
523  {
524  cout << "Time-integration : " << intTime << "s" << endl;
525  }
526 
528  {
529  cout << "-------------------------------------------" << endl
530  << "Total Implicit Stages: " << m_TotImpStages << endl
531  << "Total Newton Its : " << m_TotNewtonIts << endl
532  << "Total Linear Its : " << m_TotLinIts << endl
533  << "-------------------------------------------" << endl;
534  }
535  }
536 
537  // If homogeneous, transform back into physical space if necessary.
539  {
540  for (i = 0; i < nfields; i++)
541  {
542  if (m_fields[i]->GetWaveSpace())
543  {
544  m_fields[i]->SetWaveSpace(false);
545  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
546  m_fields[i]->UpdatePhys());
547  m_fields[i]->SetPhysState(true);
548  }
549  }
550  }
551  else
552  {
553  for (i = 0; i < nvariables; ++i)
554  {
555  m_fields[m_intVariables[i]]->SetPhys(fields[i]);
556  m_fields[m_intVariables[i]]->SetPhysState(true);
557  }
558  }
559 
560  // Finalise filters
561  for (auto &x : m_filters)
562  {
563  x.second->Finalise(m_fields, m_time);
564  }
565 
566  // Print for 1D problems
567  if (m_spacedim == 1)
568  {
569  v_AppendOutput1D(fields);
570  }
571 }
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.
SOLVER_UTILS_EXPORT NekDouble GetTimeStep()
NekDouble m_timestepMax
Time step size.
enum HomogeneousType m_HomogeneousType
int m_steps
Number of steps to take.
int m_checksteps
Number of steps between checkpoints.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
NekDouble m_CFLGrowth
CFL growth rate.
bool m_CalcPhysicalAV
flag to update artificial viscosity
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_RequireFwdTrans()
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate(int step)
NekDouble m_CFLEnd
maximun cfl in cfl growth
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate(int step)
int m_steadyStateSteps
Check for steady state at step interval.
virtual SOLVER_UTILS_EXPORT bool UpdateTimeStepCheck()
int m_filtersInfosteps
Number of time steps between outputting filters information.
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D(Array< OneD, Array< OneD, NekDouble >> &solution1D)
Print the solution at each solution point in a txt file.
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:1076

References ASSERTL0, Nektar::SolverUtils::EquationSystem::Checkpoint_Output(), CheckSteadyState(), Nektar::SolverUtils::EquationSystem::eNotHomogeneous, Nektar::SolverUtils::EquationSystem::GetTimeStep(), Nektar::NekConstants::kNekZeroTol, m_abortSteps, m_CalcPhysicalAV, 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, m_flagImplicitItsStatistics, m_flagImplicitSolver, m_flagUpdatePreconMat, Nektar::SolverUtils::EquationSystem::m_HomogeneousType, m_homoInitialFwd, 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_StagesPerStep, m_steadyStateSteps, m_steadyStateTol, Nektar::SolverUtils::EquationSystem::m_steps, Nektar::SolverUtils::EquationSystem::m_time, Nektar::SolverUtils::EquationSystem::m_TimeIncrementFactor, Nektar::SolverUtils::EquationSystem::m_timestep, Nektar::SolverUtils::EquationSystem::m_timestepMax, m_TotImpStages, m_TotLinItePerStep, m_TotLinIts, m_TotNewtonIts, Vmath::Nnan(), Nektar::LibUtilities::ReduceMax, Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), Nektar::LibUtilities::Timer::TimePerTest(), UpdateTimeStepCheck(), v_AppendOutput1D(), v_PostIntegrate(), v_PreIntegrate(), and v_RequireFwdTrans().

Referenced by Nektar::CoupledLinearNS::v_DoSolve().

◆ v_GenerateSummary()

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

Print a summary of time stepping parameters.

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

Reimplemented from Nektar::SolverUtils::EquationSystem.

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

Definition at line 588 of file UnsteadySystem.cpp.

589 {
591  AddSummaryItem(s, "Advect. advancement",
592  m_explicitAdvection ? "explicit" : "implicit");
593 
594  AddSummaryItem(s, "Diffuse. advancement",
595  m_explicitDiffusion ? "explicit" : "implicit");
596 
597  if (m_session->GetSolverInfo("EQTYPE") ==
598  "SteadyAdvectionDiffusionReaction")
599  {
600  AddSummaryItem(s, "React. advancement",
601  m_explicitReaction ? "explicit" : "implicit");
602  }
603 
604  AddSummaryItem(s, "Time Step", m_timestep);
605  AddSummaryItem(s, "No. of Steps", m_steps);
606  AddSummaryItem(s, "Checkpoints (steps)", m_checksteps);
607  AddSummaryItem(s, "Integration Type", m_intScheme->GetName());
608 }
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::IsentropicVortex::v_GenerateSummary(), Nektar::RinglebFlow::v_GenerateSummary(), Nektar::ImageWarpingSystem::v_GenerateSummary(), Nektar::VelocityCorrectionScheme::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 719 of file UnsteadySystem.cpp.

721 {
722  boost::ignore_unused(inarray);
723  NEKERROR(ErrorUtil::efatal, "Not defined for this class");
724  return 0.0;
725 }
#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)
protectedvirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::EquationSystem.

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

Definition at line 79 of file UnsteadySystem.cpp.

80 {
81  EquationSystem::v_InitObject(DeclareField);
82 
83  m_initialStep = 0;
84 
85  // Load SolverInfo parameters
86  m_session->MatchSolverInfo("DIFFUSIONADVANCEMENT", "Explicit",
87  m_explicitDiffusion, true);
88  m_session->MatchSolverInfo("ADVECTIONADVANCEMENT", "Explicit",
89  m_explicitAdvection, true);
90  m_session->MatchSolverInfo("REACTIONADVANCEMENT", "Explicit",
91  m_explicitReaction, true);
92 
93  m_session->MatchSolverInfo("FLAGIMPLICITITSSTATISTICS", "True",
95 
96  m_session->LoadParameter("CheckAbortSteps", m_abortSteps, 1);
97  // Steady state tolerance
98  m_session->LoadParameter("SteadyStateTol", m_steadyStateTol, 0.0);
99  // Frequency for checking steady state
100  m_session->LoadParameter("SteadyStateSteps", m_steadyStateSteps, 1);
101 
102  // For steady problems, we do not initialise the time integration
103  if (m_session->DefinesSolverInfo("TimeIntegrationMethod") ||
104  m_session->DefinesTimeIntScheme())
105  {
106  LibUtilities::TimeIntScheme timeInt;
107  if (m_session->DefinesTimeIntScheme())
108  {
109  timeInt = m_session->GetTimeIntScheme();
110  }
111  else
112  {
113  timeInt.method = m_session->GetSolverInfo("TimeIntegrationMethod");
114  }
115 
116  m_intScheme =
118  timeInt.method, timeInt.variant, timeInt.order,
119  timeInt.freeParams);
120 
121  // Load generic input parameters
122  m_session->LoadParameter("IO_InfoSteps", m_infosteps, 0);
123  m_session->LoadParameter("IO_FiltersInfoSteps", m_filtersInfosteps,
124  10.0 * m_infosteps);
125  m_session->LoadParameter("CFL", m_cflSafetyFactor, 0.0);
126  m_session->LoadParameter("CFLEnd", m_CFLEnd, 0.0);
127  m_session->LoadParameter("CFLGrowth", m_CFLGrowth, 1.0);
128  m_session->LoadParameter("CFLGrowth", m_CFLGrowth, 1.0);
129 
130  // Time tolerance between filter update time and time integration
131  m_session->LoadParameter("FilterTimeWarning", m_filterTimeWarning, 1);
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:144
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareFeld=true)
Initialisation object for EquationSystem.
NekDouble m_filterTimeWarning
Number of time steps between outputting status information.
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, m_filterTimeWarning, Nektar::SolverUtils::EquationSystem::m_fintime, m_flagImplicitItsStatistics, m_homoInitialFwd, 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::Bidomain::v_InitObject(), Nektar::BidomainRoth::v_InitObject(), Nektar::Monodomain::v_InitObject(), Nektar::PulseWaveSystem::v_InitObject(), Nektar::SolverUtils::AdvectionSystem::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::VelocityCorrectionScheme, Nektar::Dummy, and Nektar::SolverUtils::AdvectionSystem.

Definition at line 733 of file UnsteadySystem.cpp.

734 {
735  boost::ignore_unused(step);
736  return false;
737 }

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::IncNavierStokes, Nektar::Dummy, Nektar::UnsteadyAdvectionDiffusion, and Nektar::AcousticSystem.

Definition at line 727 of file UnsteadySystem.cpp.

728 {
729  boost::ignore_unused(step);
730  return false;
731 }

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

◆ v_RequireFwdTrans()

virtual SOLVER_UTILS_EXPORT bool Nektar::SolverUtils::UnsteadySystem::v_RequireFwdTrans ( )
inlineprotectedvirtual

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

Definition at line 178 of file UnsteadySystem.h.

179  {
180  return true;
181  }

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

882 {
883  boost::ignore_unused(step);
884  const int nPoints = GetTotPoints();
885  const int nFields = m_fields.size();
886 
887  // Holds L2 errors.
888  Array<OneD, NekDouble> RHSL2(nFields);
889  Array<OneD, NekDouble> residual(nFields);
890  Array<OneD, NekDouble> reference(nFields);
891 
892  for (int i = 0; i < nFields; ++i)
893  {
894  Array<OneD, NekDouble> tmp(nPoints);
895 
896  Vmath::Vsub(nPoints, m_fields[i]->GetPhys(), 1, m_previousSolution[i],
897  1, tmp, 1);
898  Vmath::Vmul(nPoints, tmp, 1, tmp, 1, tmp, 1);
899  residual[i] = Vmath::Vsum(nPoints, tmp, 1);
900 
901  Vmath::Vmul(nPoints, m_previousSolution[i], 1, m_previousSolution[i], 1,
902  tmp, 1);
903  reference[i] = Vmath::Vsum(nPoints, tmp, 1);
904  }
905 
906  m_comm->AllReduce(residual, LibUtilities::ReduceSum);
907  m_comm->AllReduce(reference, LibUtilities::ReduceSum);
908 
909  // L2 error
910  for (int i = 0; i < nFields; ++i)
911  {
912  reference[i] = (reference[i] == 0) ? 1 : reference[i];
913  L2[i] = sqrt(residual[i] / reference[i]);
914  }
915 }
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:895
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:419
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:291

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

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

◆ m_abortSteps

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

Number of steps between checks for abort conditions.

Definition at line 78 of file UnsteadySystem.h.

Referenced by v_DoSolve(), and v_InitObject().

◆ m_CalcPhysicalAV

bool Nektar::SolverUtils::UnsteadySystem::m_CalcPhysicalAV = true
protected

flag to update artificial viscosity

Definition at line 146 of file UnsteadySystem.h.

Referenced by v_DoSolve().

◆ m_CFLEnd

NekDouble Nektar::SolverUtils::UnsteadySystem::m_CFLEnd

maximun cfl in cfl growth

Definition at line 70 of file UnsteadySystem.h.

Referenced by v_DoSolve(), and v_InitObject().

◆ m_CFLGrowth

NekDouble Nektar::SolverUtils::UnsteadySystem::m_CFLGrowth

CFL growth rate.

Definition at line 68 of file UnsteadySystem.h.

Referenced by v_DoSolve(), and v_InitObject().

◆ m_cflNonAcoustic

NekDouble Nektar::SolverUtils::UnsteadySystem::m_cflNonAcoustic

Definition at line 66 of file UnsteadySystem.h.

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

◆ m_cflSafetyFactor

NekDouble Nektar::SolverUtils::UnsteadySystem::m_cflSafetyFactor

◆ m_epsilon

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

Definition at line 87 of file UnsteadySystem.h.

◆ m_errFile

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

Definition at line 108 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 93 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 80 of file UnsteadySystem.h.

Referenced by v_DoSolve(), and v_InitObject().

◆ m_filterTimeWarning

NekDouble Nektar::SolverUtils::UnsteadySystem::m_filterTimeWarning
protected

Number of time steps between outputting status information.

Definition at line 115 of file UnsteadySystem.h.

Referenced by v_InitObject().

◆ m_flagImplicitItsStatistics

bool Nektar::SolverUtils::UnsteadySystem::m_flagImplicitItsStatistics
protected

Definition at line 121 of file UnsteadySystem.h.

Referenced by v_DoSolve(), and v_InitObject().

◆ m_flagImplicitSolver

bool Nektar::SolverUtils::UnsteadySystem::m_flagImplicitSolver = false
protected

Definition at line 122 of file UnsteadySystem.h.

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

◆ m_flagUpdatePreconMat

bool Nektar::SolverUtils::UnsteadySystem::m_flagUpdatePreconMat
protected

Definition at line 137 of file UnsteadySystem.h.

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

◆ m_homoInitialFwd

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

◆ m_inArrayNorm

NekDouble Nektar::SolverUtils::UnsteadySystem::m_inArrayNorm = -1.0
protected

◆ m_infosteps

int Nektar::SolverUtils::UnsteadySystem::m_infosteps
protected

◆ m_intScheme

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

◆ m_intVariables

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

◆ m_locTimeStep

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

local time step(notice only for jfnk other see m_cflSafetyFactor)

Definition at line 128 of file UnsteadySystem.h.

◆ m_magnitdEstimat

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

estimate the magnitude of each conserved varibles

Definition at line 125 of file UnsteadySystem.h.

Referenced by Nektar::CFSImplicit::CalcRefValues(), and Nektar::CFSImplicit::NumCalcRiemFluxJac().

◆ m_maxLinItePerNewton

int Nektar::SolverUtils::UnsteadySystem::m_maxLinItePerNewton
protected

Definition at line 139 of file UnsteadySystem.h.

◆ m_nanSteps

int Nektar::SolverUtils::UnsteadySystem::m_nanSteps
protected

Definition at line 81 of file UnsteadySystem.h.

Referenced by Nektar::Dummy::v_InitObject().

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

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

◆ m_StagesPerStep

int Nektar::SolverUtils::UnsteadySystem::m_StagesPerStep = 1
protected

Definition at line 133 of file UnsteadySystem.h.

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

◆ m_steadyStateRes

NekDouble Nektar::SolverUtils::UnsteadySystem::m_steadyStateRes = 1.0
protected

Definition at line 102 of file UnsteadySystem.h.

Referenced by CheckSteadyState().

◆ m_steadyStateRes0

NekDouble Nektar::SolverUtils::UnsteadySystem::m_steadyStateRes0 = 1.0
protected

Definition at line 103 of file UnsteadySystem.h.

Referenced by CheckSteadyState().

◆ m_steadyStateSteps

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

Check for steady state at step interval.

Definition at line 101 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 99 of file UnsteadySystem.h.

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

◆ m_TimeIntegLambda

NekDouble Nektar::SolverUtils::UnsteadySystem::m_TimeIntegLambda = 0.0
protected

coefff of spacial derivatives(rhs or m_F in GLM) in calculating the residual of the whole equation(used in unsteady time integrations)

Definition at line 119 of file UnsteadySystem.h.

Referenced by Nektar::CFSImplicit::CalcPreconMatBRJCoeff(), Nektar::CFSImplicit::DoImplicitSolveCoeff(), Nektar::CFSImplicit::NonlinSysEvaluatorCoeff(), Nektar::CFSImplicit::PreconCoeff(), and Nektar::CFSImplicit::UpdateTimeStepCheck().

◆ m_TotImpStages

int Nektar::SolverUtils::UnsteadySystem::m_TotImpStages = 0
protected

Definition at line 143 of file UnsteadySystem.h.

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

◆ m_TotLinItePerStep

int Nektar::SolverUtils::UnsteadySystem::m_TotLinItePerStep = 0
protected

Definition at line 132 of file UnsteadySystem.h.

Referenced by v_DoSolve().

◆ m_TotLinIts

int Nektar::SolverUtils::UnsteadySystem::m_TotLinIts = 0
protected

Definition at line 142 of file UnsteadySystem.h.

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

◆ m_TotNewtonIts

int Nektar::SolverUtils::UnsteadySystem::m_TotNewtonIts = 0
protected

Definition at line 141 of file UnsteadySystem.h.

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