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

virtual ~NonlinearSWE ()
 
- Public Member Functions inherited from Nektar::ShallowWaterSystem
virtual ~ShallowWaterSystem ()
 Destructor. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep (const Array< OneD, const Array< OneD, NekDouble >> &inarray)
 Calculate the larger time-step mantaining the problem stable. More...
 
SOLVER_UTILS_EXPORT void SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void SetUpTraceNormals (void)
 
SOLVER_UTILS_EXPORT void InitObject (bool DeclareField=true)
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise ()
 Perform any initialisation necessary before solving the problem. More...
 
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem. More...
 
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space. More...
 
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space. More...
 
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve. More...
 
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. More...
 
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name. More...
 
template<class T >
std::shared_ptr< T > as ()
 
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name. More...
 
SOLVER_UTILS_EXPORT LibUtilities::SessionReaderSharedPtr GetSession ()
 Get Session name. More...
 
SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr GetPressure ()
 Get pressure field if available. More...
 
SOLVER_UTILS_EXPORT void ExtraFldOutput (std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables)
 
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics. More...
 
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda. More...
 
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction (std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
 Get a SessionFunction by name. More...
 
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields. More...
 
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the L2 error between fields and a given exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf]. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables)
 Write checkpoint file of custom data fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables)
 Write input fields to the given filename. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFldToMultiDomains (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int ndomains)
 Input field data from the given file to multiple domains. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, std::vector< std::string > &fieldStr, Array< OneD, Array< OneD, NekDouble >> &coeffs)
 Output a field. Input field data into array from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, MultiRegions::ExpListSharedPtr &pField, std::string &pFieldName)
 Output a field. Input field data into ExpList from the given file. More...
 
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary. More...
 
SOLVER_UTILS_EXPORT Array< OneD, MultiRegions::ExpListSharedPtr > & UpdateFields ()
 
SOLVER_UTILS_EXPORT LibUtilities::FieldMetaDataMapUpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated. More...
 
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time. More...
 
SOLVER_UTILS_EXPORT int GetNcoeffs ()
 
SOLVER_UTILS_EXPORT int GetNcoeffs (const int eid)
 
SOLVER_UTILS_EXPORT int GetNumExpModes ()
 
SOLVER_UTILS_EXPORT const Array< OneD, int > GetNumExpModesPerExp ()
 
SOLVER_UTILS_EXPORT int GetNvariables ()
 
SOLVER_UTILS_EXPORT const std::string GetVariable (unsigned int i)
 
SOLVER_UTILS_EXPORT int GetTraceTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTraceNpoints ()
 
SOLVER_UTILS_EXPORT int GetExpSize ()
 
SOLVER_UTILS_EXPORT int GetPhys_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetCoeff_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTotPoints (int n)
 
SOLVER_UTILS_EXPORT int GetNpoints ()
 
SOLVER_UTILS_EXPORT int GetSteps ()
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep ()
 
SOLVER_UTILS_EXPORT void CopyFromPhysField (const int i, Array< OneD, NekDouble > &output)
 
SOLVER_UTILS_EXPORT void CopyToPhysField (const int i, const Array< OneD, const NekDouble > &input)
 
SOLVER_UTILS_EXPORT void SetSteps (const int steps)
 
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
 
SOLVER_UTILS_EXPORT void FwdTransFields ()
 
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
 
SOLVER_UTILS_EXPORT int GetCheckpointNumber ()
 
SOLVER_UTILS_EXPORT void SetCheckpointNumber (int num)
 
SOLVER_UTILS_EXPORT int GetCheckpointSteps ()
 
SOLVER_UTILS_EXPORT void SetCheckpointSteps (int num)
 
SOLVER_UTILS_EXPORT int GetInfoSteps ()
 
SOLVER_UTILS_EXPORT void SetInfoSteps (int num)
 
SOLVER_UTILS_EXPORT int GetPararealIterationNumber ()
 
SOLVER_UTILS_EXPORT void SetPararealIterationNumber (int num)
 
SOLVER_UTILS_EXPORT bool GetUseInitialCondition ()
 
SOLVER_UTILS_EXPORT void SetUseInitialCondition (bool num)
 
