opm-common
Loading...
Searching...
No Matches
WetGasPvt.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_GAS_PVT_HPP
28#define OPM_WET_GAS_PVT_HPP
29
31#include <opm/common/OpmLog/OpmLog.hpp>
32
36
37#include <cstddef>
38
39namespace Opm {
40
41class EclipseState;
42class Schedule;
43class SimpleTable;
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
63 void initFromState(const EclipseState& eclState, const Schedule& schedule);
64
65private:
66 void extendPvtgTable_(unsigned regionIdx,
67 unsigned xIdx,
68 const SimpleTable& curTable,
69 const SimpleTable& masterTable);
70
71public:
72 void setNumRegions(std::size_t numRegions);
73
74 void setVapPars(const Scalar par1, const Scalar)
75 {
76 vapPar1_ = par1;
77 }
78
82 void setReferenceDensities(unsigned regionIdx,
83 Scalar rhoRefOil,
84 Scalar rhoRefGas,
85 Scalar /*rhoRefWater*/);
86
93 void setSaturatedGasOilVaporizationFactor(unsigned regionIdx,
94 const SamplingPoints& samplePoints)
95 { saturatedOilVaporizationFactorTable_[regionIdx].setContainerOfTuples(samplePoints); }
96
106 void setSaturatedGasFormationVolumeFactor(unsigned regionIdx,
107 const SamplingPoints& samplePoints);
108
121 void setInverseGasFormationVolumeFactor(unsigned regionIdx,
122 const TabulatedTwoDFunction& invBg)
123 { inverseGasB_[regionIdx] = invBg; }
124
130 void setGasViscosity(unsigned regionIdx, const TabulatedTwoDFunction& mug)
131 { gasMu_[regionIdx] = mug; }
132
140 void setSaturatedGasViscosity(unsigned regionIdx,
141 const SamplingPoints& samplePoints);
142
146 void initEnd();
147
151 unsigned numRegions() const
152 { return gasReferenceDensity_.size(); }
153
157 template <class Evaluation>
158 Evaluation internalEnergy(unsigned,
159 const Evaluation&,
160 const Evaluation&,
161 const Evaluation&,
162 const Evaluation&) const
163 {
164 throw std::runtime_error("Requested the enthalpy of gas 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& Rv,
182 const Evaluation& /*Rvw*/) const
183 {
184 const Evaluation& invBg = inverseGasB_[regionIdx].eval(pressure, Rv, /*extrapolate=*/true);
185 const Evaluation& invMugBg = inverseGasBMu_[regionIdx].eval(pressure, Rv, /*extrapolate=*/true);
186
187 return invBg / invMugBg;
188 }
189
193 template <class Evaluation>
194 Evaluation saturatedViscosity(unsigned regionIdx,
195 const Evaluation& /*temperature*/,
196 const Evaluation& pressure) const
197 {
198 const Evaluation& invBg = inverseSaturatedGasB_[regionIdx].eval(pressure, /*extrapolate=*/true);
199 const Evaluation& invMugBg = inverseSaturatedGasBMu_[regionIdx].eval(pressure, /*extrapolate=*/true);
200
201 return invBg / invMugBg;
202 }
203
207 template <class Evaluation>
208 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
209 const Evaluation& /*temperature*/,
210 const Evaluation& pressure,
211 const Evaluation& Rv,
212 const Evaluation& /*Rvw*/) const
213 { return inverseGasB_[regionIdx].eval(pressure, Rv, /*extrapolate=*/true); }
214
218 template <class FluidState, class LhsEval = typename FluidState::ValueType>
219 std::pair<LhsEval, LhsEval>
220 inverseFormationVolumeFactorAndViscosity(const FluidState& fluidState, unsigned regionIdx)
221 {
222 const LhsEval& p = decay<LhsEval>(fluidState.pressure(FluidState::gasPhaseIdx));
223 const LhsEval& Rv = decay<LhsEval>(fluidState.Rv());
224
225 const auto satSegIdx = this->saturatedOilVaporizationFactorTable_[regionIdx].findSegmentIndex(p, /*extrapolate=*/ true);
226 const auto RvSat = this->saturatedOilVaporizationFactorTable_[regionIdx].eval(p, SegmentIndex{satSegIdx});
227 const bool useSaturatedTables = (fluidState.saturation(FluidState::oilPhaseIdx) > 0.0) && (Rv >= (1.0 - 1e-10) * RvSat);
228
229 if (useSaturatedTables) {
230 const LhsEval b = this->inverseSaturatedGasB_[regionIdx].eval(p, SegmentIndex{satSegIdx});
231 const LhsEval invBMu = this->inverseSaturatedGasBMu_[regionIdx].eval(p, SegmentIndex{satSegIdx});
232 const LhsEval mu = b / invBMu;
233 return { b, mu };
234 } else {
235 unsigned ii, jj1, jj2;
236 LhsEval alpha, beta1, beta2;
237 this->inverseGasB_[regionIdx].findPoints(ii, jj1, jj2, alpha, beta1, beta2, p, Rv, /*extrapolate =*/ true);
238 const LhsEval b = this->inverseGasB_[regionIdx].eval(ii, jj1, jj2, alpha, beta1, beta2);
239 const LhsEval invBMu = this->inverseGasBMu_[regionIdx].eval(ii, jj1, jj2, alpha, beta1, beta2);
240 const LhsEval mu = b / invBMu;
241 return { b, mu };
242 }
243 }
244
248 template <class Evaluation>
249 Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
250 const Evaluation& /*temperature*/,
251 const Evaluation& pressure) const
252 { return inverseSaturatedGasB_[regionIdx].eval(pressure, /*extrapolate=*/true); }
253
257 template <class Evaluation>
258 Evaluation saturatedWaterVaporizationFactor(unsigned /*regionIdx*/,
259 const Evaluation& /*temperature*/,
260 const Evaluation& /*pressure*/) const
261 { return 0.0; /* this is non-humid gas! */ }
262
266 template <class Evaluation = Scalar>
267 Evaluation saturatedWaterVaporizationFactor(unsigned /*regionIdx*/,
268 const Evaluation& /*temperature*/,
269 const Evaluation& /*pressure*/,
270 const Evaluation& /*saltConcentration*/) const
271 { return 0.0; }
272
276 template <class Evaluation>
277 Evaluation saturatedOilVaporizationFactor(unsigned regionIdx,
278 const Evaluation& /*temperature*/,
279 const Evaluation& pressure) const
280 {
281 return saturatedOilVaporizationFactorTable_[regionIdx].eval(pressure, /*extrapolate=*/true);
282 }
283
291 template <class Evaluation>
292 Evaluation saturatedOilVaporizationFactor(unsigned regionIdx,
293 const Evaluation& /*temperature*/,
294 const Evaluation& pressure,
295 const Evaluation& oilSaturation,
296 Evaluation maxOilSaturation) const
297 {
298 Evaluation tmp =
299 saturatedOilVaporizationFactorTable_[regionIdx].eval(pressure, /*extrapolate=*/true);
300
301 // apply the vaporization parameters for the gas phase (cf. the Eclipse VAPPARS
302 // keyword)
303 maxOilSaturation = min(maxOilSaturation, Scalar(1.0));
304 if (vapPar1_ > 0.0 && maxOilSaturation > 0.01 && oilSaturation < maxOilSaturation) {
305 constexpr const Scalar eps = 0.001;
306 const Evaluation& So = max(oilSaturation, eps);
307 tmp *= max(1e-3, pow(So / maxOilSaturation, vapPar1_));
308 }
309
310 return tmp;
311 }
312
324 template <class Evaluation>
325 Evaluation saturationPressure(unsigned regionIdx,
326 const Evaluation&,
327 const Evaluation& Rv) const
328 {
329 using Toolbox = MathToolbox<Evaluation>;
330
331 const auto& RvTable = saturatedOilVaporizationFactorTable_[regionIdx];
332 constexpr const Scalar eps = std::numeric_limits<typename Toolbox::Scalar>::epsilon() * 1e6;
333
334 // use the tabulated saturation pressure function to get a pretty good initial value
335 Evaluation pSat = saturationPressure_[regionIdx].eval(Rv, /*extrapolate=*/true);
336
337 // Newton method to do the remaining work. If the initial
338 // value is good, this should only take two to three
339 // iterations...
340 bool onProbation = false;
341 for (unsigned i = 0; i < 20; ++i) {
342 const Evaluation& f = RvTable.eval(pSat, /*extrapolate=*/true) - Rv;
343 const Evaluation& fPrime = RvTable.evalDerivative(pSat, /*extrapolate=*/true);
344
345 // If the derivative is "zero" Newton will not converge,
346 // so simply return our initial guess.
347 if (std::abs(scalarValue(fPrime)) < 1.0e-30) {
348 return pSat;
349 }
350
351 const Evaluation& delta = f / fPrime;
352
353 pSat -= delta;
354
355 if (pSat < 0.0) {
356 // if the pressure is lower than 0 Pascals, we set it back to 0. if this
357 // happens twice, we give up and just return 0 Pa...
358 if (onProbation) {
359 return 0.0;
360 }
361
362 onProbation = true;
363 pSat = 0.0;
364 }
365
366 if (std::abs(scalarValue(delta)) < std::abs(scalarValue(pSat))*eps) {
367 return pSat;
368 }
369 }
370
371 const std::string msg =
372 "Finding saturation pressure did not converge: "
373 "pSat = " + std::to_string(getValue(pSat)) +
374 ", Rv = " + std::to_string(getValue(Rv));
375 OpmLog::debug("Wet gas saturation pressure", msg);
376 throw NumericalProblem(msg);
377 }
378
379 template <class Evaluation>
380 Evaluation diffusionCoefficient(const Evaluation& /*temperature*/,
381 const Evaluation& /*pressure*/,
382 unsigned /*compIdx*/) const
383 {
384 throw std::runtime_error("Not implemented: The PVT model does not provide "
385 "a diffusionCoefficient()");
386 }
387
388 Scalar gasReferenceDensity(unsigned regionIdx) const
389 { return gasReferenceDensity_[regionIdx]; }
390
391 Scalar oilReferenceDensity(unsigned regionIdx) const
392 { return oilReferenceDensity_[regionIdx]; }
393
394 const std::vector<TabulatedTwoDFunction>& inverseGasB() const
395 { return inverseGasB_; }
396
397 const std::vector<TabulatedOneDFunction>& inverseSaturatedGasB() const
398 { return inverseSaturatedGasB_; }
399
400 const std::vector<TabulatedTwoDFunction>& gasMu() const
401 { return gasMu_; }
402
403 const std::vector<TabulatedTwoDFunction>& inverseGasBMu() const
404 { return inverseGasBMu_; }
405
406 const std::vector<TabulatedOneDFunction>& inverseSaturatedGasBMu() const
407 { return inverseSaturatedGasBMu_; }
408
409 const std::vector<TabulatedOneDFunction>& saturatedOilVaporizationFactorTable() const
410 { return saturatedOilVaporizationFactorTable_; }
411
412 const std::vector<TabulatedOneDFunction>& saturationPressure() const
413 { return saturationPressure_; }
414
415 Scalar vapPar1() const
416 { return vapPar1_; }
417
418private:
419 void updateSaturationPressure_(unsigned regionIdx);
420
421 std::vector<Scalar> gasReferenceDensity_{};
422 std::vector<Scalar> oilReferenceDensity_{};
423 std::vector<TabulatedTwoDFunction> inverseGasB_{};
424 std::vector<TabulatedOneDFunction> inverseSaturatedGasB_{};
425 std::vector<TabulatedTwoDFunction> gasMu_{};
426 std::vector<TabulatedTwoDFunction> inverseGasBMu_{};
427 std::vector<TabulatedOneDFunction> inverseSaturatedGasBMu_{};
428 std::vector<TabulatedOneDFunction> saturatedOilVaporizationFactorTable_{};
429 std::vector<TabulatedOneDFunction> saturationPressure_{};
430
431 Scalar vapPar1_ = 0.0;
432};
433
434} // namespace Opm
435
436#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 phas with vaporized oil.
Definition WetGasPvt.hpp:51
Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure) const
Returns the formation volume factor [-] of oil saturated gas at a given pressure.
Definition WetGasPvt.hpp:249
Evaluation saturatedOilVaporizationFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure) const
Returns the oil vaporization factor [m^3/m^3] of the gas phase.
Definition WetGasPvt.hpp:277
void initEnd()
Finish initializing the gas phase PVT properties.
Definition WetGasPvt.cpp:327
Evaluation saturationPressure(unsigned regionIdx, const Evaluation &, const Evaluation &Rv) const
Returns the saturation pressure of the gas phase [Pa] depending on its mass fraction of the oil compo...
Definition WetGasPvt.hpp:325
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 WetGasPvt.hpp:158
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 WetGasPvt.hpp:292
void setGasViscosity(unsigned regionIdx, const TabulatedTwoDFunction &mug)
Initialize the viscosity of the gas phase.
Definition WetGasPvt.hpp:130
Evaluation saturatedWaterVaporizationFactor(unsigned, const Evaluation &, const Evaluation &) const
Returns the water vaporization factor [m^3/m^3] of the gasphase.
Definition WetGasPvt.hpp:258
Evaluation saturatedWaterVaporizationFactor(unsigned, const Evaluation &, const Evaluation &, const Evaluation &) const
Returns the water vaporization factor [m^3/m^3] of water saturated gas.
Definition WetGasPvt.hpp:267
unsigned numRegions() const
Return the number of PVT regions which are considered by this PVT-object.
Definition WetGasPvt.hpp:151
void setSaturatedGasViscosity(unsigned regionIdx, const SamplingPoints &samplePoints)
Initialize the phase viscosity for oil saturated gas.
Definition WetGasPvt.cpp:293
void setSaturatedGasFormationVolumeFactor(unsigned regionIdx, const SamplingPoints &samplePoints)
Initialize the function for the gas formation volume factor.
Definition WetGasPvt.cpp:242
Evaluation inverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure, const Evaluation &Rv, const Evaluation &) const
Returns the formation volume factor [-] of the fluid phase.
Definition WetGasPvt.hpp:208
void setSaturatedGasOilVaporizationFactor(unsigned regionIdx, const SamplingPoints &samplePoints)
Initialize the function for the oil vaporization factor .
Definition WetGasPvt.hpp:93
std::pair< LhsEval, LhsEval > inverseFormationVolumeFactorAndViscosity(const FluidState &fluidState, unsigned regionIdx)
Returns the formation volume factor [-] and viscosity [Pa s] of the fluid phase.
Definition WetGasPvt.hpp:220
void setInverseGasFormationVolumeFactor(unsigned regionIdx, const TabulatedTwoDFunction &invBg)
Initialize the function for the gas formation volume factor.
Definition WetGasPvt.hpp:121
void initFromState(const EclipseState &eclState, const Schedule &schedule)
Initialize the parameters for wet gas using an ECL deck.
Definition WetGasPvt.cpp:40
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 WetGasPvt.hpp:194
Evaluation viscosity(unsigned regionIdx, const Evaluation &, const Evaluation &pressure, const Evaluation &Rv, const Evaluation &) const
Returns the dynamic viscosity [Pa s] of the fluid phase given a set of parameters.
Definition WetGasPvt.hpp:178
void setReferenceDensities(unsigned regionIdx, Scalar rhoRefOil, Scalar rhoRefGas, Scalar)
Initialize the reference densities of all fluids for a given PVT region.
Definition WetGasPvt.cpp:231
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