RegularizedVanGenuchten.hpp
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 This file is part of the Open Porous Media project (OPM).
5
6 OPM 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 2 of the License, or
9 (at your option) any later version.
10
11 OPM 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 OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef OPM_REGULARIZED_VAN_GENUCHTEN_HPP
28#define OPM_REGULARIZED_VAN_GENUCHTEN_HPP
29
30#include "VanGenuchten.hpp"
32
34
35#include <algorithm>
36
37namespace Opm {
38
70template <class TraitsT, class ParamsT = RegularizedVanGenuchtenParams<TraitsT> >
71class RegularizedVanGenuchten : public TraitsT
72{
73 typedef ::Opm::VanGenuchten<TraitsT, ParamsT> VanGenuchten;
74
75public:
76 typedef TraitsT Traits;
77 typedef ParamsT Params;
78
79 typedef typename Traits::Scalar Scalar;
80
82 static const int numPhases = Traits::numPhases;
83 static_assert(numPhases == 2,
84 "The regularized van Genuchten capillary pressure law only "
85 "applies to the case of two fluid phases");
86
89 static const bool implementsTwoPhaseApi = true;
90
93 static const bool implementsTwoPhaseSatApi = true;
94
97 static const bool isSaturationDependent = true;
98
101 static const bool isPressureDependent = false;
102
105 static const bool isTemperatureDependent = false;
106
109 static const bool isCompositionDependent = false;
110
115 template <class Container, class FluidState>
116 static void capillaryPressures(Container& values, const Params& params, const FluidState& fs)
117 {
118 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
119
120 values[Traits::wettingPhaseIdx] = 0.0; // reference phase
121 values[Traits::nonWettingPhaseIdx] = pcnw<FluidState, Evaluation>(params, fs);
122 }
123
128 template <class Container, class FluidState>
129 static void saturations(Container& values, const Params& params, const FluidState& fs)
130 {
131 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
132
133 values[Traits::wettingPhaseIdx] = Sw<FluidState, Evaluation>(params, fs);
134 values[Traits::nonWettingPhaseIdx] = 1 - values[Traits::wettingPhaseIdx];
135 }
136
141 template <class Container, class FluidState>
142 static void relativePermeabilities(Container& values, const Params& params, const FluidState& fs)
143 {
144 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
145
146 values[Traits::wettingPhaseIdx] = krw<FluidState, Evaluation>(params, fs);
147 values[Traits::nonWettingPhaseIdx] = krn<FluidState, Evaluation>(params, fs);
148 }
149
162 template <class FluidState, class Evaluation = typename FluidState::Scalar>
163 static Evaluation pcnw(const Params& params, const FluidState& fs)
164 {
165 const auto& Sw = decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
166 return twoPhaseSatPcnw(params, Sw);
167 }
168
169 template <class Evaluation>
170 static Evaluation twoPhaseSatPcnw(const Params& params, const Evaluation& Sw)
171 {
172 // retrieve the low and the high threshold saturations for the
173 // unregularized capillary pressure curve from the parameters
174 const Scalar SwThLow = params.pcnwLowSw();
175 const Scalar SwThHigh = params.pcnwHighSw();
176
177 // make sure that the capillary pressure observes a derivative
178 // != 0 for 'illegal' saturations. This is favourable for the
179 // newton solver (if the derivative is calculated numerically)
180 // in order to get the saturation moving to the right
181 // direction if it temporarily is in an 'illegal' range.
182 if (Sw < SwThLow) {
183 return params.pcnwLow() + params.pcnwSlopeLow()*(Sw - SwThLow);
184 }
185 else if (Sw > SwThHigh)
186 {
187 Scalar yTh = params.pcnwHigh();
188 Scalar m1 = (0.0 - yTh)/(1.0 - SwThHigh)*2;
189
190 if (Sw < 1.0) {
191 // use spline between threshold Sw and 1.0
192 const Spline<Scalar>& sp = params.pcnwHighSpline();
193
194 return sp.eval(Sw);
195 }
196 else {
197 // straight line for Sw > 1.0
198 return m1*(Sw - 1.0) + 0.0;
199 }
200 }
201
202 // if the effective saturation is in an 'reasonable'
203 // range, we use the real van genuchten law...
204 return VanGenuchten::twoPhaseSatPcnw(params, Sw);
205 }
206
221 template <class FluidState, class Evaluation = typename FluidState::Scalar>
222 static Evaluation Sw(const Params& params, const FluidState& fs)
223 {
224 const Evaluation& pC =
225 decay<Evaluation>(fs.pressure(Traits::nonWettingPhaseIdx))
226 - decay<Evaluation>(fs.pressure(Traits::wettingPhaseIdx));
227 return twoPhaseSatSw(params, pC);
228 }
229
230 template <class Evaluation>
231 static Evaluation twoPhaseSatSw(const Params& params, const Evaluation& pC)
232 {
233 // retrieve the low and the high threshold saturations for the
234 // unregularized capillary pressure curve from the parameters
235 const Scalar SwThLow = params.pcnwLowSw();
236 const Scalar SwThHigh = params.pcnwHighSw();
237
238 // calculate the saturation which corrosponds to the
239 // saturation in the non-regularized verision of van
240 // Genuchten's law
241 if (pC <= 0) {
242 // invert straight line for Sw > 1.0
243 Scalar m1 = params.pcnwSlopeHigh();
244 return pC/m1 + 1.0;
245 }
246
247 Evaluation Sw = VanGenuchten::twoPhaseSatSw(params, pC);
248
249 // invert the regularization if necessary
250 if (Sw <= SwThLow) {
251 // invert the low saturation regularization of pC()
252 Scalar pC_SwLow = VanGenuchten::twoPhaseSatPcnw(params, SwThLow);
253 return (pC - pC_SwLow)/params.pcnwSlopeLow() + SwThLow;
254 }
255 else if (SwThHigh < Sw /* && Sw < 1.0*/)
256 {
257 // invert spline between threshold saturation and 1.0
258 const Spline<Scalar>& spline = params.pcnwHighSpline();
259
260 return spline.template intersectInterval<Evaluation>(/*x0=*/SwThHigh, /*x1=*/1.0,
261 /*a=*/0, /*b=*/0, /*c=*/0, /*d=*/pC);
262 }
263
264 return Sw;
265 }
266
271 template <class FluidState, class Evaluation = typename FluidState::Scalar>
272 static Evaluation Sn(const Params& params, const FluidState& fs)
273 { return 1 - Sw<FluidState, Evaluation>(params, fs); }
274
275 template <class Evaluation>
276 static Evaluation twoPhaseSatSn(const Params& params, const Evaluation& pc)
277 { return 1 - twoPhaseSatSw(params, pc); }
278
293 template <class FluidState, class Evaluation = typename FluidState::Scalar>
294 static Evaluation krw(const Params& params, const FluidState& fs)
295 {
296 const auto& Sw = decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
297 return twoPhaseSatKrw(params, Sw);
298 }
299
300 template <class Evaluation>
301 static Evaluation twoPhaseSatKrw(const Params& params, const Evaluation& Sw)
302 {
303 // regularize
304 if (Sw <= 0)
305 return 0.0;
306 else if (Sw >= 1)
307 return 1.0;
308
309 return VanGenuchten::twoPhaseSatKrw(params, Sw);
310 }
311
326 template <class FluidState, class Evaluation = typename FluidState::Scalar>
327 static Evaluation krn(const Params& params, const FluidState& fs)
328 {
329 const Evaluation& Sw =
330 1.0 - decay<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
331 return twoPhaseSatKrn(params, Sw);
332 }
333
334 template <class Evaluation>
335 static Evaluation twoPhaseSatKrn(const Params& params, const Evaluation& Sw)
336 {
337 // regularize
338 if (Sw <= 0)
339 return 1.0;
340 else if (Sw >= 1)
341 return 0.0;
342
343 return VanGenuchten::twoPhaseSatKrn(params, Sw);
344 }
345};
346
347} // namespace Opm
348
349#endif
Implementation of the regularized van Genuchten's capillary pressure / relative permeability <-> satu...
Definition: RegularizedVanGenuchten.hpp:72
TraitsT Traits
Definition: RegularizedVanGenuchten.hpp:76
static const int numPhases
The number of fluid phases.
Definition: RegularizedVanGenuchten.hpp:82
static const bool isCompositionDependent
Definition: RegularizedVanGenuchten.hpp:109
static Evaluation krn(const Params &params, const FluidState &fs)
Regularized version of the relative permeability for the non-wetting phase of the medium implied by t...
Definition: RegularizedVanGenuchten.hpp:327
static const bool implementsTwoPhaseApi
Definition: RegularizedVanGenuchten.hpp:89
static void saturations(Container &values, const Params &params, const FluidState &fs)
Calculate the saturations of the phases starting from their pressure differences.
Definition: RegularizedVanGenuchten.hpp:129
static const bool isSaturationDependent
Definition: RegularizedVanGenuchten.hpp:97
static Evaluation Sn(const Params &params, const FluidState &fs)
Calculate the non-wetting phase saturations depending on the phase pressures.
Definition: RegularizedVanGenuchten.hpp:272
static const bool isPressureDependent
Definition: RegularizedVanGenuchten.hpp:101
static Evaluation pcnw(const Params &params, const FluidState &fs)
A regularized van Genuchten capillary pressure-saturation curve.
Definition: RegularizedVanGenuchten.hpp:163
static Evaluation twoPhaseSatPcnw(const Params &params, const Evaluation &Sw)
Definition: RegularizedVanGenuchten.hpp:170
ParamsT Params
Definition: RegularizedVanGenuchten.hpp:77
Traits::Scalar Scalar
Definition: RegularizedVanGenuchten.hpp:79
static const bool isTemperatureDependent
Definition: RegularizedVanGenuchten.hpp:105
static Evaluation twoPhaseSatSn(const Params &params, const Evaluation &pc)
Definition: RegularizedVanGenuchten.hpp:276
static Evaluation twoPhaseSatKrn(const Params &params, const Evaluation &Sw)
Definition: RegularizedVanGenuchten.hpp:335
static Evaluation twoPhaseSatSw(const Params &params, const Evaluation &pC)
Definition: RegularizedVanGenuchten.hpp:231
static Evaluation Sw(const Params &params, const FluidState &fs)
A regularized van Genuchten saturation-capillary pressure curve.
Definition: RegularizedVanGenuchten.hpp:222
static Evaluation twoPhaseSatKrw(const Params &params, const Evaluation &Sw)
Definition: RegularizedVanGenuchten.hpp:301
static void relativePermeabilities(Container &values, const Params &params, const FluidState &fs)
Returns the relative permeabilities of the phases dependening on the phase saturations.
Definition: RegularizedVanGenuchten.hpp:142
static const bool implementsTwoPhaseSatApi
Definition: RegularizedVanGenuchten.hpp:93
static Evaluation krw(const Params &params, const FluidState &fs)
Regularized version of the relative permeability for the wetting phase of the medium implied by the v...
Definition: RegularizedVanGenuchten.hpp:294
static void capillaryPressures(Container &values, const Params &params, const FluidState &fs)
Calculate the pressure difference of the phases in the most generic way.
Definition: RegularizedVanGenuchten.hpp:116
Class implementing cubic splines.
Definition: Spline.hpp:91
Evaluation eval(const Evaluation &x, bool extrapolate=false) const
Evaluate the spline at a given position.
Definition: Spline.hpp:797
Implementation of the van Genuchten capillary pressure - saturation relation.
Definition: VanGenuchten.hpp:56
static Evaluation twoPhaseSatPcnw(const Params &params, const Evaluation &Sw)
The saturation-capillary pressure curve according to van Genuchten using a material law specific API.
Definition: VanGenuchten.hpp:194
static Evaluation twoPhaseSatKrn(const Params &params, Evaluation Sw)
Definition: VanGenuchten.hpp:287
static Evaluation twoPhaseSatKrw(const Params &params, const Evaluation &Sw)
Definition: VanGenuchten.hpp:260
static Evaluation twoPhaseSatSw(const Params &params, const Evaluation &pC)
Definition: VanGenuchten.hpp:221
Definition: Air_Mesitylene.hpp:34