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& rho, const NekDouble& e)
60 {
61  return GetTemperatureKernel(rho, e);
62 }
63 
65  const vec_t& rho, const vec_t& e)
66 {
67  return GetTemperatureKernel(rho, e);
68 }
69 
71  const NekDouble& rho, const NekDouble& e)
72 {
73  return GetPressureKernel(rho, e);
74 }
75 
77  const vec_t& rho, const vec_t& e)
78 {
79  return GetPressureKernel(rho, e);
80 }
81 
83  const NekDouble &e)
84 {
85  NekDouble T = GetTemperature(rho, e);
86  NekDouble logTerm = LogTerm(rho);
87  // Entropy for an ideal gas
88  NekDouble sIg =
89  m_gasConstant / (m_gamma - 1) * log(T) - m_gasConstant * log(rho);
90 
91  NekDouble deltaS = m_gasConstant * log(1 - m_b * rho);
92  deltaS -= m_a * Alpha(T) * logTerm / (2 * m_b * T);
93 
94  return sIg + deltaS;
95 }
96 
98  const NekDouble &e)
99 {
100  NekDouble T = GetTemperature(rho, e);
101  NekDouble alpha = Alpha(T);
102  NekDouble dPde = GetDPDe_rho(rho, e);
103 
104  // Calculate dPdrho_T
105  NekDouble dPdrho_T =
106  m_gasConstant * T / ((1.0 - m_b * rho) * (1.0 - m_b * rho)) -
107  m_a * alpha * rho * (m_b * rho + 2) /
108  ((1 + m_b * rho) * (1 + m_b * rho));
109 
110  // Calculate dedrho_T
111  NekDouble dedrho_T = -3 * m_a * alpha / (2 * (1 + m_b * rho));
112 
113  // The result is dPdrho_e = dPdrho_T - dPde_rho * dedrho_T
114  return dPdrho_T - dPde * dedrho_T;
115 }
116 
118  const NekDouble &e)
119 {
120  NekDouble T = GetTemperature(rho, e);
121  NekDouble alpha = Alpha(T);
122  NekDouble logTerm = LogTerm(rho);
123 
124  // First calculate the denominator 1/rho^2 + 2*b/rho - b^2
125  // and sqrt(Alpha) = 1+f_w*(1-sqrt(Tr))
126  NekDouble denom = 1.0 / (rho * rho) + m_b / rho;
127 
128  // Compute cv = dedT_rho
129  NekDouble cv = m_gasConstant / (m_gamma - 1) +
130  3 * m_a * alpha * logTerm / (4 * m_b * T);
131 
132  // Now we obtain dPdT_rho
133  NekDouble dPdT =
134  m_gasConstant / (1.0 / rho - m_b) + m_a * alpha / (denom * 2 * T);
135 
136  // The result is dPde_rho = dPdT_rho / cv
137  return dPdT / cv;
138 }
139 
141  const NekDouble &p)
142 {
143  NekDouble logTerm = LogTerm(rho);
144  // First calculate the temperature, which can be expressed as
145  // (T^1/2)^3 + A* T^1/2 + B = 0
146  NekDouble A, B;
147 
148  A = -p * (1.0 / rho - m_b) / m_gasConstant;
149  B = -m_a * sqrt(m_Tc) * (1.0 / rho - m_b) /
150  (1.0 / (rho * rho) + m_b / rho) / m_gasConstant;
151 
152  // Use ideal gas solution as starting guess for iteration
153  NekDouble sqrtT = sqrt(p / (rho * (m_gamma - 1)));
154  // Newton-Raphson iteration to find T^(1/2)
155  NekDouble tol = 1e-6;
156  NekDouble maxIter = 100;
157  NekDouble residual = 1;
158  NekDouble f, df;
159  unsigned int cnt = 0;
160  while (abs(residual) > tol && cnt < maxIter)
161  {
162  f = sqrtT * sqrtT * sqrtT + A * sqrtT + B;
163  df = 3 * sqrtT * sqrtT + A;
164  residual = f / df;
165  sqrtT -= residual;
166  ++cnt;
167  }
168  if (cnt == maxIter)
169  {
170  cout << "Newton-Raphson in RedlichKwongEoS::v_GetEFromRhoP did not "
171  "converge in "
172  << maxIter << " iterations (residual = " << residual << ")"
173  << endl;
174  }
175 
176  // Calculate T
177  NekDouble T = sqrtT * sqrtT;
178 
179  // Calculate internal energy
180  return m_gasConstant * T / (m_gamma - 1) -
181  3 * m_a * Alpha(T) / (2 * m_b) * logTerm;
182 }
183 
185  const NekDouble &T)
186 {
187  // First solve for the compressibility factor Z using the cubic equation
188  // Z^3 + k1 * Z^2 + k2 * Z + k3 = 0
189  // for RedlichKwong:
190  // k1 = -1.0, k2 = A - B - B^2, k3 = -AB
191  // where A = a*alpha(T)*P/(RT)^2, B = bP/(RT)
192  NekDouble A = m_a * Alpha(T) * p / (m_gasConstant * m_gasConstant * T * T);
193  NekDouble B = m_b * p / (m_gasConstant * T);
194 
195  NekDouble k1 = -1.0;
196  NekDouble k2 = A - B - B * B;
197  NekDouble k3 = -A * B;
198 
199  // Use ideal gas (Z=1) as starting guess for iteration
200  NekDouble Z = 1.0;
201  // Newton-Raphson iteration to find Z
202  NekDouble tol = 1e-6;
203  NekDouble maxIter = 100;
204  NekDouble residual = 1;
205  NekDouble f, df;
206  unsigned int cnt = 0;
207  while (abs(residual) > tol && cnt < maxIter)
208  {
209  f = Z * Z * Z + k1 * Z * Z + k2 * Z + k3;
210  df = 3 * Z * Z + 2 * k1 * Z + k2;
211  residual = f / df;
212  Z -= residual;
213  ++cnt;
214  }
215  if (cnt == maxIter)
216  {
217  cout << "Newton-Raphson in RedlichKwongEoS::v_GetRhoFromPT did not "
218  "converge in "
219  << maxIter << " iterations (residual = " << residual << ")"
220  << endl;
221  }
222 
223  // Now calculate rho = p/(ZRT)
224  return p / (Z * m_gasConstant * T);
225 }
226 
227 }
Encapsulates equations of state allowing us to obtain thermodynamic properties: most relations are in...
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:200
T GetTemperatureKernel(const T &rho, const T &e)
NekDouble v_GetDPDrho_e(const NekDouble &rho, const NekDouble &e) final
NekDouble GetPressure(const NekDouble &rho, const NekDouble &e) final
Calculate the pressure.
NekDouble GetTemperature(const NekDouble &rho, const NekDouble &e) final
Calculate the temperature.
T Alpha(const T &temp)
NekDouble v_GetRhoFromPT(const NekDouble &rho, const NekDouble &p) final
NekDouble v_GetDPDe_rho(const NekDouble &rho, const NekDouble &e) final
NekDouble v_GetEntropy(const NekDouble &rho, const NekDouble &e) final
T GetPressureKernel(const T &rho, const T &e)
NekDouble v_GetEFromRhoP(const NekDouble &rho, const NekDouble &p) final
std::shared_ptr< SessionReader > SessionReaderSharedPtr
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
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:278
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:272
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267