3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
exponentialcubic.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 *****************************************************************************/
28#ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_EXPONENTIAL_CUBIC
29#define DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_EXPONENTIAL_CUBIC
30
31#include <cmath>
32#include <dune/common/exceptions.hh>
33#include <dune/common/float_cmp.hh>
35
36namespace Dumux::FluidMatrix {
37
44{
45public:
46
47 template<class Scalar>
48 struct Params
49 {
50 Params(Scalar a1 = 0, Scalar a2 = 0, Scalar a3 = 0)
51 : a1_(a1), a2_(a2), a3_(a3) {}
52
53 Scalar a1() const { return a1_; }
54 void setA1(Scalar a1) { a1_ = a1; }
55
56 Scalar a2() const { return a2_; }
57 void setA2(Scalar a2) { a2_ = a2; }
58
59 Scalar a3() const { return a3_; }
60 void setA3(Scalar a3) { a3_ = a3; }
61
62 bool operator== (const Params& p) const
63 {
64 return 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 a1_, a2_, a3_;
71 };
72
77 template<class Scalar = double>
78 static Params<Scalar> makeParams(const std::string& paramGroup)
79 {
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 {a1, a2, a3};
84 }
85
97 template<class Scalar>
98 static Scalar area(const Scalar swe, const Scalar pc, const Params<Scalar>& params)
99 {
100 // TODO think about awn surface for relative saturation
101 const Scalar a1 = params.a1();
102 const Scalar a2 = params.a2();
103 const Scalar a3 = params.a3();
104
105 using std::exp;
106 return a1 * exp( a2 * swe) + a3 * pc * pc * pc ;
107 }
108};
109
110} // namespace Dumux::FluidMatrix
111
112#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Definition: brookscorey.hh:35
Implementation of a exponential function relating specific interfacial area to wetting phase saturati...
Definition: exponentialcubic.hh:44
static Params< Scalar > makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: exponentialcubic.hh:78
static Scalar area(const Scalar swe, const Scalar pc, const Params< Scalar > &params)
The interfacial area.
Definition: exponentialcubic.hh:98
void setA2(Scalar a2)
Definition: exponentialcubic.hh:57
Params(Scalar a1=0, Scalar a2=0, Scalar a3=0)
Definition: exponentialcubic.hh:50
Scalar a1() const
Definition: exponentialcubic.hh:53
void setA1(Scalar a1)
Definition: exponentialcubic.hh:54
bool operator==(const Params &p) const
Definition: exponentialcubic.hh:62
Scalar a2() const
Definition: exponentialcubic.hh:56
Scalar a3() const
Definition: exponentialcubic.hh:59
void setA3(Scalar a3)
Definition: exponentialcubic.hh:60