Base class for high pressure gas transport models. More...
#include <HighPressureGasTransport.h>
Base class for high pressure gas transport models.
This class handles common functionality used by high pressure transport models, including critical property storage and binary diffusion coefficient correction.
Definition at line 26 of file HighPressureGasTransport.h.
Public Member Functions | |
void | getBinaryDiffCoeffs (const size_t ld, double *const d) override |
Computes the matrix of binary diffusion coefficients using the Takahashi correction factor. | |
void | getMixDiffCoeffs (double *const d) override |
Returns the mixture-averaged diffusion coefficients [m^2/s]. | |
void | getMixDiffCoeffsMole (double *const d) override |
Returns the mixture-averaged diffusion coefficients [m^2/s]. | |
void | getMixDiffCoeffsMass (double *const d) override |
Returns the mixture-averaged diffusion coefficients [m^2/s]. | |
virtual void | updateCorrectionFactors () |
Updates the matrix of species-pair Takahashi correction factors for use in computing the binary diffusion coefficients, see takahashiCorrectionFactor() | |
![]() | |
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. | |
double | thermalConductivity () override |
Returns the mixture thermal conductivity (W/m /K) | |
void | getMobilities (double *const mobil) override |
Get the Electrical mobilities (m^2/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 wrt 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. | |
![]() | |
double | viscosity () override |
Viscosity of the mixture (kg /m /s) | |
void | getSpeciesViscosities (double *const visc) override |
Get the pure-species viscosities. | |
void | getBinaryDiffCoeffs (const size_t ld, double *const d) override |
Computes the matrix of binary diffusion coefficients. | |
void | getMixDiffCoeffs (double *const d) override |
Returns the Mixture-averaged diffusion coefficients [m^2/s]. | |
void | getMixDiffCoeffsMole (double *const d) override |
Returns the mixture-averaged diffusion coefficients [m^2/s]. | |
void | getMixDiffCoeffsMass (double *const d) override |
Returns the mixture-averaged diffusion coefficients [m^2/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. | |
![]() | |
Transport ()=default | |
Constructor. | |
Transport (const Transport &)=delete | |
Transport & | operator= (const Transport &)=delete |
virtual string | transportModel () const |
Identifies the model represented by this Transport object. | |
ThermoPhase & | thermo () |
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 nSpecies(). | |
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. | |
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. | |
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. | |
virtual void | getThermalDiffCoeffs (double *const dt) |
Return a vector of Thermal diffusion coefficients [kg/m/sec]. | |
virtual void | getBinaryDiffCoeffs (const size_t ld, double *const d) |
Returns the matrix of binary diffusion coefficients [m^2/s]. | |
virtual void | getMultiDiffCoeffs (const size_t ld, double *const d) |
Return the Multicomponent diffusion coefficients. Units: [m^2/s]. | |
virtual void | getMixDiffCoeffs (double *const d) |
Returns a vector of mixture averaged diffusion coefficients. | |
virtual void | getMixDiffCoeffsMole (double *const d) |
Returns a vector of mixture averaged diffusion coefficients. | |
virtual void | getMixDiffCoeffsMass (double *const d) |
Returns a vector of mixture averaged diffusion coefficients. | |
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 in Pa-s. | |
virtual double | electricalConductivity () |
The electrical conductivity (Siemens/m). | |
Protected Member Functions | |
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'. | |
![]() | |
void | updateCond_T () |
Update the temperature dependent parts of the species thermal conductivities. | |
![]() | |
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 [30]. | |
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. | |
Protected Attributes | |
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^3/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. Size is nsp x nsp. | |
![]() | |
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. | |
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. | |
![]() | |
vector< double > | m_molefracs |
Vector of species mole fractions. | |
double | m_viscmix = 0.0 |
Internal storage for the viscosity of the mixture (kg /m /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 |
m_phi is a Viscosity Weighting Function. size = m_nsp * n_nsp | |
vector< double > | m_spwork |
work space length = m_kk | |
vector< double > | m_visc |
vector of species viscosities (kg /m /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 sqrt(kg /m /s). | |
vector< double > | m_polytempvec |
Powers of the ln temperature, up to fourth order. | |
double | m_temp = -1.0 |
Current value of the temperature at which the properties in this object are calculated (Kelvin). | |
double | m_kbt = 0.0 |
Current value of Boltzmann constant times the temperature (Joules) | |
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 nsp x 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 of each species in the phase. | |
vector< double > | m_eps |
Lennard-Jones well-depth of the species in the current phase. | |
vector< double > | m_sigma |
Lennard-Jones diameter of the species in the current phase. | |
DenseMatrix | m_reducedMass |
This is the reduced mass of the interaction between species i and j. | |
DenseMatrix | m_diam |
hard-sphere diameter for (i,j) collision | |
DenseMatrix | m_epsilon |
The effective well depth for (i,j) collisions. | |
DenseMatrix | m_dipole |
The effective dipole moment for (i,j) collisions. | |
DenseMatrix | m_delta |
Reduced dipole moment of the interaction between two species. | |
vector< double > | m_w_ac |
Pitzer acentric factor. | |
vector< double > | m_disp |
Dispersion coefficient. | |
vector< double > | m_quad_polar |
Quadrupole polarizability. | |
![]() | |
ThermoPhase * | m_thermo |
pointer to the object representing the phase | |
size_t | m_nsp = 0 |
Number of species. | |
AnyMap | m_fittingErrors |
Maximum errors associated with fitting pure species transport properties. | |
Friends | |
class | TransportFactory |
|
protecteddefault |
default constructor
|
overridevirtual |
Computes the matrix of binary diffusion coefficients using the Takahashi correction factor.
Units are m^2/s.
The matrix is dimension m_nsp x m_nsp, where m_nsp is the number of species. The matrix is stored in row-major order, so that d[ld*j + i] contains the binary diffusion coefficient of species i with respect to species j.
d[ld*j + i] = (DP)_R * m_bdiff(i,j) / p
ld | Inner stride for writing the two dimension diffusion coefficients into a one dimensional vector |
d | Diffusion coefficient matrix (must be at least m_nsp * m_nsp in length. |
Reimplemented from GasTransport.
Definition at line 219 of file HighPressureGasTransport.cpp.
|
overridevirtual |
Returns the mixture-averaged diffusion coefficients [m^2/s].
This method is the same as GasTransport::getMixDiffCoeffs() with the exception that the binary diffusion coefficients are multiplied by the Takahashi correction factor, which is described in takahashiCorrectionFactor() .
[out] | d | Vector of mixture diffusion coefficients, \( D_{km}' \) , for each species (m^2/s). length m_nsp |
Reimplemented from GasTransport.
Definition at line 244 of file HighPressureGasTransport.cpp.
|
overridevirtual |
Returns the mixture-averaged diffusion coefficients [m^2/s].
This method is the same as GasTransport::getMixDiffCoeffsMole() with the exception that the binary diffusion coefficients are multiplied by the Takahashi correction factor, which is described in takahashiCorrectionFactor() .
[out] | d | vector of mixture-averaged diffusion coefficients for each species, length m_nsp. |
Reimplemented from GasTransport.
Definition at line 276 of file HighPressureGasTransport.cpp.
|
overridevirtual |
Returns the mixture-averaged diffusion coefficients [m^2/s].
This method is the same as GasTransport::getMixDiffCoeffsMass() with the exception that the binary diffusion coefficients are multiplied by the Takahashi correction factor, which is described in takahashiCorrectionFactor() .
[out] | d | vector of mixture-averaged diffusion coefficients for each species, length m_nsp. |
Reimplemented from GasTransport.
Definition at line 307 of file HighPressureGasTransport.cpp.
|
virtual |
Updates the matrix of species-pair Takahashi correction factors for use in computing the binary diffusion coefficients, see takahashiCorrectionFactor()
Definition at line 188 of file HighPressureGasTransport.cpp.
|
overrideprotectedvirtual |
Obtain required parameters from the 'critical-parameters' species input section, and checks the critical-properties.yaml file if an acentric factor is not specified.
The way that GasTransport parses the input file is that if an acentric factor is not specified, it is quietly set to 0.0. A user may have the proper acentric factor in the critical-properties.yaml file, so that file is checked if a zero value is present.
Reimplemented from GasTransport.
Definition at line 98 of file HighPressureGasTransport.cpp.
|
protected |
Computes and stores the estimate of the critical properties for each species.
This method sets the species composition vector to unity for species i and zero for all other species, and then queries the thermo object for the critical properties and stores them. It then resets the composition vector to the original state. This method only needs to be called once, as the critical properties for the pure species do not change.
All species must have critical properties defined in the input file, either via critical properties or by specific values of the equation of state that are not zero.
Definition at line 130 of file HighPressureGasTransport.cpp.
|
protected |
Returns the stored value of the critical temperature for species 'i'.
Definition at line 168 of file HighPressureGasTransport.cpp.
|
protected |
Returns the stored value of the critical pressure for species 'i'.
Definition at line 173 of file HighPressureGasTransport.cpp.
|
protected |
Returns the stored value of the critical volume for species 'i'.
Definition at line 178 of file HighPressureGasTransport.cpp.
|
protected |
Returns the stored value of the critical compressibility for species 'i'.
Definition at line 183 of file HighPressureGasTransport.cpp.
|
friend |
Definition at line 141 of file HighPressureGasTransport.h.
|
protected |
Critical temperature [K] of each species.
Definition at line 133 of file HighPressureGasTransport.h.
|
protected |
Critical pressure [Pa] of each species.
Definition at line 134 of file HighPressureGasTransport.h.
|
protected |
Critical volume [m^3/kmol] of each species.
Definition at line 135 of file HighPressureGasTransport.h.
|
protected |
Critical compressibility [unitless] of each species.
Definition at line 136 of file HighPressureGasTransport.h.
|
protected |
Matrix of Takahashi binary diffusion coefficient corrections. Size is nsp x nsp.
Definition at line 139 of file HighPressureGasTransport.h.