3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
awnsurfacepolynomial2ndorder.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_2ND_ORDER_HH
24#define AWN_SURFACE_POLYNOMIAL_2ND_ORDER_HH
25
27
28#include <dune/common/exceptions.hh>
29
30#include <algorithm>
31#include <cmath>
32#include <assert.h>
33
34namespace Dumux {
35
41template <class ScalarT, class ParamsT = AwnSurfacePolynomial2ndOrderParams<ScalarT> >
43{
44public:
45 using Params = ParamsT;
46 using Scalar = typename Params::Scalar;
47
59 static Scalar interfacialArea(const Params & params, const Scalar Sw, const Scalar pc)
60 {
61 const Scalar a00 = params.a00();
62 const Scalar a10 = params.a10();
63 const Scalar a20 = params.a20();
64 const Scalar a11 = params.a11();
65 const Scalar a01 = params.a01();
66 const Scalar a02 = params.a02();
67
68 using std::pow;
69 const Scalar aAlphaBeta = a00 + a10 * Sw + a20 * pow(Sw,2) + a11*Sw*pc + a01*pc + a02*pow(pc,2);
70 return aAlphaBeta;
71 }
72
73
81 static Scalar dawn_dpc (const Params & params, const Scalar Sw, const Scalar pc)
82 {
83 const Scalar value = params.a11()*Sw + params.a01() + 2.0*params.a02() * pc;
84 return value;
85 }
86
94 static Scalar dawn_dsw (const Params &params, const Scalar Sw, const Scalar pc)
95 {
96 const Scalar value = params.a11()*pc + params.a10() + 2.0*params.a20()*Sw;
97 return value;
98 }
99
100};
101} // end namespace Dumux
102
103#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: awnsurfacepolynomial2ndorder.hh:43
ParamsT Params
Definition: awnsurfacepolynomial2ndorder.hh:45
static Scalar interfacialArea(const Params &params, const Scalar Sw, const Scalar pc)
The awn surface.
Definition: awnsurfacepolynomial2ndorder.hh:59
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: awnsurfacepolynomial2ndorder.hh:81
typename Params::Scalar Scalar
Definition: awnsurfacepolynomial2ndorder.hh:46
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: awnsurfacepolynomial2ndorder.hh:94