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)
 
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) override
 HLL Riemann solver. More...
 
 HLLSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
void v_PointSolve (NekDouble hL, NekDouble huL, NekDouble hvL, NekDouble hR, NekDouble huR, NekDouble hvR, NekDouble &hf, NekDouble &huf, NekDouble &hvf) override
 HLL Riemann solver for the Nonlinear Shallow Water Equations. More...
 
- Protected Member Functions inherited from Nektar::CompressibleSolver
 CompressibleSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 Session ctor. More...
 
 CompressibleSolver ()
 Programmatic ctor. More...
 
void v_Solve (const int nDim, const Array< OneD, const Array< OneD, ND > > &Fwd, const Array< OneD, const Array< OneD, ND > > &Bwd, Array< OneD, Array< OneD, ND > > &flux) override
 
virtual void v_ArraySolve (const Array< OneD, const Array< OneD, ND > > &Fwd, const Array< OneD, const Array< OneD, ND > > &Bwd, Array< OneD, Array< OneD, ND > > &flux)
 
virtual void v_PointSolve (ND rhoL, ND rhouL, ND rhovL, ND rhowL, ND EL, ND rhoR, ND rhouR, ND rhovR, ND rhowR, ND ER, ND &rhof, ND &rhouf, ND &rhovf, ND &rhowf, ND &Ef)
 
ND GetRoeSoundSpeed (ND rhoL, ND pL, ND eL, ND HL, ND srL, ND rhoR, ND pR, ND eR, ND HR, ND srR, ND HRoe, ND URoe2, ND srLR)
 
- Protected Member Functions inherited from Nektar::SolverUtils::RiemannSolver
SOLVER_UTILS_EXPORT RiemannSolver ()
 
SOLVER_UTILS_EXPORT RiemannSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual SOLVER_UTILS_EXPORT ~RiemannSolver ()
 
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)=0
 
SOLVER_UTILS_EXPORT 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...
 
virtual SOLVER_UTILS_EXPORT void v_CalcFluxJacobian (const int nDim, const Array< OneD, const Array< OneD, NekDouble > > &Fwd, const Array< OneD, const Array< OneD, NekDouble > > &Bwd, const Array< OneD, const Array< OneD, NekDouble > > &normals, DNekBlkMatSharedPtr &FJac, DNekBlkMatSharedPtr &BJac)
 
- Protected Member Functions inherited from Nektar::NonlinearSWESolver
 NonlinearSWESolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
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) override
 
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)
 

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 SetALEFlag (bool &ALE)
 
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)
 
void SetAuxVec (std::string name, RSVecFuncType fp)
 
std::map< std::string, RSScalarFuncType > & GetScalars ()
 
std::map< std::string, RSVecFuncType > & GetVectors ()
 
std::map< std::string, RSParamFuncType > & GetParams ()
 
SOLVER_UTILS_EXPORT void CalcFluxJacobian (const int nDim, const Array< OneD, const Array< OneD, NekDouble > > &Fwd, const Array< OneD, const Array< OneD, NekDouble > > &Bwd, DNekBlkMatSharedPtr &FJac, DNekBlkMatSharedPtr &BJac)
 Calculate the flux jacobian of Fwd and Bwd. More...
 
- Public Attributes inherited from Nektar::SolverUtils::RiemannSolver
int m_spacedim
 
- Protected Types inherited from Nektar::CompressibleSolver
using ND = NekDouble
 
- Protected Attributes inherited from Nektar::CompressibleSolver
bool m_pointSolve
 
EquationOfStateSharedPtr m_eos
 
bool m_idealGas
 
- 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...
 
bool m_ALESolver = false
 Flag if using the ALE formulation. More...
 
- Protected Attributes inherited from Nektar::NonlinearSWESolver
bool m_pointSolve
 

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 43 of file CompressibleFlowSolver/RiemannSolvers/HLLSolver.cpp.

44 : CompressibleSolver(pSession)
45{
46}
CompressibleSolver()
Programmatic ctor.

Referenced by create().

◆ 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.

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

References HLLSolver().

◆ create() [2/2]

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

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

47 {
48 return RiemannSolverSharedPtr(new HLLSolver(pSession));
49 }

References HLLSolver().

◆ v_PointSolve() [1/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 
)
overrideprotectedvirtual

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 67 of file CompressibleFlowSolver/RiemannSolvers/HLLSolver.cpp.

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

References Nektar::CompressibleSolver::GetRoeSoundSpeed(), Nektar::CompressibleSolver::m_eos, and tinysimd::sqrt().

◆ v_PointSolve() [2/2]

void Nektar::HLLSolver::v_PointSolve ( NekDouble  hL,
NekDouble  huL,
NekDouble  hvL,
NekDouble  hR,
NekDouble  huR,
NekDouble  hvR,
NekDouble hf,
NekDouble huf,
NekDouble hvf 
)
overrideprotectedvirtual

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 61 of file ShallowWaterSolver/RiemannSolvers/HLLSolver.cpp.

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

References tinysimd::sqrt().

Member Data Documentation

◆ solverName

static std::string Nektar::HLLSolver::solverName
static
Initial value:
=
"HLL", HLLSolver::create, "HLL Riemann solver")
static RiemannSolverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
RiemannSolverFactory & GetRiemannSolverFactory()

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