13#ifndef DUMUX_PNM_2P_LOCAL_RULES_HH
14#define DUMUX_PNM_2P_LOCAL_RULES_HH
26template<
class ScalarT>
41template<
class ScalarT>
52 template<
class SpatialParams,
class Element,
class SubControlVolume,
class ElemSol>
54 const Element& element,
55 const SubControlVolume& scv,
56 const ElemSol& elemSol)
58 shape_ = spatialParams.
gridGeometry().poreGeometry(scv.dofIndex());
67 platonicBodyParams_ = std::make_unique<PlatonicBodyParams<Scalar>>(spatialParams, element, scv, elemSol);
70 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
76 {
return *platonicBodyParams_; }
84 platonicBodyParams_ = std::make_unique<PlatonicBodyParams<Scalar>>(
platonicBodyParams);
88 std::unique_ptr<PlatonicBodyParams<Scalar>> platonicBodyParams_;
103 const std::string& paramGroup =
"")
109 localRulesForTetrahedron_ = std::make_shared<typename LocalRules::Tetrahedron>(baseParams.
platonicBodyParams(),
114 localRulesForCube_ = std::make_shared<typename LocalRules::Cube>(baseParams.
platonicBodyParams(),
119 localRulesForOctahedron_ = std::make_shared<typename LocalRules::Octahedron>(baseParams.
platonicBodyParams(),
124 localRulesForIcosahedron_ = std::make_shared<typename LocalRules::Icosahedron>(baseParams.
platonicBodyParams(),
129 localRulesForDodecahedron_ = std::make_shared<typename LocalRules::Dodecahedron>(baseParams.
platonicBodyParams(),
134 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
139 template<
class SpatialParams,
class Element,
class SubControlVolume,
class ElemSol>
141 const Element& element,
142 const SubControlVolume& scv,
143 const ElemSol& elemSol)
145 shape_ = spatialParams.
gridGeometry().poreGeometry(scv.dofIndex());
149 return localRulesForTetrahedron_->updateParams(spatialParams, element, scv, elemSol);
151 return localRulesForCube_->updateParams(spatialParams, element, scv, elemSol);
153 return localRulesForOctahedron_->updateParams(spatialParams, element, scv, elemSol);
155 return localRulesForIcosahedron_->updateParams(spatialParams, element, scv, elemSol);
157 return localRulesForDodecahedron_->updateParams(spatialParams, element, scv, elemSol);
159 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
171 return localRulesForTetrahedron_->pc(
sw);
173 return localRulesForCube_->pc(
sw);
175 return localRulesForOctahedron_->pc(
sw);
177 return localRulesForIcosahedron_->pc(
sw);
179 return localRulesForDodecahedron_->pc(
sw);
181 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
193 return localRulesForTetrahedron_->sw(
pc);
195 return localRulesForCube_->sw(
pc);
197 return localRulesForOctahedron_->sw(
pc);
199 return localRulesForIcosahedron_->sw(
pc);
201 return localRulesForDodecahedron_->sw(
pc);
203 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
215 return localRulesForTetrahedron_->dpc_dsw(
sw);
217 return localRulesForCube_->dpc_dsw(
sw);
219 return localRulesForOctahedron_->dpc_dsw(
sw);
221 return localRulesForIcosahedron_->dpc_dsw(
sw);
223 return localRulesForDodecahedron_->dpc_dsw(
sw);
225 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
237 return localRulesForTetrahedron_->dsw_dpc(
pc);
239 return localRulesForCube_->dsw_dpc(
pc);
241 return localRulesForOctahedron_->dsw_dpc(
pc);
243 return localRulesForIcosahedron_->dsw_dpc(
pc);
245 return localRulesForDodecahedron_->dsw_dpc(
pc);
247 DUNE_THROW(Dune::NotImplemented,
"Invalid shape");
280 std::shared_ptr<typename LocalRules::Tetrahedron> localRulesForTetrahedron_;
281 std::shared_ptr<typename LocalRules::Cube> localRulesForCube_;
282 std::shared_ptr<typename LocalRules::Octahedron> localRulesForOctahedron_;
283 std::shared_ptr<typename LocalRules::Icosahedron> localRulesForIcosahedron_;
284 std::shared_ptr<typename LocalRules::Dodecahedron> localRulesForDodecahedron_;
const GridGeometry & gridGeometry() const
The finite volume grid geometry.
Definition: common/fvspatialparams.hh:130
Implementation of capillary pressure curves for multiple pore body geometries.
Definition: multishapelocalrules.hh:43
static constexpr bool supportsMultipleGeometries()
Definition: multishapelocalrules.hh:92
Scalar dkrw_dsw(const Scalar sw) const
The derivative of the relative permeability for the wetting phase w.r.t. saturation.
Definition: multishapelocalrules.hh:262
MultiShapeTwoPLocalRules(const BasicParams &baseParams, const std::string ¶mGroup="")
Definition: multishapelocalrules.hh:101
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:276
Scalar sw(const Scalar pc) const
The saturation-capilllary-pressure curve.
Definition: multishapelocalrules.hh:188
Scalar dpc_dsw(const Scalar sw) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: multishapelocalrules.hh:210
Scalar pc(const Scalar sw) const
The capillary pressure-saturation curve.
Definition: multishapelocalrules.hh:166
Scalar krw(const Scalar sw) const
The relative permeability for the wetting phase.
Definition: multishapelocalrules.hh:255
static constexpr int numFluidPhases()
Return the number of fluid phases.
Definition: multishapelocalrules.hh:98
Scalar krn(const Scalar sw) const
The relative permeability for the non-wetting phase.
Definition: multishapelocalrules.hh:269
Scalar dsw_dpc(const Scalar pc) const
The partial derivative of the saturation to the capillary pressure.
Definition: multishapelocalrules.hh:232
void updateParams(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: multishapelocalrules.hh:140
ScalarT Scalar
Definition: multishapelocalrules.hh:45
Base class for all standard pore-local pc-Sw curves.
Definition: singleshapelocalrules.hh:32
typename Regularization::template Params< Scalar > RegularizationParams
Definition: singleshapelocalrules.hh:38
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....
Pore-local pc-Sw curves for for platonic bodies (tetrahedron, cube, octahedron, dodecahedron,...
Definition: localrulesforplatonicbody.hh:28
Shape
Collection of different pore-body shapes.
Definition: poreproperties.hh:23
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
LocalRulesTraits for implementation of capillary pressure curves for multiple pore body geometries.
Definition: multishapelocalrules.hh:28
Definition: multishapelocalrules.hh:49
BasicParams(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: multishapelocalrules.hh:53
BasicParams()
Definition: multishapelocalrules.hh:50
PlatonicBodyParams< Scalar > platonicBodyParams() const
Definition: multishapelocalrules.hh:75
void setParams(const PlatonicBodyParams< Scalar > &platonicBodyParams)
Definition: multishapelocalrules.hh:81
Pore::Shape poreShape() const
Definition: multishapelocalrules.hh:78
The parameter type.
Definition: localrulesforplatonicbody.hh:36