Loading [MathJax]/extensions/tex2jax.js
Cantera  3.2.0a1
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
GasTransport.cpp
Go to the documentation of this file.
1//! @file GasTransport.cpp
2
3// This file is part of Cantera. See License.txt in the top-level directory or
4// at https://cantera.org/license.txt for license and copyright information.
5
7#include "MMCollisionInt.h"
14#include "cantera/base/global.h"
15
16namespace Cantera
17{
18
19//! polynomial degree used for fitting collision integrals
20//! except in CK mode, where the degree is 6.
21#define COLL_INT_POLY_DEGREE 8
22
23GasTransport::GasTransport() :
24 m_polytempvec(5)
25{
26}
27
28void GasTransport::update_T()
29{
30 if (m_thermo->nSpecies() != m_nsp) {
31 // Rebuild data structures if number of species has changed
32 init(m_thermo, m_mode);
33 }
34
35 double T = m_thermo->temperature();
36 if (T == m_temp) {
37 return;
38 }
39
40 m_temp = T;
41 m_kbt = Boltzmann * m_temp;
42 m_logt = log(m_temp);
43 m_sqrt_t = sqrt(m_temp);
44 m_t14 = sqrt(m_sqrt_t);
45
46 // compute powers of log(T)
47 m_polytempvec[0] = 1.0;
48 m_polytempvec[1] = m_logt;
49 m_polytempvec[2] = m_logt*m_logt;
50 m_polytempvec[3] = m_logt*m_logt*m_logt;
51 m_polytempvec[4] = m_logt*m_logt*m_logt*m_logt;
52
53 // temperature has changed, so polynomial fits will need to be redone
54 m_visc_ok = false;
55 m_spvisc_ok = false;
56 m_viscwt_ok = false;
57 m_bindiff_ok = false;
58}
59
60double GasTransport::viscosity()
61{
62 update_T();
63 update_C();
64
65 if (m_visc_ok) {
66 return m_viscmix;
67 }
68
69 double vismix = 0.0;
70 // update m_visc and m_phi if necessary
71 if (!m_viscwt_ok) {
72 updateViscosity_T();
73 }
74
75 multiply(m_phi, m_molefracs.data(), m_spwork.data());
76
77 for (size_t k = 0; k < m_nsp; k++) {
78 vismix += m_molefracs[k] * m_visc[k]/m_spwork[k]; //denom;
79 }
80 m_viscmix = vismix;
81 return vismix;
82}
83
84void GasTransport::updateViscosity_T()
85{
86 if (!m_spvisc_ok) {
87 updateSpeciesViscosities();
88 }
89
90 // see Eq. (9-5.14) of Poling et al. (2001)
91 for (size_t j = 0; j < m_nsp; j++) {
92 for (size_t k = j; k < m_nsp; k++) {
93 double vratiokj = m_visc[k]/m_visc[j];
94 double wratiojk = m_mw[j]/m_mw[k];
95
96 // Note that m_wratjk(k,j) holds the square root of m_wratjk(j,k)!
97 double factor1 = 1.0 + (m_sqvisc[k]/m_sqvisc[j]) * m_wratjk(k,j);
98 m_phi(k,j) = factor1*factor1 / (sqrt(8.0) * m_wratkj1(j,k));
99 m_phi(j,k) = m_phi(k,j)/(vratiokj * wratiojk);
100 }
101 }
102 m_viscwt_ok = true;
103}
104
105void GasTransport::updateSpeciesViscosities()
106{
107 update_T();
108 if (m_mode == CK_Mode) {
109 for (size_t k = 0; k < m_nsp; k++) {
110 m_visc[k] = exp(dot4(m_polytempvec, m_visccoeffs[k]));
111 m_sqvisc[k] = sqrt(m_visc[k]);
112 }
113 } else {
114 for (size_t k = 0; k < m_nsp; k++) {
115 // the polynomial fit is done for sqrt(visc/sqrt(T))
116 m_sqvisc[k] = m_t14 * dot5(m_polytempvec, m_visccoeffs[k]);
117 m_visc[k] = (m_sqvisc[k] * m_sqvisc[k]);
118 }
119 }
120 m_spvisc_ok = true;
121}
122
123void GasTransport::updateDiff_T()
124{
125 update_T();
126 // evaluate binary diffusion coefficients at unit pressure
127 size_t ic = 0;
128 if (m_mode == CK_Mode) {
129 for (size_t i = 0; i < m_nsp; i++) {
130 for (size_t j = i; j < m_nsp; j++) {
131 m_bdiff(i,j) = exp(dot4(m_polytempvec, m_diffcoeffs[ic]));
132 m_bdiff(j,i) = m_bdiff(i,j);
133 ic++;
134 }
135 }
136 } else {
137 for (size_t i = 0; i < m_nsp; i++) {
138 for (size_t j = i; j < m_nsp; j++) {
139 m_bdiff(i,j) = m_temp * m_sqrt_t*dot5(m_polytempvec,
140 m_diffcoeffs[ic]);
141 m_bdiff(j,i) = m_bdiff(i,j);
142 ic++;
143 }
144 }
145 }
146 m_bindiff_ok = true;
147}
148
149void GasTransport::getBinaryDiffCoeffs(const size_t ld, double* const d)
150{
151 update_T();
152 // if necessary, evaluate the binary diffusion coefficients from the polynomial fits
153 if (!m_bindiff_ok) {
154 updateDiff_T();
155 }
156 if (ld < m_nsp) {
157 throw CanteraError("GasTransport::getBinaryDiffCoeffs", "ld is too small");
158 }
159 double rp = 1.0/m_thermo->pressure();
160 for (size_t i = 0; i < m_nsp; i++) {
161 for (size_t j = 0; j < m_nsp; j++) {
162 d[ld*j + i] = rp * m_bdiff(i,j);
163 }
164 }
165}
166
167void GasTransport::getMixDiffCoeffs(double* const d)
168{
169 update_T();
170 update_C();
171
172 // update the binary diffusion coefficients if necessary
173 if (!m_bindiff_ok) {
174 updateDiff_T();
175 }
176
177 double mmw = m_thermo->meanMolecularWeight();
178 double p = m_thermo->pressure();
179 if (m_nsp == 1) {
180 d[0] = m_bdiff(0,0) / p;
181 } else {
182 for (size_t k = 0; k < m_nsp; k++) {
183 double sum2 = 0.0;
184 for (size_t j = 0; j < m_nsp; j++) {
185 if (j != k) {
186 sum2 += m_molefracs[j] / m_bdiff(j,k);
187 }
188 }
189 if (sum2 <= 0.0) {
190 d[k] = m_bdiff(k,k) / p;
191 } else {
192 d[k] = (mmw - m_molefracs[k] * m_mw[k])/(p * mmw * sum2);
193 }
194 }
195 }
196}
197
198void GasTransport::getMixDiffCoeffsMole(double* const d)
199{
200 update_T();
201 update_C();
202
203 // update the binary diffusion coefficients if necessary
204 if (!m_bindiff_ok) {
205 updateDiff_T();
206 }
207
208 double p = m_thermo->pressure();
209 if (m_nsp == 1) {
210 d[0] = m_bdiff(0,0) / p;
211 } else {
212 for (size_t k = 0; k < m_nsp; k++) {
213 double sum2 = 0.0;
214 for (size_t j = 0; j < m_nsp; j++) {
215 if (j != k) {
216 sum2 += m_molefracs[j] / m_bdiff(j,k);
217 }
218 }
219 if (sum2 <= 0.0) {
220 d[k] = m_bdiff(k,k) / p;
221 } else {
222 d[k] = (1 - m_molefracs[k]) / (p * sum2);
223 }
224 }
225 }
226}
227
228void GasTransport::getMixDiffCoeffsMass(double* const d)
229{
230 update_T();
231 update_C();
232
233 // update the binary diffusion coefficients if necessary
234 if (!m_bindiff_ok) {
235 updateDiff_T();
236 }
237
238 double mmw = m_thermo->meanMolecularWeight();
239 double p = m_thermo->pressure();
240
241 if (m_nsp == 1) {
242 d[0] = m_bdiff(0,0) / p;
243 } else {
244 for (size_t k=0; k<m_nsp; k++) {
245 double sum1 = 0.0;
246 double sum2 = 0.0;
247 for (size_t i=0; i<m_nsp; i++) {
248 if (i==k) {
249 continue;
250 }
251 sum1 += m_molefracs[i] / m_bdiff(k,i);
252 sum2 += m_molefracs[i] * m_mw[i] / m_bdiff(k,i);
253 }
254 sum1 *= p;
255 sum2 *= p * m_molefracs[k] / (mmw - m_mw[k]*m_molefracs[k]);
256 d[k] = 1.0 / (sum1 + sum2);
257 }
258 }
259}
260
261void GasTransport::init(ThermoPhase* thermo, int mode)
262{
263 m_thermo = thermo;
264 m_nsp = m_thermo->nSpecies();
265 m_mode = mode;
266
267 // set up Monchick and Mason collision integrals
268 setupCollisionParameters();
269 setupCollisionIntegral();
270
271 m_molefracs.resize(m_nsp);
272 m_spwork.resize(m_nsp);
273 m_visc.resize(m_nsp);
274 m_sqvisc.resize(m_nsp);
275 m_phi.resize(m_nsp, m_nsp, 0.0);
276 m_bdiff.resize(m_nsp, m_nsp);
277
278 // make a local copy of the molecular weights
279 m_mw = m_thermo->molecularWeights();
280
281 m_wratjk.resize(m_nsp, m_nsp, 0.0);
282 m_wratkj1.resize(m_nsp, m_nsp, 0.0);
283 for (size_t j = 0; j < m_nsp; j++) {
284 for (size_t k = j; k < m_nsp; k++) {
285 m_wratjk(j,k) = sqrt(m_mw[j]/m_mw[k]);
286 m_wratjk(k,j) = sqrt(m_wratjk(j,k));
287 m_wratkj1(j,k) = sqrt(1.0 + m_mw[k]/m_mw[j]);
288 }
289 }
290}
291
292void GasTransport::setupCollisionParameters()
293{
294 m_epsilon.resize(m_nsp, m_nsp, 0.0);
295 m_delta.resize(m_nsp, m_nsp, 0.0);
296 m_reducedMass.resize(m_nsp, m_nsp, 0.0);
297 m_dipole.resize(m_nsp, m_nsp, 0.0);
298 m_diam.resize(m_nsp, m_nsp, 0.0);
299 m_crot.resize(m_nsp);
300 m_zrot.resize(m_nsp);
301 m_polar.resize(m_nsp, false);
302 m_alpha.resize(m_nsp, 0.0);
303 m_poly.resize(m_nsp);
304 m_star_poly_uses_actualT.resize(m_nsp);
305 m_sigma.resize(m_nsp);
306 m_eps.resize(m_nsp);
307 m_w_ac.resize(m_nsp);
308 m_disp.resize(m_nsp, 0.0);
309 m_quad_polar.resize(m_nsp, 0.0);
310
311 const vector<double>& mw = m_thermo->molecularWeights();
312 getTransportData();
313
314 for (size_t i = 0; i < m_nsp; i++) {
315 m_poly[i].resize(m_nsp);
316 m_star_poly_uses_actualT[i].resize(m_nsp);
317 }
318
319 double f_eps, f_sigma;
320
321 for (size_t i = 0; i < m_nsp; i++) {
322 for (size_t j = i; j < m_nsp; j++) {
323 // the reduced mass
324 m_reducedMass(i,j) = mw[i] * mw[j] / (Avogadro * (mw[i] + mw[j]));
325
326 // hard-sphere diameter for (i,j) collisions
327 m_diam(i,j) = 0.5*(m_sigma[i] + m_sigma[j]);
328
329 // the effective well depth for (i,j) collisions
330 m_epsilon(i,j) = sqrt(m_eps[i]*m_eps[j]);
331
332 // the effective dipole moment for (i,j) collisions
333 m_dipole(i,j) = sqrt(m_dipole(i,i)*m_dipole(j,j));
334
335 // reduced dipole moment delta* (nondimensional)
336 double d = m_diam(i,j);
337 m_delta(i,j) = 0.5 * m_dipole(i,j)*m_dipole(i,j)
338 / (4 * Pi * epsilon_0 * m_epsilon(i,j) * d * d * d);
339 makePolarCorrections(i, j, f_eps, f_sigma);
340 m_diam(i,j) *= f_sigma;
341 m_epsilon(i,j) *= f_eps;
342
343 // properties are symmetric
344 m_reducedMass(j,i) = m_reducedMass(i,j);
345 m_diam(j,i) = m_diam(i,j);
346 m_epsilon(j,i) = m_epsilon(i,j);
347 m_dipole(j,i) = m_dipole(i,j);
348 m_delta(j,i) = m_delta(i,j);
349 }
350 }
351}
352
353void GasTransport::invalidateCache()
354{
355 Transport::invalidateCache();
356 m_temp = Undef;
357 m_visc_ok = false;
358 m_spvisc_ok = false;
359 m_viscwt_ok = false;
360 m_bindiff_ok = false;
361}
362
363void GasTransport::setupCollisionIntegral()
364{
365 double tstar_min = 1.e8, tstar_max = 0.0;
366 for (size_t i = 0; i < m_nsp; i++) {
367 for (size_t j = i; j < m_nsp; j++) {
368 // The polynomial fits of collision integrals vs. T*
369 // will be done for the T* from tstar_min to tstar_max
370 tstar_min = std::min(tstar_min, Boltzmann * m_thermo->minTemp()/m_epsilon(i,j));
371 tstar_max = std::max(tstar_max, Boltzmann * m_thermo->maxTemp()/m_epsilon(i,j));
372 }
373 }
374 // Chemkin fits the entire T* range in the Monchick and Mason tables,
375 // so modify tstar_min and tstar_max if in Chemkin compatibility mode
376 if (m_mode == CK_Mode) {
377 tstar_min = 0.101;
378 tstar_max = 99.9;
379 }
380
381 // initialize the collision integral calculator for the desired T* range
382 MMCollisionInt integrals;
383 integrals.init(tstar_min, tstar_max);
384 fitCollisionIntegrals(integrals);
385 // make polynomial fits
386 fitProperties(integrals);
387}
388
389void GasTransport::getTransportData()
390{
391 for (size_t k = 0; k < m_thermo->nSpecies(); k++) {
392 shared_ptr<Species> s = m_thermo->species(k);
393 const GasTransportData* sptran =
394 dynamic_cast<GasTransportData*>(s->transport.get());
395 if (!sptran) {
396 throw CanteraError("GasTransport::getTransportData",
397 "Missing gas-phase transport data for species '{}'.", s->name);
398 }
399
400 if (sptran->geometry == "atom") {
401 m_crot[k] = 0.0;
402 } else if (sptran->geometry == "linear") {
403 m_crot[k] = 1.0;
404 } else if (sptran->geometry == "nonlinear") {
405 m_crot[k] = 1.5;
406 }
407
408 m_sigma[k] = sptran->diameter;
409 m_eps[k] = sptran->well_depth;
410 m_dipole(k,k) = sptran->dipole;
411 m_polar[k] = (sptran->dipole > 0);
412 m_alpha[k] = sptran->polarizability;
413 m_zrot[k] = sptran->rotational_relaxation;
414 if (s->input.hasKey("critical-parameters") &&
415 s->input["critical-parameters"].hasKey("acentric-factor"))
416 {
417 m_w_ac[k] = s->input["critical-parameters"]["acentric-factor"].asDouble();
418 } else {
419 m_w_ac[k] = sptran->acentric_factor;
420 }
421 m_disp[k] = sptran->dispersion_coefficient;
422 m_quad_polar[k] = sptran->quadrupole_polarizability;
423 }
424}
425
426void GasTransport::makePolarCorrections(size_t i, size_t j,
427 double& f_eps, double& f_sigma)
428{
429 // no correction if both are nonpolar, or both are polar
430 if (m_polar[i] == m_polar[j]) {
431 f_eps = 1.0;
432 f_sigma = 1.0;
433 return;
434 }
435
436 // corrections to the effective diameter and well depth
437 // if one is polar and one is non-polar
438 size_t kp = (m_polar[i] ? i : j); // the polar one
439 size_t knp = (i == kp ? j : i); // the nonpolar one
440 double d3np, d3p, alpha_star, mu_p_star, xi;
441 d3np = pow(m_sigma[knp],3);
442 d3p = pow(m_sigma[kp],3);
443 alpha_star = m_alpha[knp]/d3np;
444 mu_p_star = m_dipole(kp,kp)/sqrt(4 * Pi * epsilon_0 * d3p * m_eps[kp]);
445 xi = 1.0 + 0.25 * alpha_star * mu_p_star * mu_p_star *
446 sqrt(m_eps[kp]/m_eps[knp]);
447 f_sigma = pow(xi, -1.0/6.0);
448 f_eps = xi*xi;
449}
450
451void GasTransport::fitCollisionIntegrals(MMCollisionInt& integrals)
452{
453 // Chemkin fits to sixth order polynomials
454 int degree = (m_mode == CK_Mode ? 6 : COLL_INT_POLY_DEGREE);
455 vector<double> fitlist;
456 m_omega22_poly.clear();
457 m_astar_poly.clear();
458 m_bstar_poly.clear();
459 m_cstar_poly.clear();
460 for (size_t i = 0; i < m_nsp; i++) {
461 for (size_t j = i; j < m_nsp; j++) {
462 // Chemkin fits only delta* = 0
463 double dstar = (m_mode != CK_Mode) ? m_delta(i,j) : 0.0;
464
465 // if a fit has already been generated for delta* = m_delta(i,j),
466 // then use it. Otherwise, make a new fit, and add m_delta(i,j) to
467 // the list of delta* values for which fits have been done.
468
469 // 'find' returns a pointer to end() if not found
470 auto dptr = find(fitlist.begin(), fitlist.end(), dstar);
471 if (dptr == fitlist.end()) {
472 vector<double> ca(degree+1), cb(degree+1), cc(degree+1);
473 vector<double> co22(degree+1);
474 integrals.fit(degree, dstar, ca.data(), cb.data(), cc.data());
475 integrals.fit_omega22(degree, dstar, co22.data());
476 m_omega22_poly.push_back(co22);
477 m_astar_poly.push_back(ca);
478 m_bstar_poly.push_back(cb);
479 m_cstar_poly.push_back(cc);
480 m_poly[i][j] = static_cast<int>(m_astar_poly.size()) - 1;
481 m_star_poly_uses_actualT[i][j] = 0;
482 fitlist.push_back(dstar);
483 } else {
484 // delta* found in fitlist, so just point to this polynomial
485 m_poly[i][j] = static_cast<int>((dptr - fitlist.begin()));
486 m_star_poly_uses_actualT[i][j] = 0;
487 }
488 m_poly[j][i] = m_poly[i][j];
489 m_star_poly_uses_actualT[j][i] = m_star_poly_uses_actualT[i][j];
490 }
491 }
492}
493
494void GasTransport::fitProperties(MMCollisionInt& integrals)
495{
496 // number of points to use in generating fit data
497 const size_t np = 50;
498 int degree = (m_mode == CK_Mode ? 3 : 4);
499 double dt = (m_thermo->maxTemp() - m_thermo->minTemp())/(np-1);
500 vector<double> tlog(np), spvisc(np), spcond(np);
501 vector<double> w(np), w2(np);
502
503 m_visccoeffs.clear();
504 m_condcoeffs.clear();
505
506 // generate array of log(t) values
507 for (size_t n = 0; n < np; n++) {
508 double t = m_thermo->minTemp() + dt*n;
509 tlog[n] = log(t);
510 }
511
512 // vector of polynomial coefficients
513 vector<double> c(degree + 1), c2(degree + 1);
514
515 // fit the pure-species viscosity and thermal conductivity for each species
516 double visc, err, relerr,
517 mxerr = 0.0, mxrelerr = 0.0, mxerr_cond = 0.0, mxrelerr_cond = 0.0;
518 double T_save = m_thermo->temperature();
519 const vector<double>& mw = m_thermo->molecularWeights();
520 for (size_t k = 0; k < m_nsp; k++) {
521 double tstar = Boltzmann * 298.0 / m_eps[k];
522 // Scaling factor for temperature dependence of z_rot. [Kee2003] Eq.
523 // 12.112 or [Kee2017] Eq. 11.115
524 double fz_298 = 1.0 + pow(Pi, 1.5) / sqrt(tstar) * (0.5 + 1.0 / tstar) +
525 (0.25 * Pi * Pi + 2) / tstar;
526
527 for (size_t n = 0; n < np; n++) {
528 double t = m_thermo->minTemp() + dt*n;
529 m_thermo->setTemperature(t);
530 vector<double> cp_R_all(m_thermo->nSpecies());
531 m_thermo->getCp_R_ref(&cp_R_all[0]);
532 double cp_R = cp_R_all[k];
533 tstar = Boltzmann * t / m_eps[k];
534 double sqrt_T = sqrt(t);
535 double om22 = integrals.omega22(tstar, m_delta(k,k));
536 double om11 = integrals.omega11(tstar, m_delta(k,k));
537
538 // self-diffusion coefficient, without polar corrections
539 double diffcoeff = 3.0/16.0 * sqrt(2.0 * Pi/m_reducedMass(k,k)) *
540 pow((Boltzmann * t), 1.5)/
541 (Pi * m_sigma[k] * m_sigma[k] * om11);
542
543 // viscosity
544 visc = 5.0/16.0 * sqrt(Pi * mw[k] * Boltzmann * t / Avogadro) /
545 (om22 * Pi * m_sigma[k]*m_sigma[k]);
546
547 // thermal conductivity
548 double f_int = mw[k]/(GasConstant * t) * diffcoeff/visc;
549 double cv_rot = m_crot[k];
550 double A_factor = 2.5 - f_int;
551 double fz_tstar = 1.0 + pow(Pi, 1.5) / sqrt(tstar) * (0.5 + 1.0 / tstar) +
552 (0.25 * Pi * Pi + 2) / tstar;
553 double B_factor = m_zrot[k] * fz_298 / fz_tstar + 2.0/Pi * (5.0/3.0 * cv_rot + f_int);
554 double c1 = 2.0/Pi * A_factor/B_factor;
555 double cv_int = cp_R - 2.5 - cv_rot;
556 double f_rot = f_int * (1.0 + c1);
557 double f_trans = 2.5 * (1.0 - c1 * cv_rot/1.5);
558 double cond = (visc/mw[k])*GasConstant*(f_trans * 1.5
559 + f_rot * cv_rot + f_int * cv_int);
560
561 if (m_mode == CK_Mode) {
562 spvisc[n] = log(visc);
563 spcond[n] = log(cond);
564 w[n] = -1.0;
565 w2[n] = -1.0;
566 } else {
567 // the viscosity should be proportional approximately to
568 // sqrt(T); therefore, visc/sqrt(T) should have only a weak
569 // temperature dependence. And since the mixture rule requires
570 // the square root of the pure-species viscosity, fit the square
571 // root of (visc/sqrt(T)) to avoid having to compute square
572 // roots in the mixture rule.
573 spvisc[n] = sqrt(visc/sqrt_T);
574
575 // the pure-species conductivity scales approximately with
576 // sqrt(T). Unlike the viscosity, there is no reason here to fit
577 // the square root, since a different mixture rule is used.
578 spcond[n] = cond/sqrt_T;
579 w[n] = 1.0/(spvisc[n]*spvisc[n]);
580 w2[n] = 1.0/(spcond[n]*spcond[n]);
581 }
582 }
583 polyfit(np, degree, tlog.data(), spvisc.data(), w.data(), c.data());
584 polyfit(np, degree, tlog.data(), spcond.data(), w2.data(), c2.data());
585
586 // evaluate max fit errors for viscosity
587 for (size_t n = 0; n < np; n++) {
588 double val, fit;
589 if (m_mode == CK_Mode) {
590 val = exp(spvisc[n]);
591 fit = exp(poly3(tlog[n], c.data()));
592 } else {
593 double sqrt_T = exp(0.5*tlog[n]);
594 val = sqrt_T * pow(spvisc[n],2);
595 fit = sqrt_T * pow(poly4(tlog[n], c.data()),2);
596 }
597 err = fit - val;
598 relerr = err/val;
599 mxerr = std::max(mxerr, fabs(err));
600 mxrelerr = std::max(mxrelerr, fabs(relerr));
601 }
602 m_fittingErrors["viscosity-max-abs-error"] = mxerr;
603 m_fittingErrors["viscosity-max-rel-error"] = mxrelerr;
604
605 // evaluate max fit errors for conductivity
606 for (size_t n = 0; n < np; n++) {
607 double val, fit;
608 if (m_mode == CK_Mode) {
609 val = exp(spcond[n]);
610 fit = exp(poly3(tlog[n], c2.data()));
611 } else {
612 double sqrt_T = exp(0.5*tlog[n]);
613 val = sqrt_T * spcond[n];
614 fit = sqrt_T * poly4(tlog[n], c2.data());
615 }
616 err = fit - val;
617 relerr = err/val;
618 mxerr_cond = std::max(mxerr_cond, fabs(err));
619 mxrelerr_cond = std::max(mxrelerr_cond, fabs(relerr));
620 }
621 m_visccoeffs.push_back(c);
622 m_condcoeffs.push_back(c2);
623 m_fittingErrors["conductivity-max-abs-error"] = mxerr_cond;
624 m_fittingErrors["conductivity-max-rel-error"] = mxrelerr_cond;
625 }
626 m_thermo->setTemperature(T_save);
627 fitDiffCoeffs(integrals);
628}
629
630void GasTransport::fitDiffCoeffs(MMCollisionInt& integrals)
631{
632 // number of points to use in generating fit data
633 const size_t np = 50;
634 int degree = (m_mode == CK_Mode ? 3 : 4);
635 double dt = (m_thermo->maxTemp() - m_thermo->minTemp())/(np-1);
636 vector<double> tlog(np);
637 vector<double> w(np), w2(np);
638
639 // generate array of log(t) values
640 for (size_t n = 0; n < np; n++) {
641 double t = m_thermo->minTemp() + dt*n;
642 tlog[n] = log(t);
643 }
644
645 // vector of polynomial coefficients
646 vector<double> c(degree + 1), c2(degree + 1);
647 double err, relerr, mxerr = 0.0, mxrelerr = 0.0;
648
649 vector<double> diff(np + 1);
650 m_diffcoeffs.clear();
651 for (size_t k = 0; k < m_nsp; k++) {
652 for (size_t j = k; j < m_nsp; j++) {
653 for (size_t n = 0; n < np; n++) {
654 double t = m_thermo->minTemp() + dt*n;
655 double eps = m_epsilon(j,k);
656 double tstar = Boltzmann * t/eps;
657 double sigma = m_diam(j,k);
658 double om11 = integrals.omega11(tstar, m_delta(j,k));
659 double diffcoeff = 3.0/16.0 * sqrt(2.0 * Pi/m_reducedMass(k,j))
660 * pow(Boltzmann * t, 1.5) / (Pi * sigma * sigma * om11);
661
662 // 2nd order correction
663 // NOTE: THIS CORRECTION IS NOT APPLIED
664 double fkj, fjk;
665 getBinDiffCorrection(t, integrals, k, j, 1.0, 1.0, fkj, fjk);
666
667 if (m_mode == CK_Mode) {
668 diff[n] = log(diffcoeff);
669 w[n] = -1.0;
670 } else {
671 diff[n] = diffcoeff/pow(t, 1.5);
672 w[n] = 1.0/(diff[n]*diff[n]);
673 }
674 }
675 polyfit(np, degree, tlog.data(), diff.data(), w.data(), c.data());
676
677 for (size_t n = 0; n < np; n++) {
678 double val, fit;
679 if (m_mode == CK_Mode) {
680 val = exp(diff[n]);
681 fit = exp(poly3(tlog[n], c.data()));
682 } else {
683 double t = exp(tlog[n]);
684 double pre = pow(t, 1.5);
685 val = pre * diff[n];
686 fit = pre * poly4(tlog[n], c.data());
687 }
688 err = fit - val;
689 relerr = err/val;
690 mxerr = std::max(mxerr, fabs(err));
691 mxrelerr = std::max(mxrelerr, fabs(relerr));
692 }
693 m_diffcoeffs.push_back(c);
694 }
695 }
696
697 m_fittingErrors["diff-coeff-max-abs-error"] = mxerr;
698 m_fittingErrors["diff-coeff-max-rel-error"] = mxrelerr;
699}
700
701void GasTransport::getBinDiffCorrection(double t, MMCollisionInt& integrals,
702 size_t k, size_t j, double xk, double xj, double& fkj, double& fjk)
703{
704 double w1 = m_thermo->molecularWeight(k);
705 double w2 = m_thermo->molecularWeight(j);
706 double wsum = w1 + w2;
707 double wmwp = (w1 - w2)/wsum;
708 double sqw12 = sqrt(w1*w2);
709 double sig1 = m_sigma[k];
710 double sig2 = m_sigma[j];
711 double sig12 = 0.5*(m_sigma[k] + m_sigma[j]);
712 double sigratio = sig1*sig1/(sig2*sig2);
713 double sigratio2 = sig1*sig1/(sig12*sig12);
714 double sigratio3 = sig2*sig2/(sig12*sig12);
715 double tstar1 = Boltzmann * t / m_eps[k];
716 double tstar2 = Boltzmann * t / m_eps[j];
717 double tstar12 = Boltzmann * t / sqrt(m_eps[k] * m_eps[j]);
718 double om22_1 = integrals.omega22(tstar1, m_delta(k,k));
719 double om22_2 = integrals.omega22(tstar2, m_delta(j,j));
720 double om11_12 = integrals.omega11(tstar12, m_delta(k,j));
721 double astar_12 = integrals.astar(tstar12, m_delta(k,j));
722 double bstar_12 = integrals.bstar(tstar12, m_delta(k,j));
723 double cstar_12 = integrals.cstar(tstar12, m_delta(k,j));
724
725 double cnst = sigratio * sqrt(2.0*w2/wsum) * 2.0 * w1*w1/(wsum * w2);
726 double p1 = cnst * om22_1 / om11_12;
727
728 cnst = (1.0/sigratio) * sqrt(2.0*w1/wsum) * 2.0*w2*w2/(wsum*w1);
729 double p2 = cnst * om22_2 / om11_12;
730 double p12 = 15.0 * wmwp*wmwp + 8.0*w1*w2*astar_12/(wsum*wsum);
731
732 cnst = (2.0/(w2*wsum))*sqrt(2.0*w2/wsum)*sigratio2;
733 double q1 = cnst*((2.5 - 1.2*bstar_12)*w1*w1 + 3.0*w2*w2
734 + 1.6*w1*w2*astar_12);
735
736 cnst = (2.0/(w1*wsum))*sqrt(2.0*w1/wsum)*sigratio3;
737 double q2 = cnst*((2.5 - 1.2*bstar_12)*w2*w2 + 3.0*w1*w1
738 + 1.6*w1*w2*astar_12);
739 double q12 = wmwp*wmwp*15.0*(2.5 - 1.2*bstar_12)
740 + 4.0*w1*w2*astar_12*(11.0 - 2.4*bstar_12)/(wsum*wsum)
741 + 1.6*wsum*om22_1*om22_2/(om11_12*om11_12*sqw12)
742 * sigratio2 * sigratio3;
743
744 cnst = 6.0*cstar_12 - 5.0;
745 fkj = 1.0 + 0.1*cnst*cnst *
746 (p1*xk*xk + p2*xj*xj + p12*xk*xj)/
747 (q1*xk*xk + q2*xj*xj + q12*xk*xj);
748 fjk = 1.0 + 0.1*cnst*cnst *
749 (p2*xk*xk + p1*xj*xj + p12*xk*xj)/
750 (q2*xk*xk + q1*xj*xj + q12*xk*xj);
751}
752
753void GasTransport::getViscosityPolynomial(size_t i, double* coeffs) const
754{
755 checkSpeciesIndex(i);
756 for (int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
757 coeffs[k] = m_visccoeffs[i][k];
758 }
759}
760
761void GasTransport::getConductivityPolynomial(size_t i, double* coeffs) const
762{
763 checkSpeciesIndex(i);
764 for (int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
765 coeffs[k] = m_condcoeffs[i][k];
766 }
767}
768
769void GasTransport::getBinDiffusivityPolynomial(size_t i, size_t j, double* coeffs) const
770{
771 checkSpeciesIndex(i);
772 checkSpeciesIndex(j);
773 size_t mi = (j >= i? i : j);
774 size_t mj = (j >= i? j : i);
775 size_t ic = 0;
776 for (size_t ii = 0; ii < mi; ii++) {
777 ic += m_nsp - ii;
778 }
779 ic += mj - mi;
780
781 for (int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
782 coeffs[k] = m_diffcoeffs[ic][k];
783 }
784}
785
786void GasTransport::getCollisionIntegralPolynomial(size_t i, size_t j,
787 double* astar_coeffs,
788 double* bstar_coeffs,
789 double* cstar_coeffs) const
790{
791 checkSpeciesIndex(i);
792 checkSpeciesIndex(j);
793 for (int k = 0; k < (m_mode == CK_Mode ? 6 : COLL_INT_POLY_DEGREE) + 1; k++) {
794 astar_coeffs[k] = m_astar_poly[m_poly[i][j]][k];
795 bstar_coeffs[k] = m_bstar_poly[m_poly[i][j]][k];
796 cstar_coeffs[k] = m_cstar_poly[m_poly[i][j]][k];
797 }
798}
799
800void GasTransport::setViscosityPolynomial(size_t i, double* coeffs)
801{
802 checkSpeciesIndex(i);
803 for (int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
804 m_visccoeffs[i][k] = coeffs[k];
805 }
806 invalidateCache();
807}
808
809void GasTransport::setConductivityPolynomial(size_t i, double* coeffs)
810{
811 checkSpeciesIndex(i);
812 for (int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
813 m_condcoeffs[i][k] = coeffs[k];
814 }
815 invalidateCache();
816}
817
818void GasTransport::setBinDiffusivityPolynomial(size_t i, size_t j, double* coeffs)
819{
820 checkSpeciesIndex(i);
821 checkSpeciesIndex(j);
822 size_t mi = (j >= i? i : j);
823 size_t mj = (j >= i? j : i);
824 size_t ic = 0;
825 for (size_t ii = 0; ii < mi; ii++) {
826 ic += m_nsp - ii;
827 }
828 ic += mj - mi;
829
830 for (int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
831 m_diffcoeffs[ic][k] = coeffs[k];
832 }
833 invalidateCache();
834}
835
836void GasTransport::setCollisionIntegralPolynomial(size_t i, size_t j,
837 double* astar_coeffs,
838 double* bstar_coeffs,
839 double* cstar_coeffs, bool actualT)
840{
841 checkSpeciesIndex(i);
842 checkSpeciesIndex(j);
843 size_t degree = (m_mode == CK_Mode ? 6 : COLL_INT_POLY_DEGREE);
844 vector<double> ca(degree+1), cb(degree+1), cc(degree+1);
845
846 for (size_t k = 0; k < degree+1; k++) {
847 ca[k] = astar_coeffs[k];
848 cb[k] = bstar_coeffs[k];
849 cc[k] = cstar_coeffs[k];
850 }
851
852 m_astar_poly.push_back(ca);
853 m_bstar_poly.push_back(cb);
854 m_cstar_poly.push_back(cc);
855 m_poly[j][i] = m_poly[i][j] = static_cast<int>(m_astar_poly.size()) - 1;
856 m_star_poly_uses_actualT[j][i] = m_star_poly_uses_actualT[i][j] = (actualT) ? 1 : 0;
857 invalidateCache();
858}
859
860}
#define COLL_INT_POLY_DEGREE
polynomial degree used for fitting collision integrals except in CK mode, where the degree is 6.
Monchick and Mason collision integrals.
Declaration for class Cantera::Species.
Header file for class ThermoPhase, the base class for phases with thermodynamic properties,...
Base class for exceptions thrown by Cantera classes.
Transport data for a single gas-phase species which can be used in mixture-averaged or multicomponent...
double polarizability
The polarizability of the molecule [m^3]. Default 0.0.
double diameter
The Lennard-Jones collision diameter [m].
double acentric_factor
Pitzer's acentric factor [dimensionless]. Default 0.0.
double quadrupole_polarizability
quadrupole. Default 0.0.
double rotational_relaxation
The rotational relaxation number (the number of collisions it takes to equilibrate the rotational deg...
double dispersion_coefficient
dispersion normalized by e^2. [m^5] Default 0.0.
double dipole
The permanent dipole moment of the molecule [Coulomb-m]. Default 0.0.
double well_depth
The Lennard-Jones well depth [J].
string geometry
A string specifying the molecular geometry.
Calculation of Collision integrals.
void init(double tsmin, double tsmax)
Initialize the object for calculation.
size_t nSpecies() const
Returns the number of species in the phase.
Definition Phase.h:232
Base class for a phase with thermodynamic properties.
This file contains definitions for utility functions and text for modules, inputfiles and logging,...
double dot5(const V &x, const V &y)
Templated Inner product of two vectors of length 5.
Definition utilities.h:55
R poly4(D x, R *c)
Evaluates a polynomial of order 4.
Definition utilities.h:153
R poly3(D x, R *c)
Templated evaluation of a polynomial of order 3.
Definition utilities.h:165
double dot4(const V &x, const V &y)
Templated Inner product of two vectors of length 4.
Definition utilities.h:40
double polyfit(size_t n, size_t deg, const double *xp, const double *yp, const double *wp, double *pp)
Fits a polynomial function to a set of data points.
Definition polyfit.cpp:14
const double Boltzmann
Boltzmann constant [J/K].
Definition ct_defs.h:84
const double Avogadro
Avogadro's Number [number/kmol].
Definition ct_defs.h:81
const double epsilon_0
Permittivity of free space [F/m].
Definition ct_defs.h:137
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition ct_defs.h:120
const double Pi
Pi.
Definition ct_defs.h:68
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
const double Undef
Fairly random number to be used to initialize variables against to see if they are subsequently defin...
Definition ct_defs.h:164
void multiply(const DenseMatrix &A, const double *const b, double *const prod)
Multiply A*b and return the result in prod. Uses BLAS routine DGEMV.
Contains declarations for string manipulation functions within Cantera.
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...