3.1-git
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
32namespace Dumux {
33
59template <class EffLawT, class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params> >
60class EffToAbsLaw
61{
62 using EffLaw = EffLawT;
63
64public:
65 using Params = AbsParamsT;
66 using Scalar = typename EffLaw::Scalar;
67
78 static Scalar pc(const Params &params, const Scalar sw)
79 {
80 return EffLaw::pc(params, swToSwe(params, sw));
81 }
82
88 static Scalar pcgw(const Params &params, const Scalar sw)
89 {
90 return EffLaw::pcgw(params, swToSwe(params, sw));
91 }
92
98 static Scalar pcnw(const Params &params, const Scalar sw)
99 {
100 return EffLaw::pcnw(params, swToSwe(params, sw));
101 }
102
108 static Scalar pcgn(const Params &params, const Scalar st)
109 {
110 return EffLaw::pcgn(params, stToSte(params, st));
111 }
112
118 static Scalar pcAlpha(const Params &params, const Scalar sn)
119 {
120 return EffLaw::pcAlpha(params, sn);
121 }
122
134 static Scalar sw(const Params &params, const Scalar pc)
135 {
136 return EffLaw::sw(params, pc);
137 }
138
154 static Scalar dpc_dsw(const Params &params, const Scalar sw)
155 {
156 return EffLaw::dpc_dswe(params, pc);
157 }
158
176 static Scalar dsw_dpc(const Params &params, const Scalar pc)
177 {
178 return EffLaw::dsw_dpc(params, pc);
179 }
180
195 static Scalar krw(const Params &params, Scalar sw, const Scalar sn)
196 {
197 return EffLaw::krw(params, swToSwe(params, sw));
198 }
199
213 static Scalar krn(const Params &params, const Scalar sw, const Scalar sn)
214 {
215 const Scalar st = sw+sn;
216 return EffLaw::krn(params, swToSwe(params, sw), sn, stToSte(params, st));
217 }
218
232 static Scalar krg(const Params &params, const Scalar sw, const Scalar sn)
233 {
234 const Scalar st = sw+sn;
235 return EffLaw::krg(params, stToSte(params, st));
236 }
237
246 static Scalar kr(const Params &params, const int phaseIdx, const Scalar sw, const Scalar sn, const Scalar sg)
247 {
248 const Scalar st = sw+sn;
249 return EffLaw::kr(params, phaseIdx, swToSwe(params, sw), sn, stToSte(params, st));
250 }
251
257 {
258 return EffLaw::bulkDensTimesAdsorpCoeff(params);
259 }
260
270 static Scalar swToSwe(const Params &params, const Scalar sw)
271 {
272 return (sw-params.swr())/(1.-params.swr());
273 }
274
284 static Scalar snToSne(const Params &params, const Scalar sn)
285 {
286 return sn; // sne equals sn
287 }
288
298 static Scalar stToSte(const Params &params, const Scalar st)
299 {
300 return (st-params.swr()) / (1-params.swr());
301 }
302
312 static Scalar sgToSge(const Params &params, Scalar sg)
313 {
314 DUNE_THROW(Dune::NotImplemented, "sgTosge for three phases not implemented!");
315 }
316
317//private:
327 static Scalar sweToSw_(const Params &params, Scalar swe)
328 {
329 DUNE_THROW(Dune::NotImplemented, "sweTosw for three phases not implemented!");
330 }
331
332 static Scalar sneToSn_(const Params &params, Scalar sne)
333 {
334 DUNE_THROW(Dune::NotImplemented, "sneTosn for three phases not implemented!");
335 }
336
337 static Scalar sgeToSg_(const Params &params, Scalar sge)
338 {
339 DUNE_THROW(Dune::NotImplemented, "sgeTosg for three phases not implemented!");
340 }
349 static Scalar dswe_dsw_(const Params &params)
350 {
351 DUNE_THROW(Dune::NotImplemented, "dswe_dsw for three phases not implemented!");
352 }
353
362 static Scalar dsw_dswe_(const Params &params)
363 {
364 DUNE_THROW(Dune::NotImplemented, "dsw_dswe for three phases not implemented!");
365 }
366};
367} // end namespace Dumux
368
369#endif
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
static Scalar pcgn(const Params &params, const Scalar st)
The capillary pressure-saturation curve for the gas and non-wetting phase.
Definition: 3p/efftoabslaw.hh:108
static Scalar pcgw(const Params &params, const Scalar sw)
The capillary pressure-saturation curve for the gas and wetting phase.
Definition: 3p/efftoabslaw.hh:88
static Scalar krw(const Params &params, Scalar sw, const Scalar sn)
The relative permeability for the wetting phase.
Definition: 3p/efftoabslaw.hh:195
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:154
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:118
static Scalar dswe_dsw_(const Params &params)
Derivative of the effective saturation w.r.t. the absolute saturation.
Definition: 3p/efftoabslaw.hh:349
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:176
static Scalar sw(const Params &params, Scalar pc)
The saturation-capillary pressure curve.
Definition: 2p/efftoabslaw.hh:92
static Scalar snToSne(const Params &params, const Scalar sn)
Convert an absolute non-wetting saturation to an effective one.
Definition: 3p/efftoabslaw.hh:284
static Scalar pcnw(const Params &params, const Scalar sw)
The capillary pressure-saturation curve the non-wetting and wetting phase.
Definition: 3p/efftoabslaw.hh:98
static Scalar krg(const Params &params, const Scalar sw, const Scalar sn)
The relative permeability for the gas phase.
Definition: 3p/efftoabslaw.hh:232
static Scalar krn(const Params &params, const Scalar sw, const Scalar sn)
The relative permeability for the non-wetting phase.
Definition: 3p/efftoabslaw.hh:213
static Scalar sgToSge(const Params &params, Scalar sg)
Convert an absolute gas saturation to an effective one.
Definition: 3p/efftoabslaw.hh:312
static Scalar sneToSn_(const Params &params, Scalar sne)
Definition: 3p/efftoabslaw.hh:332
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 sw(const Params &params, const Scalar pc)
The saturation-capillary pressure curve.
Definition: 3p/efftoabslaw.hh:134
static Scalar bulkDensTimesAdsorpCoeff(const Params &params)
the basis for calculating adsorbed NAPL in storage term
Definition: 3p/efftoabslaw.hh:256
static Scalar dsw_dswe_(const Params &params)
Derivative of the absolute saturation w.r.t. the effective saturation.
Definition: 3p/efftoabslaw.hh:362
static Scalar swToSwe(const Params &params, const Scalar sw)
Convert an absolute wetting saturation to an effective one.
Definition: 3p/efftoabslaw.hh:270
static Scalar sgeToSg_(const Params &params, Scalar sge)
Definition: 3p/efftoabslaw.hh:337
static Scalar sweToSw_(const Params &params, Scalar swe)
Convert an effective wetting saturation to an absolute one.
Definition: 3p/efftoabslaw.hh:327
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:246
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 stToSte(const Params &params, const Scalar st)
Convert an absolute total liquid saturation to an effective one.
Definition: 3p/efftoabslaw.hh:298
static Scalar pc(const Params &params, const Scalar sw)
The capillary pressure-saturation curve.
Definition: 3p/efftoabslaw.hh:78
A default implementation of the parameters for the adapter class to convert material laws from effect...