3.3.0
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#warning "This header is deprecated. Removal after 3.3. Use new material laws."
30
31#include "efftoabslawparams.hh"
32
33namespace Dumux {
34
60template <class EffLawT, class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params> >
62{
63 using EffLaw = EffLawT;
64
65public:
66 using Params = AbsParamsT;
67 using Scalar = typename EffLaw::Scalar;
68
79 static Scalar pc(const Params &params, Scalar sw)
80 {
81 return EffLaw::pc(params, swToSwe(params, sw));
82 }
83
94 static Scalar sw(const Params &params, Scalar pc)
95 {
96 return sweToSw(params, EffLaw::sw(params, pc));
97 }
98
106 static Scalar endPointPc(const Params &params)
107 { return EffLaw::endPointPc(params); }
108
125 static Scalar dpc_dsw(const Params &params, Scalar sw)
126 {
127 return EffLaw::dpc_dswe(params, swToSwe(params, sw) )*dswe_dsw(params);
128 }
129
148 static Scalar dsw_dpc(const Params &params, Scalar pc)
149 {
150 return EffLaw::dswe_dpc(params, pc)*dsw_dswe(params);
151 }
152
164 static Scalar krw(const Params &params, Scalar sw)
165 {
166 return EffLaw::krw(params, swToSwe(params, sw));
167 }
168
177 static Scalar dkrw_dsw(const Params &params, Scalar sw)
178 {
179 return EffLaw::dkrw_dswe(params, swToSwe(params, sw))*dswe_dsw(params);
180 }
181
193 static Scalar krn(const Params &params, Scalar sw)
194 {
195 return EffLaw::krn(params, swToSwe(params, sw));
196 }
197
206 static Scalar dkrn_dsw(const Params &params, Scalar sw)
207 {
208 return EffLaw::dkrn_dswe(params, swToSwe(params, sw))*dswe_dsw(params);
209 }
210
220 static Scalar swToSwe(const Params &params, Scalar sw)
221 {
222 return (sw - params.swr())/(1. - params.swr() - params.snr());
223 }
224
234 static Scalar snToSne(const Params &params, Scalar sn)
235 {
236 return (sn - params.snr())/(1. - params.swr() - params.snr());
237 }
238
248 static Scalar sweToSw(const Params &params, Scalar swe)
249 { return swe*(1. - params.swr() - params.snr()) + params.swr(); }
250
260 static Scalar sneToSn(const Params &params, Scalar sne)
261 { return sne*(1. - params.swr() - params.snr()) + params.snr(); }
262
271 static Scalar dswe_dsw(const Params &params)
272 { return 1.0/(1. - params.swr() - params.snr()); }
273
282 static Scalar dsw_dswe(const Params &params)
283 { return 1. - params.swr() - params.snr(); }
284};
285} // end namespace Dumux
286
287#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:62
static Scalar dkrn_dsw(const Params &params, Scalar sw)
Returns the partial derivative of the relative permeability of the nonwetting phase with respect to t...
Definition: 2p/efftoabslaw.hh:206
static Scalar krn(const Params &params, Scalar sw)
The relative permeability for the nonwetting phase.
Definition: 2p/efftoabslaw.hh:193
static Scalar snToSne(const Params &params, Scalar sn)
Convert an absolute nonwetting saturation to an effective one.
Definition: 2p/efftoabslaw.hh:234
static Scalar endPointPc(const Params &params)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: 2p/efftoabslaw.hh:106
static Scalar sweToSw(const Params &params, Scalar swe)
Convert an effective wetting saturation to an absolute one.
Definition: 2p/efftoabslaw.hh:248
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:125
static Scalar dswe_dsw(const Params &params)
Derivative of the effective saturation w.r.t. the absolute saturation.
Definition: 2p/efftoabslaw.hh:271
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:148
static Scalar sw(const Params &params, Scalar pc)
The saturation-capillary pressure curve.
Definition: 2p/efftoabslaw.hh:94
AbsParamsT Params
Definition: 2p/efftoabslaw.hh:66
static Scalar swToSwe(const Params &params, Scalar sw)
Convert an absolute wetting saturation to an effective one.
Definition: 2p/efftoabslaw.hh:220
static Scalar dsw_dswe(const Params &params)
Derivative of the absolute saturation w.r.t. the effective saturation.
Definition: 2p/efftoabslaw.hh:282
static Scalar pc(const Params &params, Scalar sw)
The capillary pressure-saturation curve.
Definition: 2p/efftoabslaw.hh:79
typename EffLaw::Scalar Scalar
Definition: 2p/efftoabslaw.hh:67
static Scalar sneToSn(const Params &params, Scalar sne)
Convert an effective nonwetting saturation to an absolute one.
Definition: 2p/efftoabslaw.hh:260
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:177
static Scalar krw(const Params &params, Scalar sw)
The relative permeability for the wetting phase.
Definition: 2p/efftoabslaw.hh:164
A default implementation of the parameters for the adapter class to convert material laws from effect...