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

#include <LinearSWE.h>

Inheritance diagram for Nektar::LinearSWE:
[legend]

Public Member Functions

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

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 Member Functions inherited from Nektar::ShallowWaterSystem
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::ShallowWaterSystem
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

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

Private Attributes

Array< OneD, NekDoublem_dFwd
 Still water depth traces. More...
 
Array< OneD, NekDoublem_dBwd
 

Friends

class MemoryManager< LinearSWE >
 

Additional Inherited Members

- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D , eHomogeneous2D , eHomogeneous3D , eNotHomogeneous }
 Parameter for homogeneous expansions. More...
 
- Protected Attributes inherited from Nektar::ShallowWaterSystem
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
 
SolverUtils::AdvectionSharedPtr m_advection
 
SolverUtils::DiffusionSharedPtr m_diffusion
 
int m_TotNewtonIts = 0
 
int m_TotLinIts = 0
 
int m_TotImpStages = 0
 
NekDouble m_jacobiFreeEps = 5.0E-08
 
NekDouble m_bndEvaluateTime = 0.0
 
NekDouble m_TimeIntegLambda = 0.0
 
NekDouble m_inArrayNorm = -1.0
 
LibUtilities::NekNonlinSysIterSharedPtr m_nonlinsol
 
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
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::vector< int > m_intVariables
 
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 Maximun cfl in cfl growth. More...
 
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
std::ofstream m_errFile
 
NekDouble m_epsilon
 Diffusion coefficient. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
bool m_verbose
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh. More...
 
std::string m_sessionName
 Name of the session. More...
 
NekDouble m_time
 Current time of simulation. More...
 
int m_initialStep
 Number of the step where the simulation should begin. More...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
NekDouble m_lastCheckTime
 
NekDouble m_TimeIncrementFactor
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_infosteps
 Number of time steps between outputting status information. More...
 
int m_iterPIT = 0
 Number of parallel-in-time time iteration. More...
 
int m_windowPIT = 0
 Index of windows for parallel-in-time time iteration. More...
 
int m_spacedim
 Spatial dimension (>= expansion dim). More...
 
int m_expdim
 Expansion dimension. More...
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used. More...
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used. More...
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used. More...
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform. More...
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations. More...
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation. More...
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction. More...
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity. More...
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields. More...
 
Array< OneD, NekDoublem_movingFrameData
 Moving reference frame status in the inertial frame X, Y, Z, Theta_x, Theta_y, Theta_z, U, V, W, Omega_x, Omega_y, Omega_z, A_x, A_y, A_z, DOmega_x, DOmega_y, DOmega_z, pivot_x, pivot_y, pivot_z. More...
 
std::vector< std::string > m_strFrameData
 variable name in m_movingFrameData 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...
 
- Protected Attributes inherited from Nektar::SolverUtils::ALEHelper
Array< OneD, MultiRegions::ExpListSharedPtrm_fieldsALE
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocity
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
 
std::vector< ALEBaseShPtrm_ALEs
 
bool m_ALESolver = false
 
bool m_ImplicitALESolver = false
 
NekDouble m_prevStageTime = 0.0
 
int m_spaceDim
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 
static std::string projectionTypeLookupIds []
 

Detailed Description

Definition at line 43 of file LinearSWE.h.

Constructor & Destructor Documentation

◆ ~LinearSWE()

Nektar::LinearSWE::~LinearSWE ( )
overridedefault

◆ LinearSWE()

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

Definition at line 48 of file LinearSWE.cpp.

50 : ShallowWaterSystem(pSession, pGraph)
51{
52}
ShallowWaterSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.

Member Function Documentation

◆ CopyBoundaryTrace()

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

Definition at line 286 of file LinearSWE.cpp.

288{
289 int cnt = 0;
290
291 // loop over Boundary Regions
292 for (int bcRegion = 0; bcRegion < m_fields[0]->GetBndConditions().size();
293 ++bcRegion)
294 {
295 if (m_fields[0]
296 ->GetBndConditions()[bcRegion]
297 ->GetBoundaryConditionType() == SpatialDomains::ePeriodic)
298 {
299 continue;
300 }
301
302 // Copy the forward trace of the field to the backward trace
303 int id2, npts;
304
305 for (int e = 0;
306 e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
307 ++e)
308 {
309 npts = m_fields[0]
310 ->GetBndCondExpansions()[bcRegion]
311 ->GetExp(e)
312 ->GetTotPoints();
313 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
314 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt +
315 e));
316
317 Vmath::Vcopy(npts, &Fwd[id2], 1, &Bwd[id2], 1);
318 }
319
320 cnt += m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
321 }
322}
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.hpp:825

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

