3.2-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 DUMUX_VAN_GENUCHTEN_HH
26#define DUMUX_VAN_GENUCHTEN_HH
27
28#include "vangenuchtenparams.hh"
29
30#include <algorithm>
31#include <cmath>
32#include <cassert>
33
34namespace Dumux {
35
50template <class ScalarT, class ParamsT = VanGenuchtenParams<ScalarT> >
52{
53public:
54 using Params = ParamsT;
55 using Scalar = typename Params::Scalar;
56
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 const Scalar pc = pow(pow(swe, -1.0/params.vgm()) - 1, 1.0/params.vgn())/params.vgAlpha();
81 return pc;
82 }
83
100 static Scalar sw(const Params &params, Scalar pc)
101 {
102 using std::pow;
103 using std::max;
104
105 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
106
107 const Scalar sw = pow(pow(params.vgAlpha()*pc, params.vgn()) + 1, -params.vgm());
108 return sw;
109 }
110
119 static Scalar endPointPc(const Params &params)
120 { return 0.0; }
121
142 static Scalar dpc_dswe(const Params &params, Scalar swe)
143 {
144 using std::pow;
145 using std::min;
146 using std::max;
147
148 swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
149
150 const Scalar powSwe = pow(swe, -1/params.vgm());
151 return - 1.0/params.vgAlpha() * pow(powSwe - 1, 1.0/params.vgn() - 1)/params.vgn()
152 * powSwe/swe/params.vgm();
153 }
154
168 static Scalar dswe_dpc(const Params &params, Scalar pc)
169 {
170 using std::pow;
171 using std::max;
172
173 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
174
175 const Scalar powAlphaPc = pow(params.vgAlpha()*pc, params.vgn());
176 return -pow(powAlphaPc + 1, -params.vgm()-1)*params.vgm()*powAlphaPc/pc*params.vgn();
177 }
178
193 static Scalar krw(const Params &params, Scalar swe)
194 {
195 using std::pow;
196 using std::min;
197 using std::max;
198
199 swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
200
201 const Scalar r = 1.0 - pow(1.0 - pow(swe, 1.0/params.vgm()), params.vgm());
202 return pow(swe, params.vgl())*r*r;
203 }
204
219 static Scalar dkrw_dswe(const Params &params, Scalar swe)
220 {
221 using std::pow;
222 using std::min;
223 using std::max;
224
225 swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
226
227 const Scalar x = 1.0 - pow(swe, 1.0/params.vgm());
228 const Scalar xToM = pow(x, params.vgm());
229 return (1.0 - xToM)*pow(swe, params.vgl()-1) * ( (1.0 - xToM)*params.vgl() + 2*xToM*(1.0-x)/x );
230 }
231
247 static Scalar krn(const Params &params, Scalar swe)
248 {
249 using std::pow;
250 using std::min;
251 using std::max;
252
253 swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
254
255 return pow(1 - swe, params.vgl()) * pow(1 - pow(swe, 1.0/params.vgm()), 2*params.vgm());
256 }
257
273 static Scalar dkrn_dswe(const Params &params, Scalar swe)
274 {
275 using std::pow;
276 using std::min;
277 using std::max;
278
279 swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
280
281 const auto sne = 1.0 - swe;
282 const auto x = 1.0 - pow(swe, 1.0/params.vgm());
283 return -pow(sne, params.vgl()-1.0) * pow(x, 2*params.vgm() - 1.0) * ( params.vgl()*x + 2.0*sne/swe*(1.0 - x) );
284 }
285
286};
287
288} // end namespace Dumux
289
290#endif
Specification of the material parameters for the van Genuchten-Mualem constitutive relations.
Definition: adapt.hh:29
Implementation of the van Genuchten capillary pressure <-> saturation relation. This class bundles th...
Definition: vangenuchten.hh:52
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:219
static Scalar sw(const Params &params, Scalar pc)
The saturation-capillary pressure curve according to van Genuchten.
Definition: vangenuchten.hh:100
typename Params::Scalar Scalar
Definition: vangenuchten.hh:55
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:142
static Scalar krw(const Params &params, Scalar swe)
The relative permeability for the wetting phase of the medium implied by van Genuchten / Mualem param...
Definition: vangenuchten.hh:193
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:168
static Scalar endPointPc(const Params &params)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: vangenuchten.hh:119
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:247
ParamsT Params
Definition: vangenuchten.hh:54
static Scalar pc(const Params &params, Scalar swe)
The capillary pressure-saturation curve according to van Genuchten.
Definition: vangenuchten.hh:72
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:273