21#define COLL_INT_POLY_DEGREE 8
23GasTransport::GasTransport(ThermoPhase* thermo) :
29void GasTransport::update_T()
31 if (m_thermo->nSpecies() != m_nsp) {
33 init(m_thermo, m_mode, m_log_level);
36 double T = m_thermo->temperature();
44 m_sqrt_t = sqrt(m_temp);
45 m_t14 = sqrt(m_sqrt_t);
48 m_polytempvec[0] = 1.0;
49 m_polytempvec[1] = m_logt;
50 m_polytempvec[2] = m_logt*m_logt;
51 m_polytempvec[3] = m_logt*m_logt*m_logt;
52 m_polytempvec[4] = m_logt*m_logt*m_logt*m_logt;
61double GasTransport::viscosity()
76 multiply(m_phi, m_molefracs.data(), m_spwork.data());
78 for (
size_t k = 0; k < m_nsp; k++) {
79 vismix += m_molefracs[k] * m_visc[k]/m_spwork[k];
85void GasTransport::updateViscosity_T()
88 updateSpeciesViscosities();
92 for (
size_t j = 0; j < m_nsp; j++) {
93 for (
size_t k = j; k < m_nsp; k++) {
94 double vratiokj = m_visc[k]/m_visc[j];
95 double wratiojk = m_mw[j]/m_mw[k];
98 double factor1 = 1.0 + (m_sqvisc[k]/m_sqvisc[j]) * m_wratjk(k,j);
99 m_phi(k,j) = factor1*factor1 / (sqrt(8.0) * m_wratkj1(j,k));
100 m_phi(j,k) = m_phi(k,j)/(vratiokj * wratiojk);
106void GasTransport::updateSpeciesViscosities()
109 if (m_mode == CK_Mode) {
110 for (
size_t k = 0; k < m_nsp; k++) {
111 m_visc[k] = exp(
dot4(m_polytempvec, m_visccoeffs[k]));
112 m_sqvisc[k] = sqrt(m_visc[k]);
115 for (
size_t k = 0; k < m_nsp; k++) {
117 m_sqvisc[k] = m_t14 *
dot5(m_polytempvec, m_visccoeffs[k]);
118 m_visc[k] = (m_sqvisc[k] * m_sqvisc[k]);
124void GasTransport::updateDiff_T()
129 if (m_mode == CK_Mode) {
130 for (
size_t i = 0; i < m_nsp; i++) {
131 for (
size_t j = i; j < m_nsp; j++) {
132 m_bdiff(i,j) = exp(
dot4(m_polytempvec, m_diffcoeffs[ic]));
133 m_bdiff(j,i) = m_bdiff(i,j);
138 for (
size_t i = 0; i < m_nsp; i++) {
139 for (
size_t j = i; j < m_nsp; j++) {
140 m_bdiff(i,j) = m_temp * m_sqrt_t*
dot5(m_polytempvec,
142 m_bdiff(j,i) = m_bdiff(i,j);
150void GasTransport::getBinaryDiffCoeffs(
const size_t ld,
double*
const d)
158 throw CanteraError(
"GasTransport::getBinaryDiffCoeffs",
"ld is too small");
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(
double*
const d)
178 double mmw = m_thermo->meanMolecularWeight();
179 double p = m_thermo->pressure();
181 d[0] = m_bdiff(0,0) / p;
183 for (
size_t k = 0; k < m_nsp; k++) {
185 for (
size_t j = 0; j < m_nsp; j++) {
187 sum2 += m_molefracs[j] / m_bdiff(j,k);
191 d[k] = m_bdiff(k,k) / p;
193 d[k] = (mmw - m_molefracs[k] * m_mw[k])/(p * mmw * sum2);
199void GasTransport::getMixDiffCoeffsMole(
double*
const d)
209 double p = m_thermo->pressure();
211 d[0] = m_bdiff(0,0) / p;
213 for (
size_t k = 0; k < m_nsp; k++) {
215 for (
size_t j = 0; j < m_nsp; j++) {
217 sum2 += m_molefracs[j] / m_bdiff(j,k);
221 d[k] = m_bdiff(k,k) / p;
223 d[k] = (1 - m_molefracs[k]) / (p * sum2);
229void GasTransport::getMixDiffCoeffsMass(
double*
const d)
239 double mmw = m_thermo->meanMolecularWeight();
240 double p = m_thermo->pressure();
243 d[0] = m_bdiff(0,0) / p;
245 for (
size_t k=0; k<m_nsp; k++) {
248 for (
size_t i=0; i<m_nsp; i++) {
252 sum1 += m_molefracs[i] / m_bdiff(k,i);
253 sum2 += m_molefracs[i] * m_mw[i] / m_bdiff(k,i);
256 sum2 *= p * m_molefracs[k] / (mmw - m_mw[k]*m_molefracs[k]);
257 d[k] = 1.0 / (sum1 + sum2);
262void GasTransport::init(
ThermoPhase* thermo,
int mode,
int log_level)
267 m_log_level = log_level;
270 setupCollisionParameters();
271 setupCollisionIntegral();
273 m_molefracs.resize(m_nsp);
274 m_spwork.resize(m_nsp);
275 m_visc.resize(m_nsp);
276 m_sqvisc.resize(m_nsp);
277 m_phi.resize(m_nsp, m_nsp, 0.0);
278 m_bdiff.resize(m_nsp, m_nsp);
281 m_mw = m_thermo->molecularWeights();
283 m_wratjk.resize(m_nsp, m_nsp, 0.0);
284 m_wratkj1.resize(m_nsp, m_nsp, 0.0);
285 for (
size_t j = 0; j < m_nsp; j++) {
286 for (
size_t k = j; k < m_nsp; k++) {
287 m_wratjk(j,k) = sqrt(m_mw[j]/m_mw[k]);
288 m_wratjk(k,j) = sqrt(m_wratjk(j,k));
289 m_wratkj1(j,k) = sqrt(1.0 + m_mw[k]/m_mw[j]);
294void GasTransport::setupCollisionParameters()
296 m_epsilon.resize(m_nsp, m_nsp, 0.0);
297 m_delta.resize(m_nsp, m_nsp, 0.0);
298 m_reducedMass.resize(m_nsp, m_nsp, 0.0);
299 m_dipole.resize(m_nsp, m_nsp, 0.0);
300 m_diam.resize(m_nsp, m_nsp, 0.0);
301 m_crot.resize(m_nsp);
302 m_zrot.resize(m_nsp);
303 m_polar.resize(m_nsp,
false);
304 m_alpha.resize(m_nsp, 0.0);
305 m_poly.resize(m_nsp);
306 m_star_poly_uses_actualT.resize(m_nsp);
307 m_sigma.resize(m_nsp);
309 m_w_ac.resize(m_nsp);
310 m_disp.resize(m_nsp, 0.0);
311 m_quad_polar.resize(m_nsp, 0.0);
313 const vector<double>& mw = m_thermo->molecularWeights();
316 for (
size_t i = 0; i < m_nsp; i++) {
317 m_poly[i].resize(m_nsp);
318 m_star_poly_uses_actualT[i].resize(m_nsp);
321 double f_eps, f_sigma;
323 for (
size_t i = 0; i < m_nsp; i++) {
324 for (
size_t j = i; j < m_nsp; j++) {
326 m_reducedMass(i,j) = mw[i] * mw[j] / (
Avogadro * (mw[i] + mw[j]));
329 m_diam(i,j) = 0.5*(m_sigma[i] + m_sigma[j]);
332 m_epsilon(i,j) = sqrt(m_eps[i]*m_eps[j]);
335 m_dipole(i,j) = sqrt(m_dipole(i,i)*m_dipole(j,j));
338 double d = m_diam(i,j);
339 m_delta(i,j) = 0.5 * m_dipole(i,j)*m_dipole(i,j)
340 / (4 *
Pi *
epsilon_0 * m_epsilon(i,j) * d * d * d);
341 makePolarCorrections(i, j, f_eps, f_sigma);
342 m_diam(i,j) *= f_sigma;
343 m_epsilon(i,j) *= f_eps;
346 m_reducedMass(j,i) = m_reducedMass(i,j);
347 m_diam(j,i) = m_diam(i,j);
348 m_epsilon(j,i) = m_epsilon(i,j);
349 m_dipole(j,i) = m_dipole(i,j);
350 m_delta(j,i) = m_delta(i,j);
355void GasTransport::setupCollisionIntegral()
357 double tstar_min = 1.e8, tstar_max = 0.0;
358 for (
size_t i = 0; i < m_nsp; i++) {
359 for (
size_t j = i; j < m_nsp; j++) {
362 tstar_min = std::min(tstar_min,
Boltzmann * m_thermo->minTemp()/m_epsilon(i,j));
363 tstar_max = std::max(tstar_max,
Boltzmann * m_thermo->maxTemp()/m_epsilon(i,j));
368 if (m_mode == CK_Mode) {
374 debuglog(
"*** collision_integrals ***\n", m_log_level);
376 integrals.
init(tstar_min, tstar_max, m_log_level);
377 fitCollisionIntegrals(integrals);
378 debuglog(
"*** end of collision_integrals ***\n", m_log_level);
380 debuglog(
"*** property fits ***\n", m_log_level);
381 fitProperties(integrals);
382 debuglog(
"*** end of property fits ***\n", m_log_level);
385void GasTransport::getTransportData()
387 for (
size_t k = 0; k < m_thermo->nSpecies(); k++) {
388 shared_ptr<Species> s = m_thermo->species(k);
393 "Missing gas-phase transport data for species '{}'.", s->name);
398 }
else if (sptran->
geometry ==
"linear") {
400 }
else if (sptran->
geometry ==
"nonlinear") {
406 m_dipole(k,k) = sptran->
dipole;
407 m_polar[k] = (sptran->
dipole > 0);
410 if (s->input.hasKey(
"critical-parameters") &&
411 s->input[
"critical-parameters"].hasKey(
"acentric-factor"))
413 m_w_ac[k] = s->input[
"critical-parameters"][
"acentric-factor"].asDouble();
422void GasTransport::makePolarCorrections(
size_t i,
size_t j,
423 double& f_eps,
double& f_sigma)
426 if (m_polar[i] == m_polar[j]) {
434 size_t kp = (m_polar[i] ? i : j);
435 size_t knp = (i == kp ? j : i);
436 double d3np, d3p, alpha_star, mu_p_star, xi;
437 d3np = pow(m_sigma[knp],3);
438 d3p = pow(m_sigma[kp],3);
439 alpha_star = m_alpha[knp]/d3np;
440 mu_p_star = m_dipole(kp,kp)/sqrt(4 *
Pi *
epsilon_0 * d3p * m_eps[kp]);
441 xi = 1.0 + 0.25 * alpha_star * mu_p_star * mu_p_star *
442 sqrt(m_eps[kp]/m_eps[knp]);
443 f_sigma = pow(xi, -1.0/6.0);
453 "fits to A*, B*, and C* vs. log(T*).\n"
454 "These are done only for the required dstar(j,k) values.\n\n");
455 if (m_log_level < 3) {
456 writelog(
"*** polynomial coefficients not printed (log_level < 3) ***\n");
459 vector<double> fitlist;
460 m_omega22_poly.clear();
461 m_astar_poly.clear();
462 m_bstar_poly.clear();
463 m_cstar_poly.clear();
464 for (
size_t i = 0; i < m_nsp; i++) {
465 for (
size_t j = i; j < m_nsp; j++) {
467 double dstar = (m_mode != CK_Mode) ? m_delta(i,j) : 0.0;
474 auto dptr = find(fitlist.begin(), fitlist.end(), dstar);
475 if (dptr == fitlist.end()) {
476 vector<double> ca(degree+1), cb(degree+1), cc(degree+1);
477 vector<double> co22(degree+1);
478 integrals.fit(degree, dstar, ca.data(), cb.data(), cc.data());
479 integrals.fit_omega22(degree, dstar, co22.data());
480 m_omega22_poly.push_back(co22);
481 m_astar_poly.push_back(ca);
482 m_bstar_poly.push_back(cb);
483 m_cstar_poly.push_back(cc);
484 m_poly[i][j] =
static_cast<int>(m_astar_poly.size()) - 1;
485 m_star_poly_uses_actualT[i][j] = 0;
486 fitlist.push_back(dstar);
489 m_poly[i][j] =
static_cast<int>((dptr - fitlist.begin()));
490 m_star_poly_uses_actualT[i][j] = 0;
492 m_poly[j][i] = m_poly[i][j];
493 m_star_poly_uses_actualT[j][i] = m_star_poly_uses_actualT[i][j];
501 const size_t np = 50;
502 int degree = (m_mode == CK_Mode ? 3 : 4);
503 double dt = (m_thermo->maxTemp() - m_thermo->minTemp())/(np-1);
504 vector<double> tlog(np), spvisc(np), spcond(np);
505 vector<double> w(np), w2(np);
507 m_visccoeffs.clear();
508 m_condcoeffs.clear();
511 for (
size_t n = 0; n < np; n++) {
512 double t = m_thermo->minTemp() + dt*n;
517 vector<double> c(degree + 1), c2(degree + 1);
520 if (m_log_level && m_log_level < 2) {
521 writelog(
"*** polynomial coefficients not printed (log_level < 2) ***\n");
523 double visc, err, relerr,
524 mxerr = 0.0, mxrelerr = 0.0, mxerr_cond = 0.0, mxrelerr_cond = 0.0;
527 writelog(
"Polynomial fits for viscosity:\n");
528 if (m_mode == CK_Mode) {
529 writelog(
"log(viscosity) fit to cubic polynomial in log(T)\n");
531 writelogf(
"viscosity/sqrt(T) fit to polynomial of degree "
532 "%d in log(T)", degree);
536 double T_save = m_thermo->temperature();
537 const vector<double>& mw = m_thermo->molecularWeights();
538 for (
size_t k = 0; k < m_nsp; k++) {
539 double tstar =
Boltzmann * 298.0 / m_eps[k];
542 double fz_298 = 1.0 + pow(
Pi, 1.5) / sqrt(tstar) * (0.5 + 1.0 / tstar) +
543 (0.25 *
Pi *
Pi + 2) / tstar;
545 for (
size_t n = 0; n < np; n++) {
546 double t = m_thermo->minTemp() + dt*n;
547 m_thermo->setTemperature(t);
548 vector<double> cp_R_all(m_thermo->nSpecies());
549 m_thermo->getCp_R_ref(&cp_R_all[0]);
550 double cp_R = cp_R_all[k];
552 double sqrt_T = sqrt(t);
553 double om22 = integrals.omega22(tstar, m_delta(k,k));
554 double om11 = integrals.omega11(tstar, m_delta(k,k));
557 double diffcoeff = 3.0/16.0 * sqrt(2.0 *
Pi/m_reducedMass(k,k)) *
559 (
Pi * m_sigma[k] * m_sigma[k] * om11);
563 (om22 *
Pi * m_sigma[k]*m_sigma[k]);
566 double f_int = mw[k]/(
GasConstant * t) * diffcoeff/visc;
567 double cv_rot = m_crot[k];
568 double A_factor = 2.5 - f_int;
569 double fz_tstar = 1.0 + pow(
Pi, 1.5) / sqrt(tstar) * (0.5 + 1.0 / tstar) +
570 (0.25 *
Pi *
Pi + 2) / tstar;
571 double B_factor = m_zrot[k] * fz_298 / fz_tstar + 2.0/
Pi * (5.0/3.0 * cv_rot + f_int);
572 double c1 = 2.0/
Pi * A_factor/B_factor;
573 double cv_int = cp_R - 2.5 - cv_rot;
574 double f_rot = f_int * (1.0 + c1);
575 double f_trans = 2.5 * (1.0 - c1 * cv_rot/1.5);
576 double cond = (visc/mw[k])*
GasConstant*(f_trans * 1.5
577 + f_rot * cv_rot + f_int * cv_int);
579 if (m_mode == CK_Mode) {
580 spvisc[n] = log(visc);
581 spcond[n] = log(cond);
591 spvisc[n] = sqrt(visc/sqrt_T);
596 spcond[n] = cond/sqrt_T;
597 w[n] = 1.0/(spvisc[n]*spvisc[n]);
598 w2[n] = 1.0/(spcond[n]*spcond[n]);
601 polyfit(np, degree, tlog.data(), spvisc.data(), w.data(), c.data());
602 polyfit(np, degree, tlog.data(), spcond.data(), w2.data(), c2.data());
605 for (
size_t n = 0; n < np; n++) {
607 if (m_mode == CK_Mode) {
608 val = exp(spvisc[n]);
609 fit = exp(
poly3(tlog[n], c.data()));
611 double sqrt_T = exp(0.5*tlog[n]);
612 val = sqrt_T * pow(spvisc[n],2);
613 fit = sqrt_T * pow(
poly4(tlog[n], c.data()),2);
617 mxerr = std::max(mxerr, fabs(err));
618 mxrelerr = std::max(mxrelerr, fabs(relerr));
622 for (
size_t n = 0; n < np; n++) {
624 if (m_mode == CK_Mode) {
625 val = exp(spcond[n]);
626 fit = exp(
poly3(tlog[n], c2.data()));
628 double sqrt_T = exp(0.5*tlog[n]);
629 val = sqrt_T * spcond[n];
630 fit = sqrt_T *
poly4(tlog[n], c2.data());
634 mxerr_cond = std::max(mxerr_cond, fabs(err));
635 mxrelerr_cond = std::max(mxrelerr_cond, fabs(relerr));
637 m_visccoeffs.push_back(c);
638 m_condcoeffs.push_back(c2);
640 if (m_log_level >= 2) {
644 m_thermo->setTemperature(T_save);
647 writelogf(
"Maximum viscosity absolute error: %12.6g\n", mxerr);
648 writelogf(
"Maximum viscosity relative error: %12.6g\n", mxrelerr);
649 writelog(
"\nPolynomial fits for conductivity:\n");
650 if (m_mode == CK_Mode) {
651 writelog(
"log(conductivity) fit to cubic polynomial in log(T)");
654 "polynomial of degree %d in log(T)", degree);
656 if (m_log_level >= 2) {
657 for (
size_t k = 0; k < m_nsp; k++) {
658 writelog(m_thermo->speciesName(k) +
": [" +
659 vec2str(m_condcoeffs[k]) +
"]\n");
662 writelogf(
"Maximum conductivity absolute error: %12.6g\n", mxerr_cond);
663 writelogf(
"Maximum conductivity relative error: %12.6g\n", mxrelerr_cond);
666 writelogf(
"\nbinary diffusion coefficients:\n");
667 if (m_mode == CK_Mode) {
668 writelog(
"log(D) fit to cubic polynomial in log(T)");
670 writelogf(
"D/T**(3/2) fit to polynomial of degree %d in log(T)",degree);
674 fitDiffCoeffs(integrals);
680 const size_t np = 50;
681 int degree = (m_mode == CK_Mode ? 3 : 4);
682 double dt = (m_thermo->maxTemp() - m_thermo->minTemp())/(np-1);
683 vector<double> tlog(np);
684 vector<double> w(np), w2(np);
687 for (
size_t n = 0; n < np; n++) {
688 double t = m_thermo->minTemp() + dt*n;
693 vector<double> c(degree + 1), c2(degree + 1);
695 mxerr = 0.0, mxrelerr = 0.0;
697 vector<double> diff(np + 1);
698 m_diffcoeffs.clear();
699 for (
size_t k = 0; k < m_nsp; k++) {
700 for (
size_t j = k; j < m_nsp; j++) {
701 for (
size_t n = 0; n < np; n++) {
702 double t = m_thermo->minTemp() + dt*n;
703 double eps = m_epsilon(j,k);
705 double sigma = m_diam(j,k);
706 double om11 = integrals.omega11(tstar, m_delta(j,k));
707 double diffcoeff = 3.0/16.0 * sqrt(2.0 *
Pi/m_reducedMass(k,j))
708 * pow(
Boltzmann * t, 1.5) / (
Pi * sigma * sigma * om11);
713 getBinDiffCorrection(t, integrals, k, j, 1.0, 1.0, fkj, fjk);
715 if (m_mode == CK_Mode) {
716 diff[n] = log(diffcoeff);
719 diff[n] = diffcoeff/pow(t, 1.5);
720 w[n] = 1.0/(diff[n]*diff[n]);
723 polyfit(np, degree, tlog.data(), diff.data(), w.data(), c.data());
725 for (
size_t n = 0; n < np; n++) {
727 if (m_mode == CK_Mode) {
729 fit = exp(
poly3(tlog[n], c.data()));
731 double t = exp(tlog[n]);
732 double pre = pow(t, 1.5);
734 fit = pre *
poly4(tlog[n], c.data());
738 mxerr = std::max(mxerr, fabs(err));
739 mxrelerr = std::max(mxrelerr, fabs(relerr));
741 m_diffcoeffs.push_back(c);
742 if (m_log_level >= 2) {
743 writelog(m_thermo->speciesName(k) +
"__" +
744 m_thermo->speciesName(j) +
": [" +
vec2str(c) +
"]\n");
749 writelogf(
"Maximum binary diffusion coefficient absolute error:"
751 writelogf(
"Maximum binary diffusion coefficient relative error:"
757 size_t k,
size_t j,
double xk,
double xj,
double& fkj,
double& fjk)
759 double w1 = m_thermo->molecularWeight(k);
760 double w2 = m_thermo->molecularWeight(j);
761 double wsum = w1 + w2;
762 double wmwp = (w1 - w2)/wsum;
763 double sqw12 = sqrt(w1*w2);
764 double sig1 = m_sigma[k];
765 double sig2 = m_sigma[j];
766 double sig12 = 0.5*(m_sigma[k] + m_sigma[j]);
767 double sigratio = sig1*sig1/(sig2*sig2);
768 double sigratio2 = sig1*sig1/(sig12*sig12);
769 double sigratio3 = sig2*sig2/(sig12*sig12);
770 double tstar1 =
Boltzmann * t / m_eps[k];
771 double tstar2 =
Boltzmann * t / m_eps[j];
772 double tstar12 =
Boltzmann * t / sqrt(m_eps[k] * m_eps[j]);
773 double om22_1 = integrals.omega22(tstar1, m_delta(k,k));
774 double om22_2 = integrals.omega22(tstar2, m_delta(j,j));
775 double om11_12 = integrals.omega11(tstar12, m_delta(k,j));
776 double astar_12 = integrals.astar(tstar12, m_delta(k,j));
777 double bstar_12 = integrals.bstar(tstar12, m_delta(k,j));
778 double cstar_12 = integrals.cstar(tstar12, m_delta(k,j));
780 double cnst = sigratio * sqrt(2.0*w2/wsum) * 2.0 * w1*w1/(wsum * w2);
781 double p1 = cnst * om22_1 / om11_12;
783 cnst = (1.0/sigratio) * sqrt(2.0*w1/wsum) * 2.0*w2*w2/(wsum*w1);
784 double p2 = cnst * om22_2 / om11_12;
785 double p12 = 15.0 * wmwp*wmwp + 8.0*w1*w2*astar_12/(wsum*wsum);
787 cnst = (2.0/(w2*wsum))*sqrt(2.0*w2/wsum)*sigratio2;
788 double q1 = cnst*((2.5 - 1.2*bstar_12)*w1*w1 + 3.0*w2*w2
789 + 1.6*w1*w2*astar_12);
791 cnst = (2.0/(w1*wsum))*sqrt(2.0*w1/wsum)*sigratio3;
792 double q2 = cnst*((2.5 - 1.2*bstar_12)*w2*w2 + 3.0*w1*w1
793 + 1.6*w1*w2*astar_12);
794 double q12 = wmwp*wmwp*15.0*(2.5 - 1.2*bstar_12)
795 + 4.0*w1*w2*astar_12*(11.0 - 2.4*bstar_12)/(wsum*wsum)
796 + 1.6*wsum*om22_1*om22_2/(om11_12*om11_12*sqw12)
797 * sigratio2 * sigratio3;
799 cnst = 6.0*cstar_12 - 5.0;
800 fkj = 1.0 + 0.1*cnst*cnst *
801 (p1*xk*xk + p2*xj*xj + p12*xk*xj)/
802 (q1*xk*xk + q2*xj*xj + q12*xk*xj);
803 fjk = 1.0 + 0.1*cnst*cnst *
804 (p2*xk*xk + p1*xj*xj + p12*xk*xj)/
805 (q2*xk*xk + q1*xj*xj + q12*xk*xj);
808void GasTransport::getViscosityPolynomial(
size_t i,
double* coeffs)
const
810 for (
int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
811 coeffs[k] = m_visccoeffs[i][k];
815void GasTransport::getConductivityPolynomial(
size_t i,
double* coeffs)
const
817 for (
int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
818 coeffs[k] = m_condcoeffs[i][k];
822void GasTransport::getBinDiffusivityPolynomial(
size_t i,
size_t j,
double* coeffs)
const
824 size_t mi = (j >= i? i : j);
825 size_t mj = (j >= i? j : i);
827 for (
size_t ii = 0; ii < mi; ii++) {
832 for (
int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
833 coeffs[k] = m_diffcoeffs[ic][k];
837void GasTransport::getCollisionIntegralPolynomial(
size_t i,
size_t j,
838 double* astar_coeffs,
839 double* bstar_coeffs,
840 double* cstar_coeffs)
const
843 astar_coeffs[k] = m_astar_poly[m_poly[i][j]][k];
844 bstar_coeffs[k] = m_bstar_poly[m_poly[i][j]][k];
845 cstar_coeffs[k] = m_cstar_poly[m_poly[i][j]][k];
849void GasTransport::setViscosityPolynomial(
size_t i,
double* coeffs)
851 for (
int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
852 m_visccoeffs[i][k] = coeffs[k];
858 m_bindiff_ok =
false;
862void GasTransport::setConductivityPolynomial(
size_t i,
double* coeffs)
864 for (
int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
865 m_condcoeffs[i][k] = coeffs[k];
871 m_bindiff_ok =
false;
875void GasTransport::setBinDiffusivityPolynomial(
size_t i,
size_t j,
double* coeffs)
877 size_t mi = (j >= i? i : j);
878 size_t mj = (j >= i? j : i);
880 for (
size_t ii = 0; ii < mi; ii++) {
885 for (
int k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
886 m_diffcoeffs[ic][k] = coeffs[k];
892 m_bindiff_ok =
false;
896void GasTransport::setCollisionIntegralPolynomial(
size_t i,
size_t j,
897 double* astar_coeffs,
898 double* bstar_coeffs,
899 double* cstar_coeffs,
bool actualT)
902 vector<double> ca(degree+1), cb(degree+1), cc(degree+1);
904 for (
size_t k = 0; k < degree+1; k++) {
905 ca[k] = astar_coeffs[k];
906 cb[k] = bstar_coeffs[k];
907 cc[k] = cstar_coeffs[k];
910 m_astar_poly.push_back(ca);
911 m_bstar_poly.push_back(cb);
912 m_cstar_poly.push_back(cc);
913 m_poly[i][j] =
static_cast<int>(m_astar_poly.size()) - 1;
914 m_poly[j][i] = m_poly[i][j];
916 m_star_poly_uses_actualT[i][j] = 1;
917 m_star_poly_uses_actualT[j][i] = m_star_poly_uses_actualT[i][j];
923 m_bindiff_ok =
false;
#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=0)
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.
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...