Cantera  3.1.0a1
IonGasTransport Class Reference

Class IonGasTransport implements Stockmayer-(n,6,4) model for transport of ions. More...

#include <IonGasTransport.h>

Inheritance diagram for IonGasTransport:
[legend]

Detailed Description

Class IonGasTransport implements Stockmayer-(n,6,4) model for transport of ions.

As implemented here, only binary transport between neutrals and ions is considered for calculating mixture-average diffusion coefficients and mobilities. When polarizability is not provide for an ion, LJ model is used instead of n64 model. Only neutral species are considered for thermal conductivity and viscosity.

References for Stockmayer-(n,6,4) model: Selle and Riedel [36], [37]; Han et al. [11]; Chiflikian [4]; and Viehland et al. [46].

Stockmayer-(n,6,4) model is not suitable for collision between O2/O2- due to resonant charge transfer. Therefore, an experimental collision data is used instead.

Data taken from [34].

Definition at line 34 of file IonGasTransport.h.

Public Member Functions

string transportModel () const override
 Identifies the model represented by this Transport object. More...
 
void init (ThermoPhase *thermo, int mode, int log_level) override
 Initialize a transport manager. More...
 
double viscosity () override
 Viscosity of the mixture (kg/m/s). More...
 
double thermalConductivity () override
 Returns the mixture thermal conductivity (W/m/K). More...
 
void getMobilities (double *const mobi) override
 The mobilities for ions in gas. More...
 
void getMixDiffCoeffs (double *const d) override
 The mixture transport for ionized gas. More...
 
double electricalConductivity () override
 The electrical conductivity (Siemens/m). More...
 
- Public Member Functions inherited from MixTransport
 MixTransport ()=default
 Default constructor. More...
 
string transportModel () const override
 Identifies the model represented by this Transport object. More...
 
void getThermalDiffCoeffs (double *const dt) override
 Return the thermal diffusion coefficients. More...
 
double thermalConductivity () override
 Returns the mixture thermal conductivity (W/m /K) More...
 
void getMobilities (double *const mobil) override
 Get the Electrical mobilities (m^2/V/s). More...
 
void update_T () override
 Update the internal parameters whenever the temperature has changed. More...
 
void update_C () override
 Update the internal parameters whenever the concentrations have changed. More...
 
void getSpeciesFluxes (size_t ndim, const double *const grad_T, size_t ldx, const double *const grad_X, size_t ldf, double *const fluxes) override
 Get the species diffusive mass fluxes wrt to the mass averaged velocity, given the gradients in mole fraction and temperature. More...
 
void init (ThermoPhase *thermo, int mode=0, int log_level=0) override
 Initialize a transport manager. More...
 
- Public Member Functions inherited from GasTransport
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 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...
 
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
 
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 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 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...
 

Protected Member Functions

void setupN64 ()
 setup parameters for n64 model More...
 
void fitDiffCoeffs (MMCollisionInt &integrals) override
 Generate polynomial fits to the binary diffusion coefficients. More...
 
double omega11_n64 (const double tstar, const double gamma)
 Collision integral of omega11 of n64 collision model. More...
 
