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