Nektar++
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Member Functions | Friends | List of all members
Nektar::PengRobinsonEoS Class Reference

Peng-Robinson equation of state: p = RT/(1/rho - b) - a*Alpha(T/Tc) / (1/rho^2 + 2*b/rho - b^2) with a = 0.45724 * (R*Tc)^2 / Pc b = 0.0778 * (R*Tc) / Pc Alpha(T/Tc) = [1 + fw * (1 - sqrt(T/ Tc))]^2 fw = 0.37464 + 1.54226*omega - 0.2699*omega*omega. More...

#include <PengRobinsonEoS.h>

Inheritance diagram for Nektar::PengRobinsonEoS:
[legend]

Static Public Member Functions

static EquationOfStateSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of the class. More...
 

Protected Member Functions

virtual NekDouble v_GetTemperature (const NekDouble &rho, const NekDouble &e)
 
virtual NekDouble v_GetPressure (const NekDouble &rho, const NekDouble &e)
 
virtual NekDouble v_GetEntropy (const NekDouble &rho, const NekDouble &e)
 
virtual NekDouble v_GetDPDrho_e (const NekDouble &rho, const NekDouble &e)
 
virtual NekDouble v_GetDPDe_rho (const NekDouble &rho, const NekDouble &e)
 
virtual NekDouble v_GetEFromRhoP (const NekDouble &rho, const NekDouble &p)
 
virtual NekDouble v_GetRhoFromPT (const NekDouble &rho, const NekDouble &p)
 
- Protected Member Functions inherited from Nektar::EquationOfState
 EquationOfState (const LibUtilities::SessionReaderSharedPtr &pSession)
 Constructor. More...
 
virtual NekDouble v_GetSoundSpeed (const NekDouble &rho, const NekDouble &e)
 

Protected Attributes

NekDouble m_a
 
NekDouble m_b
 
NekDouble m_Tc
 
NekDouble m_Pc
 
NekDouble m_omega
 
NekDouble m_fw
 
- Protected Attributes inherited from Nektar::EquationOfState
NekDouble m_gamma
 
NekDouble m_gasConstant
 

Private Member Functions

 PengRobinsonEoS (const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual ~PengRobinsonEoS (void)
 
NekDouble Alpha (const NekDouble &T)
 
NekDouble LogTerm (const NekDouble &rho)
 

Friends

class MemoryManager< PengRobinsonEoS >
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::EquationOfState
virtual ~EquationOfState ()
 
NekDouble GetTemperature (const NekDouble &rho, const NekDouble &e)
 Calculate the temperature. More...
 
NekDouble GetPressure (const NekDouble &rho, const NekDouble &e)
 Calculate the pressure. More...
 
NekDouble GetSoundSpeed (const NekDouble &rho, const NekDouble &e)
 Calculate the sound speed. More...
 
NekDouble GetEntropy (const NekDouble &rho, const NekDouble &e)
 Calculate the entropy. More...
 
NekDouble GetDPDrho_e (const NekDouble &rho, const NekDouble &e)
 Calculate the partial derivative of P(rho,e) with respect to rho. More...
 
NekDouble GetDPDe_rho (const NekDouble &rho, const NekDouble &e)
 Calculate the partial derivative of P(rho,e) with respect to e. More...
 
NekDouble GetEFromRhoP (const NekDouble &rho, const NekDouble &p)
 Obtain the internal energy from rho and P. More...
 
NekDouble GetRhoFromPT (const NekDouble &p, const NekDouble &T)
 Obtain the density from P and T. More...
 

Detailed Description

Peng-Robinson equation of state: p = RT/(1/rho - b) - a*Alpha(T/Tc) / (1/rho^2 + 2*b/rho - b^2) with a = 0.45724 * (R*Tc)^2 / Pc b = 0.0778 * (R*Tc) / Pc Alpha(T/Tc) = [1 + fw * (1 - sqrt(T/ Tc))]^2 fw = 0.37464 + 1.54226*omega - 0.2699*omega*omega.

Definition at line 51 of file PengRobinsonEoS.h.

Constructor & Destructor Documentation

◆ PengRobinsonEoS()

Nektar::PengRobinsonEoS::PengRobinsonEoS ( const LibUtilities::SessionReaderSharedPtr pSession)
private

Definition at line 47 of file PengRobinsonEoS.cpp.

References m_a, m_b, m_fw, Nektar::EquationOfState::m_gasConstant, m_omega, m_Pc, and m_Tc.

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 }
EquationOfState(const LibUtilities::SessionReaderSharedPtr &pSession)
Constructor.

