Line | Branch | Exec | Source |
---|---|---|---|
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 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup Fluidmatrixinteractions | ||
10 | * \brief Pc- and Kr-sw curves based on monotone splines through given data points | ||
11 | */ | ||
12 | #ifndef DUMUX_MATERIAL_FLUIDMATRIX_DATA_SPLINE_MATERIAL_HH | ||
13 | #define DUMUX_MATERIAL_FLUIDMATRIX_DATA_SPLINE_MATERIAL_HH | ||
14 | |||
15 | #include <algorithm> | ||
16 | #include <dumux/common/parameters.hh> | ||
17 | #include <dumux/common/monotonecubicspline.hh> | ||
18 | #include <dumux/material/fluidmatrixinteractions/fluidmatrixinteraction.hh> | ||
19 | #include <dumux/material/fluidmatrixinteractions/2p/materiallaw.hh> | ||
20 | |||
21 | namespace Dumux::FluidMatrix { | ||
22 | |||
23 | /*! | ||
24 | * \ingroup Fluidmatrixinteractions | ||
25 | * \brief Pc- and Kr-sw curves based on monotone splines through given data points | ||
26 | * \tparam S the type for scalar numbers | ||
27 | * \tparam approximatePcSwInverse if this is set true, the | ||
28 | * spline approximates sw(pc) and evaluating pc(sw) needs spline inversion. | ||
29 | * if this is false, the spline approximates pc(sw) and evaluating | ||
30 | * sw(pc) needs spline inversion. Spline inversion is rather expensive | ||
31 | * since it has to be done numerically. | ||
32 | */ | ||
33 | template <class S, bool approximatePcSwInverse = false> | ||
34 | class DataSplineTwoPMaterialLaw | ||
35 | : public Adapter<DataSplineTwoPMaterialLaw<S, approximatePcSwInverse>, PcKrSw> | ||
36 | { | ||
37 | public: | ||
38 | |||
39 | using Scalar = S; | ||
40 | |||
41 | /*! | ||
42 | * \brief Return the number of fluid phases | ||
43 | */ | ||
44 | static constexpr int numFluidPhases() | ||
45 | { return 2; } | ||
46 | |||
47 | static constexpr bool isRegularized() | ||
48 | { return false; } | ||
49 | |||
50 | /*! | ||
51 | * \brief Deleted default constructor (so we are never in an undefined state) | ||
52 | * \note store owning pointers to laws instead if you need default-constructible objects | ||
53 | */ | ||
54 | DataSplineTwoPMaterialLaw() = delete; | ||
55 | |||
56 | /*! | ||
57 | * \brief Construct from a subgroup from the global parameter tree | ||
58 | */ | ||
59 | 1 | explicit DataSplineTwoPMaterialLaw(const std::string& paramGroup) | |
60 | 3 | { | |
61 | using V = std::vector<Scalar>; | ||
62 |
3/12✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 11 taken 1 times.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
|
2 | const std::string swPcGroup = paramGroup.empty() ? "Pc" : paramGroup + ".Pc"; |
63 |
2/6✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
|
2 | const auto swPc = getParamFromGroup<V>(swPcGroup, "SwData"); |
64 |
4/14✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 11 taken 1 times.
✗ Branch 12 not taken.
✓ Branch 13 taken 1 times.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
|
3 | const std::string krwPcGroup = paramGroup.empty() ? "Krw" : paramGroup + ".Krw"; |
65 |
2/6✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
|
2 | const auto swKrw = getParamFromGroup<V>(krwPcGroup, "SwData"); |
66 |
4/14✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 11 taken 1 times.
✗ Branch 12 not taken.
✓ Branch 13 taken 1 times.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
|
3 | const std::string krnPcGroup = paramGroup.empty() ? "Krn" : paramGroup + ".Krn"; |
67 |
2/6✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
|
2 | const auto swKrn = getParamFromGroup<V>(krnPcGroup, "SwData"); |
68 | |||
69 |
2/6✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
|
2 | const auto pc = getParamFromGroup<V>(paramGroup, "PcData"); |
70 |
2/6✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
|
2 | const auto krw = getParamFromGroup<V>(paramGroup, "KrwData"); |
71 |
3/8✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
|
3 | const auto krn = getParamFromGroup<V>(paramGroup, "KrnData"); |
72 | |||
73 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
1 | updateData_(swPc, pc, swKrw, krw, swKrn, krn); |
74 | 1 | } | |
75 | |||
76 | /*! | ||
77 | * \brief Construct from user data | ||
78 | */ | ||
79 | DataSplineTwoPMaterialLaw(const std::vector<Scalar>& swPc, | ||
80 | const std::vector<Scalar>& pc, | ||
81 | const std::vector<Scalar>& swKrw, | ||
82 | const std::vector<Scalar>& krw, | ||
83 | const std::vector<Scalar>& swKrn, | ||
84 | const std::vector<Scalar>& krn) | ||
85 | { | ||
86 | updateData_(swPc, pc, swKrw, krw, swKrn, krn); | ||
87 | } | ||
88 | |||
89 | /*! | ||
90 | * \brief The capillary pressure-saturation curve | ||
91 | */ | ||
92 | Scalar pc(const Scalar sw) const | ||
93 | { | ||
94 | if constexpr (approximatePcSwInverse) | ||
95 | return pcSpline_.evalInverse(sw); | ||
96 | else | ||
97 | ✗ | return pcSpline_.eval(sw); | |
98 | } | ||
99 | |||
100 | /*! | ||
101 | * \brief The partial derivative of the capillary pressure w.r.t. the saturation | ||
102 | */ | ||
103 | Scalar dpc_dsw(const Scalar sw) const | ||
104 | { | ||
105 | if constexpr (approximatePcSwInverse) | ||
106 | return 1.0/pcSpline_.evalDerivative(pcSpline_.evalInverse(sw)); | ||
107 | else | ||
108 | ✗ | return pcSpline_.evalDerivative(sw); | |
109 | } | ||
110 | |||
111 | /*! | ||
112 | * \brief The saturation-capillary pressure curve | ||
113 | */ | ||
114 | Scalar sw(const Scalar pc) const | ||
115 | { | ||
116 | if constexpr (approximatePcSwInverse) | ||
117 | return pcSpline_.eval(pc); | ||
118 | else | ||
119 | ✗ | return pcSpline_.evalInverse(pc); | |
120 | } | ||
121 | |||
122 | /*! | ||
123 | * \brief The partial derivative of the saturation to the capillary pressure | ||
124 | */ | ||
125 | ✗ | Scalar dsw_dpc(const Scalar pc) const | |
126 | { | ||
127 | if constexpr (approximatePcSwInverse) | ||
128 | return pcSpline_.evalDerivative(pc); | ||
129 | else | ||
130 | ✗ | return 1.0/pcSpline_.evalDerivative(pcSpline_.evalInverse(pc)); | |
131 | } | ||
132 | |||
133 | /*! | ||
134 | * \brief The relative permeability for the wetting phase | ||
135 | */ | ||
136 | Scalar krw(const Scalar sw) const | ||
137 | { | ||
138 | ✗ | return krwSpline_.eval(sw); | |
139 | } | ||
140 | |||
141 | /*! | ||
142 | * \brief The derivative of the relative permeability for the wetting phase w.r.t. saturation | ||
143 | */ | ||
144 | Scalar dkrw_dsw(const Scalar sw) const | ||
145 | { | ||
146 | ✗ | return krwSpline_.evalDerivative(sw); | |
147 | } | ||
148 | |||
149 | /*! | ||
150 | * \brief The relative permeability for the non-wetting phase | ||
151 | */ | ||
152 | Scalar krn(const Scalar sw) const | ||
153 | { | ||
154 | ✗ | return krnSpline_.eval(sw); | |
155 | } | ||
156 | |||
157 | /*! | ||
158 | * \brief The derivative of the relative permeability for the non-wetting phase w.r.t. saturation | ||
159 | */ | ||
160 | Scalar dkrn_dsw(const Scalar sw) const | ||
161 | { | ||
162 | ✗ | return krnSpline_.evalDerivative(sw); | |
163 | } | ||
164 | |||
165 | private: | ||
166 | 1 | void updateData_(const std::vector<Scalar>& swPc, | |
167 | const std::vector<Scalar>& pc, | ||
168 | const std::vector<Scalar>& swKrw, | ||
169 | const std::vector<Scalar>& krw, | ||
170 | const std::vector<Scalar>& swKrn, | ||
171 | const std::vector<Scalar>& krn) | ||
172 | { | ||
173 | if constexpr (approximatePcSwInverse) | ||
174 | pcSpline_.updatePoints(pc, swPc); | ||
175 | else | ||
176 | 1 | pcSpline_.updatePoints(swPc, pc); | |
177 | |||
178 | 1 | krwSpline_.updatePoints(swKrw, krw); | |
179 | 1 | krnSpline_.updatePoints(swKrn, krn); | |
180 | 1 | } | |
181 | |||
182 | MonotoneCubicSpline<Scalar> pcSpline_; | ||
183 | MonotoneCubicSpline<Scalar> krwSpline_; | ||
184 | MonotoneCubicSpline<Scalar> krnSpline_; | ||
185 | }; | ||
186 | |||
187 | } // end namespace Dumux::FluidMatrix | ||
188 | |||
189 | #endif | ||
190 |