Cantera  4.0.0a1
Loading...
Searching...
No Matches
MultiTransport.cpp
Go to the documentation of this file.
1/**
2 * @file MultiTransport.cpp
3 * Implementation file for class MultiTransport
4 */
5
6// This file is part of Cantera. See License.txt in the top-level directory or
7// at https://cantera.org/license.txt for license and copyright information.
8
13
14using namespace std;
15
16namespace Cantera
17{
18
19///////////////////// helper functions /////////////////////////
20
21/**
22 * The Parker temperature correction to the rotational collision number.
23 *
24 * @param tr Reduced temperature @f$ \epsilon/kT @f$
25 * @param sqtr square root of tr.
26 */
27double Frot(double tr, double sqtr)
28{
29 const double c1 = 0.5*sqrt(Pi)*Pi;
30 const double c2 = 0.25*Pi*Pi + 2.0;
31 const double c3 = sqrt(Pi)*Pi;
32 return 1.0 + c1*sqtr + c2*tr + c3*sqtr*tr;
33}
34
35//////////////////// class MultiTransport methods //////////////
36
37void MultiTransport::init(shared_ptr<ThermoPhase> thermo, int mode)
38{
40
41 // the L matrix
42 m_Lmatrix.resize(3*m_nsp, 3*m_nsp);
43 m_a.resize(3*m_nsp, 1.0);
44 m_b.resize(3*m_nsp, 0.0);
45 m_aa.resize(m_nsp, m_nsp, 0.0);
46 m_molefracs_last.resize(m_nsp, -1.0);
47 m_frot_298.resize(m_nsp);
48 m_rotrelax.resize(m_nsp);
49 m_cinternal.resize(m_nsp);
53
54 // set flags all false
55 m_l0000_ok = false;
56 m_lmatrix_soln_ok = false;
57 m_thermal_tlast = Undef;
58
59 // some work space
60 m_spwork1.resize(m_nsp);
61 m_spwork2.resize(m_nsp);
62 m_spwork3.resize(m_nsp);
63
64 // precompute and store log(epsilon_ij/k_B)
65 m_log_eps_k.resize(m_nsp, m_nsp);
66 for (size_t i = 0; i < m_nsp; i++) {
67 for (size_t j = i; j < m_nsp; j++) {
68 m_log_eps_k(i,j) = log(m_epsilon(i,j)/Boltzmann);
69 m_log_eps_k(j,i) = m_log_eps_k(i,j);
70 }
71 }
72
73 // precompute and store constant parts of the Parker rotational
74 // collision number temperature correction
75 const double sq298 = sqrt(298.0);
76 const double kb298 = Boltzmann * 298.0;
77 m_sqrt_eps_k.resize(m_nsp);
78 for (size_t k = 0; k < m_nsp; k++) {
79 m_sqrt_eps_k[k] = sqrt(m_eps[k]/Boltzmann);
80 m_frot_298[k] = Frot(m_eps[k]/kb298, m_sqrt_eps_k[k]/sq298);
81 }
82}
83
85{
87 m_thermal_tlast = Undef;
88 m_l0000_ok = false;
89 m_lmatrix_soln_ok = false;
90 m_molefracs_last[0] += 1.23e-4;
91}
92
94{
95 solveLMatrixEquation();
96 double sum = 0.0;
97 for (size_t k = 0; k < 2*m_nsp; k++) {
98 sum += m_b[k + m_nsp] * m_a[k + m_nsp];
99 }
100 return -4.0*sum;
101}
102
104{
105 checkArraySize("MultiTransport::getThermalDiffCoeffs", dt.size(), m_nsp);
106 solveLMatrixEquation();
107 const double c = 1.6/GasConstant;
108 for (size_t k = 0; k < m_nsp; k++) {
109 dt[k] = c * m_mw[k] * m_molefracs[k] * m_a[k];
110 }
111}
112
113double MultiTransport::pressure_ig()
114{
115 return m_thermo->molarDensity() * GasConstant * m_thermo->temperature();
116}
117
118void MultiTransport::solveLMatrixEquation()
119{
120 // if T has changed, update the temperature-dependent properties.
122 update_C();
123 if (m_lmatrix_soln_ok) {
124 return;
125 }
126
127 // Copy the mole fractions twice into the last two blocks of the right-hand-
128 // side vector m_b. The first block of m_b was set to zero when it was
129 // created, and is not modified so doesn't need to be reset to zero.
130 for (size_t k = 0; k < m_nsp; k++) {
131 m_b[k] = 0.0;
132 m_b[k + m_nsp] = m_molefracs[k];
133 m_b[k + 2*m_nsp] = m_molefracs[k];
134 }
135
136 // Set the right-hand side vector to zero in the 3rd block for all species
137 // with no internal energy modes. The corresponding third-block rows and
138 // columns will be set to zero, except on the diagonal of L01,01, where they
139 // are set to 1.0. This has the effect of eliminating these equations from
140 // the system, since the equation becomes: m_a[2*m_nsp + k] = 0.0.
141
142 // Note that this differs from the Chemkin procedure, where all *monatomic*
143 // species are excluded. Since monatomic radicals can have non-zero internal
144 // heat capacities due to electronic excitation, they should be retained.
145 for (size_t k = 0; k < m_nsp; k++) {
146 if (!hasInternalModes(k)) {
147 m_b[2*m_nsp + k] = 0.0;
148 }
149 }
150
151 // evaluate the submatrices of the L matrix
152 m_Lmatrix.resize(3*m_nsp, 3*m_nsp, 0.0);
153
154 //! Evaluate the upper-left block of the L matrix.
157 eval_L0001();
158 eval_L1000();
159 eval_L1010(m_molefracs);
160 eval_L1001(m_molefracs);
161 eval_L0100();
162 eval_L0110();
163 eval_L0101(m_molefracs);
164
165 // Solve it using GMRES or LU decomposition. The last solution in m_a should
166 // provide a good starting guess, so convergence should be fast.
167 m_a = m_b;
168 solve(m_Lmatrix, m_a);
169 m_lmatrix_soln_ok = true;
171 // L matrix is overwritten with LU decomposition
172 m_l0000_ok = false;
173}
174
175void MultiTransport::getSpeciesFluxes(size_t ndim, span<const double> grad_T,
176 size_t ldx, span<const double> grad_X,
177 size_t ldf, span<double> fluxes)
178{
179 if (ldx < m_nsp) {
180 throw CanteraError("MultiTransport::getSpeciesFluxes", "ldx is too small");
181 }
182 if (ldf < m_nsp) {
183 throw CanteraError("MultiTransport::getSpeciesFluxes", "ldf is too small");
184 }
185 checkArraySize("MultiTransport::getSpeciesFluxes: grad_T", grad_T.size(), ndim);
186 checkArraySize("MultiTransport::getSpeciesFluxes: grad_X", grad_X.size(), ldx * ndim);
187 checkArraySize("MultiTransport::getSpeciesFluxes: fluxes", fluxes.size(), ldf * ndim);
188 // update the binary diffusion coefficients if necessary
189 update_T();
190 updateDiff_T();
191
192 // If any component of grad_T is non-zero, then get the
193 // thermal diffusion coefficients
194 bool addThermalDiffusion = false;
195 for (size_t i = 0; i < ndim; i++) {
196 if (grad_T[i] != 0.0) {
197 addThermalDiffusion = true;
198 }
199 }
200 if (addThermalDiffusion) {
202 }
203
204 auto y = m_thermo->massFractions();
205 double rho = m_thermo->density();
206
207 for (size_t i = 0; i < m_nsp; i++) {
208 double sum = 0.0;
209 for (size_t j = 0; j < m_nsp; j++) {
210 m_aa(i,j) = m_molefracs[j]*m_molefracs[i]/m_bdiff(i,j);
211 sum += m_aa(i,j);
212 }
213 m_aa(i,i) -= sum;
214 }
215
216 // enforce the condition \sum Y_k V_k = 0. This is done by replacing
217 // the flux equation with the largest gradx component in the first
218 // coordinate direction with the flux balance condition.
219 size_t jmax = 0;
220 double gradmax = -1.0;
221 for (size_t j = 0; j < m_nsp; j++) {
222 if (fabs(grad_X[j]) > gradmax) {
223 gradmax = fabs(grad_X[j]);
224 jmax = j;
225 }
226 }
227
228 // set the matrix elements in this row to the mass fractions,
229 // and set the entry in gradx to zero
230 for (size_t j = 0; j < m_nsp; j++) {
231 m_aa(jmax,j) = y[j];
232 }
233 vector<double> gsave(ndim), grx(ldx*m_nsp);
234 for (size_t n = 0; n < ldx*ndim; n++) {
235 grx[n] = grad_X[n];
236 }
237
238 // copy grad_X to fluxes
239 for (size_t n = 0; n < ndim; n++) {
240 auto gx = grad_X.subspan(ldx * n, m_nsp);
241 auto fx = fluxes.subspan(ldf * n, m_nsp);
242 std::copy(gx.begin(), gx.end(), fx.begin());
243 fx[jmax] = 0.0;
244 }
245
246 // solve the equations
247 solve(m_aa, fluxes, ndim, ldf);
248 double pp = pressure_ig();
249
250 // multiply diffusion velocities by rho * V to create mass fluxes, and
251 // restore the gradx elements that were modified
252 for (size_t n = 0; n < ndim; n++) {
253 size_t offset = n*ldf;
254 for (size_t i = 0; i < m_nsp; i++) {
255 fluxes[i + offset] *= rho * y[i] / pp;
256 }
257 }
258
259 // thermal diffusion
260 if (addThermalDiffusion) {
261 for (size_t n = 0; n < ndim; n++) {
262 size_t offset = n*ldf;
263 double grad_logt = grad_T[n]/m_temp;
264 for (size_t i = 0; i < m_nsp; i++) {
265 fluxes[i + offset] -= m_spwork[i]*grad_logt;
266 }
267 }
268 }
269}
270
271void MultiTransport::getMassFluxes(span<const double> state1,
272 span<const double> state2,
273 double delta, span<double> fluxes)
274{
275 checkArraySize("MultiTransport::getMassFluxes: state1", state1.size(), m_nsp + 2);
276 checkArraySize("MultiTransport::getMassFluxes: state2", state2.size(), m_nsp + 2);
277 checkArraySize("MultiTransport::getMassFluxes: fluxes", fluxes.size(), m_nsp);
278 span<double> x1(m_spwork1);
279 span<double> x2(m_spwork2);
280 span<double> x3(m_spwork3);
281 size_t nsp = m_thermo->nSpecies();
282 m_thermo->restoreState(state1);
283 double p1 = m_thermo->pressure();
284 double t1 = state1[0];
285 m_thermo->getMoleFractions(x1);
286
287 m_thermo->restoreState(state2);
288 double p2 = m_thermo->pressure();
289 double t2 = state2[0];
290 m_thermo->getMoleFractions(x2);
291
292 double p = 0.5*(p1 + p2);
293 double t = 0.5*(state1[0] + state2[0]);
294
295 for (size_t n = 0; n < nsp; n++) {
296 x3[n] = 0.5*(x1[n] + x2[n]);
297 }
298 m_thermo->setState_TPX(t, p, x3);
299 m_thermo->getMoleFractions(m_molefracs);
300
301 // update the binary diffusion coefficients if necessary
302 update_T();
303 updateDiff_T();
304
305 // If there is a temperature gradient, then get the
306 // thermal diffusion coefficients
307 bool addThermalDiffusion = false;
308 if (state1[0] != state2[0]) {
309 addThermalDiffusion = true;
311 }
312
313 auto y = m_thermo->massFractions();
314 double rho = m_thermo->density();
315 for (size_t i = 0; i < m_nsp; i++) {
316 double sum = 0.0;
317 for (size_t j = 0; j < m_nsp; j++) {
318 m_aa(i,j) = m_molefracs[j]*m_molefracs[i]/m_bdiff(i,j);
319 sum += m_aa(i,j);
320 }
321 m_aa(i,i) -= sum;
322 }
323
324 // enforce the condition \sum Y_k V_k = 0. This is done by replacing the
325 // flux equation with the largest gradx component with the flux balance
326 // condition.
327 size_t jmax = 0;
328 double gradmax = -1.0;
329 for (size_t j = 0; j < m_nsp; j++) {
330 if (fabs(x2[j] - x1[j]) > gradmax) {
331 gradmax = fabs(x1[j] - x2[j]) / delta;
332 jmax = j;
333 }
334 }
335
336 // set the matrix elements in this row to the mass fractions,
337 // and set the entry in gradx to zero
338 for (size_t j = 0; j < m_nsp; j++) {
339 m_aa(jmax,j) = y[j];
340 fluxes[j] = (x2[j] - x1[j]) / delta;
341 }
342 fluxes[jmax] = 0.0;
343
344 // Solve the equations
345 solve(m_aa, fluxes);
346
347 double pp = pressure_ig();
348 // multiply diffusion velocities by rho * Y_k to create
349 // mass fluxes, and divide by pressure
350 for (size_t i = 0; i < m_nsp; i++) {
351 fluxes[i] *= rho * y[i] / pp;
352 }
353
354 // thermal diffusion
355 if (addThermalDiffusion) {
356 double grad_logt = (t2 - t1) / m_temp / delta;
357 for (size_t i = 0; i < m_nsp; i++) {
358 fluxes[i] -= m_spwork[i]*grad_logt;
359 }
360 }
361}
362
363void MultiTransport::getMolarFluxes(span<const double> state1,
364 span<const double> state2, const double delta, span<double> fluxes)
365{
366 getMassFluxes(state1, state2, delta, fluxes);
367 for (size_t k = 0; k < m_thermo->nSpecies(); k++) {
368 fluxes[k] /= m_mw[k];
369 }
370}
371
372void MultiTransport::getMultiDiffCoeffs(const size_t ld, span<double> d)
373{
374 if (ld < m_nsp) {
375 throw CanteraError("MultiTransport::getMultiDiffCoeffs", "ld is too small");
376 }
377 checkArraySize("MultiTransport::getMultiDiffCoeffs", d.size(), ld * m_nsp);
378 double p = pressure_ig();
379
380 // update the mole fractions
381 update_C();
382
383 // update the binary diffusion coefficients
384 update_T();
386
387 // evaluate L0000 if the temperature or concentrations have
388 // changed since it was last evaluated.
389 if (!m_l0000_ok) {
391 }
392
393 // invert L00,00
394 invert(m_Lmatrix, m_nsp);
395 m_l0000_ok = false; // matrix is overwritten by inverse
396 m_lmatrix_soln_ok = false;
397
398 double prefactor = 16.0 * m_temp
399 * m_thermo->meanMolecularWeight()/(25.0 * p);
400 for (size_t i = 0; i < m_nsp; i++) {
401 for (size_t j = 0; j < m_nsp; j++) {
402 double c = prefactor/m_mw[j];
403 d[ld*j + i] = c*m_molefracs[i]*
404 (m_Lmatrix(i,j) - m_Lmatrix(i,i));
405 }
406 }
407
408}
409
411{
412 if (m_temp == m_thermo->temperature() && m_nsp == m_thermo->nSpecies()) {
413 return;
414 }
415 GasTransport::update_T();
416 // temperature has changed, so polynomial fits will need to be
417 // redone, and the L matrix reevaluated.
418 m_lmatrix_soln_ok = false;
419 m_l0000_ok = false;
420}
421
423{
424 // Update the local mole fraction array
425 m_thermo->getMoleFractions(m_molefracs);
426
427 for (size_t k = 0; k < m_nsp; k++) {
428 // add an offset to avoid a pure species condition
429 m_molefracs[k] = std::max(Tiny, m_molefracs[k]);
430 if (m_molefracs[k] != m_molefracs_last[k]) {
431 // If any mole fractions have changed, signal that concentration-
432 // dependent quantities will need to be recomputed before use.
433 m_l0000_ok = false;
434 m_lmatrix_soln_ok = false;
435 }
436 }
437}
438
440{
441 if (m_thermal_tlast == m_thermo->temperature()) {
442 return;
443 }
444 // we need species viscosities and binary diffusion coefficients
446 updateDiff_T();
447
448 // evaluate polynomial fits for A*, B*, C*
449 for (size_t i = 0; i < m_nsp; i++) {
450 for (size_t j = i; j < m_nsp; j++) {
451 double z = m_star_poly_uses_actualT[i][j] == 1 ? m_logt : m_logt - m_log_eps_k(i,j);
452 int ipoly = m_poly[i][j];
453 if (m_mode == CK_Mode) {
454 m_astar(i,j) = poly6(z, m_astar_poly[ipoly]);
455 m_bstar(i,j) = poly6(z, m_bstar_poly[ipoly]);
456 m_cstar(i,j) = poly6(z, m_cstar_poly[ipoly]);
457 } else {
458 m_astar(i,j) = poly8(z, m_astar_poly[ipoly]);
459 m_bstar(i,j) = poly8(z, m_bstar_poly[ipoly]);
460 m_cstar(i,j) = poly8(z, m_cstar_poly[ipoly]);
461 }
462 m_astar(j,i) = m_astar(i,j);
463 m_bstar(j,i) = m_bstar(i,j);
464 m_cstar(j,i) = m_cstar(i,j);
465 }
466 }
467
468 // evaluate the temperature-dependent rotational relaxation rate
469 for (size_t k = 0; k < m_nsp; k++) {
470 double tr = m_eps[k]/ m_kbt;
471 double sqtr = m_sqrt_eps_k[k] / m_sqrt_t;
472 m_rotrelax[k] = std::max(1.0,m_zrot[k]) * m_frot_298[k]/Frot(tr, sqtr);
473 }
474
475 double c = 1.2*GasConstant*m_temp;
476 for (size_t k = 0; k < m_nsp; k++) {
477 m_bdiff(k,k) = c * m_visc[k] * m_astar(k,k)/m_mw[k];
478 }
479
480 // Calculate the internal heat capacities by subtracting off the translational contributions
481 /*
482 * HKM Exploratory comment:
483 * The translational component is 1.5
484 * The rotational component is 1.0 for a linear molecule and 1.5 for a nonlinear molecule
485 * and zero for a monatomic.
486 * Chemkin has traditionally subtracted 1.5 here (SAND86-8246).
487 * The original Dixon-Lewis paper subtracted 1.5 here.
488 */
489 vector<double> cp(m_thermo->nSpecies());
490 m_thermo->getCp_R_ref(cp);
491 for (size_t k = 0; k < m_nsp; k++) {
492 m_cinternal[k] = cp[k] - 2.5;
493 }
494 m_thermal_tlast = m_thermo->temperature();
495}
496
497//! Constant to compare dimensionless heat capacities against zero
498static const double Min_C_Internal = 0.001;
499
500bool MultiTransport::hasInternalModes(size_t j)
501{
502 return (m_cinternal[j] > Min_C_Internal);
503}
504
505void MultiTransport::eval_L0000(span<const double> x)
506{
507 double prefactor = 16.0*m_temp/25.0;
508 double sum;
509 for (size_t i = 0; i < m_nsp; i++) {
510 // subtract-off the k=i term to account for the first delta
511 // function in Eq. (12.121)
512 sum = -x[i]/m_bdiff(i,i);
513 for (size_t k = 0; k < m_nsp; k++) {
514 sum += x[k]/m_bdiff(i,k);
515 }
516
517 sum /= m_mw[i];
518 for (size_t j = 0; j != m_nsp; ++j) {
519 m_Lmatrix(i,j) = prefactor * x[j]
520 * (m_mw[j] * sum + x[i]/m_bdiff(i,j));
521 }
522 // diagonal term is zero
523 m_Lmatrix(i,i) = 0.0;
524 }
525}
526
527void MultiTransport::eval_L0010(span<const double> x)
528{
529 double prefactor = 1.6*m_temp;
530 for (size_t j = 0; j < m_nsp; j++) {
531 double xj = x[j];
532 double wj = m_mw[j];
533 double sum = 0.0;
534 for (size_t i = 0; i < m_nsp; i++) {
535 m_Lmatrix(i,j + m_nsp) = - prefactor * x[i] * xj * m_mw[i] *
536 (1.2 * m_cstar(j,i) - 1.0) /
537 ((wj + m_mw[i]) * m_bdiff(j,i));
538
539 // the next term is independent of "j";
540 // need to do it for the "j,j" term
541 sum -= m_Lmatrix(i,j+m_nsp);
542 }
543 m_Lmatrix(j,j+m_nsp) += sum;
544 }
545}
546
548{
549 for (size_t j = 0; j < m_nsp; j++) {
550 for (size_t i = 0; i < m_nsp; i++) {
551 m_Lmatrix(i+m_nsp,j) = m_Lmatrix(j,i+m_nsp);
552 }
553 }
554}
555
556void MultiTransport::eval_L1010(span<const double> x)
557{
558 const double fiveover3pi = 5.0/(3.0*Pi);
559 double prefactor = (16.0*m_temp)/25.0;
560
561 for (size_t j = 0; j < m_nsp; j++) {
562 // get constant terms that depend on just species "j"
563 double constant1 = prefactor*x[j];
564 double wjsq = m_mw[j]*m_mw[j];
565 double constant2 = 13.75*wjsq;
566 double constant3 = m_crot[j]/m_rotrelax[j];
567 double constant4 = 7.5*wjsq;
568 double fourmj = 4.0*m_mw[j];
569 double threemjsq = 3.0*m_mw[j]*m_mw[j];
570 double sum = 0.0;
571 for (size_t i = 0; i < m_nsp; i++) {
572 double sumwij = m_mw[i] + m_mw[j];
573 double term1 = m_bdiff(i,j) * sumwij*sumwij;
574 double term2 = fourmj*m_astar(i,j)*(1.0 + fiveover3pi*
575 (constant3 + (m_crot[i]/m_rotrelax[i]))); // see Eq. (12.125)
576
577 m_Lmatrix(i+m_nsp,j+m_nsp) = constant1*x[i]*m_mw[i] /(m_mw[j]*term1) *
578 (constant2 - threemjsq*m_bstar(i,j)
579 - term2*m_mw[j]);
580
581 sum += x[i] /(term1) *
582 (constant4 + m_mw[i]*m_mw[i]*
583 (6.25 - 3.0*m_bstar(i,j)) + term2*m_mw[i]);
584 }
585
586 m_Lmatrix(j+m_nsp,j+m_nsp) -= sum*constant1;
587 }
588}
589
590void MultiTransport::eval_L1001(span<const double> x)
591{
592 double prefactor = 32.00*m_temp/(5.00*Pi);
593 for (size_t j = 0; j < m_nsp; j++) {
594 // collect terms that depend only on "j"
595 if (hasInternalModes(j)) {
596 double constant = prefactor*m_mw[j]*x[j]*m_crot[j]/(m_cinternal[j]*m_rotrelax[j]);
597 double sum = 0.0;
598 for (size_t i = 0; i < m_nsp; i++) {
599 // see Eq. (12.127)
600 m_Lmatrix(i+m_nsp,j+2*m_nsp) = constant * m_astar(j,i) * x[i] /
601 ((m_mw[j] + m_mw[i]) * m_bdiff(j,i));
602 sum += m_Lmatrix(i+m_nsp,j+2*m_nsp);
603 }
604 m_Lmatrix(j+m_nsp,j+2*m_nsp) += sum;
605 } else {
606 for (size_t i = 0; i < m_nsp; i++) {
607 m_Lmatrix(i+m_nsp,j+2*m_nsp) = 0.0;
608 }
609 }
610 }
611}
612
613void MultiTransport::eval_L0001()
614{
615 for (size_t j = 0; j < m_nsp; j++) {
616 for (size_t i = 0; i < m_nsp; i++) {
617 m_Lmatrix(i,j+2*m_nsp) = 0.0;
618 }
619 }
620}
621
622void MultiTransport::eval_L0100()
623{
624 for (size_t j = 0; j < m_nsp; j++) {
625 for (size_t i = 0; i < m_nsp; i++) {
626 m_Lmatrix(i+2*m_nsp,j) = 0.0; // see Eq. (12.123)
627 }
628 }
629}
630
631void MultiTransport::eval_L0110()
632{
633 for (size_t j = 0; j < m_nsp; j++) {
634 for (size_t i = 0; i < m_nsp; i++) {
635 m_Lmatrix(i+2*m_nsp,j+m_nsp) = m_Lmatrix(j+m_nsp,i+2*m_nsp); // see Eq. (12.123)
636 }
637 }
638}
639
640void MultiTransport::eval_L0101(span<const double> x)
641{
642 for (size_t i = 0; i < m_nsp; i++) {
643 if (hasInternalModes(i)) {
644 // collect terms that depend only on "i"
645 double constant1 = 4*m_temp*x[i]/m_cinternal[i];
646 double constant2 = 12*m_mw[i]*m_crot[i] /
647 (5*Pi*m_cinternal[i]*m_rotrelax[i]);
648 double sum = 0.0;
649 for (size_t k = 0; k < m_nsp; k++) {
650 // see Eq. (12.131)
651 double diff_int = m_bdiff(i,k);
652 m_Lmatrix(k+2*m_nsp,i+2*m_nsp) = 0.0;
653 sum += x[k]/diff_int;
654 if (k != i) {
655 sum += x[k]*m_astar(i,k)*constant2 / (m_mw[k]*diff_int);
656 }
657 }
658 // see Eq. (12.130)
659 m_Lmatrix(i+2*m_nsp,i+2*m_nsp) =
660 - 8/Pi*m_mw[i]*x[i]*x[i]*m_crot[i] /
661 (m_cinternal[i]*m_cinternal[i]*GasConstant*m_visc[i]*m_rotrelax[i])
662 - constant1*sum;
663 } else {
664 for (size_t k = 0; k < m_nsp; k++) {
665 m_Lmatrix(i+2*m_nsp,i+2*m_nsp) = 1.0;
666 }
667 }
668 }
669}
670
671}
ThermoPhase object for the ideal gas equation of state - workhorse for Cantera (see Thermodynamic Pro...
Interface for class MultiTransport.
Base class for exceptions thrown by Cantera classes.
void resize(size_t n, size_t m, double v=0.0) override
Resize the matrix.
vector< double > m_mw
Local copy of the species molecular weights.
vector< double > m_molefracs
Vector of species mole fractions.
double m_temp
Current value of the temperature [K] at which the properties in this object are calculated.
vector< double > m_eps
Lennard-Jones well-depth [J] of the species in the current phase.
virtual void updateDiff_T()
Update the binary diffusion coefficients.
DenseMatrix m_epsilon
The effective well depth [J] for (i,j) collisions.
vector< double > m_spwork
work space length = m_nsp
vector< double > m_zrot
Rotational relaxation number for each species.
int m_mode
Type of the polynomial fits to temperature.
double m_logt
Current value of the log of the temperature.
vector< double > m_visc
vector of species viscosities [Pa·s].
virtual void updateSpeciesViscosities()
Update the pure-species viscosities.
double m_sqrt_t
current value of temperature to 1/2 power
DenseMatrix m_bdiff
Matrix of binary diffusion coefficients at the reference pressure and the current temperature Size is...
vector< double > m_crot
Dimensionless rotational heat capacity of each species.
void invalidateCache() override
Invalidate any cached values which are normally updated only when a change in state is detected.
double m_kbt
Current value of Boltzmann constant times the temperature [J].
vector< vector< int > > m_star_poly_uses_actualT
Flag to indicate for which (i,j) interaction pairs the actual temperature is used instead of the redu...
vector< vector< double > > m_bstar_poly
Fit for bstar collision integral.
vector< vector< double > > m_astar_poly
Fit for astar collision integral.
vector< vector< int > > m_poly
Indices for the (i,j) interaction in collision integral fits.
void init(shared_ptr< ThermoPhase > thermo, int mode=0) override
Initialize a transport manager.
vector< vector< double > > m_cstar_poly
Fit for cstar collision integral.
bool m_l0000_ok
Boolean indicating viscosity is up to date.
void update_T() override
Update basic temperature-dependent quantities if the temperature has changed.
double thermalConductivity() override
Get the mixture thermal conductivity [W/m/K].
DenseMatrix m_astar
Dense matrix for astar.
void updateThermal_T()
Update the temperature-dependent terms needed to compute the thermal conductivity and thermal diffusi...
void eval_L0010(span< const double > x)
Evaluate the L0010 matrices.
vector< double > m_molefracs_last
Mole fraction vector from last L-matrix evaluation.
DenseMatrix m_cstar
Dense matrix for cstar.
void update_C() override
Update basic concentration-dependent quantities if the concentrations have changed.
void invalidateCache() override
Invalidate any cached values which are normally updated only when a change in state is detected.
void getThermalDiffCoeffs(span< double > dt) override
Return the thermal diffusion coefficients [kg/m/s].
void getMolarFluxes(span< const double > state1, span< const double > state2, const double delta, span< double > fluxes) override
Get the molar diffusional fluxes [kmol/m²/s] of the species, given the thermodynamic state at two nea...
void init(shared_ptr< ThermoPhase > thermo, int mode=0) override
Initialize a transport manager.
void getSpeciesFluxes(size_t ndim, span< const double > grad_T, size_t ldx, span< const double > grad_X, size_t ldf, span< double > fluxes) override
Get the species diffusive mass fluxes [kg/m²/s] with respect to the mass averaged velocity,...
void eval_L0000(span< const double > x)
Evaluate the L0000 matrices.
void getMassFluxes(span< const double > state1, span< const double > state2, double delta, span< double > fluxes) override
Get the mass diffusional fluxes [kg/m²/s] of the species, given the thermodynamic state at two nearby...
DenseMatrix m_bstar
Dense matrix for bstar.
void eval_L1000()
Evaluate the L1000 matrices.
void getMultiDiffCoeffs(const size_t ld, span< double > d) override
Return the multicomponent diffusion coefficients [m²/s].
shared_ptr< ThermoPhase > m_thermo
pointer to the object representing the phase
Definition Transport.h:432
size_t m_nsp
Number of species in the phase.
Definition Transport.h:435
ThermoPhase & thermo()
Phase object.
Definition Transport.h:111
auto poly8(D x, const C &c)
Templated evaluation of a polynomial of order 8.
Definition utilities.h:147
auto poly6(D x, const C &c)
Templated evaluation of a polynomial of order 6.
Definition utilities.h:131
const double Boltzmann
Boltzmann constant [J/K].
Definition ct_defs.h:87
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition ct_defs.h:123
const double Pi
Pi.
Definition ct_defs.h:71
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
const double Tiny
Small number to compare differences of mole fractions against.
Definition ct_defs.h:176
const double Undef
Fairly random number to be used to initialize variables against to see if they are subsequently defin...
Definition ct_defs.h:167
void solve(DenseMatrix &A, span< double > b, size_t nrhs, size_t ldb)
Solve Ax = b. Array b is overwritten on exit with x.
offset
Offsets of solution components in the 1D solution array.
Definition Flow1D.h:25
static const double Min_C_Internal
Constant to compare dimensionless heat capacities against zero.
void invert(DenseMatrix &A, size_t nn)
invert A. A is overwritten with A^-1.
double Frot(double tr, double sqtr)
The Parker temperature correction to the rotational collision number.
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...