25#ifndef DUMUX_PNM_2P_SINGLE_SHAPE_LOCAL_RULES_HH
26#define DUMUX_PNM_2P_SINGLE_SHAPE_LOCAL_RULES_HH
40template<
class ScalarType,
49 using BasicParams =
typename BaseLaw::template Params<Scalar>;
55 template<
class SpatialParams,
class Element,
class SubControlVolume,
class ElemSol>
57 const Element& element,
58 const SubControlVolume& scv,
59 const ElemSol& elemSol)
61 basicParams_.update(spatialParams, element, scv, elemSol);
62 regularization_.updateParams(spatialParams, element, scv, elemSol);
75 {
return !std::is_same<Regularization, Dumux::FluidMatrix::NoRegularization>::value; }
83 const std::string& paramGroup =
"")
84 : basicParams_(baseParams)
88 regularization_.init(
this, regParams, paramGroup);
94 template<
bool enableRegularization = isRegularized()>
97 if constexpr (enableRegularization)
99 const auto regularized = regularization_.pc(
sw);
101 return regularized.value();
104 return BaseLaw::pc(
sw, basicParams_);
110 template<
bool enableRegularization = isRegularized()>
113 if constexpr (enableRegularization)
115 const auto regularized = regularization_.dpc_dsw(
sw);
117 return regularized.value();
120 return BaseLaw::dpc_dsw(
sw, basicParams_);
126 template<
bool enableRegularization = isRegularized()>
129 if constexpr (enableRegularization)
131 const auto regularized = regularization_.sw(
pc);
133 return regularized.value();
136 return BaseLaw::sw(
pc, basicParams_);
142 template<
bool enableRegularization = isRegularized()>
145 if constexpr (enableRegularization)
147 const auto regularized = regularization_.dsw_dpc(
pc);
149 return regularized.value();
152 return BaseLaw::dsw_dpc(
pc, basicParams_);
159 template<
bool enableRegularization = isRegularized()>
167 template<
bool enableRegularization = isRegularized()>
175 template<
bool enableRegularization = isRegularized()>
183 template<
bool enableRegularization = isRegularized()>
192 return basicParams_ == o.basicParams_
193 && regularization_ == o.regularization_;
200 {
return basicParams_; }
205 Regularization regularization_;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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....
This file contains functions related to calculate pore-body properties.
Definition: localrulesforplatonicbody.hh:40
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:44
void updateParams(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: singleshapelocalrules.hh:56
Scalar dpc_dsw(const Scalar sw) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: singleshapelocalrules.hh:111
Scalar krw(const Scalar sw) const
The relative permeability for the wetting phase.
Definition: singleshapelocalrules.hh:160
static constexpr bool supportsMultipleGeometries()
Definition: singleshapelocalrules.hh:52
ScalarType Scalar
Definition: singleshapelocalrules.hh:47
Scalar dsw_dpc(const Scalar pc) const
The partial derivative of the saturation to the capillary pressure.
Definition: singleshapelocalrules.hh:143
Scalar sw(const Scalar pc) const
The saturation-capillary-pressure curve.
Definition: singleshapelocalrules.hh:127
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:184
Scalar pc(const Scalar sw) const
The capillary pressure-saturation curve.
Definition: singleshapelocalrules.hh:95
Scalar dkrw_dsw(const Scalar sw) const
The derivative of the relative permeability for the wetting phase w.r.t. saturation.
Definition: singleshapelocalrules.hh:168
typename BaseLaw::template Params< Scalar > BasicParams
Definition: singleshapelocalrules.hh:49
typename Regularization::template Params< Scalar > RegularizationParams
Definition: singleshapelocalrules.hh:50
static constexpr bool isRegularized()
Return whether this law is regularized.
Definition: singleshapelocalrules.hh:74
const BasicParams & basicParams() const
Return the base law's parameters.
Definition: singleshapelocalrules.hh:199
bool operator==(const SingleShapeTwoPLocalRules &o) const
Equality comparison with another instance.
Definition: singleshapelocalrules.hh:190
SingleShapeTwoPLocalRules(const BasicParams &baseParams={}, const RegularizationParams ®Params={}, const std::string ¶mGroup="")
Construct from parameter structs.
Definition: singleshapelocalrules.hh:81
Scalar krn(const Scalar sw) const
The relative permeability for the non-wetting phase.
Definition: singleshapelocalrules.hh:176
static constexpr int numFluidPhases()
Return the number of fluid phases.
Definition: singleshapelocalrules.hh:68
The base class for spatial parameters for pore-network models.
Definition: porenetwork/common/spatialparams.hh:45