3.3.0
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#warning "This header is deprecated. Removal after 3.3. Use new material laws."
29
30#include "brookscorey.hh"
32
34
35namespace Dumux {
36
62template <class ScalarT, class ParamsT = RegularizedBrooksCoreyParams<ScalarT> >
64{
66
67public:
68 using Params = ParamsT;
69 using Scalar = typename Params::Scalar;
70
84 static Scalar pc(const Params &params, Scalar swe)
85 {
86 const Scalar sThres = params.thresholdSw();
87
88 // make sure that the capilary pressure observes a
89 // derivative != 0 for 'illegal' saturations. This is
90 // required for example by newton solvers (if the
91 // derivative is calculated numerically) in order to get the
92 // saturation moving to the right direction if it
93 // temporarily is in an 'illegal' range.
94 if (swe <= sThres) {
95 Scalar m = BrooksCorey::dpc_dswe(params, sThres);
96 Scalar pcsweLow = BrooksCorey::pc(params, sThres);
97 return pcsweLow + m*(swe - sThres);
98 }
99 else if (swe > 1.0) {
100 Scalar m = BrooksCorey::dpc_dswe(params, 1.0);
101 Scalar pcsweHigh = BrooksCorey::pc(params, 1.0);
102 return pcsweHigh + m*(swe - 1.0);
103 }
104
105 // if the effective saturation is in an 'reasonable'
106 // range, we use the real Brooks-Corey law...
107 return BrooksCorey::pc(params, swe);
108 }
109
124 static Scalar sw(const Params &params, Scalar pc)
125 {
126 const Scalar sThres = params.thresholdSw();
127
128 // calculate the saturation which corrosponds to the
129 // saturation in the non-regularized version of
130 // the Brooks-Corey law
131 Scalar swe = BrooksCorey::sw(params, pc);
132
133 // make sure that the capilary pressure observes a
134 // derivative != 0 for 'illegal' saturations. This is
135 // required for example by newton solvers (if the
136 // derivative calculated numerically) in order to get the
137 // saturation moving to the right direction if it
138 // temporarily is in an 'illegal' range.
139 if (swe <= sThres) {
140 // invert the low saturation regularization of pc()
141 Scalar m = BrooksCorey::dpc_dswe(params, sThres);
142 Scalar pcsweLow = BrooksCorey::pc(params, sThres);
143 return sThres + (pc - pcsweLow)/m;
144 }
145 else if (swe > 1.0) {
146 Scalar m = BrooksCorey::dpc_dswe(params, 1.0);
147 Scalar pcsweHigh = BrooksCorey::pc(params, 1.0);
148 return 1.0 + (pc - pcsweHigh)/m;
149 }
150
151 return BrooksCorey::sw(params, pc);
152 }
153
161 static Scalar endPointPc(const Params &params)
162 { return params.pe(); }
163
178 static Scalar dpc_dswe(const Params &params, Scalar swe)
179 {
180 const Scalar sThres = params.thresholdSw();
181
182 // derivative of the regualarization
183 if (swe <= sThres) {
184 // calculate the slope of the straight line used in pc()
185 Scalar m = BrooksCorey::dpc_dswe(params, sThres);
186 return m;
187 }
188 else if (swe > 1.0) {
189 // calculate the slope of the straight line used in pc()
190 Scalar m = BrooksCorey::dpc_dswe(params, 1.0);
191 return m;
192 }
193
194 return BrooksCorey::dpc_dswe(params, swe);
195 }
196
211 static Scalar dswe_dpc(const Params &params, Scalar pc)
212 {
213 const Scalar sThres = params.thresholdSw();
214
215 //instead of return value = inf, return a very large number
216 if (params.pe() == 0.0)
217 {
218 return 1e100;
219 }
220
221 // calculate the saturation which corresponds to the
222 // saturation in the non-regularized version of the
223 // Brooks-Corey law
224 Scalar swe;
225 if (pc < 0)
226 swe = 1.5; // make sure we regularize below
227 else
228 swe = BrooksCorey::sw(params, pc);
229
230 // derivative of the regularization
231 if (swe <= sThres) {
232 // calculate the slope of the straight line used in pc()
233 Scalar m = BrooksCorey::dpc_dswe(params, sThres);
234 return 1/m;
235 }
236 else if (swe > 1.0) {
237 // calculate the slope of the straight line used in pc()
238 Scalar m = BrooksCorey::dpc_dswe(params, 1.0);
239 return 1/m;
240 }
241 return 1.0/BrooksCorey::dpc_dswe(params, swe);
242 }
243
258 static Scalar krw(const Params &params, Scalar swe)
259 {
260 if (swe <= 0.0)
261 return 0.0;
262 else if (swe >= 1.0)
263 return 1.0;
264
265 return BrooksCorey::krw(params, swe);
266 }
267
275 static Scalar dkrw_dswe(const Params &params, Scalar swe)
276 {
277 // derivative of the regularization
278 // the slope is zero below sw=0.0 and above sw=1.0
279 if (swe <= 0.0)
280 return 0.0;
281 else if (swe >= 1.0)
282 return 0.0;
283
284 return BrooksCorey::dkrw_dswe(params, swe);
285 }
286
301 static Scalar krn(const Params &params, Scalar swe)
302 {
303 if (swe >= 1.0)
304 return 0.0;
305 else if (swe <= 0.0)
306 return 1.0;
307
308 return BrooksCorey::krn(params, swe);
309 }
310
318 static Scalar dkrn_dswe(const Params &params, Scalar swe)
319 {
320 // derivative of the regularization
321 // the slope is zero below sw=0.0 and above sw=1.0
322 if (swe <= 0)
323 return 0.0;
324 else if (swe >= 1)
325 return 0.0;
326
327 return BrooksCorey::dkrn_dswe(params, swe);
328 }
329};
330} // end namespace Dumux
331
332#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.
Definition: adapt.hh:29
Implementation of the Brooks-Corey capillary pressure <-> saturation relation. This class bundles the...
Definition: brookscorey.hh:50
static Scalar dkrw_dswe(const Params &params, Scalar swe)
The derivative of the relative permeability for the wetting phase with regard to the wetting saturati...
Definition: brookscorey.hh:208
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:183
static Scalar pc(const Params &params, Scalar swe)
The capillary pressure-saturation curve according to Brooks & Corey.
Definition: brookscorey.hh:73
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 nonwetting phase of the medium as implied by the Brooks-Corey param...
Definition: brookscorey.hh:232
static Scalar dkrn_dswe(const Params &params, Scalar swe)
The derivative of the relative permeability for the nonwetting phase in regard to the wetting saturat...
Definition: brookscorey.hh:260
Implementation of the regularized Brooks-Corey capillary pressure / relative permeability <-> saturat...
Definition: regularizedbrookscorey.hh:64
static Scalar dkrn_dswe(const Params &params, Scalar swe)
A regularized version of the derivative of the relative permeability for the nonwetting phase in rega...
Definition: regularizedbrookscorey.hh:318
static Scalar pc(const Params &params, Scalar swe)
A regularized Brooks-Corey capillary pressure-saturation curve.
Definition: regularizedbrookscorey.hh:84
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:258
static Scalar krn(const Params &params, Scalar swe)
Regularized version of the relative permeability for the nonwetting phase of the medium implied by th...
Definition: regularizedbrookscorey.hh:301
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:178
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:211
static Scalar sw(const Params &params, Scalar pc)
A regularized Brooks-Corey saturation-capillary pressure curve.
Definition: regularizedbrookscorey.hh:124
ParamsT Params
Definition: regularizedbrookscorey.hh:68
typename Params::Scalar Scalar
Definition: regularizedbrookscorey.hh:69
static Scalar dkrw_dswe(const Params &params, Scalar swe)
A regularized version of the derivative of the relative permeability for the wetting phase in regard ...
Definition: regularizedbrookscorey.hh:275
static Scalar endPointPc(const Params &params)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: regularizedbrookscorey.hh:161