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.data(), m_spwork.data());
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,
double*
const d)
157 throw CanteraError(
"GasTransport::getBinaryDiffCoeffs",
"ld is too small");
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);
167void GasTransport::getMixDiffCoeffs(
double*
const d)
177 double mmw = m_thermo->meanMolecularWeight();
178 double p = m_thermo->pressure();
180 d[0] = m_bdiff(0,0) / p;
182 for (
size_t k = 0; k < m_nsp; k++) {
184 for (
size_t j = 0; j < m_nsp; j++) {
186 sum2 += m_molefracs[j] / m_bdiff(j,k);
190 d[k] = m_bdiff(k,k) / p;
192 d[k] = (mmw - m_molefracs[k] * m_mw[k])/(p * mmw * sum2);
198void GasTransport::getMixDiffCoeffsMole(
double*
const d)
208 double p = m_thermo->pressure();
210 d[0] = m_bdiff(0,0) / p;
212 for (
size_t k = 0; k < m_nsp; k++) {
214 for (
size_t j = 0; j < m_nsp; j++) {
216 sum2 += m_molefracs[j] / m_bdiff(j,k);
220 d[k] = m_bdiff(k,k) / p;
222 d[k] = (1 - m_molefracs[k]) / (p * sum2);
228void GasTransport::getMixDiffCoeffsMass(
double*
const d)
238 double mmw = m_thermo->meanMolecularWeight();
239 double p = m_thermo->pressure();
242 d[0] = m_bdiff(0,0) / p;
244 for (
size_t k=0; k<m_nsp; k++) {
247 for (
size_t i=0; i<m_nsp; i++) {
251 sum1 += m_molefracs[i] / m_bdiff(k,i);
252 sum2 += m_molefracs[i] * m_mw[i] / m_bdiff(k,i);
255 sum2 *= p * m_molefracs[k] / (mmw - m_mw[k]*m_molefracs[k]);
256 d[k] = 1.0 / (sum1 + sum2);
268 setupCollisionParameters();
269 setupCollisionIntegral();
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);
279 m_mw = m_thermo->molecularWeights();
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]);
292void GasTransport::setupCollisionParameters()
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);
307 m_w_ac.resize(m_nsp);
308 m_disp.resize(m_nsp, 0.0);
309 m_quad_polar.resize(m_nsp, 0.0);
311 const vector<double>& mw = m_thermo->molecularWeights();
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);
319 double f_eps, f_sigma;
321 for (
size_t i = 0; i < m_nsp; i++) {
322 for (
size_t j = i; j < m_nsp; j++) {
324 m_reducedMass(i,j) = mw[i] * mw[j] / (
Avogadro * (mw[i] + mw[j]));
327 m_diam(i,j) = 0.5*(m_sigma[i] + m_sigma[j]);
330 m_epsilon(i,j) = sqrt(m_eps[i]*m_eps[j]);
333 m_dipole(i,j) = sqrt(m_dipole(i,i)*m_dipole(j,j));
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;
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);
353void GasTransport::invalidateCache()
355 Transport::invalidateCache();
360 m_bindiff_ok =
false;
363void GasTransport::setupCollisionIntegral()
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++) {
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));
376 if (m_mode == CK_Mode) {
383 integrals.
init(tstar_min, tstar_max);
384 fitCollisionIntegrals(integrals);
386 fitProperties(integrals);
389void GasTransport::getTransportData()
391 for (
size_t k = 0; k < m_thermo->nSpecies(); k++) {
392 shared_ptr<Species> s = m_thermo->species(k);
397 "Missing gas-phase transport data for species '{}'.", s->name);
402 }
else if (sptran->
geometry ==
"linear") {
404 }
else if (sptran->
geometry ==
"nonlinear") {
410 m_dipole(k,k) = sptran->
dipole;
411 m_polar[k] = (sptran->
dipole > 0);
414 if (s->input.hasKey(
"critical-parameters") &&
415 s->input[
"critical-parameters"].hasKey(
"acentric-factor"))
417 m_w_ac[k] = s->input[
"critical-parameters"][
"acentric-factor"].asDouble();
426void GasTransport::makePolarCorrections(
size_t i,
size_t j,
427 double& f_eps,
double& f_sigma)
430 if (m_polar[i] == m_polar[j]) {
438 size_t kp = (m_polar[i] ? i : j);
439 size_t knp = (i == kp ? j : i);
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);
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++) {
463 double dstar = (m_mode != CK_Mode) ? m_delta(i,j) : 0.0;
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);
485 m_poly[i][j] =
static_cast<int>((dptr - fitlist.begin()));
486 m_star_poly_uses_actualT[i][j] = 0;
488 m_poly[j][i] = m_poly[i][j];
489 m_star_poly_uses_actualT[j][i] = m_star_poly_uses_actualT[i][j];
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);
503 m_visccoeffs.clear();
504 m_condcoeffs.clear();
507 for (
size_t n = 0; n < np; n++) {
508 double t = m_thermo->minTemp() + dt*n;
513 vector<double> c(degree + 1), c2(degree + 1);
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];
524 double fz_298 = 1.0 + pow(
Pi, 1.5) / sqrt(tstar) * (0.5 + 1.0 / tstar) +
525 (0.25 *
Pi *
Pi + 2) / tstar;
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];
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));
539 double diffcoeff = 3.0/16.0 * sqrt(2.0 *
Pi/m_reducedMass(k,k)) *
541 (
Pi * m_sigma[k] * m_sigma[k] * om11);
545 (om22 *
Pi * m_sigma[k]*m_sigma[k]);
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);
561 if (m_mode == CK_Mode) {
562 spvisc[n] = log(visc);
563 spcond[n] = log(cond);
573 spvisc[n] = sqrt(visc/sqrt_T);
578 spcond[n] = cond/sqrt_T;
579 w[n] = 1.0/(spvisc[n]*spvisc[n]);
580 w2[n] = 1.0/(spcond[n]*spcond[n]);
583 polyfit(np, degree, tlog.data(), spvisc.data(), w.data(), c.data());
584 polyfit(np, degree, tlog.data(), spcond.data(), w2.data(), c2.data());
587 for (
size_t n = 0; n < np; n++) {
589 if (m_mode == CK_Mode) {
590 val = exp(spvisc[n]);
591 fit = exp(
poly3(tlog[n], c.data()));
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);
599 mxerr = std::max(mxerr, fabs(err));
600 mxrelerr = std::max(mxrelerr, fabs(relerr));
602 m_fittingErrors[
"viscosity-max-abs-error"] = mxerr;
603 m_fittingErrors[
"viscosity-max-rel-error"] = mxrelerr;
606 for (
size_t n = 0; n < np; n++) {
608 if (m_mode == CK_Mode) {
609 val = exp(spcond[n]);
610 fit = exp(
poly3(tlog[n], c2.data()));
612 double sqrt_T = exp(0.5*tlog[n]);
613 val = sqrt_T * spcond[n];
614 fit = sqrt_T *
poly4(tlog[n], c2.data());
618 mxerr_cond = std::max(mxerr_cond, fabs(err));
619 mxrelerr_cond = std::max(mxrelerr_cond, fabs(relerr));
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;
626 m_thermo->setTemperature(T_save);
627 fitDiffCoeffs(integrals);
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);
640 for (
size_t n = 0; n < np; n++) {
641 double t = m_thermo->minTemp() + dt*n;
646 vector<double> c(degree + 1), c2(degree + 1);
647 double err, relerr, mxerr = 0.0, mxrelerr = 0.0;
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);
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);
665 getBinDiffCorrection(t, integrals, k, j, 1.0, 1.0, fkj, fjk);
667 if (m_mode == CK_Mode) {
668 diff[n] = log(diffcoeff);
671 diff[n] = diffcoeff/pow(t, 1.5);
672 w[n] = 1.0/(diff[n]*diff[n]);
675 polyfit(np, degree, tlog.data(), diff.data(), w.data(), c.data());
677 for (
size_t n = 0; n < np; n++) {
679 if (m_mode == CK_Mode) {
681 fit = exp(
poly3(tlog[n], c.data()));
683 double t = exp(tlog[n]);
684 double pre = pow(t, 1.5);
686 fit = pre *
poly4(tlog[n], c.data());
690 mxerr = std::max(mxerr, fabs(err));
691 mxrelerr = std::max(mxrelerr, fabs(relerr));
693 m_diffcoeffs.push_back(c);
697 m_fittingErrors[
"diff-coeff-max-abs-error"] = mxerr;
698 m_fittingErrors[
"diff-coeff-max-rel-error"] = mxrelerr;
702 size_t k,
size_t j,
double xk,
double xj,
double& fkj,
double& fjk)
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));
725 double cnst = sigratio * sqrt(2.0*w2/wsum) * 2.0 * w1*w1/(wsum * w2);
726 double p1 = cnst * om22_1 / om11_12;
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);
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);
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;
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);
753void GasTransport::getViscosityPolynomial(
size_t i,
double* coeffs)
const
755 checkSpeciesIndex(i);
756 for (
int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
757 coeffs[k] = m_visccoeffs[i][k];
761void GasTransport::getConductivityPolynomial(
size_t i,
double* coeffs)
const
763 checkSpeciesIndex(i);
764 for (
int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
765 coeffs[k] = m_condcoeffs[i][k];
769void GasTransport::getBinDiffusivityPolynomial(
size_t i,
size_t j,
double* coeffs)
const
771 checkSpeciesIndex(i);
772 checkSpeciesIndex(j);
773 size_t mi = (j >= i? i : j);
774 size_t mj = (j >= i? j : i);
776 for (
size_t ii = 0; ii < mi; ii++) {
781 for (
int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
782 coeffs[k] = m_diffcoeffs[ic][k];
786void GasTransport::getCollisionIntegralPolynomial(
size_t i,
size_t j,
787 double* astar_coeffs,
788 double* bstar_coeffs,
789 double* cstar_coeffs)
const
791 checkSpeciesIndex(i);
792 checkSpeciesIndex(j);
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];
800void GasTransport::setViscosityPolynomial(
size_t i,
double* coeffs)
802 checkSpeciesIndex(i);
803 for (
int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
804 m_visccoeffs[i][k] = coeffs[k];
809void GasTransport::setConductivityPolynomial(
size_t i,
double* coeffs)
811 checkSpeciesIndex(i);
812 for (
int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
813 m_condcoeffs[i][k] = coeffs[k];
818void GasTransport::setBinDiffusivityPolynomial(
size_t i,
size_t j,
double* coeffs)
820 checkSpeciesIndex(i);
821 checkSpeciesIndex(j);
822 size_t mi = (j >= i? i : j);
823 size_t mj = (j >= i? j : i);
825 for (
size_t ii = 0; ii < mi; ii++) {
830 for (
int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
831 m_diffcoeffs[ic][k] = coeffs[k];
836void GasTransport::setCollisionIntegralPolynomial(
size_t i,
size_t j,
837 double* astar_coeffs,
838 double* bstar_coeffs,
839 double* cstar_coeffs,
bool actualT)
841 checkSpeciesIndex(i);
842 checkSpeciesIndex(j);
844 vector<double> ca(degree+1), cb(degree+1), cc(degree+1);
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];
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;
#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.
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.
R poly4(D x, R *c)
Evaluates a polynomial of order 4.
R poly3(D x, R *c)
Templated evaluation of a polynomial of order 3.
double dot4(const V &x, const V &y)
Templated Inner product of two vectors of length 4.
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.
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, 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...