DryHumidGasPvt.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_DRY_HUMID_GAS_PVT_HPP
28#define OPM_DRY_HUMID_GAS_PVT_HPP
29
34
35#if HAVE_ECL_INPUT
36#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
37#include <opm/input/eclipse/EclipseState/Tables/TableManager.hpp>
38
39#endif
40
41#if HAVE_OPM_COMMON
42#include <opm/common/OpmLog/OpmLog.hpp>
43#endif
44
45namespace Opm {
46
51template <class Scalar>
53{
54 using SamplingPoints = std::vector<std::pair<Scalar, Scalar>>;
55
56public:
59
61 {
62 vapPar1_ = 0.0;
63 }
64
65 DryHumidGasPvt(const std::vector<Scalar>& gasReferenceDensity,
66 const std::vector<Scalar>& waterReferenceDensity,
67 const std::vector<TabulatedTwoDFunction>& inverseGasB,
68 const std::vector<TabulatedOneDFunction>& inverseSaturatedGasB,
69 const std::vector<TabulatedTwoDFunction>& gasMu,
70 const std::vector<TabulatedTwoDFunction>& inverseGasBMu,
71 const std::vector<TabulatedOneDFunction>& inverseSaturatedGasBMu,
72 const std::vector<TabulatedOneDFunction>& saturatedWaterVaporizationFactorTable,
73 const std::vector<TabulatedOneDFunction>& saturationPressure,
74 Scalar vapPar1)
75 : gasReferenceDensity_(gasReferenceDensity)
76 , waterReferenceDensity_(waterReferenceDensity)
77 , inverseGasB_(inverseGasB)
78 , inverseSaturatedGasB_(inverseSaturatedGasB)
79 , gasMu_(gasMu)
80 , inverseGasBMu_(inverseGasBMu)
81 , inverseSaturatedGasBMu_(inverseSaturatedGasBMu)
82 , saturatedWaterVaporizationFactorTable_(saturatedWaterVaporizationFactorTable)
83 , saturationPressure_(saturationPressure)
84 , vapPar1_(vapPar1)
85 {
86 }
87
88
89#if HAVE_ECL_INPUT
95 void initFromState(const EclipseState& eclState, const Schedule&)
96 {
97 const auto& pvtgwTables = eclState.getTableManager().getPvtgwTables();
98 const auto& densityTable = eclState.getTableManager().getDensityTable();
99
100 assert(pvtgwTables.size() == densityTable.size());
101
102 size_t numRegions = pvtgwTables.size();
104
105 for (unsigned regionIdx = 0; regionIdx < numRegions; ++ regionIdx) {
106 Scalar rhoRefO = densityTable[regionIdx].oil;
107 Scalar rhoRefG = densityTable[regionIdx].gas;
108 Scalar rhoRefW = densityTable[regionIdx].water;
109
110 setReferenceDensities(regionIdx, rhoRefO, rhoRefG, rhoRefW);
111 }
112
113 enableRwgSalt_ = !eclState.getTableManager().getRwgSaltTables().empty();
114 if (enableRwgSalt_)
115 {
116 const auto& rwgsaltTables = eclState.getTableManager().getRwgSaltTables();
117
118 for (unsigned regionIdx = 0; regionIdx < numRegions; ++ regionIdx) {
119 const auto& rwgsaltTable = rwgsaltTables[regionIdx];
120 const auto& saturatedTable = rwgsaltTable.getSaturatedTable();
121 assert(saturatedTable.numRows() > 1);
122
123 auto& waterVaporizationFac = saturatedWaterVaporizationSaltFactorTable_[regionIdx];
124 for (unsigned outerIdx = 0; outerIdx < saturatedTable.numRows(); ++ outerIdx) {
125 const auto& underSaturatedTable = rwgsaltTable.getUnderSaturatedTable(outerIdx);
126 Scalar pg = saturatedTable.get("PG" , outerIdx);
127 waterVaporizationFac.appendXPos(pg);
128
129 size_t numRows = underSaturatedTable.numRows();
130 for (size_t innerIdx = 0; innerIdx < numRows; ++ innerIdx) {
131 Scalar saltConcentration = underSaturatedTable.get("C_SALT" , innerIdx);
132 Scalar rvwSat= underSaturatedTable.get("RVW" , innerIdx);
133
134 waterVaporizationFac.appendSamplePoint(outerIdx, saltConcentration, rvwSat);
135 }
136 }
137 }
138 }
139
140 for (unsigned regionIdx = 0; regionIdx < numRegions; ++ regionIdx) {
141 const auto& pvtgwTable = pvtgwTables[regionIdx];
142
143 const auto& saturatedTable = pvtgwTable.getSaturatedTable();
144 assert(saturatedTable.numRows() > 1);
145
146 auto& gasMu = gasMu_[regionIdx];
147 auto& invGasB = inverseGasB_[regionIdx];
148 auto& invSatGasB = inverseSaturatedGasB_[regionIdx];
149 auto& invSatGasBMu = inverseSaturatedGasBMu_[regionIdx];
150 auto& waterVaporizationFac = saturatedWaterVaporizationFactorTable_[regionIdx];
151
152 waterVaporizationFac.setXYArrays(saturatedTable.numRows(),
153 saturatedTable.getColumn("PG"),
154 saturatedTable.getColumn("RW"));
155
156 std::vector<Scalar> invSatGasBArray;
157 std::vector<Scalar> invSatGasBMuArray;
158
159 // extract the table for the gas dissolution and the oil formation volume factors
160 for (unsigned outerIdx = 0; outerIdx < saturatedTable.numRows(); ++ outerIdx) {
161 Scalar pg = saturatedTable.get("PG" , outerIdx);
162 Scalar B = saturatedTable.get("BG" , outerIdx);
163 Scalar mu = saturatedTable.get("MUG" , outerIdx);
164
165 invGasB.appendXPos(pg);
166 gasMu.appendXPos(pg);
167
168 invSatGasBArray.push_back(1.0/B);
169 invSatGasBMuArray.push_back(1.0/(mu*B));
170
171 assert(invGasB.numX() == outerIdx + 1);
172 assert(gasMu.numX() == outerIdx + 1);
173
174 const auto& underSaturatedTable = pvtgwTable.getUnderSaturatedTable(outerIdx);
175 size_t numRows = underSaturatedTable.numRows();
176 for (size_t innerIdx = 0; innerIdx < numRows; ++ innerIdx) {
177 Scalar Rw = underSaturatedTable.get("RW" , innerIdx);
178 Scalar Bg = underSaturatedTable.get("BG" , innerIdx);
179 Scalar mug = underSaturatedTable.get("MUG" , innerIdx);
180
181 invGasB.appendSamplePoint(outerIdx, Rw, 1.0/Bg);
182 gasMu.appendSamplePoint(outerIdx, Rw, mug);
183 }
184 }
185
186 {
187 std::vector<double> tmpPressure = saturatedTable.getColumn("PG").vectorCopy( );
188
189 invSatGasB.setXYContainers(tmpPressure, invSatGasBArray);
190 invSatGasBMu.setXYContainers(tmpPressure, invSatGasBMuArray);
191 }
192
193 // make sure to have at least two sample points per gas pressure value
194 for (unsigned xIdx = 0; xIdx < invGasB.numX(); ++xIdx) {
195 // a single sample point is definitely needed
196 assert(invGasB.numY(xIdx) > 0);
197
198 // everything is fine if the current table has two or more sampling points
199 // for a given mole fraction
200 if (invGasB.numY(xIdx) > 1)
201 continue;
202
203 // find the master table which will be used as a template to extend the
204 // current line. We define master table as the first table which has values
205 // for undersaturated gas...
206 size_t masterTableIdx = xIdx + 1;
207 for (; masterTableIdx < saturatedTable.numRows(); ++masterTableIdx)
208 {
209 if (pvtgwTable.getUnderSaturatedTable(masterTableIdx).numRows() > 1)
210 break;
211 }
212
213 if (masterTableIdx >= saturatedTable.numRows())
214 throw std::runtime_error("PVTGW tables are invalid: The last table must exhibit at least one "
215 "entry for undersaturated gas!");
216
217
218 // extend the current table using the master table.
219 extendPvtgwTable_(regionIdx,
220 xIdx,
221 pvtgwTable.getUnderSaturatedTable(xIdx),
222 pvtgwTable.getUnderSaturatedTable(masterTableIdx));
223 }
224 }
225
226 vapPar1_ = 0.0;
227
228 initEnd();
229 }
230
231private:
232 void extendPvtgwTable_(unsigned regionIdx,
233 unsigned xIdx,
234 const SimpleTable& curTable,
235 const SimpleTable& masterTable)
236 {
237 std::vector<double> RwArray = curTable.getColumn("RW").vectorCopy();
238 std::vector<double> gasBArray = curTable.getColumn("BG").vectorCopy();
239 std::vector<double> gasMuArray = curTable.getColumn("MUG").vectorCopy();
240
241 auto& invGasB = inverseGasB_[regionIdx];
242 auto& gasMu = gasMu_[regionIdx];
243
244 for (size_t newRowIdx = 1; newRowIdx < masterTable.numRows(); ++ newRowIdx) {
245 const auto& RWColumn = masterTable.getColumn("RW");
246 const auto& BGColumn = masterTable.getColumn("BG");
247 const auto& viscosityColumn = masterTable.getColumn("MUG");
248
249 // compute the vaporized water factor Rw for the new entry
250 Scalar diffRw = RWColumn[newRowIdx] - RWColumn[newRowIdx - 1];
251 Scalar newRw = RwArray.back() + diffRw;
252
253 // calculate the compressibility of the master table
254 Scalar B1 = BGColumn[newRowIdx];
255 Scalar B2 = BGColumn[newRowIdx - 1];
256 Scalar x = (B1 - B2)/( (B1 + B2)/2.0 );
257
258 // calculate the gas formation volume factor which exhibits the same
259 // "compressibility" for the new value of Rw
260 Scalar newBg = gasBArray.back()*(1.0 + x/2.0)/(1.0 - x/2.0);
261
262 // calculate the "viscosibility" of the master table
263 Scalar mu1 = viscosityColumn[newRowIdx];
264 Scalar mu2 = viscosityColumn[newRowIdx - 1];
265 Scalar xMu = (mu1 - mu2)/( (mu1 + mu2)/2.0 );
266
267 // calculate the viscosity which exhibits the same
268 // "viscosibility" for the new Rw value
269 Scalar newMug = gasMuArray.back()*(1.0 + xMu/2)/(1.0 - xMu/2.0);
270
271 // append the new values to the arrays which we use to compute the additional
272 // values ...
273 RwArray.push_back(newRw);
274 gasBArray.push_back(newBg);
275 gasMuArray.push_back(newMug);
276
277 // ... and register them with the internal table objects
278 invGasB.appendSamplePoint(xIdx, newRw, 1.0/newBg);
279 gasMu.appendSamplePoint(xIdx, newRw, newMug);
280 }
281 }
282
283public:
284#endif // HAVE_ECL_INPUT
285
287 {
288 waterReferenceDensity_.resize(numRegions);
289 gasReferenceDensity_.resize(numRegions);
290 inverseGasB_.resize(numRegions, TabulatedTwoDFunction{TabulatedTwoDFunction::InterpolationPolicy::RightExtreme});
291 inverseGasBMu_.resize(numRegions, TabulatedTwoDFunction{TabulatedTwoDFunction::InterpolationPolicy::RightExtreme});
292 inverseSaturatedGasB_.resize(numRegions);
293 inverseSaturatedGasBMu_.resize(numRegions);
294 gasMu_.resize(numRegions, TabulatedTwoDFunction{TabulatedTwoDFunction::InterpolationPolicy::RightExtreme});
295 saturatedWaterVaporizationFactorTable_.resize(numRegions);
296 saturationPressure_.resize(numRegions);
297 }
298
302 void setReferenceDensities(unsigned regionIdx,
303 Scalar /*rhoRefOil*/,
304 Scalar rhoRefGas,
305 Scalar rhoRefWater)
306 {
307 waterReferenceDensity_[regionIdx] = rhoRefWater;
308 gasReferenceDensity_[regionIdx] = rhoRefGas;
309 }
310
316 void setSaturatedGasWaterVaporizationFactor(unsigned regionIdx, const SamplingPoints& samplePoints)
317 { saturatedWaterVaporizationFactorTable_[regionIdx].setContainerOfTuples(samplePoints); }
318
331 void setInverseGasFormationVolumeFactor(unsigned regionIdx, const TabulatedTwoDFunction& invBg)
332 { inverseGasB_[regionIdx] = invBg; }
333
339 void setGasViscosity(unsigned regionIdx, const TabulatedTwoDFunction& mug)
340 { gasMu_[regionIdx] = mug; }
341
349 void setSaturatedGasViscosity(unsigned regionIdx, const SamplingPoints& samplePoints )
350 {
351 auto& waterVaporizationFac = saturatedWaterVaporizationFactorTable_[regionIdx];
352
353 constexpr const Scalar RwMin = 0.0;
354 Scalar RwMax = waterVaporizationFac.eval(saturatedWaterVaporizationFactorTable_[regionIdx].xMax(), /*extrapolate=*/true);
355
356 Scalar poMin = samplePoints.front().first;
357 Scalar poMax = samplePoints.back().first;
358
359 constexpr const size_t nRw = 20;
360 size_t nP = samplePoints.size()*2;
361
362 TabulatedOneDFunction mugTable;
363 mugTable.setContainerOfTuples(samplePoints);
364
365 // calculate a table of estimated densities depending on pressure and gas mass
366 // fraction
367 for (size_t RwIdx = 0; RwIdx < nRw; ++RwIdx) {
368 Scalar Rw = RwMin + (RwMax - RwMin)*RwIdx/nRw;
369
370 gasMu_[regionIdx].appendXPos(Rw);
371
372 for (size_t pIdx = 0; pIdx < nP; ++pIdx) {
373 Scalar pg = poMin + (poMax - poMin)*pIdx/nP;
374 Scalar mug = mugTable.eval(pg, /*extrapolate=*/true);
375
376 gasMu_[regionIdx].appendSamplePoint(RwIdx, pg, mug);
377 }
378 }
379 }
380
384 void initEnd()
385 {
386 // calculate the final 2D functions which are used for interpolation.
387 size_t numRegions = gasMu_.size();
388 for (unsigned regionIdx = 0; regionIdx < numRegions; ++ regionIdx) {
389 // calculate the table which stores the inverse of the product of the gas
390 // formation volume factor and the gas viscosity
391 const auto& gasMu = gasMu_[regionIdx];
392 const auto& invGasB = inverseGasB_[regionIdx];
393 assert(gasMu.numX() == invGasB.numX());
394
395 auto& invGasBMu = inverseGasBMu_[regionIdx];
396 auto& invSatGasB = inverseSaturatedGasB_[regionIdx];
397 auto& invSatGasBMu = inverseSaturatedGasBMu_[regionIdx];
398
399 std::vector<Scalar> satPressuresArray;
400 std::vector<Scalar> invSatGasBArray;
401 std::vector<Scalar> invSatGasBMuArray;
402 for (size_t pIdx = 0; pIdx < gasMu.numX(); ++pIdx) {
403 invGasBMu.appendXPos(gasMu.xAt(pIdx));
404
405 assert(gasMu.numY(pIdx) == invGasB.numY(pIdx));
406
407 size_t numRw = gasMu.numY(pIdx);
408 for (size_t RwIdx = 0; RwIdx < numRw; ++RwIdx)
409 invGasBMu.appendSamplePoint(pIdx,
410 gasMu.yAt(pIdx, RwIdx),
411 invGasB.valueAt(pIdx, RwIdx)
412 / gasMu.valueAt(pIdx, RwIdx));
413
414 // the sampling points in UniformXTabulated2DFunction are always sorted
415 // in ascending order. Thus, the value for saturated gas is the last one
416 // (i.e., the one with the largest Rw value)
417 satPressuresArray.push_back(gasMu.xAt(pIdx));
418 invSatGasBArray.push_back(invGasB.valueAt(pIdx, numRw - 1));
419 invSatGasBMuArray.push_back(invGasBMu.valueAt(pIdx, numRw - 1));
420 }
421
422 invSatGasB.setXYContainers(satPressuresArray, invSatGasBArray);
423 invSatGasBMu.setXYContainers(satPressuresArray, invSatGasBMuArray);
424
425 updateSaturationPressure_(regionIdx);
426 }
427 }
428
432 unsigned numRegions() const
433 { return gasReferenceDensity_.size(); }
434
438 template <class Evaluation>
439 Evaluation internalEnergy(unsigned,
440 const Evaluation&,
441 const Evaluation&,
442 const Evaluation&) const
443 {
444 throw std::runtime_error("Requested the enthalpy of gas but the thermal option is not enabled");
445 }
446
450 template <class Evaluation>
451 Evaluation viscosity(unsigned regionIdx,
452 const Evaluation& /*temperature*/,
453 const Evaluation& pressure,
454 const Evaluation& /*Rv*/,
455 const Evaluation& Rvw) const
456 {
457 const Evaluation& invBg = inverseGasB_[regionIdx].eval(pressure, Rvw, /*extrapolate=*/true);
458 const Evaluation& invMugBg = inverseGasBMu_[regionIdx].eval(pressure, Rvw, /*extrapolate=*/true);
459
460 return invBg/invMugBg;
461 }
462
466 template <class Evaluation>
467 Evaluation saturatedViscosity(unsigned regionIdx,
468 const Evaluation& /*temperature*/,
469 const Evaluation& pressure) const
470 {
471 const Evaluation& invBg = inverseSaturatedGasB_[regionIdx].eval(pressure, /*extrapolate=*/true);
472 const Evaluation& invMugBg = inverseSaturatedGasBMu_[regionIdx].eval(pressure, /*extrapolate=*/true);
473
474 return invBg/invMugBg;
475 }
476
480 template <class Evaluation>
481 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
482 const Evaluation& /*temperature*/,
483 const Evaluation& pressure,
484 const Evaluation& /*Rv*/,
485 const Evaluation& Rvw) const
486 { return inverseGasB_[regionIdx].eval(pressure, Rvw, /*extrapolate=*/true); }
487
491 template <class Evaluation>
492 Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
493 const Evaluation& /*temperature*/,
494 const Evaluation& pressure) const
495 { return inverseSaturatedGasB_[regionIdx].eval(pressure, /*extrapolate=*/true); }
496
500 template <class Evaluation>
501 Evaluation saturatedWaterVaporizationFactor(unsigned regionIdx,
502 const Evaluation& /*temperature*/,
503 const Evaluation& pressure) const
504 {
505 return saturatedWaterVaporizationFactorTable_[regionIdx].eval(pressure, /*extrapolate=*/true);
506 }
507
511 template <class Evaluation>
512 Evaluation saturatedWaterVaporizationFactor(unsigned regionIdx,
513 const Evaluation& /*temperature*/,
514 const Evaluation& pressure,
515 const Evaluation& saltConcentration) const
516 {
517 if (enableRwgSalt_)
518 return saturatedWaterVaporizationSaltFactorTable_[regionIdx].eval(pressure, saltConcentration, /*extrapolate=*/true);
519 else {
520 return saturatedWaterVaporizationFactorTable_[regionIdx].eval(pressure, /*extrapolate=*/true);
521 }
522 }
523
527 template <class Evaluation>
528 Evaluation saturatedOilVaporizationFactor(unsigned /*regionIdx*/,
529 const Evaluation& /*temperature*/,
530 const Evaluation& /*pressure*/,
531 const Evaluation& /*oilSaturation*/,
532 const Evaluation& /*maxOilSaturation*/) const
533 { return 0.0; /* this is dry humid gas! */ }
534
538 template <class Evaluation>
539 Evaluation saturatedOilVaporizationFactor(unsigned /*regionIdx*/,
540 const Evaluation& /*temperature*/,
541 const Evaluation& /*pressure*/) const
542 { return 0.0; /* this is dry humid gas! */ }
543
551 template <class Evaluation>
552 Evaluation saturationPressure(unsigned regionIdx,
553 const Evaluation&,
554 const Evaluation& Rw) const
555 {
556 typedef MathToolbox<Evaluation> Toolbox;
557
558 const auto& RwTable = saturatedWaterVaporizationFactorTable_[regionIdx];
559 const Scalar eps = std::numeric_limits<typename Toolbox::Scalar>::epsilon()*1e6;
560
561 // use the tabulated saturation pressure function to get a pretty good initial value
562 Evaluation pSat = saturationPressure_[regionIdx].eval(Rw, /*extrapolate=*/true);
563
564 // Newton method to do the remaining work. If the initial
565 // value is good, this should only take two to three
566 // iterations...
567 bool onProbation = false;
568 for (unsigned i = 0; i < 20; ++i) {
569 const Evaluation& f = RwTable.eval(pSat, /*extrapolate=*/true) - Rw;
570 const Evaluation& fPrime = RwTable.evalDerivative(pSat, /*extrapolate=*/true);
571
572 // If the derivative is "zero" Newton will not converge,
573 // so simply return our initial guess.
574 if (std::abs(scalarValue(fPrime)) < 1.0e-30) {
575 return pSat;
576 }
577
578 const Evaluation& delta = f/fPrime;
579
580 pSat -= delta;
581
582 if (pSat < 0.0) {
583 // if the pressure is lower than 0 Pascals, we set it back to 0. if this
584 // happens twice, we give up and just return 0 Pa...
585 if (onProbation)
586 return 0.0;
587
588 onProbation = true;
589 pSat = 0.0;
590 }
591
592 if (std::abs(scalarValue(delta)) < std::abs(scalarValue(pSat))*eps)
593 return pSat;
594 }
595
596 std::stringstream errlog;
597 errlog << "Finding saturation pressure did not converge:"
598 << " pSat = " << pSat
599 << ", Rw = " << Rw;
600#if HAVE_OPM_COMMON
601 OpmLog::debug("Wet gas saturation pressure", errlog.str());
602#endif
603 throw NumericalIssue(errlog.str());
604 }
605
606 template <class Evaluation>
607 Evaluation diffusionCoefficient(const Evaluation& /*temperature*/,
608 const Evaluation& /*pressure*/,
609 unsigned /*compIdx*/) const
610 {
611 throw std::runtime_error("Not implemented: The PVT model does not provide a diffusionCoefficient()");
612 }
613
614 const Scalar gasReferenceDensity(unsigned regionIdx) const
615 { return gasReferenceDensity_[regionIdx]; }
616
617 const Scalar waterReferenceDensity(unsigned regionIdx) const
618 { return waterReferenceDensity_[regionIdx]; }
619
620 const std::vector<TabulatedTwoDFunction>& inverseGasB() const {
621 return inverseGasB_;
622 }
623
624 const std::vector<TabulatedOneDFunction>& inverseSaturatedGasB() const {
625 return inverseSaturatedGasB_;
626 }
627
628 const std::vector<TabulatedTwoDFunction>& gasMu() const {
629 return gasMu_;
630 }
631
632 const std::vector<TabulatedTwoDFunction>& inverseGasBMu() const {
633 return inverseGasBMu_;
634 }
635
636 const std::vector<TabulatedOneDFunction>& inverseSaturatedGasBMu() const {
637 return inverseSaturatedGasBMu_;
638 }
639
640 const std::vector<TabulatedOneDFunction>& saturatedWaterVaporizationFactorTable() const {
641 return saturatedWaterVaporizationFactorTable_;
642 }
643
644 const std::vector<TabulatedTwoDFunction>& saturatedWaterVaporizationSaltFactorTable() const {
645 return saturatedWaterVaporizationSaltFactorTable_;
646 }
647
648 const std::vector<TabulatedOneDFunction>& saturationPressure() const {
649 return saturationPressure_;
650 }
651
652 Scalar vapPar1() const {
653 return vapPar1_;
654 }
655
656 bool operator==(const DryHumidGasPvt<Scalar>& data) const
657 {
658 return this->gasReferenceDensity_ == data.gasReferenceDensity_ &&
659 this->waterReferenceDensity_ == data.waterReferenceDensity_ &&
660 this->inverseGasB() == data.inverseGasB() &&
661 this->inverseSaturatedGasB() == data.inverseSaturatedGasB() &&
662 this->gasMu() == data.gasMu() &&
663 this->inverseGasBMu() == data.inverseGasBMu() &&
666 this->saturationPressure() == data.saturationPressure() &&
667 this->vapPar1() == data.vapPar1();
668 }
669
670private:
671 void updateSaturationPressure_(unsigned regionIdx)
672 {
673 typedef std::pair<Scalar, Scalar> Pair;
674 const auto& waterVaporizationFac = saturatedWaterVaporizationFactorTable_[regionIdx];
675
676 // create the taublated function representing saturation pressure depending of
677 // Rw
678 size_t n = waterVaporizationFac.numSamples();
679 Scalar delta = (waterVaporizationFac.xMax() - waterVaporizationFac.xMin())/Scalar(n + 1);
680
681 SamplingPoints pSatSamplePoints;
682 Scalar Rw = 0;
683 for (size_t i = 0; i <= n; ++ i) {
684 Scalar pSat = waterVaporizationFac.xMin() + Scalar(i)*delta;
685 Rw = saturatedWaterVaporizationFactor(regionIdx, /*temperature=*/Scalar(1e30), pSat);
686
687 Pair val(Rw, pSat);
688 pSatSamplePoints.push_back(val);
689 }
690
691 //Prune duplicate Rv values (can occur, and will cause problems in further interpolation)
692 auto x_coord_comparator = [](const Pair& a, const Pair& b) { return a.first == b.first; };
693 auto last = std::unique(pSatSamplePoints.begin(), pSatSamplePoints.end(), x_coord_comparator);
694 if (std::distance(pSatSamplePoints.begin(), last) > 1) // only remove them if there are more than two points
695 pSatSamplePoints.erase(last, pSatSamplePoints.end());
696
697 saturationPressure_[regionIdx].setContainerOfTuples(pSatSamplePoints);
698 }
699
700 std::vector<Scalar> gasReferenceDensity_;
701 std::vector<Scalar> waterReferenceDensity_;
702 std::vector<TabulatedTwoDFunction> inverseGasB_;
703 std::vector<TabulatedOneDFunction> inverseSaturatedGasB_;
704 std::vector<TabulatedTwoDFunction> gasMu_;
705 std::vector<TabulatedTwoDFunction> inverseGasBMu_;
706 std::vector<TabulatedOneDFunction> inverseSaturatedGasBMu_;
707 std::vector<TabulatedOneDFunction> saturatedWaterVaporizationFactorTable_;
708 std::vector<TabulatedTwoDFunction> saturatedWaterVaporizationSaltFactorTable_;
709 std::vector<TabulatedOneDFunction> saturationPressure_;
710
711 bool enableRwgSalt_;
712 Scalar vapPar1_;
713};
714
715} // namespace Opm
716
717#endif
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
This class represents the Pressure-Volume-Temperature relations of the gas phase with vaporized water...
Definition: DryHumidGasPvt.hpp:53
Scalar vapPar1() const
Definition: DryHumidGasPvt.hpp:652
const Scalar gasReferenceDensity(unsigned regionIdx) const
Definition: DryHumidGasPvt.hpp:614
Evaluation inverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure, const Evaluation &, const Evaluation &Rvw) const
Returns the formation volume factor [-] of the fluid phase.
Definition: DryHumidGasPvt.hpp:481
Evaluation saturatedOilVaporizationFactor(unsigned, const Evaluation &, const Evaluation &, const Evaluation &, const Evaluation &) const
Returns the oil vaporization factor [m^3/m^3] of the oil phase.
Definition: DryHumidGasPvt.hpp:528
DryHumidGasPvt()
Definition: DryHumidGasPvt.hpp:60
const std::vector< TabulatedOneDFunction > & inverseSaturatedGasB() const
Definition: DryHumidGasPvt.hpp:624
void setGasViscosity(unsigned regionIdx, const TabulatedTwoDFunction &mug)
Initialize the viscosity of the gas phase.
Definition: DryHumidGasPvt.hpp:339
Evaluation saturationPressure(unsigned regionIdx, const Evaluation &, const Evaluation &Rw) const
Returns the saturation pressure of the gas phase [Pa] depending on its mass fraction of the water com...
Definition: DryHumidGasPvt.hpp:552
DryHumidGasPvt(const std::vector< Scalar > &gasReferenceDensity, const std::vector< Scalar > &waterReferenceDensity, const std::vector< TabulatedTwoDFunction > &inverseGasB, const std::vector< TabulatedOneDFunction > &inverseSaturatedGasB, const std::vector< TabulatedTwoDFunction > &gasMu, const std::vector< TabulatedTwoDFunction > &inverseGasBMu, const std::vector< TabulatedOneDFunction > &inverseSaturatedGasBMu, const std::vector< TabulatedOneDFunction > &saturatedWaterVaporizationFactorTable, const std::vector< TabulatedOneDFunction > &saturationPressure, Scalar vapPar1)
Definition: DryHumidGasPvt.hpp:65
bool operator==(const DryHumidGasPvt< Scalar > &data) const
Definition: DryHumidGasPvt.hpp:656
void setInverseGasFormationVolumeFactor(unsigned regionIdx, const TabulatedTwoDFunction &invBg)
Initialize the function for the gas formation volume factor.
Definition: DryHumidGasPvt.hpp:331
const std::vector< TabulatedTwoDFunction > & saturatedWaterVaporizationSaltFactorTable() const
Definition: DryHumidGasPvt.hpp:644
void initEnd()
Finish initializing the gas phase PVT properties.
Definition: DryHumidGasPvt.hpp:384
unsigned numRegions() const
Return the number of PVT regions which are considered by this PVT-object.
Definition: DryHumidGasPvt.hpp:432
void setReferenceDensities(unsigned regionIdx, Scalar, Scalar rhoRefGas, Scalar rhoRefWater)
Initialize the reference densities of all fluids for a given PVT region.
Definition: DryHumidGasPvt.hpp:302
const std::vector< TabulatedTwoDFunction > & inverseGasB() const
Definition: DryHumidGasPvt.hpp:620
const std::vector< TabulatedTwoDFunction > & gasMu() const
Definition: DryHumidGasPvt.hpp:628
Evaluation saturatedWaterVaporizationFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure) const
Returns the water vaporization factor [m^3/m^3] of the water phase.
Definition: DryHumidGasPvt.hpp:501
const std::vector< TabulatedOneDFunction > & inverseSaturatedGasBMu() const
Definition: DryHumidGasPvt.hpp:636
void setSaturatedGasWaterVaporizationFactor(unsigned regionIdx, const SamplingPoints &samplePoints)
Initialize the function for the oil vaporization factor .
Definition: DryHumidGasPvt.hpp:316
const std::vector< TabulatedTwoDFunction > & inverseGasBMu() const
Definition: DryHumidGasPvt.hpp:632
Evaluation viscosity(unsigned regionIdx, const Evaluation &, const Evaluation &pressure, const Evaluation &, const Evaluation &Rvw) const
Returns the dynamic viscosity [Pa s] of the fluid phase given a set of parameters.
Definition: DryHumidGasPvt.hpp:451
void setSaturatedGasViscosity(unsigned regionIdx, const SamplingPoints &samplePoints)
Initialize the phase viscosity for oil saturated gas.
Definition: DryHumidGasPvt.hpp:349
const Scalar waterReferenceDensity(unsigned regionIdx) const
Definition: DryHumidGasPvt.hpp:617
Evaluation internalEnergy(unsigned, const Evaluation &, const Evaluation &, const Evaluation &) const
Returns the specific enthalpy [J/kg] of gas given a set of parameters.
Definition: DryHumidGasPvt.hpp:439
Evaluation saturatedOilVaporizationFactor(unsigned, const Evaluation &, const Evaluation &) const
Returns the oil vaporization factor [m^3/m^3] of the oil phase.
Definition: DryHumidGasPvt.hpp:539
Evaluation diffusionCoefficient(const Evaluation &, const Evaluation &, unsigned) const
Definition: DryHumidGasPvt.hpp:607
Evaluation saturatedWaterVaporizationFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure, const Evaluation &saltConcentration) const
Returns the water vaporization factor [m^3/m^3] of the water phase.
Definition: DryHumidGasPvt.hpp:512
const std::vector< TabulatedOneDFunction > & saturatedWaterVaporizationFactorTable() const
Definition: DryHumidGasPvt.hpp:640
void setNumRegions(size_t numRegions)
Definition: DryHumidGasPvt.hpp:286
const std::vector< TabulatedOneDFunction > & saturationPressure() const
Definition: DryHumidGasPvt.hpp:648
Evaluation saturatedViscosity(unsigned regionIdx, const Evaluation &, const Evaluation &pressure) const
Returns the dynamic viscosity [Pa s] of oil saturated gas at a given pressure.
Definition: DryHumidGasPvt.hpp:467
Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure) const
Returns the formation volume factor [-] of water saturated gas at a given pressure.
Definition: DryHumidGasPvt.hpp:492
Definition: Exceptions.hpp:46
Implements a linearly interpolated scalar function that depends on one variable.
Definition: Tabulated1DFunction.hpp:48
Evaluation eval(const Evaluation &x, bool extrapolate=false) const
Evaluate the spline at a given position.
Definition: Tabulated1DFunction.hpp:259
void setContainerOfTuples(const XYContainer &points, bool sortInputs=true)
Set the sampling points of the piecewise linear function using a STL-compatible container of tuple-li...
Definition: Tabulated1DFunction.hpp:185
Implements a scalar function that depends on two variables and which is sampled uniformly in the X di...
Definition: UniformXTabulated2DFunction.hpp:54
Definition: Air_Mesitylene.hpp:34
auto scalarValue(const Evaluation &val) -> decltype(MathToolbox< Evaluation >::scalarValue(val))
Definition: MathToolbox.hpp:335
Evaluation abs(const Evaluation &value)
Definition: MathToolbox.hpp:350
Definition: MathToolbox.hpp:50