3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
3p/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 <dune/common/exceptions.hh>
30#include "efftoabslawparams.hh"
31
32#warning "This header is deprecated. Removal after 3.3. Use new material laws."
33
34namespace Dumux {
35
61template <class EffLawT, class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params> >
62class EffToAbsLaw
63{
64 using EffLaw = EffLawT;
65
66public:
67 using Params = AbsParamsT;
68 using Scalar = typename EffLaw::Scalar;
69
80 static Scalar pc(const Params &params, const Scalar sw)
81 {
82 return EffLaw::pc(params, swToSwe(params, sw));
83 }
84
90 static Scalar pcgw(const Params &params, const Scalar sw)
91 {
92 return EffLaw::pcgw(params, swToSwe(params, sw));
93 }
94
100 static Scalar pcnw(const Params &params, const Scalar sw)
101 {
102 return EffLaw::pcnw(params, swToSwe(params, sw));
103 }
104
110 static Scalar pcgn(const Params &params, const Scalar st)
111 {
112 return EffLaw::pcgn(params, stToSte(params, st));
113 }
114
120 static Scalar pcAlpha(const Params &params, const Scalar sn)
121 {
122 return EffLaw::pcAlpha(params, sn);
123 }
124
136 static Scalar sw(const Params &params, const Scalar pc)
137 {
138 return EffLaw::sw(params, pc);
139 }
140
156 static Scalar dpc_dsw(const Params &params, const Scalar sw)
157 {
158 return EffLaw::dpc_dswe(params, pc);
159 }
160
178 static Scalar dsw_dpc(const Params &params, const Scalar pc)
179 {
180 return EffLaw::dsw_dpc(params, pc);
181 }
182
197 static Scalar krw(const Params &params, Scalar sw, const Scalar sn)
198 {
199 return EffLaw::krw(params, swToSwe(params, sw));
200 }
201
215 static Scalar krn(const Params &params, const Scalar sw, const Scalar sn)
216 {
217 const Scalar st = sw+sn;
218 return EffLaw::krn(params, swToSwe(params, sw), sn, stToSte(params, st));
219 }
220
234 static Scalar krg(const Params &params, const Scalar sw, const Scalar sn)
235 {
236 const Scalar st = sw+sn;
237 return EffLaw::krg(params, stToSte(params, st));
238 }
239
248 static Scalar kr(const Params &params, const int phaseIdx, const Scalar sw, const Scalar sn, const Scalar sg)
249 {
250 const Scalar st = sw+sn;
251 return EffLaw::kr(params, phaseIdx, swToSwe(params, sw), sn, stToSte(params, st));
252 }
253
259 {
260 return EffLaw::bulkDensTimesAdsorpCoeff(params);
261 }
262
272 static Scalar swToSwe(const Params &params, const Scalar sw)
273 {
274 return (sw-params.swr())/(1.-params.swr());
275 }
276
286 static Scalar snToSne(const Params &params, const Scalar sn)
287 {
288 return sn; // sne equals sn
289 }
290
300 static Scalar stToSte(const Params &params, const Scalar st)
301 {
302 return (st-params.swr()) / (1-params.swr());
303 }
304
314 static Scalar sgToSge(const Params &params, Scalar sg)
315 {
316 DUNE_THROW(Dune::NotImplemented, "sgTosge for three phases not implemented!");
317 }
318
319//private:
329 static Scalar sweToSw_(const Params &params, Scalar swe)
330 {
331 DUNE_THROW(Dune::NotImplemented, "sweTosw for three phases not implemented!");
332 }
333
334 static Scalar sneToSn_(const Params &params, Scalar sne)
335 {
336 DUNE_THROW(Dune::NotImplemented, "sneTosn for three phases not implemented!");
337 }
338
339 static Scalar sgeToSg_(const Params &params, Scalar sge)
340 {
341 DUNE_THROW(Dune::NotImplemented, "sgeTosg for three phases not implemented!");
342 }
351 static Scalar dswe_dsw_(const Params &params)
352 {
353 DUNE_THROW(Dune::NotImplemented, "dswe_dsw for three phases not implemented!");
354 }
355
364 static Scalar dsw_dswe_(const Params &params)
365 {
366 DUNE_THROW(Dune::NotImplemented, "dsw_dswe for three phases not implemented!");
367 }
368};
369} // end namespace Dumux
370
371#endif
Definition: adapt.hh:29
static Scalar pcgn(const Params &params, const Scalar st)
The capillary pressure-saturation curve for the gas and nonwetting phase.
Definition: 3p/efftoabslaw.hh:110
static Scalar pcgw(const Params &params, const Scalar sw)
The capillary pressure-saturation curve for the gas and wetting phase.
Definition: 3p/efftoabslaw.hh:90
static Scalar krw(const Params &params, Scalar sw, const Scalar sn)
The relative permeability for the wetting phase.
Definition: 3p/efftoabslaw.hh:197
static Scalar dpc_dsw(const Params &params, const Scalar sw)
Returns the partial derivative of the capillary pressure w.r.t the absolute saturation....
Definition: 3p/efftoabslaw.hh:156
static Scalar pcAlpha(const Params &params, const Scalar sn)
This function ensures a continuous transition from 2 to 3 phases and vice versa.
Definition: 3p/efftoabslaw.hh:120
static Scalar dswe_dsw_(const Params &params)
Derivative of the effective saturation w.r.t. the absolute saturation.
Definition: 3p/efftoabslaw.hh:351
static Scalar dsw_dpc(const Params &params, const Scalar pc)
Returns the partial derivative of the absolute saturation w.r.t. the capillary pressure....
Definition: 3p/efftoabslaw.hh:178
static Scalar sw(const Params &params, Scalar pc)
The saturation-capillary pressure curve.
Definition: 2p/efftoabslaw.hh:94
static Scalar snToSne(const Params &params, const Scalar sn)
Convert an absolute nonwetting saturation to an effective one.
Definition: 3p/efftoabslaw.hh:286
static Scalar pcnw(const Params &params, const Scalar sw)
The capillary pressure-saturation curve the nonwetting and wetting phase.
Definition: 3p/efftoabslaw.hh:100
static Scalar krg(const Params &params, const Scalar sw, const Scalar sn)
The relative permeability for the gas phase.
Definition: 3p/efftoabslaw.hh:234
static Scalar krn(const Params &params, const Scalar sw, const Scalar sn)
The relative permeability for the nonwetting phase.
Definition: 3p/efftoabslaw.hh:215
static Scalar sgToSge(const Params &params, Scalar sg)
Convert an absolute gas saturation to an effective one.
Definition: 3p/efftoabslaw.hh:314
static Scalar sneToSn_(const Params &params, Scalar sne)
Definition: 3p/efftoabslaw.hh:334
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 sw(const Params &params, const Scalar pc)
The saturation-capillary pressure curve.
Definition: 3p/efftoabslaw.hh:136
static Scalar bulkDensTimesAdsorpCoeff(const Params &params)
the basis for calculating adsorbed NAPL in storage term
Definition: 3p/efftoabslaw.hh:258
static Scalar dsw_dswe_(const Params &params)
Derivative of the absolute saturation w.r.t. the effective saturation.
Definition: 3p/efftoabslaw.hh:364
static Scalar swToSwe(const Params &params, const Scalar sw)
Convert an absolute wetting saturation to an effective one.
Definition: 3p/efftoabslaw.hh:272
static Scalar sgeToSg_(const Params &params, Scalar sge)
Definition: 3p/efftoabslaw.hh:339
static Scalar sweToSw_(const Params &params, Scalar swe)
Convert an effective wetting saturation to an absolute one.
Definition: 3p/efftoabslaw.hh:329
static Scalar kr(const Params &params, const int phaseIdx, const Scalar sw, const Scalar sn, const Scalar sg)
The relative permeability for a phase.
Definition: 3p/efftoabslaw.hh:248
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 stToSte(const Params &params, const Scalar st)
Convert an absolute total liquid saturation to an effective one.
Definition: 3p/efftoabslaw.hh:300
static Scalar pc(const Params &params, const Scalar sw)
The capillary pressure-saturation curve.
Definition: 3p/efftoabslaw.hh:80
A default implementation of the parameters for the adapter class to convert material laws from effect...