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:
Inheritance graph
[legend]
Collaboration diagram for Nektar::NonlinearSWE:
Collaboration graph
[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...
 
- 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 >
boost::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 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 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. More...
 
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. More...
 
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. More...
 
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. More...
 
SOLVER_UTILS_EXPORT void InitialiseBaseFlow (Array< OneD, Array< OneD, NekDouble > > &base)
 Perform initialisation of the base flow. 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 WeakAdvectionGreensDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\nabla \phi \cdot F) $. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\phi, \nabla \cdot F) $. More...
 
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) $. More...
 
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. More...
 
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. More...
 
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. 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 ScanForHistoryPoints ()
 Builds map of which element holds each history point. More...
 
SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out)
 Probe each history point and write to 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 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. 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)
 Creates an instance of this class. More...
 
- Static Public Member Functions inherited from Nektar::ShallowWaterSystem
static SolverUtils::EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession)
 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

 NonlinearSWE (const LibUtilities::SessionReaderSharedPtr &pSession)
 
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 ()
 
- Protected Member Functions inherited from Nektar::ShallowWaterSystem
 ShallowWaterSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 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)
 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 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. 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_SteadyStateCheck (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. More...
 
int nocase_cmp (const string &s1, const string &s2)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. 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)
 
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)
 

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 > > &physarray)
 
