OGS
PengRobinson.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
10
11namespace MaterialPropertyLib
12{
52class PengRobinson final : public Property
53{
54public:
55 explicit PengRobinson(const double Tc, const double pc, const double omega)
56 : Tc_(Tc), pc_(pc), omega_(omega)
57 {
58 const double gas_constant =
60 a_ = 0.457235 * gas_constant * gas_constant * Tc_ * Tc_ / pc_;
61 b_ = 0.077796 * gas_constant * Tc_ / pc_;
62 }
63
65 MaterialPropertyLib::VariableArray const& variable_array,
66 ParameterLib::SpatialPosition const& pos, double const t,
67 double const dt) const override;
69 MaterialPropertyLib::VariableArray const& variable_array,
70 Variable const variable, ParameterLib::SpatialPosition const& pos,
71 double const t, double const dt) const override;
72
73private:
75 const double Tc_;
77 const double pc_;
79 const double omega_;
81 double a_;
83 double b_;
84};
85} // namespace MaterialPropertyLib
const double pc_
Critical pressure.
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.
PengRobinson(const double Tc, const double pc, const double omega)
const double omega_
Acentric factor.
double b_
Parameter b (covolume)
double a_
Parameter a (cohesion pressure)
virtual PropertyDataType value() const
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