multishapelocalrules.hh 12.8 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
 *   See the file COPYING for full copying permissions.                      *
 *                                                                           *
 *   This program is free software: you can redistribute it and/or modify    *
 *   it under the terms of the GNU General Public License as published by    *
 *   the Free Software Foundation, either version 3 of the License, or       *
 *   (at your option) any later version.                                     *
 *                                                                           *
 *   This program is distributed in the hope that it will be useful,         *
 *   but WITHOUT ANY WARRANTY; without even the implied warranty of          *
 *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the            *
 *   GNU General Public License for more details.                            *
 *                                                                           *
 *   You should have received a copy of the GNU General Public License       *
 *   along with this program.  If not, see <http://www.gnu.org/licenses/>.   *
 *****************************************************************************/
/*!
 * \file
 *
 * \brief Implementation of capillary pressure curves for multiple pore body geometries
 *
 */
#ifndef DUMUX_PNM_2P_LOCAL_RULES_HH
#define DUMUX_PNM_2P_LOCAL_RULES_HH

#include <dumux/material/fluidmatrixinteractions/fluidmatrixinteraction.hh>
#include <dumux/porenetwork/common/poreproperties.hh>
#include "localrulesforplatonicbody.hh"

namespace Dumux::PoreNetwork::FluidMatrix {

template<class ScalarT>
struct LocalRulesTraits
{
37
38
39
40
41
    using Tetrahedron = TwoPLocalRulesPlatonicBodyDefault<Pore::Shape::tetrahedron, ScalarT>;
    using Cube = TwoPLocalRulesPlatonicBodyDefault<Pore::Shape::cube, ScalarT>;
    using Octahedron = TwoPLocalRulesPlatonicBodyDefault<Pore::Shape::octahedron, ScalarT>;
    using Icosahedron = TwoPLocalRulesPlatonicBodyDefault<Pore::Shape::icosahedron, ScalarT>;
    using Dodecahedron = TwoPLocalRulesPlatonicBodyDefault<Pore::Shape::dodecahedron, ScalarT>;
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
};

template<class ScalarT>
class MultiShapeTwoPLocalRules : public Dumux::FluidMatrix::Adapter<MultiShapeTwoPLocalRules<ScalarT>, Dumux::FluidMatrix::PcKrSw>
{
public:
    using Scalar = ScalarT;
    using LocalRules = LocalRulesTraits<Scalar>;

    struct BasicParams
    {
        BasicParams() {}

        template<class SpatialParams, class Element, class SubControlVolume, class ElemSol>
        BasicParams(const SpatialParams& spatialParams,
                    const Element& element,
                    const SubControlVolume& scv,
                    const ElemSol& elemSol)
        {
            shape_ = spatialParams.gridGeometry().poreGeometry(scv.dofIndex());

            switch (shape_)
            {
                case Pore::Shape::tetrahedron:
                case Pore::Shape::cube:
                case Pore::Shape::octahedron:
                case Pore::Shape::icosahedron:
                case Pore::Shape::dodecahedron:
                    platonicBodyParams_ = std::make_unique<PlatonicBodyParams<Scalar>>(spatialParams, element, scv, elemSol);
                    break;
                default:
                    DUNE_THROW(Dune::NotImplemented, "Invalid shape");
            }
        }

         // we use cube specialization here, but all platonic bodies have the same params
        PlatonicBodyParams<Scalar> platonicBodyParams() const
        { return *platonicBodyParams_; }

        Pore::Shape poreShape() const
        { return shape_; }

        void setParams(const PlatonicBodyParams<Scalar>& platonicBodyParams)
        {
            shape_ = platonicBodyParams.poreShape();
            platonicBodyParams_ = std::make_unique<PlatonicBodyParams<Scalar>>(platonicBodyParams);
        }

    private:
        std::unique_ptr<PlatonicBodyParams<Scalar>> platonicBodyParams_;
        Pore::Shape shape_;
    };

    static constexpr bool supportsMultipleGeometries()
    { return true; }

