Cantera  3.2.0a2
Loading...
Searching...
No Matches
ChungHighPressureGasTransport Class Reference

Transport properties for high pressure gas mixtures using the Chung method for viscosity and thermal conductivity. More...

#include <HighPressureGasTransport.h>

Inheritance diagram for ChungHighPressureGasTransport:
[legend]

Detailed Description

Transport properties for high pressure gas mixtures using the Chung method for viscosity and thermal conductivity.

The implementation employs a method of corresponding states, using the Takahashi [48] approach for binary diffusion coefficients (using mixture averaging rules for the mixture properties), and the Chung method for the viscosity and thermal conductivity of a high-pressure gas mixture. All methods are described in Poling et al. [38] (viscosity in Ch. 9, thermal conductivity in Ch. 10, and diffusion coefficients in Ch. 11).

Note
All equations that are cited in this implementation are from the 5th edition of the book "The Properties of Gases and Liquids" by Poling, Prausnitz, and O'Connell.

Definition at line 826 of file HighPressureGasTransport.h.

Public Member Functions

string transportModel () const override
 Identifies the model represented by this Transport object.
 
void init (ThermoPhase *thermo, int mode=0) override
 Initialize a transport manager.
 
double viscosity () override
 Returns the high-pressure mixture viscosity [Pa·s] using the Chung method.
 
double thermalConductivity () override
 Calculates the high-pressure mixture thermal conductivity using the Chung method.
 
- Public Member Functions inherited from HighPressureGasTransportBase
void getBinaryDiffCoeffs (const size_t ld, double *const d) override
 Computes the matrix of binary diffusion coefficients [m²/s] using the Takahashi correction factor.
 
void getMixDiffCoeffs (double *const d) override
 Returns the mixture-averaged diffusion coefficients [m²/s].
 
void getMixDiffCoeffsMole (double *const d) override
 Returns the mixture-averaged diffusion coefficients [m²/s].
 
void getMixDiffCoeffsMass (double *const d) override
 Returns the mixture-averaged diffusion coefficients [m²/s].
 
virtual void updateCorrectionFactors ()
 Updates the matrix of species-pair Takahashi correction factors for use in computing the binary diffusion coefficients, see takahashiCorrectionFactor()
 
- Public Member Functions inherited from MixTransport
 MixTransport ()=default
 Default constructor.
 
string transportModel () const override
 Identifies the model represented by this Transport object.
 
void getThermalDiffCoeffs (double *const dt) override
 Return the thermal diffusion coefficients [kg/m/s].
 
double thermalConductivity () override
 Returns the mixture thermal conductivity [W/m/K].
 
void getMobilities (double *const mobil) override
 Get the electrical mobilities [m²/V/s].
 
void update_T () override
 Update the internal parameters whenever the temperature has changed.
 
void update_C () override
 Update the internal parameters whenever the concentrations have changed.
 
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 [kg/m²/s] with respect to the mass averaged velocity, given the gradients in mole fraction and temperature.
 
void init (ThermoPhase *thermo, int mode=0) override
 Initialize a transport manager.
 
- Public Member Functions inherited from GasTransport
double viscosity () override
 Get the viscosity [Pa·s] of the mixture.
 
void getSpeciesViscosities (double *const visc) override
 Get the pure-species viscosities [Pa·s].
 
void getBinaryDiffCoeffs (const size_t ld, double *const d) override
 Returns the matrix of binary diffusion coefficients [m²/s].
 
void getMixDiffCoeffs (double *const d) override
 Returns the Mixture-averaged diffusion coefficients [m²/s].
 
void getMixDiffCoeffsMole (double *const d) override
 Returns the mixture-averaged diffusion coefficients [m²/s].
 
void getMixDiffCoeffsMass (double *const d) override
 Returns the mixture-averaged diffusion coefficients [m²/s].
 
void getViscosityPolynomial (size_t i, double *coeffs) const override
 Return the polynomial fits to the viscosity of species i.
 
void getConductivityPolynomial (size_t i, double *coeffs) const override
 Return the temperature fits of the heat conductivity of species i.
 
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)
 
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)
 
void setViscosityPolynomial (size_t i, double *coeffs) override
 Modify the polynomial fits to the viscosity of species i
 
void setConductivityPolynomial (size_t i, double *coeffs) override
 Modify the temperature fits of the heat conductivity of species i
 
