Cantera  4.0.0a1
Loading...
Searching...
No Matches
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, m_spwork);
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, span<double> 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 checkArraySize("GasTransport::getBinaryDiffCoeffs", d.size(), ld * m_nsp);
160 double rp = 1.0/m_thermo->pressure();
161 for (size_t i = 0; i < m_nsp; i++) {
162 for (size_t j = 0; j < m_nsp; j++) {
163 d[ld*j + i] = rp * m_bdiff(i,j);
164 }
165 }
166}
167
168void GasTransport::getMixDiffCoeffs(span<double> d)
169{
170 checkArraySize("GasTransport::getMixDiffCoeffs", d.size(), m_nsp);
171 update_T();
172 update_C();
173
174 // update the binary diffusion coefficients if necessary
175 if (!m_bindiff_ok) {
176 updateDiff_T();
177 }
178
179 double mmw = m_thermo->meanMolecularWeight();
180 double p = m_thermo->pressure();
181 if (m_nsp == 1) {
182 d[0] = m_bdiff(0,0) / p;
183 } else {
184 for (size_t k = 0; k < m_nsp; k++) {
185 double sum2 = 0.0;
186 for (size_t j = 0; j < m_nsp; j++) {
187 if (j != k) {
188 sum2 += m_molefracs[j] / m_bdiff(j,k);
189 }
190 }
191 if (sum2 <= 0.0) {
192 d[k] = m_bdiff(k,k) / p;
193 } else {
194 d[k] = (mmw - m_molefracs[k] * m_mw[k])/(p * mmw * sum2);
195 }
196 }
197 }
198}
199
200void GasTransport::getMixDiffCoeffsMole(span<double> d)
201{
202 checkArraySize("GasTransport::getMixDiffCoeffsMole", d.size(), m_nsp);
203 update_T();
204 update_C();
205
206 // update the binary diffusion coefficients if necessary
207 if (!m_bindiff_ok) {
208 updateDiff_T();
209 }
210
211 double p = m_thermo->pressure();
212 if (m_nsp == 1) {
213 d[0] = m_bdiff(0,0) / p;
214 } else {
215 for (size_t k = 0; k < m_nsp; k++) {
216 double sum2 = 0.0;
217 for (size_t j = 0; j < m_nsp; j++) {
218 if (j != k) {
219 sum2 += m_molefracs[j] / m_bdiff(j,k);
220 }
221 }
222 if (sum2 <= 0.0) {
223 d[k] = m_bdiff(k,k) / p;
224 } else {
225 d[k] = (1 - m_molefracs[k]) / (p * sum2);
226 }
227 }
228 }
229}
230
231void GasTransport::getMixDiffCoeffsMass(span<double> d)
232{
233 checkArraySize("GasTransport::getMixDiffCoeffsMass", d.size(), m_nsp);
234 update_T();
235 update_C();
236
237 // update the binary diffusion coefficients if necessary
238 if (!m_bindiff_ok) {
239 updateDiff_T();
240 }
241
242 double mmw = m_thermo->meanMolecularWeight();
243 double p = m_thermo->pressure();
244
245 if (m_nsp == 1) {
246 d[0] = m_bdiff(0,0) / p;
247 } else {
248 for (size_t k=0; k<m_nsp; k++) {
249 double sum1 = 0.0;
250 double sum2 = 0.0;
251 for (size_t i=0; i<m_nsp; i++) {
252 if (i==k) {
253 continue;
254 }
255 sum1 += m_molefracs[i] / m_bdiff(k,i);
256 sum2 += m_molefracs[i] * m_mw[i] / m_bdiff(k,i);
257 }
258 sum1 *= p;
259 sum2 *= p * m_molefracs[k] / (mmw - m_mw[k]*m_molefracs[k]);
260 d[k] = 1.0 / (sum1 + sum2);
261 }
262 }
263}
264
265void GasTransport::init(shared_ptr<ThermoPhase> thermo, int mode)
266{
267 m_thermo = thermo;
268 m_nsp = m_thermo->nSpecies();
269 m_mode = mode;
270
271 // set up Monchick and Mason collision integrals
272 setupCollisionParameters();
273 setupCollisionIntegral();
274
275 m_molefracs.resize(m_nsp);
276 m_spwork.resize(m_nsp);
277 m_visc.resize(m_nsp);
278 m_sqvisc.resize(m_nsp);
279 m_phi.resize(m_nsp, m_nsp, 0.0);
280 m_bdiff.resize(m_nsp, m_nsp);
281
282 // make a local copy of the molecular weights
283 m_mw.resize(m_nsp);
284 m_thermo->getMolecularWeights(m_mw);
285
286 m_wratjk.resize(m_nsp, m_nsp, 0.0);
287 m_wratkj1.resize(m_nsp, m_nsp, 0.0);
288 for (size_t j = 0; j < m_nsp; j++) {
289 for (size_t k = j; k < m_nsp; k++) {
290 m_wratjk(j,k) = sqrt(m_mw[j]/m_mw[k]);
291 m_wratjk(k,j) = sqrt(m_wratjk(j,k));
292 m_wratkj1(j,k) = sqrt(1.0 + m_mw[k]/m_mw[j]);
293 }
294 }
295}
296
297void GasTransport::setupCollisionParameters()
298{
299 m_epsilon.resize(m_nsp, m_nsp, 0.0);
300 m_delta.resize(m_nsp, m_nsp, 0.0);
301 m_reducedMass.resize(m_nsp, m_nsp, 0.0);
302 m_dipole.resize(m_nsp, m_nsp, 0.0);
303 m_diam.resize(m_nsp, m_nsp, 0.0);
304 m_crot.resize(m_nsp);
305 m_zrot.resize(m_nsp);
306 m_polar.resize(m_nsp, false);
307 m_alpha.resize(m_nsp, 0.0);
308 m_poly.resize(m_nsp);
309 m_star_poly_uses_actualT.resize(m_nsp);
310 m_sigma.resize(m_nsp);
311 m_eps.resize(m_nsp);
312 m_w_ac.resize(m_nsp);
313 m_disp.resize(m_nsp, 0.0);
314 m_quad_polar.resize(m_nsp, 0.0);
315
316 auto mw = m_thermo->molecularWeights();
317 getTransportData();
318
319 for (size_t i = 0; i < m_nsp; i++) {
320 m_poly[i].resize(m_nsp);
321 m_star_poly_uses_actualT[i].resize(m_nsp);
322 }
323
324 double f_eps, f_sigma;
325
326 for (size_t i = 0; i < m_nsp; i++) {
327 for (size_t j = i; j < m_nsp; j++) {
328 // the reduced mass
329 m_reducedMass(i,j) = mw[i] * mw[j] / (Avogadro * (mw[i] + mw[j]));
330
331 // hard-sphere diameter for (i,j) collisions
332 m_diam(i,j) = 0.5*(m_sigma[i] + m_sigma[j]);
333
334 // the effective well depth for (i,j) collisions
335 m_epsilon(i,j) = sqrt(m_eps[i]*m_eps[j]);
336
337 // the effective dipole moment for (i,j) collisions
338 m_dipole(i,j) = sqrt(m_dipole(i,i)*m_dipole(j,j));
339
340 // reduced dipole moment delta* (nondimensional)
341 double d = m_diam(i,j);
342 m_delta(i,j) = 0.5 * m_dipole(i,j)*m_dipole(i,j)
343 / (4 * Pi * epsilon_0 * m_epsilon(i,j) * d * d * d);
344 makePolarCorrections(i, j, f_eps, f_sigma);
345 m_diam(i,j) *= f_sigma;
346 m_epsilon(i,j) *= f_eps;
347
348 // properties are symmetric
349 m_reducedMass(j,i) = m_reducedMass(i,j);
350 m_diam(j,i) = m_diam(i,j);
351 m_epsilon(j,i) = m_epsilon(i,j);
352 m_dipole(j,i) = m_dipole(i,j);
353 m_delta(j,i) = m_delta(i,j);
354 }
355 }
356}
357
358void GasTransport::invalidateCache()
359{
360 Transport::invalidateCache();
361 m_temp = Undef;
362 m_visc_ok = false;
363 m_spvisc_ok = false;
364 m_viscwt_ok = false;
365 m_bindiff_ok = false;
366}
367
368void GasTransport::setupCollisionIntegral()
369{
370 double tstar_min = 1.e8, tstar_max = 0.0;
371 for (size_t i = 0; i < m_nsp; i++) {
372 for (size_t j = i; j < m_nsp; j++) {
373 // The polynomial fits of collision integrals vs. T*
374 // will be done for the T* from tstar_min to tstar_max
375 tstar_min = std::min(tstar_min, Boltzmann * m_thermo->minTemp()/m_epsilon(i,j));
376 tstar_max = std::max(tstar_max, Boltzmann * m_thermo->maxTemp()/m_epsilon(i,j));
377 }
378 }
379 // Chemkin fits the entire T* range in the Monchick and Mason tables,
380 // so modify tstar_min and tstar_max if in Chemkin compatibility mode
381 if (m_mode == CK_Mode) {
382 tstar_min = 0.101;
383 tstar_max = 99.9;
384 }
385
386 // initialize the collision integral calculator for the desired T* range
387 MMCollisionInt integrals;
388 integrals.init(tstar_min, tstar_max);
389 fitCollisionIntegrals(integrals);
390 // make polynomial fits
391 fitProperties(integrals);
392}
393
394void GasTransport::getTransportData()
395{
396 for (size_t k = 0; k < m_thermo->nSpecies(); k++) {
397 shared_ptr<Species> s = m_thermo->species(k);
398 const GasTransportData* sptran =
399 dynamic_cast<GasTransportData*>(s->transport.get());
400 if (!sptran) {
401 throw CanteraError("GasTransport::getTransportData",
402 "Missing gas-phase transport data for species '{}'.", s->name);
403 }
404
405 if (sptran->geometry == "atom") {
406 m_crot[k] = 0.0;
407 } else if (sptran->geometry == "linear") {
408 m_crot[k] = 1.0;
409 } else if (sptran->geometry == "nonlinear") {
410 m_crot[k] = 1.5;
411 }
412
413 m_sigma[k] = sptran->diameter;
414 m_eps[k] = sptran->well_depth;
415 m_dipole(k,k) = sptran->dipole;
416 m_polar[k] = (sptran->dipole > 0);
417 m_alpha[k] = sptran->polarizability;
418 m_zrot[k] = sptran->rotational_relaxation;
419 if (s->input.hasKey("critical-parameters") &&
420 s->input["critical-parameters"].hasKey("acentric-factor"))
421 {
422 m_w_ac[k] = s->input["critical-parameters"]["acentric-factor"].asDouble();
423 } else {
424 m_w_ac[k] = sptran->acentric_factor;
425 }
426 m_disp[k] = sptran->dispersion_coefficient;
427 m_quad_polar[k] = sptran->quadrupole_polarizability;
428 }
429}
430
431void GasTransport::makePolarCorrections(size_t i, size_t j,
432 double& f_eps, double& f_sigma)
433{
434 // no correction if both are nonpolar, or both are polar
435 if (m_polar[i] == m_polar[j]) {
436 f_eps = 1.0;
437 f_sigma = 1.0;
438 return;
439 }
440
441 // corrections to the effective diameter and well depth
442 // if one is polar and one is non-polar
443 size_t kp = (m_polar[i] ? i : j); // the polar one
444 size_t knp = (i == kp ? j : i); // the nonpolar one
445 double d3np, d3p, alpha_star, mu_p_star, xi;
446 d3np = pow(m_sigma[knp],3);
447 d3p = pow(m_sigma[kp],3);
448 alpha_star = m_alpha[knp]/d3np;
449 mu_p_star = m_dipole(kp,kp)/sqrt(4 * Pi * epsilon_0 * d3p * m_eps[kp]);
450 xi = 1.0 + 0.25 * alpha_star * mu_p_star * mu_p_star *
451 sqrt(m_eps[kp]/m_eps[knp]);
452 f_sigma = pow(xi, -1.0/6.0);
453 f_eps = xi*xi;
454}
455
456void GasTransport::fitCollisionIntegrals(MMCollisionInt& integrals)
457{
458 // Chemkin fits to sixth order polynomials
459 int degree = (m_mode == CK_Mode ? 6 : COLL_INT_POLY_DEGREE);
460 vector<double> fitlist;
461 m_omega22_poly.clear();
462 m_astar_poly.clear();
463 m_bstar_poly.clear();
464 m_cstar_poly.clear();
465 for (size_t i = 0; i < m_nsp; i++) {
466 for (size_t j = i; j < m_nsp; j++) {
467 // Chemkin fits only delta* = 0
468 double dstar = (m_mode != CK_Mode) ? m_delta(i,j) : 0.0;
469
470 // if a fit has already been generated for delta* = m_delta(i,j),
471 // then use it. Otherwise, make a new fit, and add m_delta(i,j) to
472 // the list of delta* values for which fits have been done.
473
474 // 'find' returns a pointer to end() if not found
475 auto dptr = find(fitlist.begin(), fitlist.end(), dstar);
476 if (dptr == fitlist.end()) {
477 vector<double> ca(degree+1), cb(degree+1), cc(degree+1);
478 vector<double> co22(degree+1);
479 integrals.fit(degree, dstar, ca, cb, cc);
480 integrals.fit_omega22(degree, dstar, co22);
481 m_omega22_poly.push_back(co22);
482 m_astar_poly.push_back(ca);
483 m_bstar_poly.push_back(cb);
484 m_cstar_poly.push_back(cc);
485 m_poly[i][j] = static_cast<int>(m_astar_poly.size()) - 1;
486 m_star_poly_uses_actualT[i][j] = 0;
487 fitlist.push_back(dstar);
488 } else {
489 // delta* found in fitlist, so just point to this polynomial
490 m_poly[i][j] = static_cast<int>((dptr - fitlist.begin()));
491 m_star_poly_uses_actualT[i][j] = 0;
492 }
493 m_poly[j][i] = m_poly[i][j];
494 m_star_poly_uses_actualT[j][i] = m_star_poly_uses_actualT[i][j];
495 }
496 }
497}
498
499void GasTransport::fitProperties(MMCollisionInt& integrals)
500{
501 // number of points to use in generating fit data
502 const size_t np = 50;
503 int degree = (m_mode == CK_Mode ? 3 : 4);
504 double dt = (m_thermo->maxTemp() - m_thermo->minTemp())/(np-1);
505 vector<double> tlog(np), spvisc(np), spcond(np);
506 vector<double> w(np), w2(np);
507
508 m_visccoeffs.clear();
509 m_condcoeffs.clear();
510
511 // generate array of log(t) values
512 for (size_t n = 0; n < np; n++) {
513 double t = m_thermo->minTemp() + dt*n;
514 tlog[n] = log(t);
515 }
516
517 // vector of polynomial coefficients
518 vector<double> c(degree + 1), c2(degree + 1);
519
520 // fit the pure-species viscosity and thermal conductivity for each species
521 double visc, err, relerr,
522 mxerr = 0.0, mxrelerr = 0.0, mxerr_cond = 0.0, mxrelerr_cond = 0.0;
523 double T_save = m_thermo->temperature();
524 auto mw = m_thermo->molecularWeights();
525 for (size_t k = 0; k < m_nsp; k++) {
526 double tstar = Boltzmann * 298.0 / m_eps[k];
527 // Scaling factor for temperature dependence of z_rot. [Kee2003] Eq.
528 // 12.112 or [Kee2017] Eq. 11.115
529 double fz_298 = 1.0 + pow(Pi, 1.5) / sqrt(tstar) * (0.5 + 1.0 / tstar) +
530 (0.25 * Pi * Pi + 2) / tstar;
531
532 for (size_t n = 0; n < np; n++) {
533 double t = m_thermo->minTemp() + dt*n;
534 m_thermo->setTemperature(t);
535 vector<double> cp_R_all(m_thermo->nSpecies());
536 m_thermo->getCp_R_ref(cp_R_all);
537 double cp_R = cp_R_all[k];
538 tstar = Boltzmann * t / m_eps[k];
539 double sqrt_T = sqrt(t);
540 double om22 = integrals.omega22(tstar, m_delta(k,k));
541 double om11 = integrals.omega11(tstar, m_delta(k,k));
542
543 // self-diffusion coefficient, without polar corrections
544 double diffcoeff = 3.0/16.0 * sqrt(2.0 * Pi/m_reducedMass(k,k)) *
545 pow((Boltzmann * t), 1.5)/
546 (Pi * m_sigma[k] * m_sigma[k] * om11);
547
548 // viscosity
549 visc = 5.0/16.0 * sqrt(Pi * mw[k] * Boltzmann * t / Avogadro) /
550 (om22 * Pi * m_sigma[k]*m_sigma[k]);
551
552 // thermal conductivity
553 double f_int = mw[k]/(GasConstant * t) * diffcoeff/visc;
554 double cv_rot = m_crot[k];
555 double A_factor = 2.5 - f_int;
556 double fz_tstar = 1.0 + pow(Pi, 1.5) / sqrt(tstar) * (0.5 + 1.0 / tstar) +
557 (0.25 * Pi * Pi + 2) / tstar;
558 double B_factor = m_zrot[k] * fz_298 / fz_tstar + 2.0/Pi * (5.0/3.0 * cv_rot + f_int);
559 double c1 = 2.0/Pi * A_factor/B_factor;
560 double cv_int = cp_R - 2.5 - cv_rot;
561 double f_rot = f_int * (1.0 + c1);
562 double f_trans = 2.5 * (1.0 - c1 * cv_rot/1.5);
563 double cond = (visc/mw[k])*GasConstant*(f_trans * 1.5
564 + f_rot * cv_rot + f_int * cv_int);
565
566 if (m_mode == CK_Mode) {
567 spvisc[n] = log(visc);
568 spcond[n] = log(cond);
569 w[n] = -1.0;
570 w2[n] = -1.0;
571 } else {
572 // the viscosity should be proportional approximately to
573 // sqrt(T); therefore, visc/sqrt(T) should have only a weak
574 // temperature dependence. And since the mixture rule requires
575 // the square root of the pure-species viscosity, fit the square
576 // root of (visc/sqrt(T)) to avoid having to compute square
577 // roots in the mixture rule.
578 spvisc[n] = sqrt(visc/sqrt_T);
579
580 // the pure-species conductivity scales approximately with
581 // sqrt(T). Unlike the viscosity, there is no reason here to fit
582 // the square root, since a different mixture rule is used.
583 spcond[n] = cond/sqrt_T;
584 w[n] = 1.0/(spvisc[n]*spvisc[n]);
585 w2[n] = 1.0/(spcond[n]*spcond[n]);
586 }
587 }
588 polyfit(degree, tlog, spvisc, w, c);
589 polyfit(degree, tlog, spcond, w2, c2);
590
591 // evaluate max fit errors for viscosity
592 for (size_t n = 0; n < np; n++) {
593 double val, fit;
594 if (m_mode == CK_Mode) {
595 val = exp(spvisc[n]);
596 fit = exp(poly3(tlog[n], c));
597 } else {
598 double sqrt_T = exp(0.5*tlog[n]);
599 val = sqrt_T * pow(spvisc[n],2);
600 fit = sqrt_T * pow(poly4(tlog[n], c),2);
601 }
602 err = fit - val;
603 relerr = err/val;
604 mxerr = std::max(mxerr, fabs(err));
605 mxrelerr = std::max(mxrelerr, fabs(relerr));
606 }
607 m_fittingErrors["viscosity-max-abs-error"] = mxerr;
608 m_fittingErrors["viscosity-max-rel-error"] = mxrelerr;
609
610 // evaluate max fit errors for conductivity
611 for (size_t n = 0; n < np; n++) {
612 double val, fit;
613 if (m_mode == CK_Mode) {
614 val = exp(spcond[n]);
615 fit = exp(poly3(tlog[n], c2));
616 } else {
617 double sqrt_T = exp(0.5*tlog[n]);
618 val = sqrt_T * spcond[n];
619 fit = sqrt_T * poly4(tlog[n], c2);
620 }
621 err = fit - val;
622 relerr = err/val;
623 mxerr_cond = std::max(mxerr_cond, fabs(err));
624 mxrelerr_cond = std::max(mxrelerr_cond, fabs(relerr));
625 }
626 m_visccoeffs.push_back(c);
627 m_condcoeffs.push_back(c2);
628 m_fittingErrors["conductivity-max-abs-error"] = mxerr_cond;
629 m_fittingErrors["conductivity-max-rel-error"] = mxrelerr_cond;
630 }
631 m_thermo->setTemperature(T_save);
632 fitDiffCoeffs(integrals);
633}
634
635void GasTransport::fitDiffCoeffs(MMCollisionInt& integrals)
636{
637 // number of points to use in generating fit data
638 const size_t np = 50;
639 int degree = (m_mode == CK_Mode ? 3 : 4);
640 double dt = (m_thermo->maxTemp() - m_thermo->minTemp())/(np-1);
641 vector<double> tlog(np);
642 vector<double> w(np), w2(np);
643
644 // generate array of log(t) values
645 for (size_t n = 0; n < np; n++) {
646 double t = m_thermo->minTemp() + dt*n;
647 tlog[n] = log(t);
648 }
649
650 // vector of polynomial coefficients
651 vector<double> c(degree + 1), c2(degree + 1);
652 double err, relerr, mxerr = 0.0, mxrelerr = 0.0;
653
654 vector<double> diff(np + 1);
655 m_diffcoeffs.clear();
656 for (size_t k = 0; k < m_nsp; k++) {
657 for (size_t j = k; j < m_nsp; j++) {
658 for (size_t n = 0; n < np; n++) {
659 double t = m_thermo->minTemp() + dt*n;
660 double eps = m_epsilon(j,k);
661 double tstar = Boltzmann * t/eps;
662 double sigma = m_diam(j,k);
663 double om11 = integrals.omega11(tstar, m_delta(j,k));
664 double diffcoeff = 3.0/16.0 * sqrt(2.0 * Pi/m_reducedMass(k,j))
665 * pow(Boltzmann * t, 1.5) / (Pi * sigma * sigma * om11);
666
667 // 2nd order correction
668 // NOTE: THIS CORRECTION IS NOT APPLIED
669 double fkj, fjk;
670 getBinDiffCorrection(t, integrals, k, j, 1.0, 1.0, fkj, fjk);
671
672 if (m_mode == CK_Mode) {
673 diff[n] = log(diffcoeff);
674 w[n] = -1.0;
675 } else {
676 diff[n] = diffcoeff/pow(t, 1.5);
677 w[n] = 1.0/(diff[n]*diff[n]);
678 }
679 }
680 polyfit(degree, tlog, diff, w, c);
681
682 for (size_t n = 0; n < np; n++) {
683 double val, fit;
684 if (m_mode == CK_Mode) {
685 val = exp(diff[n]);
686 fit = exp(poly3(tlog[n], c));
687 } else {
688 double t = exp(tlog[n]);
689 double pre = pow(t, 1.5);
690 val = pre * diff[n];
691 fit = pre * poly4(tlog[n], c);
692 }
693 err = fit - val;
694 relerr = err/val;
695 mxerr = std::max(mxerr, fabs(err));
696 mxrelerr = std::max(mxrelerr, fabs(relerr));
697 }
698 m_diffcoeffs.push_back(c);
699 }
700 }
701
702 m_fittingErrors["diff-coeff-max-abs-error"] = mxerr;
703 m_fittingErrors["diff-coeff-max-rel-error"] = mxrelerr;
704}
705
706void GasTransport::getBinDiffCorrection(double t, MMCollisionInt& integrals,
707 size_t k, size_t j, double xk, double xj, double& fkj, double& fjk)
708{
709 double w1 = m_thermo->molecularWeight(k);
710 double w2 = m_thermo->molecularWeight(j);
711 double wsum = w1 + w2;
712 double wmwp = (w1 - w2)/wsum;
713 double sqw12 = sqrt(w1*w2);
714 double sig1 = m_sigma[k];
715 double sig2 = m_sigma[j];
716 double sig12 = 0.5*(m_sigma[k] + m_sigma[j]);
717 double sigratio = sig1*sig1/(sig2*sig2);
718 double sigratio2 = sig1*sig1/(sig12*sig12);
719 double sigratio3 = sig2*sig2/(sig12*sig12);
720 double tstar1 = Boltzmann * t / m_eps[k];
721 double tstar2 = Boltzmann * t / m_eps[j];
722 double tstar12 = Boltzmann * t / sqrt(m_eps[k] * m_eps[j]);
723 double om22_1 = integrals.omega22(tstar1, m_delta(k,k));
724 double om22_2 = integrals.omega22(tstar2, m_delta(j,j));
725 double om11_12 = integrals.omega11(tstar12, m_delta(k,j));
726 double astar_12 = integrals.astar(tstar12, m_delta(k,j));
727 double bstar_12 = integrals.bstar(tstar12, m_delta(k,j));
728 double cstar_12 = integrals.cstar(tstar12, m_delta(k,j));
729
730 double cnst = sigratio * sqrt(2.0*w2/wsum) * 2.0 * w1*w1/(wsum * w2);
731 double p1 = cnst * om22_1 / om11_12;
732
733 cnst = (1.0/sigratio) * sqrt(2.0*w1/wsum) * 2.0*w2*w2/(wsum*w1);
734 double p2 = cnst * om22_2 / om11_12;
735 double p12 = 15.0 * wmwp*wmwp + 8.0*w1*w2*astar_12/(wsum*wsum);
736
737 cnst = (2.0/(w2*wsum))*sqrt(2.0*w2/wsum)*sigratio2;
738 double q1 = cnst*((2.5 - 1.2*bstar_12)*w1*w1 + 3.0*w2*w2
739 + 1.6*w1*w2*astar_12);
740
741 cnst = (2.0/(w1*wsum))*sqrt(2.0*w1/wsum)*sigratio3;
742 double q2 = cnst*((2.5 - 1.2*bstar_12)*w2*w2 + 3.0*w1*w1
743 + 1.6*w1*w2*astar_12);
744 double q12 = wmwp*wmwp*15.0*(2.5 - 1.2*bstar_12)
745 + 4.0*w1*w2*astar_12*(11.0 - 2.4*bstar_12)/(wsum*wsum)
746 + 1.6*wsum*om22_1*om22_2/(om11_12*om11_12*sqw12)
747 * sigratio2 * sigratio3;
748
749 cnst = 6.0*cstar_12 - 5.0;
750 fkj = 1.0 + 0.1*cnst*cnst *
751 (p1*xk*xk + p2*xj*xj + p12*xk*xj)/
752 (q1*xk*xk + q2*xj*xj + q12*xk*xj);
753 fjk = 1.0 + 0.1*cnst*cnst *
754 (p2*xk*xk + p1*xj*xj + p12*xk*xj)/
755 (q2*xk*xk + q1*xj*xj + q12*xk*xj);
756}
757
758void GasTransport::getViscosityPolynomial(size_t i, span<double> coeffs) const
759{
760 checkSpeciesIndex(i);
761 size_t N = (m_mode == CK_Mode) ? 4 : 5;
762 checkArraySize("GasTransport::getViscosityPolynomial", coeffs.size(), N);
763 for (size_t k = 0; k < N; k++) {
764 coeffs[k] = m_visccoeffs[i][k];
765 }
766}
767
768void GasTransport::getConductivityPolynomial(size_t i, span<double> coeffs) const
769{
770 checkSpeciesIndex(i);
771 size_t N = (m_mode == CK_Mode) ? 4 : 5;
772 checkArraySize("GasTransport::getConductivityPolynomial", coeffs.size(), N);
773 for (size_t k = 0; k < N; k++) {
774 coeffs[k] = m_condcoeffs[i][k];
775 }
776}
777
778void GasTransport::getBinDiffusivityPolynomial(size_t i, size_t j, span<double> coeffs) const
779{
780 checkSpeciesIndex(i);
781 checkSpeciesIndex(j);
782 size_t N = (m_mode == CK_Mode) ? 4 : 5;
783 checkArraySize("GasTransport::getBinDiffusivityPolynomial", coeffs.size(), N);
784 size_t mi = (j >= i? i : j);
785 size_t mj = (j >= i? j : i);
786 size_t ic = 0;
787 for (size_t ii = 0; ii < mi; ii++) {
788 ic += m_nsp - ii;
789 }
790 ic += mj - mi;
791
792 for (size_t k = 0; k < N; k++) {
793 coeffs[k] = m_diffcoeffs[ic][k];
794 }
795}
796
797void GasTransport::getCollisionIntegralPolynomial(size_t i, size_t j,
798 span<double> astar_coeffs,
799 span<double> bstar_coeffs,
800 span<double> cstar_coeffs) const
801{
802 checkSpeciesIndex(i);
803 checkSpeciesIndex(j);
804 size_t N = (m_mode == CK_Mode) ? 7 : COLL_INT_POLY_DEGREE + 1;
805 checkArraySize("GasTransport::getCollisionIntegralPolynomial[astar]",
806 astar_coeffs.size(), N);
807 checkArraySize("GasTransport::getCollisionIntegralPolynomial[bstar]",
808 bstar_coeffs.size(), N);
809 checkArraySize("GasTransport::getCollisionIntegralPolynomial[cstar]",
810 cstar_coeffs.size(), N);
811 for (size_t k = 0; k < N; k++) {
812 astar_coeffs[k] = m_astar_poly[m_poly[i][j]][k];
813 bstar_coeffs[k] = m_bstar_poly[m_poly[i][j]][k];
814 cstar_coeffs[k] = m_cstar_poly[m_poly[i][j]][k];
815 }
816}
817
818void GasTransport::setViscosityPolynomial(size_t i, span<double> coeffs)
819{
820 checkSpeciesIndex(i);
821 size_t N = (m_mode == CK_Mode) ? 4 : 5;
822 checkArraySize("GasTransport::setViscosityPolynomial", coeffs.size(), N);
823 for (size_t k = 0; k < N; k++) {
824 m_visccoeffs[i][k] = coeffs[k];
825 }
826 invalidateCache();
827}
828
829void GasTransport::setConductivityPolynomial(size_t i, span<double> coeffs)
830{
831 checkSpeciesIndex(i);
832 size_t N = (m_mode == CK_Mode) ? 4 : 5;
833 checkArraySize("GasTransport::setConductivityPolynomial", coeffs.size(), N);
834 for (size_t k = 0; k < N; k++) {
835 m_condcoeffs[i][k] = coeffs[k];
836 }
837 invalidateCache();
838}
839
840void GasTransport::setBinDiffusivityPolynomial(size_t i, size_t j, span<double> coeffs)
841{
842 checkSpeciesIndex(i);
843 checkSpeciesIndex(j);
844 size_t N = (m_mode == CK_Mode) ? 4 : 5;
845 checkArraySize("GasTransport::setBinDiffusivityPolynomial", coeffs.size(), N);
846 size_t mi = (j >= i? i : j);
847 size_t mj = (j >= i? j : i);
848 size_t ic = 0;
849 for (size_t ii = 0; ii < mi; ii++) {
850 ic += m_nsp - ii;
851 }
852 ic += mj - mi;
853
854 for (size_t k = 0; k < N; k++) {
855 m_diffcoeffs[ic][k] = coeffs[k];
856 }
857 invalidateCache();
858}
859
860void GasTransport::setCollisionIntegralPolynomial(size_t i, size_t j,
861 span<double> astar_coeffs, span<double> bstar_coeffs, span<double> cstar_coeffs,
862 bool actualT)
863{
864 checkSpeciesIndex(i);
865 checkSpeciesIndex(j);
866 size_t N = (m_mode == CK_Mode) ? 7 : COLL_INT_POLY_DEGREE + 1;
867 checkArraySize("GasTransport::setCollisionIntegralPolynomial[astar]",
868 astar_coeffs.size(), N);
869 checkArraySize("GasTransport::setCollisionIntegralPolynomial[bstar]",
870 bstar_coeffs.size(), N);
871 checkArraySize("GasTransport::setCollisionIntegralPolynomial[cstar]",
872 cstar_coeffs.size(), N);
873 vector<double> ca(N), cb(N), cc(N);
874
875 for (size_t k = 0; k < N; k++) {
876 ca[k] = astar_coeffs[k];
877 cb[k] = bstar_coeffs[k];
878 cc[k] = cstar_coeffs[k];
879 }
880
881 m_astar_poly.push_back(ca);
882 m_bstar_poly.push_back(cb);
883 m_cstar_poly.push_back(cc);
884 m_poly[j][i] = m_poly[i][j] = static_cast<int>(m_astar_poly.size()) - 1;
885 m_star_poly_uses_actualT[j][i] = m_star_poly_uses_actualT[i][j] = (actualT) ? 1 : 0;
886 invalidateCache();
887}
888
889}
#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³]. 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 the square of the elementary charge. [m⁵] 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.
This file contains definitions for utility functions and text for modules, inputfiles and logging,...
auto dot5(const X &x, const Y &y)
Templated Inner product of two vectors of length 5.
Definition utilities.h:63
auto poly4(D x, const C &c)
Evaluates a polynomial of order 4.
Definition utilities.h:179
auto poly3(D x, const C &c)
Templated evaluation of a polynomial of order 3.
Definition utilities.h:195
auto dot4(const X &x, const Y &y)
Templated Inner product of two vectors of length 4.
Definition utilities.h:42
double polyfit(size_t deg, span< const double > x, span< const double > y, span< const double > w, span< double > p)
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:87
const double Avogadro
Avogadro's Number [number/kmol].
Definition ct_defs.h:84
const double epsilon_0
Permittivity of free space [F/m].
Definition ct_defs.h:140
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition ct_defs.h:123
const double Pi
Pi.
Definition ct_defs.h:71
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:167
void multiply(const DenseMatrix &A, span< const double > b, span< double > prod)
Multiply A*b and return the result in prod. Uses BLAS routine DGEMV.
void checkArraySize(const char *procedure, size_t available, size_t required)
Wrapper for throwing ArraySizeError.
Contains declarations for string manipulation functions within Cantera.
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...