21#define COLL_INT_POLY_DEGREE 8
23GasTransport::GasTransport() :
28void GasTransport::update_T()
30 if (m_thermo->nSpecies() != m_nsp) {
32 init(m_thermo, m_mode);
35 double T = m_thermo->temperature();
43 m_sqrt_t = sqrt(m_temp);
44 m_t14 = sqrt(m_sqrt_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;
60double GasTransport::viscosity()
75 multiply(m_phi, m_molefracs, m_spwork);
77 for (
size_t k = 0; k < m_nsp; k++) {
78 vismix += m_molefracs[k] * m_visc[k]/m_spwork[k];
84void GasTransport::updateViscosity_T()
87 updateSpeciesViscosities();
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];
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);
105void GasTransport::updateSpeciesViscosities()
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]);
114 for (
size_t k = 0; k < m_nsp; k++) {
116 m_sqvisc[k] = m_t14 *
dot5(m_polytempvec, m_visccoeffs[k]);
117 m_visc[k] = (m_sqvisc[k] * m_sqvisc[k]);
123void GasTransport::updateDiff_T()
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);
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,
141 m_bdiff(j,i) = m_bdiff(i,j);
149void GasTransport::getBinaryDiffCoeffs(
const size_t ld, span<double> d)
157 throw CanteraError(
"GasTransport::getBinaryDiffCoeffs",
"ld is too small");
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);
168void GasTransport::getMixDiffCoeffs(span<double> d)
170 checkArraySize(
"GasTransport::getMixDiffCoeffs", d.size(), m_nsp);
179 double mmw = m_thermo->meanMolecularWeight();
180 double p = m_thermo->pressure();
182 d[0] = m_bdiff(0,0) / p;
184 for (
size_t k = 0; k < m_nsp; k++) {
186 for (
size_t j = 0; j < m_nsp; j++) {
188 sum2 += m_molefracs[j] / m_bdiff(j,k);
192 d[k] = m_bdiff(k,k) / p;
194 d[k] = (mmw - m_molefracs[k] * m_mw[k])/(p * mmw * sum2);
200void GasTransport::getMixDiffCoeffsMole(span<double> d)
202 checkArraySize(
"GasTransport::getMixDiffCoeffsMole", d.size(), m_nsp);
211 double p = m_thermo->pressure();
213 d[0] = m_bdiff(0,0) / p;
215 for (
size_t k = 0; k < m_nsp; k++) {
217 for (
size_t j = 0; j < m_nsp; j++) {
219 sum2 += m_molefracs[j] / m_bdiff(j,k);
223 d[k] = m_bdiff(k,k) / p;
225 d[k] = (1 - m_molefracs[k]) / (p * sum2);
231void GasTransport::getMixDiffCoeffsMass(span<double> d)
233 checkArraySize(
"GasTransport::getMixDiffCoeffsMass", d.size(), m_nsp);
242 double mmw = m_thermo->meanMolecularWeight();
243 double p = m_thermo->pressure();
246 d[0] = m_bdiff(0,0) / p;
248 for (
size_t k=0; k<m_nsp; k++) {
251 for (
size_t i=0; i<m_nsp; i++) {
255 sum1 += m_molefracs[i] / m_bdiff(k,i);
256 sum2 += m_molefracs[i] * m_mw[i] / m_bdiff(k,i);
259 sum2 *= p * m_molefracs[k] / (mmw - m_mw[k]*m_molefracs[k]);
260 d[k] = 1.0 / (sum1 + sum2);
265void GasTransport::init(shared_ptr<ThermoPhase> thermo,
int mode)
268 m_nsp = m_thermo->nSpecies();
272 setupCollisionParameters();
273 setupCollisionIntegral();
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);
284 m_thermo->getMolecularWeights(m_mw);
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]);
297void GasTransport::setupCollisionParameters()
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);
312 m_w_ac.resize(m_nsp);
313 m_disp.resize(m_nsp, 0.0);
314 m_quad_polar.resize(m_nsp, 0.0);
316 auto mw = m_thermo->molecularWeights();
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);
324 double f_eps, f_sigma;
326 for (
size_t i = 0; i < m_nsp; i++) {
327 for (
size_t j = i; j < m_nsp; j++) {
329 m_reducedMass(i,j) = mw[i] * mw[j] / (
Avogadro * (mw[i] + mw[j]));
332 m_diam(i,j) = 0.5*(m_sigma[i] + m_sigma[j]);
335 m_epsilon(i,j) = sqrt(m_eps[i]*m_eps[j]);
338 m_dipole(i,j) = sqrt(m_dipole(i,i)*m_dipole(j,j));
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;
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);
358void GasTransport::invalidateCache()
360 Transport::invalidateCache();
365 m_bindiff_ok =
false;
368void GasTransport::setupCollisionIntegral()
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++) {
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));
381 if (m_mode == CK_Mode) {
388 integrals.
init(tstar_min, tstar_max);
389 fitCollisionIntegrals(integrals);
391 fitProperties(integrals);
394void GasTransport::getTransportData()
396 for (
size_t k = 0; k < m_thermo->nSpecies(); k++) {
397 shared_ptr<Species> s = m_thermo->species(k);
402 "Missing gas-phase transport data for species '{}'.", s->name);
407 }
else if (sptran->
geometry ==
"linear") {
409 }
else if (sptran->
geometry ==
"nonlinear") {
415 m_dipole(k,k) = sptran->
dipole;
416 m_polar[k] = (sptran->
dipole > 0);
419 if (s->input.hasKey(
"critical-parameters") &&
420 s->input[
"critical-parameters"].hasKey(
"acentric-factor"))
422 m_w_ac[k] = s->input[
"critical-parameters"][
"acentric-factor"].asDouble();
431void GasTransport::makePolarCorrections(
size_t i,
size_t j,
432 double& f_eps,
double& f_sigma)
435 if (m_polar[i] == m_polar[j]) {
443 size_t kp = (m_polar[i] ? i : j);
444 size_t knp = (i == kp ? j : i);
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);
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++) {
468 double dstar = (m_mode != CK_Mode) ? m_delta(i,j) : 0.0;
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);
490 m_poly[i][j] =
static_cast<int>((dptr - fitlist.begin()));
491 m_star_poly_uses_actualT[i][j] = 0;
493 m_poly[j][i] = m_poly[i][j];
494 m_star_poly_uses_actualT[j][i] = m_star_poly_uses_actualT[i][j];
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);
508 m_visccoeffs.clear();
509 m_condcoeffs.clear();
512 for (
size_t n = 0; n < np; n++) {
513 double t = m_thermo->minTemp() + dt*n;
518 vector<double> c(degree + 1), c2(degree + 1);
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];
529 double fz_298 = 1.0 + pow(
Pi, 1.5) / sqrt(tstar) * (0.5 + 1.0 / tstar) +
530 (0.25 *
Pi *
Pi + 2) / tstar;
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];
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));
544 double diffcoeff = 3.0/16.0 * sqrt(2.0 *
Pi/m_reducedMass(k,k)) *
546 (
Pi * m_sigma[k] * m_sigma[k] * om11);
550 (om22 *
Pi * m_sigma[k]*m_sigma[k]);
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);
566 if (m_mode == CK_Mode) {
567 spvisc[n] = log(visc);
568 spcond[n] = log(cond);
578 spvisc[n] = sqrt(visc/sqrt_T);
583 spcond[n] = cond/sqrt_T;
584 w[n] = 1.0/(spvisc[n]*spvisc[n]);
585 w2[n] = 1.0/(spcond[n]*spcond[n]);
588 polyfit(degree, tlog, spvisc, w, c);
589 polyfit(degree, tlog, spcond, w2, c2);
592 for (
size_t n = 0; n < np; n++) {
594 if (m_mode == CK_Mode) {
595 val = exp(spvisc[n]);
596 fit = exp(
poly3(tlog[n], c));
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);
604 mxerr = std::max(mxerr, fabs(err));
605 mxrelerr = std::max(mxrelerr, fabs(relerr));
607 m_fittingErrors[
"viscosity-max-abs-error"] = mxerr;
608 m_fittingErrors[
"viscosity-max-rel-error"] = mxrelerr;
611 for (
size_t n = 0; n < np; n++) {
613 if (m_mode == CK_Mode) {
614 val = exp(spcond[n]);
615 fit = exp(
poly3(tlog[n], c2));
617 double sqrt_T = exp(0.5*tlog[n]);
618 val = sqrt_T * spcond[n];
619 fit = sqrt_T *
poly4(tlog[n], c2);
623 mxerr_cond = std::max(mxerr_cond, fabs(err));
624 mxrelerr_cond = std::max(mxrelerr_cond, fabs(relerr));
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;
631 m_thermo->setTemperature(T_save);
632 fitDiffCoeffs(integrals);
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);
645 for (
size_t n = 0; n < np; n++) {
646 double t = m_thermo->minTemp() + dt*n;
651 vector<double> c(degree + 1), c2(degree + 1);
652 double err, relerr, mxerr = 0.0, mxrelerr = 0.0;
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);
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);
670 getBinDiffCorrection(t, integrals, k, j, 1.0, 1.0, fkj, fjk);
672 if (m_mode == CK_Mode) {
673 diff[n] = log(diffcoeff);
676 diff[n] = diffcoeff/pow(t, 1.5);
677 w[n] = 1.0/(diff[n]*diff[n]);
680 polyfit(degree, tlog, diff, w, c);
682 for (
size_t n = 0; n < np; n++) {
684 if (m_mode == CK_Mode) {
686 fit = exp(
poly3(tlog[n], c));
688 double t = exp(tlog[n]);
689 double pre = pow(t, 1.5);
691 fit = pre *
poly4(tlog[n], c);
695 mxerr = std::max(mxerr, fabs(err));
696 mxrelerr = std::max(mxrelerr, fabs(relerr));
698 m_diffcoeffs.push_back(c);
702 m_fittingErrors[
"diff-coeff-max-abs-error"] = mxerr;
703 m_fittingErrors[
"diff-coeff-max-rel-error"] = mxrelerr;
707 size_t k,
size_t j,
double xk,
double xj,
double& fkj,
double& fjk)
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));
730 double cnst = sigratio * sqrt(2.0*w2/wsum) * 2.0 * w1*w1/(wsum * w2);
731 double p1 = cnst * om22_1 / om11_12;
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);
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);
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;
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);
758void GasTransport::getViscosityPolynomial(
size_t i, span<double> coeffs)
const
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];
768void GasTransport::getConductivityPolynomial(
size_t i, span<double> coeffs)
const
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];
778void GasTransport::getBinDiffusivityPolynomial(
size_t i,
size_t j, span<double> coeffs)
const
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);
787 for (
size_t ii = 0; ii < mi; ii++) {
792 for (
size_t k = 0; k < N; k++) {
793 coeffs[k] = m_diffcoeffs[ic][k];
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
802 checkSpeciesIndex(i);
803 checkSpeciesIndex(j);
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];
818void GasTransport::setViscosityPolynomial(
size_t i, span<double> coeffs)
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];
829void GasTransport::setConductivityPolynomial(
size_t i, span<double> coeffs)
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];
840void GasTransport::setBinDiffusivityPolynomial(
size_t i,
size_t j, span<double> coeffs)
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);
849 for (
size_t ii = 0; ii < mi; ii++) {
854 for (
size_t k = 0; k < N; k++) {
855 m_diffcoeffs[ic][k] = coeffs[k];
860void GasTransport::setCollisionIntegralPolynomial(
size_t i,
size_t j,
861 span<double> astar_coeffs, span<double> bstar_coeffs, span<double> cstar_coeffs,
864 checkSpeciesIndex(i);
865 checkSpeciesIndex(j);
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);
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];
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;
#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.
auto poly4(D x, const C &c)
Evaluates a polynomial of order 4.
auto poly3(D x, const C &c)
Templated evaluation of a polynomial of order 3.
auto dot4(const X &x, const Y &y)
Templated Inner product of two vectors of length 4.
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.
const double Boltzmann
Boltzmann constant [J/K].
const double Avogadro
Avogadro's Number [number/kmol].
const double epsilon_0
Permittivity of free space [F/m].
const double GasConstant
Universal Gas Constant [J/kmol/K].
Namespace for the Cantera kernel.
const double Undef
Fairly random number to be used to initialize variables against to see if they are subsequently defin...
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...