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

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

Protected Attributes

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

Private Member Functions

void EvaluateWaterDepth (void)
 
void EvaluateCoriolis (void)
 

Friends

class MemoryManager< ShallowWaterSystem >
 

Additional Inherited Members

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

Detailed Description

Base class for unsteady solvers.

Provides the underlying timestepping framework for 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 46 of file ShallowWaterSystem.h.

Constructor & Destructor Documentation

◆ ~ShallowWaterSystem()

Nektar::ShallowWaterSystem::~ShallowWaterSystem ( )
overridedefault

Destructor.

◆ 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

◆ AddCoriolis()

void Nektar::ShallowWaterSystem::AddCoriolis ( const Array< OneD, const Array< OneD, NekDouble > > &  physarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protected

Definition at line 369 of file ShallowWaterSystem.cpp.

372{
373
374 int ncoeffs = GetNcoeffs();
375 int nq = GetTotPoints();
376
377 Array<OneD, NekDouble> tmp(nq);
378 Array<OneD, NekDouble> mod(ncoeffs);
379
380 switch (m_projectionType)
381 {
383 {
384 // add to u equation
385 Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
386 m_fields[0]->IProductWRTBase(tmp, mod);
387 m_fields[0]->MultiplyByElmtInvMass(mod, mod);
388 m_fields[0]->BwdTrans(mod, tmp);
389 Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
390
391 // add to v equation
392 Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
393 Vmath::Neg(nq, tmp, 1);
394 m_fields[0]->IProductWRTBase(tmp, mod);
395 m_fields[0]->MultiplyByElmtInvMass(mod, mod);
396 m_fields[0]->BwdTrans(mod, tmp);
397 Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
398 }
399 break;
402 {
403 // add to u equation
404 Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
405 Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
406
407 // add to v equation
408 Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
409 Vmath::Neg(nq, tmp, 1);
410 Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
411 }
412 break;
413 default:
414 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
415 break;
416 }
417}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
Array< OneD, NekDouble > m_coriolis
Coriolis force.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetNcoeffs()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT int GetTotPoints()
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 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 ASSERTL0, Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::SolverUtils::EquationSystem::GetNcoeffs(), Nektar::SolverUtils::EquationSystem::GetTotPoints(), m_coriolis, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_projectionType, Vmath::Neg(), Vmath::Vadd(), and Vmath::Vmul().

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

◆ ConservativeToPrimitive()

void Nektar::ShallowWaterSystem::ConservativeToPrimitive ( void  )
protected

Definition at line 419 of file ShallowWaterSystem.cpp.

420{
421 int nq = GetTotPoints();
422
423 // u = hu/h
424 Vmath::Vdiv(nq, m_fields[1]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
425 m_fields[1]->UpdatePhys(), 1);
426
427 // v = hv/ v
428 Vmath::Vdiv(nq, m_fields[2]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
429 m_fields[2]->UpdatePhys(), 1);
430
431 // \eta = h - d
432 Vmath::Vsub(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
433 m_fields[0]->UpdatePhys(), 1);
434}
Array< OneD, NekDouble > m_depth
Still water depth.
void Vdiv(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:126
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.hpp:220

References Nektar::SolverUtils::EquationSystem::GetTotPoints(), m_depth, Nektar::SolverUtils::EquationSystem::m_fields, Vmath::Vdiv(), and Vmath::Vsub().

◆ CopyBoundaryTrace()

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

Definition at line 463 of file ShallowWaterSystem.cpp.

465{
466
467 int cnt = 0;
468 // loop over Boundary Regions
469 for (int bcRegion = 0; bcRegion < m_fields[0]->GetBndConditions().size();
470 ++bcRegion)
471 {
472 if (m_fields[0]
473 ->GetBndConditions()[bcRegion]
474 ->GetBoundaryConditionType() == SpatialDomains::ePeriodic)
475 {
476 continue;
477 }
478
479 // Copy the forward trace of the field to the backward trace
480 int e, id2, npts;
481
482 for (e = 0;
483 e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
484 ++e)
485 {
486 npts = m_fields[0]
487 ->GetBndCondExpansions()[bcRegion]
488 ->GetExp(e)
489 ->GetTotPoints();
490 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
491 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt +
492 e));
493
494 Vmath::Vcopy(npts, &Fwd[id2], 1, &Bwd[id2], 1);
495 }
496
497 cnt += m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
498 }
499}
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 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 52 of file ShallowWaterSystem.h.

