45 "Peng-Robinson equation of state.");
51 pSession->LoadParameter(
"Tcrit",
m_Tc);
52 pSession->LoadParameter(
"Pcrit",
m_Pc);
53 pSession->LoadParameter(
"AcentricFactor",
m_omega);
96 deltaS +=
m_a * sqrtA *
m_fw * logTerm * (sqrtTr / T) / (
m_b *
sqrt(8));
115 2 *
m_a * alpha * rho * (1 +
m_b * rho) /
116 ((denom * rho * rho) * (denom * rho * rho));
120 -
m_a *
sqrt(alpha) * (1.0 +
m_fw) / (denom * rho * rho);
123 return dPdrho_T - dPde * dedrho_T;
175 (alpha +
sqrt(alpha) *
m_fw * sqrtTr);
200 unsigned int cnt = 0;
201 while ((fabs(residual) > tol) && cnt < maxIter)
203 f = Z * Z * Z + k1 * Z * Z + k2 * Z + k3;
204 df = 3 * Z * Z + 2 * k1 * Z + k2;
211 cout <<
"Newton-Raphson in PengRobinsonEoS::v_GetRhoFromPT did not "
213 << maxIter <<
" iterations (residual = " << residual <<
")"
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.
PengRobinsonEoS(const LibUtilities::SessionReaderSharedPtr &pSession)
static std::string className
Name of the class.
static EquationOfStateSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
Creates an instance of this class.
T GetTemperatureKernel(const T &rho, const T &e)
T GetPressureKernel(const T &rho, const T &e)
NekDouble v_GetPressure(const NekDouble &rho, const NekDouble &e) final
NekDouble v_GetEFromRhoP(const NekDouble &rho, const NekDouble &p) final
NekDouble v_GetDPDe_rho(const NekDouble &rho, const NekDouble &e) final
NekDouble v_GetTemperature(const NekDouble &rho, const NekDouble &e) final
NekDouble v_GetRhoFromPT(const NekDouble &rho, const NekDouble &p) final
NekDouble v_GetDPDrho_e(const NekDouble &rho, const NekDouble &e) final
NekDouble v_GetEntropy(const NekDouble &rho, const NekDouble &e) final
std::shared_ptr< SessionReader > SessionReaderSharedPtr
EquationOfStateFactory & GetEquationOfStateFactory()
Declaration of the equation of state factory singleton.
tinysimd::simd< NekDouble > vec_t
scalarT< T > log(scalarT< T > in)
scalarT< T > sqrt(scalarT< T > in)