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

virtual ~LinearSWE ()
 
- 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...
 
- 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 ()
 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, Array< OneD, NekDouble > &output)
 
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 void SetTime (const NekDouble time)
 
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...
 

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

Protected Member Functions

 LinearSWE (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
virtual void v_InitObject ()
 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)
 Print a summary of time stepping parameters. More...
 
virtual void v_PrimitiveToConservative ()
 
virtual void v_ConservativeToPrimitive ()
 
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 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 ()
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise ()
 Sets up initial conditions. More...
 
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D (Array< OneD, Array< OneD, NekDouble >> &solution1D)
 Print the solution at each solution point in a txt file. 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 ()
 
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...
 
- 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)
 

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::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_infosteps
 Number of time steps between outputting status information. More...
 
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::TimeIntegrationWrapperSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
 
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...
 
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...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
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...
 
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_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...
 
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

- Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
- 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 []
 

Detailed Description

Definition at line 49 of file LinearSWE.h.

Constructor & Destructor Documentation

◆ ~LinearSWE()

Nektar::LinearSWE::~LinearSWE ( )
virtual

Definition at line 168 of file LinearSWE.cpp.

169  {
170 
171  }

◆ LinearSWE()

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

Definition at line 51 of file LinearSWE.cpp.

54  : ShallowWaterSystem(pSession, pGraph)
55  {
56  }
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 174 of file LinearSWE.cpp.

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

176  {
177 
178  int ncoeffs = GetNcoeffs();
179  int nq = GetTotPoints();
180 
181  Array<OneD, NekDouble> tmp(nq);
182  Array<OneD, NekDouble> mod(ncoeffs);
183 
184  switch(m_projectionType)
185  {
187  {
188  // add to u equation
189  Vmath::Vmul(nq,m_coriolis,1,physarray[2],1,tmp,1);
190  m_fields[0]->IProductWRTBase(tmp,mod);
191  m_fields[0]->MultiplyByElmtInvMass(mod,mod);
192  m_fields[0]->BwdTrans(mod,tmp);
193  Vmath::Vadd(nq,tmp,1,outarray[1],1,outarray[1],1);
194 
195  // add to v equation
196  Vmath::Vmul(nq,m_coriolis,1,physarray[1],1,tmp,1);
197  Vmath::Neg(nq,tmp,1);
198  m_fields[0]->IProductWRTBase(tmp,mod);
199  m_fields[0]->MultiplyByElmtInvMass(mod,mod);
200  m_fields[0]->BwdTrans(mod,tmp);
201  Vmath::Vadd(nq,tmp,1,outarray[2],1,outarray[2],1);
202  }
203  break;
206  {
207  // add to u equation
208  Vmath::Vmul(nq,m_coriolis,1,physarray[2],1,tmp,1);
209  Vmath::Vadd(nq,tmp,1,outarray[1],1,outarray[1],1);
210 
211  // add to v equation
212  Vmath::Vmul(nq,m_coriolis,1,physarray[1],1,tmp,1);
213  Vmath::Neg(nq,tmp,1);
214  Vmath::Vadd(nq,tmp,1,outarray[2],1,outarray[2],1);
215  }
216  break;
217  default:
218  ASSERTL0(false,"Unknown projection scheme for the NonlinearSWE");
219  break;
220  }
221 
222 
223  }
Array< OneD, NekDouble > m_coriolis
Coriolis force.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT int GetTotPoints()
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:399
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetNcoeffs()
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:302
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:186

◆ ConservativeToPrimitive()

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

Definition at line 583 of file LinearSWE.cpp.

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

585  {
586  int nq = GetTotPoints();
587 
588  if(physin.get() == physout.get())
589  {
590  // copy indata and work with tmp array
591  Array<OneD, Array<OneD, NekDouble> >tmp(3);
592  for (int i = 0; i < 3; ++i)
593  {
594  // deep copy
595  tmp[i] = Array<OneD, NekDouble>(nq);
596  Vmath::Vcopy(nq,physin[i],1,tmp[i],1);
597  }
598 
599  // \eta = h - d
600  Vmath::Vsub(nq,tmp[0],1,m_depth,1,physout[0],1);
601 
602  // u = hu/h
603  Vmath::Vdiv(nq,tmp[1],1,tmp[0],1,physout[1],1);
604 
605  // v = hv/ v
606  Vmath::Vdiv(nq,tmp[2],1,tmp[0],1,physout[2],1);
607  }
608  else
609  {
610  // \eta = h - d
611  Vmath::Vsub(nq,physin[0],1,m_depth,1,physout[0],1);
612 
613  // u = hu/h
614  Vmath::Vdiv(nq,physin[1],1,physin[0],1,physout[1],1);
615 
616  // v = hv/ v
617  Vmath::Vdiv(nq,physin[2],1,physin[0],1,physout[2],1);
618  }
619  }
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:244
SOLVER_UTILS_EXPORT int GetTotPoints()
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:346
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

