3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
xylene.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_XYLENE_HH
25#define DUMUX_XYLENE_HH
26
27#include <cmath>
30
34
35namespace Dumux {
36namespace Components {
37
44template <class Scalar>
45class Xylene
46: public Components::Base<Scalar, Xylene<Scalar> >
47, public Components::Liquid<Scalar, Xylene<Scalar> >
48, public Components::Gas<Scalar, Xylene<Scalar> >
49{
52
53public:
57 static std::string name()
58 { return "xylene"; }
59
63 constexpr static Scalar molarMass()
64 { return 0.106; }
65
69 constexpr static Scalar criticalTemperature()
70 { return 617.1; }
71
75 constexpr static Scalar criticalPressure()
76 { return 35.4e5; }
77
81 constexpr static Scalar boilingTemperature()
82 { return 412.3; }
83
88 {
89 DUNE_THROW(Dune::NotImplemented, "tripleTemperature for xylene");
90 }
91
96 {
97 DUNE_THROW(Dune::NotImplemented, "triplePressure for xylene");
98 }
99
107 {
108 const Scalar A = 7.00909;
109 const Scalar B = 1462.266;
110 const Scalar C = 215.110;
111
112 Scalar T = temperature - 273.15;
113
114 using std::pow;
115 Scalar psat = 1.334*pow(10.0, (A - (B/(T + C)))); // in [mbar]
116 psat *= 100.0; // in [Pa] (0.001*1.E5)
117
118 return psat;
119 }
120
130 {
131 Scalar CH3,C6H5,H;
132 // after Reid et al. : Missenard group contrib. method (s. example 5-8) \cite reid1987 <BR>
133 // Xylene: C9H12 : 3* CH3 ; 1* C6H5 (phenyl-ring) ; -2* H (this was too much!)
134 // linear interpolation between table values [J/(mol K)]
135
136 if(temp < 298.0){ // take care: extrapolation for Temp<273
137 H = 13.4 + 1.2*(temp - 273.0)/25.0; // 13.4 + 1.2 = 14.6 = H(T=298K) i.e. interpolation of table values 273<T<298
138 CH3 = 40.0 + 1.6*(temp - 273.0)/25.0; // 40 + 1.6 = 41.6 = CH3(T=298K)
139 C6H5 = 113.0 + 4.2*(temp - 273.0)/25.0; // 113 + 4.2 = 117.2 = C6H5(T=298K)
140 }
141 else if(temp < 323.0){
142 H = 14.6 + 0.9*(temp - 298.0)/25.0; // i.e. interpolation of table values 298<T<323
143 CH3 = 41.6 + 1.9*(temp - 298.0)/25.0;
144 C6H5 = 117.2 + 6.2*(temp - 298.0)/25.0;
145 }
146 else if(temp < 348.0){
147 H = 15.5 + 1.2*(temp - 323.0)/25.0; // i.e. interpolation of table values 323<T<348
148 CH3 = 43.5 + 2.3*(temp - 323.0)/25.0;
149 C6H5 = 123.4 + 6.3*(temp - 323.0)/25.0;
150 }
151 else {
152 H = 16.7 + 2.1*(temp - 348.0)/25.0; // i.e. interpolation of table values 348<T<373
153 CH3 = 45.8 + 2.5*(temp - 348.0)/25.0; // take care: extrapolation for Temp>373
154 C6H5 = 129.7 + 6.3*(temp - 348.0)/25.0; // most likely leads to underestimation
155 }
156
157 return (C6H5 + 2*CH3 - H)/molarMass();// J/(mol K) -> J/(kg K)
158 }
159
160
168 const Scalar pressure)
169 {
170 // Gauss quadrature rule:
171 // Interval: [0K; temperature (K)]
172 // Gauss-Legendre-Integration with variable transformation:
173 // \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 )
174 // with: n=2, legendre -> x_i = +/- \sqrt(1/3), \apha_i=1
175 // here: a=273.15K, b=actual temperature in Kelvin
176 // \leadsto h(T) = \int_273.15^T c_p(T) dT
177 // \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))
178
179 // Enthalpy may have arbitrary reference state, but the empirical/fitted heatCapacity function needs Kelvin as input and is
180 // fit over a certain temperature range. This suggests choosing an interval of integration being in the actual fit range.
181 // I.e. choosing T=273.15K as reference point for liquid enthalpy.
182 using std::sqrt;
183 const Scalar sqrt1over3 = sqrt(1./3.);
184 // evaluation points according to Gauss-Legendre integration
185 const Scalar TEval1 = 0.5*(temperature-273.15)* sqrt1over3 + 0.5*(273.15+temperature);
186 // evaluation points according to Gauss-Legendre integration
187 const Scalar TEval2 = 0.5*(temperature-273.15)* (-1)* sqrt1over3 + 0.5*(273.15+temperature);
188
189 const Scalar h_n = 0.5 * (temperature-273.15) * ( liquidHeatCapacity(TEval1, pressure) + liquidHeatCapacity(TEval2, pressure) );
190
191 return h_n;
192 }
193
203 const Scalar pressure)
204 {
205 using std::min;
206 using std::max;
207 temperature = min(temperature, criticalTemperature()); // regularization
208 temperature = max(temperature, 0.0); // regularization
209
210 constexpr Scalar T_crit = criticalTemperature();
212 constexpr Scalar p_crit = criticalPressure();
213
214 // Chen method, eq. 7-11.4 (at boiling)
215 using std::log;
216 const Scalar DH_v_boil = Consts::R * T_crit * Tr1
217 * (3.978 * Tr1 - 3.958 + 1.555*log(p_crit * 1e-5 /*Pa->bar*/ ) )
218 / (1.07 - Tr1); /* [J/mol] */
219
220 /* Variation with temp according to Watson relation eq 7-12.1*/
221 using std::pow;
223 const Scalar n = 0.375;
224 const Scalar DH_vap = DH_v_boil * pow(((1.0 - Tr2)/(1.0 - Tr1)), n);
225
226 return (DH_vap/molarMass()); // we need [J/kg]
227 }
228
239 {
241 }
242
250 {
253 pressure);
254 }
255
264
275 {
276 // saturated molar volume according to Lide, CRC Handbook of
277 // Thermophysical and Thermochemical Data, CRC Press, 1994
278 // valid for 245 < Temp < 600
279 using std::min;
280 using std::max;
281 temp = min(temp, 500.0); // regularization
282 temp = max(temp, 250.0); // regularization
283
284 using std::pow;
285 const Scalar A1 = 0.25919; // from table
286 const Scalar A2 = 0.0014569; // from table
287 const Scalar expo = 1.0 + pow((1.0 - temp/criticalTemperature()), (2.0/7.0));
288 const Scalar V = A2*pow(A1, expo); // liquid molar volume [m^3/mol]
289
290 return 1.0/V; // molar density [mol/m^3]
291 }
292
300 {
302 }
303
307 static constexpr bool gasIsCompressible()
308 { return true; }
309
313 static constexpr bool gasIsIdeal()
314 { return true; }
315
319 static constexpr bool liquidIsCompressible()
320 { return false; }
321
329 {
330 using std::min;
331 using std::max;
332 temp = min(temp, 500.0); // regularization
333 temp = max(temp, 250.0); // regularization
334
335 using std::pow;
336 using std::exp;
337 const Scalar Tr = max(temp/criticalTemperature(), 1e-10);
338 const Scalar Fp0 = 1.0;
339 const Scalar xi = 0.004623;
340 const Scalar eta_xi = Fp0*(0.807*pow(Tr, 0.618)
341 - 0.357*exp(-0.449*Tr)
342 + 0.34*exp(-4.058*Tr)
343 + 0.018);
344 Scalar r = eta_xi/xi; // [1e-6 P]
345 r /= 1.0e7; // [Pa s]
346
347 return r;
348 }
349
357 {
358 using std::min;
359 using std::max;
360 temp = min(temp, 500.0); // regularization
361 temp = max(temp, 250.0); // regularization
362
363 const Scalar A = -3.82;
364 const Scalar B = 1027.0;
365 const Scalar C = -6.38e-4;
366 const Scalar D = 4.52e-7;
367
368 using std::exp;
369 Scalar r = exp(A + B/temp + C*temp + D*temp*temp); // in [cP]
370 r *= 1.0e-3; // in [Pa s]
371
372 return r; // [Pa s]
373 }
374
385 {
386 return 0.13;
387 }
388};
389
390} // end namespace Components
391
392} // end namespace Dumux
393
394#endif
Interface for components that have a gas state.
Interface for components that have a liquid state.
Relations valid for an ideal gas.
A central place for various physical constants occurring in some equations.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
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
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
Interface for components that have a liquid state.
Definition: liquid.hh:41
Properties of xylene.
Definition: xylene.hh:49
static Scalar gasViscosity(Scalar temp, Scalar pressure)
The dynamic viscosity of xylene vapor.
Definition: xylene.hh:328
static std::string name()
A human readable name for the xylene.
Definition: xylene.hh:57
static Scalar liquidDensity(Scalar temperature, Scalar pressure)
The density of pure xylene at a given pressure and temperature .
Definition: xylene.hh:299
static constexpr Scalar criticalPressure()
Returns the critical pressure of xylene.
Definition: xylene.hh:75
static Scalar liquidEnthalpy(const Scalar temperature, const Scalar pressure)
Specific enthalpy of liquid xylene .
Definition: xylene.hh:167
static Scalar gasMolarDensity(Scalar temperature, Scalar pressure)
The molar gas density of xylene gas at a given pressure and temperature.
Definition: xylene.hh:262
static constexpr bool gasIsIdeal()
Returns true if the gas phase is assumed to be ideal.
Definition: xylene.hh:313
static Scalar vaporPressure(Scalar temperature)
The saturation vapor pressure in of pure xylene at a given temperature according to Antoine after Be...
Definition: xylene.hh:106
static Scalar tripleTemperature()
Returns the temperature at xylene's triple point.
Definition: xylene.hh:87
static constexpr Scalar boilingTemperature()
Returns the temperature at xylene's boiling point (1 atm).
Definition: xylene.hh:81
static Scalar liquidThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of xylene.
Definition: xylene.hh:384
static constexpr bool liquidIsCompressible()
Returns true if the liquid phase is assumed to be compressible.
Definition: xylene.hh:319
static Scalar gasEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of xylene vapor .
Definition: xylene.hh:238
static constexpr Scalar criticalTemperature()
Returns the critical temperature of xylene.
Definition: xylene.hh:69
static Scalar triplePressure()
Returns the pressure at xylene's triple point.
Definition: xylene.hh:95
static Scalar heatVap(Scalar temperature, const Scalar pressure)
Latent heat of vaporization for xylene .
Definition: xylene.hh:202
static Scalar liquidHeatCapacity(Scalar temp, Scalar pressure)
Specific heat cap of liquid xylene .
Definition: xylene.hh:129
static constexpr Scalar molarMass()
The molar mass in of xylene.
Definition: xylene.hh:63
static Scalar gasDensity(Scalar temperature, Scalar pressure)
The density of xylene gas at a given pressure and temperature.
Definition: xylene.hh:249
static Scalar liquidMolarDensity(Scalar temp, Scalar pressure)
The molar liquid density of pure xylene at a given pressure and temperature .
Definition: xylene.hh:274
static constexpr bool gasIsCompressible()
Returns true if the gas phase is assumed to be compressible.
Definition: xylene.hh:307
static Scalar liquidViscosity(Scalar temp, Scalar pressure)
The dynamic viscosity of pure xylene.
Definition: xylene.hh:356
A central place for various physical constants occurring 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,...