25#ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_EXPONENTIAL
26#define DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_EXPONENTIAL
29#include <dune/common/exceptions.hh>
30#include <dune/common/float_cmp.hh>
45 template<
class Scalar>
51 Scalar
swr()
const {
return swr_; }
54 Scalar
a1()
const {
return a1_; }
57 Scalar
a2()
const {
return a2_; }
60 Scalar
a3()
const {
return a3_; }
65 return Dune::FloatCmp::eq(
swr(), p.
swr(), 1e-6)
66 && Dune::FloatCmp::eq(
a1(), p.
a1(), 1e-6)
67 && Dune::FloatCmp::eq(
a2(), p.
a2(), 1e-6)
68 && Dune::FloatCmp::eq(
a3(), p.
a3(), 1e-6);
72 Scalar swr_, a1_, a2_, a3_;
79 template<
class Scalar =
double>
82 const auto swr = getParamFromGroup<Scalar>(paramGroup,
"Swr", 0.0);
83 const auto a1 = getParamFromGroup<Scalar>(paramGroup,
"A1");
84 const auto a2 = getParamFromGroup<Scalar>(paramGroup,
"A2");
85 const auto a3 = getParamFromGroup<Scalar>(paramGroup,
"A3");
86 return {swr, a1, a2, a3};
101 template<
class Scalar>
104 const Scalar a1 = params.
a1();
105 const Scalar a2 = params.
a2();
106 const Scalar a3 = params.
a3();
107 const Scalar swr = params.
swr();
109 return a1 * (swr-swe) * (1-swe) + a2 * (swr-swe) * (1-swe) * exp( a3 * pc) ;
118 template<
class Scalar>
121 const Scalar a2 = params.
a2();
122 const Scalar a3 = params.
a3();
123 const Scalar swr = params.
swr();
125 return a2 * a3 * (swr-swe) * (1-swe) * exp(a3*pc);
134 template<
class Scalar>
137 Scalar a1 = params.
a1();
138 Scalar a2 = params.
a2();
139 Scalar a3 = params.
a3();
140 Scalar swr = params.
swr();
142 return - a1 *(1+swr-2*swe) - a2 * exp(a3*pc) * (1+swr-2*swe);
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Definition: brookscorey.hh:35
Implementation of the exponential function relating specific interfacial area to wetting phase satura...
Definition: exponential.hh:42
static Scalar darea_dsw(const Scalar swe, const Scalar pc, const Params< Scalar > ¶ms)
the derivative of specific interfacial area function w.r.t. saturation
Definition: exponential.hh:135
static Params< Scalar > makeParams(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: exponential.hh:80
static Scalar darea_dpc(const Scalar swe, const Scalar pc, const Params< Scalar > ¶ms)
the derivative of specific interfacial area function w.r.t. capillary pressure
Definition: exponential.hh:119
static Scalar area(const Scalar swe, const Scalar pc, const Params< Scalar > ¶ms)
The interfacial area.
Definition: exponential.hh:102
Definition: exponential.hh:47
void setSwr(Scalar swr)
Definition: exponential.hh:52
Scalar a3() const
Definition: exponential.hh:60
Params(Scalar swr=0, Scalar a1=0, Scalar a2=0, Scalar a3=0)
Definition: exponential.hh:48
Scalar a2() const
Definition: exponential.hh:57
void setA1(Scalar a1)
Definition: exponential.hh:55
void setA2(Scalar a2)
Definition: exponential.hh:58
Scalar a1() const
Definition: exponential.hh:54
Scalar swr() const
Definition: exponential.hh:51
bool operator==(const Params &p) const
Definition: exponential.hh:63
void setA3(Scalar a3)
Definition: exponential.hh:61