◆ 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 55 of file LinearSWE.h.

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

58  {
60  ::AllocateSharedPtr(pSession, pGraph);
61  p->InitObject();
62  return p;
63  }
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.

◆ 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 328 of file LinearSWE.cpp.

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

331  {
332  int i;
333  int nvariables = inarray.num_elements();
334 
335 
336  switch(m_projectionType)
337  {
339  {
340 
341  // Just copy over array
342  int npoints = GetNpoints();
343 
344  for(i = 0; i < nvariables; ++i)
345  {
346  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
347  }
348  SetBoundaryConditions(outarray, time);
349  break;
350  }
353  {
354 
356  Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs());
357 
358  for(i = 0; i < nvariables; ++i)
359  {
360  m_fields[i]->FwdTrans(inarray[i],coeffs);
361  m_fields[i]->BwdTrans_IterPerExp(coeffs,outarray[i]);
362  }
363  break;
364  }
365  default:
366  ASSERTL0(false,"Unknown projection scheme");
367  break;
368  }
369  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
void SetBoundaryConditions(Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
Definition: LinearSWE.cpp:373
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetNcoeffs()
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

◆ 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 225 of file LinearSWE.cpp.

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

228  {
229  int i, j;
230  int ndim = m_spacedim;
231  int nvariables = inarray.num_elements();
232  int nq = GetTotPoints();
233 
234 
235  switch(m_projectionType)
236  {
238  {
239 
240  //-------------------------------------------------------
241  // Compute the DG advection including the numerical flux
242  // by using SolverUtils/Advection
243  // Input and output in physical space
244  Array<OneD, Array<OneD, NekDouble> > advVel;
245 
246  m_advection->Advect(nvariables, m_fields, advVel, inarray,
247  outarray, time);
248  //-------------------------------------------------------
249 
250 
251  //-------------------------------------------------------
252  // negate the outarray since moving terms to the rhs
253  for(i = 0; i < nvariables; ++i)
254  {
255  Vmath::Neg(nq,outarray[i],1);
256  }
257  //-------------------------------------------------------
258 
259 
260  //-------------------------------------------------
261  // Add "source terms"
262  // Input and output in physical space
263 
264  // Coriolis forcing
265  if (m_coriolis.num_elements() != 0)
266  {
267  AddCoriolis(inarray,outarray);
268  }
269  //-------------------------------------------------
270 
271  }
272  break;
275  {
276 
277  //-------------------------------------------------------
278  // Compute the fluxvector in physical space
279  Array<OneD, Array<OneD, Array<OneD, NekDouble> > >
280  fluxvector(nvariables);
281 
282  for (i = 0; i < nvariables; ++i)
283  {
284  fluxvector[i] = Array<OneD, Array<OneD, NekDouble> >(ndim);
285  for(j = 0; j < ndim; ++j)
286  {
287  fluxvector[i][j] = Array<OneD, NekDouble>(nq);
288  }
289  }
290 
291  LinearSWE::GetFluxVector(inarray, fluxvector);
292  //-------------------------------------------------------
293 
294 
295  //-------------------------------------------------------
296  // Take the derivative of the flux terms
297  // and negate the outarray since moving terms to the rhs
298  Array<OneD,NekDouble> tmp(nq);
299  Array<OneD, NekDouble>tmp1(nq);
300 
301  for(i = 0; i < nvariables; ++i)
302  {
303  m_fields[i]->PhysDeriv(MultiRegions::DirCartesianMap[0],fluxvector[i][0],tmp);
304  m_fields[i]->PhysDeriv(MultiRegions::DirCartesianMap[1],fluxvector[i][1],tmp1);
305  Vmath::Vadd(nq,tmp,1,tmp1,1,outarray[i],1);
306  Vmath::Neg(nq,outarray[i],1);
307  }
308 
309  //-------------------------------------------------
310  // Add "source terms"
311  // Input and output in physical space
312 
313  // Coriolis forcing
314  if (m_coriolis.num_elements() != 0)
315  {
316  AddCoriolis(inarray,outarray);
317  }
318  //-------------------------------------------------
319  }
320  break;
321  default:
322  ASSERTL0(false,"Unknown projection scheme for the NonlinearSWE");
323  break;
324  }
325  }
Array< OneD, NekDouble > m_coriolis
Coriolis force.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Definition: LinearSWE.cpp:549
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
void AddCoriolis(const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: LinearSWE.cpp:174
SolverUtils::AdvectionSharedPtr m_advection
SOLVER_UTILS_EXPORT int GetTotPoints()
int m_spacedim
Spatial dimension (>= expansion dim).
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:399
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:88
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
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:302

◆ GetDepthBwd()

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

◆ GetDepthFwd()

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

Definition at line 98 of file LinearSWE.h.

References m_dFwd.

Referenced by v_InitObject().

99  {
100  return m_dFwd;
101  }
Array< OneD, NekDouble > m_dFwd
Still water depth traces.
Definition: LinearSWE.h:77

◆ 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 549 of file LinearSWE.cpp.

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

552  {
553  int i, j;
554  int nq = m_fields[0]->GetTotPoints();
555 
556  NekDouble g = m_g;
557 
558  // Flux vector for the mass equation
559  for (i = 0; i < m_spacedim; ++i)
560  {
561  Vmath::Vmul(nq, m_depth, 1, physfield[i+1], 1, flux[0][i], 1);
562  }
563 
564  // Put (g eta) in tmp
565  Array<OneD, NekDouble> tmp(nq);
566  Vmath::Smul(nq, g, physfield[0], 1, tmp, 1);
567 
568  // Flux vector for the momentum equations
569  for (i = 0; i < m_spacedim; ++i)
570  {
571  for (j = 0; j < m_spacedim; ++j)
572  {
573  // must zero fluxes as not initialised to zero in AdvectionWeakDG ...
574  Vmath::Zero(nq, flux[i+1][j], 1);
575  }
576 
577  // Add (g eta) to appropriate field
578  Vmath::Vadd(nq, flux[i+1][i], 1, tmp, 1, flux[i+1][i], 1);
579  }
580 
581  }
Array< OneD, NekDouble > m_depth
Still water depth.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216
int m_spacedim
Spatial dimension (>= expansion dim).
double NekDouble
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376
NekDouble m_g
Acceleration of gravity.
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:302
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:186

◆ 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 700 of file LinearSWE.cpp.

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

Referenced by GetDepthBwd().

703  {
704  const int npts = physfield[0].num_elements();
705 
706  for (int i = 0; i < m_spacedim; ++i)
707  {
708  Vmath::Vcopy(npts, physfield[1+i], 1, velocity[i], 1);
709  }
710  }
int m_spacedim
Spatial dimension (>= expansion dim).
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

◆ NumericalFlux1D()

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

Referenced by GetDepthBwd().

◆ NumericalFlux2D()

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

Referenced by GetDepthBwd().

◆ PrimitiveToConservative()

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

Definition at line 636 of file LinearSWE.cpp.

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

638  {
639 
640  int nq = GetTotPoints();
641 
642  if(physin.get() == physout.get())
643  {
644  // copy indata and work with tmp array
645  Array<OneD, Array<OneD, NekDouble> >tmp(3);
646  for (int i = 0; i < 3; ++i)
647  {
648  // deep copy
649  tmp[i] = Array<OneD, NekDouble>(nq);
650  Vmath::Vcopy(nq,physin[i],1,tmp[i],1);
651  }
652 
653  // h = \eta + d
654  Vmath::Vadd(nq,tmp[0],1,m_depth,1,physout[0],1);
655 
656  // hu = h * u
657  Vmath::Vmul(nq,physout[0],1,tmp[1],1,physout[1],1);
658 
659  // hv = h * v
660  Vmath::Vmul(nq,physout[0],1,tmp[2],1,physout[2],1);
661 
662  }
663  else
664  {
665  // h = \eta + d
666  Vmath::Vadd(nq,physin[0],1,m_depth,1,physout[0],1);
667 
668  // hu = h * u
669  Vmath::Vmul(nq,physout[0],1,physin[1],1,physout[1],1);
670 
671  // hv = h * v
672  Vmath::Vmul(nq,physout[0],1,physin[2],1,physout[2],1);
673 
674  }
675 
676  }
Array< OneD, NekDouble > m_depth
Still water depth.
SOLVER_UTILS_EXPORT int GetTotPoints()
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
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:302
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:186

◆ SetBoundaryConditions()

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

Definition at line 373 of file LinearSWE.cpp.

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

Referenced by DoOdeProjection(), and GetDepthBwd().

376  {
377  std::string varName;
378  int nvariables = m_fields.num_elements();
379  int cnt = 0;
380  int nTracePts = GetTraceTotPoints();
381 
382  // Extract trace for boundaries. Needs to be done on all processors to avoid
383  // deadlock.
384  Array<OneD, Array<OneD, NekDouble> > Fwd(nvariables);
385  for (int i = 0; i < nvariables; ++i)
386  {
387  Fwd[i] = Array<OneD, NekDouble>(nTracePts);
388  m_fields[i]->ExtractTracePhys(inarray[i], Fwd[i]);
389  }
390 
391  // loop over Boundary Regions
392  for(int n = 0; n < m_fields[0]->GetBndConditions().num_elements(); ++n)
393  {
394  if (m_fields[0]->GetBndConditions()[n]->GetBoundaryConditionType()
396  {
397  continue;
398  }
399 
400  // Wall Boundary Condition
401  if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),"Wall"))
402  {
403  WallBoundary2D(n, cnt, Fwd, inarray);
404  }
405 
406  // Time Dependent Boundary Condition (specified in meshfile)
407  if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
408  {
409  for (int i = 0; i < nvariables; ++i)
410  {
411  varName = m_session->GetVariable(i);
412  m_fields[i]->EvaluateBoundaryConditions(time, varName);
413  }
414  }
415  cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
416  }
417  }
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
void WallBoundary2D(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray)
Definition: LinearSWE.cpp:484

