Cantera  2.1.2
SimpleTransport.h
Go to the documentation of this file.
1 /**
2  * @file SimpleTransport.h
3  * Header file for the class SimpleTransport which provides simple
4  * transport properties for liquids and solids
5  * (see \ref tranprops and \link Cantera::SimpleTransport SimpleTransport \endlink) .
6  */
7 #ifndef CT_SIMPLETRAN_H
8 #define CT_SIMPLETRAN_H
9 
10 #include "TransportBase.h"
12 #include "TransportParams.h"
13 #include "LiquidTransportParams.h"
14 
15 namespace Cantera
16 {
17 
18 class LiquidTransportParams;
19 
20 
21 //! Class SimpleTransport implements mixture-averaged transport
22 //! properties for liquid phases.
23 /*!
24  * The model is based on that described by Newman, Electrochemical Systems
25  *
26  * The velocity of species i may be described by the
27  * following equation p. 297 (12.1)
28  *
29  * \f[
30  * c_i \nabla \mu_i = R T \sum_j \frac{c_i c_j}{c_T D_{ij}}
31  * (\mathbf{v}_j - \mathbf{v}_i)
32  * \f]
33  *
34  * This as written is degenerate by 1 dof.
35  *
36  * To fix this we must add in the definition of the mass averaged
37  * velocity of the solution. We will call the simple bold-faced
38  * \f$\mathbf{v} \f$
39  * symbol the mass-averaged velocity. Then, the relation
40  * between \f$\mathbf{v}\f$ and the individual species velocities is
41  * \f$\mathbf{v}_i\f$
42  *
43  * \f[
44  * \rho_i \mathbf{v}_i = \rho_i \mathbf{v} + \mathbf{j}_i
45  * \f]
46  * where \f$\mathbf{j}_i\f$ are the diffusional fluxes of species i
47  * with respect to the mass averaged velocity and
48  *
49  * \f[
50  * \sum_i \mathbf{j}_i = 0
51  * \f]
52  *
53  * and
54  *
55  * \f[
56  * \sum_i \rho_i \mathbf{v}_i = \rho \mathbf{v}
57  * \f]
58  *
59  * Using these definitions, we can write
60  *
61  * \f[
62  * \mathbf{v}_i = \mathbf{v} + \frac{\mathbf{j}_i}{\rho_i}
63  * \f]
64  *
65  * \f[
66  * c_i \nabla \mu_i = R T \sum_j \frac{c_i c_j}{c_T D_{ij}}
67  * (\frac{\mathbf{j}_j}{\rho_j} - \frac{\mathbf{j}_i}{\rho_i})
68  * = R T \sum_j \frac{1}{D_{ij}}
69  * (\frac{x_i \mathbf{j}_j}{M_j} - \frac{x_j \mathbf{j}_i}{M_i})
70  * \f]
71  *
72  * The equations that we actually solve are
73  *
74  * \f[
75  * c_i \nabla \mu_i =
76  * = R T \sum_j \frac{1}{D_{ij}}
77  * (\frac{x_i \mathbf{j}_j}{M_j} - \frac{x_j \mathbf{j}_i}{M_i})
78  * \f]
79  * and we replace the 0th equation with the following:
80  *
81  * \f[
82  * \sum_i \mathbf{j}_i = 0
83  * \f]
84  *
85  * When there are charged species, we replace the rhs with the
86  * gradient of the electrochemical potential to obtain the
87  * modified equation
88  *
89  * \f[
90  * c_i \nabla \mu_i + c_i F z_i \nabla \Phi
91  * = R T \sum_j \frac{1}{D_{ij}}
92  * (\frac{x_i \mathbf{j}_j}{M_j} - \frac{x_j \mathbf{j}_i}{M_i})
93  * \f]
94  *
95  * With this formulation we may solve for the diffusion velocities,
96  * without having to worry about what the mass averaged velocity
97  * is.
98  *
99  * <H2> Viscosity Calculation </H2>
100  *
101  * The viscosity calculation may be broken down into two parts.
102  * In the first part, the viscosity of the pure species are calculated
103  * In the second part, a mixing rule is applied. There are two mixing rules.
104  * Solvent-only and mixture-averaged.
105  *
106  * For the solvent-only mixing rule, we use the pure species viscosity calculated for
107  * the solvent as the viscosity of the entire mixture. For the mixture averaged rule
108  * we do a mole fraction based average of the pure species viscosities:
109  *
110  * Solvent-only:
111  * \f[
112  * \mu = \mu_0
113  * \f]
114  * Mixture-average:
115  * \f[
116  * \mu = \sum_k {\mu_k X_k}
117  * \f]
118  *
119  * <H2> Calculate of the Binary Diffusion Coefficients </H2>
120  *
121  * The binary diffusion coefficients are obtained from the pure species diffusion coefficients
122  * using an additive process
123  *
124  * \f[
125  * D_{i,j} = \frac{1}{2} \left( D^0_i(T) + D^0_j(T) \right)
126  * \f]
127  *
128  * <H2> Electrical Mobilities </H2>
129  *
130  * The mobility \f$ \mu^e_k \f$ is calculated from the diffusion coefficient using the Einstein relation.
131  *
132  * \f[
133  * \mu^e_k = \frac{F D_k}{R T}
134  * \f]
135  *
136  * The diffusion coefficients, \f$ D_k \f$ , is calculated from a call to the mixture diffusion
137  * coefficient routine.
138  *
139  * <H2> Species Diffusive Fluxes </H2>
140  *
141  * The diffusive mass flux of species \e k is computed from the following
142  * formula
143  *
144  * Usually the specified solution average velocity is the mass averaged velocity.
145  * This is changed in some subclasses, however.
146  *
147  * \f[
148  * j_k = - c^T M_k D_k \nabla X_k - \rho Y_k V_c
149  * \f]
150  *
151  * where V_c is the correction velocity
152  *
153  * \f[
154  * \rho V_c = - \sum_j {c^T M_j D_j \nabla X_j}
155  * \f]
156  *
157  * In the above equation, \f$ D_k \f$ is the mixture diffusivity for species k calculated for the current
158  * conditions, which may depend on T, P, and X_k. \f$ C^T \f$ is the total concentration of the phase.
159  *
160  * When this is electrical migration, the formulas above are enhanced to
161  *
162  * \f[
163  * j_k = - C^T M_k D_k \nabla X_k + F C^T M_k \frac{D_k}{ R T } X_k z_k \nabla V - \rho Y_k V_c
164  * \f]
165  *
166  * where V_c is the correction velocity
167  *
168  * \f[
169  * \rho V_c = - \sum_j {c^T M_j D_j \nabla X_j} + \sum_j F C^T M_j \frac{D_j}{ R T } X_j z_j \nabla V
170  * \f]
171  *
172  * <H2> Species Diffusional Velocities </H2>
173  *
174  * Species diffusional velocities are calculated from the species diffusional fluxes, within this object,
175  * using the following formula for the diffusional velocity of the kth species, \f$ V_k^d \f$
176  *
177  * \f[
178  * j_k = \rho Y_k V_k^d
179  * \f]
180  *
181  * TODO
182  * This object has to be made compatible with different types of reference velocities. Right now, elements
183  * of the formulas are only compatible with the mass-averaged velocity.
184  *
185  * @ingroup tranprops
186  */
188 {
189 public:
190  //! Default constructor.
191  /*!
192  * This requires call to initLiquid(LiquidTransportParams& tr)
193  * after filling LiquidTransportParams to complete instantiation.
194  * The filling of LiquidTransportParams is currently carried out
195  * in the TransportFactory class, but might be moved at some point.
196  *
197  * @param thermo ThermoPhase object holding species information.
198  * @param ndim Number of spatial dimensions.
199  */
200  SimpleTransport(thermo_t* thermo = 0, int ndim = 1);
201 
202  SimpleTransport(const SimpleTransport& right);
203  SimpleTransport& operator=(const SimpleTransport& right);
204  virtual Transport* duplMyselfAsTransport() const;
205  virtual ~SimpleTransport();
206 
207  //! Initialize the transport object
208  /*!
209  * Here we change all of the internal dimensions to be sufficient.
210  * We get the object ready to do property evaluations.
211  *
212  * @param tr Transport parameters for all of the species in the phase.
213  */
214  virtual bool initLiquid(LiquidTransportParams& tr);
215 
216  virtual int model() const {
217  return cSimpleTransport;
218  }
219 
220  //! Returns the mixture viscosity of the solution
221  /*!
222  * The viscosity is computed using the general mixture rules
223  * specified in the variable compositionDepType_.
224  *
225  * Solvent-only:
226  * \f[
227  * \mu = \mu_0
228  * \f]
229  * Mixture-average:
230  * \f[
231  * \mu = \sum_k {\mu_k X_k}
232  * \f]
233  *
234  * Here \f$ \mu_k \f$ is the viscosity of pure species \e k.
235  *
236  * units are Pa s or kg/m/s
237  *
238  * @see updateViscosity_T();
239  */
240  virtual doublereal viscosity();
241 
242  //! Returns the pure species viscosities
243  /*!
244  * The pure species viscosities are to be given in an Arrhenius
245  * form in accordance with activated-jump-process dominated transport.
246  *
247  * units are Pa s or kg/m/s
248  *
249  * @param visc Return the species viscosities as a vector of length m_nsp
250  */
251  virtual void getSpeciesViscosities(doublereal* const visc);
252 
253  //! Returns the binary diffusion coefficients
254  /*!
255  * @param ld
256  * @param d
257  */
258  virtual void getBinaryDiffCoeffs(const size_t ld, doublereal* const d);
259 
260  //! Get the Mixture diffusion coefficients
261  /*!
262  * @param d vector of mixture diffusion coefficients
263  * units = m2 s-1. length = number of species
264  */
265  virtual void getMixDiffCoeffs(doublereal* const d);
266 
267  //! Return the thermal diffusion coefficients
268  /*!
269  * These are all zero for this simple implementation
270  *
271  * @param dt thermal diffusion coefficients
272  */
273  virtual void getThermalDiffCoeffs(doublereal* const dt);
274 
275  //! Returns the mixture thermal conductivity of the solution
276  /*!
277  * The thermal is computed using the general mixture rules
278  * specified in the variable compositionDepType_.
279  *
280  * Controlling update boolean = m_condmix_ok
281  *
282  * Units are in W/m/K or equivalently kg m / s3 / K
283  *
284  * Solvent-only:
285  * \f[
286  * \lambda = \lambda_0
287  *
288  * \f]
289  * Mixture-average:
290  * \f[
291  * \lambda = \sum_k {\lambda_k X_k}
292  * \f]
293  *
294  * Here \f$ \lambda_k \f$ is the thermal conductivity of pure species \e k.
295  *
296  * @see updateCond_T();
297  */
298  virtual doublereal thermalConductivity();
299 
300  virtual void getMobilities(doublereal* const mobil_e);
301 
302  virtual void getFluidMobilities(doublereal* const mobil_f);
303 
304  //! Specify the value of the gradient of the voltage
305  /*!
306  * @param grad_V Gradient of the voltage (length num dimensions);
307  */
308  virtual void set_Grad_V(const doublereal* const grad_V);
309 
310  //! Specify the value of the gradient of the temperature
311  /*!
312  * @param grad_T Gradient of the temperature (length num dimensions);
313  */
314  virtual void set_Grad_T(const doublereal* const grad_T);
315 
316  //! Specify the value of the gradient of the MoleFractions
317  /*!
318  * @param grad_X Gradient of the mole fractions(length nsp * num dimensions);
319  */
320  virtual void set_Grad_X(const doublereal* const grad_X);
321 
322  //! Get the species diffusive velocities wrt to the averaged velocity,
323  //! given the gradients in mole fraction and temperature
324  /*!
325  * The average velocity can be computed on a mole-weighted
326  * or mass-weighted basis, or the diffusion velocities may
327  * be specified as relative to a specific species (i.e. a
328  * solvent) all according to the velocityBasis input parameter.
329  *
330  * Units for the returned velocities are m s-1.
331  *
332  * @param ndim Number of dimensions in the flux expressions
333  * @param grad_T Gradient of the temperature
334  * (length = ndim)
335  * @param ldx Leading dimension of the grad_X array
336  * (usually equal to m_nsp but not always)
337  * @param grad_X Gradients of the mole fraction
338  * Flat vector with the m_nsp in the inner loop.
339  * length = ldx * ndim
340  * @param ldf Leading dimension of the fluxes array
341  * (usually equal to m_nsp but not always)
342  * @param Vdiff Output of the diffusive velocities.
343  * Flat vector with the m_nsp in the inner loop.
344  * length = ldx * ndim
345  */
346  virtual void getSpeciesVdiff(size_t ndim,
347  const doublereal* grad_T,
348  int ldx,
349  const doublereal* grad_X,
350  int ldf,
351  doublereal* Vdiff);
352 
353  //! Get the species diffusive velocities wrt to the averaged velocity,
354  //! given the gradients in mole fraction, temperature and electrostatic potential.
355  /*!
356  * The average velocity can be computed on a mole-weighted
357  * or mass-weighted basis, or the diffusion velocities may
358  * be specified as relative to a specific species (i.e. a
359  * solvent) all according to the velocityBasis input parameter.
360  *
361  * Units for the returned velocities are m s-1.
362  *
363  * @param ndim Number of dimensions in the flux expressions
364  * @param grad_T Gradient of the temperature
365  * (length = ndim)
366  * @param ldx Leading dimension of the grad_X array
367  * (usually equal to m_nsp but not always)
368  * @param grad_X Gradients of the mole fraction
369  * Flat vector with the m_nsp in the inner loop.
370  * length = ldx * ndim
371  * @param ldf Leading dimension of the fluxes array
372  * (usually equal to m_nsp but not always)
373  * @param grad_Phi Gradients of the electrostatic potential
374  * (length = ndim)
375  * @param Vdiff Output of the species diffusion velocities
376  * Flat vector with the m_nsp in the inner loop.
377  * length = ldx * ndim
378  */
379  virtual void getSpeciesVdiffES(size_t ndim, const doublereal* grad_T,
380  int ldx, const doublereal* grad_X,
381  int ldf, const doublereal* grad_Phi,
382  doublereal* Vdiff);
383 
384  //! Get the species diffusive mass fluxes wrt to the specified solution averaged velocity,
385  //! given the gradients in mole fraction and temperature
386  /*!
387  * units = kg/m2/s
388  *
389  * The diffusive mass flux of species \e k is computed from the following
390  * formula
391  *
392  * Usually the specified solution average velocity is the mass averaged velocity.
393  * This is changed in some subclasses, however.
394  *
395  * \f[
396  * j_k = - \rho M_k D_k \nabla X_k - Y_k V_c
397  * \f]
398  *
399  * where V_c is the correction velocity
400  *
401  * \f[
402  * V_c = - \sum_j {\rho M_j D_j \nabla X_j}
403  * \f]
404  *
405  * @param ndim The number of spatial dimensions (1, 2, or 3).
406  * @param grad_T The temperature gradient (ignored in this model).
407  * @param ldx Leading dimension of the grad_X array.
408  * @param grad_X Gradient of the mole fractions(length nsp * num dimensions);
409  * @param ldf Leading dimension of the fluxes array.
410  * @param fluxes Output fluxes of species.
411  */
412  virtual void getSpeciesFluxes(size_t ndim, const doublereal* const grad_T,
413  size_t ldx, const doublereal* const grad_X,
414  size_t ldf, doublereal* const fluxes);
415 
416  //! Return the species diffusive mass fluxes wrt to
417  //! the mass averaged velocity,
418  /*!
419  * units = kg/m2/s
420  *
421  * Internally, gradients in the in mole fraction, temperature
422  * and electrostatic potential contribute to the diffusive flux
423  *
424  * The diffusive mass flux of species \e k is computed from the following
425  * formula
426  *
427  * \f[
428  * j_k = - \rho M_k D_k \nabla X_k - Y_k V_c
429  * \f]
430  *
431  * where V_c is the correction velocity
432  *
433  * \f[
434  * V_c = - \sum_j {\rho M_j D_j \nabla X_j}
435  * \f]
436  *
437  * @param ldf stride of the fluxes array. Must be equal to
438  * or greater than the number of species.
439  * @param fluxes Vector of calculated fluxes
440  */
441  virtual void getSpeciesFluxesExt(size_t ldf, doublereal* fluxes);
442 
443 protected:
444  //! Handles the effects of changes in the Temperature, internally
445  //! within the object.
446  /*!
447  * This is called whenever a transport property is requested.
448  * The first task is to check whether the temperature has changed
449  * since the last call to update_T().
450  * If it hasn't then an immediate return is carried out.
451  *
452  * @return Returns true if the temperature has changed, and false otherwise
453  */
454  virtual bool update_T();
455 
456  //! Handles the effects of changes in the mixture concentration
457  /*!
458  * This is called for every interface call to check whether
459  * the concentrations have changed. Concentrations change
460  * whenever the pressure or the mole fraction has changed.
461  * If it has changed, the recalculations should be done.
462  *
463  * Note this should be a lightweight function since it's
464  * part of all of the interfaces.
465  */
466  virtual bool update_C();
467 
468  //! Update the temperature-dependent viscosity terms.
469  //! Updates the array of pure species viscosities, and the
470  //! weighting functions in the viscosity mixture rule.
471  /*!
472  * The flag m_visc_temp_ok is set to true.
473  */
474  void updateViscosity_T();
475 
476  //! Update the temperature-dependent parts of the mixture-averaged
477  //! thermal conductivity.
478  void updateCond_T();
479 
480  //! Update the concentration parts of the viscosities
481  /*!
482  * Internal routine is run whenever the update_boolean is false. This
483  * routine will calculate internal values for the species viscosities.
484  */
485  void updateViscosities_C();
486 
487  //! Update the binary diffusion coefficients wrt T.
488  /*!
489  * These are evaluated from the polynomial fits at unit pressure (1 Pa).
490  */
491  void updateDiff_T();
492 
493 private:
494  //! Temperature dependence type
495  /*!
496  * The following coefficients are allowed to have simple
497  * temperature dependencies:
498  * mixture viscosity
499  * mixture thermal conductivity
500  * diffusitivy
501  *
502  * Types of temperature dependencies:
503  * 0 - Independent of temperature (only one implemented so far)
504  * 1 - extended arrhenius form
505  * 2 - polynomial in temperature form
506  */
508 
509  //! Composition dependence of the transport properties
510  /*!
511  * The following coefficients are allowed to have simple
512  * composition dependencies
513  * mixture viscosity
514  * mixture thermal conductivity
515  *
516  * Types of composition dependencies
517  * 0 - Solvent values (i.e., species 0) contributes only
518  * 1 - linear combination of mole fractions;
519  */
521 
522  //! Boolean indicating whether to use the hydrodynamic radius formulation
523  /*!
524  * If true, then the diffusion coefficient is calculated from the
525  * hydrodynamic radius.
526  */
528 
529  //! Boolean indicating whether electro-migration term should be added
531 
532  //! Local Copy of the molecular weights of the species
533  /*!
534  * Length is Equal to the number of species in the mechanism.
535  */
537 
538  //! Pure species viscosities in Arrhenius temperature-dependent form.
539  std::vector<LTPspecies*> m_coeffVisc_Ns;
540 
541  //! Pure species thermal conductivities in Arrhenius temperature-dependent form.
542  std::vector<LTPspecies*> m_coeffLambda_Ns;
543 
544  //! Pure species viscosities in Arrhenius temperature-dependent form.
545  std::vector<LTPspecies*> m_coeffDiff_Ns;
546 
547  //! Hydrodynamic radius in LTPspecies form
548  std::vector<LTPspecies*> m_coeffHydroRadius_Ns;
549 
550  //! Internal value of the gradient of the mole fraction vector
551  /*!
552  * Note, this is the only gradient value that can and perhaps
553  * should reflect the true state of the mole fractions in the
554  * application solution vector. In other words no cropping or
555  * massaging of the values to make sure they are above zero
556  * should occur. - developing ....
557  *
558  * m_nsp is the number of species in the fluid
559  * k is the species index
560  * n is the dimensional index (x, y, or z). It has a length
561  * equal to m_nDim
562  *
563  * m_Grad_X[n*m_nsp + k]
564  */
566 
567  //! Internal value of the gradient of the Temperature vector
568  /*!
569  * Generally, if a transport property needs this in its evaluation it
570  * will look to this place to get it.
571  *
572  * No internal property is precalculated based on gradients. Gradients
573  * are assumed to be freshly updated before every property call.
574  */
576 
577  //! Internal value of the gradient of the Pressure vector
578  /*!
579  * Generally, if a transport property needs this in its evaluation it
580  * will look to this place to get it.
581  *
582  * No internal property is precalculated based on gradients. Gradients
583  * are assumed to be freshly updated before every property call.
584  */
586 
587  //! Internal value of the gradient of the Electric Voltage
588  /*!
589  * Generally, if a transport property needs this in its evaluation it
590  * will look to this place to get it.
591  *
592  * No internal property is precalculated based on gradients. Gradients
593  * are assumed to be freshly updated before every property call.
594  */
596 
597  // property values
598 
599  //! Vector of Species Diffusivities
600  /*!
601  * Depends on the temperature. We have set the pressure dependence
602  * to zero for this liquid phase constituitve model
603  *
604  * units m2/s
605  */
607 
608  //! Species viscosities
609  /*!
610  * Viscosity of the species
611  * Length = number of species
612  *
613  * Depends on the temperature. We have set the pressure dependence
614  * to zero for this model
615  *
616  * controlling update boolean -> m_visc_temp_ok
617  */
619 
620  //! Internal value of the species individual thermal conductivities
621  /*!
622  * Then a mixture rule is applied to get the solution conductivities
623  *
624  * Depends on the temperature and perhaps pressure, but
625  * not the species concentrations
626  *
627  * controlling update boolean -> m_cond_temp_ok
628  */
630 
631  //! State of the mole fraction vector.
633 
634  //! Local copy of the mole fractions of the species in the phase
635  /*!
636  * The mole fractions here are assumed to be bounded by 0.0 and 1.0
637  * and they are assumed to add up to one exactly. This mole
638  * fraction vector comes from the ThermoPhase object. Derivative
639  * quantities from this are referred to as bounded.
640  *
641  * Update info?
642  * length = m_nsp
643  */
645 
646  //! Local copy of the concentrations of the species in the phase
647  /*!
648  * The concentrations are consistent with the m_molefracs
649  * vector which is bounded and sums to one.
650  *
651  * Update info?
652  * length = m_nsp
653  */
655 
656  //! Local copy of the total concentration.
657  /*!
658  * This is consistent with the m_concentrations[] and
659  * m_molefracs[] vector.
660  */
661  doublereal concTot_;
662 
663  //! Mean molecular weight
665 
666  //! Density
667  doublereal dens_;
668 
669  //! Local copy of the charge of each species
670  /*!
671  * Contains the charge of each species (length m_nsp)
672  */
674 
675  //! Current Temperature -> locally stored
676  /*!
677  * This is used to test whether new temperature computations
678  * should be performed.
679  */
680  doublereal m_temp;
681 
682  //! Current value of the pressure
683  doublereal m_press;
684 
685  //! Saved value of the mixture thermal conductivity
686  doublereal m_lambda;
687 
688  //! Saved value of the mixture viscosity
689  doublereal m_viscmix;
690 
691  //! work space
692  /*!
693  * Length is equal to m_nsp
694  */
696 
697  vector_fp m_fluxes;
698 
699 private:
700  //! Boolean indicating that the top-level mixture viscosity is current
701  /*!
702  * This is turned false for every change in T, P, or C.
703  */
705 
706  //! Boolean indicating that weight factors wrt viscosity is current
708 
709  //! Boolean indicating that mixture diffusion coeffs are current
711 
712  //! Boolean indicating that binary diffusion coeffs are current
714 
715  //! Flag to indicate that the pure species conductivities
716  //! are current wrt the temperature
718 
719  //! Boolean indicating that mixture conductivity is current
721 
722  //! Number of dimensions
723  /*!
724  * Either 1, 2, or 3
725  */
726  size_t m_nDim;
727 
728  //! Temporary variable that stores the rho Vc value
729  double rhoVc[3];
730 
731 private:
732 
733  //! Throw an exception if this method is invoked.
734  /*!
735  * This probably indicates something is not yet implemented.
736  *
737  * @param msg Indicates the member function which is not implemented
738  */
739  doublereal err(const std::string& msg) const;
740 
741 };
742 }
743 #endif
bool m_visc_temp_ok
Boolean indicating that weight factors wrt viscosity is current.
vector_fp m_Grad_T
Internal value of the gradient of the Temperature vector.
vector_fp m_Grad_V
Internal value of the gradient of the Electric Voltage.
vector_fp m_Grad_P
Internal value of the gradient of the Pressure vector.
virtual int model() const
Transport model.
virtual bool initLiquid(LiquidTransportParams &tr)
Initialize the transport object.
vector_fp m_spwork
work space
bool m_visc_mix_ok
Boolean indicating that the top-level mixture viscosity is current.
virtual void getBinaryDiffCoeffs(const size_t ld, doublereal *const d)
Returns the binary diffusion coefficients.
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 specified solution averaged velocity, given the gradients in mole fraction and temperature.
Headers for the Transport object, which is the virtual base class for all transport property evaluato...
doublereal m_viscmix
Saved value of the mixture viscosity.
double rhoVc[3]
Temporary variable that stores the rho Vc value.
virtual bool update_T()
Handles the effects of changes in the Temperature, internally within the object.
void updateViscosities_C()
Update the concentration parts of the viscosities.
bool m_cond_mix_ok
Boolean indicating that mixture conductivity is current.
Base class for transport property managers.
void updateDiff_T()
Update the binary diffusion coefficients wrt T.
virtual void getMobilities(doublereal *const mobil_e)
Get the Electrical mobilities (m^2/V/s).
virtual bool update_C()
Handles the effects of changes in the mixture concentration.
doublereal m_lambda
Saved value of the mixture thermal conductivity.
vector_fp m_Grad_X
Internal value of the gradient of the mole fraction vector.
void updateCond_T()
Update the temperature-dependent parts of the mixture-averaged thermal conductivity.
int m_iStateMF
State of the mole fraction vector.
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:101
Class LiquidTransportParams holds transport model parameters relevant to transport in mixtures...
Header file defining class LiquidTransportParams.
doublereal meanMolecularWeight_
Mean molecular weight.
bool m_diff_temp_ok
Boolean indicating that binary diffusion coeffs are current.
virtual void getFluidMobilities(doublereal *const mobil_f)
Get the fluid mobilities (s kmol/kg).
std::vector< LTPspecies * > m_coeffVisc_Ns
Pure species viscosities in Arrhenius temperature-dependent form.
virtual void getMixDiffCoeffs(doublereal *const d)
Get the Mixture diffusion coefficients.
virtual void getSpeciesVdiffES(size_t ndim, const doublereal *grad_T, int ldx, const doublereal *grad_X, int ldf, const doublereal *grad_Phi, doublereal *Vdiff)
Get the species diffusive velocities wrt to the averaged velocity, given the gradients in mole fracti...
vector_fp m_diffSpecies
Vector of Species Diffusivities.
vector_fp m_viscSpecies
Species viscosities.
int tempDepType_
Temperature dependence type.
std::vector< LTPspecies * > m_coeffLambda_Ns
Pure species thermal conductivities in Arrhenius temperature-dependent form.
SimpleTransport(thermo_t *thermo=0, int ndim=1)
Default constructor.
doublereal m_press
Current value of the pressure.
std::vector< LTPspecies * > m_coeffHydroRadius_Ns
Hydrodynamic radius in LTPspecies form.
virtual void set_Grad_X(const doublereal *const grad_X)
Specify the value of the gradient of the MoleFractions.
virtual void getSpeciesVdiff(size_t ndim, const doublereal *grad_T, int ldx, const doublereal *grad_X, int ldf, doublereal *Vdiff)
Get the species diffusive velocities wrt to the averaged velocity, given the gradients in mole fracti...
virtual doublereal thermalConductivity()
Returns the mixture thermal conductivity of the solution.
thermo_t & thermo()
virtual void set_Grad_V(const doublereal *const grad_V)
Specify the value of the gradient of the voltage.
doublereal concTot_
Local copy of the total concentration.
doublereal err(const std::string &msg) const
Throw an exception if this method is invoked.
vector_fp m_concentrations
Local copy of the concentrations of the species in the phase.
void updateViscosity_T()
Update the temperature-dependent viscosity terms.
virtual void getThermalDiffCoeffs(doublereal *const dt)
Return the thermal diffusion coefficients.
doublereal dens_
Density.
virtual doublereal viscosity()
Returns the mixture viscosity of the solution.
doublereal m_temp
Current Temperature -> locally stored.
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
vector_fp m_condSpecies
Internal value of the species individual thermal conductivities.
bool doMigration_
Boolean indicating whether electro-migration term should be added.
bool useHydroRadius_
Boolean indicating whether to use the hydrodynamic radius formulation.
virtual void getSpeciesFluxesExt(size_t ldf, doublereal *fluxes)
Return the species diffusive mass fluxes wrt to the mass averaged velocity,.
Headers for the DenseMatrix object, which deals with dense rectangular matrices and description of th...
virtual Transport * duplMyselfAsTransport() const
Duplication routine for objects which inherit from Transport.
int compositionDepType_
Composition dependence of the transport properties.
size_t m_nDim
Number of dimensions.
vector_fp m_chargeSpecies
Local copy of the charge of each species.
bool m_cond_temp_ok
Flag to indicate that the pure species conductivities are current wrt the temperature.
vector_fp m_molefracs
Local copy of the mole fractions of the species in the phase.
virtual void getSpeciesViscosities(doublereal *const visc)
Returns the pure species viscosities.
virtual void set_Grad_T(const doublereal *const grad_T)
Specify the value of the gradient of the temperature.
std::vector< LTPspecies * > m_coeffDiff_Ns
Pure species viscosities in Arrhenius temperature-dependent form.
Class that holds the data that is read in from the xml file, and which is used for processing of the ...
Class SimpleTransport implements mixture-averaged transport properties for liquid phases...
bool m_diff_mix_ok
Boolean indicating that mixture diffusion coeffs are current.
vector_fp m_mw
Local Copy of the molecular weights of the species.