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::LinearSWE Class Reference

#include <LinearSWE.h>

Inheritance diagram for Nektar::LinearSWE:
[legend]

Public Member Functions

 ~LinearSWE () 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

 LinearSWE (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
 
const Array< OneD, NekDouble > & GetDepthFwd ()
 
const Array< OneD, NekDouble > & GetDepthBwd ()
 
- Protected Member Functions inherited from Nektar::ShallowWaterSystem
 ShallowWaterSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
void v_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)
 

Protected Attributes

Array< OneD, NekDoublem_dFwd
 Still water depth traces. More...
 
Array< OneD, NekDoublem_dBwd
 
- 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...
 

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 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< LinearSWE >
 

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

Definition at line 48 of file LinearSWE.h.

Constructor & Destructor Documentation

◆ ~LinearSWE()

Nektar::LinearSWE::~LinearSWE ( )
override

Definition at line 133 of file LinearSWE.cpp.

134{
135}

◆ LinearSWE()

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

Definition at line 51 of file LinearSWE.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::LinearSWE::AddCoriolis ( const Array< OneD, const Array< OneD, NekDouble > > &  physarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
private

Definition at line 138 of file LinearSWE.cpp.

141{
142
143 int ncoeffs = GetNcoeffs();
144 int nq = GetTotPoints();
145
146 Array<OneD, NekDouble> tmp(nq);
147 Array<OneD, NekDouble> mod(ncoeffs);
148
149 switch (m_projectionType)
150 {
152 {
153 // add to u equation
154 Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
155 m_fields[0]->IProductWRTBase(tmp, mod);
156 m_fields[0]->MultiplyByElmtInvMass(mod, mod);
157 m_fields[0]->BwdTrans(mod, tmp);
158 Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
159
160 // add to v equation
161 Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
162 Vmath::Neg(nq, tmp, 1);
163 m_fields[0]->IProductWRTBase(tmp, mod);
164 m_fields[0]->MultiplyByElmtInvMass(mod, mod);
165 m_fields[0]->BwdTrans(mod, tmp);
166 Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
167 }
168 break;
171 {
172 // add to u equation
173 Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
174 Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
175
176 // add to v equation
177 Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
178 Vmath::Neg(nq, tmp, 1);
179 Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
180 }
181 break;
182 default:
183 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
184 break;
185 }
186}
#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().

◆ ConservativeToPrimitive()

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

Definition at line 553 of file LinearSWE.cpp.

556{
557 int nq = GetTotPoints();
558
559 if (physin.get() == physout.get())
560 {
561 // copy indata and work with tmp array
562 Array<OneD, Array<OneD, NekDouble>> tmp(3);
563 for (int i = 0; i < 3; ++i)
564 {
565 // deep copy
566 tmp[i] = Array<OneD, NekDouble>(nq);
567 Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
568 }
569
570 // \eta = h - d
571 Vmath::Vsub(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
572
573 // u = hu/h
574 Vmath::Vdiv(nq, tmp[1], 1, tmp[0], 1, physout[1], 1);
575
576 // v = hv/ v
577 Vmath::Vdiv(nq, tmp[2], 1, tmp[0], 1, physout[2], 1);
578 }
579 else
580 {
581 // \eta = h - d
582 Vmath::Vsub(nq, physin[0], 1, m_depth, 1, physout[0], 1);
583
584 // u = hu/h
585 Vmath::Vdiv(nq, physin[1], 1, physin[0], 1, physout[1], 1);
586
587 // v = hv/ v
588 Vmath::Vdiv(nq, physin[2], 1, physin[0], 1, physout[2], 1);
589 }
590}
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::LinearSWE::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
inlinestatic

Creates an instance of this class.

Definition at line 54 of file LinearSWE.h.

57 {
60 p->InitObject();
61 return p;
62 }
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::LinearSWE::DoOdeProjection ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Definition at line 287 of file LinearSWE.cpp.

290{
291 int i;
292 int nvariables = inarray.size();
293
294 switch (m_projectionType)
295 {
297 {
298
299 // Just copy over array
300 if (inarray != outarray)
301 {
302 int npoints = GetNpoints();
303
304 for (i = 0; i < nvariables; ++i)
305 {
306 Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
307 }
308 }
309 SetBoundaryConditions(outarray, time);
310 break;
311 }
314 {
315
317 Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs(), 0.0);
318
319 for (i = 0; i < nvariables; ++i)
320 {
321 m_fields[i]->FwdTrans(inarray[i], coeffs);
322 m_fields[i]->BwdTrans(coeffs, outarray[i]);
323 }
324 break;
325 }
326 default:
327 ASSERTL0(false, "Unknown projection scheme");
328 break;
329 }
330}
void SetBoundaryConditions(Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
Definition: LinearSWE.cpp:333
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::LinearSWE::DoOdeRhs ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Definition at line 188 of file LinearSWE.cpp.

191{
192 int i, j;
193 int ndim = m_spacedim;
194 int nvariables = inarray.size();
195 int nq = GetTotPoints();
196
197 switch (m_projectionType)
198 {
200 {
201
202 //-------------------------------------------------------
203 // Compute the DG advection including the numerical flux
204 // by using SolverUtils/Advection
205 // Input and output in physical space
206 Array<OneD, Array<OneD, NekDouble>> advVel;
207
208 m_advection->Advect(nvariables, m_fields, advVel, inarray, outarray,
209 time);
210 //-------------------------------------------------------
211
212 //-------------------------------------------------------
213 // negate the outarray since moving terms to the rhs
214 for (i = 0; i < nvariables; ++i)
215 {
216 Vmath::Neg(nq, outarray[i], 1);
217 }
218 //-------------------------------------------------------
219
220 //-------------------------------------------------
221 // Add "source terms"
222 // Input and output in physical space
223
224 // Coriolis forcing
225 if (m_coriolis.size() != 0)
226 {
227 AddCoriolis(inarray, outarray);
228 }
229 //-------------------------------------------------
230 }
231 break;
234 {
235
236 //-------------------------------------------------------
237 // Compute the fluxvector in physical space
238 Array<OneD, Array<OneD, Array<OneD, NekDouble>>> fluxvector(
239 nvariables);
240
241 for (i = 0; i < nvariables; ++i)
242 {
243 fluxvector[i] = Array<OneD, Array<OneD, NekDouble>>(ndim);
244 for (j = 0; j < ndim; ++j)
245 {
246 fluxvector[i][j] = Array<OneD, NekDouble>(nq);
247 }
248 }
249
250 LinearSWE::GetFluxVector(inarray, fluxvector);
251 //-------------------------------------------------------
252
253 //-------------------------------------------------------
254 // Take the derivative of the flux terms
255 // and negate the outarray since moving terms to the rhs
256 Array<OneD, NekDouble> tmp(nq);
257 Array<OneD, NekDouble> tmp1(nq);
258
259 for (i = 0; i < nvariables; ++i)
260 {
262 fluxvector[i][0], tmp);
264 fluxvector[i][1], tmp1);
265 Vmath::Vadd(nq, tmp, 1, tmp1, 1, outarray[i], 1);
266 Vmath::Neg(nq, outarray[i], 1);
267 }
268
269 //-------------------------------------------------
270 // Add "source terms"
271 // Input and output in physical space
272
273 // Coriolis forcing
274 if (m_coriolis.size() != 0)
275 {
276 AddCoriolis(inarray, outarray);
277 }
278 //-------------------------------------------------
279 }
280 break;
281 default:
282 ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
283 break;
284 }
285}
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Definition: LinearSWE.cpp:519
void AddCoriolis(const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: LinearSWE.cpp:138
SolverUtils::AdvectionSharedPtr m_advection
int m_spacedim
Spatial dimension (>= expansion dim).
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86

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

◆ GetDepthBwd()

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

Definition at line 100 of file LinearSWE.h.

101 {
102 return m_dBwd;
103 }
Array< OneD, NekDouble > m_dBwd
Definition: LinearSWE.h:76

References m_dBwd.

Referenced by v_InitObject().

◆ GetDepthFwd()

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

Definition at line 96 of file LinearSWE.h.

97 {
98 return m_dFwd;
99 }
Array< OneD, NekDouble > m_dFwd
Still water depth traces.
Definition: LinearSWE.h:75

References m_dFwd.

Referenced by v_InitObject().

◆ GetFluxVector()

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

Definition at line 519 of file LinearSWE.cpp.

522{
523 int i, j;
524 int nq = m_fields[0]->GetTotPoints();
525
526 NekDouble g = m_g;
527
528 // Flux vector for the mass equation
529 for (i = 0; i < m_spacedim; ++i)
530 {
531 Vmath::Vmul(nq, m_depth, 1, physfield[i + 1], 1, flux[0][i], 1);
532 }
533
534 // Put (g eta) in tmp
535 Array<OneD, NekDouble> tmp(nq);
536 Vmath::Smul(nq, g, physfield[0], 1, tmp, 1);
537
538 // Flux vector for the momentum equations
539 for (i = 0; i < m_spacedim; ++i)
540 {
541 for (j = 0; j < m_spacedim; ++j)
542 {
543 // must zero fluxes as not initialised to zero in AdvectionWeakDG
544 // ...
545 Vmath::Zero(nq, flux[i + 1][j], 1);
546 }
547
548 // Add (g eta) to appropriate field
549 Vmath::Vadd(nq, flux[i + 1][i], 1, tmp, 1, flux[i + 1][i], 1);
550 }
551}
NekDouble m_g
Acceleration of gravity.
double NekDouble
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273

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

Referenced by DoOdeRhs(), and v_InitObject().

◆ GetVelocityVector()

void Nektar::LinearSWE::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
physfieldVelocity field.
velocityVelocity field.

Definition at line 673 of file LinearSWE.cpp.

676{
677 const int npts = physfield[0].size();
678
679 for (int i = 0; i < m_spacedim; ++i)
680 {
681 Vmath::Vcopy(npts, physfield[1 + i], 1, velocity[i], 1);
682 }
683}
const std::vector< NekDouble > velocity

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

◆ NumericalFlux1D()

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

◆ NumericalFlux2D()

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

◆ PrimitiveToConservative()

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

Definition at line 609 of file LinearSWE.cpp.

612{
613
614 int nq = GetTotPoints();
615
616 if (physin.get() == physout.get())
617 {
618 // copy indata and work with tmp array
619 Array<OneD, Array<OneD, NekDouble>> tmp(3);
620 for (int i = 0; i < 3; ++i)
621 {
622 // deep copy
623 tmp[i] = Array<OneD, NekDouble>(nq);
624 Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
625 }
626
627 // h = \eta + d
628 Vmath::Vadd(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
629
630 // hu = h * u
631 Vmath::Vmul(nq, physout[0], 1, tmp[1], 1, physout[1], 1);
632
633 // hv = h * v
634 Vmath::Vmul(nq, physout[0], 1, tmp[2], 1, physout[2], 1);
635 }
636 else
637 {
638 // h = \eta + d
639 Vmath::Vadd(nq, physin[0], 1, m_depth, 1, physout[0], 1);
640
641 // hu = h * u
642 Vmath::Vmul(nq, physout[0], 1, physin[1], 1, physout[1], 1);
643
644 // hv = h * v
645 Vmath::Vmul(nq, physout[0], 1, physin[2], 1, physout[2], 1);
646 }
647}

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

◆ SetBoundaryConditions()

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

Definition at line 333 of file LinearSWE.cpp.

335{
336 std::string varName;
337 int nvariables = m_fields.size();
338 int cnt = 0;
339 int nTracePts = GetTraceTotPoints();
340
341 // Extract trace for boundaries. Needs to be done on all processors to avoid
342 // deadlock.
343 Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
344 for (int i = 0; i < nvariables; ++i)
345 {
346 Fwd[i] = Array<OneD, NekDouble>(nTracePts);
347 m_fields[i]->ExtractTracePhys(inarray[i], Fwd[i]);
348 }
349
350 // loop over Boundary Regions
351 for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
352 {
353 if (m_fields[0]->GetBndConditions()[n]->GetBoundaryConditionType() ==
355 {
356 continue;
357 }
358
359 // Wall Boundary Condition
360 if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),
361 "Wall"))
362 {
363 WallBoundary2D(n, cnt, Fwd, inarray);
364 }
365
366 // Time Dependent Boundary Condition (specified in meshfile)
367 if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
368 {
369 for (int i = 0; i < nvariables; ++i)
370 {
371 varName = m_session->GetVariable(i);
372 m_fields[i]->EvaluateBoundaryConditions(time, varName);
373 }
374 }
375 cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
376 }
377}
void WallBoundary2D(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray)
Definition: LinearSWE.cpp:437
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::LinearSWE::v_ConservativeToPrimitive ( )
overrideprotectedvirtual

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 592 of file LinearSWE.cpp.

