3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
regularizedparkervangen3p.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 REGULARIZED_PARKERVANGEN_3P_HH
26#define REGULARIZED_PARKERVANGEN_3P_HH
27
28#include "parkervangen3p.hh"
30
32
33namespace Dumux {
34
60template <class ScalarT, class ParamsT = RegularizedParkerVanGen3PParams<ScalarT> >
62{
64
65public:
66 using Params = ParamsT;
67 using Scalar = typename Params::Scalar;
68
82 static Scalar pc(const Params &params, Scalar sw)
83 {
84 return ParkerVanGen3P::pc(params, sw);
85 }
86
93 static Scalar pcgw(const Params &params, Scalar swe)
94 {
95 //if specified, a constant value is used for regularization
96 if(params.constRegularization())
97 {
98 if(swe < 0.0)
99 swe = 0.0;
100 if(swe > 1.0)
101 swe = 1.0;
102 }
103
104 // retrieve the low and the high threshold saturations for the
105 // unregularized capillary pressure curve from the parameters
106 const Scalar swThLow = params.pcLowS();
107 const Scalar swThHigh = params.pcHighS();
108
109 // make sure that the capillary pressure observes a derivative
110 // != 0 for 'illegal' saturations. This is favourable for the
111 // newton solver (if the derivative is calculated numerically)
112 // in order to get the saturation moving to the right
113 // direction if it temporarily is in an 'illegal' range.
114 if (swe < swThLow)
115 {
116 const Scalar mLow = ParkerVanGen3P::dpcgw_dswe(params, swThLow);
117 return ParkerVanGen3P::pcgw(params, swThLow) + mLow*(swe - swThLow);
118 }
119 else if (swe > swThHigh)
120 {
121 Scalar yTh = ParkerVanGen3P::pcgw(params, swThHigh);
122 Scalar m1 = (0.0 - yTh)/(1.0 - swThHigh)*2;
123
124 if (swe < 1.0) {
125 // use spline between threshold swe and 1.0
126 Scalar mTh = ParkerVanGen3P::dpcgw_dswe(params, swThHigh);
127 Spline<Scalar> sp(swThHigh, 1.0, // x0, x1
128 yTh, 0, // y0, y1
129 mTh, m1); // m0, m1
130 return sp.eval(swe);
131 }
132 else {
133 // straight line for swe > 1.0
134 return m1*(swe - 1.0) + 0.0;
135 }
136 }
137
138 // if the effective saturation is in an 'reasonable'
139 // range, we use the real van genuchten law...
140 return ParkerVanGen3P::pcgw(params, swe);
141 }
142
148 static Scalar pcnw(const Params &params, Scalar swe)
149 {
150 //if specified, a constant value is used for regularization
151 if(params.constRegularization())
152 {
153 if(swe < 0.0)
154 swe = 0.0;
155 if(swe > 1.0)
156 swe = 1.0;
157 }
158
159 // retrieve the low and the high threshold saturations for the
160 // unregularized capillary pressure curve from the parameters
161 const Scalar swThLow = params.pcLowS();
162 const Scalar swThHigh = params.pcHighS();
163
164 // make sure that the capillary pressure observes a derivative
165 // != 0 for 'illegal' saturations. This is favourable for the
166 // newton solver (if the derivative is calculated numerically)
167 // in order to get the saturation moving to the right
168 // direction if it temporarily is in an 'illegal' range.
169 if (swe < swThLow)
170 {
171 const Scalar mLow = ParkerVanGen3P::dpcnw_dswe(params, swThLow);
172 return ParkerVanGen3P::pcnw(params, swThLow) + mLow*(swe - swThLow);
173 }
174 else if (swe > swThHigh)
175 {
176 Scalar yTh = ParkerVanGen3P::pcnw(params, swThHigh);
177 Scalar m1 = (0.0 - yTh)/(1.0 - swThHigh)*2;
178
179 if (swe < 1.0) {
180 // use spline between threshold swe and 1.0
181 Scalar mTh = ParkerVanGen3P::dpcnw_dswe(params, swThHigh);
182 Spline<Scalar> sp(swThHigh, 1.0, // x0, x1
183 yTh, 0, // y0, y1
184 mTh, m1); // m0, m1
185 return sp.eval(swe);
186 }
187 else {
188 // straight line for swe > 1.0
189 return m1*(swe - 1.0) + 0.0;
190 }
191 }
192
193 // if the effective saturation is in an 'reasonable'
194 // range, we use the real van genuchten law...
195 return ParkerVanGen3P::pcnw(params, swe);
196 }
197
203 static Scalar pcgn(const Params &params, Scalar ste)
204 {
205 //if specified, a constant value is used for regularization
206 if(params.constRegularization())
207 {
208 if(ste < 0.0)
209 ste = 0.0;
210 if(ste > 1.0)
211 ste = 1.0;
212 }
213
214 // retrieve the low and the high threshold saturations for the
215 // unregularized capillary pressure curve from the parameters
216 const Scalar swThLow = params.pcLowS();
217 const Scalar swThHigh = params.pcHighS();
218
219 // make sure that the capillary pressure observes a derivative
220 // != 0 for 'illegal' saturations. This is favourable for the
221 // newton solver (if the derivative is calculated numerically)
222 // in order to get the saturation moving to the right
223 // direction if it temporarily is in an 'illegal' range.
224 if (ste < swThLow)
225 {
226 const Scalar mLow = ParkerVanGen3P::dpcgn_dste(params, swThLow);
227 return ParkerVanGen3P::pcgn(params, swThLow) + mLow*(ste - swThLow);
228 }
229 else if (ste > swThHigh)
230 {
231 Scalar yTh = ParkerVanGen3P::pcgn(params, swThHigh);
232 Scalar m1 = (0.0 - yTh)/(1.0 - swThHigh)*2;
233
234 if (ste < 1.0) {
235 // use spline between threshold swe and 1.0
236 Scalar mTh = ParkerVanGen3P::dpcgn_dste(params, swThHigh);
237 Spline<Scalar> sp(swThHigh, 1.0, // x0, x1
238 yTh, 0, // y0, y1
239 mTh, m1); // m0, m1
240 return sp.eval(ste);
241 }
242 else {
243 // straight line for swe > 1.0
244 return m1*(ste - 1.0) + 0.0;
245 }
246 }
247
248 // if the effective saturation is in an 'reasonable'
249 // range, we use the real van genuchten law...
250 return ParkerVanGen3P::pcgn(params, ste);
251 }
252
258 static Scalar pcAlpha(const Params &params, Scalar sne)
259 {
260 return ParkerVanGen3P::pcAlpha(params, sne);
261 }
262
269 static Scalar sw(const Params &params, Scalar pc)
270 {
271 return ParkerVanGen3P::sw(params, pc);
272 }
273
280 static Scalar dpc_dswe(const Params &params, Scalar swe)
281 {
282 return ParkerVanGen3P::dpc_dswe(params, swe);
283 }
284
291 static Scalar dswe_dpc(const Params &params, Scalar pc)
292 {
293 return ParkerVanGen3P::dswe_dpc(params, pc);
294 }
295
308 static Scalar krw(const Params &params, const Scalar swe)
309 {
310 //use regularization
311 if(swe > 1.0) return 1.0;
312 if(swe < 0.0) return 0.0;
313
314 //or use actual material law
315 return ParkerVanGen3P::krw(params, swe);
316 }
317
334 static Scalar krn(const Params &params, Scalar swe, Scalar sn, Scalar ste)
335 {
336 using std::max;
337 using std::min;
338 swe = max(min(swe, 1.0), 0.0);
339 ste = max(min(ste, 1.0), 0.0);
340
341 if(ste - swe <= 0.0)
342 return 0.0;
343
344 //or use actual material law
345 return ParkerVanGen3P::krn(params, swe, sn, ste);
346 }
347
348
361 static Scalar krg(const Params &params, const Scalar ste)
362 {
363 //return 0 if there is no gas
364 if(ste > 1.0)
365 return 0.0;
366
367 // use linear regularization for very high gas saturations
368 // to avoid a kink in the curve and to maintain a slope for
369 // the Newton solver
370 const Scalar threshold = 1e-3;
371 if(ste <= threshold)
372 {
373 const Scalar mSwr = ParkerVanGen3P::dkrg_dste(params, threshold);
374 const Scalar ySwr = ParkerVanGen3P::krg(params, threshold);
375 return ySwr + mSwr*(ste - threshold);
376 }
377
378 // For very low gas saturations:
379 // We use a scaling factor that decreases the gas phase permeability quite fast a very low gas phase
380 // saturations, thus making that phase virtually immobile.
381 // This prevents numerical issues related to the degeneration of the gas phase mass balance for the 3p3c model
382 // at very low gas phase saturations.
383
384 //get the absolute gas phase saturation
385 const Scalar st = ste*(1 - params.swr()) + params.swr();
386 const Scalar sg = 1.0 - st;
387 using std::max;
388 const Scalar scalFact = (sg > 0.1) ? 1.0 : max(0.0,
389 (sg - params.sgr())/(0.1 - params.sgr()));
390
391 //or use actual material law
392 return scalFact * ParkerVanGen3P::krg(params, ste);
393 }
394
403 static Scalar kr(const Params &params, const int phaseIdx, const Scalar swe, const Scalar sn, const Scalar ste)
404 {
405 switch (phaseIdx)
406 {
407 case 0:
408 return krw(params, swe);
409 case 1:
410 return krn(params, swe, sn, ste);
411 case 2:
412 return krg(params, ste);
413 }
414 DUNE_THROW(Dune::InvalidStateException,
415 "Invalid phase index ");
416 }
417
423 {
425 }
426
427};
428} // end namespace Dumux
429
430#endif
Provides 3rd order polynomial splines.
Parameters that are necessary for the regularization of the Parker - Van Genuchten capillary pressure...
Implementation of van Genuchten's capillary pressure-saturation relation for three phases.
Definition: adapt.hh:29
A 3rd order polynomial spline.
Definition: spline.hh:55
Scalar eval(Scalar x, bool extrapolate=false) const
Evaluate the spline at a given position.
Definition: splinecommon_.hh:142
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 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
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
Implementation of the regularized van Genuchten's capillary pressure <-> saturation relation....
Definition: regularizedparkervangen3p.hh:62
typename Params::Scalar Scalar
Definition: regularizedparkervangen3p.hh:67
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: regularizedparkervangen3p.hh:361
static Scalar pcAlpha(const Params &params, Scalar sne)
This function ensures a continuous transition from 2 to 3 phases and vice versa.
Definition: regularizedparkervangen3p.hh:258
static Scalar pcgw(const Params &params, Scalar swe)
The capillary pressure-saturation curve for the gas and wetting phase.
Definition: regularizedparkervangen3p.hh:93
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: regularizedparkervangen3p.hh:403
static Scalar dpc_dswe(const Params &params, Scalar swe)
Returns the partial derivative of the capillary pressure to the effective saturation.
Definition: regularizedparkervangen3p.hh:280
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: regularizedparkervangen3p.hh:308
static Scalar krn(const Params &params, Scalar swe, Scalar sn, Scalar ste)
The relative permeability for the non-wetting phase after the Model of Parker et al....
Definition: regularizedparkervangen3p.hh:334
ParamsT Params
Definition: regularizedparkervangen3p.hh:66
static Scalar pcgn(const Params &params, Scalar ste)
The capillary pressure-saturation curve for the gas and non-wetting phase.
Definition: regularizedparkervangen3p.hh:203
static Scalar sw(const Params &params, Scalar pc)
The saturation-capillary pressure curve.
Definition: regularizedparkervangen3p.hh:269
static Scalar bulkDensTimesAdsorpCoeff(const Params &params)
the basis for calculating adsorbed NAPL in storage term
Definition: regularizedparkervangen3p.hh:422
static Scalar pcnw(const Params &params, Scalar swe)
The capillary pressure-saturation curve for the non-wettigng and wetting phase.
Definition: regularizedparkervangen3p.hh:148
static Scalar pc(const Params &params, Scalar sw)
A regularized Parker- van Genuchten capillary pressure-saturation curve.
Definition: regularizedparkervangen3p.hh:82
static Scalar dswe_dpc(const Params &params, Scalar pc)
Returns the partial derivative of the effective saturation to the capillary pressure.
Definition: regularizedparkervangen3p.hh:291