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

#include <NonlinearSWESolver.h>

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

Protected Member Functions

 NonlinearSWESolver ()
virtual void v_Solve (const int nDim, const Array< OneD, const Array< OneD, NekDouble > > &Fwd, const Array< OneD, const Array< OneD, NekDouble > > &Bwd, Array< OneD, Array< OneD, NekDouble > > &flux)
virtual void v_ArraySolve (const Array< OneD, const Array< OneD, NekDouble > > &Fwd, const Array< OneD, const Array< OneD, NekDouble > > &Bwd, Array< OneD, Array< OneD, NekDouble > > &flux)
virtual void v_PointSolve (NekDouble hL, NekDouble huL, NekDouble hvL, NekDouble hR, NekDouble huR, NekDouble hvR, NekDouble &hf, NekDouble &huf, NekDouble &hvf)
- Protected Member Functions inherited from Nektar::SolverUtils::RiemannSolver
SOLVER_UTILS_EXPORT RiemannSolver ()
void GenerateRotationMatrices (const Array< OneD, const Array< OneD, NekDouble > > &normals)
 Generate rotation matrices for 3D expansions.
void FromToRotation (Array< OneD, const NekDouble > &from, Array< OneD, const NekDouble > &to, NekDouble *mat)
 A function for creating a rotation matrix that rotates a vector from into another vector to.
SOLVER_UTILS_EXPORT void rotateToNormal (const Array< OneD, const Array< OneD, NekDouble > > &inarray, const Array< OneD, const Array< OneD, NekDouble > > &normals, const Array< OneD, const Array< OneD, NekDouble > > &vecLocs, Array< OneD, Array< OneD, NekDouble > > &outarray)
 Rotate a vector field to trace normal.
SOLVER_UTILS_EXPORT void rotateFromNormal (const Array< OneD, const Array< OneD, NekDouble > > &inarray, const Array< OneD, const Array< OneD, NekDouble > > &normals, const Array< OneD, const Array< OneD, NekDouble > > &vecLocs, Array< OneD, Array< OneD, NekDouble > > &outarray)
 Rotate a vector field from trace normal.
bool CheckScalars (std::string name)
 Determine whether a scalar has been defined in m_scalars.
bool CheckVectors (std::string name)
 Determine whether a vector has been defined in m_vectors.
bool CheckParams (std::string name)
 Determine whether a parameter has been defined in m_params.
bool CheckAuxScal (std::string name)
 Determine whether a scalar has been defined in m_auxScal.
bool CheckAuxVec (std::string name)
 Determine whether a vector has been defined in m_auxVec.

Protected Attributes

bool m_pointSolve
- Protected Attributes inherited from Nektar::SolverUtils::RiemannSolver
bool m_requiresRotation
 Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields.
std::map< std::string,
RSScalarFuncType
m_scalars
 Map of scalar function types.
std::map< std::string,
RSVecFuncType
m_vectors
 Map of vector function types.
std::map< std::string,
RSParamFuncType
m_params
 Map of parameter function types.
std::map< std::string,
RSScalarFuncType
m_auxScal
 Map of auxiliary scalar function types.
std::map< std::string,
RSVecFuncType
m_auxVec
 Map of auxiliary vector function types.
Array< OneD, Array< OneD,
NekDouble > > 
m_rotMat
 Rotation matrices for each trace quadrature point.
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_rotStorage
 Rotation storage.

Additional Inherited Members

- Public Member Functions inherited from Nektar::SolverUtils::RiemannSolver
SOLVER_UTILS_EXPORT void Solve (const int nDim, const Array< OneD, const Array< OneD, NekDouble > > &Fwd, const Array< OneD, const Array< OneD, NekDouble > > &Bwd, Array< OneD, Array< OneD, NekDouble > > &flux)
 Perform the Riemann solve given the forwards and backwards spaces.