- 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
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 [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...
 
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_speciesCharge
 electrical properties More...
 
vector< size_t > m_kIon
 index of ions (exclude electron.) More...
 
vector< size_t > m_kNeutral
 index of neutral species More...
 
size_t m_kElectron = npos
 index of electron More...
 
DenseMatrix m_gamma
 parameter of omega11 of n64 More...
 
vector< double > m_om11_O2
 polynomial of the collision integral for O2/O2- More...
 
- Protected Attributes inherited from MixTransport
vector< double > m_cond
 vector of species thermal conductivities (W/m /K) More...
 
double m_lambda = 0.0
 Internal storage for the calculated mixture thermal conductivity. More...
 
bool m_spcond_ok = false
 Update boolean for the species thermal conductivities. More...
 
bool m_condmix_ok = false
 Update boolean for the mixture rule for the mixture thermal conductivity. More...
 
- Protected Attributes inherited from GasTransport
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

◆ transportModel()

string transportModel ( ) const
inlineoverridevirtual

Identifies the model represented by this Transport object.

Each derived class should override this method to return a meaningful identifier.

Since
New in Cantera 3.0. The name returned by this method corresponds to the canonical name used in the YAML input format.

Reimplemented from Transport.

Definition at line 39 of file IonGasTransport.h.

◆ init()

void init ( ThermoPhase thermo,
int  mode,
int  log_level 
)
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 GasTransport.

Definition at line 17 of file IonGasTransport.cpp.

◆ viscosity()

double viscosity ( )
overridevirtual

Viscosity of the mixture (kg/m/s).

Only Neutral species contribute to Viscosity.

Reimplemented from GasTransport.

Definition at line 96 of file IonGasTransport.cpp.

◆ thermalConductivity()

double thermalConductivity ( )
overridevirtual

Returns the mixture thermal conductivity (W/m/K).

Only Neutral species contribute to thermal conductivity.

Reimplemented from Transport.

Definition at line 120 of file IonGasTransport.cpp.

◆ getMobilities()

void getMobilities ( double *const  mobi)
overridevirtual

The mobilities for ions in gas.

The ion mobilities are calculated by Blanc's law.

Reimplemented from Transport.

Definition at line 382 of file IonGasTransport.cpp.

◆ getMixDiffCoeffs()

void getMixDiffCoeffs ( double *const  d)
overridevirtual

The mixture transport for ionized gas.

The binary transport between two charged species is neglected.

Reimplemented from GasTransport.

Definition at line 347 of file IonGasTransport.cpp.

◆ electricalConductivity()

double electricalConductivity ( )
overridevirtual

The electrical conductivity (Siemens/m).

\[ \sigma = \sum_k{\left|C_k\right| \mu_k \frac{X_k P}{k_b T}} \]

Reimplemented from Transport.

Definition at line 139 of file IonGasTransport.cpp.

◆ setupN64()

void setupN64 ( )
protected

setup parameters for n64 model

Definition at line 247 of file IonGasTransport.cpp.

◆ fitDiffCoeffs()

void fitDiffCoeffs ( MMCollisionInt integrals)
overrideprotectedvirtual

Generate polynomial fits to the binary diffusion coefficients.

Use Stockmayer-(n,6,4) model for collision between charged and neutral species.

Reimplemented from GasTransport.

Definition at line 156 of file IonGasTransport.cpp.

◆ omega11_n64()

double omega11_n64 ( const double  tstar,
const double  gamma 
)
protected

Collision integral of omega11 of n64 collision model.

The collision integral was fitted by Han et al. using the table by Viehlan et al. Note: Han release the range to 1000, but Selle suggested that a high temperature model is needed for T* > 10.

Definition at line 315 of file IonGasTransport.cpp.

Member Data Documentation

◆ m_speciesCharge

vector<double> m_speciesCharge
protected

electrical properties

Definition at line 87 of file IonGasTransport.h.

◆ m_kIon

vector<size_t> m_kIon
protected

index of ions (exclude electron.)

Definition at line 90 of file IonGasTransport.h.

◆ m_kNeutral

vector<size_t> m_kNeutral
protected

index of neutral species

Definition at line 93 of file IonGasTransport.h.

◆ m_kElectron

size_t m_kElectron = npos
protected

index of electron

Definition at line 96 of file IonGasTransport.h.

◆ m_gamma

DenseMatrix m_gamma
protected

parameter of omega11 of n64

Definition at line 99 of file IonGasTransport.h.

◆ m_om11_O2

vector<double> m_om11_O2
protected

polynomial of the collision integral for O2/O2-

Definition at line 102 of file IonGasTransport.h.


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