OGS
TwoPhaseFlowWithPrhoMaterialProperties.cpp
Go to the documentation of this file.
1
12
13#include <Eigen/LU>
14#include <utility>
15
16#include "BaseLib/Logging.h"
23#include "MeshLib/Mesh.h"
28
32namespace ProcessLib
33{
34namespace TwoPhaseFlowWithPrho
35{
37 MeshLib::PropertyVector<int> const* const material_ids,
38 std::vector<std::unique_ptr<
40 capillary_pressure_models,
41 std::vector<
42 std::unique_ptr<MaterialLib::PorousMedium::RelativePermeability>>&&
43 relative_permeability_models)
44 : _material_ids(material_ids),
45 _capillary_pressure_models(std::move(capillary_pressure_models)),
46 _relative_permeability_models(std::move(relative_permeability_models))
47{
48 DBUG("Create material properties for Two-Phase flow with P-RHO model.");
49}
50
52 const std::size_t element_id)
53{
54 if (!_material_ids)
55 {
56 return 0;
57 }
58
59 assert(element_id < _material_ids->size());
60 return (*_material_ids)[element_id];
61}
62
64 const double /*t*/, const ParameterLib::SpatialPosition& /*pos*/,
65 const double /*p*/, const double /*T*/, const double saturation) const
66{
67 return _relative_permeability_models[0]->getValue(saturation);
68}
69
71 const double /*t*/, const ParameterLib::SpatialPosition& /*pos*/,
72 const double /*p*/, const double /*T*/, const double saturation) const
73{
74 return _relative_permeability_models[1]->getValue(saturation);
75}
76
78 const int material_id, const double /*t*/,
79 const ParameterLib::SpatialPosition& /*pos*/, const double /*p*/,
80 const double /*T*/, const double saturation) const
81{
82 return _capillary_pressure_models[material_id]->getCapillaryPressure(
83 saturation);
84}
85
87 const int material_id, const double /*t*/,
88 const ParameterLib::SpatialPosition& /*pos*/, const double /*p*/,
89 const double /*T*/, const double saturation) const
90{
91 return _capillary_pressure_models[material_id]->getdPcdS(saturation);
92}
93
95 double const t,
97 int const material_id,
98 double const pg,
99 double const X,
100 double const T,
101 double& Sw,
102 double& X_m,
103 double& dsw_dpg,
104 double& dsw_dX,
105 double& dxm_dpg,
106 double& dxm_dX)
107{
108 { // Local Newton solver
109 using LocalJacobianMatrix =
110 Eigen::Matrix<double, 2, 2, Eigen::RowMajor>;
111 using LocalResidualVector = Eigen::Matrix<double, 2, 1>;
112 using LocalUnknownVector = Eigen::Matrix<double, 2, 1>;
113 LocalJacobianMatrix J_loc;
114
115 Eigen::PartialPivLU<LocalJacobianMatrix> linear_solver(2);
116 auto const update_residual = [&](LocalResidualVector& residual)
117 { calculateResidual(material_id, pg, X, T, Sw, X_m, residual); };
118
119 auto const update_jacobian = [&](LocalJacobianMatrix& jacobian)
120 {
121 calculateJacobian(material_id, t, x, pg, X, T, jacobian, Sw,
122 X_m); // for solution dependent Jacobians
123 };
124
125 auto const update_solution = [&](LocalUnknownVector const& increment)
126 {
127 // increment solution vectors
128 Sw += increment[0];
129 X_m += increment[1];
130 };
131
132 // TODO Make the following choice of maximum iterations and convergence
133 // criteria available from the input file configuration. See Ehlers
134 // material model implementation for the example.
135 const int maximum_iterations(20);
136 const double residuum_tolerance(1.e-14);
137 const double increment_tolerance(0);
138
139 auto newton_solver = NumLib::NewtonRaphson(
140 linear_solver, update_jacobian, update_residual, update_solution,
141 {maximum_iterations, residuum_tolerance, increment_tolerance});
142
143 auto const success_iterations = newton_solver.solve(J_loc);
144
145 if (!success_iterations)
146 {
147 return false;
148 }
149 }
150 dsw_dpg = calculatedSwdP(pg, Sw, X_m, T, material_id);
151 dsw_dX = calculatedSwdX(pg, X, Sw, X_m, T, material_id);
152 dxm_dpg = calculatedXmdP(pg, Sw, X_m, dsw_dpg, material_id);
153 dxm_dX = calculatedXmdX(pg, Sw, X_m, dsw_dX, material_id);
154 return true;
155}
157 const int material_id, double const pl, double const X, double const T,
158 double const Sw, double const rho_h2_wet, ResidualVector& res)
159{
160 const double pg =
161 pl + _capillary_pressure_models[material_id]->getCapillaryPressure(Sw);
162 const double rho_h2_nonwet = pg * H2 / IdealGasConstant / T;
163
164 // calculating residual
165 res[0] = calculateEquilibiumRhoWetLight(pg, Sw, rho_h2_wet);
166 res[1] = calculateSaturation(pl, X, Sw, rho_h2_wet, rho_h2_nonwet, T);
167}
168
170 const int material_id, double const /*t*/,
171 ParameterLib::SpatialPosition const& /*x*/, double const pl,
172 double const /*X*/, double const T, JacobianMatrix& Jac, double const Sw,
173 double const rho_h2_wet)
174{
175 const double pg =
176 pl + _capillary_pressure_models[material_id]->getCapillaryPressure(Sw);
177 const double rho_h2_nonwet = pg * H2 / IdealGasConstant / T;
178 double const rho_equili_h2_wet = pg * HenryConstantH2 * H2;
179 double const dPC_dSw =
180 _capillary_pressure_models[material_id]->getdPcdS(Sw);
181 double const drhoh2wet_dpg = HenryConstantH2 * H2;
182 Jac.setZero();
183 if ((1 - Sw) < (rho_equili_h2_wet - rho_h2_wet))
184 {
185 Jac(0, 0) = -1;
186 }
187 else
188 {
189 Jac(0, 0) = drhoh2wet_dpg * dPC_dSw;
190 Jac(0, 1) = -1;
191 }
192
193 Jac(1, 0) = rho_h2_nonwet - rho_h2_wet;
194 Jac(1, 1) = -Sw;
195}
196
201 double const pg, double const Sw, double const rho_wet_h2)
202{
203 double const rho_equilibrium_wet_h2 = pg * HenryConstantH2 * H2;
204 return std::min(1 - Sw, rho_equilibrium_wet_h2 - rho_wet_h2);
205}
206
211 double /*PL*/, double X, double Sw, double rho_wet_h2, double rho_nonwet_h2,
212 double /*T*/)
213{
214 return X - (Sw * rho_wet_h2 + (1 - Sw) * rho_nonwet_h2);
215}
216
221 double pl, double S, double rho_wet_h2, double const T,
222 int current_material_id) const
223{
224 const double pg =
225 pl +
226 _capillary_pressure_models[current_material_id]->getCapillaryPressure(
227 S);
228 double const rho_equilibrium_wet_h2 = pg * HenryConstantH2 * H2;
229 if ((1 - S) < (rho_equilibrium_wet_h2 - rho_wet_h2))
230 {
231 return 0.0;
232 }
233 double const drhoh2wet_dpg = HenryConstantH2 * H2;
234 double const drhoh2nonwet_dpg = H2 / IdealGasConstant / T;
235 double const alpha =
236 ((drhoh2nonwet_dpg - drhoh2wet_dpg) * (1 - S) + drhoh2wet_dpg);
237 double const beta = (drhoh2nonwet_dpg - drhoh2wet_dpg) *
238 pg; // NOTE here should be PG^h, but we ignore vapor
239 double const dPC_dSw =
240 _capillary_pressure_models[current_material_id]->getdPcdS(S);
241 return alpha / (beta - alpha * dPC_dSw);
242}
247 double const pl, const double /*X*/, const double S,
248 const double rho_wet_h2, double const T, int current_material_id) const
249{
250 const double pg =
251 pl +
252 _capillary_pressure_models[current_material_id]->getCapillaryPressure(
253 S);
254 double const rho_equilibrium_wet_h2 = pg * HenryConstantH2 * H2;
255 if ((1 - S) < (rho_equilibrium_wet_h2 - rho_wet_h2))
256 {
257 return 0.0;
258 }
259 double const drhoh2wet_dpg = HenryConstantH2 * H2;
260 double const drhoh2nonwet_dpg = H2 / IdealGasConstant / T;
261 double const alpha =
262 ((drhoh2nonwet_dpg - drhoh2wet_dpg) * (1 - S) + drhoh2wet_dpg);
263 double const beta = (drhoh2nonwet_dpg - drhoh2wet_dpg) *
264 pg; // NOTE here should be PG^h, but we ignore vapor
265 double const dPC_dSw =
266 _capillary_pressure_models[current_material_id]->getdPcdS(S);
267 return -1 / (beta - alpha * dPC_dSw);
268}
273 double pl, double Sw, double rho_wet_h2, double dSwdX,
274 int current_material_id) const
275{
276 const double pg =
277 pl +
278 _capillary_pressure_models[current_material_id]->getCapillaryPressure(
279 Sw);
280 double const rho_equilibrium_wet_h2 = pg * HenryConstantH2 * H2;
281 double const dPC_dSw =
282 _capillary_pressure_models[current_material_id]->getdPcdS(Sw);
283 if ((1 - Sw) < (rho_equilibrium_wet_h2 - rho_wet_h2))
284 {
285 return 1.0;
286 }
287 return HenryConstantH2 * H2 * dPC_dSw * dSwdX;
288}
293 double pl, double Sw, double rho_wet_h2, double dSwdP,
294 int current_material_id) const
295{
296 const double pg =
297 pl +
298 _capillary_pressure_models[current_material_id]->getCapillaryPressure(
299 Sw);
300 double const rho_equilibrium_wet_h2 = pg * HenryConstantH2 * H2;
301 double const dPC_dSw =
302 _capillary_pressure_models[current_material_id]->getdPcdS(Sw);
303 if ((1 - Sw) < (rho_equilibrium_wet_h2 - rho_wet_h2))
304 {
305 return 0.0;
306 }
307 return HenryConstantH2 * H2 * (1 + dPC_dSw * dSwdP);
308}
309} // namespace TwoPhaseFlowWithPrho
310} // namespace ProcessLib
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
Definition of the Mesh class.
Definition of the PiecewiseLinearInterpolation class.
std::optional< int > solve(JacobianMatrix &jacobian) const
static double calculateSaturation(double, double X, double Sw, double rho_wet_h2, double rho_nonwet_h2, double)
double calculatedXmdP(double pl, double Sw, double rho_wet_h2, double dSwdP, int current_material_id) const
std::vector< std::unique_ptr< MaterialLib::PorousMedium::CapillaryPressureSaturation > > _capillary_pressure_models
double calculatedSwdX(double const pl, const double, const double S, const double rho_wet_h2, double const T, int current_material_id) const
double getNonwetRelativePermeability(const double t, const ParameterLib::SpatialPosition &pos, const double p, const double T, const double saturation) const
void calculateResidual(const int material_id, double const pl, double const X, double const T, double Sw, double rho_h2_wet, ResidualVector &res)
double getWetRelativePermeability(const double t, const ParameterLib::SpatialPosition &pos, const double p, const double T, const double saturation) const
void calculateJacobian(const int material_id, double const t, ParameterLib::SpatialPosition const &x, double const pl, double const X, double const T, JacobianMatrix &Jac, double Sw, double rho_h2_wet)
TwoPhaseFlowWithPrhoMaterialProperties(MeshLib::PropertyVector< int > const *material_ids, std::vector< std::unique_ptr< MaterialLib::PorousMedium::CapillaryPressureSaturation > > &&capillary_pressure_models, std::vector< std::unique_ptr< MaterialLib::PorousMedium::RelativePermeability > > &&relative_permeability_models)
static double calculateEquilibiumRhoWetLight(double const pg, double const Sw, double const rho_wet_h2)
double getCapillaryPressureDerivative(const int material_id, const double t, const ParameterLib::SpatialPosition &pos, const double p, const double T, const double saturation) const
Eigen::Matrix< double, jacobian_residual_size, jacobian_residual_size, Eigen::RowMajor > JacobianMatrix
double calculatedSwdP(double pl, double S, double rho_wet_h2, double const T, int current_material_id) const
double getCapillaryPressure(const int material_id, const double t, const ParameterLib::SpatialPosition &pos, const double p, const double T, const double saturation) const
bool computeConstitutiveRelation(double const t, ParameterLib::SpatialPosition const &x, const int material_id, double const pg, double const X, double const T, double &Sw, double &X_m, double &dsw_dpg, double &dsw_dX, double &dxm_dpg, double &dxm_dX)
double calculatedXmdX(double pl, double Sw, double rho_wet_h2, double dSwdX, int current_material_id) const
std::vector< std::unique_ptr< MaterialLib::PorousMedium::RelativePermeability > > _relative_permeability_models