3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
2p/efftoabslaw.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 *****************************************************************************/
26#ifndef DUMUX_EFF_TO_ABS_LAW_HH
27#define DUMUX_EFF_TO_ABS_LAW_HH
28
29#include "efftoabslawparams.hh"
30
31namespace Dumux {
32
58template <class EffLawT, class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params> >
60{
61 using EffLaw = EffLawT;
62
63public:
64 using Params = AbsParamsT;
65 using Scalar = typename EffLaw::Scalar;
66
77 static Scalar pc(const Params &params, Scalar sw)
78 {
79 return EffLaw::pc(params, swToSwe(params, sw));
80 }
81
92 static Scalar sw(const Params &params, Scalar pc)
93 {
94 return sweToSw(params, EffLaw::sw(params, pc));
95 }
96
104 static Scalar endPointPc(const Params &params)
105 { return EffLaw::endPointPc(params); }
106
123 static Scalar dpc_dsw(const Params &params, Scalar sw)
124 {
125 return EffLaw::dpc_dswe(params, swToSwe(params, sw) )*dswe_dsw(params);
126 }
127
146 static Scalar dsw_dpc(const Params &params, Scalar pc)
147 {
148 return EffLaw::dswe_dpc(params, pc)*dsw_dswe(params);
149 }
150
162 static Scalar krw(const Params &params, Scalar sw)
163 {
164 return EffLaw::krw(params, swToSwe(params, sw));
165 }
166
175 static Scalar dkrw_dsw(const Params &params, Scalar sw)
176 {
177 return EffLaw::dkrw_dswe(params, swToSwe(params, sw))*dswe_dsw(params);
178 }
179
191 static Scalar krn(const Params &params, Scalar sw)
192 {
193 return EffLaw::krn(params, swToSwe(params, sw));
194 }
195
204 static Scalar dkrn_dsw(const Params &params, Scalar sw)
205 {
206 return EffLaw::dkrn_dswe(params, swToSwe(params, sw))*dswe_dsw(params);
207 }
208
218 static Scalar swToSwe(const Params &params, Scalar sw)
219 {
220 return (sw - params.swr())/(1. - params.swr() - params.snr());
221 }
222
232 static Scalar snToSne(const Params &params, Scalar sn)
233 {
234 return (sn - params.snr())/(1. - params.swr() - params.snr());
235 }
236
246 static Scalar sweToSw(const Params &params, Scalar swe)
247 { return swe*(1. - params.swr() - params.snr()) + params.swr(); }
248
249 [[deprecated("Will be removed after 3.2. Use sweToSw (without underscore suffix) instead!")]]
250 static Scalar sweToSw_(const Params &params, Scalar swe)
251 { return sweToSw(params, swe); }
252
262 static Scalar sneToSn(const Params &params, Scalar sne)
263 {
264 return sne*(1. - params.swr() - params.snr()) + params.snr();
265 }
266
275 static Scalar dswe_dsw(const Params &params)
276 { return 1.0/(1. - params.swr() - params.snr()); }
277
278 [[deprecated("Will be removed after 3.2. Use dswe_dsw (without underscore suffix) instead!")]]
279 static Scalar dswe_dsw_(const Params &params)
280 { return dswe_dsw(params); }
281
290 static Scalar dsw_dswe(const Params &params)
291 { return 1. - params.swr() - params.snr(); }
292
293 [[deprecated("Will be removed after 3.2. Use dsw_dswe (without underscore suffix) instead!")]]
294 static Scalar dsw_dswe_(const Params &params)
295 { return dsw_dswe(params); }
296};
297} // end namespace Dumux
298
299#endif
Definition: adapt.hh:29
This material law takes a material law defined for effective saturations and converts it to a materia...
Definition: 2p/efftoabslaw.hh:60
static Scalar dkrn_dsw(const Params &params, Scalar sw)
Returns the partial derivative of the relative permeability of the non-wetting phase with respect to ...
Definition: 2p/efftoabslaw.hh:204
static Scalar krn(const Params &params, Scalar sw)
The relative permeability for the non-wetting phase.
Definition: 2p/efftoabslaw.hh:191
static Scalar snToSne(const Params &params, Scalar sn)
Convert an absolute non-wetting saturation to an effective one.
Definition: 2p/efftoabslaw.hh:232
static Scalar endPointPc(const Params &params)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: 2p/efftoabslaw.hh:104
static Scalar sweToSw(const Params &params, Scalar swe)
Convert an effective wetting saturation to an absolute one.
Definition: 2p/efftoabslaw.hh:246
static Scalar dpc_dsw(const Params &params, Scalar sw)
Returns the partial derivative of the capillary pressure w.r.t the absolute saturation.
Definition: 2p/efftoabslaw.hh:123
static Scalar dswe_dsw(const Params &params)
Derivative of the effective saturation w.r.t. the absolute saturation.
Definition: 2p/efftoabslaw.hh:275
static Scalar dswe_dsw_(const Params &params)
Definition: 2p/efftoabslaw.hh:279
static Scalar dsw_dpc(const Params &params, Scalar pc)
Returns the partial derivative of the absolute saturation w.r.t. the capillary pressure.
Definition: 2p/efftoabslaw.hh:146
static Scalar sw(const Params &params, Scalar pc)
The saturation-capillary pressure curve.
Definition: 2p/efftoabslaw.hh:92
AbsParamsT Params
Definition: 2p/efftoabslaw.hh:64
static Scalar swToSwe(const Params &params, Scalar sw)
Convert an absolute wetting saturation to an effective one.
Definition: 2p/efftoabslaw.hh:218
static Scalar dsw_dswe_(const Params &params)
Definition: 2p/efftoabslaw.hh:294
static Scalar dsw_dswe(const Params &params)
Derivative of the absolute saturation w.r.t. the effective saturation.
Definition: 2p/efftoabslaw.hh:290
static Scalar sweToSw_(const Params &params, Scalar swe)
Definition: 2p/efftoabslaw.hh:250
static Scalar pc(const Params &params, Scalar sw)
The capillary pressure-saturation curve.
Definition: 2p/efftoabslaw.hh:77
typename EffLaw::Scalar Scalar
Definition: 2p/efftoabslaw.hh:65
static Scalar sneToSn(const Params &params, Scalar sne)
Convert an effective non-wetting saturation to an absolute one.
Definition: 2p/efftoabslaw.hh:262
static Scalar dkrw_dsw(const Params &params, Scalar sw)
Returns the partial derivative of the relative permeability of the wetting phase with respect to the ...
Definition: 2p/efftoabslaw.hh:175
static Scalar krw(const Params &params, Scalar sw)
The relative permeability for the wetting phase.
Definition: 2p/efftoabslaw.hh:162
A default implementation of the parameters for the adapter class to convert material laws from effect...