13#ifndef DUMUX_PNM_2P_SINGLE_SHAPE_LOCAL_RULES_HH
14#define DUMUX_PNM_2P_SINGLE_SHAPE_LOCAL_RULES_HH
28template<
class ScalarType,
37 using BasicParams =
typename BaseLaw::template Params<Scalar>;
43 template<
class SpatialParams,
class Element,
class SubControlVolume,
class ElemSol>
45 const Element& element,
46 const SubControlVolume& scv,
47 const ElemSol& elemSol)
49 basicParams_.update(spatialParams, element, scv, elemSol);
50 regularization_.updateParams(spatialParams, element, scv, elemSol);
63 {
return !std::is_same<Regularization, Dumux::FluidMatrix::NoRegularization>::value; }
71 const std::string& paramGroup =
"")
72 : basicParams_(baseParams)
76 regularization_.init(
this, regParams, paramGroup);
82 template<
bool enableRegularization = isRegularized()>
85 if constexpr (enableRegularization)
87 const auto regularized = regularization_.pc(
sw);
89 return regularized.value();
92 return BaseLaw::pc(
sw, basicParams_);
98 template<
bool enableRegularization = isRegularized()>
101 if constexpr (enableRegularization)
103 const auto regularized = regularization_.dpc_dsw(
sw);
105 return regularized.value();
108 return BaseLaw::dpc_dsw(
sw, basicParams_);
114 template<
bool enableRegularization = isRegularized()>
117 if constexpr (enableRegularization)
119 const auto regularized = regularization_.sw(
pc);
121 return regularized.value();
124 return BaseLaw::sw(
pc, basicParams_);
130 template<
bool enableRegularization = isRegularized()>
133 if constexpr (enableRegularization)
135 const auto regularized = regularization_.dsw_dpc(
pc);
137 return regularized.value();
140 return BaseLaw::dsw_dpc(
pc, basicParams_);
147 template<
bool enableRegularization = isRegularized()>
155 template<
bool enableRegularization = isRegularized()>
163 template<
bool enableRegularization = isRegularized()>
171 template<
bool enableRegularization = isRegularized()>
180 return basicParams_ == o.basicParams_
181 && regularization_ == o.regularization_;
188 {
return basicParams_; }
193 Regularization regularization_;
Base class for all standard pore-local pc-Sw curves.
Definition: singleshapelocalrules.hh:32
void updateParams(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: singleshapelocalrules.hh:44
Scalar dpc_dsw(const Scalar sw) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: singleshapelocalrules.hh:99
Scalar krw(const Scalar sw) const
The relative permeability for the wetting phase.
Definition: singleshapelocalrules.hh:148
static constexpr bool supportsMultipleGeometries()
Definition: singleshapelocalrules.hh:40
ScalarType Scalar
Definition: singleshapelocalrules.hh:35
Scalar dsw_dpc(const Scalar pc) const
The partial derivative of the saturation to the capillary pressure.
Definition: singleshapelocalrules.hh:131
Scalar sw(const Scalar pc) const
The saturation-capillary-pressure curve.
Definition: singleshapelocalrules.hh:115
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:172
Scalar pc(const Scalar sw) const
The capillary pressure-saturation curve.
Definition: singleshapelocalrules.hh:83
Scalar dkrw_dsw(const Scalar sw) const
The derivative of the relative permeability for the wetting phase w.r.t. saturation.
Definition: singleshapelocalrules.hh:156
typename BaseLaw::template Params< Scalar > BasicParams
Definition: singleshapelocalrules.hh:37
typename Regularization::template Params< Scalar > RegularizationParams
Definition: singleshapelocalrules.hh:38
static constexpr bool isRegularized()
Return whether this law is regularized.
Definition: singleshapelocalrules.hh:62
const BasicParams & basicParams() const
Return the base law's parameters.
Definition: singleshapelocalrules.hh:187
bool operator==(const SingleShapeTwoPLocalRules &o) const
Equality comparison with another instance.
Definition: singleshapelocalrules.hh:178
SingleShapeTwoPLocalRules(const BasicParams &baseParams={}, const RegularizationParams ®Params={}, const std::string ¶mGroup="")
Construct from parameter structs.
Definition: singleshapelocalrules.hh:69
Scalar krn(const Scalar sw) const
The relative permeability for the non-wetting phase.
Definition: singleshapelocalrules.hh:164
static constexpr int numFluidPhases()
Return the number of fluid phases.
Definition: singleshapelocalrules.hh:56
The base class for spatial parameters for pore-network models.
Definition: porenetwork/common/spatialparams.hh:33
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
Definition: localrulesforplatonicbody.hh:28
A tag to turn off regularization and it's overhead.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
This file contains functions related to calculate pore-body properties.
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:55
A tag to turn off regularization and it's overhead.
Definition: noregularization.hh:22