3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
vangenuchten.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 VAN_GENUCHTEN_HH
26#define VAN_GENUCHTEN_HH
27
28#include "vangenuchtenparams.hh"
29
30#include <algorithm>
31#include <cmath>
32#include <cassert>
33
34namespace Dumux {
35
47template <class ScalarT, class ParamsT = VanGenuchtenParams<ScalarT> >
49{
50public:
51 using Params = ParamsT;
52 using Scalar = typename Params::Scalar;
53
69 static Scalar pc(const Params &params, Scalar swe)
70 {
71 using std::pow;
72 using std::min;
73 using std::max;
74
75 swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
76
77 const Scalar pc = pow(pow(swe, -1.0/params.vgm()) - 1, 1.0/params.vgn())/params.vgAlpha();
78 return pc;
79 }
80
97 static Scalar sw(const Params &params, Scalar pc)
98 {
99 using std::pow;
100 using std::max;
101
102 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
103
104 const Scalar sw = pow(pow(params.vgAlpha()*pc, params.vgn()) + 1, -params.vgm());
105 return sw;
106 }
107
115 static Scalar endPointPc(const Params &params)
116 { return 0.0; }
117
137 static Scalar dpc_dswe(const Params &params, Scalar swe)
138 {
139 using std::pow;
140 using std::min;
141 using std::max;
142
143 swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
144
145 const Scalar powSwe = pow(swe, -1/params.vgm());
146 return - 1.0/params.vgAlpha() * pow(powSwe - 1, 1.0/params.vgn() - 1)/params.vgn()
147 * powSwe/swe/params.vgm();
148 }
149
162 static Scalar dswe_dpc(const Params &params, Scalar pc)
163 {
164 using std::pow;
165 using std::max;
166
167 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
168
169 const Scalar powAlphaPc = pow(params.vgAlpha()*pc, params.vgn());
170 return -pow(powAlphaPc + 1, -params.vgm()-1)*params.vgm()*powAlphaPc/pc*params.vgn();
171 }
172
186 static Scalar krw(const Params &params, Scalar swe)
187 {
188 using std::pow;
189 using std::sqrt;
190 using std::min;
191 using std::max;
192
193 swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
194
195 const Scalar r = 1.0 - pow(1.0 - pow(swe, 1.0/params.vgm()), params.vgm());
196 return sqrt(swe)*r*r;
197 }
198
212 static Scalar dkrw_dswe(const Params &params, Scalar swe)
213 {
214 using std::pow;
215 using std::sqrt;
216 using std::min;
217 using std::max;
218
219 swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
220
221 const Scalar x = 1.0 - pow(swe, 1.0/params.vgm());
222 const Scalar xToM = pow(x, params.vgm());
223 return (1.0 - xToM)/sqrt(swe) * ( (1.0 - xToM)/2 + 2*xToM*(1.0-x)/x );
224 }
225
239 static Scalar krn(const Params &params, Scalar swe)
240 {
241 using std::pow;
242 using std::min;
243 using std::max;
244
245 swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
246
247 return pow(1 - swe, 1.0/3) * pow(1 - pow(swe, 1.0/params.vgm()), 2*params.vgm());
248 }
249
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 const Scalar x = pow(swe, 1.0/params.vgm());
273 return -pow(1.0 - x, 2*params.vgm()) * pow(1.0 - swe, -2.0/3) * (1.0/3 + 2*x/swe);
274 }
275
276};
277
278} // end namespace Dumux
279
280#endif
Specification of the material parameters for the van Genuchten constitutive relations.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Implementation of the van Genuchten capillary pressure <-> saturation relation. This class bundles th...
Definition: vangenuchten.hh:49
static Scalar dkrw_dswe(const Params &params, Scalar swe)
The derivative of the relative permeability for the wetting phase in regard to the wetting saturation...
Definition: vangenuchten.hh:212
static Scalar sw(const Params &params, Scalar pc)
The saturation-capillary pressure curve according to van Genuchten.
Definition: vangenuchten.hh:97
typename Params::Scalar Scalar
Definition: vangenuchten.hh:52
static Scalar dpc_dswe(const Params &params, Scalar swe)
The partial derivative of the capillary pressure w.r.t. the effective saturation according to van Gen...
Definition: vangenuchten.hh:137
static Scalar krw(const Params &params, Scalar swe)
The relative permeability for the wetting phase of the medium implied by van Genuchten's parameteriza...
Definition: vangenuchten.hh:186
static Scalar dswe_dpc(const Params &params, Scalar pc)
The partial derivative of the effective saturation to the capillary pressure according to van Genucht...
Definition: vangenuchten.hh:162
static Scalar endPointPc(const Params &params)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: vangenuchten.hh:115
static Scalar krn(const Params &params, Scalar swe)
The relative permeability for the non-wetting phase of the medium implied by van Genuchten's paramete...
Definition: vangenuchten.hh:239
ParamsT Params
Definition: vangenuchten.hh:51
static Scalar pc(const Params &params, Scalar swe)
The capillary pressure-saturation curve according to van Genuchten.
Definition: vangenuchten.hh:69
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: vangenuchten.hh:264