void setBinDiffusivityPolynomial (size_t i, size_t j, double *coeffs) override
 Modify the polynomial fits to the binary diffusivity of species pair (i, j)
 
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)
 
void init (ThermoPhase *thermo, int mode=0) override
 Initialize a transport manager.
 
bool CKMode () const override
 Boolean indicating the form of the transport properties polynomial fits.
 
void invalidateCache () override
 Invalidate any cached values which are normally updated only when a change in state is detected.
 
- Public Member Functions inherited from Transport
 Transport ()=default
 Constructor.
 
 Transport (const Transport &)=delete
 
Transportoperator= (const Transport &)=delete
 
virtual string transportModel () const
 Identifies the model represented by this Transport object.
 
ThermoPhasethermo ()
 Phase object.
 
void checkSpeciesIndex (size_t k) const
 Check that the specified species index is in range.
 
void checkSpeciesArraySize (size_t kk) const
 Check that an array size is at least m_nsp.
 
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 [kg/m²/s] with respect to the specified solution averaged velocity, given the mole fraction and temperature gradients.
 
virtual void getMolarFluxes (const double *const state1, const double *const state2, const double delta, double *const cfluxes)
 Get the molar fluxes [kmol/m²/s], given the thermodynamic state at two nearby points.
 
virtual void getMassFluxes (const double *state1, const double *state2, double delta, double *mfluxes)
 Get the mass fluxes [kg/m²/s], given the thermodynamic state at two nearby points.
 
virtual void getThermalDiffCoeffs (double *const dt)
 Return a vector of thermal diffusion coefficients [kg/m/s].
 
virtual void getBinaryDiffCoeffs (const size_t ld, double *const d)
 Returns the matrix of binary diffusion coefficients [m²/s].
 
virtual void getMultiDiffCoeffs (const size_t ld, double *const d)
 Return the multicomponent diffusion coefficients [m²/s].
 
virtual void getMixDiffCoeffs (double *const d)
 Return a vector of mixture averaged diffusion coefficients [m²/s].
 
virtual void getMixDiffCoeffsMole (double *const d)
 Returns a vector of mixture averaged diffusion coefficients [m²/s].
 
virtual void getMixDiffCoeffsMass (double *const d)
 Returns a vector of mixture averaged diffusion coefficients [m²/s].
 
virtual void getViscosityPolynomial (size_t i, double *coeffs) const
 Return the polynomial fits to the viscosity of species i.
 
virtual void getConductivityPolynomial (size_t i, double *coeffs) const
 Return the temperature fits of the heat conductivity of species i.
 
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)
 
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)
 
virtual void setViscosityPolynomial (size_t i, double *coeffs)
 Modify the polynomial fits to the viscosity of species i
 
virtual void setConductivityPolynomial (size_t i, double *coeffs)
 Modify the temperature fits of the heat conductivity of species i
 
virtual void setBinDiffusivityPolynomial (size_t i, size_t j, double *coeffs)
 Modify the polynomial fits to the binary diffusivity of species pair (i, j)
 
virtual void setCollisionIntegralPolynomial (size_t i, size_t j, double *astar_coeffs, double *bstar_coeffs, double *cstar_coeffs, bool flag)
 Modify the polynomial fits to the collision integral of species pair (i, j)
 
AnyMap parameters () const
 Return the parameters for a phase definition which are needed to reconstruct an identical object using the newTransport() function.
 
AnyMap fittingErrors () const
 Get error metrics about any functional fits calculated for pure species transport properties.
 
virtual void invalidateCache ()
 Invalidate any cached values which are normally updated only when a change in state is detected.
 
virtual double bulkViscosity ()
 The bulk viscosity [Pa·s].
 
virtual double electricalConductivity ()
 Get the electrical conductivity [siemens/m].
 

Protected Member Functions

 ChungHighPressureGasTransport ()=default
 default constructor
 
void initializePureFluidProperties ()
 Computes and stores pure-fluid specific properties each species.
 
void computeMixtureParameters ()
 Computes the composition-dependent values of the parameters needed for the Chung viscosity model.
 
double lowPressureViscosity (double T, double T_star, double MW, double acentric_factor, double mu_r, double sigma, double kappa)
 Returns the low-pressure mixture viscosity [Pa·s] using the Chung method.
 
