Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Member Functions | Friends | List of all members
Nektar::ShallowWaterSystem Class Reference

Base class for unsteady solvers. More...

#include <ShallowWaterSystem.h>

Inheritance diagram for Nektar::ShallowWaterSystem:
Inheritance graph
[legend]
Collaboration diagram for Nektar::ShallowWaterSystem:
Collaboration graph
[legend]

Public Member Functions

virtual ~ShallowWaterSystem ()
 Destructor.
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem ()
 Destructor.
SOLVER_UTILS_EXPORT NekDouble GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Calculate the larger time-step mantaining the problem stable.
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor.
SOLVER_UTILS_EXPORT void SetUpTraceNormals (void)
SOLVER_UTILS_EXPORT void InitObject ()
 Initialises the members of this object.
SOLVER_UTILS_EXPORT void DoInitialise ()
 Perform any initialisation necessary before solving the problem.
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem.
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space.
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space.
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve.
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation.
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name.
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name.
SOLVER_UTILS_EXPORT
LibUtilities::SessionReaderSharedPtr 
GetSession ()
 Get Session name.
SOLVER_UTILS_EXPORT
MultiRegions::ExpListSharedPtr 
GetPressure ()
 Get pressure field if available.
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics.
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda.
SOLVER_UTILS_EXPORT void EvaluateFunction (Array< OneD, Array< OneD, NekDouble > > &pArray, std::string pFunctionName, const NekDouble pTime=0.0, const int domain=0)
 Evaluates a function as specified in the session file.
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, Array< OneD, NekDouble > > &pFields, const std::string &pName, const int domain=0)
 Populate given fields with the function from session.
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const std::string &pName, const int domain=0)
 Populate given fields with the function from session.
SOLVER_UTILS_EXPORT void EvaluateFunction (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
SOLVER_UTILS_EXPORT std::string DescribeFunction (std::string pFieldName, const std::string &pFunctionName, const int domain)
 Provide a description of a function for a given field name.
SOLVER_UTILS_EXPORT void InitialiseBaseFlow (Array< OneD, Array< OneD, NekDouble > > &base)
 Perform initialisation of the base flow.
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields.
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution.
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.
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields.
SOLVER_UTILS_EXPORT Array
< OneD, NekDouble
ErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf].
SOLVER_UTILS_EXPORT void WeakAdvectionGreensDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\nabla \phi \cdot F) $.
SOLVER_UTILS_EXPORT void WeakAdvectionDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\phi, \nabla \cdot F) $.
SOLVER_UTILS_EXPORT void WeakAdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, bool UseContCoeffs=false)
 Compute the inner product $ (\phi, V\cdot \nabla u) $.
f SOLVER_UTILS_EXPORT void AdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray)
 Compute the non-conservative advection.
SOLVER_UTILS_EXPORT void WeakDGAdvection (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false, int nvariables=0)
 Calculate the weak discontinuous Galerkin advection.
SOLVER_UTILS_EXPORT void WeakDGDiffusion (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false)
 Calculate weak DG Diffusion in the LDG form.
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields.
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.
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename.
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.
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file.
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.
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.
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.
SOLVER_UTILS_EXPORT void ScanForHistoryPoints ()
 Builds map of which element holds each history point.
SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out)
 Probe each history point and write to file.
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary.
SOLVER_UTILS_EXPORT Array
< OneD,
MultiRegions::ExpListSharedPtr > & 
UpdateFields ()
SOLVER_UTILS_EXPORT
LibUtilities::FieldMetaDataMap
UpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated.
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time.
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 GetNumElmVelocity ()
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 SetStepsToOne ()
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
SOLVER_UTILS_EXPORT void FwdTransFields ()
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &fluxX, Array< OneD, Array< OneD, NekDouble > > &fluxY)
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, const int j, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
SOLVER_UTILS_EXPORT void NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
SOLVER_UTILS_EXPORT void NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
SOLVER_UTILS_EXPORT int NoCaseStringCompare (const string &s1, const string &s2)
 Perform a case-insensitive string comparison.

Static Public Member Functions

static
SolverUtils::EquationSystemSharedPtr 
create (const LibUtilities::SessionReaderSharedPtr &pSession)
 Creates an instance of this class.

Static Public Attributes

static std::string className
 Name of class.

Protected Member Functions

 ShallowWaterSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises UnsteadySystem class members.
virtual void v_InitObject ()
 Init object for UnsteadySystem class.
virtual void v_GenerateSummary (SolverUtils::SummaryList &s)
 Print a summary of time stepping parameters.
