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

#include <LinearHLLSolver.h>

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

Static Public Member Functions

static RiemannSolverSharedPtr create ()

Static Public Attributes

static std::string solverName

Protected Member Functions

 LinearHLLSolver ()
virtual void v_PointSolve (double etaL, double uL, double vL, double dL, double etaR, double uR, double vR, double dR, double &etaf, double &uf, double &vf)
 HLL Riemann solver for the linear Shallow Water Equations.
- Protected Member Functions inherited from Nektar::LinearSWESolver
 LinearSWESolver ()
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)
- 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.
SOLVER_UTILS_EXPORT bool CheckScalars (std::string name)
 Determine whether a scalar has been defined in m_scalars.
SOLVER_UTILS_EXPORT bool CheckVectors (std::string name)
 Determine whether a vector has been defined in m_vectors.
SOLVER_UTILS_EXPORT bool CheckParams (std::string name)
 Determine whether a parameter has been defined in m_params.
SOLVER_UTILS_EXPORT bool CheckAuxScal (std::string name)
 Determine whether a scalar has been defined in m_auxScal.
SOLVER_UTILS_EXPORT bool CheckAuxVec (std::string name)
 Determine whether a vector has been defined in m_auxVec.

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
- Protected Attributes inherited from Nektar::LinearSWESolver
bool m_pointSolve

Detailed Description

Definition at line 43 of file LinearHLLSolver.h.

Constructor & Destructor Documentation

Nektar::LinearHLLSolver::LinearHLLSolver ( )
protected

Definition at line 47 of file LinearHLLSolver.cpp.

Referenced by create().

{
}

Member Function Documentation

static RiemannSolverSharedPtr Nektar::LinearHLLSolver::create ( )
inlinestatic

Definition at line 46 of file LinearHLLSolver.h.

References LinearHLLSolver().

void Nektar::LinearHLLSolver::v_PointSolve ( double  etaL,
double  uL,
double  vL,
double  dL,
double  etaR,
double  uR,
double  vR,
double  dR,
double &  etaf,
double &  uf,
double &  vf 
)
protectedvirtual

HLL Riemann solver for the linear Shallow Water Equations.

Parameters
etaLFree surface elevation left state.
etaRFree surface elevation right state.
uLx-velocity left state.
uRx-velocity right state.
vLy-velocity left state.
vRy-velocity right state.
dLstill water depth component left state.
dRstill water depth component right state.
etafComputed Riemann flux for continuity.
ufComputed Riemann flux for x-momentum component
vfComputed Riemann flux for y-momentum component

Reimplemented from Nektar::LinearSWESolver.

Definition at line 67 of file LinearHLLSolver.cpp.

References ASSERTL0, and Nektar::SolverUtils::RiemannSolver::m_params.

{
static NekDouble g = m_params["gravity"]();
if (dL != dR)
{
ASSERTL0(false,"not constant depth in LinearHLL");
}
// Left and right wave speeds
NekDouble cL = sqrt(g * dL);
NekDouble cR = sqrt(g * dR);
NekDouble SL = uL - cL;
NekDouble SR = uR + cR;
if (SL >= 0)
{
etaf = dL * uL;
uf = g * etaL;
vf = 0.0;
}
else if (SR <= 0)
{
etaf = dR * uR;
uf = g * etaR;
vf = 0.0;
}
else
{
etaf = (SR * dL * uL - SL * dR * uR +
SL * SR * (etaR - etaL)) / (SR - SL);
uf = (SR * g * etaL - SL * g * etaR +
SL * SR * (uR - uL)) / (SR - SL);
vf = (SL * SR * (vR - vL)) / (SR - SL);
}
}

Member Data Documentation

std::string Nektar::LinearHLLSolver::solverName
static
Initial value:

Definition at line 52 of file LinearHLLSolver.h.