3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
regularizedlinearmaterial.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 DUMUX_REGULARIZED_LINEAR_MATERIAL_HH
26#define DUMUX_REGULARIZED_LINEAR_MATERIAL_HH
27
28#include "linearmaterial.hh"
30
32
33namespace Dumux {
34
56template <class ScalarT, class ParamsT = RegularizedLinearMaterialParams<ScalarT> >
58{
60
61public:
62 using Params = ParamsT;
63 using Scalar = typename Params::Scalar;
64
78 static Scalar pc(const Params &params, Scalar swe)
79 {
80 return LinearMaterial::pc(params, swe);
81 }
82
97 static Scalar sw(const Params &params, Scalar pc)
98 {
99 return LinearMaterial::sw(params, pc);
100 }
101
109 static Scalar endPointPc(const Params &params)
110 { return params.entryPc(); }
111
126 static Scalar dpc_dswe(const Params &params, Scalar swe)
127 {
128 return LinearMaterial::dpc_dswe(params, swe);
129 }
130
140 static Scalar dswe_dpc(const Params &params, Scalar pc)
141 {
142 return LinearMaterial::dswe_dpc(params, pc);
143 }
144
153 static Scalar krw(const Params &params, Scalar swe)
154 {
155 return relperm_(params, swe);
156 }
157
166 static Scalar krn(const Params &params, Scalar swe)
167 {
168 Scalar sne = 1 - swe;
169 return relperm_(params, sne);
170 }
171
172private:
173 static Scalar relperm_(const Params &params, Scalar S)
174 {
175 const Scalar lowS = params.krLowS();
176 const Scalar highS = params.krHighS();
177
178
179 const Scalar m = (1 - ((1 - highS) + lowS)/2 ) / (1 - (1 - highS) - lowS);
180
181 // check whether the saturation is unpyhsical
182 if (S >= 1.0)
183 return 1.0;
184 else if (S <= 0.0)
185 return 0;
186 // check wether the permeability needs to be regularized
187 else if (S < lowS) {
189 Spline sp(0, lowS,
190 0, lowS/2,
191 0, m);
192 return sp.eval(S);
193 }
194 else if (S > highS) {
195 using Spline = Dumux::Spline<Scalar>;
196 Spline sp(highS, 1,
197 1 - (1 - highS)/2, 1,
198 m, 0);
199 return sp.eval(S);
200 }
201
202 // straight line for S \in [lowS, highS]
203 return lowS/2 + m*(S - lowS);
204 }
205};
206} // end namespace Dumux
207
208#endif
Provides 3rd order polynomial splines.
Linear capillary pressure and relative permeability <-> saturation relations.
Parameters that are necessary for the regularization of the linear constitutive relations.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
A 3rd order polynomial spline.
Definition: spline.hh:55
Linear capillary pressure and relative permeability <-> saturation relations.
Definition: linearmaterial.hh:48
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 pc(const Params &params, Scalar swe)
The linear capillary pressure-saturation curve.
Definition: linearmaterial.hh:67
Implements a linear saturation-capillary pressure relation.
Definition: regularizedlinearmaterial.hh:58
static Scalar endPointPc(const Params &params)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: regularizedlinearmaterial.hh:109
ParamsT Params
Definition: regularizedlinearmaterial.hh:62
static Scalar dswe_dpc(const Params &params, Scalar pc)
Returns the partial derivative of the effective saturation to the capillary pressure.
Definition: regularizedlinearmaterial.hh:140
static Scalar sw(const Params &params, Scalar pc)
The saturation-capillary pressure curve.
Definition: regularizedlinearmaterial.hh:97
static Scalar dpc_dswe(const Params &params, Scalar swe)
Returns the partial derivative of the capillary pressure to the effective saturation.
Definition: regularizedlinearmaterial.hh:126
static Scalar pc(const Params &params, Scalar swe)
The linear capillary pressure-saturation curve.
Definition: regularizedlinearmaterial.hh:78
static Scalar krn(const Params &params, Scalar swe)
The relative permeability for the non-wetting phase.
Definition: regularizedlinearmaterial.hh:166
typename Params::Scalar Scalar
Definition: regularizedlinearmaterial.hh:63
static Scalar krw(const Params &params, Scalar swe)
The relative permeability for the wetting phase.
Definition: regularizedlinearmaterial.hh:153