◆ v_ConservativeToPrimitive()

void Nektar::LinearSWE::v_ConservativeToPrimitive ( )
protectedvirtual

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 622 of file LinearSWE.cpp.

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

623  {
624  int nq = GetTotPoints();
625 
626  // u = hu/h
627  Vmath::Vdiv(nq,m_fields[1]->GetPhys(),1,m_fields[0]->GetPhys(),1,m_fields[1]->UpdatePhys(),1);
628 
629  // v = hv/ v
630  Vmath::Vdiv(nq,m_fields[2]->GetPhys(),1,m_fields[0]->GetPhys(),1,m_fields[2]->UpdatePhys(),1);
631 
632  // \eta = h - d
633  Vmath::Vsub(nq,m_fields[0]->GetPhys(),1,m_depth,1,m_fields[0]->UpdatePhys(),1);
634  }
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:244
SOLVER_UTILS_EXPORT int GetTotPoints()
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:346
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.

◆ v_GenerateSummary()

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

Print a summary of time stepping parameters.

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 713 of file LinearSWE.cpp.

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

714  {
716  if (m_session->DefinesSolverInfo("UpwindType"))
717  {
718  std::string UpwindType;
719  UpwindType = m_session->GetSolverInfo("UpwindType");
720  if (UpwindType == "LinearAverage")
721  {
722  SolverUtils::AddSummaryItem(s, "Riemann Solver", "Linear Average");
723  }
724  if (UpwindType == "LinearHLL")
725  {
726  SolverUtils::AddSummaryItem(s, "Riemann Solver", "Linear HLL");
727  }
728  }
729  SolverUtils::AddSummaryItem(s, "Variables", "eta should be in field[0]");
730  SolverUtils::AddSummaryItem(s, "", "u should be in field[1]");
731  SolverUtils::AddSummaryItem(s, "", "v should be in field[2]");
732  }
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
virtual void v_GenerateSummary(SolverUtils::SummaryList &s)
Print a summary of time stepping parameters.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.