SOLVER_UTILS_EXPORT Array< OneD, const Array< OneD, NekDouble > > GetTraceNormals ()
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetTimeStep (const NekDouble timestep)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to identify if operator is negated in DoSolve. More...
 
SOLVER_UTILS_EXPORT bool ParallelInTime ()
 Check if solver use Parallel-in-Time. More...
 

Static Public Member Functions

static SolverUtils::EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class. More...
 
- Static Public 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)
 
virtual 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)
 
virtual void v_GenerateSummary (SolverUtils::SummaryList &s) override
 Print a summary of time stepping parameters. More...
 
virtual void v_PrimitiveToConservative () override
 
virtual 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 PrimitiveToConservative ()
 
void ConservativeToPrimitive ()
 
NekDouble GetGravity ()
 
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs ()
 
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals ()
 
const Array< OneD, NekDouble > & GetDepth ()
 
bool IsConstantDepth ()
 
void CopyBoundaryTrace (const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
 
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve () override
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise () override
 Sets up initial conditions. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble >> &inarray)
 Return the timestep to be used for the next step in the time-marching loop. More...
 
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_RequireFwdTrans ()
 
virtual SOLVER_UTILS_EXPORT void v_SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time, int &nchk)
 
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff (const Array< OneD, Array< OneD, NekDouble >> vel, StdRegions::VarCoeffMap &varCoeffMap)
 Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h time velocity. More...
 
virtual SOLVER_UTILS_EXPORT bool v_UpdateTimeStepCheck ()
 
SOLVER_UTILS_EXPORT void DoDummyProjection (const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
 Perform dummy projection. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises EquationSystem class members. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space. More...
 
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
virtual SOLVER_UTILS_EXPORT void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure (void)
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables)
 

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

- Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
NekDouble m_cflNonAcoustic
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 maximun cfl in cfl growth More...
 
- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D , eHomogeneous2D , eHomogeneous3D , eNotHomogeneous }
 Parameter for homogeneous expansions. More...
 
- Protected Attributes inherited from Nektar::ShallowWaterSystem
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
 
SolverUtils::RiemannSolverSharedPtr m_riemannSolverLDG
 
SolverUtils::AdvectionSharedPtr m_advection
 
SolverUtils::DiffusionSharedPtr m_diffusion
 
bool m_primitive
 Indicates if variables are primitive or conservative. More...
 
bool m_constantDepth
 Indicates if constant depth case. More...
 
NekDouble m_g
 Acceleration of gravity. More...
 
Array< OneD, NekDoublem_depth
 Still water depth. More...
 
Array< OneD, Array< OneD, NekDouble > > m_bottomSlope
 
Array< OneD, NekDoublem_coriolis
 Coriolis force. More...
 
Array< OneD, Array< OneD, NekDouble > > m_vecLocs
 
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
int m_nanSteps
 
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
NekDouble m_epsilon
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateRes = 1.0
 
NekDouble m_steadyStateRes0 = 1.0
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::ofstream m_errFile
 
std::vector< int > m_intVariables
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
NekDouble m_filterTimeWarning
 Number of time steps between outputting status information. More...
 
NekDouble m_TimeIntegLambda = 0.0
 coefff of spacial derivatives(rhs or m_F in GLM) in calculating the residual of the whole equation(used in unsteady time integrations) More...
 
bool m_flagImplicitItsStatistics
 
bool m_flagImplicitSolver = false
 
Array< OneD, NekDoublem_magnitdEstimat
 estimate the magnitude of each conserved varibles More...
 
Array< OneD, NekDoublem_locTimeStep
 local time step(notice only for jfnk other see m_cflSafetyFactor) More...
 
NekDouble m_inArrayNorm = -1.0
 
int m_TotLinItePerStep = 0
 
int m_StagesPerStep = 1
 
bool m_flagUpdatePreconMat
 
int m_maxLinItePerNewton
 
int m_TotNewtonIts = 0
 
int m_TotLinIts = 0
 
int m_TotImpStages = 0
 
bool m_CalcPhysicalAV = true
 flag to update artificial viscosity More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
bool m_verbose
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh. More...
 
