My Project
Loading...
Searching...
No Matches
TwoPhaseLETCurves.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 3 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_TWO_PHASE_LET_CURVES_HPP
28#define OPM_TWO_PHASE_LET_CURVES_HPP
29
31
33
35
36namespace Opm {
48template <class TraitsT, class ParamsT = TwoPhaseLETCurvesParams<TraitsT> >
49class TwoPhaseLETCurves : public TraitsT
50{
51public:
52 using Traits = TraitsT;
53 using Params = ParamsT;
54 using Scalar = typename Traits::Scalar;
55
56 static_assert(Traits::numPhases == 2,
57 "The number of fluid phases must be two if you want to use "
58 "this material law!");
59
60 static constexpr Scalar eps = 1.0e-10; //tolerance
61
63 static constexpr int numPhases = Traits::numPhases;
64
67 static constexpr bool implementsTwoPhaseApi = true;
68
71 static constexpr bool implementsTwoPhaseSatApi = true;
72
75 static constexpr bool isSaturationDependent = true;
76
79 static constexpr bool isPressureDependent = false;
80
83 static constexpr bool isTemperatureDependent = false;
84
87 static constexpr bool isCompositionDependent = false;
88
92 template <class Container, class FluidState>
93 static void capillaryPressures(Container& /* values */, const Params& /* params */, const FluidState& /* fs */)
94 {
95 throw std::logic_error("The capillaryPressures(fs) method is not yet implemented");
96 }
97
102 template <class Container, class FluidState>
103 static void saturations(Container& /* pc */, const Params& /* params */, const FluidState& /* fs */)
104 {
105 throw std::logic_error("The saturations(fs) method is not yet implemented");
106 }
107
118 template <class Container, class FluidState>
119 static void relativePermeabilities(Container& /* pc */, const Params& /* params */, const FluidState& /* fs */)
120 {
121 throw std::logic_error("The relativePermeabilities(fs) method is not yet implemented");
122 }
123
129 template <class FluidState, class Evaluation = typename FluidState::Scalar>
130 static Evaluation pcnw(const Params& /* params */, const FluidState& /* fs */)
131 {
132 throw std::logic_error("TwoPhaseLETCurves::pcnw"
133 " not implemented!");
134
135 }
136
137 template <class Evaluation>
138 static Evaluation twoPhaseSatPcnw(const Params& params, const Evaluation& Sw)
139 {
140 Evaluation Ss = (Sw-params.Sminpc())/params.dSpc();
141 if (Ss < 0.0) {
142 Ss -= (Opm::decay<Scalar>(Ss));
143 } else if (Ss > 1.0) {
144 Ss -= (Opm::decay<Scalar>(Ss-1.0));
145 }
146
147 const Evaluation powS = Opm::pow(Ss,params.Tpc());
148 const Evaluation pow1mS = Opm::pow(1.0-Ss,params.Lpc());
149
150 const Evaluation F = pow1mS/(pow1mS+powS*params.Epc());
151 Evaluation tmp = params.Pct()+(params.Pcir()-params.Pct())*F;
152
153 return tmp;
154 }
155
156 template <class Evaluation>
157 static Evaluation twoPhaseSatPcnwInv(const Params& /* params */, const Evaluation&)
158 {
159 throw std::logic_error("TwoPhaseLETCurves::twoPhaseSatPcnwInv"
160 " not implemented!");
161 }
162
163 template <class FluidState, class Evaluation = typename FluidState::Scalar>
164 static Evaluation Sw(const Params& /* params */, const FluidState& /* fs */)
165 {
166 throw std::logic_error("The Sw(fs) method is not yet implemented");
167 }
168
169 template <class Evaluation>
170 static Evaluation twoPhaseSatSw(const Params& /* params */, const Evaluation& /* pc */)
171 {
172 throw std::logic_error("The twoPhaseSatSw(fs) method is not yet implemented");
173 }
174
175 template <class FluidState, class Evaluation = typename FluidState::Scalar>
176 static Evaluation Sn(const Params& /* params */, const FluidState& /* fs */)
177 {
178 throw std::logic_error("The Sn(fs) method is not yet implemented");
179 }
180
181 template <class Evaluation>
182 static Evaluation twoPhaseSatSn(const Params& /* params */, const Evaluation& /* pc */)
183 {
184 throw std::logic_error("The twoPhaseSatSn(fs) method is not yet implemented");
185 }
192 template <class FluidState, class Evaluation = typename FluidState::Scalar>
193 static Evaluation krw(const Params& /* params */, const FluidState& /* fs */)
194 {
195 throw std::logic_error("TwoPhaseLETCurves::krw"
196 " not implemented!");
197 }
198
199 template <class Evaluation>
200 static Evaluation twoPhaseSatKrw(const Params& params, const Evaluation& Sw)
201 {
202 return twoPhaseSatKrLET(Params::wIdx, params, Sw);
203 }
204
205 template <class Evaluation>
206 static Evaluation twoPhaseSatKrLET(const unsigned phaseIdx, const Params& params, const Evaluation& S)
207 {
208 Evaluation Ss = (S-params.Smin(phaseIdx))/params.dS(phaseIdx);
209 if (Ss < 0.0) {
210 Ss -= (Opm::decay<Scalar>(Ss));
211 } else if (Ss > 1.0) {
212 Ss -= (Opm::decay<Scalar>(Ss-1.0));
213 }
214
215 const Evaluation powS = Opm::pow(Ss,params.L(phaseIdx));
216 const Evaluation pow1mS = Opm::pow(1.0-Ss,params.T(phaseIdx));
217
218 const Evaluation tmp = params.Krt(phaseIdx)*powS/(powS+pow1mS*params.E(phaseIdx));
219
220 return tmp;
221 }
222
223 template <class Evaluation>
224 static Evaluation twoPhaseSatKrwInv(const Params& /* params */, const Evaluation& /* krw */)
225 {
226 throw std::logic_error("TwoPhaseLETCurves::twoPhaseSatKrwInv"
227 " not implemented!");
228 }
229
236 template <class FluidState, class Evaluation = typename FluidState::Scalar>
237 static Evaluation krn(const Params& /* params */, const FluidState& /* fs */)
238 {
239 throw std::logic_error("TwoPhaseLETCurves::krn"
240 " not implemented!");
241 }
242
243 template <class Evaluation>
244 static Evaluation twoPhaseSatKrn(const Params& params, const Evaluation& Sw)
245 {
246 const Evaluation Sn = 1.0 - Sw;
247
248 return twoPhaseSatKrLET(Params::nwIdx, params, Sn);
249 }
250
251 template <class Evaluation>
252 static Evaluation twoPhaseSatKrnInv(const Params& params, const Evaluation& krn)
253 {
254 // since inverting the formula for krn is hard to do analytically, we use the
255 // Newton-Raphson method
256 Evaluation Sw = 0.5;
257 //Scalar eps = 1e-10;
258 for (int i = 0; i < 20; ++i) {
259 Evaluation f = krn - twoPhaseSatKrn(params, Sw);
260 if (Opm::abs(f) < 1e-10)
261 return Sw;
262 Evaluation fStar = krn - twoPhaseSatKrn(params, Sw + eps);
263 Evaluation fPrime = (fStar - f)/eps;
264 Evaluation delta = f/fPrime;
265
266 Sw -= delta;
267 if (Sw < 0)
268 Sw = 0.0;
269 if (Sw > 1.0)
270 Sw = 1.0;
271 if (Opm::abs(delta) < 1e-10)
272 return Sw;
273 }
274
275 // Fallback to simple bisection
276 Evaluation SL = 0.0;
277 Evaluation fL = krn - twoPhaseSatKrn(params, SL);
278 if (Opm::abs(fL) < eps)
279 return SL;
280 Evaluation SR = 1.0;
281 Evaluation fR = krn - twoPhaseSatKrn(params, SR);
282 if (Opm::abs(fR) < eps)
283 return SR;
284 if (fL*fR < 0.0) {
285 for (int i = 0; i < 50; ++i) {
286 Sw = 0.5*(SL+SR);
287 if (abs(SR-SL) < eps)
288 return Sw;
289 Evaluation fw = krn - twoPhaseSatKrn(params, Sw);
290 if (Opm::abs(fw) < eps)
291 return Sw;
292 if (fw * fR > 0) {
293 SR = Sw;
294 fR = fw;
295 } else if (fw * fL > 0) {
296 SL = Sw;
297 fL = fw;
298 }
299 }
300
301 }
302
303 throw NumericalProblem("Couldn't invert the TwoPhaseLETCurves non-wetting phase"
304 " relperm within 20 newton iterations and 50 bisection iterations");
305 }
306};
307
308} // namespace Opm
309
310#endif // OPM_TWO_PHASE_LET_CURVES_HPP
Provides the OPM specific exception classes.
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
Implementation of the LET curve saturation functions.
Definition TwoPhaseLETCurves.hpp:50
static Evaluation krn(const Params &, const FluidState &)
The relative permeability for the non-wetting phase of the medium as implied by the LET parameterizat...
Definition TwoPhaseLETCurves.hpp:237
static Evaluation pcnw(const Params &, const FluidState &)
The capillary pressure-saturation curve.
Definition TwoPhaseLETCurves.hpp:130
static constexpr bool implementsTwoPhaseSatApi
Specify whether this material law implements the two-phase convenience API which only depends on the ...
Definition TwoPhaseLETCurves.hpp:71
static Evaluation krw(const Params &, const FluidState &)
The relative permeability for the wetting phase of the medium implied by the LET parameterization.
Definition TwoPhaseLETCurves.hpp:193
static constexpr bool isSaturationDependent
Specify whether the quantities defined by this material law are saturation dependent.
Definition TwoPhaseLETCurves.hpp:75
static constexpr int numPhases
The number of fluid phases to which this material law applies.
Definition TwoPhaseLETCurves.hpp:63
static constexpr bool isCompositionDependent
Specify whether the quantities defined by this material law are dependent on the phase composition.
Definition TwoPhaseLETCurves.hpp:87
static constexpr bool isTemperatureDependent
Specify whether the quantities defined by this material law are temperature dependent.
Definition TwoPhaseLETCurves.hpp:83
static void saturations(Container &, const Params &, const FluidState &)
Calculate the saturations of the phases starting from their pressure differences.
Definition TwoPhaseLETCurves.hpp:103
static constexpr bool isPressureDependent
Specify whether the quantities defined by this material law are dependent on the absolute pressure.
Definition TwoPhaseLETCurves.hpp:79
static void capillaryPressures(Container &, const Params &, const FluidState &)
The capillary pressure-saturation curves.
Definition TwoPhaseLETCurves.hpp:93
static constexpr bool implementsTwoPhaseApi
Specify whether this material law implements the two-phase convenience API.
Definition TwoPhaseLETCurves.hpp:67
static void relativePermeabilities(Container &, const Params &, const FluidState &)
The relative permeability-saturation curves.
Definition TwoPhaseLETCurves.hpp:119
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30