void WallBoundary (int bcRegion, int cnt, 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...
 
- 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_infosteps
 Number of time steps between outputting status information. More...
 
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...
 
std::vector< int > m_intVariables
 
std::vector< FilterSharedPtrm_filters
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
map< std::string, Array< OneD, Array< OneD, float > > > m_interpWeights
 Map of the interpolation weights for a specific filename. More...
 
map< std::string, Array< OneD, Array< OneD, unsigned int > > > m_interpInds
 Map of the interpolation indices for a specific filename. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_base
 Base fields. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_derivedfields
 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...
 
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_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, Array< OneD, Array< OneD, NekDouble > > > m_gradtan
 1 x nvariable x nq More...
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_tanbasis
 2 x m_spacedim x nq 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...
 
int m_NumMode
 Mode to use in case of single mode analysis. More...
 

Detailed Description

Definition at line 50 of file NonlinearSWE.h.

Constructor & Destructor Documentation

Nektar::NonlinearSWE::~NonlinearSWE ( )
virtual

Definition at line 152 of file NonlinearSWE.cpp.

153  {
154 
155  }
Nektar::NonlinearSWE::NonlinearSWE ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Definition at line 50 of file NonlinearSWE.cpp.

52  : ShallowWaterSystem(pSession)
53  {
54  }
ShallowWaterSystem(const LibUtilities::SessionReaderSharedPtr &pSession)
Initialises UnsteadySystem class members.

Member Function Documentation

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

Definition at line 158 of file NonlinearSWE.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().

160  {
161 
162  int ncoeffs = GetNcoeffs();
163  int nq = GetTotPoints();
164 
165  Array<OneD, NekDouble> tmp(nq);
166  Array<OneD, NekDouble> mod(ncoeffs);
167 
168  switch(m_projectionType)
169  {
171  {
172  // add to hu equation
173  Vmath::Vmul(nq,m_coriolis,1,physarray[2],1,tmp,1);
174  m_fields[0]->IProductWRTBase(tmp,mod);
175  m_fields[0]->MultiplyByElmtInvMass(mod,mod);
176  m_fields[0]->BwdTrans(mod,tmp);
177  Vmath::Vadd(nq,tmp,1,outarray[1],1,outarray[1],1);
178 
179  // add to hv equation
180  Vmath::Vmul(nq,m_coriolis,1,physarray[1],1,tmp,1);
181  Vmath::Neg(nq,tmp,1);
182  m_fields[0]->IProductWRTBase(tmp,mod);
183  m_fields[0]->MultiplyByElmtInvMass(mod,mod);
184  m_fields[0]->BwdTrans(mod,tmp);
185  Vmath::Vadd(nq,tmp,1,outarray[2],1,outarray[2],1);
186  }
187  break;
190  {
191  // add to hu equation
192  Vmath::Vmul(nq,m_coriolis,1,physarray[2],1,tmp,1);
193  Vmath::Vadd(nq,tmp,1,outarray[1],1,outarray[1],1);
194 
195  // add to hv equation
196  Vmath::Vmul(nq,m_coriolis,1,physarray[1],1,tmp,1);
197  Vmath::Neg(nq,tmp,1);
198  Vmath::Vadd(nq,tmp,1,outarray[2],1,outarray[2],1);
199  }
200  break;
201  default:
202  ASSERTL0(false,"Unknown projection scheme for the NonlinearSWE");
203  break;
204  }
205 
206 
207  }
Array< OneD, NekDouble > m_coriolis
Coriolis force.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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:382
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:285
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:169
void Nektar::NonlinearSWE::AddVariableDepth ( const Array< OneD, const Array< OneD, NekDouble > > &  physarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
private

Definition at line 211 of file NonlinearSWE.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_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().

213  {
214 
215  int ncoeffs = GetNcoeffs();
216  int nq = GetTotPoints();
217 
218  Array<OneD, NekDouble> tmp(nq);
219  Array<OneD, NekDouble> mod(ncoeffs);
220 
221  switch(m_projectionType)
222  {
224  {
225  for (int i = 0; i < m_spacedim; ++i)
226  {
227  Vmath::Vmul(nq,m_bottomSlope[i],1,physarray[0],1,tmp,1);
228  Vmath::Smul(nq,m_g,tmp,1,tmp,1);
229  m_fields[0]->IProductWRTBase(tmp,mod);
230  m_fields[0]->MultiplyByElmtInvMass(mod,mod);
231  m_fields[0]->BwdTrans(mod,tmp);
232  Vmath::Vadd(nq,tmp,1,outarray[i+1],1,outarray[i+1],1);
233  }
234  }
235  break;
238  {
239  for (int i = 0; i < m_spacedim; ++i)
240  {
241  Vmath::Vmul(nq,m_bottomSlope[i],1,physarray[0],1,tmp,1);
242  Vmath::Smul(nq,m_g,tmp,1,tmp,1);
243  Vmath::Vadd(nq,tmp,1,outarray[i+1],1,outarray[i+1],1);
244  }
245  }
246  break;
247  default:
248  ASSERTL0(false,"Unknown projection scheme for the NonlinearSWE");
249  break;
250  }
251 
252 
253  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT int GetTotPoints()
Array< OneD, Array< OneD, NekDouble > > m_bottomSlope
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:199
int m_spacedim
Spatial dimension (>= expansion dim).
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetNcoeffs()
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:285
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:169
void Nektar::NonlinearSWE::ConservativeToPrimitive ( const Array< OneD, const Array< OneD, NekDouble > > &  physin,
Array< OneD, Array< OneD, NekDouble > > &  physout 
)
private

Definition at line 634 of file NonlinearSWE.cpp.

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

636  {
637  int nq = GetTotPoints();
638 
639  if(physin.get() == physout.get())
640  {
641  // copy indata and work with tmp array
643  for (int i = 0; i < 3; ++i)
644  {
645  // deep copy
646  tmp[i] = Array<OneD, NekDouble>(nq);
647  Vmath::Vcopy(nq,physin[i],1,tmp[i],1);
648  }
649 
650  // \eta = h - d
651  Vmath::Vsub(nq,tmp[0],1,m_depth,1,physout[0],1);
652 
653  // u = hu/h
654  Vmath::Vdiv(nq,tmp[1],1,tmp[0],1,physout[1],1);
655 
656  // v = hv/ v
657  Vmath::Vdiv(nq,tmp[2],1,tmp[0],1,physout[2],1);
658  }
659  else
660  {
661  // \eta = h - d
662  Vmath::Vsub(nq,physin[0],1,m_depth,1,physout[0],1);
663 
664  // u = hu/h
665  Vmath::Vdiv(nq,physin[1],1,physin[0],1,physout[1],1);
666 
667  // v = hv/ v
668  Vmath::Vdiv(nq,physin[2],1,physin[0],1,physout[2],1);
669  }
670  }
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:227
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:329
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1038
static SolverUtils::EquationSystemSharedPtr Nektar::NonlinearSWE::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Creates an instance of this class.

Definition at line 56 of file NonlinearSWE.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

58  {
60  p->InitObject();
61  return p;
62  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.
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 372 of file NonlinearSWE.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(), and Vmath::Vcopy().

Referenced by v_InitObject().

375  {
376  int i;
377  int nvariables = inarray.num_elements();
378 
379 
380  switch(m_projectionType)
381  {
383  {
384 
385  // Just copy over array
386  int npoints = GetNpoints();
387 
388  for(i = 0; i < nvariables; ++i)
389  {
390  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
391  }
392  SetBoundaryConditions(outarray, time);
393  break;
394  }
397  {
398 
399  EquationSystem::SetBoundaryConditions(time);
401 
402  for(i = 0; i < nvariables; ++i)
403  {
404  m_fields[i]->FwdTrans(inarray[i],coeffs);
405  m_fields[i]->BwdTrans_IterPerExp(coeffs,outarray[i]);
406  }
407  break;
408  }
409  default:
410  ASSERTL0(false,"Unknown projection scheme");
411  break;
412  }
413  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
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:1038
void SetBoundaryConditions(Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
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 255 of file NonlinearSWE.cpp.

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

258  {
259  int i, j;
260  int ndim = m_spacedim;
261  int nvariables = inarray.num_elements();
262  int nq = GetTotPoints();
263 
264 
265  switch(m_projectionType)
266  {
268  {
269 
270  //-------------------------------------------------------
271  // Compute the DG advection including the numerical flux
272  // by using SolverUtils/Advection
273  // Input and output in physical space
275 
276  m_advection->Advect(nvariables, m_fields, advVel, inarray,
277  outarray, time);
278  //-------------------------------------------------------
279 
280 
281  //-------------------------------------------------------
282  // negate the outarray since moving terms to the rhs
283  for(i = 0; i < nvariables; ++i)
284  {
285  Vmath::Neg(nq,outarray[i],1);
286  }
287  //-------------------------------------------------------
288 
289 
290  //-------------------------------------------------
291  // Add "source terms"
292  // Input and output in physical space
293 
294  // Coriolis forcing
295  if (m_coriolis.num_elements() != 0)
296  {
297  AddCoriolis(inarray,outarray);
298  }
299 
300  // Variable Depth
301  if (m_constantDepth != true)
302  {
303  AddVariableDepth(inarray,outarray);
304  }
305  //-------------------------------------------------
306 
307  }
308  break;
311  {
312 
313  //-------------------------------------------------------
314  // Compute the fluxvector in physical space
316  fluxvector(nvariables);
317 
318  for (i = 0; i < nvariables; ++i)
319  {
320  fluxvector[i] = Array<OneD, Array<OneD, NekDouble> >(ndim);
321  for(j = 0; j < ndim; ++j)
322  {
323  fluxvector[i][j] = Array<OneD, NekDouble>(nq);
324  }
325  }
326 
327  NonlinearSWE::GetFluxVector(inarray, fluxvector);
328  //-------------------------------------------------------
329 
330 
331 
332  //-------------------------------------------------------
333  // Take the derivative of the flux terms
334  // and negate the outarray since moving terms to the rhs
335  Array<OneD,NekDouble> tmp(nq);
336  Array<OneD, NekDouble>tmp1(nq);
337 
338  for(i = 0; i < nvariables; ++i)
339  {
340  m_fields[i]->PhysDeriv(MultiRegions::DirCartesianMap[0],fluxvector[i][0],tmp);
341  m_fields[i]->PhysDeriv(MultiRegions::DirCartesianMap[1],fluxvector[i][1],tmp1);
342  Vmath::Vadd(nq,tmp,1,tmp1,1,outarray[i],1);
343  Vmath::Neg(nq,outarray[i],1);
344  }
345 
346 
347  //-------------------------------------------------
348  // Add "source terms"
349  // Input and output in physical space
350 
351  // Coriolis forcing
352  if (m_coriolis.num_elements() != 0)
353  {
354  AddCoriolis(inarray,outarray);
355  }
356 
357  // Variable Depth
358  if (m_constantDepth != true)
359  {
360  AddVariableDepth(inarray,outarray);
361  }
362  //-------------------------------------------------
363  }
364  break;
365  default:
366  ASSERTL0(false,"Unknown projection scheme for the NonlinearSWE");
367  break;
368  }
369  }
Array< OneD, NekDouble > m_coriolis
Coriolis force.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SolverUtils::AdvectionSharedPtr m_advection
void AddVariableDepth(const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
SOLVER_UTILS_EXPORT int GetTotPoints()
void AddCoriolis(const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
int m_spacedim
Spatial dimension (>= expansion dim).
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:382
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
bool m_constantDepth
Indicates if constant depth case.
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:285
void Nektar::NonlinearSWE::GetFluxVector ( const Array< OneD, const Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  flux 
)
protected

Definition at line 595 of file NonlinearSWE.cpp.

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

598  {
599  int i, j;
600  int nq = m_fields[0]->GetTotPoints();
601 
602  NekDouble g = m_g;
604 
605  // Flux vector for the mass equation
606  for (i = 0; i < m_spacedim; ++i)
607  {
608  velocity[i] = Array<OneD, NekDouble>(nq);
609  Vmath::Vcopy(nq, physfield[i+1], 1, flux[0][i], 1);
610  }
611 
612  GetVelocityVector(physfield, velocity);
613 
614  // Put (0.5 g h h) in tmp
615  Array<OneD, NekDouble> tmp(nq);
616  Vmath::Vmul(nq, physfield[0], 1, physfield[0], 1, tmp, 1);
617  Vmath::Smul(nq, 0.5*g, tmp, 1, tmp, 1);
618 
619  // Flux vector for the momentum equations
620  for (i = 0; i < m_spacedim; ++i)
621  {
622  for (j = 0; j < m_spacedim; ++j)
623  {
624  Vmath::Vmul(nq, velocity[j], 1, physfield[i+1], 1,
625  flux[i+1][j], 1);
626  }
627 
628  // Add (0.5 g h h) to appropriate field
629  Vmath::Vadd(nq, flux[i+1][i], 1, tmp, 1, flux[i+1][i], 1);
630  }
631 
632  }
void GetVelocityVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)
Compute the velocity field given the momentum .
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:199
int m_spacedim
Spatial dimension (>= expansion dim).
double NekDouble
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
NekDouble m_g
Acceleration of gravity.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1038
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:285
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:169
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 751 of file NonlinearSWE.cpp.

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

Referenced by GetFluxVector().

754  {
755  const int npts = physfield[0].num_elements();
756 
757  for (int i = 0; i < m_spacedim; ++i)
758  {
759  Vmath::Vdiv(npts, physfield[1+i], 1, physfield[0], 1,
760  velocity[i], 1);
761  }
762  }
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:227
static std::string npts
Definition: InputFld.cpp:43
int m_spacedim
Spatial dimension (>= expansion dim).
void Nektar::NonlinearSWE::NumericalFlux1D ( Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  numfluxX 
)
private
void Nektar::NonlinearSWE::NumericalFlux2D ( Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  numfluxX,
Array< OneD, Array< OneD, NekDouble > > &  numfluxY 
)
private
void Nektar::NonlinearSWE::PrimitiveToConservative ( const Array< OneD, const Array< OneD, NekDouble > > &  physin,
Array< OneD, Array< OneD, NekDouble > > &  physout 
)
private

Definition at line 687 of file NonlinearSWE.cpp.

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

689  {
690 
691  int nq = GetTotPoints();
692 
693  if(physin.get() == physout.get())
694  {
695  // copy indata and work with tmp array
697  for (int i = 0; i < 3; ++i)
698  {
699  // deep copy
700  tmp[i] = Array<OneD, NekDouble>(nq);
701  Vmath::Vcopy(nq,physin[i],1,tmp[i],1);
702  }
703 
704  // h = \eta + d
705  Vmath::Vadd(nq,tmp[0],1,m_depth,1,physout[0],1);
706 
707  // hu = h * u
708  Vmath::Vmul(nq,physout[0],1,tmp[1],1,physout[1],1);
709 
710  // hv = h * v
711  Vmath::Vmul(nq,physout[0],1,tmp[2],1,physout[2],1);
712 
713  }
714  else
715  {
716  // h = \eta + d
717  Vmath::Vadd(nq,physin[0],1,m_depth,1,physout[0],1);
718 
719  // hu = h * u
720  Vmath::Vmul(nq,physout[0],1,physin[1],1,physout[1],1);
721 
722  // hv = h * v
723  Vmath::Vmul(nq,physout[0],1,physin[2],1,physout[2],1);
724 
725  }
726 
727  }
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:1038
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:285
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:169
void Nektar::NonlinearSWE::SetBoundaryConditions ( Array< OneD, Array< OneD, NekDouble > > &  physarray,
NekDouble  time 
)
private

Definition at line 417 of file NonlinearSWE.cpp.

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

Referenced by DoOdeProjection().

420  {
421  std::string varName;
422  int nvariables = m_fields.num_elements();
423  int cnt = 0;
424 
425  // Loop over Boundary Regions
426  for (int n = 0; n < m_fields[0]->GetBndConditions().num_elements(); ++n)
427  {
428 
429  // Wall Boundary Condition
430  if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),"Wall"))
431  {
432  WallBoundary2D(n, cnt, inarray);
433  }
434 
435  // Time Dependent Boundary Condition (specified in meshfile)
436  if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
437  {
438  for (int i = 0; i < nvariables; ++i)
439  {
440  varName = m_session->GetVariable(i);
441  m_fields[i]->EvaluateBoundaryConditions(time, varName);
442  }
443  }
444  cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
445  }
446  }
void WallBoundary2D(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &physarray)
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
void Nektar::NonlinearSWE::v_ConservativeToPrimitive ( )
protectedvirtual

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 673 of file NonlinearSWE.cpp.

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

674  {
675  int nq = GetTotPoints();
676 
677  // u = hu/h
678  Vmath::Vdiv(nq,m_fields[1]->GetPhys(),1,m_fields[0]->GetPhys(),1,m_fields[1]->UpdatePhys(),1);
679 
680  // v = hv/ v
681  Vmath::Vdiv(nq,m_fields[2]->GetPhys(),1,m_fields[0]->GetPhys(),1,m_fields[2]->UpdatePhys(),1);
682 
683  // \eta = h - d
684  Vmath::Vsub(nq,m_fields[0]->GetPhys(),1,m_depth,1,m_fields[0]->UpdatePhys(),1);
685  }
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:227
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:329
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Nektar::NonlinearSWE::v_GenerateSummary ( SolverUtils::SummaryList s)
protectedvirtual

Print a summary of time stepping parameters.

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 765 of file NonlinearSWE.cpp.

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

766  {
768  SolverUtils::AddSummaryItem(s, "Variables", "h should be in field[0]");
769  SolverUtils::AddSummaryItem(s, "", "hu should be in field[1]");
770  SolverUtils::AddSummaryItem(s, "", "hv should be in field[2]");
771  }
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:50
virtual void v_GenerateSummary(SolverUtils::SummaryList &s)
Print a summary of time stepping parameters.
void Nektar::NonlinearSWE::v_InitObject ( )
protectedvirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 56 of file NonlinearSWE.cpp.

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, >::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().

57  {
59 
61  {
64  }
65  else
66  {
67  ASSERTL0(false, "Implicit SWE not set up.");
68  }
69 
70  // Type of advection class to be used
71  switch(m_projectionType)
72  {
73  // Continuous field
75  {
76  // Do nothing
77  break;
78  }
79  // Discontinuous field
81  {
82  string advName;
83  string diffName;
84  string riemName;
85 
86  //---------------------------------------------------------------
87  // Setting up advection and diffusion operators
88  // NB: diffusion not set up for SWE at the moment
89  // but kept here for future use ...
90  m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
91  // m_session->LoadSolverInfo("DiffusionType", diffName, "LDGEddy");
93  .CreateInstance(advName, advName);
94  // m_diffusion = SolverUtils::GetDiffusionFactory()
95  // .CreateInstance(diffName, diffName);
96 
97  m_advection->SetFluxVector(&NonlinearSWE::GetFluxVector, this);
98  // m_diffusion->SetFluxVectorNS(&ShallowWaterSystem::
99  // GetEddyViscosityFluxVector, this);
100 
101  // Setting up Riemann solver for advection operator
102  m_session->LoadSolverInfo("UpwindType", riemName, "Average");
104  .CreateInstance(riemName);
105 
106  // Setting up upwind solver for diffusion operator
107  // m_riemannSolverLDG = SolverUtils::GetRiemannSolverFactory()
108  // .CreateInstance("UpwindLDG");
109 
110  // Setting up parameters for advection operator Riemann solver
111  m_riemannSolver->SetParam (
112  "gravity",
113  &NonlinearSWE::GetGravity, this);
114  m_riemannSolver->SetAuxVec(
115  "vecLocs",
116  &NonlinearSWE::GetVecLocs, this);
117  m_riemannSolver->SetVector(
118  "N",
119  &NonlinearSWE::GetNormals, this);
120  m_riemannSolver->SetScalar(
121  "depth",
122  &NonlinearSWE::GetDepth, this);
123 
124  // Setting up parameters for diffusion operator Riemann solver
125  // m_riemannSolverLDG->AddParam (
126  // "gravity",
127  // &NonlinearSWE::GetGravity, this);
128  // m_riemannSolverLDG->SetAuxVec(
129  // "vecLocs",
130  // &NonlinearSWE::GetVecLocs, this);
131  // m_riemannSolverLDG->AddVector(
132  // "N",
133  // &NonlinearSWE::GetNormals, this);
134 
135  // Concluding initialisation of advection / diffusion operators
136  m_advection->SetRiemannSolver (m_riemannSolver);
137  //m_diffusion->SetRiemannSolver (m_riemannSolverLDG);
138  m_advection->InitObject (m_session, m_fields);
139  //m_diffusion->InitObject (m_session, m_fields);
140  break;
141  }
142  default:
143  {
144  ASSERTL0(false, "Unsupported projection type.");
145  break;
146  }
147  }
148 
149 
150  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
tBaseSharedPtr CreateInstance(tKey idKey BOOST_PP_COMMA_IF(MAX_PARAM) BOOST_PP_ENUM_BINARY_PARAMS(MAX_PARAM, tParam, x))
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:162
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
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()
const Array< OneD, NekDouble > & GetDepth()
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
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
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:46
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void Nektar::NonlinearSWE::v_PrimitiveToConservative ( )
protectedvirtual

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 729 of file NonlinearSWE.cpp.

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

730  {
731  int nq = GetTotPoints();
732 
733  // h = \eta + d
734  Vmath::Vadd(nq,m_fields[0]->GetPhys(),1,m_depth,1,m_fields[0]->UpdatePhys(),1);
735 
736  // hu = h * u
737  Vmath::Vmul(nq,m_fields[0]->GetPhys(),1,m_fields[1]->GetPhys(),1,m_fields[1]->UpdatePhys(),1);
738 
739  // hv = h * v
740  Vmath::Vmul(nq,m_fields[0]->GetPhys(),1,m_fields[2]->GetPhys(),1,m_fields[2]->UpdatePhys(),1);
741  }
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:285
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:169
void Nektar::NonlinearSWE::WallBoundary ( int  bcRegion,
int  cnt,
Array< OneD, Array< OneD, NekDouble > > &  physarray 
)
private

Wall boundary condition.

Definition at line 452 of file NonlinearSWE.cpp.

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

456  {
457  int i;
458  int nTracePts = GetTraceTotPoints();
459  int nvariables = physarray.num_elements();
460 
461  // get physical values of the forward trace
462  Array<OneD, Array<OneD, NekDouble> > Fwd(nvariables);
463  for (i = 0; i < nvariables; ++i)
464  {
465  Fwd[i] = Array<OneD, NekDouble>(nTracePts);
466  m_fields[i]->ExtractTracePhys(physarray[i], Fwd[i]);
467  }
468 
469  // Adjust the physical values of the trace to take
470  // user defined boundaries into account
471  int e, id1, id2, npts;
472 
473  for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]
474  ->GetExpSize(); ++e)
475  {
476  npts = m_fields[0]->GetBndCondExpansions()[bcRegion]->
477  GetExp(e)->GetTotPoints();
478  id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->
479  GetPhys_Offset(e);
480  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
481  m_fields[0]->GetTraceMap()->
482  GetBndCondCoeffsToGlobalCoeffsMap(cnt+e));
483 
484  // For 2D/3D, define: v* = v - 2(v.n)n
485  Array<OneD, NekDouble> tmp(npts, 0.0);
486 
487  // Calculate (v.n)
488  for (i = 0; i < m_spacedim; ++i)
489  {
490  Vmath::Vvtvp(npts,
491  &Fwd[1+i][id2], 1,
492  &m_traceNormals[i][id2], 1,
493  &tmp[0], 1,
494  &tmp[0], 1);
495  }
496 
497  // Calculate 2.0(v.n)
498  Vmath::Smul(npts, -2.0, &tmp[0], 1, &tmp[0], 1);
499 
500  // Calculate v* = v - 2.0(v.n)n
501  for (i = 0; i < m_spacedim; ++i)
502  {
503  Vmath::Vvtvp(npts,
504  &tmp[0], 1,
505  &m_traceNormals[i][id2], 1,
506  &Fwd[1+i][id2], 1,
507  &Fwd[1+i][id2], 1);
508  }
509 
510  // copy boundary adjusted values into the boundary expansion
511  for (i = 0; i < nvariables; ++i)
512  {
513  Vmath::Vcopy(npts, &Fwd[i][id2], 1,
514  &(m_fields[i]->GetBndCondExpansions()[bcRegion]->
515  UpdatePhys())[id1], 1);
516  }
517  }
518  }
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:428
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:199
static std::string npts
Definition: InputFld.cpp:43
int m_spacedim
Spatial dimension (>= expansion dim).
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
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:1038
void Nektar::NonlinearSWE::WallBoundary2D ( int  bcRegion,
int  cnt,
Array< OneD, Array< OneD, NekDouble > > &  physarray 
)
private

Definition at line 521 of file NonlinearSWE.cpp.

References ASSERTL0, Nektar::SolverUtils::EquationSystem::GetTraceTotPoints(), Nektar::SolverUtils::EquationSystem::m_expdim, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_traceNormals, Vmath::Neg(), npts, Vmath::Vcopy(), Vmath::Vmul(), Vmath::Vvtvm(), and Vmath::Vvtvp().

Referenced by SetBoundaryConditions().

522  {
523 
524  int i;
525  int nTraceNumPoints = GetTraceTotPoints();
526  int nvariables = physarray.num_elements();
527 
528  // get physical values of the forward trace
529  Array<OneD, Array<OneD, NekDouble> > Fwd(nvariables);
530  for (i = 0; i < nvariables; ++i)
531  {
532  Fwd[i] = Array<OneD, NekDouble>(nTraceNumPoints);
533  m_fields[i]->ExtractTracePhys(physarray[i],Fwd[i]);
534  }
535 
536  // Adjust the physical values of the trace to take
537  // user defined boundaries into account
538  int e, id1, id2, npts;
539 
540  for(e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize(); ++e)
541  {
542  npts = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExp(e)->GetNumPoints(0);
543  id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e) ;
544  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(m_fields[0]->GetTraceMap()->GetBndCondCoeffsToGlobalCoeffsMap(cnt+e));
545 
546  switch(m_expdim)
547  {
548  case 1:
549  {
550  // negate the forward flux
551  Vmath::Neg(npts,&Fwd[1][id2],1);
552  }
553  break;
554  case 2:
555  {
556  Array<OneD, NekDouble> tmp_n(npts);
557  Array<OneD, NekDouble> tmp_t(npts);
558 
559  Vmath::Vmul(npts,&Fwd[1][id2],1,&m_traceNormals[0][id2],1,&tmp_n[0],1);
560  Vmath::Vvtvp(npts,&Fwd[2][id2],1,&m_traceNormals[1][id2],1,&tmp_n[0],1,&tmp_n[0],1);
561 
562  Vmath::Vmul(npts,&Fwd[1][id2],1,&m_traceNormals[1][id2],1,&tmp_t[0],1);
563  Vmath::Vvtvm(npts,&Fwd[2][id2],1,&m_traceNormals[0][id2],1,&tmp_t[0],1,&tmp_t[0],1);
564 
565  // negate the normal flux
566  Vmath::Neg(npts,tmp_n,1);
567 
568  // rotate back to Cartesian
569  Vmath::Vmul(npts,&tmp_t[0],1,&m_traceNormals[1][id2],1,&Fwd[1][id2],1);
570  Vmath::Vvtvm(npts,&tmp_n[0],1,&m_traceNormals[0][id2],1,&Fwd[1][id2],1,&Fwd[1][id2],1);
571 
572  Vmath::Vmul(npts,&tmp_t[0],1,&m_traceNormals[0][id2],1,&Fwd[2][id2],1);
573  Vmath::Vvtvp(npts,&tmp_n[0],1,&m_traceNormals[1][id2],1,&Fwd[2][id2],1,&Fwd[2][id2],1);
574  }
575  break;
576  case 3:
577  ASSERTL0(false,"3D not implemented for Shallow Water Equations");
578  break;
579  default:
580  ASSERTL0(false,"Illegal expansion dimension");
581  }
582 
583 
584 
585  // copy boundary adjusted values into the boundary expansion
586  for (i = 0; i < nvariables; ++i)
587  {
588  Vmath::Vcopy(npts,&Fwd[i][id2], 1,&(m_fields[i]->GetBndCondExpansions()[bcRegion]->UpdatePhys())[id1],1);
589  }
590  }
591  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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:428
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
static std::string npts
Definition: InputFld.cpp:43
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:382
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
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:451
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1038
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:169

Friends And Related Function Documentation

friend class MemoryManager< NonlinearSWE >
friend

Definition at line 53 of file NonlinearSWE.h.

Member Data Documentation

string Nektar::NonlinearSWE::className
static
Initial value:
=
"NonlinearSWE", NonlinearSWE::create,
"Nonlinear shallow water equation in conservative variables.")

Name of class.

Definition at line 64 of file NonlinearSWE.h.