Nektar++
PengRobinsonEoS.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: PengRobinsonEoS.h
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#ifndef NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_MISC_PENGROBINSONEOS
36#define NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_MISC_PENGROBINSONEOS
37
38#include "EquationOfState.h"
39
40namespace Nektar
41{
42
43/**
44 * @brief Peng-Robinson equation of state:
45 * p = RT/(1/rho - b) - a*Alpha(T/Tc) / (1/rho^2 + 2*b/rho - b^2)
46 * with a = 0.45724 * (R*Tc)^2 / Pc
47 * b = 0.0778 * (R*Tc) / Pc
48 * Alpha(T/Tc) = [1 + fw * (1 - sqrt(T/ Tc))]^2
49 * fw = 0.37464 + 1.54226*omega - 0.2699*omega*omega
50 */
52{
53public:
54 friend class MemoryManager<PengRobinsonEoS>;
55
56 /// Creates an instance of this class
59 {
62 return p;
63 }
64
65 /// Name of the class
66 static std::string className;
67
68protected:
75
76 NekDouble v_GetTemperature(const NekDouble &rho, const NekDouble &e) final;
77
78 vec_t v_GetTemperature(const vec_t &rho, const vec_t &e) final;
79
80 NekDouble v_GetPressure(const NekDouble &rho, const NekDouble &e) final;
81
82 vec_t v_GetPressure(const vec_t &rho, const vec_t &e) final;
83
84 NekDouble v_GetEntropy(const NekDouble &rho, const NekDouble &e) final;
85
86 NekDouble v_GetDPDrho_e(const NekDouble &rho, const NekDouble &e) final;
87
88 NekDouble v_GetDPDe_rho(const NekDouble &rho, const NekDouble &e) final;
89
90 NekDouble v_GetEFromRhoP(const NekDouble &rho, const NekDouble &p) final;
91
92 NekDouble v_GetRhoFromPT(const NekDouble &rho, const NekDouble &p) final;
93
94private:
96
97 ~PengRobinsonEoS(void) override{};
98
99 // Alpha term of Peng-Robinson EoS
100 template <class T, typename = typename std::enable_if<
101 std::is_floating_point<T>::value ||
103 inline T Alpha(const T &temp)
104 {
105 T sqrtAlpha = 1.0 + m_fw * (1.0 - sqrt(temp / m_Tc));
106 return sqrtAlpha * sqrtAlpha;
107 }
108
109 // Log term term of Peng-Robinson EoS
110 template <class T, typename = typename std::enable_if<
111 std::is_floating_point<T>::value ||
113 inline T LogTerm(const T &rho)
114 {
115 return log((1.0 / rho + m_b - m_b * sqrt(2)) /
116 (1.0 / rho + m_b + m_b * sqrt(2)));
117 }
118
119 template <class T, typename = typename std::enable_if<
120 std::is_floating_point<T>::value ||
122 inline T GetTemperatureKernel(const T &rho, const T &e)
123 {
124 // First we need to evaluate the log term
125 // ln[(1/rho + b - b*sqrt(2)) / (1/rho + b + b*sqrt(2))]
126 NekDouble sqrt2 = sqrt(2.0);
127 T logTerm = LogTerm(rho);
128
129 // The temperature can be expressed as an equation in the form
130 // A * (T^1/2)^2 + B * T^1/2 + C = 0
131
132 NekDouble A = m_gasConstant / (m_gamma - 1.0);
133 NekDouble f1 = m_a / (m_b * 2.0 * sqrt2);
134 NekDouble f2 = (1.0 + m_fw);
135 T B = -f1 * logTerm / sqrt(m_Tc) * m_fw * f2;
136 T C = f1 * logTerm * f2 * f2 - e;
137
138 // Solve for T^1/2 (positive root)
139 T sqrtT = (sqrt(B * B - 4 * A * C) - B) / (2 * A);
140 // Calculate the temperature
141 return sqrtT * sqrtT;
142 }
143
144 template <class T, typename = typename std::enable_if<
145 std::is_floating_point<T>::value ||
147 inline T GetPressureKernel(const T &rho, const T &e)
148 {
149 T temp = GetTemperatureKernel(rho, e);
150 T oneOrho = 1.0 / rho;
151 T p = m_gasConstant * temp / (oneOrho - m_b) -
152 m_a * Alpha(temp) /
153 (oneOrho * oneOrho + 2.0 * m_b * oneOrho - m_b * m_b);
154
155 return p;
156 }
157};
158} // namespace Nektar
159
160#endif
Encapsulates equations of state allowing us to obtain thermodynamic properties: most relations are in...
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Peng-Robinson equation of state: p = RT/(1/rho - b) - a*Alpha(T/Tc) / (1/rho^2 + 2*b/rho - b^2) with ...
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
~PengRobinsonEoS(void) override
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
T Alpha(const T &temp)
NekDouble v_GetEntropy(const NekDouble &rho, const NekDouble &e) final
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< EquationOfState > EquationOfStateSharedPtr
A shared pointer to an equation of state object.
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