593{
594 int nq = GetTotPoints();
595
596 // u = hu/h
597 Vmath::Vdiv(nq, m_fields[1]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
598 m_fields[1]->UpdatePhys(), 1);
599
600 // v = hv/ v
601 Vmath::Vdiv(nq, m_fields[2]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
602 m_fields[2]->UpdatePhys(), 1);
603
604 // \eta = h - d
605 Vmath::Vsub(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
606 m_fields[0]->UpdatePhys(), 1);
607}

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

◆ v_GenerateSummary()

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

Print a summary of time stepping parameters.

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

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 685 of file LinearSWE.cpp.

686{
688 if (m_session->DefinesSolverInfo("UpwindType"))
689 {
690 std::string UpwindType;
691 UpwindType = m_session->GetSolverInfo("UpwindType");
692 if (UpwindType == "LinearAverage")
693 {
694 SolverUtils::AddSummaryItem(s, "Riemann Solver", "Linear Average");
695 }
696 if (UpwindType == "LinearHLL")
697 {
698 SolverUtils::AddSummaryItem(s, "Riemann Solver", "Linear HLL");
699 }
700 }
701 SolverUtils::AddSummaryItem(s, "Variables", "eta should be in field[0]");
702 SolverUtils::AddSummaryItem(s, "", "u should be in field[1]");
703 SolverUtils::AddSummaryItem(s, "", "v should be in field[2]");
704}
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(), Nektar::SolverUtils::EquationSystem::m_session, and Nektar::ShallowWaterSystem::v_GenerateSummary().