◆ v_InitObject()

void Nektar::LinearSWE::v_InitObject ( )
protectedvirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 58 of file LinearSWE.cpp.

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

59  {
61 
63  {
66  }
67  else
68  {
69  ASSERTL0(false, "Implicit SWE not set up.");
70  }
71 
72  // Type of advection class to be used
73  switch(m_projectionType)
74  {
75  // Continuous field
77  {
78  // Do nothing
79  break;
80  }
81  // Discontinuous field
83  {
84  string advName;
85  string diffName;
86  string riemName;
87 
88  //---------------------------------------------------------------
89  // Setting up advection and diffusion operators
90  // NB: diffusion not set up for SWE at the moment
91  // but kept here for future use ...
92  m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
93  // m_session->LoadSolverInfo("DiffusionType", diffName, "LDGEddy");
95  .CreateInstance(advName, advName);
96  // m_diffusion = SolverUtils::GetDiffusionFactory()
97  // .CreateInstance(diffName, diffName);
98 
99  m_advection->SetFluxVector(&LinearSWE::GetFluxVector, this);
100  // m_diffusion->SetFluxVectorNS(&ShallowWaterSystem::
101  // GetEddyViscosityFluxVector, this);
102 
103  // Setting up Riemann solver for advection operator
104  m_session->LoadSolverInfo("UpwindType", riemName, "NoSolver");
105  if ((riemName == "LinearHLL") && (m_constantDepth != true))
106  {
107  ASSERTL0(false,"LinearHLL only valid for constant depth");
108  }
110  .CreateInstance(riemName, m_session);
111 
112  // Setting up upwind solver for diffusion operator
113  // m_riemannSolverLDG = SolverUtils::GetRiemannSolverFactory()
114  // .CreateInstance("UpwindLDG");
115 
116  // Setting up parameters for advection operator Riemann solver
117  m_riemannSolver->SetParam (
118  "gravity",
119  &LinearSWE::GetGravity, this);
120  m_riemannSolver->SetAuxVec(
121  "vecLocs",
122  &LinearSWE::GetVecLocs, this);
123  m_riemannSolver->SetVector(
124  "N",
125  &LinearSWE::GetNormals, this);
126 
127  // The numerical flux for linear SWE requires depth information
128  int nTracePointsTot = m_fields[0]->GetTrace()->GetTotPoints();
129  m_dFwd = Array<OneD, NekDouble>(nTracePointsTot);
130  m_dBwd = Array<OneD, NekDouble>(nTracePointsTot);
131  m_fields[0]->GetFwdBwdTracePhys(m_depth, m_dFwd, m_dBwd);
132  CopyBoundaryTrace(m_dFwd,m_dBwd);
133  m_riemannSolver->SetScalar(
134  "depthFwd",
135  &LinearSWE::GetDepthFwd, this);
136  m_riemannSolver->SetScalar(
137  "depthBwd",
138  &LinearSWE::GetDepthBwd, this);
139 
140  // Setting up parameters for diffusion operator Riemann solver
141  // m_riemannSolverLDG->AddParam (
142  // "gravity",
143  // &NonlinearSWE::GetGravity, this);
144  // m_riemannSolverLDG->SetAuxVec(
145  // "vecLocs",
146  // &NonlinearSWE::GetVecLocs, this);
147  // m_riemannSolverLDG->AddVector(
148  // "N",
149  // &NonlinearSWE::GetNormals, this);
150 
151  // Concluding initialisation of advection / diffusion operators
152  m_advection->SetRiemannSolver (m_riemannSolver);
153  //m_diffusion->SetRiemannSolver (m_riemannSolverLDG);
154  m_advection->InitObject (m_session, m_fields);
155  //m_diffusion->InitObject (m_session, m_fields);
156  break;
157  }
158  default:
159  {
160  ASSERTL0(false, "Unsupported projection type.");
161  break;
162  }
163  }
164 
165 
166  }
const Array< OneD, NekDouble > & GetDepthFwd()
Definition: LinearSWE.h:98
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
Array< OneD, NekDouble > m_depth
Still water depth.
Array< OneD, NekDouble > m_dFwd
Still water depth traces.
Definition: LinearSWE.h:77
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Definition: LinearSWE.cpp:549
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs()
SolverUtils::AdvectionSharedPtr m_advection
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
virtual void v_InitObject()
Init object for UnsteadySystem class.
RiemannSolverFactory & GetRiemannSolverFactory()
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
void CopyBoundaryTrace(const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
Array< OneD, NekDouble > m_dBwd
Definition: LinearSWE.h:78
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Definition: LinearSWE.cpp:328
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
bool m_constantDepth
Indicates if constant depth case.
const Array< OneD, NekDouble > & GetDepthBwd()
Definition: LinearSWE.h:102
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Definition: LinearSWE.cpp:225

◆ v_PrimitiveToConservative()

void Nektar::LinearSWE::v_PrimitiveToConservative ( )
protectedvirtual

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 678 of file LinearSWE.cpp.

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

679  {
680  int nq = GetTotPoints();
681 
682  // h = \eta + d
683  Vmath::Vadd(nq,m_fields[0]->GetPhys(),1,m_depth,1,m_fields[0]->UpdatePhys(),1);
684 
685  // hu = h * u
686  Vmath::Vmul(nq,m_fields[0]->GetPhys(),1,m_fields[1]->GetPhys(),1,m_fields[1]->UpdatePhys(),1);
687 
688  // hv = h * v
689  Vmath::Vmul(nq,m_fields[0]->GetPhys(),1,m_fields[2]->GetPhys(),1,m_fields[2]->UpdatePhys(),1);
690  }
Array< OneD, NekDouble > m_depth
Still water depth.
SOLVER_UTILS_EXPORT int GetTotPoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
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:302
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:186

◆ 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 423 of file LinearSWE.cpp.

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

Referenced by GetDepthBwd().

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

◆ 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 484 of file LinearSWE.cpp.

References ASSERTL0, 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 GetDepthBwd(), and SetBoundaryConditions().

485  {
486 
487  int i;
488  int nvariables = physarray.num_elements();
489 
490  // Adjust the physical values of the trace to take
491  // user defined boundaries into account
492  int e, id1, id2, npts;
493 
494  for(e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize(); ++e)
495  {
496  npts = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExp(e)->GetNumPoints(0);
497  id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e) ;
498  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(m_fields[0]->GetTraceMap()->GetBndCondCoeffsToGlobalCoeffsMap(cnt+e));
499 
500  switch(m_expdim)
501  {
502  case 1:
503  {
504  // negate the forward flux
505  Vmath::Neg(npts,&Fwd[1][id2],1);
506  }
507  break;
508  case 2:
509  {
510  Array<OneD, NekDouble> tmp_n(npts);
511  Array<OneD, NekDouble> tmp_t(npts);
512 
513  Vmath::Vmul(npts,&Fwd[1][id2],1,&m_traceNormals[0][id2],1,&tmp_n[0],1);
514  Vmath::Vvtvp(npts,&Fwd[2][id2],1,&m_traceNormals[1][id2],1,&tmp_n[0],1,&tmp_n[0],1);
515 
516  Vmath::Vmul(npts,&Fwd[1][id2],1,&m_traceNormals[1][id2],1,&tmp_t[0],1);
517  Vmath::Vvtvm(npts,&Fwd[2][id2],1,&m_traceNormals[0][id2],1,&tmp_t[0],1,&tmp_t[0],1);
518 
519  // negate the normal flux
520  Vmath::Neg(npts,tmp_n,1);
521 
522  // rotate back to Cartesian
523  Vmath::Vmul(npts,&tmp_t[0],1,&m_traceNormals[1][id2],1,&Fwd[1][id2],1);
524  Vmath::Vvtvm(npts,&tmp_n[0],1,&m_traceNormals[0][id2],1,&Fwd[1][id2],1,&Fwd[1][id2],1);
525 
526  Vmath::Vmul(npts,&tmp_t[0],1,&m_traceNormals[0][id2],1,&Fwd[2][id2],1);
527  Vmath::Vvtvp(npts,&tmp_n[0],1,&m_traceNormals[1][id2],1,&Fwd[2][id2],1,&Fwd[2][id2],1);
528  }
529  break;
530  case 3:
531  ASSERTL0(false,"3D not implemented for Shallow Water Equations");
532  break;
533  default:
534  ASSERTL0(false,"Illegal expansion dimension");
535  }
536 
537 
538 
539  // copy boundary adjusted values into the boundary expansion
540  for (i = 0; i < nvariables; ++i)
541  {
542  Vmath::Vcopy(npts,&Fwd[i][id2], 1,&(m_fields[i]->GetBndCondExpansions()[bcRegion]->UpdatePhys())[id1],1);
543  }
544  }
545  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
int m_expdim
Expansion dimension.
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:445
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:399
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
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 plus vector): z = w*x - y
Definition: Vmath.cpp:468
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
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:186

Friends And Related Function Documentation

◆ MemoryManager< LinearSWE >

friend class MemoryManager< LinearSWE >
friend

Definition at line 52 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.")

Name of class.

Definition at line 65 of file LinearSWE.h.

◆ m_dBwd

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

Definition at line 78 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 77 of file LinearSWE.h.

Referenced by GetDepthFwd(), and v_InitObject().