27 #ifndef OPM_LIVE_OIL_PVT_HPP
28 #define OPM_LIVE_OIL_PVT_HPP
35 #include <opm/input/eclipse/EclipseState/EclipseState.hpp>
36 #include <opm/input/eclipse/Schedule/OilVaporizationProperties.hpp>
37 #include <opm/input/eclipse/Schedule/Schedule.hpp>
38 #include <opm/input/eclipse/EclipseState/Tables/TableManager.hpp>
42 #include <opm/common/OpmLog/OpmLog.hpp>
50 template <
class Scalar>
53 using SamplingPoints = std::vector<std::pair<Scalar, Scalar>>;
64 LiveOilPvt(
const std::vector<Scalar>& gasReferenceDensity,
65 const std::vector<Scalar>& oilReferenceDensity,
66 const std::vector<TabulatedTwoDFunction>& inverseOilBTable,
67 const std::vector<TabulatedTwoDFunction>& oilMuTable,
68 const std::vector<TabulatedTwoDFunction>& inverseOilBMuTable,
69 const std::vector<TabulatedOneDFunction>& saturatedOilMuTable,
70 const std::vector<TabulatedOneDFunction>& inverseSaturatedOilBTable,
71 const std::vector<TabulatedOneDFunction>& inverseSaturatedOilBMuTable,
72 const std::vector<TabulatedOneDFunction>& saturatedGasDissolutionFactorTable,
75 : gasReferenceDensity_(gasReferenceDensity)
76 , oilReferenceDensity_(oilReferenceDensity)
77 , inverseOilBTable_(inverseOilBTable)
78 , oilMuTable_(oilMuTable)
79 , inverseOilBMuTable_(inverseOilBMuTable)
80 , saturatedOilMuTable_(saturatedOilMuTable)
81 , inverseSaturatedOilBTable_(inverseSaturatedOilBTable)
82 , inverseSaturatedOilBMuTable_(inverseSaturatedOilBMuTable)
83 , saturatedGasDissolutionFactorTable_(saturatedGasDissolutionFactorTable)
84 , saturationPressure_(saturationPressure)
92 void initFromState(
const EclipseState& eclState,
const Schedule& schedule)
94 const auto& pvtoTables = eclState.getTableManager().getPvtoTables();
95 const auto& densityTable = eclState.getTableManager().getDensityTable();
97 assert(pvtoTables.size() == densityTable.size());
102 for (
unsigned regionIdx = 0; regionIdx <
numRegions; ++ regionIdx) {
103 Scalar rhoRefO = densityTable[regionIdx].oil;
104 Scalar rhoRefG = densityTable[regionIdx].gas;
105 Scalar rhoRefW = densityTable[regionIdx].water;
111 for (
unsigned regionIdx = 0; regionIdx <
numRegions; ++ regionIdx) {
112 const auto& pvtoTable = pvtoTables[regionIdx];
114 const auto& saturatedTable = pvtoTable.getSaturatedTable();
115 assert(saturatedTable.numRows() > 1);
117 auto& oilMu = oilMuTable_[regionIdx];
118 auto& satOilMu = saturatedOilMuTable_[regionIdx];
119 auto& invOilB = inverseOilBTable_[regionIdx];
120 auto& invSatOilB = inverseSaturatedOilBTable_[regionIdx];
121 auto& gasDissolutionFac = saturatedGasDissolutionFactorTable_[regionIdx];
122 std::vector<Scalar> invSatOilBArray;
123 std::vector<Scalar> satOilMuArray;
126 for (
unsigned outerIdx = 0; outerIdx < saturatedTable.numRows(); ++ outerIdx) {
127 Scalar Rs = saturatedTable.get(
"RS", outerIdx);
128 Scalar BoSat = saturatedTable.get(
"BO", outerIdx);
129 Scalar muoSat = saturatedTable.get(
"MU", outerIdx);
131 satOilMuArray.push_back(muoSat);
132 invSatOilBArray.push_back(1.0/BoSat);
134 invOilB.appendXPos(Rs);
135 oilMu.appendXPos(Rs);
137 assert(invOilB.numX() == outerIdx + 1);
138 assert(oilMu.numX() == outerIdx + 1);
140 const auto& underSaturatedTable = pvtoTable.getUnderSaturatedTable(outerIdx);
141 size_t numRows = underSaturatedTable.numRows();
142 for (
unsigned innerIdx = 0; innerIdx < numRows; ++ innerIdx) {
143 Scalar po = underSaturatedTable.get(
"P", innerIdx);
144 Scalar Bo = underSaturatedTable.get(
"BO", innerIdx);
145 Scalar muo = underSaturatedTable.get(
"MU", innerIdx);
147 invOilB.appendSamplePoint(outerIdx, po, 1.0/Bo);
148 oilMu.appendSamplePoint(outerIdx, po, muo);
155 const auto& tmpPressureColumn = saturatedTable.getColumn(
"P");
156 const auto& tmpGasSolubilityColumn = saturatedTable.getColumn(
"RS");
158 invSatOilB.setXYContainers(tmpPressureColumn, invSatOilBArray);
159 satOilMu.setXYContainers(tmpPressureColumn, satOilMuArray);
160 gasDissolutionFac.setXYContainers(tmpPressureColumn, tmpGasSolubilityColumn);
163 updateSaturationPressure_(regionIdx);
165 for (
unsigned xIdx = 0; xIdx < invOilB.numX(); ++xIdx) {
167 assert(invOilB.numY(xIdx) > 0);
171 if (invOilB.numY(xIdx) > 1)
177 size_t masterTableIdx = xIdx + 1;
178 for (; masterTableIdx < saturatedTable.numRows(); ++masterTableIdx)
180 if (pvtoTable.getUnderSaturatedTable(masterTableIdx).numRows() > 1)
184 if (masterTableIdx >= saturatedTable.numRows())
185 throw std::runtime_error(
"PVTO tables are invalid: The last table must exhibit at least one "
186 "entry for undersaturated oil!");
189 extendPvtoTable_(regionIdx,
191 pvtoTable.getUnderSaturatedTable(xIdx),
192 pvtoTable.getUnderSaturatedTable(masterTableIdx));
197 const auto& oilVap = schedule[0].oilvap();
198 if (oilVap.getType() == OilVaporizationProperties::OilVaporization::VAPPARS) {
199 vapPar2_ = oilVap.vap2();
206 void extendPvtoTable_(
unsigned regionIdx,
208 const SimpleTable& curTable,
209 const SimpleTable& masterTable)
211 std::vector<double> pressuresArray = curTable.getColumn(
"P").vectorCopy();
212 std::vector<double> oilBArray = curTable.getColumn(
"BO").vectorCopy();
213 std::vector<double> oilMuArray = curTable.getColumn(
"MU").vectorCopy();
215 auto& invOilB = inverseOilBTable_[regionIdx];
216 auto& oilMu = oilMuTable_[regionIdx];
218 for (
unsigned newRowIdx = 1; newRowIdx < masterTable.numRows(); ++ newRowIdx) {
219 const auto& pressureColumn = masterTable.getColumn(
"P");
220 const auto& BOColumn = masterTable.getColumn(
"BO");
221 const auto& viscosityColumn = masterTable.getColumn(
"MU");
224 Scalar diffPo = pressureColumn[newRowIdx] - pressureColumn[newRowIdx - 1];
225 Scalar newPo = pressuresArray.back() + diffPo;
228 Scalar B1 = BOColumn[newRowIdx];
229 Scalar B2 = BOColumn[newRowIdx - 1];
230 Scalar x = (B1 - B2)/( (B1 + B2)/2.0 );
234 Scalar newBo = oilBArray.back()*(1.0 + x/2.0)/(1.0 - x/2.0);
237 Scalar mu1 = viscosityColumn[newRowIdx];
238 Scalar mu2 = viscosityColumn[newRowIdx - 1];
239 Scalar xMu = (mu1 - mu2)/( (mu1 + mu2)/2.0 );
243 Scalar newMuo = oilMuArray.back()*(1.0 + xMu/2)/(1.0 - xMu/2.0);
247 pressuresArray.push_back(newPo);
248 oilBArray.push_back(newBo);
249 oilMuArray.push_back(newMuo);
252 invOilB.appendSamplePoint(xIdx, newPo, 1.0/newBo);
253 oilMu.appendSamplePoint(xIdx, newPo, newMuo);
266 inverseSaturatedOilBTable_.resize(
numRegions);
267 inverseSaturatedOilBMuTable_.resize(
numRegions);
270 saturatedGasDissolutionFactorTable_.resize(
numRegions);
282 oilReferenceDensity_[regionIdx] = rhoRefOil;
283 gasReferenceDensity_[regionIdx] = rhoRefGas;
292 { saturatedGasDissolutionFactorTable_[regionIdx].setContainerOfTuples(samplePoints); }
305 constexpr
const Scalar T = 273.15 + 15.56;
306 auto& invOilB = inverseOilBTable_[regionIdx];
308 updateSaturationPressure_(regionIdx);
311 for (
size_t pIdx = 0; pIdx < samplePoints.size(); ++pIdx) {
312 Scalar p1 = std::get<0>(samplePoints[pIdx]);
313 Scalar p2 = p1 * 2.0;
315 Scalar Bo1 = std::get<1>(samplePoints[pIdx]);
316 Scalar drhoo_dp = (1.1200 - 1.1189)/((5000 - 4000)*6894.76);
317 Scalar Bo2 = Bo1/(1.0 + (p2 - p1)*drhoo_dp);
321 invOilB.appendXPos(Rs);
322 invOilB.appendSamplePoint(pIdx, p1, 1.0/Bo1);
323 invOilB.appendSamplePoint(pIdx, p2, 1.0/Bo2);
340 { inverseOilBTable_[regionIdx] = invBo; }
348 { oilMuTable_[regionIdx] = muo; }
359 constexpr
const Scalar T = 273.15 + 15.56;
362 saturatedOilMuTable_[regionIdx].setContainerOfTuples(samplePoints);
366 for (
size_t pIdx = 0; pIdx < samplePoints.size(); ++pIdx) {
367 Scalar p1 = std::get<0>(samplePoints[pIdx]);
368 Scalar p2 = p1 * 2.0;
371 Scalar mu1 = std::get<1>(samplePoints[pIdx]);
376 oilMuTable_[regionIdx].appendXPos(Rs);
377 oilMuTable_[regionIdx].appendSamplePoint(pIdx, p1, mu1);
378 oilMuTable_[regionIdx].appendSamplePoint(pIdx, p2, mu2);
389 for (
unsigned regionIdx = 0; regionIdx <
numRegions; ++ regionIdx) {
392 const auto& oilMu = oilMuTable_[regionIdx];
393 const auto& satOilMu = saturatedOilMuTable_[regionIdx];
394 const auto& invOilB = inverseOilBTable_[regionIdx];
395 assert(oilMu.numX() == invOilB.numX());
397 auto& invOilBMu = inverseOilBMuTable_[regionIdx];
398 auto& invSatOilB = inverseSaturatedOilBTable_[regionIdx];
399 auto& invSatOilBMu = inverseSaturatedOilBMuTable_[regionIdx];
401 std::vector<Scalar> satPressuresArray;
402 std::vector<Scalar> invSatOilBArray;
403 std::vector<Scalar> invSatOilBMuArray;
404 for (
unsigned rsIdx = 0; rsIdx < oilMu.numX(); ++rsIdx) {
405 invOilBMu.appendXPos(oilMu.xAt(rsIdx));
407 assert(oilMu.numY(rsIdx) == invOilB.numY(rsIdx));
409 size_t numPressures = oilMu.numY(rsIdx);
410 for (
unsigned pIdx = 0; pIdx < numPressures; ++pIdx)
411 invOilBMu.appendSamplePoint(rsIdx,
412 oilMu.yAt(rsIdx, pIdx),
413 invOilB.valueAt(rsIdx, pIdx)
414 / oilMu.valueAt(rsIdx, pIdx));
419 satPressuresArray.push_back(oilMu.yAt(rsIdx, 0));
420 invSatOilBArray.push_back(invOilB.valueAt(rsIdx, 0));
421 invSatOilBMuArray.push_back(invSatOilBArray.back()/satOilMu.valueAt(rsIdx));
424 invSatOilB.setXYContainers(satPressuresArray, invSatOilBArray);
425 invSatOilBMu.setXYContainers(satPressuresArray, invSatOilBMuArray);
427 updateSaturationPressure_(regionIdx);
435 {
return inverseOilBMuTable_.size(); }
440 template <
class Evaluation>
444 const Evaluation&)
const
446 throw std::runtime_error(
"Requested the enthalpy of oil but the thermal option is not enabled");
452 template <
class Evaluation>
455 const Evaluation& pressure,
456 const Evaluation& Rs)
const
459 const Evaluation& invBo = inverseOilBTable_[regionIdx].eval(Rs, pressure,
true);
460 const Evaluation& invMuoBo = inverseOilBMuTable_[regionIdx].eval(Rs, pressure,
true);
462 return invBo/invMuoBo;
468 template <
class Evaluation>
471 const Evaluation& pressure)
const
474 const Evaluation& invBo = inverseSaturatedOilBTable_[regionIdx].eval(pressure,
true);
475 const Evaluation& invMuoBo = inverseSaturatedOilBMuTable_[regionIdx].eval(pressure,
true);
477 return invBo/invMuoBo;
483 template <
class Evaluation>
486 const Evaluation& pressure,
487 const Evaluation& Rs)
const
490 return inverseOilBTable_[regionIdx].eval(Rs, pressure,
true);
496 template <
class Evaluation>
499 const Evaluation& pressure)
const
502 return inverseSaturatedOilBTable_[regionIdx].eval(pressure,
true);
508 template <
class Evaluation>
511 const Evaluation& pressure)
const
512 {
return saturatedGasDissolutionFactorTable_[regionIdx].eval(pressure,
true); }
521 template <
class Evaluation>
524 const Evaluation& pressure,
525 const Evaluation& oilSaturation,
526 Evaluation maxOilSaturation)
const
529 saturatedGasDissolutionFactorTable_[regionIdx].eval(pressure,
true);
533 maxOilSaturation = min(maxOilSaturation, Scalar(1.0));
534 if (vapPar2_ > 0.0 && maxOilSaturation > 0.01 && oilSaturation < maxOilSaturation) {
535 constexpr
const Scalar eps = 0.001;
536 const Evaluation& So = max(oilSaturation, eps);
537 tmp *= max(1e-3, pow(So/maxOilSaturation, vapPar2_));
549 template <
class Evaluation>
552 const Evaluation& Rs)
const
556 const auto& RsTable = saturatedGasDissolutionFactorTable_[regionIdx];
557 constexpr
const Scalar eps = std::numeric_limits<typename Toolbox::Scalar>::epsilon()*1e6;
560 Evaluation pSat = saturationPressure_[regionIdx].eval(Rs,
true);
565 bool onProbation =
false;
566 for (
int i = 0; i < 20; ++i) {
567 const Evaluation& f = RsTable.eval(pSat,
true) - Rs;
568 const Evaluation& fPrime = RsTable.evalDerivative(pSat,
true);
572 if (std::abs(scalarValue(fPrime)) < 1.0e-30) {
576 const Evaluation& delta = f/fPrime;
590 if (std::abs(scalarValue(delta)) < std::abs(scalarValue(pSat))*eps)
594 std::stringstream errlog;
595 errlog <<
"Finding saturation pressure did not converge:"
596 <<
" pSat = " << pSat
599 OpmLog::debug(
"Live oil saturation pressure", errlog.str());
604 template <
class Evaluation>
605 Evaluation diffusionCoefficient(
const Evaluation& ,
609 throw std::runtime_error(
"Not implemented: The PVT model does not provide a diffusionCoefficient()");
611 const Scalar gasReferenceDensity(
unsigned regionIdx)
const
612 {
return gasReferenceDensity_[regionIdx]; }
614 const Scalar oilReferenceDensity(
unsigned regionIdx)
const
615 {
return oilReferenceDensity_[regionIdx]; }
617 const std::vector<TabulatedTwoDFunction>& inverseOilBTable()
const
618 {
return inverseOilBTable_; }
620 const std::vector<TabulatedTwoDFunction>& oilMuTable()
const
621 {
return oilMuTable_; }
623 const std::vector<TabulatedTwoDFunction>& inverseOilBMuTable()
const
624 {
return inverseOilBMuTable_; }
626 const std::vector<TabulatedOneDFunction>& saturatedOilMuTable()
const
627 {
return saturatedOilMuTable_; }
629 const std::vector<TabulatedOneDFunction>& inverseSaturatedOilBTable()
const
630 {
return inverseSaturatedOilBTable_; }
632 const std::vector<TabulatedOneDFunction>& inverseSaturatedOilBMuTable()
const
633 {
return inverseSaturatedOilBMuTable_; }
635 const std::vector<TabulatedOneDFunction>& saturatedGasDissolutionFactorTable()
const
636 {
return saturatedGasDissolutionFactorTable_; }
638 const std::vector<TabulatedOneDFunction>& saturationPressure()
const
639 {
return saturationPressure_; }
641 Scalar vapPar2()
const
644 bool operator==(
const LiveOilPvt<Scalar>& data)
const
646 return this->gasReferenceDensity_ == data.gasReferenceDensity_ &&
647 this->oilReferenceDensity_ == data.oilReferenceDensity_ &&
648 this->inverseOilBTable() == data.inverseOilBTable() &&
649 this->oilMuTable() == data.oilMuTable() &&
650 this->inverseOilBMuTable() == data.inverseOilBMuTable() &&
651 this->saturatedOilMuTable() == data.saturatedOilMuTable() &&
652 this->inverseSaturatedOilBTable() == data.inverseSaturatedOilBTable() &&
653 this->inverseSaturatedOilBMuTable() == data.inverseSaturatedOilBMuTable() &&
654 this->saturatedGasDissolutionFactorTable() == data.saturatedGasDissolutionFactorTable() &&
655 this->vapPar2() == data.vapPar2();
659 void updateSaturationPressure_(
unsigned regionIdx)
661 typedef std::pair<Scalar, Scalar> Pair;
662 const auto& gasDissolutionFac = saturatedGasDissolutionFactorTable_[regionIdx];
666 size_t n = gasDissolutionFac.numSamples();
667 Scalar delta = (gasDissolutionFac.xMax() - gasDissolutionFac.xMin())/Scalar(n + 1);
669 SamplingPoints pSatSamplePoints;
671 for (
size_t i=0; i <= n; ++ i) {
672 Scalar pSat = gasDissolutionFac.xMin() + Scalar(i)*delta;
678 pSatSamplePoints.push_back(val);
682 auto x_coord_comparator = [](
const Pair& a,
const Pair& b) {
return a.first == b.first; };
683 auto last = std::unique(pSatSamplePoints.begin(), pSatSamplePoints.end(), x_coord_comparator);
684 if (std::distance(pSatSamplePoints.begin(), last) > 1)
685 pSatSamplePoints.erase(last, pSatSamplePoints.end());
687 saturationPressure_[regionIdx].setContainerOfTuples(pSatSamplePoints);
690 std::vector<Scalar> gasReferenceDensity_;
691 std::vector<Scalar> oilReferenceDensity_;
692 std::vector<TabulatedTwoDFunction> inverseOilBTable_;
693 std::vector<TabulatedTwoDFunction> oilMuTable_;
694 std::vector<TabulatedTwoDFunction> inverseOilBMuTable_;
695 std::vector<TabulatedOneDFunction> saturatedOilMuTable_;
696 std::vector<TabulatedOneDFunction> inverseSaturatedOilBTable_;
697 std::vector<TabulatedOneDFunction> inverseSaturatedOilBMuTable_;
698 std::vector<TabulatedOneDFunction> saturatedGasDissolutionFactorTable_;
699 std::vector<TabulatedOneDFunction> saturationPressure_;
Implements a linearly interpolated scalar function that depends on one variable.
This class represents the Pressure-Volume-Temperature relations of the oil phas with dissolved gas.
Definition: LiveOilPvt.hpp:52
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:441
void setSaturatedOilFormationVolumeFactor(unsigned regionIdx, const SamplingPoints &samplePoints)
Initialize the function for the oil formation volume factor.
Definition: LiveOilPvt.hpp:303
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:484
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:522
Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure) const
Returns the formation volume factor [-] of the fluid phase.
Definition: LiveOilPvt.hpp:497
void setSaturatedOilViscosity(unsigned regionIdx, const SamplingPoints &samplePoints)
Initialize the phase viscosity for gas saturated oil.
Definition: LiveOilPvt.hpp:357
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:550
void setSaturatedOilGasDissolutionFactor(unsigned regionIdx, const SamplingPoints &samplePoints)
Initialize the function for the gas dissolution factor .
Definition: LiveOilPvt.hpp:291
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:453
void setReferenceDensities(unsigned regionIdx, Scalar rhoRefOil, Scalar rhoRefGas, Scalar)
Initialize the reference densities of all fluids for a given PVT region.
Definition: LiveOilPvt.hpp:277
void setOilViscosity(unsigned regionIdx, const TabulatedTwoDFunction &muo)
Initialize the viscosity of the oil phase.
Definition: LiveOilPvt.hpp:347
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:509
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:469
unsigned numRegions() const
Return the number of PVT regions which are considered by this PVT-object.
Definition: LiveOilPvt.hpp:434
void initEnd()
Finish initializing the oil phase PVT properties.
Definition: LiveOilPvt.hpp:385
void setInverseOilFormationVolumeFactor(unsigned regionIdx, const TabulatedTwoDFunction &invBo)
Initialize the function for the oil formation volume factor.
Definition: LiveOilPvt.hpp:339
Definition: Exceptions.hpp:46
Implements a linearly interpolated scalar function that depends on one variable.
Definition: Tabulated1DFunction.hpp:48