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 <cmath>
38 #include "EquationOfState.h"
39 
41 
42 using namespace std;
43 
44 namespace 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 {
55 public:
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 
70 protected:
75 
76  NekDouble GetTemperature(const NekDouble& rho, const NekDouble& e) final;
77 
78  vec_t GetTemperature(const vec_t& rho, const vec_t& e) final;
79 
80  NekDouble GetPressure(const NekDouble& rho, const NekDouble& e) final;
81 
82  vec_t 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 
94 private:
96 
98 
99  // Alpha term of Redlich-Kwong EoS ( 1.0/sqrt(Tr))
100  template <class T, typename = typename std::enable_if
101  <
102  std::is_floating_point<T>::value ||
104  >::type
105  >
106  inline T Alpha(const T& temp)
107  {
108  return 1.0 / sqrt(temp / m_Tc);
109  }
110 
111  // Log term term of Peng-Robinson EoS
112  template <class T, typename = typename std::enable_if
113  <
114  std::is_floating_point<T>::value ||
116  >::type
117  >
118  inline T LogTerm(const T& rho)
119  {
120  return log(1.0 + m_b * rho);
121  }
122 
123  template <class T, typename = typename std::enable_if
124  <
125  std::is_floating_point<T>::value ||
127  >::type
128  >
129  inline T GetTemperatureKernel(const T& rho, const T& e)
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  }
167 
168  template <class T, typename = typename std::enable_if
169  <
170  std::is_floating_point<T>::value ||
172  >::type
173  >
174  inline T GetPressureKernel(const T& rho, const T& e)
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  }
182 
183 
184 
185 
186 };
187 }
188 #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...
Redlich-Kwong equation of state: p = RT/(1/rho - b) - a/( sqrt(T / Tc) * (1/rho^2 + b/rho) with a = 0...
T GetTemperatureKernel(const T &rho, const T &e)
static std::string className
Name of the class.
static EquationOfStateSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
Creates an instance of this class.
T Alpha(const T &temp)
T GetPressureKernel(const T &rho, const T &e)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
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:278
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:272
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267