24#ifndef PARKERVANGEN_3P_HH
25#define PARKERVANGEN_3P_HH
42template <
class ScalarT,
class ParamsT = ParkerVanGen3PParams<ScalarT> >
48 using Scalar =
typename Params::Scalar;
57 DUNE_THROW(Dune::NotImplemented,
"Capillary pressures for three phases is not so simple! Use pcgn, pcnw, and pcgw");
67 assert(0 <= swe && swe <= 1);
68 return pc_(params, swe);
78 assert(0 <= swe && swe <= 1);
79 return pc_(params, swe)/params.betaNw();
89 assert(0 <= ste && ste <= 1);
90 return pc_(params,ste)/params.betaGn();
103 if (sne<=0.001) sne=0.0;
104 if (sne>=1.0) sne=1.0;
106 if (sne>params.snr())
110 if (params.snr()>=0.001)
111 alpha = sne/params.snr();
126 DUNE_THROW(Dune::NotImplemented,
"sw(pc) for three phases not implemented! Do it yourself!");
137 DUNE_THROW(Dune::NotImplemented,
"dpc/dswe for three phases not implemented! Do it yourself!");
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();
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();
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();
190 DUNE_THROW(Dune::NotImplemented,
"dswe/dpc for three phases not implemented! Do it yourself!");
209 const Scalar r = 1.0 - pow(1 - pow(swe, 1/params.vgm()), params.vgm());
210 return sqrt(swe)*r*r;
234 krn = pow(1 - pow(swe, 1/params.vgm()), params.vgm());
235 krn -= pow(1 - pow(ste, 1/params.vgm()), params.vgm());
241 if (params.krRegardsSnr())
244 Scalar resIncluded = max(min((sn - params.snr()/ (1-params.swr())), 1.0), 0.0);
245 krn *= sqrt(resIncluded );
248 krn *= sqrt(sn / (1 - params.swr()));
267 assert(0 <= ste && ste <= 1);
270 return cbrt(1 - ste) * pow(1 - pow(ste, 1/params.vgm()), 2*params.vgm());
286 assert(0 < ste && ste <= 1);
289 const Scalar x = pow(ste, 1.0/params.vgm());
291 -pow(1.0 - x, 2*params.vgm())
292 *pow(1.0 - ste, -2.0/3)
309 return krw(params, swe);
311 return krn(params, swe, sn, ste);
313 return krg(params, ste);
315 DUNE_THROW(Dune::InvalidStateException,
316 "Invalid phase index ");
325 return params.rhoBulk() * params.KdNAPL();
339 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:44
static Scalar pc(const Params ¶ms, const Scalar sw)
The capillary pressure-saturation curve.
Definition: parkervangen3p.hh:55
static Scalar krg(const Params ¶ms, 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 ¶ms, 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 ¶ms, const Scalar pc)
The saturation-capillary pressure curve.
Definition: parkervangen3p.hh:124
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:174
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:304
static Scalar krn(const Params ¶ms, 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 ¶ms, const Scalar pc)
Returns the partial derivative of the effective saturation to the capillary pressure.
Definition: parkervangen3p.hh:188
static Scalar bulkDensTimesAdsorpCoeff(const Params ¶ms)
the basis for calculating adsorbed NAPL in storage term
Definition: parkervangen3p.hh:323
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:205
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:135
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:160
static Scalar pcgn(const Params ¶ms, const Scalar ste)
The capillary pressure-saturation curve for the gas and non-wetting phase.
Definition: parkervangen3p.hh:87
static Scalar pcgw(const Params ¶ms, 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 ¶ms, 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 ¶ms, 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 ¶ms, 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