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

NekDouble GetTemperature (const NekDouble &rho, const NekDouble &e) final
 Calculate the temperature. More...
 
vec_t GetTemperature (const vec_t &rho, const vec_t &e) final
 
NekDouble GetPressure (const NekDouble &rho, const NekDouble &e) final
 Calculate the pressure. More...
 
vec_t GetPressure (const vec_t &rho, const vec_t &e) final
 
NekDouble v_GetEntropy (const NekDouble &rho, const NekDouble &e) final
 
NekDouble v_GetDPDrho_e (const NekDouble &rho, const NekDouble &e) final
 
NekDouble v_GetDPDe_rho (const NekDouble &rho, const NekDouble &e) final
 
NekDouble v_GetEFromRhoP (const NekDouble &rho, const NekDouble &p) final
 
NekDouble v_GetRhoFromPT (const NekDouble &rho, const NekDouble &p) final
 
- Protected Member Functions inherited from Nektar::EquationOfState
 EquationOfState (const LibUtilities::SessionReaderSharedPtr &pSession)
 Constructor. More...
 
 EquationOfState (const NekDouble &gamma, const NekDouble &gasConstant)
 Programmatic 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
 
NekDouble m_gammaMone
 
NekDouble m_gammaMoneOgasConst
 

Private Member Functions

 RedlichKwongEoS (const LibUtilities::SessionReaderSharedPtr &pSession)
 
 ~RedlichKwongEoS (void)
 
template<class T , typename = typename std::enable_if < std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value >::type>
Alpha (const T &temp)
 
template<class T , typename = typename std::enable_if < std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value >::type>
LogTerm (const T &rho)
 
template<class T , typename = typename std::enable_if < std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value >::type>
GetTemperatureKernel (const T &rho, const T &e)
 
template<class T , typename = typename std::enable_if < std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value >::type>
GetPressureKernel (const T &rho, const T &e)
 

Friends

class MemoryManager< RedlichKwongEoS >
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::EquationOfState
virtual ~EquationOfState ()
 
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 53 of file RedlichKwongEoS.h.

Constructor & Destructor Documentation

◆ RedlichKwongEoS()

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

Definition at line 47 of file RedlichKwongEoS.cpp.

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.

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

◆ ~RedlichKwongEoS()

Nektar::RedlichKwongEoS::~RedlichKwongEoS ( void  )
inlineprivate

Definition at line 97 of file RedlichKwongEoS.h.

97 {};

Member Function Documentation

◆ Alpha()

template<class T , typename = typename std::enable_if < std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value >::type>
T Nektar::RedlichKwongEoS::Alpha ( const T &  temp)
inlineprivate

Definition at line 106 of file RedlichKwongEoS.h.

107  {
108  return 1.0 / sqrt(temp / m_Tc);
109  }
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267

References tinysimd::sqrt().

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

◆ create()

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

Creates an instance of this class.

Definition at line 59 of file RedlichKwongEoS.h.

61  {
64  return p;
65  }
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.

References CellMLToNektar.cellml_metadata::p.

◆ GetPressure() [1/2]

NekDouble Nektar::RedlichKwongEoS::GetPressure ( const NekDouble rho,
const NekDouble e 
)
finalprotectedvirtual

Calculate the pressure.

Implements Nektar::EquationOfState.

Definition at line 70 of file RedlichKwongEoS.cpp.

72 {
73  return GetPressureKernel(rho, e);
74 }
T GetPressureKernel(const T &rho, const T &e)

References GetPressureKernel().

◆ GetPressure() [2/2]

vec_t Nektar::RedlichKwongEoS::GetPressure ( const vec_t rho,
const vec_t e 
)
finalprotectedvirtual

Implements Nektar::EquationOfState.

Definition at line 76 of file RedlichKwongEoS.cpp.

78 {
79  return GetPressureKernel(rho, e);
80 }

References GetPressureKernel().

◆ GetPressureKernel()

template<class T , typename = typename std::enable_if < std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value >::type>
T Nektar::RedlichKwongEoS::GetPressureKernel ( const T &  rho,
const T &  e 
)
inlineprivate

Definition at line 174 of file RedlichKwongEoS.h.

175  {
176  T temp = GetTemperatureKernel(rho, e);
177  T oneOrho = 1.0 / rho;
178  T p = m_gasConstant * temp / (oneOrho - m_b) - m_a * Alpha(temp) /
179  (oneOrho * (oneOrho + m_b));
180  return p;
181  }
T GetTemperatureKernel(const T &rho, const T &e)
T Alpha(const T &temp)

References CellMLToNektar.cellml_metadata::p.

Referenced by GetPressure().

◆ GetTemperature() [1/2]

NekDouble Nektar::RedlichKwongEoS::GetTemperature ( const NekDouble rho,
const NekDouble e 
)
finalprotectedvirtual

Calculate the temperature.

Implements Nektar::EquationOfState.

Definition at line 58 of file RedlichKwongEoS.cpp.

60 {
61  return GetTemperatureKernel(rho, e);
62 }

References GetTemperatureKernel().