void PrimitiveToConservative ()
virtual void v_PrimitiveToConservative ()
void ConservativeToPrimitive ()
virtual void v_ConservativeToPrimitive ()
NekDouble GetGravity ()
const Array< OneD, const Array
< OneD, NekDouble > > & 
GetVecLocs ()
const Array< OneD, const Array
< OneD, NekDouble > > & 
GetNormals ()
const Array< OneD, NekDouble > & GetDepth ()
bool IsConstantDepth ()
void CopyBoundaryTrace (const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises UnsteadySystem class members.
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control.
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Solves an unsteady problem.
virtual SOLVER_UTILS_EXPORT void v_DoInitialise ()
 Sets up initial conditions.
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D (Array< OneD, Array< OneD, NekDouble > > &solution1D)
 Print the solution at each solution point in a txt file.
virtual SOLVER_UTILS_EXPORT void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
virtual SOLVER_UTILS_EXPORT void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
virtual SOLVER_UTILS_EXPORT void v_NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
virtual SOLVER_UTILS_EXPORT void v_NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
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.
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate (int step)
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time)
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises EquationSystem class members.
int nocase_cmp (const string &s1, const string &s2)
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time.
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.
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.
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space.
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space.
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)
SOLVER_UTILS_EXPORT void SetUpBaseFields (SpatialDomains::MeshGraphSharedPtr &mesh)
SOLVER_UTILS_EXPORT void ImportFldBase (std::string pInfile, SpatialDomains::MeshGraphSharedPtr pGraph)
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

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.
bool m_constantDepth
 Indicates if constant depth case.
NekDouble m_g
 Acceleration of gravity.
Array< OneD, NekDoublem_depth
 Still water depth.
Array< OneD, Array< OneD,
NekDouble > > 
m_bottomSlope
Array< OneD, NekDoublem_coriolis
 Coriolis force.
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.
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
 Wrapper to the time integration scheme.
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use.
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
NekDouble m_epsilon
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used.
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used.
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used.
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state.
std::vector< int > m_intVariables
std::vector< FilterSharedPtrm_filters
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator.
LibUtilities::SessionReaderSharedPtr m_session
 The session reader.
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output.
Array< OneD,
MultiRegions::ExpListSharedPtr
m_fields
 Array holding all dependent variables.
Array< OneD,
MultiRegions::ExpListSharedPtr
m_base
 Base fields.
Array< OneD,
MultiRegions::ExpListSharedPtr
m_derivedfields
 Array holding all dependent variables.
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object.
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh.
std::string m_filename
 Filename.
std::string m_sessionName
 Name of the session.
NekDouble m_time
 Current time of simulation.
NekDouble m_fintime
 Finish time of the simulation.
NekDouble m_timestep
 Time step size.
NekDouble m_lambda
 Lambda constant in real system if one required.
NekDouble m_checktime
 Time between checkpoints.
int m_steps
 Number of steps to take.
int m_checksteps
 Number of steps between checkpoints.
int m_spacedim
 Spatial dimension (>= expansion dim).
int m_expdim
 Expansion dimension.
bool m_SingleMode
 Flag to determine if single homogeneous mode is used.
bool m_HalfMode
 Flag to determine if half homogeneous mode is used.
bool m_MultipleModes
 Flag to determine if use multiple homogenenous modes are used.
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform.
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations.
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous.
Array< OneD, Array< OneD,
NekDouble > > 
m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction.
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_gradtan
 1 x nvariable x nq
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_tanbasis
 2 x m_spacedim x nq
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity.
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields.
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error.
enum HomogeneousType m_HomogeneousType
NekDouble m_LhomX
 physical length in X direction (if homogeneous)
NekDouble m_LhomY
 physical length in Y direction (if homogeneous)
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous)
int m_npointsX
 number of points in X direction (if homogeneous)
int m_npointsY
 number of points in Y direction (if homogeneous)
int m_npointsZ
 number of points in Z direction (if homogeneous)
int m_HomoDirec
 number of homogenous directions
int m_NumMode
 Mode to use in case of single mode analysis.

Private Member Functions

void EvaluateWaterDepth (void)
void EvaluateCoriolis (void)

Friends

class MemoryManager< ShallowWaterSystem >

Additional Inherited Members

- Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1).
- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D, eHomogeneous2D, eHomogeneous3D, eNotHomogeneous }
 Parameter for homogeneous expansions. More...

Detailed Description

Base class for unsteady solvers.

Provides the underlying timestepping framework for shallow water flow solvers including the general timestepping routines. This class is not intended to be directly instantiated, but rather is a base class on which to define shallow water solvers, e.g. SWE, Boussinesq, linear and nonlinear versions.

For details on implementing unsteady solvers see sectionADRSolverModuleImplementation

Definition at line 48 of file ShallowWaterSystem.h.

Constructor & Destructor Documentation

Nektar::ShallowWaterSystem::~ShallowWaterSystem ( )
virtual

