3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
awnsurfaceexpswpcto3.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_SW_PC_TO_3
24#define AWN_SURFACE_EXP_SW_PC_TO_3
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 =AwnSurfaceExpSwPcTo3Params<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 // TODO think about awn surface for relative saturation
62 const Scalar a1 = params.a1();
63 const Scalar a2 = params.a2();
64 const Scalar a3 = params.a3();
65
66 using std::exp;
67 const Scalar aAlphaBeta = a1 * exp( a2 * Sw) + a3 * pc * pc * pc ;
68 return aAlphaBeta;
69 }
70
78 static Scalar dawn_dpc (const Params &params, const Scalar Sw, const Scalar pc)
79 {
80 DUNE_THROW(Dune::NotImplemented, __FILE__ << " dawndpc()");
81 }
82
90 static Scalar dawn_dsw (const Params &params, const Scalar Sw, const Scalar pc)
91 {
92 DUNE_THROW(Dune::NotImplemented, __FILE__ << " dawndSw()");
93 }
94
95};
96} // end namespace Dumux
97
98#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 a exponential function relating specific interfacial area to wetting phase saturati...
Definition: awnsurfaceexpswpcto3.hh:43
ParamsT Params
Definition: awnsurfaceexpswpcto3.hh:45
typename Params::Scalar Scalar
Definition: awnsurfaceexpswpcto3.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: awnsurfaceexpswpcto3.hh:90
static Scalar interfacialArea(const Params &params, const Scalar Sw, const Scalar pc)
The awn surface.
Definition: awnsurfaceexpswpcto3.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: awnsurfaceexpswpcto3.hh:78