opm-common
Loading...
Searching...
No Matches
LiveOilPvt.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_LIVE_OIL_PVT_HPP
28#define OPM_LIVE_OIL_PVT_HPP
29
31#include <opm/common/OpmLog/OpmLog.hpp>
32
36
37namespace Opm {
38
39#if HAVE_ECL_INPUT
40class EclipseState;
41class Schedule;
42class SimpleTable;
43#endif
44
49template <class Scalar>
51{
52 using SamplingPoints = std::vector<std::pair<Scalar, Scalar>>;
53
54public:
55 using TabulatedTwoDFunction = UniformXTabulated2DFunction<Scalar>;
56 using TabulatedOneDFunction = Tabulated1DFunction<Scalar>;
57
58#if HAVE_ECL_INPUT
62 void initFromState(const EclipseState& eclState, const Schedule& schedule);
63
64private:
65 void extendPvtoTable_(unsigned regionIdx,
66 unsigned xIdx,
67 const SimpleTable& curTable,
68 const SimpleTable& masterTable);
69
70public:
71#endif // HAVE_ECL_INPUT
72
73 void setNumRegions(size_t numRegions);
74
75 void setVapPars(const Scalar, const Scalar par2)
76 {
77 vapPar2_ = par2;
78 }
79
83 void setReferenceDensities(unsigned regionIdx,
84 Scalar rhoRefOil,
85 Scalar rhoRefGas,
86 Scalar /*rhoRefWater*/);
87
94 void setSaturatedOilGasDissolutionFactor(unsigned regionIdx,
95 const SamplingPoints& samplePoints)
96 { saturatedGasDissolutionFactorTable_[regionIdx].setContainerOfTuples(samplePoints); }
97
107 void setSaturatedOilFormationVolumeFactor(unsigned regionIdx,
108 const SamplingPoints& samplePoints);
109
122 void setInverseOilFormationVolumeFactor(unsigned regionIdx,
123 const TabulatedTwoDFunction& invBo)
124 { inverseOilBTable_[regionIdx] = invBo; }
125
131 void setOilViscosity(unsigned regionIdx, const TabulatedTwoDFunction& muo)
132 { oilMuTable_[regionIdx] = muo; }
133
141 void setSaturatedOilViscosity(unsigned regionIdx,
142 const SamplingPoints& samplePoints);
143
147 void initEnd();
148
152 unsigned numRegions() const
153 { return inverseOilBMuTable_.size(); }
154
158 template <class Evaluation>
159 Evaluation internalEnergy(unsigned,
160 const Evaluation&,
161 const Evaluation&,
162 const Evaluation&) const
163 {
164 throw std::runtime_error("Requested the enthalpy of oil but the thermal "
165 "option is not enabled");
166 }
167
168 Scalar hVap(unsigned) const
169 {
170 throw std::runtime_error("Requested the hvap of oil but the thermal "
171 "option is not enabled");
172 }
173
177 template <class Evaluation>
178 Evaluation viscosity(unsigned regionIdx,
179 const Evaluation& /*temperature*/,
180 const Evaluation& pressure,
181 const Evaluation& Rs) const
182 {
183 // ATTENTION: Rs is the first axis!
184 const Evaluation& invBo = inverseOilBTable_[regionIdx].eval(Rs, pressure, /*extrapolate=*/true);
185 const Evaluation& invMuoBo = inverseOilBMuTable_[regionIdx].eval(Rs, pressure, /*extrapolate=*/true);
186
187 return invBo / invMuoBo;
188 }
189
193 template <class Evaluation>
194 Evaluation saturatedViscosity(unsigned regionIdx,
195 const Evaluation& /*temperature*/,
196 const Evaluation& pressure) const
197 {
198 // ATTENTION: Rs is the first axis!
199 const Evaluation& invBo = inverseSaturatedOilBTable_[regionIdx].eval(pressure, /*extrapolate=*/true);
200 const Evaluation& invMuoBo = inverseSaturatedOilBMuTable_[regionIdx].eval(pressure, /*extrapolate=*/true);
201
202 return invBo / invMuoBo;
203 }
204
208 template <class Evaluation>
209 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
210 const Evaluation& /*temperature*/,
211 const Evaluation& pressure,
212 const Evaluation& Rs) const
213 {
214 // ATTENTION: Rs is represented by the _first_ axis!
215 return inverseOilBTable_[regionIdx].eval(Rs, pressure, /*extrapolate=*/true);
216 }
217
221 template <class FluidState, class LhsEval = typename FluidState::Scalar>
222 std::pair<LhsEval, LhsEval>
223 inverseFormationVolumeFactorAndViscosity(const FluidState& fluidState, unsigned regionIdx)
224 {
225 const LhsEval& p = decay<LhsEval>(fluidState.pressure(FluidState::oilPhaseIdx));
226 const LhsEval& Rs = decay<LhsEval>(fluidState.Rs());
227
228 const auto satSegIdx = this->saturatedGasDissolutionFactorTable_[regionIdx].findSegmentIndex(p, /*extrapolate=*/ true);
229 const auto RsSat = this->saturatedGasDissolutionFactorTable_[regionIdx].eval(p, SegmentIndex{satSegIdx});
230 const bool useSaturatedTables = (fluidState.saturation(FluidState::gasPhaseIdx) > 0.0) && (Rs >= (1.0 - 1e-10) * RsSat);
231
232 if (useSaturatedTables) {
233 const LhsEval b = this->inverseSaturatedOilBTable_[regionIdx].eval(p, SegmentIndex{satSegIdx});
234 const LhsEval invBMu = this->inverseSaturatedOilBMuTable_[regionIdx].eval(p, SegmentIndex{satSegIdx});
235 const LhsEval mu = b / invBMu;
236 return { b, mu };
237 } else {
238 unsigned ii, jj1, jj2;
239 LhsEval alpha, beta1, beta2;
240 this->inverseOilBTable_[regionIdx].findPoints(ii, jj1, jj2, alpha, beta1, beta2, Rs, p, /*extrapolate =*/ true);
241 const LhsEval b = this->inverseOilBTable_[regionIdx].eval(ii, jj1, jj2, alpha, beta1, beta2);
242 const LhsEval invBMu = this->inverseOilBMuTable_[regionIdx].eval(ii, jj1, jj2, alpha, beta1, beta2);
243 const LhsEval mu = b / invBMu;
244 return { b, mu };
245 }
246 }
247
251 template <class Evaluation>
252 Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
253 const Evaluation& /*temperature*/,
254 const Evaluation& pressure) const
255 {
256 // ATTENTION: Rs is represented by the _first_ axis!
257 return inverseSaturatedOilBTable_[regionIdx].eval(pressure, /*extrapolate=*/true);
258 }
259
263 template <class Evaluation>
264 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
265 const Evaluation& /*temperature*/,
266 const Evaluation& pressure) const
267 { return saturatedGasDissolutionFactorTable_[regionIdx].eval(pressure, /*extrapolate=*/true); }
268
276 template <class Evaluation>
277 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
278 const Evaluation& /*temperature*/,
279 const Evaluation& pressure,
280 const Evaluation& oilSaturation,
281 Evaluation maxOilSaturation) const
282 {
283 Evaluation tmp =
284 saturatedGasDissolutionFactorTable_[regionIdx].eval(pressure, /*extrapolate=*/true);
285
286 // apply the vaporization parameters for the gas phase (cf. the Eclipse VAPPARS
287 // keyword)
288 maxOilSaturation = min(maxOilSaturation, Scalar(1.0));
289 if (vapPar2_ > 0.0 && maxOilSaturation > 0.01 && oilSaturation < maxOilSaturation) {
290 constexpr const Scalar eps = 0.001;
291 const Evaluation& So = max(oilSaturation, eps);
292 tmp *= max(1e-3, pow(So / maxOilSaturation, vapPar2_));
293 }
294
295 return tmp;
296 }
297
305 template <class Evaluation>
306 Evaluation saturationPressure(unsigned regionIdx,
307 const Evaluation&,
308 const Evaluation& Rs) const
309 {
310 using Toolbox = MathToolbox<Evaluation>;
311
312 const auto& RsTable = saturatedGasDissolutionFactorTable_[regionIdx];
313 constexpr const Scalar eps = std::numeric_limits<typename Toolbox::Scalar>::epsilon() * 1e6;
314
315 // use the saturation pressure function to get a pretty good initial value
316 Evaluation pSat = saturationPressure_[regionIdx].eval(Rs, /*extrapolate=*/true);
317
318 // Newton method to do the remaining work. If the initial
319 // value is good, this should only take two to three
320 // iterations...
321 bool onProbation = false;
322 for (int i = 0; i < 20; ++i) {
323 const Evaluation& f = RsTable.eval(pSat, /*extrapolate=*/true) - Rs;
324 const Evaluation& fPrime = RsTable.evalDerivative(pSat, /*extrapolate=*/true);
325
326 // If the derivative is "zero" Newton will not converge,
327 // so simply return our initial guess.
328 if (std::abs(scalarValue(fPrime)) < 1.0e-30) {
329 return pSat;
330 }
331
332 const Evaluation& delta = f / fPrime;
333
334 pSat -= delta;
335
336 if (pSat < 0.0) {
337 // if the pressure is lower than 0 Pascals, we set it back to 0. if this
338 // happens twice, we give up and just return 0 Pa...
339 if (onProbation) {
340 return 0.0;
341 }
342
343 onProbation = true;
344 pSat = 0.0;
345 }
346
347 if (std::abs(scalarValue(delta)) < std::abs(scalarValue(pSat))*eps) {
348 return pSat;
349 }
350 }
351
352 const std::string msg =
353 "Finding saturation pressure did not converge: "
354 " pSat = " + std::to_string(getValue(pSat)) +
355 ", Rs = " + std::to_string(getValue(Rs));
356 OpmLog::debug("Live oil saturation pressure", msg);
357 throw NumericalProblem(msg);
358 }
359
360 template <class Evaluation>
361 Evaluation diffusionCoefficient(const Evaluation& /*temperature*/,
362 const Evaluation& /*pressure*/,
363 unsigned /*compIdx*/) const
364 {
365 throw std::runtime_error("Not implemented: The PVT model does not provide "
366 "a diffusionCoefficient()");
367 }
368
369 Scalar gasReferenceDensity(unsigned regionIdx) const
370 { return gasReferenceDensity_[regionIdx]; }
371
372 Scalar oilReferenceDensity(unsigned regionIdx) const
373 { return oilReferenceDensity_[regionIdx]; }
374
375 const std::vector<TabulatedTwoDFunction>& inverseOilBTable() const
376 { return inverseOilBTable_; }
377
378 const std::vector<TabulatedTwoDFunction>& oilMuTable() const
379 { return oilMuTable_; }
380
381 const std::vector<TabulatedTwoDFunction>& inverseOilBMuTable() const
382 { return inverseOilBMuTable_; }
383
384 const std::vector<TabulatedOneDFunction>& saturatedOilMuTable() const
385 { return saturatedOilMuTable_; }
386
387 const std::vector<TabulatedOneDFunction>& inverseSaturatedOilBTable() const
388 { return inverseSaturatedOilBTable_; }
389
390 const std::vector<TabulatedOneDFunction>& inverseSaturatedOilBMuTable() const
391 { return inverseSaturatedOilBMuTable_; }
392
393 const std::vector<TabulatedOneDFunction>& saturatedGasDissolutionFactorTable() const
394 { return saturatedGasDissolutionFactorTable_; }
395
396 const std::vector<TabulatedOneDFunction>& saturationPressure() const
397 { return saturationPressure_; }
398
399 Scalar vapPar2() const
400 { return vapPar2_; }
401
402private:
403 void updateSaturationPressure_(unsigned regionIdx);
404
405 std::vector<Scalar> gasReferenceDensity_{};
406 std::vector<Scalar> oilReferenceDensity_{};
407 std::vector<TabulatedTwoDFunction> inverseOilBTable_{};
408 std::vector<TabulatedTwoDFunction> oilMuTable_{};
409 std::vector<TabulatedTwoDFunction> inverseOilBMuTable_{};
410 std::vector<TabulatedOneDFunction> saturatedOilMuTable_{};
411 std::vector<TabulatedOneDFunction> inverseSaturatedOilBTable_{};
412 std::vector<TabulatedOneDFunction> inverseSaturatedOilBMuTable_{};
413 std::vector<TabulatedOneDFunction> saturatedGasDissolutionFactorTable_{};
414 std::vector<TabulatedOneDFunction> saturationPressure_{};
415
416 Scalar vapPar2_ = 0.0;
417};
418
419} // namespace Opm
420
421#endif
Provides the OPM specific exception classes.
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
Implements a linearly interpolated scalar function that depends on one variable.
Implements a scalar function that depends on two variables and which is sampled uniformly in the X di...
Definition EclipseState.hpp:62
This class represents the Pressure-Volume-Temperature relations of the oil phas with dissolved gas.
Definition LiveOilPvt.hpp:51
Evaluation internalEnergy(unsigned, const Evaluation &, const Evaluation &, const Evaluation &) const
Returns the specific enthalpy [J/kg] of oil given a set of parameters.
Definition LiveOilPvt.hpp:159
void setSaturatedOilFormationVolumeFactor(unsigned regionIdx, const SamplingPoints &samplePoints)
Initialize the function for the oil formation volume factor.
Definition LiveOilPvt.cpp:246
Evaluation inverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure, const Evaluation &Rs) const
Returns the formation volume factor [-] of the fluid phase.
Definition LiveOilPvt.hpp:209
Evaluation saturatedGasDissolutionFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure, const Evaluation &oilSaturation, Evaluation maxOilSaturation) const
Returns the gas dissolution factor [m^3/m^3] of the oil phase.
Definition LiveOilPvt.hpp:277
Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure) const
Returns the formation volume factor [-] of the fluid phase.
Definition LiveOilPvt.hpp:252
void setSaturatedOilViscosity(unsigned regionIdx, const SamplingPoints &samplePoints)
Initialize the phase viscosity for gas saturated oil.
Definition LiveOilPvt.cpp:273
Evaluation saturationPressure(unsigned regionIdx, const Evaluation &, const Evaluation &Rs) const
Returns the saturation pressure of the oil phase [Pa] depending on its mass fraction of the gas compo...
Definition LiveOilPvt.hpp:306
void setSaturatedOilGasDissolutionFactor(unsigned regionIdx, const SamplingPoints &samplePoints)
Initialize the function for the gas dissolution factor .
Definition LiveOilPvt.hpp:94
Evaluation viscosity(unsigned regionIdx, const Evaluation &, const Evaluation &pressure, const Evaluation &Rs) const
Returns the dynamic viscosity [Pa s] of the fluid phase given a set of parameters.
Definition LiveOilPvt.hpp:178
void setReferenceDensities(unsigned regionIdx, Scalar rhoRefOil, Scalar rhoRefGas, Scalar)
Initialize the reference densities of all fluids for a given PVT region.
Definition LiveOilPvt.cpp:235
void setOilViscosity(unsigned regionIdx, const TabulatedTwoDFunction &muo)
Initialize the viscosity of the oil phase.
Definition LiveOilPvt.hpp:131
Evaluation saturatedGasDissolutionFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure) const
Returns the gas dissolution factor [m^3/m^3] of the oil phase.
Definition LiveOilPvt.hpp:264
Evaluation saturatedViscosity(unsigned regionIdx, const Evaluation &, const Evaluation &pressure) const
Returns the dynamic viscosity [Pa s] of the fluid phase given a set of parameters.
Definition LiveOilPvt.hpp:194
unsigned numRegions() const
Return the number of PVT regions which are considered by this PVT-object.
Definition LiveOilPvt.hpp:152
std::pair< LhsEval, LhsEval > inverseFormationVolumeFactorAndViscosity(const FluidState &fluidState, unsigned regionIdx)
Returns the formation volume factor [-] and viscosity [Pa s] of the fluid phase.
Definition LiveOilPvt.hpp:223
void initEnd()
Finish initializing the oil phase PVT properties.
Definition LiveOilPvt.cpp:300
void setInverseOilFormationVolumeFactor(unsigned regionIdx, const TabulatedTwoDFunction &invBo)
Initialize the function for the oil formation volume factor.
Definition LiveOilPvt.hpp:122
Definition Exceptions.hpp:40
Definition Schedule.hpp:101
Definition SimpleTable.hpp:35
Implements a linearly interpolated scalar function that depends on one variable.
Definition Tabulated1DFunction.hpp:51
Implements a scalar function that depends on two variables and which is sampled uniformly in the X di...
Definition UniformXTabulated2DFunction.hpp:55
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Definition MathToolbox.hpp:51
Definition Tabulated1DFunction.hpp:41