18 "scaling", 1.0,
"scaling>=0",
"Relative permeability is multiplied by this factor");
22 "s_res >= 0 & s_res < 1",
23 "The residual saturation of the phase j. Must be between 0 and 1");
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");
39 _scaling(getParam<
Real>(
"scaling")),
42 ? getGenericMaterialProperty<
std::vector<
Real>, is_ad>(
"PorousFlow_saturation_nodal")
43 : getGenericMaterialProperty<
std::vector<
Real>, is_ad>(
"PorousFlow_saturation_qp")),
44 _relative_permeability(
46 ? declareGenericProperty<
Real, is_ad>(
"PorousFlow_relative_permeability_nodal" + _phase)
47 : declareGenericProperty<
Real, is_ad>(
"PorousFlow_relative_permeability_qp" + _phase)),
48 _drelative_permeability_ds(
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))
60 mooseError(
"Sum of residual saturations sum_s_res cannot be smaller than s_res in ",
name());
78 else if (seff >= 0.0 && seff <= 1)
90 _relative_permeability[_qp] = relperm * _scaling;
93 (*_drelative_permeability_ds)[_qp] = drelperm * _dseff_ds * _scaling;
101 return (saturation - _s_res) / (1.0 - _sum_s_res);
static InputParameters validParams()
static InputParameters validParams()
void mooseError(Args &&... args)
const Real _sum_s_res
Sum of residual saturations over all phases.
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.
virtual GenericReal< is_ad > effectiveSaturation(GenericReal< is_ad > saturation) const
Effective saturation of fluid phase.
virtual void computeQpProperties() override
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 ¶meters)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
typename Moose::GenericType< Real, is_ad > GenericReal