24#ifndef PARKERVANGEN_3P_HH
25#define PARKERVANGEN_3P_HH
31#warning "This header is deprecated. Removal after 3.3. Use new material laws."
44template <
class ScalarT,
class ParamsT = ParkerVanGen3PParams<ScalarT> >
50 using Scalar =
typename Params::Scalar;
59 DUNE_THROW(Dune::NotImplemented,
"Capillary pressures for three phases is not so simple! Use pcgn, pcnw, and pcgw");
69 assert(0 <= swe && swe <= 1);
70 return pc_(params, swe);
80 assert(0 <= swe && swe <= 1);
81 return pc_(params, swe)/params.betaNw();
91 assert(0 <= ste && ste <= 1);
92 return pc_(params,ste)/params.betaGn();
105 if (sne<=0.001) sne=0.0;
106 if (sne>=1.0) sne=1.0;
108 if (sne>params.snr())
112 if (params.snr()>=0.001)
113 alpha = sne/params.snr();
128 DUNE_THROW(Dune::NotImplemented,
"sw(pc) for three phases not implemented! Do it yourself!");
139 DUNE_THROW(Dune::NotImplemented,
"dpc/dswe for three phases not implemented! Do it yourself!");
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();
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();
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();
192 DUNE_THROW(Dune::NotImplemented,
"dswe/dpc for three phases not implemented! Do it yourself!");
211 const Scalar r = 1.0 - pow(1 - pow(swe, 1/params.vgm()), params.vgm());
212 return sqrt(swe)*r*r;
236 krn = pow(1 - pow(swe, 1/params.vgm()), params.vgm());
237 krn -= pow(1 - pow(ste, 1/params.vgm()), params.vgm());
242 if (params.krRegardsSnr())
245 const Scalar resIncluded = clamp(sn - params.snr()/ (1-params.swr()), 0.0, 1.0);
246 krn *= sqrt(resIncluded );
249 krn *= sqrt(sn / (1 - params.swr()));
268 assert(0 <= ste && ste <= 1);
271 return cbrt(1 - ste) * pow(1 - pow(ste, 1/params.vgm()), 2*params.vgm());
287 assert(0 < ste && ste <= 1);
290 const Scalar x = pow(ste, 1.0/params.vgm());
292 -pow(1.0 - x, 2*params.vgm())
293 *pow(1.0 - ste, -2.0/3)
310 return krw(params, swe);
312 return krn(params, swe, sn, ste);
314 return krg(params, ste);
316 DUNE_THROW(Dune::InvalidStateException,
317 "Invalid phase index ");
326 return params.rhoBulk() * params.KdNAPL();
340 return pow(pow(se, -1/params.vgm()) - 1, 1/params.vgn())/params.vgAlpha();
Specification of the material params for the van Genuchten capillary pressure model.
Implementation of van Genuchten's capillary pressure <-> saturation relation. This class bundles the ...
Definition: parkervangen3p.hh:46
static Scalar pc(const Params ¶ms, const Scalar sw)
The capillary pressure-saturation curve.
Definition: parkervangen3p.hh:57
static Scalar krg(const Params ¶ms, 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 ¶ms, 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 ¶ms, const Scalar pc)
The saturation-capillary pressure curve.
Definition: parkervangen3p.hh:126
static Scalar dpcgn_dste(const Params ¶ms, const Scalar seRegu)
Returns the partial derivative of the capillary pressure to the effective saturation.
Definition: parkervangen3p.hh:176
static Scalar kr(const Params ¶ms, 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 ¶ms, 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 ¶ms, const Scalar pc)
Returns the partial derivative of the effective saturation to the capillary pressure.
Definition: parkervangen3p.hh:190
static Scalar bulkDensTimesAdsorpCoeff(const Params ¶ms)
the basis for calculating adsorbed NAPL in storage term
Definition: parkervangen3p.hh:324
static Scalar krw(const Params ¶ms, 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 ¶ms, 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 ¶ms, const Scalar seRegu)
Returns the partial derivative of the capillary pressure to the effective saturation.
Definition: parkervangen3p.hh:162
static Scalar pcgn(const Params ¶ms, const Scalar ste)
The capillary pressure-saturation curve for the gas and nonwetting phase.
Definition: parkervangen3p.hh:89
static Scalar pcgw(const Params ¶ms, 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 ¶ms, 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 ¶ms, 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 ¶ms, 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