    /*!
     * \brief Return the number of fluid phases
     */
    static constexpr int numFluidPhases()
    { return 2; }

    template<class SpatialParams, class Element, class SubControlVolume, class ElemSol>
    static BasicParams makeParams(const SpatialParams& spatialParams,
                                  const Element& element,
                                  const SubControlVolume& scv,
                                  const ElemSol& elemSol)
    {
        return BasicParams(spatialParams, element, scv, elemSol);
    }

    MultiShapeTwoPLocalRules(const BasicParams& baseParams,
                             //const RegularizationParams& regParams = {}, TODO
                             const std::string& paramGroup = "")
    {
        shape_ = baseParams.poreShape();
        switch (shape_)
        {
            case Pore::Shape::tetrahedron:
                localRulesForTetrahedron_ = std::make_shared<typename LocalRules::Tetrahedron>(baseParams.platonicBodyParams(),
                                                                                               typename LocalRules::Tetrahedron::RegularizationParams{},
                                                                                               paramGroup);
                break;
            case Pore::Shape::cube:
                localRulesForCube_ = std::make_shared<typename LocalRules::Cube>(baseParams.platonicBodyParams(),
                                                                                 typename LocalRules::Cube::RegularizationParams{},
                                                                                 paramGroup);
                break;
            case Pore::Shape::octahedron:
                localRulesForOctahedron_ = std::make_shared<typename LocalRules::Octahedron>(baseParams.platonicBodyParams(),
                                                                                             typename LocalRules::Octahedron::RegularizationParams{},
                                                                                             paramGroup);
                break;
            case Pore::Shape::icosahedron:
                localRulesForIcosahedron_ = std::make_shared<typename LocalRules::Icosahedron>(baseParams.platonicBodyParams(),
                                                                                               typename LocalRules::Icosahedron::RegularizationParams{},
                                                                                               paramGroup);
                break;
            case Pore::Shape::dodecahedron:
                localRulesForDodecahedron_ = std::make_shared<typename LocalRules::Dodecahedron>(baseParams.platonicBodyParams(),
                                                                                                 typename LocalRules::Dodecahedron::RegularizationParams{},
                                                                                                 paramGroup);
                break;
            default:
                DUNE_THROW(Dune::NotImplemented, "Invalid shape");
        }
    }


    template<class SpatialParams, class Element, class SubControlVolume, class ElemSol>
    void updateParams(const SpatialParams& spatialParams,
                      const Element& element,
                      const SubControlVolume& scv,
                      const ElemSol& elemSol)
    {
        shape_ = spatialParams.gridGeometry().poreGeometry(scv.dofIndex());
        switch (shape_)
        {
            case Pore::Shape::tetrahedron:
                return localRulesForTetrahedron_->updateParams(spatialParams, element, scv, elemSol);
            case Pore::Shape::cube:
                return localRulesForCube_->updateParams(spatialParams, element, scv, elemSol);
            case Pore::Shape::octahedron:
                return localRulesForOctahedron_->updateParams(spatialParams, element, scv, elemSol);
            case Pore::Shape::icosahedron:
                return localRulesForIcosahedron_->updateParams(spatialParams, element, scv, elemSol);
            case Pore::Shape::dodecahedron:
                return localRulesForDodecahedron_->updateParams(spatialParams, element, scv, elemSol);
            default:
                DUNE_THROW(Dune::NotImplemented, "Invalid shape");
        }
    }

    /*!
     * \brief The capillary pressure-saturation curve
     */
    Scalar pc(const Scalar sw) const
    {
        switch (shape_)
        {
            case Pore::Shape::tetrahedron:
                return localRulesForTetrahedron_->pc(sw);
            case Pore::Shape::cube:
                return localRulesForCube_->pc(sw);
            case Pore::Shape::octahedron:
                return localRulesForOctahedron_->pc(sw);
            case Pore::Shape::icosahedron:
                return localRulesForIcosahedron_->pc(sw);
            case Pore::Shape::dodecahedron:
                return localRulesForDodecahedron_->pc(sw);
            default:
                DUNE_THROW(Dune::NotImplemented, "Invalid shape");
        }
    }

