version 3.10-dev
multishapelocalrules.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_PNM_2P_LOCAL_RULES_HH
14#define DUMUX_PNM_2P_LOCAL_RULES_HH
15
19
26template<class ScalarT>
28{
34};
35
41template<class ScalarT>
42class MultiShapeTwoPLocalRules : public Dumux::FluidMatrix::Adapter<MultiShapeTwoPLocalRules<ScalarT>, Dumux::FluidMatrix::PcKrSw>
43{
44public:
45 using Scalar = ScalarT;
47
49 {
51
52 template<class SpatialParams, class Element, class SubControlVolume, class ElemSol>
53 BasicParams(const SpatialParams& spatialParams,
54 const Element& element,
55 const SubControlVolume& scv,
56 const ElemSol& elemSol)
57 {
58 shape_ = spatialParams.gridGeometry().poreGeometry(scv.dofIndex());
59
60 switch (shape_)
61 {
67 platonicBodyParams_ = std::make_unique<PlatonicBodyParams<Scalar>>(spatialParams, element, scv, elemSol);
68 break;
69 default:
70 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
71 }
72 }
73
74 // we use cube specialization here, but all platonic bodies have the same params
76 { return *platonicBodyParams_; }
77
79 { return shape_; }
80
82 {
83 shape_ = platonicBodyParams.poreShape();
84 platonicBodyParams_ = std::make_unique<PlatonicBodyParams<Scalar>>(platonicBodyParams);
85 }
86
87 private:
88 std::unique_ptr<PlatonicBodyParams<Scalar>> platonicBodyParams_;
89 Pore::Shape shape_;
90 };
91
92 static constexpr bool supportsMultipleGeometries()
93 { return true; }
94
98 static constexpr int numFluidPhases()
99 { return 2; }
100
102 //const RegularizationParams& regParams = {}, TODO
103 const std::string& paramGroup = "")
104 {
105 shape_ = baseParams.poreShape();
106 switch (shape_)
107 {
109 localRulesForTetrahedron_ = std::make_shared<typename LocalRules::Tetrahedron>(baseParams.platonicBodyParams(),
111 paramGroup);
112 break;
114 localRulesForCube_ = std::make_shared<typename LocalRules::Cube>(baseParams.platonicBodyParams(),
116 paramGroup);
117 break;
119 localRulesForOctahedron_ = std::make_shared<typename LocalRules::Octahedron>(baseParams.platonicBodyParams(),
121 paramGroup);
122 break;
124 localRulesForIcosahedron_ = std::make_shared<typename LocalRules::Icosahedron>(baseParams.platonicBodyParams(),
126 paramGroup);
127 break;
129 localRulesForDodecahedron_ = std::make_shared<typename LocalRules::Dodecahedron>(baseParams.platonicBodyParams(),
131 paramGroup);
132 break;
133 default:
134 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
135 }
136 }
137
138
139 template<class SpatialParams, class Element, class SubControlVolume, class ElemSol>
140 void updateParams(const SpatialParams& spatialParams,
141 const Element& element,
142 const SubControlVolume& scv,
143 const ElemSol& elemSol)
144 {
145 shape_ = spatialParams.gridGeometry().poreGeometry(scv.dofIndex());
146 switch (shape_)
147 {
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);
158 default:
159 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
160 }
161 }
162
166 Scalar pc(const Scalar sw) const
167 {
168 switch (shape_)
169 {
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);
180 default:
181 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
182 }
183 }
184
188 Scalar sw(const Scalar pc) const
189 {
190 switch (shape_)
191 {
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);
202 default:
203 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
204 }
205 }
206
210 Scalar dpc_dsw(const Scalar sw) const
211 {
212 switch (shape_)
213 {
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);
224 default:
225 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
226 }
227 }
228
232 Scalar dsw_dpc(const Scalar pc) const
233 {
234 switch (shape_)
235 {
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);
246 default:
247 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
248 }
249 }
250
255 Scalar krw(const Scalar sw) const
256 { return 1.0; }
257
263 { return 0; }
264
269 Scalar krn(const Scalar sw) const
270 { return 1.0; }
271
277 { return 0.0; }
278
279private:
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_;
285
286 // TODO add more shapes here
287
288 Pore::Shape shape_;
289};
290
291} // end namespace Dumux::FluidMatrix
292
293#endif
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 &paramGroup="")
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
BasicParams(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: multishapelocalrules.hh:53
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