Destructor.

Definition at line 135 of file ShallowWaterSystem.cpp.

{
}
Nektar::ShallowWaterSystem::ShallowWaterSystem ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Initialises UnsteadySystem class members.

Definition at line 67 of file ShallowWaterSystem.cpp.

: UnsteadySystem(pSession)
{
}

Member Function Documentation

void Nektar::ShallowWaterSystem::ConservativeToPrimitive ( )
inlineprotected

Definition at line 100 of file ShallowWaterSystem.h.

References v_ConservativeToPrimitive().

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

Definition at line 176 of file ShallowWaterSystem.cpp.

References Nektar::SolverUtils::EquationSystem::GetPhys_Offset(), Nektar::SolverUtils::EquationSystem::m_fields, npts, and Vmath::Vcopy().

Referenced by Nektar::LinearSWE::v_InitObject().

{
int cnt = 0;
// loop over Boundary Regions
for(int bcRegion = 0; bcRegion < m_fields[0]->GetBndConditions().num_elements(); ++bcRegion)
{
// Copy the forward trace of the field to the backward trace
int e, id1, id2, npts;
for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]
->GetExpSize(); ++e)
{
npts = m_fields[0]->GetBndCondExpansions()[bcRegion]->
GetExp(e)->GetTotPoints();
id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->
id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
m_fields[0]->GetTraceMap()->
GetBndCondCoeffsToGlobalCoeffsMap(cnt+e));
Vmath::Vcopy(npts, &Fwd[id2], 1, &Bwd[id2], 1);
}
cnt +=m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
}
}
static SolverUtils::EquationSystemSharedPtr Nektar::ShallowWaterSystem::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Creates an instance of this class.

Reimplemented in Nektar::LinearSWE, and Nektar::NonlinearSWE.

Definition at line 54 of file ShallowWaterSystem.h.

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

Definition at line 171 of file ShallowWaterSystem.cpp.

References Nektar::SolverUtils::EquationSystem::EvaluateFunction(), and m_coriolis.

Referenced by v_InitObject().

{
EvaluateFunction("f",m_coriolis,"Coriolis");
}
void Nektar::ShallowWaterSystem::EvaluateWaterDepth ( void  )
private

Definition at line 165 of file ShallowWaterSystem.cpp.

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

Referenced by v_InitObject().

{
EvaluateFunction("d",m_depth,"WaterDepth");
}
const Array<OneD, NekDouble>& Nektar::ShallowWaterSystem::GetDepth ( )
inlineprotected

Definition at line 122 of file ShallowWaterSystem.h.

References m_depth.

Referenced by Nektar::NonlinearSWE::v_InitObject().

{
return m_depth;
}
NekDouble Nektar::ShallowWaterSystem::GetGravity ( )
inlineprotected

Definition at line 107 of file ShallowWaterSystem.h.

References m_g.

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

{
return m_g;
}
const Array<OneD, const Array<OneD, NekDouble> >& Nektar::ShallowWaterSystem::GetNormals ( )
inlineprotected
const Array<OneD, const Array<OneD, NekDouble> >& Nektar::ShallowWaterSystem::GetVecLocs ( )
inlineprotected

Definition at line 112 of file ShallowWaterSystem.h.

References m_vecLocs.

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

{
return m_vecLocs;
}
bool Nektar::ShallowWaterSystem::IsConstantDepth ( )
inlineprotected

Definition at line 127 of file ShallowWaterSystem.h.

References m_constantDepth.

{
}
void Nektar::ShallowWaterSystem::PrimitiveToConservative ( )
inlineprotected

Definition at line 94 of file ShallowWaterSystem.h.

References v_PrimitiveToConservative().

void Nektar::ShallowWaterSystem::v_ConservativeToPrimitive ( )
protectedvirtual

Reimplemented in Nektar::LinearSWE, and Nektar::NonlinearSWE.

Definition at line 160 of file ShallowWaterSystem.cpp.

References ASSERTL0.

Referenced by ConservativeToPrimitive().

{
ASSERTL0(false, "This function is not implemented for this equation.");
}
void Nektar::ShallowWaterSystem::v_GenerateSummary ( SolverUtils::SummaryList s)
protectedvirtual

Print a summary of time stepping parameters.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Reimplemented in Nektar::LinearSWE, and Nektar::NonlinearSWE.

Definition at line 140 of file ShallowWaterSystem.cpp.

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

{
if (m_constantDepth == true)
{
SolverUtils::AddSummaryItem(s, "Depth", "constant");
}
else
{
SolverUtils::AddSummaryItem(s, "Depth", "variable");
}
}
void Nektar::ShallowWaterSystem::v_InitObject ( )
protectedvirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Reimplemented in Nektar::LinearSWE, and Nektar::NonlinearSWE.

