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#include <boost/core/ignore_unused.hpp>
38
39namespace Nektar
40{
43 : RiemannSolver(pSession), m_pointSolve(true)
44{
45 m_requiresRotation = true;
46
47 // Create equation of state object
48 std::string eosType;
49 pSession->LoadSolverInfo("EquationOfState", eosType, "IdealGas");
50 m_eos = GetEquationOfStateFactory().CreateInstance(eosType, pSession);
51 // Check if using ideal gas
52 m_idealGas = boost::iequals(eosType, "IdealGas");
53}
54
56 : RiemannSolver(), m_pointSolve(true), m_idealGas(true)
57{
58 m_requiresRotation = true;
59}
60
61/**
62 *
63 */
65 const int nDim, const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
66 const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
68{
69 if (m_pointSolve)
70 {
71 size_t expDim = nDim;
72 NekDouble rhouf{}, rhovf{};
73
74 if (expDim == 1)
75 {
76 for (size_t i = 0; i < Fwd[0].size(); ++i)
77 {
78 v_PointSolve(Fwd[0][i], Fwd[1][i], 0.0, 0.0, Fwd[2][i],
79 Bwd[0][i], Bwd[1][i], 0.0, 0.0, Bwd[2][i],
80 flux[0][i], flux[1][i], rhouf, rhovf, flux[2][i]);
81 }
82 }
83 else if (expDim == 2)
84 {
85 for (size_t i = 0; i < Fwd[0].size(); ++i)
86 {
87 v_PointSolve(Fwd[0][i], Fwd[1][i], Fwd[2][i], 0.0, Fwd[3][i],
88 Bwd[0][i], Bwd[1][i], Bwd[2][i], 0.0, Bwd[3][i],
89 flux[0][i], flux[1][i], flux[2][i], rhovf,
90 flux[3][i]);
91 }
92 }
93 else if (expDim == 3)
94 {
95 for (size_t i = 0; i < Fwd[0].size(); ++i)
96 {
97 v_PointSolve(Fwd[0][i], Fwd[1][i], Fwd[2][i], Fwd[3][i],
98 Fwd[4][i], Bwd[0][i], Bwd[1][i], Bwd[2][i],
99 Bwd[3][i], Bwd[4][i], flux[0][i], flux[1][i],
100 flux[2][i], flux[3][i], flux[4][i]);
101 }
102 }
103 }
104 else
105 {
106 v_ArraySolve(Fwd, Bwd, flux);
107 }
108}
109
110/**
111 *
112 */
114 NekDouble rhoL, NekDouble pL, NekDouble eL, NekDouble HL, NekDouble srL,
115 NekDouble rhoR, NekDouble pR, NekDouble eR, NekDouble HR, NekDouble srR,
116 NekDouble HRoe, NekDouble URoe2, NekDouble srLR)
117{
118 boost::ignore_unused(HL, srL, HR, srR, srLR);
119
120 static NekDouble gamma = m_params["gamma"]();
121 NekDouble cRoe;
122 if (m_idealGas)
123 {
124 cRoe = sqrt((gamma - 1.0) * (HRoe - 0.5 * URoe2));
125 }
126 else
127 {
128 // Calculate static enthalpy of left and right states
129 NekDouble hL = eL + pL / rhoL;
130 NekDouble hR = eR + pR / rhoR;
131
132 // Get partial derivatives of P(rho,e)
133 NekDouble dpdeL = m_eos->GetDPDe_rho(rhoL, eL);
134 NekDouble dpdeR = m_eos->GetDPDe_rho(rhoR, eR);
135 NekDouble dpdrhoL = m_eos->GetDPDrho_e(rhoL, eL);
136 NekDouble dpdrhoR = m_eos->GetDPDrho_e(rhoR, eR);
137
138 // Evaluate chi and kappa parameters
139 NekDouble chiL = dpdrhoL - eL / rhoL * dpdeL;
140 NekDouble kappaL = dpdeL / rhoL;
141 NekDouble chiR = dpdrhoR - eR / rhoR * dpdeR;
142 NekDouble kappaR = dpdeR / rhoR;
143
144 //
145 // Calculate interface speed of sound using procedure from
146 // Vinokur, M.; Montagné, J.-L. "Generalized Flux-Vector
147 // Splitting and Roe Average for an Equilibrium Real Gas",
148 // JCP (1990).
149 //
150
151 // Calculate averages
152 NekDouble avgChi = 0.5 * (chiL + chiR);
153 NekDouble avgKappa = 0.5 * (kappaL + kappaR);
154 NekDouble avgKappaH = 0.5 * (kappaL * hL + kappaR * hR);
155
156 // Calculate jumps
157 NekDouble deltaP = pR - pL;
158 NekDouble deltaRho = rhoR - rhoL;
159 NekDouble deltaRhoe = rhoR * eR - rhoL * eL;
160
161 // Evaluate dP: equation (64) from Vinokur-Montagné
162 NekDouble dP = deltaP - avgChi * deltaRho - avgKappa * deltaRhoe;
163 // s (eq 66)
164 NekDouble s = avgChi + avgKappaH;
165 // D (eq 65)
166 NekDouble D = (s * deltaRho) * (s * deltaRho) + deltaP * deltaP;
167 // chiRoe and kappaRoe (eq 66)
168 NekDouble chiRoe, kappaRoe;
169 NekDouble fac = D - deltaP * deltaRho;
171 {
172 chiRoe = (D * avgChi + s * s * deltaRho * dP) / fac;
173 kappaRoe = D * avgKappa / fac;
174 }
175 else
176 {
177 chiRoe = avgChi;
178 kappaRoe = avgKappa;
179 }
180 // Speed of sound (eq 53)
181 cRoe = sqrt(chiRoe + kappaRoe * (HRoe - 0.5 * URoe2));
182 }
183 return cRoe;
184}
185
186} // 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)
virtual 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.
Definition: NekFactory.hpp:144
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
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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