Nektar++
CompressibleFlowSolver/RiemannSolvers/HLLCSolver.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: HLLCSolver.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: HLLC Riemann solver.
32//
33///////////////////////////////////////////////////////////////////////////////
34
36
37namespace Nektar
38{
39std::string HLLCSolver::solverName =
41 "HLLC", HLLCSolver::create, "HLLC Riemann solver");
42
44 : CompressibleSolver(pSession)
45{
46}
47
48/**
49 * @brief HLLC Riemann solver
50 *
51 * @param rhoL Density left state.
52 * @param rhoR Density right state.
53 * @param rhouL x-momentum component left state.
54 * @param rhouR x-momentum component right state.
55 * @param rhovL y-momentum component left state.
56 * @param rhovR y-momentum component right state.
57 * @param rhowL z-momentum component left state.
58 * @param rhowR z-momentum component right state.
59 * @param EL Energy left state.
60 * @param ER Energy right state.
61 * @param rhof Computed Riemann flux for density.
62 * @param rhouf Computed Riemann flux for x-momentum component
63 * @param rhovf Computed Riemann flux for y-momentum component
64 * @param rhowf Computed Riemann flux for z-momentum component
65 * @param Ef Computed Riemann flux for energy.
66 */
68 NekDouble rhowL, NekDouble EL, NekDouble rhoR,
69 NekDouble rhouR, NekDouble rhovR, NekDouble rhowR,
70 NekDouble ER, NekDouble &rhof, NekDouble &rhouf,
71 NekDouble &rhovf, NekDouble &rhowf, NekDouble &Ef)
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 // HLLC 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 // HLLC 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 // HLLC Riemann fluxes (general case (SL < 0 | SR > 0)
132 else
133 {
134 NekDouble SM = (pR - pL + rhouL * (SL - uL) - rhouR * (SR - uR)) /
135 (rhoL * (SL - uL) - rhoR * (SR - uR));
136 NekDouble rhoML = rhoL * (SL - uL) / (SL - SM);
137 NekDouble rhouML = rhoML * SM;
138 NekDouble rhovML = rhoML * vL;
139 NekDouble rhowML = rhoML * wL;
140 NekDouble EML =
141 rhoML * (EL / rhoL + (SM - uL) * (SM + pL / (rhoL * (SL - uL))));
142
143 NekDouble rhoMR = rhoR * (SR - uR) / (SR - SM);
144 NekDouble rhouMR = rhoMR * SM;
145 NekDouble rhovMR = rhoMR * vR;
146 NekDouble rhowMR = rhoMR * wR;
147 NekDouble EMR =
148 rhoMR * (ER / rhoR + (SM - uR) * (SM + pR / (rhoR * (SR - uR))));
149
150 if (SL < 0.0 && SM >= 0.0)
151 {
152 rhof = rhouL + SL * (rhoML - rhoL);
153 rhouf = rhouL * uL + pL + SL * (rhouML - rhouL);
154 rhovf = rhouL * vL + SL * (rhovML - rhovL);
155 rhowf = rhouL * wL + SL * (rhowML - rhowL);
156 Ef = uL * (EL + pL) + SL * (EML - EL);
157 }
158 else if (SM < 0.0 && SR > 0.0)
159 {
160 rhof = rhouR + SR * (rhoMR - rhoR);
161 rhouf = rhouR * uR + pR + SR * (rhouMR - rhouR);
162 rhovf = rhouR * vR + SR * (rhovMR - rhovR);
163 rhowf = rhouR * wR + SR * (rhowMR - rhowR);
164 Ef = uR * (ER + pR) + SR * (EMR - ER);
165 }
166 }
167}
168} // namespace Nektar
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
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) override
HLLC Riemann solver.
static RiemannSolverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
HLLCSolver(const LibUtilities::SessionReaderSharedPtr &pSession)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
std::shared_ptr< SessionReader > SessionReaderSharedPtr
RiemannSolverFactory & GetRiemannSolverFactory()
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294