3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
parkervangen3pparams.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 *****************************************************************************/
28#ifndef PARKERVANGEN_PARAMS_3P_HH
29#define PARKERVANGEN_PARAMS_3P_HH
30
31#include <dune/common/fvector.hh>
32#include <iostream>
33
34namespace Dumux {
35
40template<class ScalarT>
42{
43public:
44 using Scalar = ScalarT;
45
47 {betaGw_ = betaNw_ = betaGn_ = 1.0;}
48
50 Dune::FieldVector<Scalar, 4> residualSaturation, Scalar betaNw = 1.0,
51 Scalar betaGn = 1.0, Scalar betaGw = 1.0, bool regardSnr=false)
52 {
54 setVgn(vgn);
55 setSwr(residualSaturation[0]);
56 setSnr(residualSaturation[1]);
57 setSgr(residualSaturation[2]);
58 setKrRegardsSnr(regardSnr);
64 }
65
71 { return vgAlpha_; }
72
79 { vgAlpha_ = v; }
80
85 Scalar vgm() const
86 { return vgm_; }
87
95 void setVgm(Scalar m)
96 { vgm_ = m; vgn_ = 1/(1 - vgm_); }
97
102 Scalar vgn() const
103 { return vgn_; }
104
113 { vgn_ = n; vgm_ = 1 - 1/vgn_; }
114
119 Scalar satResidual(int phaseIdx) const
120 {
121 switch (phaseIdx)
122 {
123 case 0:
124 return swr_;
125 case 1:
126 return snr_;
127 case 2:
128 DUNE_THROW(Dune::NotImplemented, "sgr for three phases not required and therefore not implemented");
129 }
130 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
131 }
132
136 void setResiduals(Dune::FieldVector<Scalar, 3> residualSaturation)
137 {
138 setSwr(residualSaturation[0]);
139 setSnr(residualSaturation[1]);
140 setSgr(residualSaturation[2]);
141 }
142
143
147 Scalar swr() const
148 { return swr_; }
149
154 void setSwr(Scalar input)
155 { swr_ = input; }
156
160 Scalar snr() const
161 { return snr_; }
162
167 void setSnr(Scalar input)
168 { snr_ = input; }
169
173 Scalar sgr() const
174 {
175 return sgr_;
176 }
177
182 void setSgr(Scalar input)
183 {
184 sgr_ = input;
185 }
186
190 Scalar swrx() const
191 {
192 std::cerr << "swrx for three phases not implemented anymore. Equals swr" << std::endl;
193 return swr_;
194 }
195
201 {
202 std::cerr << "swrx for three phases not implemented anymore. Equals swr" << std::endl;
203 }
207 void setBetaNw(Scalar input)
208 { betaNw_ = input; }
209
210 void setBetaGn(Scalar input)
211 { betaGn_ = input; }
212
213 void setBetaGw(Scalar input)
214 { betaGw_ = input; }
215
220 { return betaNw_; }
221
223 { return betaGn_; }
224
226 { return betaGw_; }
227
232 void setKrRegardsSnr(bool input)
233 { krRegardsSnr_ = input; }
234
238 bool krRegardsSnr() const
239 { return krRegardsSnr_; }
240
241
246 { return rhoBulk_; }
247
252 void setRhoBulk(Scalar input)
253 { rhoBulk_ = input; }
254
259 { return KdNAPL_; }
260
265 void setKdNAPL(Scalar input)
266 { KdNAPL_ = input; }
267
268
269private:
270 Scalar vgAlpha_;
271 Scalar vgm_;
272 Scalar vgn_;
273 Scalar swr_;
274 Scalar snr_;
275 Scalar sgr_;
276
277 Scalar KdNAPL_;
278 Scalar rhoBulk_;
279
280 Scalar betaNw_;
281 Scalar betaGn_;
282 Scalar betaGw_;
283
284 bool krRegardsSnr_ ;
285};
286} // namespace Dumux
287
288#endif
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Reference implementation of a van Genuchten params.
Definition: parkervangen3pparams.hh:42
ScalarT Scalar
Definition: parkervangen3pparams.hh:44
Scalar vgm() const
Return the shape parameter of van Genuchten's curve.
Definition: parkervangen3pparams.hh:85
Scalar satResidual(int phaseIdx) const
Return the residual saturation.
Definition: parkervangen3pparams.hh:119
void setKdNAPL(Scalar input)
Set the adsorption coefficient.
Definition: parkervangen3pparams.hh:265
Scalar snr() const
Return the residual non-wetting saturation.
Definition: parkervangen3pparams.hh:160
Scalar rhoBulk() const
Return the bulk density of the porous medium in .
Definition: parkervangen3pparams.hh:245
void setResiduals(Dune::FieldVector< Scalar, 3 > residualSaturation)
Set all residual saturations.
Definition: parkervangen3pparams.hh:136
Scalar sgr() const
Return the residual gas saturation.
Definition: parkervangen3pparams.hh:173
bool krRegardsSnr() const
Calls if residual n-phase saturation should be regarded in its relative permeability.
Definition: parkervangen3pparams.hh:238
void setBetaGn(Scalar input)
Definition: parkervangen3pparams.hh:210
void setVgAlpha(Scalar v)
Set the shape parameter of van Genuchten's curve.
Definition: parkervangen3pparams.hh:78
void setSnr(Scalar input)
Set the residual non-wetting saturation.
Definition: parkervangen3pparams.hh:167
void setRhoBulk(Scalar input)
Set the bulk density of the porous medium.
Definition: parkervangen3pparams.hh:252
void setVgn(Scalar n)
Set the shape parameter of van Genuchten's curve.
Definition: parkervangen3pparams.hh:112
Scalar betaGn() const
Definition: parkervangen3pparams.hh:222
Scalar KdNAPL() const
Return the adsorption coefficient.
Definition: parkervangen3pparams.hh:258
Scalar betaNw() const
Return the values for the beta scaling parameters of capillary pressure between the phases.
Definition: parkervangen3pparams.hh:219
void setVgm(Scalar m)
Set the shape parameter of van Genuchten's curve.
Definition: parkervangen3pparams.hh:95
void setSwr(Scalar input)
Set the residual wetting saturation.
Definition: parkervangen3pparams.hh:154
void setKrRegardsSnr(bool input)
defines if residual n-phase saturation should be regarded in its relative permeability.
Definition: parkervangen3pparams.hh:232
void setSwrx(Scalar v)
Set the residual total liquid saturation.
Definition: parkervangen3pparams.hh:200
void setSgr(Scalar input)
Set the residual gas saturation.
Definition: parkervangen3pparams.hh:182
Scalar swr() const
Return the residual wetting saturation.
Definition: parkervangen3pparams.hh:147
Scalar betaGw() const
Definition: parkervangen3pparams.hh:225
Scalar swrx() const
Set the residual total liquid saturation.
Definition: parkervangen3pparams.hh:190
Scalar vgAlpha() const
Return the shape parameter of van Genuchten's curve.
Definition: parkervangen3pparams.hh:70
Scalar vgn() const
Return the shape parameter of van Genuchten's curve.
Definition: parkervangen3pparams.hh:102
ParkerVanGen3PParams()
Definition: parkervangen3pparams.hh:46
void setBetaGw(Scalar input)
Definition: parkervangen3pparams.hh:213
void setBetaNw(Scalar input)
defines the scaling parameters of capillary pressure between the phases (=1 for Gas-Water)
Definition: parkervangen3pparams.hh:207
ParkerVanGen3PParams(Scalar vgAlpha, Scalar vgn, Scalar KdNAPL, Scalar rhoBulk, Dune::FieldVector< Scalar, 4 > residualSaturation, Scalar betaNw=1.0, Scalar betaGn=1.0, Scalar betaGw=1.0, bool regardSnr=false)
Definition: parkervangen3pparams.hh:49