25#ifndef DUMUX_PNM_2P_LOCAL_RULES_HH
26#define DUMUX_PNM_2P_LOCAL_RULES_HH
38template<
class ScalarT>
53template<
class ScalarT>
64 template<
class SpatialParams,
class Element,
class SubControlVolume,
class ElemSol>
66 const Element& element,
67 const SubControlVolume& scv,
68 const ElemSol& elemSol)
70 shape_ = spatialParams.
gridGeometry().poreGeometry(scv.dofIndex());
79 platonicBodyParams_ = std::make_unique<PlatonicBodyParams<Scalar>>(spatialParams, element, scv, elemSol);
82 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
88 {
return *platonicBodyParams_; }
96 platonicBodyParams_ = std::make_unique<PlatonicBodyParams<Scalar>>(
platonicBodyParams);
100 std::unique_ptr<PlatonicBodyParams<Scalar>> platonicBodyParams_;
115 const std::string& paramGroup =
"")
121 localRulesForTetrahedron_ = std::make_shared<typename LocalRules::Tetrahedron>(baseParams.
platonicBodyParams(),
126 localRulesForCube_ = std::make_shared<typename LocalRules::Cube>(baseParams.
platonicBodyParams(),
131 localRulesForOctahedron_ = std::make_shared<typename LocalRules::Octahedron>(baseParams.
platonicBodyParams(),
136 localRulesForIcosahedron_ = std::make_shared<typename LocalRules::Icosahedron>(baseParams.
platonicBodyParams(),
141 localRulesForDodecahedron_ = std::make_shared<typename LocalRules::Dodecahedron>(baseParams.
platonicBodyParams(),
146 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
151 template<
class SpatialParams,
class Element,
class SubControlVolume,
class ElemSol>
153 const Element& element,
154 const SubControlVolume& scv,
155 const ElemSol& elemSol)
157 shape_ = spatialParams.
gridGeometry().poreGeometry(scv.dofIndex());
161 return localRulesForTetrahedron_->updateParams(spatialParams, element, scv, elemSol);
163 return localRulesForCube_->updateParams(spatialParams, element, scv, elemSol);
165 return localRulesForOctahedron_->updateParams(spatialParams, element, scv, elemSol);
167 return localRulesForIcosahedron_->updateParams(spatialParams, element, scv, elemSol);
169 return localRulesForDodecahedron_->updateParams(spatialParams, element, scv, elemSol);
171 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
183 return localRulesForTetrahedron_->pc(
sw);
185 return localRulesForCube_->pc(
sw);
187 return localRulesForOctahedron_->pc(
sw);
189 return localRulesForIcosahedron_->pc(
sw);
191 return localRulesForDodecahedron_->pc(
sw);
193 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
205 return localRulesForTetrahedron_->sw(
pc);
207 return localRulesForCube_->sw(
pc);
209 return localRulesForOctahedron_->sw(
pc);
211 return localRulesForIcosahedron_->sw(
pc);
213 return localRulesForDodecahedron_->sw(
pc);
215 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
227 return localRulesForTetrahedron_->dpc_dsw(
sw);
229 return localRulesForCube_->dpc_dsw(
sw);
231 return localRulesForOctahedron_->dpc_dsw(
sw);
233 return localRulesForIcosahedron_->dpc_dsw(
sw);
235 return localRulesForDodecahedron_->dpc_dsw(
sw);
237 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
249 return localRulesForTetrahedron_->dsw_dpc(
pc);
251 return localRulesForCube_->dsw_dpc(
pc);
253 return localRulesForOctahedron_->dsw_dpc(
pc);
255 return localRulesForIcosahedron_->dsw_dpc(
pc);
257 return localRulesForDodecahedron_->dsw_dpc(
pc);
259 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
292 std::shared_ptr<typename LocalRules::Tetrahedron> localRulesForTetrahedron_;
293 std::shared_ptr<typename LocalRules::Cube> localRulesForCube_;
294 std::shared_ptr<typename LocalRules::Octahedron> localRulesForOctahedron_;
295 std::shared_ptr<typename LocalRules::Icosahedron> localRulesForIcosahedron_;
296 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:40
Shape
Collection of different pore-body shapes.
Definition: poreproperties.hh:35
const GridGeometry & gridGeometry() const
The finite volume grid geometry.
Definition: common/fvspatialparams.hh:142
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:67
The parameter type.
Definition: localrulesforplatonicbody.hh:48
LocalRulesTraits for implementation of capillary pressure curves for multiple pore body geometries.
Definition: multishapelocalrules.hh:40
Implementation of capillary pressure curves for multiple pore body geometries.
Definition: multishapelocalrules.hh:55
static constexpr bool supportsMultipleGeometries()
Definition: multishapelocalrules.hh:104
Scalar dkrw_dsw(const Scalar sw) const
The derivative of the relative permeability for the wetting phase w.r.t. saturation.
Definition: multishapelocalrules.hh:274
MultiShapeTwoPLocalRules(const BasicParams &baseParams, const std::string ¶mGroup="")
Definition: multishapelocalrules.hh:113
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:288
Scalar sw(const Scalar pc) const
The saturation-capilllary-pressure curve.
Definition: multishapelocalrules.hh:200
Scalar dpc_dsw(const Scalar sw) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: multishapelocalrules.hh:222
Scalar pc(const Scalar sw) const
The capillary pressure-saturation curve.
Definition: multishapelocalrules.hh:178
Scalar krw(const Scalar sw) const
The relative permeability for the wetting phase.
Definition: multishapelocalrules.hh:267
static constexpr int numFluidPhases()
Return the number of fluid phases.
Definition: multishapelocalrules.hh:110
Scalar krn(const Scalar sw) const
The relative permeability for the non-wetting phase.
Definition: multishapelocalrules.hh:281
Scalar dsw_dpc(const Scalar pc) const
The partial derivative of the saturation to the capillary pressure.
Definition: multishapelocalrules.hh:244
void updateParams(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: multishapelocalrules.hh:152
ScalarT Scalar
Definition: multishapelocalrules.hh:57
Definition: multishapelocalrules.hh:61
BasicParams(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: multishapelocalrules.hh:65
BasicParams()
Definition: multishapelocalrules.hh:62
PlatonicBodyParams< Scalar > platonicBodyParams() const
Definition: multishapelocalrules.hh:87
void setParams(const PlatonicBodyParams< Scalar > &platonicBodyParams)
Definition: multishapelocalrules.hh:93
Pore::Shape poreShape() const
Definition: multishapelocalrules.hh:90
Base class for all standard pore-local pc-Sw curves.
Definition: singleshapelocalrules.hh:44
typename Regularization::template Params< Scalar > RegularizationParams
Definition: singleshapelocalrules.hh:50
The base class for spatial parameters for pore-network models.
Definition: porenetwork/common/spatialparams.hh:45