27#ifndef OPM_PARKER_LENHARD_HPP
28#define OPM_PARKER_LENHARD_HPP
46template <
class ScalarT>
50 using Scalar = ScalarT;
166 return this->
Sw() < SwReversal && SwReversal < prev_->Sw();
173 return prev_->Sw() < SwReversal && SwReversal < this->
Sw();
245template <
class TraitsT,
class ParamsT = ParkerLenhardParams<TraitsT> >
249 using Traits = TraitsT;
250 using Params = ParamsT;
251 using Scalar =
typename Traits::Scalar;
256 "The Parker-Lenhard capillary pressure law only "
257 "applies to the case of two fluid phases");
283 static_assert(Traits::numPhases == 2,
284 "The number of fluid phases must be two if you want to use "
285 "this material law!");
288 typedef typename ParamsT::VanGenuchten VanGenuchten;
299 params.setMdc(
new ScanningCurve(params.SwrPc()));
300 params.setCsc(params.mdc());
301 params.setPisc(
nullptr);
302 params.setCurrentSnr(0.0);
309 template <
class Flu
idState>
310 static void update(Params& params,
const FluidState& fs)
312 const Scalar sw = scalarValue(fs.saturation(Traits::wettingPhaseIdx));
323 ScanningCurve* curve = findScanningCurve_(params, sw);
329 Scalar Sw_mic = VanGenuchten::twoPhaseSatSw(params.micParams(), pc);
330 Scalar Sw_mdc = VanGenuchten::twoPhaseSatSw(params.mdcParams(), pc);
332 curve->
setNext(sw, pc, Sw_mic, Sw_mdc);
336 params.setCsc(curve);
339 if (params.csc() == params.mdc()) {
340 params.setPisc(params.mdc()->next());
341 params.setCurrentSnr(computeCurrentSnr_(params, sw));
349 template <
class Container,
class Flu
idState>
352 using Evaluation = std::remove_reference_t<
decltype(values[0])>;
354 values[Traits::wettingPhaseIdx] = 0.0;
362 template <
class Container,
class Flu
idState>
363 static void saturations(Container& ,
const Params& ,
const FluidState& )
364 {
throw std::logic_error(
"Not implemented: ParkerLenhard::saturations()"); }
370 template <
class Container,
class Flu
idState>
373 using Evaluation = std::remove_reference_t<
decltype(values[0])>;
383 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
384 static Evaluation
pcnw(
const Params& params,
const FluidState& fs)
386 const Evaluation& sw =
387 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
389 return twoPhaseSatPcnw(params, sw);
392 template <
class Evaluation>
393 static Evaluation twoPhaseSatPcnw(
const Params& params,
const Evaluation&
Sw)
396 ScanningCurve* sc = findScanningCurve_(params, scalarValue(
Sw));
411 const Evaluation& pos = (Sw_app - SwAppCurSC)/(SwAppPrevSC - SwAppCurSC);
413 const Evaluation& SwMic =
414 pos * (sc->prev()->SwMic() - sc->SwMic()) + sc->SwMic();
419 const Evaluation SwMdc =
420 pos * (sc->prev()->SwMdc() - sc->SwMdc()) + sc->SwMdc();
430 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
431 static Evaluation
Sw(
const Params& ,
const FluidState& )
432 {
throw std::logic_error(
"Not implemented: ParkerLenhard::Sw()"); }
434 template <
class Evaluation>
435 static Evaluation twoPhaseSatSw(
const Params& ,
const Evaluation& )
436 {
throw std::logic_error(
"Not implemented: ParkerLenhard::twoPhaseSatSw()"); }
442 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
443 static Evaluation
Sn(
const Params& params,
const FluidState& fs)
446 template <
class Evaluation>
447 static Evaluation twoPhaseSatSn(
const Params& ,
const Evaluation& )
448 {
throw std::logic_error(
"Not implemented: ParkerLenhard::twoPhaseSatSn()"); }
454 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
455 static Evaluation
krw(
const Params& params,
const FluidState& fs)
457 const Evaluation& sw =
458 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
460 return twoPhaseSatKrw(params, sw);
463 template <
class Evaluation>
464 static Evaluation twoPhaseSatKrw(
const Params& params,
const Evaluation&
Sw)
469 return VanGenuchten::twoPhaseSatKrw(params.mdcParams(), Sw_app);
476 template <
class Flu
idState,
class Evaluation =
typename Flu
idState::Scalar>
477 static Evaluation
krn(
const Params& params,
const FluidState& fs)
479 const Evaluation& sw =
480 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
482 return twoPhaseSatKrn(params, sw);
485 template <
class Evaluation>
486 static Evaluation twoPhaseSatKrn(
const Params& params,
const Evaluation&
Sw)
491 return VanGenuchten::twoPhaseSatKrn(params.mdcParams(), Sw_app);
497 template <
class Evaluation>
500 return effectiveToApparentSw_(params, absoluteToEffectiveSw_(params,
Sw));
513 template <
class Evaluation>
514 static Evaluation absoluteToEffectiveSw_(
const Params& params,
const Evaluation&
Sw)
515 {
return (
Sw - params.SwrPc()) / (1 - params.SwrPc()); }
526 template <
class Evaluation>
527 static Evaluation effectiveToAbsoluteSw_(
const Params& params,
const Evaluation& Swe)
528 {
return Swe * (1 - params.SwrPc()) + params.SwrPc(); }
532 template <
class Evaluation>
533 static Evaluation computeCurrentSnr_(
const Params& params,
const Evaluation&
Sw)
536 if (
Sw > 1 - params.Snr()) {
539 if (
Sw < params.SwrPc()) {
543 if (params.Snr() == 0.0) {
548 const Scalar R = 1.0 / params.Snr() - 1;
549 const Evaluation& curSnr = (1 -
Sw) / (1 + R * (1 -
Sw));
553 assert(curSnr <= params.Snr());
560 template <
class Evaluation>
561 static Evaluation trappedEffectiveSn_(
const Params& params,
const Evaluation&
Sw)
563 const Evaluation& Swe = absoluteToEffectiveSw_(params,
Sw);
564 const Scalar SwePisc = absoluteToEffectiveSw_(params, params.pisc()->Sw());
566 const Scalar Snre = absoluteToEffectiveSw_(params, params.currentSnr());
567 return Snre * (Swe - SwePisc) / (1 - Snre - SwePisc);
572 template <
class Evaluation>
573 static Evaluation effectiveToApparentSw_(
const Params& params,
const Evaluation& Swe)
575 if (params.pisc() ==
nullptr ||
576 Swe <= absoluteToEffectiveSw_(params, params.pisc()->Sw()))
587 return Swe + trappedEffectiveSn_(params, Swe);
591 template <
class Evaluation>
592 static Evaluation apparentToEffectiveSw_(
const Params& params,
const Evaluation& Swapp)
594 const Scalar SwePisc = absoluteToEffectiveSw_(params, params.pisc()->Sw());
595 if (params.pisc() ==
nullptr || Swapp <= SwePisc) {
602 const Scalar Snre = absoluteToEffectiveSw_(params.currentSnr());
604 (Swapp * (1 - Snre - SwePisc) + Snre * SwePisc)
611 static ScanningCurve* findScanningCurve_(
const Params& params, Scalar
Sw)
613 if (params.pisc() ==
nullptr ||
Sw <= params.pisc()->Sw()) {
624 if (params.pisc()->next() ==
nullptr ||
625 params.pisc()->next()->Sw() <
Sw)
627 return params.pisc();
630 ScanningCurve* curve = params.csc()->next();
632 assert(curve != params.mdc()->prev());
633 if (curve->isValidAt_Sw(
Sw)) {
636 curve = curve->prev();
Default parameter class for the Parker-Lenhard hysteresis model.
Implementation of the van Genuchten capillary pressure - saturation relation.
Represents a scanning curve in the Parker-Lenhard hysteresis model.
Definition ParkerLenhard.hpp:48
bool isDrain() const
Returns true iff the scanning curve is a drainage curve.
Definition ParkerLenhard.hpp:188
bool isValidAt_Sw(Scalar SwReversal) const
Returns true iff the given effective saturation Swei is within the scope of the curve,...
Definition ParkerLenhard.hpp:159
PLScanningCurve * next() const
Return the next scanning curve, i.e.
Definition ParkerLenhard.hpp:124
~PLScanningCurve()
Destructor.
Definition ParkerLenhard.hpp:103
bool isImbib() const
Returns true iff the scanning curve is a imbibition curve.
Definition ParkerLenhard.hpp:181
void setNext(Scalar SwReversal, Scalar pcnwReversal, Scalar SwMiCurve, Scalar SwMdCurve)
Set the next scanning curve.
Definition ParkerLenhard.hpp:135
Scalar Sw() const
Absolute wetting-phase saturation at the scanning curve's reversal point.
Definition ParkerLenhard.hpp:203
Scalar pcnw() const
Capillary pressure at the last reversal point.
Definition ParkerLenhard.hpp:209
Scalar SwMic() const
Apparent saturation of the last reversal point on the pressure MIC.
Definition ParkerLenhard.hpp:216
Scalar SwMdc() const
Apparent saturation of the last reversal point on the pressure MDC.
Definition ParkerLenhard.hpp:223
PLScanningCurve * prev() const
Return the previous scanning curve, i.e.
Definition ParkerLenhard.hpp:117
PLScanningCurve(const Scalar Swr)
Constructs main imbibition curve.
Definition ParkerLenhard.hpp:58
int loopNum() const
The loop number of the scanning curve.
Definition ParkerLenhard.hpp:196
Implements the Parker-Lenhard twophase p_c-Sw hysteresis model.
Definition ParkerLenhard.hpp:247
static void reset(Params ¶ms)
Resets the hysteresis model to the initial parameters on the main drainage curve.
Definition ParkerLenhard.hpp:296
static Evaluation Sw(const Params &, const FluidState &)
Calculate the wetting phase saturations depending on the phase pressures.
Definition ParkerLenhard.hpp:431
static Evaluation krn(const Params ¶ms, const FluidState &fs)
The relative permeability for the non-wetting phase of the params.
Definition ParkerLenhard.hpp:477
static Evaluation krw(const Params ¶ms, const FluidState &fs)
The relative permeability for the wetting phase of the medium.
Definition ParkerLenhard.hpp:455
static void update(Params ¶ms, const FluidState &fs)
Set the current absolute saturation for the current timestep.
Definition ParkerLenhard.hpp:310
static void relativePermeabilities(Container &values, const Params ¶ms, const FluidState &fs)
Returns the relative permeabilities of the phases dependening on the phase saturations.
Definition ParkerLenhard.hpp:371
static constexpr bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition ParkerLenhard.hpp:265
static Evaluation absoluteToApparentSw_(const Params ¶ms, const Evaluation &Sw)
Convert an absolute wetting saturation to an apparent one.
Definition ParkerLenhard.hpp:498
static constexpr bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure.
Definition ParkerLenhard.hpp:273
static constexpr bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition ParkerLenhard.hpp:269
static constexpr bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition ParkerLenhard.hpp:277
static Evaluation Sn(const Params ¶ms, const FluidState &fs)
Calculate the non-wetting phase saturations depending on the phase pressures.
Definition ParkerLenhard.hpp:443
static Evaluation pcnw(const Params ¶ms, const FluidState &fs)
Returns the capillary pressure dependend on the phase saturations.
Definition ParkerLenhard.hpp:384
static constexpr int numPhases
The number of fluid phases.
Definition ParkerLenhard.hpp:254
static void capillaryPressures(Container &values, const Params ¶ms, const FluidState &fs)
Returns the capillary pressure dependening on the phase saturations.
Definition ParkerLenhard.hpp:350
static void saturations(Container &, const Params &, const FluidState &)
Returns the capillary pressure dependening on the phase saturations.
Definition ParkerLenhard.hpp:363
static constexpr bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition ParkerLenhard.hpp:261
static constexpr bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition.
Definition ParkerLenhard.hpp:281
static Evaluation twoPhaseSatPcnw(const Params ¶ms, const Evaluation &Sw)
The saturation-capillary pressure curve according to van Genuchten using a material law specific API.
Definition VanGenuchten.hpp:194
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30