vangenuchten.hh 12.8 KB
Newer Older
1
2
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// 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
12
 *   This program is distributed in the hope that it will be useful,         *
 *   but WITHOUT ANY WARRANTY; without even the implied warranty of          *
13
 *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the            *
14
15
16
17
 *   GNU General Public License for more details.                            *
 *                                                                           *
 *   You should have received a copy of the GNU General Public License       *
 *   along with this program.  If not, see <http://www.gnu.org/licenses/>.   *
18
19
 *****************************************************************************/
/*!
Bernd Flemisch's avatar
Bernd Flemisch committed
20
 * \file
21
 * \ingroup Fluidmatrixinteractions
22
23
 * \brief   Implementation of the capillary pressure and
 *          relative permeability <-> saturation relations according to van Genuchten.
24
 */
Timo Koch's avatar
Timo Koch committed
25
26
#ifndef DUMUX_VAN_GENUCHTEN_HH
#define DUMUX_VAN_GENUCHTEN_HH
27
28
29
30

#include "vangenuchtenparams.hh"

#include <algorithm>
31
32
#include <cmath>
#include <cassert>
33

34
35
namespace Dumux {

36
/*!
37
 * \ingroup Fluidmatrixinteractions
38
 * \brief Implementation of the van Genuchten capillary pressure <->
39
40
 *        saturation relation. This class bundles the "raw" curves
 *        as static members and doesn't concern itself converting
41
 *        absolute to effective saturations and vice versa.
42
 *
43
44
 * For general info: EffToAbsLaw
 *
45
46
47
 * \note Capillary pressure model from van Genuchten (1980),
 *       relative permeability model from Mualem (1976)
 *
48
 * \see VanGenuchtenParams
49
50
51
52
53
 */
template <class ScalarT, class ParamsT = VanGenuchtenParams<ScalarT> >
class VanGenuchten
{
public:
54
55
    using Params = ParamsT;
    using Scalar = typename Params::Scalar;
56
57

    /*!
58
     * \brief The capillary pressure-saturation curve according to van Genuchten.
59
60
61
     *
     * Van Genuchten's empirical capillary pressure <-> saturation
     * function is given by
Alexander Kissinger's avatar
Alexander Kissinger committed
62
     * \f$\mathrm{
63
     p_C = (\overline{S}_w^{-1/m} - 1)^{1/n}/\alpha
Alexander Kissinger's avatar
Alexander Kissinger committed
64
65
66
     }\f$
     * \param swe Effective saturation of the wetting phase \f$\mathrm{\overline{S}_w}\f$
     * \param params A container object that is populated with the appropriate coefficients for the respective law.
67
     *                  Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container
68
     *                  is constructed accordingly. Afterwards the values are set there, too.
69
70
     * \note Instead of undefined behaviour if swe is not in the valid range, we return a valid number,
     *       by clamping the input. Note that for pc(swe = 0.0) = inf, have a look at RegularizedVanGenuchten if this is a problem.
71
     */
72
    static Scalar pc(const Params &params, Scalar swe)
73
    {
74
75
76
77
78
79
80
81
        using std::pow;
        using std::min;
        using std::max;

        swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0

        const Scalar pc = pow(pow(swe, -1.0/params.vgm()) - 1, 1.0/params.vgn())/params.vgAlpha();
        return pc;
82
83
84
    }

