opm-common
Loading...
Searching...
No Matches
ParkerLenhard.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 2 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_PARKER_LENHARD_HPP
28#define OPM_PARKER_LENHARD_HPP
29
32
33#include <cassert>
34#include <stdexcept>
35#include <type_traits>
36
37namespace Opm {
38
46template <class ScalarT>
48{
49public:
50 using Scalar = ScalarT;
51
58 explicit PLScanningCurve(const Scalar Swr)
59 {
60 loopNum_ = 0;
61 prev_ = new PLScanningCurve(nullptr, // prev
62 this, // next
63 -1, // loop number
64 Swr, // Sw
65 1e12, // pcnw
66 Swr, // SwMic
67 Swr); // SwMdc
68 next_ = nullptr;
69
70 Sw_ = 1.0;
71 pcnw_ = 0.0;
72 SwMic_ = 1.0;
73 SwMdc_ = 1.0;
74 }
75
76 PLScanningCurve& operator=(const PLScanningCurve&) = delete;
77 PLScanningCurve(const PLScanningCurve&) = delete;
78
79protected:
81 PLScanningCurve* nextSC,
82 int loopN,
83 Scalar SwReversal,
84 Scalar pcnwReversal,
85 Scalar SwMiCurve,
86 Scalar SwMdCurve)
87 : prev_(prevSC)
88 , next_(nextSC)
89 , loopNum_(loopN)
90 , Sw_(SwReversal)
91 , pcnw_(pcnwReversal)
92 , SwMdc_(SwMdCurve)
93 , SwMic_(SwMiCurve)
94 {
95 }
96
97public:
104 {
105 if (loopNum_ == 0) {
106 delete prev_;
107 }
108 if (loopNum_ >= 0) {
109 delete next_;
110 }
111 }
112
118 { return prev_; }
119
125 { return next_; }
126
135 void setNext(Scalar SwReversal,
136 Scalar pcnwReversal,
137 Scalar SwMiCurve,
138 Scalar SwMdCurve)
139 {
140 // if next_ is nullptr, delete does nothing, so
141 // this is valid!!
142 delete next_;
143
144 next_ = new PLScanningCurve(this, // prev
145 nullptr, // next
146 loopNum() + 1,
147 SwReversal,
148 pcnwReversal,
149 SwMiCurve,
150 SwMdCurve);
151 }
152
159 bool isValidAt_Sw(Scalar SwReversal) const
160 {
161 if (isImbib()) {
162 // for inbibition the given saturation
163 // must be between the start of the
164 // current imbibition and the the start
165 // of the last drainage
166 return this->Sw() < SwReversal && SwReversal < prev_->Sw();
167 }
168 else {
169 // for drainage the given saturation
170 // must be between the start of the
171 // last imbibition and the start
172 // of the current drainage
173 return prev_->Sw() < SwReversal && SwReversal < this->Sw();
174 }
175 }
176
181 bool isImbib() const
182 { return loopNum() % 2 == 1; }
183
188 bool isDrain() const
189 { return !isImbib(); }
190
196 int loopNum() const
197 { return loopNum_; }
198
203 Scalar Sw() const
204 { return Sw_; }
205
209 Scalar pcnw() const
210 { return pcnw_; }
211
216 Scalar SwMic() const
217 { return SwMic_; }
218
223 Scalar SwMdc() const
224 { return SwMdc_; }
225
226private:
227 PLScanningCurve* prev_;
228 PLScanningCurve* next_;
229
230 int loopNum_;
231
232 Scalar Sw_;
233 Scalar pcnw_;
234
235 Scalar SwMdc_;
236 Scalar SwMic_;
237};
238
245template <class TraitsT, class ParamsT = ParkerLenhardParams<TraitsT> >
246class ParkerLenhard : public TraitsT
247{
248public:
249 using Traits = TraitsT;
250 using Params = ParamsT;
251 using Scalar = typename Traits::Scalar;
252
254 static constexpr int numPhases = Traits::numPhases;
255 static_assert(numPhases == 2,
256 "The Parker-Lenhard capillary pressure law only "
257 "applies to the case of two fluid phases");
258
261 static constexpr bool implementsTwoPhaseApi = true;
262
265 static constexpr bool implementsTwoPhaseSatApi = true;
266
269 static constexpr bool isSaturationDependent = true;
270
273 static constexpr bool isPressureDependent = false;
274
277 static constexpr bool isTemperatureDependent = false;
278
281 static constexpr bool isCompositionDependent = false;
282
283 static_assert(Traits::numPhases == 2,
284 "The number of fluid phases must be two if you want to use "
285 "this material law!");
286
287private:
288 typedef typename ParamsT::VanGenuchten VanGenuchten;
289 typedef PLScanningCurve<Scalar> ScanningCurve;
290
291public:
296 static void reset(Params& params)
297 {
298 delete params.mdc(); // this will work even if mdc_ == nullptr!
299 params.setMdc(new ScanningCurve(params.SwrPc()));
300 params.setCsc(params.mdc());
301 params.setPisc(nullptr);
302 params.setCurrentSnr(0.0);
303 }
304
309 template <class FluidState>
310 static void update(Params& params, const FluidState& fs)
311 {
312 const Scalar sw = scalarValue(fs.saturation(Traits::wettingPhaseIdx));
313
314 if (sw > 1 - 1e-5) {
315 // if the absolute saturation is almost 1,
316 // it means that we're back to the beginning
317 reset(params);
318 return;
319 }
320
321 // find the loop number which corrosponds to the
322 // given effective saturation
323 ScanningCurve* curve = findScanningCurve_(params, sw);
324
325 // calculate the apparent saturation on the MIC and MDC
326 // which yield the same capillary pressure as the
327 // Sw at the current scanning curve
328 Scalar pc = pcnw<FluidState, Scalar>(params, fs);
329 Scalar Sw_mic = VanGenuchten::twoPhaseSatSw(params.micParams(), pc);
330 Scalar Sw_mdc = VanGenuchten::twoPhaseSatSw(params.mdcParams(), pc);
331
332 curve->setNext(sw, pc, Sw_mic, Sw_mdc);
333 if (!curve->next())
334 return;
335
336 params.setCsc(curve);
337
338 // if we're back on the MDC, we also have a new PISC!
339 if (params.csc() == params.mdc()) {
340 params.setPisc(params.mdc()->next());
341 params.setCurrentSnr(computeCurrentSnr_(params, sw));
342 }
343 }
344
349 template <class Container, class FluidState>
350 static void capillaryPressures(Container& values, const Params& params, const FluidState& fs)
351 {
352 using Evaluation = std::remove_reference_t<decltype(values[0])>;
353
354 values[Traits::wettingPhaseIdx] = 0.0; // reference phase
355 values[Traits::nonWettingPhaseIdx] = pcnw<FluidState, Evaluation>(params, fs);
356 }
357
362 template <class Container, class FluidState>
363 static void saturations(Container& /*values*/, const Params& /*params*/, const FluidState& /*fs*/)
364 { throw std::logic_error("Not implemented: ParkerLenhard::saturations()"); }
365
370 template <class Container, class FluidState>
371 static void relativePermeabilities(Container& values, const Params& params, const FluidState& fs)
372 {
373 using Evaluation = std::remove_reference_t<decltype(values[0])>;
374
375 values[Traits::wettingPhaseIdx] = krw<FluidState, Evaluation>(params, fs);
376 values[Traits::nonWettingPhaseIdx] = krn<FluidState, Evaluation>(params, fs);
377 }
378
383 template <class FluidState, class Evaluation = typename FluidState::ValueType>
384 static Evaluation pcnw(const Params& params, const FluidState& fs)
385 {
386 const Evaluation& sw =
387 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
388
389 return twoPhaseSatPcnw(params, sw);
390 }
391
392 template <class Evaluation>
393 static Evaluation twoPhaseSatPcnw(const Params& params, const Evaluation& Sw)
394 {
395 // calculate the current apparent saturation
396 ScanningCurve* sc = findScanningCurve_(params, scalarValue(Sw));
397
398 // calculate the apparant saturation
399 const Evaluation& Sw_app = absoluteToApparentSw_(params, Sw);
400
401 // if the apparent saturation exceeds the 'legal' limits,
402 // we also the underlying material law decide what to do.
403 if (Sw_app > 1) {
404 return 0.0; // VanGenuchten::pcnw(params.mdcParams(), Sw_app);
405 }
406
407 // put the apparent saturation into the main imbibition or
408 // drainage curve
409 Scalar SwAppCurSC = absoluteToApparentSw_(params, sc->Sw());
410 Scalar SwAppPrevSC = absoluteToApparentSw_(params, sc->prev()->Sw());
411 const Evaluation& pos = (Sw_app - SwAppCurSC)/(SwAppPrevSC - SwAppCurSC);
412 if (sc->isImbib()) {
413 const Evaluation& SwMic =
414 pos * (sc->prev()->SwMic() - sc->SwMic()) + sc->SwMic();
415
416 return VanGenuchten::twoPhaseSatPcnw(params.micParams(), SwMic);
417 }
418 else { // sc->isDrain()
419 const Evaluation SwMdc =
420 pos * (sc->prev()->SwMdc() - sc->SwMdc()) + sc->SwMdc();
421
422 return VanGenuchten::twoPhaseSatPcnw(params.mdcParams(), SwMdc);
423 }
424 }
425
430 template <class FluidState, class Evaluation = typename FluidState::ValueType>
431 static Evaluation Sw(const Params& /*params*/, const FluidState& /*fs*/)
432 { throw std::logic_error("Not implemented: ParkerLenhard::Sw()"); }
433
434 template <class Evaluation>
435 static Evaluation twoPhaseSatSw(const Params& /*params*/, const Evaluation& /*pc*/)
436 { throw std::logic_error("Not implemented: ParkerLenhard::twoPhaseSatSw()"); }
437
442 template <class FluidState, class Evaluation = typename FluidState::ValueType>
443 static Evaluation Sn(const Params& params, const FluidState& fs)
444 { return 1 - Sw<FluidState, Evaluation>(params, fs); }
445
446 template <class Evaluation>
447 static Evaluation twoPhaseSatSn(const Params& /*params*/, const Evaluation& /*pc*/)
448 { throw std::logic_error("Not implemented: ParkerLenhard::twoPhaseSatSn()"); }
449
454 template <class FluidState, class Evaluation = typename FluidState::ValueType>
455 static Evaluation krw(const Params& params, const FluidState& fs)
456 {
457 const Evaluation& sw =
458 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
459
460 return twoPhaseSatKrw(params, sw);
461 }
462
463 template <class Evaluation>
464 static Evaluation twoPhaseSatKrw(const Params& params, const Evaluation& Sw)
465 {
466 // for the effective permeability we only use Land's law and
467 // the relative permeability of the main drainage curve.
468 const Evaluation& Sw_app = absoluteToApparentSw_(params, Sw);
469 return VanGenuchten::twoPhaseSatKrw(params.mdcParams(), Sw_app);
470 }
471
476 template <class FluidState, class Evaluation = typename FluidState::ValueType>
477 static Evaluation krn(const Params& params, const FluidState& fs)
478 {
479 const Evaluation& sw =
480 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
481
482 return twoPhaseSatKrn(params, sw);
483 }
484
485 template <class Evaluation>
486 static Evaluation twoPhaseSatKrn(const Params& params, const Evaluation& Sw)
487 {
488 // for the effective permeability we only use Land's law and
489 // the relative permeability of the main drainage curve.
490 const Evaluation& Sw_app = absoluteToApparentSw_(params, Sw);
491 return VanGenuchten::twoPhaseSatKrn(params.mdcParams(), Sw_app);
492 }
493
497 template <class Evaluation>
498 static Evaluation absoluteToApparentSw_(const Params& params, const Evaluation& Sw)
499 {
500 return effectiveToApparentSw_(params, absoluteToEffectiveSw_(params, Sw));
501 }
502
503private:
513 template <class Evaluation>
514 static Evaluation absoluteToEffectiveSw_(const Params& params, const Evaluation& Sw)
515 { return (Sw - params.SwrPc()) / (1 - params.SwrPc()); }
516
526 template <class Evaluation>
527 static Evaluation effectiveToAbsoluteSw_(const Params& params, const Evaluation& Swe)
528 { return Swe * (1 - params.SwrPc()) + params.SwrPc(); }
529
530 // return the effctive residual non-wetting saturation, given an
531 // effective wetting saturation
532 template <class Evaluation>
533 static Evaluation computeCurrentSnr_(const Params& params, const Evaluation& Sw)
534 {
535 // regularize
536 if (Sw > 1 - params.Snr()) {
537 return 0.0;
538 }
539 if (Sw < params.SwrPc()) {
540 return params.Snr();
541 }
542
543 if (params.Snr() == 0.0) {
544 return 0.0;
545 }
546
547 // use Land's law
548 const Scalar R = 1.0 / params.Snr() - 1;
549 const Evaluation& curSnr = (1 - Sw) / (1 + R * (1 - Sw));
550
551 // the current effective residual non-wetting saturation must
552 // be smaller than the residual non-wetting saturation
553 assert(curSnr <= params.Snr());
554
555 return curSnr;
556 }
557
558 // returns the trapped effective non-wetting saturation for a
559 // given wetting phase saturation
560 template <class Evaluation>
561 static Evaluation trappedEffectiveSn_(const Params& params, const Evaluation& Sw)
562 {
563 const Evaluation& Swe = absoluteToEffectiveSw_(params, Sw);
564 const Scalar SwePisc = absoluteToEffectiveSw_(params, params.pisc()->Sw());
565
566 const Scalar Snre = absoluteToEffectiveSw_(params, params.currentSnr());
567 return Snre * (Swe - SwePisc) / (1 - Snre - SwePisc);
568 }
569
570 // returns the apparent saturation of the wetting phase depending
571 // on the effective saturation
572 template <class Evaluation>
573 static Evaluation effectiveToApparentSw_(const Params& params, const Evaluation& Swe)
574 {
575 if (params.pisc() == nullptr ||
576 Swe <= absoluteToEffectiveSw_(params, params.pisc()->Sw()))
577 {
578 // we are on the main drainage curve, i.e. no non-wetting
579 // fluid is trapped -> apparent saturation == effective
580 // saturation
581 return Swe;
582 }
583
584 // we are on a imbibition or drainage curve which is not the
585 // main drainage curve -> apparent saturation == effective
586 // saturation + trapped effective saturation
587 return Swe + trappedEffectiveSn_(params, Swe);
588 }
589
590 // Returns the effective saturation to a given apparent one
591 template <class Evaluation>
592 static Evaluation apparentToEffectiveSw_(const Params& params, const Evaluation& Swapp)
593 {
594 const Scalar SwePisc = absoluteToEffectiveSw_(params, params.pisc()->Sw());
595 if (params.pisc() == nullptr || Swapp <= SwePisc) {
596 // we are on the main drainage curve, i.e.
597 // no non-wetting fluid is trapped
598 // -> apparent saturation == effective saturation
599 return Swapp;
600 }
601
602 const Scalar Snre = absoluteToEffectiveSw_(params.currentSnr());
603 return
604 (Swapp * (1 - Snre - SwePisc) + Snre * SwePisc)
605 / (1 - SwePisc);
606 }
607
608
609 // find the loop on which the an effective
610 // saturation has to be
611 static ScanningCurve* findScanningCurve_(const Params& params, Scalar Sw)
612 {
613 if (params.pisc() == nullptr || Sw <= params.pisc()->Sw()) {
614 // we don't have a PISC yet, or the effective
615 // saturation is smaller than the saturation where the
616 // PISC begins. In this case are on the MDC
617 return params.mdc();
618 }
619
620 // if we have a primary imbibition curve, and our current
621 // effective saturation is higher than the beginning of
622 // the secondary drainage curve. this means we are on the
623 // PISC again.
624 if (params.pisc()->next() == nullptr ||
625 params.pisc()->next()->Sw() < Sw)
626 {
627 return params.pisc();
628 }
629
630 ScanningCurve* curve = params.csc()->next();
631 while (true) {
632 assert(curve != params.mdc()->prev());
633 if (curve->isValidAt_Sw(Sw)) {
634 return curve;
635 }
636 curve = curve->prev();
637 }
638 }
639};
640
641} // namespace Opm
642
643#endif // PARKER_LENHARD_HPP
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 &params)
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 &params, const FluidState &fs)
The relative permeability for the non-wetting phase of the params.
Definition ParkerLenhard.hpp:477
static Evaluation krw(const Params &params, const FluidState &fs)
The relative permeability for the wetting phase of the medium.
Definition ParkerLenhard.hpp:455
static void update(Params &params, const FluidState &fs)
Set the current absolute saturation for the current timestep.
Definition ParkerLenhard.hpp:310
static void relativePermeabilities(Container &values, const Params &params, 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 &params, 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 &params, const FluidState &fs)
Calculate the non-wetting phase saturations depending on the phase pressures.
Definition ParkerLenhard.hpp:443
static Evaluation pcnw(const Params &params, 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 &params, 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 &params, 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