25#ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWOP_MATERIAL_LAW_HH
26#define DUMUX_MATERIAL_FLUIDMATRIX_TWOP_MATERIAL_LAW_HH
52template<
class ScalarType,
62 using BasicParams =
typename BaseLaw::template Params<Scalar>;
78 {
return !std::is_same<Regularization, NoRegularization>::value; }
95 regularization_.init(
this, paramGroup);
105 : basicParams_(baseParams)
109 regularization_.init(
this, baseParams,
effToAbsParams, regParams);
115 template<
bool enableRegularization = isRegularized()>
119 if constexpr (enableRegularization)
121 const auto regularized = regularization_.pc(swe);
123 return regularized.value();
126 return BaseLaw::pc(swe, basicParams_);
132 template<
bool enableRegularization = isRegularized()>
136 if constexpr (enableRegularization)
138 const auto regularized = regularization_.dpc_dswe(swe);
151 return BaseLaw::endPointPc(basicParams_);
157 template<
bool enableRegularization = isRegularized()>
160 if constexpr (enableRegularization)
162 const auto regularized = regularization_.swe(
pc);
173 template<
bool enableRegularization = isRegularized()>
176 if constexpr (enableRegularization)
178 const auto regularized = regularization_.dswe_dpc(
pc);
189 template<
bool enableRegularization = isRegularized()>
193 if constexpr (enableRegularization)
195 const auto regularized = regularization_.krw(swe);
197 return regularized.value();
200 return BaseLaw::krw(swe, basicParams_);
206 template<
bool enableRegularization = isRegularized()>
210 if constexpr (enableRegularization)
212 const auto regularized = regularization_.dkrw_dswe(swe);
223 template<
bool enableRegularization = isRegularized()>
227 if constexpr (enableRegularization)
229 const auto regularized = regularization_.krn(swe);
231 return regularized.value();
234 return BaseLaw::krn(swe, basicParams_);
240 template<
bool enableRegularization = isRegularized()>
244 if constexpr (enableRegularization)
246 const auto regularized = regularization_.dkrn_dswe(swe);
259 return basicParams_ == o.basicParams_
260 && effToAbsParams_ == o.effToAbsParams_
261 && regularization_ == o.regularization_;
270 return BaseLaw::template makeParams<Scalar>(paramGroup);
277 {
return basicParams_; }
285 return EffToAbs::template makeParams<Scalar>(paramGroup);
292 {
return effToAbsParams_; }
297 Regularization regularization_;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
This is a policy for 2p material laws how to convert absolute to relative saturations and vice versa.
A tag to turn off regularization and it's overhead.
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
Definition brookscorey.hh:35
This is a policy for 2p material laws how to convert absolute to relative saturations and vice versa.
Definition efftoabsdefaultpolicy.hh:46
static Scalar swToSwe(const Scalar sw, const Params< Scalar > ¶ms)
Convert an absolute wetting saturation to an effective one.
Definition efftoabsdefaultpolicy.hh:118
static Scalar sweToSw(const Scalar swe, const Params< Scalar > ¶ms)
Convert an effective wetting saturation to an absolute one.
Definition efftoabsdefaultpolicy.hh:133
static Scalar dswe_dsw(const Params< Scalar > ¶ms)
Derivative of the effective saturation w.r.t. the absolute saturation.
Definition efftoabsdefaultpolicy.hh:147
static Scalar dsw_dswe(const Params< Scalar > ¶ms)
Derivative of the absolute saturation w.r.t. the effective saturation.
Definition efftoabsdefaultpolicy.hh:161
TwoPMaterialLaw()=delete
Deleted default constructor (so we are never in an undefined state).
static constexpr int numFluidPhases()
Return the number of fluid phases.
Definition materiallaw.hh:71
Scalar endPointPc() const
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition materiallaw.hh:149
Scalar dsw_dpc(const Scalar pc) const
The partial derivative of the saturation to the capillary pressure.
Definition materiallaw.hh:174
TwoPMaterialLaw(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition materiallaw.hh:90
TwoPMaterialLaw(const BasicParams &baseParams, const EffToAbsParams &effToAbsParams={}, const RegularizationParams ®Params={})
Construct from parameter structs.
Definition materiallaw.hh:102
static constexpr bool isRegularized()
Return whether this law is regularized.
Definition materiallaw.hh:77
Scalar dkrn_dsw(const Scalar sw) const
The derivative of the relative permeability for the non-wetting phase w.r.t. saturation.
Definition materiallaw.hh:241
static EffToAbsParams makeEffToAbsParams(const std::string ¶mGroup)
Definition materiallaw.hh:283
Scalar Scalar
Definition materiallaw.hh:60
Scalar krn(const Scalar sw) const
The relative permeability for the non-wetting phase.
Definition materiallaw.hh:224
const BasicParams & basicParams() const
Return the base law's parameters.
Definition materiallaw.hh:276
static BasicParams makeBasicParams(const std::string ¶mGroup)
Definition materiallaw.hh:268
Scalar dkrw_dsw(const Scalar sw) const
The derivative of the relative permeability for the wetting phase w.r.t. saturation.
Definition materiallaw.hh:207
typename BrooksCoreyRegularization< Scalar >::template Params< Scalar > RegularizationParams
Definition materiallaw.hh:64
const EffToAbsParams & effToAbsParams() const
Definition materiallaw.hh:291
bool operator==(const TwoPMaterialLaw &o) const
Equality comparison with another instance.
Definition materiallaw.hh:257
Scalar dpc_dsw(const Scalar sw) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition materiallaw.hh:133
TwoPEffToAbsDefaultPolicy EffToAbs
Definition materiallaw.hh:66
Scalar sw(const Scalar pc) const
Definition materiallaw.hh:158
typename TwoPEffToAbsDefaultPolicy::template Params< Scalar > EffToAbsParams
Definition materiallaw.hh:63
typename BrooksCorey::template Params< Scalar > BasicParams
Definition materiallaw.hh:62
Scalar pc(const Scalar sw) const
The capillary pressure-saturation curve.
Definition materiallaw.hh:116
Scalar krw(const Scalar sw) const
The relative permeability for the wetting phase.
Definition materiallaw.hh:190
A tag to turn off regularization and it's overhead.
Definition noregularization.hh:34
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition fluidmatrixinteraction.hh:67