ParkerLenhardParams.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) 2012-2013 by Andreas Lauser
5 
6  This file is part of the Open Porous Media project (OPM).
7 
8  OPM is free software: you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation, either version 2 of the License, or
11  (at your option) any later version.
12 
13  OPM is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with OPM. If not, see <http://www.gnu.org/licenses/>.
20 */
25 #ifndef OPM_PARKER_LENHARD_PARAMS_HPP
26 #define OPM_PARKER_LENHARD_PARAMS_HPP
27 
29 
30 #include <cassert>
31 
32 namespace Opm
33 {
34 // forward declaration
35 template <class ScalarT>
36 class PLScanningCurve;
37 
42 template <class TraitsT>
44 {
45 public:
46  typedef typename TraitsT::Scalar Scalar;
50 
52  {
53  currentSnr_ = 0;
54  mdc_ = new ScanningCurve(/*Swr=*/0);
55  pisc_ = csc_ = NULL;
56 
57 #ifndef NDEBUG
58  finalized_ = false;
59 #endif
60  }
61 
63  {
64  currentSnr_ = 0;
65  SwrPc_ = p.SwrPc_;
66  mdc_ = new ScanningCurve(SwrPc_);
67  pisc_ = csc_ = NULL;
68 
69 #ifndef NDEBUG
70  finalized_ = p.finalized_;
71 #endif
72  }
73 
75  { delete mdc_; }
76 
81  void finalize()
82  {
83 #ifndef NDEBUG
84  finalized_ = true;
85 #endif
86  }
87 
92  const VanGenuchtenParams &micParams() const
93  { assertFinalized_(); return *micParams_; }
94 
99  void setMicParams(const VanGenuchtenParams *val)
100  { micParams_ = val; }
101 
106  const VanGenuchtenParams &mdcParams() const
107  { assertFinalized_(); return *mdcParams_; }
108 
113  void setMdcParams(const VanGenuchtenParams *val)
114  { mdcParams_ = val; }
115 
119  Scalar Snr() const
120  { assertFinalized_(); return Snr_; }
121 
125  void setSnr(Scalar val)
126  { Snr_ = val; }
127 
131  Scalar SwrPc() const
132  { assertFinalized_(); return SwrPc_; }
133 
137  Scalar SwrKr() const
138  { assertFinalized_(); return SwrKr_; }
139 
144  void setSwr(Scalar pcSwr, Scalar krSwr = -1)
145  {
146  SwrPc_ = pcSwr;
147  if (krSwr < 0)
148  SwrKr_ = pcSwr;
149  else
150  SwrKr_ = krSwr;
151  }
152 
156  Scalar currentSnr() const
157  { assertFinalized_(); return currentSnr_; }
158 
162  void setCurrentSnr(Scalar val)
163  { currentSnr_ = val; }
164 
168  ScanningCurve *mdc() const
169  { assertFinalized_(); return mdc_; }
170 
174  void setMdc(ScanningCurve *val)
175  { mdc_ = val; }
176 
180  ScanningCurve *pisc() const
181  { assertFinalized_(); return pisc_; }
182 
186  void setPisc(ScanningCurve *val)
187  { pisc_ = val; }
188 
192  ScanningCurve *csc() const
193  { assertFinalized_(); return csc_; }
194 
198  void setCsc(ScanningCurve *val)
199  { csc_ = val; }
200 
201 private:
202 #ifndef NDEBUG
203  void assertFinalized_() const
204  { assert(finalized_); }
205 
206  bool finalized_;
207 #else
208  void assertFinalized_() const
209  { }
210 #endif
211 
212  const VanGenuchtenParams *micParams_;
213  const VanGenuchtenParams *mdcParams_;
214  Scalar SwrPc_;
215  Scalar SwrKr_;
216  Scalar Snr_;
217  Scalar currentSnr_;
218  mutable ScanningCurve *mdc_;
219  mutable ScanningCurve *pisc_;
220  mutable ScanningCurve *csc_;
221 };
222 } // namespace Opm
223 
224 #endif
ParkerLenhardParams()
Definition: ParkerLenhardParams.hpp:51
void setMicParams(const VanGenuchtenParams *val)
Sets the parameters of the main imbibition curve (which uses the van Genuchten capillary pressure mod...
Definition: ParkerLenhardParams.hpp:99
VanGenuchten::Params VanGenuchtenParams
Definition: ParkerLenhardParams.hpp:48
Definition: Air_Mesitylene.hpp:31
void setPisc(ScanningCurve *val)
Set the primary imbibition scanning curve.
Definition: ParkerLenhardParams.hpp:186
const VanGenuchtenParams & mdcParams() const
Returns the parameters of the main drainage curve (which uses the van Genuchten capillary pressure mo...
Definition: ParkerLenhardParams.hpp:106
const VanGenuchtenParams & micParams() const
Returns the parameters of the main imbibition curve (which uses the van Genuchten capillary pressure ...
Definition: ParkerLenhardParams.hpp:92
ParamsT Params
Definition: RegularizedVanGenuchten.hpp:77
Scalar SwrKr() const
Returns wetting phase residual saturation for the residual saturation curves.
Definition: ParkerLenhardParams.hpp:137
~ParkerLenhardParams()
Definition: ParkerLenhardParams.hpp:74
void finalize()
Calculate all dependent quantities once the independent quantities of the parameter object have been ...
Definition: ParkerLenhardParams.hpp:81
Implementation of the regularized van Genuchten's capillary pressure / relative permeability <-> satu...
Definition: RegularizedVanGenuchten.hpp:71
Opm::RegularizedVanGenuchten< TraitsT > VanGenuchten
Definition: ParkerLenhardParams.hpp:47
ScanningCurve * csc() const
Returns the current scanning curve.
Definition: ParkerLenhardParams.hpp:192
PLScanningCurve< Scalar > ScanningCurve
Definition: ParkerLenhardParams.hpp:49
ScanningCurve * mdc() const
Returns the main drainage curve.
Definition: ParkerLenhardParams.hpp:168
Scalar SwrPc() const
Returns wetting phase residual saturation for the capillary pressure curve.
Definition: ParkerLenhardParams.hpp:131
Default parameter class for the Parker-Lenhard hysteresis model.
Definition: ParkerLenhardParams.hpp:43
Represents a scanning curve in the Parker-Lenhard hysteresis model.
Definition: ParkerLenhard.hpp:46
TraitsT::Scalar Scalar
Definition: ParkerLenhardParams.hpp:46
void setMdcParams(const VanGenuchtenParams *val)
Sets the parameters of the main drainage curve (which uses the van Genuchten capillary pressure model...
Definition: ParkerLenhardParams.hpp:113
Scalar Snr() const
Returns non-wetting phase residual saturation.
Definition: ParkerLenhardParams.hpp:119
Implementation of the regularized van Genuchten's capillary pressure / relative permeability <-> satu...
ScanningCurve * pisc() const
Returns the primary imbibition scanning curve.
Definition: ParkerLenhardParams.hpp:180
void setCurrentSnr(Scalar val)
Set the current effective residual saturation.
Definition: ParkerLenhardParams.hpp:162
void setSnr(Scalar val)
Set the non-wetting phase residual saturation.
Definition: ParkerLenhardParams.hpp:125
void setMdc(ScanningCurve *val)
Set the main drainage curve.
Definition: ParkerLenhardParams.hpp:174
ParkerLenhardParams(const ParkerLenhardParams &p)
Definition: ParkerLenhardParams.hpp:62
void setCsc(ScanningCurve *val)
Set the current scanning curve.
Definition: ParkerLenhardParams.hpp:198
void setSwr(Scalar pcSwr, Scalar krSwr=-1)
Set the wetting phase residual saturation for the capillary pressure and the relative permeabilities...
Definition: ParkerLenhardParams.hpp:144
Scalar currentSnr() const
Returns the current effective residual saturation.
Definition: ParkerLenhardParams.hpp:156