Nektar++
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. More...
 
- 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. More...
 
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. More...
 
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. More...
 
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. More...
 
SOLVER_UTILS_EXPORT bool CheckScalars (std::string name)
 Determine whether a scalar has been defined in m_scalars. More...
 
SOLVER_UTILS_EXPORT bool CheckVectors (std::string name)
 Determine whether a vector has been defined in m_vectors. More...
 
SOLVER_UTILS_EXPORT bool CheckParams (std::string name)
 Determine whether a parameter has been defined in m_params. More...
 
SOLVER_UTILS_EXPORT bool CheckAuxScal (std::string name)
 Determine whether a scalar has been defined in m_auxScal. More...
 
SOLVER_UTILS_EXPORT bool CheckAuxVec (std::string name)
 Determine whether a vector has been defined in m_auxVec. More...
 

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. More...
 
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
 
- 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. More...
 
std::map< std::string, RSScalarFuncTypem_scalars
 Map of scalar function types. More...
 
std::map< std::string, RSVecFuncTypem_vectors
 Map of vector function types. More...
 
std::map< std::string, RSParamFuncTypem_params
 Map of parameter function types. More...
 
std::map< std::string, RSScalarFuncTypem_auxScal
 Map of auxiliary scalar function types. More...
 
std::map< std::string, RSVecFuncTypem_auxVec
 Map of auxiliary vector function types. More...
 
Array< OneD, Array< OneD, NekDouble > > m_rotMat
 Rotation matrices for each trace quadrature point. More...
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_rotStorage
 Rotation storage. More...
 

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

47  : LinearSWESolver()
48  {
49 
50  }

Member Function Documentation

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

Definition at line 46 of file LinearHLLSolver.h.

References LinearHLLSolver().

47  {
49  new LinearHLLSolver());
50  }
boost::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.
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.

71  {
72  static NekDouble g = m_params["gravity"]();
73 
74  if (dL != dR)
75  {
76  ASSERTL0(false,"not constant depth in LinearHLL");
77  }
78 
79  // Left and right wave speeds
80  NekDouble cL = sqrt(g * dL);
81  NekDouble cR = sqrt(g * dR);
82 
83  NekDouble SL = uL - cL;
84  NekDouble SR = uR + cR;
85 
86  if (SL >= 0)
87  {
88  etaf = dL * uL;
89  uf = g * etaL;
90  vf = 0.0;
91  }
92  else if (SR <= 0)
93  {
94  etaf = dR * uR;
95  uf = g * etaR;
96  vf = 0.0;
97  }
98  else
99  {
100  etaf = (SR * dL * uL - SL * dR * uR +
101  SL * SR * (etaR - etaL)) / (SR - SL);
102  uf = (SR * g * etaL - SL * g * etaR +
103  SL * SR * (uR - uL)) / (SR - SL);
104  vf = (SL * SR * (vR - vL)) / (SR - SL);
105  }
106  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
double NekDouble
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.

Member Data Documentation

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

Definition at line 52 of file LinearHLLSolver.h.