opm-common
Loading...
Searching...
No Matches
OilPvtThermal.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_OIL_PVT_THERMAL_HPP
28#define OPM_OIL_PVT_THERMAL_HPP
29
31
32#include <cstddef>
33
34namespace Opm {
35
36#if HAVE_ECL_INPUT
37class EclipseState;
38class Schedule;
39#endif
40
41template <class Scalar, bool enableThermal>
43
50template <class Scalar>
51class OilPvtThermal
52{
53public:
54 using TabulatedOneDFunction = Tabulated1DFunction<Scalar>;
55 using IsothermalPvt = OilPvtMultiplexer<Scalar, /*enableThermal=*/false>;
56
57 OilPvtThermal() = default;
58
59 OilPvtThermal(IsothermalPvt* isothermalPvt,
60 const std::vector<TabulatedOneDFunction>& oilvisctCurves,
61 const std::vector<Scalar>& viscrefPress,
62 const std::vector<Scalar>& viscrefRs,
63 const std::vector<Scalar>& viscRef,
64 const std::vector<Scalar>& oildentRefTemp,
65 const std::vector<Scalar>& oildentCT1,
66 const std::vector<Scalar>& oildentCT2,
67 const std::vector<Scalar>& oilJTRefPres,
68 const std::vector<Scalar>& oilJTC,
69 const std::vector<TabulatedOneDFunction>& internalEnergyCurves,
73 bool enableInternalEnergy)
74 : isothermalPvt_(isothermalPvt)
75 , oilvisctCurves_(oilvisctCurves)
76 , viscrefPress_(viscrefPress)
77 , viscrefRs_(viscrefRs)
78 , viscRef_(viscRef)
79 , oildentRefTemp_(oildentRefTemp)
80 , oildentCT1_(oildentCT1)
81 , oildentCT2_(oildentCT2)
82 , oilJTRefPres_(oilJTRefPres)
83 , oilJTC_(oilJTC)
84 , internalEnergyCurves_(internalEnergyCurves)
85 , enableThermalDensity_(enableThermalDensity)
86 , enableJouleThomson_(enableJouleThomson)
87 , enableThermalViscosity_(enableThermalViscosity)
88 , enableInternalEnergy_(enableInternalEnergy)
89 { }
90
91 OilPvtThermal(const OilPvtThermal& data)
92 { *this = data; }
93
94 ~OilPvtThermal()
95 { delete isothermalPvt_; }
96
97#if HAVE_ECL_INPUT
101 void initFromState(const EclipseState& eclState, const Schedule& schedule);
102#endif // HAVE_ECL_INPUT
103
107 void setNumRegions(std::size_t numRegions);
108
109 void setVapPars(const Scalar par1, const Scalar par2)
110 {
111 isothermalPvt_->setVapPars(par1, par2);
112 }
113
117 void initEnd()
118 { }
119
124 { return enableThermalDensity_; }
125
130 { return enableJouleThomson_; }
131
136 { return enableThermalViscosity_; }
137
138 std::size_t numRegions() const
139 { return viscrefRs_.size(); }
140
144 template <class Evaluation>
145 Evaluation internalEnergy(unsigned regionIdx,
146 const Evaluation& temperature,
147 const Evaluation& pressure,
148 const Evaluation& Rs) const
149 {
150 if (!enableInternalEnergy_) {
151 throw std::runtime_error("Requested the internal energy of oil but it is disabled");
152 }
153
154 if (!enableJouleThomson_) {
155 // compute the specific internal energy for the specified tempature. We use linear
156 // interpolation here despite the fact that the underlying heat capacities are
157 // piecewise linear (which leads to a quadratic function)
158 return internalEnergyCurves_[regionIdx].eval(temperature, /*extrapolate=*/true);
159 }
160 else {
161 OpmLog::warning("Experimental code for jouleThomson: simulation will be slower");
162 Evaluation Tref = oildentRefTemp_[regionIdx];
163 Evaluation Pref = oilJTRefPres_[regionIdx];
164 Scalar JTC = oilJTC_[regionIdx]; // if JTC is default then JTC is calculated
165
166 Evaluation invB = inverseFormationVolumeFactor(regionIdx, temperature, pressure, Rs);
167 Evaluation Cp = internalEnergyCurves_[regionIdx].eval(temperature, /*extrapolate=*/true)/temperature;
168 Evaluation density = invB * (oilReferenceDensity(regionIdx) + Rs * rhoRefG_[regionIdx]);
169
170 Evaluation enthalpyPres;
171 if (JTC != 0) {
172 enthalpyPres = -Cp * JTC * (pressure - Pref);
173 }
174 else if (enableThermalDensity_) {
175 Scalar c1T = oildentCT1_[regionIdx];
176 Scalar c2T = oildentCT2_[regionIdx];
177
178 Evaluation alpha = (c1T + 2 * c2T * (temperature - Tref)) /
179 (1 + c1T *(temperature - Tref) + c2T * (temperature - Tref) * (temperature - Tref));
180
181 const int N = 100; // value is experimental
182 Evaluation deltaP = (pressure - Pref) / N;
183 Evaluation enthalpyPresPrev = 0;
184 for (std::size_t i = 0; i < N; ++i) {
185 Evaluation Pnew = Pref + i * deltaP;
186 Evaluation rho = inverseFormationVolumeFactor(regionIdx, temperature, Pnew, Rs) *
187 (oilReferenceDensity(regionIdx) + Rs * rhoRefG_[regionIdx]) ;
188 // see e.g.https://en.wikipedia.org/wiki/Joule-Thomson_effect for a derivation of the Joule-Thomson coeff.
189 Evaluation jouleThomsonCoefficient = -(1.0/Cp) * (1.0 - alpha * temperature)/rho;
190 Evaluation deltaEnthalpyPres = -Cp * jouleThomsonCoefficient * deltaP;
191 enthalpyPres = enthalpyPresPrev + deltaEnthalpyPres;
192 enthalpyPresPrev = enthalpyPres;
193 }
194 }
195 else {
196 throw std::runtime_error("Requested Joule-thomson calculation but thermal oil"
197 "density (OILDENT) is not provided");
198 }
199
200 Evaluation enthalpy = Cp * (temperature - Tref) + enthalpyPres;
201
202 return enthalpy - pressure/density;
203 }
204 }
205
209 template <class Evaluation>
210 Evaluation viscosity(unsigned regionIdx,
211 const Evaluation& temperature,
212 const Evaluation& pressure,
213 const Evaluation& Rs) const
214 {
215 const auto& isothermalMu = isothermalPvt_->viscosity(regionIdx, temperature, pressure, Rs);
216 if (!enableThermalViscosity()) {
217 return isothermalMu;
218 }
219
220 // compute the viscosity deviation due to temperature
221 const auto& muOilvisct = oilvisctCurves_[regionIdx].eval(temperature, /*extrapolate=*/true);
222 return muOilvisct / viscRef_[regionIdx] * isothermalMu;
223 }
224
228 template <class Evaluation>
229 Evaluation saturatedViscosity(unsigned regionIdx,
230 const Evaluation& temperature,
231 const Evaluation& pressure) const
232 {
233 const auto& isothermalMu = isothermalPvt_->saturatedViscosity(regionIdx, temperature, pressure);
234 if (!enableThermalViscosity()) {
235 return isothermalMu;
236 }
237
238 // compute the viscosity deviation due to temperature
239 const auto& muOilvisct = oilvisctCurves_[regionIdx].eval(temperature, /*extrapolate=*/true);
240 return muOilvisct / viscRef_[regionIdx] * isothermalMu;
241 }
242
246 template <class Evaluation>
247 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
248 const Evaluation& temperature,
249 const Evaluation& pressure,
250 const Evaluation& Rs) const
251 {
252 const auto& b =
253 isothermalPvt_->inverseFormationVolumeFactor(regionIdx, temperature, pressure, Rs);
254
255 if (!enableThermalDensity()) {
256 return b;
257 }
258
259 // we use the same approach as for the for water here, but with the OPM-specific
260 // OILDENT keyword.
261 Scalar TRef = oildentRefTemp_[regionIdx];
262 Scalar cT1 = oildentCT1_[regionIdx];
263 Scalar cT2 = oildentCT2_[regionIdx];
264 const Evaluation& Y = temperature - TRef;
265
266 return b / (1 + (cT1 + cT2 * Y) * Y);
267 }
268
272 template <class FluidState, class LhsEval = typename FluidState::Scalar>
273 std::pair<LhsEval, LhsEval>
274 inverseFormationVolumeFactorAndViscosity(const FluidState& fluidState, unsigned regionIdx)
275 {
276 auto [b, mu] = isothermalPvt_->inverseFormationVolumeFactorAndViscosity(fluidState, regionIdx);
277 const LhsEval& temperature = decay<LhsEval>(fluidState.temperature(FluidState::oilPhaseIdx));
278 if (enableThermalDensity()) {
279 // we use the same approach as for the for water here, but with the OPM-specific
280 // OILDENT keyword.
281 Scalar TRef = oildentRefTemp_[regionIdx];
282 Scalar cT1 = oildentCT1_[regionIdx];
283 Scalar cT2 = oildentCT2_[regionIdx];
284 const LhsEval Y = temperature - TRef;
285 b /= (1.0 + (cT1 + cT2 * Y) * Y);
286 }
288 // compute the viscosity deviation due to temperature
289 const auto muOilvisct = oilvisctCurves_[regionIdx].eval(temperature, /*extrapolate=*/true);
290 mu *= (muOilvisct / viscRef_[regionIdx]);
291 }
292 return { b, mu };
293 }
294
298 template <class Evaluation>
299 Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
300 const Evaluation& temperature,
301 const Evaluation& pressure) const
302 {
303 const auto& b =
304 isothermalPvt_->saturatedInverseFormationVolumeFactor(regionIdx, temperature, pressure);
305
306 if (!enableThermalDensity()) {
307 return b;
308 }
309
310 // we use the same approach as for the for water here, but with the OPM-specific
311 // OILDENT keyword.
312 Scalar TRef = oildentRefTemp_[regionIdx];
313 Scalar cT1 = oildentCT1_[regionIdx];
314 Scalar cT2 = oildentCT2_[regionIdx];
315 const Evaluation& Y = temperature - TRef;
316
317 return b / (1 + (cT1 + cT2 * Y) * Y);
318 }
319
327 template <class Evaluation>
328 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
329 const Evaluation& temperature,
330 const Evaluation& pressure) const
331 { return isothermalPvt_->saturatedGasDissolutionFactor(regionIdx, temperature, pressure); }
332
340 template <class Evaluation>
341 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
342 const Evaluation& temperature,
343 const Evaluation& pressure,
344 const Evaluation& oilSaturation,
345 const Evaluation& maxOilSaturation) const
346 { return isothermalPvt_->saturatedGasDissolutionFactor(regionIdx, temperature, pressure, oilSaturation, maxOilSaturation); }
347
355 template <class Evaluation>
356 Evaluation saturationPressure(unsigned regionIdx,
357 const Evaluation& temperature,
358 const Evaluation& pressure) const
359 { return isothermalPvt_->saturationPressure(regionIdx, temperature, pressure); }
360
361 template <class Evaluation>
362 Evaluation diffusionCoefficient(const Evaluation& temperature,
363 const Evaluation& pressure,
364 unsigned compIdx) const
365 {
366 return isothermalPvt_->diffusionCoefficient(temperature, pressure, compIdx);
367 }
368
369 const IsothermalPvt* isoThermalPvt() const
370 { return isothermalPvt_; }
371
372 Scalar oilReferenceDensity(unsigned regionIdx) const
373 { return isothermalPvt_->oilReferenceDensity(regionIdx); }
374
375 Scalar hVap(unsigned regionIdx) const
376 { return this->hVap_[regionIdx]; }
377
378 const std::vector<TabulatedOneDFunction>& oilvisctCurves() const
379 { return oilvisctCurves_; }
380
381 const std::vector<Scalar>& viscrefPress() const
382 { return viscrefPress_; }
383
384 const std::vector<Scalar>& viscrefRs() const
385 { return viscrefRs_; }
386
387 const std::vector<Scalar>& viscRef() const
388 { return viscRef_; }
389
390 const std::vector<Scalar>& oildentRefTemp() const
391 { return oildentRefTemp_; }
392
393 const std::vector<Scalar>& oildentCT1() const
394 { return oildentCT1_; }
395
396 const std::vector<Scalar>& oildentCT2() const
397 { return oildentCT2_; }
398
399 const std::vector<TabulatedOneDFunction>& internalEnergyCurves() const
400 { return internalEnergyCurves_; }
401
402 bool enableInternalEnergy() const
403 { return enableInternalEnergy_; }
404
405 const std::vector<Scalar>& oilJTRefPres() const
406 { return oilJTRefPres_; }
407
408 const std::vector<Scalar>& oilJTC() const
409 { return oilJTC_; }
410
411 bool operator==(const OilPvtThermal<Scalar>& data) const;
412
413 OilPvtThermal<Scalar>& operator=(const OilPvtThermal<Scalar>& data);
414
415private:
416 IsothermalPvt* isothermalPvt_{nullptr};
417
418 // The PVT properties needed for temperature dependence of the viscosity. We need
419 // to store one value per PVT region.
420 std::vector<TabulatedOneDFunction> oilvisctCurves_{};
421 std::vector<Scalar> viscrefPress_{};
422 std::vector<Scalar> viscrefRs_{};
423 std::vector<Scalar> viscRef_{};
424
425 // The PVT properties needed for temperature dependence of the density.
426 std::vector<Scalar> oildentRefTemp_{};
427 std::vector<Scalar> oildentCT1_{};
428 std::vector<Scalar> oildentCT2_{};
429
430 std::vector<Scalar> oilJTRefPres_{};
431 std::vector<Scalar> oilJTC_{};
432
433 std::vector<Scalar> rhoRefG_{};
434 std::vector<Scalar> hVap_{};
435
436 // piecewise linear curve representing the internal energy of oil
437 std::vector<TabulatedOneDFunction> internalEnergyCurves_{};
438
439 bool enableThermalDensity_{false};
440 bool enableJouleThomson_{false};
441 bool enableThermalViscosity_{false};
442 bool enableInternalEnergy_{false};
443};
444
445} // namespace Opm
446
447#endif
Implements a linearly interpolated scalar function that depends on one variable.
Definition EclipseState.hpp:62
This class represents the Pressure-Volume-Temperature relations of the oil phase in the black-oil mod...
Definition OilPvtMultiplexer.hpp:105
Scalar oilReferenceDensity(unsigned regionIdx) const
Return the reference density which are considered by this PVT-object.
Definition OilPvtMultiplexer.cpp:116
Evaluation diffusionCoefficient(const Evaluation &temperature, const Evaluation &pressure, unsigned compIdx) const
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition OilPvtMultiplexer.hpp:253
std::pair< LhsEval, LhsEval > inverseFormationVolumeFactorAndViscosity(const FluidState &fluidState, unsigned regionIdx)
Returns the formation volume factor [-] and viscosity [Pa s] of the fluid phase.
Definition OilPvtThermal.hpp:274
void initEnd()
Finish initializing the thermal part of the oil phase PVT properties.
Definition OilPvtThermal.hpp:117
bool enableThermalDensity() const
Returns true iff the density of the oil phase is temperature dependent.
Definition OilPvtThermal.hpp:123
Evaluation saturatedGasDissolutionFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns the gas dissolution factor [m^3/m^3] of the oil phase.
Definition OilPvtThermal.hpp:328
Evaluation inverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &Rs) const
Returns the formation volume factor [-] of the fluid phase.
Definition OilPvtThermal.hpp:247
Evaluation saturatedViscosity(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns the dynamic viscosity [Pa s] of the fluid phase given a set of parameters.
Definition OilPvtThermal.hpp:229
Evaluation internalEnergy(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &Rs) const
Returns the specific internal energy [J/kg] of oil given a set of parameters.
Definition OilPvtThermal.hpp:145
Evaluation viscosity(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &Rs) const
Returns the dynamic viscosity [Pa s] of the fluid phase given a set of parameters.
Definition OilPvtThermal.hpp:210
void setNumRegions(std::size_t numRegions)
Set the number of PVT-regions considered by this object.
Definition OilPvtThermal.cpp:183
bool enableThermalViscosity() const
Returns true iff the viscosity of the oil phase is temperature dependent.
Definition OilPvtThermal.hpp:135
Evaluation saturationPressure(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns the saturation pressure of the oil phase [Pa].
Definition OilPvtThermal.hpp:356
Evaluation saturatedGasDissolutionFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &oilSaturation, const Evaluation &maxOilSaturation) const
Returns the gas dissolution factor [m^3/m^3] of the oil phase.
Definition OilPvtThermal.hpp:341
Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns the formation volume factor [-] of gas-saturated oil phase.
Definition OilPvtThermal.hpp:299
bool enableJouleThomson() const
Returns true iff Joule-Thomson effect for the oil phase is active.
Definition OilPvtThermal.hpp:129
Definition Schedule.hpp:101
Implements a linearly interpolated scalar function that depends on one variable.
Definition Tabulated1DFunction.hpp:51
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30