3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
linearmaterial.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef LINEAR_MATERIAL_HH
26#define LINEAR_MATERIAL_HH
27
29
30#include <algorithm>
31
32namespace Dumux {
33
46template <class ScalarT, class ParamsT = LinearMaterialParams<ScalarT> >
48{
49public:
50 using Params = ParamsT;
51 using Scalar = typename Params::Scalar;
52
67 static Scalar pc(const Params &params, Scalar swe)
68 {
69 return (1 - swe)*(params.maxPc() - params.entryPc()) + params.entryPc();
70 }
71
86 static Scalar sw(const Params &params, Scalar pc)
87 {
88 return 1 - (pc - params.entryPc())/(params.maxPc() - params.entryPc());
89 }
90
98 static Scalar endPointPc(const Params &params)
99 { return params.entryPc(); }
100
116 static Scalar dpc_dswe(const Params &params, Scalar swe)
117 {
118 return - (params.maxPc() - params.entryPc());
119 }
120
131 static Scalar dsw_dpc(const Params &params, Scalar pc)
132 {
133 return - 1/(params.maxPc() - params.entryPc());
134 }
135
145 static Scalar krw(const Params &params, Scalar swe)
146 {
147 using std::max;
148 using std::min;
149 return max(min(swe,1.0),0.0);
150 }
151
161 static Scalar krn(const Params &params, Scalar swe)
162 {
163 Scalar sne = 1 - swe;
164 using std::max;
165 using std::min;
166 return max(min(sne,1.0),0.0);
167 }
168};
169} // end namespace Dumux
170
171#endif
Parameters for the linear capillary pressure and relative permeability <-> saturation relations.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Linear capillary pressure and relative permeability <-> saturation relations.
Definition: linearmaterial.hh:48
static Scalar krw(const Params &params, Scalar swe)
The relative permeability for the wetting phase.
Definition: linearmaterial.hh:145
static Scalar sw(const Params &params, Scalar pc)
The saturation-capillary pressure curve.
Definition: linearmaterial.hh:86
static Scalar dpc_dswe(const Params &params, Scalar swe)
Returns the partial derivative of the capillary pressure w.r.t. the effective saturation.
Definition: linearmaterial.hh:116
static Scalar endPointPc(const Params &params)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: linearmaterial.hh:98
static Scalar pc(const Params &params, Scalar swe)
The linear capillary pressure-saturation curve.
Definition: linearmaterial.hh:67
typename Params::Scalar Scalar
Definition: linearmaterial.hh:51
static Scalar dsw_dpc(const Params &params, Scalar pc)
Returns the partial derivative of the effective saturation to the capillary pressure.
Definition: linearmaterial.hh:131
static Scalar krn(const Params &params, Scalar swe)
The relative permeability for the non-wetting phase.
Definition: linearmaterial.hh:161
ParamsT Params
Definition: linearmaterial.hh:50