3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
regularizedparkervangen3pparams.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 *****************************************************************************/
25#ifndef DUMUX_REGULARIZED_PARKERVANGEN_3P_PARAMS_HH
26#define DUMUX_REGULARIZED_PARKERVANGEN_3P_PARAMS_HH
27
29
30namespace Dumux {
31
37template <class ScalarT>
39{
41
42public:
43 using Scalar = ScalarT;
44
46 : ParkerVanGen3PParams(), constRegularization_(false)
47 {
48 pcLowS_ = 1e-2;
49 pcHighS_ = 99e-2;
50 }
51
53 Dune::FieldVector<Scalar, 4> residualSaturation, Scalar betaNw = 1.,
54 Scalar betaGn = 1., Scalar betaGw = 1., bool regardSnr=false)
56 residualSaturation, betaNw,
57 betaGn , betaGw, regardSnr), constRegularization_(false)
58 {
59 pcLowS_ = 1e-2;
60 pcHighS_ = 99e-2;
61 }
62
70 Scalar pcLowS() const
71 {
72 // Most problems are very sensitive to this value
73 // (e.g. making it smaller might result in negative
74 // pressures)
75 //
76 // If you want to use a different regularization threshold,
77 // overload this class and supply the new class as second
78 // template parameter for the RegularizedVanGenuchten law!
79 return pcLowS_;
80 }
81
90 {
91 // Most problems are very sensitive to this value
92 // (e.g. making it smaller might result in negative
93 // pressures)
94 //
95 // If you want to use a different regularization threshold,
96 // overload this class and supply the new class as second
97 // template parameter for the RegularizedVanGenuchten law!
98 return pcHighS_;
99 }
100
105 void setPcLowS(const Scalar input)
106 { pcLowS_ = input; }
107
112 void setPcHighS(const Scalar input)
113 { pcHighS_ = input; }
114
120 void useConstRegularization(const bool input)
121 {
122 constRegularization_ = input;
123 }
124
130 {
131 return constRegularization_;
132 }
133
134private:
135 Scalar pcLowS_;
136 Scalar pcHighS_;
137 bool constRegularization_;
138
139};
140} // namespace Dumux
141
142#endif
Specification of the material params for the van Genuchten capillary pressure model.
Definition: adapt.hh:29
Reference implementation of a van Genuchten params.
Definition: parkervangen3pparams.hh:42
ScalarT Scalar
Definition: parkervangen3pparams.hh:44
Scalar rhoBulk() const
Return the bulk density of the porous medium in .
Definition: parkervangen3pparams.hh:245
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
Scalar betaGw() const
Definition: parkervangen3pparams.hh:225
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
Parameters that are necessary for the regularization of the Parker - van Genuchten capillary pressure...
Definition: regularizedparkervangen3pparams.hh:39
void setPcHighS(const Scalar input)
Set the upper saturation threshold value.
Definition: regularizedparkervangen3pparams.hh:112
RegularizedParkerVanGen3PParams()
Definition: regularizedparkervangen3pparams.hh:45
Scalar pcLowS() const
Threshold saturation below which the capillary pressure is regularized.
Definition: regularizedparkervangen3pparams.hh:70
Scalar pcHighS() const
Threshold saturation above which the capillary pressure is regularized.
Definition: regularizedparkervangen3pparams.hh:89
ScalarT Scalar
Definition: regularizedparkervangen3pparams.hh:43
void useConstRegularization(const bool input)
Choose whether to use a constant value for regularization of the pc-S curves or not.
Definition: regularizedparkervangen3pparams.hh:120
RegularizedParkerVanGen3PParams(Scalar vgAlpha, Scalar vgn, Scalar KdNAPL, Scalar rhoBulk, Dune::FieldVector< Scalar, 4 > residualSaturation, Scalar betaNw=1., Scalar betaGn=1., Scalar betaGw=1., bool regardSnr=false)
Definition: regularizedparkervangen3pparams.hh:52
void setPcLowS(const Scalar input)
Set the lower saturation threshold value.
Definition: regularizedparkervangen3pparams.hh:105
bool constRegularization() const
Returns whether to use a constant value for regularization of the pc-S curves or not.
Definition: regularizedparkervangen3pparams.hh:129