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 
39 namespace 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 
62  const int nDim, const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
63  const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
65 {
66  if (m_pointSolve)
67  {
68  size_t expDim = nDim;
69  NekDouble rhouf{}, rhovf{};
70 
71  if (expDim == 1)
72  {
73  for (size_t i = 0; i < Fwd[0].size(); ++i)
74  {
75  v_PointSolve(Fwd[0][i], Fwd[1][i], 0.0, 0.0, Fwd[2][i],
76  Bwd[0][i], Bwd[1][i], 0.0, 0.0, Bwd[2][i],
77  flux[0][i], flux[1][i], rhouf, rhovf, flux[2][i]);
78  }
79  }
80  else if (expDim == 2)
81  {
82  for (size_t i = 0; i < Fwd[0].size(); ++i)
83  {
84  v_PointSolve(Fwd[0][i], Fwd[1][i], Fwd[2][i], 0.0, Fwd[3][i],
85  Bwd[0][i], Bwd[1][i], Bwd[2][i], 0.0, Bwd[3][i],
86  flux[0][i], flux[1][i], flux[2][i], rhovf,
87  flux[3][i]);
88  }
89  }
90  else if (expDim == 3)
91  {
92  for (size_t i = 0; i < Fwd[0].size(); ++i)
93  {
94  v_PointSolve(Fwd[0][i], Fwd[1][i], Fwd[2][i], Fwd[3][i],
95  Fwd[4][i], Bwd[0][i], Bwd[1][i], Bwd[2][i],
96  Bwd[3][i], Bwd[4][i], flux[0][i], flux[1][i],
97  flux[2][i], flux[3][i], flux[4][i]);
98  }
99  }
100  }
101  else
102  {
103  v_ArraySolve(Fwd, Bwd, flux);
104  }
105 }
106 
108  NekDouble rhoL, NekDouble pL, NekDouble eL, NekDouble HL, NekDouble srL,
109  NekDouble rhoR, NekDouble pR, NekDouble eR, NekDouble HR, NekDouble srR,
110  NekDouble HRoe, NekDouble URoe2, NekDouble srLR)
111 {
112  boost::ignore_unused(HL, srL, HR, srR, srLR);
113 
114  static NekDouble gamma = m_params["gamma"]();
115  NekDouble cRoe;
116  if (m_idealGas)
117  {
118  cRoe = sqrt((gamma - 1.0) * (HRoe - 0.5 * URoe2));
119  }
120  else
121  {
122  // Calculate static enthalpy of left and right states
123  NekDouble hL = eL + pL / rhoL;
124  NekDouble hR = eR + pR / rhoR;
125 
126  // Get partial derivatives of P(rho,e)
127  NekDouble dpdeL = m_eos->GetDPDe_rho(rhoL, eL);
128  NekDouble dpdeR = m_eos->GetDPDe_rho(rhoR, eR);
129  NekDouble dpdrhoL = m_eos->GetDPDrho_e(rhoL, eL);
130  NekDouble dpdrhoR = m_eos->GetDPDrho_e(rhoR, eR);
131 
132  // Evaluate chi and kappa parameters
133  NekDouble chiL = dpdrhoL - eL / rhoL * dpdeL;
134  NekDouble kappaL = dpdeL / rhoL;
135  NekDouble chiR = dpdrhoR - eR / rhoR * dpdeR;
136  NekDouble kappaR = dpdeR / rhoR;
137 
138  //
139  // Calculate interface speed of sound using procedure from
140  // Vinokur, M.; Montagné, J.-L. "Generalized Flux-Vector
141  // Splitting and Roe Average for an Equilibrium Real Gas",
142  // JCP (1990).
143  //
144 
145  // Calculate averages
146  NekDouble avgChi = 0.5 * (chiL + chiR);
147  NekDouble avgKappa = 0.5 * (kappaL + kappaR);
148  NekDouble avgKappaH = 0.5 * (kappaL * hL + kappaR * hR);
149 
150  // Calculate jumps
151  NekDouble deltaP = pR - pL;
152  NekDouble deltaRho = rhoR - rhoL;
153  NekDouble deltaRhoe = rhoR * eR - rhoL * eL;
154 
155  // Evaluate dP: equation (64) from Vinokur-Montagné
156  NekDouble dP = deltaP - avgChi * deltaRho - avgKappa * deltaRhoe;
157  // s (eq 66)
158  NekDouble s = avgChi + avgKappaH;
159  // D (eq 65)
160  NekDouble D = (s * deltaRho) * (s * deltaRho) + deltaP * deltaP;
161  // chiRoe and kappaRoe (eq 66)
162  NekDouble chiRoe, kappaRoe;
163  NekDouble fac = D - deltaP * deltaRho;
165  {
166  chiRoe = (D * avgChi + s * s * deltaRho * dP) / fac;
167  kappaRoe = D * avgKappa / fac;
168  }
169  else
170  {
171  chiRoe = avgChi;
172  kappaRoe = avgKappa;
173  }
174  // Speed of sound (eq 53)
175  cRoe = sqrt(chiRoe + kappaRoe * (HRoe - 0.5 * URoe2));
176  }
177  return cRoe;
178 }
179 
180 } // 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)
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)
CompressibleSolver()
Programmatic ctor.
EquationOfStateSharedPtr m_eos
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