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

#include <HLLCSolver.h>

Inheritance diagram for Nektar::HLLCSolver:
[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

 HLLCSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual void v_PointSolve (NekDouble rhoL, NekDouble rhouL, NekDouble rhovL, NekDouble rhowL, NekDouble EL, NekDouble rhoR, NekDouble rhouR, NekDouble rhovR, NekDouble rhowR, NekDouble ER, NekDouble &rhof, NekDouble &rhouf, NekDouble &rhovf, NekDouble &rhowf, NekDouble &Ef)
 HLLC Riemann solver. More...
 
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)
 
 HLLCSolver (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)
 HLLC 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)
 
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/HLLCSolver.h.

Constructor & Destructor Documentation

◆ HLLCSolver() [1/2]

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

Definition at line 43 of file CompressibleFlowSolver/RiemannSolvers/HLLCSolver.cpp.

Referenced by create().

44  : CompressibleSolver(pSession)
45  {
46 
47  }
CompressibleSolver(const LibUtilities::SessionReaderSharedPtr &pSession)

◆ HLLCSolver() [2/2]

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

Member Function Documentation

◆ create() [1/2]

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

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

References HLLCSolver().

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

◆ create() [2/2]

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

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

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

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

◆ v_PointSolve() [1/2]

void Nektar::HLLCSolver::v_PointSolve ( NekDouble  rhoL,
NekDouble  rhouL,
NekDouble  rhovL,
NekDouble  rhowL,
NekDouble  EL,
NekDouble  rhoR,
NekDouble  rhouR,
NekDouble  rhovR,
NekDouble  rhowR,
NekDouble  ER,
NekDouble rhof,
NekDouble rhouf,
NekDouble rhovf,
NekDouble rhowf,
NekDouble Ef 
)
protectedvirtual

HLLC 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 68 of file CompressibleFlowSolver/RiemannSolvers/HLLCSolver.cpp.

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

Referenced by create().

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 =
83  (EL - 0.5 * (rhouL * uL + rhovL * vL + rhowL * wL)) / rhoL;
84  NekDouble eR =
85  (ER - 0.5 * (rhouR * uR + rhovR * vR + rhowR * wR)) / rhoR;
86  // Pressure
87  NekDouble pL = m_eos->GetPressure(rhoL, eL);
88  NekDouble pR = m_eos->GetPressure(rhoR, eR);
89  // Speed of sound
90  NekDouble cL = m_eos->GetSoundSpeed(rhoL, eL);
91  NekDouble cR = m_eos->GetSoundSpeed(rhoR, eR);
92 
93  // Left and right total enthalpy
94  NekDouble HL = (EL + pL) / rhoL;
95  NekDouble HR = (ER + pR) / rhoR;
96 
97  // Square root of rhoL and rhoR.
98  NekDouble srL = sqrt(rhoL);
99  NekDouble srR = sqrt(rhoR);
100  NekDouble srLR = srL + srR;
101 
102  // Roe average state
103  NekDouble uRoe = (srL * uL + srR * uR) / srLR;
104  NekDouble vRoe = (srL * vL + srR * vR) / srLR;
105  NekDouble wRoe = (srL * wL + srR * wR) / srLR;
106  NekDouble URoe2 = uRoe*uRoe + vRoe*vRoe + wRoe*wRoe;
107  NekDouble HRoe = (srL * HL + srR * HR) / srLR;
109  rhoL, pL, eL, HL, srL,
110  rhoR, pR, eR, HR, srR,
111  HRoe, URoe2, srLR);
112 
113  // Maximum wave speeds
114  NekDouble SL = std::min(uL-cL, uRoe-cRoe);
115  NekDouble SR = std::max(uR+cR, uRoe+cRoe);
116 
117  // HLLC Riemann fluxes (positive case)
118  if (SL >= 0)
119  {
120  rhof = rhouL;
121  rhouf = rhouL * uL + pL;
122  rhovf = rhouL * vL;
123  rhowf = rhouL * wL;
124  Ef = uL * (EL + pL);
125  }
126  // HLLC Riemann fluxes (negative case)
127  else if (SR <= 0)
128  {
129  rhof = rhouR;
130  rhouf = rhouR * uR + pR;
131  rhovf = rhouR * vR;
132  rhowf = rhouR * wR;
133  Ef = uR * (ER + pR);
134  }
135  // HLLC Riemann fluxes (general case (SL < 0 | SR > 0)
136  else
137  {
138  NekDouble SM = (pR - pL + rhouL * (SL - uL) - rhouR * (SR - uR)) /
139  (rhoL * (SL - uL) - rhoR * (SR - uR));
140  NekDouble rhoML = rhoL * (SL - uL) / (SL - SM);
141  NekDouble rhouML = rhoML * SM;
142  NekDouble rhovML = rhoML * vL;
143  NekDouble rhowML = rhoML * wL;
144  NekDouble EML = rhoML * (EL / rhoL +
145  (SM - uL) * (SM + pL / (rhoL * (SL - uL))));
146 
147  NekDouble rhoMR = rhoR * (SR - uR) / (SR - SM);
148  NekDouble rhouMR = rhoMR * SM;
149  NekDouble rhovMR = rhoMR * vR;
150  NekDouble rhowMR = rhoMR * wR;
151  NekDouble EMR = rhoMR * (ER / rhoR +
152  (SM - uR) * (SM + pR / (rhoR * (SR - uR))));
153 
154  if (SL < 0.0 && SM >= 0.0)
155  {
156  rhof = rhouL + SL * (rhoML - rhoL);
157  rhouf = rhouL * uL + pL + SL * (rhouML - rhouL);
158  rhovf = rhouL * vL + SL * (rhovML - rhovL);
159  rhowf = rhouL * wL + SL * (rhowML - rhowL);
160  Ef = uL * (EL + pL) + SL * (EML - EL);
161  }
162  else if(SM < 0.0 && SR > 0.0)
163  {
164  rhof = rhouR + SR * (rhoMR - rhoR);
165  rhouf = rhouR * uR + pR + SR * (rhouMR - rhouR);
166  rhovf = rhouR * vR + SR * (rhovMR - rhovR);
167  rhowf = rhouR * wR + SR * (rhowMR - rhowR);
168  Ef = uR * (ER + pR) + SR * (EMR - ER);
169  }
170  }
171  }
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

