3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
awnsurfaceexpfct.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 AWN_SURFACE_EXP_FCT_HH
24#define AWN_SURFACE_EXP_FCT_HH
25
27
28#include <dune/common/exceptions.hh>
29
30#include <algorithm>
31#include <cmath>
32#include <assert.h>
33
34namespace Dumux {
35
42template <class ParamsT>
44{
45public:
46 using Params = ParamsT;
47 using Scalar = typename Params::Scalar;
48
60 static Scalar interfacialArea(const Params & params, const Scalar Sw, const Scalar pc)
61 {
62 const Scalar a1 = params.a1();
63 const Scalar a2 = params.a2();
64 const Scalar a3 = params.a3();
65 const Scalar Swr = params.Swr();
66 using std::exp;
67 const Scalar aAlphaBeta = a1 * (Swr-Sw) * (1-Sw) + a2 * (Swr-Sw) * (1-Sw) * exp( a3 * pc) ;
68 return aAlphaBeta;
69 }
70
77 static Scalar dawn_dpc (const Params & params, const Scalar Sw, const Scalar pc)
78 {
79 const Scalar a2 = params.a2();
80 const Scalar a3 = params.a3();
81 const Scalar Swr = params.Swr();
82 using std::exp;
83 const Scalar value = a2 * a3 * (Swr-Sw) * (1-Sw) * exp(a3*pc);
84 return value;
85 }
86
93 static Scalar dawn_dsw (const Params & params, const Scalar Sw, const Scalar pc)
94 {
95 Scalar value;
96 Scalar a1 = params.a1();
97 Scalar a2 = params.a2();
98 Scalar a3 = params.a3();
99 Scalar Swr = params.Swr();
100 using std::exp;
101 value = - a1 *( 1+Swr-2*Sw ) - a2 * exp(a3*pc) * ( 1+Swr-2*Sw );
102 return value;
103 }
104
105};
106} // namespace Dumux
107
108#endif
Specification of the parameters for a function relating volume specific interfacial area to capillary...
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Implementation of the exponential function relating specific interfacial area to wetting phase satura...
Definition: awnsurfaceexpfct.hh:44
static Scalar dawn_dpc(const Params &params, const Scalar Sw, const Scalar pc)
the derivative of specific interfacial area function w.r.t. capillary pressure
Definition: awnsurfaceexpfct.hh:77
ParamsT Params
Definition: awnsurfaceexpfct.hh:46
typename Params::Scalar Scalar
Definition: awnsurfaceexpfct.hh:47
static Scalar interfacialArea(const Params &params, const Scalar Sw, const Scalar pc)
The interfacial area surface.
Definition: awnsurfaceexpfct.hh:60
static Scalar dawn_dsw(const Params &params, const Scalar Sw, const Scalar pc)
the derivative of specific interfacial area function w.r.t. saturation
Definition: awnsurfaceexpfct.hh:93