3.3.0
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
31#warning "This header is deprecated. Removal after 3.3. Use new material laws."
32
33namespace Dumux {
34
44template <class ScalarT, class ParamsT = ParkerVanGen3PParams<ScalarT> >
46{
47
48public:
49 using Params = ParamsT;
50 using Scalar = typename Params::Scalar;
51
57 static Scalar pc(const Params &params, const Scalar sw)
58 {
59 DUNE_THROW(Dune::NotImplemented, "Capillary pressures for three phases is not so simple! Use pcgn, pcnw, and pcgw");
60 }
67 static Scalar pcgw(const Params &params, const Scalar swe)
68 {
69 assert(0 <= swe && swe <= 1);
70 return pc_(params, swe);
71 }
72
78 static Scalar pcnw(const Params &params, const Scalar swe)
79 {
80 assert(0 <= swe && swe <= 1);
81 return pc_(params, swe)/params.betaNw();
82 }
83
89 static Scalar pcgn(const Params &params, const Scalar ste)
90 {
91 assert(0 <= ste && ste <= 1);
92 return pc_(params,ste)/params.betaGn();
93 }
94
100 static Scalar pcAlpha(const Params &params, Scalar sne)
101 {
102 Scalar alpha;
103
104 /* regularization */
105 if (sne<=0.001) sne=0.0;
106 if (sne>=1.0) sne=1.0;
107
108 if (sne>params.snr())
109 alpha = 1.0;
110 else
111 {
112 if (params.snr()>=0.001)
113 alpha = sne/params.snr();
114 else
115 alpha = 0.0;
116 }
117 return alpha;
118 }
119
126 static Scalar sw(const Params &params, const Scalar pc)
127 {
128 DUNE_THROW(Dune::NotImplemented, "sw(pc) for three phases not implemented! Do it yourself!");
129 }
130
137 static Scalar dpc_dswe(const Params &params, const Scalar swe)
138 {
139 DUNE_THROW(Dune::NotImplemented, "dpc/dswe for three phases not implemented! Do it yourself!");
140 }
141
148 static Scalar dpcgw_dswe(const Params &params, const Scalar seRegu)
149 {
150 using std::pow;
151 const Scalar powSeRegu = pow(seRegu, -1/params.vgm());
152 return - 1.0/params.vgAlpha() * pow(powSeRegu - 1, 1.0/params.vgn() - 1)/params.vgn()
153 * powSeRegu/seRegu/params.vgm()/params.betaGw();
154 }
155
162 static Scalar dpcnw_dswe(const Params &params, const Scalar seRegu)
163 {
164 using std::pow;
165 const Scalar powSeRegu = pow(seRegu, -1/params.vgm());
166 return - 1.0/params.vgAlpha() * pow(powSeRegu - 1, 1.0/params.vgn() - 1)/params.vgn()
167 * powSeRegu/seRegu/params.vgm()/params.betaNw();
168 }
169
176 static Scalar dpcgn_dste(const Params &params, const Scalar seRegu)
177 {
178 using std::pow;
179 const Scalar powSeRegu = pow(seRegu, -1/params.vgm());
180 return - 1.0/params.vgAlpha() * pow(powSeRegu - 1, 1.0/params.vgn() - 1)/params.vgn()
181 * powSeRegu/seRegu/params.vgm()/params.betaGn();
182 }
183
190 static Scalar dswe_dpc(const Params &params, const Scalar pc)
191 {
192 DUNE_THROW(Dune::NotImplemented, "dswe/dpc for three phases not implemented! Do it yourself!");
193 }
194
207 static Scalar krw(const Params &params, const Scalar swe)
208 {
209 using std::pow;
210 using std::sqrt;
211 const Scalar r = 1.0 - pow(1 - pow(swe, 1/params.vgm()), params.vgm());
212 return sqrt(swe)*r*r;
213 }
214
232 static Scalar krn(const Params &params, const Scalar swe, const Scalar sn, const Scalar ste)
233 {
234 Scalar krn;
235 using std::pow;
236 krn = pow(1 - pow(swe, 1/params.vgm()), params.vgm());
237 krn -= pow(1 - pow(ste, 1/params.vgm()), params.vgm());
238 krn *= krn;
239
240 using std::clamp;
241 using std::sqrt;
242 if (params.krRegardsSnr())
243 {
244 // regard Snr in the permeability of the n-phase, see Helmig1997
245 const Scalar resIncluded = clamp(sn - params.snr()/ (1-params.swr()), 0.0, 1.0);
246 krn *= sqrt(resIncluded );
247 }
248 else
249 krn *= sqrt(sn / (1 - params.swr())); // Hint: (ste - swe) = sn / (1-Swr)
250
251 return krn;
252 }
253
266 static Scalar krg(const Params &params, const Scalar ste)
267 {
268 assert(0 <= ste && ste <= 1);
269 using std::cbrt;
270 using std::pow;
271 return cbrt(1 - ste) * pow(1 - pow(ste, 1/params.vgm()), 2*params.vgm());
272 }
273
285 static Scalar dkrg_dste(const Params &params, Scalar ste)
286 {
287 assert(0 < ste && ste <= 1);
288
289 using std::pow;
290 const Scalar x = pow(ste, 1.0/params.vgm());
291 return
292 -pow(1.0 - x, 2*params.vgm())
293 *pow(1.0 - ste, -2.0/3)
294 *(1.0/3 + 2*x/ste);
295 }
296
305 static Scalar kr(const Params &params, const int phaseIdx, const Scalar swe, const Scalar sn, const Scalar ste)
306 {
307 switch (phaseIdx)
308 {
309 case 0:
310 return krw(params, swe);
311 case 1:
312 return krn(params, swe, sn, ste);
313 case 2:
314 return krg(params, ste);
315 }
316 DUNE_THROW(Dune::InvalidStateException,
317 "Invalid phase index ");
318 }
319
325 {
326 return params.rhoBulk() * params.KdNAPL();
327 }
328
329private:
330
337 const static Scalar pc_(const Params &params, const Scalar se)
338 {
339 using std::pow;
340 return pow(pow(se, -1/params.vgm()) - 1, 1/params.vgn())/params.vgAlpha();
341 }
342
343};
344} // end namespace Dumux
345
346#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:46
static Scalar pc(const Params &params, const Scalar sw)
The capillary pressure-saturation curve.
Definition: parkervangen3p.hh:57
static Scalar krg(const Params &params, const Scalar ste)
The relative permeability for the nonwetting phase of the medium implied by van Genuchten's parameter...
Definition: parkervangen3p.hh:266
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:100
static Scalar sw(const Params &params, const Scalar pc)
The saturation-capillary pressure curve.
Definition: parkervangen3p.hh:126
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:176
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:305
static Scalar krn(const Params &params, const Scalar swe, const Scalar sn, const Scalar ste)
The relative permeability for the nonwetting phase after the Model of Parker et al....
Definition: parkervangen3p.hh:232
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:190
static Scalar bulkDensTimesAdsorpCoeff(const Params &params)
the basis for calculating adsorbed NAPL in storage term
Definition: parkervangen3p.hh:324
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:207
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:137
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:162
static Scalar pcgn(const Params &params, const Scalar ste)
The capillary pressure-saturation curve for the gas and nonwetting phase.
Definition: parkervangen3p.hh:89
static Scalar pcgw(const Params &params, const Scalar swe)
The capillary pressure-saturation curve for the gas and wetting phase.
Definition: parkervangen3p.hh:67
typename Params::Scalar Scalar
Definition: parkervangen3p.hh:50
static Scalar pcnw(const Params &params, const Scalar swe)
The capillary pressure-saturation curve for the non-wettigng and wetting phase.
Definition: parkervangen3p.hh:78
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:148
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:285
ParamsT Params
Definition: parkervangen3p.hh:49