Nektar++
Static Public Member Functions | Static Public Attributes | Protected Member Functions | List of all members
Nektar::HLLSolver Class Reference

#include <HLLSolver.h>

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

Static Public Member Functions

static RiemannSolverSharedPtr create ()
 
static RiemannSolverSharedPtr create ()
 

Static Public Attributes

static std::string solverName
 

Protected Member Functions

 HLLSolver ()
 
virtual void v_PointSolve (double rhoL, double rhouL, double rhovL, double rhowL, double EL, double rhoR, double rhouR, double rhovR, double rhowR, double ER, double &rhof, double &rhouf, double &rhovf, double &rhowf, double &Ef)
 HLL Riemann solver. More...
 
 HLLSolver ()
 
virtual void v_PointSolve (double hL, double huL, double hvL, double hR, double huR, double hvR, double &hf, double &huf, double &hvf)
 HLL Riemann solver for the Nonlinear Shallow Water Equations. More...
 
- Protected Member Functions inherited from Nektar::NonlinearSWESolver
 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)
 
- 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...
 
- Protected Member Functions inherited from Nektar::CompressibleSolver
 CompressibleSolver ()
 
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_PointSolveVisc (NekDouble rhoL, NekDouble rhouL, NekDouble rhovL, NekDouble rhowL, NekDouble EL, NekDouble EpsL, NekDouble rhoR, NekDouble rhouR, NekDouble rhovR, NekDouble rhowR, NekDouble ER, NekDouble EpsR, NekDouble &rhof, NekDouble &rhouf, NekDouble &rhovf, NekDouble &rhowf, NekDouble &Ef, NekDouble &Epsf)
 

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

Detailed Description

Definition at line 43 of file CompressibleFlowSolver/RiemannSolvers/HLLSolver.h.

Constructor & Destructor Documentation

Nektar::HLLSolver::HLLSolver ( )
protected

Definition at line 46 of file CompressibleFlowSolver/RiemannSolvers/HLLSolver.cpp.

Referenced by create().

47  {
48 
49  }
Nektar::HLLSolver::HLLSolver ( )
protected

Member Function Documentation

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

Definition at line 46 of file CompressibleFlowSolver/RiemannSolvers/HLLSolver.h.

References HLLSolver().

47  {
49  new HLLSolver());
50  }
boost::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.
static RiemannSolverSharedPtr Nektar::HLLSolver::create ( )
inlinestatic

Definition at line 46 of file ShallowWaterSolver/RiemannSolvers/HLLSolver.h.

References HLLSolver().

47  {
49  new HLLSolver());
50  }
boost::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.
void Nektar::HLLSolver::v_PointSolve ( double  hL,
double  huL,
double  hvL,
double  hR,
double  huR,
double  hvR,
double &  hf,
double &  huf,
double &  hvf 
)
protectedvirtual

HLL Riemann solver for the Nonlinear Shallow Water Equations.

Parameters
hLWater depth left state.
hRWater depth right state.
huLx-momentum component left state.
huRx-momentum component right state.
hvLy-momentum component left state.
hvRy-momentum component right state.
hfComputed Riemann flux for density.
hufComputed Riemann flux for x-momentum component
hvfComputed Riemann flux for y-momentum component

Reimplemented from Nektar::NonlinearSWESolver.

Definition at line 64 of file ShallowWaterSolver/RiemannSolvers/HLLSolver.cpp.

References Nektar::SolverUtils::RiemannSolver::m_params.

68  {
69  static NekDouble g = m_params["gravity"]();
70 
71  // Left and Right velocities
72  NekDouble uL = huL / hL;
73  NekDouble vL = hvL / hL;
74  NekDouble uR = huR / hR;
75  NekDouble vR = hvR / hR;
76 
77 
78  // Left and right wave speeds
79  NekDouble cL = sqrt(g * hL);
80  NekDouble cR = sqrt(g * hR);
81 
82  // the two-rarefaction wave assumption
83  NekDouble hstar,fL,fR;
84  hstar = 0.5 * (cL + cR) + 0.25 * (uL - uR);
85  hstar *= hstar;
86  hstar *= (1.0/g);
87 
88 
89  // Compute SL
90  NekDouble SL;
91  if (hstar > hL)
92  SL = uL - cL * sqrt(0.5*((hstar*hstar + hstar*hL)/(hL*hL)));
93  else
94  SL = uL - cL;
95 
96  // Compute SR
97  NekDouble SR;
98  if (hstar > hR)
99  SR = uR + cR * sqrt(0.5*((hstar*hstar + hstar*hR)/(hR*hR)));
100  else
101  SR = uR + cR;
102 
103  if (SL >= 0)
104  {
105  hf = hL * uL;
106  huf = uL * uL * hL + 0.5 * g * hL * hL;
107  hvf = hL * uL * vL;
108  }
109  else if (SR <= 0)
110  {
111  hf = hR * uR;
112  huf = uR * uR * hR + 0.5 * g * hR * hR;
113  hvf = hR * uR *vR;
114  }
115  else
116  {
117  hf = (SR * hL * uL - SL * hR * uR +
118  SL * SR * (hR - hL)) / (SR - SL);
119  fL = uL * uL * hL + 0.5 * g * hL * hL;
120  fR = uR * uR * hR + 0.5 * g * hR * hR;
121  huf =(SR * fL - SL * fR +
122  SL * SR * (hR * uR - hL * uL)) / (SR - SL);
123  fL = uL * vL * hL;
124  fR = uR * vR * hR;
125  hvf =(SR * fL - SL * fR +
126  SL * SR * (hR * vR - hL * vL)) / (SR - SL);
127  }
128  }
double NekDouble
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
void Nektar::HLLSolver::v_PointSolve ( double  rhoL,
double  rhouL,
double  rhovL,
double  rhowL,
double  EL,
double  rhoR,
double  rhouR,
double  rhovR,
double  rhowR,
double  ER,
double &  rhof,
double &  rhouf,
double &  rhovf,
double &  rhowf,
double &  Ef 
)
protectedvirtual

