26#ifndef DUMUX_EFF_TO_ABS_LAW_HH
27#define DUMUX_EFF_TO_ABS_LAW_HH
29#include <dune/common/exceptions.hh>
59template <
class EffLawT,
class AbsParamsT = EffToAbsLawParams<
typename EffLawT::Params> >
62 using EffLaw = EffLawT;
66 using Scalar =
typename EffLaw::Scalar;
80 return EffLaw::pc(params,
swToSwe(params,
sw));
90 return EffLaw::pcgw(params,
swToSwe(params,
sw));
100 return EffLaw::pcnw(params,
swToSwe(params,
sw));
110 return EffLaw::pcgn(params,
stToSte(params, st));
120 return EffLaw::pcAlpha(params, sn);
136 return EffLaw::sw(params,
pc);
156 return EffLaw::dpc_dswe(params,
pc);
178 return EffLaw::dsw_dpc(params,
pc);
197 return EffLaw::krw(params,
swToSwe(params,
sw));
235 return EffLaw::krg(params,
stToSte(params, st));
249 return EffLaw::kr(params, phaseIdx,
swToSwe(params,
sw), sn,
stToSte(params, st));
258 return EffLaw::bulkDensTimesAdsorpCoeff(params);
272 return (
sw-params.swr())/(1.-params.swr());
300 return (st-params.swr()) / (1-params.swr());
314 DUNE_THROW(Dune::NotImplemented,
"sgTosge for three phases not implemented!");
329 DUNE_THROW(Dune::NotImplemented,
"sweTosw for three phases not implemented!");
334 DUNE_THROW(Dune::NotImplemented,
"sneTosn for three phases not implemented!");
339 DUNE_THROW(Dune::NotImplemented,
"sgeTosg for three phases not implemented!");
351 DUNE_THROW(Dune::NotImplemented,
"dswe_dsw for three phases not implemented!");
364 DUNE_THROW(Dune::NotImplemented,
"dsw_dswe for three phases not implemented!");
static Scalar pcgn(const Params ¶ms, const Scalar st)
The capillary pressure-saturation curve for the gas and non-wetting phase.
Definition: 3p/efftoabslaw.hh:108
static Scalar pcgw(const Params ¶ms, const Scalar sw)
The capillary pressure-saturation curve for the gas and wetting phase.
Definition: 3p/efftoabslaw.hh:88
static Scalar krw(const Params ¶ms, Scalar sw, const Scalar sn)
The relative permeability for the wetting phase.
Definition: 3p/efftoabslaw.hh:195
static Scalar dpc_dsw(const Params ¶ms, const Scalar sw)
Returns the partial derivative of the capillary pressure w.r.t the absolute saturation....
Definition: 3p/efftoabslaw.hh:154
static Scalar pcAlpha(const Params ¶ms, const Scalar sn)
This function ensures a continuous transition from 2 to 3 phases and vice versa.
Definition: 3p/efftoabslaw.hh:118
static Scalar dswe_dsw_(const Params ¶ms)
Derivative of the effective saturation w.r.t. the absolute saturation.
Definition: 3p/efftoabslaw.hh:349
static Scalar dsw_dpc(const Params ¶ms, const Scalar pc)
Returns the partial derivative of the absolute saturation w.r.t. the capillary pressure....
Definition: 3p/efftoabslaw.hh:176
static Scalar sw(const Params ¶ms, Scalar pc)
The saturation-capillary pressure curve.
Definition: 2p/efftoabslaw.hh:92
static Scalar snToSne(const Params ¶ms, const Scalar sn)
Convert an absolute non-wetting saturation to an effective one.
Definition: 3p/efftoabslaw.hh:284
static Scalar pcnw(const Params ¶ms, const Scalar sw)
The capillary pressure-saturation curve the non-wetting and wetting phase.
Definition: 3p/efftoabslaw.hh:98
static Scalar krg(const Params ¶ms, const Scalar sw, const Scalar sn)
The relative permeability for the gas phase.
Definition: 3p/efftoabslaw.hh:232
static Scalar krn(const Params ¶ms, const Scalar sw, const Scalar sn)
The relative permeability for the non-wetting phase.
Definition: 3p/efftoabslaw.hh:213
static Scalar sgToSge(const Params ¶ms, Scalar sg)
Convert an absolute gas saturation to an effective one.
Definition: 3p/efftoabslaw.hh:312
static Scalar sneToSn_(const Params ¶ms, Scalar sne)
Definition: 3p/efftoabslaw.hh:332
AbsParamsT Params
Definition: 2p/efftoabslaw.hh:64
static Scalar swToSwe(const Params ¶ms, Scalar sw)
Convert an absolute wetting saturation to an effective one.
Definition: 2p/efftoabslaw.hh:218
static Scalar sw(const Params ¶ms, const Scalar pc)
The saturation-capillary pressure curve.
Definition: 3p/efftoabslaw.hh:134
static Scalar bulkDensTimesAdsorpCoeff(const Params ¶ms)
the basis for calculating adsorbed NAPL in storage term
Definition: 3p/efftoabslaw.hh:256
static Scalar dsw_dswe_(const Params ¶ms)
Derivative of the absolute saturation w.r.t. the effective saturation.
Definition: 3p/efftoabslaw.hh:362
static Scalar swToSwe(const Params ¶ms, const Scalar sw)
Convert an absolute wetting saturation to an effective one.
Definition: 3p/efftoabslaw.hh:270
static Scalar sgeToSg_(const Params ¶ms, Scalar sge)
Definition: 3p/efftoabslaw.hh:337
static Scalar sweToSw_(const Params ¶ms, Scalar swe)
Convert an effective wetting saturation to an absolute one.
Definition: 3p/efftoabslaw.hh:327
static Scalar kr(const Params ¶ms, const int phaseIdx, const Scalar sw, const Scalar sn, const Scalar sg)
The relative permeability for a phase.
Definition: 3p/efftoabslaw.hh:246
static Scalar pc(const Params ¶ms, Scalar sw)
The capillary pressure-saturation curve.
Definition: 2p/efftoabslaw.hh:77
typename EffLaw::Scalar Scalar
Definition: 2p/efftoabslaw.hh:65
static Scalar stToSte(const Params ¶ms, const Scalar st)
Convert an absolute total liquid saturation to an effective one.
Definition: 3p/efftoabslaw.hh:298
static Scalar pc(const Params ¶ms, const Scalar sw)
The capillary pressure-saturation curve.
Definition: 3p/efftoabslaw.hh:78
A default implementation of the parameters for the adapter class to convert material laws from effect...