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  Copyright (C) 2008-2013 by Andreas Lauser
5  Copyright (C) 2010 by Philipp Nuske
6 
7  This file is part of the Open Porous Media project (OPM).
8 
9  OPM is free software: you can redistribute it and/or modify
10  it under the terms of the GNU General Public License as published by
11  the Free Software Foundation, either version 2 of the License, or
12  (at your option) any later version.
13 
14  OPM is distributed in the hope that it will be useful,
15  but WITHOUT ANY WARRANTY; without even the implied warranty of
16  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  GNU General Public License for more details.
18 
19  You should have received a copy of the GNU General Public License
20  along with OPM. If not, see <http://www.gnu.org/licenses/>.
21 */
26 #ifndef OPM_PARKER_LENHARD_HPP
27 #define OPM_PARKER_LENHARD_HPP
28 
29 #include "ParkerLenhardParams.hpp"
30 
32 
33 #include <algorithm>
34 #include <cassert>
35 
36 namespace Opm {
37 
45 template <class ScalarT>
47 {
48 public:
49  typedef ScalarT Scalar;
50 
57  PLScanningCurve(Scalar Swr)
58  {
59  loopNum_ = 0;
60  prev_ = new PLScanningCurve(NULL, // prev
61  this, // next
62  -1, // loop number
63  Swr, // Sw
64  1e12, // pcnw
65  Swr, // SwMic
66  Swr); // SwMdc
67  next_ = NULL;
68 
69  Sw_ = 1.0;
70  pcnw_ = 0.0;
71  SwMic_ = 1.0;
72  SwMdc_ = 1.0;
73  }
74 
75 protected:
78  int loopN,
79  Scalar Sw,
80  Scalar pcnw,
81  Scalar SwMic,
82  Scalar SwMdc)
83  {
84  prev_ = prev;
85  next_ = next;
86  loopNum_ = loopN;
87  Sw_ = Sw;
88  pcnw_ = pcnw;
89  SwMic_ = SwMic;
90  SwMdc_ = SwMdc;
91  }
92 
93 public:
100  {
101  if (loopNum_ == 0)
102  delete prev_;
103  if (loopNum_ >= 0)
104  delete next_;
105  }
106 
112  { return prev_; }
113 
119  { return next_; }
120 
129  void setNext(Scalar Sw,
130  Scalar pcnw,
131  Scalar SwMic,
132  Scalar SwMdc)
133  {
134  // if next_ is NULL, delete does nothing, so
135  // this is valid!!
136  delete next_;
137 
138  next_ = new PLScanningCurve(this, // prev
139  NULL, // next
140  loopNum() + 1,
141  Sw,
142  pcnw,
143  SwMic,
144  SwMdc);
145  }
146 
153  bool isValidAt_Sw(Scalar Sw)
154  {
155  if (isImbib())
156  // for inbibition the given saturation
157  // must be between the start of the
158  // current imbibition and the the start
159  // of the last drainage
160  return this->Sw() < Sw && Sw < prev_->Sw();
161  else
162  // for drainage the given saturation
163  // must be between the start of the
164  // last imbibition and the start
165  // of the current drainage
166  return prev_->Sw() < Sw && Sw < this->Sw();
167  }
168 
173  bool isImbib()
174  { return loopNum()%2 == 1; }
175 
180  bool isDrain()
181  { return !isImbib(); }
182 
188  int loopNum()
189  { return loopNum_; }
190 
195  Scalar Sw() const
196  { return Sw_; }
197 
201  Scalar pcnw() const
202  { return pcnw_; }
203 
208  Scalar SwMic()
209  { return SwMic_; }
210 
215  Scalar SwMdc()
216  { return SwMdc_; }
217 
218 private:
219  PLScanningCurve *prev_;
220  PLScanningCurve *next_;
221 
222  int loopNum_;
223 
224  Scalar Sw_;
225  Scalar pcnw_;
226 
227  Scalar SwMdc_;
228  Scalar SwMic_;
229 };
230 
237 template <class TraitsT, class ParamsT = ParkerLenhardParams<TraitsT> >
238 class ParkerLenhard : public TraitsT
239 {
240 public:
241  typedef TraitsT Traits;
242  typedef ParamsT Params;
243  typedef typename Traits::Scalar Scalar;
244 
246  static const int numPhases = Traits::numPhases;
247  static_assert(numPhases == 2,
248  "The Parker-Lenhard capillary pressure law only "
249  "applies to the case of two fluid phases");
250 
253  static const bool implementsTwoPhaseApi = true;
254 
257  static const bool implementsTwoPhaseSatApi = true;
258 
261  static const bool isSaturationDependent = true;
262 
265  static const bool isPressureDependent = false;
266 
269  static const bool isTemperatureDependent = false;
270 
273  static const bool isCompositionDependent = false;
274 
275  static_assert(Traits::numPhases == 2,
276  "The number of fluid phases must be two if you want to use "
277  "this material law!");
278 
279 private:
280  typedef typename ParamsT::VanGenuchten VanGenuchten;
281  typedef Opm::PLScanningCurve<Scalar> ScanningCurve;
282 
283 public:
288  static void reset(Params &params)
289  {
290  delete params.mdc(); // this will work even if mdc_ == NULL!
291  params.setMdc(new ScanningCurve(params.SwrPc()));
292  params.setCsc(params.mdc());
293  params.setPisc(NULL);
294  params.setCurrentSnr(0.0);
295  }
296 
301  template <class FluidState>
302  static void update(Params &params, const FluidState &fs)
303  {
305 
306  Scalar Sw = FsToolbox::value(fs.saturation(Traits::wettingPhaseIdx));
307 
308  if (Sw > 1 - 1e-5) {
309  // if the absolute saturation is almost 1,
310  // it means that we're back to the beginning
311  reset(params);
312  return;
313  }
314 
315  // find the loop number which corrosponds to the
316  // given effective saturation
317  ScanningCurve *curve = findScanningCurve_(params, Sw);
318 
319  // calculate the apparent saturation on the MIC and MDC
320  // which yield the same capillary pressure as the
321  // Sw at the current scanning curve
322  Scalar pc = pcnw<FluidState, Scalar>(params, fs);
323  Scalar Sw_mic = VanGenuchten::twoPhaseSatSw(params.micParams(), pc);
324  Scalar Sw_mdc = VanGenuchten::twoPhaseSatSw(params.mdcParams(), pc);
325 
326  curve->setNext(Sw, pc, Sw_mic, Sw_mdc);
327  if (!curve->next())
328  return;
329 
330  params.setCsc(curve);
331 
332  // if we're back on the MDC, we also have a new PISC!
333  if (params.csc() == params.mdc()) {
334  params.setPisc(params.mdc()->next());
335  params.setCurrentSnr(computeCurrentSnr_(params, Sw));
336  }
337  }
338 
343  template <class Container, class FluidState>
344  static void capillaryPressures(Container &values, const Params &params, const FluidState &fs)
345  {
346  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
347 
348  values[Traits::wettingPhaseIdx] = 0.0; // reference phase
349  values[Traits::nonWettingPhaseIdx] = pcnw<FluidState, Evaluation>(params, fs);
350  }
351 
356  template <class Container, class FluidState>
357  static void saturations(Container &/*values*/, const Params &/*params*/, const FluidState &/*fs*/)
358  { OPM_THROW(std::logic_error, "Not implemented: ParkerLenhard::saturations()"); }
359 
364  template <class Container, class FluidState>
365  static void relativePermeabilities(Container &values, const Params &params, const FluidState &fs)
366  {
367  typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
368 
369  values[Traits::wettingPhaseIdx] = krw<FluidState, Evaluation>(params, fs);
370  values[Traits::nonWettingPhaseIdx] = krn<FluidState, Evaluation>(params, fs);
371  }
372 
377  template <class FluidState, class Evaluation = typename FluidState::Scalar>
378  static Evaluation pcnw(const Params &params, const FluidState &fs)
379  {
380  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
381 
382  const Evaluation& Sw =
383  FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
384 
385  return twoPhaseSatPcnw(params, Sw);
386  }
387 
388  template <class Evaluation>
389  static Evaluation twoPhaseSatPcnw(const Params &params, const Evaluation& Sw)
390  {
391  typedef MathToolbox<Evaluation> Toolbox;
392 
393  // calculate the current apparent saturation
394  ScanningCurve *sc = findScanningCurve_(params, Toolbox::value(Sw));
395 
396  // calculate the apparant saturation
397  const Evaluation& Sw_app = absoluteToApparentSw_(params, Sw);
398 
399  // if the apparent saturation exceeds the 'legal' limits,
400  // we also the underlying material law decide what to do.
401  if (Sw_app > 1) {
402  return 0.0; // VanGenuchten::pcnw(params.mdcParams(), Sw_app);
403  }
404 
405  // put the apparent saturation into the main imbibition or
406  // drainage curve
407  Scalar SwAppCurSC = absoluteToApparentSw_(params, sc->Sw());
408  Scalar SwAppPrevSC = absoluteToApparentSw_(params, sc->prev()->Sw());
409  const Evaluation& pos = (Sw_app - SwAppCurSC)/(SwAppPrevSC - SwAppCurSC);
410  if (sc->isImbib()) {
411  const Evaluation& SwMic =
412  pos * (sc->prev()->SwMic() - sc->SwMic()) + sc->SwMic();
413 
414  return VanGenuchten::twoPhaseSatPcnw(params.micParams(), SwMic);
415  }
416  else { // sc->isDrain()
417  const Evaluation& SwMdc =
418  pos*(sc->prev()->SwMdc() - sc->SwMdc()) + sc->SwMdc();
419 
420  return VanGenuchten::twoPhaseSatPcnw(params.mdcParams(), SwMdc);
421  }
422  }
423 
428  template <class FluidState, class Evaluation = typename FluidState::Scalar>
429  static Evaluation Sw(const Params &/*params*/, const FluidState &/*fs*/)
430  { OPM_THROW(std::logic_error, "Not implemented: ParkerLenhard::Sw()"); }
431 
432  template <class Evaluation>
433  static Evaluation twoPhaseSatSw(const Params &/*params*/, const Evaluation& /*pc*/)
434  { OPM_THROW(std::logic_error, "Not implemented: ParkerLenhard::twoPhaseSatSw()"); }
435 
440  template <class FluidState, class Evaluation = typename FluidState::Scalar>
441  static Evaluation Sn(const Params &params, const FluidState &fs)
442  { return 1 - Sw<FluidState, Evaluation>(params, fs); }
443 
444  template <class Evaluation>
445  static Evaluation twoPhaseSatSn(const Params &/*params*/, const Evaluation& /*pc*/)
446  { OPM_THROW(std::logic_error, "Not implemented: ParkerLenhard::twoPhaseSatSn()"); }
447 
452  template <class FluidState, class Evaluation = typename FluidState::Scalar>
453  static Evaluation krw(const Params &params, const FluidState &fs)
454  {
455  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
456 
457  const Evaluation& Sw =
458  FsToolbox::template toLhs<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::Scalar>
477  static Evaluation krn(const Params &params, const FluidState &fs)
478  {
479  typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
480 
481  const Evaluation& Sw =
482  FsToolbox::template toLhs<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
483 
484  return twoPhaseSatKrn(params, Sw);
485  }
486 
487  template <class Evaluation>
488  static Evaluation twoPhaseSatKrn(const Params &params, const Evaluation& Sw)
489  {
490  // for the effective permeability we only use Land's law and
491  // the relative permeability of the main drainage curve.
492  const Evaluation& Sw_app = absoluteToApparentSw_(params, Sw);
493  return VanGenuchten::twoPhaseSatKrn(params.mdcParams(), Sw_app);
494  }
495 
499  template <class Evaluation>
500  static Evaluation absoluteToApparentSw_(const Params &params, const Evaluation& Sw)
501  {
502  return effectiveToApparentSw_(params, absoluteToEffectiveSw_(params, Sw));
503  }
504 
505 private:
515  template <class Evaluation>
516  static Evaluation absoluteToEffectiveSw_(const Params &params, const Evaluation& Sw)
517  { return (Sw - params.SwrPc())/(1 - params.SwrPc()); }
518 
528  template <class Evaluation>
529  static Evaluation effectiveToAbsoluteSw_(const Params &params, const Evaluation& Swe)
530  { return Swe*(1 - params.SwrPc()) + params.SwrPc(); }
531 
532  // return the effctive residual non-wetting saturation, given an
533  // effective wetting saturation
534  template <class Evaluation>
535  static Evaluation computeCurrentSnr_(const Params &params, const Evaluation& Sw)
536  {
537  // regularize
538  if (Sw > 1 - params.Snr())
539  return 0.0;
540  if (Sw < params.SwrPc())
541  return params.Snr();
542 
543  if (params.Snr() == 0.0)
544  return 0.0;
545 
546  // use Land's law
547  Scalar R = 1.0/params.Snr() - 1;
548  const Evaluation& curSnr = (1 - Sw)/(1 + R*(1 - Sw));
549 
550  // the current effective residual non-wetting saturation must
551  // be smaller than the residual non-wetting saturation
552  assert(curSnr <= params.Snr());
553 
554  return curSnr;
555  }
556 
557  // returns the trapped effective non-wetting saturation for a
558  // given wetting phase saturation
559  template <class Evaluation>
560  static Evaluation trappedEffectiveSn_(const Params &params, const Evaluation& Sw)
561  {
562  const Evaluation& Swe = absoluteToEffectiveSw_(params, Sw);
563  Scalar SwePisc = absoluteToEffectiveSw_(params, params.pisc()->Sw());
564 
565  Scalar Snre = absoluteToEffectiveSw_(params, params.currentSnr());
566  return Snre*(Swe - SwePisc) / (1 - Snre - SwePisc);
567  }
568 
569  // returns the apparent saturation of the wetting phase depending
570  // on the effective saturation
571  template <class Evaluation>
572  static Evaluation effectiveToApparentSw_(const Params &params, const Evaluation& Swe)
573  {
574  if (params.pisc() == NULL ||
575  Swe <= absoluteToEffectiveSw_(params, params.pisc()->Sw()))
576  {
577  // we are on the main drainage curve, i.e. no non-wetting
578  // fluid is trapped -> apparent saturation == effective
579  // saturation
580  return Swe;
581  }
582 
583  // we are on a imbibition or drainage curve which is not the
584  // main drainage curve -> apparent saturation == effective
585  // saturation + trapped effective saturation
586  return Swe + trappedEffectiveSn_(params, Swe);
587  }
588 
589  // Returns the effective saturation to a given apparent one
590  template <class Evaluation>
591  static Evaluation apparentToEffectiveSw_(const Params &params, const Evaluation& Swapp)
592  {
593  Scalar SwePisc = absoluteToEffectiveSw_(params, params.pisc()->Sw());
594  if (params.pisc() == NULL || Swapp <= SwePisc) {
595  // we are on the main drainage curve, i.e.
596  // no non-wetting fluid is trapped
597  // -> apparent saturation == effective saturation
598  return Swapp;
599  }
600 
601  Scalar Snre = absoluteToEffectiveSw_(params.currentSnr());
602  return
603  (Swapp*(1 - Snre - SwePisc) + Snre*SwePisc)
604  /(1 - SwePisc);
605  }
606 
607 
608  // find the loop on which the an effective
609  // saturation has to be
610  static ScanningCurve *findScanningCurve_(const Params &params, Scalar Sw)
611  {
612  if (params.pisc() == NULL || Sw <= params.pisc()->Sw()) {
613  // we don't have a PISC yet, or the effective
614  // saturation is smaller than the saturation where the
615  // PISC begins. In this case are on the MDC
616  return params.mdc();
617  }
618 
619  // if we have a primary imbibition curve, and our current
620  // effective saturation is higher than the beginning of
621  // the secondary drainage curve. this means we are on the
622  // PISC again.
623  if (params.pisc()->next() == NULL ||
624  params.pisc()->next()->Sw() < Sw)
625  {
626  return params.pisc();
627  }
628 
629  ScanningCurve *curve = params.csc()->next();
630  while (true) {
631  assert(curve != params.mdc()->prev());
632  if (curve->isValidAt_Sw(Sw)) {
633  return curve;
634  }
635  curve = curve->prev();
636  }
637  }
638 };
639 
640 } // namespace Opm
641 
642 #endif // PARKER_LENHARD_HPP
Scalar SwMic()
Apparent saturation of the last reversal point on the pressure MIC.
Definition: ParkerLenhard.hpp:208
static const int numPhases
The number of fluid phases.
Definition: ParkerLenhard.hpp:246
Definition: MathToolbox.hpp:39
Definition: Air_Mesitylene.hpp:31
Implements the Parker-Lenhard twophase p_c-Sw hysteresis model. This class adheres to the twophase ca...
Definition: ParkerLenhard.hpp:238
Implementation of the van Genuchten capillary pressure - saturation relation.
Definition: VanGenuchten.hpp:54
Implementation of the van Genuchten capillary pressure - saturation relation.
void setNext(Scalar Sw, Scalar pcnw, Scalar SwMic, Scalar SwMdc)
Set the next scanning curve.
Definition: ParkerLenhard.hpp:129
Traits::Scalar Scalar
Definition: ParkerLenhard.hpp:243
PLScanningCurve(Scalar Swr)
Constructs main imbibition curve.
Definition: ParkerLenhard.hpp:57
ScalarT Scalar
Definition: ParkerLenhard.hpp:49
Scalar pcnw() const
Capillary pressure at the last reversal point.
Definition: ParkerLenhard.hpp:201
bool isDrain()
Returns true iff the scanning curve is a drainage curve.
Definition: ParkerLenhard.hpp:180
int loopNum()
The loop number of the scanning curve.
Definition: ParkerLenhard.hpp:188
Default parameter class for the Parker-Lenhard hysteresis model.
~PLScanningCurve()
Destructor. After it was called all references to the next() curve are invalid!
Definition: ParkerLenhard.hpp:99
bool isValidAt_Sw(Scalar Sw)
Returns true iff the given effective saturation Swei is within the scope of the curve, i.e. whether Swei is part of the curve's domain and the curve thus applies to Swi.
Definition: ParkerLenhard.hpp:153
Scalar Sw() const
Absolute wetting-phase saturation at the scanning curve's reversal point.
Definition: ParkerLenhard.hpp:195
Represents a scanning curve in the Parker-Lenhard hysteresis model.
Definition: ParkerLenhard.hpp:46
PLScanningCurve * prev() const
Return the previous scanning curve, i.e. the curve with one less reversal than the current one...
Definition: ParkerLenhard.hpp:111
Scalar SwMdc()
Apparent saturation of the last reversal point on the pressure MDC.
Definition: ParkerLenhard.hpp:215
PLScanningCurve(PLScanningCurve *prev, PLScanningCurve *next, int loopN, Scalar Sw, Scalar pcnw, Scalar SwMic, Scalar SwMdc)
Definition: ParkerLenhard.hpp:76
ParamsT Params
Definition: ParkerLenhard.hpp:242
PLScanningCurve * next() const
Return the next scanning curve, i.e. the curve with one more reversal than the current one...
Definition: ParkerLenhard.hpp:118
TraitsT Traits
Definition: ParkerLenhard.hpp:241
bool isImbib()
Returns true iff the scanning curve is a imbibition curve.
Definition: ParkerLenhard.hpp:173