Cantera 2.6.0
|
Class UnityLewisTransport implements the unity Lewis number approximation for the mixture-averaged species diffusion coefficients. More...
#include <UnityLewisTransport.h>
Public Member Functions | |
virtual std::string | transportType () const |
Identifies the Transport object type. More... | |
virtual void | getMixDiffCoeffs (double *const d) |
Returns the unity Lewis number approximation based diffusion coefficients [m^2/s]. More... | |
virtual void | getMixDiffCoeffsMole (double *const d) |
Not implemented for unity Lewis number approximation. More... | |
virtual void | getMixDiffCoeffsMass (double *const d) |
Returns the unity Lewis number approximation based diffusion coefficients [m^2/s]. More... | |
Public Member Functions inherited from MixTransport | |
MixTransport () | |
Default constructor. More... | |
virtual void | getThermalDiffCoeffs (doublereal *const dt) |
Return the thermal diffusion coefficients. More... | |
virtual doublereal | thermalConductivity () |
Returns the mixture thermal conductivity (W/m /K) More... | |
virtual void | getMobilities (doublereal *const mobil) |
Get the Electrical mobilities (m^2/V/s). More... | |
virtual void | update_T () |
Update the internal parameters whenever the temperature has changed. More... | |
virtual void | update_C () |
Update the internal parameters whenever the concentrations have changed. More... | |
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 fraction and temperature. More... | |
virtual void | init (ThermoPhase *thermo, int mode=0, int log_level=0) |
Initialize a transport manager. More... | |
Public Member Functions inherited from GasTransport | |
virtual doublereal | viscosity () |
Viscosity of the mixture (kg /m /s) More... | |
virtual void | getSpeciesViscosities (doublereal *const visc) |
Get the pure-species viscosities. More... | |
virtual void | getBinaryDiffCoeffs (const size_t ld, doublereal *const d) |
Returns the matrix of binary diffusion coefficients. More... | |
virtual void | getViscosityPolynomial (size_t i, double *coeffs) const |
Return the polynomial fits to the viscosity of species i. More... | |
virtual void | getConductivityPolynomial (size_t i, double *coeffs) const |
Return the temperature fits of the heat conductivity of species i. More... | |
virtual void | getBinDiffusivityPolynomial (size_t i, size_t j, double *coeffs) const |
Return the polynomial fits to the binary diffusivity of species pair (i, j) More... | |
virtual void | getCollisionIntegralPolynomial (size_t i, size_t j, double *astar_coeffs, double *bstar_coeffs, double *cstar_coeffs) const |
Return the polynomial fits to the collision integral of species pair (i, j) More... | |
virtual void | setViscosityPolynomial (size_t i, double *coeffs) |
Modify the polynomial fits to the viscosity of species i. More... | |
virtual void | setConductivityPolynomial (size_t i, double *coeffs) |
Modify the temperature fits of the heat conductivity of species i. More... | |
virtual void | setBinDiffusivityPolynomial (size_t i, size_t j, double *coeffs) |
Modify the polynomial fits to the binary diffusivity of species pair (i, j) More... | |
virtual void | setCollisionIntegralPolynomial (size_t i, size_t j, double *astar_coeffs, double *bstar_coeffs, double *cstar_coeffs, bool actualT) |
Modify the polynomial fits to the collision integral of species pair (i, j) More... | |
bool | CKMode () const |
Boolean indicating the form of the transport properties polynomial fits. More... | |
Public Member Functions inherited from Transport | |
Transport (ThermoPhase *thermo=0, size_t ndim=1) | |
Constructor. More... | |
Transport (const Transport &)=delete | |
Transport & | operator= (const Transport &)=delete |
ThermoPhase & | thermo () |
bool | ready () |
void | setNDim (const int ndim) |
Set the number of dimensions to be expected in flux expressions. More... | |
size_t | nDim () const |
Return the number of dimensions in flux expressions. More... | |
void | checkSpeciesIndex (size_t k) const |
Check that the specified species index is in range. More... | |
void | checkSpeciesArraySize (size_t kk) const |
Check that an array size is at least nSpecies(). More... | |
virtual doublereal | getElectricConduct () |
Compute the mixture electrical conductivity (S m-1) at the current conditions of the phase (Siemens m-1) More... | |
virtual void | getElectricCurrent (int ndim, const doublereal *grad_T, int ldx, const doublereal *grad_X, int ldf, const doublereal *grad_V, doublereal *current) |
Compute the electric current density in A/m^2. More... | |
virtual void | getSpeciesFluxesES (size_t ndim, const doublereal *grad_T, size_t ldx, const doublereal *grad_X, size_t ldf, const doublereal *grad_Phi, doublereal *fluxes) |
Get the species diffusive mass fluxes wrt to the mass averaged velocity, given the gradients in mole fraction, temperature and electrostatic potential. More... | |
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 mass averaged velocity, given the gradients in mole fraction and temperature. More... | |
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 mass averaged velocity, given the gradients in mole fraction, temperature, and electrostatic potential. More... | |
virtual void | getMolarFluxes (const doublereal *const state1, const doublereal *const state2, const doublereal delta, doublereal *const cfluxes) |
Get the molar fluxes [kmol/m^2/s], given the thermodynamic state at two nearby points. More... | |
virtual void | getMassFluxes (const doublereal *state1, const doublereal *state2, doublereal delta, doublereal *mfluxes) |
Get the mass fluxes [kg/m^2/s], given the thermodynamic state at two nearby points. More... | |
virtual void | getMultiDiffCoeffs (const size_t ld, doublereal *const d) |
Return the Multicomponent diffusion coefficients. Units: [m^2/s]. More... | |
virtual void | setParameters (const int type, const int k, const doublereal *const p) |
Set model parameters for derived classes. More... | |
AnyMap | parameters () const |
Return the parameters for a phase definition which are needed to reconstruct an identical object using the newTransport function. More... | |
void | setVelocityBasis (VelocityBasis ivb) |
Sets the velocity basis. More... | |
VelocityBasis | getVelocityBasis () const |
Gets the velocity basis. More... | |
virtual doublereal | bulkViscosity () |
The bulk viscosity in Pa-s. More... | |
virtual doublereal | ionConductivity () |
The ionic conductivity in 1/ohm/m. More... | |
virtual void | getSpeciesIonConductivity (doublereal *const ionCond) |
Returns the pure species ionic conductivity. More... | |
virtual void | mobilityRatio (double *mobRat) |
Returns the pointer to the mobility ratios of the species in the phase. More... | |
virtual void | getSpeciesMobilityRatio (double **mobRat) |
Returns the pure species limit of the mobility ratios. More... | |
virtual doublereal | electricalConductivity () |
The electrical conductivity (Siemens/m). More... | |
virtual void | getFluidMobilities (doublereal *const mobil_f) |
Get the fluid mobilities (s kmol/kg). More... | |
virtual void | setThermo (ThermoPhase &thermo) |
Specifies the ThermoPhase object. More... | |
virtual void | setRoot (std::shared_ptr< Solution > root) |
Set root Solution holding all phase information. More... | |
Additional Inherited Members | |
Protected Member Functions inherited from MixTransport | |
void | updateCond_T () |
Update the temperature dependent parts of the species thermal conductivities. More... | |
Protected Member Functions inherited from GasTransport | |
GasTransport (ThermoPhase *thermo=0) | |
virtual void | updateViscosity_T () |
Update the temperature-dependent viscosity terms. More... | |
virtual void | updateSpeciesViscosities () |
Update the pure-species viscosities. More... | |
virtual void | updateDiff_T () |
Update the binary diffusion coefficients. More... | |
virtual void | setupCollisionParameters () |
Setup parameters for a new kinetic-theory-based transport manager for low-density gases. More... | |
void | setupCollisionIntegral () |
Setup range for polynomial fits to collision integrals of Monchick & Mason. More... | |
void | getTransportData () |
Read the transport database. More... | |
void | makePolarCorrections (size_t i, size_t j, doublereal &f_eps, doublereal &f_sigma) |
Corrections for polar-nonpolar binary diffusion coefficients. More... | |
void | fitCollisionIntegrals (MMCollisionInt &integrals) |
Generate polynomial fits to collision integrals. More... | |
virtual void | fitProperties (MMCollisionInt &integrals) |
Generate polynomial fits to the viscosity and conductivity. More... | |
virtual void | fitDiffCoeffs (MMCollisionInt &integrals) |
Generate polynomial fits to the binary diffusion coefficients. More... | |
void | getBinDiffCorrection (doublereal t, MMCollisionInt &integrals, size_t k, size_t j, doublereal xk, doublereal xj, doublereal &fkj, doublereal &fjk) |
Second-order correction to the binary diffusion coefficients. More... | |
Protected Member Functions inherited from Transport | |
void | finalize () |
Enable the transport object for use. More... | |
Protected Attributes inherited from MixTransport | |
vector_fp | m_cond |
vector of species thermal conductivities (W/m /K) More... | |
doublereal | m_lambda |
Internal storage for the calculated mixture thermal conductivity. More... | |
bool | m_spcond_ok |
Update boolean for the species thermal conductivities. More... | |
bool | m_condmix_ok |
Update boolean for the mixture rule for the mixture thermal conductivity. More... | |
Protected Attributes inherited from GasTransport | |
vector_fp | m_molefracs |
Vector of species mole fractions. More... | |
doublereal | m_viscmix |
Internal storage for the viscosity of the mixture (kg /m /s) More... | |
bool | m_visc_ok |
Update boolean for mixture rule for the mixture viscosity. More... | |
bool | m_viscwt_ok |
Update boolean for the weighting factors for the mixture viscosity. More... | |
bool | m_spvisc_ok |
Update boolean for the species viscosities. More... | |
bool | m_bindiff_ok |
Update boolean for the binary diffusivities at unit pressure. More... | |
int | m_mode |
Type of the polynomial fits to temperature. More... | |
DenseMatrix | m_phi |
m_phi is a Viscosity Weighting Function. size = m_nsp * n_nsp More... | |
vector_fp | m_spwork |
work space length = m_kk More... | |
vector_fp | m_visc |
vector of species viscosities (kg /m /s). More... | |
std::vector< vector_fp > | m_visccoeffs |
Polynomial fits to the viscosity of each species. More... | |
vector_fp | m_mw |
Local copy of the species molecular weights. More... | |
DenseMatrix | m_wratjk |
Holds square roots of molecular weight ratios. More... | |
DenseMatrix | m_wratkj1 |
Holds square roots of molecular weight ratios. More... | |
vector_fp | m_sqvisc |
vector of square root of species viscosities sqrt(kg /m /s). More... | |
vector_fp | m_polytempvec |
Powers of the ln temperature, up to fourth order. More... | |
doublereal | m_temp |
Current value of the temperature at which the properties in this object are calculated (Kelvin). More... | |
doublereal | m_kbt |
Current value of Boltzmann constant times the temperature (Joules) More... | |
doublereal | m_sqrt_kbt |
current value of Boltzmann constant times the temperature. More... | |
doublereal | m_sqrt_t |
current value of temperature to 1/2 power More... | |
doublereal | m_logt |
Current value of the log of the temperature. More... | |
doublereal | m_t14 |
Current value of temperature to 1/4 power. More... | |
doublereal | m_t32 |
Current value of temperature to the 3/2 power. More... | |
std::vector< vector_fp > | m_diffcoeffs |
Polynomial fits to the binary diffusivity of each species. More... | |
DenseMatrix | m_bdiff |
Matrix of binary diffusion coefficients at the reference pressure and the current temperature Size is nsp x nsp. More... | |
std::vector< vector_fp > | m_condcoeffs |
temperature fits of the heat conduction More... | |
std::vector< vector_int > | m_poly |
Indices for the (i,j) interaction in collision integral fits. More... | |
std::vector< vector_fp > | m_omega22_poly |
Fit for omega22 collision integral. More... | |
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 reduced temperature. More... | |
std::vector< vector_fp > | m_astar_poly |
Fit for astar collision integral. More... | |
std::vector< vector_fp > | m_bstar_poly |
Fit for bstar collision integral. More... | |
std::vector< vector_fp > | m_cstar_poly |
Fit for cstar collision integral. More... | |
vector_fp | m_zrot |
Rotational relaxation number for each species. More... | |
vector_fp | m_crot |
Dimensionless rotational heat capacity of each species. More... | |
std::vector< bool > | m_polar |
Vector of booleans indicating whether a species is a polar molecule. More... | |
vector_fp | m_alpha |
Polarizability of each species in the phase. More... | |
vector_fp | m_eps |
Lennard-Jones well-depth of the species in the current phase. More... | |
vector_fp | m_sigma |
Lennard-Jones diameter of the species in the current phase. More... | |
DenseMatrix | m_reducedMass |
This is the reduced mass of the interaction between species i and j. More... | |
DenseMatrix | m_diam |
hard-sphere diameter for (i,j) collision More... | |
DenseMatrix | m_epsilon |
The effective well depth for (i,j) collisions. More... | |
DenseMatrix | m_dipole |
The effective dipole moment for (i,j) collisions. More... | |
DenseMatrix | m_delta |
Reduced dipole moment of the interaction between two species. More... | |
vector_fp | m_w_ac |
Pitzer acentric factor. More... | |
vector_fp | m_disp |
Dispersion coefficient. More... | |
vector_fp | m_quad_polar |
Quadrupole polarizability. More... | |
int | m_log_level |
Level of verbose printing during initialization. More... | |
Protected Attributes inherited from Transport | |
ThermoPhase * | m_thermo |
pointer to the object representing the phase More... | |
bool | m_ready |
true if finalize has been called More... | |
size_t | m_nsp |
Number of species. More... | |
size_t | m_nDim |
Number of dimensions used in flux expressions. More... | |
int | m_velocityBasis |
Velocity basis from which diffusion velocities are computed. More... | |
std::weak_ptr< Solution > | m_root |
reference to Solution More... | |
Class UnityLewisTransport implements the unity Lewis number approximation for the mixture-averaged species diffusion coefficients.
Mixture-averaged transport properties for viscosity and thermal conductivity are inherited from the MixTransport class.
Definition at line 25 of file UnityLewisTransport.h.
|
inlinevirtual |
Identifies the Transport object type.
Each derived class should override this method to return a meaningful identifier.
Reimplemented from MixTransport.
Definition at line 30 of file UnityLewisTransport.h.
|
inlinevirtual |
Returns the unity Lewis number approximation based diffusion coefficients [m^2/s].
Returns the unity Lewis number approximation based diffusion coefficients for a gas, appropriate for calculating the mass averaged diffusive flux with respect to the mass averaged velocity using gradients of the mole fraction.
\[ D^\prime_{km} = \frac{\lambda}{\rho c_p} \]
In order to obtain the expected behavior from a unity Lewis number model, this formulation requires that the correction velocity be computed as
\[ V_c = \sum \frac{W_k}{\overline{W}} D^\prime_{km} \nabla X_k \]
[out] | d | Vector of diffusion coefficients for each species (m^2/s). length m_nsp. |
Reimplemented from GasTransport.
Definition at line 56 of file UnityLewisTransport.h.
References ThermoPhase::cp_mass(), Phase::density(), Transport::m_nsp, Transport::m_thermo, and MixTransport::thermalConductivity().
|
inlinevirtual |
Not implemented for unity Lewis number approximation.
Reimplemented from GasTransport.
Definition at line 64 of file UnityLewisTransport.h.
|
inlinevirtual |
Returns the unity Lewis number approximation based diffusion coefficients [m^2/s].
These are the coefficients for calculating the diffusive mass fluxes from the species mass fraction gradients, computed as
\[ D_{km} = \frac{\lambda}{\rho c_p} \]
[out] | d | Vector of diffusion coefficients for each species (m^2/s). length m_nsp. |
Reimplemented from GasTransport.
Definition at line 81 of file UnityLewisTransport.h.
References ThermoPhase::cp_mass(), Phase::density(), Transport::m_nsp, Transport::m_thermo, and MixTransport::thermalConductivity().