std::string m_sessionName
 Name of the session. More...
 
NekDouble m_time
 Current time of simulation. More...
 
int m_initialStep
 Number of the step where the simulation should begin. More...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_timestepMax = -1.0
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
NekDouble m_lastCheckTime
 
NekDouble m_TimeIncrementFactor
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_infosteps
 Number of time steps between outputting status information. More...
 
int m_pararealIter
 Number of parareal time iteration. More...
 
int m_spacedim
 Spatial dimension (>= expansion dim). More...
 
int m_expdim
 Expansion dimension. More...
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used. More...
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used. More...
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used. More...
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform. More...
 
bool m_useInitialCondition
 Flag to determine if IC are used. More...
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations. More...
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation. More...
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction. More...
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity. More...
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields. More...
 
Array< OneD, NekDoublem_movingFrameVelsxyz
 Moving frame of reference velocities. More...
 
Array< OneD, NekDoublem_movingFrameTheta
 Moving frame of reference angles with respect to the. More...
 
boost::numeric::ublas::matrix< NekDoublem_movingFrameProjMat
 Projection matrix for transformation between inertial and moving. More...
 
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error. More...
 
enum HomogeneousType m_HomogeneousType
 
NekDouble m_LhomX
 physical length in X direction (if homogeneous) More...
 
NekDouble m_LhomY
 physical length in Y direction (if homogeneous) More...
 
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous) More...
 
int m_npointsX
 number of points in X direction (if homogeneous) More...
 
int m_npointsY
 number of points in Y direction (if homogeneous) More...
 
int m_npointsZ
 number of points in Z direction (if homogeneous) More...
 
int m_HomoDirec
 number of homogenous directions More...
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 

Detailed Description

Definition at line 48 of file NonlinearSWE.h.

Constructor & Destructor Documentation

◆ ~NonlinearSWE()

Nektar::NonlinearSWE::~NonlinearSWE ( )
virtual

Definition at line 147 of file NonlinearSWE.cpp.

148 {
149 }

◆ 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 152 of file NonlinearSWE.cpp.