55 {
58 pGraph);
59 p->InitObject();
60 return p;
61 }
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.

◆ DoOdeProjection()

void Nektar::ShallowWaterSystem::DoOdeProjection ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Definition at line 142 of file ShallowWaterSystem.cpp.

145{
146 int i;
147 int nvariables = inarray.size();
148
149 switch (m_projectionType)
150 {
152 {
153
154 // Just copy over array
155 if (inarray != outarray)
156 {
157 int npoints = GetNpoints();
158
159 for (i = 0; i < nvariables; ++i)
160 {
161 Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
162 }
163 }
164
165 SetBoundaryConditions(outarray, time);
166 break;
167 }
170 {
171
173 Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs(), 0.0);
174
175 for (i = 0; i < nvariables; ++i)
176 {
177 m_fields[i]->FwdTrans(inarray[i], coeffs);
178 m_fields[i]->BwdTrans(coeffs, outarray[i]);
179 }
180 break;
181 }
182 default:
183 ASSERTL0(false, "Unknown projection scheme");
184 break;
185 }
186}
void SetBoundaryConditions(Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
SOLVER_UTILS_EXPORT int GetNpoints()
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.

References ASSERTL0, Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::SolverUtils::EquationSystem::GetNcoeffs(), Nektar::SolverUtils::EquationSystem::GetNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_projectionType, SetBoundaryConditions(), Nektar::SolverUtils::EquationSystem::SetBoundaryConditions(), and Vmath::Vcopy().

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

◆ EvaluateCoriolis()

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

Definition at line 458 of file ShallowWaterSystem.cpp.

459{
460 GetFunction("Coriolis")->Evaluate("f", m_coriolis);
461}
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 453 of file ShallowWaterSystem.cpp.

454{
455 GetFunction("WaterDepth")->Evaluate("d", m_depth);
456}

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

Referenced by v_InitObject().

◆ GetDepth()

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

Definition at line 133 of file ShallowWaterSystem.h.

134 {
135 return m_depth;
136 }

References m_depth.

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

◆ GetGravity()

NekDouble Nektar::ShallowWaterSystem::GetGravity ( )
inlineprotected

Definition at line 118 of file ShallowWaterSystem.h.

119 {
120 return m_g;
121 }
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 128 of file ShallowWaterSystem.h.

129 {
130 return m_traceNormals;
131 }
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 123 of file ShallowWaterSystem.h.

124 {
125 return m_vecLocs;
126 }
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 138 of file ShallowWaterSystem.h.

139 {
140 return m_constantDepth;
141 }
bool m_constantDepth
Indicates if constant depth case.

References m_constantDepth.

◆ PrimitiveToConservative()

void Nektar::ShallowWaterSystem::PrimitiveToConservative ( void  )
protected

Definition at line 436 of file ShallowWaterSystem.cpp.

437{
438 int nq = GetTotPoints();
439
440 // h = \eta + d
441 Vmath::Vadd(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
442 m_fields[0]->UpdatePhys(), 1);
443
444 // hu = h * u
445 Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[1]->GetPhys(), 1,
446 m_fields[1]->UpdatePhys(), 1);
447
448 // hv = h * v
449 Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[2]->GetPhys(), 1,
450 m_fields[2]->UpdatePhys(), 1);
451}

References Nektar::SolverUtils::EquationSystem::GetTotPoints(), m_depth, Nektar::SolverUtils::EquationSystem::m_fields, Vmath::Vadd(), and Vmath::Vmul().

◆ SetBoundaryConditions()

void Nektar::ShallowWaterSystem::SetBoundaryConditions ( Array< OneD, Array< OneD, NekDouble > > &  physarray,
NekDouble  time 
)
protected

Definition at line 189 of file ShallowWaterSystem.cpp.

