vangenuchten.hh 9.54 KB
Newer Older
1
2
3
4
5
6
/*****************************************************************************
 *   Copyright (C) 2008 by Andreas Lauser, Bernd Flemisch                    *
 *   Institute of Hydraulic Engineering                                      *
 *   University of Stuttgart, Germany                                        *
 *   email: <givenname>.<name>@iws.uni-stuttgart.de                          *
 *                                                                           *
7
 *   This program is free software: you can redistribute it and/or modify    *
8
 *   it under the terms of the GNU General Public License as published by    *
9
10
 *   the Free Software Foundation, either version 2 of the License, or       *
 *   (at your option) any later version.                                     *
11
 *                                                                           *
12
13
14
15
16
17
18
 *   This program is distributed in the hope that it will be useful,         *
 *   but WITHOUT ANY WARRANTY; without even the implied warranty of          *
 *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the           *
 *   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/>.   *
19
20
 *****************************************************************************/
/*!
Bernd Flemisch's avatar
Bernd Flemisch committed
21
22
 * \file
 *
23
24
 * \brief   Implementation of the capillary pressure and
 *          relative permeability <-> saturation relations according to van Genuchten.
25
26
27
28
29
30
31
32
 */
#ifndef VAN_GENUCHTEN_HH
#define VAN_GENUCHTEN_HH

#include "vangenuchtenparams.hh"

#include <algorithm>

33
34
#include <math.h>
#include <assert.h>
35
36
37
38

namespace Dumux
{
/*!
39
 * \ingroup fluidmatrixinteractionslaws
40
 *
41
 * \brief Implementation of the van Genuchten capillary pressure <->
42
43
 *        saturation relation. This class bundles the "raw" curves
 *        as static members and doesn't concern itself converting
44
 *        absolute to effective saturations and vice versa.
45
 *
46
47
48
 * For general info: EffToAbsLaw
 *
 * \see VanGenuchtenParams
49
50
51
52
53
 */
template <class ScalarT, class ParamsT = VanGenuchtenParams<ScalarT> >
class VanGenuchten
{
public:
54
55
    typedef ParamsT     Params;
    typedef typename    Params::Scalar Scalar;
56
57

    /*!
58
     * \brief The capillary pressure-saturation curve according to van Genuchten.
59
60
61
62
63
64
     *
     * Van Genuchten's empirical capillary pressure <-> saturation
     * function is given by
     * \f[
     p_C = (\overline{S}_w^{-1/m} - 1)^{1/n}/\alpha
     \f]
65
66
67
68
     * \param Swe       Effective saturation of the wetting phase \f$\overline{S}_w\f$
     * \param params    A container object that is populated with the appropriate coefficients for the respective law.
     *                  Therefore, in the (problem specific) spatialParameters  first, the material law is chosen, and then the params container
     *                  is constructed accordingly. Afterwards the values are set there, too.
69
     */
Bernd Flemisch's avatar
Bernd Flemisch committed
70
    static Scalar pC(const Params &params, Scalar Swe)
71
    {
Bernd Flemisch's avatar
Bernd Flemisch committed
72
73
        assert(0 <= Swe && Swe <= 1);
        return pow(pow(Swe, -1.0/params.vgM()) - 1, 1.0/params.vgN())/params.vgAlpha();
74
75
76
    }

    /*!
77
     * \brief The saturation-capillary pressure curve according to van Genuchten.
78
79
80
81
82
83
     *
     * This is the inverse of the capillary pressure-saturation curve:
     * \f[
     \overline{S}_w = {p_C}^{-1} = ((\alpha p_C)^n + 1)^{-m}
     \f]
     *
84
85
86
87
88
     * \param pC        Capillary pressure \f$p_C\f$
     * \param params    A container object that is populated with the appropriate coefficients for the respective law.
     *                  Therefore, in the (problem specific) spatialParameters  first, the material law is chosen, and then the params container
     *                  is constructed accordingly. Afterwards the values are set there, too.
     * \return          The effective saturation of the wetting phase \f$\overline{S}_w\f$
89
90
91
92
93
94
95
96
97
     */
    static Scalar Sw(const Params &params, Scalar pC)
    {
        assert(pC >= 0);

        return pow(pow(params.vgAlpha()*pC, params.vgN()) + 1, -params.vgM());
    }

    /*!
98
99
     * \brief The partial derivative of the capillary
     *        pressure w.r.t. the effective saturation according to van Genuchten.
100
101
102
103
104
105
106
     *
     * This is equivalent to
     * \f[
     \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
     \f]
Bernd Flemisch's avatar
Bernd Flemisch committed
107
     *
108
109
110
111
     * \param Swe       Effective saturation of the wetting phase \f$\overline{S}_w\f$
     * \param params    A container object that is populated with the appropriate coefficients for the respective law.
     *                  Therefore, in the (problem specific) spatialParameters  first, the material law is chosen, and then the params container
     *                  is constructed accordingly. Afterwards the values are set there, too.
112
    */
Bernd Flemisch's avatar
Bernd Flemisch committed
113
    static Scalar dpC_dSw(const Params &params, Scalar Swe)
114
    {
Bernd Flemisch's avatar
Bernd Flemisch committed
115
        assert(0 <= Swe && Swe <= 1);
116

117
118
119
        Scalar powSwe = pow(Swe, -1/params.vgM());
        return - 1/params.vgAlpha() * pow(powSwe - 1, 1/params.vgN() - 1)/params.vgN()
            * powSwe/Swe/params.vgM();
120
121
122
    }

