3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
pcmax.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 *****************************************************************************/
26#ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_PC_MAX
27#define DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_PC_MAX
28
29#include <cmath>
30#include <dune/common/exceptions.hh>
31#include <dune/common/float_cmp.hh>
33
34namespace Dumux::FluidMatrix {
35
42{
43public:
44
45 template<class Scalar>
46 struct Params
47 {
48 Params(Scalar pcMax = 0, Scalar a1 = 0, Scalar a2 = 0, Scalar a3 = 0)
49 : pcMax_(pcMax), a1_(a1), a2_(a2), a3_(a3) {}
50
51 Scalar pcMax() const { return pcMax_; }
52 void setPcMax(Scalar pcMax) { pcMax_ = pcMax; }
53
54 Scalar a1() const { return a1_; }
55 void setA1(Scalar a1) { a1_ = a1; }
56
57 Scalar a2() const { return a2_; }
58 void setA2(Scalar a2) { a2_ = a2; }
59
60 Scalar a3() const { return a3_; }
61 void setA3(Scalar a3) { a3_ = a3; }
62
63 bool operator== (const Params& p) const
64 {
65 return Dune::FloatCmp::eq(pcMax(), p.pcMax(), 1e-6)
66 && Dune::FloatCmp::eq(a1(), p.a1(), 1e-6)
67 && Dune::FloatCmp::eq(a2(), p.a2(), 1e-6)
68 && Dune::FloatCmp::eq(a3(), p.a3(), 1e-6);
69 }
70
71 private:
72 Scalar pcMax_, a1_, a2_, a3_;
73 };
74
79 template<class Scalar = double>
80 static Params<Scalar> makeParams(const std::string& paramGroup)
81 {
82 const auto pcMax = getParamFromGroup<Scalar>(paramGroup, "PcMax");
83 const auto a1 = getParamFromGroup<Scalar>(paramGroup, "A1");
84 const auto a2 = getParamFromGroup<Scalar>(paramGroup, "A2");
85 const auto a3 = getParamFromGroup<Scalar>(paramGroup, "A3");
86 return {pcMax, a1, a2, a3};
87 }
88
100 template<class Scalar>
101 static Scalar area(const Scalar swe, const Scalar pc, const Params<Scalar>& params)
102 {
103 const Scalar a1 = params.a1();
104 const Scalar a2 = params.a2();
105 const Scalar a3 = params.a3();
106 const Scalar pcMax = params.pcMax();
107 return a1 * (pcMax-pc) * (1.0-swe) + a2*(pcMax-pc)*(pcMax-pc) * (1.-swe) + a3 * (pcMax-pc)*(1-swe)*(1-swe);
108 }
109};
110
111} // namespace Dumux::FluidMatrix
112
113#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: pcmax.hh:42
static Scalar area(const Scalar swe, const Scalar pc, const Params< Scalar > &params)
The interfacial area.
Definition: pcmax.hh:101
static Params< Scalar > makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: pcmax.hh:80
Scalar a2() const
Definition: pcmax.hh:57
void setA1(Scalar a1)
Definition: pcmax.hh:55
Params(Scalar pcMax=0, Scalar a1=0, Scalar a2=0, Scalar a3=0)
Definition: pcmax.hh:48
void setA3(Scalar a3)
Definition: pcmax.hh:61
Scalar a3() const
Definition: pcmax.hh:60
bool operator==(const Params &p) const
Definition: pcmax.hh:63
Scalar a1() const
Definition: pcmax.hh:54
void setA2(Scalar a2)
Definition: pcmax.hh:58
void setPcMax(Scalar pcMax)
Definition: pcmax.hh:52
Scalar pcMax() const
Definition: pcmax.hh:51