Referenced by v_GetDPDe_rho(), v_GetDPDrho_e(), and v_GetEntropy().

◆ GetTemperature() [2/2]

vec_t Nektar::RedlichKwongEoS::GetTemperature ( const vec_t rho,
const vec_t e 
)
finalprotectedvirtual

Implements Nektar::EquationOfState.

Definition at line 64 of file RedlichKwongEoS.cpp.

66 {
67  return GetTemperatureKernel(rho, e);
68 }

References GetTemperatureKernel().

◆ GetTemperatureKernel()

template<class T , typename = typename std::enable_if < std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value >::type>
T Nektar::RedlichKwongEoS::GetTemperatureKernel ( const T &  rho,
const T &  e 
)
inlineprivate

Definition at line 129 of file RedlichKwongEoS.h.

130  {
131  // First we need to evaluate the log term
132  // ln[1 + b*rho]
133  T logTerm = LogTerm(rho);
134 
135  // The temperature can be expressed as an equation in the form
136  // (T^1/2)^3 + A* T^1/2 + B = 0, which we solve iteratively
137  T A = e * (1.0 - m_gamma) / m_gasConstant;
138  T B = -3.0 * m_a / (2.0 * m_b * m_gasConstant) * (m_gamma - 1) * sqrt(m_Tc)
139  * logTerm;
140 
141  // Use ideal gas solution as starting guess for iteration
142  T sqrtT = sqrt(e * (m_gamma - 1) / m_gasConstant);
143  // Newton-Raphson iteration to find T^(1/2)
144  T tol = 1e-6;
145  T residual = 1;
146  unsigned int maxIter = 100;
147  unsigned int cnt = 0;
148  while (abs(residual) > tol && cnt < maxIter)
149  {
150  T f = sqrtT * sqrtT * sqrtT + A * sqrtT + B;
151  T df = 3 * sqrtT * sqrtT + A;
152  residual = f / df;
153  sqrtT -= residual;
154  ++cnt;
155  }
156  if (cnt == maxIter)
157  {
158  std::cout << "Newton-Raphson in RedlichKwongEoS::v_GetTemperature did not "
159  "converge in "
160  << maxIter << " iterations (residual = " << residual << ")"
161  << std::endl;
162  }
163 
164  // Calculate the temperature
165  return sqrtT * sqrtT;
166  }
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:272

References tinysimd::abs(), and tinysimd::sqrt().

Referenced by GetTemperature().

◆ LogTerm()

template<class T , typename = typename std::enable_if < std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value >::type>
T Nektar::RedlichKwongEoS::LogTerm ( const T &  rho)
inlineprivate

Definition at line 118 of file RedlichKwongEoS.h.

119  {
120  return log(1.0 + m_b * rho);
121  }
scalarT< T > log(scalarT< T > in)
Definition: scalar.hpp:278

References tinysimd::log().

Referenced by v_GetDPDe_rho(), v_GetEFromRhoP(), and v_GetEntropy().

◆ v_GetDPDe_rho()

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

Implements Nektar::EquationOfState.

Definition at line 117 of file RedlichKwongEoS.cpp.

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 }
NekDouble GetTemperature(const NekDouble &rho, const NekDouble &e) final
Calculate the temperature.
double NekDouble

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

◆ v_GetDPDrho_e()

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

Implements Nektar::EquationOfState.

Definition at line 97 of file RedlichKwongEoS.cpp.

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 }
NekDouble GetDPDe_rho(const NekDouble &rho, const NekDouble &e)
Calculate the partial derivative of P(rho,e) with respect to e.

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

◆ v_GetEFromRhoP()

NekDouble Nektar::RedlichKwongEoS::v_GetEFromRhoP ( const NekDouble rho,
const NekDouble p 
)
finalprotectedvirtual

Implements Nektar::EquationOfState.

Definition at line 140 of file RedlichKwongEoS.cpp.

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 }

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

◆ v_GetEntropy()

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

Implements Nektar::EquationOfState.

Definition at line 82 of file RedlichKwongEoS.cpp.

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 }

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

◆ v_GetRhoFromPT()

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

Implements Nektar::EquationOfState.

Definition at line 184 of file RedlichKwongEoS.cpp.

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 }

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

Friends And Related Function Documentation

◆ MemoryManager< RedlichKwongEoS >

friend class MemoryManager< RedlichKwongEoS >
friend

Definition at line 1 of file RedlichKwongEoS.h.

Member Data Documentation

◆ className

std::string Nektar::RedlichKwongEoS::className
static
Initial value:
=
"RedlichKwong", RedlichKwongEoS::create,
"Redlich-Kwong equation of state.")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
static EquationOfStateSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
Creates an instance of this class.
EquationOfStateFactory & GetEquationOfStateFactory()
Declaration of the equation of state factory singleton.

Name of the class.

Definition at line 68 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 74 of file RedlichKwongEoS.h.

Referenced by RedlichKwongEoS().

◆ m_Tc

NekDouble Nektar::RedlichKwongEoS::m_Tc
protected

Definition at line 73 of file RedlichKwongEoS.h.

Referenced by RedlichKwongEoS(), and v_GetEFromRhoP().