3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
brookscorey.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 DUMUX_BROOKS_COREY_HH
26#define DUMUX_BROOKS_COREY_HH
27
28#include "brookscoreyparams.hh"
29
30#include <algorithm>
31#include <cmath>
32
33namespace Dumux {
34
47template <class ScalarT, class ParamsT = BrooksCoreyParams<ScalarT> >
49{
50public:
51 using Params = ParamsT;
52 using Scalar = typename Params::Scalar;
53
72 static Scalar pc(const Params &params, Scalar swe)
73 {
74 using std::pow;
75 using std::min;
76 using std::max;
77
78 swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
79
80 return params.pe()*pow(swe, -1.0/params.lambda());
81 }
82
98 static Scalar sw(const Params &params, Scalar pc)
99 {
100 using std::pow;
101 using std::max;
102
103 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
104
105 return pow(pc/params.pe(), -params.lambda());
106 }
107
115 static Scalar endPointPc(const Params &params)
116 { return params.pe(); }
117
136 static Scalar dpc_dswe(const Params &params, Scalar swe)
137 {
138 using std::pow;
139 using std::min;
140 using std::max;
141
142 swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
143
144 return - params.pe()/params.lambda() * pow(swe, -1/params.lambda() - 1);
145 }
146
160 static Scalar dswe_dpc(const Params &params, Scalar pc)
161 {
162 using std::pow;
163 using std::max;
164
165 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
166
167 return -params.lambda()/params.pe() * pow(pc/params.pe(), - params.lambda() - 1);
168 }
169
184 static Scalar krw(const Params &params, Scalar swe)
185 {
186 using std::pow;
187 using std::min;
188 using std::max;
189
190 swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
191
192 return pow(swe, 2.0/params.lambda() + 3);
193 }
194
210 static Scalar dkrw_dswe(const Params &params, Scalar swe)
211 {
212 using std::pow;
213 using std::min;
214 using std::max;
215
216 swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
217
218 return (2.0/params.lambda() + 3)*pow(swe, 2.0/params.lambda() + 2);
219 }
220
235 static Scalar krn(const Params &params, Scalar swe)
236 {
237 using std::pow;
238 using std::min;
239 using std::max;
240
241 swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
242
243 const Scalar exponent = 2.0/params.lambda() + 1;
244 const Scalar tmp = 1.0 - swe;
245 return tmp*tmp*(1.0 - pow(swe, exponent));
246 }
247
264 static Scalar dkrn_dswe(const Params &params, Scalar swe)
265 {
266 using std::pow;
267 using std::min;
268 using std::max;
269
270 swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
271
272 return 2.0*(swe - 1)*(1 + pow(swe, 2.0/params.lambda())*(1.0/params.lambda() + 1.0/2
273 - swe*(1.0/params.lambda() + 1.0/2)
274 )
275 );
276 }
277
278};
279
280} // end namespace Dumux
281
282#endif // BROOKS_COREY_HH
Specification of the material parameters for the Brooks Corey constitutive relations.
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 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:210
static Scalar dswe_dpc(const Params &params, Scalar pc)
The partial derivative of the effective saturation w.r.t. the capillary pressure according to Brooks ...
Definition: brookscorey.hh:160
ParamsT Params
Definition: brookscorey.hh:51
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
typename Params::Scalar Scalar
Definition: brookscorey.hh:52
static Scalar endPointPc(const Params &params)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: brookscorey.hh:115
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
static Scalar dkrn_dswe(const Params &params, Scalar swe)
The derivative of the relative permeability for the non-wetting phase in regard to the wetting satura...
Definition: brookscorey.hh:264