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-FileCopyrightText: 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 Implementation helper for capillary pressure and | ||
11 | * relative permeability <-> saturation relations for two-phase models | ||
12 | */ | ||
13 | #ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWOP_MATERIAL_LAW_HH | ||
14 | #define DUMUX_MATERIAL_FLUIDMATRIX_TWOP_MATERIAL_LAW_HH | ||
15 | |||
16 | #include <dumux/common/parameters.hh> | ||
17 | #include <dumux/material/fluidmatrixinteractions/fluidmatrixinteraction.hh> | ||
18 | #include <dumux/material/fluidmatrixinteractions/2p/efftoabsdefaultpolicy.hh> | ||
19 | #include <dumux/material/fluidmatrixinteractions/2p/noregularization.hh> | ||
20 | |||
21 | namespace Dumux::FluidMatrix { | ||
22 | |||
23 | /*! | ||
24 | * \ingroup Fluidmatrixinteractions | ||
25 | * \brief Wrapper class to implement regularized material laws (pc-sw, kr-sw) | ||
26 | * with a conversion policy between absolution and effective saturations | ||
27 | * \note See vangenuchten.hh / brookscorey.hh for default configurations using this class | ||
28 | * \tparam ScalarType the scalar type | ||
29 | * \tparam BaseLaw the base law (e.g. VanGenuchten, BrooksCorey, Linear, ...) | ||
30 | * \tparam Regularization the regularization type (set to NoRegularization to turn it off) | ||
31 | * \tparam EffToAbsPolicy the policy how to convert effective <-> absolute saturations | ||
32 | * | ||
33 | * \note The regularization interface is expected to return Dumux::OptionalScalars which | ||
34 | * are wrappers around a Scalar type that provide a boolean operator to | ||
35 | * check whether the result is valid. If the regularization returns | ||
36 | * a non-valid value, it means that the given parameter | ||
37 | * range is outside the regularized region. | ||
38 | * For that case we forward to the call to the standard law. | ||
39 | */ | ||
40 | template<class ScalarType, | ||
41 | class BaseLaw, | ||
42 | class Regularization = NoRegularization, | ||
43 | class EffToAbsPolicy = TwoPEffToAbsDefaultPolicy> | ||
44 | class TwoPMaterialLaw : public Adapter<TwoPMaterialLaw<ScalarType, BaseLaw, Regularization, EffToAbsPolicy>, PcKrSw> | ||
45 | { | ||
46 | public: | ||
47 | |||
48 | using Scalar = ScalarType; | ||
49 | |||
50 | using BasicParams = typename BaseLaw::template Params<Scalar>; | ||
51 | using EffToAbsParams = typename EffToAbsPolicy::template Params<Scalar>; | ||
52 | using RegularizationParams = typename Regularization::template Params<Scalar>; | ||
53 | |||
54 | using EffToAbs = EffToAbsPolicy; | ||
55 | |||
56 | /*! | ||
57 | * \brief Return the number of fluid phases | ||
58 | */ | ||
59 | static constexpr int numFluidPhases() | ||
60 | { return 2; } | ||
61 | |||
62 | /*! | ||
63 | * \brief Return whether this law is regularized | ||
64 | */ | ||
65 | static constexpr bool isRegularized() | ||
66 | { return !std::is_same<Regularization, NoRegularization>::value; } | ||
67 | |||
68 | /*! | ||
69 | * \brief Deleted default constructor (so we are never in an undefined state) | ||
70 | * \note store owning pointers to laws instead if you need default-constructible objects | ||
71 | */ | ||
72 | TwoPMaterialLaw() = delete; | ||
73 | |||
74 | /*! | ||
75 | * \brief Construct from a subgroup from the global parameter tree | ||
76 | * \note This will give you nice error messages if a mandatory parameter is missing | ||
77 | */ | ||
78 | 219 | explicit TwoPMaterialLaw(const std::string& paramGroup) | |
79 | 219 | : basicParams_(makeBasicParams(paramGroup)) | |
80 | 219 | , effToAbsParams_(makeEffToAbsParams(paramGroup)) | |
81 | { | ||
82 | if constexpr (isRegularized()) | ||
83 | 196 | regularization_.init(this, paramGroup); | |
84 | 219 | } | |
85 | |||
86 | /*! | ||
87 | * \brief Construct from parameter structs | ||
88 | * \note More efficient constructor but you need to ensure all parameters are initialized | ||
89 | */ | ||
90 | 15 | explicit TwoPMaterialLaw(const BasicParams& baseParams, | |
91 | const EffToAbsParams& effToAbsParams = {}, | ||
92 | const RegularizationParams& regParams = {}) | ||
93 | 15 | : basicParams_(baseParams) | |
94 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
|
14 | , effToAbsParams_(effToAbsParams) |
95 | { | ||
96 | if constexpr (isRegularized()) | ||
97 | 12 | regularization_.init(this, baseParams, effToAbsParams, regParams); | |
98 | 1 | } | |
99 | |||
100 | /*! | ||
101 | * \brief The capillary pressure-saturation curve | ||
102 | */ | ||
103 | template<bool enableRegularization = isRegularized()> | ||
104 | 170558547 | Scalar pc(const Scalar sw) const | |
105 | { | ||
106 |
3/4✗ Branch 0 not taken.
✓ Branch 1 taken 892926 times.
✓ Branch 2 taken 130400335 times.
✓ Branch 3 taken 635 times.
|
170558306 | const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_); |
107 | if constexpr (enableRegularization) | ||
108 | { | ||
109 |
3/3✓ Branch 1 taken 156649219 times.
✓ Branch 2 taken 9646071 times.
✓ Branch 0 taken 892992 times.
|
167188282 | const auto regularized = regularization_.pc(swe); |
110 |
2/2✓ Branch 0 taken 135375739 times.
✓ Branch 1 taken 31812543 times.
|
167188282 | if (regularized) |
111 | 135375739 | return regularized.value(); | |
112 | } | ||
113 | |||
114 |
0/4✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
35182567 | return BaseLaw::pc(swe, basicParams_); |
115 | } | ||
116 | |||
117 | /*! | ||
118 | * \brief The partial derivative of the capillary pressure w.r.t. the saturation | ||
119 | */ | ||
120 | template<bool enableRegularization = isRegularized()> | ||
121 | 13310452 | Scalar dpc_dsw(const Scalar sw) const | |
122 | { | ||
123 |
2/2✓ Branch 0 taken 22 times.
✓ Branch 1 taken 178 times.
|
13310251 | const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_); |
124 | if constexpr (enableRegularization) | ||
125 | { | ||
126 |
3/3✓ Branch 0 taken 22 times.
✓ Branch 1 taken 1651878 times.
✓ Branch 2 taken 92116 times.
|
1744016 | const auto regularized = regularization_.dpc_dswe(swe); |
127 |
2/2✓ Branch 0 taken 1651742 times.
✓ Branch 1 taken 92274 times.
|
1744016 | if (regularized) |
128 | 1651742 | return regularized.value()*EffToAbs::dswe_dsw(effToAbsParams_); | |
129 | } | ||
130 | |||
131 | 11658509 | return BaseLaw::dpc_dswe(swe, basicParams_)*EffToAbs::dswe_dsw(effToAbsParams_); | |
132 | } | ||
133 | |||
134 | /*! | ||
135 | * \brief The capillary pressure at Swe = 1.0 also called end point capillary pressure | ||
136 | */ | ||
137 | 10109616 | Scalar endPointPc() const | |
138 | { | ||
139 |
4/4✓ Branch 0 taken 4363952 times.
✓ Branch 1 taken 4420067 times.
✓ Branch 2 taken 8332 times.
✓ Branch 3 taken 564597 times.
|
91765595 | return BaseLaw::endPointPc(basicParams_); |
140 | } | ||
141 | |||
142 | /*! | ||
143 | * \brief The saturation-capillary pressure curve | ||
144 | */ | ||
145 | template<bool enableRegularization = isRegularized()> | ||
146 | 790206031 | Scalar sw(const Scalar pc) const | |
147 | { | ||
148 | if constexpr (enableRegularization) | ||
149 | { | ||
150 |
2/2✓ Branch 0 taken 40 times.
✓ Branch 1 taken 39122 times.
|
810672894 | const auto regularized = regularization_.swe(pc); |
151 |
2/2✓ Branch 0 taken 23156005 times.
✓ Branch 1 taken 787516889 times.
|
810672894 | if (regularized) |
152 | 23156005 | return EffToAbs::sweToSw(regularized.value(), effToAbsParams_); | |
153 | } | ||
154 | |||
155 |
1/3✗ Branch 0 not taken.
✓ Branch 1 taken 30006 times.
✗ Branch 2 not taken.
|
800602442 | return EffToAbs::sweToSw(BaseLaw::swe(pc, basicParams_), effToAbsParams_); |
156 | } | ||
157 | |||
158 | /*! | ||
159 | * \brief The partial derivative of the saturation to the capillary pressure | ||
160 | */ | ||
161 | template<bool enableRegularization = isRegularized()> | ||
162 | 42607504 | Scalar dsw_dpc(const Scalar pc) const | |
163 | { | ||
164 | if constexpr (enableRegularization) | ||
165 | { | ||
166 |
2/2✓ Branch 0 taken 20 times.
✓ Branch 1 taken 180 times.
|
41607502 | const auto regularized = regularization_.dswe_dpc(pc); |
167 |
2/2✓ Branch 0 taken 4246933 times.
✓ Branch 1 taken 37360569 times.
|
41607502 | if (regularized) |
168 | 4246933 | return regularized.value()*EffToAbs::dsw_dswe(effToAbsParams_); | |
169 | } | ||
170 | |||
171 | 38360571 | return BaseLaw::dswe_dpc(pc, basicParams_)*EffToAbs::dsw_dswe(effToAbsParams_); | |
172 | } | ||
173 | |||
174 | /*! | ||
175 | * \brief The relative permeability for the wetting phase | ||
176 | */ | ||
177 | template<bool enableRegularization = isRegularized()> | ||
178 | 1006845415 | Scalar krw(const Scalar sw) const | |
179 | { | ||
180 |
2/2✓ Branch 0 taken 127880716 times.
✓ Branch 1 taken 7575680 times.
|
1006845415 | const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_); |
181 | if constexpr (enableRegularization) | ||
182 | { | ||
183 |
3/3✓ Branch 1 taken 50956109 times.
✓ Branch 2 taken 793425394 times.
✓ Branch 0 taken 127880716 times.
|
972262219 | const auto regularized = regularization_.krw(swe); |
184 |
2/2✓ Branch 0 taken 158001175 times.
✓ Branch 1 taken 814261044 times.
|
972262219 | if (regularized) |
185 | 158001175 | return regularized.value(); | |
186 | } | ||
187 | |||
188 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 7122944 times.
|
848844240 | return BaseLaw::krw(swe, basicParams_); |
189 | } | ||
190 | |||
191 | /*! | ||
192 | * \brief The derivative of the relative permeability for the wetting phase w.r.t. saturation | ||
193 | */ | ||
194 | template<bool enableRegularization = isRegularized()> | ||
195 | 35308052 | Scalar dkrw_dsw(const Scalar sw) const | |
196 | { | ||
197 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 90 times.
|
35308052 | const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_); |
198 | if constexpr (enableRegularization) | ||
199 | { | ||
200 |
3/3✓ Branch 0 taken 10 times.
✓ Branch 1 taken 5740388 times.
✓ Branch 2 taken 28567532 times.
|
34307930 | const auto regularized = regularization_.dkrw_dswe(swe); |
201 |
2/2✓ Branch 0 taken 5740318 times.
✓ Branch 1 taken 28567612 times.
|
34307930 | if (regularized) |
202 | 5740318 | return regularized.value()*EffToAbs::dswe_dsw(effToAbsParams_); | |
203 | } | ||
204 | |||
205 | 29567734 | return BaseLaw::dkrw_dswe(swe, basicParams_)*EffToAbs::dswe_dsw(effToAbsParams_); | |
206 | } | ||
207 | |||
208 | /*! | ||
209 | * \brief The relative permeability for the non-wetting phase | ||
210 | */ | ||
211 | template<bool enableRegularization = isRegularized()> | ||
212 | 163544120 | Scalar krn(const Scalar sw) const | |
213 | { | ||
214 |
2/2✓ Branch 0 taken 127880716 times.
✓ Branch 1 taken 452736 times.
|
163544120 | const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_); |
215 | if constexpr (enableRegularization) | ||
216 | { | ||
217 |
3/3✓ Branch 1 taken 18559461 times.
✓ Branch 2 taken 15734209 times.
✓ Branch 0 taken 127880716 times.
|
162174386 | const auto regularized = regularization_.krn(swe); |
218 |
2/2✓ Branch 0 taken 125604527 times.
✓ Branch 1 taken 36569859 times.
|
162174386 | if (regularized) |
219 | 125604527 | return regularized.value(); | |
220 | } | ||
221 | |||
222 |
0/2✗ Branch 0 not taken.
✗ Branch 1 not taken.
|
37939593 | return BaseLaw::krn(swe, basicParams_); |
223 | } | ||
224 | |||
225 | /*! | ||
226 | * \brief The derivative of the relative permeability for the non-wetting phase w.r.t. saturation | ||
227 | */ | ||
228 | template<bool enableRegularization = isRegularized()> | ||
229 | 2190978 | Scalar dkrn_dsw(const Scalar sw) const | |
230 | { | ||
231 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 90 times.
|
2190978 | const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_); |
232 | if constexpr (enableRegularization) | ||
233 | { | ||
234 |
3/3✓ Branch 0 taken 10 times.
✓ Branch 1 taken 1067304 times.
✓ Branch 2 taken 123542 times.
|
1190856 | const auto regularized = regularization_.dkrn_dswe(swe); |
235 |
2/2✓ Branch 0 taken 1067234 times.
✓ Branch 1 taken 123622 times.
|
1190856 | if (regularized) |
236 | 1067234 | return regularized.value()*EffToAbs::dswe_dsw(effToAbsParams_); | |
237 | } | ||
238 | |||
239 |
2/3✓ Branch 1 taken 111 times.
✓ Branch 2 taken 9 times.
✗ Branch 3 not taken.
|
1123744 | return BaseLaw::dkrn_dswe(swe, basicParams_)*EffToAbs::dswe_dsw(effToAbsParams_); |
240 | } | ||
241 | |||
242 | /*! | ||
243 | * \brief Equality comparison with another instance | ||
244 | */ | ||
245 | 5038371 | bool operator== (const TwoPMaterialLaw& o) const | |
246 | { | ||
247 | 5038371 | return basicParams_ == o.basicParams_ | |
248 |
1/2✓ Branch 0 taken 1991433 times.
✗ Branch 1 not taken.
|
1991433 | && effToAbsParams_ == o.effToAbsParams_ |
249 |
2/2✓ Branch 0 taken 1991433 times.
✓ Branch 1 taken 3046938 times.
|
7029804 | && regularization_ == o.regularization_; |
250 | } | ||
251 | |||
252 | /*! | ||
253 | * \brief Create the base law's parameters using | ||
254 | * input file parameters | ||
255 | */ | ||
256 | 207 | static BasicParams makeBasicParams(const std::string& paramGroup) | |
257 | { | ||
258 |
1/2✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
|
207 | return BaseLaw::template makeParams<Scalar>(paramGroup); |
259 | } | ||
260 | |||
261 | /*! | ||
262 | * \brief Return the base law's parameters | ||
263 | */ | ||
264 | const BasicParams& basicParams() const | ||
265 | { return basicParams_; } | ||
266 | |||
267 | /*! | ||
268 | * \brief Create the parameters of the EffToAbs policy using | ||
269 | * input file parameters | ||
270 | */ | ||
271 | 207 | static EffToAbsParams makeEffToAbsParams(const std::string& paramGroup) | |
272 | { | ||
273 |
1/2✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
|
207 | return EffToAbs::template makeParams<Scalar>(paramGroup); |
274 | } | ||
275 | |||
276 | /*! | ||
277 | * \brief Return the parameters of the EffToAbs policy | ||
278 | */ | ||
279 | 1 | const EffToAbsParams& effToAbsParams() const | |
280 | 326 | { return effToAbsParams_; } | |
281 | |||
282 | private: | ||
283 | BasicParams basicParams_; | ||
284 | EffToAbsParams effToAbsParams_; | ||
285 | Regularization regularization_; | ||
286 | }; | ||
287 | |||
288 | } // end namespace Dumux::FluidMatrix | ||
289 | |||
290 | #endif | ||
291 |