3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_PNM_2P_LOCAL_RULES_HH
26#define DUMUX_PNM_2P_LOCAL_RULES_HH
27
31
33
34template<class ScalarT>
36{
42};
43
44template<class ScalarT>
45class MultiShapeTwoPLocalRules : public Dumux::FluidMatrix::Adapter<MultiShapeTwoPLocalRules<ScalarT>, Dumux::FluidMatrix::PcKrSw>
46{
47public:
48 using Scalar = ScalarT;
50
52 {
54
55 template<class SpatialParams, class Element, class SubControlVolume, class ElemSol>
56 BasicParams(const SpatialParams& spatialParams,
57 const Element& element,
58 const SubControlVolume& scv,
59 const ElemSol& elemSol)
60 {
61 shape_ = spatialParams.gridGeometry().poreGeometry(scv.dofIndex());
62
63 switch (shape_)
64 {
70 platonicBodyParams_ = std::make_unique<PlatonicBodyParams<Scalar>>(spatialParams, element, scv, elemSol);
71 break;
72 default:
73 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
74 }
75 }
76
77 // we use cube specialization here, but all platonic bodies have the same params
79 { return *platonicBodyParams_; }
80
82 { return shape_; }
83
85 {
86 shape_ = platonicBodyParams.poreShape();
87 platonicBodyParams_ = std::make_unique<PlatonicBodyParams<Scalar>>(platonicBodyParams);
88 }
89
90 private:
91 std::unique_ptr<PlatonicBodyParams<Scalar>> platonicBodyParams_;
92 Pore::Shape shape_;
93 };
94
95 static constexpr bool supportsMultipleGeometries()
96 { return true; }
97
101 static constexpr int numFluidPhases()
102 { return 2; }
103
104 template<class SpatialParams, class Element, class SubControlVolume, class ElemSol>
105 [[deprecated("Use constructor of BasicParams directly. Will be removed before release 3.4")]]
106 static BasicParams makeParams(const SpatialParams& spatialParams,
107 const Element& element,
108 const SubControlVolume& scv,
109 const ElemSol& elemSol)
110 {
111 return BasicParams(spatialParams, element, scv, elemSol);
112 }
113
115 //const RegularizationParams& regParams = {}, TODO
116 const std::string& paramGroup = "")
117 {
118 shape_ = baseParams.poreShape();
119 switch (shape_)
120 {
122 localRulesForTetrahedron_ = std::make_shared<typename LocalRules::Tetrahedron>(baseParams.platonicBodyParams(),
124 paramGroup);
125 break;
127 localRulesForCube_ = std::make_shared<typename LocalRules::Cube>(baseParams.platonicBodyParams(),
129 paramGroup);
130 break;
132 localRulesForOctahedron_ = std::make_shared<typename LocalRules::Octahedron>(baseParams.platonicBodyParams(),
134 paramGroup);
135 break;
137 localRulesForIcosahedron_ = std::make_shared<typename LocalRules::Icosahedron>(baseParams.platonicBodyParams(),
139 paramGroup);
140 break;
142 localRulesForDodecahedron_ = std::make_shared<typename LocalRules::Dodecahedron>(baseParams.platonicBodyParams(),
144 paramGroup);
145 break;
146 default:
147 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
148 }
149 }
150
151
152 template<class SpatialParams, class Element, class SubControlVolume, class ElemSol>
153 void updateParams(const SpatialParams& spatialParams,
154 const Element& element,
155 const SubControlVolume& scv,
156 const ElemSol& elemSol)
157 {
158 shape_ = spatialParams.gridGeometry().poreGeometry(scv.dofIndex());
159 switch (shape_)
160 {
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);
171 default:
172 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
173 }
174 }
175
179 Scalar pc(const Scalar sw) const
180 {
181 switch (shape_)
182 {
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);
193 default:
194 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
195 }
196 }
197
201 Scalar sw(const Scalar pc) const
202 {
203 switch (shape_)
204 {
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);
215 default:
216 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
217 }
218 }
219
223 Scalar dpc_dsw(const Scalar sw) const
224 {
225 switch (shape_)
226 {
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);
237 default:
238 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
239 }
240 }
241
245 Scalar dsw_dpc(const Scalar pc) const
246 {
247 switch (shape_)
248 {
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);
259 default:
260 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
261 }
262 }
263
268 Scalar krw(const Scalar sw) const
269 { return 1.0; }
270
276 { return 0; }
277
282 Scalar krn(const Scalar sw) const
283 { return 1.0; }
284
290 { return 0.0; }
291
292private:
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_;
298
299 // TODO add more shapes here
300
301 Pore::Shape shape_;
302};
303
304} // end namespace Dumux::FluidMatrix
305
306#endif
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 &paramGroup="")
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
BasicParams(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: multishapelocalrules.hh:56
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