◆ ~PengRobinsonEoS()

virtual Nektar::PengRobinsonEoS::~PengRobinsonEoS ( void  )
inlineprivatevirtual

Definition at line 94 of file PengRobinsonEoS.h.

References Alpha(), and LogTerm().

94 {};

Member Function Documentation

◆ Alpha()

NekDouble Nektar::PengRobinsonEoS::Alpha ( const NekDouble T)
private

Definition at line 233 of file PengRobinsonEoS.cpp.

References m_fw, and m_Tc.

Referenced by v_GetDPDe_rho(), v_GetDPDrho_e(), v_GetEFromRhoP(), v_GetEntropy(), v_GetPressure(), v_GetRhoFromPT(), and ~PengRobinsonEoS().

234 {
235  NekDouble sqrtAlpha = 1.0 + m_fw * (1.0 - sqrt(T / m_Tc));
236  return sqrtAlpha * sqrtAlpha;
237 }
double NekDouble

◆ create()

static EquationOfStateSharedPtr Nektar::PengRobinsonEoS::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Creates an instance of this class.

Definition at line 57 of file PengRobinsonEoS.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

59  {
62  return p;
63  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< EquationOfState > EquationOfStateSharedPtr
A shared pointer to an equation of state object.

◆ LogTerm()

NekDouble Nektar::PengRobinsonEoS::LogTerm ( const NekDouble rho)
private

Definition at line 239 of file PengRobinsonEoS.cpp.

References m_b.

Referenced by v_GetDPDe_rho(), v_GetEFromRhoP(), v_GetEntropy(), v_GetTemperature(), and ~PengRobinsonEoS().

240 {
241  return log((1.0 / rho + m_b - m_b * sqrt(2)) /
242  (1.0 / rho + m_b + m_b * sqrt(2)));
243 }

◆ v_GetDPDe_rho()

NekDouble Nektar::PengRobinsonEoS::v_GetDPDe_rho ( const NekDouble rho,
const NekDouble e 
)
protectedvirtual

Implements Nektar::EquationOfState.

Definition at line 138 of file PengRobinsonEoS.cpp.

References Alpha(), Nektar::EquationOfState::GetTemperature(), LogTerm(), m_a, m_b, m_fw, Nektar::EquationOfState::m_gamma, Nektar::EquationOfState::m_gasConstant, and m_Tc.

140 {
141  NekDouble T = GetTemperature(rho, e);
142  NekDouble logTerm = LogTerm(rho);
143 
144  // First calculate the denominator 1/rho^2 + 2*b/rho - b^2
145  // and sqrt(Alpha) = 1+f_w*(1-sqrt(Tr))
146  NekDouble denom = 1.0 / (rho * rho) + 2.0 * m_b / rho - m_b * m_b;
147  NekDouble sqrtA = sqrt(Alpha(T));
148 
149  // Compute cv = dedT_rho
150  NekDouble cv = m_gasConstant / (m_gamma - 1) -
151  m_a / (2 * m_b * sqrt(8)) * logTerm * (m_fw * (1 + m_fw)) /
152  sqrt(T * m_Tc);
153 
154  // Now we obtain dPdT_rho
155  NekDouble dPdT = m_gasConstant / (1.0 / rho - m_b) +
156  m_a / sqrt(T * m_Tc) * m_fw * sqrtA / denom;
157 
158  // The result is dPde_rho = dPdT_rho / cv
159  return dPdT / cv;
160 }
NekDouble LogTerm(const NekDouble &rho)
NekDouble Alpha(const NekDouble &T)
double NekDouble
NekDouble GetTemperature(const NekDouble &rho, const NekDouble &e)
Calculate the temperature.

◆ v_GetDPDrho_e()

NekDouble Nektar::PengRobinsonEoS::v_GetDPDrho_e ( const NekDouble rho,
const NekDouble e 
)
protectedvirtual

Implements Nektar::EquationOfState.

Definition at line 113 of file PengRobinsonEoS.cpp.

References Alpha(), Nektar::EquationOfState::GetDPDe_rho(), Nektar::EquationOfState::GetTemperature(), m_a, m_b, m_fw, and Nektar::EquationOfState::m_gasConstant.

115 {
116  NekDouble T = GetTemperature(rho, e);
117  NekDouble dPde = GetDPDe_rho(rho, e);
118 
119  // First calculate the denominator 1/rho^2 + 2*b/rho - b^2
120  // and alpha = [1+f_w*(1-sqrt(Tr))]^2
121  NekDouble denom = 1.0 / (rho * rho) + 2.0 * m_b / rho - m_b * m_b;
122  NekDouble alpha = Alpha(T);
123 
124  // Calculate dPdrho_T
125  NekDouble dPdrho_T =
126  m_gasConstant * T / ((1.0 - m_b * rho) * (1.0 - m_b * rho)) -
127  2 * m_a * alpha * rho * (1 + m_b * rho) /
128  ((denom * rho * rho) * (denom * rho * rho));
129 
130  // Calculate dedrho_T
131  NekDouble dedrho_T =
132  -m_a * sqrt(alpha) * (1.0 + m_fw) / (denom * rho * rho);
133 
134  // The result is dPdrho_e = dPdrho_T - dPde_rho * dedrho_T
135  return dPdrho_T - dPde * dedrho_T;
136 }
NekDouble GetDPDe_rho(const NekDouble &rho, const NekDouble &e)
Calculate the partial derivative of P(rho,e) with respect to e.
NekDouble Alpha(const NekDouble &T)
double NekDouble
NekDouble GetTemperature(const NekDouble &rho, const NekDouble &e)
Calculate the temperature.

◆ v_GetEFromRhoP()

NekDouble Nektar::PengRobinsonEoS::v_GetEFromRhoP ( const NekDouble rho,
const NekDouble p 
)
protectedvirtual

Implements Nektar::EquationOfState.

Definition at line 162 of file PengRobinsonEoS.cpp.

References Alpha(), LogTerm(), m_a, m_b, m_fw, Nektar::EquationOfState::m_gamma, Nektar::EquationOfState::m_gasConstant, m_Tc, and CellMLToNektar.cellml_metadata::p.

164 {
165  NekDouble denom = 1.0 / (rho * rho) + 2.0 * m_b / rho - m_b * m_b;
166  NekDouble logTerm = LogTerm(rho);
167  // First we solve for the temperature, which can be expressed as
168  // A * (T^1/2)^2 + B * T^1/2 + C = 0
169  NekDouble A, B, C;
170 
171  A = m_gasConstant / (1.0 / rho - m_b) -
172  (m_a * m_fw * m_fw) / (denom * m_Tc);
173  B = 2 * m_a / denom * m_fw * (1.0 + m_fw) / sqrt(m_Tc);
174  C = -m_a * (1.0 + m_fw) * (1 + m_fw) / denom - p;
175 
176  // Solve for T^1/2 (positive root)
177  NekDouble T = (-B + sqrt(B * B - 4 * A * C)) / (2 * A);
178  T = T * T;
179 
180  // Calculate alpha(T))
181  NekDouble alpha = Alpha(T);
182  // sqrt(Tr)
183  NekDouble sqrtTr = sqrt(T / m_Tc);
184  // Calculate internal energy
185  return m_gasConstant * T / (m_gamma - 1) +
186  m_a / (m_b * sqrt(8)) * logTerm *
187  (alpha + sqrt(alpha) * m_fw * sqrtTr);
188 }
NekDouble LogTerm(const NekDouble &rho)
NekDouble Alpha(const NekDouble &T)
double NekDouble

◆ v_GetEntropy()

NekDouble Nektar::PengRobinsonEoS::v_GetEntropy ( const NekDouble rho,
const NekDouble e 
)
protectedvirtual

Implements Nektar::EquationOfState.

Definition at line 94 of file PengRobinsonEoS.cpp.

References Alpha(), Nektar::EquationOfState::GetTemperature(), LogTerm(), m_a, m_b, m_fw, Nektar::EquationOfState::m_gamma, Nektar::EquationOfState::m_gasConstant, and m_Tc.

96 {
97  NekDouble T = GetTemperature(rho, e);
98  NekDouble logTerm = LogTerm(rho);
99  // Entropy for an ideal gas
100  NekDouble sIg =
101  m_gasConstant / (m_gamma - 1) * log(T) - m_gasConstant * log(rho);
102 
103  // First sqrt(Alpha) = 1+_w*(1-sqrt(Tr)) and sqrt(Tr)
104  NekDouble sqrtA = sqrt(Alpha(T));
105  NekDouble sqrtTr = sqrt(T / m_Tc);
106 
107  NekDouble deltaS = m_gasConstant * log(1 - m_b * rho);
108  deltaS += m_a * sqrtA * m_fw * logTerm * (sqrtTr / T) / (m_b * sqrt(8));
109 
110  return sIg + deltaS;
111 }
NekDouble LogTerm(const NekDouble &rho)
NekDouble Alpha(const NekDouble &T)
double NekDouble
NekDouble GetTemperature(const NekDouble &rho, const NekDouble &e)
Calculate the temperature.

◆ v_GetPressure()

NekDouble Nektar::PengRobinsonEoS::v_GetPressure ( const NekDouble rho,
const NekDouble e 
)
protectedvirtual

Implements Nektar::EquationOfState.

Definition at line 82 of file PengRobinsonEoS.cpp.

References Alpha(), Nektar::EquationOfState::GetTemperature(), m_a, m_b, Nektar::EquationOfState::m_gasConstant, and CellMLToNektar.cellml_metadata::p.

84 {
85  NekDouble T = GetTemperature(rho, e);
86 
87  NekDouble p =
88  m_gasConstant * T / (1.0 / rho - m_b) -
89  m_a * Alpha(T) / (1.0 / (rho * rho) + 2.0 * m_b / rho - m_b * m_b);
90 
91  return p;
92 }
NekDouble Alpha(const NekDouble &T)
double NekDouble
NekDouble GetTemperature(const NekDouble &rho, const NekDouble &e)
Calculate the temperature.

◆ v_GetRhoFromPT()

NekDouble Nektar::PengRobinsonEoS::v_GetRhoFromPT ( const NekDouble rho,
const NekDouble p 
)
protectedvirtual

Implements Nektar::EquationOfState.

Definition at line 190 of file PengRobinsonEoS.cpp.

References Alpha(), m_a, m_b, and Nektar::EquationOfState::m_gasConstant.

192 {
193  // First solve for the compressibility factor Z using the cubic equation
194  // Z^3 + k1 * Z^2 + k2 * Z + k3 = 0
195  // for PengRobinson:
196  // k1 = B-1, k2 = A - 2*B - 3*B^2, k3 = - AB + B^2 + B^3
197  // where A = a*alpha(T)*P/(RT)^2, B = bP/(RT)
198  NekDouble A = m_a * Alpha(T) * p / (m_gasConstant * m_gasConstant * T * T);
199  NekDouble B = m_b * p / (m_gasConstant * T);
200 
201  NekDouble k1 = B - 1.0;
202  NekDouble k2 = A - 2 * B - 3 * B * B;
203  NekDouble k3 = -A * B + B * B + B * B * B;
204 
205  // Use ideal gas (Z=1) as starting guess for iteration
206  NekDouble Z = 1.0;
207  // Newton-Raphson iteration to find Z
208  NekDouble tol = 1e-6;
209  NekDouble maxIter = 100;
210  NekDouble residual = 1;
211  NekDouble f, df;
212  unsigned int cnt = 0;
213  while (abs(residual) > tol && cnt < maxIter)
214  {
215  f = Z * Z * Z + k1 * Z * Z + k2 * Z + k3;
216  df = 3 * Z * Z + 2 * k1 * Z + k2;
217  residual = f / df;
218  Z -= residual;
219  ++cnt;
220  }
221  if (cnt == maxIter)
222  {
223  cout << "Newton-Raphson in PengRobinsonEoS::v_GetRhoFromPT did not "
224  "converge in "
225  << maxIter << " iterations (residual = " << residual << ")"
226  << endl;
227  }
228 
229  // Now calculate rho = p/(ZRT)
230  return p / (Z * m_gasConstant * T);
231 }
NekDouble Alpha(const NekDouble &T)
double NekDouble

◆ v_GetTemperature()

NekDouble Nektar::PengRobinsonEoS::v_GetTemperature ( const NekDouble rho,
const NekDouble e 
)
protectedvirtual

Implements Nektar::EquationOfState.

Definition at line 60 of file PengRobinsonEoS.cpp.

References LogTerm(), m_a, m_b, m_fw, Nektar::EquationOfState::m_gamma, Nektar::EquationOfState::m_gasConstant, and m_Tc.

62 {
63  // First we need to evaluate the log term
64  // ln[(1/rho + b - b*sqrt(2)) / (1/rho + b + b*sqrt(2))]
65  NekDouble sqrt2 = sqrt(2.0);
66  NekDouble logTerm = LogTerm(rho);
67 
68  // The temperature can be expressed as an equation in the form
69  // A * (T^1/2)^2 + B * T^1/2 + C = 0
70  NekDouble A, B, C;
71 
72  A = m_gasConstant / (m_gamma - 1);
73  B = -m_a / (m_b * 2 * sqrt2) * logTerm / sqrt(m_Tc) * m_fw * (1 + m_fw);
74  C = m_a / (m_b * 2 * sqrt2) * logTerm * (1 + m_fw) * (1 + m_fw) - e;
75 
76  // Solve for T^1/2 (positive root)
77  NekDouble sqrtT = (-B + sqrt(B * B - 4 * A * C)) / (2 * A);
78  // Calculate the temperature
79  return sqrtT * sqrtT;
80 }
NekDouble LogTerm(const NekDouble &rho)
double NekDouble

Friends And Related Function Documentation

◆ MemoryManager< PengRobinsonEoS >

friend class MemoryManager< PengRobinsonEoS >
friend

Definition at line 54 of file PengRobinsonEoS.h.

Member Data Documentation

◆ className

std::string Nektar::PengRobinsonEoS::className
static
Initial value:
=
"PengRobinson", PengRobinsonEoS::create,
"Peng-Robinson equation of state.")

Name of the class.

Definition at line 66 of file PengRobinsonEoS.h.

◆ m_a

NekDouble Nektar::PengRobinsonEoS::m_a
protected

◆ m_b

NekDouble Nektar::PengRobinsonEoS::m_b
protected

◆ m_fw

NekDouble Nektar::PengRobinsonEoS::m_fw
protected

◆ m_omega

NekDouble Nektar::PengRobinsonEoS::m_omega
protected

Definition at line 73 of file PengRobinsonEoS.h.

Referenced by PengRobinsonEoS().

◆ m_Pc

NekDouble Nektar::PengRobinsonEoS::m_Pc
protected

Definition at line 72 of file PengRobinsonEoS.h.

Referenced by PengRobinsonEoS().

◆ m_Tc

NekDouble Nektar::PengRobinsonEoS::m_Tc
protected