3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
polynomial2ndorder.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_POLYNOMIAL_SECOND_ORDER_HH
26#define DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_POLYNOMIAL_SECOND_ORDER_HH
27
28#include <cmath>
29#include <dune/common/exceptions.hh>
30#include <dune/common/float_cmp.hh>
32
33namespace Dumux::FluidMatrix {
34
41{
42public:
43
44 template<class Scalar>
45 struct Params
46 {
47 Params(Scalar a00 = 0, Scalar a01 = 0, Scalar a02 = 0, Scalar a10 = 0, Scalar a11 = 0, Scalar a20 = 0)
48 : a00_(a00), a01_(a01), a02_(a02), a10_(a10), a11_(a11), a20_(a20) {}
49
50 Scalar a00() const { return a00_; }
51 void setA00(Scalar a00) { a00_ = a00; }
52
53 Scalar a01() const { return a01_; }
54 void setA01(Scalar a01) { a01_ = a01; }
55
56 Scalar a02() const { return a02_; }
57 void setA02(Scalar a02) { a02_ = a02; }
58
59 Scalar a10() const { return a10_; }
60 void setA10(Scalar a10) { a10_ = a10; }
61
62 Scalar a11() const { return a11_; }
63 void setA11(Scalar a11) { a11_ = a11; }
64
65 Scalar a20() const { return a20_; }
66 void setA20(Scalar a20) { a20_ = a20; }
67
68 bool operator== (const Params& p) const
69 {
70 return Dune::FloatCmp::eq(a00(), p.a00(), 1e-6)
71 && Dune::FloatCmp::eq(a01(), p.a01(), 1e-6)
72 && Dune::FloatCmp::eq(a02(), p.a02(), 1e-6)
73 && Dune::FloatCmp::eq(a10(), p.a10(), 1e-6)
74 && Dune::FloatCmp::eq(a11(), p.a11(), 1e-6)
75 && Dune::FloatCmp::eq(a20(), p.a20(), 1e-6);
76 }
77
78 private:
79 Scalar a00_, a01_, a02_, a10_, a11_, a20_;
80 };
81
86 template<class Scalar = double>
87 static Params<Scalar> makeParams(const std::string& paramGroup)
88 {
89 const auto a00 = getParamFromGroup<Scalar>(paramGroup, "A00", 0.0);
90 const auto a01 = getParamFromGroup<Scalar>(paramGroup, "A01");
91 const auto a02 = getParamFromGroup<Scalar>(paramGroup, "A02");
92 const auto a10 = getParamFromGroup<Scalar>(paramGroup, "A10");
93 const auto a11 = getParamFromGroup<Scalar>(paramGroup, "A11");
94 const auto a20 = getParamFromGroup<Scalar>(paramGroup, "A20");
95 return {a00, a01, a02, a10, a11, a20};
96 }
97
109 template<class Scalar>
110 static Scalar area(const Scalar swe, const Scalar pc, const Params<Scalar>& params)
111 {
112 const Scalar a00 = params.a00();
113 const Scalar a10 = params.a10();
114 const Scalar a20 = params.a20();
115 const Scalar a11 = params.a11();
116 const Scalar a01 = params.a01();
117 const Scalar a02 = params.a02();
118
119 return a00 + a10 * swe + a20 * swe*swe + a11*swe*pc + a01*pc + a02*pc*pc;
120 }
121
129 template<class Scalar>
130 static Scalar darea_dpc(const Scalar swe, const Scalar pc, const Params<Scalar>& params)
131 {
132 return params.a11()*swe + params.a01() + 2.0*params.a02() * pc;
133 }
134
142 template<class Scalar>
143 static Scalar darea_dsw(const Scalar swe, const Scalar pc, const Params<Scalar>& params)
144 {
145 return params.a11()*pc + params.a10() + 2.0*params.a20()*swe;
146 }
147};
148
149} // end namespace Dumux::FluidMatrix
150
151#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: polynomial2ndorder.hh:41
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:130
static Params< Scalar > makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: polynomial2ndorder.hh:87
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:143
static Scalar area(const Scalar swe, const Scalar pc, const Params< Scalar > &params)
The interfacial area.
Definition: polynomial2ndorder.hh:110
Params(Scalar a00=0, Scalar a01=0, Scalar a02=0, Scalar a10=0, Scalar a11=0, Scalar a20=0)
Definition: polynomial2ndorder.hh:47
void setA20(Scalar a20)
Definition: polynomial2ndorder.hh:66
void setA02(Scalar a02)
Definition: polynomial2ndorder.hh:57
bool operator==(const Params &p) const
Definition: polynomial2ndorder.hh:68
Scalar a20() const
Definition: polynomial2ndorder.hh:65
Scalar a00() const
Definition: polynomial2ndorder.hh:50
Scalar a11() const
Definition: polynomial2ndorder.hh:62
void setA01(Scalar a01)
Definition: polynomial2ndorder.hh:54
Scalar a10() const
Definition: polynomial2ndorder.hh:59
Scalar a02() const
Definition: polynomial2ndorder.hh:56
void setA00(Scalar a00)
Definition: polynomial2ndorder.hh:51
void setA11(Scalar a11)
Definition: polynomial2ndorder.hh:63
Scalar a01() const
Definition: polynomial2ndorder.hh:53
void setA10(Scalar a10)
Definition: polynomial2ndorder.hh:60