3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
polynomial2ndorder.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_POLYNOMIAL_SECOND_ORDER_HH
24#define DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_POLYNOMIAL_SECOND_ORDER_HH
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 a00 = 0, Scalar a01 = 0, Scalar a02 = 0, Scalar a10 = 0, Scalar a11 = 0, Scalar a20 = 0)
46 : a00_(a00), a01_(a01), a02_(a02), a10_(a10), a11_(a11), a20_(a20) {}
47
48 Scalar a00() const { return a00_; }
49 void setA00(Scalar a00) { a00_ = a00; }
50
51 Scalar a01() const { return a01_; }
52 void setA01(Scalar a01) { a01_ = a01; }
53
54 Scalar a02() const { return a02_; }
55 void setA02(Scalar a02) { a02_ = a02; }
56
57 Scalar a10() const { return a10_; }
58 void setA10(Scalar a10) { a10_ = a10; }
59
60 Scalar a11() const { return a11_; }
61 void setA11(Scalar a11) { a11_ = a11; }
62
63 Scalar a20() const { return a20_; }
64 void setA20(Scalar a20) { a20_ = a20; }
65
66 bool operator== (const Params& p) const
67 {
68 return Dune::FloatCmp::eq(a00(), p.a00(), 1e-6)
69 && Dune::FloatCmp::eq(a01(), p.a01(), 1e-6)
70 && Dune::FloatCmp::eq(a02(), p.a02(), 1e-6)
71 && Dune::FloatCmp::eq(a10(), p.a10(), 1e-6)
72 && Dune::FloatCmp::eq(a11(), p.a11(), 1e-6)
73 && Dune::FloatCmp::eq(a20(), p.a20(), 1e-6);
74 }
75
76 private:
77 Scalar a00_, a01_, a02_, a10_, a11_, a20_;
78 };
79
84 template<class Scalar = double>
85 static Params<Scalar> makeParams(const std::string& paramGroup)
86 {
87 const auto a00 = getParamFromGroup<Scalar>(paramGroup, "A00", 0.0);
88 const auto a01 = getParamFromGroup<Scalar>(paramGroup, "A01");
89 const auto a02 = getParamFromGroup<Scalar>(paramGroup, "A02");
90 const auto a10 = getParamFromGroup<Scalar>(paramGroup, "A10");
91 const auto a11 = getParamFromGroup<Scalar>(paramGroup, "A11");
92 const auto a20 = getParamFromGroup<Scalar>(paramGroup, "A20");
93 return {a00, a01, a02, a10, a11, a20};
94 }
95
107 template<class Scalar>
108 static Scalar area(const Scalar swe, const Scalar pc, const Params<Scalar>& params)
109 {
110 const Scalar a00 = params.a00();
111 const Scalar a10 = params.a10();
112 const Scalar a20 = params.a20();
113 const Scalar a11 = params.a11();
114 const Scalar a01 = params.a01();
115 const Scalar a02 = params.a02();
116
117 return a00 + a10 * swe + a20 * swe*swe + a11*swe*pc + a01*pc + a02*pc*pc;
118 }
119
127 template<class Scalar>
128 static Scalar darea_dpc(const Scalar swe, const Scalar pc, const Params<Scalar>& params)
129 {
130 return params.a11()*swe + params.a01() + 2.0*params.a02() * pc;
131 }
132
140 template<class Scalar>
141 static Scalar darea_dsw(const Scalar swe, const Scalar pc, const Params<Scalar>& params)
142 {
143 return params.a11()*pc + params.a10() + 2.0*params.a20()*swe;
144 }
145};
146
147} // end namespace Dumux::FluidMatrix
148
149#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Definition: brookscorey.hh:286
Implementation of the polynomial of second order relating specific interfacial area to wetting phase ...
Definition: polynomial2ndorder.hh:39
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: polynomial2ndorder.hh:128
static Params< Scalar > makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: polynomial2ndorder.hh:85
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: polynomial2ndorder.hh:141
static Scalar area(const Scalar swe, const Scalar pc, const Params< Scalar > &params)
The interfacial area.
Definition: polynomial2ndorder.hh:108
Params(Scalar a00=0, Scalar a01=0, Scalar a02=0, Scalar a10=0, Scalar a11=0, Scalar a20=0)
Definition: polynomial2ndorder.hh:45
void setA20(Scalar a20)
Definition: polynomial2ndorder.hh:64
void setA02(Scalar a02)
Definition: polynomial2ndorder.hh:55
bool operator==(const Params &p) const
Definition: polynomial2ndorder.hh:66
Scalar a20() const
Definition: polynomial2ndorder.hh:63
Scalar a00() const
Definition: polynomial2ndorder.hh:48
Scalar a11() const
Definition: polynomial2ndorder.hh:60
void setA01(Scalar a01)
Definition: polynomial2ndorder.hh:52
Scalar a10() const
Definition: polynomial2ndorder.hh:57
Scalar a02() const
Definition: polynomial2ndorder.hh:54
void setA00(Scalar a00)
Definition: polynomial2ndorder.hh:49
void setA11(Scalar a11)
Definition: polynomial2ndorder.hh:61
Scalar a01() const
Definition: polynomial2ndorder.hh:51
void setA10(Scalar a10)
Definition: polynomial2ndorder.hh:58