opm-common
Loading...
Searching...
No Matches
SatCurveMultiplexer.hpp
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 This file is part of the Open Porous Media project (OPM).
5
6 OPM 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 OPM 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 OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef OPM_SAT_CURVE_MULTIPLEXER_HPP
28#define OPM_SAT_CURVE_MULTIPLEXER_HPP
29
31
32#include <stdexcept>
33
34
35// The static_assert does not compile with gcc 12 and earlier when placed in the multiplexer calls below.
36#if defined(__GNUC__) && (__GNUC__ < 13)
37 #define STATIC_ASSERT_SATCURVE_MULTIPLEXER_UNLESS_GCC_LT_13 throw std::logic_error("Unhandled SatCurveMultiplexerApproach")
38#else
39 #define STATIC_ASSERT_SATCURVE_MULTIPLEXER_UNLESS_GCC_LT_13 static_assert(false, "Unhandled SatCurveMultiplexerApproach")
40#endif
41
42namespace Opm {
50template <class TraitsT, class ParamsT = SatCurveMultiplexerParams<TraitsT> >
51class SatCurveMultiplexer : public TraitsT
52{
53public:
54 using Traits = TraitsT;
55 using Params = ParamsT;
56 using Scalar = typename Traits::Scalar;
57
58 using LETTwoPhaseLaw = TwoPhaseLETCurves<Traits>;
59 using PLTwoPhaseLaw = PiecewiseLinearTwoPhaseMaterial<Traits>;
60
62 using Traits::numPhases;
63 static_assert(numPhases == 2,
64 "The Brooks-Corey capillary pressure law only applies "
65 "to the case of two fluid phases");
66
69 static constexpr bool implementsTwoPhaseApi = true;
70
73 static constexpr bool implementsTwoPhaseSatApi = true;
74
77 static constexpr bool isSaturationDependent = true;
78
81 static constexpr bool isPressureDependent = false;
82
85 static constexpr bool isTemperatureDependent = false;
86
89 static constexpr bool isCompositionDependent = false;
90
91 static_assert(Traits::numPhases == 2,
92 "The number of fluid phases must be two if you want to use "
93 "this material law!");
94
98 template <class Container, class FluidState>
99 static void capillaryPressures(Container& values, const Params& params, const FluidState& fluidState)
100 {
101 switch (params.approach()) {
102 case SatCurveMultiplexerApproach::LET:
104 params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
105 fluidState);
106 break;
107
108 case SatCurveMultiplexerApproach::PiecewiseLinear:
110 params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
111 fluidState);
112 break;
113 }
114 }
115
120 template <class Container, class FluidState>
121 static void saturations(Container& values, const Params& params, const FluidState& fluidState)
122 {
123 switch (params.approach()) {
124 case SatCurveMultiplexerApproach::LET:
126 params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
127 fluidState);
128 break;
129
130 case SatCurveMultiplexerApproach::PiecewiseLinear:
132 params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
133 fluidState);
134 break;
135 }
136 }
137
148 template <class Container, class FluidState>
149 static void relativePermeabilities(Container& values, const Params& params, const FluidState& fluidState)
150 {
151 switch (params.approach()) {
152 case SatCurveMultiplexerApproach::LET:
154 params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
155 fluidState);
156 break;
157
158 case SatCurveMultiplexerApproach::PiecewiseLinear:
160 params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
161 fluidState);
162 break;
163 }
164 }
165
169 template <class FluidState, class Evaluation = typename FluidState::Scalar>
170 static Evaluation pcnw(const Params& params, const FluidState& fluidState)
171 {
172 switch (params.approach()) {
173 case SatCurveMultiplexerApproach::LET:
174 return LETTwoPhaseLaw::pcnw(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
175 fluidState);
176 break;
177
178 case SatCurveMultiplexerApproach::PiecewiseLinear:
179 return PLTwoPhaseLaw::pcnw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
180 fluidState);
181 break;
182 }
183
184 return 0.0;
185 }
186
187 template <class Evaluation, class ...Args>
188 static Evaluation twoPhaseSatPcnw(const Params& params, const Evaluation& Sw)
189 {
190 if constexpr (FrontIsSatCurveMultiplexerDispatchV<Args...>) {
191 return twoPhaseSatPcnwT<Evaluation, Args...>(params, Sw);
192 }
193
194 switch (params.approach()) {
195 case SatCurveMultiplexerApproach::LET:
196 return LETTwoPhaseLaw::twoPhaseSatPcnw(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
197 Sw);
198 break;
199
200 case SatCurveMultiplexerApproach::PiecewiseLinear:
201 return PLTwoPhaseLaw::twoPhaseSatPcnw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
202 Sw);
203 break;
204 }
205
206 return 0.0;
207 }
208
209 template <class Evaluation, class Head, class ...Args>
210 static Evaluation twoPhaseSatPcnwT(const Params& params, const Evaluation& Sw)
211 {
212 if constexpr (Head::approach == SatCurveMultiplexerApproach::LET) {
213 return LETTwoPhaseLaw::twoPhaseSatPcnw(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
214 Sw);
215 } else if constexpr (Head::approach == SatCurveMultiplexerApproach::PiecewiseLinear) {
216 return PLTwoPhaseLaw::twoPhaseSatPcnw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
217 Sw);
218 } else {
219 STATIC_ASSERT_SATCURVE_MULTIPLEXER_UNLESS_GCC_LT_13;
220 }
221 }
222
223 template <class Evaluation>
224 static Evaluation twoPhaseSatPcnwInv(const Params&, const Evaluation&)
225 {
226 throw std::logic_error("SatCurveMultiplexer::twoPhaseSatPcnwInv"
227 " not implemented!");
228 }
229
233 template <class FluidState, class Evaluation = typename FluidState::Scalar>
234 static Evaluation Sw(const Params& params, const FluidState& fluidstate)
235 {
236 switch (params.approach()) {
237 case SatCurveMultiplexerApproach::LET:
238 return LETTwoPhaseLaw::Sw(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
239 fluidstate);
240 break;
241
242 case SatCurveMultiplexerApproach::PiecewiseLinear:
243 return PLTwoPhaseLaw::Sw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
244 fluidstate);
245 break;
246 }
247
248 return 0.0;
249 }
250
251 template <class Evaluation>
252 static Evaluation twoPhaseSatSw(const Params&, const Evaluation&)
253 {
254 throw std::logic_error("SatCurveMultiplexer::twoPhaseSatSw"
255 " not implemented!");
256 }
257
262 template <class FluidState, class Evaluation = typename FluidState::Scalar>
263 static Evaluation Sn(const Params& params, const FluidState& fluidstate)
264 {
265 switch (params.approach()) {
266 case SatCurveMultiplexerApproach::LET:
267 return LETTwoPhaseLaw::Sn(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
268 fluidstate);
269 break;
270
271 case SatCurveMultiplexerApproach::PiecewiseLinear:
272 return PLTwoPhaseLaw::Sn(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
273 fluidstate);
274 break;
275 }
276
277 return 0.0;
278 }
279
280 template <class Evaluation>
281 static Evaluation twoPhaseSatSn(const Params& params, const Evaluation& pc)
282 {
283 switch (params.approach()) {
284 case SatCurveMultiplexerApproach::LET:
285 return LETTwoPhaseLaw::twoPhaseSatSn(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
286 pc);
287 break;
288
289 case SatCurveMultiplexerApproach::PiecewiseLinear:
290 return PLTwoPhaseLaw::twoPhaseSatSn(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
291 pc);
292 break;
293 }
294
295 return 0.0;
296 }
297
302 template <class FluidState, class Evaluation = typename FluidState::Scalar>
303 static Evaluation krw(const Params& params, const FluidState& fluidstate)
304 {
305 switch (params.approach()) {
306 case SatCurveMultiplexerApproach::LET:
307 return LETTwoPhaseLaw::krw(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
308 fluidstate);
309 break;
310
311 case SatCurveMultiplexerApproach::PiecewiseLinear:
312 return PLTwoPhaseLaw::krw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
313 fluidstate);
314 break;
315 }
316
317 return 0.0;
318 }
319
320 template <class Evaluation, class ...Args>
321 static Evaluation twoPhaseSatKrw(const Params& params, const Evaluation& Sw)
322 {
323 if constexpr (FrontIsSatCurveMultiplexerDispatchV<Args...>) {
324 return twoPhaseSatKrwT<Evaluation, Args...>(params, Sw);
325 }
326
327 switch (params.approach()) {
328 case SatCurveMultiplexerApproach::LET:
329 return LETTwoPhaseLaw::twoPhaseSatKrw(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
330 Sw);
331 break;
332
333 case SatCurveMultiplexerApproach::PiecewiseLinear:
334 return PLTwoPhaseLaw::twoPhaseSatKrw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
335 Sw);
336 break;
337 }
338
339 return 0.0;
340 }
341
342 template <class Evaluation, class Head, class ...Args>
343 static Evaluation twoPhaseSatKrwT(const Params& params, const Evaluation& Sw)
344 {
345 if constexpr (Head::approach == SatCurveMultiplexerApproach::LET) {
346 return LETTwoPhaseLaw::twoPhaseSatKrw(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
347 Sw);
348 } else if constexpr (Head::approach == SatCurveMultiplexerApproach::PiecewiseLinear) {
349 return PLTwoPhaseLaw::twoPhaseSatKrw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
350 Sw);
351 } else {
352 STATIC_ASSERT_SATCURVE_MULTIPLEXER_UNLESS_GCC_LT_13;
353 }
354 }
355
356 template <class Evaluation>
357 static Evaluation twoPhaseSatKrwInv(const Params& params, const Evaluation& krw)
358 {
359 switch (params.approach()) {
360 case SatCurveMultiplexerApproach::LET:
361 return LETTwoPhaseLaw::twoPhaseSatKrwInv(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
362 krw);
363 break;
364
365 case SatCurveMultiplexerApproach::PiecewiseLinear:
366 return PLTwoPhaseLaw::twoPhaseSatKrwInv(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
367 krw);
368 break;
369 }
370
371 return 0.0;
372 }
373
378 template <class FluidState, class Evaluation = typename FluidState::Scalar>
379 static Evaluation krn(const Params& params, const FluidState& fluidstate)
380 {
381 switch (params.approach()) {
382 case SatCurveMultiplexerApproach::LET:
383 return LETTwoPhaseLaw::krn(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
384 fluidstate);
385 break;
386
387 case SatCurveMultiplexerApproach::PiecewiseLinear:
388 return PLTwoPhaseLaw::krn(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
389 fluidstate);
390 break;
391 }
392
393 return 0.0;
394 }
395
396 template <class Evaluation, class ...Args>
397 static Evaluation twoPhaseSatKrn(const Params& params, const Evaluation& Sw)
398 {
399 if constexpr (FrontIsSatCurveMultiplexerDispatchV<Args...>) {
400 return twoPhaseSatKrnT<Evaluation, Args...>(params, Sw);
401 }
402
403 switch (params.approach()) {
404 case SatCurveMultiplexerApproach::LET:
405 return LETTwoPhaseLaw::twoPhaseSatKrn(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
406 Sw);
407 break;
408
409 case SatCurveMultiplexerApproach::PiecewiseLinear:
410 return PLTwoPhaseLaw::twoPhaseSatKrn(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
411 Sw);
412 break;
413 }
414
415 return 0.0;
416 }
417
418 template <class Evaluation, class Head, class ...Args>
419 static Evaluation twoPhaseSatKrnT(const Params& params, const Evaluation& Sw)
420 {
421 if constexpr (Head::approach == SatCurveMultiplexerApproach::LET) {
422 return LETTwoPhaseLaw::twoPhaseSatKrn(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
423 Sw);
424 } else if constexpr (Head::approach == SatCurveMultiplexerApproach::PiecewiseLinear) {
425 return PLTwoPhaseLaw::twoPhaseSatKrn(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
426 Sw);
427 } else {
428 STATIC_ASSERT_SATCURVE_MULTIPLEXER_UNLESS_GCC_LT_13;
429 }
430 }
431
432 template <class Evaluation>
433 static Evaluation twoPhaseSatKrnInv(const Params& params, const Evaluation& krn)
434 {
435 switch (params.approach()) {
436 case SatCurveMultiplexerApproach::LET:
437 return LETTwoPhaseLaw::twoPhaseSatKrnInv(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
438 krn);
439 break;
440
441 case SatCurveMultiplexerApproach::PiecewiseLinear:
442 return PLTwoPhaseLaw::twoPhaseSatKrnInv(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
443 krn);
444 break;
445 }
446
447 return 0.0;
448 }
449};
450
451} // namespace Opm
452
453#endif // OPM_SAT_CURVE_MULTIPLEXER_HPP
Specification of the material parameters for the saturation function multiplexer.
Implementation of a tabulated, piecewise linear capillary pressure law.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:53
static OPM_HOST_DEVICE void capillaryPressures(Container &values, const Params &params, const FluidState &fs)
Definition PiecewiseLinearTwoPhaseMaterial.hpp:100
static OPM_HOST_DEVICE Evaluation Sn(const Params &params, const FluidState &fs)
Definition PiecewiseLinearTwoPhaseMaterial.hpp:173
static OPM_HOST_DEVICE void relativePermeabilities(Container &values, const Params &params, const FluidState &fs)
Definition PiecewiseLinearTwoPhaseMaterial.hpp:121
static OPM_HOST_DEVICE Evaluation twoPhaseSatPcnw(const Params &params, const Evaluation &Sw)
Definition PiecewiseLinearTwoPhaseMaterial.hpp:147
static OPM_HOST_DEVICE Evaluation krw(const Params &params, const FluidState &fs)
Definition PiecewiseLinearTwoPhaseMaterial.hpp:185
static OPM_HOST_DEVICE void saturations(Container &, const Params &, const FluidState &)
Definition PiecewiseLinearTwoPhaseMaterial.hpp:114
static OPM_HOST_DEVICE Evaluation krn(const Params &params, const FluidState &fs)
Definition PiecewiseLinearTwoPhaseMaterial.hpp:213
static OPM_HOST_DEVICE Evaluation Sw(const Params &, const FluidState &)
Definition PiecewiseLinearTwoPhaseMaterial.hpp:161
static OPM_HOST_DEVICE Evaluation pcnw(const Params &params, const FluidState &fs)
Definition PiecewiseLinearTwoPhaseMaterial.hpp:134
Implements a multiplexer class that provides LET curves and piecewise linear saturation functions.
Definition SatCurveMultiplexer.hpp:52
static Evaluation Sw(const Params &params, const FluidState &fluidstate)
The saturation-capillary pressure curve.
Definition SatCurveMultiplexer.hpp:234
static Evaluation krn(const Params &params, const FluidState &fluidstate)
The relative permeability for the non-wetting phase of the medium.
Definition SatCurveMultiplexer.hpp:379
static Evaluation krw(const Params &params, const FluidState &fluidstate)
The relative permeability for the wetting phase of the medium.
Definition SatCurveMultiplexer.hpp:303
static constexpr bool isPressureDependent
Definition SatCurveMultiplexer.hpp:81
static void saturations(Container &values, const Params &params, const FluidState &fluidState)
Calculate the saturations of the phases starting from their pressure differences.
Definition SatCurveMultiplexer.hpp:121
static constexpr bool implementsTwoPhaseSatApi
Definition SatCurveMultiplexer.hpp:73
static constexpr bool isTemperatureDependent
Definition SatCurveMultiplexer.hpp:85
static constexpr bool isCompositionDependent
Definition SatCurveMultiplexer.hpp:89
static void capillaryPressures(Container &values, const Params &params, const FluidState &fluidState)
The capillary pressure-saturation curves.
Definition SatCurveMultiplexer.hpp:99
static Evaluation Sn(const Params &params, const FluidState &fluidstate)
Calculate the non-wetting phase saturations depending on the phase pressures.
Definition SatCurveMultiplexer.hpp:263
static void relativePermeabilities(Container &values, const Params &params, const FluidState &fluidState)
The relative permeability-saturation curves.
Definition SatCurveMultiplexer.hpp:149
static constexpr bool implementsTwoPhaseApi
Definition SatCurveMultiplexer.hpp:69
static constexpr bool isSaturationDependent
Definition SatCurveMultiplexer.hpp:77
static Evaluation pcnw(const Params &params, const FluidState &fluidState)
The capillary pressure-saturation curve.
Definition SatCurveMultiplexer.hpp:170
Implementation of the LET curve saturation functions.
Definition TwoPhaseLETCurves.hpp:50
static Evaluation krn(const Params &, const FluidState &)
Definition TwoPhaseLETCurves.hpp:237
static Evaluation pcnw(const Params &, const FluidState &)
Definition TwoPhaseLETCurves.hpp:130
static Evaluation krw(const Params &, const FluidState &)
Definition TwoPhaseLETCurves.hpp:193
static void saturations(Container &, const Params &, const FluidState &)
Definition TwoPhaseLETCurves.hpp:103
static void capillaryPressures(Container &, const Params &, const FluidState &)
Definition TwoPhaseLETCurves.hpp:93
static void relativePermeabilities(Container &, const Params &, const FluidState &)
Definition TwoPhaseLETCurves.hpp:119
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30