Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Member Functions | Friends | List of all members
Nektar::ShallowWaterSystem Class Reference

Base class for unsteady solvers. More...

#include <ShallowWaterSystem.h>

Inheritance diagram for Nektar::ShallowWaterSystem:
[legend]

Public Member Functions

virtual ~ShallowWaterSystem ()
 Destructor. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
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, 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 int GetPararealIterationNumber ()
 
SOLVER_UTILS_EXPORT void SetPararealIterationNumber (int num)
 
SOLVER_UTILS_EXPORT bool GetUseInitialCondition ()
 
SOLVER_UTILS_EXPORT void SetUseInitialCondition (bool 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...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to identify if operator is negated in DoSolve. More...
 
SOLVER_UTILS_EXPORT bool ParallelInTime ()
 Check if solver use Parallel-in-Time. More...
 

Static Public Member Functions

static SolverUtils::EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of class. More...
 
- Static Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
static std::string cmdSetStartTime
 
static std::string cmdSetStartChkNum
 

Protected Member Functions

 ShallowWaterSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
virtual void v_InitObject (bool DeclareFields=true) override
 Init object for UnsteadySystem class. More...
 
virtual void v_GenerateSummary (SolverUtils::SummaryList &s) override
 Print a summary of time stepping parameters. More...
 
void PrimitiveToConservative ()
 
virtual void v_PrimitiveToConservative ()
 
void ConservativeToPrimitive ()
 
virtual void v_ConservativeToPrimitive ()
 
NekDouble GetGravity ()
 
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs ()
 
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals ()
 
const Array< OneD, NekDouble > & GetDepth ()
 
bool IsConstantDepth ()
 
void CopyBoundaryTrace (const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
 
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve () override
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise () override
 Sets up initial conditions. 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 v_UpdateTimeStepCheck ()
 
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 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

SolverUtils::RiemannSolverSharedPtr m_riemannSolver
 
SolverUtils::RiemannSolverSharedPtr m_riemannSolverLDG
 
SolverUtils::AdvectionSharedPtr m_advection
 
SolverUtils::DiffusionSharedPtr m_diffusion
 
bool m_primitive
 Indicates if variables are primitive or conservative. More...
 
bool m_constantDepth
 Indicates if constant depth case. More...
 
NekDouble m_g
 Acceleration of gravity. More...
 
Array< OneD, NekDoublem_depth
 Still water depth. More...
 
Array< OneD, Array< OneD, NekDouble > > m_bottomSlope
 
Array< OneD, NekDoublem_coriolis
 Coriolis force. More...
 
Array< OneD, Array< OneD, NekDouble > > m_vecLocs
 
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
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
 
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_infosteps
 Number of time steps between outputting status information. More...
 
int m_pararealIter
 Number of parareal 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_useInitialCondition
 Flag to determine if IC are used. 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 EvaluateWaterDepth (void)
 
void EvaluateCoriolis (void)
 

Friends

class MemoryManager< ShallowWaterSystem >
 

Additional Inherited Members

- Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
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...
 
- 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 shallow water flow 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 shallow water solvers, e.g. SWE, Boussinesq, linear and nonlinear versions.

For details on implementing unsteady solvers see sectionADRSolverModuleImplementation

Definition at line 47 of file ShallowWaterSystem.h.

Constructor & Destructor Documentation

◆ ~ShallowWaterSystem()

Nektar::ShallowWaterSystem::~ShallowWaterSystem ( )
virtual

Destructor.

Definition at line 135 of file ShallowWaterSystem.cpp.

136 {
137 }

◆ ShallowWaterSystem()

Nektar::ShallowWaterSystem::ShallowWaterSystem ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
protected

Initialises UnsteadySystem class members.

Definition at line 67 of file ShallowWaterSystem.cpp.

70  : UnsteadySystem(pSession, pGraph)
71 {
72 }
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.

Member Function Documentation

◆ ConservativeToPrimitive()

void Nektar::ShallowWaterSystem::ConservativeToPrimitive ( void  )
inlineprotected

Definition at line 102 of file ShallowWaterSystem.h.

103  {
105  }

References v_ConservativeToPrimitive().

◆ CopyBoundaryTrace()

void Nektar::ShallowWaterSystem::CopyBoundaryTrace ( const Array< OneD, NekDouble > &  Fwd,
Array< OneD, NekDouble > &  Bwd 
)
protected

Definition at line 172 of file ShallowWaterSystem.cpp.

174 {
175 
176  int cnt = 0;
177  // loop over Boundary Regions
178  for (int bcRegion = 0; bcRegion < m_fields[0]->GetBndConditions().size();
179  ++bcRegion)
180  {
181  if (m_fields[0]
182  ->GetBndConditions()[bcRegion]
183  ->GetBoundaryConditionType() == SpatialDomains::ePeriodic)
184  {
185  continue;
186  }
187 
188  // Copy the forward trace of the field to the backward trace
189  int e, id2, npts;
190 
191  for (e = 0;
192  e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
193  ++e)
194  {
195  npts = m_fields[0]
196  ->GetBndCondExpansions()[bcRegion]
197  ->GetExp(e)
198  ->GetTotPoints();
199  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
200  m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt +
201  e));
202 
203  Vmath::Vcopy(npts, &Fwd[id2], 1, &Bwd[id2], 1);
204  }
205 
206  cnt += m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
207  }
208 }
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetExpSize()
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

