www.mooseframework.org
PorousFlowRelativePermeabilityBase.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
11 
12 template <bool is_ad>
15 {
17  params.addRangeCheckedParam<Real>(
18  "scaling", 1.0, "scaling>=0", "Relative permeability is multiplied by this factor");
19  params.addRangeCheckedParam<Real>(
20  "s_res",
21  0,
22  "s_res >= 0 & s_res < 1",
23  "The residual saturation of the phase j. Must be between 0 and 1");
24  params.addRangeCheckedParam<Real>(
25  "sum_s_res",
26  0,
27  "sum_s_res >= 0 & sum_s_res < 1",
28  "Sum of residual saturations over all phases. Must be between 0 and 1");
29  params.addPrivateParam<std::string>("pf_material_type", "relative_permeability");
30  params.addPrivateParam<bool>("is_ad", is_ad);
31  params.addClassDescription("Base class for PorousFlow relative permeability materials");
32  return params;
33 }
34 
35 template <bool is_ad>
37  const InputParameters & parameters)
38  : PorousFlowMaterialBase(parameters),
39  _scaling(getParam<Real>("scaling")),
40  _saturation(
41  _nodal_material
42  ? getGenericMaterialProperty<std::vector<Real>, is_ad>("PorousFlow_saturation_nodal")
43  : getGenericMaterialProperty<std::vector<Real>, is_ad>("PorousFlow_saturation_qp")),
44  _relative_permeability(
45  _nodal_material
46  ? declareGenericProperty<Real, is_ad>("PorousFlow_relative_permeability_nodal" + _phase)
47  : declareGenericProperty<Real, is_ad>("PorousFlow_relative_permeability_qp" + _phase)),
48  _drelative_permeability_ds(
49  is_ad ? nullptr
50  : _nodal_material
51  ? &declarePropertyDerivative<Real>("PorousFlow_relative_permeability_nodal" + _phase,
52  _saturation_variable_name)
53  : &declarePropertyDerivative<Real>("PorousFlow_relative_permeability_qp" + _phase,
54  _saturation_variable_name)),
55  _s_res(getParam<Real>("s_res")),
56  _sum_s_res(getParam<Real>("sum_s_res")),
57  _dseff_ds(1.0 / (1.0 - _sum_s_res))
58 {
59  if (_sum_s_res < _s_res)
60  mooseError("Sum of residual saturations sum_s_res cannot be smaller than s_res in ", name());
61 }
62 
63 template <bool is_ad>
64 void
66 {
67  // Effective saturation
68  GenericReal<is_ad> seff = effectiveSaturation(_saturation[_qp][_phase_num]);
69  GenericReal<is_ad> relperm;
70  Real drelperm;
71 
72  if (seff < 0.0)
73  {
74  // Relative permeability is 0 for saturation less than residual
75  relperm = 0.0;
76  drelperm = 0.0;
77  }
78  else if (seff >= 0.0 && seff <= 1)
79  {
80  relperm = relativePermeability(seff);
82  }
83  else // seff > 1
84  {
85  // Relative permeability is 1 when fully saturated
86  relperm = 1.0;
87  drelperm = 0.0;
88  }
89 
90  _relative_permeability[_qp] = relperm * _scaling;
91 
92  if (!is_ad)
93  (*_drelative_permeability_ds)[_qp] = drelperm * _dseff_ds * _scaling;
94 }
95 
96 template <bool is_ad>
99  GenericReal<is_ad> saturation) const
100 {
101  return (saturation - _s_res) / (1.0 - _sum_s_res);
102 }
103 
static InputParameters validParams()
void addPrivateParam(const std::string &name, const T &value)
void mooseError(Args &&... args)
const Real _sum_s_res
Sum of residual saturations over all phases.
auto raw_value(const Eigen::Map< T > &in)
T relativePermeability(const T &s, Real c, Real sn, Real ss, Real kn, Real ks)
Relative permeability as a function of saturation.
Real effectiveSaturation(Real pc, Real c, Real sn, Real ss, Real las)
Effective saturation as a function of capillary pressure If pc>=0 this will yield 1...
const Real _s_res
Residual saturation of specified phase.
Base class for all PorousFlow materials that provide phase-dependent properties.
Base class for PorousFlow relative permeability materials.
const std::string name
Definition: Setup.h:20
virtual GenericReal< is_ad > effectiveSaturation(GenericReal< is_ad > saturation) const
Effective saturation of fluid phase.
Real dRelativePermeability(Real s, Real c, Real sn, Real ss, Real kn, Real ks)
Derivative of relative permeability with respect to saturation.
PorousFlowRelativePermeabilityBaseTempl(const InputParameters &parameters)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void addClassDescription(const std::string &doc_string)
void addRangeCheckedParam(const std::string &name, const T &value, const std::string &parsed_function, const std::string &doc_string)
typename Moose::GenericType< Real, is_ad > GenericReal