OGS
TwoPhaseFlowWithPrhoMaterialProperties.cpp
Go to the documentation of this file.
1 
12 
13 #include <utility>
14 
15 #include "BaseLib/Logging.h"
24 #include "MeshLib/Mesh.h"
25 #include "MeshLib/PropertyVector.h"
26 #include "NumLib/NewtonRaphson.h"
27 #include "ParameterLib/Parameter.h"
29 
33 namespace ProcessLib
34 {
35 namespace TwoPhaseFlowWithPrho
36 {
38  MeshLib::PropertyVector<int> const* const material_ids,
39  std::unique_ptr<MaterialLib::Fluid::FluidProperty>
40  liquid_density,
41  std::unique_ptr<MaterialLib::Fluid::FluidProperty>
42  viscosity,
43  std::unique_ptr<MaterialLib::Fluid::FluidProperty>
44  gas_density,
45  std::unique_ptr<MaterialLib::Fluid::FluidProperty>
46  gas_viscosity,
47  std::vector<std::unique_ptr<MaterialLib::PorousMedium::Permeability>>&&
48  intrinsic_permeability_models,
49  std::vector<std::unique_ptr<MaterialLib::PorousMedium::Porosity>>&&
50  porosity_models,
51  std::vector<std::unique_ptr<MaterialLib::PorousMedium::Storage>>&&
52  storage_models,
53  std::vector<std::unique_ptr<
55  capillary_pressure_models,
56  std::vector<
57  std::unique_ptr<MaterialLib::PorousMedium::RelativePermeability>>&&
58  relative_permeability_models)
59  : _liquid_density(std::move(liquid_density)),
60  _viscosity(std::move(viscosity)),
61  _gas_density(std::move(gas_density)),
62  _gas_viscosity(std::move(gas_viscosity)),
63  _material_ids(material_ids),
64  _intrinsic_permeability_models(std::move(intrinsic_permeability_models)),
65  _porosity_models(std::move(porosity_models)),
66  _storage_models(std::move(storage_models)),
67  _capillary_pressure_models(std::move(capillary_pressure_models)),
68  _relative_permeability_models(std::move(relative_permeability_models))
69 {
70  DBUG("Create material properties for Two-Phase flow with P-RHO model.");
71 }
72 
74  const std::size_t element_id)
75 {
76  if (!_material_ids)
77  {
78  return 0;
79  }
80 
81  assert(element_id < _material_ids->size());
82  return (*_material_ids)[element_id];
83 }
84 
86  const double p, const double T) const
87 {
88  ArrayType vars;
89  vars[static_cast<int>(MaterialLib::Fluid::PropertyVariableType::T)] = T;
90  vars[static_cast<int>(MaterialLib::Fluid::PropertyVariableType::p)] = p;
91  return _liquid_density->getValue(vars);
92 }
93 
95  const double p, const double T) const
96 {
97  ArrayType vars;
98  vars[static_cast<int>(MaterialLib::Fluid::PropertyVariableType::T)] = T;
99  vars[static_cast<int>(MaterialLib::Fluid::PropertyVariableType::p)] = p;
100  return _gas_density->getValue(vars);
101 }
102 
104  const double p, const double T) const
105 {
106  ArrayType vars;
107  vars[static_cast<int>(MaterialLib::Fluid::PropertyVariableType::T)] = T;
108  vars[static_cast<int>(MaterialLib::Fluid::PropertyVariableType::p)] = p;
109  return _viscosity->getValue(vars);
110 }
111 
113  const double p, const double T) const
114 {
115  ArrayType vars;
116  vars[static_cast<int>(MaterialLib::Fluid::PropertyVariableType::T)] = T;
117  vars[static_cast<int>(MaterialLib::Fluid::PropertyVariableType::p)] = p;
118  return _gas_viscosity->getValue(vars);
119 }
120 
122  const int material_id, const double t,
123  const ParameterLib::SpatialPosition& pos, const int /*dim*/) const
124 {
125  return _intrinsic_permeability_models[material_id]->getValue(t, pos, 0.0,
126  0.0);
127 }
128 
130  const int material_id, const double t,
131  const ParameterLib::SpatialPosition& pos, const double /*p*/,
132  const double T, const double porosity_variable) const
133 {
134  return _porosity_models[material_id]->getValue(t, pos, porosity_variable,
135  T);
136 }
137 
139  const double /*t*/, const ParameterLib::SpatialPosition& /*pos*/,
140  const double /*p*/, const double /*T*/, const double saturation) const
141 {
142  return _relative_permeability_models[0]->getValue(saturation);
143 }
144 
146  const double /*t*/, const ParameterLib::SpatialPosition& /*pos*/,
147  const double /*p*/, const double /*T*/, const double saturation) const
148 {
149  return _relative_permeability_models[1]->getValue(saturation);
150 }
151 
153  const int material_id, const double /*t*/,
154  const ParameterLib::SpatialPosition& /*pos*/, const double /*p*/,
155  const double /*T*/, const double saturation) const
156 {
157  return _capillary_pressure_models[material_id]->getCapillaryPressure(
158  saturation);
159 }
160 
162  const int material_id, const double /*t*/,
163  const ParameterLib::SpatialPosition& /*pos*/, const double /*p*/,
164  const double /*T*/, const double saturation) const
165 {
166  return _capillary_pressure_models[material_id]->getdPcdS(saturation);
167 }
168 
170  double const t,
172  int const material_id,
173  double const pg,
174  double const X,
175  double const T,
176  double& Sw,
177  double& X_m,
178  double& dsw_dpg,
179  double& dsw_dX,
180  double& dxm_dpg,
181  double& dxm_dX)
182 {
183  { // Local Newton solver
184  using LocalJacobianMatrix =
185  Eigen::Matrix<double, 2, 2, Eigen::RowMajor>;
186  using LocalResidualVector = Eigen::Matrix<double, 2, 1>;
187  using LocalUnknownVector = Eigen::Matrix<double, 2, 1>;
188  LocalJacobianMatrix J_loc;
189 
190  Eigen::PartialPivLU<LocalJacobianMatrix> linear_solver(2);
191  auto const update_residual = [&](LocalResidualVector& residual)
192  { calculateResidual(material_id, pg, X, T, Sw, X_m, residual); };
193 
194  auto const update_jacobian = [&](LocalJacobianMatrix& jacobian)
195  {
196  calculateJacobian(material_id, t, x, pg, X, T, jacobian, Sw,
197  X_m); // for solution dependent Jacobians
198  };
199 
200  auto const update_solution = [&](LocalUnknownVector const& increment)
201  {
202  // increment solution vectors
203  Sw += increment[0];
204  X_m += increment[1];
205  };
206 
207  // TODO Make the following choice of maximum iterations and convergence
208  // criteria available from the input file configuration. See Ehlers
209  // material model implementation for the example.
210  const int maximum_iterations(20);
211  const double residuum_tolerance(1.e-14);
212  const double increment_tolerance(0);
213 
214  auto newton_solver = NumLib::NewtonRaphson<
215  decltype(linear_solver), LocalJacobianMatrix,
216  decltype(update_jacobian), LocalResidualVector,
217  decltype(update_residual), decltype(update_solution)>(
218  linear_solver, update_jacobian, update_residual, update_solution,
219  {maximum_iterations, residuum_tolerance, increment_tolerance});
220 
221  auto const success_iterations = newton_solver.solve(J_loc);
222 
223  if (!success_iterations)
224  {
225  return false;
226  }
227  }
228  dsw_dpg = calculatedSwdP(pg, Sw, X_m, T, material_id);
229  dsw_dX = calculatedSwdX(pg, X, Sw, X_m, T, material_id);
230  dxm_dpg = calculatedXmdP(pg, Sw, X_m, dsw_dpg, material_id);
231  dxm_dX = calculatedXmdX(pg, Sw, X_m, dsw_dX, material_id);
232  return true;
233 }
235  const int material_id, double const pl, double const X, double const T,
236  double const Sw, double const rho_h2_wet, ResidualVector& res)
237 {
238  const double pg =
239  pl + _capillary_pressure_models[material_id]->getCapillaryPressure(Sw);
240  const double rho_h2_nonwet = pg * H2 / IdealGasConstant / T;
241 
242  // calculating residual
243  res(0) = calculateEquilibiumRhoWetLight(pg, Sw, rho_h2_wet);
244  res(1) = calculateSaturation(pl, X, Sw, rho_h2_wet, rho_h2_nonwet, T);
245 }
246 
248  const int material_id, double const /*t*/,
249  ParameterLib::SpatialPosition const& /*x*/, double const pl,
250  double const /*X*/, double const T, JacobianMatrix& Jac, double const Sw,
251  double const rho_h2_wet)
252 {
253  const double pg =
254  pl + _capillary_pressure_models[material_id]->getCapillaryPressure(Sw);
255  const double rho_h2_nonwet = pg * H2 / IdealGasConstant / T;
256  double const rho_equili_h2_wet = pg * HenryConstantH2 * H2;
257  double const dPC_dSw =
258  _capillary_pressure_models[material_id]->getdPcdS(Sw);
259  double const drhoh2wet_dpg = HenryConstantH2 * H2;
260  Jac.setZero();
261  if ((1 - Sw) < (rho_equili_h2_wet - rho_h2_wet))
262  {
263  Jac(0, 0) = -1;
264  }
265  else
266  {
267  Jac(0, 0) = drhoh2wet_dpg * dPC_dSw;
268  Jac(0, 1) = -1;
269  }
270 
271  Jac(1, 0) = rho_h2_nonwet - rho_h2_wet;
272  Jac(1, 1) = -Sw;
273 }
274 
279  double const pg, double const Sw, double const rho_wet_h2)
280 {
281  double const rho_equilibrium_wet_h2 = pg * HenryConstantH2 * H2;
282  return std::min(1 - Sw, rho_equilibrium_wet_h2 - rho_wet_h2);
283 }
284 
289  double /*PL*/, double X, double Sw, double rho_wet_h2, double rho_nonwet_h2,
290  double /*T*/)
291 {
292  return X - (Sw * rho_wet_h2 + (1 - Sw) * rho_nonwet_h2);
293 }
294 
299  double pl, double S, double rho_wet_h2, double const T,
300  int current_material_id) const
301 {
302  const double pg =
303  pl +
304  _capillary_pressure_models[current_material_id]->getCapillaryPressure(
305  S);
306  double const rho_equilibrium_wet_h2 = pg * HenryConstantH2 * H2;
307  if ((1 - S) < (rho_equilibrium_wet_h2 - rho_wet_h2))
308  {
309  return 0.0;
310  }
311  double const drhoh2wet_dpg = HenryConstantH2 * H2;
312  double const drhoh2nonwet_dpg = H2 / IdealGasConstant / T;
313  double const alpha =
314  ((drhoh2nonwet_dpg - drhoh2wet_dpg) * (1 - S) + drhoh2wet_dpg);
315  double const beta = (drhoh2nonwet_dpg - drhoh2wet_dpg) *
316  pg; // NOTE here should be PG^h, but we ignore vapor
317  double const dPC_dSw =
318  _capillary_pressure_models[current_material_id]->getdPcdS(S);
319  return alpha / (beta - alpha * dPC_dSw);
320 }
325  double const pl, const double /*X*/, const double S,
326  const double rho_wet_h2, double const T, int current_material_id) const
327 {
328  const double pg =
329  pl +
330  _capillary_pressure_models[current_material_id]->getCapillaryPressure(
331  S);
332  double const rho_equilibrium_wet_h2 = pg * HenryConstantH2 * H2;
333  if ((1 - S) < (rho_equilibrium_wet_h2 - rho_wet_h2))
334  {
335  return 0.0;
336  }
337  double const drhoh2wet_dpg = HenryConstantH2 * H2;
338  double const drhoh2nonwet_dpg = H2 / IdealGasConstant / T;
339  double const alpha =
340  ((drhoh2nonwet_dpg - drhoh2wet_dpg) * (1 - S) + drhoh2wet_dpg);
341  double const beta = (drhoh2nonwet_dpg - drhoh2wet_dpg) *
342  pg; // NOTE here should be PG^h, but we ignore vapor
343  double const dPC_dSw =
344  _capillary_pressure_models[current_material_id]->getdPcdS(S);
345  return -1 / (beta - alpha * dPC_dSw);
346 }
351  double pl, double Sw, double rho_wet_h2, double dSwdX,
352  int current_material_id) const
353 {
354  const double pg =
355  pl +
356  _capillary_pressure_models[current_material_id]->getCapillaryPressure(
357  Sw);
358  double const rho_equilibrium_wet_h2 = pg * HenryConstantH2 * H2;
359  double const dPC_dSw =
360  _capillary_pressure_models[current_material_id]->getdPcdS(Sw);
361  if ((1 - Sw) < (rho_equilibrium_wet_h2 - rho_wet_h2))
362  {
363  return 1.0;
364  }
365  return HenryConstantH2 * H2 * dPC_dSw * dSwdX;
366 }
371  double pl, double Sw, double rho_wet_h2, double dSwdP,
372  int current_material_id) const
373 {
374  const double pg =
375  pl +
376  _capillary_pressure_models[current_material_id]->getCapillaryPressure(
377  Sw);
378  double const rho_equilibrium_wet_h2 = pg * HenryConstantH2 * H2;
379  double const dPC_dSw =
380  _capillary_pressure_models[current_material_id]->getdPcdS(Sw);
381  if ((1 - Sw) < (rho_equilibrium_wet_h2 - rho_wet_h2))
382  {
383  return 0.0;
384  }
385  return HenryConstantH2 * H2 * (1 + dPC_dSw * dSwdP);
386 }
387 } // namespace TwoPhaseFlowWithPrho
388 } // namespace ProcessLib
void DBUG(char const *fmt, Args const &... args)
Definition: Logging.h:27
Definition of the Mesh class.
Definition of the PiecewiseLinearInterpolation class.
std::optional< int > solve(JacobianMatrix &jacobian) const
Definition: NewtonRaphson.h:56
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
TwoPhaseFlowWithPrhoMaterialProperties(MeshLib::PropertyVector< int > const *material_ids, std::unique_ptr< MaterialLib::Fluid::FluidProperty > liquid_density, std::unique_ptr< MaterialLib::Fluid::FluidProperty > viscosity, std::unique_ptr< MaterialLib::Fluid::FluidProperty > gas_density, std::unique_ptr< MaterialLib::Fluid::FluidProperty > gas_viscosity, std::vector< std::unique_ptr< MaterialLib::PorousMedium::Permeability >> &&intrinsic_permeability_models, std::vector< std::unique_ptr< MaterialLib::PorousMedium::Porosity >> &&porosity_models, std::vector< std::unique_ptr< MaterialLib::PorousMedium::Storage >> &&storage_models, std::vector< std::unique_ptr< MaterialLib::PorousMedium::CapillaryPressureSaturation >> &&capillary_pressure_models, std::vector< std::unique_ptr< MaterialLib::PorousMedium::RelativePermeability >> &&relative_permeability_models)
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
std::vector< std::unique_ptr< MaterialLib::PorousMedium::Permeability > > _intrinsic_permeability_models
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 getPorosity(const int material_id, const double t, const ParameterLib::SpatialPosition &pos, const double p, const double T, const double porosity_variable) const
bool computeConstitutiveRelation(double const t, ParameterLib::SpatialPosition const &x_position, 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 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)
static double calculateEquilibiumRhoWetLight(double const pg, double const Sw, double const rho_wet_h2)
std::vector< std::unique_ptr< MaterialLib::PorousMedium::Porosity > > _porosity_models
double getCapillaryPressureDerivative(const int material_id, const double t, const ParameterLib::SpatialPosition &pos, const double p, const double T, const double saturation) const
double calculatedSwdP(double pl, double S, double rho_wet_h2, double const T, int current_material_id) const
Eigen::Matrix< double, jacobian_residual_size, jacobian_residual_size, Eigen::RowMajor > JacobianMatrix
double getCapillaryPressure(const int material_id, const double t, const ParameterLib::SpatialPosition &pos, const double p, const double T, const double saturation) const
Eigen::MatrixXd getPermeability(const int material_id, const double t, const ParameterLib::SpatialPosition &pos, const int dim) const
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