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 &e)
62 {
63  return GetTemperatureKernel(rho, e);
64 }
65 
67 {
68  return GetTemperatureKernel(rho, e);
69 }
70 
72  const NekDouble &e)
73 {
74  return GetPressureKernel(rho, e);
75 }
76 
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  // First sqrt(Alpha) = 1+_w*(1-sqrt(Tr)) and sqrt(Tr)
92  NekDouble sqrtA = sqrt(Alpha(T));
93  NekDouble sqrtTr = sqrt(T / m_Tc);
94 
95  NekDouble deltaS = m_gasConstant * log(1 - m_b * rho);
96  deltaS += m_a * sqrtA * m_fw * logTerm * (sqrtTr / T) / (m_b * sqrt(8));
97 
98  return sIg + deltaS;
99 }
100 
102  const NekDouble &e)
103 {
104  NekDouble T = GetTemperature(rho, e);
105  NekDouble dPde = GetDPDe_rho(rho, e);
106 
107  // First calculate the denominator 1/rho^2 + 2*b/rho - b^2
108  // and alpha = [1+f_w*(1-sqrt(Tr))]^2
109  NekDouble denom = 1.0 / (rho * rho) + 2.0 * m_b / rho - m_b * m_b;
110  NekDouble alpha = Alpha(T);
111 
112  // Calculate dPdrho_T
113  NekDouble dPdrho_T =
114  m_gasConstant * T / ((1.0 - m_b * rho) * (1.0 - m_b * rho)) -
115  2 * m_a * alpha * rho * (1 + m_b * rho) /
116  ((denom * rho * rho) * (denom * rho * rho));
117 
118  // Calculate dedrho_T
119  NekDouble dedrho_T =
120  -m_a * sqrt(alpha) * (1.0 + m_fw) / (denom * rho * rho);
121 
122  // The result is dPdrho_e = dPdrho_T - dPde_rho * dedrho_T
123  return dPdrho_T - dPde * dedrho_T;
124 }
125 
127  const NekDouble &e)
128 {
129  NekDouble T = GetTemperature(rho, e);
130  NekDouble logTerm = LogTerm(rho);
131 
132  // First calculate the denominator 1/rho^2 + 2*b/rho - b^2
133  // and sqrt(Alpha) = 1+f_w*(1-sqrt(Tr))
134  NekDouble denom = 1.0 / (rho * rho) + 2.0 * m_b / rho - m_b * m_b;
135  NekDouble sqrtA = sqrt(Alpha(T));
136 
137  // Compute cv = dedT_rho
138  NekDouble cv = m_gasConstant / (m_gamma - 1) -
139  m_a / (2 * m_b * sqrt(8)) * logTerm * (m_fw * (1 + m_fw)) /
140  sqrt(T * m_Tc);
141 
142  // Now we obtain dPdT_rho
143  NekDouble dPdT = m_gasConstant / (1.0 / rho - m_b) +
144  m_a / sqrt(T * m_Tc) * m_fw * sqrtA / denom;
145 
146  // The result is dPde_rho = dPdT_rho / cv
147  return dPdT / cv;
148 }
149 
151  const NekDouble &p)
152 {
153  NekDouble denom = 1.0 / (rho * rho) + 2.0 * m_b / rho - m_b * m_b;
154  NekDouble logTerm = LogTerm(rho);
155  // First we solve for the temperature, which can be expressed as
156  // A * (T^1/2)^2 + B * T^1/2 + C = 0
157  NekDouble A, B, C;
158 
159  A = m_gasConstant / (1.0 / rho - m_b) -
160  (m_a * m_fw * m_fw) / (denom * m_Tc);
161  B = 2 * m_a / denom * m_fw * (1.0 + m_fw) / sqrt(m_Tc);
162  C = -m_a * (1.0 + m_fw) * (1 + m_fw) / denom - p;
163 
164  // Solve for T^1/2 (positive root)
165  NekDouble T = (-B + sqrt(B * B - 4 * A * C)) / (2 * A);
166  T = T * T;
167 
168  // Calculate alpha(T))
169  NekDouble alpha = Alpha(T);
170  // sqrt(Tr)
171  NekDouble sqrtTr = sqrt(T / m_Tc);
172  // Calculate internal energy
173  return m_gasConstant * T / (m_gamma - 1) +
174  m_a / (m_b * sqrt(8)) * logTerm *
175  (alpha + sqrt(alpha) * m_fw * sqrtTr);
176 }
177 
179  const NekDouble &T)
180 {
181  // First solve for the compressibility factor Z using the cubic equation
182  // Z^3 + k1 * Z^2 + k2 * Z + k3 = 0
183  // for PengRobinson:
184  // k1 = B-1, k2 = A - 2*B - 3*B^2, k3 = - AB + B^2 + B^3
185  // where A = a*alpha(T)*P/(RT)^2, B = bP/(RT)
186  NekDouble A = m_a * Alpha(T) * p / (m_gasConstant * m_gasConstant * T * T);
187  NekDouble B = m_b * p / (m_gasConstant * T);
188 
189  NekDouble k1 = B - 1.0;
190  NekDouble k2 = A - 2 * B - 3 * B * B;
191  NekDouble k3 = -A * B + B * B + B * B * B;
192 
193  // Use ideal gas (Z=1) as starting guess for iteration
194  NekDouble Z = 1.0;
195  // Newton-Raphson iteration to find Z
196  NekDouble tol = 1e-6;
197  NekDouble maxIter = 100;
198  NekDouble residual = 1;
199  NekDouble f, df;
200  unsigned int cnt = 0;
201  while ((fabs(residual) > tol) && cnt < maxIter)
202  {
203  f = Z * Z * Z + k1 * Z * Z + k2 * Z + k3;
204  df = 3 * Z * Z + 2 * k1 * Z + k2;
205  residual = f / df;
206  Z -= residual;
207  ++cnt;
208  }
209  if (cnt == maxIter)
210  {
211  cout << "Newton-Raphson in PengRobinsonEoS::v_GetRhoFromPT did not "
212  "converge in "
213  << maxIter << " iterations (residual = " << residual << ")"
214  << endl;
215  }
216 
217  // Now calculate rho = p/(ZRT)
218  return p / (Z * m_gasConstant * T);
219 }
220 
221 } // 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
virtual NekDouble v_GetDPDe_rho(const NekDouble &rho, const NekDouble &e) override final
T GetTemperatureKernel(const T &rho, const T &e)
virtual NekDouble v_GetPressure(const NekDouble &rho, const NekDouble &e) override final
T GetPressureKernel(const T &rho, const T &e)
virtual NekDouble v_GetDPDrho_e(const NekDouble &rho, const NekDouble &e) override final
virtual NekDouble v_GetEntropy(const NekDouble &rho, const NekDouble &e) override final
virtual NekDouble v_GetTemperature(const NekDouble &rho, const NekDouble &e) override final
virtual NekDouble v_GetRhoFromPT(const NekDouble &rho, const NekDouble &p) override final
virtual NekDouble v_GetEFromRhoP(const NekDouble &rho, const NekDouble &p) override final
T Alpha(const T &temp)
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 > sqrt(scalarT< T > in)
Definition: scalar.hpp:294