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 
36 #include <boost/core/ignore_unused.hpp>
37 #include <boost/algorithm/string/predicate.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",
50  eosType, "IdealGas");
52  .CreateInstance(eosType, pSession);
53  // Check if using ideal gas
54  m_idealGas = boost::iequals(eosType,"IdealGas");
55  }
56 
58  const int nDim,
59  const Array<OneD, const Array<OneD, NekDouble> > &Fwd,
60  const Array<OneD, const Array<OneD, NekDouble> > &Bwd,
62  {
63  if (m_pointSolve)
64  {
65  int expDim = nDim;
66  int nvariables = Fwd.num_elements();
67 
68  NekDouble rhouf, rhovf;
69 
70  // Check if PDE-based SC is used
71  if (expDim == 1)
72  {
73  for (int i = 0; i < Fwd[0].num_elements(); ++i)
74  {
76  Fwd [0][i], Fwd [1][i], 0.0, 0.0, Fwd [2][i],
77  Bwd [0][i], Bwd [1][i], 0.0, 0.0, Bwd [2][i],
78  flux[0][i], flux[1][i], rhouf, rhovf, flux[2][i]);
79  }
80  }
81  else if (expDim == 2)
82  {
83  if (nvariables == expDim+2)
84  {
85  for (int i = 0; i < Fwd[0].num_elements(); ++i)
86  {
88  Fwd [0][i], Fwd [1][i], Fwd [2][i], 0.0, Fwd [3][i],
89  Bwd [0][i], Bwd [1][i], Bwd [2][i], 0.0, Bwd [3][i],
90  flux[0][i], flux[1][i], flux[2][i], rhovf, flux[3][i]);
91  }
92  }
93 
94  if (nvariables > expDim+2)
95  {
96  for (int i = 0; i < Fwd[0].num_elements(); ++i)
97  {
99  Fwd [0][i], Fwd [1][i], Fwd [2][i], 0.0, Fwd [3][i], Fwd [4][i],
100  Bwd [0][i], Bwd [1][i], Bwd [2][i], 0.0, Bwd [3][i], Bwd [4][i],
101  flux[0][i], flux[1][i], flux[2][i], rhovf, flux[3][i], flux[4][i]);
102  }
103  }
104 
105  }
106  else if (expDim == 3)
107  {
108  for (int i = 0; i < Fwd[0].num_elements(); ++i)
109  {
110  v_PointSolve(
111  Fwd [0][i], Fwd [1][i], Fwd [2][i], Fwd [3][i], Fwd [4][i],
112  Bwd [0][i], Bwd [1][i], Bwd [2][i], Bwd [3][i], Bwd [4][i],
113  flux[0][i], flux[1][i], flux[2][i], flux[3][i], flux[4][i]);
114  }
115  if (nvariables > expDim+2)
116  {
117  for (int i = 0; i < Fwd[0].num_elements(); ++i)
118  {
120  Fwd [0][i], Fwd [1][i], Fwd [2][i], Fwd [3][i], Fwd [4][i], Fwd [5][i],
121  Bwd [0][i], Bwd [1][i], Bwd [2][i], Bwd [3][i], Bwd [4][i], Bwd [5][i],
122  flux[0][i], flux[1][i], flux[2][i], flux[3][i], flux[4][i], flux[5][i]);
123  }
124  }
125  }
126  }
127  else
128  {
129  v_ArraySolve(Fwd, Bwd, flux);
130  }
131  }
132 
134  NekDouble rhoL, NekDouble pL, NekDouble eL, NekDouble HL, NekDouble srL,
135  NekDouble rhoR, NekDouble pR, NekDouble eR, NekDouble HR, NekDouble srR,
136  NekDouble HRoe, NekDouble URoe2, NekDouble srLR)
137  {
138  boost::ignore_unused(HL, srL, HR, srR, srLR);
139 
140  static NekDouble gamma = m_params["gamma"]();
141  NekDouble cRoe;
142  if(m_idealGas)
143  {
144  cRoe = sqrt((gamma - 1.0)*(HRoe - 0.5 * URoe2));
145  }
146  else
147  {
148  // Calculate static enthalpy of left and right states
149  NekDouble hL = eL + pL/rhoL;
150  NekDouble hR = eR + pR/rhoR;
151 
152  // Get partial derivatives of P(rho,e)
153  NekDouble dpdeL = m_eos->GetDPDe_rho(rhoL,eL);
154  NekDouble dpdeR = m_eos->GetDPDe_rho(rhoR,eR);
155  NekDouble dpdrhoL = m_eos->GetDPDrho_e(rhoL,eL);
156  NekDouble dpdrhoR = m_eos->GetDPDrho_e(rhoR,eR);
157 
158  // Evaluate chi and kappa parameters
159  NekDouble chiL = dpdrhoL - eL / rhoL * dpdeL;
160  NekDouble kappaL = dpdeL / rhoL;
161  NekDouble chiR = dpdrhoR - eR / rhoR * dpdeR;
162  NekDouble kappaR = dpdeR / rhoR;
163 
164  //
165  // Calculate interface speed of sound using procedure from
166  // Vinokur, M.; Montagné, J.-L. "Generalized Flux-Vector
167  // Splitting and Roe Average for an Equilibrium Real Gas",
168  // JCP (1990).
169  //
170 
171  // Calculate averages
172  NekDouble avgChi = 0.5 * (chiL + chiR);
173  NekDouble avgKappa = 0.5 * (kappaL + kappaR);
174  NekDouble avgKappaH = 0.5 * (kappaL*hL + kappaR*hR);
175 
176  // Calculate jumps
177  NekDouble deltaP = pR - pL;
178  NekDouble deltaRho = rhoR - rhoL;
179  NekDouble deltaRhoe = rhoR*eR - rhoL*eL;
180 
181  // Evaluate dP: equation (64) from Vinokur-Montagné
182  NekDouble dP = deltaP - avgChi * deltaRho - avgKappa * deltaRhoe;
183  // s (eq 66)
184  NekDouble s = avgChi + avgKappaH;
185  // D (eq 65)
186  NekDouble D = (s*deltaRho)*(s*deltaRho) + deltaP*deltaP;
187  // chiRoe and kappaRoe (eq 66)
188  NekDouble chiRoe, kappaRoe;
189  NekDouble fac = D - deltaP*deltaRho;
190  if( std::abs(fac) > NekConstants::kNekZeroTol)
191  {
192  chiRoe = (D*avgChi + s*s*deltaRho*dP) / fac;
193  kappaRoe = D*avgKappa / fac;
194  }
195  else
196  {
197  chiRoe = avgChi;
198  kappaRoe = avgKappa;
199  }
200  // Speed of sound (eq 53)
201  cRoe = sqrt( chiRoe + kappaRoe*(HRoe - 0.5 * URoe2));
202  }
203  return cRoe;
204  }
205 
206 }
CompressibleSolver(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)
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
static const NekDouble kNekZeroTol
EquationOfStateFactory & GetEquationOfStateFactory()
Declaration of the equation of state factory singleton.
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)
double NekDouble
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)
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)
bool m_requiresRotation
Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields...
The RiemannSolver class provides an abstract interface under which solvers for various Riemann proble...
Definition: RiemannSolver.h:59
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
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
std::shared_ptr< SessionReader > SessionReaderSharedPtr