Cantera  3.1.0a1
GasTransport Class Referenceabstract

Class GasTransport implements some functions and properties that are shared by the MixTransport and MultiTransport classes. More...

#include <GasTransport.h>

Inheritance diagram for GasTransport:
[legend]

Detailed Description

Class GasTransport implements some functions and properties that are shared by the MixTransport and MultiTransport classes.

For details, see Kee, et al. [15] and [16].

Definition at line 25 of file GasTransport.h.

Public Member Functions

double viscosity () override
 Viscosity of the mixture (kg /m /s) More...
 
void getSpeciesViscosities (double *const visc) override
 Get the pure-species viscosities. More...
 
void getBinaryDiffCoeffs (const size_t ld, double *const d) override
 Returns the matrix of binary diffusion coefficients. More...
 
void getMixDiffCoeffs (double *const d) override
 Returns the Mixture-averaged diffusion coefficients [m^2/s]. More...
 
void getMixDiffCoeffsMole (double *const d) override
 Returns the mixture-averaged diffusion coefficients [m^2/s]. More...
 
void getMixDiffCoeffsMass (double *const d) override
 Returns the mixture-averaged diffusion coefficients [m^2/s]. More...
 
void getViscosityPolynomial (size_t i, double *coeffs) const override
 Return the polynomial fits to the viscosity of species i. More...
 
void getConductivityPolynomial (size_t i, double *coeffs) const override
 Return the temperature fits of the heat conductivity of species i. More...
 
void getBinDiffusivityPolynomial (size_t i, size_t j, double *coeffs) const override
 Return the polynomial fits to the binary diffusivity of species pair (i, j) More...
 
void getCollisionIntegralPolynomial (size_t i, size_t j, double *astar_coeffs, double *bstar_coeffs, double *cstar_coeffs) const override
 Return the polynomial fits to the collision integral of species pair (i, j) More...
 
void setViscosityPolynomial (size_t i, double *coeffs) override
 Modify the polynomial fits to the viscosity of species i. More...
 
void setConductivityPolynomial (size_t i, double *coeffs) override
 Modify the temperature fits of the heat conductivity of species i. More...
 
void setBinDiffusivityPolynomial (size_t i, size_t j, double *coeffs) override
 Modify the polynomial fits to the binary diffusivity of species pair (i, j) More...
 
void setCollisionIntegralPolynomial (size_t i, size_t j, double *astar_coeffs, double *bstar_coeffs, double *cstar_coeffs, bool actualT) override
 Modify the polynomial fits to the collision integral of species pair (i, j) More...
 
void init (ThermoPhase *thermo, int mode=0, int log_level=0) override
 Initialize a transport manager. More...
 
bool CKMode () const override
 Boolean indicating the form of the transport properties polynomial fits. More...
 
- Public Member Functions inherited from Transport
 Transport ()=default
 Constructor. More...
 
 Transport (const Transport &)=delete
 
Transportoperator= (const Transport &)=delete
 
virtual string transportModel () const
 Identifies the model represented by this Transport object. More...
 
ThermoPhasethermo ()
 Phase object. 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 void getSpeciesFluxes (size_t ndim, const double *const grad_T, size_t ldx, const double *const grad_X, size_t ldf, double *const fluxes)
 Get the species diffusive mass fluxes wrt to the specified solution averaged velocity, given the gradients in mole fraction and temperature. More...
 
virtual void getMolarFluxes (const double *const state1, const double *const state2, const double delta, double *const cfluxes)
 Get the molar fluxes [kmol/m^2/s], given the thermodynamic state at two nearby points. More...
 
virtual void getMassFluxes (const double *state1, const double *state2, double delta, double *mfluxes)
 Get the mass fluxes [kg/m^2/s], given the thermodynamic state at two nearby points. More...
 
virtual void getThermalDiffCoeffs (double *const dt)
 Return a vector of Thermal diffusion coefficients [kg/m/sec]. More...
 
virtual void getMultiDiffCoeffs (const size_t ld, double *const d)
 Return the Multicomponent diffusion coefficients. Units: [m^2/s]. More...
 