155 {
156 
157  int ncoeffs = GetNcoeffs();
158  int nq = GetTotPoints();
159 
160  Array<OneD, NekDouble> tmp(nq);
161  Array<OneD, NekDouble> mod(ncoeffs);
162 
163  switch (m_projectionType)
164  {
166  {
167  // add to hu equation
168  Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
169  m_fields[0]->IProductWRTBase(tmp, mod);
170  m_fields[0]->MultiplyByElmtInvMass(mod, mod);
171  m_fields[0]->BwdTrans(mod, tmp);
172  Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
173 
174  // add to hv equation
175  Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
176  Vmath::Neg(nq, tmp, 1);
177  m_fields[0]->IProductWRTBase(tmp, mod);
178  m_fields[0]->MultiplyByElmtInvMass(mod, mod);
179  m_fields[0]->BwdTrans(mod, tmp);
180  Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
181  }
182  break;
185  {
186  // add to hu equation
187  Vmath::Vmul(nq, m_coriolis, 1, physarray[2], 1, tmp, 1);
188  Vmath::Vadd(nq, tmp, 1, outarray[1], 1, outarray[1], 1);
189 
190  // add to hv equation
191  Vmath::Vmul(nq, m_coriolis, 1, physarray[1], 1, tmp, 1);
192  Vmath::Neg(nq, tmp, 1);
193  Vmath::Vadd(nq, tmp, 1, outarray[2], 1, outarray[2], 1);
194  }
195  break;
196  default:
197  ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
198  break;
199  }
200 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
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.cpp:209
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518
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.cpp:359

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 203 of file NonlinearSWE.cpp.

206 {
207 
208  int ncoeffs = GetNcoeffs();
209  int nq = GetTotPoints();
210 
211  Array<OneD, NekDouble> tmp(nq);
212  Array<OneD, NekDouble> mod(ncoeffs);
213 
214  switch (m_projectionType)
215  {
217  {
218  for (int i = 0; i < m_spacedim; ++i)
219  {
220  Vmath::Vmul(nq, m_bottomSlope[i], 1, physarray[0], 1, tmp, 1);
221  Vmath::Smul(nq, m_g, tmp, 1, tmp, 1);
222  m_fields[0]->IProductWRTBase(tmp, mod);
223  m_fields[0]->MultiplyByElmtInvMass(mod, mod);
224  m_fields[0]->BwdTrans(mod, tmp);
225  Vmath::Vadd(nq, tmp, 1, outarray[i + 1], 1, outarray[i + 1], 1);
226  }
227  }
228  break;
231  {
232  for (int i = 0; i < m_spacedim; ++i)
233  {
234  Vmath::Vmul(nq, m_bottomSlope[i], 1, physarray[0], 1, tmp, 1);
235  Vmath::Smul(nq, m_g, tmp, 1, tmp, 1);
236  Vmath::Vadd(nq, tmp, 1, outarray[i + 1], 1, outarray[i + 1], 1);
237  }
238  }
239  break;
240  default:
241  ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
242  break;
243  }
244 }
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.cpp:248

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 625 of file NonlinearSWE.cpp.

628 {
629  int nq = GetTotPoints();
630 
631  if (physin.get() == physout.get())
632  {
633  // copy indata and work with tmp array
634  Array<OneD, Array<OneD, NekDouble>> tmp(3);
635  for (int i = 0; i < 3; ++i)
636  {
637  // deep copy
638  tmp[i] = Array<OneD, NekDouble>(nq);
639  Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
640  }
641 
642  // \eta = h - d
643  Vmath::Vsub(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
644 
645  // u = hu/h
646  Vmath::Vdiv(nq, tmp[1], 1, tmp[0], 1, physout[1], 1);
647 
648  // v = hv/ v
649  Vmath::Vdiv(nq, tmp[2], 1, tmp[0], 1, physout[2], 1);
650  }
651  else
652  {
653  // \eta = h - d
654  Vmath::Vsub(nq, physin[0], 1, m_depth, 1, physout[0], 1);
655 
656  // u = hu/h
657  Vmath::Vdiv(nq, physin[1], 1, physin[0], 1, physout[1], 1);
658 
659  // v = hv/ v
660  Vmath::Vdiv(nq, physin[2], 1, physin[0], 1, physout[2], 1);
661  }
662 }
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.cpp:284
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255
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.cpp:419

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 54 of file NonlinearSWE.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::NonlinearSWE::DoOdeProjection ( const Array< OneD, const Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray,
const NekDouble  time 
)
protected

Definition at line 357 of file NonlinearSWE.cpp.

360 {
361  int i;
362  int nvariables = inarray.size();
363 
364  switch (m_projectionType)
365  {
367  {
368 
369  // Just copy over array
370  int npoints = GetNpoints();
371 
372  for (i = 0; i < nvariables; ++i)
373  {
374  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
375  }
376  SetBoundaryConditions(outarray, time);
377  break;
378  }
381  {
382 
384  Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs(), 0.0);
385 
386  for (i = 0; i < nvariables; ++i)
387  {
388  m_fields[i]->FwdTrans(inarray[i], coeffs);
389  m_fields[i]->BwdTrans(coeffs, outarray[i]);
390  }
391  break;
392  }
393  default:
394  ASSERTL0(false, "Unknown projection scheme");
395  break;
396  }
397 }
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 246 of file NonlinearSWE.cpp.

249 {
250  int i, j;
251  int ndim = m_spacedim;
252  int nvariables = inarray.size();
253  int nq = GetTotPoints();
254 
255  switch (m_projectionType)
256  {
258  {
259 
260  //-------------------------------------------------------
261  // Compute the DG advection including the numerical flux
262  // by using SolverUtils/Advection
263  // Input and output in physical space
264  Array<OneD, Array<OneD, NekDouble>> advVel;
265 
266  m_advection->Advect(nvariables, m_fields, advVel, inarray, outarray,
267  time);
268  //-------------------------------------------------------
269 
270  //-------------------------------------------------------
271  // negate the outarray since moving terms to the rhs
272  for (i = 0; i < nvariables; ++i)
273  {
274  Vmath::Neg(nq, outarray[i], 1);
275  }
276  //-------------------------------------------------------
277 
278  //-------------------------------------------------
279  // Add "source terms"
280  // Input and output in physical space
281 
282  // Coriolis forcing
283  if (m_coriolis.size() != 0)
284  {
285  AddCoriolis(inarray, outarray);
286  }
287 
288  // Variable Depth
289  if (m_constantDepth != true)
290  {
291  AddVariableDepth(inarray, outarray);
292  }
293  //-------------------------------------------------
294  }
295  break;
298  {
299 
300  //-------------------------------------------------------
301  // Compute the fluxvector in physical space
302  Array<OneD, Array<OneD, Array<OneD, NekDouble>>> fluxvector(
303  nvariables);
304 
305  for (i = 0; i < nvariables; ++i)
306  {
307  fluxvector[i] = Array<OneD, Array<OneD, NekDouble>>(ndim);
308  for (j = 0; j < ndim; ++j)
309  {
310  fluxvector[i][j] = Array<OneD, NekDouble>(nq);
311  }
312  }
313 
314  NonlinearSWE::GetFluxVector(inarray, fluxvector);
315  //-------------------------------------------------------
316 
317  //-------------------------------------------------------
318  // Take the derivative of the flux terms
319  // and negate the outarray since moving terms to the rhs
320  Array<OneD, NekDouble> tmp(nq);
321  Array<OneD, NekDouble> tmp1(nq);
322 
323  for (i = 0; i < nvariables; ++i)
324  {
325  m_fields[i]->PhysDeriv(MultiRegions::DirCartesianMap[0],
326  fluxvector[i][0], tmp);
327  m_fields[i]->PhysDeriv(MultiRegions::DirCartesianMap[1],
328  fluxvector[i][1], tmp1);
329  Vmath::Vadd(nq, tmp, 1, tmp1, 1, outarray[i], 1);
330  Vmath::Neg(nq, outarray[i], 1);
331  }
332 
333  //-------------------------------------------------
334  // Add "source terms"
335  // Input and output in physical space
336 
337  // Coriolis forcing
338  if (m_coriolis.size() != 0)
339  {
340  AddCoriolis(inarray, outarray);
341  }
342 
343  // Variable Depth
344  if (m_constantDepth != true)
345  {
346  AddVariableDepth(inarray, outarray);
347  }
348  //-------------------------------------------------
349  }
350  break;
351  default:
352  ASSERTL0(false, "Unknown projection scheme for the NonlinearSWE");
353  break;
354  }
355 }
void AddCoriolis(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)
void AddVariableDepth(const Array< OneD, const Array< OneD, NekDouble >> &physarray, Array< OneD, Array< OneD, NekDouble >> &outarray)
SolverUtils::AdvectionSharedPtr m_advection
bool m_constantDepth
Indicates if constant depth case.
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:91

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 587 of file NonlinearSWE.cpp.

590 {
591  int i, j;
592  int nq = m_fields[0]->GetTotPoints();
593 
594  NekDouble g = m_g;
595  Array<OneD, Array<OneD, NekDouble>> velocity(m_spacedim);
596 
597  // Flux vector for the mass equation
598  for (i = 0; i < m_spacedim; ++i)
599  {
600  velocity[i] = Array<OneD, NekDouble>(nq);
601  Vmath::Vcopy(nq, physfield[i + 1], 1, flux[0][i], 1);
602  }
603 
604  GetVelocityVector(physfield, velocity);
605 
606  // Put (0.5 g h h) in tmp
607  Array<OneD, NekDouble> tmp(nq);
608  Vmath::Vmul(nq, physfield[0], 1, physfield[0], 1, tmp, 1);
609  Vmath::Smul(nq, 0.5 * g, tmp, 1, tmp, 1);
610 
611  // Flux vector for the momentum equations
612  for (i = 0; i < m_spacedim; ++i)
613  {
614  for (j = 0; j < m_spacedim; ++j)
615  {
616  Vmath::Vmul(nq, velocity[j], 1, physfield[i + 1], 1, flux[i + 1][j],
617  1);
618  }
619 
620  // Add (0.5 g h h) to appropriate field
621  Vmath::Vadd(nq, flux[i + 1][i], 1, tmp, 1, flux[i + 1][i], 1);
622  }
623 }
void GetVelocityVector(const Array< OneD, Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, NekDouble >> &velocity)
Compute the velocity field given the momentum .
double NekDouble

References GetVelocityVector(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::ShallowWaterSystem::m_g, Nektar::SolverUtils::EquationSystem::m_spacedim, Vmath::Smul(), Vmath::Vadd(), Vmath::Vcopy(), 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 745 of file NonlinearSWE.cpp.

748 {
749  const int npts = physfield[0].size();
750 
751  for (int i = 0; i < m_spacedim; ++i)
752  {
753  Vmath::Vdiv(npts, physfield[1 + i], 1, physfield[0], 1, velocity[i], 1);
754  }
755 }

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

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 681 of file NonlinearSWE.cpp.

684 {
685 
686  int nq = GetTotPoints();
687 
688  if (physin.get() == physout.get())
689  {
690  // copy indata and work with tmp array
691  Array<OneD, Array<OneD, NekDouble>> tmp(3);
692  for (int i = 0; i < 3; ++i)
693  {
694  // deep copy
695  tmp[i] = Array<OneD, NekDouble>(nq);
696  Vmath::Vcopy(nq, physin[i], 1, tmp[i], 1);
697  }
698 
699  // h = \eta + d
700  Vmath::Vadd(nq, tmp[0], 1, m_depth, 1, physout[0], 1);
701 
702  // hu = h * u
703  Vmath::Vmul(nq, physout[0], 1, tmp[1], 1, physout[1], 1);
704 
705  // hv = h * v
706  Vmath::Vmul(nq, physout[0], 1, tmp[2], 1, physout[2], 1);
707  }
708  else
709  {
710  // h = \eta + d
711  Vmath::Vadd(nq, physin[0], 1, m_depth, 1, physout[0], 1);
712 
713  // hu = h * u
714  Vmath::Vmul(nq, physout[0], 1, physin[1], 1, physout[1], 1);
715 
716  // hv = h * v
717  Vmath::Vmul(nq, physout[0], 1, physin[2], 1, physout[2], 1);
718  }
719 }

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 400 of file NonlinearSWE.cpp.

402 {
403  std::string varName;
404  int nvariables = m_fields.size();
405  int cnt = 0;
406  int nTracePts = GetTraceTotPoints();
407 
408  // Extract trace for boundaries. Needs to be done on all processors to avoid
409  // deadlock.
410  Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
411  for (int i = 0; i < nvariables; ++i)
412  {
413  Fwd[i] = Array<OneD, NekDouble>(nTracePts);
414  m_fields[i]->ExtractTracePhys(inarray[i], Fwd[i]);
415  }
416 
417  // Loop over Boundary Regions
418  for (int n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
419  {
420 
421  if (m_fields[0]->GetBndConditions()[n]->GetBoundaryConditionType() ==
423  {
424  continue;
425  }
426 
427  // Wall Boundary Condition
428  if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),
429  "Wall"))
430  {
431  WallBoundary2D(n, cnt, Fwd, inarray);
432  }
433 
434  // Time Dependent Boundary Condition (specified in meshfile)
435  if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
436  {
437  for (int i = 0; i < nvariables; ++i)
438  {
439  varName = m_session->GetVariable(i);
440  m_fields[i]->EvaluateBoundaryConditions(time, varName);
441  }
442  }
443  cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
444  }
445 }
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 664 of file NonlinearSWE.cpp.

665 {
666  int nq = GetTotPoints();
667 
668  // u = hu/h
669  Vmath::Vdiv(nq, m_fields[1]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
670  m_fields[1]->UpdatePhys(), 1);
671 
672  // v = hv/ v
673  Vmath::Vdiv(nq, m_fields[2]->GetPhys(), 1, m_fields[0]->GetPhys(), 1,
674  m_fields[2]->UpdatePhys(), 1);
675 
676  // \eta = h - d
677  Vmath::Vsub(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
678  m_fields[0]->UpdatePhys(), 1);
679 }

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.

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 757 of file NonlinearSWE.cpp.

758 {
760  SolverUtils::AddSummaryItem(s, "Variables", "h should be in field[0]");
761  SolverUtils::AddSummaryItem(s, "", "hu should be in field[1]");
762  SolverUtils::AddSummaryItem(s, "", "hv should be in field[2]");
763 }
virtual 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:49

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::ShallowWaterSystem.

Definition at line 57 of file NonlinearSWE.cpp.

58 {
59  ShallowWaterSystem::v_InitObject(DeclareFields);
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  // NB: diffusion not set up for SWE at the moment
90  // but kept here for future use ...
91  m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
92  // m_session->LoadSolverInfo("DiffusionType", diffName, "LDGEddy");
94  advName, advName);
95  // m_diffusion = SolverUtils::GetDiffusionFactory()
96  // .CreateInstance(diffName, diffName);
97 
98  m_advection->SetFluxVector(&NonlinearSWE::GetFluxVector, this);
99  // m_diffusion->SetFluxVectorNS(&ShallowWaterSystem::
100  // GetEddyViscosityFluxVector,
101  // this);
102 
103  // Setting up Riemann solver for advection operator
104  m_session->LoadSolverInfo("UpwindType", riemName, "Average");
107  riemName, m_session);
108 
109  // Setting up upwind solver for diffusion operator
110  // m_riemannSolverLDG = SolverUtils::GetRiemannSolverFactory()
111  // .CreateInstance("UpwindLDG");
112 
113  // Setting up parameters for advection operator Riemann solver
114  m_riemannSolver->SetParam("gravity", &NonlinearSWE::GetGravity,
115  this);
116  m_riemannSolver->SetAuxVec("vecLocs", &NonlinearSWE::GetVecLocs,
117  this);
118  m_riemannSolver->SetVector("N", &NonlinearSWE::GetNormals, this);
119  m_riemannSolver->SetScalar("depth", &NonlinearSWE::GetDepth, this);
120 
121  // Setting up parameters for diffusion operator Riemann solver
122  // m_riemannSolverLDG->AddParam (
123  // "gravity",
124  // &NonlinearSWE::GetGravity, this);
125  // m_riemannSolverLDG->SetAuxVec(
126  // "vecLocs",
127  // &NonlinearSWE::GetVecLocs, this);
128  // m_riemannSolverLDG->AddVector(
129  // "N",
130  // &NonlinearSWE::GetNormals, this);
131 
132  // Concluding initialisation of advection / diffusion operators
133  m_advection->SetRiemannSolver(m_riemannSolver);
134  // m_diffusion->SetRiemannSolver (m_riemannSolverLDG);
135  m_advection->InitObject(m_session, m_fields);
136  // m_diffusion->InitObject (m_session, m_fields);
137  break;
138  }
139  default:
140  {
141  ASSERTL0(false, "Unsupported projection type.");
142  break;
143  }
144  }
145 }
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
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()
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
virtual void v_InitObject(bool DeclareFields=true) override
Init object for UnsteadySystem class.
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:47
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 721 of file NonlinearSWE.cpp.

