3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
heavyoil.hh
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program 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 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_HEAVYOIL_HH
25#define DUMUX_HEAVYOIL_HH
26
27#include <algorithm>
28
29#include <dune/common/math.hh>
30
36
37namespace Dumux::Components {
38
45template <class Scalar>
47: public Components::Base<Scalar, HeavyOil<Scalar> >
48, public Components::Liquid<Scalar, HeavyOil<Scalar> >
49, public Components::Gas<Scalar, HeavyOil<Scalar> >
50{
53public:
57 static std::string name()
58 { return "heavyoil"; }
59
63 constexpr static Scalar molarMass()
64 { return .350; }
65
70 { return .400; }
71
76 constexpr static Scalar molecularWeight()
77 { return .350; }
78
82 constexpr static Scalar specificGravity()
83 { return 0.91; }
84
89 {
90 DUNE_THROW(Dune::NotImplemented, "tripleTemperature for heavyoil");
91 }
92
97 {
98 DUNE_THROW(Dune::NotImplemented, "triplePressure for heavyoil");
99 }
100
102 {
103 constexpr Scalar A = 0.83;
104 constexpr Scalar B = 89.9513;
105 constexpr Scalar C = 139.6612;
106 constexpr Scalar D = 3.2033;
107 constexpr Scalar E = 1.0564;
108
109 const Scalar mW = refComponentMolecularWeight() *1000. ; // in [g/mol];
110
111 using std::pow;
112 return A+(B/mW)-(C/pow((mW+D),E));
113 }
114
116 {
117 constexpr Scalar A = -7.4120e-2; //All factors for 1 atm / 101325 pascals [760 mmHg]
118 constexpr Scalar B = -7.5041e-3;
119 constexpr Scalar C = -2.6031;
120 constexpr Scalar D = 9.0180e-2;
121 constexpr Scalar E = -1.0482;
122
123 using std::log;
124 const Scalar deltaSpecificGravity = log(refComponentSpecificGravity()/specificGravity());
125 const Scalar deltaMolecularWeight = log(refComponentMolecularWeight()/molecularWeight());
126
127 using Dune::power;
128 return A*power(deltaSpecificGravity,2) + B*deltaSpecificGravity + C*power(deltaMolecularWeight,2) + D*deltaMolecularWeight
129 + E*deltaSpecificGravity*deltaMolecularWeight;
130 }
131
133 {
134 constexpr Scalar A = -6.1294e-2;
135 constexpr Scalar B = -7.0862e-2;
136 constexpr Scalar C = 6.1976e-1;
137 constexpr Scalar D = -5.7090e-2;
138 constexpr Scalar E = -8.4583e-2;
139
140 using std::log;
141 const Scalar deltaSpecificGravity = log(refComponentSpecificGravity()/specificGravity());
142 const Scalar deltaMolecularWeight = log(refComponentMolecularWeight()/molecularWeight());
143
144 using Dune::power;
145 return A*power(deltaSpecificGravity,2) + B*deltaSpecificGravity + C*power(deltaMolecularWeight,2) + D*deltaMolecularWeight
146 + E*deltaSpecificGravity*deltaMolecularWeight;
147 }
148
150 {
151 constexpr Scalar A = 1.8270e-1;
152 constexpr Scalar B = -2.4864e-1;
153 constexpr Scalar C = 8.3611;
154 constexpr Scalar D = -2.2389e-1;
155 constexpr Scalar E = 2.6984;
156
157 using std::log;
158 const Scalar deltaSpecificGravity = log(refComponentSpecificGravity()/specificGravity());
159 const Scalar deltaMolecularWeight = log(refComponentMolecularWeight()/molecularWeight());
160
161 using Dune::power;
162 return A*power(deltaSpecificGravity,2) + B*deltaSpecificGravity + C*power(deltaMolecularWeight,2) + D*deltaMolecularWeight
163 + E*deltaSpecificGravity*deltaMolecularWeight;
164 }
165
167 {
168 constexpr Scalar A = 477.63; //All factors for 1 atm / 101325 pascals [760 mmHg]
169 constexpr Scalar B = 88.51;
170 constexpr Scalar C = 1007;
171 constexpr Scalar D = 1214.40;
172
173 using std::log;
174 return A*log((1000.*refComponentMolecularWeight() + B)/(1000.*refComponentMolecularWeight()+C)) + D;
175 }
176
178 {
179 constexpr Scalar A = 226.50;
180 constexpr Scalar B = 6.78;
181 constexpr Scalar C = 1.282e6;
182 constexpr Scalar D = 2668;
183
184 using std::log;
185 return A*log((1000.*refComponentMolecularWeight() + B)/(1000.*refComponentMolecularWeight()+C)) + D ;
186 }
187
189 {
190 constexpr Scalar A = 141.20;
191 constexpr Scalar B = 45.66e-2;
192 constexpr Scalar C = 16.59e-3;
193 constexpr Scalar D = 2.19;
194
195 using std::pow;
196 return (A*1000.*molecularWeight())/(pow(B + (C*1000.*molecularWeight()),D)) ;
197 }
198
203 {
204 using Dune::power;
206 }
207
212 {
213 using Dune::power;
215 }
216
221 {
222 using Dune::power;
224 }
225
232 {
233 constexpr Scalar A = 8.25990;
234 constexpr Scalar B = 2830.065;
235 constexpr Scalar C = 42.95101;
236
237 const Scalar T = temperature - 273.15;
238
239 using std::pow;
240 return 100*1.334*pow(10.0, (A - (B/(T + C)))); // in [Pa]
241 }
242
244 {
245 constexpr Scalar A = 8.25990;
246 constexpr Scalar B = 2830.065;
247 constexpr Scalar C = 42.95101;
248
249 const Scalar P = pressure;
250
251 using std::log10;
252 return Scalar ((B/(A-log10(P/100*1.334)))-C);
253 }
254
262 const Scalar pressure)
263 {
264 // Gauss quadrature rule:
265 // Interval: [0K; temperature (K)]
266 // Gauss-Legendre-Integration with variable transformation:
267 // \int_a^b f(T) dT \approx (b-a)/2 \sum_i=1^n \alpha_i f( (b-a)/2 x_i + (a+b)/2 )
268 // with: n=2, legendre -> x_i = +/- \sqrt(1/3), \apha_i=1
269 // here: a=273.15K, b=actual temperature in Kelvin
270 // \leadsto h(T) = \int_273.15^T c_p(T) dT
271 // \approx 0.5 (T-273.15) * (cp( 0.5(temperature-273.15)sqrt(1/3) ) + cp(0.5(temperature-273.15)(-1)sqrt(1/3))
272
273 // Enthalpy may have arbitrary reference state, but the empirical/fitted heatCapacity function needs Kelvin as input and is
274 // fit over a certain temperature range. This suggests choosing an interval of integration being in the actual fit range.
275 // I.e. choosing T=273.15K as reference point for liquid enthalpy.
276 using std::sqrt;
277 const Scalar sqrt1over3 = sqrt(1./3.);
278 const Scalar TEval1 = 0.5*(temperature-273.15)* sqrt1over3 + 0.5*(273.15+temperature) ; // evaluation points according to Gauss-Legendre integration
279 const Scalar TEval2 = 0.5*(temperature-273.15)* (-1)* sqrt1over3 + 0.5*(273.15+temperature) ; // evaluation points according to Gauss-Legendre integration
280
281 const Scalar h_n = 0.5 * (temperature-273.15) * ( liquidHeatCapacity(TEval1, pressure) + liquidHeatCapacity(TEval2, pressure) ) ;
282
283 return h_n;
284 }
285
295 const Scalar pressure)
296 {
297 using std::clamp;
298 temperature = clamp(temperature, 0.0, criticalTemperature()); // regularization
299
300 const Scalar T_crit = criticalTemperature();
302 const Scalar p_crit = criticalPressure();
303
304 // Chen method, eq. 7-11.4 (at boiling)
305 using std::log;
306 const Scalar DH_v_boil = Consts::R * T_crit * Tr1
307 * (3.978 * Tr1 - 3.958 + 1.555*log(p_crit * 1e-5 /*Pa->bar*/ ) )
308 / (1.07 - Tr1); /* [J/mol] */
309
310 /* Variation with temp according to Watson relation eq 7-12.1*/
311 using std::pow;
313 const Scalar n = 0.375;
314 const Scalar DH_vap = DH_v_boil * pow(((1.0 - Tr2)/(1.0 - Tr1)), n);
315
316 return (DH_vap/molarMass()); // we need [J/kg]
317 }
318
319
330 {
332 }
333
341 {
343 }
344
353
361 {
362
363 /* according to Lashanizadegan et al (2008) in Chemical Engineering Communications: */
364 /* Simultaneous Heat and Fluid Flow in Porous Media: Case Study: Steam Injection for Tertiary Oil Recovery */
365 constexpr Scalar rhoReference = 906.; // [kg/m^3] at reference pressure and temperature
366 constexpr Scalar compressCoeff = 1.e-8; // just a value without justification
367 constexpr Scalar expansCoeff = 1.e-7; // also just a value
368 return rhoReference * (1. + (pressure - 1.e5)*compressCoeff) * (1. - (temperature - 293.)*expansCoeff);
369 }
370
380
384 static constexpr bool gasIsCompressible()
385 { return true; }
386
390 static constexpr bool gasIsIdeal()
391 { return true; }
392
396 static constexpr bool liquidIsCompressible()
397 { return true; }
398
406 static Scalar gasViscosity(Scalar temperature, Scalar pressure, bool regularize=true)
407 {
408 using std::clamp;
409 temperature = clamp(temperature, 250.0, 500.0); // regularization
410
411 // reduced temperature
413
414 using std::pow;
415 using std::exp;
416 constexpr Scalar Fp0 = 1.0;
417 constexpr Scalar xi = 0.00474;
418 const Scalar eta_xi =
419 Fp0*(0.807*pow(Tr,0.618)
420 - 0.357*exp(-0.449*Tr)
421 + 0.34*exp(-4.058*Tr)
422 + 0.018);
423
424 return eta_xi/xi/1e7; // [Pa s]
425 }
426
436 {
437
438 /* according to Lashanizadegan et al (2008) in Chemical Engineering Communications: */
439 /* Simultaneous Heat and Fluid Flow in Porous Media: Case Study: Steam Injection for Tertiary Oil Recovery */
440
441 //using std::exp;
442 //return 1027919.422*exp(-0.04862*temperature); // [Pa s]
443
444 //according to http://www.ecltechnology.com/subsur/reports/pvt_tgb.pdf [Page 10]
445 const Scalar temperatureFahrenheit = (9/5)*(temperature-273.15)+32;
446 constexpr Scalar API = 9;
447
448 using std::pow;
449 using Dune::power;
450 return ((pow(10,0.10231*power(API,2)-3.9464*API+46.5037))*(pow(temperatureFahrenheit,-0.04542*power(API,2)+1.70405*API-19.18)))*0.001;
451
452 }
463 const Scalar pressure)
464 {
465 return 618.; // J/(kg K)
466 }
467
477 {
478 return 0.127;
479 }
480};
481
482} // end namespace Dumux::Components
483
484#endif
Interface for components that have a gas state.
Interface for components that have a liquid state.
A central place for various physical constants occuring in some equations.
Relations valid for an ideal gas.
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
Definition: air.hh:35
Base class for all components Components provide the thermodynamic relations for the liquid,...
Definition: components/base.hh:59
Scalar Scalar
export the scalar type used by the component
Definition: components/base.hh:63
Interface for components that have a gas state.
Definition: gas.hh:41
Properties of the component heavyoil.
Definition: heavyoil.hh:50
static Scalar liquidViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of pure heavyoil.
Definition: heavyoil.hh:435
static constexpr Scalar molarMass()
The molar mass in of heavyoil.
Definition: heavyoil.hh:63
static Scalar vaporPressure(Scalar temperature)
The saturation vapor pressure in of.
Definition: heavyoil.hh:231
static constexpr Scalar refComponentMolecularWeight()
The MolecularWeight in of refComponent.
Definition: heavyoil.hh:69
static Scalar perbutationFactorCriticalTemperature()
Definition: heavyoil.hh:132
static Scalar liquidHeatCapacity(const Scalar temperature, const Scalar pressure)
Specific heat cap of liquid heavyoil .
Definition: heavyoil.hh:462
static constexpr Scalar molecularWeight()
The molar mass in of heavyoil.
Definition: heavyoil.hh:76
static Scalar refComponentCriticalPressure()
Definition: heavyoil.hh:188
static constexpr bool gasIsIdeal()
Returns true if the gas phase is assumed to be ideal.
Definition: heavyoil.hh:390
static Scalar liquidMolarDensity(Scalar temperature, Scalar pressure)
The molar density of pure heavyoil in at a given pressure and temperature.
Definition: heavyoil.hh:378
static Scalar perbutationFactorBoilingTemperature()
Definition: heavyoil.hh:115
static Scalar tripleTemperature()
Returns the temperature at heavyoil's triple point.
Definition: heavyoil.hh:88
static Scalar criticalPressure()
Returns the critical pressure of heavyoil.
Definition: heavyoil.hh:220
static Scalar gasMolarDensity(Scalar temperature, Scalar pressure)
The molar density of pure heavyoil in , depending on pressure and temperature.
Definition: heavyoil.hh:351
static Scalar vaporTemperature(Scalar pressure)
Definition: heavyoil.hh:243
static constexpr bool gasIsCompressible()
Returns true if the gas phase is assumed to be compressible.
Definition: heavyoil.hh:384
static Scalar liquidThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of heavy oil.
Definition: heavyoil.hh:476
static constexpr Scalar specificGravity()
The Specific Gravity of heavyoil.
Definition: heavyoil.hh:82
static Scalar gasViscosity(Scalar temperature, Scalar pressure, bool regularize=true)
The dynamic viscosity of heavyoil vapor.
Definition: heavyoil.hh:406
static Scalar criticalTemperature()
Returns the critical temperature of heavyoil.
Definition: heavyoil.hh:211
static Scalar gasDensity(Scalar temperature, Scalar pressure)
The (ideal) gas density of heavyoil vapor at a given temperature and pressure .
Definition: heavyoil.hh:340
static Scalar triplePressure()
Returns the pressure at heavyoil's triple point.
Definition: heavyoil.hh:96
static Scalar refComponentCriticalTemperature()
Definition: heavyoil.hh:177
static constexpr bool liquidIsCompressible()
Returns true if the liquid phase is assumed to be compressible.
Definition: heavyoil.hh:396
static Scalar liquidEnthalpy(const Scalar temperature, const Scalar pressure)
Specific enthalpy of liquid heavyoil .
Definition: heavyoil.hh:261
static Scalar refComponentBoilingTemperature()
Definition: heavyoil.hh:166
static std::string name()
A human readable name for heavyoil.
Definition: heavyoil.hh:57
static Scalar perbutationFactorCriticalPressure()
Definition: heavyoil.hh:149
static Scalar refComponentSpecificGravity()
Definition: heavyoil.hh:101
static Scalar heatVap(Scalar temperature, const Scalar pressure)
Latent heat of vaporization for heavyoil .
Definition: heavyoil.hh:294
static Scalar boilingTemperature()
Returns the temperature at heavyoil's boiling point (1 atm)
Definition: heavyoil.hh:202
static Scalar gasEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of heavyoil vapor .
Definition: heavyoil.hh:329
static Scalar liquidDensity(Scalar temperature, Scalar pressure)
The density of pure heavyoil at a given pressure and temperature .
Definition: heavyoil.hh:360
Interface for components that have a liquid state.
Definition: liquid.hh:41
A central place for various physical constants occuring in some equations.
Definition: constants.hh:39
static constexpr Scalar R
The ideal gas constant .
Definition: constants.hh:44
Relations valid for an ideal gas.
Definition: idealgas.hh:37
static constexpr Scalar density(Scalar avgMolarMass, Scalar temperature, Scalar pressure)
The density of the gas in , depending on pressure, temperature and average molar mass of the gas.
Definition: idealgas.hh:49
static constexpr Scalar molarDensity(Scalar temperature, Scalar pressure)
The molar density of the gas , depending on pressure and temperature.
Definition: idealgas.hh:70
Base class for all components Components provide the thermodynamic relations for the liquid,...