◆ v_InitObject()

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

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 57 of file LinearSWE.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(&LinearSWE::GetFluxVector, this);
93
94 // Setting up Riemann solver for advection operator
95 m_session->LoadSolverInfo("UpwindType", riemName, "NoSolver");
96 if ((riemName == "LinearHLL") && (m_constantDepth != true))
97 {
98 ASSERTL0(false, "LinearHLL only valid for constant depth");
99 }
102 riemName, m_session);
103
104 // Setting up parameters for advection operator Riemann solver
105 m_riemannSolver->SetParam("gravity", &LinearSWE::GetGravity, this);
106 m_riemannSolver->SetAuxVec("vecLocs", &LinearSWE::GetVecLocs, this);
107 m_riemannSolver->SetVector("N", &LinearSWE::GetNormals, this);
108
109 // The numerical flux for linear SWE requires depth information
110 int nTracePointsTot = m_fields[0]->GetTrace()->GetTotPoints();
111 m_dFwd = Array<OneD, NekDouble>(nTracePointsTot);
112 m_dBwd = Array<OneD, NekDouble>(nTracePointsTot);
113 m_fields[0]->GetFwdBwdTracePhys(m_depth, m_dFwd, m_dBwd);
115 m_riemannSolver->SetScalar("depthFwd", &LinearSWE::GetDepthFwd,
116 this);
117 m_riemannSolver->SetScalar("depthBwd", &LinearSWE::GetDepthBwd,
118 this);
119
120 // Concluding initialisation of advection operators
121 m_advection->SetRiemannSolver(m_riemannSolver);
122 m_advection->InitObject(m_session, m_fields);
123 break;
124 }
125 default:
126 {
127 ASSERTL0(false, "Unsupported projection type.");
128 break;
129 }
130 }
131}
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 DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Definition: LinearSWE.cpp:188
const Array< OneD, NekDouble > & GetDepthBwd()
Definition: LinearSWE.h:100
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Definition: LinearSWE.cpp:287
const Array< OneD, NekDouble > & GetDepthFwd()
Definition: LinearSWE.h:96
void CopyBoundaryTrace(const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
bool m_constantDepth
Indicates if constant depth case.
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::ShallowWaterSystem::CopyBoundaryTrace(), 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(), GetDepthBwd(), GetDepthFwd(), GetFluxVector(), Nektar::ShallowWaterSystem::GetGravity(), Nektar::ShallowWaterSystem::GetNormals(), Nektar::SolverUtils::GetRiemannSolverFactory(), Nektar::ShallowWaterSystem::GetVecLocs(), Nektar::ShallowWaterSystem::m_advection, Nektar::ShallowWaterSystem::m_constantDepth, m_dBwd, Nektar::ShallowWaterSystem::m_depth, m_dFwd, Nektar::SolverUtils::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::LinearSWE::v_PrimitiveToConservative ( )
overrideprotectedvirtual

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 649 of file LinearSWE.cpp.

650{
651 int nq = GetTotPoints();
652
653 // h = \eta + d
654 Vmath::Vadd(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
655 m_fields[0]->UpdatePhys(), 1);
656
657 // hu = h * u
658 Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[1]->GetPhys(), 1,
659 m_fields[1]->UpdatePhys(), 1);
660
661 // hv = h * v
662 Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[2]->GetPhys(), 1,
663 m_fields[2]->UpdatePhys(), 1);
664}

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

