3.3.0
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
34#warning "This header is deprecated. Removal after 3.3. Use new material laws."
35
36namespace Dumux {
37
44template <class ParamsT>
45class [[deprecated("Use new material laws and FluidMatrix::InterfacialAreaExponential instead!")]] AwnSurfaceExpFct
46{
47public:
48 using Params = ParamsT;
49 using Scalar = typename Params::Scalar;
50
62 static Scalar interfacialArea(const Params & params, const Scalar Sw, const Scalar pc)
63 {
64 const Scalar a1 = params.a1();
65 const Scalar a2 = params.a2();
66 const Scalar a3 = params.a3();
67 const Scalar Swr = params.Swr();
68 using std::exp;
69 const Scalar aAlphaBeta = a1 * (Swr-Sw) * (1-Sw) + a2 * (Swr-Sw) * (1-Sw) * exp( a3 * pc) ;
70 return aAlphaBeta;
71 }
72
79 static Scalar dawn_dpc (const Params & params, const Scalar Sw, const Scalar pc)
80 {
81 const Scalar a2 = params.a2();
82 const Scalar a3 = params.a3();
83 const Scalar Swr = params.Swr();
84 using std::exp;
85 const Scalar value = a2 * a3 * (Swr-Sw) * (1-Sw) * exp(a3*pc);
86 return value;
87 }
88
95 static Scalar dawn_dsw (const Params & params, const Scalar Sw, const Scalar pc)
96 {
97 Scalar value;
98 Scalar a1 = params.a1();
99 Scalar a2 = params.a2();
100 Scalar a3 = params.a3();
101 Scalar Swr = params.Swr();
102 using std::exp;
103 value = - a1 *( 1+Swr-2*Sw ) - a2 * exp(a3*pc) * ( 1+Swr-2*Sw );
104 return value;
105 }
106
107};
108} // namespace Dumux
109
110#endif
Specification of the parameters for a function relating volume specific interfacial area to capillary...
Definition: adapt.hh:29
Implementation of the exponential function relating specific interfacial area to wetting phase satura...
Definition: awnsurfaceexpfct.hh:46
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:79
ParamsT Params
Definition: awnsurfaceexpfct.hh:48
typename Params::Scalar Scalar
Definition: awnsurfaceexpfct.hh:49
static Scalar interfacialArea(const Params &params, const Scalar Sw, const Scalar pc)
The interfacial area surface.
Definition: awnsurfaceexpfct.hh:62
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:95