Nektar++
RedlichKwongEoS.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: RedlichKwongEoS.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: Redlich-Kwong equation of state
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include "RedlichKwongEoS.h"
36 
37 using namespace std;
38 
39 namespace Nektar
40 {
41 
42 std::string RedlichKwongEoS::className =
44  "RedlichKwong", RedlichKwongEoS::create,
45  "Redlich-Kwong equation of state.");
46 
47 RedlichKwongEoS::RedlichKwongEoS(
49  : EquationOfState(pSession)
50 {
51  pSession->LoadParameter("Tcrit", m_Tc);
52  pSession->LoadParameter("Pcrit", m_Pc);
53 
54  m_a = 0.42748 * m_gasConstant * m_gasConstant * m_Tc * m_Tc / m_Pc;
55  m_b = 0.08664 * m_gasConstant * m_Tc / m_Pc;
56 }
57 
59  const NekDouble &e)
60 {
61  return GetTemperatureKernel(rho, e);
62 }
63 
65 {
66  return GetTemperatureKernel(rho, e);
67 }
68 
70  const NekDouble &e)
71 {
72  return GetPressureKernel(rho, e);
73 }
74 
76 {
77  return GetPressureKernel(rho, e);
78 }
79 
81  const NekDouble &e)
82 {
83  NekDouble T = GetTemperature(rho, e);
84  NekDouble logTerm = LogTerm(rho);
85  // Entropy for an ideal gas
86  NekDouble sIg =
87  m_gasConstant / (m_gamma - 1) * log(T) - m_gasConstant * log(rho);
88 
89  NekDouble deltaS = m_gasConstant * log(1 - m_b * rho);
90  deltaS -= m_a * Alpha(T) * logTerm / (2 * m_b * T);
91 
92  return sIg + deltaS;
93 }
94 
96  const NekDouble &e)
97 {
98  NekDouble T = GetTemperature(rho, e);
99  NekDouble alpha = Alpha(T);
100  NekDouble dPde = GetDPDe_rho(rho, e);
101 
102  // Calculate dPdrho_T
103  NekDouble dPdrho_T =
104  m_gasConstant * T / ((1.0 - m_b * rho) * (1.0 - m_b * rho)) -
105  m_a * alpha * rho * (m_b * rho + 2) /
106  ((1 + m_b * rho) * (1 + m_b * rho));
107 
108  // Calculate dedrho_T
109  NekDouble dedrho_T = -3 * m_a * alpha / (2 * (1 + m_b * rho));
110 
111  // The result is dPdrho_e = dPdrho_T - dPde_rho * dedrho_T
112  return dPdrho_T - dPde * dedrho_T;
113 }
114 
116  const NekDouble &e)
117 {
118  NekDouble T = GetTemperature(rho, e);
119  NekDouble alpha = Alpha(T);
120  NekDouble logTerm = LogTerm(rho);
121 
122  // First calculate the denominator 1/rho^2 + 2*b/rho - b^2
123  // and sqrt(Alpha) = 1+f_w*(1-sqrt(Tr))
124  NekDouble denom = 1.0 / (rho * rho) + m_b / rho;
125 
126  // Compute cv = dedT_rho
127  NekDouble cv = m_gasConstant / (m_gamma - 1) +
128  3 * m_a * alpha * logTerm / (4 * m_b * T);
129 
130  // Now we obtain dPdT_rho
131  NekDouble dPdT =
132  m_gasConstant / (1.0 / rho - m_b) + m_a * alpha / (denom * 2 * T);
133 
134  // The result is dPde_rho = dPdT_rho / cv
135  return dPdT / cv;
136 }
137 
139  const NekDouble &p)
140 {
141  NekDouble logTerm = LogTerm(rho);
142  // First calculate the temperature, which can be expressed as
143  // (T^1/2)^3 + A* T^1/2 + B = 0
144  NekDouble A, B;
145 
146  A = -p * (1.0 / rho - m_b) / m_gasConstant;
147  B = -m_a * sqrt(m_Tc) * (1.0 / rho - m_b) /
148  (1.0 / (rho * rho) + m_b / rho) / m_gasConstant;
149 
150  // Use ideal gas solution as starting guess for iteration
151  NekDouble sqrtT = sqrt(p / (rho * (m_gamma - 1)));
152  // Newton-Raphson iteration to find T^(1/2)
153  NekDouble tol = 1e-6;
154  NekDouble maxIter = 100;
155  NekDouble residual = 1;
156  NekDouble f, df;
157  unsigned int cnt = 0;
158  while (abs(residual) > tol && cnt < maxIter)
159  {
160  f = sqrtT * sqrtT * sqrtT + A * sqrtT + B;
161  df = 3 * sqrtT * sqrtT + A;
162  residual = f / df;
163  sqrtT -= residual;
164  ++cnt;
165  }
166  if (cnt == maxIter)
167  {
168  cout << "Newton-Raphson in RedlichKwongEoS::v_GetEFromRhoP did not "
169  "converge in "
170  << maxIter << " iterations (residual = " << residual << ")"
171  << endl;
172  }
173 
174  // Calculate T
175  NekDouble T = sqrtT * sqrtT;
176 
177  // Calculate internal energy
178  return m_gasConstant * T / (m_gamma - 1) -
179  3 * m_a * Alpha(T) / (2 * m_b) * logTerm;
180 }
181 
183  const NekDouble &T)
184 {
185  // First solve for the compressibility factor Z using the cubic equation
186  // Z^3 + k1 * Z^2 + k2 * Z + k3 = 0
187  // for RedlichKwong:
188  // k1 = -1.0, k2 = A - B - B^2, k3 = -AB
189  // where A = a*alpha(T)*P/(RT)^2, B = bP/(RT)
190  NekDouble A = m_a * Alpha(T) * p / (m_gasConstant * m_gasConstant * T * T);
191  NekDouble B = m_b * p / (m_gasConstant * T);
192 
193  NekDouble k1 = -1.0;
194  NekDouble k2 = A - B - B * B;
195  NekDouble k3 = -A * B;
196 
197  // Use ideal gas (Z=1) as starting guess for iteration
198  NekDouble Z = 1.0;
199  // Newton-Raphson iteration to find Z
200  NekDouble tol = 1e-6;
201  NekDouble maxIter = 100;
202  NekDouble residual = 1;
203  NekDouble f, df;
204  unsigned int cnt = 0;
205  while (abs(residual) > tol && cnt < maxIter)
206  {
207  f = Z * Z * Z + k1 * Z * Z + k2 * Z + k3;
208  df = 3 * Z * Z + 2 * k1 * Z + k2;
209  residual = f / df;
210  Z -= residual;
211  ++cnt;
212  }
213  if (cnt == maxIter)
214  {
215  cout << "Newton-Raphson in RedlichKwongEoS::v_GetRhoFromPT did not "
216  "converge in "
217  << maxIter << " iterations (residual = " << residual << ")"
218  << endl;
219  }
220 
221  // Now calculate rho = p/(ZRT)
222  return p / (Z * m_gasConstant * T);
223 }
224 
225 } // namespace Nektar
Encapsulates equations of state allowing us to obtain thermodynamic properties: most relations are in...
NekDouble GetTemperature(const NekDouble &rho, const NekDouble &e)
Calculate the temperature.
NekDouble GetDPDe_rho(const NekDouble &rho, const NekDouble &e)
Calculate the partial derivative of P(rho,e) with respect to e.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
T GetTemperatureKernel(const T &rho, const T &e)
virtual NekDouble v_GetRhoFromPT(const NekDouble &rho, const NekDouble &p) override final
virtual NekDouble v_GetDPDe_rho(const NekDouble &rho, const NekDouble &e) override final
virtual NekDouble v_GetDPDrho_e(const NekDouble &rho, const NekDouble &e) override final
virtual NekDouble v_GetEFromRhoP(const NekDouble &rho, const NekDouble &p) override final
T Alpha(const T &temp)
virtual NekDouble v_GetEntropy(const NekDouble &rho, const NekDouble &e) override final
T GetPressureKernel(const T &rho, const T &e)
virtual NekDouble v_GetPressure(const NekDouble &rho, const NekDouble &e) override final
virtual NekDouble v_GetTemperature(const NekDouble &rho, const NekDouble &e) override final
std::shared_ptr< SessionReader > SessionReaderSharedPtr
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.
tinysimd::simd< NekDouble > vec_t
double NekDouble
scalarT< T > log(scalarT< T > in)
Definition: scalar.hpp:303
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:298
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294