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

Redlich-Kwong equation of state: p = RT/(1/rho - b) - a/( sqrt(T / Tc) * (1/rho^2 + b/rho) with a = 0.42748 * (R*Tc)^2 / Pc b = 0.08664 * (R*Tc) / Pc. More...

#include <RedlichKwongEoS.h>

Inheritance diagram for Nektar::RedlichKwongEoS:
[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
 
- Protected Attributes inherited from Nektar::EquationOfState
NekDouble m_gamma
 
NekDouble m_gasConstant
 

Private Member Functions

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

Friends

class MemoryManager< RedlichKwongEoS >
 

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

Redlich-Kwong equation of state: p = RT/(1/rho - b) - a/( sqrt(T / Tc) * (1/rho^2 + b/rho) with a = 0.42748 * (R*Tc)^2 / Pc b = 0.08664 * (R*Tc) / Pc.

Definition at line 49 of file RedlichKwongEoS.h.

Constructor & Destructor Documentation

◆ RedlichKwongEoS()

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

Definition at line 47 of file RedlichKwongEoS.cpp.

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

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

◆ ~RedlichKwongEoS()

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

Definition at line 90 of file RedlichKwongEoS.h.

References Alpha(), and LogTerm().

90 {};

Member Function Documentation

◆ Alpha()

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

Definition at line 257 of file RedlichKwongEoS.cpp.

References m_Tc.

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

258 {
259  return 1.0 / sqrt(T / m_Tc);
260 }

◆ create()

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

Creates an instance of this class.

Definition at line 55 of file RedlichKwongEoS.h.

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

57  {
60  return p;
61  }
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::RedlichKwongEoS::LogTerm ( const NekDouble rho)
private

Definition at line 262 of file RedlichKwongEoS.cpp.

References m_b.

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

263 {
264  return log(1 + m_b * rho);
265 }

◆ v_GetDPDe_rho()

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

Implements Nektar::EquationOfState.

Definition at line 147 of file RedlichKwongEoS.cpp.

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

149 {
150  NekDouble T = GetTemperature(rho, e);
151  NekDouble alpha = Alpha(T);
152  NekDouble logTerm = LogTerm(rho);
153 
154  // First calculate the denominator 1/rho^2 + 2*b/rho - b^2
155  // and sqrt(Alpha) = 1+f_w*(1-sqrt(Tr))
156  NekDouble denom = 1.0 / (rho * rho) + m_b / rho;
157 
158  // Compute cv = dedT_rho
159  NekDouble cv = m_gasConstant / (m_gamma - 1) +
160  3 * m_a * alpha * logTerm / (4 * m_b * T);
161 
162  // Now we obtain dPdT_rho
163  NekDouble dPdT =
164  m_gasConstant / (1.0 / rho - m_b) + m_a * alpha / (denom * 2 * T);
165 
166  // The result is dPde_rho = dPdT_rho / cv
167  return dPdT / cv;
168 }
NekDouble Alpha(const NekDouble &T)
double NekDouble
NekDouble GetTemperature(const NekDouble &rho, const NekDouble &e)
Calculate the temperature.
NekDouble LogTerm(const NekDouble &rho)

◆ v_GetDPDrho_e()

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

Implements Nektar::EquationOfState.

Definition at line 127 of file RedlichKwongEoS.cpp.

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

129 {
130  NekDouble T = GetTemperature(rho, e);
131  NekDouble alpha = Alpha(T);
132  NekDouble dPde = GetDPDe_rho(rho, e);
133 
134  // Calculate dPdrho_T
135  NekDouble dPdrho_T =
136  m_gasConstant * T / ((1.0 - m_b * rho) * (1.0 - m_b * rho)) -
137  m_a * alpha * rho * (m_b * rho + 2) /
138  ((1 + m_b * rho) * (1 + m_b * rho));
139 
140  // Calculate dedrho_T
141  NekDouble dedrho_T = -3 * m_a * alpha / (2 * (1 + m_b * rho));
142 
143  // The result is dPdrho_e = dPdrho_T - dPde_rho * dedrho_T
144  return dPdrho_T - dPde * dedrho_T;
145 }
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::RedlichKwongEoS::v_GetEFromRhoP ( const NekDouble rho,
const NekDouble p 
)
protectedvirtual

Implements Nektar::EquationOfState.

Definition at line 170 of file RedlichKwongEoS.cpp.

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

172 {
173  NekDouble logTerm = LogTerm(rho);
174  // First calculate the temperature, which can be expressed as
175  // (T^1/2)^3 + A* T^1/2 + B = 0
176  NekDouble A, B;
177 
178  A = -p * (1.0 / rho - m_b) / m_gasConstant;
179  B = -m_a * sqrt(m_Tc) * (1.0 / rho - m_b) /
180  (1.0 / (rho * rho) + m_b / rho) / m_gasConstant;
181 
182  // Use ideal gas solution as starting guess for iteration
183  NekDouble sqrtT = sqrt(p / (rho * (m_gamma - 1)));
184  // Newton-Raphson iteration to find T^(1/2)
185  NekDouble tol = 1e-6;
186  NekDouble maxIter = 100;
187  NekDouble residual = 1;
188  NekDouble f, df;
189  unsigned int cnt = 0;
190  while (abs(residual) > tol && cnt < maxIter)
191  {
192  f = sqrtT * sqrtT * sqrtT + A * sqrtT + B;
193  df = 3 * sqrtT * sqrtT + A;
194  residual = f / df;
195  sqrtT -= residual;
196  ++cnt;
197  }
198  if (cnt == maxIter)
199  {
200  cout << "Newton-Raphson in RedlichKwongEoS::v_GetEFromRhoP did not "
201  "converge in "
202  << maxIter << " iterations (residual = " << residual << ")"
203  << endl;
204  }
205 
206  // Calculate T
207  NekDouble T = sqrtT * sqrtT;
208 
209  // Calculate internal energy
210  return m_gasConstant * T / (m_gamma - 1) -
211  3 * m_a * Alpha(T) / (2 * m_b) * logTerm;
212 }
NekDouble Alpha(const NekDouble &T)
double NekDouble
NekDouble LogTerm(const NekDouble &rho)

◆ v_GetEntropy()

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

Implements Nektar::EquationOfState.

Definition at line 112 of file RedlichKwongEoS.cpp.

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

114 {
115  NekDouble T = GetTemperature(rho, e);
116  NekDouble logTerm = LogTerm(rho);
117  // Entropy for an ideal gas
118  NekDouble sIg =
119  m_gasConstant / (m_gamma - 1) * log(T) - m_gasConstant * log(rho);
120 
121  NekDouble deltaS = m_gasConstant * log(1 - m_b * rho);
122  deltaS -= m_a * Alpha(T) * logTerm / (2 * m_b * T);
123 
124  return sIg + deltaS;
125 }
NekDouble Alpha(const NekDouble &T)
double NekDouble
NekDouble GetTemperature(const NekDouble &rho, const NekDouble &e)
Calculate the temperature.
NekDouble LogTerm(const NekDouble &rho)

◆ v_GetPressure()

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

Implements Nektar::EquationOfState.

Definition at line 101 of file RedlichKwongEoS.cpp.

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

103 {
104  NekDouble T = GetTemperature(rho, e);
105 
106  NekDouble p = m_gasConstant * T / (1.0 / rho - m_b) -
107  m_a * Alpha(T) / (1.0 / (rho * rho) + m_b / rho);
108 
109  return p;
110 }
NekDouble Alpha(const NekDouble &T)
double NekDouble
NekDouble GetTemperature(const NekDouble &rho, const NekDouble &e)
Calculate the temperature.

◆ v_GetRhoFromPT()

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

Implements Nektar::EquationOfState.

Definition at line 214 of file RedlichKwongEoS.cpp.

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

216 {
217  // First solve for the compressibility factor Z using the cubic equation
218  // Z^3 + k1 * Z^2 + k2 * Z + k3 = 0
219  // for RedlichKwong:
220  // k1 = -1.0, k2 = A - B - B^2, k3 = -AB
221  // where A = a*alpha(T)*P/(RT)^2, B = bP/(RT)
222  NekDouble A = m_a * Alpha(T) * p / (m_gasConstant * m_gasConstant * T * T);
223  NekDouble B = m_b * p / (m_gasConstant * T);
224 
225  NekDouble k1 = -1.0;
226  NekDouble k2 = A - B - B * B;
227  NekDouble k3 = -A * B;
228 
229  // Use ideal gas (Z=1) as starting guess for iteration
230  NekDouble Z = 1.0;
231  // Newton-Raphson iteration to find Z
232  NekDouble tol = 1e-6;
233  NekDouble maxIter = 100;
234  NekDouble residual = 1;
235  NekDouble f, df;
236  unsigned int cnt = 0;
237  while (abs(residual) > tol && cnt < maxIter)
238  {
239  f = Z * Z * Z + k1 * Z * Z + k2 * Z + k3;
240  df = 3 * Z * Z + 2 * k1 * Z + k2;
241  residual = f / df;
242  Z -= residual;
243  ++cnt;
244  }
245  if (cnt == maxIter)
246  {
247  cout << "Newton-Raphson in RedlichKwongEoS::v_GetRhoFromPT did not "
248  "converge in "
249  << maxIter << " iterations (residual = " << residual << ")"
250  << endl;
251  }
252 
253  // Now calculate rho = p/(ZRT)
254  return p / (Z * m_gasConstant * T);
255 }
NekDouble Alpha(const NekDouble &T)
double NekDouble

◆ v_GetTemperature()

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

Implements Nektar::EquationOfState.

Definition at line 58 of file RedlichKwongEoS.cpp.

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

60 {
61  // First we need to evaluate the log term
62  // ln[1 + b*rho]
63  NekDouble logTerm = LogTerm(rho);
64 
65  // The temperature can be expressed as an equation in the form
66  // (T^1/2)^3 + A* T^1/2 + B = 0, which we solve iteratively
67  NekDouble A, B;
68 
69  A = -e * (m_gamma - 1) / m_gasConstant;
70  B = -3.0 * m_a / (2.0 * m_b * m_gasConstant) * (m_gamma - 1) * sqrt(m_Tc) *
71  logTerm;
72 
73  // Use ideal gas solution as starting guess for iteration
74  NekDouble sqrtT = sqrt(e * (m_gamma - 1) / m_gasConstant);
75  // Newton-Raphson iteration to find T^(1/2)
76  NekDouble tol = 1e-6;
77  NekDouble maxIter = 100;
78  NekDouble residual = 1;
79  NekDouble f, df;
80  unsigned int cnt = 0;
81  while (abs(residual) > tol && cnt < maxIter)
82  {
83  f = sqrtT * sqrtT * sqrtT + A * sqrtT + B;
84  df = 3 * sqrtT * sqrtT + A;
85  residual = f / df;
86  sqrtT -= residual;
87  ++cnt;
88  }
89  if (cnt == maxIter)
90  {
91  cout << "Newton-Raphson in RedlichKwongEoS::v_GetTemperature did not "
92  "converge in "
93  << maxIter << " iterations (residual = " << residual << ")"
94  << endl;
95  }
96 
97  // Calculate the temperature
98  return sqrtT * sqrtT;
99 }
double NekDouble
NekDouble LogTerm(const NekDouble &rho)

Friends And Related Function Documentation

◆ MemoryManager< RedlichKwongEoS >

friend class MemoryManager< RedlichKwongEoS >
friend

Definition at line 52 of file RedlichKwongEoS.h.

Member Data Documentation

◆ className

std::string Nektar::RedlichKwongEoS::className
static
Initial value:
=
"RedlichKwong", RedlichKwongEoS::create,
"Redlich-Kwong equation of state.")

Name of the class.

Definition at line 64 of file RedlichKwongEoS.h.

◆ m_a

NekDouble Nektar::RedlichKwongEoS::m_a
protected

◆ m_b

NekDouble Nektar::RedlichKwongEoS::m_b
protected

◆ m_Pc

NekDouble Nektar::RedlichKwongEoS::m_Pc
protected

Definition at line 70 of file RedlichKwongEoS.h.

Referenced by RedlichKwongEoS().

◆ m_Tc

NekDouble Nektar::RedlichKwongEoS::m_Tc
protected

Definition at line 69 of file RedlichKwongEoS.h.

Referenced by Alpha(), RedlichKwongEoS(), v_GetEFromRhoP(), and v_GetTemperature().