Methods in the py_KineticGas class
KineticGas - Latest (beta)
The py_KineticGas
class, found in pykingas/py_KineticGas.py
, is the core of the KineticGas Python interface. All Revised Enskog Theory models inherit from py_KineticGas
. This is the class that contains the interface to all practical calculations that can be done from the python-side of KineticGas. Derived classes only implement specific functions for parameter handling etc.
Table of contents
- The constructor
- TV-property interfaces
- Tp-property interfaces
- Frame of Reference transformations
- Interfaces to C++ methods
- Utility methods
- Deprecated methods
- Internal methods
The constructor
The constructor
Table of contents
__init__(self, comps, mole_weights=None, N=3, is_idealgas=False)
Args:
comps (str):
Comma-separated list of components, following ThermoPack-convention
mole_weights (1d array) :
Mole weights [g/mol]. Will be used instead of database values if provided
is_idealgas (bool) :
If true, radial distribution function is unity, if false use radial distribution function of modelIn addition, several density-dependent factors are set to zero, to ensure consistency with the ideal gas law / Gibbs-Duhem for an ideal gas.
TV-property interfaces
Computing properties as a function of temperature and volume.
Table of contents
bulk_viscosity(self, T, Vm, x, N=None)
Not implemented
Raises: NotImplementedError
conductivity_matrix(self, T, Vm, x, N=2, formulation='T-psi', frame_of_reference='CoM', use_thermal_conductivity=None)
Compute the conductivity matrix $L$, for use in NET calculations. The Flux/Force formulation used in the NET
model is selected using the formulation
kwarg. Currently implemented formulations are:
———————————————————————————————-
'T-psi'
:
\[J_q = L_{qq} \nabla (1 / T) - \sum_{i=1}^{N_c-1}(1 / T) L_{qi} \nabla_T \Psi_{i; p, x}\] \[J_i = L_{iq} \nabla (1 / T) - \sum_{j=1}^{N_c-1}(1 / T) L_{ij} \nabla_T \Psi_{j; p, x}\]
Where
\[\Psi_i = \mu_i - \mu_{N_c}\]and $N_c$ denotes the number of components. The last component is used as the dependent component. The fluxes in this formulation are on a mass basis, and the formulation is implemented in the centre of mass frame of reference. The formulation is only implemented for ideal gases. ———————————————————————————————-
Args:
T (float) :
Temperature [K]
Vm (float) :
Molar volume [m3 / mol]
x (1darray) :
Molar composition [-]
N (int, optional) :
Enskog approximation order (must be >= 2 for thermal effects)
formulation (str, optional) :
The NET formulation for which to compute the conductivity matrix.
frame_of_reference (str, optional) :
The frame of reference (‘CoM’, ‘CoN’, ‘CoV’, or ‘solvent’). Default: ‘CoM’.
use_thermal_conductivity (callable, optional) :
External thermal conductivity model. Assumed to have the signatureuse_thermal_conductivity(T, Vm, x), returning the thermal conductivity in [W / m K]. Defaults to None. If no model is supplied, KineticGas is used to compute thermal conductivity.
Returns:
ndarray :
The conductivity matrix, contents will vary depending on the formulation
kwarg.
interdiffusion(self, T, Vm, x, N=None, use_independent=True, dependent_idx=-1, frame_of_reference='CoN', use_binary=True, solvent_idx=-1)
Compute the interdiffusion coefficients [m^2 / s]. Default definition is
\[J_i^{(n, n)} = - \sum_{j \neq l} D_{ij} \nabla n_j, \nabla T = \nabla p = F_k = 0 \forall k\]where the flux, $J_i^{(n, n)}$ is on a molar basis, in the molar frame of reference, and $j \neq l$ is an
independent set of forces with $l=$ dependent_idx
.
For fluxes in other frames of reference, use the frame_of_reference
kwarg.
For the diffusion coefficients describing the fluxes’ response to all forces (not independent) the definition
is:
use the use_independent
and dependent_idx
kwargs to switch between these definitions.
See: Eq. (17-20) in RET for Mie fluids (https://doi.org/10.1063/5.0149865)
Args:
T (float) :
Temperature (K)
Vm (float) :
Molar volume (m3 / mol)
x (array_like) :
composition (mole fractions)
N (int, optional) :
Enskog approximation order. Default set on model initialisation.
use_binary (bool, optional) :
If the mixture is binary, and an independent set of fluxes and forces is considered, i.e.use_independent=True
, the diffusion coefficients will be exactly equal with opposite sign. Setting use_binary=True
will in that case only return the coefficient describing the independent flux-force relation.
use_independent (bool, optional) :
Return diffusion coefficients for independent set of forces.
dependent_idx (int, optional) :
Index of the dependent molar density gradient (only if use_independent=True
, the default behaviour).Defaults to last component, except when frame_of_reference='solvent'
, in which case default is equal to solvent_idx
.
frame_of_reference (str, optional) :
Which frame of reference the diffusion coefficients apply to. Defaultis 'CoN'
. Can be 'CoN'
(molar FoR), 'CoM'
(barycentric FoR), 'solvent'
(solvent FoR), 'zarate'
(See Memo on definitions of the diffusion coefficient), 'zarate_x'
($D^{(x)}$ as defined by Ortiz de Zárate, doi 10.1140/epje/i2019-11803-2) or 'zarate_w'
($D^{(w)}$ as defined by Ortiz de Zárate).
solvent_idx (int, optional) :
Index of component identified as solvent (only when using frame_of_reference='solvent'
)
Returns:
(ndarray or float) :
Diffusion coefficients, shape varies based on options and number of components. Unit [m^2 / s]
interdiffusion_general(self, T, Vm, x, N=None)
Compute the ‘Kinetic CoM diffusion coefficients’, defined by
\[J_i^{(n, m)} = - \sum_j D_{ij} \nabla n_j, \nabla T = \nabla p = F_k = 0 \forall k\]For end-users, see interdiffusion
See Eq. (19) in RET for Mie fluids (https://doi.org/10.1063/5.0149865)
Args:
T (float) :
Temperature (K)
Vm (float) :
Molar volume (m3 / mol)
x (array_like) :
composition (mole fractions)
N (int, optional) :
Enskog approximation order
Returns:
(2D array) :
Array of the (not independent) $D^{(K, m)}$ diffusion coefficients. Unit [m^2 / s]
kinematic_viscosity(self, T, Vm, x, N=None)
Compute the kinematic viscosity
Args:
T (float) :
Temperature [K]
Vm (float) :
Molar volume [m3 / mol]
x (array_like) :
Molar composition [-]
N (int, optional) :
Enskog approximation order
Returns:
(float) :
The kinematic viscosity of the mixture [m2 / s]
resistivity_matrix(self, T, Vm, x, N=2, formulation='T-psi', frame_of_reference='CoM', use_thermal_conductivity=None)
Compute the resistivity matrix $R = L^{-1}$, for use in NET calculations. The Flux/Force formulation used in the NET
model is selected using the formulation
kwarg. Currently implemented formulations are:
———————————————————————————————-
'T-psi'
:
\[J_q = L_{qq} \nabla (1 / T) - \sum_{i=1}^{N_c-1}(1 / T) L_{qi} \nabla_T \Psi_i\] \[J_i = L_{iq} \nabla (1 / T) - \sum_{j=1}^{N_c-1}(1 / T) L_{ij} \nabla_T \Psi_j\]
Where
\[\Psi_i = \mu_i - \mu_{N_c}\]and $N_c$ denotes the number of components. The last component is used as the dependent component. The fluxes in this formulation are on a mass basis, and the formulation is implemented in the centre of mass frame of reference. The formulation is only implemented for ideal gases. ———————————————————————————————-
Args:
T (float) :
Temperature [K]
Vm (float) :
Molar volume [m3 / mol]
x (1darray) :
Molar composition [-]
N (int, optional) :
Enskog approximation order (must be >= 2 for thermal effects)
formulation (str, optional) :
The NET formulation for which to compute the conductivity matrix.
frame_of_reference (str, optional) :
The frame of reference (‘CoM’, ‘CoN’, ‘CoV’, or ‘solvent’). Default: ‘CoM’.
use_thermal_conductivity (callable, optional) :
External thermal conductivity model. Assumed to have the signaturethermal_conductivity(T, Vm, x), returning the thermal conductivity in [W / m K]. Defaults to None. If no model is supplied, KineticGas is used to compute the thermal conductivity.
Returns:
ndarray :
The resistivity matrix, contents will vary depending on the formulation
kwarg.
soret_coefficient(self, T, Vm, x, N=None, use_zarate=True, dependent_idx=-1)
Compute the Soret coefficients, $S_{T, ij}$. If use_zarate=False
, the Soret coefficient is defined by
where $\alpha_{ij}$ are the thermal diffusion factors. If use_zarate=True
, uses the definition proposed by
Ortiz de Zarate in (doi 10.1140/epje/i2019-11803-2), i.e.
where $(S_T)$ and $(D_T)$ indicate the vectors of Soret coefficients and thermal diffusion coefficients. Or, following the notation in the memo on definitions of diffusion and thermal diffusion coefficients,
\[D^{(z)} (S_T) = (D_T).\]The Soret coefficients defined this way satisfy
\[X (S_T) \nabla T = - (\nabla x)\]and
\[W (S_T) \nabla T = - (\nabla w)\]which in a binary mixture reduces to the same as the definition used when use_zarate=False
, i.e.
if species 2 is the dependent species.
Args:
T (float) :
Temperature [K]
Vm (float) :
Molar volume [m3 / mol]
x (array_like) :
Molar composition
N (int, optional) :
Enskog approximation order (>= 2)
use_zarate (bool, optional) :
Use Ortiz de Zarate formulation of the Soret coefficient, as given in themethod description and (doi 10.1140/epje/i2019-11803-2). Defaults to True
.
dependent_idx (int) :
Only applicable when use_zarate=True
(default behaviour). The index of the dependentspecies. Defaults to the last species.
Returns:
(ndarray) :
Soret coefficient matrix (K$^{-1}$)
thermal_conductivity(self, T, Vm, x, N=None, idealgas=None, contributions='all')
Compute the thermal conductivity, $\lambda$. For models initialized with is_idealgas=True
, the thermal
conductivity is not a function of density (i.e. $d \lambda / d V_m = 0$).
See Eq. (13) in RET for Mie fluids (https://doi.org/10.1063/5.0149865)
Args:
T (float) :
Temperature [K]
Vm (float) :
Molar volume [m3 / mol]
x (array_like) :
Molar composition [-]
N (int, optional) :
Enskog approximation order (>= 2)
idealgas (bool, optional) :
Return infinite dilution value? Defaults to model default (set on init).
contributions (str, optional) :
Return only specific contributions, can be (‘all’ (default), ‘(i)nternal’,’(t)ranslational’, ‘(d)ensity’, or several of the above, such as ‘tid’ or ‘td’. If several contributions are selected, these are returned in an array of contributions in the same order as indicated in the supplied flag.
Returns:
(float) :
The thermal conductivity of the mixture [W / m K].
thermal_diffusion_coeff(self, T, Vm, x, N=None, use_independent=False, dependent_idx=-1, frame_of_reference='CoN', solvent_idx=-1)
Compute thermal diffusion coefficients, $D_{T,i}$ [mol / m^2 s] Default definition is
\[J_i^{(n, n)} = D_{T,i} \nabla \ln T - \sum_j D_{ij} \nabla n_j, \nabla p = F_k = 0 \forall k\]where the flux, $J_i^{(n, n)}$ is on a molar basis, in the molar frame of reference. For fluxes in other frames of reference, use the ‘frame_of_reference’ kwarg. For the diffusion coefficients corresponding to an independent set of forces, defined by
\[J_i^{(n, n)} = D_{T,i} \nabla \ln T - \sum_{j \neq l} D_{ij} \nabla n_j, \nabla p = F_k = 0 \forall k\]where $l$ is the index of the dependent molar density gradient, use the ‘use_independent’ and ‘dependent_idx’ kwargs. See Eq. (23) in RET for Mie fluids (https://doi.org/10.1063/5.0149865)
Args:
T (float) :
Temperature [K]
Vm (float) :
Molar volume [m^3 / mol]
x (array_like) :
Mole fractions [-]
N (int, optional) :
Enskog approximation order (>=2). Defaults to 2.
use_independent (bool, optional) :
Return diffusion coefficients for independent set of forces.
dependent_idx (int, optional) :
Index of the dependent molar density gradient (only if use_dependent=True).Defaults to last component.
frame_of_reference (str, optional) :
What frame of reference the coefficients apply to. Valid options are'CoM'
(centre of mass / barycentric), 'CoN'
(centre of moles), 'CoV'
(centre of volume) 'solvent'
(together with solvent_idx
) or 'zarate'
, for the coefficients as defined by Ortiz de Zarate (doi 10.1140/epje/i2019-11803-2).
solvent_idx (int, optional) :
Index of component identified as solvent (only when using frame_of_reference='solvent'
)
Returns:
(1D array) :
Thermal diffusion coefficients. Unit [mol m^2 / s]
thermal_diffusion_factor(self, T, Vm, x, N=None)
Compute the thermal diffusion factors $\alpha_{ij}$, defined by
\[\alpha_{ij} = k_{T, i} - k_{T, j}\]where $k_{T,i}$ are the thermal diffusion ratios. This definition implies that
\[\nabla \ln (n_i / n_j) = - \alpha_{ij} \nabla \ln T\]when the mass fluxes vanish. The thermal diffusion factors are independent of the frame of reference. See Eq. (29) in RET for Mie fluids (https://doi.org/10.1063/5.0149865)
Args:
T (float) :
Temperature [K]
Vm (float) :
Molar volume [m3 / mol]
x (array_like) :
Molar composition [-]
N (int, optional) :
Enskog approximation order (>= 2)
Returns:
(2D array) :
The thermal diffusion factors of the mixture. Dimensionless.
thermal_diffusion_ratio(self, T, Vm, x, N=None)
Calculate the “independent” thermal diffusion ratios, $k_{T, i}$ defined by
\[J_i^{(n, n)} = - \sum_{j \neq l} D_{ij}^{(I, n)} ( \nabla n_j + n_j k_{T, j} D_{T, j}^{(I, n)} \nabla \ln T )\]and
\[\sum_i x_i \sum_j [\delta_{ij} + (4 \pi / 3) \sigma_{ij}^3 n_j M_{ij} \chi_{ij} - (n_j k_{T,j} / k_B T) ( d \mu_i / n_j )_{T,n_{l \neq j}}] = 0\]This definition implies that
\[\nabla n_j = -n_j k_{T,j} \nabla \ln T\ \forall j\]when all mass fluxes vanish. For models initialised with is_idealgas=True
, the second equation is replaced with
The thermal diffusion ratios are independent of the frame of reference. See Eq. (26-27) in RET for Mie fluids (https://doi.org/10.1063/5.0149865)
Args:
T (float) :
Temperature [K]
Vm (float) :
Molar volume [m3 / mol]
x (array_like) :
Molar composition [-]
N (int, optional) :
Enskog approximation order (>= 2)
Returns:
(1D array) :
The thermal diffusion ratio of each component. Unit Dimensionless.
thermal_diffusivity(self, T, Vm, x, N=None)
Compute the thermal diffusivity
Args:
T (float) :
Temperature [K]
Vm (float) :
Molar volume [m3 / mol]
x (array_like) :
Molar composition [-]
N (int, optional) :
Enskog approximation order
Returns:
(float) :
Thermal diffusivity of the mixture [m2 / s].
viscosity(self, T, Vm, x, N=None, idealgas=None)
Compute the shear viscosity, $\eta$. For models initialized with is_idealgas=True
, the shear viscosity
is not a function of density (i.e. $d \eta / d V_m = 0). See Eq. (12) in RET for Mie fluids (https://doi.org/10.1063/5.0149865)
Args:
T (float) :
Temperature [K]
Vm (float) :
Molar volume [m3 / mol]
x (array_like) :
Molar composition [-]
N (int, optional) :
Enskog approximation order
idealgas (bool, optional) :
Use infinite dilution value? Defaults to model default value (set on init)
Returns:
(float) :
The shear viscosity of the mixture [Pa s].
Tp-property interfaces
Computing properties as a function of temperature and pressure. Simply forwards calls to the appropriate TV-property methods after computing molar volume using the internal equation of state object.
Table of contents
interdiffusion_tp(self, T, p, x, N=None, use_independent=True, dependent_idx=None, frame_of_reference='CoN', use_binary=True, solvent_idx=None)
Compute molar volume using the internal equation of state (self.eos
), assuming vapour, and pass the call to
self.interdiffusion
. See self.interdiffusion
for documentation.
kinematic_viscosity_tp(self, T, p, x, N=None)
Compute molar volume using the internal equation of state (self.eos
), assuming vapour, and pass the call to
self.kinematic_viscosity
. See self.kinematic_viscosity
for documentation.
thermal_conductivity_tp(self, T, p, x, N=None, contributions='all')
Compute molar volume using the internal equation of state (self.eos
), assuming vapour, and pass the call to
self.thermal_conductivity
. See self.thermal_conductivity
for documentation.
thermal_diffusion_coeff_tp(self, T, p, x, N=None, use_independent=False, dependent_idx=None, frame_of_reference='CoN', solvent_idx=None)
Compute molar volume using the internal equation of state (self.eos
), assuming vapour, and pass the call to
self.thermal_diffusion_coeff
. See self.thermal_diffusion_coeff
for documentation.
thermal_diffusion_factor_tp(self, T, p, x, N=None)
Compute molar volume using the internal equation of state (self.eos
), assuming vapour, and pass the call to
self.thermal_diffusion_factor
. See self.thermal_diffusion_factor
for documentation.
thermal_diffusivity_tp(self, T, p, x, N=None)
Compute molar volume using the internal equation of state (self.eos
), assuming vapour, and pass the call to
self.thermal_diffusivity
. See self.thermal_diffusivity
for documentation.
viscosity_tp(self, T, p, x, N=None)
Compute molar volume using the internal equation of state (self.eos
), assuming vapour, and pass the call to
self.viscosity
. See self.viscosity
for documentation.
Frame of Reference transformations
Generate matrices for Frame of Reference transformations. See the supportingmaterial of Revised Enskog Theory for Mie fluids for details.
Table of contents
get_com_2_con_matr(self, x)
Get transformation matrix from centre of mass (CoM) to centre of moles (CoN).
Args:
x (array_like) :
Molar composition [-]
Returns:
(2d array) :
Transformation matrix $\Psi^{n \leftmapsto m}$
get_com_2_cov_matr(self, T, Vm, x)
Get centre of mass (CoM) to centre of volume (CoV) transformation matrix
Args:
T (float) :
Temperature [K]
Vm (float) :
Molar volume [m3 / mol]
x (array_like) :
Molar composition [-]
Returns:
2d array :
The transformation matrix $\Psi^{V \leftmapsto m}$.
get_com_2_for_matr(self, T, Vm, x, FoR, **kwargs)
Dispatcher to get a specific ‘change of frame of reference’ matrix for transformation from centre of mass to ‘FoR’. Returns the appropriate matrix for the transformations derived in Appendix A of … to transform a flux from the centre of mass frame of reference to the ‘FoR’ frame of reference by the transformation
\[J^{(n, FoR)} = \psi @ J^{(n, m)}\]where $\psi$ is the matrix returned by this method, and $J$ is the vector of (all) molar fluxes, with the subscript indicating the frame of reference.
Args:
T (float) :
Temperature [K]
Vm (float) :
Molar volume [m3 / mol]
x (array_like) :
Molar composition [-]
Returns:
(2Darray) :
The $N$ x $N$ transformation matrix to transform the fluxes, with $N$ being the number ofcomponents.
get_com_2_solv_matr(self, x, solvent_idx)
Get transformation matrix from centre of mass (CoM) to solvent (solvent) frame of reference
Args:
x (array_like) :
Molar composition [-]
solvent_idx (int) :
The component index of the solvent.
Returns:
2darray :
The transformation matrix $\Psi^{n_k \leftmapsto m}$, where $k$ is the solvent_idx
.
get_solv_2_solv_matr(self, x, prev_solv_idx, new_solv_idx)
Get solvent-to-solvent frame of reference transformation matrix
Args:
x (array_like) :
Molar composition [-]
prev_solv_idx (int) :
Component index of the old (current) solvent
new_solv_idx (int) :
Component index of the new solvent
Returns:
2d array :
The transformation matrix $\Psi^{n_k \leftmapsto n_l}$, where $k$ is new_solv_idx
and $l$ is prev_solv_idx
.
get_zarate_W_matr(self, x, dependent_idx)
Compute the matrix $W$ as defined by Zárate. See (Definition of frame-invariant thermodiffusion and Soret coefficients for ternary mixtures) and memo on diffusion coefficient definitions.
Args:
x (array_like) :
Molar composition [-]
dependent_idx (int) :
Index of the dependent species
Returns:
2d array :
The transformation matrix $W$
get_zarate_X_matr(self, x, dependent_idx)
Compute the matrix $X$ as defined by Zárate. See (Definition of frame-invariant thermodiffusion and Soret coefficients for ternary mixtures) and memo on diffusion coefficient definitions.
Args:
x (array_like) :
Molar composition [-]
dependent_idx (int) :
Index of the dependent species
Returns:
2d array :
The transformation matrix $X$
Interfaces to C++ methods
Lightweight wrappers for the most commonly used C++ side methods.
Table of contents
get_conductivity_matrix(self, particle_density, T, mole_fracs, N=None)
Compute the elements of the matrix corresponding to the set of equations that must be solved for the thermal response function sonine polynomial expansion coefficients: Eq. (6) in RET for Mie fluids (https://doi.org/10.1063/5.0149865)
Args:
particle_density (float) :
Particle density (not molar!) [1 / m3]
T (float) :
Temperature [K]
mole_fracs (list[float]) :
Molar composition [-]
N (Optional, int) :
Enskog approximation order.
Returns:
2Darray :
($N N_c$ x $N N_c$) matrix, where $N$ is the Enskog approximation order and $N_c$ isthe number of components.
get_conductivity_vector(self, particle_density, T, mole_fracs, N)
Compute the right-hand side vector to the set of equations that must be solved for the thermal response function Sonine polynomial expansion coefficients: Eq. (6) in RET for Mie fluids (https://doi.org/10.1063/5.0149865)
Args:
particle_density (float) :
Particle density (not molar!) [1 / m3]
T (float) :
Temperature [K]
**mole_fracs (list
Molar composition [-]
N (Optional, int) :
Enskog approximation order.
Returns:
(1Darray) :
($N N_c$,) vector, where $N$ is the Enskog approximation order and $N_c$ isthe number of components.
get_diffusion_vector(self, particle_density, T, mole_fracs, N=None)
Compute the right-hand side vector to the set of equations that must be solved for the diffusive response function Sonine polynomial expansion coefficients. Eq. (10) in RET for Mie fluids (https://doi.org/10.1063/5.0149865)
Args:
particle_density (float) :
Particle density (not molar!) [1 / m3]
T (float) :
Temperature [K]
mole_fracs (list[float]) :
Molar composition [-]
N (Optional, int) :
Enskog approximation order.
Returns:
(1Darray) :
($N N_c^2$,) vector, where $N$ is the Enskog approximation order and $N_c$ isthe number of components.
get_etl(self, particle_density, T, x)
Compute energy transfer lengths
Args:
particle_density (float) :
Particle density (not molar!) [1 / m3]
T (float) :
Temperature [K]
x (list[float]) :
Molar composition [-]
Returns:
2d array :
Collision diameters [m], indexed by component pair.
get_mtl(self, particle_density, T, x)
Compute energy transfer lengths
Args:
particle_density (float) :
Particle density (not molar!) [1 / m3]
T (float) :
Temperature [K]
x (list[float]) :
Molar composition [-]
Returns:
2d array :
Collision diameters [m], indexed by component pair.
get_rdf(self, particle_density, T, x)
Compute the radial distribution function at contact
Args:
particle_density (float) :
Particle density (not molar!) [1 / m3]
T (float) :
Temperature [K]
x (list[float]) :
Molar composition [-]
Returns:
2d array :
RDF at contact, indexed by component pair.
Utility methods
Methods for various helper computations
Table of contents
get_reducing_units(self, ci=0, cj=None)
Get reducing units for this model, as a Units
struct. See units.py
.
Args:
ci (int, optional) :
Which component to use for reducing units, defaults to first component
cj (int, optional) :
Other species to use, if cross-interactions are used for reducing units. Defaults to None
.
Returns:
Units :
Struct holding the reducing units
get_tl_model(self)
Get the currently active transfer length model
Use get_valid_tl_models
for an overview of valid models.
Returns:
(tuple[int, str]) :
The index (id) of the current model, and a description
get_valid_tl_models(self)
Get valid transfer length model identifiers
Returns:
(dict) :
Keys are the model identifiers (used with set_tl_model
), values are model descriptions.
set_tl_model(self, model)
Set the transfer length model
Use get_valid_tl_models
for an overview of valid models.
Args:
model (int) :
The model id
Raises:
RuntimeError :
If the model id is invalid.
Deprecated methods
Deprecated methods are not maintained, and may be removed in the future.
Table of contents
thermal_coductivity_tp(self, T, p, x, N=None)
Slightly embarrasing typo in method name… Keeping alive for a while because some code out there uses this one.
Internal methods
Internal methods are not intended for use by end-users.
Table of contents
check_valid_composition(self, x)
Check that enough mole fractions are supplied for the initialised model. Also check that they sum to unity. for single-component mixtures, returns [0.5, 0.5], such that [1] can be supplied, even though single-component fluids are treated as binaries internally.
Args:
x (array_like) :
Molar composition
Raises:
IndexError :
If wrong number of mole fractions is supplied.
RuntimeWarning :
If mole fractions do not sum to unity.
Returns:
x (1d array) :
Sanitized composition
compute_cond_vector(self, particle_density, T, mole_fracs, N=None)
Compute the thermal response function Sonine polynomial expansion coefficients by solving the set of equations
\[\Lambda \ell = \lambda\]Corresponding to Eq. (6) in RET for Mie fluids (https://doi.org/10.1063/5.0149865)
Where $\Lambda$ is the matrix returned by the c++ method get_conductivity_matrix
, and $\lambda$ is the vector
returned by the c++ method get_conductivity_vector
.
Args:
particle_density (float) :
Particle density (not molar!) [1 / m3]
T (float) :
Temperature [K]
mole_fracs (list[float]) :
Molar composition [-]
N (Optional, int) :
Enskog approximation order.
Returns:
1Darray :
($N N_c$,) vector, where $N$ is the Enskog approximation order and $N_c$ isthe number of components. The vector is ordered as
l[:
N_c] = [$\ell_1^{(0)}$, $\ell_2^{(0)}$, …, $\ell_{N_c}^{(0)}$]
l[N_c :
2 * N_c] = [$\ell_1^{(1)}$, $\ell_2^{(1)}$, …, $\ell_{N_c}^{(1)}$]… etc … Where subscripts indicate component indices, and superscripts are Enskog approximation order summation indices.
compute_diffusion_coeff_vector(self, particle_density, T, mole_fracs, N=None)
Compute the diffusive response function Sonine polynomial expansion coefficients by solving the set of equations
\[D d = \delta\]Corresponding to Eq. (10) in RET for Mie fluids (https://doi.org/10.1063/5.0149865)
Where $D$ is the matrix returned by the c++ method get_diffusion_matrix
, and $\delta$ is the vector
returned by the c++ method get_diffusion_vector
.
Args:
particle_density (float) :
Particle density (not molar!) [1 / m3]
T (float) :
Temperature [K]
mole_fracs (list[float]) :
Molar composition [-]
N (Optional, int) :
Enskog approximation order.
Returns:
1Darray :
($N N_c^2$,) vector, where $N$ is the Enskog approximation order and $N_c$ isthe number of components, containing the diffusive response function sonine polynomial expansion coefficients ($d_{i, j}^{(q)}$). See reshape_diffusive_coeff_vector
for help on practical usage.
compute_dth_vector(self, particle_density, T, mole_fracs, N=None)
Compute the coefficients $d_i^{(J = 0)}$, by solving the set of equations
\(\sum_j d_{i,j}^{(0)} d_j^{\vec{J} = 0} = \ell_i^{(0)}\),
i.e. Eq. (15) in RET for Mie fluids (https://doi.org/10.1063/5.0149865)
Args:
particle_density (float) :
Particle density (not molar!) [1 / m3]
T (float) :
Temperature [K]
**mole_fracs (list
Molar composition [-]
N (Optional, int) :
Enskog approximation order.
Returns:
1Darray :
Vector of the $d_i^{(J = 0)}$ coefficients.
compute_visc_vector(self, T, particle_density, mole_fracs, N=None)
Compute the viscous response function Sonine polynomial expansion coefficients by solving the set of equations
\[\Beta b = \beta\]Corresponding to Eq. (8) in RET for Mie fluids (https://doi.org/10.1063/5.0149865)
Where $\Beta$ is the matrix returned by the c++ method get_viscosity_matrix
, and $\beta$ is the vector
returned by the c++ method get_viscosity_vector
.
Args:
particle_density (float) :
Particle density (not molar!) [1 / m3]
T (float) :
Temperature [K]
**mole_fracs (list
Molar composition [-]
N (Optional, int) :
Enskog approximation order.
Returns:
1D array :
($N N_c$,) vector, where $N$ is the Enskog approximation order and $N_c$ isthe number of components, ordered as
b[:
N_c] = [b_1^{(0)}, b_2^{(0)}, …, b_{N_c}^{(0)}]
b[N_c:
2 * N_c] = [b_1^{(1)}, b_2^{(1)}, …, b_{N_c}^{(1)}]… etc … where subscripts denote component indices and superscripts denote Enskog approximation summation indices.
get_Eij(self, Vm, T, x)
Compute the factors
\[( n_i / k_B T ) (d \mu_i / d n_j)_{T, n_{k \neq j}},\]where $n_i$ is the molar density of species $i$.
Args:
Vm (float) :
Molar volume [m3 / mol]
T (float) :
Temperature [K]
x (array_like) :
Molar composition
Returns:
(2D array) :
The factors E[i][j] = $ ( n_i / k_B T ) (d \mu_i / d n_j){T, n{k \neq j}}$, where $n_i$is the molar density of species $i$. Unit [1 / mol]
get_P_factors(self, Vm, T, x)
Compute the factors $\Xi_i = \sum_j E_{ij}$, where $E_{ij}$ are the factors computed by get_Eij
.
Args:
Vm (float) :
Molar volume [m3 / mol]
T (float) :
Temperature [K]
x (array_like) :
Molar composition
Returns:
(1D array) :
The factors $\Xi_i$, Unit [1 / mol]
reshape_diffusion_coeff_vector(self, d)
The vector returned by compute_diffusion_coeff_vector
contains the diffusive response function sonine
polynomial expansion coefficients (eg. $d_{i, j}^{(q)}$ ).
To more easily access the correct coefficients, this method reshapes the vector to a matrix indiced
as d[i][q][j]
where i and j are component indices, and q refferes to the approximation order summation index.
Args:
d (1D array) :
The array returned by get_diffusion_coeff_vector
Returns:
3D array :
The matrix of $d_{i, j}^{(q)}$ coefficients ordered as d[i][q][j]