opm-common
Loading...
Searching...
No Matches
BrooksCorey.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_BROOKS_COREY_HPP
28#define OPM_BROOKS_COREY_HPP
29
30#include "BrooksCoreyParams.hpp"
31
33
35
36#include <algorithm>
37#include <cassert>
38#include <cmath>
39
40namespace Opm {
53template <class TraitsT, class ParamsT = BrooksCoreyParams<TraitsT> >
54class BrooksCorey : public TraitsT
55{
56public:
57 typedef TraitsT Traits;
58 typedef ParamsT Params;
59 typedef typename Traits::Scalar Scalar;
60
62 static const int numPhases = Traits::numPhases;
63 static_assert(numPhases == 2,
64 "The Brooks-Corey capillary pressure law only applies "
65 "to the case of two fluid phases");
66
69 static const bool implementsTwoPhaseApi = true;
70
73 static const bool implementsTwoPhaseSatApi = true;
74
77 static const bool isSaturationDependent = true;
78
81 static const bool isPressureDependent = false;
82
85 static const bool isTemperatureDependent = false;
86
89 static const bool isCompositionDependent = false;
90
91 static_assert(Traits::numPhases == 2,
92 "The number of fluid phases must be two if you want to use "
93 "this material law!");
94
98 template <class Container, class FluidState>
99 static void capillaryPressures(Container& values, const Params& params, const FluidState& fs)
100 {
101 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
102
103 values[Traits::wettingPhaseIdx] = 0.0; // reference phase
104 values[Traits::nonWettingPhaseIdx] = pcnw<FluidState, Evaluation>(params, fs);
105 }
106
111 template <class Container, class FluidState>
112 static void saturations(Container& values, const Params& params, const FluidState& fs)
113 {
114 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
115
116 values[Traits::wettingPhaseIdx] = Sw<FluidState, Evaluation>(params, fs);
117 values[Traits::nonWettingPhaseIdx] = 1.0 - values[Traits::wettingPhaseIdx];
118 }
119
130 template <class Container, class FluidState>
131 static void relativePermeabilities(Container& values, const Params& params, const FluidState& fs)
132 {
133 typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
134
135 values[Traits::wettingPhaseIdx] = krw<FluidState, Evaluation>(params, fs);
136 values[Traits::nonWettingPhaseIdx] = krn<FluidState, Evaluation>(params, fs);
137 }
138
153 template <class FluidState, class Evaluation = typename FluidState::ValueType>
154 static Evaluation pcnw(const Params& params, const FluidState& fs)
155 {
156 const Evaluation& Sw =
157 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
158
159 assert(0.0 <= Sw && Sw <= 1.0);
160
161 return twoPhaseSatPcnw(params, Sw);
162 }
163
164 template <class Evaluation>
165 static Evaluation twoPhaseSatPcnw(const Params& params, const Evaluation& Sw)
166 {
167 assert(0.0 <= Sw && Sw <= 1.0);
168
169 return params.entryPressure()*pow(Sw, -1/params.lambda());
170 }
171
172 template <class Evaluation>
173 static Evaluation twoPhaseSatPcnwInv(const Params& params, const Evaluation& pcnw)
174 {
175 assert(pcnw > 0.0);
176
177 return pow(params.entryPressure()/pcnw, -params.lambda());
178 }
179
193 template <class FluidState, class Evaluation = typename FluidState::ValueType>
194 static Evaluation Sw(const Params& params, const FluidState& fs)
195 {
196 Evaluation pC =
197 decay<Evaluation>(fs.pressure(Traits::nonWettingPhaseIdx))
198 - decay<Evaluation>(fs.pressure(Traits::wettingPhaseIdx));
199 return twoPhaseSatSw(params, pC);
200 }
201
202 template <class Evaluation>
203 static Evaluation twoPhaseSatSw(const Params& params, const Evaluation& pc)
204 {
205 assert(pc > 0.0); // if we don't assume that, std::pow will screw up!
206
207 return pow(pc/params.entryPressure(), -params.lambda());
208 }
209
214 template <class FluidState, class Evaluation = typename FluidState::ValueType>
215 static Evaluation Sn(const Params& params, const FluidState& fs)
216 { return 1.0 - Sw<FluidState, Evaluation>(params, fs); }
217
218 template <class Evaluation>
219 static Evaluation twoPhaseSatSn(const Params& params, const Evaluation& pc)
220 { return 1.0 - twoPhaseSatSw(params, pc); }
230 template <class FluidState, class Evaluation = typename FluidState::ValueType>
231 static Evaluation krw(const Params& params, const FluidState& fs)
232 {
233 const auto& sw =
234 decay<Evaluation>(fs.saturation(Traits::wettingPhaseIdx));
235
236 return twoPhaseSatKrw(params, sw);
237 }
238
239 template <class Evaluation>
240 static Evaluation twoPhaseSatKrw(const Params& params, const Evaluation& Sw)
241 {
242 assert(0.0 <= Sw && Sw <= 1.0);
243
244 return pow(Sw, 2.0/params.lambda() + 3.0);
245 }
246
247 template <class Evaluation>
248 static Evaluation twoPhaseSatKrwInv(const Params& params, const Evaluation& krw)
249 {
250 return pow(krw, 1.0/(2.0/params.lambda() + 3.0));
251 }
252
262 template <class FluidState, class Evaluation = typename FluidState::ValueType>
263 static Evaluation krn(const Params& params, const FluidState& fs)
264 {
265 const Evaluation& sw =
266 1.0 - decay<Evaluation>(fs.saturation(Traits::nonWettingPhaseIdx));
267
268 return twoPhaseSatKrn(params, sw);
269 }
270
271 template <class Evaluation>
272 static Evaluation twoPhaseSatKrn(const Params& params, const Evaluation& Sw)
273 {
274 assert(0.0 <= Sw && Sw <= 1.0);
275
276 Scalar exponent = 2.0/params.lambda() + 1.0;
277 const Evaluation sn = 1.0 - Sw;
278 return sn*sn*(1. - pow(Sw, exponent));
279 }
280
281 template <class Evaluation>
282 static Evaluation twoPhaseSatKrnInv(const Params& params, const Evaluation& krn)
283 {
284 // since inverting the formula for krn is hard to do analytically, we use the
285 // Newton-Raphson method
286 Evaluation sw = 0.5;
287 Scalar eps = 1e-10;
288 for (int i = 0; i < 20; ++i) {
289 Evaluation f = krn - twoPhaseSatKrn(params, sw);
290 Evaluation fStar = krn - twoPhaseSatKrn(params, sw + eps);
291 Evaluation fPrime = (fStar - f) / eps;
292
293 Evaluation delta = f / fPrime;
294 sw -= delta;
295 if (sw < 0) {
296 sw = 0.0;
297 }
298 if (abs(delta) < 1e-10) {
299 return sw;
300 }
301 }
302
303 throw NumericalProblem("Couldn't invert the Brooks-Corey non-wetting phase"
304 " relperm within 20 iterations");
305 }
306};
307} // namespace Opm
308
309#endif // BROOKS_COREY_HPP
Specification of the material parameters for the Brooks-Corey constitutive relations.
Provides the OPM specific exception classes.
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
Implementation of the Brooks-Corey capillary pressure <-> saturation relation.
Definition BrooksCorey.hpp:55
static void saturations(Container &values, const Params &params, const FluidState &fs)
Calculate the saturations of the phases starting from their pressure differences.
Definition BrooksCorey.hpp:112
static const bool isPressureDependent
Definition BrooksCorey.hpp:81
static const bool isSaturationDependent
Definition BrooksCorey.hpp:77
static const int numPhases
Definition BrooksCorey.hpp:62
static void capillaryPressures(Container &values, const Params &params, const FluidState &fs)
The capillary pressure-saturation curves.
Definition BrooksCorey.hpp:99
static void relativePermeabilities(Container &values, const Params &params, const FluidState &fs)
The relative permeability-saturation curves.
Definition BrooksCorey.hpp:131
static const bool implementsTwoPhaseApi
Definition BrooksCorey.hpp:69
static Evaluation Sn(const Params &params, const FluidState &fs)
Calculate the non-wetting phase saturations depending on the phase pressures.
Definition BrooksCorey.hpp:215
static const bool isTemperatureDependent
Definition BrooksCorey.hpp:85
static Evaluation Sw(const Params &params, const FluidState &fs)
The saturation-capillary pressure curve according to Brooks & Corey.
Definition BrooksCorey.hpp:194
static Evaluation krw(const Params &params, const FluidState &fs)
The relative permeability for the wetting phase of the medium implied by the Brooks-Corey parameteriz...
Definition BrooksCorey.hpp:231
static Evaluation krn(const Params &params, const FluidState &fs)
The relative permeability for the non-wetting phase of the medium as implied by the Brooks-Corey para...
Definition BrooksCorey.hpp:263
static const bool isCompositionDependent
Definition BrooksCorey.hpp:89
static Evaluation pcnw(const Params &params, const FluidState &fs)
The capillary pressure-saturation curve according to Brooks and Corey.
Definition BrooksCorey.hpp:154
static const bool implementsTwoPhaseSatApi
Definition BrooksCorey.hpp:73
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30