AnyMap parameters () const
 Return the parameters for a phase definition which are needed to reconstruct an identical object using the newTransport function. More...
 
virtual double bulkViscosity ()
 The bulk viscosity in Pa-s. More...
 
virtual double thermalConductivity ()
 Returns the mixture thermal conductivity in W/m/K. More...
 
virtual double electricalConductivity ()
 The electrical conductivity (Siemens/m). More...
 
virtual void getMobilities (double *const mobil_e)
 Get the Electrical mobilities (m^2/V/s). More...
 

Protected Member Functions

virtual void update_T ()
 
virtual void update_C ()=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...
 
Initialization
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 [27]. More...
 
void getTransportData ()
 Read the transport database. More...
 
void makePolarCorrections (size_t i, size_t j, double &f_eps, double &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 \( \eta \) and conductivity \( \lambda \). More...
 
virtual void fitDiffCoeffs (MMCollisionInt &integrals)
 Generate polynomial fits to the binary diffusion coefficients. More...
 
void getBinDiffCorrection (double t, MMCollisionInt &integrals, size_t k, size_t j, double xk, double xj, double &fkj, double &fjk)
 Second-order correction to the binary diffusion coefficients. More...
 

Protected Attributes

vector< double > m_molefracs
 Vector of species mole fractions. More...
 
double m_viscmix = 0.0
 Internal storage for the viscosity of the mixture (kg /m /s) More...
 
bool m_visc_ok = false
 Update boolean for mixture rule for the mixture viscosity. More...
 
bool m_viscwt_ok = false
 Update boolean for the weighting factors for the mixture viscosity. More...
 
bool m_spvisc_ok = false
 Update boolean for the species viscosities. More...
 
bool m_bindiff_ok = false
 Update boolean for the binary diffusivities at unit pressure. More...
 
int m_mode = 0
 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< double > m_spwork
 work space length = m_kk More...
 
vector< double > m_visc
 vector of species viscosities (kg /m /s). More...
 
vector< vector< double > > m_visccoeffs
 Polynomial fits to the viscosity of each species. More...
 
vector< double > 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< double > m_sqvisc
 vector of square root of species viscosities sqrt(kg /m /s). More...
 
vector< double > m_polytempvec
 Powers of the ln temperature, up to fourth order. More...
 
double m_temp = -1.0
 Current value of the temperature at which the properties in this object are calculated (Kelvin). More...
 
double m_kbt = 0.0
 Current value of Boltzmann constant times the temperature (Joules) More...
 
double m_sqrt_t = 0.0
 current value of temperature to 1/2 power More...
 
double m_logt = 0.0
 Current value of the log of the temperature. More...
 
double m_t14 = 0.0
 Current value of temperature to 1/4 power. More...
 
vector< vector< double > > 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...
 
vector< vector< double > > m_condcoeffs
 temperature fits of the heat conduction More...
 
vector< vector< int > > m_poly
 Indices for the (i,j) interaction in collision integral fits. More...
 
vector< vector< double > > m_omega22_poly
 Fit for omega22 collision integral. More...
 
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...
 
vector< vector< double > > m_astar_poly
 Fit for astar collision integral. More...
 
vector< vector< double > > m_bstar_poly
 Fit for bstar collision integral. More...
 
vector< vector< double > > m_cstar_poly
 Fit for cstar collision integral. More...
 
vector< double > m_zrot
 Rotational relaxation number for each species. More...
 
vector< double > m_crot
 Dimensionless rotational heat capacity of each species. More...
 
vector< bool > m_polar
 Vector of booleans indicating whether a species is a polar molecule. More...
 
vector< double > m_alpha
 Polarizability of each species in the phase. More...
 
vector< double > m_eps
 Lennard-Jones well-depth of the species in the current phase. More...
 
vector< double > 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< double > m_w_ac
 Pitzer acentric factor. More...
 
vector< double > m_disp
 Dispersion coefficient. More...
 
vector< double > m_quad_polar
 Quadrupole polarizability. More...
 
int m_log_level = 0
 Level of verbose printing during initialization. More...
 
- Protected Attributes inherited from Transport
ThermoPhasem_thermo
 pointer to the object representing the phase More...
 
size_t m_nsp = 0
 Number of species. More...
 

Member Function Documentation

◆ viscosity()

double viscosity ( )
overridevirtual

Viscosity of the mixture (kg /m /s)

The viscosity is computed using the Wilke mixture rule (kg /m /s)

\[ \mu = \sum_k \frac{\mu_k X_k}{\sum_j \Phi_{k,j} X_j}. \]

Here \( \mu_k \) is the viscosity of pure species k, and

\[ \Phi_{k,j} = \frac{\left[1 + \sqrt{\left(\frac{\mu_k}{\mu_j}\sqrt{\frac{M_j}{M_k}}\right)}\right]^2} {\sqrt{8}\sqrt{1 + M_k/M_j}} \]

Returns
the viscosity of the mixture (units = Pa s = kg /m /s)
See also
updateViscosity_T()

Reimplemented from Transport.

Reimplemented in IonGasTransport, and HighPressureGasTransport.

Definition at line 60 of file GasTransport.cpp.

◆ getSpeciesViscosities()

void getSpeciesViscosities ( double *const  visc)
inlineoverridevirtual

Get the pure-species viscosities.

Reimplemented from Transport.

Definition at line 51 of file GasTransport.h.

◆ getBinaryDiffCoeffs()

void getBinaryDiffCoeffs ( const size_t  ld,
double *const  d 
)
overridevirtual

Returns the matrix of binary diffusion coefficients.

d[ld*j + i] = rp * m_bdiff(i,j);

Parameters
ldoffset of rows in the storage
doutput vector of diffusion coefficients. Units of m**2 / s

Reimplemented from Transport.

Reimplemented in HighPressureGasTransport.

Definition at line 149 of file GasTransport.cpp.

◆ getMixDiffCoeffs()

void getMixDiffCoeffs ( double *const  d)
overridevirtual

Returns the Mixture-averaged diffusion coefficients [m^2/s].

Returns the mixture averaged 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. Note, for the single species case or the pure fluid case the routine returns the self-diffusion coefficient. This is needed to avoid a Nan result in the formula below.

This is Eqn. 12.180 from "Chemically Reacting Flow"

\[ D_{km}' = \frac{\left( \bar{M} - X_k M_k \right)}{ \bar{\qquad M \qquad } } {\left( \sum_{j \ne k} \frac{X_j}{D_{kj}} \right) }^{-1} \]

Parameters
[out]dVector of mixture diffusion coefficients, \( D_{km}' \) , for each species (m^2/s). length m_nsp

Reimplemented from Transport.

Reimplemented in UnityLewisTransport, and IonGasTransport.

Definition at line 167 of file GasTransport.cpp.

◆ getMixDiffCoeffsMole()

void getMixDiffCoeffsMole ( double *const  d)
overridevirtual

Returns the mixture-averaged diffusion coefficients [m^2/s].

These are the coefficients for calculating the molar diffusive fluxes from the species mole fraction gradients, computed according to Eq. 12.176 in "Chemically Reacting Flow":

\[ D_{km}^* = \frac{1-X_k}{\sum_{j \ne k}^K X_j/\mathcal{D}_{kj}} \]

Parameters
[out]dvector of mixture-averaged diffusion coefficients for each species, length m_nsp.

Reimplemented from Transport.

Reimplemented in UnityLewisTransport.

Definition at line 198 of file GasTransport.cpp.

◆ getMixDiffCoeffsMass()

void getMixDiffCoeffsMass ( double *const  d)
overridevirtual

Returns the mixture-averaged diffusion coefficients [m^2/s].

These are the coefficients for calculating the diffusive mass fluxes from the species mass fraction gradients, computed according to Eq. 12.178 in "Chemically Reacting Flow":

\[ \frac{1}{D_{km}} = \sum_{j \ne k}^K \frac{X_j}{\mathcal{D}_{kj}} + \frac{X_k}{1-Y_k} \sum_{j \ne k}^K \frac{Y_j}{\mathcal{D}_{kj}} \]

Parameters
[out]dvector of mixture-averaged diffusion coefficients for each species, length m_nsp.

Reimplemented from Transport.

Reimplemented in UnityLewisTransport.

Definition at line 228 of file GasTransport.cpp.

◆ getViscosityPolynomial()

void getViscosityPolynomial ( size_t  i,
double *  coeffs 
) const
overridevirtual

Return the polynomial fits to the viscosity of species i.

See also
fitProperties()

Reimplemented from Transport.

Definition at line 807 of file GasTransport.cpp.

◆ getConductivityPolynomial()

void getConductivityPolynomial ( size_t  i,
double *  coeffs 
) const
overridevirtual

Return the temperature fits of the heat conductivity of species i.

See also
fitProperties()

Reimplemented from Transport.

Definition at line 814 of file GasTransport.cpp.

◆ getBinDiffusivityPolynomial()

void getBinDiffusivityPolynomial ( size_t  i,
size_t  j,
double *  coeffs 
) const
overridevirtual

Return the polynomial fits to the binary diffusivity of species pair (i, j)

See also
fitDiffCoeffs()

Reimplemented from Transport.

Definition at line 821 of file GasTransport.cpp.

◆ getCollisionIntegralPolynomial()

void getCollisionIntegralPolynomial ( size_t  i,
size_t  j,
double *  astar_coeffs,
double *  bstar_coeffs,
double *  cstar_coeffs 
) const
overridevirtual

Return the polynomial fits to the collision integral of species pair (i, j)

See also
fitCollisionIntegrals()

Reimplemented from Transport.

Definition at line 836 of file GasTransport.cpp.

◆ setViscosityPolynomial()

void setViscosityPolynomial ( size_t  i,
double *  coeffs 
)
overridevirtual

Modify the polynomial fits to the viscosity of species i.

See also
fitProperties()

Reimplemented from Transport.

Definition at line 848 of file GasTransport.cpp.

◆ setConductivityPolynomial()

void setConductivityPolynomial ( size_t  i,
double *  coeffs 
)
overridevirtual

Modify the temperature fits of the heat conductivity of species i.

See also
fitProperties()

Reimplemented from Transport.

Definition at line 861 of file GasTransport.cpp.

◆ setBinDiffusivityPolynomial()

void setBinDiffusivityPolynomial ( size_t  i,
size_t  j,
double *  coeffs 
)
overridevirtual

Modify the polynomial fits to the binary diffusivity of species pair (i, j)

See also
fitDiffCoeffs()

Reimplemented from Transport.

Definition at line 874 of file GasTransport.cpp.

◆ setCollisionIntegralPolynomial()

void setCollisionIntegralPolynomial ( size_t  i,
size_t  j,
double *  astar_coeffs,
double *  bstar_coeffs,
double *  cstar_coeffs,
bool  actualT 
)
overridevirtual

Modify the polynomial fits to the collision integral of species pair (i, j)

See also
fitCollisionIntegrals()

Reimplemented from Transport.

Definition at line 895 of file GasTransport.cpp.

◆ init()

void init ( ThermoPhase thermo,
int  mode = 0,
int  log_level = 0 
)
overridevirtual

Initialize a transport manager.

This routine sets up a transport manager. It calculates the collision integrals and populates species-dependent data structures.

Parameters
thermoPointer to the ThermoPhase object
modeChemkin compatible mode or not. This alters the specification of the collision integrals. defaults to no.
log_levelDefaults to zero, no logging

Reimplemented from Transport.

Reimplemented in MultiTransport, MixTransport, and IonGasTransport.

Definition at line 261 of file GasTransport.cpp.

◆ CKMode()

bool CKMode ( ) const
inlineoverridevirtual

Boolean indicating the form of the transport properties polynomial fits.

Returns true if the Chemkin form is used.

Reimplemented from Transport.

Definition at line 153 of file GasTransport.h.

◆ updateViscosity_T()

void updateViscosity_T ( )
protectedvirtual

Update the temperature-dependent viscosity terms.

Updates the array of pure species viscosities, and the weighting functions in the viscosity mixture rule. The flag m_visc_ok is set to true.

The formula for the weighting function is from Poling et al. [33], Eq. (9-5.14):

\[ \phi_{ij} = \frac{ \left[ 1 + \left( \mu_i / \mu_j \right)^{1/2} \left( M_j / M_i \right)^{1/4} \right]^2 } {\left[ 8 \left( 1 + M_i / M_j \right) \right]^{1/2}} \]

Definition at line 84 of file GasTransport.cpp.

◆ updateSpeciesViscosities()

void updateSpeciesViscosities ( )
protectedvirtual

Update the pure-species viscosities.

These are evaluated from the polynomial fits of the temperature and are assumed to be independent of pressure.

Definition at line 105 of file GasTransport.cpp.

◆ updateDiff_T()

void updateDiff_T ( )
protectedvirtual

Update the binary diffusion coefficients.

These are evaluated from the polynomial fits of the temperature at the unit pressure of 1 Pa.

Definition at line 123 of file GasTransport.cpp.

◆ setupCollisionParameters()

void setupCollisionParameters ( )
protectedvirtual

Setup parameters for a new kinetic-theory-based transport manager for low-density gases.

Definition at line 293 of file GasTransport.cpp.

◆ setupCollisionIntegral()

void setupCollisionIntegral ( )
protected

Setup range for polynomial fits to collision integrals of Monchick & Mason [27].

Definition at line 354 of file GasTransport.cpp.

◆ getTransportData()

void getTransportData ( )
protected

Read the transport database.

Read transport property data from a file for a list of species. Given the name of a file containing transport property parameters and a list of species names.

Definition at line 384 of file GasTransport.cpp.

◆ makePolarCorrections()

void makePolarCorrections ( size_t  i,
size_t  j,
double &  f_eps,
double &  f_sigma 
)
protected

Corrections for polar-nonpolar binary diffusion coefficients.

Calculate corrections to the well depth parameter and the diameter for use in computing the binary diffusion coefficient of polar-nonpolar pairs. For more information about this correction, see Dixon-Lewis [6].

Parameters
iSpecies one - this is a bimolecular correction routine
jspecies two - this is a bimolecular correction routine
f_epsMultiplicative correction factor to be applied to epsilon(i,j)
f_sigmaMultiplicative correction factor to be applied to diam(i,j)

Definition at line 421 of file GasTransport.cpp.

◆ fitCollisionIntegrals()

void fitCollisionIntegrals ( MMCollisionInt integrals)
protected

Generate polynomial fits to collision integrals.

Parameters
integralsinterpolator for the collision integrals

Definition at line 446 of file GasTransport.cpp.

◆ fitProperties()

void fitProperties ( MMCollisionInt integrals)
protectedvirtual

Generate polynomial fits to the viscosity \( \eta \) and conductivity \( \lambda \).

If CK_mode, then the fits are of the form

\[ \ln \eta(i) = \sum_{n=0}^3 a_n(i) \, (\ln T)^n \]

and

\[ \ln \lambda(i) = \sum_{n=0}^3 b_n(i) \, (\ln T)^n \]

Otherwise the fits are of the form

\[ \left(\eta(i)\right)^{1/2} = T^{1/4} \sum_{n=0}^4 a_n(i) \, (\ln T)^n \]

and

\[ \lambda(i) = T^{1/2} \sum_{n=0}^4 b_n(i) \, (\ln T)^n \]

Parameters
integralsinterpolator for the collision integrals

Definition at line 497 of file GasTransport.cpp.

◆ fitDiffCoeffs()

void fitDiffCoeffs ( MMCollisionInt integrals)
protectedvirtual

Generate polynomial fits to the binary diffusion coefficients.

If CK_mode, then the fits are of the form

\[ \ln D(i,j) = \sum_{n=0}^3 c_n(i,j) \, (\ln T)^n \]

Otherwise the fits are of the form

\[ D(i,j) = T^{3/2} \sum_{n=0}^4 c_n(i,j) \, (\ln T)^n \]

Parameters
integralsinterpolator for the collision integrals

Reimplemented in IonGasTransport.

Definition at line 676 of file GasTransport.cpp.

◆ getBinDiffCorrection()

void getBinDiffCorrection ( double  t,
MMCollisionInt integrals,
size_t  k,
size_t  j,
double  xk,
double  xj,
double &  fkj,
double &  fjk 
)
protected

Second-order correction to the binary diffusion coefficients.

Calculate second-order corrections to binary diffusion coefficient pair (dkj, djk). At first order, the binary diffusion coefficients are independent of composition, and d(k,j) = d(j,k). But at second order, there is a weak dependence on composition, with the result that d(k,j) != d(j,k). This method computes the multiplier by which the first-order binary diffusion coefficient should be multiplied to produce the value correct to second order. The expressions here are taken from Marerro and Mason [24].

Parameters
tTemperature (K)
integralsinterpolator for the collision integrals
kindex of first species
jindex of second species
xkMole fraction of species k
xjMole fraction of species j
fkjmultiplier for d(k,j)
fjkmultiplier for d(j,k)
Note
This method is not used currently.

Definition at line 755 of file GasTransport.cpp.

Member Data Documentation

◆ m_molefracs

vector<double> m_molefracs
protected

Vector of species mole fractions.

These are processed so that all mole fractions are >= Tiny. Length = m_kk.

Definition at line 297 of file GasTransport.h.

◆ m_viscmix

double m_viscmix = 0.0
protected

Internal storage for the viscosity of the mixture (kg /m /s)

Definition at line 300 of file GasTransport.h.

◆ m_visc_ok

bool m_visc_ok = false
protected

Update boolean for mixture rule for the mixture viscosity.

Definition at line 303 of file GasTransport.h.

◆ m_viscwt_ok

bool m_viscwt_ok = false
protected

Update boolean for the weighting factors for the mixture viscosity.

Definition at line 306 of file GasTransport.h.

◆ m_spvisc_ok

bool m_spvisc_ok = false
protected

Update boolean for the species viscosities.

Definition at line 309 of file GasTransport.h.

◆ m_bindiff_ok

bool m_bindiff_ok = false
protected

Update boolean for the binary diffusivities at unit pressure.

Definition at line 312 of file GasTransport.h.

◆ m_mode

int m_mode = 0
protected

Type of the polynomial fits to temperature.

CK_Mode means Chemkin mode. Any other value means to use Cantera's preferred fitting functions.

Definition at line 316 of file GasTransport.h.

◆ m_phi

DenseMatrix m_phi
protected

m_phi is a Viscosity Weighting Function. size = m_nsp * n_nsp

Definition at line 319 of file GasTransport.h.

◆ m_spwork

vector<double> m_spwork
protected

work space length = m_kk

Definition at line 322 of file GasTransport.h.

◆ m_visc

vector<double> m_visc
protected

vector of species viscosities (kg /m /s).

These are used in Wilke's rule to calculate the viscosity of the solution. length = m_kk.

Definition at line 326 of file GasTransport.h.

◆ m_visccoeffs

vector<vector<double> > m_visccoeffs
protected

Polynomial fits to the viscosity of each species.

m_visccoeffs[k] is the vector of polynomial coefficients for species k that fits the viscosity as a function of temperature.

Definition at line 331 of file GasTransport.h.

◆ m_mw

vector<double> m_mw
protected

Local copy of the species molecular weights.

Definition at line 334 of file GasTransport.h.

◆ m_wratjk

DenseMatrix m_wratjk
protected

Holds square roots of molecular weight ratios.

m_wratjk(j,k) = sqrt(mw[j]/mw[k]) j < k
m_wratjk(k,j) = sqrt(sqrt(mw[j]/mw[k])) j < k
DenseMatrix m_wratjk
Holds square roots of molecular weight ratios.
Definition: GasTransport.h:343

Definition at line 343 of file GasTransport.h.

◆ m_wratkj1

DenseMatrix m_wratkj1
protected

Holds square roots of molecular weight ratios.

m_wratjk1(j,k) = sqrt(1.0 + mw[k]/mw[j]) j < k

Definition at line 349 of file GasTransport.h.

◆ m_sqvisc

vector<double> m_sqvisc
protected

vector of square root of species viscosities sqrt(kg /m /s).

These are used in Wilke's rule to calculate the viscosity of the solution. length = m_kk.

Definition at line 354 of file GasTransport.h.

◆ m_polytempvec

vector<double> m_polytempvec
protected

Powers of the ln temperature, up to fourth order.

Definition at line 357 of file GasTransport.h.

◆ m_temp

double m_temp = -1.0
protected

Current value of the temperature at which the properties in this object are calculated (Kelvin).

Definition at line 361 of file GasTransport.h.

◆ m_kbt

double m_kbt = 0.0
protected

Current value of Boltzmann constant times the temperature (Joules)

Definition at line 364 of file GasTransport.h.

◆ m_sqrt_t

double m_sqrt_t = 0.0
protected

current value of temperature to 1/2 power

Definition at line 367 of file GasTransport.h.

◆ m_logt

double m_logt = 0.0
protected

Current value of the log of the temperature.

Definition at line 370 of file GasTransport.h.

◆ m_t14

double m_t14 = 0.0
protected

Current value of temperature to 1/4 power.

Definition at line 373 of file GasTransport.h.

◆ m_diffcoeffs

vector<vector<double> > m_diffcoeffs
protected

Polynomial fits to the binary diffusivity of each species.

m_diffcoeff[ic] is vector of polynomial coefficients for species i species j that fits the binary diffusion coefficient. The relationship between i j and ic is determined from the following algorithm:

 int ic = 0;
 for (i = 0; i < m_nsp; i++) {
    for (j = i; j < m_nsp; j++) {
      ic++;
    }
 }

Definition at line 388 of file GasTransport.h.

◆ m_bdiff

DenseMatrix m_bdiff
protected

Matrix of binary diffusion coefficients at the reference pressure and the current temperature Size is nsp x nsp.

Definition at line 392 of file GasTransport.h.

◆ m_condcoeffs

vector<vector<double> > m_condcoeffs
protected

temperature fits of the heat conduction

Dimensions are number of species (nsp) polynomial order of the collision integral fit (degree+1).

Definition at line 399 of file GasTransport.h.

◆ m_poly

vector<vector<int> > m_poly
protected

Indices for the (i,j) interaction in collision integral fits.

m_poly[i][j] contains the index for (i,j) interactions in m_omega22_poly, m_astar_poly, m_bstar_poly, and m_cstar_poly.

Definition at line 406 of file GasTransport.h.

◆ m_omega22_poly

vector<vector<double> > m_omega22_poly
protected

Fit for omega22 collision integral.

m_omega22_poly[m_poly[i][j]] is the vector of polynomial coefficients (length degree+1) for the collision integral fit for the species pair (i,j).

Definition at line 414 of file GasTransport.h.

◆ m_star_poly_uses_actualT

vector<vector<int> > m_star_poly_uses_actualT
protected

Flag to indicate for which (i,j) interaction pairs the actual temperature is used instead of the reduced temperature.

Definition at line 418 of file GasTransport.h.

◆ m_astar_poly

vector<vector<double> > m_astar_poly
protected

Fit for astar collision integral.

m_astar_poly[m_poly[i][j]] is the vector of polynomial coefficients (length degree+1) for the collision integral fit for the species pair (i,j).

Definition at line 426 of file GasTransport.h.

◆ m_bstar_poly

vector<vector<double> > m_bstar_poly
protected

Fit for bstar collision integral.

m_bstar_poly[m_poly[i][j]] is the vector of polynomial coefficients (length degree+1) for the collision integral fit for the species pair (i,j).

Definition at line 434 of file GasTransport.h.

◆ m_cstar_poly

vector<vector<double> > m_cstar_poly
protected

Fit for cstar collision integral.

m_bstar_poly[m_poly[i][j]] is the vector of polynomial coefficients (length degree+1) for the collision integral fit for the species pair (i,j).

Definition at line 442 of file GasTransport.h.

◆ m_zrot

vector<double> m_zrot
protected

Rotational relaxation number for each species.

length is the number of species in the phase. units are dimensionless

Definition at line 448 of file GasTransport.h.

◆ m_crot

vector<double> m_crot
protected

Dimensionless rotational heat capacity of each species.

These values are 0, 1 and 1.5 for single-molecule, linear, and nonlinear species respectively length is the number of species in the phase. Dimensionless (Cr / R)

Definition at line 456 of file GasTransport.h.

◆ m_polar

vector<bool> m_polar
protected

Vector of booleans indicating whether a species is a polar molecule.

Length is nsp

Definition at line 462 of file GasTransport.h.

◆ m_alpha

vector<double> m_alpha
protected

Polarizability of each species in the phase.

Length = nsp. Units = m^3

Definition at line 468 of file GasTransport.h.

◆ m_eps

vector<double> m_eps
protected

Lennard-Jones well-depth of the species in the current phase.

length is the number of species in the phase. Units are Joules (Note this is not Joules/kmol) (note, no kmol -> this is a per molecule amount)

Definition at line 475 of file GasTransport.h.

◆ m_sigma

vector<double> m_sigma
protected

Lennard-Jones diameter of the species in the current phase.

length is the number of species in the phase. units are in meters.

Definition at line 481 of file GasTransport.h.

◆ m_reducedMass

DenseMatrix m_reducedMass
protected

This is the reduced mass of the interaction between species i and j.

reducedMass(i,j) = mw[i] * mw[j] / (Avogadro * (mw[i] + mw[j]));

Units are kg (note, no kmol -> this is a per molecule amount)

Length nsp * nsp. This is a symmetric matrix

Definition at line 491 of file GasTransport.h.

◆ m_diam

DenseMatrix m_diam
protected

hard-sphere diameter for (i,j) collision

diam(i,j) = 0.5*(sigma[i] + sigma[j]); Units are m (note, no kmol -> this is a per molecule amount)

Length nsp * nsp. This is a symmetric matrix.

Definition at line 500 of file GasTransport.h.

◆ m_epsilon

DenseMatrix m_epsilon
protected

The effective well depth for (i,j) collisions.

epsilon(i,j) = sqrt(eps[i]*eps[j]);
Units are Joules (note, no kmol -> this is a per molecule amount)

Length nsp * nsp. This is a symmetric matrix.

Definition at line 509 of file GasTransport.h.

◆ m_dipole

DenseMatrix m_dipole
protected

The effective dipole moment for (i,j) collisions.

Given dipoleMoment in Debye (a Debye is 3.335e-30 C-m):

dipole(i,i) = 1.e-21 / lightSpeed * dipoleMoment; dipole(i,j) = sqrt(dipole(i,i) * dipole(j,j)); (note, no kmol -> this is a per molecule amount)

Length nsp * nsp. This is a symmetric matrix.

Definition at line 521 of file GasTransport.h.

◆ m_delta

DenseMatrix m_delta
protected

Reduced dipole moment of the interaction between two species.

This is the reduced dipole moment of the interaction between two species 0.5 * dipole(i,j)^2 / (4 * Pi * epsilon_0 * epsilon(i,j) * d^3);

Length nsp * nsp .This is a symmetric matrix

Definition at line 530 of file GasTransport.h.

◆ m_w_ac

vector<double> m_w_ac
protected

Pitzer acentric factor.

Length is the number of species in the phase. Dimensionless.

Definition at line 536 of file GasTransport.h.

◆ m_disp

vector<double> m_disp
protected

Dispersion coefficient.

Definition at line 539 of file GasTransport.h.

◆ m_quad_polar

vector<double> m_quad_polar
protected

Quadrupole polarizability.

Definition at line 542 of file GasTransport.h.

◆ m_log_level

int m_log_level = 0
protected

Level of verbose printing during initialization.

Definition at line 545 of file GasTransport.h.


The documentation for this class was generated from the following files: