24#ifndef DUMUX_PNM_2P_SINGLE_SHAPE_LOCAL_RULES_HH
25#define DUMUX_PNM_2P_SINGLE_SHAPE_LOCAL_RULES_HH
38template<
class ScalarType,
47 using BasicParams =
typename BaseLaw::template Params<Scalar>;
53 template<
class SpatialParams,
class Element,
class SubControlVolume,
class ElemSol>
55 const Element& element,
56 const SubControlVolume& scv,
57 const ElemSol& elemSol)
59 basicParams_.update(spatialParams, element, scv, elemSol);
60 regularization_.updateParams(spatialParams, element, scv, elemSol);
73 {
return !std::is_same<Regularization, Dumux::FluidMatrix::NoRegularization>::value; }
81 const std::string& paramGroup =
"")
82 : basicParams_(baseParams)
86 regularization_.init(
this, regParams, paramGroup);
92 template<
bool enableRegularization = isRegularized()>
95 if constexpr (enableRegularization)
97 const auto regularized = regularization_.pc(
sw);
99 return regularized.value();
102 return BaseLaw::pc(
sw, basicParams_);
108 template<
bool enableRegularization = isRegularized()>
111 if constexpr (enableRegularization)
113 const auto regularized = regularization_.dpc_dsw(
sw);
115 return regularized.value();
118 return BaseLaw::dpc_dsw(
sw, basicParams_);
124 template<
bool enableRegularization = isRegularized()>
127 if constexpr (enableRegularization)
129 const auto regularized = regularization_.sw(
pc);
131 return regularized.value();
134 return BaseLaw::sw(
pc, basicParams_);
140 template<
bool enableRegularization = isRegularized()>
143 if constexpr (enableRegularization)
145 const auto regularized = regularization_.dsw_dpc(
pc);
147 return regularized.value();
150 return BaseLaw::dsw_dpc(
pc, basicParams_);
157 template<
bool enableRegularization = isRegularized()>
165 template<
bool enableRegularization = isRegularized()>
173 template<
bool enableRegularization = isRegularized()>
181 template<
bool enableRegularization = isRegularized()>
190 return basicParams_ == o.basicParams_
191 && regularization_ == o.regularization_;
198 {
return basicParams_; }
203 Regularization regularization_;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
A tag to turn off regularization and it's overhead.
This file contains functions related to calculate pore-body properties.
Definition: localrulesforplatonicbody.hh:38
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
Base class for all standard pore-local pc-Sw curves.
Definition: singleshapelocalrules.hh:42
void updateParams(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: singleshapelocalrules.hh:54
Scalar dpc_dsw(const Scalar sw) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: singleshapelocalrules.hh:109
Scalar krw(const Scalar sw) const
The relative permeability for the wetting phase.
Definition: singleshapelocalrules.hh:158
static constexpr bool supportsMultipleGeometries()
Definition: singleshapelocalrules.hh:50
ScalarType Scalar
Definition: singleshapelocalrules.hh:45
Scalar dsw_dpc(const Scalar pc) const
The partial derivative of the saturation to the capillary pressure.
Definition: singleshapelocalrules.hh:141
Scalar sw(const Scalar pc) const
The saturation-capillary-pressure curve.
Definition: singleshapelocalrules.hh:125
Scalar dkrn_dsw(const Scalar sw) const
The derivative of the relative permeability for the non-wetting phase w.r.t. saturation.
Definition: singleshapelocalrules.hh:182
Scalar pc(const Scalar sw) const
The capillary pressure-saturation curve.
Definition: singleshapelocalrules.hh:93
Scalar dkrw_dsw(const Scalar sw) const
The derivative of the relative permeability for the wetting phase w.r.t. saturation.
Definition: singleshapelocalrules.hh:166
typename BaseLaw::template Params< Scalar > BasicParams
Definition: singleshapelocalrules.hh:47
typename Regularization::template Params< Scalar > RegularizationParams
Definition: singleshapelocalrules.hh:48
static constexpr bool isRegularized()
Return whether this law is regularized.
Definition: singleshapelocalrules.hh:72
const BasicParams & basicParams() const
Return the base law's parameters.
Definition: singleshapelocalrules.hh:197
bool operator==(const SingleShapeTwoPLocalRules &o) const
Equality comparison with another instance.
Definition: singleshapelocalrules.hh:188
SingleShapeTwoPLocalRules(const BasicParams &baseParams={}, const RegularizationParams ®Params={}, const std::string ¶mGroup="")
Construct from parameter structs.
Definition: singleshapelocalrules.hh:79
Scalar krn(const Scalar sw) const
The relative permeability for the non-wetting phase.
Definition: singleshapelocalrules.hh:174
static constexpr int numFluidPhases()
Return the number of fluid phases.
Definition: singleshapelocalrules.hh:66