191{
192 std::string varName;
193 int nvariables = 3;
194 int cnt = 0;
195 int nTracePts = GetTraceTotPoints();
196
197 // Extract trace for boundaries. Needs to be done on all processors to avoid
198 // deadlock.
199 Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
200 for (int i = 0; i < nvariables; ++i)
201 {
202 Fwd[i] = Array<OneD, NekDouble>(nTracePts);
203 m_fields[i]->ExtractTracePhys(inarray[i], Fwd[i]);
204 }
205
206 // Loop over Boundary Regions
207 for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
208 {
209 if (m_fields[0]->GetBndConditions()[n]->GetBoundaryConditionType() ==
211 {
212 continue;
213 }
214
215 // Wall Boundary Condition
216 if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),
217 "Wall"))
218 {
219 WallBoundary2D(n, cnt, Fwd);
220 }
221
222 // Time Dependent Boundary Condition (specified in meshfile)
223 if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
224 {
225 for (int i = 0; i < nvariables; ++i)
226 {
227 varName = m_session->GetVariable(i);
228 m_fields[i]->EvaluateBoundaryConditions(time, varName);
229 }
230 }
231 cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
232 }
233}
void WallBoundary2D(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &Fwd)
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
LibUtilities::SessionReaderSharedPtr m_session
The session reader.

References Nektar::SpatialDomains::ePeriodic, Nektar::SolverUtils::EquationSystem::GetTraceTotPoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_session, and WallBoundary2D().

Referenced by DoOdeProjection().

◆ v_GenerateSummary()

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

Print a summary of time stepping parameters.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 129 of file ShallowWaterSystem.cpp.

130{
132 if (m_constantDepth == true)
133 {
134 SolverUtils::AddSummaryItem(s, "Depth", "constant");
135 }
136 else
137 {
138 SolverUtils::AddSummaryItem(s, "Depth", "variable");
139 }
140}
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:47

References Nektar::SolverUtils::AddSummaryItem(), m_constantDepth, and Nektar::SolverUtils::UnsteadySystem::v_GenerateSummary().

Referenced by Nektar::LinearSWE::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.

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 acceleration of gravity
94 m_session->LoadParameter("Gravity", m_g, 9.81);
95
96 // input/output in primitive variables
97 m_primitive = true;
98
100
101 m_constantDepth = true;
102 NekDouble depth = m_depth[0];
103 for (int i = 0; i < GetTotPoints(); ++i)
104 {
105 if (m_depth[i] != depth)
106 {
107 m_constantDepth = false;
108 break;
109 }
110 }
111
112 // Compute the bottom slopes
113 int nq = GetTotPoints();
114 if (m_constantDepth != true)
115 {
116 m_bottomSlope = Array<OneD, Array<OneD, NekDouble>>(m_spacedim);
117 for (int i = 0; i < m_spacedim; ++i)
118 {
119 m_bottomSlope[i] = Array<OneD, NekDouble>(nq);
121 m_bottomSlope[i]);
122 Vmath::Neg(nq, m_bottomSlope[i], 1);
123 }
124 }
125
127}
Array< OneD, Array< OneD, NekDouble > > m_bottomSlope
bool m_primitive
Indicates if variables are primitive or conservative.
int m_spacedim
Spatial dimension (>= expansion dim).
SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:87
double NekDouble

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

◆ WallBoundary()

void Nektar::ShallowWaterSystem::WallBoundary ( int  bcRegion,
int  cnt,
Array< OneD, Array< OneD, NekDouble > > &  Fwd,
Array< OneD, Array< OneD, NekDouble > > &  physarray 
)
protected

Definition at line 314 of file ShallowWaterSystem.cpp.

317{
318 int i;
319 int nvariables = physarray.size();
320
321 // Adjust the physical values of the trace to take
322 // user defined boundaries into account
323 int e, id1, id2, npts;
324
325 for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
326 ++e)
327 {
328 npts = m_fields[0]
329 ->GetBndCondExpansions()[bcRegion]
330 ->GetExp(e)
331 ->GetTotPoints();
332 id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
333 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
334 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
335
336 // For 2D/3D, define: v* = v - 2(v.n)n
337 Array<OneD, NekDouble> tmp(npts, 0.0);
338
339 // Calculate (v.n)
340 for (i = 0; i < m_spacedim; ++i)
341 {
342 Vmath::Vvtvp(npts, &Fwd[1 + i][id2], 1, &m_traceNormals[i][id2], 1,
343 &tmp[0], 1, &tmp[0], 1);
344 }
345
346 // Calculate 2.0(v.n)
347 Vmath::Smul(npts, -2.0, &tmp[0], 1, &tmp[0], 1);
348
349 // Calculate v* = v - 2.0(v.n)n
350 for (i = 0; i < m_spacedim; ++i)
351 {
352 Vmath::Vvtvp(npts, &tmp[0], 1, &m_traceNormals[i][id2], 1,
353 &Fwd[1 + i][id2], 1, &Fwd[1 + i][id2], 1);
354 }
355
356 // copy boundary adjusted values into the boundary expansion
357 for (i = 0; i < nvariables; ++i)
358 {
359 Vmath::Vcopy(npts, &Fwd[i][id2], 1,
360 &(m_fields[i]
361 ->GetBndCondExpansions()[bcRegion]
362 ->UpdatePhys())[id1],
363 1);
364 }
365 }
366}
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100