References Nektar::SpatialDomains::ePeriodic, Nektar::SolverUtils::EquationSystem::GetExpSize(), Nektar::SolverUtils::EquationSystem::m_fields, and Vmath::Vcopy().

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

◆ create()

static SolverUtils::EquationSystemSharedPtr Nektar::ShallowWaterSystem::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
inlinestatic

Creates an instance of this class.

Definition at line 53 of file ShallowWaterSystem.h.

56  {
58  pGraph);
59  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

◆ EvaluateCoriolis()

void Nektar::ShallowWaterSystem::EvaluateCoriolis ( void  )
private

Definition at line 167 of file ShallowWaterSystem.cpp.

168 {
169  GetFunction("Coriolis")->Evaluate("f", m_coriolis);
170 }
Array< OneD, NekDouble > m_coriolis
Coriolis force.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.

References Nektar::SolverUtils::EquationSystem::GetFunction(), and m_coriolis.

Referenced by v_InitObject().

◆ EvaluateWaterDepth()

void Nektar::ShallowWaterSystem::EvaluateWaterDepth ( void  )
private

Definition at line 162 of file ShallowWaterSystem.cpp.

163 {
164  GetFunction("WaterDepth")->Evaluate("d", m_depth);
165 }
Array< OneD, NekDouble > m_depth
Still water depth.

References Nektar::SolverUtils::EquationSystem::GetFunction(), and m_depth.

Referenced by v_InitObject().

◆ GetDepth()

const Array<OneD, NekDouble>& Nektar::ShallowWaterSystem::GetDepth ( )
inlineprotected

Definition at line 123 of file ShallowWaterSystem.h.

124  {
125  return m_depth;
126  }

References m_depth.

Referenced by Nektar::NonlinearPeregrine::v_InitObject(), and Nektar::NonlinearSWE::v_InitObject().

◆ GetGravity()

NekDouble Nektar::ShallowWaterSystem::GetGravity ( )
inlineprotected

Definition at line 108 of file ShallowWaterSystem.h.

109  {
110  return m_g;
111  }
NekDouble m_g
Acceleration of gravity.

References m_g.

Referenced by Nektar::LinearSWE::v_InitObject(), Nektar::NonlinearPeregrine::v_InitObject(), and Nektar::NonlinearSWE::v_InitObject().

◆ GetNormals()

const Array<OneD, const Array<OneD, NekDouble> >& Nektar::ShallowWaterSystem::GetNormals ( )
inlineprotected

Definition at line 118 of file ShallowWaterSystem.h.

119  {
120  return m_traceNormals;
121  }
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.

References Nektar::SolverUtils::EquationSystem::m_traceNormals.

Referenced by Nektar::LinearSWE::v_InitObject(), Nektar::NonlinearPeregrine::v_InitObject(), and Nektar::NonlinearSWE::v_InitObject().

◆ GetVecLocs()

const Array<OneD, const Array<OneD, NekDouble> >& Nektar::ShallowWaterSystem::GetVecLocs ( )
inlineprotected

