Cantera  2.1.2
MultiTransport.cpp
Go to the documentation of this file.
1 /**
2  * @file MultiTransport.cpp
3  * Implementation file for class MultiTransport
4  */
5 /*
6  * Copyright 2001 California Institute of Technology
7  * See file License.txt for licensing information
8  */
9 
12 #include "cantera/base/utilities.h"
13 #include "L_matrix.h"
18 
19 using namespace std;
20 
21 namespace Cantera
22 {
23 
24 ///////////////////// helper functions /////////////////////////
25 
26 /**
27  * The Parker temperature correction to the rotational collision number.
28  *
29  * @param tr Reduced temperature \f$ \epsilon/kT \f$
30  * @param sqtr square root of tr.
31  */
32 inline doublereal Frot(doublereal tr, doublereal sqtr)
33 {
34  const doublereal c1 = 0.5*SqrtPi*Pi;
35  const doublereal c2 = 0.25*Pi*Pi + 2.0;
36  const doublereal c3 = SqrtPi*Pi;
37  return 1.0 + c1*sqtr + c2*tr + c3*sqtr*tr;
38 }
39 
40 //////////////////// class MultiTransport methods //////////////
41 
42 MultiTransport::MultiTransport(thermo_t* thermo)
43  : GasTransport(thermo)
44 {
45 }
46 
48 {
50 
51  // copy polynomials and parameters into local storage
52  m_poly = tr.poly;
53  m_astar_poly = tr.astar_poly;
54  m_bstar_poly = tr.bstar_poly;
55  m_cstar_poly = tr.cstar_poly;
56  m_om22_poly = tr.omega22_poly;
57  m_zrot = tr.zrot;
58  m_crot = tr.crot;
59  m_eps = tr.eps;
60  m_sigma = tr.sigma;
61  m_alpha = tr.alpha;
62  m_dipole = tr.dipole;
63  m_zrot = tr.zrot;
64 
65  // the L matrix
66  m_Lmatrix.resize(3*m_nsp, 3*m_nsp);
67  m_a.resize(3*m_nsp, 1.0);
68  m_b.resize(3*m_nsp, 0.0);
69  m_aa.resize(m_nsp, m_nsp, 0.0);
70  m_molefracs_last.resize(m_nsp, -1.0);
71 
72  m_frot_298.resize(m_nsp);
73  m_rotrelax.resize(m_nsp);
74 
75  m_cinternal.resize(m_nsp);
76 
81 
82  // set flags all false
83  m_abc_ok = false;
84  m_l0000_ok = false;
85  m_lmatrix_soln_ok = false;
86 
87  m_thermal_tlast = 0.0;
88 
89  // some work space
90  m_spwork1.resize(m_nsp);
91  m_spwork2.resize(m_nsp);
92  m_spwork3.resize(m_nsp);
93 
94  // precompute and store log(epsilon_ij/k_B)
95  m_log_eps_k.resize(m_nsp, m_nsp);
96  // int j;
97  for (size_t i = 0; i < m_nsp; i++) {
98  for (size_t j = i; j < m_nsp; j++) {
99  m_log_eps_k(i,j) = log(tr.epsilon(i,j)/Boltzmann);
100  m_log_eps_k(j,i) = m_log_eps_k(i,j);
101  }
102  }
103 
104  // precompute and store constant parts of the Parker rotational
105  // collision number temperature correction
106  const doublereal sq298 = sqrt(298.0);
107  const doublereal kb298 = Boltzmann * 298.0;
108  m_sqrt_eps_k.resize(m_nsp);
109  for (size_t k = 0; k < m_nsp; k++) {
110  m_sqrt_eps_k[k] = sqrt(tr.eps[k]/Boltzmann);
111  m_frot_298[k] = Frot(tr.eps[k]/kb298,
112  m_sqrt_eps_k[k]/sq298);
113  }
114 
115  return true;
116 }
117 
119 {
120  solveLMatrixEquation();
121  doublereal sum = 0.0;
122  for (size_t k = 0; k < 2*m_nsp; k++) {
123  sum += m_b[k + m_nsp] * m_a[k + m_nsp];
124  }
125  return -4.0*sum;
126 }
127 
128 void MultiTransport::getThermalDiffCoeffs(doublereal* const dt)
129 {
130  solveLMatrixEquation();
131  const doublereal c = 1.6/GasConstant;
132  for (size_t k = 0; k < m_nsp; k++) {
133  dt[k] = c * m_mw[k] * m_molefracs[k] * m_a[k];
134  }
135 }
136 
137 void MultiTransport::solveLMatrixEquation()
138 {
139  // if T has changed, update the temperature-dependent properties.
140  updateThermal_T();
141  update_C();
142  if (m_lmatrix_soln_ok) {
143  return;
144  }
145 
146  // Copy the mole fractions twice into the last two blocks of
147  // the right-hand-side vector m_b. The first block of m_b was
148  // set to zero when it was created, and is not modified so
149  // doesn't need to be reset to zero.
150  for (size_t k = 0; k < m_nsp; k++) {
151  m_b[k] = 0.0;
152  m_b[k + m_nsp] = m_molefracs[k];
153  m_b[k + 2*m_nsp] = m_molefracs[k];
154  }
155 
156  // Set the right-hand side vector to zero in the 3rd block for
157  // all species with no internal energy modes. The
158  // corresponding third-block rows and columns will be set to
159  // zero, except on the diagonal of L01,01, where they are set
160  // to 1.0. This has the effect of eliminating these equations
161  // from the system, since the equation becomes: m_a[2*m_nsp +
162  // k] = 0.0.
163 
164  // Note that this differs from the Chemkin procedure, where
165  // all *monatomic* species are excluded. Since monatomic
166  // radicals can have non-zero internal heat capacities due to
167  // electronic excitation, they should be retained.
168 
169  for (size_t k = 0; k < m_nsp; k++) {
170  if (!hasInternalModes(k)) {
171  m_b[2*m_nsp + k] = 0.0;
172  }
173  }
174 
175  // evaluate the submatrices of the L matrix
176  m_Lmatrix.resize(3*m_nsp, 3*m_nsp, 0.0);
177 
178  //! Evaluate the upper-left block of the L matrix.
181  eval_L0001();
182  eval_L1000();
183  eval_L1010(DATA_PTR(m_molefracs));
184  eval_L1001(DATA_PTR(m_molefracs));
185  eval_L0100();
186  eval_L0110();
187  eval_L0101(DATA_PTR(m_molefracs));
188 
189  // Solve it using GMRES or LU decomposition. The last solution
190  // in m_a should provide a good starting guess, so convergence
191  // should be fast.
192 
193  //if (m_gmres) {
194  // gmres(m_mgmres, 3*m_nsp, m_Lmatrix, m_b.begin(),
195  // m_a.begin(), m_eps_gmres);
196  // m_lmatrix_soln_ok = true;
197  // m_l0000_ok = true; // L matrix not modified by GMRES
198  //}
199  //else {
200  copy(m_b.begin(), m_b.end(), m_a.begin());
201  try {
202  solve(m_Lmatrix, DATA_PTR(m_a));
203  } catch (CanteraError& err) {
204  err.save();
205  //if (info != 0) {
206  throw CanteraError("MultiTransport::solveLMatrixEquation",
207  "error in solving L matrix.");
208  }
209  m_lmatrix_soln_ok = true;
211  // L matrix is overwritten with LU decomposition
212  m_l0000_ok = false;
213 }
214 
215 void MultiTransport::getSpeciesFluxes(size_t ndim, const doublereal* const grad_T,
216  size_t ldx, const doublereal* const grad_X,
217  size_t ldf, doublereal* const fluxes)
218 {
219  // update the binary diffusion coefficients if necessary
220  update_T();
221  updateDiff_T();
222 
223  // If any component of grad_T is non-zero, then get the
224  // thermal diffusion coefficients
225  bool addThermalDiffusion = false;
226  for (size_t i = 0; i < ndim; i++) {
227  if (grad_T[i] != 0.0) {
228  addThermalDiffusion = true;
229  }
230  }
231  if (addThermalDiffusion) {
233  }
234 
235  const doublereal* y = m_thermo->massFractions();
236  doublereal rho = m_thermo->density();
237 
238  for (size_t i = 0; i < m_nsp; i++) {
239  double sum = 0.0;
240  for (size_t j = 0; j < m_nsp; j++) {
241  m_aa(i,j) = m_molefracs[j]*m_molefracs[i]/m_bdiff(i,j);
242  sum += m_aa(i,j);
243  }
244  m_aa(i,i) -= sum;
245  }
246 
247  // enforce the condition \sum Y_k V_k = 0. This is done by replacing
248  // the flux equation with the largest gradx component in the first
249  // coordinate direction with the flux balance condition.
250  size_t jmax = 0;
251  doublereal gradmax = -1.0;
252  for (size_t j = 0; j < m_nsp; j++) {
253  if (fabs(grad_X[j]) > gradmax) {
254  gradmax = fabs(grad_X[j]);
255  jmax = j;
256  }
257  }
258 
259  // set the matrix elements in this row to the mass fractions,
260  // and set the entry in gradx to zero
261  for (size_t j = 0; j < m_nsp; j++) {
262  m_aa(jmax,j) = y[j];
263  }
264  vector_fp gsave(ndim), grx(ldx*m_nsp);
265  for (size_t n = 0; n < ldx*ndim; n++) {
266  grx[n] = grad_X[n];
267  }
268  //for (n = 0; n < ndim; n++) {
269  // gsave[n] = grad_X[jmax + n*ldx]; // save the input mole frac gradient
270  //grad_X[jmax + n*ldx] = 0.0;
271  // grx[jmax + n*ldx] = 0.0;
272  // }
273 
274  // copy grad_X to fluxes
275  const doublereal* gx;
276  for (size_t n = 0; n < ndim; n++) {
277  gx = grad_X + ldx*n;
278  copy(gx, gx + m_nsp, fluxes + ldf*n);
279  fluxes[jmax + n*ldf] = 0.0;
280  }
281 
282  // use LAPACK to solve the equations
283  int info=0;
284  ct_dgetrf(static_cast<int>(m_aa.nRows()),
285  static_cast<int>(m_aa.nColumns()), m_aa.ptrColumn(0),
286  static_cast<int>(m_aa.nRows()),
287  &m_aa.ipiv()[0], info);
288  if (info == 0) {
289  ct_dgetrs(ctlapack::NoTranspose,
290  static_cast<int>(m_aa.nRows()), ndim,
291  m_aa.ptrColumn(0), static_cast<int>(m_aa.nRows()),
292  &m_aa.ipiv()[0], fluxes, ldf, info);
293  if (info != 0) {
294  info += 100;
295  }
296  } else
297  throw CanteraError("MultiTransport::getSpeciesFluxes",
298  "Error in DGETRF");
299  if (info > 50)
300  throw CanteraError("MultiTransport::getSpeciesFluxes",
301  "Error in DGETRS");
302 
303  size_t offset;
304  doublereal pp = pressure_ig();
305 
306  // multiply diffusion velocities by rho * V to create
307  // mass fluxes, and restore the gradx elements that were
308  // modified
309  for (size_t n = 0; n < ndim; n++) {
310  offset = n*ldf;
311  for (size_t i = 0; i < m_nsp; i++) {
312  fluxes[i + offset] *= rho * y[i] / pp;
313  }
314  //grad_X[jmax + n*ldx] = gsave[n];
315  }
316 
317  // thermal diffusion
318  if (addThermalDiffusion) {
319  for (size_t n = 0; n < ndim; n++) {
320  offset = n*ldf;
321  doublereal grad_logt = grad_T[n]/m_temp;
322  for (size_t i = 0; i < m_nsp; i++) {
323  fluxes[i + offset] -= m_spwork[i]*grad_logt;
324  }
325  }
326  }
327 }
328 
329 void MultiTransport::getMassFluxes(const doublereal* state1, const doublereal* state2, doublereal delta,
330  doublereal* fluxes)
331 {
332  double* x1 = DATA_PTR(m_spwork1);
333  double* x2 = DATA_PTR(m_spwork2);
334  double* x3 = DATA_PTR(m_spwork3);
335  size_t n, nsp = m_thermo->nSpecies();
336  m_thermo->restoreState(nsp+2, state1);
337  double p1 = m_thermo->pressure();
338  double t1 = state1[0];
340 
341  m_thermo->restoreState(nsp+2, state2);
342  double p2 = m_thermo->pressure();
343  double t2 = state2[0];
345 
346  double p = 0.5*(p1 + p2);
347  double t = 0.5*(state1[0] + state2[0]);
348 
349  for (n = 0; n < nsp; n++) {
350  x3[n] = 0.5*(x1[n] + x2[n]);
351  }
352  m_thermo->setState_TPX(t, p, x3);
354 
355  // update the binary diffusion coefficients if necessary
356  update_T();
357  updateDiff_T();
358 
359  // If there is a temperature gradient, then get the
360  // thermal diffusion coefficients
361 
362  bool addThermalDiffusion = false;
363  if (state1[0] != state2[0]) {
364  addThermalDiffusion = true;
366  }
367 
368  const doublereal* y = m_thermo->massFractions();
369  doublereal rho = m_thermo->density();
370 
371  for (size_t i = 0; i < m_nsp; i++) {
372  doublereal sum = 0.0;
373  for (size_t j = 0; j < m_nsp; j++) {
374  m_aa(i,j) = m_molefracs[j]*m_molefracs[i]/m_bdiff(i,j);
375  sum += m_aa(i,j);
376  }
377  m_aa(i,i) -= sum;
378  }
379 
380  // enforce the condition \sum Y_k V_k = 0. This is done by
381  // replacing the flux equation with the largest gradx
382  // component with the flux balance condition.
383  size_t jmax = 0;
384  doublereal gradmax = -1.0;
385  for (size_t j = 0; j < m_nsp; j++) {
386  if (fabs(x2[j] - x1[j]) > gradmax) {
387  gradmax = fabs(x1[j] - x2[j]);
388  jmax = j;
389  }
390  }
391 
392  // set the matrix elements in this row to the mass fractions,
393  // and set the entry in gradx to zero
394 
395  for (size_t j = 0; j < m_nsp; j++) {
396  m_aa(jmax,j) = y[j];
397  fluxes[j] = x2[j] - x1[j];
398  }
399  fluxes[jmax] = 0.0;
400 
401  // use LAPACK to solve the equations
402  int info=0;
403  size_t nr = m_aa.nRows();
404  size_t nc = m_aa.nColumns();
405 
406  ct_dgetrf(nr, nc, m_aa.ptrColumn(0), nr, &m_aa.ipiv()[0], info);
407  if (info == 0) {
408  int ndim = 1;
409  ct_dgetrs(ctlapack::NoTranspose, nr, ndim,
410  m_aa.ptrColumn(0), nr, &m_aa.ipiv()[0], fluxes, nr, info);
411  if (info != 0)
412  throw CanteraError("MultiTransport::getMassFluxes",
413  "Error in DGETRS. Info = "+int2str(info));
414  } else
415  throw CanteraError("MultiTransport::getMassFluxes",
416  "Error in DGETRF. Info = "+int2str(info));
417 
418  doublereal pp = pressure_ig();
419 
420  // multiply diffusion velocities by rho * Y_k to create
421  // mass fluxes, and divide by pressure
422  for (size_t i = 0; i < m_nsp; i++) {
423  fluxes[i] *= rho * y[i] / pp;
424  }
425 
426  // thermal diffusion
427  if (addThermalDiffusion) {
428  doublereal grad_logt = (t2 - t1)/m_temp;
429  for (size_t i = 0; i < m_nsp; i++) {
430  fluxes[i] -= m_spwork[i]*grad_logt;
431  }
432  }
433 }
434 
435 void MultiTransport::getMolarFluxes(const doublereal* const state1,
436  const doublereal* const state2,
437  const doublereal delta,
438  doublereal* const fluxes)
439 {
440  getMassFluxes(state1, state2, delta, fluxes);
441  for (size_t k = 0; k < m_thermo->nSpecies(); k++) {
442  fluxes[k] /= m_mw[k];
443  }
444 }
445 
446 void MultiTransport::getMultiDiffCoeffs(const size_t ld, doublereal* const d)
447 {
448  doublereal p = pressure_ig();
449 
450  // update the mole fractions
451  update_C();
452 
453  // update the binary diffusion coefficients
454  update_T();
455  updateThermal_T();
456 
457  // evaluate L0000 if the temperature or concentrations have
458  // changed since it was last evaluated.
459  if (!m_l0000_ok) {
461  }
462 
463  // invert L00,00
464  int ierr = invert(m_Lmatrix, m_nsp);
465  if (ierr != 0) {
466  throw CanteraError("MultiTransport::getMultiDiffCoeffs",
467  string(" invert returned ierr = ")+int2str(ierr));
468  }
469  m_l0000_ok = false; // matrix is overwritten by inverse
470  m_lmatrix_soln_ok = false;
471 
472  //doublereal pres = m_thermo->pressure();
473  doublereal prefactor = 16.0 * m_temp
474  * m_thermo->meanMolecularWeight()/(25.0 * p);
475  doublereal c;
476 
477  for (size_t i = 0; i < m_nsp; i++) {
478  for (size_t j = 0; j < m_nsp; j++) {
479  c = prefactor/m_mw[j];
480  d[ld*j + i] = c*m_molefracs[i]*
481  (m_Lmatrix(i,j) - m_Lmatrix(i,i));
482  }
483  }
484 }
485 
487 {
488  if (m_temp == m_thermo->temperature()) {
489  return;
490  }
491 
492  GasTransport::update_T();
493 
494  // temperature has changed, so polynomial fits will need to be
495  // redone, and the L matrix reevaluated.
496  m_abc_ok = false;
497  m_lmatrix_soln_ok = false;
498  m_l0000_ok = false;
499 }
500 
502 {
503  // Update the local mole fraction array
505 
506  for (size_t k = 0; k < m_nsp; k++) {
507  // add an offset to avoid a pure species condition
508  m_molefracs[k] = std::max(Tiny, m_molefracs[k]);
509  if (m_molefracs[k] != m_molefracs_last[k]) {
510  // If any mole fractions have changed, signal that concentration-
511  // dependent quantities will need to be recomputed before use.
512  m_l0000_ok = false;
513  m_lmatrix_soln_ok = false;
514  }
515  }
516 }
517 
519 {
520  if (m_thermal_tlast == m_thermo->temperature()) {
521  return;
522  }
523  // we need species viscosities and binary diffusion coefficients
525  updateDiff_T();
526 
527  // evaluate polynomial fits for A*, B*, C*
528  doublereal z;
529  int ipoly;
530  for (size_t i = 0; i < m_nsp; i++) {
531  for (size_t j = i; j < m_nsp; j++) {
532  z = m_logt - m_log_eps_k(i,j);
533  ipoly = m_poly[i][j];
534  if (m_mode == CK_Mode) {
535  m_om22(i,j) = poly6(z, DATA_PTR(m_om22_poly[ipoly]));
536  m_astar(i,j) = poly6(z, DATA_PTR(m_astar_poly[ipoly]));
537  m_bstar(i,j) = poly6(z, DATA_PTR(m_bstar_poly[ipoly]));
538  m_cstar(i,j) = poly6(z, DATA_PTR(m_cstar_poly[ipoly]));
539  } else {
540  m_om22(i,j) = poly8(z, DATA_PTR(m_om22_poly[ipoly]));
541  m_astar(i,j) = poly8(z, DATA_PTR(m_astar_poly[ipoly]));
542  m_bstar(i,j) = poly8(z, DATA_PTR(m_bstar_poly[ipoly]));
543  m_cstar(i,j) = poly8(z, DATA_PTR(m_cstar_poly[ipoly]));
544  }
545  m_om22(j,i) = m_om22(i,j);
546  m_astar(j,i) = m_astar(i,j);
547  m_bstar(j,i) = m_bstar(i,j);
548  m_cstar(j,i) = m_cstar(i,j);
549  }
550  }
551  m_abc_ok = true;
552 
553  // evaluate the temperature-dependent rotational relaxation rate
554  doublereal tr, sqtr;
555  for (size_t k = 0; k < m_nsp; k++) {
556  tr = m_eps[k]/ m_kbt;
557  sqtr = m_sqrt_eps_k[k] / m_sqrt_t;
558  m_rotrelax[k] = std::max(1.0,m_zrot[k]) * m_frot_298[k]/Frot(tr, sqtr);
559  }
560 
561  doublereal d;
562  doublereal c = 1.2*GasConstant*m_temp;
563  for (size_t k = 0; k < m_nsp; k++) {
564  d = c * m_visc[k] * m_astar(k,k)/m_mw[k];
565  m_bdiff(k,k) = d;
566  }
567 
568  // Calculate the internal heat capacities by subtracting off the translational contributions
569  /*
570  * HKM Exploratory comment:
571  * The translational component is 1.5
572  * The rotational component is 1.0 for a linear molecule and 1.5 for a nonlinear molecule
573  * and zero for a monatomic.
574  * Chemkin has traditionally subtracted 1.5 here (SAND86-8246).
575  * The original Dixon-Lewis paper subtracted 1.5 here.
576  */
577  const vector_fp& cp = ((IdealGasPhase*)m_thermo)->cp_R_ref();
578  for (size_t k = 0; k < m_nsp; k++) {
579  m_cinternal[k] = cp[k] - 2.5;
580  }
581 
582  // m_thermo->update_T(m_update_thermal_T);
583  m_thermal_tlast = m_thermo->temperature();
584 }
585 
586 }
virtual doublereal thermalConductivity()
Returns the mixture thermal conductivity in W/m/K.
void eval_L1000()
Evaluate the L1000 matrices.
Definition: L_matrix.h:69
size_t nRows() const
Number of rows.
Definition: Array.h:317
virtual void updateSpeciesViscosities()
Update the pure-species viscosities.
virtual doublereal density() const
Density (kg/m^3).
Definition: Phase.h:534
std::string int2str(const int n, const std::string &fmt)
Convert an int to a string using a format converter.
Definition: stringUtils.cpp:40
DenseMatrix m_bstar
Dense matrix for bstar.
int invert(DenseMatrix &A, size_t nn)
invert A. A is overwritten with A^-1.
void restoreState(const vector_fp &state)
Restore a state saved on a previous call to saveState.
Definition: Phase.cpp:289
This structure holds transport model parameters relevant to transport in ideal gases with a kinetic t...
vector_fp alpha
Polarizability of each species in the phase.
virtual void getThermalDiffCoeffs(doublereal *const dt)
Return the thermal diffusion coefficients (kg/m/s)
std::vector< vector_fp > omega22_poly
This is vector of vectors containing the astar fit.
Various templated functions that carry out common vector operations (see Templated Utility Functions)...
std::vector< vector_fp > astar_poly
This is vector of vectors containing the astar fit.
virtual void setState_TPX(doublereal t, doublereal p, const doublereal *x)
Set the temperature (K), pressure (Pa), and mole fractions.
Class IdealGasPhase represents low-density gases that obey the ideal gas equation of state...
Interface for class MultiTransport.
DenseMatrix m_astar
Dense matrix for astar.
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...
thermo_t * m_thermo
pointer to the object representing the phase
DenseMatrix m_cstar
Dense matrix for cstar.
void update_T()
Update basic temperature-dependent quantities if the temperature has changed.
vector_fp sigma
Lennard-Jones diameter of the species in the current phase.
int solve(DenseMatrix &A, double *b)
Solve Ax = b. Array b is overwritten on exit with x.
DenseMatrix m_bdiff
Matrix of binary diffusion coefficients at the reference pressure and the current temperature Size is...
Definition: GasTransport.h:243
const doublereal Pi
Pi.
Definition: ct_defs.h:51
doublereal * ptrColumn(size_t j)
Return a pointer to the top of column j, columns are contiguous in memory.
Definition: Array.h:363
void getMoleFractions(doublereal *const x) const
Get the species mole fraction vector.
Definition: Phase.cpp:519
doublereal m_temp
Current value of the temperature at which the properties in this object are calculated (Kelvin)...
Definition: GasTransport.h:205
std::vector< vector_fp > cstar_poly
This is vector of vectors containing the astar fit.
Header file defining class TransportFactory (see TransportFactory)
doublereal Frot(doublereal tr, doublereal sqtr)
The Parker temperature correction to the rotational collision number.
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:101
Functions to evaluate portions of the L matrix needed for multicomponent transport properties...
int m_mode
Type of the polynomial fits to temperature.
Definition: GasTransport.h:160
doublereal err(const std::string &msg) const
Error routine.
vector_fp zrot
Rotational relaxation number for the species in the current phase.
R poly8(D x, R *c)
Templated evaluation of a polynomial of order 8.
Definition: utilities.h:601
void resize(size_t n, size_t m, doublereal v=0.0)
Resize the matrix.
Definition: DenseMatrix.cpp:71
vector_fp eps
Lennard-Jones well-depth of the species in the current phase.
ThermoPhase object for the ideal gas equation of state - workhorse for Cantera (see Thermodynamic Pro...
R poly6(D x, R *c)
Templated evaluation of a polynomial of order 6.
Definition: utilities.h:589
vector_fp m_molefracs_last
Mole fraction vector from last L-matrix evaluation.
virtual void updateDiff_T()
Update the binary diffusion coefficients.
vector_fp m_visc
vector of species viscosities (kg /m /s).
Definition: GasTransport.h:170
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_spwork
work space length = m_kk
Definition: GasTransport.h:166
size_t nColumns() const
Number of columns.
Definition: Array.h:322
void eval_L0000(const doublereal *const x)
Evaluate the L0000 matrices.
Definition: L_matrix.h:23
const doublereal * massFractions() const
Return a const pointer to the mass fraction array.
Definition: Phase.h:454
virtual void getMultiDiffCoeffs(const size_t ld, doublereal *const d)
Return the Multicomponent diffusion coefficients. Units: [m^2/s].
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:68
virtual bool initGas(GasTransportParams &tr)
Initialize the transport operator with parameters from GasTransportParams object. ...
void updateThermal_T()
Update the temperature-dependent terms needed to compute the thermal conductivity and thermal diffusi...
vector_fp crot
Dimensionless rotational heat capacity of the species in the current phase.
DenseMatrix epsilon
The effective well depth for (i,j) collisions.
virtual bool initGas(GasTransportParams &tr)
Called by TransportFactory to set parameters.
std::vector< vector_fp > bstar_poly
This is vector of vectors containing the astar fit.
size_t nSpecies() const
Returns the number of species in the phase.
Definition: Phase.h:252
std::vector< std::vector< int > > poly
This is vector of vectors containing the integer lookup value for the (i,j) interaction.
virtual doublereal pressure() const
Return the thermodynamic pressure (Pa).
Definition: ThermoPhase.h:314
doublereal temperature() const
Temperature (K).
Definition: Phase.h:528
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:165
DenseMatrix dipole
The effective dipole moment for (i,j) collisions.
doublereal m_sqrt_t
current value of temperature to 1/2 power
Definition: GasTransport.h:215
const doublereal Tiny
Small number to compare differences of mole fractions against.
Definition: ct_defs.h:155
const doublereal SqrtPi
sqrt(Pi)
Definition: ct_defs.h:53
doublereal meanMolecularWeight() const
The mean molecular weight. Units: (kg/kmol)
Definition: Phase.h:588
const doublereal GasConstant
Universal Gas Constant. [J/kmol/K].
Definition: ct_defs.h:66
Contains declarations for string manipulation functions within Cantera.
#define DATA_PTR(vec)
Creates a pointer to the start of the raw data for a vector.
Definition: ct_defs.h:36
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 ...
size_t m_nsp
Number of species.
doublereal m_logt
Current value of the log of the temperature.
Definition: GasTransport.h:218
bool m_abc_ok
Boolean indicating viscosity is up to date.
doublereal m_kbt
Current value of Boltzman's constant times the temperature (Joules)
Definition: GasTransport.h:208
void eval_L0010(const doublereal *const x)
Evaluate the L0010 matrices.
Definition: L_matrix.h:46
Class that holds the data that is read in from the xml file, and which is used for processing of the ...
const doublereal Boltzmann
Boltzmann's constant [J/K].
Definition: ct_defs.h:78
DenseMatrix m_om22
Dense matrix for omega22.
vector_int & ipiv()
Return a changeable value of the pivot vector.
vector_fp m_molefracs
Vector of species mole fractions.
Definition: GasTransport.h:141
Class GasTransport implements some functions and properties that are shared by the MixTransport and M...
Definition: GasTransport.h:16
void update_C()
Update basic concentration-dependent quantities if the concentrations have changed.
vector_fp m_mw
Local copy of the species molecular weights.
Definition: GasTransport.h:178