Definition at line 73 of file ShallowWaterSystem.cpp.

References ASSERTL0, Nektar::MultiRegions::DirCartesianMap, Nektar::MultiRegions::eDiscontinuous, EvaluateCoriolis(), EvaluateWaterDepth(), Nektar::SolverUtils::EquationSystem::GetTotPoints(), m_bottomSlope, m_constantDepth, m_depth, Nektar::SolverUtils::EquationSystem::m_fields, m_g, Nektar::SolverUtils::UnsteadySystem::m_infosteps, m_primitive, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_spacedim, m_vecLocs, and Vmath::Neg().

{
// if discontinuous Galerkin determine numerical flux to use
{
ASSERTL0(m_session->DefinesSolverInfo("UPWINDTYPE"),
"No UPWINDTYPE defined in session.");
}
// Set up locations of velocity vector.
m_vecLocs = Array<OneD, Array<OneD, NekDouble> >(1);
m_vecLocs[0] = Array<OneD, NekDouble>(m_spacedim);
for (int i = 0; i < m_spacedim; ++i)
{
m_vecLocs[0][i] = 1 + i;
}
// Load generic input parameters
m_session->LoadParameter("IO_InfoSteps", m_infosteps, 0);
// Load acceleration of gravity
m_session->LoadParameter("Gravity", m_g, 9.81);
// input/output in primitive variables
m_primitive = true;
NekDouble depth = m_depth[0];
for (int i = 0; i < GetTotPoints(); ++i)
{
if (m_depth[i] != depth)
{
m_constantDepth = false;
break;
}
}
// Compute the bottom slopes
int nq = GetTotPoints();
if (m_constantDepth != true)
{
m_bottomSlope = Array<OneD, Array<OneD, NekDouble> >(m_spacedim);
for (int i = 0; i < m_spacedim; ++i)
{
m_bottomSlope[i] = Array<OneD, NekDouble>(nq);
}
}
}
void Nektar::ShallowWaterSystem::v_PrimitiveToConservative ( )
protectedvirtual

Reimplemented in Nektar::LinearSWE, and Nektar::NonlinearSWE.

Definition at line 155 of file ShallowWaterSystem.cpp.

References ASSERTL0.

Referenced by PrimitiveToConservative().

{
ASSERTL0(false, "This function is not implemented for this equation.");
}

Friends And Related Function Documentation

friend class MemoryManager< ShallowWaterSystem >
friend

Definition at line 51 of file ShallowWaterSystem.h.

Member Data Documentation

string Nektar::ShallowWaterSystem::className
static
Initial value:
"ShallowWaterSystem",
"Auxiliary functions for the shallow water system.")

Name of class.

Processes SolverInfo parameters from the session file and sets up timestepping-specific code.

Parameters
pSessionSession object to read parameters from.

Definition at line 61 of file ShallowWaterSystem.h.

SolverUtils::AdvectionSharedPtr Nektar::ShallowWaterSystem::m_advection
protected
Array<OneD, Array<OneD, NekDouble> > Nektar::ShallowWaterSystem::m_bottomSlope
protected

Definition at line 81 of file ShallowWaterSystem.h.

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

bool Nektar::ShallowWaterSystem::m_constantDepth
protected

Indicates if constant depth case.

Definition at line 75 of file ShallowWaterSystem.h.

Referenced by Nektar::NonlinearSWE::DoOdeRhs(), IsConstantDepth(), v_GenerateSummary(), Nektar::LinearSWE::v_InitObject(), and v_InitObject().

Array<OneD, NekDouble> Nektar::ShallowWaterSystem::m_coriolis
protected
Array<OneD, NekDouble> Nektar::ShallowWaterSystem::m_depth
protected
SolverUtils::DiffusionSharedPtr Nektar::ShallowWaterSystem::m_diffusion
protected

Definition at line 70 of file ShallowWaterSystem.h.

NekDouble Nektar::ShallowWaterSystem::m_g
protected
bool Nektar::ShallowWaterSystem::m_primitive
protected

Indicates if variables are primitive or conservative.

Definition at line 73 of file ShallowWaterSystem.h.

Referenced by v_InitObject().

SolverUtils::RiemannSolverSharedPtr Nektar::ShallowWaterSystem::m_riemannSolver
protected
SolverUtils::RiemannSolverSharedPtr Nektar::ShallowWaterSystem::m_riemannSolverLDG
protected

Definition at line 68 of file ShallowWaterSystem.h.

Array<OneD, Array<OneD, NekDouble> > Nektar::ShallowWaterSystem::m_vecLocs
protected

Definition at line 85 of file ShallowWaterSystem.h.

Referenced by GetVecLocs(), and v_InitObject().