3.5-git
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
38template<class ScalarT>
40{
46};
47
53template<class ScalarT>
54class MultiShapeTwoPLocalRules : public Dumux::FluidMatrix::Adapter<MultiShapeTwoPLocalRules<ScalarT>, Dumux::FluidMatrix::PcKrSw>
55{
56public:
57 using Scalar = ScalarT;
59
61 {
63
64 template<class SpatialParams, class Element, class SubControlVolume, class ElemSol>
65 BasicParams(const SpatialParams& spatialParams,
66 const Element& element,
67 const SubControlVolume& scv,
68 const ElemSol& elemSol)
69 {
70 shape_ = spatialParams.gridGeometry().poreGeometry(scv.dofIndex());
71
72 switch (shape_)
73 {
79 platonicBodyParams_ = std::make_unique<PlatonicBodyParams<Scalar>>(spatialParams, element, scv, elemSol);
80 break;
81 default:
82 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
83 }
84 }
85
86 // we use cube specialization here, but all platonic bodies have the same params
88 { return *platonicBodyParams_; }
89
91 { return shape_; }
92
94 {
95 shape_ = platonicBodyParams.poreShape();
96 platonicBodyParams_ = std::make_unique<PlatonicBodyParams<Scalar>>(platonicBodyParams);
97 }
98
99 private:
100 std::unique_ptr<PlatonicBodyParams<Scalar>> platonicBodyParams_;
101 Pore::Shape shape_;
102 };
103
104 static constexpr bool supportsMultipleGeometries()
105 { return true; }
106
110 static constexpr int numFluidPhases()
111 { return 2; }
112
114 //const RegularizationParams& regParams = {}, TODO
115 const std::string& paramGroup = "")
116 {
117 shape_ = baseParams.poreShape();
118 switch (shape_)
119 {
121 localRulesForTetrahedron_ = std::make_shared<typename LocalRules::Tetrahedron>(baseParams.platonicBodyParams(),
123 paramGroup);
124 break;
126 localRulesForCube_ = std::make_shared<typename LocalRules::Cube>(baseParams.platonicBodyParams(),
128 paramGroup);
129 break;
131 localRulesForOctahedron_ = std::make_shared<typename LocalRules::Octahedron>(baseParams.platonicBodyParams(),
133 paramGroup);
134 break;
136 localRulesForIcosahedron_ = std::make_shared<typename LocalRules::Icosahedron>(baseParams.platonicBodyParams(),
138 paramGroup);
139 break;
141 localRulesForDodecahedron_ = std::make_shared<typename LocalRules::Dodecahedron>(baseParams.platonicBodyParams(),
143 paramGroup);
144 break;
145 default:
146 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
147 }
148 }
149
150
151 template<class SpatialParams, class Element, class SubControlVolume, class ElemSol>
152 void updateParams(const SpatialParams& spatialParams,
153 const Element& element,
154 const SubControlVolume& scv,
155 const ElemSol& elemSol)
156 {
157 shape_ = spatialParams.gridGeometry().poreGeometry(scv.dofIndex());
158 switch (shape_)
159 {
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);
170 default:
171 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
172 }
173 }
174
178 Scalar pc(const Scalar sw) const
179 {
180 switch (shape_)
181 {
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);
192 default:
193 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
194 }
195 }
196
200 Scalar sw(const Scalar pc) const
201 {
202 switch (shape_)
203 {
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);
214 default:
215 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
216 }
217 }
218
222 Scalar dpc_dsw(const Scalar sw) const
223 {
224 switch (shape_)
225 {
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);
236 default:
237 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
238 }
239 }
240
244 Scalar dsw_dpc(const Scalar pc) const
245 {
246 switch (shape_)
247 {
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);
258 default:
259 DUNE_THROW(Dune::NotImplemented, "Invalid shape");
260 }
261 }
262
267 Scalar krw(const Scalar sw) const
268 { return 1.0; }
269
275 { return 0; }
276
281 Scalar krn(const Scalar sw) const
282 { return 1.0; }
283
289 { return 0.0; }
290
291private:
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_;
297
298 // TODO add more shapes here
299
300 Pore::Shape shape_;
301};
302
303} // end namespace Dumux::FluidMatrix
304
305#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: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 &paramGroup="")
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
BasicParams(const SpatialParams &spatialParams, const Element &element, const SubControlVolume &scv, const ElemSol &elemSol)
Definition: multishapelocalrules.hh:65
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