Referenced by v_InitObject().

◆ create()

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

Creates an instance of this class.

Definition at line 49 of file LinearSWE.h.

52 {
55 p->InitObject();
56 return p;
57 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ GetDepthBwd()

const Array< OneD, NekDouble > & Nektar::LinearSWE::GetDepthBwd ( )
inlineprotected

Definition at line 92 of file LinearSWE.h.

93 {
94 return m_dBwd;
95 }
Array< OneD, NekDouble > m_dBwd
Definition: LinearSWE.h:100

References m_dBwd.

Referenced by v_InitObject().

◆ GetDepthFwd()

const Array< OneD, NekDouble > & Nektar::LinearSWE::GetDepthFwd ( )
inlineprotected

Definition at line 87 of file LinearSWE.h.

88 {
89 return m_dFwd;
90 }
Array< OneD, NekDouble > m_dFwd
Still water depth traces.
Definition: LinearSWE.h:99

References m_dFwd.

Referenced by v_InitObject().

◆ GetFluxVector()

void Nektar::LinearSWE::GetFluxVector ( const Array< OneD, const Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  flux 
)
protected

Definition at line 239 of file LinearSWE.cpp.

242{
243 int nq = m_fields[0]->GetTotPoints();
244
245 NekDouble g = m_g;
246
247 // Flux vector for the mass equation
248 for (int i = 0; i < m_spacedim; ++i)
249 {
250 Vmath::Vmul(nq, m_depth, 1, physfield[i + 1], 1, flux[0][i], 1);
251 }
252
253 // Flux vector for the momentum equations
254 for (int i = 0; i < m_spacedim; ++i)
255 {
256 for (int j = 0; j < m_spacedim; ++j)
257 {
258 // must zero fluxes as not initialised to zero in AdvectionWeakDG
259 Vmath::Zero(nq, flux[i + 1][j], 1);
260 }
261
262 // Add (g eta) to appropriate field
263 Vmath::Smul(nq, g, physfield[0], 1, flux[i + 1][i], 1);
264 }
265}
NekDouble m_g
Acceleration of gravity.
Array< OneD, NekDouble > m_depth
Still water depth.
int m_spacedim
Spatial dimension (>= expansion dim).
double NekDouble
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273

References Nektar::ShallowWaterSystem::m_depth, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::ShallowWaterSystem::m_g, Nektar::SolverUtils::EquationSystem::m_spacedim, Vmath::Smul(), Vmath::Vmul(), and Vmath::Zero().

Referenced by v_DoOdeRhs(), and v_InitObject().

◆ GetVelocityVector()

void Nektar::LinearSWE::GetVelocityVector ( const Array< OneD, const Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  velocity 
)
protected

Compute the velocity field \( \mathbf{v} \) given the momentum \( h\mathbf{v} \).

Parameters
physfieldVelocity field.
velocityVelocity field.

Definition at line 274 of file LinearSWE.cpp.

277{
278 const int npts = physfield[0].size();
279
280 for (int i = 0; i < m_spacedim; ++i)
281 {
282 Vmath::Vcopy(npts, physfield[1 + i], 1, velocity[i], 1);
283 }
284}

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

◆ v_DoOdeRhs()

void Nektar::LinearSWE::v_DoOdeRhs ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
overrideprotectedvirtual

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 124 of file LinearSWE.cpp.

127{
128 int ndim = m_spacedim;
129 int nvariables = inarray.size();
130 int nq = GetTotPoints();
131
132 switch (m_projectionType)
133 {
135 {
136 //-------------------------------------------------------
137 // Compute the DG advection including the numerical flux
138 // by using SolverUtils/Advection
139 // Input and output in physical space
141 inarray, outarray, time);
142 //-------------------------------------------------------
143
144 //-------------------------------------------------------
145 // negate the outarray since moving terms to the rhs
146 for (int i = 0; i < nvariables; ++i)
147 {
148 Vmath::Neg(nq, outarray[i], 1);
149 }
150 //-------------------------------------------------------
151
152 //-------------------------------------------------
153 // Add "source terms"
154 // Input and output in physical space
155
156 // Coriolis forcing
157 if (m_coriolis.size() != 0)
158 {
159 AddCoriolis(inarray, outarray);
160 }
161 //-------------------------------------------------
162 }
163 break;
165 {
166 //-------------------------------------------------------
167 // Compute the fluxvector in physical space
168 Array<OneD, Array<OneD, Array<OneD, NekDouble>>> fluxvector(
169 nvariables);
170
171 for (int i = 0; i < nvariables; ++i)
172 {
173 fluxvector[i] = Array<OneD, Array<OneD, NekDouble>>(ndim);
174 for (int j = 0; j < ndim; ++j)
175 {
176 fluxvector[i][j] = Array<OneD, NekDouble>(nq);
177 }
178 }
179
180 LinearSWE::GetFluxVector(inarray, fluxvector);
181 //-------------------------------------------------------
182
183 //-------------------------------------------------------
184 // Take the derivative of the flux terms
185 // and negate the outarray since moving terms to the rhs
186 Array<OneD, NekDouble> tmp0(nq);
187 Array<OneD, NekDouble> tmp1(nq);
188
189 for (int i = 0; i < nvariables; ++i)
190 {
192 fluxvector[i][0], tmp0);
194 fluxvector[i][1], tmp1);
195 Vmath::Vadd(nq, tmp0, 1, tmp1, 1, outarray[i], 1);
196 Vmath::Neg(nq, outarray[i], 1);
197 }
198
199 //-------------------------------------------------
200 // Add "source terms"
201 // Input and output in physical space
202
203 // Coriolis forcing
204 if (m_coriolis.size() != 0)
205 {
206 AddCoriolis(inarray, outarray);
207 }
208 //-------------------------------------------------
209 }
210 break;
211 default:
212 ASSERTL0(false, "Unknown projection scheme for the LinearSWE");
213 break;
214 }
215}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Definition: LinearSWE.cpp:239
void AddCoriolis(const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
SolverUtils::AdvectionSharedPtr m_advection
Array< OneD, NekDouble > m_coriolis
Coriolis force.
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT int GetTotPoints()
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:87
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180

References Nektar::ShallowWaterSystem::AddCoriolis(), ASSERTL0, Nektar::MultiRegions::DirCartesianMap, Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, GetFluxVector(), Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::ShallowWaterSystem::m_advection, Nektar::ShallowWaterSystem::m_coriolis, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::SolverUtils::EquationSystem::m_spacedim, Vmath::Neg(), Nektar::NullNekDoubleArrayOfArray, and Vmath::Vadd().

Referenced by v_InitObject().

◆ v_GenerateSummary()

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

Print a summary of time stepping parameters.

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

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 217 of file LinearSWE.cpp.

218{
220 if (m_session->DefinesSolverInfo("UpwindType"))
221 {
222 std::string UpwindType;
223 UpwindType = m_session->GetSolverInfo("UpwindType");
224 if (UpwindType == "LinearAverage")
225 {
226 SolverUtils::AddSummaryItem(s, "Riemann Solver", "Linear Average");
227 }
228 else if (UpwindType == "LinearHLL")
229 {
230 SolverUtils::AddSummaryItem(s, "Riemann Solver", "Linear HLL");
231 }
232 }
233 SolverUtils::AddSummaryItem(s, "Variables", "eta should be in field[0]");
234 SolverUtils::AddSummaryItem(s, "", "u should be in field[1]");
235 SolverUtils::AddSummaryItem(s, "", "v should be in field[2]");
236}
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:47

References Nektar::SolverUtils::AddSummaryItem(), Nektar::SolverUtils::EquationSystem::m_session, and Nektar::ShallowWaterSystem::v_GenerateSummary().

◆ v_InitObject()

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

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 54 of file LinearSWE.cpp.

55{
57
58 // Type of advection class to be used
59 switch (m_projectionType)
60 {
61 // Continuous field
63 {
64 // Do nothing
65 break;
66 }
67 // Discontinuous field
69 {
70 std::string advName;
71 std::string diffName;
72 std::string riemName;
73
74 //---------------------------------------------------------------
75 // Setting up advection and diffusion operators
76 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
78 advName, advName);
79 m_advection->SetFluxVector(&LinearSWE::GetFluxVector, this);
80
81 // Setting up Riemann solver for advection operator
82 m_session->LoadSolverInfo("UpwindType", riemName, "NoSolver");
83 if ((riemName == "LinearHLL") && (m_constantDepth != true))
84 {
85 ASSERTL0(false, "LinearHLL only valid for constant depth");
86 }
89 riemName, m_session);
90
91 // Setting up parameters for advection operator Riemann solver
92 m_riemannSolver->SetParam("gravity", &LinearSWE::GetGravity, this);
93 m_riemannSolver->SetAuxVec("vecLocs", &LinearSWE::GetVecLocs, this);
94 m_riemannSolver->SetVector("N", &LinearSWE::GetNormals, this);
95
96 // The numerical flux for linear SWE requires depth information
97 int nTracePointsTot = m_fields[0]->GetTrace()->GetTotPoints();
98 m_dFwd = Array<OneD, NekDouble>(nTracePointsTot);
99 m_dBwd = Array<OneD, NekDouble>(nTracePointsTot);
100 m_fields[0]->GetFwdBwdTracePhys(m_depth, m_dFwd, m_dBwd);
102 m_riemannSolver->SetScalar("depthFwd", &LinearSWE::GetDepthFwd,
103 this);
104 m_riemannSolver->SetScalar("depthBwd", &LinearSWE::GetDepthBwd,
105 this);
106
107 // Concluding initialisation of advection operators
108 m_advection->SetRiemannSolver(m_riemannSolver);
109 m_advection->InitObject(m_session, m_fields);
110 break;
111 }
112 default:
113 {
114 ASSERTL0(false, "Unsupported projection type.");
115 break;
116 }
117 }
118
122}
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
const Array< OneD, NekDouble > & GetDepthBwd()
Definition: LinearSWE.h:92
const Array< OneD, NekDouble > & GetDepthFwd()
Definition: LinearSWE.h:87
void v_DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time) override
Definition: LinearSWE.cpp:124
void CopyBoundaryTrace(const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: LinearSWE.cpp:286
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
bool m_constantDepth
Indicates if constant depth case.
void DoImplicitSolve(const Array< OneD, const Array< OneD, NekDouble > > &inpnts, Array< OneD, Array< OneD, NekDouble > > &outpnt, const NekDouble time, const NekDouble lambda)
void v_InitObject(bool DeclareFields=true) override
Init object for UnsteadySystem class.
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs()
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:43
RiemannSolverFactory & GetRiemannSolverFactory()

References ASSERTL0, CopyBoundaryTrace(), Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineImplicitSolve(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), Nektar::ShallowWaterSystem::DoImplicitSolve(), Nektar::ShallowWaterSystem::DoOdeProjection(), Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::SolverUtils::GetAdvectionFactory(), GetDepthBwd(), GetDepthFwd(), GetFluxVector(), Nektar::ShallowWaterSystem::GetGravity(), Nektar::ShallowWaterSystem::GetNormals(), Nektar::SolverUtils::GetRiemannSolverFactory(), Nektar::ShallowWaterSystem::GetVecLocs(), Nektar::ShallowWaterSystem::m_advection, Nektar::ShallowWaterSystem::m_constantDepth, m_dBwd, Nektar::ShallowWaterSystem::m_depth, m_dFwd, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::UnsteadySystem::m_ode, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::ShallowWaterSystem::m_riemannSolver, Nektar::SolverUtils::EquationSystem::m_session, v_DoOdeRhs(), and Nektar::ShallowWaterSystem::v_InitObject().

Friends And Related Function Documentation

◆ MemoryManager< LinearSWE >

friend class MemoryManager< LinearSWE >
friend

Definition at line 1 of file LinearSWE.h.

Member Data Documentation

◆ className

std::string Nektar::LinearSWE::className
static
Initial value:
=
"LinearSWE", LinearSWE::create,
"Linear shallow water equation in primitive variables.")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
Definition: LinearSWE.h:49
EquationSystemFactory & GetEquationSystemFactory()

Name of class.

Definition at line 60 of file LinearSWE.h.

◆ m_dBwd

Array<OneD, NekDouble> Nektar::LinearSWE::m_dBwd
private

Definition at line 100 of file LinearSWE.h.

Referenced by GetDepthBwd(), and v_InitObject().

◆ m_dFwd

Array<OneD, NekDouble> Nektar::LinearSWE::m_dFwd
private

Still water depth traces.

Definition at line 99 of file LinearSWE.h.

Referenced by GetDepthFwd(), and v_InitObject().