My Project
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
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()) {
94 case SatCurveMultiplexerApproach::LET:
96 params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
97 fluidState);
98 break;
99
100 case SatCurveMultiplexerApproach::PiecewiseLinear:
102 params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
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()) {
116 case SatCurveMultiplexerApproach::LET:
118 params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
119 fluidState);
120 break;
121
122 case SatCurveMultiplexerApproach::PiecewiseLinear:
124 params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
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()) {
144 case SatCurveMultiplexerApproach::LET:
146 params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
147 fluidState);
148 break;
149
150 case SatCurveMultiplexerApproach::PiecewiseLinear:
152 params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
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()) {
165 case SatCurveMultiplexerApproach::LET:
166 return LETTwoPhaseLaw::pcnw(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
167 fluidState);
168 break;
169
170 case SatCurveMultiplexerApproach::PiecewiseLinear:
171 return PLTwoPhaseLaw::pcnw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
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()) {
183 case SatCurveMultiplexerApproach::LET:
184 return LETTwoPhaseLaw::twoPhaseSatPcnw(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
185 Sw);
186 break;
187
188 case SatCurveMultiplexerApproach::PiecewiseLinear:
189 return PLTwoPhaseLaw::twoPhaseSatPcnw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
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()) {
211 case SatCurveMultiplexerApproach::LET:
212 return LETTwoPhaseLaw::Sw(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
213 fluidstate);
214 break;
215
216 case SatCurveMultiplexerApproach::PiecewiseLinear:
217 return PLTwoPhaseLaw::Sw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
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()) {
240 case SatCurveMultiplexerApproach::LET:
241 return LETTwoPhaseLaw::Sn(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
242 fluidstate);
243 break;
244
245 case SatCurveMultiplexerApproach::PiecewiseLinear:
246 return PLTwoPhaseLaw::Sn(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
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()) {
258 case SatCurveMultiplexerApproach::LET:
259 return LETTwoPhaseLaw::twoPhaseSatSn(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
260 pc);
261 break;
262
263 case SatCurveMultiplexerApproach::PiecewiseLinear:
264 return PLTwoPhaseLaw::twoPhaseSatSn(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
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()) {
280 case SatCurveMultiplexerApproach::LET:
281 return LETTwoPhaseLaw::krw(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
282 fluidstate);
283 break;
284
285 case SatCurveMultiplexerApproach::PiecewiseLinear:
286 return PLTwoPhaseLaw::krw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
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()) {
298 case SatCurveMultiplexerApproach::LET:
299 return LETTwoPhaseLaw::twoPhaseSatKrw(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
300 Sw);
301 break;
302
303 case SatCurveMultiplexerApproach::PiecewiseLinear:
304 return PLTwoPhaseLaw::twoPhaseSatKrw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
305 Sw);
306 break;
307 }
308
309 return 0.0;
310 }
311
312 template <class Evaluation>
313 static Evaluation twoPhaseSatKrwInv(const Params& params, const Evaluation& krw)
314 {
315 switch (params.approach()) {
316 case SatCurveMultiplexerApproach::LET:
317 return LETTwoPhaseLaw::twoPhaseSatKrwInv(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
318 krw);
319 break;
320
321 case SatCurveMultiplexerApproach::PiecewiseLinear:
322 return PLTwoPhaseLaw::twoPhaseSatKrwInv(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
323 krw);
324 break;
325 }
326
327 return 0.0;
328 }
329
334 template <class FluidState, class Evaluation = typename FluidState::Scalar>
335 static Evaluation krn(const Params& params, const FluidState& fluidstate)
336 {
337 switch (params.approach()) {
338 case SatCurveMultiplexerApproach::LET:
339 return LETTwoPhaseLaw::krn(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
340 fluidstate);
341 break;
342
343 case SatCurveMultiplexerApproach::PiecewiseLinear:
344 return PLTwoPhaseLaw::krn(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
345 fluidstate);
346 break;
347 }
348
349 return 0.0;
350 }
351
352 template <class Evaluation>
353 static Evaluation twoPhaseSatKrn(const Params& params, const Evaluation& Sw)
354 {
355 switch (params.approach()) {
356 case SatCurveMultiplexerApproach::LET:
357 return LETTwoPhaseLaw::twoPhaseSatKrn(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
358 Sw);
359 break;
360
361 case SatCurveMultiplexerApproach::PiecewiseLinear:
362 return PLTwoPhaseLaw::twoPhaseSatKrn(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
363 Sw);
364 break;
365 }
366
367 return 0.0;
368 }
369
370 template <class Evaluation>
371 static Evaluation twoPhaseSatKrnInv(const Params& params, const Evaluation& krn)
372 {
373 switch (params.approach()) {
374 case SatCurveMultiplexerApproach::LET:
375 return LETTwoPhaseLaw::twoPhaseSatKrnInv(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
376 krn);
377 break;
378
379 case SatCurveMultiplexerApproach::PiecewiseLinear:
380 return PLTwoPhaseLaw::twoPhaseSatKrnInv(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
381 krn);
382 break;
383 }
384
385 return 0.0;
386 }
387};
388
389} // namespace Opm
390
391#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: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:210
static Evaluation twoPhaseSatPcnw(const Params &params, const Evaluation &Sw)
The saturation-capillary pressure curve.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:144
static Evaluation krw(const Params &params, const FluidState &fs)
The relative permeability for the wetting phase of the porous medium.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:182
static void relativePermeabilities(Container &values, const Params &params, const FluidState &fs)
The relative permeabilities.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:118
static Evaluation Sn(const Params &params, const FluidState &fs)
Calculate the non-wetting phase saturations depending on the phase pressures.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:170
static void saturations(Container &, const Params &, const FluidState &)
The saturations of the fluid phases starting from their pressure differences.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:111
static Evaluation pcnw(const Params &params, const FluidState &fs)
The capillary pressure-saturation curve.
Definition PiecewiseLinearTwoPhaseMaterial.hpp:131
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:158
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 krn(const Params &params, const FluidState &fluidstate)
The relative permeability for the non-wetting phase of the medium.
Definition SatCurveMultiplexer.hpp:335
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
Specify whether the quantities defined by this material law are dependent on the absolute pressure.
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 constexpr bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition SatCurveMultiplexer.hpp:65
static constexpr bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition SatCurveMultiplexer.hpp:77
static constexpr bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition.
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 Sn(const Params &params, const FluidState &fluidstate)
Calculate the non-wetting phase saturations depending on the phase pressures.
Definition SatCurveMultiplexer.hpp:237
static void relativePermeabilities(Container &values, const Params &params, const FluidState &fluidState)
The relative permeability-saturation curves.
Definition SatCurveMultiplexer.hpp:141
static constexpr bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition SatCurveMultiplexer.hpp:61
static constexpr bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
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:50
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:237
static Evaluation pcnw(const Params &, const FluidState &)
The capillary pressure-saturation curve.
Definition TwoPhaseLETCurves.hpp:130
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:193
static void saturations(Container &, const Params &, const FluidState &)
Calculate the saturations of the phases starting from their pressure differences.
Definition TwoPhaseLETCurves.hpp:103
static void capillaryPressures(Container &, const Params &, const FluidState &)
The capillary pressure-saturation curves.
Definition TwoPhaseLETCurves.hpp:93
static void relativePermeabilities(Container &, const Params &, const FluidState &)
The relative permeability-saturation curves.
Definition TwoPhaseLETCurves.hpp:119
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30