Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties 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::LinearSWE Class Reference

#include <LinearSWE.h>

Inheritance diagram for Nektar::LinearSWE:
Inheritance graph
[legend]
Collaboration diagram for Nektar::LinearSWE:
Collaboration graph
[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 >
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 NekDouble &pTime=0.0, 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 NekDouble &pTime=0.0, 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, NekDouble
ErrorExtraPoints (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::FieldMetaDataMap
UpdateFieldMetaDataMap ()
 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 SetSteps (const int steps)
 
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 std::string &s1, const std::string &s2)
 Perform a case-insensitive string comparison. More...
 
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)
 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

 LinearSWE (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 ()
 
const Array< OneD, NekDouble > & GetDepthFwd ()
 
const Array< OneD, NekDouble > & GetDepthBwd ()
 
- 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)
 
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)
 Initialises EquationSystem class members. More...
 
int nocase_cmp (const std::string &s1, const std::string &s2)
 
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 EvaluateFunctionExp (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
 
SOLVER_UTILS_EXPORT void EvaluateFunctionFld (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
 
SOLVER_UTILS_EXPORT void EvaluateFunctionPts (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
 
SOLVER_UTILS_EXPORT void LoadPts (std::string funcFilename, std::string filename, Nektar::LibUtilities::PtsFieldSharedPtr &outPts)
 
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

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_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...
 
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...
 
std::map< std::string,
FieldUtils::Interpolator
m_interpolators
 Map of interpolator objects. More...
 
std::map< std::string,
std::pair< std::string,
LibUtilities::PtsFieldSharedPtr > > 
m_loadedPtsFields
 pts fields we already read from disk: {funcFilename: (filename, ptsfield)} More...
 
std::map< std::string,
std::pair< std::string,
loadedFldField > > 
m_loadedFldFields
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_fields
 Array holding all dependent variables. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_base
 Base fields. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_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...
 
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, 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...
 

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

Detailed Description

Definition at line 50 of file LinearSWE.h.

Constructor & Destructor Documentation

Nektar::LinearSWE::~LinearSWE ( )
virtual

Definition at line 168 of file LinearSWE.cpp.

169  {
170 
171  }
Nektar::LinearSWE::LinearSWE ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Definition at line 52 of file LinearSWE.cpp.

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

Member Function Documentation

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

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:198
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:396
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:299
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:183
void Nektar::LinearSWE::ConservativeToPrimitive ( const Array< OneD, const Array< OneD, NekDouble > > &  physin,
Array< OneD, Array< OneD, NekDouble > > &  physout 
)
private

Definition at line 577 of file LinearSWE.cpp.

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

579  {
580  int nq = GetTotPoints();
581 
582  if(physin.get() == physout.get())
583  {
584  // copy indata and work with tmp array
585  Array<OneD, Array<OneD, NekDouble> >tmp(3);
586  for (int i = 0; i < 3; ++i)
587  {
588  // deep copy
589  tmp[i] = Array<OneD, NekDouble>(nq);
590  Vmath::Vcopy(nq,physin[i],1,tmp[i],1);
591  }
592 
593  // \eta = h - d
594  Vmath::Vsub(nq,tmp[0],1,m_depth,1,physout[0],1);
595 
596  // u = hu/h
597  Vmath::Vdiv(nq,tmp[1],1,tmp[0],1,physout[1],1);
598 
599  // v = hv/ v
600  Vmath::Vdiv(nq,tmp[2],1,tmp[0],1,physout[2],1);
601  }
602  else
603  {
604  // \eta = h - d
605  Vmath::Vsub(nq,physin[0],1,m_depth,1,physout[0],1);
606 
607  // u = hu/h
608  Vmath::Vdiv(nq,physin[1],1,physin[0],1,physout[1],1);
609 
610  // v = hv/ v
611  Vmath::Vdiv(nq,physin[2],1,physin[0],1,physout[2],1);
612  }
613  }
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:241
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:343
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
static SolverUtils::EquationSystemSharedPtr Nektar::LinearSWE::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Creates an instance of this class.

Definition at line 56 of file LinearSWE.h.

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

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::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:198
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:1061
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:198
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Definition: LinearSWE.cpp:543
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:396
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
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:299
const Array<OneD, NekDouble>& Nektar::LinearSWE::GetDepthBwd ( )
inlineprotected

Definition at line 100 of file LinearSWE.h.

References m_dBwd.

Referenced by v_InitObject().

101  {
102  return m_dBwd;
103  }
Array< OneD, NekDouble > m_dBwd
Definition: LinearSWE.h:76
const Array<OneD, NekDouble>& Nektar::LinearSWE::GetDepthFwd ( )
inlineprotected

Definition at line 96 of file LinearSWE.h.

References m_dFwd.

Referenced by v_InitObject().

97  {
98  return m_dFwd;
99  }
Array< OneD, NekDouble > m_dFwd
Still water depth traces.
Definition: LinearSWE.h:75
void Nektar::LinearSWE::GetFluxVector ( const Array< OneD, const Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  flux 
)
protected

Definition at line 543 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().

546  {
547  int i, j;
548  int nq = m_fields[0]->GetTotPoints();
549 
550  NekDouble g = m_g;
551 
552  // Flux vector for the mass equation
553  for (i = 0; i < m_spacedim; ++i)
554  {
555  Vmath::Vmul(nq, m_depth, 1, physfield[i+1], 1, flux[0][i], 1);
556  }
557 
558  // Put (g eta) in tmp
559  Array<OneD, NekDouble> tmp(nq);
560  Vmath::Smul(nq, g, physfield[0], 1, tmp, 1);
561 
562  // Flux vector for the momentum equations
563  for (i = 0; i < m_spacedim; ++i)
564  {
565  for (j = 0; j < m_spacedim; ++j)
566  {
567  // must zero fluxes as not initialised to zero in AdvectionWeakDG ...
568  Vmath::Zero(nq, flux[i+1][j], 1);
569  }
570 
571  // Add (g eta) to appropriate field
572  Vmath::Vadd(nq, flux[i+1][i], 1, tmp, 1, flux[i+1][i], 1);
573  }
574 
575  }
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:213
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:373
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:299
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:183
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 694 of file LinearSWE.cpp.

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

697  {
698  const int npts = physfield[0].num_elements();
699 
700  for (int i = 0; i < m_spacedim; ++i)
701  {
702  Vmath::Vcopy(npts, physfield[1+i], 1, velocity[i], 1);
703  }
704  }
static std::string npts
Definition: InputFld.cpp:43
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:1061
void Nektar::LinearSWE::NumericalFlux1D ( Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  numfluxX 
)
private
void Nektar::LinearSWE::NumericalFlux2D ( Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  numfluxX,
Array< OneD, Array< OneD, NekDouble > > &  numfluxY 
)
private
void Nektar::LinearSWE::PrimitiveToConservative ( const Array< OneD, const Array< OneD, NekDouble > > &  physin,
Array< OneD, Array< OneD, NekDouble > > &  physout 
)
private

Definition at line 630 of file LinearSWE.cpp.

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

632  {
633 
634  int nq = GetTotPoints();
635 
636  if(physin.get() == physout.get())
637  {
638  // copy indata and work with tmp array
639  Array<OneD, Array<OneD, NekDouble> >tmp(3);
640  for (int i = 0; i < 3; ++i)
641  {
642  // deep copy
643  tmp[i] = Array<OneD, NekDouble>(nq);
644  Vmath::Vcopy(nq,physin[i],1,tmp[i],1);
645  }
646 
647  // h = \eta + d
648  Vmath::Vadd(nq,tmp[0],1,m_depth,1,physout[0],1);
649 
650  // hu = h * u
651  Vmath::Vmul(nq,physout[0],1,tmp[1],1,physout[1],1);
652 
653  // hv = h * v
654  Vmath::Vmul(nq,physout[0],1,tmp[2],1,physout[2],1);
655 
656  }
657  else
658  {
659  // h = \eta + d
660  Vmath::Vadd(nq,physin[0],1,m_depth,1,physout[0],1);
661 
662  // hu = h * u
663  Vmath::Vmul(nq,physout[0],1,physin[1],1,physout[1],1);
664 
665  // hv = h * v
666  Vmath::Vmul(nq,physout[0],1,physin[2],1,physout[2],1);
667 
668  }
669 
670  }
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:1061
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:299
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:183
void Nektar::LinearSWE::SetBoundaryConditions ( Array< OneD, Array< OneD, NekDouble > > &  physarray,
NekDouble  time 
)
private

Definition at line 373 of file LinearSWE.cpp.

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

Referenced by DoOdeProjection().

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  // Wall Boundary Condition
395  if (boost::iequals(m_fields[0]->GetBndConditions()[n]->GetUserDefined(),"Wall"))
396  {
397  WallBoundary2D(n, cnt, Fwd, inarray);
398  }
399 
400  // Time Dependent Boundary Condition (specified in meshfile)
401  if (m_fields[0]->GetBndConditions()[n]->IsTimeDependent())
402  {
403  for (int i = 0; i < nvariables; ++i)
404  {
405  varName = m_session->GetVariable(i);
406  m_fields[i]->EvaluateBoundaryConditions(time, varName);
407  }
408  }
409  cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
410  }
411  }
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:478
void Nektar::LinearSWE::v_ConservativeToPrimitive ( )
protectedvirtual

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 616 of file LinearSWE.cpp.

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

617  {
618  int nq = GetTotPoints();
619 
620  // u = hu/h
621  Vmath::Vdiv(nq,m_fields[1]->GetPhys(),1,m_fields[0]->GetPhys(),1,m_fields[1]->UpdatePhys(),1);
622 
623  // v = hv/ v
624  Vmath::Vdiv(nq,m_fields[2]->GetPhys(),1,m_fields[0]->GetPhys(),1,m_fields[2]->UpdatePhys(),1);
625 
626  // \eta = h - d
627  Vmath::Vsub(nq,m_fields[0]->GetPhys(),1,m_depth,1,m_fields[0]->UpdatePhys(),1);
628  }
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:241
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:343
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Nektar::LinearSWE::v_GenerateSummary ( SolverUtils::SummaryList s)
protectedvirtual

Print a summary of time stepping parameters.

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 707 of file LinearSWE.cpp.

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

708  {
710  if (m_session->DefinesSolverInfo("UpwindType"))
711  {
712  std::string UpwindType;
713  UpwindType = m_session->GetSolverInfo("UpwindType");
714  if (UpwindType == "LinearAverage")
715  {
716  SolverUtils::AddSummaryItem(s, "Riemann Solver", "Linear Average");
717  }
718  if (UpwindType == "LinearHLL")
719  {
720  SolverUtils::AddSummaryItem(s, "Riemann Solver", "Linear HLL");
721  }
722  }
723  SolverUtils::AddSummaryItem(s, "Variables", "eta should be in field[0]");
724  SolverUtils::AddSummaryItem(s, "", "u should be in field[1]");
725  SolverUtils::AddSummaryItem(s, "", "v should be in field[2]");
726  }
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.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
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, >::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);
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:96
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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
Array< OneD, NekDouble > m_depth
Still water depth.
Array< OneD, NekDouble > m_dFwd
Still water depth traces.
Definition: LinearSWE.h:75
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Definition: LinearSWE.cpp:543
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)
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:46
Array< OneD, NekDouble > m_dBwd
Definition: LinearSWE.h:76
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:100
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Definition: LinearSWE.cpp:225
void Nektar::LinearSWE::v_PrimitiveToConservative ( )
protectedvirtual