Definition at line 113 of file ShallowWaterSystem.h.

114  {
115  return m_vecLocs;
116  }
Array< OneD, Array< OneD, NekDouble > > m_vecLocs

References m_vecLocs.

Referenced by Nektar::LinearSWE::v_InitObject(), Nektar::NonlinearPeregrine::v_InitObject(), and Nektar::NonlinearSWE::v_InitObject().

◆ IsConstantDepth()

bool Nektar::ShallowWaterSystem::IsConstantDepth ( )
inlineprotected

Definition at line 128 of file ShallowWaterSystem.h.

129  {
130  return m_constantDepth;
131  }
bool m_constantDepth
Indicates if constant depth case.

References m_constantDepth.

◆ PrimitiveToConservative()

void Nektar::ShallowWaterSystem::PrimitiveToConservative ( void  )
inlineprotected

Definition at line 96 of file ShallowWaterSystem.h.

97  {
99  }

References v_PrimitiveToConservative().

◆ v_ConservativeToPrimitive()

void Nektar::ShallowWaterSystem::v_ConservativeToPrimitive ( )
protectedvirtual

Reimplemented in Nektar::NonlinearSWE, Nektar::NonlinearPeregrine, and Nektar::LinearSWE.

Definition at line 157 of file ShallowWaterSystem.cpp.

158 {
159  ASSERTL0(false, "This function is not implemented for this equation.");
160 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215

References ASSERTL0.

Referenced by ConservativeToPrimitive().

◆ v_GenerateSummary()

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

Print a summary of time stepping parameters.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Reimplemented in Nektar::NonlinearSWE, Nektar::NonlinearPeregrine, and Nektar::LinearSWE.

Definition at line 139 of file ShallowWaterSystem.cpp.

140 {
142  if (m_constantDepth == true)
143  {
144  SolverUtils::AddSummaryItem(s, "Depth", "constant");
145  }
146  else
147  {
148  SolverUtils::AddSummaryItem(s, "Depth", "variable");
149  }
150 }
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.
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(), m_constantDepth, and Nektar::SolverUtils::UnsteadySystem::v_GenerateSummary().

Referenced by Nektar::LinearSWE::v_GenerateSummary(), Nektar::NonlinearPeregrine::v_GenerateSummary(), and Nektar::NonlinearSWE::v_GenerateSummary().

◆ v_InitObject()

void Nektar::ShallowWaterSystem::v_InitObject ( bool  DeclareField = true)
overrideprotectedvirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Reimplemented in Nektar::NonlinearSWE, Nektar::NonlinearPeregrine, and Nektar::LinearSWE.

Definition at line 74 of file ShallowWaterSystem.cpp.

75 {
76  UnsteadySystem::v_InitObject(DeclareFields);
77 
78  // if discontinuous Galerkin determine numerical flux to use
80  {
81  ASSERTL0(m_session->DefinesSolverInfo("UPWINDTYPE"),
82  "No UPWINDTYPE defined in session.");
83  }
84 
85  // Set up locations of velocity vector.
86  m_vecLocs = Array<OneD, Array<OneD, NekDouble>>(1);
87  m_vecLocs[0] = Array<OneD, NekDouble>(m_spacedim);
88  for (int i = 0; i < m_spacedim; ++i)
89  {
90  m_vecLocs[0][i] = 1 + i;
91  }
92 
93  // Load generic input parameters
94  m_session->LoadParameter("IO_InfoSteps", m_infosteps, 0);
95 
96  // Load acceleration of gravity
97  m_session->LoadParameter("Gravity", m_g, 9.81);
98 
99  // input/output in primitive variables
100  m_primitive = true;
101 
103 
104  m_constantDepth = true;
105  NekDouble depth = m_depth[0];
106  for (int i = 0; i < GetTotPoints(); ++i)
107  {
108  if (m_depth[i] != depth)
109  {
110  m_constantDepth = false;
111  break;
112  }
113  }
114 
115  // Compute the bottom slopes
116  int nq = GetTotPoints();
117  if (m_constantDepth != true)
118  {
119  m_bottomSlope = Array<OneD, Array<OneD, NekDouble>>(m_spacedim);
120  for (int i = 0; i < m_spacedim; ++i)
121  {
122  m_bottomSlope[i] = Array<OneD, NekDouble>(nq);
124  m_bottomSlope[i]);
125  Vmath::Neg(nq, m_bottomSlope[i], 1);
126  }
127  }
128 
130 }
Array< OneD, Array< OneD, NekDouble > > m_bottomSlope
bool m_primitive
Indicates if variables are primitive or conservative.
int m_spacedim
Spatial dimension (>= expansion dim).
int m_infosteps
Number of time steps between outputting status information.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT int GetTotPoints()
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:91
double NekDouble
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518

