25#ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_EXPONENTIAL_CUBIC
26#define DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_EXPONENTIAL_CUBIC
29#include <dune/common/exceptions.hh>
30#include <dune/common/float_cmp.hh>
44 template<
class Scalar>
50 Scalar
swr()
const {
return swr_; }
53 Scalar
a1()
const {
return a1_; }
56 Scalar
a2()
const {
return a2_; }
59 Scalar
a3()
const {
return a3_; }
64 return Dune::FloatCmp::eq(
swr(), p.
swr(), 1e-6)
65 && Dune::FloatCmp::eq(
a1(), p.
a1(), 1e-6)
66 && Dune::FloatCmp::eq(
a2(), p.
a2(), 1e-6)
67 && Dune::FloatCmp::eq(
a3(), p.
a3(), 1e-6);
71 Scalar swr_, a1_, a2_, a3_;
78 template<
class Scalar =
double>
81 const auto swr = getParamFromGroup<Scalar>(paramGroup,
"Swr", 0.0);
82 const auto a1 = getParamFromGroup<Scalar>(paramGroup,
"A1");
83 const auto a2 = getParamFromGroup<Scalar>(paramGroup,
"A2");
84 const auto a3 = getParamFromGroup<Scalar>(paramGroup,
"A3");
85 return {swr, a1, a2, a3};
98 template<
class Scalar>
101 const Scalar a1 = params.
a1();
102 const Scalar a2 = params.
a2();
103 const Scalar a3 = params.
a3();
104 const Scalar swr = params.
swr();
105 const Scalar factor = (swr-sw)*(1.0-sw) ;
106 return a1*factor + a2*factor*pc + a3*factor*pc*pc;
116 template<
class Scalar>
119 const Scalar swr = params.
swr();
120 const Scalar a1 = params.
a1();
121 const Scalar a2 = params.
a2();
122 const Scalar a3 = params.
a3();
123 return a2*(swr-sw)*(1-sw) + a3*(swr-sw)*(1-sw)*pc;
133 template<
class Scalar>
136 const Scalar swr = params.
swr();
137 const Scalar a1 = params.
a1();
138 const Scalar a2 = params.
a2();
139 const Scalar a3 = params.
a3();
140 const Scalar derivativeFactor = (sw-1.0)+(sw-swr);
141 return a1 * derivativeFactor + a2 * derivativeFactor * pc + a3 * derivativeFactor * pc*pc ;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Definition: brookscorey.hh:35
Implementation of the polynomial of second order relating specific interfacial area to wetting phase ...
Definition: polynomialedgezero2ndorder.hh:41
static Params< Scalar > makeParams(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: polynomialedgezero2ndorder.hh:79
static Scalar darea_dpc(const Scalar sw, const Scalar pc, const Params< Scalar > ¶ms)
the derivative of specific interfacial area function w.r.t. capillary pressure
Definition: polynomialedgezero2ndorder.hh:117
static Scalar darea_dsw(const Scalar sw, const Scalar pc, const Params< Scalar > ¶ms)
the derivative of specific interfacial area function w.r.t. saturation
Definition: polynomialedgezero2ndorder.hh:134
static Scalar area(const Scalar sw, const Scalar pc, const Params< Scalar > ¶ms)
The interfacial area the suggested (as estimated from pore network models) awn surface: ].
Definition: polynomialedgezero2ndorder.hh:99
Definition: polynomialedgezero2ndorder.hh:46
void setA2(Scalar a2)
Definition: polynomialedgezero2ndorder.hh:57
Scalar a1() const
Definition: polynomialedgezero2ndorder.hh:53
void setSwr(Scalar swr)
Definition: polynomialedgezero2ndorder.hh:51
Scalar a3() const
Definition: polynomialedgezero2ndorder.hh:59
Params(Scalar swr=0, Scalar a1=0, Scalar a2=0, Scalar a3=0)
Definition: polynomialedgezero2ndorder.hh:47
bool operator==(const Params &p) const
Definition: polynomialedgezero2ndorder.hh:62
void setA3(Scalar a3)
Definition: polynomialedgezero2ndorder.hh:60
void setA1(Scalar a1)
Definition: polynomialedgezero2ndorder.hh:54
Scalar swr() const
Definition: polynomialedgezero2ndorder.hh:50
Scalar a2() const
Definition: polynomialedgezero2ndorder.hh:56