template<typename FuncPointerT , typename ObjectPointerT >
void SetScalar (std::string name, FuncPointerT func, ObjectPointerT obj)
void SetScalar (std::string name, RSScalarFuncType fp)
template<typename FuncPointerT , typename ObjectPointerT >
void SetVector (std::string name, FuncPointerT func, ObjectPointerT obj)
void SetVector (std::string name, RSVecFuncType fp)
template<typename FuncPointerT , typename ObjectPointerT >
void SetParam (std::string name, FuncPointerT func, ObjectPointerT obj)
void SetParam (std::string name, RSParamFuncType fp)
template<typename FuncPointerT , typename ObjectPointerT >
void SetAuxScal (std::string name, FuncPointerT func, ObjectPointerT obj)
template<typename FuncPointerT , typename ObjectPointerT >
void SetAuxVec (std::string name, FuncPointerT func, ObjectPointerT obj)
std::map< std::string,
RSScalarFuncType > & 
GetScalars ()
std::map< std::string,
RSVecFuncType > & 
GetVectors ()
std::map< std::string,
RSParamFuncType > & 
GetParams ()
- Public Attributes inherited from Nektar::SolverUtils::RiemannSolver
int m_spacedim

Detailed Description

Definition at line 45 of file NonlinearSWESolver.h.

Constructor & Destructor Documentation

Nektar::NonlinearSWESolver::NonlinearSWESolver ( )
protected

Member Function Documentation

virtual void Nektar::NonlinearSWESolver::v_ArraySolve ( const Array< OneD, const Array< OneD, NekDouble > > &  Fwd,
const Array< OneD, const Array< OneD, NekDouble > > &  Bwd,
Array< OneD, Array< OneD, NekDouble > > &  flux 
)
inlineprotectedvirtual

Reimplemented in Nektar::AverageSolver.

Definition at line 58 of file NonlinearSWESolver.h.

References ASSERTL0.

Referenced by v_Solve().

{
ASSERTL0(false, "This function should be defined by subclasses.");
}
virtual void Nektar::NonlinearSWESolver::v_PointSolve ( NekDouble  hL,
NekDouble  huL,
NekDouble  hvL,
NekDouble  hR,
NekDouble  huR,
NekDouble  hvR,
NekDouble hf,
NekDouble huf,
NekDouble hvf 
)
inlineprotectedvirtual

Reimplemented in Nektar::AverageSolver, Nektar::HLLCSolver, Nektar::HLLSolver, and Nektar::LaxFriedrichsSolver.

Definition at line 66 of file NonlinearSWESolver.h.

References ASSERTL0.

Referenced by v_Solve().

{
ASSERTL0(false, "This function should be defined by subclasses.");
}
void Nektar::NonlinearSWESolver::v_Solve ( const int  nDim,
const Array< OneD, const Array< OneD, NekDouble > > &  Fwd,
const Array< OneD, const Array< OneD, NekDouble > > &  Bwd,
Array< OneD, Array< OneD, NekDouble > > &  flux 
)
protectedvirtual

Implements Nektar::SolverUtils::RiemannSolver.

Definition at line 46 of file NonlinearSWESolver.cpp.

References ASSERTL0, m_pointSolve, v_ArraySolve(), and v_PointSolve().

{
{
int expDim = Fwd.num_elements()-1;
NekDouble hvf;
if (expDim == 1)
{
for (int i = 0; i < Fwd[0].num_elements(); ++i)
{
Fwd [0][i], Fwd [1][i], 0.0,
Bwd [0][i], Bwd [1][i], 0.0,
flux[0][i], flux[1][i], hvf);
}
}
else if (expDim == 2)
{
for (int i = 0; i < Fwd[0].num_elements(); ++i)
{
Fwd [0][i], Fwd [1][i], Fwd [2][i],
Bwd [0][i], Bwd [1][i], Bwd [2][i],
flux[0][i], flux[1][i], flux[2][i]);
}
}
else if (expDim == 3)
{
ASSERTL0(false, "No 3D Shallow water supported.");
}
}
else
{
v_ArraySolve(Fwd, Bwd, flux);
}
}

Member Data Documentation

bool Nektar::NonlinearSWESolver::m_pointSolve
protected

Definition at line 48 of file NonlinearSWESolver.h.

Referenced by v_Solve().