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:
[legend]

Static Public Member Functions

static RiemannSolverSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession)
 
static RiemannSolverSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession)
 

Static Public Attributes

static std::string solverName
 

Protected Member Functions

 HLLSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
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 (const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual void v_PointSolve (NekDouble hL, NekDouble huL, NekDouble hvL, NekDouble hR, NekDouble huR, NekDouble hvR, NekDouble &hf, NekDouble &huf, NekDouble &hvf)
 HLL Riemann solver for the Nonlinear Shallow Water Equations. More...
 
- Protected Member Functions inherited from Nektar::NonlinearSWESolver
 NonlinearSWESolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
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 (const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual 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 (const LibUtilities::SessionReaderSharedPtr &pSession)
 
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)
 
NekDouble GetRoeSoundSpeed (NekDouble rhoL, NekDouble pL, NekDouble eL, NekDouble HL, NekDouble srL, NekDouble rhoR, NekDouble pR, NekDouble eR, NekDouble HR, NekDouble srR, NekDouble HRoe, NekDouble URoe2, NekDouble srLR)
 

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
 
EquationOfStateSharedPtr m_eos
 
bool m_idealGas
 

Detailed Description

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

Constructor & Destructor Documentation

◆ HLLSolver() [1/2]

Nektar::HLLSolver::HLLSolver ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

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

Referenced by create().

46  : CompressibleSolver(pSession)
47  {
48 
49  }
CompressibleSolver(const LibUtilities::SessionReaderSharedPtr &pSession)

◆ HLLSolver() [2/2]

Nektar::HLLSolver::HLLSolver ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Member Function Documentation

◆ create() [1/2]

static RiemannSolverSharedPtr Nektar::HLLSolver::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

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

References HLLSolver().

47  {
49  new HLLSolver(pSession));
50  }
HLLSolver(const LibUtilities::SessionReaderSharedPtr &pSession)
std::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.

◆ create() [2/2]

static RiemannSolverSharedPtr Nektar::HLLSolver::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

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

References HLLSolver(), solverName, and v_PointSolve().

47  {
49  new HLLSolver(pSession));
50  }
HLLSolver(const LibUtilities::SessionReaderSharedPtr &pSession)
std::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.

◆ v_PointSolve() [1/2]

void Nektar::HLLSolver::v_PointSolve ( NekDouble  hL,
NekDouble  huL,
NekDouble  hvL,
NekDouble  hR,
NekDouble  huR,
NekDouble  hvR,
NekDouble hf,
NekDouble huf,
NekDouble 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 65 of file ShallowWaterSolver/RiemannSolvers/HLLSolver.cpp.

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

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

◆ v_PointSolve() [2/2]

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::CompressibleSolver::GetRoeSoundSpeed(), and Nektar::CompressibleSolver::m_eos.

Referenced by create().

74  {
75  // Left and Right velocities
76  NekDouble uL = rhouL / rhoL;
77  NekDouble vL = rhovL / rhoL;
78  NekDouble wL = rhowL / rhoL;
79  NekDouble uR = rhouR / rhoR;
80  NekDouble vR = rhovR / rhoR;
81  NekDouble wR = rhowR / rhoR;
82 
83  // Internal energy (per unit mass)
84  NekDouble eL =
85  (EL - 0.5 * (rhouL * uL + rhovL * vL + rhowL * wL)) / rhoL;
86  NekDouble eR =
87  (ER - 0.5 * (rhouR * uR + rhovR * vR + rhowR * wR)) / rhoR;
88  // Pressure
89  NekDouble pL = m_eos->GetPressure(rhoL, eL);
90  NekDouble pR = m_eos->GetPressure(rhoR, eR);
91  // Speed of sound
92  NekDouble cL = m_eos->GetSoundSpeed(rhoL, eL);
93  NekDouble cR = m_eos->GetSoundSpeed(rhoR, eR);
94 
95  // Left and right total enthalpy
96  NekDouble HL = (EL + pL) / rhoL;
97  NekDouble HR = (ER + pR) / rhoR;
98 
99  // Square root of rhoL and rhoR.
100  NekDouble srL = sqrt(rhoL);
101  NekDouble srR = sqrt(rhoR);
102  NekDouble srLR = srL + srR;
103 
104  // Roe average state
105  NekDouble uRoe = (srL * uL + srR * uR) / srLR;
106  NekDouble vRoe = (srL * vL + srR * vR) / srLR;
107  NekDouble wRoe = (srL * wL + srR * wR) / srLR;
108  NekDouble URoe2 = uRoe*uRoe + vRoe*vRoe + wRoe*wRoe;
109  NekDouble HRoe = (srL * HL + srR * HR) / srLR;
111  rhoL, pL, eL, HL, srL,
112  rhoR, pR, eR, HR, srR,
113  HRoe, URoe2, srLR);
114 
115  // Maximum wave speeds
116  NekDouble SL = std::min(uL-cL, uRoe-cRoe);
117  NekDouble SR = std::max(uR+cR, uRoe+cRoe);
118 
119  // HLL Riemann fluxes (positive case)
120  if (SL >= 0)
121  {
122  rhof = rhouL;
123  rhouf = rhouL * uL + pL;
124  rhovf = rhouL * vL;
125  rhowf = rhouL * wL;
126  Ef = uL * (EL + pL);
127  }
128  // HLL Riemann fluxes (negative case)
129  else if (SR <= 0)
130  {
131  rhof = rhouR;
132  rhouf = rhouR * uR + pR;
133  rhovf = rhouR * vR;
134  rhowf = rhouR * wR;
135  Ef = uR * (ER + pR);
136  }
137  // HLL Riemann fluxes (general case (SL < 0 | SR > 0)
138  else
139  {
140  NekDouble tmp1 = 1.0 / (SR - SL);
141  NekDouble tmp2 = SR * SL;
142  rhof = (SR * rhouL - SL * rhouR + tmp2 * (rhoR - rhoL)) * tmp1;
143  rhouf = (SR * (rhouL * uL + pL) - SL * (rhouR * uR + pR) +
144  tmp2 * (rhouR - rhouL)) * tmp1;
145  rhovf = (SR * rhouL * vL - SL * rhouR * vR +
146  tmp2 * (rhovR - rhovL)) * tmp1;
147  rhowf = (SR * rhouL * wL - SL * rhouR * wR +
148  tmp2 * (rhowR - rhowL)) * tmp1;
149  Ef = (SR * uL * (EL + pL) - SL * uR * (ER + pR) +
150  tmp2 * (ER - EL)) * tmp1;
151  }
152  }
double NekDouble
NekDouble GetRoeSoundSpeed(NekDouble rhoL, NekDouble pL, NekDouble eL, NekDouble HL, NekDouble srL, NekDouble rhoR, NekDouble pR, NekDouble eR, NekDouble HR, NekDouble srR, NekDouble HRoe, NekDouble URoe2, NekDouble srLR)
EquationOfStateSharedPtr m_eos

Member Data Documentation

◆ solverName

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

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

Referenced by create().