3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
parkervangen3p.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 *****************************************************************************/
24#ifndef PARKERVANGEN_3P_HH
25#define PARKERVANGEN_3P_HH
26
27#include <algorithm>
28
30
31namespace Dumux {
32
42template <class ScalarT, class ParamsT = ParkerVanGen3PParams<ScalarT> >
44{
45
46public:
47 using Params = ParamsT;
48 using Scalar = typename Params::Scalar;
49
55 static Scalar pc(const Params &params, const Scalar sw)
56 {
57 DUNE_THROW(Dune::NotImplemented, "Capillary pressures for three phases is not so simple! Use pcgn, pcnw, and pcgw");
58 }
65 static Scalar pcgw(const Params &params, const Scalar swe)
66 {
67 assert(0 <= swe && swe <= 1);
68 return pc_(params, swe);
69 }
70
76 static Scalar pcnw(const Params &params, const Scalar swe)
77 {
78 assert(0 <= swe && swe <= 1);
79 return pc_(params, swe)/params.betaNw();
80 }
81
87 static Scalar pcgn(const Params &params, const Scalar ste)
88 {
89 assert(0 <= ste && ste <= 1);
90 return pc_(params,ste)/params.betaGn();
91 }
92
98 static Scalar pcAlpha(const Params &params, Scalar sne)
99 {
100 Scalar alpha;
101
102 /* regularization */
103 if (sne<=0.001) sne=0.0;
104 if (sne>=1.0) sne=1.0;
105
106 if (sne>params.snr())
107 alpha = 1.0;
108 else
109 {
110 if (params.snr()>=0.001)
111 alpha = sne/params.snr();
112 else
113 alpha = 0.0;
114 }
115 return alpha;
116 }
117
124 static Scalar sw(const Params &params, const Scalar pc)
125 {
126 DUNE_THROW(Dune::NotImplemented, "sw(pc) for three phases not implemented! Do it yourself!");
127 }
128
135 static Scalar dpc_dswe(const Params &params, const Scalar swe)
136 {
137 DUNE_THROW(Dune::NotImplemented, "dpc/dswe for three phases not implemented! Do it yourself!");
138 }
139
146 static Scalar dpcgw_dswe(const Params &params, const Scalar seRegu)
147 {
148 using std::pow;
149 const Scalar powSeRegu = pow(seRegu, -1/params.vgm());
150 return - 1.0/params.vgAlpha() * pow(powSeRegu - 1, 1.0/params.vgn() - 1)/params.vgn()
151 * powSeRegu/seRegu/params.vgm()/params.betaGw();
152 }
153
160 static Scalar dpcnw_dswe(const Params &params, const Scalar seRegu)
161 {
162 using std::pow;
163 const Scalar powSeRegu = pow(seRegu, -1/params.vgm());
164 return - 1.0/params.vgAlpha() * pow(powSeRegu - 1, 1.0/params.vgn() - 1)/params.vgn()
165 * powSeRegu/seRegu/params.vgm()/params.betaNw();
166 }
167
174 static Scalar dpcgn_dste(const Params &params, const Scalar seRegu)
175 {
176 using std::pow;
177 const Scalar powSeRegu = pow(seRegu, -1/params.vgm());
178 return - 1.0/params.vgAlpha() * pow(powSeRegu - 1, 1.0/params.vgn() - 1)/params.vgn()
179 * powSeRegu/seRegu/params.vgm()/params.betaGn();
180 }
181
188 static Scalar dswe_dpc(const Params &params, const Scalar pc)
189 {
190 DUNE_THROW(Dune::NotImplemented, "dswe/dpc for three phases not implemented! Do it yourself!");
191 }
192
205 static Scalar krw(const Params &params, const Scalar swe)
206 {
207 using std::pow;
208 using std::sqrt;
209 const Scalar r = 1.0 - pow(1 - pow(swe, 1/params.vgm()), params.vgm());
210 return sqrt(swe)*r*r;
211 }
212
230 static Scalar krn(const Params &params, const Scalar swe, const Scalar sn, const Scalar ste)
231 {
232 Scalar krn;
233 using std::pow;
234 krn = pow(1 - pow(swe, 1/params.vgm()), params.vgm());
235 krn -= pow(1 - pow(ste, 1/params.vgm()), params.vgm());
236 krn *= krn;
237
238 using std::max;
239 using std::min;
240 using std::sqrt;
241 if (params.krRegardsSnr())
242 {
243 // regard Snr in the permeability of the n-phase, see Helmig1997
244 Scalar resIncluded = max(min((sn - params.snr()/ (1-params.swr())), 1.0), 0.0);
245 krn *= sqrt(resIncluded );
246 }
247 else
248 krn *= sqrt(sn / (1 - params.swr())); // Hint: (ste - swe) = sn / (1-Swr)
249
250 return krn;
251 }
252
265 static Scalar krg(const Params &params, const Scalar ste)
266 {
267 assert(0 <= ste && ste <= 1);
268 using std::cbrt;
269 using std::pow;
270 return cbrt(1 - ste) * pow(1 - pow(ste, 1/params.vgm()), 2*params.vgm());
271 }
272
284 static Scalar dkrg_dste(const Params &params, Scalar ste)
285 {
286 assert(0 < ste && ste <= 1);
287
288 using std::pow;
289 const Scalar x = pow(ste, 1.0/params.vgm());
290 return
291 -pow(1.0 - x, 2*params.vgm())
292 *pow(1.0 - ste, -2.0/3)
293 *(1.0/3 + 2*x/ste);
294 }
295
304 static Scalar kr(const Params &params, const int phaseIdx, const Scalar swe, const Scalar sn, const Scalar ste)
305 {
306 switch (phaseIdx)
307 {
308 case 0:
309 return krw(params, swe);
310 case 1:
311 return krn(params, swe, sn, ste);
312 case 2:
313 return krg(params, ste);
314 }
315 DUNE_THROW(Dune::InvalidStateException,
316 "Invalid phase index ");
317 }
318
324 {
325 return params.rhoBulk() * params.KdNAPL();
326 }
327
328private:
329
336 const static Scalar pc_(const Params &params, const Scalar se)
337 {
338 using std::pow;
339 return pow(pow(se, -1/params.vgm()) - 1, 1/params.vgn())/params.vgAlpha();
340 }
341
342};
343} // end namespace Dumux
344
345#endif
Specification of the material params for the van Genuchten capillary pressure model.
Definition: adapt.hh:29
Implementation of van Genuchten's capillary pressure <-> saturation relation. This class bundles the ...
Definition: parkervangen3p.hh:44
static Scalar pc(const Params &params, const Scalar sw)
The capillary pressure-saturation curve.
Definition: parkervangen3p.hh:55
static Scalar krg(const Params &params, const Scalar ste)
The relative permeability for the non-wetting phase of the medium implied by van Genuchten's paramete...
Definition: parkervangen3p.hh:265
static Scalar pcAlpha(const Params &params, Scalar sne)
This function ensures a continuous transition from 2 to 3 phases and vice versa.
Definition: parkervangen3p.hh:98
static Scalar sw(const Params &params, const Scalar pc)
The saturation-capillary pressure curve.
Definition: parkervangen3p.hh:124
static Scalar dpcgn_dste(const Params &params, const Scalar seRegu)
Returns the partial derivative of the capillary pressure to the effective saturation.
Definition: parkervangen3p.hh:174
static Scalar kr(const Params &params, const int phaseIdx, const Scalar swe, const Scalar sn, const Scalar ste)
The relative permeability for a phase.
Definition: parkervangen3p.hh:304
static Scalar krn(const Params &params, const Scalar swe, const Scalar sn, const Scalar ste)
The relative permeability for the non-wetting phase after the Model of Parker et al....
Definition: parkervangen3p.hh:230
static Scalar dswe_dpc(const Params &params, const Scalar pc)
Returns the partial derivative of the effective saturation to the capillary pressure.
Definition: parkervangen3p.hh:188
static Scalar bulkDensTimesAdsorpCoeff(const Params &params)
the basis for calculating adsorbed NAPL in storage term
Definition: parkervangen3p.hh:323
static Scalar krw(const Params &params, const Scalar swe)
The relative permeability for the wetting phase of the medium implied by van Genuchten's parameteriza...
Definition: parkervangen3p.hh:205
static Scalar dpc_dswe(const Params &params, const Scalar swe)
Returns the partial derivative of the capillary pressure to the effective saturation.
Definition: parkervangen3p.hh:135
static Scalar dpcnw_dswe(const Params &params, const Scalar seRegu)
Returns the partial derivative of the capillary pressure to the effective saturation.
Definition: parkervangen3p.hh:160
static Scalar pcgn(const Params &params, const Scalar ste)
The capillary pressure-saturation curve for the gas and non-wetting phase.
Definition: parkervangen3p.hh:87
static Scalar pcgw(const Params &params, const Scalar swe)
The capillary pressure-saturation curve for the gas and wetting phase.
Definition: parkervangen3p.hh:65
typename Params::Scalar Scalar
Definition: parkervangen3p.hh:48
static Scalar pcnw(const Params &params, const Scalar swe)
The capillary pressure-saturation curve for the non-wettigng and wetting phase.
Definition: parkervangen3p.hh:76
static Scalar dpcgw_dswe(const Params &params, const Scalar seRegu)
Returns the partial derivative of the capillary pressure to the effective saturation.
Definition: parkervangen3p.hh:146
static Scalar dkrg_dste(const Params &params, Scalar ste)
The derivative of the relative permeability for the gas phase in regard to the total liquid saturatio...
Definition: parkervangen3p.hh:284
ParamsT Params
Definition: parkervangen3p.hh:47