    /*!
     * \brief The saturation-capilllary-pressure curve
     */
    Scalar sw(const Scalar pc) const
    {
        switch (shape_)
        {
            case Pore::Shape::tetrahedron:
                return localRulesForTetrahedron_->sw(pc);
            case Pore::Shape::cube:
                return localRulesForCube_->sw(pc);
            case Pore::Shape::octahedron:
                return localRulesForOctahedron_->sw(pc);
            case Pore::Shape::icosahedron:
                return localRulesForIcosahedron_->sw(pc);
            case Pore::Shape::dodecahedron:
                return localRulesForDodecahedron_->sw(pc);
            default:
                DUNE_THROW(Dune::NotImplemented, "Invalid shape");
        }
    }

    /*!
     * \brief The partial derivative of the capillary pressure w.r.t. the saturation
     */
    Scalar dpc_dsw(const Scalar sw) const
    {
        switch (shape_)
        {
            case Pore::Shape::tetrahedron:
                return localRulesForTetrahedron_->dpc_dsw(sw);
            case Pore::Shape::cube:
                return localRulesForCube_->dpc_dsw(sw);
            case Pore::Shape::octahedron:
                return localRulesForOctahedron_->dpc_dsw(sw);
            case Pore::Shape::icosahedron:
                return localRulesForIcosahedron_->dpc_dsw(sw);
            case Pore::Shape::dodecahedron:
                return localRulesForDodecahedron_->dpc_dsw(sw);
            default:
                DUNE_THROW(Dune::NotImplemented, "Invalid shape");
        }
    }

    /*!
     * \brief The partial derivative of the saturation to the capillary pressure
     */
    Scalar dsw_dpc(const Scalar pc) const
    {
        switch (shape_)
        {
            case Pore::Shape::tetrahedron:
                return localRulesForTetrahedron_->dsw_dpc(pc);
            case Pore::Shape::cube:
                return localRulesForCube_->dsw_dpc(pc);
            case Pore::Shape::octahedron:
                return localRulesForOctahedron_->dsw_dpc(pc);
            case Pore::Shape::icosahedron:
                return localRulesForIcosahedron_->dsw_dpc(pc);
            case Pore::Shape::dodecahedron:
                return localRulesForDodecahedron_->dsw_dpc(pc);
            default:
                DUNE_THROW(Dune::NotImplemented, "Invalid shape");
        }
    }

    /*!
     * \brief The relative permeability for the wetting phase.
     * \note This is only for compatibility. Will not be used.
     */
    Scalar krw(const Scalar sw) const
    { return 1.0; }

    /*!
     * \brief The derivative of the relative permeability for the wetting phase w.r.t. saturation
     * \note This is only for compatibility. Will not be used.
     */
    Scalar dkrw_dsw(const Scalar sw) const
    { return 0; }

    /*!
     * \brief The relative permeability for the non-wetting phase
     * \note This is only for compatibility. Will not be used.
     */
    Scalar krn(const Scalar sw) const
    { return 1.0; }

    /*!
     * \brief The derivative of the relative permeability for the non-wetting phase w.r.t. saturation
     * \note This is only for compatibility. Will not be used.
     */
    Scalar dkrn_dsw(const Scalar sw) const
    { return 0.0; }

private:
    std::shared_ptr<typename LocalRules::Tetrahedron> localRulesForTetrahedron_;
    std::shared_ptr<typename LocalRules::Cube> localRulesForCube_;
    std::shared_ptr<typename LocalRules::Octahedron> localRulesForOctahedron_;
    std::shared_ptr<typename LocalRules::Icosahedron> localRulesForIcosahedron_;
    std::shared_ptr<typename LocalRules::Dodecahedron> localRulesForDodecahedron_;

    // TODO add more shapes here

    Pore::Shape shape_;
};

} // end namespace Dumux::FluidMatrix

#endif