double highPressureViscosity (double T_star, double MW, double rho, double Vc, double Tc, double acentric_factor, double mu_r, double kappa)
 Returns the high-pressure mixture viscosity [micropoise] using the Chung method.
 
double highPressureThermalConductivity (double T, double T_star, double MW, double rho, double Cv, double Vc, double Tc, double sigma, double acentric_factor, double mu_r, double kappa)
 Computes the high-pressure thermal conductivity [W/m/K] using the Chung method.
 
- Protected Member Functions inherited from HighPressureGasTransportBase
 HighPressureGasTransportBase ()=default
 default constructor
 
void getTransportData () override
 Obtain required parameters from the 'critical-parameters' species input section, and checks the critical-properties.yaml file if an acentric factor is not specified.
 
void initializeCriticalProperties ()
 Computes and stores the estimate of the critical properties for each species.
 
double Tcrit_i (size_t i)
 Returns the stored value of the critical temperature for species 'i'.
 
double Pcrit_i (size_t i)
 Returns the stored value of the critical pressure for species 'i'.
 
double Vcrit_i (size_t i)
 Returns the stored value of the critical volume for species 'i'.
 
double Zcrit_i (size_t i)
 Returns the stored value of the critical compressibility for species 'i'.
 
- Protected Member Functions inherited from MixTransport
void updateCond_T ()
 Update the temperature dependent parts of the species thermal conductivities.
 
- Protected Member Functions inherited from GasTransport
virtual void update_T ()
 
virtual void update_C ()=0
 
virtual void updateViscosity_T ()
 Update the temperature-dependent viscosity terms.
 
virtual void updateSpeciesViscosities ()
 Update the pure-species viscosities.
 
virtual void updateDiff_T ()
 Update the binary diffusion coefficients.
 
virtual void setupCollisionParameters ()
 Setup parameters for a new kinetic-theory-based transport manager for low-density gases.
 
void setupCollisionIntegral ()
 Setup range for polynomial fits to collision integrals of Monchick & Mason [31].
 
void makePolarCorrections (size_t i, size_t j, double &f_eps, double &f_sigma)
 Corrections for polar-nonpolar binary diffusion coefficients.
 
void fitCollisionIntegrals (MMCollisionInt &integrals)
 Generate polynomial fits to collision integrals.
 
virtual void fitProperties (MMCollisionInt &integrals)
 Generate polynomial fits to the viscosity \( \eta \) and conductivity \( \lambda \).
 
virtual void fitDiffCoeffs (MMCollisionInt &integrals)
 Generate polynomial fits to the binary diffusion coefficients.
 
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.
 

Private Attributes

Pure fluid properties

These are the pure fluid properties of each species that are needed to compute the mixture properties for the Chung viscosity and thermal conductivity models.

vector< double > m_sigma_i
 Effective molecular diameter [Angstroms].
 
vector< double > m_epsilon_over_k_i
 Characteristic temperature [K].
 
vector< double > m_MW_i
 Molecular weight [kg/kmol].
 
vector< double > m_acentric_factor_i
 Acentric factor [unitless].
 
vector< double > m_kappa_i
 Association factor [unitless].
 
Chung mixture parameters

These are the parameters that are needed to calculate the viscosity and thermal conductivity using the Chung method.

double m_Vc_mix = 0
 Mixture critical volume [m³/kmol].
 
double m_Tc_mix = 0
 Mixture critical temperature [K].
 
double m_sigma_mix = 0
 Effective mixture molecular diameter [Angstroms].
 
double m_epsilon_over_k_mix = 0
 Mixture characteristic temperature [K].
 
double m_MW_mix = 0
 Effective mixture molecular weight [kg/kmol].
 
double m_mu_mix = 0
 Mixture dipole moment [Debye].
 
double m_mu_r_mix = 0
 Mixture reduced dipole moment [unitless].
 
double m_acentric_factor_mix = 0
 Mixture acentric factor [unitless].
 
double m_kappa_mix = 0
 Mixture association factor [unitless].
 

Friends

class TransportFactory
 

Additional Inherited Members

- Protected Attributes inherited from HighPressureGasTransportBase
vector< double > m_Tcrit
 Critical temperature [K] of each species.
 
vector< double > m_Pcrit
 Critical pressure [Pa] of each species.
 
vector< double > m_Vcrit
 Critical volume [m³/kmol] of each species.
 
