RegularizedBrooksCoreyParams.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 
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_REGULARIZED_BROOKS_COREY_PARAMS_HPP
26 #define OPM_REGULARIZED_BROOKS_COREY_PARAMS_HPP
27 
28 #include "BrooksCorey.hpp"
29 #include "BrooksCoreyParams.hpp"
30 
31 #include <dune/common/deprecated.hh>
32 
33 #include <cassert>
34 
35 namespace Opm {
42 template <class TraitsT>
44 {
47  typedef typename TraitsT::Scalar Scalar;
48 
49 public:
50  typedef TraitsT Traits;
51 
53  : BrooksCoreyParams()
54  , pcnwLowSw_(0.01)
55  {
56 #ifndef NDEBUG
57  finalized_ = false;
58 #endif
59  }
60 
62  : BrooksCoreyParams(entryPressure, lambda)
63  , pcnwLowSw_(0.01)
64  { finalize(); }
65 
70  void finalize()
71  {
73 
74  pcnwLow_ = BrooksCorey::twoPhaseSatPcnw(*this, pcnwLowSw_);
75  pcnwSlopeLow_ = dPcnw_dSw_(pcnwLowSw_);
76  pcnwHigh_ = BrooksCorey::twoPhaseSatPcnw(*this, 1.0);
77  pcnwSlopeHigh_ = dPcnw_dSw_(1.0);
78 
79 #ifndef NDEBUG
80  finalized_ = true;
81 #endif
82  }
83 
88  Scalar pcnwLowSw() const
89  { assertFinalized_(); return pcnwLowSw_; }
90 
95  Scalar pcnwLow() const
96  { assertFinalized_(); return pcnwLow_; }
97 
104  Scalar pcnwSlopeLow() const
105  { assertFinalized_(); return pcnwSlopeLow_; }
106 
111  void setPcLowSw(Scalar value)
112  { pcnwLowSw_ = value; }
113 
114  void setThresholdSw(Scalar value) DUNE_DEPRECATED_MSG("this method has been renamed to setPcLowSw()")
115  { pcnwLowSw_ = value; }
116 
121  Scalar pcnwHigh() const
122  { assertFinalized_(); return pcnwHigh_; }
123 
130  Scalar pcnwSlopeHigh() const
131  { assertFinalized_(); return pcnwSlopeHigh_; }
132 
133 private:
134 #ifndef NDEBUG
135  void assertFinalized_() const
136  { assert(finalized_); }
137 
138  bool finalized_;
139 #else
140  void assertFinalized_() const
141  { }
142 #endif
143 
144  Scalar dPcnw_dSw_(Scalar Sw) const
145  {
146  // use finite differences to calculate the derivative w.r.t. Sw of the
147  // unregularized curve's capillary pressure.
148  const Scalar eps = 1e-7;
149  Scalar delta = 0.0;
150  Scalar pc1;
151  Scalar pc2;
152  if (Sw + eps < 1.0) {
153  pc2 = BrooksCorey::twoPhaseSatPcnw(*this, Sw + eps);
154  delta += eps;
155  }
156  else
157  pc2 = BrooksCorey::twoPhaseSatPcnw(*this, Sw);
158 
159  if (Sw - eps > 0.0) {
160  pc1 = BrooksCorey::twoPhaseSatPcnw(*this, Sw - eps);
161  delta += eps;
162  }
163  else
164  pc1 = BrooksCorey::twoPhaseSatPcnw(*this, Sw);
165 
166  return (pc2 - pc1)/delta;
167  }
168 
169  Scalar pcnwLowSw_;
170  Scalar pcnwLow_;
171  Scalar pcnwSlopeLow_;
172  Scalar pcnwHigh_;
173  Scalar pcnwSlopeHigh_;
174 };
175 } // namespace Opm
176 
177 #endif
TraitsT Traits
Definition: RegularizedBrooksCoreyParams.hpp:50
Implementation of the Brooks-Corey capillary pressure <-> saturation relation.
Definition: Air_Mesitylene.hpp:31
Scalar pcnwSlopeHigh() const
Return the slope capillary pressure curve if Sw is larger or equal to 1.
Definition: RegularizedBrooksCoreyParams.hpp:130
Scalar lambda() const
Returns the lambda shape parameter.
Definition: BrooksCoreyParams.hpp:92
void setThresholdSw(Scalar value) DUNE_DEPRECATED_MSG("this method has been renamed to setPcLowSw()")
Definition: RegularizedBrooksCoreyParams.hpp:114
Scalar pcnwSlopeLow() const
Return the slope capillary pressure curve if Sw is smaller or equal to the low threshold saturation...
Definition: RegularizedBrooksCoreyParams.hpp:104
RegularizedBrooksCoreyParams(Scalar entryPressure, Scalar lambda)
Definition: RegularizedBrooksCoreyParams.hpp:61
RegularizedBrooksCoreyParams()
Definition: RegularizedBrooksCoreyParams.hpp:52
Specification of the material parameters for the Brooks-Corey constitutive relations.
Implementation of the Brooks-Corey capillary pressure <-> saturation relation.
Definition: BrooksCorey.hpp:54
Scalar pcnwLow() const
Return the capillary pressure at the low threshold saturation of the wetting phase.
Definition: RegularizedBrooksCoreyParams.hpp:95
void setPcLowSw(Scalar value)
Set the threshold saturation below which the capillary pressure is regularized.
Definition: RegularizedBrooksCoreyParams.hpp:111
Parameters that are necessary for the regularization of the Brooks-Corey capillary pressure model...
Definition: RegularizedBrooksCoreyParams.hpp:43
Scalar pcnwHigh() const
Return the capillary pressure at the high threshold saturation of the wetting phase.
Definition: RegularizedBrooksCoreyParams.hpp:121
void finalize()
Calculate all dependent quantities once the independent quantities of the parameter object have been ...
Definition: BrooksCoreyParams.hpp:69
Specification of the material parameters for the Brooks-Corey constitutive relations.
Definition: BrooksCoreyParams.hpp:44
void finalize()
Calculate all dependent quantities once the independent quantities of the parameter object have been ...
Definition: RegularizedBrooksCoreyParams.hpp:70
BrooksCoreyParams()
Definition: BrooksCoreyParams.hpp:51
Scalar entryPressure() const
Returns the entry pressure [Pa].
Definition: BrooksCoreyParams.hpp:79
Scalar pcnwLowSw() const
Return the threshold saturation below which the capillary pressure is regularized.
Definition: RegularizedBrooksCoreyParams.hpp:88