◆ WallBoundary()

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

Wall boundary condition.

Definition at line 383 of file LinearSWE.cpp.

386{
387 int i;
388 int nvariables = physarray.size();
389
390 // Adjust the physical values of the trace to take
391 // user defined boundaries into account
392 int e, id1, id2, npts;
393
394 for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
395 ++e)
396 {
397 npts = m_fields[0]
398 ->GetBndCondExpansions()[bcRegion]
399 ->GetExp(e)
400 ->GetTotPoints();
401 id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
402 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
403 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
404
405 // For 2D/3D, define: v* = v - 2(v.n)n
406 Array<OneD, NekDouble> tmp(npts, 0.0);
407
408 // Calculate (v.n)
409 for (i = 0; i < m_spacedim; ++i)
410 {
411 Vmath::Vvtvp(npts, &Fwd[1 + i][id2], 1, &m_traceNormals[i][id2], 1,
412 &tmp[0], 1, &tmp[0], 1);
413 }
414
415 // Calculate 2.0(v.n)
416 Vmath::Smul(npts, -2.0, &tmp[0], 1, &tmp[0], 1);
417
418 // Calculate v* = v - 2.0(v.n)n
419 for (i = 0; i < m_spacedim; ++i)
420 {
421 Vmath::Vvtvp(npts, &tmp[0], 1, &m_traceNormals[i][id2], 1,
422 &Fwd[1 + i][id2], 1, &Fwd[1 + i][id2], 1);
423 }
424
425 // copy boundary adjusted values into the boundary expansion
426 for (i = 0; i < nvariables; ++i)
427 {
428 Vmath::Vcopy(npts, &Fwd[i][id2], 1,
429 &(m_fields[i]
430 ->GetBndCondExpansions()[bcRegion]
431 ->UpdatePhys())[id1],
432 1);
433 }
434 }
435}
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::LinearSWE::WallBoundary2D ( int  bcRegion,
int  cnt,
Array< OneD, Array< OneD, NekDouble > > &  Fwd,
Array< OneD, Array< OneD, NekDouble > > &  physarray 
)
private