vector< double > m_Zcrit
 Critical compressibility [unitless] of each species.
 
DenseMatrix m_P_corr_ij
 Matrix of Takahashi binary diffusion coefficient corrections.
 
- Protected Attributes inherited from MixTransport
vector< double > m_cond
 vector of species thermal conductivities [W/m/K]
 
double m_lambda = 0.0
 Internal storage for the calculated mixture thermal conductivity [W/m/K].
 
bool m_spcond_ok = false
 Update boolean for the species thermal conductivities.
 
bool m_condmix_ok = false
 Update boolean for the mixture rule for the mixture thermal conductivity.
 
- Protected Attributes inherited from GasTransport
vector< double > m_molefracs
 Vector of species mole fractions.
 
double m_viscmix = 0.0
 Internal storage for the viscosity of the mixture [Pa·s].
 
bool m_visc_ok = false
 Update boolean for mixture rule for the mixture viscosity.
 
bool m_viscwt_ok = false
 Update boolean for the weighting factors for the mixture viscosity.
 
bool m_spvisc_ok = false
 Update boolean for the species viscosities.
 
bool m_bindiff_ok = false
 Update boolean for the binary diffusivities at unit pressure.
 
int m_mode = 0
 Type of the polynomial fits to temperature.
 
DenseMatrix m_phi
 Viscosity weighting function. size = m_nsp * m_nsp.
 
vector< double > m_spwork
 work space length = m_nsp
 
vector< double > m_visc
 vector of species viscosities [Pa·s].
 
vector< vector< double > > m_visccoeffs
 Polynomial fits to the viscosity of each species.
 
vector< double > m_mw
 Local copy of the species molecular weights.
 
DenseMatrix m_wratjk
 Holds square roots of molecular weight ratios.
 
DenseMatrix m_wratkj1
 Holds square roots of molecular weight ratios.
 
vector< double > m_sqvisc
 vector of square root of species viscosities.
 
vector< double > m_polytempvec
 Powers of the ln temperature, up to fourth order.
 
double m_temp = -1.0
 Current value of the temperature [K] at which the properties in this object are calculated.
 
double m_kbt = 0.0
 Current value of Boltzmann constant times the temperature [J].
 
double m_sqrt_t = 0.0
 current value of temperature to 1/2 power
 
double m_logt = 0.0
 Current value of the log of the temperature.
 
double m_t14 = 0.0
 Current value of temperature to 1/4 power.
 
vector< vector< double > > m_diffcoeffs
 Polynomial fits to the binary diffusivity of each species.
 
DenseMatrix m_bdiff
 Matrix of binary diffusion coefficients at the reference pressure and the current temperature Size is m_nsp x m_nsp.
 
vector< vector< double > > m_condcoeffs
 temperature fits of the heat conduction
 
vector< vector< int > > m_poly
 Indices for the (i,j) interaction in collision integral fits.
 
vector< vector< double > > m_omega22_poly
 Fit for omega22 collision integral.
 
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.
 
vector< vector< double > > m_astar_poly
 Fit for astar collision integral.
 
vector< vector< double > > m_bstar_poly
 Fit for bstar collision integral.
 
vector< vector< double > > m_cstar_poly
 Fit for cstar collision integral.
 
vector< double > m_zrot
 Rotational relaxation number for each species.
 
vector< double > m_crot
 Dimensionless rotational heat capacity of each species.
 
vector< bool > m_polar
 Vector of booleans indicating whether a species is a polar molecule.
 
vector< double > m_alpha
 Polarizability [m³] of each species in the phase.
 
vector< double > m_eps
 Lennard-Jones well-depth [J] of the species in the current phase.
 
vector< double > m_sigma
 Lennard-Jones diameter [m] of the species in the current phase.
 
DenseMatrix m_reducedMass
 This is the reduced mass [kg] of the interaction between species i and j.
 
DenseMatrix m_diam
 hard-sphere diameter [m] for (i,j) collision
 
DenseMatrix m_epsilon
 The effective well depth [J] for (i,j) collisions.
 
DenseMatrix m_dipole
 The effective dipole moment [Coulomb·m] for (i,j) collisions.
 
DenseMatrix m_delta
 Reduced dipole moment of the interaction between two species.
 
vector< double > m_w_ac
 Pitzer acentric factor [dimensionless].
 
