3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
regularizedbrookscorey.hh
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program 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 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef REGULARIZED_BROOKS_COREY_HH
26#define REGULARIZED_BROOKS_COREY_HH
27
28#include "brookscorey.hh"
30
32
33namespace Dumux {
34
60template <class ScalarT, class ParamsT = RegularizedBrooksCoreyParams<ScalarT> >
62{
64
65public:
66 using Params = ParamsT;
67 using Scalar = typename Params::Scalar;
68
82 static Scalar pc(const Params &params, Scalar swe)
83 {
84 const Scalar sThres = params.thresholdSw();
85
86 // make sure that the capilary pressure observes a
87 // derivative != 0 for 'illegal' saturations. This is
88 // required for example by newton solvers (if the
89 // derivative is calculated numerically) in order to get the
90 // saturation moving to the right direction if it
91 // temporarily is in an 'illegal' range.
92 if (swe <= sThres) {
93 Scalar m = BrooksCorey::dpc_dswe(params, sThres);
94 Scalar pcsweLow = BrooksCorey::pc(params, sThres);
95 return pcsweLow + m*(swe - sThres);
96 }
97 else if (swe > 1.0) {
98 Scalar m = BrooksCorey::dpc_dswe(params, 1.0);
99 Scalar pcsweHigh = BrooksCorey::pc(params, 1.0);
100 return pcsweHigh + m*(swe - 1.0);
101 }
102
103 // if the effective saturation is in an 'reasonable'
104 // range, we use the real Brooks-Corey law...
105 return BrooksCorey::pc(params, swe);
106 }
107
122 static Scalar sw(const Params &params, Scalar pc)
123 {
124 const Scalar sThres = params.thresholdSw();
125
126 // calculate the saturation which corrosponds to the
127 // saturation in the non-regularized version of
128 // the Brooks-Corey law
129 Scalar swe = BrooksCorey::sw(params, pc);
130
131 // make sure that the capilary pressure observes a
132 // derivative != 0 for 'illegal' saturations. This is
133 // required for example by newton solvers (if the
134 // derivative calculated numerically) in order to get the
135 // saturation moving to the right direction if it
136 // temporarily is in an 'illegal' range.
137 if (swe <= sThres) {
138 // invert the low saturation regularization of pc()
139 Scalar m = BrooksCorey::dpc_dswe(params, sThres);
140 Scalar pcsweLow = BrooksCorey::pc(params, sThres);
141 return sThres + (pc - pcsweLow)/m;
142 }
143 else if (swe > 1.0) {
144 Scalar m = BrooksCorey::dpc_dswe(params, 1.0);
145 Scalar pcsweHigh = BrooksCorey::pc(params, 1.0);
146 return 1.0 + (pc - pcsweHigh)/m;
147 }
148
149 return BrooksCorey::sw(params, pc);
150 }
151
159 static Scalar endPointPc(const Params &params)
160 { return params.pe(); }
161
176 static Scalar dpc_dswe(const Params &params, Scalar swe)
177 {
178 const Scalar sThres = params.thresholdSw();
179
180 // derivative of the regualarization
181 if (swe <= sThres) {
182 // calculate the slope of the straight line used in pc()
183 Scalar m = BrooksCorey::dpc_dswe(params, sThres);
184 return m;
185 }
186 else if (swe > 1.0) {
187 // calculate the slope of the straight line used in pc()
188 Scalar m = BrooksCorey::dpc_dswe(params, 1.0);
189 return m;
190 }
191
192 return BrooksCorey::dpc_dswe(params, swe);
193 }
194
209 static Scalar dswe_dpc(const Params &params, Scalar pc)
210 {
211 const Scalar sThres = params.thresholdSw();
212
213 //instead of return value = inf, return a very large number
214 if (params.pe() == 0.0)
215 {
216 return 1e100;
217 }
218
219 // calculate the saturation which corresponds to the
220 // saturation in the non-regularized version of the
221 // Brooks-Corey law
222 Scalar swe;
223 if (pc < 0)
224 swe = 1.5; // make sure we regularize below
225 else
226 swe = BrooksCorey::sw(params, pc);
227
228 // derivative of the regularization
229 if (swe <= sThres) {
230 // calculate the slope of the straight line used in pc()
231 Scalar m = BrooksCorey::dpc_dswe(params, sThres);
232 return 1/m;
233 }
234 else if (swe > 1.0) {
235 // calculate the slope of the straight line used in pc()
236 Scalar m = BrooksCorey::dpc_dswe(params, 1.0);
237 return 1/m;
238 }
239 return 1.0/BrooksCorey::dpc_dswe(params, swe);
240 }
241
256 static Scalar krw(const Params &params, Scalar swe)
257 {
258 if (swe <= 0.0)
259 return 0.0;
260 else if (swe >= 1.0)
261 return 1.0;
262
263 return BrooksCorey::krw(params, swe);
264 }
265
280 static Scalar krn(const Params &params, Scalar swe)
281 {
282 if (swe >= 1.0)
283 return 0.0;
284 else if (swe <= 0.0)
285 return 1.0;
286
287 return BrooksCorey::krn(params, swe);
288 }
289};
290} // end namespace Dumux
291
292#endif
Provides 3rd order polynomial splines.
Implementation of the capillary pressure and relative permeability <-> saturation relations according...
Parameters that are necessary for the regularization of the Brooks-Corey capillary pressure model.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Implementation of the Brooks-Corey capillary pressure <-> saturation relation. This class bundles the...
Definition: brookscorey.hh:49
static Scalar dpc_dswe(const Params &params, Scalar swe)
The partial derivative of the capillary pressure w.r.t. the effective saturation according to Brooks ...
Definition: brookscorey.hh:136
static Scalar krw(const Params &params, Scalar swe)
The relative permeability for the wetting phase of the medium implied by the Brooks-Corey parameteriz...
Definition: brookscorey.hh:184
static Scalar pc(const Params &params, Scalar swe)
The capillary pressure-saturation curve according to Brooks & Corey.
Definition: brookscorey.hh:72
static Scalar sw(const Params &params, Scalar pc)
The saturation-capillary pressure curve according to Brooks & Corey.
Definition: brookscorey.hh:98
static Scalar krn(const Params &params, Scalar swe)
The relative permeability for the non-wetting phase of the medium as implied by the Brooks-Corey para...
Definition: brookscorey.hh:235
Implementation of the regularized Brooks-Corey capillary pressure / relative permeability <-> saturat...
Definition: regularizedbrookscorey.hh:62
static Scalar pc(const Params &params, Scalar swe)
A regularized Brooks-Corey capillary pressure-saturation curve.
Definition: regularizedbrookscorey.hh:82
static Scalar krw(const Params &params, Scalar swe)
Regularized version of the relative permeability for the wetting phase of the medium implied by the B...
Definition: regularizedbrookscorey.hh:256
static Scalar krn(const Params &params, Scalar swe)
Regularized version of the relative permeability for the non-wetting phase of the medium implied by t...
Definition: regularizedbrookscorey.hh:280
static Scalar dpc_dswe(const Params &params, Scalar swe)
A regularized version of the partial derivative of the w.r.t. effective saturation according to Broo...
Definition: regularizedbrookscorey.hh:176
static Scalar dswe_dpc(const Params &params, Scalar pc)
A regularized version of the partial derivative of the w.r.t. cap.pressure according to Brooks & Cor...
Definition: regularizedbrookscorey.hh:209
static Scalar sw(const Params &params, Scalar pc)
A regularized Brooks-Corey saturation-capillary pressure curve.
Definition: regularizedbrookscorey.hh:122
ParamsT Params
Definition: regularizedbrookscorey.hh:66
typename Params::Scalar Scalar
Definition: regularizedbrookscorey.hh:67
static Scalar endPointPc(const Params &params)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: regularizedbrookscorey.hh:159