Nektar++
RedlichKwongEoS.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: RedlichKwongEoS.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: Redlich-Kwong equation of state
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#ifndef NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_MISC_REDLICHKWONGEOS
36#define NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_MISC_REDLICHKWONGEOS
37#include "EquationOfState.h"
38#include <cmath>
39
41
42using namespace std;
43
44namespace Nektar
45{
46
47/**
48 * @brief Redlich-Kwong equation of state:
49 * p = RT/(1/rho - b) - a/( sqrt(T / Tc) * (1/rho^2 + b/rho)
50 * with a = 0.42748 * (R*Tc)^2 / Pc
51 * b = 0.08664 * (R*Tc) / Pc
52 */
54{
55public:
56 friend class MemoryManager<RedlichKwongEoS>;
57
58 /// Creates an instance of this class
61 {
64 return p;
65 }
66
67 /// Name of the class
68 static std::string className;
69
70protected:
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 ~RedlichKwongEoS(void) override{};
98
99 // Alpha term of Redlich-Kwong EoS ( 1.0/sqrt(Tr))
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 return 1.0 / sqrt(temp / m_Tc);
106 }
107
108 // Log term term of Peng-Robinson EoS
109 template <class T, typename = typename std::enable_if<
110 std::is_floating_point<T>::value ||
112 inline T LogTerm(const T &rho)
113 {
114 return log(1.0 + m_b * rho);
115 }
116
117 template <class T, typename = typename std::enable_if<
118 std::is_floating_point<T>::value ||
120 inline T GetTemperatureKernel(const T &rho, const T &e)
121 {
122 // First we need to evaluate the log term
123 // ln[1 + b*rho]
124 T logTerm = LogTerm(rho);
125
126 // The temperature can be expressed as an equation in the form
127 // (T^1/2)^3 + A* T^1/2 + B = 0, which we solve iteratively
128 T A = e * (1.0 - m_gamma) / m_gasConstant;
129 T B = -3.0 * m_a / (2.0 * m_b * m_gasConstant) * (m_gamma - 1) *
130 sqrt(m_Tc) * logTerm;
131
132 // Use ideal gas solution as starting guess for iteration
133 T sqrtT = sqrt(e * (m_gamma - 1) / m_gasConstant);
134 // Newton-Raphson iteration to find T^(1/2)
135 T tol = 1e-6;
136 T residual = 1;
137 unsigned int maxIter = 100;
138 unsigned int cnt = 0;
139 while (abs(residual) > tol && cnt < maxIter)
140 {
141 T f = sqrtT * sqrtT * sqrtT + A * sqrtT + B;
142 T df = 3 * sqrtT * sqrtT + A;
143 residual = f / df;
144 sqrtT -= residual;
145 ++cnt;
146 }
147 if (cnt == maxIter)
148 {
149 std::cout << "Newton-Raphson in RedlichKwongEoS::v_GetTemperature "
150 "did not "
151 "converge in "
152 << maxIter << " iterations (residual = " << residual
153 << ")" << std::endl;
154 }
155
156 // Calculate the temperature
157 return sqrtT * sqrtT;
158 }
159
160 template <class T, typename = typename std::enable_if<
161 std::is_floating_point<T>::value ||
163 inline T GetPressureKernel(const T &rho, const T &e)
164 {
165 T temp = GetTemperatureKernel(rho, e);
166 T oneOrho = 1.0 / rho;
167 T p = m_gasConstant * temp / (oneOrho - m_b) -
168 m_a * Alpha(temp) / (oneOrho * (oneOrho + m_b));
169 return p;
170 }
171};
172} // namespace Nektar
173#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.
Redlich-Kwong equation of state: p = RT/(1/rho - b) - a/( sqrt(T / Tc) * (1/rho^2 + b/rho) with a = 0...
~RedlichKwongEoS(void) override
T GetTemperatureKernel(const T &rho, const T &e)
static std::string className
Name of the class.
NekDouble v_GetDPDrho_e(const NekDouble &rho, const NekDouble &e) final
RedlichKwongEoS(const LibUtilities::SessionReaderSharedPtr &pSession)
static EquationOfStateSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
Creates an instance of this class.
T Alpha(const T &temp)
NekDouble v_GetRhoFromPT(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_GetEntropy(const NekDouble &rho, const NekDouble &e) final
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
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 > abs(scalarT< T > in)
Definition: scalar.hpp:298
scalarT< T > log(scalarT< T > in)
Definition: scalar.hpp:303
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294