vector< double > m_disp
 Dispersion coefficient normalized by the square of the elementary charge [m⁵].
 
vector< double > m_quad_polar
 Quadrupole polarizability.
 
- Protected Attributes inherited from Transport
ThermoPhasem_thermo
 pointer to the object representing the phase
 
size_t m_nsp = 0
 Number of species in the phase.
 
AnyMap m_fittingErrors
 Maximum errors associated with fitting pure species transport properties.
 

Constructor & Destructor Documentation

◆ ChungHighPressureGasTransport()

ChungHighPressureGasTransport ( )
protecteddefault

default constructor

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 833 of file HighPressureGasTransport.h.

◆ init()

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

Reimplemented from Transport.

Definition at line 776 of file HighPressureGasTransport.cpp.

◆ viscosity()

double viscosity ( )
overridevirtual

Returns the high-pressure mixture viscosity [Pa·s] using the Chung method.

Based on the high-pressure gas mixture viscosity model of Chung described in chapter 9-7 of Poling. This method uses the pure species high-pressure viscosity relation of Chung with the difference being that mixture values of the model are computed using a set of mixing rules given by Chung (see computeMixtureParameters() ). The mixing rules are defined in section 9-5 of [38].

Because this method is using the high-pressure viscosity model with mixture parameters, see highPressureViscosity() for details on the model.

Reimplemented from Transport.

Definition at line 900 of file HighPressureGasTransport.cpp.

◆ thermalConductivity()

double thermalConductivity ( )
overridevirtual

Calculates the high-pressure mixture thermal conductivity using the Chung method.

This method obtains appropriate mixture values of the parameters needed for the Chung model and then calls the highPressureThermalConductivity() method to obtain the mixture thermal conductivity.

The mixture values of the pseudo-critical temperature and other model parameters are calculated using the Chung mixing rules defined on page 9.25 (see computeMixtureParameters() ).

The mixture value of the specific heat is computed using equation 10-6.6, which is the mole fraction weighted sum of the pure species specific heats. This value is not directly computed by the computeMixtureParameters() method.

\[ C_{v,m} = \sum_i X_i C_{v,i} \quad \text{(Equation 10-6.6)} \]

Where \( C_{v,i} \) is the specific heat of species i, and \( X_i \) is the mole fraction of species i.

Reimplemented from Transport.

Definition at line 812 of file HighPressureGasTransport.cpp.

◆ initializePureFluidProperties()

void initializePureFluidProperties ( )
protected

Computes and stores pure-fluid specific properties each species.

Definition at line 784 of file HighPressureGasTransport.cpp.

◆ computeMixtureParameters()

void computeMixtureParameters ( )
protected

Computes the composition-dependent values of the parameters needed for the Chung viscosity model.

The equations for the mixing rules defined on page 9.25 for the Chung method's composition dependent parameters. The primary mixing rules are defined below.

\[ \sigma_m^3 = \sum_{i} \sum_{j} X_i X_j \sigma_{ij}^3 \quad \text{(Equation 9-5.25)} \]

Where \( \sigma_{ij} \) is the molecular diameter

\[ T_m^* = \frac{T}{\left( \frac{\epsilon}{k} \right )_m} \quad \text{(Equation 9-5.26)} \]

Where \( k \) is the Boltzmann constant and \( \epsilon \) is the minimum of the pair-potential energy. In these equations, we do not need to worry about what the values of \( \epsilon \) and \( k \) are.

\[ \left( \frac{\epsilon}{k} \right)_m = \frac{\sum_{i} \sum_{j} X_i X_j \left( \frac{\epsilon_{ij}}{k} \right) \sigma_{ij}^3}{\sigma_m^3} \quad \text{(Equation 9-5.27)} \]

\[ MW_m = \left[ \frac{\sum_{i} \sum_{j} X_i X_j \left( \frac{\epsilon_{ij}}{k} \right) \sigma_{ij}^2 MW_{ij}^{\frac{1}{2}}}{\left( \frac{\epsilon}{k} \right)_m \sigma_m^2} \right]^2 \quad \text{(Equation 9-5.28)} \]

Where MW is the molecular weight.

\[ \omega_m = \frac{\sum_{i} \sum_{j} X_i X_j \omega_{ij} \sigma{ij}^3}{\sigma_m^3} \quad \text{(Equation 9-5.29)} \]

