27#ifndef OPM_SAT_CURVE_MULTIPLEXER_HPP
28#define OPM_SAT_CURVE_MULTIPLEXER_HPP
36#if defined(__GNUC__) && (__GNUC__ < 13)
37 #define STATIC_ASSERT_SATCURVE_MULTIPLEXER_UNLESS_GCC_LT_13 throw std::logic_error("Unhandled SatCurveMultiplexerApproach")
39 #define STATIC_ASSERT_SATCURVE_MULTIPLEXER_UNLESS_GCC_LT_13 static_assert(false, "Unhandled SatCurveMultiplexerApproach")
50template <
class TraitsT,
class ParamsT = SatCurveMultiplexerParams<TraitsT> >
54 using Traits = TraitsT;
55 using Params = ParamsT;
56 using Scalar =
typename Traits::Scalar;
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");
91 static_assert(Traits::numPhases == 2,
92 "The number of fluid phases must be two if you want to use "
93 "this material law!");
98 template <
class Container,
class Flu
idState>
99 static void capillaryPressures(Container& values,
const Params& params,
const FluidState& fluidState)
101 switch (params.approach()) {
102 case SatCurveMultiplexerApproach::LET:
104 params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
108 case SatCurveMultiplexerApproach::PiecewiseLinear:
110 params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
120 template <
class Container,
class Flu
idState>
121 static void saturations(Container& values,
const Params& params,
const FluidState& fluidState)
123 switch (params.approach()) {
124 case SatCurveMultiplexerApproach::LET:
126 params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
130 case SatCurveMultiplexerApproach::PiecewiseLinear:
132 params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
148 template <
class Container,
class Flu
idState>
151 switch (params.approach()) {
152 case SatCurveMultiplexerApproach::LET:
154 params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
158 case SatCurveMultiplexerApproach::PiecewiseLinear:
160 params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
169 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
170 static Evaluation
pcnw(
const Params& params,
const FluidState& fluidState)
172 switch (params.approach()) {
173 case SatCurveMultiplexerApproach::LET:
178 case SatCurveMultiplexerApproach::PiecewiseLinear:
179 return PLTwoPhaseLaw::pcnw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
187 template <
class Evaluation,
class ...Args>
188 static Evaluation twoPhaseSatPcnw(
const Params& params,
const Evaluation&
Sw)
190 if constexpr (FrontIsSatCurveMultiplexerDispatchV<Args...>) {
191 return twoPhaseSatPcnwT<Evaluation, Args...>(params,
Sw);
194 switch (params.approach()) {
195 case SatCurveMultiplexerApproach::LET:
196 return LETTwoPhaseLaw::twoPhaseSatPcnw(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
200 case SatCurveMultiplexerApproach::PiecewiseLinear:
209 template <
class Evaluation,
class Head,
class ...Args>
210 static Evaluation twoPhaseSatPcnwT(
const Params& params,
const Evaluation&
Sw)
212 if constexpr (Head::approach == SatCurveMultiplexerApproach::LET) {
213 return LETTwoPhaseLaw::twoPhaseSatPcnw(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
215 }
else if constexpr (Head::approach == SatCurveMultiplexerApproach::PiecewiseLinear) {
219 STATIC_ASSERT_SATCURVE_MULTIPLEXER_UNLESS_GCC_LT_13;
223 template <
class Evaluation>
224 static Evaluation twoPhaseSatPcnwInv(
const Params&,
const Evaluation&)
226 throw std::logic_error(
"SatCurveMultiplexer::twoPhaseSatPcnwInv"
227 " not implemented!");
233 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
234 static Evaluation
Sw(
const Params& params,
const FluidState& fluidstate)
236 switch (params.approach()) {
237 case SatCurveMultiplexerApproach::LET:
238 return LETTwoPhaseLaw::Sw(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
242 case SatCurveMultiplexerApproach::PiecewiseLinear:
243 return PLTwoPhaseLaw::Sw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
251 template <
class Evaluation>
252 static Evaluation twoPhaseSatSw(
const Params&,
const Evaluation&)
254 throw std::logic_error(
"SatCurveMultiplexer::twoPhaseSatSw"
255 " not implemented!");
262 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
263 static Evaluation
Sn(
const Params& params,
const FluidState& fluidstate)
265 switch (params.approach()) {
266 case SatCurveMultiplexerApproach::LET:
267 return LETTwoPhaseLaw::Sn(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
271 case SatCurveMultiplexerApproach::PiecewiseLinear:
272 return PLTwoPhaseLaw::Sn(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
280 template <
class Evaluation>
281 static Evaluation twoPhaseSatSn(
const Params& params,
const Evaluation& pc)
283 switch (params.approach()) {
284 case SatCurveMultiplexerApproach::LET:
285 return LETTwoPhaseLaw::twoPhaseSatSn(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
289 case SatCurveMultiplexerApproach::PiecewiseLinear:
290 return PLTwoPhaseLaw::twoPhaseSatSn(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
302 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
303 static Evaluation
krw(
const Params& params,
const FluidState& fluidstate)
305 switch (params.approach()) {
306 case SatCurveMultiplexerApproach::LET:
307 return LETTwoPhaseLaw::krw(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
311 case SatCurveMultiplexerApproach::PiecewiseLinear:
312 return PLTwoPhaseLaw::krw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
320 template <
class Evaluation,
class ...Args>
321 static Evaluation twoPhaseSatKrw(
const Params& params,
const Evaluation&
Sw)
323 if constexpr (FrontIsSatCurveMultiplexerDispatchV<Args...>) {
324 return twoPhaseSatKrwT<Evaluation, Args...>(params,
Sw);
327 switch (params.approach()) {
328 case SatCurveMultiplexerApproach::LET:
329 return LETTwoPhaseLaw::twoPhaseSatKrw(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
333 case SatCurveMultiplexerApproach::PiecewiseLinear:
334 return PLTwoPhaseLaw::twoPhaseSatKrw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
342 template <
class Evaluation,
class Head,
class ...Args>
343 static Evaluation twoPhaseSatKrwT(
const Params& params,
const Evaluation&
Sw)
345 if constexpr (Head::approach == SatCurveMultiplexerApproach::LET) {
346 return LETTwoPhaseLaw::twoPhaseSatKrw(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
348 }
else if constexpr (Head::approach == SatCurveMultiplexerApproach::PiecewiseLinear) {
349 return PLTwoPhaseLaw::twoPhaseSatKrw(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
352 STATIC_ASSERT_SATCURVE_MULTIPLEXER_UNLESS_GCC_LT_13;
356 template <
class Evaluation>
357 static Evaluation twoPhaseSatKrwInv(
const Params& params,
const Evaluation&
krw)
359 switch (params.approach()) {
360 case SatCurveMultiplexerApproach::LET:
361 return LETTwoPhaseLaw::twoPhaseSatKrwInv(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
365 case SatCurveMultiplexerApproach::PiecewiseLinear:
366 return PLTwoPhaseLaw::twoPhaseSatKrwInv(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
378 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
379 static Evaluation
krn(
const Params& params,
const FluidState& fluidstate)
381 switch (params.approach()) {
382 case SatCurveMultiplexerApproach::LET:
383 return LETTwoPhaseLaw::krn(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
387 case SatCurveMultiplexerApproach::PiecewiseLinear:
388 return PLTwoPhaseLaw::krn(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
396 template <
class Evaluation,
class ...Args>
397 static Evaluation twoPhaseSatKrn(
const Params& params,
const Evaluation&
Sw)
399 if constexpr (FrontIsSatCurveMultiplexerDispatchV<Args...>) {
400 return twoPhaseSatKrnT<Evaluation, Args...>(params,
Sw);
403 switch (params.approach()) {
404 case SatCurveMultiplexerApproach::LET:
405 return LETTwoPhaseLaw::twoPhaseSatKrn(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
409 case SatCurveMultiplexerApproach::PiecewiseLinear:
410 return PLTwoPhaseLaw::twoPhaseSatKrn(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
418 template <
class Evaluation,
class Head,
class ...Args>
419 static Evaluation twoPhaseSatKrnT(
const Params& params,
const Evaluation&
Sw)
421 if constexpr (Head::approach == SatCurveMultiplexerApproach::LET) {
422 return LETTwoPhaseLaw::twoPhaseSatKrn(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
424 }
else if constexpr (Head::approach == SatCurveMultiplexerApproach::PiecewiseLinear) {
425 return PLTwoPhaseLaw::twoPhaseSatKrn(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
428 STATIC_ASSERT_SATCURVE_MULTIPLEXER_UNLESS_GCC_LT_13;
432 template <
class Evaluation>
433 static Evaluation twoPhaseSatKrnInv(
const Params& params,
const Evaluation&
krn)
435 switch (params.approach()) {
436 case SatCurveMultiplexerApproach::LET:
437 return LETTwoPhaseLaw::twoPhaseSatKrnInv(params.template getRealParams<SatCurveMultiplexerApproach::LET>(),
441 case SatCurveMultiplexerApproach::PiecewiseLinear:
442 return PLTwoPhaseLaw::twoPhaseSatKrnInv(params.template getRealParams<SatCurveMultiplexerApproach::PiecewiseLinear>(),
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 ¶ms, const FluidState &fs)
Definition PiecewiseLinearTwoPhaseMaterial.hpp:100
static OPM_HOST_DEVICE Evaluation Sn(const Params ¶ms, const FluidState &fs)
Definition PiecewiseLinearTwoPhaseMaterial.hpp:173
static OPM_HOST_DEVICE void relativePermeabilities(Container &values, const Params ¶ms, const FluidState &fs)
Definition PiecewiseLinearTwoPhaseMaterial.hpp:121
static OPM_HOST_DEVICE Evaluation twoPhaseSatPcnw(const Params ¶ms, const Evaluation &Sw)
Definition PiecewiseLinearTwoPhaseMaterial.hpp:147
static OPM_HOST_DEVICE Evaluation krw(const Params ¶ms, 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 ¶ms, 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 ¶ms, 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 ¶ms, const FluidState &fluidstate)
The saturation-capillary pressure curve.
Definition SatCurveMultiplexer.hpp:234
static Evaluation krn(const Params ¶ms, const FluidState &fluidstate)
The relative permeability for the non-wetting phase of the medium.
Definition SatCurveMultiplexer.hpp:379
static Evaluation krw(const Params ¶ms, 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 ¶ms, 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 ¶ms, const FluidState &fluidState)
The capillary pressure-saturation curves.
Definition SatCurveMultiplexer.hpp:99
static Evaluation Sn(const Params ¶ms, 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 ¶ms, 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 ¶ms, 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