722 {
723  int nq = GetTotPoints();
724 
725  // h = \eta + d
726  Vmath::Vadd(nq, m_fields[0]->GetPhys(), 1, m_depth, 1,
727  m_fields[0]->UpdatePhys(), 1);
728 
729  // hu = h * u
730  Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[1]->GetPhys(), 1,
731  m_fields[1]->UpdatePhys(), 1);
732 
733  // hv = h * v
734  Vmath::Vmul(nq, m_fields[0]->GetPhys(), 1, m_fields[2]->GetPhys(), 1,
735  m_fields[2]->UpdatePhys(), 1);
736 }

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 451 of file NonlinearSWE.cpp.

454 {
455  int i;
456  int nvariables = physarray.size();
457 
458  // Adjust the physical values of the trace to take
459  // user defined boundaries into account
460  int e, id1, id2, npts;
461 
462  for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
463  ++e)
464  {
465  npts = m_fields[0]
466  ->GetBndCondExpansions()[bcRegion]
467  ->GetExp(e)
468  ->GetTotPoints();
469  id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
470  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
471  m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
472 
473  // For 2D/3D, define: v* = v - 2(v.n)n
474  Array<OneD, NekDouble> tmp(npts, 0.0);
475 
476  // Calculate (v.n)
477  for (i = 0; i < m_spacedim; ++i)
478  {
479  Vmath::Vvtvp(npts, &Fwd[1 + i][id2], 1, &m_traceNormals[i][id2], 1,
480  &tmp[0], 1, &tmp[0], 1);
481  }
482 
483  // Calculate 2.0(v.n)
484  Vmath::Smul(npts, -2.0, &tmp[0], 1, &tmp[0], 1);
485 
486  // Calculate v* = v - 2.0(v.n)n
487  for (i = 0; i < m_spacedim; ++i)
488  {
489  Vmath::Vvtvp(npts, &tmp[0], 1, &m_traceNormals[i][id2], 1,
490  &Fwd[1 + i][id2], 1, &Fwd[1 + i][id2], 1);
491  }
492 
493  // copy boundary adjusted values into the boundary expansion
494  for (i = 0; i < nvariables; ++i)
495  {
496  Vmath::Vcopy(npts, &Fwd[i][id2], 1,
497  &(m_fields[i]
498  ->GetBndCondExpansions()[bcRegion]
499  ->UpdatePhys())[id1],
500  1);
501  }
502  }
503 }
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.cpp:574

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 505 of file NonlinearSWE.cpp.

