3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
exponential.hh
Go to the documentation of this file.
1/*****************************************************************************
2 * See the file COPYING for full copying permissions. *
3 * *
4 * This program is free software: you can redistribute it and/or modify *
5 * it under the terms of the GNU General Public License as published by *
6 * the Free Software Foundation, either version 3 of the License, or *
7 * (at your option) any later version. *
8 * *
9 * This program is distributed in the hope that it will be useful, *
10 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
11 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
12 * GNU General Public License for more details. *
13 * *
14 * You should have received a copy of the GNU General Public License *
15 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
16 *****************************************************************************/
23#ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_EXPONENTIAL
24#define DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_EXPONENTIAL
25
26#include <cmath>
27#include <dune/common/exceptions.hh>
28#include <dune/common/float_cmp.hh>
30
31namespace Dumux::FluidMatrix {
32
40{
41public:
42
43 template<class Scalar>
44 struct Params
45 {
46 Params(Scalar swr = 0, Scalar a1 = 0, Scalar a2 = 0, Scalar a3 = 0)
47 : swr_(swr), a1_(a1), a2_(a2), a3_(a3) {}
48
49 Scalar swr() const { return swr_; }
50 void setSwr(Scalar swr) { swr_ = swr; }
51
52 Scalar a1() const { return a1_; }
53 void setA1(Scalar a1) { a1_ = a1; }
54
55 Scalar a2() const { return a2_; }
56 void setA2(Scalar a2) { a2_ = a2; }
57
58 Scalar a3() const { return a3_; }
59 void setA3(Scalar a3) { a3_ = a3; }
60
61 bool operator== (const Params& p) const
62 {
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);
67 }
68
69 private:
70 Scalar swr_, a1_, a2_, a3_;
71 };
72
77 template<class Scalar = double>
78 static Params<Scalar> makeParams(const std::string& paramGroup)
79 {
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};
85 }
86
87
99 template<class Scalar>
100 static Scalar area(const Scalar swe, const Scalar pc, const Params<Scalar>& params)
101 {
102 const Scalar a1 = params.a1();
103 const Scalar a2 = params.a2();
104 const Scalar a3 = params.a3();
105 const Scalar swr = params.swr();
106 using std::exp;
107 return a1 * (swr-swe) * (1-swe) + a2 * (swr-swe) * (1-swe) * exp( a3 * pc) ;
108 }
109
116 template<class Scalar>
117 static Scalar darea_dpc(const Scalar swe, const Scalar pc, const Params<Scalar>& params)
118 {
119 const Scalar a2 = params.a2();
120 const Scalar a3 = params.a3();
121 const Scalar swr = params.swr();
122 using std::exp;
123 return a2 * a3 * (swr-swe) * (1-swe) * exp(a3*pc);
124 }
125
132 template<class Scalar>
133 static Scalar darea_dsw(const Scalar swe, const Scalar pc, const Params<Scalar>& params)
134 {
135 Scalar a1 = params.a1();
136 Scalar a2 = params.a2();
137 Scalar a3 = params.a3();
138 Scalar swr = params.swr();
139 using std::exp;
140 return - a1 *(1+swr-2*swe) - a2 * exp(a3*pc) * (1+swr-2*swe);
141 }
142
143};
144
145} // namespace Dumux::FluidMatrix
146
147#endif
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:40
static Scalar darea_dsw(const Scalar swe, const Scalar pc, const Params< Scalar > &params)
the derivative of specific interfacial area function w.r.t. saturation
Definition: exponential.hh:133
static Params< Scalar > makeParams(const std::string &paramGroup)
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 > &params)
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 > &params)
The interfacial area.
Definition: exponential.hh:100
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