◆ v_PointSolve() [2/2]

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

HLLC 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/HLLCSolver.cpp.

References Nektar::ErrorUtil::efatal, Nektar::SolverUtils::RiemannSolver::m_params, and NEKERROR.

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 hC,huC,hvC,SL,SR,hstar,ustar,Sstar;
84  hstar = 0.5*(cL + cR) + 0.25*(uL - uR);
85  hstar *= hstar;
86  hstar *= (1.0/g);
87  ustar = 0.5*(uL + uR) + cL - cR;
88 
89 
90  // Compute 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  if (hstar > hR)
98  SR = uR + cR * sqrt(0.5*((hstar*hstar + hstar*hR)/(hR*hR)));
99  else
100  SR = uR + cR;
101 
102  if (fabs(hR*(uR-SR)-hL*(uL-SL)) <= 1.0e-10)
103  Sstar = ustar;
104  else
105  Sstar = (SL*hR*(uR-SR)-SR*hL*(uL-SL))/(hR*(uR-SR)-hL*(uL-SL));
106 
107  if (SL >= 0)
108  {
109  hf = hL * uL;
110  huf = uL * uL * hL + 0.5 * g * hL * hL;
111  hvf = hL * uL * vL;
112  }
113  else if (SR <= 0)
114  {
115  hf = hR * uR;
116  huf = uR * uR * hR + 0.5 * g * hR * hR;
117  hvf = hR * uR *vR;
118  }
119  else if ((SL < 0) && (Sstar >= 0))
120  {
121  hC = hL * ((SL - uL) / (SL - Sstar));
122  huC = hC * Sstar;
123  hvC = hC * vL;
124 
125  hf = hL*uL + SL * (hC - hL);
126  huf = (uL*uL*hL+0.5*g*hL*hL) + SL * (huC - hL*uL);
127  hvf = (uL*vL*hL) + SL * (hvC - hL*vL);
128  }
129  else if ((SR > 0) && (Sstar <= 0))
130  {
131  hC = hR * ((SR - uR) / (SR - Sstar));
132  huC = hC * Sstar;
133  hvC = hC * vR;
134 
135  hf = hR*uR + SR * (hC - hR);
136  huf = (uR*uR*hR+0.5*g*hR*hR) + SR * (huC - hR*uR);
137  hvf = (uR*vR*hR) + SR * (hvC - hR*vR);
138  }
139  else
140  {
142  "Error in HLLC solver -- non physical combination of "
143  "SR, SL and Sstar");
144  }
145  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:209
double NekDouble
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.

◆ v_PointSolveVisc()

void Nektar::HLLCSolver::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 
)
protectedvirtual

Reimplemented from Nektar::CompressibleSolver.

Definition at line 173 of file CompressibleFlowSolver/RiemannSolvers/HLLCSolver.cpp.

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

