3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
polynomialedgezero2ndorder.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_CUBIC
24#define DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_EXPONENTIAL_CUBIC
25
26#include <cmath>
27#include <dune/common/exceptions.hh>
28#include <dune/common/float_cmp.hh>
30
31namespace Dumux::FluidMatrix {
32
39{
40public:
41
42 template<class Scalar>
43 struct Params
44 {
45 Params(Scalar swr = 0, Scalar a1 = 0, Scalar a2 = 0, Scalar a3 = 0)
46 : swr_(swr), a1_(a1), a2_(a2), a3_(a3) {}
47
48 Scalar swr() const { return swr_; }
49 void setSwr(Scalar swr) { swr_ = swr; }
50
51 Scalar a1() const { return a1_; }
52 void setA1(Scalar a1) { a1_ = a1; }
53
54 Scalar a2() const { return a2_; }
55 void setA2(Scalar a2) { a2_ = a2; }
56
57 Scalar a3() const { return a3_; }
58 void setA3(Scalar a3) { a3_ = a3; }
59
60 bool operator== (const Params& p) const
61 {
62 return Dune::FloatCmp::eq(swr(), p.swr(), 1e-6)
63 && Dune::FloatCmp::eq(a1(), p.a1(), 1e-6)
64 && Dune::FloatCmp::eq(a2(), p.a2(), 1e-6)
65 && Dune::FloatCmp::eq(a3(), p.a3(), 1e-6);
66 }
67
68 private:
69 Scalar swr_, a1_, a2_, a3_;
70 };
71
76 template<class Scalar = double>
77 static Params<Scalar> makeParams(const std::string& paramGroup)
78 {
79 const auto swr = getParamFromGroup<Scalar>(paramGroup, "Swr", 0.0);
80 const auto a1 = getParamFromGroup<Scalar>(paramGroup, "A1");
81 const auto a2 = getParamFromGroup<Scalar>(paramGroup, "A2");
82 const auto a3 = getParamFromGroup<Scalar>(paramGroup, "A3");
83 return {swr, a1, a2, a3};
84 }
85
96 template<class Scalar>
97 static Scalar area(const Scalar sw, const Scalar pc, const Params<Scalar>& params)
98 {
99 const Scalar a1 = params.a1();
100 const Scalar a2 = params.a2();
101 const Scalar a3 = params.a3();
102 const Scalar swr = params.swr();
103 const Scalar factor = (swr-sw)*(1.0-sw) ;
104 return a1*factor + a2*factor*pc + a3*factor*pc*pc;
105 }
106
114 template<class Scalar>
115 static Scalar darea_dpc(const Scalar sw, const Scalar pc, const Params<Scalar>& params)
116 {
117 const Scalar swr = params.swr();
118 const Scalar a1 = params.a1();
119 const Scalar a2 = params.a2();
120 const Scalar a3 = params.a3();
121 return a2*(swr-sw)*(1-sw) + a3*(swr-sw)*(1-sw)*pc;
122 }
123
131 template<class Scalar>
132 static Scalar darea_dsw(const Scalar sw, const Scalar pc, const Params<Scalar>& params)
133 {
134 const Scalar swr = params.swr();
135 const Scalar a1 = params.a1();
136 const Scalar a2 = params.a2();
137 const Scalar a3 = params.a3();
138 const Scalar derivativeFactor = (sw-1.0)+(sw-swr);
139 return a1 * derivativeFactor + a2 * derivativeFactor * pc + a3 * derivativeFactor * pc*pc ;
140 }
141};
142
143} // namespace Dumux::FluidMatrix
144
145#endif
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:39
static Params< Scalar > makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: polynomialedgezero2ndorder.hh:77
static Scalar darea_dpc(const Scalar sw, const Scalar pc, const Params< Scalar > &params)
the derivative of specific interfacial area function w.r.t. capillary pressure
Definition: polynomialedgezero2ndorder.hh:115
static Scalar darea_dsw(const Scalar sw, const Scalar pc, const Params< Scalar > &params)
the derivative of specific interfacial area function w.r.t. saturation
Definition: polynomialedgezero2ndorder.hh:132
static Scalar area(const Scalar sw, const Scalar pc, const Params< Scalar > &params)
The interfacial area the suggested (as estimated from pore network models) awn surface: ].
Definition: polynomialedgezero2ndorder.hh:97
void setA2(Scalar a2)
Definition: polynomialedgezero2ndorder.hh:55
Scalar a1() const
Definition: polynomialedgezero2ndorder.hh:51
void setSwr(Scalar swr)
Definition: polynomialedgezero2ndorder.hh:49
Scalar a3() const
Definition: polynomialedgezero2ndorder.hh:57
Params(Scalar swr=0, Scalar a1=0, Scalar a2=0, Scalar a3=0)
Definition: polynomialedgezero2ndorder.hh:45
bool operator==(const Params &p) const
Definition: polynomialedgezero2ndorder.hh:60
void setA3(Scalar a3)
Definition: polynomialedgezero2ndorder.hh:58
void setA1(Scalar a1)
Definition: polynomialedgezero2ndorder.hh:52
Scalar swr() const
Definition: polynomialedgezero2ndorder.hh:48
Scalar a2() const
Definition: polynomialedgezero2ndorder.hh:54