Where \( \omega \) is the acentric factor.

\[ \mu_m^4 = \sigma_m^3 \sum_{i} \sum_{j} \left( \frac{X_i X_j \mu_i^2 \mu_j^2} {\sigma_{ij}^3} \right) \quad \text{(Equation 9-5.30)} \]

Where \( \mu \) is the dipole moment.

\[ \kappa_m = \sum_{i} \sum_{j} X_i X_j \kappa_{ij} \quad \text{(Equation 9-5.31)} \]

Where \( \kappa \) is the association factor, which is used for highly polar molecules. In this work, the value is assumed to be zero for all species.

The combining rules for species-species values (subscripted with ij in the above equations) are defined below.

\[ \sigma_{i} = 0.809 V_{c,i}^{1/3} \quad \text{(Equation 9-5.32)} \]

Where \( V_{c,i} \) is the critical volume of species i.

\[ \sigma_{ij} = \xi_{ij} \left( \sigma_{i} \sigma_{j} \right)^{1/2} \quad \text{(Equation 9-5.33)} \]

Where \( \xi \) is a binary interaction parameter.

\[ \left( \frac{\epsilon_i}{k} \right) = \frac{T_{c,i}}{1.2593} \quad \text{(Equation 9-5.34)} \]

\[ \frac{\epsilon_{ij}}{k} = \zeta_{ij} \left( \right)^{\frac{1}{2}} \quad \text{(Equation 9-5.35)} \]

Where \( \zeta \) is a binary interaction parameter.

\[ \omega_{ij} = \frac{\omega_i + \omega_j}{2} \quad \text{(Equation 9-5.37)} \]

Where \( \omega \) is the acentric factor.

\[ \kappa_{ij} = \left( \kappa_i \kappa_j \right)^{1/2} \quad \text{(Equation 9-5.39)} \]

Where \( \kappa \) is the association factor.

\[ MW_{ij} = \frac{2 MW_i MW_j}{MW_i + MW_j} \quad \text{(Equation 9-5.40)} \]

\( \xi \) and \( \zeta \) are the binary interaction parameters, and are assumed to be unity in this implementation, in keeping with the Chung method.

The Chung viscosity correction factor is defined as:

\[ F_{c,m} = 1 - 0.275 \omega_m + 0.059035 \mu_{r,m}^4 + \kappa_m \quad \text{(Equation 9-5.41)} \]

Where \( \omega_m \) is the mixture acentric factor, \( \mu_{r,m} \) is the mixture reduced dipole moment, and \( \kappa_m \) is the mixture association factor.

The mixture reduced dipole moment is computed using:

\[ \mu_{r,m} = \frac{131.3 \mu_m}{( V_{c,m} T_{c,m})^{\frac{1}{2}}} \quad \text{(Equation 9-5.42)} \]

Where \( V_{c,m} \) and \( T_{c,m} \) are computed using the following equations.

\[ V_{c,m} = \left( \frac{\sigma_m}{0.809} \right)^3 \quad \text{(Equation 9-5.43)} \]

\[ T_{c,m} = 1.2593 \left( \frac{\epsilon}{k} \right)_m \quad \text{(Equation 9-5.44)} \]

In the equations, \( T_c \) must be in units of K, \( V_c \) must be in units of cm^3/mol, and \( \mu \) must be in units of Debye.

Definition at line 925 of file HighPressureGasTransport.cpp.

◆ lowPressureViscosity()

double lowPressureViscosity ( double  T,
double  T_star,
double  MW,
double  acentric_factor,
double  mu_r,
double  sigma,
double  kappa 
)
protected

Returns the low-pressure mixture viscosity [Pa·s] using the Chung method.

\[ \eta = 26.69 F_c \frac{(MW*T)^(1/2)}{\sigma^2 \Omega} \quad \text{(Equation 9-4.10)} \]

T must be in units of K, MW must be units of kg/kmol, and \( \sigma \) must be in units of Angstroms. The viscosity is computed in micropoise, but the return value is in standard SI units (Pa*s).

This function is structured such that it can be used for pure species or mixtures, with the only difference being the values that are passed to the function (pure values versus mixture values).

Parameters
TTemperature [K]
T_starReduced temperature [unitless]
MWMolecular weight [kg/kmol]
acentric_factorAcentric factor [unitless]
mu_rDipole moment [Debye]
sigmaLennard-Jones collision diameter [Angstroms]
kappaPolar correction factor [unitless]

