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

#include <NonlinearSWE.h>

Inheritance diagram for Nektar::NonlinearSWE:
[legend]

Public Member Functions

 ~NonlinearSWE () override
 
- Public Member Functions inherited from Nektar::ShallowWaterSystem
 ~ShallowWaterSystem () override
 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...
 

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

 NonlinearSWE (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
void v_InitObject (bool DeclareFields=true) override
 Init object for UnsteadySystem class. More...
 
void DoOdeRhs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 
void DoOdeProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 
void GetFluxVector (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
 
void v_GenerateSummary (SolverUtils::SummaryList &s) override
 Print a summary of time stepping parameters. More...
 
void v_PrimitiveToConservative () override
 
void v_ConservativeToPrimitive () override
 
- Protected Member Functions inherited from Nektar::ShallowWaterSystem
 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 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 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 Member Functions

void NumericalFlux1D (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX)
 
void NumericalFlux2D (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
 
void SetBoundaryConditions (Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
 
void WallBoundary2D (int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray)
 
void WallBoundary (int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray)
 Wall boundary condition. More...
 
void AddCoriolis (const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
void AddVariableDepth (const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
void ConservativeToPrimitive (const Array< OneD, const Array< OneD, NekDouble > > &physin, Array< OneD, Array< OneD, NekDouble > > &physout)
 
void PrimitiveToConservative (const Array< OneD, const Array< OneD, NekDouble > > &physin, Array< OneD, Array< OneD, NekDouble > > &physout)
 
void GetVelocityVector (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)
 Compute the velocity field \( \mathbf{v} \) given the momentum \( h\mathbf{v} \). More...
 

Friends

class MemoryManager< NonlinearSWE >
 

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
 
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...
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 
static std::string projectionTypeLookupIds []
 

Detailed Description

Definition at line 47 of file NonlinearSWE.h.

Constructor & Destructor Documentation

◆ ~NonlinearSWE()

Nektar::NonlinearSWE::~NonlinearSWE ( )
override

Definition at line 121 of file NonlinearSWE.cpp.

122{
123}

◆ NonlinearSWE()

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

Definition at line 51 of file NonlinearSWE.cpp.

53 : ShallowWaterSystem(pSession, pGraph)
54{
55}
ShallowWaterSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.

Member Function Documentation

◆ AddCoriolis()

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

Definition at line 126 of file NonlinearSWE.cpp.

129{
130
131 int ncoeffs = GetNcoeffs();
132 int nq = GetTotPoints();
133
134 Array<OneD, NekDouble> tmp(nq);
135 Array<OneD, NekDouble> mod(ncoeffs);
136
137 switch (m_projectionType)
138 {
140 {
141 // add to hu equation
142 Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
143 m_fields[0]->IProductWRTBase(tmp, mod);
144 m_fields[0]->MultiplyByElmtInvMass(mod, mod);
145 m_fields[0]->BwdTrans(mod, tmp);
146 Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
147
148 // add to hv equation
149 Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
150 Vmath::Neg(nq, tmp, 1);
151 m_fields[0]->IProductWRTBase(tmp, mod);
152 m_fields[0]->MultiplyByElmtInvMass(mod, mod);
153 m_fields[0]->BwdTrans(mod, tmp);
154 Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
155 }
156 break;
159 {
160 // add to hu equation
161 Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
162 Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
163
164 // add to hv equation
165 Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
166 Vmath::Neg(nq, tmp, 1);
167 Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
168 }
169 break;
170 default:
171 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
172 break;
173 }
174}
#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(), Nektar::ShallowWaterSystem::m_coriolis, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_projectionType, Vmath::Neg(), Vmath::Vadd(), and Vmath::Vmul().

Referenced by DoOdeRhs().

◆ AddVariableDepth()

void Nektar::NonlinearSWE::AddVariableDepth ( const Array< OneD, const Array< OneD, NekDouble > > &  physarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
private

Definition at line 177 of file NonlinearSWE.cpp.

180{
181
182 int ncoeffs = GetNcoeffs();
183 int nq = GetTotPoints();
184
185 Array<OneD, NekDouble> tmp(nq);
186 Array<OneD, NekDouble> mod(ncoeffs);
187
188 switch (m_projectionType)
189 {
191 {
192 for (int i = 0; i < m_spacedim; ++i)
193 {
194 Vmath::Vmul(nq, m_bottomSlope[i], 1, physarray[0], 1, tmp, 1);
195 Vmath::Smul(nq, m_g, tmp, 1, tmp, 1);
196 m_fields[0]->IProductWRTBase(tmp, mod);
197 m_fields[0]->MultiplyByElmtInvMass(mod, mod);
198 m_fields[0]->BwdTrans(mod, tmp);
199 Vmath::Vadd(nq, tmp, 1, outarray[i + 1], 1, outarray[i + 1], 1);
200 }
201 }
202 break;
205 {
206 for (int i = 0; i < m_spacedim; ++i)
207 {
208 Vmath::Vmul(nq, m_bottomSlope[i], 1, physarray[0], 1, tmp, 1);
209 Vmath::Smul(nq, m_g, tmp, 1, tmp, 1);
210 Vmath::Vadd(nq, tmp, 1, outarray[i + 1], 1, outarray[i + 1], 1);
211 }
212 }
213 break;
214 default:
215 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
216 break;
217 }
218}
NekDouble m_g
Acceleration of gravity.
Array< OneD, Array< OneD, NekDouble > > m_bottomSlope
int m_spacedim
Spatial dimension (>= expansion dim).
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 ASSERTL0, Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::SolverUtils::EquationSystem::GetNcoeffs(), Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::ShallowWaterSystem::m_bottomSlope, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::ShallowWaterSystem::m_g, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::SolverUtils::EquationSystem::m_spacedim, Vmath::Smul(), Vmath::Vadd(), and Vmath::Vmul().

Referenced by DoOdeRhs().

◆ ConservativeToPrimitive()

void Nektar::NonlinearSWE::ConservativeToPrimitive ( const Array< OneD, const Array< OneD, NekDouble > > &  physin,
Array< OneD, Array< OneD, NekDouble > > &  physout 
)
private

Definition at line 602 of file NonlinearSWE.cpp.

605{
606 int nq = GetTotPoints();
607
608 if (physin.get() == physout.get())
609 {
610 // copy indata and work with tmp array
611 Array<OneD, Array<OneD, NekDouble>> tmp(3);
612 for (int i = 0; i < 3; ++i)
613 {
614 // deep copy
615 tmp[i] = Array<OneD, NekDouble>(nq);
616 Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
617 }
618
619 // \eta = h - d
620 Vmath::Vsub(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
621
622 // u = hu/h
623 Vmath::Vdiv(nq, tmp[1], 1, tmp[0], 1, physout[1], 1);
624
625 // v = hv/ v
626 Vmath::Vdiv(nq, tmp[2], 1, tmp[0], 1, physout[2], 1);
627 }
628 else
629 {
630 // \eta = h - d
631 Vmath::Vsub(nq, physin[0], 1, m_depth, 1, physout[0], 1);
632
633 // u = hu/h
634 Vmath::Vdiv(nq, physin[1], 1, physin[0], 1, physout[1], 1);
635
636 // v = hv/ v
637 Vmath::Vdiv(nq, physin[2], 1, physin[0], 1, physout[2], 1);
638 }
639}
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 Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
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(), Nektar::ShallowWaterSystem::m_depth, Vmath::Vcopy(), Vmath::Vdiv(), and Vmath::Vsub().

◆ create()

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

Creates an instance of this class.

Definition at line 53 of file NonlinearSWE.h.

56 {
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::NonlinearSWE::DoOdeProjection ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Definition at line 331 of file NonlinearSWE.cpp.

334{
335 int i;
336 int nvariables = inarray.size();
337
338 switch (m_projectionType)
339 {
341 {
342
343 // Just copy over array
344 if (inarray != outarray)
345 {
346 int npoints = GetNpoints();
347
348 for (i = 0; i < nvariables; ++i)
349 {
350 Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
351 }
352 }
353 SetBoundaryConditions(outarray, time);
354 break;
355 }
358 {
359
361 Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs(), 0.0);
362
363 for (i = 0; i < nvariables; ++i)
364 {
365 m_fields[i]->FwdTrans(inarray[i], coeffs);
366 m_fields[i]->BwdTrans(coeffs, outarray[i]);
367 }
368 break;
369 }
370 default:
371 ASSERTL0(false, "Unknown projection scheme");
372 break;
373 }
374}
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 v_InitObject().

◆ DoOdeRhs()

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

Definition at line 220 of file NonlinearSWE.cpp.

223{
224 int i, j;
225 int ndim = m_spacedim;
226 int nvariables = inarray.size();
227 int nq = GetTotPoints();
228
229 switch (m_projectionType)
230 {
232 {
233
234 //-------------------------------------------------------
235 // Compute the DG advection including the numerical flux
236 // by using SolverUtils/Advection
237 // Input and output in physical space
238 Array<OneD, Array<OneD, NekDouble>> advVel;
239
240 m_advection->Advect(nvariables, m_fields, advVel, inarray, outarray,
241 time);
242 //-------------------------------------------------------
243
244 //-------------------------------------------------------
245 // negate the outarray since moving terms to the rhs
246 for (i = 0; i < nvariables; ++i)
247 {
248 Vmath::Neg(nq, outarray[i], 1);
249 }
250 //-------------------------------------------------------
251
252 //-------------------------------------------------
253 // Add "source terms"
254 // Input and output in physical space
255
256 // Coriolis forcing
257 if (m_coriolis.size() != 0)
258 {
259 AddCoriolis(inarray, outarray);
260 }
261
262 // Variable Depth
263 if (m_constantDepth != true)
264 {
265 AddVariableDepth(inarray, outarray);
266 }
267 //-------------------------------------------------
268 }
269 break;
272 {
273
274 //-------------------------------------------------------
275 // Compute the fluxvector in physical space
276 Array<OneD, Array<OneD, Array<OneD, NekDouble>>> fluxvector(
277 nvariables);
278
279 for (i = 0; i < nvariables; ++i)
280 {
281 fluxvector[i] = Array<OneD, Array<OneD, NekDouble>>(ndim);
282 for (j = 0; j < ndim; ++j)
283 {
284 fluxvector[i][j] = Array<OneD, NekDouble>(nq);
285 }
286 }
287
288 NonlinearSWE::GetFluxVector(inarray, fluxvector);
289 //-------------------------------------------------------
290
291 //-------------------------------------------------------
292 // Take the derivative of the flux terms
293 // and negate the outarray since moving terms to the rhs
294 Array<OneD, NekDouble> tmp(nq);
295 Array<OneD, NekDouble> tmp1(nq);
296
297 for (i = 0; i < nvariables; ++i)
298 {
300 fluxvector[i][0], tmp);
302 fluxvector[i][1], tmp1);
303 Vmath::Vadd(nq, tmp, 1, tmp1, 1, outarray[i], 1);
304 Vmath::Neg(nq, outarray[i], 1);
305 }
306
307 //-------------------------------------------------
308 // Add "source terms"
309 // Input and output in physical space
310
311 // Coriolis forcing
312 if (m_coriolis.size() != 0)
313 {
314 AddCoriolis(inarray, outarray);
315 }
316
317 // Variable Depth
318 if (m_constantDepth != true)
319 {
320 AddVariableDepth(inarray, outarray);
321 }
322 //-------------------------------------------------
323 }
324 break;
325 default:
326 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
327 break;
328 }
329}
void AddCoriolis(const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
void AddVariableDepth(const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
SolverUtils::AdvectionSharedPtr m_advection
bool m_constantDepth
Indicates if constant depth case.
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86

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

Referenced by v_InitObject().

◆ GetFluxVector()

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

Definition at line 564 of file NonlinearSWE.cpp.

567{
568 int i, j;
569 int nq = m_fields[0]->GetTotPoints();
570
571 NekDouble g = m_g;
572 Array<OneD, Array<OneD, NekDouble>> velocity(m_spacedim);
573
574 // Flux vector for the mass equation
575 for (i = 0; i < m_spacedim; ++i)
576 {
577 velocity[i] = Array<OneD, NekDouble>(nq);
578 Vmath::Vcopy(nq, physfield[i + 1], 1, flux[0][i], 1);
579 }
580
581 GetVelocityVector(physfield, velocity);
582
583 // Put (0.5 g h h) in tmp
584 Array<OneD, NekDouble> tmp(nq);
585 Vmath::Vmul(nq, physfield[0], 1, physfield[0], 1, tmp, 1);
586 Vmath::Smul(nq, 0.5 * g, tmp, 1, tmp, 1);
587
588 // Flux vector for the momentum equations
589 for (i = 0; i < m_spacedim; ++i)
590 {
591 for (j = 0; j < m_spacedim; ++j)
592 {
593 Vmath::Vmul(nq, velocity[j], 1, physfield[i + 1], 1, flux[i + 1][j],
594 1);
595 }
596
597 // Add (0.5 g h h) to appropriate field
598 Vmath::Vadd(nq, flux[i + 1][i], 1, tmp, 1, flux[i + 1][i], 1);
599 }
600}
void GetVelocityVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)
Compute the velocity field given the momentum .
const std::vector< NekDouble > velocity
double NekDouble

References GetVelocityVector(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::ShallowWaterSystem::m_g, Nektar::SolverUtils::EquationSystem::m_spacedim, Vmath::Smul(), Vmath::Vadd(), Vmath::Vcopy(), Nektar::MovementTests::velocity, and Vmath::Vmul().

Referenced by DoOdeRhs(), and v_InitObject().

◆ GetVelocityVector()

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

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

Parameters
physfieldMomentum field.
velocityVelocity field.

Definition at line 721 of file NonlinearSWE.cpp.

724{
725 const int npts = physfield[0].size();
726
727 for (int i = 0; i < m_spacedim; ++i)
728 {
729 Vmath::Vdiv(npts, physfield[1 + i], 1, physfield[0], 1, velocity[i], 1);
730 }
731}

References Nektar::SolverUtils::EquationSystem::m_spacedim, Vmath::Vdiv(), and Nektar::MovementTests::velocity.

Referenced by GetFluxVector().

◆ NumericalFlux1D()

void Nektar::NonlinearSWE::NumericalFlux1D ( Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  numfluxX 
)
private

◆ NumericalFlux2D()

void Nektar::NonlinearSWE::NumericalFlux2D ( Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  numfluxX,
Array< OneD, Array< OneD, NekDouble > > &  numfluxY 
)
private

◆ PrimitiveToConservative()

void Nektar::NonlinearSWE::PrimitiveToConservative ( const Array< OneD, const Array< OneD, NekDouble > > &  physin,
Array< OneD, Array< OneD, NekDouble > > &  physout 
)
private

Definition at line 658 of file NonlinearSWE.cpp.

661{
662 int nq = GetTotPoints();
663
664 if (physin.get() == physout.get())
665 {
666 // copy indata and work with tmp array
667 Array<OneD, Array<OneD, NekDouble>> tmp(3);
668 for (int i = 0; i < 3; ++i)
669 {
670 // deep copy
671 tmp[i] = Array<OneD, NekDouble>(nq);
672 Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
673 }
674
675 // h = \eta + d
676 Vmath::Vadd(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
677
678 // hu = h * u
679 Vmath::Vmul(nq, physout[0], 1, tmp[1], 1, physout[1], 1);
680
681 // hv = h * v
682 Vmath::Vmul(nq, physout[0], 1, tmp[2], 1, physout[2], 1);
683 }
684 else
685 {
686 // h = \eta + d
687 Vmath::Vadd(nq, physin[0], 1, m_depth, 1, physout[0], 1);
688
689 // hu = h * u
690 Vmath::Vmul(nq, physout[0], 1, physin[1], 1, physout[1], 1);
691
692 // hv = h * v
693 Vmath::Vmul(nq, physout[0], 1, physin[2], 1, physout[2], 1);
694 }
695}

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

◆ SetBoundaryConditions()

void Nektar::NonlinearSWE::SetBoundaryConditions ( Array< OneD, Array< OneD, NekDouble > > &  physarray,
NekDouble  time 
)
private

Definition at line 377 of file NonlinearSWE.cpp.

379{
380 std::string varName;
381 int nvariables = m_fields.size();
382 int cnt = 0;
383 int nTracePts = GetTraceTotPoints();
384
385 // Extract trace for boundaries. Needs to be done on all processors to avoid
386 // deadlock.
387 Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
388 for (int i = 0; i < nvariables; ++i)
389 {
390 Fwd[i] = Array<OneD, NekDouble>(nTracePts);
391 m_fields[i]->ExtractTracePhys(inarray[i], Fwd[i]);
392 }
393
394 // Loop over Boundary Regions
395 for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
396 {
397
398 if (m_fields[0]->GetBndConditions()[n]->GetBoundaryConditionType() ==
400 {
401 continue;
402 }
403
404 // Wall Boundary Condition
405 if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),
406 "Wall"))
407 {
408 WallBoundary2D(n, cnt, Fwd, inarray);
409 }
410
411 // Time Dependent Boundary Condition (specified in meshfile)
412 if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
413 {
414 for (int i = 0; i < nvariables; ++i)
415 {
416 varName = m_session->GetVariable(i);
417 m_fields[i]->EvaluateBoundaryConditions(time, varName);
418 }
419 }
420 cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
421 }
422}
void WallBoundary2D(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray)
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_ConservativeToPrimitive()

void Nektar::NonlinearSWE::v_ConservativeToPrimitive ( )
overrideprotectedvirtual

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 641 of file NonlinearSWE.cpp.

642{
643 int nq = GetTotPoints();
644
645 // u = hu/h
646 Vmath::Vdiv(nq, m_fields[1]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
647 m_fields[1]->UpdatePhys(), 1);
648
649 // v = hv/ v
650 Vmath::Vdiv(nq, m_fields[2]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
651 m_fields[2]->UpdatePhys(), 1);
652
653 // \eta = h - d
654 Vmath::Vsub(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
655 m_fields[0]->UpdatePhys(), 1);
656}

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

◆ v_GenerateSummary()

void Nektar::NonlinearSWE::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 733 of file NonlinearSWE.cpp.

734{
736 SolverUtils::AddSummaryItem(s, "Variables", "h should be in field[0]");
737 SolverUtils::AddSummaryItem(s, "", "hu should be in field[1]");
738 SolverUtils::AddSummaryItem(s, "", "hv should be in field[2]");
739}
void v_GenerateSummary(SolverUtils::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(), and Nektar::ShallowWaterSystem::v_GenerateSummary().

◆ v_InitObject()

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

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 57 of file NonlinearSWE.cpp.

58{
60
62 {
65 }
66 else
67 {
68 ASSERTL0(false, "Implicit SWE not set up.");
69 }
70
71 // Type of advection class to be used
72 switch (m_projectionType)
73 {
74 // Continuous field
76 {
77 // Do nothing
78 break;
79 }
80 // Discontinuous field
82 {
83 string advName;
84 string diffName;
85 string riemName;
86
87 //---------------------------------------------------------------
88 // Setting up advection and diffusion operators
89 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
91 advName, advName);
92 m_advection->SetFluxVector(&NonlinearSWE::GetFluxVector, this);
93
94 // Setting up Riemann solver for advection operator
95 m_session->LoadSolverInfo("UpwindType", riemName, "Average");
98 riemName, m_session);
99
100 // Setting up parameters for advection operator Riemann solver
101 m_riemannSolver->SetParam("gravity", &NonlinearSWE::GetGravity,
102 this);
103 m_riemannSolver->SetAuxVec("vecLocs", &NonlinearSWE::GetVecLocs,
104 this);
105 m_riemannSolver->SetVector("N", &NonlinearSWE::GetNormals, this);
106 m_riemannSolver->SetScalar("depth", &NonlinearSWE::GetDepth, this);
107
108 // Concluding initialisation of advection operators
109 m_advection->SetRiemannSolver(m_riemannSolver);
110 m_advection->InitObject(m_session, m_fields);
111 break;
112 }
113 default:
114 {
115 ASSERTL0(false, "Unsupported projection type.");
116 break;
117 }
118 }
119}
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
const Array< OneD, NekDouble > & GetDepth()
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.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:43
RiemannSolverFactory & GetRiemannSolverFactory()

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), DoOdeProjection(), DoOdeRhs(), Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::SolverUtils::GetAdvectionFactory(), Nektar::ShallowWaterSystem::GetDepth(), GetFluxVector(), Nektar::ShallowWaterSystem::GetGravity(), Nektar::ShallowWaterSystem::GetNormals(), Nektar::SolverUtils::GetRiemannSolverFactory(), Nektar::ShallowWaterSystem::GetVecLocs(), Nektar::ShallowWaterSystem::m_advection, Nektar::SolverUtils::UnsteadySystem::m_explicitAdvection, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::UnsteadySystem::m_ode, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::ShallowWaterSystem::m_riemannSolver, Nektar::SolverUtils::EquationSystem::m_session, and Nektar::ShallowWaterSystem::v_InitObject().

◆ v_PrimitiveToConservative()

void Nektar::NonlinearSWE::v_PrimitiveToConservative ( )
overrideprotectedvirtual

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 697 of file NonlinearSWE.cpp.

698{
699 int nq = GetTotPoints();
700
701 // h = \eta + d
702 Vmath::Vadd(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
703 m_fields[0]->UpdatePhys(), 1);
704
705 // hu = h * u
706 Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[1]->GetPhys(), 1,
707 m_fields[1]->UpdatePhys(), 1);
708
709 // hv = h * v
710 Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[2]->GetPhys(), 1,
711 m_fields[2]->UpdatePhys(), 1);
712}

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

◆ WallBoundary()

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

Wall boundary condition.

Definition at line 428 of file NonlinearSWE.cpp.

431{
432 int i;
433 int nvariables = physarray.size();
434
435 // Adjust the physical values of the trace to take
436 // user defined boundaries into account
437 int e, id1, id2, npts;
438
439 for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
440 ++e)
441 {
442 npts = m_fields[0]
443 ->GetBndCondExpansions()[bcRegion]
444 ->GetExp(e)
445 ->GetTotPoints();
446 id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
447 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
448 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
449
450 // For 2D/3D, define: v* = v - 2(v.n)n
451 Array<OneD, NekDouble> tmp(npts, 0.0);
452
453 // Calculate (v.n)
454 for (i = 0; i < m_spacedim; ++i)
455 {
456 Vmath::Vvtvp(npts, &Fwd[1 + i][id2], 1, &m_traceNormals[i][id2], 1,
457 &tmp[0], 1, &tmp[0], 1);
458 }
459
460 // Calculate 2.0(v.n)
461 Vmath::Smul(npts, -2.0, &tmp[0], 1, &tmp[0], 1);
462
463 // Calculate v* = v - 2.0(v.n)n
464 for (i = 0; i < m_spacedim; ++i)
465 {
466 Vmath::Vvtvp(npts, &tmp[0], 1, &m_traceNormals[i][id2], 1,
467 &Fwd[1 + i][id2], 1, &Fwd[1 + i][id2], 1);
468 }
469
470 // copy boundary adjusted values into the boundary expansion
471 for (i = 0; i < nvariables; ++i)
472 {
473 Vmath::Vcopy(npts, &Fwd[i][id2], 1,
474 &(m_fields[i]
475 ->GetBndCondExpansions()[bcRegion]
476 ->UpdatePhys())[id1],
477 1);
478 }
479 }
480}
SOLVER_UTILS_EXPORT int GetExpSize()
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
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

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::NonlinearSWE::WallBoundary2D ( int  bcRegion,
int  cnt,
Array< OneD, Array< OneD, NekDouble > > &  Fwd,
Array< OneD, Array< OneD, NekDouble > > &  physarray 
)
private

Definition at line 482 of file NonlinearSWE.cpp.

485{
486
487 int i;
488 int nvariables = physarray.size();
489
490 // Adjust the physical values of the trace to take
491 // user defined boundaries into account
492 int e, id1, id2, npts;
493
494 for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
495 ++e)
496 {
497 npts = m_fields[0]
498 ->GetBndCondExpansions()[bcRegion]
499 ->GetExp(e)
500 ->GetNumPoints(0);
501 id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
502 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
503 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
504
505 switch (m_expdim)
506 {
507 case 1:
508 {
509 // negate the forward flux
510 Vmath::Neg(npts, &Fwd[1][id2], 1);
511 }
512 break;
513 case 2:
514 {
515 Array<OneD, NekDouble> tmp_n(npts);
516 Array<OneD, NekDouble> tmp_t(npts);
517
518 Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[0][id2], 1,
519 &tmp_n[0], 1);
520 Vmath::Vvtvp(npts, &Fwd[2][id2], 1, &m_traceNormals[1][id2], 1,
521 &tmp_n[0], 1, &tmp_n[0], 1);
522
523 Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[1][id2], 1,
524 &tmp_t[0], 1);
525 Vmath::Vvtvm(npts, &Fwd[2][id2], 1, &m_traceNormals[0][id2], 1,
526 &tmp_t[0], 1, &tmp_t[0], 1);
527
528 // negate the normal flux
529 Vmath::Neg(npts, tmp_n, 1);
530
531 // rotate back to Cartesian
532 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[1][id2], 1,
533 &Fwd[1][id2], 1);
534 Vmath::Vvtvm(npts, &tmp_n[0], 1, &m_traceNormals[0][id2], 1,
535 &Fwd[1][id2], 1, &Fwd[1][id2], 1);
536
537 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[0][id2], 1,
538 &Fwd[2][id2], 1);
539 Vmath::Vvtvp(npts, &tmp_n[0], 1, &m_traceNormals[1][id2], 1,
540 &Fwd[2][id2], 1, &Fwd[2][id2], 1);
541 }
542 break;
543 case 3:
544 ASSERTL0(false,
545 "3D not implemented for Shallow Water Equations");
546 break;
547 default:
548 ASSERTL0(false, "Illegal expansion dimension");
549 }
550
551 // copy boundary adjusted values into the boundary expansion
552 for (i = 0; i < nvariables; ++i)
553 {
554 Vmath::Vcopy(npts, &Fwd[i][id2], 1,
555 &(m_fields[i]
556 ->GetBndCondExpansions()[bcRegion]
557 ->UpdatePhys())[id1],
558 1);
559 }
560 }
561}
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< NonlinearSWE >

friend class MemoryManager< NonlinearSWE >
friend

Definition at line 1 of file NonlinearSWE.h.

Member Data Documentation

◆ className

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

Name of class.

Definition at line 63 of file NonlinearSWE.h.