508 {
509 
510  int i;
511  int nvariables = physarray.size();
512 
513  // Adjust the physical values of the trace to take
514  // user defined boundaries into account
515  int e, id1, id2, npts;
516 
517  for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
518  ++e)
519  {
520  npts = m_fields[0]
521  ->GetBndCondExpansions()[bcRegion]
522  ->GetExp(e)
523  ->GetNumPoints(0);
524  id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
525  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
526  m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt + e));
527 
528  switch (m_expdim)
529  {
530  case 1:
531  {
532  // negate the forward flux
533  Vmath::Neg(npts, &Fwd[1][id2], 1);
534  }
535  break;
536  case 2:
537  {
538  Array<OneD, NekDouble> tmp_n(npts);
539  Array<OneD, NekDouble> tmp_t(npts);
540 
541  Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[0][id2], 1,
542  &tmp_n[0], 1);
543  Vmath::Vvtvp(npts, &Fwd[2][id2], 1, &m_traceNormals[1][id2], 1,
544  &tmp_n[0], 1, &tmp_n[0], 1);
545 
546  Vmath::Vmul(npts, &Fwd[1][id2], 1, &m_traceNormals[1][id2], 1,
547  &tmp_t[0], 1);
548  Vmath::Vvtvm(npts, &Fwd[2][id2], 1, &m_traceNormals[0][id2], 1,
549  &tmp_t[0], 1, &tmp_t[0], 1);
550 
551  // negate the normal flux
552  Vmath::Neg(npts, tmp_n, 1);
553 
554  // rotate back to Cartesian
555  Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[1][id2], 1,
556  &Fwd[1][id2], 1);
557  Vmath::Vvtvm(npts, &tmp_n[0], 1, &m_traceNormals[0][id2], 1,
558  &Fwd[1][id2], 1, &Fwd[1][id2], 1);
559 
560  Vmath::Vmul(npts, &tmp_t[0], 1, &m_traceNormals[0][id2], 1,
561  &Fwd[2][id2], 1);
562  Vmath::Vvtvp(npts, &tmp_n[0], 1, &m_traceNormals[1][id2], 1,
563  &Fwd[2][id2], 1, &Fwd[2][id2], 1);
564  }
565  break;
566  case 3:
567  ASSERTL0(false,
568  "3D not implemented for Shallow Water Equations");
569  break;
570  default:
571  ASSERTL0(false, "Illegal expansion dimension");
572  }
573 
574  // copy boundary adjusted values into the boundary expansion
575  for (i = 0; i < nvariables; ++i)
576  {
577  Vmath::Vcopy(npts, &Fwd[i][id2], 1,
578  &(m_fields[i]
579  ->GetBndCondExpansions()[bcRegion]
580  ->UpdatePhys())[id1],
581  1);
582  }
583  }
584 }
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.cpp:598

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:198
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
Definition: NonlinearSWE.h:54
EquationSystemFactory & GetEquationSystemFactory()

Name of class.

Definition at line 64 of file NonlinearSWE.h.