    /*!
85
     * \brief The saturation-capillary pressure curve according to van Genuchten.
86
87
     *
     * This is the inverse of the capillary pressure-saturation curve:
Alexander Kissinger's avatar
Alexander Kissinger committed
88
     * \f$\mathrm{
89
     \overline{S}_w = {p_C}^{-1} = ((\alpha p_C)^n + 1)^{-m}
Alexander Kissinger's avatar
Alexander Kissinger committed
90
     }\f$
91
     *
Alexander Kissinger's avatar
Alexander Kissinger committed
92
93
     * \param pc Capillary pressure \f$\mathrm{p_C}\f$ in \f$\mathrm{[Pa]}\f$
     * \param params A container object that is populated with the appropriate coefficients for the respective law.
94
     *                  Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container
95
     *                  is constructed accordingly. Afterwards the values are set there, too.
Alexander Kissinger's avatar
Alexander Kissinger committed
96
     * \return          The effective saturation of the wetting phase \f$\mathrm{\overline{S}_w}\f$
97
98
     * \note Instead of undefined behaviour if pc is not in the valid range, we return a valid number,
     *       i.e. sw(pc < 0.0) = 0.0, by clamping the input to the physical bounds.
99
     */
100
    static Scalar sw(const Params &params, Scalar pc)
101
    {
102
103
104
105
        using std::pow;
        using std::max;

        pc = max(pc, 0.0); // the equation below is undefined for negative pcs
106

107
108
        const Scalar sw = pow(pow(params.vgAlpha()*pc, params.vgn()) + 1, -params.vgm());
        return sw;
109
110
    }

111
112
113
114
    /*!
     * \brief The capillary pressure at Swe = 1.0 also called end point capillary pressure
     *
     * \param params A container object that is populated with the appropriate coefficients for the respective law.
115
116
     *                  Therefore, in the (problem specific) spatialParameters first,
     *                  the material law is chosen, and then the params container
117
118
119
120
121
     *                  is constructed accordingly. Afterwards the values are set there, too.
     */
    static Scalar endPointPc(const Params &params)
    { return 0.0; }

122
    /*!
123
124
     * \brief The partial derivative of the capillary
     *        pressure w.r.t. the effective saturation according to van Genuchten.
125
126
     *
     * This is equivalent to
Alexander Kissinger's avatar
Alexander Kissinger committed
127
     * \f$\mathrm{
128
129
130
     \frac{\partial p_C}{\partial \overline{S}_w} =
     -\frac{1}{\alpha} (\overline{S}_w^{-1/m} - 1)^{1/n - }
     \overline{S}_w^{-1/m} / \overline{S}_w / m
Alexander Kissinger's avatar
Alexander Kissinger committed
131
     }\f$
Bernd Flemisch's avatar
Bernd Flemisch committed
132
     *
Alexander Kissinger's avatar
Alexander Kissinger committed
133
134
     * \param swe Effective saturation of the wetting phase \f$\mathrm{\overline{S}_w}\f$
     * \param params A container object that is populated with the appropriate coefficients for the respective law.
135
136
     *                  Therefore, in the (problem specific) spatialParameters
     *                  first, the material law is chosen, and then the params container
137
     *                  is constructed accordingly. Afterwards the values are set there, too.
138
139
140
141
     *
     * \note Instead of undefined behaviour if swe is not in the valid range, we return a valid number,
     *       by clamping the input.
     */
142
    static Scalar dpc_dswe(const Params &params, Scalar swe)
143
    {
144
145
146
        using std::pow;
        using std::min;
        using std::max;
147

148
149
150
        swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0

        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
    }

155
    /*!
156
157
     * \brief The partial derivative of the effective
     *        saturation to the capillary pressure according to van Genuchten.
Bernd Flemisch's avatar
Bernd Flemisch committed
158
     *
Alexander Kissinger's avatar
Alexander Kissinger committed
159
160
     * \param pc Capillary pressure \f$\mathrm{p_C}\f$ in \f$\mathrm{[Pa]}\f$
     * \param params A container object that is populated with the appropriate coefficients for the respective law.
161
162
     *                  Therefore, in the (problem specific) spatialParameters
     *                  first, the material law is chosen, and then the params container
163
     *                  is constructed accordingly. Afterwards the values are set there, too.
164
165
166
     *
     * \note Instead of undefined behaviour if pc is not in the valid range, we return a valid number,
     *       by clamping the input.
167
     */
168
    static Scalar dswe_dpc(const Params &params, Scalar pc)
169
    {
170
171
172
173
        using std::pow;
        using std::max;

        pc = max(pc, 0.0); // the equation below is undefined for negative pcs
174

175
176
        const Scalar powAlphaPc = pow(params.vgAlpha()*pc, params.vgn());
        return -pow(powAlphaPc + 1, -params.vgm()-1)*params.vgm()*powAlphaPc/pc*params.vgn();
177
178
    }

179
180
    /*!
     * \brief The relative permeability for the wetting phase of
181
     *        the medium implied by van Genuchten / Mualem
182
183
     *        parameterization.
     *
Alexander Kissinger's avatar
Alexander Kissinger committed
184
185
     * \param swe The mobile saturation of the wetting phase.
     * \param params A container object that is populated with the appropriate coefficients for the respective law.
186
187
     *                  Therefore, in the (problem specific) spatialParameters
     *                  first, the material law is chosen, and then the params container
188
189
190
191
192
     *                  is constructed accordingly. Afterwards the values are set there, too.
     *
     * \note Instead of undefined behaviour if pc is not in the valid range, we return a valid number,
     *       by clamping the input.
     */
193
    static Scalar krw(const Params &params, Scalar swe)
194
    {
195
196
197
198
199
        using std::pow;
        using std::min;
        using std::max;

        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
205
206
207

    /*!
     * \brief The derivative of the relative permeability for the
     *        wetting phase in regard to the wetting saturation of the
208
     *        medium implied by the van Genuchten / Mualem parameterization.
209
     *
Alexander Kissinger's avatar
Alexander Kissinger committed
210
211
     * \param swe The mobile saturation of the wetting phase.
     * \param params A container object that is populated with the appropriate coefficients for the respective law.
212
213
     *                  Therefore, in the (problem specific) spatialParameters
     *                  first, the material law is chosen, and then the params container
214
     *                  is constructed accordingly. Afterwards the values are set there, too.
215
216
217
     *
     * \note Instead of undefined behaviour if pc is not in the valid range, we return a valid number,
     *       by clamping the input.
218
     */
219
    static Scalar dkrw_dswe(const Params &params, Scalar swe)
220
    {
221
222
223
        using std::pow;
        using std::min;
        using std::max;
224

225
226
227
228
        swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0

        const Scalar x = 1.0 - pow(swe, 1.0/params.vgm());
        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
232
233
234
235
236

    /*!
     * \brief The relative permeability for the non-wetting phase
     *        of the medium implied by van Genuchten's
     *        parameterization.
     *
Alexander Kissinger's avatar
Alexander Kissinger committed
237
238
     * \param swe The mobile saturation of the wetting phase.
     * \param params A container object that is populated with the appropriate coefficients for the respective law.
239
240
     *                  Therefore, in the (problem specific) spatialParameters
     *                  first, the material law is chosen, and then the params container
241
     *                  is constructed accordingly. Afterwards the values are set there, too.
242
243
244
     *
     * \note Instead of undefined behaviour if pc is not in the valid range, we return a valid number,
     *       by clamping the input.
245
     * \note See e.g. Dury, Fischer, Schulin (1999) for application of Mualem model to non-wetting rel. perm. 
246
     */
247
    static Scalar krn(const Params &params, Scalar swe)
248
    {
249
250
251
        using std::pow;
        using std::min;
        using std::max;
252

253
254
        swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0

255
        return pow(1 - swe, params.vgl()) * pow(1 - pow(swe, 1.0/params.vgm()), 2*params.vgm());
256
    }
257
258
259
260
261
262
263

    /*!
     * \brief The derivative of the relative permeability for the
     *        non-wetting phase in regard to the wetting saturation of
     *        the medium as implied by the van Genuchten
     *        parameterization.
     *
Alexander Kissinger's avatar
Alexander Kissinger committed
264
265
     * \param swe The mobile saturation of the wetting phase.
     * \param params A container object that is populated with the appropriate coefficients for the respective law.
266
267
     *                  Therefore, in the (problem specific) spatialParameters
     *                  first, the material law is chosen, and then the params container
268
     *                  is constructed accordingly. Afterwards the values are set there, too.
269
270
271
     *
     * \note Instead of undefined behaviour if pc is not in the valid range, we return a valid number,
     *       by clamping the input.
272
     */
273
    static Scalar dkrn_dswe(const Params &params, Scalar swe)
274
    {
275
276
277
        using std::pow;
        using std::min;
        using std::max;
278

279
280
        swe = min(max(swe, 0.0), 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0

281
282
        const auto sne = 1.0 - swe;
        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