Go to the documentation of this file.
26 #ifndef OPM_IAPWS_REGION2_HPP
27 #define OPM_IAPWS_REGION2_HPP
49 template < class Scalar>
60 template < class Evaluation>
61 static bool isValid( const Evaluation& temperature, const Evaluation& pressure)
64 temperature <= 623.15 && pressure <= 100e6;
80 template < class Evaluation>
81 static Evaluation tau( const Evaluation& temperature)
82 { return 540.0 / temperature; }
90 template < class Evaluation>
91 static Evaluation dtau_dT( const Evaluation& temperature)
92 { return - 540.0 / (temperature*temperature); }
99 template < class Evaluation>
100 static Evaluation pi( const Evaluation& pressure)
101 { return pressure / 1e6; }
109 template < class Evaluation>
111 { return 1.0 / 1e6; }
119 template < class Evaluation>
120 static Evaluation dp_dpi( const Evaluation& )
134 template < class Evaluation>
135 static Evaluation gamma( const Evaluation& temperature, const Evaluation& pressure)
139 const Evaluation& tau_ = tau(temperature);
140 const Evaluation& pi_ = pi(pressure);
145 result = Toolbox::ln(pi_);
146 for ( int i = 0; i < 9; ++i)
150 for ( int i = 0; i < 43; ++i)
170 template < class Evaluation>
171 static Evaluation dgamma_dtau( const Evaluation& temperature, const Evaluation& pressure)
175 const Evaluation& tau_ = tau(temperature);
176 const Evaluation& pi_ = pi(pressure);
179 Evaluation result = Toolbox::createConstant(0.0);
180 for ( int i = 0; i < 9; i++) {
188 for ( int i = 0; i < 43; i++) {
193 Toolbox::pow(tau_ - 0.5, static_cast<Scalar>(J_r(i) - 1));
211 template < class Evaluation>
212 static Evaluation dgamma_dpi( const Evaluation& temperature, const Evaluation& pressure)
216 const Evaluation& tau_ = tau(temperature);
217 const Evaluation& pi_ = pi(pressure);
220 Evaluation result = 1/pi_;
223 for ( int i = 0; i < 43; i++) {
246 template < class Evaluation>
247 static Evaluation ddgamma_dtaudpi( const Evaluation& temperature, const Evaluation& pressure)
251 const Evaluation& tau_ = tau(temperature);
252 const Evaluation& pi_ = pi(pressure);
255 Evaluation result = Toolbox::createConstant(0.0);
258 for ( int i = 0; i < 43; i++) {
264 Toolbox::pow(tau_ - 0.5, static_cast<Scalar>(J_r(i) - 1));
282 template < class Evaluation>
283 static Evaluation ddgamma_ddpi( const Evaluation& temperature, const Evaluation& pressure)
287 const Evaluation& tau_ = tau(temperature);
288 const Evaluation& pi_ = pi(pressure);
291 Evaluation result = -1/(pi_*pi_);
294 for ( int i = 0; i < 43; i++) {
318 template < class Evaluation>
319 static Evaluation ddgamma_ddtau( const Evaluation& temperature, const Evaluation& pressure)
323 const Evaluation& tau_ = tau(temperature);
324 const Evaluation& pi_ = pi(pressure);
327 Evaluation result = Toolbox::createConstant(0.0);
328 for ( int i = 0; i < 9; i++) {
337 for ( int i = 0; i < 43; i++) {
343 Toolbox::pow(tau_ - 0.5, static_cast<Scalar>(J_r(i) - 2));
351 static Scalar n_g( int i)
353 static const Scalar n[9] = {
354 -0.96927686500217e1, 0.10086655968018e2, -0.56087911283020e-2,
355 0.71452738081455e-1, -0.40710498223928, 0.14240819171444e1,
356 -0.43839511319450e1, -0.28408632460772, 0.21268463753307e-1
361 static Scalar n_r( int i)
363 static const Scalar n[43] = {
364 -0.17731742473213e-2, -0.17834862292358e-1, -0.45996013696365e-1,
365 -0.57581259083432e-1, -0.50325278727930e-1, -0.33032641670203e-4,
366 -0.18948987516315e-3, -0.39392777243355e-2, -0.43797295650573e-1,
367 -0.26674547914087e-4, 0.20481737692309e-7, 0.43870667284435e-6,
368 -0.32277677238570e-4, -0.15033924542148e-2, -0.40668253562649e-1,
369 -0.78847309559367e-9, 0.12790717852285e-7, 0.48225372718507e-6,
370 0.22922076337661e-5, -0.16714766451061e-10, -0.21171472321355e-2,
371 -0.23895741934104e2, -0.59059564324270e-17, -0.12621808899101e-5,
372 -0.38946842435739e-1, 0.11256211360459e-10, -0.82311340897998e1,
373 0.19809712802088e-7, 0.10406965210174e-18, -0.10234747095929e-12,
374 -0.10018179379511e-8, -0.80882908646985e-10, 0.10693031879409,
375 -0.33662250574171, 0.89185845355421e-24, 0.30629316876232e-12,
376 -0.42002467698208e-5, -0.59056029685639e-25, 0.37826947613457e-5,
377 -0.12768608934681e-14, 0.73087610595061e-28, 0.55414715350778e-16,
383 static Scalar I_r( int i)
385 static const short int I[43] = {
405 static Scalar J_g( int i)
407 static const short int J[9] = {
415 static Scalar J_r( int i)
417 static const short int J[43] = {
static Evaluation gamma(const Evaluation &temperature, const Evaluation &pressure) The Gibbs free energy for IAPWS region 2 (i.e. sub-critical steam) (dimensionless). Definition: Region2.hpp:135
static Evaluation ddgamma_dtaudpi(const Evaluation &temperature, const Evaluation &pressure) The partial derivative of the Gibbs free energy to the normalized pressure and to the normalized temp... Definition: Region2.hpp:247
Definition: Air_Mesitylene.hpp:31
static Evaluation ddgamma_ddpi(const Evaluation &temperature, const Evaluation &pressure) The second partial derivative of the Gibbs free energy to the normalized pressure for IAPWS region 2 ... Definition: Region2.hpp:283
static Evaluation dtau_dT(const Evaluation &temperature) Returns the derivative of the reduced temperature to the temperature for IAPWS region 2... Definition: Region2.hpp:91
static Evaluation tau(const Evaluation &temperature) Returns the reduced temperature (dimensionless) for IAPWS region 2. Definition: Region2.hpp:81
static Scalar dpi_dp(const Evaluation &) Returns the derivative of the reduced pressure to the pressure for IAPWS region 2 in ... Definition: Region2.hpp:110
static bool isValid(const Evaluation &temperature, const Evaluation &pressure) Returns true if IAPWS region 2 applies for a (temperature, pressure) pair. Definition: Region2.hpp:61
static Evaluation dp_dpi(const Evaluation &) Returns the derivative of the pressure to the reduced pressure for IAPWS region 2 (dimensionless)... Definition: Region2.hpp:120
static Evaluation dgamma_dpi(const Evaluation &temperature, const Evaluation &pressure) The partial derivative of the Gibbs free energy to the normalized pressure for IAPWS region 2 (i... Definition: Region2.hpp:212
static Evaluation dgamma_dtau(const Evaluation &temperature, const Evaluation &pressure) The partial derivative of the Gibbs free energy to the normalized temperature for IAPWS region 2 (i... Definition: Region2.hpp:171
Evaluation< Scalar, VarSetTag, numVars > pow(const Evaluation< Scalar, VarSetTag, numVars > &base, Scalar exp) Definition: Math.hpp:312
static Evaluation pi(const Evaluation &pressure) Returns the reduced pressure (dimensionless) for IAPWS region 2. Definition: Region2.hpp:100
Implements the equations for region 2 of the IAPWS '97 formulation. Definition: Region2.hpp:50
static Evaluation ddgamma_ddtau(const Evaluation &temperature, const Evaluation &pressure) The second partial derivative of the Gibbs free energy to the normalized temperature for IAPWS region... Definition: Region2.hpp:319
|