21#define COLL_INT_POLY_DEGREE 8
23GasTransport::GasTransport(ThermoPhase* thermo) :
43void GasTransport::update_T()
45 if (m_thermo->nSpecies() != m_nsp) {
47 init(m_thermo, m_mode, m_log_level);
50 double T = m_thermo->temperature();
59 m_sqrt_t = sqrt(m_temp);
60 m_t14 = sqrt(m_sqrt_t);
61 m_t32 = m_temp * m_sqrt_t;
64 m_polytempvec[0] = 1.0;
65 m_polytempvec[1] = m_logt;
66 m_polytempvec[2] = m_logt*m_logt;
67 m_polytempvec[3] = m_logt*m_logt*m_logt;
68 m_polytempvec[4] = m_logt*m_logt*m_logt*m_logt;
77doublereal GasTransport::viscosity()
86 doublereal vismix = 0.0;
92 multiply(m_phi, m_molefracs.data(), m_spwork.data());
94 for (
size_t k = 0; k < m_nsp; k++) {
95 vismix += m_molefracs[k] * m_visc[k]/m_spwork[k];
101void GasTransport::updateViscosity_T()
104 updateSpeciesViscosities();
108 for (
size_t j = 0; j < m_nsp; j++) {
109 for (
size_t k = j; k < m_nsp; k++) {
110 double vratiokj = m_visc[k]/m_visc[j];
111 double wratiojk = m_mw[j]/m_mw[k];
114 double factor1 = 1.0 + (m_sqvisc[k]/m_sqvisc[j]) * m_wratjk(k,j);
115 m_phi(k,j) = factor1*factor1 / (sqrt(8.0) * m_wratkj1(j,k));
116 m_phi(j,k) = m_phi(k,j)/(vratiokj * wratiojk);
122void GasTransport::updateSpeciesViscosities()
125 if (m_mode == CK_Mode) {
126 for (
size_t k = 0; k < m_nsp; k++) {
127 m_visc[k] = exp(
dot4(m_polytempvec, m_visccoeffs[k]));
128 m_sqvisc[k] = sqrt(m_visc[k]);
131 for (
size_t k = 0; k < m_nsp; k++) {
133 m_sqvisc[k] = m_t14 *
dot5(m_polytempvec, m_visccoeffs[k]);
134 m_visc[k] = (m_sqvisc[k] * m_sqvisc[k]);
140void GasTransport::updateDiff_T()
145 if (m_mode == CK_Mode) {
146 for (
size_t i = 0; i < m_nsp; i++) {
147 for (
size_t j = i; j < m_nsp; j++) {
148 m_bdiff(i,j) = exp(
dot4(m_polytempvec, m_diffcoeffs[ic]));
149 m_bdiff(j,i) = m_bdiff(i,j);
154 for (
size_t i = 0; i < m_nsp; i++) {
155 for (
size_t j = i; j < m_nsp; j++) {
156 m_bdiff(i,j) = m_temp * m_sqrt_t*
dot5(m_polytempvec,
158 m_bdiff(j,i) = m_bdiff(i,j);
166void GasTransport::getBinaryDiffCoeffs(
const size_t ld, doublereal*
const d)
174 throw CanteraError(
"GasTransport::getBinaryDiffCoeffs",
"ld is too small");
176 doublereal rp = 1.0/m_thermo->pressure();
177 for (
size_t i = 0; i < m_nsp; i++) {
178 for (
size_t j = 0; j < m_nsp; j++) {
179 d[ld*j + i] = rp * m_bdiff(i,j);
184void GasTransport::getMixDiffCoeffs(doublereal*
const d)
194 doublereal mmw = m_thermo->meanMolecularWeight();
195 doublereal p = m_thermo->pressure();
197 d[0] = m_bdiff(0,0) / p;
199 for (
size_t k = 0; k < m_nsp; k++) {
201 for (
size_t j = 0; j < m_nsp; j++) {
203 sum2 += m_molefracs[j] / m_bdiff(j,k);
207 d[k] = m_bdiff(k,k) / p;
209 d[k] = (mmw - m_molefracs[k] * m_mw[k])/(p * mmw * sum2);
215void GasTransport::getMixDiffCoeffsMole(doublereal*
const d)
225 doublereal p = m_thermo->pressure();
227 d[0] = m_bdiff(0,0) / p;
229 for (
size_t k = 0; k < m_nsp; k++) {
231 for (
size_t j = 0; j < m_nsp; j++) {
233 sum2 += m_molefracs[j] / m_bdiff(j,k);
237 d[k] = m_bdiff(k,k) / p;
239 d[k] = (1 - m_molefracs[k]) / (p * sum2);
245void GasTransport::getMixDiffCoeffsMass(doublereal*
const d)
255 doublereal mmw = m_thermo->meanMolecularWeight();
256 doublereal p = m_thermo->pressure();
259 d[0] = m_bdiff(0,0) / p;
261 for (
size_t k=0; k<m_nsp; k++) {
264 for (
size_t i=0; i<m_nsp; i++) {
268 sum1 += m_molefracs[i] / m_bdiff(k,i);
269 sum2 += m_molefracs[i] * m_mw[i] / m_bdiff(k,i);
272 sum2 *= p * m_molefracs[k] / (mmw - m_mw[k]*m_molefracs[k]);
273 d[k] = 1.0 / (sum1 + sum2);
278void GasTransport::init(
ThermoPhase* thermo,
int mode,
int log_level)
283 m_log_level = log_level;
286 setupCollisionParameters();
287 setupCollisionIntegral();
289 m_molefracs.resize(m_nsp);
290 m_spwork.resize(m_nsp);
291 m_visc.resize(m_nsp);
292 m_sqvisc.resize(m_nsp);
293 m_phi.resize(m_nsp, m_nsp, 0.0);
294 m_bdiff.resize(m_nsp, m_nsp);
297 m_mw = m_thermo->molecularWeights();
299 m_wratjk.resize(m_nsp, m_nsp, 0.0);
300 m_wratkj1.resize(m_nsp, m_nsp, 0.0);
301 for (
size_t j = 0; j < m_nsp; j++) {
302 for (
size_t k = j; k < m_nsp; k++) {
303 m_wratjk(j,k) = sqrt(m_mw[j]/m_mw[k]);
304 m_wratjk(k,j) = sqrt(m_wratjk(j,k));
305 m_wratkj1(j,k) = sqrt(1.0 + m_mw[k]/m_mw[j]);
310void GasTransport::setupCollisionParameters()
312 m_epsilon.resize(m_nsp, m_nsp, 0.0);
313 m_delta.resize(m_nsp, m_nsp, 0.0);
314 m_reducedMass.resize(m_nsp, m_nsp, 0.0);
315 m_dipole.resize(m_nsp, m_nsp, 0.0);
316 m_diam.resize(m_nsp, m_nsp, 0.0);
317 m_crot.resize(m_nsp);
318 m_zrot.resize(m_nsp);
319 m_polar.resize(m_nsp,
false);
320 m_alpha.resize(m_nsp, 0.0);
321 m_poly.resize(m_nsp);
322 m_star_poly_uses_actualT.resize(m_nsp);
323 m_sigma.resize(m_nsp);
325 m_w_ac.resize(m_nsp);
326 m_disp.resize(m_nsp, 0.0);
327 m_quad_polar.resize(m_nsp, 0.0);
329 const vector_fp& mw = m_thermo->molecularWeights();
332 for (
size_t i = 0; i < m_nsp; i++) {
333 m_poly[i].resize(m_nsp);
334 m_star_poly_uses_actualT[i].resize(m_nsp);
337 double f_eps, f_sigma;
339 for (
size_t i = 0; i < m_nsp; i++) {
340 for (
size_t j = i; j < m_nsp; j++) {
342 m_reducedMass(i,j) = mw[i] * mw[j] / (
Avogadro * (mw[i] + mw[j]));
345 m_diam(i,j) = 0.5*(m_sigma[i] + m_sigma[j]);
348 m_epsilon(i,j) = sqrt(m_eps[i]*m_eps[j]);
351 m_dipole(i,j) = sqrt(m_dipole(i,i)*m_dipole(j,j));
354 double d = m_diam(i,j);
355 m_delta(i,j) = 0.5 * m_dipole(i,j)*m_dipole(i,j)
356 / (4 *
Pi *
epsilon_0 * m_epsilon(i,j) * d * d * d);
357 makePolarCorrections(i, j, f_eps, f_sigma);
358 m_diam(i,j) *= f_sigma;
359 m_epsilon(i,j) *= f_eps;
362 m_reducedMass(j,i) = m_reducedMass(i,j);
363 m_diam(j,i) = m_diam(i,j);
364 m_epsilon(j,i) = m_epsilon(i,j);
365 m_dipole(j,i) = m_dipole(i,j);
366 m_delta(j,i) = m_delta(i,j);
371void GasTransport::setupCollisionIntegral()
373 double tstar_min = 1.e8, tstar_max = 0.0;
374 for (
size_t i = 0; i < m_nsp; i++) {
375 for (
size_t j = i; j < m_nsp; j++) {
378 tstar_min = std::min(tstar_min,
Boltzmann * m_thermo->minTemp()/m_epsilon(i,j));
379 tstar_max = std::max(tstar_max,
Boltzmann * m_thermo->maxTemp()/m_epsilon(i,j));
384 if (m_mode == CK_Mode) {
390 debuglog(
"*** collision_integrals ***\n", m_log_level);
392 integrals.
init(tstar_min, tstar_max, m_log_level);
393 fitCollisionIntegrals(integrals);
394 debuglog(
"*** end of collision_integrals ***\n", m_log_level);
396 debuglog(
"*** property fits ***\n", m_log_level);
397 fitProperties(integrals);
398 debuglog(
"*** end of property fits ***\n", m_log_level);
401void GasTransport::getTransportData()
403 for (
size_t k = 0; k < m_thermo->nSpecies(); k++) {
404 shared_ptr<Species> s = m_thermo->species(k);
409 "Missing gas-phase transport data for species '{}'.", s->name);
414 }
else if (sptran->
geometry ==
"linear") {
416 }
else if (sptran->
geometry ==
"nonlinear") {
422 m_dipole(k,k) = sptran->
dipole;
423 m_polar[k] = (sptran->
dipole > 0);
426 if (s->input.hasKey(
"critical-parameters") &&
427 s->input[
"critical-parameters"].hasKey(
"acentric-factor"))
429 m_w_ac[k] = s->input[
"critical-parameters"][
"acentric-factor"].asDouble();
438void GasTransport::makePolarCorrections(
size_t i,
size_t j,
439 doublereal& f_eps, doublereal& f_sigma)
442 if (m_polar[i] == m_polar[j]) {
450 size_t kp = (m_polar[i] ? i : j);
451 size_t knp = (i == kp ? j : i);
452 double d3np, d3p, alpha_star, mu_p_star, xi;
453 d3np = pow(m_sigma[knp],3);
454 d3p = pow(m_sigma[kp],3);
455 alpha_star = m_alpha[knp]/d3np;
456 mu_p_star = m_dipole(kp,kp)/sqrt(4 *
Pi *
epsilon_0 * d3p * m_eps[kp]);
457 xi = 1.0 + 0.25 * alpha_star * mu_p_star * mu_p_star *
458 sqrt(m_eps[kp]/m_eps[knp]);
459 f_sigma = pow(xi, -1.0/6.0);
469 "fits to A*, B*, and C* vs. log(T*).\n"
470 "These are done only for the required dstar(j,k) values.\n\n");
471 if (m_log_level < 3) {
472 writelog(
"*** polynomial coefficients not printed (log_level < 3) ***\n");
476 m_omega22_poly.clear();
477 m_astar_poly.clear();
478 m_bstar_poly.clear();
479 m_cstar_poly.clear();
480 for (
size_t i = 0; i < m_nsp; i++) {
481 for (
size_t j = i; j < m_nsp; j++) {
483 double dstar = (m_mode != CK_Mode) ? m_delta(i,j) : 0.0;
490 auto dptr = find(fitlist.begin(), fitlist.end(), dstar);
491 if (dptr == fitlist.end()) {
492 vector_fp ca(degree+1), cb(degree+1), cc(degree+1);
494 integrals.fit(degree, dstar, ca.data(), cb.data(), cc.data());
495 integrals.fit_omega22(degree, dstar, co22.data());
496 m_omega22_poly.push_back(co22);
497 m_astar_poly.push_back(ca);
498 m_bstar_poly.push_back(cb);
499 m_cstar_poly.push_back(cc);
500 m_poly[i][j] =
static_cast<int>(m_astar_poly.size()) - 1;
501 m_star_poly_uses_actualT[i][j] = 0;
502 fitlist.push_back(dstar);
505 m_poly[i][j] =
static_cast<int>((dptr - fitlist.begin()));
506 m_star_poly_uses_actualT[i][j] = 0;
508 m_poly[j][i] = m_poly[i][j];
509 m_star_poly_uses_actualT[j][i] = m_star_poly_uses_actualT[i][j];
517 const size_t np = 50;
518 int degree = (m_mode == CK_Mode ? 3 : 4);
519 double dt = (m_thermo->maxTemp() - m_thermo->minTemp())/(np-1);
520 vector_fp tlog(np), spvisc(np), spcond(np);
523 m_visccoeffs.clear();
524 m_condcoeffs.clear();
527 for (
size_t n = 0; n < np; n++) {
528 double t = m_thermo->minTemp() + dt*n;
536 if (m_log_level && m_log_level < 2) {
537 writelog(
"*** polynomial coefficients not printed (log_level < 2) ***\n");
539 double visc, err, relerr,
540 mxerr = 0.0, mxrelerr = 0.0, mxerr_cond = 0.0, mxrelerr_cond = 0.0;
543 writelog(
"Polynomial fits for viscosity:\n");
544 if (m_mode == CK_Mode) {
545 writelog(
"log(viscosity) fit to cubic polynomial in log(T)\n");
547 writelogf(
"viscosity/sqrt(T) fit to polynomial of degree "
548 "%d in log(T)", degree);
552 double T_save = m_thermo->temperature();
553 const vector_fp& mw = m_thermo->molecularWeights();
554 for (
size_t k = 0; k < m_nsp; k++) {
555 double tstar =
Boltzmann * 298.0 / m_eps[k];
558 double fz_298 = 1.0 + pow(
Pi, 1.5) / sqrt(tstar) * (0.5 + 1.0 / tstar) +
559 (0.25 *
Pi *
Pi + 2) / tstar;
561 for (
size_t n = 0; n < np; n++) {
562 double t = m_thermo->minTemp() + dt*n;
563 m_thermo->setTemperature(t);
564 vector_fp cp_R_all(m_thermo->nSpecies());
565 m_thermo->getCp_R_ref(&cp_R_all[0]);
566 double cp_R = cp_R_all[k];
568 double sqrt_T = sqrt(t);
569 double om22 = integrals.omega22(tstar, m_delta(k,k));
570 double om11 = integrals.omega11(tstar, m_delta(k,k));
573 double diffcoeff = 3.0/16.0 * sqrt(2.0 *
Pi/m_reducedMass(k,k)) *
575 (
Pi * m_sigma[k] * m_sigma[k] * om11);
579 (om22 *
Pi * m_sigma[k]*m_sigma[k]);
582 double f_int = mw[k]/(
GasConstant * t) * diffcoeff/visc;
583 double cv_rot = m_crot[k];
584 double A_factor = 2.5 - f_int;
585 double fz_tstar = 1.0 + pow(
Pi, 1.5) / sqrt(tstar) * (0.5 + 1.0 / tstar) +
586 (0.25 *
Pi *
Pi + 2) / tstar;
587 double B_factor = m_zrot[k] * fz_298 / fz_tstar + 2.0/
Pi * (5.0/3.0 * cv_rot + f_int);
588 double c1 = 2.0/
Pi * A_factor/B_factor;
589 double cv_int = cp_R - 2.5 - cv_rot;
590 double f_rot = f_int * (1.0 + c1);
591 double f_trans = 2.5 * (1.0 - c1 * cv_rot/1.5);
592 double cond = (visc/mw[k])*
GasConstant*(f_trans * 1.5
593 + f_rot * cv_rot + f_int * cv_int);
595 if (m_mode == CK_Mode) {
596 spvisc[n] = log(visc);
597 spcond[n] = log(cond);
607 spvisc[n] = sqrt(visc/sqrt_T);
612 spcond[n] = cond/sqrt_T;
613 w[n] = 1.0/(spvisc[n]*spvisc[n]);
614 w2[n] = 1.0/(spcond[n]*spcond[n]);
617 polyfit(np, degree, tlog.data(), spvisc.data(), w.data(), c.data());
618 polyfit(np, degree, tlog.data(), spcond.data(), w2.data(), c2.data());
621 for (
size_t n = 0; n < np; n++) {
623 if (m_mode == CK_Mode) {
624 val = exp(spvisc[n]);
625 fit = exp(
poly3(tlog[n], c.data()));
627 double sqrt_T = exp(0.5*tlog[n]);
628 val = sqrt_T * pow(spvisc[n],2);
629 fit = sqrt_T * pow(
poly4(tlog[n], c.data()),2);
633 mxerr = std::max(mxerr, fabs(err));
634 mxrelerr = std::max(mxrelerr, fabs(relerr));
638 for (
size_t n = 0; n < np; n++) {
640 if (m_mode == CK_Mode) {
641 val = exp(spcond[n]);
642 fit = exp(
poly3(tlog[n], c2.data()));
644 double sqrt_T = exp(0.5*tlog[n]);
645 val = sqrt_T * spcond[n];
646 fit = sqrt_T *
poly4(tlog[n], c2.data());
650 mxerr_cond = std::max(mxerr_cond, fabs(err));
651 mxrelerr_cond = std::max(mxrelerr_cond, fabs(relerr));
653 m_visccoeffs.push_back(c);
654 m_condcoeffs.push_back(c2);
656 if (m_log_level >= 2) {
660 m_thermo->setTemperature(T_save);
663 writelogf(
"Maximum viscosity absolute error: %12.6g\n", mxerr);
664 writelogf(
"Maximum viscosity relative error: %12.6g\n", mxrelerr);
665 writelog(
"\nPolynomial fits for conductivity:\n");
666 if (m_mode == CK_Mode) {
667 writelog(
"log(conductivity) fit to cubic polynomial in log(T)");
670 "polynomial of degree %d in log(T)", degree);
672 if (m_log_level >= 2) {
673 for (
size_t k = 0; k < m_nsp; k++) {
674 writelog(m_thermo->speciesName(k) +
": [" +
675 vec2str(m_condcoeffs[k]) +
"]\n");
678 writelogf(
"Maximum conductivity absolute error: %12.6g\n", mxerr_cond);
679 writelogf(
"Maximum conductivity relative error: %12.6g\n", mxrelerr_cond);
682 writelogf(
"\nbinary diffusion coefficients:\n");
683 if (m_mode == CK_Mode) {
684 writelog(
"log(D) fit to cubic polynomial in log(T)");
686 writelogf(
"D/T**(3/2) fit to polynomial of degree %d in log(T)",degree);
690 fitDiffCoeffs(integrals);
696 const size_t np = 50;
697 int degree = (m_mode == CK_Mode ? 3 : 4);
698 double dt = (m_thermo->maxTemp() - m_thermo->minTemp())/(np-1);
703 for (
size_t n = 0; n < np; n++) {
704 double t = m_thermo->minTemp() + dt*n;
711 mxerr = 0.0, mxrelerr = 0.0;
714 m_diffcoeffs.clear();
715 for (
size_t k = 0; k < m_nsp; k++) {
716 for (
size_t j = k; j < m_nsp; j++) {
717 for (
size_t n = 0; n < np; n++) {
718 double t = m_thermo->minTemp() + dt*n;
719 double eps = m_epsilon(j,k);
721 double sigma = m_diam(j,k);
722 double om11 = integrals.omega11(tstar, m_delta(j,k));
723 double diffcoeff = 3.0/16.0 * sqrt(2.0 *
Pi/m_reducedMass(k,j))
724 * pow(
Boltzmann * t, 1.5) / (
Pi * sigma * sigma * om11);
729 getBinDiffCorrection(t, integrals, k, j, 1.0, 1.0, fkj, fjk);
731 if (m_mode == CK_Mode) {
732 diff[n] = log(diffcoeff);
735 diff[n] = diffcoeff/pow(t, 1.5);
736 w[n] = 1.0/(diff[n]*diff[n]);
739 polyfit(np, degree, tlog.data(), diff.data(), w.data(), c.data());
741 for (
size_t n = 0; n < np; n++) {
743 if (m_mode == CK_Mode) {
745 fit = exp(
poly3(tlog[n], c.data()));
747 double t = exp(tlog[n]);
748 double pre = pow(t, 1.5);
750 fit = pre *
poly4(tlog[n], c.data());
754 mxerr = std::max(mxerr, fabs(err));
755 mxrelerr = std::max(mxrelerr, fabs(relerr));
757 m_diffcoeffs.push_back(c);
758 if (m_log_level >= 2) {
759 writelog(m_thermo->speciesName(k) +
"__" +
760 m_thermo->speciesName(j) +
": [" +
vec2str(c) +
"]\n");
765 writelogf(
"Maximum binary diffusion coefficient absolute error:"
767 writelogf(
"Maximum binary diffusion coefficient relative error:"
773 size_t k,
size_t j,
double xk,
double xj,
double& fkj,
double& fjk)
775 double w1 = m_thermo->molecularWeight(k);
776 double w2 = m_thermo->molecularWeight(j);
777 double wsum = w1 + w2;
778 double wmwp = (w1 - w2)/wsum;
779 double sqw12 = sqrt(w1*w2);
780 double sig1 = m_sigma[k];
781 double sig2 = m_sigma[j];
782 double sig12 = 0.5*(m_sigma[k] + m_sigma[j]);
783 double sigratio = sig1*sig1/(sig2*sig2);
784 double sigratio2 = sig1*sig1/(sig12*sig12);
785 double sigratio3 = sig2*sig2/(sig12*sig12);
786 double tstar1 =
Boltzmann * t / m_eps[k];
787 double tstar2 =
Boltzmann * t / m_eps[j];
788 double tstar12 =
Boltzmann * t / sqrt(m_eps[k] * m_eps[j]);
789 double om22_1 = integrals.omega22(tstar1, m_delta(k,k));
790 double om22_2 = integrals.omega22(tstar2, m_delta(j,j));
791 double om11_12 = integrals.omega11(tstar12, m_delta(k,j));
792 double astar_12 = integrals.astar(tstar12, m_delta(k,j));
793 double bstar_12 = integrals.bstar(tstar12, m_delta(k,j));
794 double cstar_12 = integrals.cstar(tstar12, m_delta(k,j));
796 double cnst = sigratio * sqrt(2.0*w2/wsum) * 2.0 * w1*w1/(wsum * w2);
797 double p1 = cnst * om22_1 / om11_12;
799 cnst = (1.0/sigratio) * sqrt(2.0*w1/wsum) * 2.0*w2*w2/(wsum*w1);
800 double p2 = cnst * om22_2 / om11_12;
801 double p12 = 15.0 * wmwp*wmwp + 8.0*w1*w2*astar_12/(wsum*wsum);
803 cnst = (2.0/(w2*wsum))*sqrt(2.0*w2/wsum)*sigratio2;
804 double q1 = cnst*((2.5 - 1.2*bstar_12)*w1*w1 + 3.0*w2*w2
805 + 1.6*w1*w2*astar_12);
807 cnst = (2.0/(w1*wsum))*sqrt(2.0*w1/wsum)*sigratio3;
808 double q2 = cnst*((2.5 - 1.2*bstar_12)*w2*w2 + 3.0*w1*w1
809 + 1.6*w1*w2*astar_12);
810 double q12 = wmwp*wmwp*15.0*(2.5 - 1.2*bstar_12)
811 + 4.0*w1*w2*astar_12*(11.0 - 2.4*bstar_12)/(wsum*wsum)
812 + 1.6*wsum*om22_1*om22_2/(om11_12*om11_12*sqw12)
813 * sigratio2 * sigratio3;
815 cnst = 6.0*cstar_12 - 5.0;
816 fkj = 1.0 + 0.1*cnst*cnst *
817 (p1*xk*xk + p2*xj*xj + p12*xk*xj)/
818 (q1*xk*xk + q2*xj*xj + q12*xk*xj);
819 fjk = 1.0 + 0.1*cnst*cnst *
820 (p2*xk*xk + p1*xj*xj + p12*xk*xj)/
821 (q2*xk*xk + q1*xj*xj + q12*xk*xj);
824void GasTransport::getViscosityPolynomial(
size_t i,
double* coeffs)
const
826 for (
size_t k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
827 coeffs[k] = m_visccoeffs[i][k];
831void GasTransport::getConductivityPolynomial(
size_t i,
double* coeffs)
const
833 for (
size_t k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
834 coeffs[k] = m_condcoeffs[i][k];
838void GasTransport::getBinDiffusivityPolynomial(
size_t i,
size_t j,
double* coeffs)
const
840 size_t mi = (j >= i? i : j);
841 size_t mj = (j >= i? j : i);
843 for (
size_t ii = 0; ii < mi; ii++) {
848 for (
size_t k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
849 coeffs[k] = m_diffcoeffs[ic][k];
853void GasTransport::getCollisionIntegralPolynomial(
size_t i,
size_t j,
854 double* astar_coeffs,
855 double* bstar_coeffs,
856 double* cstar_coeffs)
const
859 astar_coeffs[k] = m_astar_poly[m_poly[i][j]][k];
860 bstar_coeffs[k] = m_bstar_poly[m_poly[i][j]][k];
861 cstar_coeffs[k] = m_cstar_poly[m_poly[i][j]][k];
865void GasTransport::setViscosityPolynomial(
size_t i,
double* coeffs)
867 for (
size_t k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
868 m_visccoeffs[i][k] = coeffs[k];
874 m_bindiff_ok =
false;
878void GasTransport::setConductivityPolynomial(
size_t i,
double* coeffs)
880 for (
size_t k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
881 m_condcoeffs[i][k] = coeffs[k];
887 m_bindiff_ok =
false;
891void GasTransport::setBinDiffusivityPolynomial(
size_t i,
size_t j,
double* coeffs)
893 size_t mi = (j >= i? i : j);
894 size_t mj = (j >= i? j : i);
896 for (
size_t ii = 0; ii < mi; ii++) {
901 for (
size_t k = 0; k < (m_mode == CK_Mode ? 4 : 5); k++) {
902 m_diffcoeffs[ic][k] = coeffs[k];
908 m_bindiff_ok =
false;
912void GasTransport::setCollisionIntegralPolynomial(
size_t i,
size_t j,
913 double* astar_coeffs,
914 double* bstar_coeffs,
915 double* cstar_coeffs,
bool actualT)
918 vector_fp ca(degree+1), cb(degree+1), cc(degree+1);
920 for (
size_t k = 0; k < degree+1; k++) {
921 ca[k] = astar_coeffs[k];
922 cb[k] = bstar_coeffs[k];
923 cc[k] = cstar_coeffs[k];
926 m_astar_poly.push_back(ca);
927 m_bstar_poly.push_back(cb);
928 m_cstar_poly.push_back(cc);
929 m_poly[i][j] =
static_cast<int>(m_astar_poly.size()) - 1;
930 m_poly[j][i] = m_poly[i][j];
932 m_star_poly_uses_actualT[i][j] = 1;
933 m_star_poly_uses_actualT[j][i] = m_star_poly_uses_actualT[i][j];
939 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.
Monk and Monchick 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.
std::string geometry
A string specifying the molecular geometry.
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].
Calculation of Collision integrals.
void init(doublereal tsmin, doublereal 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,...
void writelog(const std::string &fmt, const Args &... args)
Write a formatted message to the screen.
void writelogf(const char *fmt, const Args &... args)
Write a formatted message to the screen.
Namespace for the Cantera kernel.
const double Boltzmann
Boltzmann constant [J/K].
const double Avogadro
Avogadro's Number [number/kmol].
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.
std::string vec2str(const vector_fp &v, const std::string &fmt="%g", const std::string &sep=", ")
Convert a vector to a string (separated by commas)
void debuglog(const std::string &msg, int loglevel)
Write a message to the log only if loglevel > 0.
doublereal dot4(const V &x, const V &y)
Templated Inner product of two vectors of length 4.
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
const double epsilon_0
Permittivity of free space [F/m].
const double GasConstant
Universal Gas Constant [J/kmol/K].
doublereal dot5(const V &x, const V &y)
Templated Inner product of two vectors of length 5.
R poly3(D x, R *c)
Templated evaluation of a polynomial of order 3.
R poly4(D x, R *c)
Evaluates a polynomial of order 4.
double polyfit(size_t n, size_t deg, const double *x, const double *y, const double *w, double *p)
Fits a polynomial function to a set of data points.
Contains declarations for string manipulation functions within Cantera.
Various templated functions that carry out common vector operations (see Templated Utility Functions)...