opm-common
Loading...
Searching...
No Matches
WetHumidGasPvt.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_WET_HUMID_GAS_PVT_HPP
28#define OPM_WET_HUMID_GAS_PVT_HPP
29
31#include <opm/common/OpmLog/OpmLog.hpp>
32
37
38#include <cstddef>
39
40namespace Opm {
41
42class EclipseState;
43class Schedule;
44class SimpleTable;
45
50template <class Scalar>
52{
53 using SamplingPoints = std::vector<std::pair<Scalar, Scalar>>;
54
55public:
56 using TabulatedTwoDFunction = UniformXTabulated2DFunction<Scalar>;
57 using TabulatedOneDFunction = Tabulated1DFunction<Scalar>;
58
64 void initFromState(const EclipseState& eclState, const Schedule& schedule);
65
66private:
67 void extendPvtgwTable_(unsigned regionIdx,
68 unsigned xIdx,
69 const SimpleTable& curTable,
70 const SimpleTable& masterTable);
71
72 void extendPvtgTable_(unsigned regionIdx,
73 unsigned xIdx,
74 const SimpleTable& curTable,
75 const SimpleTable& masterTable);
76
77public:
78 void setNumRegions(std::size_t numRegions);
79
80 void setVapPars(const Scalar par1, const Scalar)
81 {
82 vapPar1_ = par1;
83 }
84
88 void setReferenceDensities(unsigned regionIdx,
89 Scalar rhoRefOil,
90 Scalar rhoRefGas,
91 Scalar rhoRefWater);
92
100 const SamplingPoints& samplePoints)
101 { saturatedWaterVaporizationFactorTable_[regionIdx].setContainerOfTuples(samplePoints); }
102
110 const SamplingPoints& samplePoints)
111 { saturatedOilVaporizationFactorTable_[regionIdx].setContainerOfTuples(samplePoints); }
112
116 void initEnd();
117
121 unsigned numRegions() const
122 { return gasReferenceDensity_.size(); }
123
127 template <class Evaluation>
128 Evaluation internalEnergy(unsigned,
129 const Evaluation&,
130 const Evaluation&,
131 const Evaluation&,
132 const Evaluation&) const
133 {
134 throw std::runtime_error("Requested the enthalpy of gas but the thermal "
135 "option is not enabled");
136 }
137
138 Scalar hVap(unsigned) const
139 {
140 throw std::runtime_error("Requested the hvap of oil but the thermal "
141 "option is not enabled");
142 }
143
147 template <class Evaluation>
148 Evaluation viscosity(unsigned regionIdx,
149 const Evaluation& /*temperature*/,
150 const Evaluation& pressure,
151 const Evaluation& Rv,
152 const Evaluation& Rvw) const
153 {
154 const Evaluation& temperature = 1E30;
155
156 if (Rv >= (1.0 - 1e-10) * saturatedOilVaporizationFactor(regionIdx, temperature, pressure)) {
157 const Evaluation& invBg = inverseGasBRvSat_[regionIdx].eval(pressure, Rvw, /*extrapolate=*/true);
158 const Evaluation& invMugBg = inverseGasBMuRvSat_[regionIdx].eval(pressure, Rvw, /*extrapolate=*/true);
159 return invBg / invMugBg;
160 }
161 else {
162 // for Rv undersaturated viscosity is evaluated at saturated Rvw values
163 const Evaluation& invBg = inverseGasBRvwSat_[regionIdx].eval(pressure, Rv, /*extrapolate=*/true);
164 const Evaluation& invMugBg = inverseGasBMuRvwSat_[regionIdx].eval(pressure, Rv, /*extrapolate=*/true);
165 return invBg / invMugBg;
166 }
167 }
168
172 template <class Evaluation>
173 Evaluation saturatedViscosity(unsigned regionIdx,
174 const Evaluation& /*temperature*/,
175 const Evaluation& pressure) const
176 {
177 const Evaluation& invBg = inverseSaturatedGasB_[regionIdx].eval(pressure, /*extrapolate=*/true);
178 const Evaluation& invMugBg = inverseSaturatedGasBMu_[regionIdx].eval(pressure, /*extrapolate=*/true);
179
180 return invBg / invMugBg;
181 }
182
186 // template <class Evaluation>
187 // Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
188 // const Evaluation& /*temperature*/,
189 // const Evaluation& pressure,
190 // const Evaluation& Rw) const
191 // { return inverseGasB_[regionIdx].eval(pressure, Rw, /*extrapolate=*/true); }
192
193 template <class Evaluation>
194 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
195 const Evaluation& /*temperature*/,
196 const Evaluation& pressure,
197 const Evaluation& Rv,
198 const Evaluation& Rvw) const
199 {
200 const Evaluation& temperature = 1E30;
201
202 if (Rv >= (1.0 - 1e-10) * saturatedOilVaporizationFactor(regionIdx, temperature, pressure)) {
203 return inverseGasBRvSat_[regionIdx].eval(pressure, Rvw, /*extrapolate=*/true);
204 }
205 else {
206 // for Rv undersaturated Bg^-1 is evaluated at saturated Rvw values
207 return inverseGasBRvwSat_[regionIdx].eval(pressure, Rv, /*extrapolate=*/true);
208 }
209 }
210
214 template <class FluidState, class LhsEval = typename FluidState::ValueType>
215 std::pair<LhsEval, LhsEval>
216 inverseFormationVolumeFactorAndViscosity(const FluidState& fluidState, unsigned regionIdx)
217 {
218 const LhsEval& p = decay<LhsEval>(fluidState.pressure(FluidState::gasPhaseIdx));
219 const LhsEval& Rv = decay<LhsEval>(fluidState.Rv());
220 const LhsEval& Rvw = decay<LhsEval>(fluidState.Rvw());
221 const LhsEval& saltConc
222 = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
223
224 // It is not guaranteed that the oil and water vaporization
225 // factor tables, and also the saturated B and Mu tables, have
226 // the same pressure sample points. Therefore we do not bother
227 // to separately call findSegmentIndex() and eval() here.
228 const auto RvSat = this->saturatedOilVaporizationFactorTable_[regionIdx].eval(p, /*extrapolate=*/ true);
229 // TODO: check that handling of salt concentration is correct, it seems to only affect the saturation curve.
230 const auto RvwSat = enableRwgSalt_
231 ? this->saturatedWaterVaporizationSaltFactorTable_[regionIdx].eval(p, saltConc, /*extrapolate=*/ true)
232 : this->saturatedWaterVaporizationFactorTable_[regionIdx].eval(p, /*extrapolate=*/true);
233
234 const bool waterSaturated = (fluidState.saturation(FluidState::waterPhaseIdx) > 0.0) && (Rvw >= (1.0 - 1e-10) * RvwSat);
235 const bool oilSaturated = (fluidState.saturation(FluidState::oilPhaseIdx) > 0.0) && (Rv >= (1.0 - 1e-10) * RvSat);
236
237 if (waterSaturated && oilSaturated) {
238 const auto satSegIdx = this->inverseSaturatedGasB_[regionIdx].findSegmentIndex(p, /*extrapolate=*/ true);
239 const LhsEval b = this->inverseSaturatedGasB_[regionIdx].eval(p, SegmentIndex{satSegIdx});
240 const LhsEval invBMu = this->inverseSaturatedGasBMu_[regionIdx].eval(p, SegmentIndex{satSegIdx});
241 const LhsEval mu = b / invBMu;
242 return { b, mu };
243 } else if (oilSaturated) {
244 unsigned ii, jj1, jj2;
245 LhsEval alpha, beta1, beta2;
246 this->inverseGasBRvSat_[regionIdx].findPoints(ii, jj1, jj2, alpha, beta1, beta2, p, Rvw, /*extrapolate =*/ true);
247 const LhsEval b = this->inverseGasBRvSat_[regionIdx].eval(ii, jj1, jj2, alpha, beta1, beta2);
248 const LhsEval invBMu = this->inverseGasBMuRvSat_[regionIdx].eval(ii, jj1, jj2, alpha, beta1, beta2);
249 const LhsEval mu = b / invBMu;
250 return { b, mu };
251 } else {
252 // At this point, we assume waterSaturated is true, but this is not checked.
253 unsigned ii, jj1, jj2;
254 LhsEval alpha, beta1, beta2;
255 this->inverseGasBRvwSat_[regionIdx].findPoints(ii, jj1, jj2, alpha, beta1, beta2, p, Rv, /*extrapolate =*/ true);
256 const LhsEval b = this->inverseGasBRvwSat_[regionIdx].eval(ii, jj1, jj2, alpha, beta1, beta2);
257 const LhsEval invBMu = this->inverseGasBMuRvwSat_[regionIdx].eval(ii, jj1, jj2, alpha, beta1, beta2);
258 const LhsEval mu = b / invBMu;
259 return { b, mu };
260 }
261 }
262
266 template <class Evaluation>
267 Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
268 const Evaluation& /*temperature*/,
269 const Evaluation& pressure) const
270 { return inverseSaturatedGasB_[regionIdx].eval(pressure, /*extrapolate=*/true); }
271
275 template <class Evaluation>
276 Evaluation saturatedWaterVaporizationFactor(unsigned regionIdx,
277 const Evaluation& /*temperature*/,
278 const Evaluation& pressure) const
279 {
280 return saturatedWaterVaporizationFactorTable_[regionIdx].eval(pressure, /*extrapolate=*/true);
281 }
282
286 template <class Evaluation>
287 Evaluation saturatedWaterVaporizationFactor(unsigned regionIdx,
288 const Evaluation& /*temperature*/,
289 const Evaluation& pressure,
290 const Evaluation& saltConcentration) const
291 {
292 if (enableRwgSalt_) {
293 return saturatedWaterVaporizationSaltFactorTable_[regionIdx].eval(pressure, saltConcentration, /*extrapolate=*/true);
294 }
295 else {
296 return saturatedWaterVaporizationFactorTable_[regionIdx].eval(pressure, /*extrapolate=*/true);
297 }
298 }
299
300 template <class Evaluation>
301 Evaluation saturatedOilVaporizationFactor(unsigned regionIdx,
302 const Evaluation& /*temperature*/,
303 const Evaluation& pressure) const
304 {
305 return saturatedOilVaporizationFactorTable_[regionIdx].eval(pressure, /*extrapolate=*/true);
306 }
307
315 template <class Evaluation>
316 Evaluation saturatedOilVaporizationFactor(unsigned regionIdx,
317 const Evaluation& /*temperature*/,
318 const Evaluation& pressure,
319 const Evaluation& oilSaturation,
320 Evaluation maxOilSaturation) const
321 {
322 Evaluation tmp =
323 saturatedOilVaporizationFactorTable_[regionIdx].eval(pressure, /*extrapolate=*/true);
324
325 // apply the vaporization parameters for the gas phase (cf. the Eclipse VAPPARS
326 // keyword)
327 maxOilSaturation = min(maxOilSaturation, Scalar(1.0));
328 if (vapPar1_ > 0.0 && maxOilSaturation > 0.01 && oilSaturation < maxOilSaturation) {
329 constexpr const Scalar eps = 0.001;
330 const Evaluation& So = max(oilSaturation, eps);
331 tmp *= max(1e-3, pow(So / maxOilSaturation, vapPar1_));
332 }
333
334 return tmp;
335 }
336
345 //PJPE assume dependence on Rv
346 template <class Evaluation>
347 Evaluation saturationPressure(unsigned regionIdx,
348 const Evaluation&,
349 const Evaluation& Rw) const
350 {
351 using Toolbox = MathToolbox<Evaluation>;
352
353 const auto& RwTable = saturatedWaterVaporizationFactorTable_[regionIdx];
354 constexpr const Scalar eps = std::numeric_limits<typename Toolbox::Scalar>::epsilon() * 1e6;
355
356 // use the tabulated saturation pressure function to get a pretty good initial value
357 Evaluation pSat = saturationPressure_[regionIdx].eval(Rw, /*extrapolate=*/true);
358
359 // Newton method to do the remaining work. If the initial
360 // value is good, this should only take two to three
361 // iterations...
362 bool onProbation = false;
363 for (unsigned i = 0; i < 20; ++i) {
364 const Evaluation& f = RwTable.eval(pSat, /*extrapolate=*/true) - Rw;
365 const Evaluation& fPrime = RwTable.evalDerivative(pSat, /*extrapolate=*/true);
366
367 // If the derivative is "zero" Newton will not converge,
368 // so simply return our initial guess.
369 if (std::abs(scalarValue(fPrime)) < 1.0e-30) {
370 return pSat;
371 }
372
373 const Evaluation& delta = f / fPrime;
374
375 pSat -= delta;
376
377 if (pSat < 0.0) {
378 // if the pressure is lower than 0 Pascals, we set it back to 0. if this
379 // happens twice, we give up and just return 0 Pa...
380 if (onProbation) {
381 return 0.0;
382 }
383
384 onProbation = true;
385 pSat = 0.0;
386 }
387
388 if (std::abs(scalarValue(delta)) < std::abs(scalarValue(pSat))*eps) {
389 return pSat;
390 }
391 }
392
393 const std::string msg =
394 "Finding saturation pressure did not converge: "
395 "pSat = " + std::to_string(getValue(pSat)) +
396 ", Rw = " + std::to_string(getValue(Rw));
397 OpmLog::debug("Wet gas saturation pressure", msg);
398 throw NumericalProblem(msg);
399 }
400
401 template <class Evaluation>
402 Evaluation diffusionCoefficient(const Evaluation& /*temperature*/,
403 const Evaluation& /*pressure*/,
404 unsigned /*compIdx*/) const
405 {
406 throw std::runtime_error("Not implemented: The PVT model does not provide "
407 "a diffusionCoefficient()");
408 }
409
410 Scalar gasReferenceDensity(unsigned regionIdx) const
411 { return gasReferenceDensity_[regionIdx]; }
412
413 Scalar oilReferenceDensity(unsigned regionIdx) const
414 { return oilReferenceDensity_[regionIdx]; }
415
416 Scalar waterReferenceDensity(unsigned regionIdx) const
417 { return waterReferenceDensity_[regionIdx]; }
418
419 const std::vector<TabulatedTwoDFunction>& inverseGasB() const
420 { return inverseGasBRvSat_; }
421
422 const std::vector<TabulatedOneDFunction>& inverseSaturatedGasB() const
423 { return inverseSaturatedGasB_; }
424
425 const std::vector<TabulatedTwoDFunction>& gasMu() const
426 { return gasMuRvSat_; }
427
428 const std::vector<TabulatedTwoDFunction>& inverseGasBMu() const
429 { return inverseGasBMuRvSat_; }
430
431 const std::vector<TabulatedOneDFunction>& inverseSaturatedGasBMu() const
432 { return inverseSaturatedGasBMu_; }
433
434 const std::vector<TabulatedOneDFunction>& saturatedWaterVaporizationFactorTable() const
435 { return saturatedWaterVaporizationFactorTable_; }
436
437 const std::vector<TabulatedTwoDFunction>& saturatedWaterVaporizationSaltFactorTable() const
438 { return saturatedWaterVaporizationSaltFactorTable_; }
439
440 const std::vector<TabulatedOneDFunction>& saturatedOilVaporizationFactorTable() const
441 { return saturatedOilVaporizationFactorTable_; }
442
443 const std::vector<TabulatedOneDFunction>& saturationPressure() const
444 { return saturationPressure_; }
445
446 Scalar vapPar1() const
447 { return vapPar1_; }
448
449private:
450 void updateSaturationPressure_(unsigned regionIdx);
451
452 std::vector<Scalar> gasReferenceDensity_{};
453 std::vector<Scalar> oilReferenceDensity_{};
454 std::vector<Scalar> waterReferenceDensity_{};
455 std::vector<TabulatedTwoDFunction> inverseGasBRvwSat_{};
456 std::vector<TabulatedTwoDFunction> inverseGasBRvSat_{};
457 std::vector<TabulatedOneDFunction> inverseSaturatedGasB_{};
458 std::vector<TabulatedTwoDFunction> gasMuRvwSat_{};
459 std::vector<TabulatedTwoDFunction> gasMuRvSat_{};
460 std::vector<TabulatedTwoDFunction> inverseGasBMuRvwSat_{};
461 std::vector<TabulatedTwoDFunction> inverseGasBMuRvSat_{};
462 std::vector<TabulatedOneDFunction> inverseSaturatedGasBMu_{};
463 std::vector<TabulatedOneDFunction> saturatedWaterVaporizationFactorTable_{};
464 std::vector<TabulatedTwoDFunction> saturatedWaterVaporizationSaltFactorTable_{};
465 std::vector<TabulatedOneDFunction> saturatedOilVaporizationFactorTable_{};
466 std::vector<TabulatedOneDFunction> saturationPressure_{};
467
468 bool enableRwgSalt_ = false;
469 Scalar vapPar1_ = 0.0;
470};
471
472} // namespace Opm
473
474#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
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 represents the Pressure-Volume-Temperature relations of the gas phase with vaporized oil a...
Definition WetHumidGasPvt.hpp:52
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 WetHumidGasPvt.hpp:347
std::pair< LhsEval, LhsEval > inverseFormationVolumeFactorAndViscosity(const FluidState &fluidState, unsigned regionIdx)
Returns the formation volume factor [-] and viscosity [Pa s] of the fluid phase.
Definition WetHumidGasPvt.hpp:216
unsigned numRegions() const
Return the number of PVT regions which are considered by this PVT-object.
Definition WetHumidGasPvt.hpp:121
void setSaturatedGasOilVaporizationFactor(unsigned regionIdx, const SamplingPoints &samplePoints)
Initialize the function for the oil vaporization factor .
Definition WetHumidGasPvt.hpp:109
Evaluation viscosity(unsigned regionIdx, const Evaluation &, const Evaluation &pressure, const Evaluation &Rv, const Evaluation &Rvw) const
Returns the dynamic viscosity [Pa s] of the fluid phase given a set of parameters.
Definition WetHumidGasPvt.hpp:148
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 WetHumidGasPvt.hpp:287
Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure) const
Returns the formation volume factor [-] of water saturated gas at a given pressure.
Definition WetHumidGasPvt.hpp:267
Evaluation saturatedWaterVaporizationFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure) const
Returns the water vaporization factor [m^3/m^3] of the water phase.
Definition WetHumidGasPvt.hpp:276
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 WetHumidGasPvt.hpp:173
Evaluation saturatedOilVaporizationFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure, const Evaluation &oilSaturation, Evaluation maxOilSaturation) const
Returns the oil vaporization factor [m^3/m^3] of the gas phase.
Definition WetHumidGasPvt.hpp:316
void setReferenceDensities(unsigned regionIdx, Scalar rhoRefOil, Scalar rhoRefGas, Scalar rhoRefWater)
Initialize the reference densities of all fluids for a given PVT region.
Definition WetHumidGasPvt.cpp:418
Evaluation inverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure, const Evaluation &Rv, const Evaluation &Rvw) const
Returns the formation volume factor [-] of the fluid phase.
Definition WetHumidGasPvt.hpp:194
Evaluation internalEnergy(unsigned, const Evaluation &, const Evaluation &, const Evaluation &, const Evaluation &) const
Returns the specific enthalpy [J/kg] of gas given a set of parameters.
Definition WetHumidGasPvt.hpp:128
void initFromState(const EclipseState &eclState, const Schedule &schedule)
Initialize the parameters for wet gas using an ECL deck.
Definition WetHumidGasPvt.cpp:39
void initEnd()
Finish initializing the gas phase PVT properties.
Definition WetHumidGasPvt.cpp:429
void setSaturatedGasWaterVaporizationFactor(unsigned regionIdx, const SamplingPoints &samplePoints)
Initialize the function for the water vaporization factor .
Definition WetHumidGasPvt.hpp:99
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