Nektar++
PengRobinsonEoS.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: PengRobinsonEoS.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: Peng-Robinson equation of state
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include "PengRobinsonEoS.h"
36 
37 using namespace std;
38 
39 namespace Nektar
40 {
41 
42 std::string PengRobinsonEoS::className =
44  "PengRobinson", PengRobinsonEoS::create,
45  "Peng-Robinson equation of state.");
46 
47 PengRobinsonEoS::PengRobinsonEoS(
49  : EquationOfState(pSession)
50 {
51  pSession->LoadParameter("Tcrit", m_Tc);
52  pSession->LoadParameter("Pcrit", m_Pc);
53  pSession->LoadParameter("AcentricFactor", m_omega);
54 
55  m_a = 0.45724 * m_gasConstant * m_gasConstant * m_Tc * m_Tc / m_Pc;
56  m_b = 0.0778 * m_gasConstant * m_Tc / m_Pc;
57  m_fw = 0.37464 + 1.54226 * m_omega - 0.2699 * m_omega * m_omega;
58 }
59 
61  const NekDouble& rho, const NekDouble& e)
62 {
63  return GetTemperatureKernel(rho, e);
64 }
65 
67  const vec_t& rho, const vec_t& e)
68 {
69  return GetTemperatureKernel(rho, e);
70 }
71 
73  const NekDouble &rho, const NekDouble &e)
74 {
75  return GetPressureKernel(rho, e);
76 }
77 
79  const vec_t &rho, const vec_t &e)
80 {
81  return GetPressureKernel(rho, e);
82 }
83 
85  const NekDouble &e)
86 {
87  NekDouble T = GetTemperature(rho, e);
88  NekDouble logTerm = LogTerm(rho);
89  // Entropy for an ideal gas
90  NekDouble sIg =
91  m_gasConstant / (m_gamma - 1) * log(T) - m_gasConstant * log(rho);
92 
93  // First sqrt(Alpha) = 1+_w*(1-sqrt(Tr)) and sqrt(Tr)
94  NekDouble sqrtA = sqrt(Alpha(T));
95  NekDouble sqrtTr = sqrt(T / m_Tc);
96 
97  NekDouble deltaS = m_gasConstant * log(1 - m_b * rho);
98  deltaS += m_a * sqrtA * m_fw * logTerm * (sqrtTr / T) / (m_b * sqrt(8));
99 
100  return sIg + deltaS;
101 }
102 
104  const NekDouble &e)
105 {
106  NekDouble T = GetTemperature(rho, e);
107  NekDouble dPde = GetDPDe_rho(rho, e);
108 
109  // First calculate the denominator 1/rho^2 + 2*b/rho - b^2
110  // and alpha = [1+f_w*(1-sqrt(Tr))]^2
111  NekDouble denom = 1.0 / (rho * rho) + 2.0 * m_b / rho - m_b * m_b;
112  NekDouble alpha = Alpha(T);
113 
114  // Calculate dPdrho_T
115  NekDouble dPdrho_T =
116  m_gasConstant * T / ((1.0 - m_b * rho) * (1.0 - m_b * rho)) -
117  2 * m_a * alpha * rho * (1 + m_b * rho) /
118  ((denom * rho * rho) * (denom * rho * rho));
119 
120  // Calculate dedrho_T
121  NekDouble dedrho_T =
122  -m_a * sqrt(alpha) * (1.0 + m_fw) / (denom * rho * rho);
123 
124  // The result is dPdrho_e = dPdrho_T - dPde_rho * dedrho_T
125  return dPdrho_T - dPde * dedrho_T;
126 }
127 
129  const NekDouble &e)
130 {
131  NekDouble T = GetTemperature(rho, e);
132  NekDouble logTerm = LogTerm(rho);
133 
134  // First calculate the denominator 1/rho^2 + 2*b/rho - b^2
135  // and sqrt(Alpha) = 1+f_w*(1-sqrt(Tr))
136  NekDouble denom = 1.0 / (rho * rho) + 2.0 * m_b / rho - m_b * m_b;
137  NekDouble sqrtA = sqrt(Alpha(T));
138 
139  // Compute cv = dedT_rho
140  NekDouble cv = m_gasConstant / (m_gamma - 1) -
141  m_a / (2 * m_b * sqrt(8)) * logTerm * (m_fw * (1 + m_fw)) /
142  sqrt(T * m_Tc);
143 
144  // Now we obtain dPdT_rho
145  NekDouble dPdT = m_gasConstant / (1.0 / rho - m_b) +
146  m_a / sqrt(T * m_Tc) * m_fw * sqrtA / denom;
147 
148  // The result is dPde_rho = dPdT_rho / cv
149  return dPdT / cv;
150 }
151 
153  const NekDouble &p)
154 {
155  NekDouble denom = 1.0 / (rho * rho) + 2.0 * m_b / rho - m_b * m_b;
156  NekDouble logTerm = LogTerm(rho);
157  // First we solve for the temperature, which can be expressed as
158  // A * (T^1/2)^2 + B * T^1/2 + C = 0
159  NekDouble A, B, C;
160 
161  A = m_gasConstant / (1.0 / rho - m_b) -
162  (m_a * m_fw * m_fw) / (denom * m_Tc);
163  B = 2 * m_a / denom * m_fw * (1.0 + m_fw) / sqrt(m_Tc);
164  C = -m_a * (1.0 + m_fw) * (1 + m_fw) / denom - p;
165 
166  // Solve for T^1/2 (positive root)
167  NekDouble T = (-B + sqrt(B * B - 4 * A * C)) / (2 * A);
168  T = T * T;
169 
170  // Calculate alpha(T))
171  NekDouble alpha = Alpha(T);
172  // sqrt(Tr)
173  NekDouble sqrtTr = sqrt(T / m_Tc);
174  // Calculate internal energy
175  return m_gasConstant * T / (m_gamma - 1) +
176  m_a / (m_b * sqrt(8)) * logTerm *
177  (alpha + sqrt(alpha) * m_fw * sqrtTr);
178 }
179 
181  const NekDouble &T)
182 {
183  // First solve for the compressibility factor Z using the cubic equation
184  // Z^3 + k1 * Z^2 + k2 * Z + k3 = 0
185  // for PengRobinson:
186  // k1 = B-1, k2 = A - 2*B - 3*B^2, k3 = - AB + B^2 + B^3
187  // where A = a*alpha(T)*P/(RT)^2, B = bP/(RT)
188  NekDouble A = m_a * Alpha(T) * p / (m_gasConstant * m_gasConstant * T * T);
189  NekDouble B = m_b * p / (m_gasConstant * T);
190 
191  NekDouble k1 = B - 1.0;
192  NekDouble k2 = A - 2 * B - 3 * B * B;
193  NekDouble k3 = -A * B + B * B + B * B * B;
194 
195  // Use ideal gas (Z=1) as starting guess for iteration
196  NekDouble Z = 1.0;
197  // Newton-Raphson iteration to find Z
198  NekDouble tol = 1e-6;
199  NekDouble maxIter = 100;
200  NekDouble residual = 1;
201  NekDouble f, df;
202  unsigned int cnt = 0;
203  while ((fabs(residual) > tol) && cnt < maxIter)
204  {
205  f = Z * Z * Z + k1 * Z * Z + k2 * Z + k3;
206  df = 3 * Z * Z + 2 * k1 * Z + k2;
207  residual = f / df;
208  Z -= residual;
209  ++cnt;
210  }
211  if (cnt == maxIter)
212  {
213  cout << "Newton-Raphson in PengRobinsonEoS::v_GetRhoFromPT did not "
214  "converge in "
215  << maxIter << " iterations (residual = " << residual << ")"
216  << endl;
217  }
218 
219  // Now calculate rho = p/(ZRT)
220  return p / (Z * m_gasConstant * T);
221 }
222 
223 }
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)
T GetPressureKernel(const T &rho, const T &e)
NekDouble v_GetEFromRhoP(const NekDouble &rho, const NekDouble &p) final
NekDouble GetTemperature(const NekDouble &rho, const NekDouble &e) final
Calculate the temperature.
NekDouble v_GetDPDe_rho(const NekDouble &rho, const NekDouble &e) final
NekDouble GetPressure(const NekDouble &rho, const NekDouble &e) final
Calculate the pressure.
NekDouble v_GetRhoFromPT(const NekDouble &rho, const NekDouble &p) final
NekDouble v_GetDPDrho_e(const NekDouble &rho, const NekDouble &e) final
T Alpha(const T &temp)
NekDouble v_GetEntropy(const NekDouble &rho, const NekDouble &e) 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 > sqrt(scalarT< T > in)
Definition: scalar.hpp:267