177  {
178  // Left and Right velocities
179  NekDouble uL = rhouL / rhoL;
180  NekDouble vL = rhovL / rhoL;
181  NekDouble wL = rhowL / rhoL;
182  NekDouble uR = rhouR / rhoR;
183  NekDouble vR = rhovR / rhoR;
184  NekDouble wR = rhowR / rhoR;
185 
186  // Internal energy (per unit mass)
187  NekDouble eL =
188  (EL - 0.5 * (rhouL * uL + rhovL * vL + rhowL * wL)) / rhoL;
189  NekDouble eR =
190  (ER - 0.5 * (rhouR * uR + rhovR * vR + rhowR * wR)) / rhoR;
191  // Pressure
192  NekDouble pL = m_eos->GetPressure(rhoL, eL);
193  NekDouble pR = m_eos->GetPressure(rhoR, eR);
194  // Speed of sound
195  NekDouble cL = m_eos->GetSoundSpeed(rhoL, eL);
196  NekDouble cR = m_eos->GetSoundSpeed(rhoR, eR);
197 
198  // Left and right total enthalpy
199  NekDouble HL = (EL + pL) / rhoL;
200  NekDouble HR = (ER + pR) / rhoR;
201 
202  // Square root of rhoL and rhoR.
203  NekDouble srL = sqrt(rhoL);
204  NekDouble srR = sqrt(rhoR);
205  NekDouble srLR = srL + srR;
206 
207  // Roe average state
208  NekDouble uRoe = (srL * uL + srR * uR) / srLR;
209  NekDouble vRoe = (srL * vL + srR * vR) / srLR;
210  NekDouble wRoe = (srL * wL + srR * wR) / srLR;
211  NekDouble URoe2 = uRoe*uRoe + vRoe*vRoe + wRoe*wRoe;
212  NekDouble HRoe = (srL * HL + srR * HR) / srLR;
214  rhoL, pL, eL, HL, srL,
215  rhoR, pR, eR, HR, srR,
216  HRoe, URoe2, srLR);
217 
218  // Maximum wave speeds
219  NekDouble SL = std::min(uL-cL, uRoe-cRoe);
220  NekDouble SR = std::max(uR+cR, uRoe+cRoe);
221 
222  // HLLC Riemann fluxes (positive case)
223  if (SL >= 0)
224  {
225  rhof = rhouL;
226  rhouf = rhouL * uL + pL;
227  rhovf = rhouL * vL;
228  rhowf = rhouL * wL;
229  Ef = uL * (EL + pL);
230  Epsf = 0.0;
231  }
232  // HLLC Riemann fluxes (negative case)
233  else if (SR <= 0)
234  {
235  rhof = rhouR;
236  rhouf = rhouR * uR + pR;
237  rhovf = rhouR * vR;
238  rhowf = rhouR * wR;
239  Ef = uR * (ER + pR);
240  Epsf = 0.0;
241  }
242  // HLLC Riemann fluxes (general case (SL < 0 | SR > 0)
243  else
244  {
245  NekDouble SM = (pR - pL + rhouL * (SL - uL) - rhouR * (SR - uR)) /
246  (rhoL * (SL - uL) - rhoR * (SR - uR));
247  NekDouble rhoML = rhoL * (SL - uL) / (SL - SM);
248  NekDouble rhouML = rhoML * SM;
249  NekDouble rhovML = rhoML * vL;
250  NekDouble rhowML = rhoML * wL;
251  NekDouble EML = rhoML * (EL / rhoL +
252  (SM - uL) * (SM + pL / (rhoL * (SL - uL))));
253  NekDouble EpsML = EpsL * (SL - uL) / (SL - SM);
254 
255  NekDouble rhoMR = rhoR * (SR - uR) / (SR - SM);
256  NekDouble rhouMR = rhoMR * SM;
257  NekDouble rhovMR = rhoMR * vR;
258  NekDouble rhowMR = rhoMR * wR;
259  NekDouble EMR = rhoMR * (ER / rhoR +
260  (SM - uR) * (SM + pR / (rhoR * (SR - uR))));
261  NekDouble EpsMR = EpsR * (SL - uR) / (SL - SM);
262 
263  if (SL < 0.0 && SM >= 0.0)
264  {
265  rhof = rhouL + SL * (rhoML - rhoL);
266  rhouf = rhouL * uL + pL + SL * (rhouML - rhouL);
267  rhovf = rhouL * vL + SL * (rhovML - rhovL);
268  rhowf = rhouL * wL + SL * (rhowML - rhowL);
269  Ef = uL * (EL + pL) + SL * (EML - EL);
270  Epsf = 0.0 + SL * (EpsML - EpsL);
271  }
272  else if(SM < 0.0 && SR > 0.0)
273  {
274  rhof = rhouR + SR * (rhoMR - rhoR);
275  rhouf = rhouR * uR + pR + SR * (rhouMR - rhouR);
276  rhovf = rhouR * vR + SR * (rhovMR - rhovR);
277  rhowf = rhouR * wR + SR * (rhowMR - rhowR);
278  Ef = uR * (ER + pR) + SR * (EMR - ER);
279  Epsf = 0.0 + SR * (EpsMR - EpsR);
280  }
281  }
282  }
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::HLLCSolver::solverName
static
Initial value:

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

Referenced by create().