HLL Riemann solver.

Parameters
rhoLDensity left state.
rhoRDensity right state.
rhouLx-momentum component left state.
rhouRx-momentum component right state.
rhovLy-momentum component left state.
rhovRy-momentum component right state.
rhowLz-momentum component left state.
rhowRz-momentum component right state.
ELEnergy left state.
EREnergy right state.
rhofComputed Riemann flux for density.
rhoufComputed Riemann flux for x-momentum component
rhovfComputed Riemann flux for y-momentum component
rhowfComputed Riemann flux for z-momentum component
EfComputed Riemann flux for energy.

Reimplemented from Nektar::CompressibleSolver.

Definition at line 70 of file CompressibleFlowSolver/RiemannSolvers/HLLSolver.cpp.

References Nektar::SolverUtils::RiemannSolver::m_params.

74  {
75  static NekDouble gamma = m_params["gamma"]();
76 
77  // Left and Right velocities
78  NekDouble uL = rhouL / rhoL;
79  NekDouble vL = rhovL / rhoL;
80  NekDouble wL = rhowL / rhoL;
81  NekDouble uR = rhouR / rhoR;
82  NekDouble vR = rhovR / rhoR;
83  NekDouble wR = rhowR / rhoR;
84 
85  // Left and right pressures
86  NekDouble pL = (gamma - 1.0) *
87  (EL - 0.5 * (rhouL * uL + rhovL * vL + rhowL * wL));
88  NekDouble pR = (gamma - 1.0) *
89  (ER - 0.5 * (rhouR * uR + rhovR * vR + rhowR * wR));
90  NekDouble cL = sqrt(gamma * pL / rhoL);
91  NekDouble cR = sqrt(gamma * pR / rhoR);
92  NekDouble hL = (EL + pL) / rhoL;
93  NekDouble hR = (ER + pR) / rhoR;
94 
95  // Square root of rhoL and rhoR.
96  NekDouble srL = sqrt(rhoL);
97  NekDouble srR = sqrt(rhoR);
98  NekDouble srLR = srL + srR;
99 
100  // Velocity Roe averages
101  NekDouble uRoe = (srL * uL + srR * uR) / srLR;
102  NekDouble vRoe = (srL * vL + srR * vR) / srLR;
103  NekDouble wRoe = (srL * wL + srR * wR) / srLR;
104  NekDouble hRoe = (srL * hL + srR * hR) / srLR;
105  NekDouble cRoe = sqrt((gamma - 1.0)*(hRoe - 0.5 *
106  (uRoe * uRoe + vRoe * vRoe + wRoe * wRoe)));
107 
108  // Maximum wave speeds
109  NekDouble SL = std::min(uL-cL, uRoe-cRoe);
110  NekDouble SR = std::max(uR+cR, uRoe+cRoe);
111 
112  // HLL Riemann fluxes (positive case)
113  if (SL >= 0)
114  {
115  rhof = rhouL;
116  rhouf = rhouL * uL + pL;
117  rhovf = rhouL * vL;
118  rhowf = rhouL * wL;
119  Ef = uL * (EL + pL);
120  }
121  // HLL Riemann fluxes (negative case)
122  else if (SR <= 0)
123  {
124  rhof = rhouR;
125  rhouf = rhouR * uR + pR;
126  rhovf = rhouR * vR;
127  rhowf = rhouR * wR;
128  Ef = uR * (ER + pR);
129  }
130  // HLL Riemann fluxes (general case (SL < 0 | SR > 0)
131  else
132  {
133  NekDouble tmp1 = 1.0 / (SR - SL);
134  NekDouble tmp2 = SR * SL;
135  rhof = (SR * rhouL - SL * rhouR + tmp2 * (rhoR - rhoL)) * tmp1;
136  rhouf = (SR * (rhouL * uL + pL) - SL * (rhouR * uR + pR) +
137  tmp2 * (rhouR - rhouL)) * tmp1;
138  rhovf = (SR * rhouL * vL - SL * rhouR * vR +
139  tmp2 * (rhovR - rhovL)) * tmp1;
140  rhowf = (SR * rhouL * wL - SL * rhouR * wR +
141  tmp2 * (rhowR - rhowL)) * tmp1;
142  Ef = (SR * uL * (EL + pL) - SL * uR * (ER + pR) +
143  tmp2 * (ER - EL)) * tmp1;
144  }
145  }
double NekDouble
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.

Member Data Documentation

static std::string Nektar::HLLSolver::solverName
static