23#ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_EXPONENTIAL
24#define DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_EXPONENTIAL
27#include <dune/common/exceptions.hh>
28#include <dune/common/float_cmp.hh>
43 template<
class Scalar>
49 Scalar
swr()
const {
return swr_; }
52 Scalar
a1()
const {
return a1_; }
55 Scalar
a2()
const {
return a2_; }
58 Scalar
a3()
const {
return a3_; }
63 return Dune::FloatCmp::eq(
swr(), p.
swr(), 1e-6)
64 && Dune::FloatCmp::eq(
a1(), p.
a1(), 1e-6)
65 && Dune::FloatCmp::eq(
a2(), p.
a2(), 1e-6)
66 && Dune::FloatCmp::eq(
a3(), p.
a3(), 1e-6);
70 Scalar swr_, a1_, a2_, a3_;
77 template<
class Scalar =
double>
80 const auto swr = getParamFromGroup<Scalar>(paramGroup,
"Swr", 0.0);
81 const auto a1 = getParamFromGroup<Scalar>(paramGroup,
"A1");
82 const auto a2 = getParamFromGroup<Scalar>(paramGroup,
"A2");
83 const auto a3 = getParamFromGroup<Scalar>(paramGroup,
"A3");
84 return {swr, a1, a2, a3};
99 template<
class Scalar>
102 const Scalar a1 = params.
a1();
103 const Scalar a2 = params.
a2();
104 const Scalar a3 = params.
a3();
105 const Scalar swr = params.
swr();
107 return a1 * (swr-swe) * (1-swe) + a2 * (swr-swe) * (1-swe) * exp( a3 * pc) ;
116 template<
class Scalar>
119 const Scalar a2 = params.
a2();
120 const Scalar a3 = params.
a3();
121 const Scalar swr = params.
swr();
123 return a2 * a3 * (swr-swe) * (1-swe) * exp(a3*pc);
132 template<
class Scalar>
135 Scalar a1 = params.
a1();
136 Scalar a2 = params.
a2();
137 Scalar a3 = params.
a3();
138 Scalar swr = params.
swr();
140 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:286
Implementation of the exponential function relating specific interfacial area to wetting phase satura...
Definition: exponential.hh:40
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:133
static Params< Scalar > makeParams(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: exponential.hh:78
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:117
static Scalar area(const Scalar swe, const Scalar pc, const Params< Scalar > ¶ms)
The interfacial area.
Definition: exponential.hh:100
Definition: exponential.hh:45
void setSwr(Scalar swr)
Definition: exponential.hh:50
Scalar a3() const
Definition: exponential.hh:58
Params(Scalar swr=0, Scalar a1=0, Scalar a2=0, Scalar a3=0)
Definition: exponential.hh:46
Scalar a2() const
Definition: exponential.hh:55
void setA1(Scalar a1)
Definition: exponential.hh:53
void setA2(Scalar a2)
Definition: exponential.hh:56
Scalar a1() const
Definition: exponential.hh:52
Scalar swr() const
Definition: exponential.hh:49
bool operator==(const Params &p) const
Definition: exponential.hh:61
void setA3(Scalar a3)
Definition: exponential.hh:59