3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
awnsurfacepolynomialedgezero2ndorder.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_POLYNOMIAL_EDGE_ZERO_2ND_ORDER_HH
24#define AWN_SURFACE_POLYNOMIAL_EDGE_ZERO_2ND_ORDER_HH
25
27
28#include <dune/common/exceptions.hh>
29
30#include <algorithm>
31#include <math.h>
32#include <assert.h>
33
34namespace Dumux {
35
41template <class ParamsT>
43{
44public:
45 using Params = ParamsT;
46 using Scalar = typename Params::Scalar;
47
58 static Scalar interfacialArea(const Params &params, const Scalar Sw, const Scalar pc)
59 {
60 const Scalar a1 = params.a1();
61 const Scalar a2 = params.a2();
62 const Scalar a3 = params.a3();
63 const Scalar Swr = params.Swr();
64 const Scalar factor = (Swr-Sw)*(1.-Sw) ;
65
66 const Scalar aAlphaBeta = a1*factor + a2*factor*pc + a3*factor*pc*pc;
67 return aAlphaBeta;
68 }
69
70
78 static Scalar dawn_dpc (const Params &params, const Scalar Sw, const Scalar pc)
79 {
80 const Scalar Swr = params.Swr();
81 const Scalar a1 = params.a1();
82 const Scalar a2 = params.a2();
83 const Scalar a3 = params.a3();
84 const Scalar value = a2*(Swr-Sw)*(1-Sw) + a3*(Swr-Sw)*(1-Sw)*pc;
85 return value;
86 }
87
95 static Scalar dawn_dsw (const Params &params, const Scalar Sw, const Scalar pc)
96 {
97 const Scalar Swr = params.Swr();
98 const Scalar a1 = params.a1();
99 const Scalar a2 = params.a2();
100 const Scalar a3 = params.a3();
101 const Scalar derivativeFactor = (Sw-1.)+(Sw-Swr);
102 const Scalar value = a1 * derivativeFactor + a2 * derivativeFactor * pc + a3 * derivativeFactor * pc*pc ;
103 return value;
104 }
105
106};
107} // namespace Dumux
108
109#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 polynomial of second order relating specific interfacial area to wetting phase ...
Definition: awnsurfacepolynomialedgezero2ndorder.hh:43
ParamsT Params
Definition: awnsurfacepolynomialedgezero2ndorder.hh:45
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: awnsurfacepolynomialedgezero2ndorder.hh:78
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: awnsurfacepolynomialedgezero2ndorder.hh:95
typename Params::Scalar Scalar
Definition: awnsurfacepolynomialedgezero2ndorder.hh:46
static Scalar interfacialArea(const Params &params, const Scalar Sw, const Scalar pc)
The awn surface the suggested (as estimated from pore network models) awn surface: ].
Definition: awnsurfacepolynomialedgezero2ndorder.hh:58