OGS
PengRobinson.cpp
Go to the documentation of this file.
1
11#include "PengRobinson.h"
12
13#include <algorithm>
14#include <boost/math/special_functions/pow.hpp>
15#include <cmath>
16
18
19namespace MaterialPropertyLib
20{
21
23 MaterialPropertyLib::VariableArray const& variable_array,
24 ParameterLib::SpatialPosition const& /*pos*/, double const /*t*/,
25 double const /*dt*/) const
26{
27 const double gas_constant = MaterialLib::PhysicalConstant::IdealGasConstant;
28 const double pressure = variable_array.gas_phase_pressure;
29 const double temperature = variable_array.temperature;
30 const double molar_mass = variable_array.molar_mass;
31
32 const double Tr = temperature / Tc_; // reduced Temperature
33
34 const double kappa = 0.37464 + 1.5422 * omega_ - 0.26992 * omega_ * omega_;
35
36 const double alpha = boost::math::pow<2>(1 + kappa * (1 - std::sqrt(Tr)));
37
38 // EOS in the form: 0 = rho^3 + z1*rho^2 + z2*rho + z3
39 const double denominator =
40 b_ *
41 (pressure * b_ * b_ + b_ * gas_constant * temperature - a_ * alpha);
42
43 const double z1 =
44 (molar_mass * a_ * alpha - 3 * molar_mass * b_ * b_ * pressure -
45 2 * molar_mass * gas_constant * temperature * b_) /
46 denominator;
47 const auto z2 = (molar_mass * molar_mass *
48 (b_ * pressure - gas_constant * temperature)) /
49 denominator;
50 const auto z3 =
51 (molar_mass * molar_mass * molar_mass * pressure) / denominator;
52
53 MathLib::CubicSolver cubic_solver_(1., z1, z2, z3);
54 return cubic_solver_.smallestPositiveRealRoot();
55}
56
58 MaterialPropertyLib::VariableArray const& variable_array,
59 Variable const variable, ParameterLib::SpatialPosition const& pos,
60 double const t, double const dt) const
61{
62 if (variable == Variable::temperature)
63 {
64 const double temperature = variable_array.temperature;
65 const double epsilon = 1.e-6;
66
67 MaterialPropertyLib::VariableArray perturbed = variable_array;
68 // Increase temperature by +epsilon
69 perturbed.temperature = temperature + epsilon;
70 const double rho_plus = std::get<double>(value(perturbed, pos, t, dt));
71
72 // Decrease temperature by -epsilon
73 perturbed.temperature = temperature - epsilon;
74 const double rho_minus = std::get<double>(value(perturbed, pos, t, dt));
75
76 // Calculate the central difference quotient
77 return (rho_plus - rho_minus) / (2 * epsilon);
78 }
79
80 if (variable == Variable::gas_phase_pressure)
81 {
82 const double pressure = variable_array.gas_phase_pressure;
83 const double epsilon = 1.e-3;
84
85 MaterialPropertyLib::VariableArray perturbed = variable_array;
86 // Increase pressure by +epsilon
87 perturbed.gas_phase_pressure = pressure + epsilon;
88 const double rho_plus = std::get<double>(value(perturbed, pos, t, dt));
89
90 // Decrease pressure by -epsilon
91 perturbed.gas_phase_pressure = pressure - epsilon;
92 const double rho_minus = std::get<double>(value(perturbed, pos, t, dt));
93
94 // Calculate the central difference quotient
95 return (rho_plus - rho_minus) / (2 * epsilon);
96 }
97
99 "PengRobinson::dValue is implemented for derivatives with respect to "
100 "gas phase pressure or temperature only.");
101
102 return 0.;
103}
104
105} // namespace MaterialPropertyLib
#define OGS_FATAL(...)
Definition Error.h:26
PropertyDataType dValue(MaterialPropertyLib::VariableArray const &variable_array, Variable const variable, ParameterLib::SpatialPosition const &pos, double const t, double const dt) const override
const double Tc_
Critical temperature.
const double omega_
Acentric factor.
double b_
Parameter b (covolume)
double a_
Parameter a (cohesion pressure)
virtual PropertyDataType value() const
Definition Property.cpp:76
double smallestPositiveRealRoot()
Definition CubicRoots.h:55
std::variant< double, Eigen::Matrix< double, 2, 1 >, Eigen::Matrix< double, 3, 1 >, Eigen::Matrix< double, 2, 2 >, Eigen::Matrix< double, 3, 3 >, Eigen::Matrix< double, 4, 1 >, Eigen::Matrix< double, 6, 1 >, Eigen::MatrixXd > PropertyDataType
Definition Property.h:31