Loading [MathJax]/extensions/tex2jax.js
3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
exponential.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_EXPONENTIAL
26#define DUMUX_MATERIAL_FLUIDMATRIX_TWO_P_INTERFACIAL_AREA_EXPONENTIAL
27
28#include <cmath>
29#include <dune/common/exceptions.hh>
30#include <dune/common/float_cmp.hh>
32
33namespace Dumux::FluidMatrix {
34
42{
43public:
44
45 template<class Scalar>
46 struct Params
47 {
48 Params(Scalar swr = 0, Scalar a1 = 0, Scalar a2 = 0, Scalar a3 = 0)
49 : swr_(swr), a1_(a1), a2_(a2), a3_(a3) {}
50
51 Scalar swr() const { return swr_; }
52 void setSwr(Scalar swr) { swr_ = swr; }
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(swr(), p.swr(), 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 swr_, a1_, a2_, a3_;
73 };
74
79 template<class Scalar = double>
80 static Params<Scalar> makeParams(const std::string& paramGroup)
81 {
82 const auto swr = getParamFromGroup<Scalar>(paramGroup, "Swr", 0.0);
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 {swr, a1, a2, a3};
87 }
88
89
101 template<class Scalar>
102 static Scalar area(const Scalar swe, const Scalar pc, const Params<Scalar>& params)
103 {
104 const Scalar a1 = params.a1();
105 const Scalar a2 = params.a2();
106 const Scalar a3 = params.a3();
107 const Scalar swr = params.swr();
108 using std::exp;
109 return a1 * (swr-swe) * (1-swe) + a2 * (swr-swe) * (1-swe) * exp( a3 * pc) ;
110 }
111
118 template<class Scalar>
119 static Scalar darea_dpc(const Scalar swe, const Scalar pc, const Params<Scalar>& params)
120 {
121 const Scalar a2 = params.a2();
122 const Scalar a3 = params.a3();
123 const Scalar swr = params.swr();
124 using std::exp;
125 return a2 * a3 * (swr-swe) * (1-swe) * exp(a3*pc);
126 }
127
134 template<class Scalar>
135 static Scalar darea_dsw(const Scalar swe, const Scalar pc, const Params<Scalar>& params)
136 {
137 Scalar a1 = params.a1();
138 Scalar a2 = params.a2();
139 Scalar a3 = params.a3();
140 Scalar swr = params.swr();
141 using std::exp;
142 return - a1 *(1+swr-2*swe) - a2 * exp(a3*pc) * (1+swr-2*swe);
143 }
144
145};
146
147} // namespace Dumux::FluidMatrix
148
149#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Definition: brookscorey.hh:35
Implementation of the exponential function relating specific interfacial area to wetting phase satura...
Definition: exponential.hh:42
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: exponential.hh:135
static Params< Scalar > makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: exponential.hh:80
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: exponential.hh:119
static Scalar area(const Scalar swe, const Scalar pc, const Params< Scalar > &params)
The interfacial area.
Definition: exponential.hh:102
void setSwr(Scalar swr)
Definition: exponential.hh:52
Scalar a3() const
Definition: exponential.hh:60
Params(Scalar swr=0, Scalar a1=0, Scalar a2=0, Scalar a3=0)
Definition: exponential.hh:48
Scalar a2() const
Definition: exponential.hh:57
void setA1(Scalar a1)
Definition: exponential.hh:55
void setA2(Scalar a2)
Definition: exponential.hh:58
Scalar a1() const
Definition: exponential.hh:54
Scalar swr() const
Definition: exponential.hh:51
bool operator==(const Params &p) const
Definition: exponential.hh:63
void setA3(Scalar a3)
Definition: exponential.hh:61