References Nektar::SolverUtils::EquationSystem::GetExpSize(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_spacedim, Nektar::SolverUtils::EquationSystem::m_traceNormals, Vmath::Smul(), Vmath::Vcopy(), and Vmath::Vvtvp().

◆ WallBoundary2D()

void Nektar::ShallowWaterSystem::WallBoundary2D ( int  bcRegion,
int  cnt,
Array< OneD, Array< OneD, NekDouble > > &  Fwd 
)
protected

Definition at line 235 of file ShallowWaterSystem.cpp.

237{
238 int i;
239 int nvariables = 3;
240
241 // Adjust the physical values of the trace to take
242 // user defined boundaries into account
243 int e, id1, id2, npts;
244
245 for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
246 ++e)
247 {
248 npts = m_fields[0]
249 ->GetBndCondExpansions()[bcRegion]
250 ->GetExp(e)
251 ->GetNumPoints(0);
252 id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
253 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
254 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
255
256 switch (m_expdim)
257 {
258 case 1:
259 {
260 // negate the forward flux
261 Vmath::Neg(npts, &Fwd[1][id2], 1);
262 }
263 break;
264 case 2:
265 {
266 Array<OneD, NekDouble> tmp_n(npts);
267 Array<OneD, NekDouble> tmp_t(npts);
268
269 Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[0][id2], 1,
270 &tmp_n[0], 1);
271 Vmath::Vvtvp(npts, &Fwd[2][id2], 1, &m_traceNormals[1][id2], 1,
272 &tmp_n[0], 1, &tmp_n[0], 1);
273
274 Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[1][id2], 1,
275 &tmp_t[0], 1);
276 Vmath::Vvtvm(npts, &Fwd[2][id2], 1, &m_traceNormals[0][id2], 1,
277 &tmp_t[0], 1, &tmp_t[0], 1);
278
279 // negate the normal flux
280 Vmath::Neg(npts, tmp_n, 1);
281
282 // rotate back to Cartesian
283 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[1][id2], 1,
284 &Fwd[1][id2], 1);
285 Vmath::Vvtvm(npts, &tmp_n[0], 1, &m_traceNormals[0][id2], 1,
286 &Fwd[1][id2], 1, &Fwd[1][id2], 1);
287
288 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[0][id2], 1,
289 &Fwd[2][id2], 1);
290 Vmath::Vvtvp(npts, &tmp_n[0], 1, &m_traceNormals[1][id2], 1,
291 &Fwd[2][id2], 1, &Fwd[2][id2], 1);
292 }
293 break;
294 case 3:
295 ASSERTL0(false,
296 "3D not implemented for Shallow Water Equations");
297 break;
298 default:
299 ASSERTL0(false, "Illegal expansion dimension");
300 }
301
302 // copy boundary adjusted values into the boundary expansion
303 for (i = 0; i < nvariables; ++i)
304 {
305 Vmath::Vcopy(npts, &Fwd[i][id2], 1,
306 &(m_fields[i]
307 ->GetBndCondExpansions()[bcRegion]
308 ->UpdatePhys())[id1],
309 1);
310 }
311 }
312}
int m_expdim
Expansion dimension.
void Vvtvm(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)
vvtvm (vector times vector minus vector): z = w*x - y
Definition: Vmath.hpp:381

References ASSERTL0, Nektar::SolverUtils::EquationSystem::GetExpSize(), Nektar::SolverUtils::EquationSystem::m_expdim, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_traceNormals, Vmath::Neg(), Vmath::Vcopy(), Vmath::Vmul(), Vmath::Vvtvm(), and Vmath::Vvtvp().

Referenced by SetBoundaryConditions().

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.
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 63 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

Definition at line 82 of file ShallowWaterSystem.h.

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

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