References ASSERTL0, Nektar::MultiRegions::DirCartesianMap, Nektar::MultiRegions::eDiscontinuous, EvaluateCoriolis(), EvaluateWaterDepth(), Nektar::SolverUtils::EquationSystem::GetTotPoints(), m_bottomSlope, m_constantDepth, m_depth, Nektar::SolverUtils::EquationSystem::m_fields, m_g, Nektar::SolverUtils::EquationSystem::m_infosteps, m_primitive, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_spacedim, m_vecLocs, Vmath::Neg(), and Nektar::SolverUtils::UnsteadySystem::v_InitObject().

Referenced by Nektar::LinearSWE::v_InitObject(), Nektar::NonlinearPeregrine::v_InitObject(), and Nektar::NonlinearSWE::v_InitObject().

◆ v_PrimitiveToConservative()

void Nektar::ShallowWaterSystem::v_PrimitiveToConservative ( )
protectedvirtual

Reimplemented in Nektar::NonlinearSWE, Nektar::NonlinearPeregrine, and Nektar::LinearSWE.

Definition at line 152 of file ShallowWaterSystem.cpp.

153 {
154  ASSERTL0(false, "This function is not implemented for this equation.");
155 }

References ASSERTL0.

Referenced by PrimitiveToConservative().

Friends And Related Function Documentation

◆ MemoryManager< ShallowWaterSystem >

friend class MemoryManager< ShallowWaterSystem >
friend

Definition at line 1 of file ShallowWaterSystem.h.

Member Data Documentation

◆ className

string Nektar::ShallowWaterSystem::className
static
Initial value:
=
"ShallowWaterSystem", ShallowWaterSystem::create,
"Auxiliary functions for the shallow water system.")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
EquationSystemFactory & GetEquationSystemFactory()

Name of class.

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

Parameters
pSessionSession object to read parameters from.

Definition at line 62 of file ShallowWaterSystem.h.

◆ m_advection

SolverUtils::AdvectionSharedPtr Nektar::ShallowWaterSystem::m_advection
protected

◆ m_bottomSlope

Array<OneD, Array<OneD, NekDouble> > Nektar::ShallowWaterSystem::m_bottomSlope
protected

◆ m_constantDepth

bool Nektar::ShallowWaterSystem::m_constantDepth
protected

◆ m_coriolis

Array<OneD, NekDouble> Nektar::ShallowWaterSystem::m_coriolis
protected

◆ m_depth

Array<OneD, NekDouble> Nektar::ShallowWaterSystem::m_depth
protected

◆ m_diffusion

SolverUtils::DiffusionSharedPtr Nektar::ShallowWaterSystem::m_diffusion
protected

Definition at line 71 of file ShallowWaterSystem.h.

◆ m_g

NekDouble Nektar::ShallowWaterSystem::m_g
protected

◆ m_primitive

bool Nektar::ShallowWaterSystem::m_primitive
protected

Indicates if variables are primitive or conservative.

Definition at line 74 of file ShallowWaterSystem.h.

Referenced by v_InitObject().

◆ m_riemannSolver

SolverUtils::RiemannSolverSharedPtr Nektar::ShallowWaterSystem::m_riemannSolver
protected

◆ m_riemannSolverLDG

SolverUtils::RiemannSolverSharedPtr Nektar::ShallowWaterSystem::m_riemannSolverLDG
protected

Definition at line 69 of file ShallowWaterSystem.h.

◆ m_vecLocs

Array<OneD, Array<OneD, NekDouble> > Nektar::ShallowWaterSystem::m_vecLocs
protected

Definition at line 86 of file ShallowWaterSystem.h.

Referenced by GetVecLocs(), and v_InitObject().