    /*!
123
124
     * \brief The partial derivative of the effective
     *        saturation to the capillary pressure according to van Genuchten.
Bernd Flemisch's avatar
Bernd Flemisch committed
125
     *
126
127
128
129
     * \param pC        Capillary pressure \f$p_C\f$
     * \param params    A container object that is populated with the appropriate coefficients for the respective law.
     *                  Therefore, in the (problem specific) spatialParameters  first, the material law is chosen, and then the params container
     *                  is constructed accordingly. Afterwards the values are set there, too.
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
     */
    static Scalar dSw_dpC(const Params &params, Scalar pC)
    {
        assert(pC >= 0);

        Scalar powAlphaPc = pow(params.vgAlpha()*pC, params.vgN());
        return -pow(powAlphaPc + 1, -params.vgM()-1)*
            params.vgM()*powAlphaPc/pC*params.vgN();
    }

    /*!
     * \brief The relative permeability for the wetting phase of
     *        the medium implied by van Genuchten's
     *        parameterization.
     *
145
146
147
148
     * \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.
     *                  Therefore, in the (problem specific) spatialParameters  first, the material law is chosen, and then the params container
     *                  is constructed accordingly. Afterwards the values are set there, too.     */
Bernd Flemisch's avatar
Bernd Flemisch committed
149
    static Scalar krw(const Params &params, Scalar Swe)
150
    {
Bernd Flemisch's avatar
Bernd Flemisch committed
151
        assert(0 <= Swe && Swe <= 1);
152

Bernd Flemisch's avatar
Bernd Flemisch committed
153
154
        Scalar r = 1. - pow(1 - pow(Swe, 1/params.vgM()), params.vgM());
        return sqrt(Swe)*r*r;
155
156
157
158
159
160
161
    };

    /*!
     * \brief The derivative of the relative permeability for the
     *        wetting phase in regard to the wetting saturation of the
     *        medium implied by the van Genuchten parameterization.
     *
162
163
164
165
     * \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.
     *                  Therefore, in the (problem specific) spatialParameters  first, the material law is chosen, and then the params container
     *                  is constructed accordingly. Afterwards the values are set there, too.
166
     */
Bernd Flemisch's avatar
Bernd Flemisch committed
167
    static Scalar dkrw_dSw(const Params &params, Scalar Swe)
168
    {
Bernd Flemisch's avatar
Bernd Flemisch committed
169
        assert(0 <= Swe && Swe <= 1);
170

Bernd Flemisch's avatar
Bernd Flemisch committed
171
        const Scalar x = 1 - std::pow(Swe, 1.0/params.vgM());
172
        const Scalar xToM = std::pow(x, params.vgM());
173
        return (1 - xToM)/std::sqrt(Swe) * ( (1 - xToM)/2 + 2*xToM*(1-x)/x );
174
175
176
177
178
179
180
181
    };


    /*!
     * \brief The relative permeability for the non-wetting phase
     *        of the medium implied by van Genuchten's
     *        parameterization.
     *
182
183
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.
     *                  Therefore, in the (problem specific) spatialParameters  first, the material law is chosen, and then the params container
     *                  is constructed accordingly. Afterwards the values are set there, too.
186
     */
Bernd Flemisch's avatar
Bernd Flemisch committed
187
    static Scalar krn(const Params &params, Scalar Swe)
188
    {
Bernd Flemisch's avatar
Bernd Flemisch committed
189
        assert(0 <= Swe && Swe <= 1);
190
191

        return
Bernd Flemisch's avatar
Bernd Flemisch committed
192
193
            pow(1 - Swe, 1.0/3) *
            pow(1 - pow(Swe, 1/params.vgM()), 2*params.vgM());
194
195
196
197
198
199
200
201
    };

    /*!
     * \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.
     *
202
203
204
205
     * \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.
     *                  Therefore, in the (problem specific) spatialParameters  first, the material law is chosen, and then the params container
     *                  is constructed accordingly. Afterwards the values are set there, too.
206
     */
Bernd Flemisch's avatar
Bernd Flemisch committed
207
    static Scalar dkrn_dSw(const Params &params, Scalar Swe)
208
    {
Bernd Flemisch's avatar
Bernd Flemisch committed
209
        assert(0 <= Swe && Swe <= 1);
210

Bernd Flemisch's avatar
Bernd Flemisch committed
211
        const Scalar x = std::pow(Swe, 1.0/params.vgM());
212
213
        return
            -std::pow(1 - x, 2*params.vgM())
214
            *std::pow(1 - Swe, -2/3)
Bernd Flemisch's avatar
Bernd Flemisch committed
215
            *(1.0/3 + 2*x/Swe);
216
217
218
219
220
221
    }

};
}

#endif