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
39class EclipseState;
40class Schedule;
41class SimpleTable;
42
47template <class Scalar>
49{
50 using SamplingPoints = std::vector<std::pair<Scalar, Scalar>>;
51
52public:
53 using TabulatedTwoDFunction = UniformXTabulated2DFunction<Scalar>;
54 using TabulatedOneDFunction = Tabulated1DFunction<Scalar>;
55
59 void initFromState(const EclipseState& eclState, const Schedule& schedule);
60
61private:
62 void extendPvtoTable_(unsigned regionIdx,
63 unsigned xIdx,
64 const SimpleTable& curTable,
65 const SimpleTable& masterTable);
66
67public:
68 void setNumRegions(size_t numRegions);
69
70 void setVapPars(const Scalar, const Scalar par2)
71 {
72 vapPar2_ = par2;
73 }
74
78 void setReferenceDensities(unsigned regionIdx,
79 Scalar rhoRefOil,
80 Scalar rhoRefGas,
81 Scalar /*rhoRefWater*/);
82
89 void setSaturatedOilGasDissolutionFactor(unsigned regionIdx,
90 const SamplingPoints& samplePoints)
91 { saturatedGasDissolutionFactorTable_[regionIdx].setContainerOfTuples(samplePoints); }
92
102 void setSaturatedOilFormationVolumeFactor(unsigned regionIdx,
103 const SamplingPoints& samplePoints);
104
117 void setInverseOilFormationVolumeFactor(unsigned regionIdx,
118 const TabulatedTwoDFunction& invBo)
119 { inverseOilBTable_[regionIdx] = invBo; }
120
126 void setOilViscosity(unsigned regionIdx, const TabulatedTwoDFunction& muo)
127 { oilMuTable_[regionIdx] = muo; }
128
136 void setSaturatedOilViscosity(unsigned regionIdx,
137 const SamplingPoints& samplePoints);
138
142 void initEnd();
143
147 unsigned numRegions() const
148 { return inverseOilBMuTable_.size(); }
149
153 template <class Evaluation>
154 Evaluation internalEnergy(unsigned,
155 const Evaluation&,
156 const Evaluation&,
157 const Evaluation&) const
158 {
159 throw std::runtime_error("Requested the enthalpy of oil but the thermal "
160 "option is not enabled");
161 }
162
163 Scalar hVap(unsigned) const
164 {
165 throw std::runtime_error("Requested the hvap of oil but the thermal "
166 "option is not enabled");
167 }
168
172 template <class Evaluation>
173 Evaluation viscosity(unsigned regionIdx,
174 const Evaluation& /*temperature*/,
175 const Evaluation& pressure,
176 const Evaluation& Rs) const
177 {
178 // ATTENTION: Rs is the first axis!
179 const Evaluation& invBo = inverseOilBTable_[regionIdx].eval(Rs, pressure, /*extrapolate=*/true);
180 const Evaluation& invMuoBo = inverseOilBMuTable_[regionIdx].eval(Rs, pressure, /*extrapolate=*/true);
181
182 return invBo / invMuoBo;
183 }
184
188 template <class Evaluation>
189 Evaluation saturatedViscosity(unsigned regionIdx,
190 const Evaluation& /*temperature*/,
191 const Evaluation& pressure) const
192 {
193 // ATTENTION: Rs is the first axis!
194 const Evaluation& invBo = inverseSaturatedOilBTable_[regionIdx].eval(pressure, /*extrapolate=*/true);
195 const Evaluation& invMuoBo = inverseSaturatedOilBMuTable_[regionIdx].eval(pressure, /*extrapolate=*/true);
196
197 return invBo / invMuoBo;
198 }
199
203 template <class Evaluation>
204 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
205 const Evaluation& /*temperature*/,
206 const Evaluation& pressure,
207 const Evaluation& Rs) const
208 {
209 // ATTENTION: Rs is represented by the _first_ axis!
210 return inverseOilBTable_[regionIdx].eval(Rs, pressure, /*extrapolate=*/true);
211 }
212
216 template <class FluidState, class LhsEval = typename FluidState::ValueType>
217 std::pair<LhsEval, LhsEval>
218 inverseFormationVolumeFactorAndViscosity(const FluidState& fluidState, unsigned regionIdx)
219 {
220 const LhsEval& p = decay<LhsEval>(fluidState.pressure(FluidState::oilPhaseIdx));
221 const LhsEval& Rs = decay<LhsEval>(fluidState.Rs());
222
223 const auto satSegIdx = this->saturatedGasDissolutionFactorTable_[regionIdx].findSegmentIndex(p, /*extrapolate=*/ true);
224 const auto RsSat = this->saturatedGasDissolutionFactorTable_[regionIdx].eval(p, SegmentIndex{satSegIdx});
225 const bool useSaturatedTables = (fluidState.saturation(FluidState::gasPhaseIdx) > 0.0) && (Rs >= (1.0 - 1e-10) * RsSat);
226
227 if (useSaturatedTables) {
228 const LhsEval b = this->inverseSaturatedOilBTable_[regionIdx].eval(p, SegmentIndex{satSegIdx});
229 const LhsEval invBMu = this->inverseSaturatedOilBMuTable_[regionIdx].eval(p, SegmentIndex{satSegIdx});
230 const LhsEval mu = b / invBMu;
231 return { b, mu };
232 } else {
233 unsigned ii, jj1, jj2;
234 LhsEval alpha, beta1, beta2;
235 this->inverseOilBTable_[regionIdx].findPoints(ii, jj1, jj2, alpha, beta1, beta2, Rs, p, /*extrapolate =*/ true);
236 const LhsEval b = this->inverseOilBTable_[regionIdx].eval(ii, jj1, jj2, alpha, beta1, beta2);
237 const LhsEval invBMu = this->inverseOilBMuTable_[regionIdx].eval(ii, jj1, jj2, alpha, beta1, beta2);
238 const LhsEval mu = b / invBMu;
239 return { b, mu };
240 }
241 }
242
246 template <class Evaluation>
247 Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
248 const Evaluation& /*temperature*/,
249 const Evaluation& pressure) const
250 {
251 // ATTENTION: Rs is represented by the _first_ axis!
252 return inverseSaturatedOilBTable_[regionIdx].eval(pressure, /*extrapolate=*/true);
253 }
254
258 template <class Evaluation>
259 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
260 const Evaluation& /*temperature*/,
261 const Evaluation& pressure) const
262 { return saturatedGasDissolutionFactorTable_[regionIdx].eval(pressure, /*extrapolate=*/true); }
263
271 template <class Evaluation>
272 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
273 const Evaluation& /*temperature*/,
274 const Evaluation& pressure,
275 const Evaluation& oilSaturation,
276 Evaluation maxOilSaturation) const
277 {
278 Evaluation tmp =
279 saturatedGasDissolutionFactorTable_[regionIdx].eval(pressure, /*extrapolate=*/true);
280
281 // apply the vaporization parameters for the gas phase (cf. the Eclipse VAPPARS
282 // keyword)
283 maxOilSaturation = min(maxOilSaturation, Scalar(1.0));
284 if (vapPar2_ > 0.0 && maxOilSaturation > 0.01 && oilSaturation < maxOilSaturation) {
285 constexpr const Scalar eps = 0.001;
286 const Evaluation& So = max(oilSaturation, eps);
287 tmp *= max(1e-3, pow(So / maxOilSaturation, vapPar2_));
288 }
289
290 return tmp;
291 }
292
300 template <class Evaluation>
301 Evaluation saturationPressure(unsigned regionIdx,
302 const Evaluation&,
303 const Evaluation& Rs) const
304 {
305 using Toolbox = MathToolbox<Evaluation>;
306
307 const auto& RsTable = saturatedGasDissolutionFactorTable_[regionIdx];
308 constexpr const Scalar eps = std::numeric_limits<typename Toolbox::Scalar>::epsilon() * 1e6;
309
310 // use the saturation pressure function to get a pretty good initial value
311 Evaluation pSat = saturationPressure_[regionIdx].eval(Rs, /*extrapolate=*/true);
312
313 // Newton method to do the remaining work. If the initial
314 // value is good, this should only take two to three
315 // iterations...
316 bool onProbation = false;
317 for (int i = 0; i < 20; ++i) {
318 const Evaluation& f = RsTable.eval(pSat, /*extrapolate=*/true) - Rs;
319 const Evaluation& fPrime = RsTable.evalDerivative(pSat, /*extrapolate=*/true);
320
321 // If the derivative is "zero" Newton will not converge,
322 // so simply return our initial guess.
323 if (std::abs(scalarValue(fPrime)) < 1.0e-30) {
324 return pSat;
325 }
326
327 const Evaluation& delta = f / fPrime;
328
329 pSat -= delta;
330
331 if (pSat < 0.0) {
332 // if the pressure is lower than 0 Pascals, we set it back to 0. if this
333 // happens twice, we give up and just return 0 Pa...
334 if (onProbation) {
335 return 0.0;
336 }
337
338 onProbation = true;
339 pSat = 0.0;
340 }
341
342 if (std::abs(scalarValue(delta)) < std::abs(scalarValue(pSat))*eps) {
343 return pSat;
344 }
345 }
346
347 const std::string msg =
348 "Finding saturation pressure did not converge: "
349 " pSat = " + std::to_string(getValue(pSat)) +
350 ", Rs = " + std::to_string(getValue(Rs));
351 OpmLog::debug("Live oil saturation pressure", msg);
352 throw NumericalProblem(msg);
353 }
354
355 template <class Evaluation>
356 Evaluation diffusionCoefficient(const Evaluation& /*temperature*/,
357 const Evaluation& /*pressure*/,
358 unsigned /*compIdx*/) const
359 {
360 throw std::runtime_error("Not implemented: The PVT model does not provide "
361 "a diffusionCoefficient()");
362 }
363
364 Scalar gasReferenceDensity(unsigned regionIdx) const
365 { return gasReferenceDensity_[regionIdx]; }
366
367 Scalar oilReferenceDensity(unsigned regionIdx) const
368 { return oilReferenceDensity_[regionIdx]; }
369
370 const std::vector<TabulatedTwoDFunction>& inverseOilBTable() const
371 { return inverseOilBTable_; }
372
373 const std::vector<TabulatedTwoDFunction>& oilMuTable() const
374 { return oilMuTable_; }
375
376 const std::vector<TabulatedTwoDFunction>& inverseOilBMuTable() const
377 { return inverseOilBMuTable_; }
378
379 const std::vector<TabulatedOneDFunction>& saturatedOilMuTable() const
380 { return saturatedOilMuTable_; }
381
382 const std::vector<TabulatedOneDFunction>& inverseSaturatedOilBTable() const
383 { return inverseSaturatedOilBTable_; }
384
385 const std::vector<TabulatedOneDFunction>& inverseSaturatedOilBMuTable() const
386 { return inverseSaturatedOilBMuTable_; }
387
388 const std::vector<TabulatedOneDFunction>& saturatedGasDissolutionFactorTable() const
389 { return saturatedGasDissolutionFactorTable_; }
390
391 const std::vector<TabulatedOneDFunction>& saturationPressure() const
392 { return saturationPressure_; }
393
394 Scalar vapPar2() const
395 { return vapPar2_; }
396
397private:
398 void updateSaturationPressure_(unsigned regionIdx);
399
400 std::vector<Scalar> gasReferenceDensity_{};
401 std::vector<Scalar> oilReferenceDensity_{};
402 std::vector<TabulatedTwoDFunction> inverseOilBTable_{};
403 std::vector<TabulatedTwoDFunction> oilMuTable_{};
404 std::vector<TabulatedTwoDFunction> inverseOilBMuTable_{};
405 std::vector<TabulatedOneDFunction> saturatedOilMuTable_{};
406 std::vector<TabulatedOneDFunction> inverseSaturatedOilBTable_{};
407 std::vector<TabulatedOneDFunction> inverseSaturatedOilBMuTable_{};
408 std::vector<TabulatedOneDFunction> saturatedGasDissolutionFactorTable_{};
409 std::vector<TabulatedOneDFunction> saturationPressure_{};
410
411 Scalar vapPar2_ = 0.0;
412};
413
414} // namespace Opm
415
416#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:66
This class represents the Pressure-Volume-Temperature relations of the oil phas with dissolved gas.
Definition LiveOilPvt.hpp:49
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:154
void setSaturatedOilFormationVolumeFactor(unsigned regionIdx, const SamplingPoints &samplePoints)
Initialize the function for the oil formation volume factor.
Definition LiveOilPvt.cpp:242
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:204
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:272
Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure) const
Returns the formation volume factor [-] of the fluid phase.
Definition LiveOilPvt.hpp:247
void initFromState(const EclipseState &eclState, const Schedule &schedule)
Initialize the oil parameters via the data specified by the PVTO ECL keyword.
Definition LiveOilPvt.cpp:40
void setSaturatedOilViscosity(unsigned regionIdx, const SamplingPoints &samplePoints)
Initialize the phase viscosity for gas saturated oil.
Definition LiveOilPvt.cpp:269
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:301
void setSaturatedOilGasDissolutionFactor(unsigned regionIdx, const SamplingPoints &samplePoints)
Initialize the function for the gas dissolution factor .
Definition LiveOilPvt.hpp:89
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:173
void setReferenceDensities(unsigned regionIdx, Scalar rhoRefOil, Scalar rhoRefGas, Scalar)
Initialize the reference densities of all fluids for a given PVT region.
Definition LiveOilPvt.cpp:231
void setOilViscosity(unsigned regionIdx, const TabulatedTwoDFunction &muo)
Initialize the viscosity of the oil phase.
Definition LiveOilPvt.hpp:126
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:259
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:189
unsigned numRegions() const
Return the number of PVT regions which are considered by this PVT-object.
Definition LiveOilPvt.hpp:147
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:218
void initEnd()
Finish initializing the oil phase PVT properties.
Definition LiveOilPvt.cpp:296
void setInverseOilFormationVolumeFactor(unsigned regionIdx, const TabulatedTwoDFunction &invBo)
Initialize the function for the oil formation volume factor.
Definition LiveOilPvt.hpp:117
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