Definition at line 1039 of file HighPressureGasTransport.cpp.

◆ highPressureViscosity()

double highPressureViscosity ( double  T_star,
double  MW,
double  rho,
double  Vc,
double  Tc,
double  acentric_factor,
double  mu_r,
double  kappa 
)
protected

Returns the high-pressure mixture viscosity [micropoise] using the Chung method.

\[ \eta = \eta^* \frac{36.344 (M*T_c)^(1/2)}{V^{\frac{2}{3}}} \quad \text{(Equation 9-6.18)} \]

where:

\[ \begin{align*} \eta &= \text{viscosity, } \mu P \\ M &= \text{molecular weight, kg/kmol} \\ T_c &= \text{critical temperature, K} \\ V_c &= \text{critical molar volume, cm}^3 / \text{mol} \\ \end{align*} \]

and,

\[ \eta^* = \frac{(T^*)^{\frac{1}{2}}}{\Omega_v} {F_c[(G_2)^{-1} + E_6 y]} + \eta^{**} \quad \text{(Equation 9-6.19)} \]

The values of \( T^* \) and \( F_c \) are defined as follows.

\[ T^* = 1.2593 T_r \quad \text{(Equation 9-4.9)} \]

\[ F_c = 1 - 0.275 \omega + 0.059035 \mu_r^4 + \kappa \quad \text{(Equation 9-4.11)} \]

The value of \( \Omega_v \) is the viscosity collision integral evaluated at the non-dimensional reduced temperature \( T^* \). For details on the collision integral see neufeldCollisionIntegral() .

This function is structured such that it can be used for pure species or mixtures, with the only difference being the values that are passed to the function (pure values versus mixture values).

Parameters
T_starReduced temperature [unitless]
MWMolecular weight [kg/kmol]
rhoDensity [mol/cm³]
VcCritical volume [cm³/mol]
TcCritical temperature [K]
acentric_factorAcentric factor [unitless]
mu_rDipole moment [Debye]
kappaPolar correction factor [unitless]

Definition at line 1056 of file HighPressureGasTransport.cpp.

◆ highPressureThermalConductivity()

double highPressureThermalConductivity ( double  T,
double  T_star,
double  MW,
double  rho,
double  Cv,
double  Vc,
double  Tc,
double  sigma,
double  acentric_factor,
double  mu_r,
double  kappa 
)
protected

Computes the high-pressure thermal conductivity [W/m/K] using the Chung method.

The Chung method for computing high-pressure thermal conductivity is described on page 10.23 of [38] .

\[ \lambda = \frac{31.2 \eta^0 \Psi}{M'} \left( G_2^{-1} + B_6 y \right) + q B_7 y^2 T_r^{1/2} G_2 \quad \text{(Equation 10-5.5)} \]

where:

\[ \begin{align*} \lambda &= \text{thermal conductivity, W/(m·K)} \\ \eta^0 &= \text{low-pressure gas viscosity, N·s/m}^2 \\ M' &= \text{molecular weight, kg/mol} \\ \Psi &= f(C_v, \omega, T_r) \text{ [as defined under Eq. (10-3.14)]} \\ q &= 3.586 \times 10^{-3} \left( \frac{T_c}{M'} \right)^{1/2} V_c^{2/3} \\ T &= \text{temperature, K} \\ T_c &= \text{critical temperature, K} \\ T_r &= \text{reduced temperature, } \frac{T}{T_c} \\ V_c &= \text{critical molar volume, cm}^3/\text{mol} \\ \end{align*} \]

where,

\[ y = \frac{V_c}{6V} \quad \text{(Equation 10-5.6)} \]

V is the molar volume of the fluid in cm^3/mol.

\[ G_1 = \frac{1 - 0.5y}{(1-y)^3} \quad \text{(Equation 10-5.7)} \]

\[ G_2 = \frac{(B_1 / y)[1 - \text{exp}(-B_4 y)] + B_2 G_1 \text{exp}(B_5 y) + B_3 G_1}{B_1 B_4 + B_2 + B_3} \quad \text{(Equation 10-5.8)} \]

The coefficients \( B_1 \), through \( B_7 \) are functions of the acentric factor, reduced dipole moment and the association factor.