Definition at line 437 of file LinearSWE.cpp.

440{
441
442 int i;
443 int nvariables = physarray.size();
444
445 // Adjust the physical values of the trace to take
446 // user defined boundaries into account
447 int e, id1, id2, npts;
448
449 for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
450 ++e)
451 {
452 npts = m_fields[0]
453 ->GetBndCondExpansions()[bcRegion]
454 ->GetExp(e)
455 ->GetNumPoints(0);
456 id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
457 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
458 m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
459
460 switch (m_expdim)
461 {
462 case 1:
463 {
464 // negate the forward flux
465 Vmath::Neg(npts, &Fwd[1][id2], 1);
466 }
467 break;
468 case 2:
469 {
470 Array<OneD, NekDouble> tmp_n(npts);
471 Array<OneD, NekDouble> tmp_t(npts);
472
473 Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[0][id2], 1,
474 &tmp_n[0], 1);
475 Vmath::Vvtvp(npts, &Fwd[2][id2], 1, &m_traceNormals[1][id2], 1,
476 &tmp_n[0], 1, &tmp_n[0], 1);
477
478 Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[1][id2], 1,
479 &tmp_t[0], 1);
480 Vmath::Vvtvm(npts, &Fwd[2][id2], 1, &m_traceNormals[0][id2], 1,
481 &tmp_t[0], 1, &tmp_t[0], 1);
482
483 // negate the normal flux
484 Vmath::Neg(npts, tmp_n, 1);
485
486 // rotate back to Cartesian
487 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[1][id2], 1,
488 &Fwd[1][id2], 1);
489 Vmath::Vvtvm(npts, &tmp_n[0], 1, &m_traceNormals[0][id2], 1,
490 &Fwd[1][id2], 1, &Fwd[1][id2], 1);
491
492 Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[0][id2], 1,
493 &Fwd[2][id2], 1);
494 Vmath::Vvtvp(npts, &tmp_n[0], 1, &m_traceNormals[1][id2], 1,
495 &Fwd[2][id2], 1, &Fwd[2][id2], 1);
496 }
497 break;
498 case 3:
499 ASSERTL0(false,
500 "3D not implemented for Shallow Water Equations");
501 break;
502 default:
503 ASSERTL0(false, "Illegal expansion dimension");
504 }
505
506 // copy boundary adjusted values into the boundary expansion
507 for (i = 0; i < nvariables; ++i)
508 {
509 Vmath::Vcopy(npts, &Fwd[i][id2], 1,
510 &(m_fields[i]
511 ->GetBndCondExpansions()[bcRegion]
512 ->UpdatePhys())[id1],
513 1);
514 }
515 }
516}
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< LinearSWE >

friend class MemoryManager< LinearSWE >
friend

Definition at line 1 of file LinearSWE.h.

Member Data Documentation

◆ className

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

Name of class.

Definition at line 64 of file LinearSWE.h.

◆ m_dBwd

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

Definition at line 76 of file LinearSWE.h.

Referenced by GetDepthBwd(), and v_InitObject().

◆ m_dFwd

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

Still water depth traces.

Definition at line 75 of file LinearSWE.h.

Referenced by GetDepthFwd(), and v_InitObject().