Reimplemented from Nektar::ShallowWaterSystem.

Definition at line 672 of file LinearSWE.cpp.

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

673  {
674  int nq = GetTotPoints();
675 
676  // h = \eta + d
677  Vmath::Vadd(nq,m_fields[0]->GetPhys(),1,m_depth,1,m_fields[0]->UpdatePhys(),1);
678 
679  // hu = h * u
680  Vmath::Vmul(nq,m_fields[0]->GetPhys(),1,m_fields[1]->GetPhys(),1,m_fields[1]->UpdatePhys(),1);
681 
682  // hv = h * v
683  Vmath::Vmul(nq,m_fields[0]->GetPhys(),1,m_fields[2]->GetPhys(),1,m_fields[2]->UpdatePhys(),1);
684  }
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:299
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:183
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 417 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, npts, Vmath::Smul(), Vmath::Vcopy(), and Vmath::Vvtvp().

422  {
423  int i;
424  int nvariables = physarray.num_elements();
425 
426  // Adjust the physical values of the trace to take
427  // user defined boundaries into account
428  int e, id1, id2, npts;
429 
430  for (e = 0; e < m_fields[0]->GetBndCondExpansions()[bcRegion]
431  ->GetExpSize(); ++e)
432  {
433  npts = m_fields[0]->GetBndCondExpansions()[bcRegion]->
434  GetExp(e)->GetTotPoints();
435  id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->
436  GetPhys_Offset(e);
437  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
438  m_fields[0]->GetTraceMap()->
439  GetBndCondCoeffsToGlobalCoeffsMap(cnt+e));
440 
441  // For 2D/3D, define: v* = v - 2(v.n)n
442  Array<OneD, NekDouble> tmp(npts, 0.0);
443 
444  // Calculate (v.n)
445  for (i = 0; i < m_spacedim; ++i)
446  {
447  Vmath::Vvtvp(npts,
448  &Fwd[1+i][id2], 1,
449  &m_traceNormals[i][id2], 1,
450  &tmp[0], 1,
451  &tmp[0], 1);
452  }
453 
454  // Calculate 2.0(v.n)
455  Vmath::Smul(npts, -2.0, &tmp[0], 1, &tmp[0], 1);
456 
457  // Calculate v* = v - 2.0(v.n)n
458  for (i = 0; i < m_spacedim; ++i)
459  {
460  Vmath::Vvtvp(npts,
461  &tmp[0], 1,
462  &m_traceNormals[i][id2], 1,
463  &Fwd[1+i][id2], 1,
464  &Fwd[1+i][id2], 1);
465  }
466 
467  // copy boundary adjusted values into the boundary expansion
468  for (i = 0; i < nvariables; ++i)
469  {
470  Vmath::Vcopy(npts, &Fwd[i][id2], 1,
471  &(m_fields[i]->GetBndCondExpansions()[bcRegion]->
472  UpdatePhys())[id1], 1);
473  }
474  }
475  }
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:442
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:213
static std::string npts
Definition: InputFld.cpp:43
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:1061
void Nektar::LinearSWE::WallBoundary2D ( int  bcRegion,
int  cnt,
Array< OneD, Array< OneD, NekDouble > > &  Fwd,
Array< OneD, Array< OneD, NekDouble > > &  physarray 
)
private

Definition at line 478 of file LinearSWE.cpp.

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

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

Friends And Related Function Documentation

friend class MemoryManager< LinearSWE >
friend

Definition at line 53 of file LinearSWE.h.

Member Data Documentation

string Nektar::LinearSWE::className
static
Initial value:
=
"LinearSWE", LinearSWE::create,
"Linear shallow water equation in primitive variables.")

Name of class.

Definition at line 64 of file LinearSWE.h.

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

Definition at line 76 of file LinearSWE.h.

Referenced by GetDepthBwd(), and v_InitObject().

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

Still water depth traces.

Definition at line 75 of file LinearSWE.h.

Referenced by GetDepthFwd(), and v_InitObject().