\[ B_i = a_i + b_i \omega + c_i \mu_r^4 + d_i \kappa \quad \text{(Equation 10-5.9)} \]

The constants in the above equation are from table 10-3 on page 10.23 of [38].

The definition of the \( \Psi \) function is given by:

\[ \Psi = 1 + \alpha {\frac{0.215 + 0.28288\alpha - 1.061\beta + 0.26665Z}{0.6366 + \beta Z + 1.061 \alpha \beta}} \]

with,

\[ \alpha = \frac{C_v}{R} - \frac{3}{2} \]

\[ \beta = 0.7862 - 0.7109 \omega + 1.3168 \omega^2 \]

\[ Z = 2.0 + 10.5 T_r^2 \]

These functions are from page 10.12 of [38] .

This method utilizes the low-pressure Chung viscosity as that is a required parameter in the model, and thus calls the low pressure viscosity implementation. This is why it requires parameters typically associated with the viscosity calculation.

This function is structured such that it can be used for pure species or mixtures, with the only difference being the values that are passed to the function (pure values versus mixture values).

For mixtures, the mixture values of the input variables are computed using the mixing rules of Chung, see computeMixtureParameters() .

Parameters
TTemperature [K]
T_starReduced temperature [unitless]
MWMolecular weight [kg/kmol]
rhoDensity [mol/cm³]
CvSpecific heat [J/kg/K]
VcCritical volume [cm³/mol]
TcCritical temperature [K]
sigmaLennard-Jones collision diameter [Angstroms]
acentric_factorAcentric factor [unitless]
mu_rDipole moment [Debye]
kappaPolar correction factor [unitless]
Returns
High pressure thermal conductivity [W/m/K]

Definition at line 840 of file HighPressureGasTransport.cpp.

Friends And Related Symbol Documentation

◆ TransportFactory

friend class TransportFactory
friend

Definition at line 880 of file HighPressureGasTransport.h.

Member Data Documentation

◆ m_sigma_i

vector<double> m_sigma_i
private

Effective molecular diameter [Angstroms].

Definition at line 1271 of file HighPressureGasTransport.h.

◆ m_epsilon_over_k_i

vector<double> m_epsilon_over_k_i
private

Characteristic temperature [K].

Definition at line 1272 of file HighPressureGasTransport.h.

◆ m_MW_i

vector<double> m_MW_i
private

Molecular weight [kg/kmol].

Definition at line 1273 of file HighPressureGasTransport.h.

◆ m_acentric_factor_i

vector<double> m_acentric_factor_i
private

Acentric factor [unitless].

Definition at line 1274 of file HighPressureGasTransport.h.

◆ m_kappa_i

vector<double> m_kappa_i
private

Association factor [unitless].

Definition at line 1275 of file HighPressureGasTransport.h.

◆ m_Vc_mix

double m_Vc_mix = 0
private

Mixture critical volume [m³/kmol].

Definition at line 1285 of file HighPressureGasTransport.h.

◆ m_Tc_mix

double m_Tc_mix = 0
private

Mixture critical temperature [K].

Definition at line 1286 of file HighPressureGasTransport.h.

◆ m_sigma_mix

double m_sigma_mix = 0
private

Effective mixture molecular diameter [Angstroms].

Definition at line 1288 of file HighPressureGasTransport.h.

◆ m_epsilon_over_k_mix

double m_epsilon_over_k_mix = 0
private

Mixture characteristic temperature [K].

Definition at line 1289 of file HighPressureGasTransport.h.

◆ m_MW_mix

double m_MW_mix = 0
private

Effective mixture molecular weight [kg/kmol].

Definition at line 1290 of file HighPressureGasTransport.h.

◆ m_mu_mix

double m_mu_mix = 0
private

Mixture dipole moment [Debye].

Definition at line 1293 of file HighPressureGasTransport.h.

◆ m_mu_r_mix

double m_mu_r_mix = 0
private

Mixture reduced dipole moment [unitless].

Definition at line 1294 of file HighPressureGasTransport.h.

◆ m_acentric_factor_mix

double m_acentric_factor_mix = 0
private

Mixture acentric factor [unitless].

Definition at line 1295 of file HighPressureGasTransport.h.

◆ m_kappa_mix

double m_kappa_mix = 0
private

Mixture association factor [unitless].

Definition at line 1296 of file HighPressureGasTransport.h.


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