27#ifndef OPM_OIL_PVT_THERMAL_HPP
28#define OPM_OIL_PVT_THERMAL_HPP
39template <
class Scalar,
bool enableThermal>
48template <
class Scalar>
55 OilPvtThermal() =
default;
57 OilPvtThermal(IsothermalPvt* isothermalPvt,
58 const std::vector<TabulatedOneDFunction>& oilvisctCurves,
59 const std::vector<Scalar>& viscrefPress,
60 const std::vector<Scalar>& viscrefRs,
61 const std::vector<Scalar>& viscRef,
62 const std::vector<Scalar>& oildentRefTemp,
63 const std::vector<Scalar>& oildentCT1,
64 const std::vector<Scalar>& oildentCT2,
65 const std::vector<Scalar>& oilJTRefPres,
66 const std::vector<Scalar>& oilJTC,
67 const std::vector<TabulatedOneDFunction>& internalEnergyCurves,
71 bool enableInternalEnergy)
72 : isothermalPvt_(isothermalPvt)
73 , oilvisctCurves_(oilvisctCurves)
74 , viscrefPress_(viscrefPress)
75 , viscrefRs_(viscrefRs)
77 , oildentRefTemp_(oildentRefTemp)
78 , oildentCT1_(oildentCT1)
79 , oildentCT2_(oildentCT2)
80 , oilJTRefPres_(oilJTRefPres)
82 , internalEnergyCurves_(internalEnergyCurves)
86 , enableInternalEnergy_(enableInternalEnergy)
89 OilPvtThermal(
const OilPvtThermal& data)
93 {
delete isothermalPvt_; }
105 void setVapPars(
const Scalar par1,
const Scalar par2)
107 isothermalPvt_->setVapPars(par1, par2);
120 {
return enableThermalDensity_; }
126 {
return enableJouleThomson_; }
132 {
return enableThermalViscosity_; }
134 std::size_t numRegions()
const
135 {
return viscrefRs_.size(); }
140 template <
class Evaluation>
142 const Evaluation& temperature,
143 const Evaluation& pressure,
144 const Evaluation& Rs)
const
146 if (!enableInternalEnergy_) {
147 throw std::runtime_error(
"Requested the internal energy of oil but it is disabled");
150 if (!enableJouleThomson_) {
154 return internalEnergyCurves_[regionIdx].eval(temperature,
true);
157 OpmLog::warning(
"Experimental code for jouleThomson: simulation will be slower");
158 Evaluation Tref = oildentRefTemp_[regionIdx];
159 Evaluation Pref = oilJTRefPres_[regionIdx];
160 Scalar JTC = oilJTC_[regionIdx];
163 Evaluation Cp = internalEnergyCurves_[regionIdx].eval(temperature,
true)/temperature;
164 Evaluation density = invB * (oilReferenceDensity(regionIdx) + Rs * rhoRefG_[regionIdx]);
166 Evaluation enthalpyPres;
168 enthalpyPres = -Cp * JTC * (pressure - Pref);
170 else if (enableThermalDensity_) {
171 Scalar c1T = oildentCT1_[regionIdx];
172 Scalar c2T = oildentCT2_[regionIdx];
174 Evaluation alpha = (c1T + 2 * c2T * (temperature - Tref)) /
175 (1 + c1T *(temperature - Tref) + c2T * (temperature - Tref) * (temperature - Tref));
178 Evaluation deltaP = (pressure - Pref) / N;
179 Evaluation enthalpyPresPrev = 0;
180 for (std::size_t i = 0; i < N; ++i) {
181 Evaluation Pnew = Pref + i * deltaP;
183 (oilReferenceDensity(regionIdx) + Rs * rhoRefG_[regionIdx]) ;
185 Evaluation jouleThomsonCoefficient = -(1.0/Cp) * (1.0 - alpha * temperature)/rho;
186 Evaluation deltaEnthalpyPres = -Cp * jouleThomsonCoefficient * deltaP;
187 enthalpyPres = enthalpyPresPrev + deltaEnthalpyPres;
188 enthalpyPresPrev = enthalpyPres;
192 throw std::runtime_error(
"Requested Joule-thomson calculation but thermal oil"
193 "density (OILDENT) is not provided");
196 Evaluation enthalpy = Cp * (temperature - Tref) + enthalpyPres;
198 return enthalpy - pressure/density;
205 template <
class Evaluation>
207 const Evaluation& temperature,
208 const Evaluation& pressure,
209 const Evaluation& Rs)
const
211 const auto& isothermalMu = isothermalPvt_->viscosity(regionIdx, temperature, pressure, Rs);
217 const auto& muOilvisct = oilvisctCurves_[regionIdx].eval(temperature,
true);
218 return muOilvisct / viscRef_[regionIdx] * isothermalMu;
224 template <
class Evaluation>
226 const Evaluation& temperature,
227 const Evaluation& pressure)
const
229 const auto& isothermalMu = isothermalPvt_->saturatedViscosity(regionIdx, temperature, pressure);
235 const auto& muOilvisct = oilvisctCurves_[regionIdx].eval(temperature,
true);
236 return muOilvisct / viscRef_[regionIdx] * isothermalMu;
242 template <
class Evaluation>
244 const Evaluation& temperature,
245 const Evaluation& pressure,
246 const Evaluation& Rs)
const
249 isothermalPvt_->inverseFormationVolumeFactor(regionIdx, temperature, pressure, Rs);
257 Scalar TRef = oildentRefTemp_[regionIdx];
258 Scalar cT1 = oildentCT1_[regionIdx];
259 Scalar cT2 = oildentCT2_[regionIdx];
260 const Evaluation& Y = temperature - TRef;
262 return b / (1 + (cT1 + cT2 * Y) * Y);
268 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::ValueType>
269 std::pair<LhsEval, LhsEval>
272 auto [b, mu] = isothermalPvt_->inverseFormationVolumeFactorAndViscosity(fluidState, regionIdx);
273 const LhsEval& temperature = decay<LhsEval>(fluidState.temperature(FluidState::oilPhaseIdx));
277 Scalar TRef = oildentRefTemp_[regionIdx];
278 Scalar cT1 = oildentCT1_[regionIdx];
279 Scalar cT2 = oildentCT2_[regionIdx];
280 const LhsEval Y = temperature - TRef;
281 b /= (1.0 + (cT1 + cT2 * Y) * Y);
285 const auto muOilvisct = oilvisctCurves_[regionIdx].eval(temperature,
true);
286 mu *= (muOilvisct / viscRef_[regionIdx]);
294 template <
class Evaluation>
296 const Evaluation& temperature,
297 const Evaluation& pressure)
const
300 isothermalPvt_->saturatedInverseFormationVolumeFactor(regionIdx, temperature, pressure);
308 Scalar TRef = oildentRefTemp_[regionIdx];
309 Scalar cT1 = oildentCT1_[regionIdx];
310 Scalar cT2 = oildentCT2_[regionIdx];
311 const Evaluation& Y = temperature - TRef;
313 return b / (1 + (cT1 + cT2 * Y) * Y);
323 template <
class Evaluation>
325 const Evaluation& temperature,
326 const Evaluation& pressure)
const
327 {
return isothermalPvt_->saturatedGasDissolutionFactor(regionIdx, temperature, pressure); }
336 template <
class Evaluation>
338 const Evaluation& temperature,
339 const Evaluation& pressure,
340 const Evaluation& oilSaturation,
341 const Evaluation& maxOilSaturation)
const
342 {
return isothermalPvt_->saturatedGasDissolutionFactor(regionIdx, temperature, pressure, oilSaturation, maxOilSaturation); }
351 template <
class Evaluation>
353 const Evaluation& temperature,
354 const Evaluation& pressure)
const
355 {
return isothermalPvt_->saturationPressure(regionIdx, temperature, pressure); }
357 template <
class Evaluation>
358 Evaluation diffusionCoefficient(
const Evaluation& temperature,
359 const Evaluation& pressure,
360 unsigned compIdx)
const
365 const IsothermalPvt* isoThermalPvt()
const
366 {
return isothermalPvt_; }
368 Scalar oilReferenceDensity(
unsigned regionIdx)
const
371 Scalar hVap(
unsigned regionIdx)
const
372 {
return this->hVap_[regionIdx]; }
374 const std::vector<TabulatedOneDFunction>& oilvisctCurves()
const
375 {
return oilvisctCurves_; }
377 const std::vector<Scalar>& viscrefPress()
const
378 {
return viscrefPress_; }
380 const std::vector<Scalar>& viscrefRs()
const
381 {
return viscrefRs_; }
383 const std::vector<Scalar>& viscRef()
const
386 const std::vector<Scalar>& oildentRefTemp()
const
387 {
return oildentRefTemp_; }
389 const std::vector<Scalar>& oildentCT1()
const
390 {
return oildentCT1_; }
392 const std::vector<Scalar>& oildentCT2()
const
393 {
return oildentCT2_; }
395 const std::vector<TabulatedOneDFunction>& internalEnergyCurves()
const
396 {
return internalEnergyCurves_; }
398 bool enableInternalEnergy()
const
399 {
return enableInternalEnergy_; }
401 const std::vector<Scalar>& oilJTRefPres()
const
402 {
return oilJTRefPres_; }
404 const std::vector<Scalar>& oilJTC()
const
407 bool operator==(
const OilPvtThermal<Scalar>& data)
const;
409 OilPvtThermal<Scalar>& operator=(
const OilPvtThermal<Scalar>& data);
412 IsothermalPvt* isothermalPvt_{
nullptr};
416 std::vector<TabulatedOneDFunction> oilvisctCurves_{};
417 std::vector<Scalar> viscrefPress_{};
418 std::vector<Scalar> viscrefRs_{};
419 std::vector<Scalar> viscRef_{};
422 std::vector<Scalar> oildentRefTemp_{};
423 std::vector<Scalar> oildentCT1_{};
424 std::vector<Scalar> oildentCT2_{};
426 std::vector<Scalar> oilJTRefPres_{};
427 std::vector<Scalar> oilJTC_{};
429 std::vector<Scalar> rhoRefG_{};
430 std::vector<Scalar> hVap_{};
433 std::vector<TabulatedOneDFunction> internalEnergyCurves_{};
435 bool enableThermalDensity_{
false};
436 bool enableJouleThomson_{
false};
437 bool enableThermalViscosity_{
false};
438 bool enableInternalEnergy_{
false};
Implements a linearly interpolated scalar function that depends on one variable.
Definition EclipseState.hpp:66
This class represents the Pressure-Volume-Temperature relations of the oil phase in the black-oil mod...
Definition OilPvtMultiplexer.hpp:110
Scalar oilReferenceDensity(unsigned regionIdx) const
Return the reference density which are considered by this PVT-object.
Definition OilPvtMultiplexer.cpp:130
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:256
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:270
void initEnd()
Finish initializing the thermal part of the oil phase PVT properties.
Definition OilPvtThermal.hpp:113
bool enableThermalDensity() const
Returns true iff the density of the oil phase is temperature dependent.
Definition OilPvtThermal.hpp:119
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:324
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:243
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:225
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:141
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:206
void setNumRegions(std::size_t numRegions)
Set the number of PVT-regions considered by this object.
Definition OilPvtThermal.cpp:181
bool enableThermalViscosity() const
Returns true iff the viscosity of the oil phase is temperature dependent.
Definition OilPvtThermal.hpp:131
void initFromState(const EclipseState &eclState, const Schedule &schedule)
Implement the temperature part of the oil PVT properties.
Definition OilPvtThermal.cpp:41
Evaluation saturationPressure(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns the saturation pressure of the oil phase [Pa].
Definition OilPvtThermal.hpp:352
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:337
Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns the formation volume factor [-] of gas-saturated oil phase.
Definition OilPvtThermal.hpp:295
bool enableJouleThomson() const
Returns true iff Joule-Thomson effect for the oil phase is active.
Definition OilPvtThermal.hpp:125
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