3.3.0
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#warning "This header is deprecated. Removal after 3.3. Use new material laws."
32
33#include <dune/common/fvector.hh>
34#include <iostream>
35
36namespace Dumux {
37
42template<class ScalarT>
44{
45public:
46 using Scalar = ScalarT;
47
49 {betaGw_ = betaNw_ = betaGn_ = 1.0;}
50
52 Dune::FieldVector<Scalar, 4> residualSaturation, Scalar betaNw = 1.0,
53 Scalar betaGn = 1.0, Scalar betaGw = 1.0, bool regardSnr=false)
54 {
56 setVgn(vgn);
57 setSwr(residualSaturation[0]);
58 setSnr(residualSaturation[1]);
59 setSgr(residualSaturation[2]);
60 setKrRegardsSnr(regardSnr);
66 }
67
73 { return vgAlpha_; }
74
81 { vgAlpha_ = v; }
82
87 Scalar vgm() const
88 { return vgm_; }
89
97 void setVgm(Scalar m)
98 { vgm_ = m; vgn_ = 1/(1 - vgm_); }
99
104 Scalar vgn() const
105 { return vgn_; }
106
115 { vgn_ = n; vgm_ = 1 - 1/vgn_; }
116
121 Scalar satResidual(int phaseIdx) const
122 {
123 switch (phaseIdx)
124 {
125 case 0:
126 return swr_;
127 case 1:
128 return snr_;
129 case 2:
130 DUNE_THROW(Dune::NotImplemented, "sgr for three phases not required and therefore not implemented");
131 }
132 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
133 }
134
138 void setResiduals(Dune::FieldVector<Scalar, 3> residualSaturation)
139 {
140 setSwr(residualSaturation[0]);
141 setSnr(residualSaturation[1]);
142 setSgr(residualSaturation[2]);
143 }
144
145
149 Scalar swr() const
150 { return swr_; }
151
156 void setSwr(Scalar input)
157 { swr_ = input; }
158
162 Scalar snr() const
163 { return snr_; }
164
169 void setSnr(Scalar input)
170 { snr_ = input; }
171
175 Scalar sgr() const
176 {
177 return sgr_;
178 }
179
184 void setSgr(Scalar input)
185 {
186 sgr_ = input;
187 }
188
192 Scalar swrx() const
193 {
194 std::cerr << "swrx for three phases not implemented anymore. Equals swr" << std::endl;
195 return swr_;
196 }
197
203 {
204 std::cerr << "swrx for three phases not implemented anymore. Equals swr" << std::endl;
205 }
209 void setBetaNw(Scalar input)
210 { betaNw_ = input; }
211
212 void setBetaGn(Scalar input)
213 { betaGn_ = input; }
214
215 void setBetaGw(Scalar input)
216 { betaGw_ = input; }
217
222 { return betaNw_; }
223
225 { return betaGn_; }
226
228 { return betaGw_; }
229
234 void setKrRegardsSnr(bool input)
235 { krRegardsSnr_ = input; }
236
240 bool krRegardsSnr() const
241 { return krRegardsSnr_; }
242
243
248 { return rhoBulk_; }
249
254 void setRhoBulk(Scalar input)
255 { rhoBulk_ = input; }
256
261 { return KdNAPL_; }
262
267 void setKdNAPL(Scalar input)
268 { KdNAPL_ = input; }
269
270
271private:
272 Scalar vgAlpha_;
273 Scalar vgm_;
274 Scalar vgn_;
275 Scalar swr_;
276 Scalar snr_;
277 Scalar sgr_;
278
279 Scalar KdNAPL_;
280 Scalar rhoBulk_;
281
282 Scalar betaNw_;
283 Scalar betaGn_;
284 Scalar betaGw_;
285
286 bool krRegardsSnr_ ;
287};
288} // namespace Dumux
289
290#endif
Definition: adapt.hh:29
Reference implementation of a van Genuchten params.
Definition: parkervangen3pparams.hh:44
ScalarT Scalar
Definition: parkervangen3pparams.hh:46
Scalar vgm() const
Return the shape parameter of van Genuchten's curve.
Definition: parkervangen3pparams.hh:87
Scalar satResidual(int phaseIdx) const
Return the residual saturation.
Definition: parkervangen3pparams.hh:121
void setKdNAPL(Scalar input)
Set the adsorption coefficient.
Definition: parkervangen3pparams.hh:267
Scalar snr() const
Return the residual nonwetting saturation.
Definition: parkervangen3pparams.hh:162
Scalar rhoBulk() const
Return the bulk density of the porous medium in .
Definition: parkervangen3pparams.hh:247
void setResiduals(Dune::FieldVector< Scalar, 3 > residualSaturation)
Set all residual saturations.
Definition: parkervangen3pparams.hh:138
Scalar sgr() const
Return the residual gas saturation.
Definition: parkervangen3pparams.hh:175
bool krRegardsSnr() const
Calls if residual n-phase saturation should be regarded in its relative permeability.
Definition: parkervangen3pparams.hh:240
void setBetaGn(Scalar input)
Definition: parkervangen3pparams.hh:212
void setVgAlpha(Scalar v)
Set the shape parameter of van Genuchten's curve.
Definition: parkervangen3pparams.hh:80
void setSnr(Scalar input)
Set the residual nonwetting saturation.
Definition: parkervangen3pparams.hh:169
void setRhoBulk(Scalar input)
Set the bulk density of the porous medium.
Definition: parkervangen3pparams.hh:254
void setVgn(Scalar n)
Set the shape parameter of van Genuchten's curve.
Definition: parkervangen3pparams.hh:114
Scalar betaGn() const
Definition: parkervangen3pparams.hh:224
Scalar KdNAPL() const
Return the adsorption coefficient.
Definition: parkervangen3pparams.hh:260
Scalar betaNw() const
Return the values for the beta scaling parameters of capillary pressure between the phases.
Definition: parkervangen3pparams.hh:221
void setVgm(Scalar m)
Set the shape parameter of van Genuchten's curve.
Definition: parkervangen3pparams.hh:97
void setSwr(Scalar input)
Set the residual wetting saturation.
Definition: parkervangen3pparams.hh:156
void setKrRegardsSnr(bool input)
defines if residual n-phase saturation should be regarded in its relative permeability.
Definition: parkervangen3pparams.hh:234
void setSwrx(Scalar v)
Set the residual total liquid saturation.
Definition: parkervangen3pparams.hh:202
void setSgr(Scalar input)
Set the residual gas saturation.
Definition: parkervangen3pparams.hh:184
Scalar swr() const
Return the residual wetting saturation.
Definition: parkervangen3pparams.hh:149
Scalar betaGw() const
Definition: parkervangen3pparams.hh:227
Scalar swrx() const
Set the residual total liquid saturation.
Definition: parkervangen3pparams.hh:192
Scalar vgAlpha() const
Return the shape parameter of van Genuchten's curve.
Definition: parkervangen3pparams.hh:72
Scalar vgn() const
Return the shape parameter of van Genuchten's curve.
Definition: parkervangen3pparams.hh:104
ParkerVanGen3PParams()
Definition: parkervangen3pparams.hh:48
void setBetaGw(Scalar input)
Definition: parkervangen3pparams.hh:215
void setBetaNw(Scalar input)
defines the scaling parameters of capillary pressure between the phases (=1 for Gas-Water)
Definition: parkervangen3pparams.hh:209
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:51