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