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
34namespace Opm {
42template <class TraitsT, class ParamsT = SatCurveMultiplexerParams<TraitsT> >
43class SatCurveMultiplexer : public TraitsT
44{
45public:
46 using Traits = TraitsT;
47 using Params = ParamsT;
48 using Scalar = typename Traits::Scalar;
49
52
54 static constexpr int numPhases = Traits::numPhases;
55 static_assert(numPhases == 2,
56 "The Brooks-Corey capillary pressure law only applies "
57 "to the case of two fluid phases");
58
61 static constexpr bool implementsTwoPhaseApi = true;
62
65 static constexpr bool implementsTwoPhaseSatApi = true;
66
69 static constexpr bool isSaturationDependent = true;
70
73 static constexpr bool isPressureDependent = false;
74
77 static constexpr bool isTemperatureDependent = false;
78
81 static constexpr bool isCompositionDependent = false;
82
83 static_assert(Traits::numPhases == 2,
84 "The number of fluid phases must be two if you want to use "
85 "this material law!");
86
90 template <class Container, class FluidState>
91 static void capillaryPressures(Container& values, const Params& params, const FluidState& fluidState)
92 {
93 switch (params.approach()) {
96 params.template getRealParams<SatCurveMultiplexerApproach::LETApproach>(),
97 fluidState);
98 break;
99
102 params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinearApproach>(),
103 fluidState);
104 break;
105 }
106 }
107
112 template <class Container, class FluidState>
113 static void saturations(Container& values, const Params& params, const FluidState& fluidState)
114 {
115 switch (params.approach()) {
118 params.template getRealParams<SatCurveMultiplexerApproach::LETApproach>(),
119 fluidState);
120 break;
121
124 params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinearApproach>(),
125 fluidState);
126 break;
127 }
128 }
129
140 template <class Container, class FluidState>
141 static void relativePermeabilities(Container& values, const Params& params, const FluidState& fluidState)
142 {
143 switch (params.approach()) {
146 params.template getRealParams<SatCurveMultiplexerApproach::LETApproach>(),
147 fluidState);
148 break;
149
152 params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinearApproach>(),
153 fluidState);
154 break;
155 }
156 }
157
161 template <class FluidState, class Evaluation = typename FluidState::Scalar>
162 static Evaluation pcnw(const Params& params, const FluidState& fluidState)
163 {
164 switch (params.approach()) {
166 return LETTwoPhaseLaw::pcnw(params.template getRealParams<SatCurveMultiplexerApproach::LETApproach>(),
167 fluidState);
168 break;
169
171 return PLTwoPhaseLaw::pcnw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinearApproach>(),
172 fluidState);
173 break;
174 }
175
176 return 0.0;
177 }
178
179 template <class Evaluation>
180 static Evaluation twoPhaseSatPcnw(const Params& params, const Evaluation& Sw)
181 {
182 switch (params.approach()) {
184 return LETTwoPhaseLaw::twoPhaseSatPcnw(params.template getRealParams<SatCurveMultiplexerApproach::LETApproach>(),
185 Sw);
186 break;
187
189 return PLTwoPhaseLaw::twoPhaseSatPcnw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinearApproach>(),
190 Sw);
191 break;
192 }
193
194 return 0.0;
195 }
196
197 template <class Evaluation>
198 static Evaluation twoPhaseSatPcnwInv(const Params&, const Evaluation&)
199 {
200 throw std::logic_error("SatCurveMultiplexer::twoPhaseSatPcnwInv"
201 " not implemented!");
202 }
203
207 template <class FluidState, class Evaluation = typename FluidState::Scalar>
208 static Evaluation Sw(const Params& params, const FluidState& fluidstate)
209 {
210 switch (params.approach()) {
212 return LETTwoPhaseLaw::Sw(params.template getRealParams<SatCurveMultiplexerApproach::LETApproach>(),
213 fluidstate);
214 break;
215
217 return PLTwoPhaseLaw::Sw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinearApproach>(),
218 fluidstate);
219 break;
220 }
221
222 return 0.0;
223 }
224
225 template <class Evaluation>
226 static Evaluation twoPhaseSatSw(const Params&, const Evaluation&)
227 {
228 throw std::logic_error("SatCurveMultiplexer::twoPhaseSatSw"
229 " not implemented!");
230 }
231
236 template <class FluidState, class Evaluation = typename FluidState::Scalar>
237 static Evaluation Sn(const Params& params, const FluidState& fluidstate)
238 {
239 switch (params.approach()) {
241 return LETTwoPhaseLaw::Sn(params.template getRealParams<SatCurveMultiplexerApproach::LETApproach>(),
242 fluidstate);
243 break;
244
246 return PLTwoPhaseLaw::Sn(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinearApproach>(),
247 fluidstate);
248 break;
249 }
250
251 return 0.0;
252 }
253
254 template <class Evaluation>
255 static Evaluation twoPhaseSatSn(const Params& params, const Evaluation& pc)
256 {
257 switch (params.approach()) {
259 return LETTwoPhaseLaw::twoPhaseSatSn(params.template getRealParams<SatCurveMultiplexerApproach::LETApproach>(),
260 pc);
261 break;
262
264 return PLTwoPhaseLaw::twoPhaseSatSn(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinearApproach>(),
265 pc);
266 break;
267 }
268
269 return 0.0;
270 }
271
276 template <class FluidState, class Evaluation = typename FluidState::Scalar>
277 static Evaluation krw(const Params& params, const FluidState& fluidstate)
278 {
279 switch (params.approach()) {
281 return LETTwoPhaseLaw::krw(params.template getRealParams<SatCurveMultiplexerApproach::LETApproach>(),
282 fluidstate);
283 break;
284
286 return PLTwoPhaseLaw::krw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinearApproach>(),
287 fluidstate);
288 break;
289 }
290
291 return 0.0;
292 }
293
294 template <class Evaluation>
295 static Evaluation twoPhaseSatKrw(const Params& params, const Evaluation& Sw)
296 {
297 switch (params.approach()) {
299 return LETTwoPhaseLaw::twoPhaseSatKrw(params.template getRealParams<SatCurveMultiplexerApproach::LETApproach>(),
300 Sw);
301 break;
302
304 return PLTwoPhaseLaw::twoPhaseSatKrw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinearApproach>(),
305 Sw);
306 break;
307 }
308
309 return 0.0;
310 }
311
312 template <class Evaluation>
313 static Evaluation twoPhaseSatKrwInv(const Params&, const Evaluation&)
314 {
315 throw std::logic_error("Not implemented: twoPhaseSatKrwInv()");
316 }
317
322 template <class FluidState, class Evaluation = typename FluidState::Scalar>
323 static Evaluation krn(const Params& params, const FluidState& fluidstate)
324 {
325 switch (params.approach()) {
327 return LETTwoPhaseLaw::krn(params.template getRealParams<SatCurveMultiplexerApproach::LETApproach>(),
328 fluidstate);
329 break;
330
332 return PLTwoPhaseLaw::krn(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinearApproach>(),
333 fluidstate);
334 break;
335 }
336
337 return 0.0;
338 }
339
340 template <class Evaluation>
341 static Evaluation twoPhaseSatKrn(const Params& params, const Evaluation& Sw)
342 {
343 switch (params.approach()) {
345 return LETTwoPhaseLaw::twoPhaseSatKrn(params.template getRealParams<SatCurveMultiplexerApproach::LETApproach>(),
346 Sw);
347 break;
348
350 return PLTwoPhaseLaw::twoPhaseSatKrn(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinearApproach>(),
351 Sw);
352 break;
353 }
354
355 return 0.0;
356 }
357
358 template <class Evaluation>
359 static Evaluation twoPhaseSatKrnInv(const Params& params, const Evaluation& krn)
360 {
361 switch (params.approach()) {
363 return LETTwoPhaseLaw::twoPhaseSatKrnInv(params.template getRealParams<SatCurveMultiplexerApproach::LETApproach>(),
364 krn);
365 break;
366
368 return PLTwoPhaseLaw::twoPhaseSatKrnInv(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinearApproach>(),
369 krn);
370 break;
371 }
372
373 return 0.0;
374 }
375};
376
377} // namespace Opm
378
379#endif // OPM_SAT_CURVE_MULTIPLEXER_HPP
Implementation of a tabulated, piecewise linear capillary pressure law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:50
static Evaluation krn(const Params &params, const FluidState &fs)
The relative permeability for the non-wetting phase of the porous medium.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:197
static Evaluation twoPhaseSatPcnw(const Params &params, const Evaluation &Sw)
The saturation-capillary pressure curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:141
static Evaluation krw(const Params &params, const FluidState &fs)
The relative permeability for the wetting phase of the porous medium.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:176
static Evaluation twoPhaseSatKrn(const Params &params, const Evaluation &Sw)
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:206
static void relativePermeabilities(Container &values, const Params &params, const FluidState &fs)
The relative permeabilities.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:117
static Evaluation Sn(const Params &params, const FluidState &fs)
Calculate the non-wetting phase saturations depending on the phase pressures.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:164
static Evaluation twoPhaseSatKrw(const Params &params, const Evaluation &Sw)
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:185
static void saturations(Container &, const Params &, const FluidState &)
The saturations of the fluid phases starting from their pressure differences.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:110
static Evaluation pcnw(const Params &params, const FluidState &fs)
The capillary pressure-saturation curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:129
static Evaluation twoPhaseSatSn(const Params &params, const Evaluation &pC)
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:168
static Evaluation twoPhaseSatKrnInv(const Params &params, const Evaluation &krn)
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:210
static void capillaryPressures(Container &values, const Params &params, const FluidState &fs)
The capillary pressure-saturation curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:97
static Evaluation Sw(const Params &, const FluidState &)
The saturation-capillary pressure curve.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:152
Implements a multiplexer class that provides LET curves and piecewise linear saturation functions.
Definition: SatCurveMultiplexer.hpp:44
static Evaluation Sw(const Params &params, const FluidState &fluidstate)
The saturation-capillary pressure curve.
Definition: SatCurveMultiplexer.hpp:208
static Evaluation twoPhaseSatKrnInv(const Params &params, const Evaluation &krn)
Definition: SatCurveMultiplexer.hpp:359
static Evaluation krn(const Params &params, const FluidState &fluidstate)
The relative permeability for the non-wetting phase of the medium.
Definition: SatCurveMultiplexer.hpp:323
static Evaluation krw(const Params &params, const FluidState &fluidstate)
The relative permeability for the wetting phase of the medium.
Definition: SatCurveMultiplexer.hpp:277
static constexpr bool isPressureDependent
Definition: SatCurveMultiplexer.hpp:73
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:113
static Evaluation twoPhaseSatSw(const Params &, const Evaluation &)
Definition: SatCurveMultiplexer.hpp:226
static constexpr bool implementsTwoPhaseSatApi
Definition: SatCurveMultiplexer.hpp:65
static Evaluation twoPhaseSatSn(const Params &params, const Evaluation &pc)
Definition: SatCurveMultiplexer.hpp:255
TraitsT Traits
Definition: SatCurveMultiplexer.hpp:46
static Evaluation twoPhaseSatPcnw(const Params &params, const Evaluation &Sw)
Definition: SatCurveMultiplexer.hpp:180
static constexpr bool isTemperatureDependent
Definition: SatCurveMultiplexer.hpp:77
static constexpr bool isCompositionDependent
Definition: SatCurveMultiplexer.hpp:81
static void capillaryPressures(Container &values, const Params &params, const FluidState &fluidState)
The capillary pressure-saturation curves.
Definition: SatCurveMultiplexer.hpp:91
static constexpr int numPhases
The number of fluid phases to which this material law applies.
Definition: SatCurveMultiplexer.hpp:54
static Evaluation twoPhaseSatKrn(const Params &params, const Evaluation &Sw)
Definition: SatCurveMultiplexer.hpp:341
static Evaluation Sn(const Params &params, const FluidState &fluidstate)
Calculate the non-wetting phase saturations depending on the phase pressures.
Definition: SatCurveMultiplexer.hpp:237
ParamsT Params
Definition: SatCurveMultiplexer.hpp:47
static void relativePermeabilities(Container &values, const Params &params, const FluidState &fluidState)
The relative permeability-saturation curves.
Definition: SatCurveMultiplexer.hpp:141
static Evaluation twoPhaseSatKrw(const Params &params, const Evaluation &Sw)
Definition: SatCurveMultiplexer.hpp:295
static Evaluation twoPhaseSatPcnwInv(const Params &, const Evaluation &)
Definition: SatCurveMultiplexer.hpp:198
static Evaluation twoPhaseSatKrwInv(const Params &, const Evaluation &)
Definition: SatCurveMultiplexer.hpp:313
typename Traits::Scalar Scalar
Definition: SatCurveMultiplexer.hpp:48
static constexpr bool implementsTwoPhaseApi
Definition: SatCurveMultiplexer.hpp:61
static constexpr bool isSaturationDependent
Definition: SatCurveMultiplexer.hpp:69
static Evaluation pcnw(const Params &params, const FluidState &fluidState)
The capillary pressure-saturation curve.
Definition: SatCurveMultiplexer.hpp:162
Implementation of the LET curve saturation functions.
Definition: TwoPhaseLETCurves.hpp:49
static Evaluation krn(const Params &, const FluidState &)
The relative permeability for the non-wetting phase of the medium as implied by the LET parameterizat...
Definition: TwoPhaseLETCurves.hpp:236
static Evaluation pcnw(const Params &, const FluidState &)
The capillary pressure-saturation curve.
Definition: TwoPhaseLETCurves.hpp:129
static Evaluation Sw(const Params &, const FluidState &)
Definition: TwoPhaseLETCurves.hpp:163
static Evaluation twoPhaseSatKrn(const Params &params, const Evaluation &Sw)
Definition: TwoPhaseLETCurves.hpp:243
static Evaluation twoPhaseSatPcnw(const Params &params, const Evaluation &Sw)
Definition: TwoPhaseLETCurves.hpp:137
static Evaluation krw(const Params &, const FluidState &)
The relative permeability for the wetting phase of the medium implied by the LET parameterization.
Definition: TwoPhaseLETCurves.hpp:192
static Evaluation Sn(const Params &, const FluidState &)
Definition: TwoPhaseLETCurves.hpp:175
static Evaluation twoPhaseSatKrw(const Params &params, const Evaluation &Sw)
Definition: TwoPhaseLETCurves.hpp:199
static Evaluation twoPhaseSatSn(const Params &, const Evaluation &)
Definition: TwoPhaseLETCurves.hpp:181
static void saturations(Container &, const Params &, const FluidState &)
Calculate the saturations of the phases starting from their pressure differences.
Definition: TwoPhaseLETCurves.hpp:102
static void capillaryPressures(Container &, const Params &, const FluidState &)
The capillary pressure-saturation curves.
Definition: TwoPhaseLETCurves.hpp:92
static Evaluation twoPhaseSatKrnInv(const Params &params, const Evaluation &krn)
Definition: TwoPhaseLETCurves.hpp:251
static void relativePermeabilities(Container &, const Params &, const FluidState &)
The relative permeability-saturation curves.
Definition: TwoPhaseLETCurves.hpp:118
Definition: Air_Mesitylene.hpp:34