35 #ifndef NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_MISC_REDLICHKWONGEOS
36 #define NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_MISC_REDLICHKWONGEOS
100 template <
class T,
typename =
typename std::enable_if
102 std::is_floating_point<T>::value ||
108 return 1.0 /
sqrt(temp / m_Tc);
112 template <
class T,
typename =
typename std::enable_if
114 std::is_floating_point<T>::value ||
120 return log(1.0 + m_b * rho);
123 template <
class T,
typename =
typename std::enable_if
125 std::is_floating_point<T>::value ||
133 T logTerm = LogTerm(rho);
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)
142 T sqrtT =
sqrt(e * (m_gamma - 1) / m_gasConstant);
146 unsigned int maxIter = 100;
147 unsigned int cnt = 0;
148 while (
abs(residual) > tol && cnt < maxIter)
150 T f = sqrtT * sqrtT * sqrtT +
A * sqrtT + B;
151 T df = 3 * sqrtT * sqrtT +
A;
158 std::cout <<
"Newton-Raphson in RedlichKwongEoS::v_GetTemperature did not "
160 << maxIter <<
" iterations (residual = " << residual <<
")"
165 return sqrtT * sqrtT;
168 template <
class T,
typename =
typename std::enable_if
170 std::is_floating_point<T>::value ||
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));
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 GetPressureKernel(const T &rho, const T &e)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
The above copyright notice and this permission notice shall be included.
std::shared_ptr< EquationOfState > EquationOfStateSharedPtr
A shared pointer to an equation of state object.
tinysimd::simd< NekDouble > vec_t
scalarT< T > log(scalarT< T > in)
scalarT< T > abs(scalarT< T > in)
scalarT< T > sqrt(scalarT< T > in)