Nektar++
CompressibleSolver.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: CompressibleSolver.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: Compressible Riemann solver.
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include "CompressibleSolver.h"
36#include <boost/algorithm/string/predicate.hpp>
37
38namespace Nektar
39{
42 : RiemannSolver(pSession), m_pointSolve(true)
43{
44 m_requiresRotation = true;
45
46 // Create equation of state object
47 std::string eosType;
48 pSession->LoadSolverInfo("EquationOfState", eosType, "IdealGas");
49 m_eos = GetEquationOfStateFactory().CreateInstance(eosType, pSession);
50 // Check if using ideal gas
51 m_idealGas = boost::iequals(eosType, "IdealGas");
52}
53
55 : RiemannSolver(), m_pointSolve(true), m_idealGas(true)
56{
57 m_requiresRotation = true;
58}
59
60/**
61 *
62 */
64 const int nDim, const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
65 const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
67{
68 if (m_pointSolve)
69 {
70 size_t expDim = nDim;
71 NekDouble rhouf{}, rhovf{};
72
73 if (expDim == 1)
74 {
75 for (size_t i = 0; i < Fwd[0].size(); ++i)
76 {
77 v_PointSolve(Fwd[0][i], Fwd[1][i], 0.0, 0.0, Fwd[2][i],
78 Bwd[0][i], Bwd[1][i], 0.0, 0.0, Bwd[2][i],
79 flux[0][i], flux[1][i], rhouf, rhovf, flux[2][i]);
80 }
81 }
82 else if (expDim == 2)
83 {
84 for (size_t i = 0; i < Fwd[0].size(); ++i)
85 {
86 v_PointSolve(Fwd[0][i], Fwd[1][i], Fwd[2][i], 0.0, Fwd[3][i],
87 Bwd[0][i], Bwd[1][i], Bwd[2][i], 0.0, Bwd[3][i],
88 flux[0][i], flux[1][i], flux[2][i], rhovf,
89 flux[3][i]);
90 }
91 }
92 else if (expDim == 3)
93 {
94 for (size_t i = 0; i < Fwd[0].size(); ++i)
95 {
96 v_PointSolve(Fwd[0][i], Fwd[1][i], Fwd[2][i], Fwd[3][i],
97 Fwd[4][i], Bwd[0][i], Bwd[1][i], Bwd[2][i],
98 Bwd[3][i], Bwd[4][i], flux[0][i], flux[1][i],
99 flux[2][i], flux[3][i], flux[4][i]);
100 }
101 }
102 }
103 else
104 {
105 v_ArraySolve(Fwd, Bwd, flux);
106 }
107}
108
109/**
110 *
111 */
113 NekDouble rhoL, NekDouble pL, NekDouble eL, [[maybe_unused]] NekDouble HL,
114 [[maybe_unused]] NekDouble srL, NekDouble rhoR, NekDouble pR, NekDouble eR,
115 [[maybe_unused]] NekDouble HR, [[maybe_unused]] NekDouble srR,
116 NekDouble HRoe, NekDouble URoe2, [[maybe_unused]] NekDouble srLR)
117{
118 static NekDouble gamma = m_params["gamma"]();
119 NekDouble cRoe;
120 if (m_idealGas)
121 {
122 cRoe = sqrt((gamma - 1.0) * (HRoe - 0.5 * URoe2));
123 }
124 else
125 {
126 // Calculate static enthalpy of left and right states
127 NekDouble hL = eL + pL / rhoL;
128 NekDouble hR = eR + pR / rhoR;
129
130 // Get partial derivatives of P(rho,e)
131 NekDouble dpdeL = m_eos->GetDPDe_rho(rhoL, eL);
132 NekDouble dpdeR = m_eos->GetDPDe_rho(rhoR, eR);
133 NekDouble dpdrhoL = m_eos->GetDPDrho_e(rhoL, eL);
134 NekDouble dpdrhoR = m_eos->GetDPDrho_e(rhoR, eR);
135
136 // Evaluate chi and kappa parameters
137 NekDouble chiL = dpdrhoL - eL / rhoL * dpdeL;
138 NekDouble kappaL = dpdeL / rhoL;
139 NekDouble chiR = dpdrhoR - eR / rhoR * dpdeR;
140 NekDouble kappaR = dpdeR / rhoR;
141
142 //
143 // Calculate interface speed of sound using procedure from
144 // Vinokur, M.; Montagné, J.-L. "Generalized Flux-Vector
145 // Splitting and Roe Average for an Equilibrium Real Gas",
146 // JCP (1990).
147 //
148
149 // Calculate averages
150 NekDouble avgChi = 0.5 * (chiL + chiR);
151 NekDouble avgKappa = 0.5 * (kappaL + kappaR);
152 NekDouble avgKappaH = 0.5 * (kappaL * hL + kappaR * hR);
153
154 // Calculate jumps
155 NekDouble deltaP = pR - pL;
156 NekDouble deltaRho = rhoR - rhoL;
157 NekDouble deltaRhoe = rhoR * eR - rhoL * eL;
158
159 // Evaluate dP: equation (64) from Vinokur-Montagné
160 NekDouble dP = deltaP - avgChi * deltaRho - avgKappa * deltaRhoe;
161 // s (eq 66)
162 NekDouble s = avgChi + avgKappaH;
163 // D (eq 65)
164 NekDouble D = (s * deltaRho) * (s * deltaRho) + deltaP * deltaP;
165 // chiRoe and kappaRoe (eq 66)
166 NekDouble chiRoe, kappaRoe;
167 NekDouble fac = D - deltaP * deltaRho;
169 {
170 chiRoe = (D * avgChi + s * s * deltaRho * dP) / fac;
171 kappaRoe = D * avgKappa / fac;
172 }
173 else
174 {
175 chiRoe = avgChi;
176 kappaRoe = avgKappa;
177 }
178 // Speed of sound (eq 53)
179 cRoe = sqrt(chiRoe + kappaRoe * (HRoe - 0.5 * URoe2));
180 }
181 return cRoe;
182}
183
184} // namespace Nektar
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)
CompressibleSolver()
Programmatic ctor.
EquationOfStateSharedPtr m_eos
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)
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
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
The RiemannSolver class provides an abstract interface under which solvers for various Riemann proble...
Definition: RiemannSolver.h:58
bool m_requiresRotation
Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields.
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
static const NekDouble kNekZeroTol
EquationOfStateFactory & GetEquationOfStateFactory()
Declaration of the equation of state factory singleton.
double NekDouble
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:298
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294