25#ifndef DUMUX_PNM_2P_LOCAL_RULES_HH
26#define DUMUX_PNM_2P_LOCAL_RULES_HH
34template<
class ScalarT>
44template<
class ScalarT>
55 template<
class SpatialParams,
class Element,
class SubControlVolume,
class ElemSol>
57 const Element& element,
58 const SubControlVolume& scv,
59 const ElemSol& elemSol)
61 shape_ = spatialParams.gridGeometry().poreGeometry(scv.dofIndex());
70 platonicBodyParams_ = std::make_unique<PlatonicBodyParams<Scalar>>(spatialParams, element, scv, elemSol);
73 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
79 {
return *platonicBodyParams_; }
87 platonicBodyParams_ = std::make_unique<PlatonicBodyParams<Scalar>>(
platonicBodyParams);
91 std::unique_ptr<PlatonicBodyParams<Scalar>> platonicBodyParams_;
104 template<
class SpatialParams,
class Element,
class SubControlVolume,
class ElemSol>
105 [[deprecated(
"Use constructor of BasicParams directly. Will be removed before release 3.4")]]
107 const Element& element,
108 const SubControlVolume& scv,
109 const ElemSol& elemSol)
111 return BasicParams(spatialParams, element, scv, elemSol);
116 const std::string& paramGroup =
"")
122 localRulesForTetrahedron_ = std::make_shared<typename LocalRules::Tetrahedron>(baseParams.
platonicBodyParams(),
127 localRulesForCube_ = std::make_shared<typename LocalRules::Cube>(baseParams.
platonicBodyParams(),
132 localRulesForOctahedron_ = std::make_shared<typename LocalRules::Octahedron>(baseParams.
platonicBodyParams(),
137 localRulesForIcosahedron_ = std::make_shared<typename LocalRules::Icosahedron>(baseParams.
platonicBodyParams(),
142 localRulesForDodecahedron_ = std::make_shared<typename LocalRules::Dodecahedron>(baseParams.
platonicBodyParams(),
147 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
152 template<
class SpatialParams,
class Element,
class SubControlVolume,
class ElemSol>
154 const Element& element,
155 const SubControlVolume& scv,
156 const ElemSol& elemSol)
158 shape_ = spatialParams.gridGeometry().poreGeometry(scv.dofIndex());
162 return localRulesForTetrahedron_->updateParams(spatialParams, element, scv, elemSol);
164 return localRulesForCube_->updateParams(spatialParams, element, scv, elemSol);
166 return localRulesForOctahedron_->updateParams(spatialParams, element, scv, elemSol);
168 return localRulesForIcosahedron_->updateParams(spatialParams, element, scv, elemSol);
170 return localRulesForDodecahedron_->updateParams(spatialParams, element, scv, elemSol);
172 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
184 return localRulesForTetrahedron_->pc(
sw);
186 return localRulesForCube_->pc(
sw);
188 return localRulesForOctahedron_->pc(
sw);
190 return localRulesForIcosahedron_->pc(
sw);
192 return localRulesForDodecahedron_->pc(
sw);
194 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
206 return localRulesForTetrahedron_->sw(
pc);
208 return localRulesForCube_->sw(
pc);
210 return localRulesForOctahedron_->sw(
pc);
212 return localRulesForIcosahedron_->sw(
pc);
214 return localRulesForDodecahedron_->sw(
pc);
216 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
228 return localRulesForTetrahedron_->dpc_dsw(
sw);
230 return localRulesForCube_->dpc_dsw(
sw);
232 return localRulesForOctahedron_->dpc_dsw(
sw);
234 return localRulesForIcosahedron_->dpc_dsw(
sw);
236 return localRulesForDodecahedron_->dpc_dsw(
sw);
238 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
250 return localRulesForTetrahedron_->dsw_dpc(
pc);
252 return localRulesForCube_->dsw_dpc(
pc);
254 return localRulesForOctahedron_->dsw_dpc(
pc);
256 return localRulesForIcosahedron_->dsw_dpc(
pc);
258 return localRulesForDodecahedron_->dsw_dpc(
pc);
260 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
293 std::shared_ptr<typename LocalRules::Tetrahedron> localRulesForTetrahedron_;
294 std::shared_ptr<typename LocalRules::Cube> localRulesForCube_;
295 std::shared_ptr<typename LocalRules::Octahedron> localRulesForOctahedron_;
296 std::shared_ptr<typename LocalRules::Icosahedron> localRulesForIcosahedron_;
297 std::shared_ptr<typename LocalRules::Dodecahedron> localRulesForDodecahedron_;
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
Pore-local pc-Sw curves for for platonic bodies (tetrahedron, cube, octahedron, dodecahedron,...
This file contains functions related to calculate pore-body properties.
Definition: localrulesforplatonicbody.hh:38
Shape
Collection of different pore-body shapes.
Definition: poreproperties.hh:35
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:67
The parameter type.
Definition: localrulesforplatonicbody.hh:46
Definition: multishapelocalrules.hh:36
Definition: multishapelocalrules.hh:46
static constexpr bool supportsMultipleGeometries()
Definition: multishapelocalrules.hh:95
Scalar dkrw_dsw(const Scalar sw) const
The derivative of the relative permeability for the wetting phase w.r.t. saturation.
Definition: multishapelocalrules.hh:275
MultiShapeTwoPLocalRules(const BasicParams &baseParams, const std::string ¶mGroup="")
Definition: multishapelocalrules.hh:114
static BasicParams makeParams(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: multishapelocalrules.hh:106
Scalar dkrn_dsw(const Scalar sw) const
The derivative of the relative permeability for the non-wetting phase w.r.t. saturation.
Definition: multishapelocalrules.hh:289
Scalar sw(const Scalar pc) const
The saturation-capilllary-pressure curve.
Definition: multishapelocalrules.hh:201
Scalar dpc_dsw(const Scalar sw) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: multishapelocalrules.hh:223
Scalar pc(const Scalar sw) const
The capillary pressure-saturation curve.
Definition: multishapelocalrules.hh:179
Scalar krw(const Scalar sw) const
The relative permeability for the wetting phase.
Definition: multishapelocalrules.hh:268
static constexpr int numFluidPhases()
Return the number of fluid phases.
Definition: multishapelocalrules.hh:101
Scalar krn(const Scalar sw) const
The relative permeability for the non-wetting phase.
Definition: multishapelocalrules.hh:282
Scalar dsw_dpc(const Scalar pc) const
The partial derivative of the saturation to the capillary pressure.
Definition: multishapelocalrules.hh:245
void updateParams(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: multishapelocalrules.hh:153
ScalarT Scalar
Definition: multishapelocalrules.hh:48
Definition: multishapelocalrules.hh:52
BasicParams(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: multishapelocalrules.hh:56
BasicParams()
Definition: multishapelocalrules.hh:53
PlatonicBodyParams< Scalar > platonicBodyParams() const
Definition: multishapelocalrules.hh:78
void setParams(const PlatonicBodyParams< Scalar > &platonicBodyParams)
Definition: multishapelocalrules.hh:84
Pore::Shape poreShape() const
Definition: multishapelocalrules.hh:81
Base class for all standard pore-local pc-Sw curves.
Definition: singleshapelocalrules.hh:42
typename Regularization::template Params< Scalar > RegularizationParams
Definition: singleshapelocalrules.hh:48