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, m_log_level != 0 ? m_log_level : -7);
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);
261void GasTransport::init(
ThermoPhase* thermo,
int mode,
int log_level)
266 if (log_level == -7) {
270 "is deprecated and will be removed after Cantera 3.1.");
272 m_log_level = log_level;
275 setupCollisionParameters();
276 setupCollisionIntegral();
278 m_molefracs.resize(m_nsp);
279 m_spwork.resize(m_nsp);
280 m_visc.resize(m_nsp);
281 m_sqvisc.resize(m_nsp);
282 m_phi.resize(m_nsp, m_nsp, 0.0);
283 m_bdiff.resize(m_nsp, m_nsp);
286 m_mw = m_thermo->molecularWeights();
288 m_wratjk.resize(m_nsp, m_nsp, 0.0);
289 m_wratkj1.resize(m_nsp, m_nsp, 0.0);
290 for (
size_t j = 0; j < m_nsp; j++) {
291 for (
size_t k = j; k < m_nsp; k++) {
292 m_wratjk(j,k) = sqrt(m_mw[j]/m_mw[k]);
293 m_wratjk(k,j) = sqrt(m_wratjk(j,k));
294 m_wratkj1(j,k) = sqrt(1.0 + m_mw[k]/m_mw[j]);
299void GasTransport::setupCollisionParameters()
301 m_epsilon.resize(m_nsp, m_nsp, 0.0);
302 m_delta.resize(m_nsp, m_nsp, 0.0);
303 m_reducedMass.resize(m_nsp, m_nsp, 0.0);
304 m_dipole.resize(m_nsp, m_nsp, 0.0);
305 m_diam.resize(m_nsp, m_nsp, 0.0);
306 m_crot.resize(m_nsp);
307 m_zrot.resize(m_nsp);
308 m_polar.resize(m_nsp,
false);
309 m_alpha.resize(m_nsp, 0.0);
310 m_poly.resize(m_nsp);
311 m_star_poly_uses_actualT.resize(m_nsp);
312 m_sigma.resize(m_nsp);
314 m_w_ac.resize(m_nsp);
315 m_disp.resize(m_nsp, 0.0);
316 m_quad_polar.resize(m_nsp, 0.0);
318 const vector<double>& mw = m_thermo->molecularWeights();
321 for (
size_t i = 0; i < m_nsp; i++) {
322 m_poly[i].resize(m_nsp);
323 m_star_poly_uses_actualT[i].resize(m_nsp);
326 double f_eps, f_sigma;
328 for (
size_t i = 0; i < m_nsp; i++) {
329 for (
size_t j = i; j < m_nsp; j++) {
331 m_reducedMass(i,j) = mw[i] * mw[j] / (
Avogadro * (mw[i] + mw[j]));
334 m_diam(i,j) = 0.5*(m_sigma[i] + m_sigma[j]);
337 m_epsilon(i,j) = sqrt(m_eps[i]*m_eps[j]);
340 m_dipole(i,j) = sqrt(m_dipole(i,i)*m_dipole(j,j));
343 double d = m_diam(i,j);
344 m_delta(i,j) = 0.5 * m_dipole(i,j)*m_dipole(i,j)
345 / (4 *
Pi *
epsilon_0 * m_epsilon(i,j) * d * d * d);
346 makePolarCorrections(i, j, f_eps, f_sigma);
347 m_diam(i,j) *= f_sigma;
348 m_epsilon(i,j) *= f_eps;
351 m_reducedMass(j,i) = m_reducedMass(i,j);
352 m_diam(j,i) = m_diam(i,j);
353 m_epsilon(j,i) = m_epsilon(i,j);
354 m_dipole(j,i) = m_dipole(i,j);
355 m_delta(j,i) = m_delta(i,j);
360void GasTransport::invalidateCache()
362 Transport::invalidateCache();
367 m_bindiff_ok =
false;
370void GasTransport::setupCollisionIntegral()
372 double tstar_min = 1.e8, tstar_max = 0.0;
373 for (
size_t i = 0; i < m_nsp; i++) {
374 for (
size_t j = i; j < m_nsp; j++) {
377 tstar_min = std::min(tstar_min,
Boltzmann * m_thermo->minTemp()/m_epsilon(i,j));
378 tstar_max = std::max(tstar_max,
Boltzmann * m_thermo->maxTemp()/m_epsilon(i,j));
383 if (m_mode == CK_Mode) {
389 debuglog(
"*** collision_integrals ***\n", m_log_level);
391 integrals.
init(tstar_min, tstar_max, m_log_level != 0 ? m_log_level : -7);
392 fitCollisionIntegrals(integrals);
393 debuglog(
"*** end of collision_integrals ***\n", m_log_level);
395 debuglog(
"*** property fits ***\n", m_log_level);
396 fitProperties(integrals);
397 debuglog(
"*** end of property fits ***\n", m_log_level);
400void GasTransport::getTransportData()
402 for (
size_t k = 0; k < m_thermo->nSpecies(); k++) {
403 shared_ptr<Species> s = m_thermo->species(k);
408 "Missing gas-phase transport data for species '{}'.", s->name);
413 }
else if (sptran->
geometry ==
"linear") {
415 }
else if (sptran->
geometry ==
"nonlinear") {
421 m_dipole(k,k) = sptran->
dipole;
422 m_polar[k] = (sptran->
dipole > 0);
425 if (s->input.hasKey(
"critical-parameters") &&
426 s->input[
"critical-parameters"].hasKey(
"acentric-factor"))
428 m_w_ac[k] = s->input[
"critical-parameters"][
"acentric-factor"].asDouble();
437void GasTransport::makePolarCorrections(
size_t i,
size_t j,
438 double& f_eps,
double& f_sigma)
441 if (m_polar[i] == m_polar[j]) {
449 size_t kp = (m_polar[i] ? i : j);
450 size_t knp = (i == kp ? j : i);
451 double d3np, d3p, alpha_star, mu_p_star, xi;
452 d3np = pow(m_sigma[knp],3);
453 d3p = pow(m_sigma[kp],3);
454 alpha_star = m_alpha[knp]/d3np;
455 mu_p_star = m_dipole(kp,kp)/sqrt(4 *
Pi *
epsilon_0 * d3p * m_eps[kp]);
456 xi = 1.0 + 0.25 * alpha_star * mu_p_star * mu_p_star *
457 sqrt(m_eps[kp]/m_eps[knp]);
458 f_sigma = pow(xi, -1.0/6.0);
468 "fits to A*, B*, and C* vs. log(T*).\n"
469 "These are done only for the required dstar(j,k) values.\n\n");
470 if (m_log_level < 3) {
471 writelog(
"*** polynomial coefficients not printed (log_level < 3) ***\n");
474 vector<double> fitlist;
475 m_omega22_poly.clear();
476 m_astar_poly.clear();
477 m_bstar_poly.clear();
478 m_cstar_poly.clear();
479 for (
size_t i = 0; i < m_nsp; i++) {
480 for (
size_t j = i; j < m_nsp; j++) {
482 double dstar = (m_mode != CK_Mode) ? m_delta(i,j) : 0.0;
489 auto dptr = find(fitlist.begin(), fitlist.end(), dstar);
490 if (dptr == fitlist.end()) {
491 vector<double> ca(degree+1), cb(degree+1), cc(degree+1);
492 vector<double> co22(degree+1);
493 integrals.fit(degree, dstar, ca.data(), cb.data(), cc.data());
494 integrals.fit_omega22(degree, dstar, co22.data());
495 m_omega22_poly.push_back(co22);
496 m_astar_poly.push_back(ca);
497 m_bstar_poly.push_back(cb);
498 m_cstar_poly.push_back(cc);
499 m_poly[i][j] =
static_cast<int>(m_astar_poly.size()) - 1;
500 m_star_poly_uses_actualT[i][j] = 0;
501 fitlist.push_back(dstar);
504 m_poly[i][j] =
static_cast<int>((dptr - fitlist.begin()));
505 m_star_poly_uses_actualT[i][j] = 0;
507 m_poly[j][i] = m_poly[i][j];
508 m_star_poly_uses_actualT[j][i] = m_star_poly_uses_actualT[i][j];
516 const size_t np = 50;
517 int degree = (m_mode == CK_Mode ? 3 : 4);
518 double dt = (m_thermo->maxTemp() - m_thermo->minTemp())/(np-1);
519 vector<double> tlog(np), spvisc(np), spcond(np);
520 vector<double> w(np), w2(np);
522 m_visccoeffs.clear();
523 m_condcoeffs.clear();
526 for (
size_t n = 0; n < np; n++) {
527 double t = m_thermo->minTemp() + dt*n;
532 vector<double> c(degree + 1), c2(degree + 1);
535 if (m_log_level && m_log_level < 2) {
536 writelog(
"*** polynomial coefficients not printed (log_level < 2) ***\n");
538 double visc, err, relerr,
539 mxerr = 0.0, mxrelerr = 0.0, mxerr_cond = 0.0, mxrelerr_cond = 0.0;
542 writelog(
"Polynomial fits for viscosity:\n");
543 if (m_mode == CK_Mode) {
544 writelog(
"log(viscosity) fit to cubic polynomial in log(T)\n");
546 writelogf(
"viscosity/sqrt(T) fit to polynomial of degree "
547 "%d in log(T)", degree);
551 double T_save = m_thermo->temperature();
552 const vector<double>& mw = m_thermo->molecularWeights();
553 for (
size_t k = 0; k < m_nsp; k++) {
554 double tstar =
Boltzmann * 298.0 / m_eps[k];
557 double fz_298 = 1.0 + pow(
Pi, 1.5) / sqrt(tstar) * (0.5 + 1.0 / tstar) +
558 (0.25 *
Pi *
Pi + 2) / tstar;
560 for (
size_t n = 0; n < np; n++) {
561 double t = m_thermo->minTemp() + dt*n;
562 m_thermo->setTemperature(t);
563 vector<double> cp_R_all(m_thermo->nSpecies());
564 m_thermo->getCp_R_ref(&cp_R_all[0]);
565 double cp_R = cp_R_all[k];
567 double sqrt_T = sqrt(t);
568 double om22 = integrals.omega22(tstar, m_delta(k,k));
569 double om11 = integrals.omega11(tstar, m_delta(k,k));
572 double diffcoeff = 3.0/16.0 * sqrt(2.0 *
Pi/m_reducedMass(k,k)) *
574 (
Pi * m_sigma[k] * m_sigma[k] * om11);
578 (om22 *
Pi * m_sigma[k]*m_sigma[k]);
581 double f_int = mw[k]/(
GasConstant * t) * diffcoeff/visc;
582 double cv_rot = m_crot[k];
583 double A_factor = 2.5 - f_int;
584 double fz_tstar = 1.0 + pow(
Pi, 1.5) / sqrt(tstar) * (0.5 + 1.0 / tstar) +
585 (0.25 *
Pi *
Pi + 2) / tstar;
586 double B_factor = m_zrot[k] * fz_298 / fz_tstar + 2.0/
Pi * (5.0/3.0 * cv_rot + f_int);
587 double c1 = 2.0/
Pi * A_factor/B_factor;
588 double cv_int = cp_R - 2.5 - cv_rot;
589 double f_rot = f_int * (1.0 + c1);
590 double f_trans = 2.5 * (1.0 - c1 * cv_rot/1.5);
591 double cond = (visc/mw[k])*
GasConstant*(f_trans * 1.5
592 + f_rot * cv_rot + f_int * cv_int);
594 if (m_mode == CK_Mode) {
595 spvisc[n] = log(visc);
596 spcond[n] = log(cond);
606 spvisc[n] = sqrt(visc/sqrt_T);
611 spcond[n] = cond/sqrt_T;
612 w[n] = 1.0/(spvisc[n]*spvisc[n]);
613 w2[n] = 1.0/(spcond[n]*spcond[n]);
616 polyfit(np, degree, tlog.data(), spvisc.data(), w.data(), c.data());
617 polyfit(np, degree, tlog.data(), spcond.data(), w2.data(), c2.data());
620 for (
size_t n = 0; n < np; n++) {
622 if (m_mode == CK_Mode) {
623 val = exp(spvisc[n]);
624 fit = exp(
poly3(tlog[n], c.data()));
626 double sqrt_T = exp(0.5*tlog[n]);
627 val = sqrt_T * pow(spvisc[n],2);
628 fit = sqrt_T * pow(
poly4(tlog[n], c.data()),2);
632 mxerr = std::max(mxerr, fabs(err));
633 mxrelerr = std::max(mxrelerr, fabs(relerr));
635 m_fittingErrors[
"viscosity-max-abs-error"] = mxerr;
636 m_fittingErrors[
"viscosity-max-rel-error"] = mxrelerr;
639 for (
size_t n = 0; n < np; n++) {
641 if (m_mode == CK_Mode) {
642 val = exp(spcond[n]);
643 fit = exp(
poly3(tlog[n], c2.data()));
645 double sqrt_T = exp(0.5*tlog[n]);
646 val = sqrt_T * spcond[n];
647 fit = sqrt_T *
poly4(tlog[n], c2.data());
651 mxerr_cond = std::max(mxerr_cond, fabs(err));
652 mxrelerr_cond = std::max(mxrelerr_cond, fabs(relerr));
654 m_visccoeffs.push_back(c);
655 m_condcoeffs.push_back(c2);
656 m_fittingErrors[
"conductivity-max-abs-error"] = mxerr_cond;
657 m_fittingErrors[
"conductivity-max-rel-error"] = mxrelerr_cond;
659 if (m_log_level >= 2) {
663 m_thermo->setTemperature(T_save);
666 writelogf(
"Maximum viscosity absolute error: %12.6g\n", mxerr);
667 writelogf(
"Maximum viscosity relative error: %12.6g\n", mxrelerr);
668 writelog(
"\nPolynomial fits for conductivity:\n");
669 if (m_mode == CK_Mode) {
670 writelog(
"log(conductivity) fit to cubic polynomial in log(T)");
673 "polynomial of degree %d in log(T)", degree);
675 if (m_log_level >= 2) {
676 for (
size_t k = 0; k < m_nsp; k++) {
677 writelog(m_thermo->speciesName(k) +
": [" +
678 vec2str(m_condcoeffs[k]) +
"]\n");
681 writelogf(
"Maximum conductivity absolute error: %12.6g\n", mxerr_cond);
682 writelogf(
"Maximum conductivity relative error: %12.6g\n", mxrelerr_cond);
685 writelogf(
"\nbinary diffusion coefficients:\n");
686 if (m_mode == CK_Mode) {
687 writelog(
"log(D) fit to cubic polynomial in log(T)");
689 writelogf(
"D/T**(3/2) fit to polynomial of degree %d in log(T)",degree);
693 fitDiffCoeffs(integrals);
699 const size_t np = 50;
700 int degree = (m_mode == CK_Mode ? 3 : 4);
701 double dt = (m_thermo->maxTemp() - m_thermo->minTemp())/(np-1);
702 vector<double> tlog(np);
703 vector<double> w(np), w2(np);
706 for (
size_t n = 0; n < np; n++) {
707 double t = m_thermo->minTemp() + dt*n;
712 vector<double> c(degree + 1), c2(degree + 1);
713 double err, relerr, mxerr = 0.0, mxrelerr = 0.0;
715 vector<double> diff(np + 1);
716 m_diffcoeffs.clear();
717 for (
size_t k = 0; k < m_nsp; k++) {
718 for (
size_t j = k; j < m_nsp; j++) {
719 for (
size_t n = 0; n < np; n++) {
720 double t = m_thermo->minTemp() + dt*n;
721 double eps = m_epsilon(j,k);
723 double sigma = m_diam(j,k);
724 double om11 = integrals.omega11(tstar, m_delta(j,k));
725 double diffcoeff = 3.0/16.0 * sqrt(2.0 *
Pi/m_reducedMass(k,j))
726 * pow(
Boltzmann * t, 1.5) / (
Pi * sigma * sigma * om11);
731 getBinDiffCorrection(t, integrals, k, j, 1.0, 1.0, fkj, fjk);
733 if (m_mode == CK_Mode) {
734 diff[n] = log(diffcoeff);
737 diff[n] = diffcoeff/pow(t, 1.5);
738 w[n] = 1.0/(diff[n]*diff[n]);
741 polyfit(np, degree, tlog.data(), diff.data(), w.data(), c.data());
743 for (
size_t n = 0; n < np; n++) {
745 if (m_mode == CK_Mode) {
747 fit = exp(
poly3(tlog[n], c.data()));
749 double t = exp(tlog[n]);
750 double pre = pow(t, 1.5);
752 fit = pre *
poly4(tlog[n], c.data());
756 mxerr = std::max(mxerr, fabs(err));
757 mxrelerr = std::max(mxrelerr, fabs(relerr));
759 m_diffcoeffs.push_back(c);
760 if (m_log_level >= 2) {
761 writelog(m_thermo->speciesName(k) +
"__" +
762 m_thermo->speciesName(j) +
": [" +
vec2str(c) +
"]\n");
767 m_fittingErrors[
"diff-coeff-max-abs-error"] = mxerr;
768 m_fittingErrors[
"diff-coeff-max-rel-error"] = mxrelerr;
770 writelogf(
"Maximum binary diffusion coefficient absolute error:"
772 writelogf(
"Maximum binary diffusion coefficient relative error:"
778 size_t k,
size_t j,
double xk,
double xj,
double& fkj,
double& fjk)
780 double w1 = m_thermo->molecularWeight(k);
781 double w2 = m_thermo->molecularWeight(j);
782 double wsum = w1 + w2;
783 double wmwp = (w1 - w2)/wsum;
784 double sqw12 = sqrt(w1*w2);
785 double sig1 = m_sigma[k];
786 double sig2 = m_sigma[j];
787 double sig12 = 0.5*(m_sigma[k] + m_sigma[j]);
788 double sigratio = sig1*sig1/(sig2*sig2);
789 double sigratio2 = sig1*sig1/(sig12*sig12);
790 double sigratio3 = sig2*sig2/(sig12*sig12);
791 double tstar1 =
Boltzmann * t / m_eps[k];
792 double tstar2 =
Boltzmann * t / m_eps[j];
793 double tstar12 =
Boltzmann * t / sqrt(m_eps[k] * m_eps[j]);
794 double om22_1 = integrals.omega22(tstar1, m_delta(k,k));
795 double om22_2 = integrals.omega22(tstar2, m_delta(j,j));
796 double om11_12 = integrals.omega11(tstar12, m_delta(k,j));
797 double astar_12 = integrals.astar(tstar12, m_delta(k,j));
798 double bstar_12 = integrals.bstar(tstar12, m_delta(k,j));
799 double cstar_12 = integrals.cstar(tstar12, m_delta(k,j));
801 double cnst = sigratio * sqrt(2.0*w2/wsum) * 2.0 * w1*w1/(wsum * w2);
802 double p1 = cnst * om22_1 / om11_12;
804 cnst = (1.0/sigratio) * sqrt(2.0*w1/wsum) * 2.0*w2*w2/(wsum*w1);
805 double p2 = cnst * om22_2 / om11_12;
806 double p12 = 15.0 * wmwp*wmwp + 8.0*w1*w2*astar_12/(wsum*wsum);
808 cnst = (2.0/(w2*wsum))*sqrt(2.0*w2/wsum)*sigratio2;
809 double q1 = cnst*((2.5 - 1.2*bstar_12)*w1*w1 + 3.0*w2*w2
810 + 1.6*w1*w2*astar_12);
812 cnst = (2.0/(w1*wsum))*sqrt(2.0*w1/wsum)*sigratio3;
813 double q2 = cnst*((2.5 - 1.2*bstar_12)*w2*w2 + 3.0*w1*w1
814 + 1.6*w1*w2*astar_12);
815 double q12 = wmwp*wmwp*15.0*(2.5 - 1.2*bstar_12)
816 + 4.0*w1*w2*astar_12*(11.0 - 2.4*bstar_12)/(wsum*wsum)
817 + 1.6*wsum*om22_1*om22_2/(om11_12*om11_12*sqw12)
818 * sigratio2 * sigratio3;
820 cnst = 6.0*cstar_12 - 5.0;
821 fkj = 1.0 + 0.1*cnst*cnst *
822 (p1*xk*xk + p2*xj*xj + p12*xk*xj)/
823 (q1*xk*xk + q2*xj*xj + q12*xk*xj);
824 fjk = 1.0 + 0.1*cnst*cnst *
825 (p2*xk*xk + p1*xj*xj + p12*xk*xj)/
826 (q2*xk*xk + q1*xj*xj + q12*xk*xj);
829void GasTransport::getViscosityPolynomial(
size_t i,
double* coeffs)
const
831 checkSpeciesIndex(i);
832 for (
int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
833 coeffs[k] = m_visccoeffs[i][k];
837void GasTransport::getConductivityPolynomial(
size_t i,
double* coeffs)
const
839 checkSpeciesIndex(i);
840 for (
int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
841 coeffs[k] = m_condcoeffs[i][k];
845void GasTransport::getBinDiffusivityPolynomial(
size_t i,
size_t j,
double* coeffs)
const
847 checkSpeciesIndex(i);
848 checkSpeciesIndex(j);
849 size_t mi = (j >= i? i : j);
850 size_t mj = (j >= i? j : i);
852 for (
size_t ii = 0; ii < mi; ii++) {
857 for (
int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
858 coeffs[k] = m_diffcoeffs[ic][k];
862void GasTransport::getCollisionIntegralPolynomial(
size_t i,
size_t j,
863 double* astar_coeffs,
864 double* bstar_coeffs,
865 double* cstar_coeffs)
const
867 checkSpeciesIndex(i);
868 checkSpeciesIndex(j);
870 astar_coeffs[k] = m_astar_poly[m_poly[i][j]][k];
871 bstar_coeffs[k] = m_bstar_poly[m_poly[i][j]][k];
872 cstar_coeffs[k] = m_cstar_poly[m_poly[i][j]][k];
876void GasTransport::setViscosityPolynomial(
size_t i,
double* coeffs)
878 checkSpeciesIndex(i);
879 for (
int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
880 m_visccoeffs[i][k] = coeffs[k];
885void GasTransport::setConductivityPolynomial(
size_t i,
double* coeffs)
887 checkSpeciesIndex(i);
888 for (
int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
889 m_condcoeffs[i][k] = coeffs[k];
894void GasTransport::setBinDiffusivityPolynomial(
size_t i,
size_t j,
double* coeffs)
896 checkSpeciesIndex(i);
897 checkSpeciesIndex(j);
898 size_t mi = (j >= i? i : j);
899 size_t mj = (j >= i? j : i);
901 for (
size_t ii = 0; ii < mi; ii++) {
906 for (
int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
907 m_diffcoeffs[ic][k] = coeffs[k];
912void GasTransport::setCollisionIntegralPolynomial(
size_t i,
size_t j,
913 double* astar_coeffs,
914 double* bstar_coeffs,
915 double* cstar_coeffs,
bool actualT)
917 checkSpeciesIndex(i);
918 checkSpeciesIndex(j);
920 vector<double> ca(degree+1), cb(degree+1), cc(degree+1);
922 for (
size_t k = 0; k < degree+1; k++) {
923 ca[k] = astar_coeffs[k];
924 cb[k] = bstar_coeffs[k];
925 cc[k] = cstar_coeffs[k];
928 m_astar_poly.push_back(ca);
929 m_bstar_poly.push_back(cb);
930 m_cstar_poly.push_back(cc);
931 m_poly[j][i] = m_poly[i][j] =
static_cast<int>(m_astar_poly.size()) - 1;
932 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, int loglevel=-7)
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,...
string vec2str(const vector< double > &v, const string &fmt, const string &sep)
Convert a vector to a string (separated by commas)
void debuglog(const string &msg, int loglevel)
Write a message to the log only if loglevel > 0.
void writelogf(const char *fmt, const Args &... args)
Write a formatted message to the screen.
void writelog(const string &fmt, const Args &... args)
Write a formatted message to the screen.
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.
void warn_deprecated(const string &source, const AnyBase &node, const string &message)
A deprecation warning for syntax in an input file.
Contains declarations for string manipulation functions within Cantera.
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...