pub struct State<E> {
pub eos: Arc<E>,
pub temperature: Temperature,
pub volume: Volume,
pub moles: Moles<Array1<f64>>,
pub total_moles: Moles,
pub partial_density: Density<Array1<f64>>,
pub density: Density,
pub molefracs: Array1<f64>,
/* private fields */
}
Expand description
Thermodynamic state of the system.
The state is always specified by the variables of the Helmholtz energy: volume $V$, temperature $T$ and mole numbers $N_i$. Additional to these variables, the state saves properties like the density, that can be calculated directly from the basic variables. The state also contains a reference to the equation of state used to create the state. Therefore, it can be used directly to calculate all state properties.
Calculated partial derivatives are cached in the state. Therefore, the second evaluation of a property like the pressure, does not require a recalculation of the equation of state. This can be used in situations where both lower and higher order derivatives are required, as in a calculation of a derivative all lower derivatives have to be calculated internally as well. Since they are cached it is more efficient to calculate the highest derivatives first. For example during the calculation of the isochoric heat capacity $c_v$, the entropy and the Helmholtz energy are calculated as well.
State
objects are meant to be immutable. If individual fields like volume
are changed, the
calculations are wrong as the internal fields of the state are not updated.
§Contents
Fields§
§eos: Arc<E>
Equation of state
temperature: Temperature
Temperature $T$
volume: Volume
Volume $V$
moles: Moles<Array1<f64>>
Mole numbers $N_i$
total_moles: Moles
Total number of moles $N=\sum_iN_i$
partial_density: Density<Array1<f64>>
Partial densities $\rho_i=\frac{N_i}{V}$
density: Density
Total density $\rho=\frac{N}{V}=\sum_i\rho_i$
molefracs: Array1<f64>
Mole fractions $x_i=\frac{N_i}{N}=\frac{\rho_i}{\rho}$
Implementations§
Source§impl<E: Residual> State<E>
§Stability analysis
impl<E: Residual> State<E>
§Stability analysis
Sourcepub fn is_stable(&self, options: SolverOptions) -> EosResult<bool>
pub fn is_stable(&self, options: SolverOptions) -> EosResult<bool>
Determine if the state is stable, i.e. if a phase split should occur or not.
Sourcepub fn stability_analysis(
&self,
options: SolverOptions,
) -> EosResult<Vec<State<E>>>
pub fn stability_analysis( &self, options: SolverOptions, ) -> EosResult<Vec<State<E>>>
Perform a stability analysis. The result is a list of States with negative tangent plane distance (i.e. lower Gibbs energy) that can be used as initial estimates for a phase equilibrium calculation.
Source§impl<E: Residual> State<E>
§Flash calculations
impl<E: Residual> State<E>
§Flash calculations
Sourcepub fn tp_flash(
&self,
initial_state: Option<&PhaseEquilibrium<E, 2>>,
options: SolverOptions,
non_volatile_components: Option<Vec<usize>>,
) -> EosResult<PhaseEquilibrium<E, 2>>
pub fn tp_flash( &self, initial_state: Option<&PhaseEquilibrium<E, 2>>, options: SolverOptions, non_volatile_components: Option<Vec<usize>>, ) -> EosResult<PhaseEquilibrium<E, 2>>
Perform a Tp-flash calculation using the State as feed. If no initial values are given, the solution is initialized using a stability analysis.
The algorithm can be use to calculate phase equilibria of systems containing non-volatile components (e.g. ions).
pub fn tp_flash_( &self, new_vle_state: PhaseEquilibrium<E, 2>, options: SolverOptions, non_volatile_components: Option<Vec<usize>>, ) -> EosResult<PhaseEquilibrium<E, 2>>
Source§impl<E: Residual + IdealGas> State<E>
impl<E: Residual + IdealGas> State<E>
Sourcepub fn chemical_potential(
&self,
contributions: Contributions,
) -> MolarEnergy<Array1<f64>>
pub fn chemical_potential( &self, contributions: Contributions, ) -> MolarEnergy<Array1<f64>>
Chemical potential: $\mu_i=\left(\frac{\partial A}{\partial N_i}\right)_{T,V,N_j}$
Sourcepub fn dmu_dt(
&self,
contributions: Contributions,
) -> <MolarEnergy<Array1<f64>> as Div<Temperature>>::Output
pub fn dmu_dt( &self, contributions: Contributions, ) -> <MolarEnergy<Array1<f64>> as Div<Temperature>>::Output
Partial derivative of chemical potential w.r.t. temperature: $\left(\frac{\partial\mu_i}{\partial T}\right)_{V,N_i}$
Sourcepub fn molar_isochoric_heat_capacity(
&self,
contributions: Contributions,
) -> MolarEntropy
pub fn molar_isochoric_heat_capacity( &self, contributions: Contributions, ) -> MolarEntropy
Molar isochoric heat capacity: $c_v=\left(\frac{\partial u}{\partial T}\right)_{V,N_i}$
Sourcepub fn dc_v_dt(
&self,
contributions: Contributions,
) -> <MolarEntropy as Div<Temperature>>::Output
pub fn dc_v_dt( &self, contributions: Contributions, ) -> <MolarEntropy as Div<Temperature>>::Output
Partial derivative of the molar isochoric heat capacity w.r.t. temperature: $\left(\frac{\partial c_V}{\partial T}\right)_{V,N_i}$
Sourcepub fn molar_isobaric_heat_capacity(
&self,
contributions: Contributions,
) -> MolarEntropy
pub fn molar_isobaric_heat_capacity( &self, contributions: Contributions, ) -> MolarEntropy
Molar isobaric heat capacity: $c_p=\left(\frac{\partial h}{\partial T}\right)_{p,N_i}$
Sourcepub fn entropy(&self, contributions: Contributions) -> Entropy
pub fn entropy(&self, contributions: Contributions) -> Entropy
Entropy: $S=-\left(\frac{\partial A}{\partial T}\right)_{V,N_i}$
Sourcepub fn molar_entropy(&self, contributions: Contributions) -> MolarEntropy
pub fn molar_entropy(&self, contributions: Contributions) -> MolarEntropy
Molar entropy: $s=\frac{S}{N}$
Sourcepub fn partial_molar_entropy(&self) -> MolarEntropy<Array1<f64>>
pub fn partial_molar_entropy(&self) -> MolarEntropy<Array1<f64>>
Partial molar entropy: $s_i=\left(\frac{\partial S}{\partial N_i}\right)_{T,p,N_j}$
Sourcepub fn ds_dt(
&self,
contributions: Contributions,
) -> <Entropy as Div<Temperature>>::Output
pub fn ds_dt( &self, contributions: Contributions, ) -> <Entropy as Div<Temperature>>::Output
Partial derivative of the entropy w.r.t. temperature: $\left(\frac{\partial S}{\partial T}\right)_{V,N_i}$
Sourcepub fn d2s_dt2(
&self,
contributions: Contributions,
) -> <<Entropy as Div<Temperature>>::Output as Div<Temperature>>::Output
pub fn d2s_dt2( &self, contributions: Contributions, ) -> <<Entropy as Div<Temperature>>::Output as Div<Temperature>>::Output
Second partial derivative of the entropy w.r.t. temperature: $\left(\frac{\partial^2 S}{\partial T^2}\right)_{V,N_i}$
Sourcepub fn enthalpy(&self, contributions: Contributions) -> Energy
pub fn enthalpy(&self, contributions: Contributions) -> Energy
Enthalpy: $H=A+TS+pV$
Sourcepub fn molar_enthalpy(&self, contributions: Contributions) -> MolarEnergy
pub fn molar_enthalpy(&self, contributions: Contributions) -> MolarEnergy
Molar enthalpy: $h=\frac{H}{N}$
Sourcepub fn partial_molar_enthalpy(&self) -> MolarEnergy<Array1<f64>>
pub fn partial_molar_enthalpy(&self) -> MolarEnergy<Array1<f64>>
Partial molar enthalpy: $h_i=\left(\frac{\partial H}{\partial N_i}\right)_{T,p,N_j}$
Sourcepub fn helmholtz_energy(&self, contributions: Contributions) -> Energy
pub fn helmholtz_energy(&self, contributions: Contributions) -> Energy
Helmholtz energy: $A$
Sourcepub fn molar_helmholtz_energy(
&self,
contributions: Contributions,
) -> MolarEnergy
pub fn molar_helmholtz_energy( &self, contributions: Contributions, ) -> MolarEnergy
Molar Helmholtz energy: $a=\frac{A}{N}$
Sourcepub fn internal_energy(&self, contributions: Contributions) -> Energy
pub fn internal_energy(&self, contributions: Contributions) -> Energy
Internal energy: $U=A+TS$
Sourcepub fn molar_internal_energy(&self, contributions: Contributions) -> MolarEnergy
pub fn molar_internal_energy(&self, contributions: Contributions) -> MolarEnergy
Molar internal energy: $u=\frac{U}{N}$
Sourcepub fn gibbs_energy(&self, contributions: Contributions) -> Energy
pub fn gibbs_energy(&self, contributions: Contributions) -> Energy
Gibbs energy: $G=A+pV$
Sourcepub fn molar_gibbs_energy(&self, contributions: Contributions) -> MolarEnergy
pub fn molar_gibbs_energy(&self, contributions: Contributions) -> MolarEnergy
Molar Gibbs energy: $g=\frac{G}{N}$
Sourcepub fn joule_thomson(&self) -> <Temperature as Div<Pressure>>::Output
pub fn joule_thomson(&self) -> <Temperature as Div<Pressure>>::Output
Joule Thomson coefficient: $\mu_{JT}=\left(\frac{\partial T}{\partial p}\right)_{H,N_i}$
Sourcepub fn isentropic_compressibility(&self) -> <f64 as Div<Pressure>>::Output
pub fn isentropic_compressibility(&self) -> <f64 as Div<Pressure>>::Output
Isentropic compressibility: $\kappa_s=-\frac{1}{V}\left(\frac{\partial V}{\partial p}\right)_{S,N_i}$
Sourcepub fn isenthalpic_compressibility(&self) -> <f64 as Div<Pressure>>::Output
pub fn isenthalpic_compressibility(&self) -> <f64 as Div<Pressure>>::Output
Isenthalpic compressibility: $\kappa_H=-\frac{1}{V}\left(\frac{\partial V}{\partial p}\right)_{H,N_i}$
Sourcepub fn thermal_expansivity(&self) -> <f64 as Div<Temperature>>::Output
pub fn thermal_expansivity(&self) -> <f64 as Div<Temperature>>::Output
Thermal expansivity: $\alpha_p=-\frac{1}{V}\left(\frac{\partial V}{\partial T}\right)_{p,N_i}$
Sourcepub fn grueneisen_parameter(&self) -> f64
pub fn grueneisen_parameter(&self) -> f64
Grueneisen parameter: $\phi=V\left(\frac{\partial p}{\partial U}\right){V,n_i}=\frac{v}{c_v}\left(\frac{\partial p}{\partial T}\right){v,n_i}=\frac{\rho}{T}\left(\frac{\partial T}{\partial \rho}\right)_{s, n_i}$
Sourcepub fn chemical_potential_contributions(
&self,
component: usize,
contributions: Contributions,
) -> Vec<(String, MolarEnergy)>
pub fn chemical_potential_contributions( &self, component: usize, contributions: Contributions, ) -> Vec<(String, MolarEnergy)>
Chemical potential $\mu_i$ evaluated for each contribution of the equation of state.
Source§impl<E: Residual + Molarweight + IdealGas> State<E>
impl<E: Residual + Molarweight + IdealGas> State<E>
Sourcepub fn specific_isochoric_heat_capacity(
&self,
contributions: Contributions,
) -> SpecificEntropy
pub fn specific_isochoric_heat_capacity( &self, contributions: Contributions, ) -> SpecificEntropy
Specific isochoric heat capacity: $c_v^{(m)}=\frac{C_v}{m}$
Sourcepub fn specific_isobaric_heat_capacity(
&self,
contributions: Contributions,
) -> SpecificEntropy
pub fn specific_isobaric_heat_capacity( &self, contributions: Contributions, ) -> SpecificEntropy
Specific isobaric heat capacity: $c_p^{(m)}=\frac{C_p}{m}$
Sourcepub fn specific_entropy(&self, contributions: Contributions) -> SpecificEntropy
pub fn specific_entropy(&self, contributions: Contributions) -> SpecificEntropy
Specific entropy: $s^{(m)}=\frac{S}{m}$
Sourcepub fn specific_enthalpy(&self, contributions: Contributions) -> SpecificEnergy
pub fn specific_enthalpy(&self, contributions: Contributions) -> SpecificEnergy
Specific enthalpy: $h^{(m)}=\frac{H}{m}$
Sourcepub fn specific_helmholtz_energy(
&self,
contributions: Contributions,
) -> SpecificEnergy
pub fn specific_helmholtz_energy( &self, contributions: Contributions, ) -> SpecificEnergy
Specific Helmholtz energy: $a^{(m)}=\frac{A}{m}$
Sourcepub fn specific_internal_energy(
&self,
contributions: Contributions,
) -> SpecificEnergy
pub fn specific_internal_energy( &self, contributions: Contributions, ) -> SpecificEnergy
Specific internal energy: $u^{(m)}=\frac{U}{m}$
Sourcepub fn specific_gibbs_energy(
&self,
contributions: Contributions,
) -> SpecificEnergy
pub fn specific_gibbs_energy( &self, contributions: Contributions, ) -> SpecificEnergy
Specific Gibbs energy: $g^{(m)}=\frac{G}{m}$
Sourcepub fn speed_of_sound(&self) -> Velocity
pub fn speed_of_sound(&self) -> Velocity
Speed of sound: $c=\sqrt{\left(\frac{\partial p}{\partial\rho^{(m)}}\right)_{S,N_i}}$
Source§impl<E: Residual> State<E>
impl<E: Residual> State<E>
Sourcepub fn residual_helmholtz_energy(&self) -> Energy
pub fn residual_helmholtz_energy(&self) -> Energy
Residual Helmholtz energy $A^\text{res}$
Sourcepub fn residual_molar_helmholtz_energy(&self) -> MolarEnergy
pub fn residual_molar_helmholtz_energy(&self) -> MolarEnergy
Residual molar Helmholtz energy $a^\text{res}$
Sourcepub fn residual_helmholtz_energy_contributions(&self) -> Vec<(String, Energy)>
pub fn residual_helmholtz_energy_contributions(&self) -> Vec<(String, Energy)>
Residual Helmholtz energy $A^\text{res}$ evaluated for each contribution of the equation of state.
Sourcepub fn residual_entropy(&self) -> Entropy
pub fn residual_entropy(&self) -> Entropy
Residual entropy $S^\text{res}=\left(\frac{\partial A^\text{res}}{\partial T}\right)_{V,N_i}$
Sourcepub fn residual_molar_entropy(&self) -> MolarEntropy
pub fn residual_molar_entropy(&self) -> MolarEntropy
Residual entropy $s^\text{res}=\left(\frac{\partial a^\text{res}}{\partial T}\right)_{V,N_i}$
Sourcepub fn pressure(&self, contributions: Contributions) -> Pressure
pub fn pressure(&self, contributions: Contributions) -> Pressure
Pressure: $p=-\left(\frac{\partial A}{\partial V}\right)_{T,N_i}$
Sourcepub fn residual_chemical_potential(&self) -> MolarEnergy<Array1<f64>>
pub fn residual_chemical_potential(&self) -> MolarEnergy<Array1<f64>>
Residual chemical potential: $\mu_i^\text{res}=\left(\frac{\partial A^\text{res}}{\partial N_i}\right)_{T,V,N_j}$
Sourcepub fn residual_chemical_potential_contributions(
&self,
component: usize,
) -> Vec<(String, MolarEnergy)>
pub fn residual_chemical_potential_contributions( &self, component: usize, ) -> Vec<(String, MolarEnergy)>
Chemical potential $\mu_i^\text{res}$ evaluated for each contribution of the equation of state.
Sourcepub fn compressibility(&self, contributions: Contributions) -> f64
pub fn compressibility(&self, contributions: Contributions) -> f64
Compressibility factor: $Z=\frac{pV}{NRT}$
Sourcepub fn dp_dv(
&self,
contributions: Contributions,
) -> <Pressure as Div<Volume>>::Output
pub fn dp_dv( &self, contributions: Contributions, ) -> <Pressure as Div<Volume>>::Output
Partial derivative of pressure w.r.t. volume: $\left(\frac{\partial p}{\partial V}\right)_{T,N_i}$
Sourcepub fn dp_drho(
&self,
contributions: Contributions,
) -> <Pressure as Div<Density>>::Output
pub fn dp_drho( &self, contributions: Contributions, ) -> <Pressure as Div<Density>>::Output
Partial derivative of pressure w.r.t. density: $\left(\frac{\partial p}{\partial \rho}\right)_{T,N_i}$
Sourcepub fn dp_dt(
&self,
contributions: Contributions,
) -> <Pressure as Div<Temperature>>::Output
pub fn dp_dt( &self, contributions: Contributions, ) -> <Pressure as Div<Temperature>>::Output
Partial derivative of pressure w.r.t. temperature: $\left(\frac{\partial p}{\partial T}\right)_{V,N_i}$
Sourcepub fn dp_dni(
&self,
contributions: Contributions,
) -> <Pressure as Div<Moles<Array1<f64>>>>::Output
pub fn dp_dni( &self, contributions: Contributions, ) -> <Pressure as Div<Moles<Array1<f64>>>>::Output
Partial derivative of pressure w.r.t. moles: $\left(\frac{\partial p}{\partial N_i}\right)_{T,V,N_j}$
Sourcepub fn d2p_dv2(
&self,
contributions: Contributions,
) -> <<Pressure as Div<Volume>>::Output as Div<Volume>>::Output
pub fn d2p_dv2( &self, contributions: Contributions, ) -> <<Pressure as Div<Volume>>::Output as Div<Volume>>::Output
Second partial derivative of pressure w.r.t. volume: $\left(\frac{\partial^2 p}{\partial V^2}\right)_{T,N_j}$
Sourcepub fn d2p_drho2(
&self,
contributions: Contributions,
) -> <<Pressure as Div<Density>>::Output as Div<Density>>::Output
pub fn d2p_drho2( &self, contributions: Contributions, ) -> <<Pressure as Div<Density>>::Output as Div<Density>>::Output
Second partial derivative of pressure w.r.t. density: $\left(\frac{\partial^2 p}{\partial \rho^2}\right)_{T,N_j}$
Sourcepub fn structure_factor(&self) -> f64
pub fn structure_factor(&self) -> f64
Structure factor: $S(0)=k_BT\left(\frac{\partial\rho}{\partial p}\right)_{T,N_i}$
Sourcepub fn partial_molar_volume(&self) -> MolarVolume<Array1<f64>>
pub fn partial_molar_volume(&self) -> MolarVolume<Array1<f64>>
Partial molar volume: $v_i=\left(\frac{\partial V}{\partial N_i}\right)_{T,p,N_j}$
Sourcepub fn dmu_dni(
&self,
contributions: Contributions,
) -> <MolarEnergy<Array2<f64>> as Div<Moles>>::Output
pub fn dmu_dni( &self, contributions: Contributions, ) -> <MolarEnergy<Array2<f64>> as Div<Moles>>::Output
Partial derivative of chemical potential w.r.t. moles: $\left(\frac{\partial\mu_i}{\partial N_j}\right)_{T,V,N_k}$
Sourcepub fn isothermal_compressibility(&self) -> <f64 as Div<Pressure>>::Output
pub fn isothermal_compressibility(&self) -> <f64 as Div<Pressure>>::Output
Isothermal compressibility: $\kappa_T=-\frac{1}{V}\left(\frac{\partial V}{\partial p}\right)_{T,N_i}$
Sourcepub fn pressure_contributions(&self) -> Vec<(String, Pressure)>
pub fn pressure_contributions(&self) -> Vec<(String, Pressure)>
Pressure $p$ evaluated for each contribution of the equation of state.
Sourcepub fn ds_res_dt(&self) -> <Entropy as Div<Temperature>>::Output
pub fn ds_res_dt(&self) -> <Entropy as Div<Temperature>>::Output
Partial derivative of the residual entropy w.r.t. temperature: $\left(\frac{\partial S^\text{res}}{\partial T}\right)_{V,N_i}$
Sourcepub fn d2s_res_dt2(
&self,
) -> <<Entropy as Div<Temperature>>::Output as Div<Temperature>>::Output
pub fn d2s_res_dt2( &self, ) -> <<Entropy as Div<Temperature>>::Output as Div<Temperature>>::Output
Second partial derivative of the residual entropy w.r.t. temperature: $\left(\frac{\partial^2S^\text{res}}{\partial T^2}\right)_{V,N_i}$
Sourcepub fn dmu_res_dt(
&self,
) -> <MolarEnergy<Array1<f64>> as Div<Temperature>>::Output
pub fn dmu_res_dt( &self, ) -> <MolarEnergy<Array1<f64>> as Div<Temperature>>::Output
Partial derivative of chemical potential w.r.t. temperature: $\left(\frac{\partial\mu_i}{\partial T}\right)_{V,N_i}$
Sourcepub fn ln_phi(&self) -> Array1<f64>
pub fn ln_phi(&self) -> Array1<f64>
Logarithm of the fugacity coefficient: $\ln\varphi_i=\beta\mu_i^\mathrm{res}\left(T,p,\lbrace N_i\rbrace\right)$
Sourcepub fn ln_phi_pure_liquid(&self) -> EosResult<Array1<f64>>
pub fn ln_phi_pure_liquid(&self) -> EosResult<Array1<f64>>
Logarithm of the fugacity coefficient of all components treated as pure substance at mixture temperature and pressure.
Sourcepub fn ln_symmetric_activity_coefficient(&self) -> EosResult<Array1<f64>>
pub fn ln_symmetric_activity_coefficient(&self) -> EosResult<Array1<f64>>
Activity coefficient $\ln \gamma_i = \ln \varphi_i(T, p, \mathbf{N}) - \ln \varphi_i^\mathrm{pure}(T, p)$
Sourcepub fn henrys_law_constant(
eos: &Arc<E>,
temperature: Temperature,
molefracs: &Array1<f64>,
) -> EosResult<Pressure<Array1<f64>>>
pub fn henrys_law_constant( eos: &Arc<E>, temperature: Temperature, molefracs: &Array1<f64>, ) -> EosResult<Pressure<Array1<f64>>>
Henry’s law constant $H_{i,s}=\lim_{x_i\to 0}\frac{y_ip}{x_i}=p_s^\mathrm{sat}\frac{\varphi_i^{\infty,\mathrm{L}}}{\varphi_i^{\infty,\mathrm{V}}}$
The composition of the (possibly mixed) solvent is determined by the molefracs. All components for which the composition is 0 are treated as solutes.
Sourcepub fn henrys_law_constant_binary(
eos: &Arc<E>,
temperature: Temperature,
) -> EosResult<Pressure>
pub fn henrys_law_constant_binary( eos: &Arc<E>, temperature: Temperature, ) -> EosResult<Pressure>
Henry’s law constant $H_{i,s}=\lim_{x_i\to 0}\frac{y_ip}{x_i}=p_s^\mathrm{sat}\frac{\varphi_i^{\infty,\mathrm{L}}}{\varphi_i^{\infty,\mathrm{V}}}$ for a binary system
The solute (i) is the first component and the solvent (s) the second component.
Sourcepub fn dln_phi_dt(&self) -> <f64 as Div<Temperature<Array1<f64>>>>::Output
pub fn dln_phi_dt(&self) -> <f64 as Div<Temperature<Array1<f64>>>>::Output
Partial derivative of the logarithm of the fugacity coefficient w.r.t. temperature: $\left(\frac{\partial\ln\varphi_i}{\partial T}\right)_{p,N_i}$
Sourcepub fn dln_phi_dp(&self) -> <f64 as Div<Pressure<Array1<f64>>>>::Output
pub fn dln_phi_dp(&self) -> <f64 as Div<Pressure<Array1<f64>>>>::Output
Partial derivative of the logarithm of the fugacity coefficient w.r.t. pressure: $\left(\frac{\partial\ln\varphi_i}{\partial p}\right)_{T,N_i}$
Sourcepub fn dln_phi_dnj(&self) -> <f64 as Div<Moles<Array2<f64>>>>::Output
pub fn dln_phi_dnj(&self) -> <f64 as Div<Moles<Array2<f64>>>>::Output
Partial derivative of the logarithm of the fugacity coefficient w.r.t. moles: $\left(\frac{\partial\ln\varphi_i}{\partial N_j}\right)_{T,p,N_k}$
Sourcepub fn thermodynamic_factor(&self) -> Array2<f64>
pub fn thermodynamic_factor(&self) -> Array2<f64>
Thermodynamic factor: $\Gamma_{ij}=\delta_{ij}+x_i\left(\frac{\partial\ln\varphi_i}{\partial x_j}\right)_{T,p,\Sigma}$
Sourcepub fn residual_molar_isochoric_heat_capacity(&self) -> MolarEntropy
pub fn residual_molar_isochoric_heat_capacity(&self) -> MolarEntropy
Residual molar isochoric heat capacity: $c_v^\text{res}=\left(\frac{\partial u^\text{res}}{\partial T}\right)_{V,N_i}$
Sourcepub fn dc_v_res_dt(&self) -> <MolarEntropy as Div<Temperature>>::Output
pub fn dc_v_res_dt(&self) -> <MolarEntropy as Div<Temperature>>::Output
Partial derivative of the residual molar isochoric heat capacity w.r.t. temperature: $\left(\frac{\partial c_V^\text{res}}{\partial T}\right)_{V,N_i}$
Sourcepub fn residual_molar_isobaric_heat_capacity(&self) -> MolarEntropy
pub fn residual_molar_isobaric_heat_capacity(&self) -> MolarEntropy
Residual molar isobaric heat capacity: $c_p^\text{res}=\left(\frac{\partial h^\text{res}}{\partial T}\right)_{p,N_i}$
Sourcepub fn residual_enthalpy(&self) -> Energy
pub fn residual_enthalpy(&self) -> Energy
Residual enthalpy: $H^\text{res}(T,p,\mathbf{n})=A^\text{res}+TS^\text{res}+p^\text{res}V$
Sourcepub fn residual_molar_enthalpy(&self) -> MolarEnergy
pub fn residual_molar_enthalpy(&self) -> MolarEnergy
Residual molar enthalpy: $h^\text{res}(T,p,\mathbf{n})=a^\text{res}+Ts^\text{res}+p^\text{res}v$
Sourcepub fn residual_internal_energy(&self) -> Energy
pub fn residual_internal_energy(&self) -> Energy
Residual internal energy: $U^\text{res}(T,V,\mathbf{n})=A^\text{res}+TS^\text{res}$
Sourcepub fn residual_molar_internal_energy(&self) -> MolarEnergy
pub fn residual_molar_internal_energy(&self) -> MolarEnergy
Residual molar internal energy: $u^\text{res}(T,V,\mathbf{n})=a^\text{res}+Ts^\text{res}$
Sourcepub fn residual_gibbs_energy(&self) -> Energy
pub fn residual_gibbs_energy(&self) -> Energy
Residual Gibbs energy: $G^\text{res}(T,p,\mathbf{n})=A^\text{res}+p^\text{res}V-NRT \ln Z$
Sourcepub fn residual_molar_gibbs_energy(&self) -> MolarEnergy
pub fn residual_molar_gibbs_energy(&self) -> MolarEnergy
Residual Gibbs energy: $g^\text{res}(T,p,\mathbf{n})=a^\text{res}+p^\text{res}v-RT \ln Z$
Source§impl<E: Residual + Molarweight> State<E>
impl<E: Residual + Molarweight> State<E>
Sourcepub fn total_molar_weight(&self) -> MolarWeight
pub fn total_molar_weight(&self) -> MolarWeight
Total molar weight: $MW=\sum_ix_iMW_i$
Sourcepub fn total_mass(&self) -> Mass
pub fn total_mass(&self) -> Mass
Total mass: $m=\sum_im_i=nMW$
Sourcepub fn mass_density(&self) -> MassDensity
pub fn mass_density(&self) -> MassDensity
Mass density: $\rho^{(m)}=\frac{m}{V}$
Source§impl<E: Residual + EntropyScaling> State<E>
§Transport properties
These properties are available for equations of state
that implement the EntropyScaling trait.
impl<E: Residual + EntropyScaling> State<E>
§Transport properties
These properties are available for equations of state that implement the EntropyScaling trait.
Sourcepub fn ln_viscosity_reduced(&self) -> EosResult<f64>
pub fn ln_viscosity_reduced(&self) -> EosResult<f64>
Return the logarithm of the reduced viscosity.
This term equals the viscosity correlation function that is used for entropy scaling.
Sourcepub fn viscosity_reference(&self) -> EosResult<Viscosity>
pub fn viscosity_reference(&self) -> EosResult<Viscosity>
Return the viscosity reference as used in entropy scaling.
Sourcepub fn diffusion(&self) -> EosResult<Diffusivity>
pub fn diffusion(&self) -> EosResult<Diffusivity>
Return the diffusion via entropy scaling.
Sourcepub fn ln_diffusion_reduced(&self) -> EosResult<f64>
pub fn ln_diffusion_reduced(&self) -> EosResult<f64>
Return the logarithm of the reduced diffusion.
This term equals the diffusion correlation function that is used for entropy scaling.
Sourcepub fn diffusion_reference(&self) -> EosResult<Diffusivity>
pub fn diffusion_reference(&self) -> EosResult<Diffusivity>
Return the diffusion reference as used in entropy scaling.
Sourcepub fn thermal_conductivity(&self) -> EosResult<ThermalConductivity>
pub fn thermal_conductivity(&self) -> EosResult<ThermalConductivity>
Return the thermal conductivity via entropy scaling.
Sourcepub fn ln_thermal_conductivity_reduced(&self) -> EosResult<f64>
pub fn ln_thermal_conductivity_reduced(&self) -> EosResult<f64>
Return the logarithm of the reduced thermal conductivity.
This term equals the thermal conductivity correlation function that is used for entropy scaling.
Sourcepub fn thermal_conductivity_reference(&self) -> EosResult<ThermalConductivity>
pub fn thermal_conductivity_reference(&self) -> EosResult<ThermalConductivity>
Return the thermal conductivity reference as used in entropy scaling.
Source§impl<R: Residual> State<R>
§Critical points
impl<R: Residual> State<R>
§Critical points
Sourcepub fn critical_point_pure(
eos: &Arc<R>,
initial_temperature: Option<Temperature>,
options: SolverOptions,
) -> EosResult<Vec<Self>>
pub fn critical_point_pure( eos: &Arc<R>, initial_temperature: Option<Temperature>, options: SolverOptions, ) -> EosResult<Vec<Self>>
Calculate the pure component critical point of all components.
pub fn critical_point_binary<TP: TemperatureOrPressure>( eos: &Arc<R>, temperature_or_pressure: TP, initial_temperature: Option<Temperature>, initial_molefracs: Option<[f64; 2]>, options: SolverOptions, ) -> EosResult<Self>
Sourcepub fn critical_point(
eos: &Arc<R>,
moles: Option<&Moles<Array1<f64>>>,
initial_temperature: Option<Temperature>,
options: SolverOptions,
) -> EosResult<Self>
pub fn critical_point( eos: &Arc<R>, moles: Option<&Moles<Array1<f64>>>, initial_temperature: Option<Temperature>, options: SolverOptions, ) -> EosResult<Self>
Calculate the critical point of a system for given moles.
pub fn spinodal( eos: &Arc<R>, temperature: Temperature, moles: Option<&Moles<Array1<f64>>>, options: SolverOptions, ) -> EosResult<[Self; 2]>
Source§impl<E: Residual> State<E>
§State constructors
impl<E: Residual> State<E>
§State constructors
Sourcepub fn new_nvt(
eos: &Arc<E>,
temperature: Temperature,
volume: Volume,
moles: &Moles<Array1<f64>>,
) -> EosResult<Self>
pub fn new_nvt( eos: &Arc<E>, temperature: Temperature, volume: Volume, moles: &Moles<Array1<f64>>, ) -> EosResult<Self>
Return a new State
given a temperature, an array of mole numbers and a volume.
This function will perform a validation of the given properties, i.e. test for signs and if values are finite. It will not validate physics, i.e. if the resulting densities are below the maximum packing fraction.
Sourcepub fn new_pure(
eos: &Arc<E>,
temperature: Temperature,
density: Density,
) -> EosResult<Self>
pub fn new_pure( eos: &Arc<E>, temperature: Temperature, density: Density, ) -> EosResult<Self>
Return a new State
for a pure component given a temperature and a density. The moles
are set to the reference value for each component.
This function will perform a validation of the given properties, i.e. test for signs and if values are finite. It will not validate physics, i.e. if the resulting densities are below the maximum packing fraction.
Sourcepub fn new(
eos: &Arc<E>,
temperature: Option<Temperature>,
volume: Option<Volume>,
density: Option<Density>,
partial_density: Option<&Density<Array1<f64>>>,
total_moles: Option<Moles>,
moles: Option<&Moles<Array1<f64>>>,
molefracs: Option<&Array1<f64>>,
pressure: Option<Pressure>,
density_initialization: DensityInitialization,
) -> EosResult<Self>
pub fn new( eos: &Arc<E>, temperature: Option<Temperature>, volume: Option<Volume>, density: Option<Density>, partial_density: Option<&Density<Array1<f64>>>, total_moles: Option<Moles>, moles: Option<&Moles<Array1<f64>>>, molefracs: Option<&Array1<f64>>, pressure: Option<Pressure>, density_initialization: DensityInitialization, ) -> EosResult<Self>
Return a new State
for the combination of inputs.
The function attempts to create a new state using the given input values. If the state is overdetermined, it will choose a method based on the following hierarchy.
- Create a state non-iteratively from the set of $T$, $V$, $\rho$, $\rho_i$, $N$, $N_i$ and $x_i$.
- Use a density iteration for a given pressure.
The StateBuilder provides a convenient way of calling this function without the need to provide all the optional input values.
§Errors
When the state cannot be created using the combination of inputs.
Sourcepub fn new_npt(
eos: &Arc<E>,
temperature: Temperature,
pressure: Pressure,
moles: &Moles<Array1<f64>>,
density_initialization: DensityInitialization,
) -> EosResult<Self>
pub fn new_npt( eos: &Arc<E>, temperature: Temperature, pressure: Pressure, moles: &Moles<Array1<f64>>, density_initialization: DensityInitialization, ) -> EosResult<Self>
Return a new State
using a density iteration. DensityInitialization is used to
influence the calculation with respect to the possible solutions.
Sourcepub fn new_npvx(
eos: &Arc<E>,
temperature: Temperature,
pressure: Pressure,
volume: Volume,
molefracs: &Array1<f64>,
density_initialization: DensityInitialization,
) -> EosResult<Self>
pub fn new_npvx( eos: &Arc<E>, temperature: Temperature, pressure: Pressure, volume: Volume, molefracs: &Array1<f64>, density_initialization: DensityInitialization, ) -> EosResult<Self>
Return a new State
for given pressure $p$, volume $V$, temperature $T$ and composition $x_i$.
Source§impl<E: Residual + IdealGas> State<E>
impl<E: Residual + IdealGas> State<E>
Sourcepub fn new_full(
eos: &Arc<E>,
temperature: Option<Temperature>,
volume: Option<Volume>,
density: Option<Density>,
partial_density: Option<&Density<Array1<f64>>>,
total_moles: Option<Moles>,
moles: Option<&Moles<Array1<f64>>>,
molefracs: Option<&Array1<f64>>,
pressure: Option<Pressure>,
molar_enthalpy: Option<MolarEnergy>,
molar_entropy: Option<MolarEntropy>,
molar_internal_energy: Option<MolarEnergy>,
density_initialization: DensityInitialization,
initial_temperature: Option<Temperature>,
) -> EosResult<Self>
pub fn new_full( eos: &Arc<E>, temperature: Option<Temperature>, volume: Option<Volume>, density: Option<Density>, partial_density: Option<&Density<Array1<f64>>>, total_moles: Option<Moles>, moles: Option<&Moles<Array1<f64>>>, molefracs: Option<&Array1<f64>>, pressure: Option<Pressure>, molar_enthalpy: Option<MolarEnergy>, molar_entropy: Option<MolarEntropy>, molar_internal_energy: Option<MolarEnergy>, density_initialization: DensityInitialization, initial_temperature: Option<Temperature>, ) -> EosResult<Self>
Return a new State
for the combination of inputs.
The function attempts to create a new state using the given input values. If the state is overdetermined, it will choose a method based on the following hierarchy.
- Create a state non-iteratively from the set of $T$, $V$, $\rho$, $\rho_i$, $N$, $N_i$ and $x_i$.
- Use a density iteration for a given pressure.
- Determine the state using a Newton iteration from (in this order): $(p, h)$, $(p, s)$, $(T, h)$, $(T, s)$, $(V, u)$
The StateBuilder provides a convenient way of calling this function without the need to provide all the optional input values.
§Errors
When the state cannot be created using the combination of inputs.
Sourcepub fn new_nph(
eos: &Arc<E>,
pressure: Pressure,
molar_enthalpy: MolarEnergy,
moles: &Moles<Array1<f64>>,
density_initialization: DensityInitialization,
initial_temperature: Option<Temperature>,
) -> EosResult<Self>
pub fn new_nph( eos: &Arc<E>, pressure: Pressure, molar_enthalpy: MolarEnergy, moles: &Moles<Array1<f64>>, density_initialization: DensityInitialization, initial_temperature: Option<Temperature>, ) -> EosResult<Self>
Return a new State
for given pressure $p$ and molar enthalpy $h$.
Sourcepub fn new_nth(
eos: &Arc<E>,
temperature: Temperature,
molar_enthalpy: MolarEnergy,
moles: &Moles<Array1<f64>>,
density_initialization: DensityInitialization,
) -> EosResult<Self>
pub fn new_nth( eos: &Arc<E>, temperature: Temperature, molar_enthalpy: MolarEnergy, moles: &Moles<Array1<f64>>, density_initialization: DensityInitialization, ) -> EosResult<Self>
Return a new State
for given temperature $T$ and molar enthalpy $h$.
Sourcepub fn new_nts(
eos: &Arc<E>,
temperature: Temperature,
molar_entropy: MolarEntropy,
moles: &Moles<Array1<f64>>,
density_initialization: DensityInitialization,
) -> EosResult<Self>
pub fn new_nts( eos: &Arc<E>, temperature: Temperature, molar_entropy: MolarEntropy, moles: &Moles<Array1<f64>>, density_initialization: DensityInitialization, ) -> EosResult<Self>
Return a new State
for given temperature $T$ and molar entropy $s$.
Sourcepub fn new_nps(
eos: &Arc<E>,
pressure: Pressure,
molar_entropy: MolarEntropy,
moles: &Moles<Array1<f64>>,
density_initialization: DensityInitialization,
initial_temperature: Option<Temperature>,
) -> EosResult<Self>
pub fn new_nps( eos: &Arc<E>, pressure: Pressure, molar_entropy: MolarEntropy, moles: &Moles<Array1<f64>>, density_initialization: DensityInitialization, initial_temperature: Option<Temperature>, ) -> EosResult<Self>
Return a new State
for given pressure $p$ and molar entropy $s$.
Source§impl<E: Residual> State<E>
impl<E: Residual> State<E>
Sourcepub fn update_temperature(&self, temperature: Temperature) -> EosResult<Self>
pub fn update_temperature(&self, temperature: Temperature) -> EosResult<Self>
Update the state with the given temperature
Sourcepub fn derive1(&self, derivative: Derivative) -> StateHD<Dual64>
pub fn derive1(&self, derivative: Derivative) -> StateHD<Dual64>
Creates a StateHD taking the first derivative.
Sourcepub fn derive2(&self, derivative: Derivative) -> StateHD<Dual2_64>
pub fn derive2(&self, derivative: Derivative) -> StateHD<Dual2_64>
Creates a StateHD taking the first and second (partial) derivatives.
Sourcepub fn derive2_mixed(
&self,
derivative1: Derivative,
derivative2: Derivative,
) -> StateHD<HyperDual64>
pub fn derive2_mixed( &self, derivative1: Derivative, derivative2: Derivative, ) -> StateHD<HyperDual64>
Creates a StateHD taking the first and second (partial) derivatives.
Trait Implementations§
Auto Trait Implementations§
impl<E> !Freeze for State<E>
impl<E> RefUnwindSafe for State<E>where
E: RefUnwindSafe,
impl<E> Send for State<E>
impl<E> Sync for State<E>
impl<E> Unpin for State<E>
impl<E> UnwindSafe for State<E>where
E: RefUnwindSafe,
Blanket Implementations§
Source§impl<Src, Scheme> ApproxFrom<Src, Scheme> for Srcwhere
Scheme: ApproxScheme,
impl<Src, Scheme> ApproxFrom<Src, Scheme> for Srcwhere
Scheme: ApproxScheme,
Source§fn approx_from(src: Src) -> Result<Src, <Src as ApproxFrom<Src, Scheme>>::Err>
fn approx_from(src: Src) -> Result<Src, <Src as ApproxFrom<Src, Scheme>>::Err>
Source§impl<Dst, Src, Scheme> ApproxInto<Dst, Scheme> for Srcwhere
Dst: ApproxFrom<Src, Scheme>,
Scheme: ApproxScheme,
impl<Dst, Src, Scheme> ApproxInto<Dst, Scheme> for Srcwhere
Dst: ApproxFrom<Src, Scheme>,
Scheme: ApproxScheme,
Source§type Err = <Dst as ApproxFrom<Src, Scheme>>::Err
type Err = <Dst as ApproxFrom<Src, Scheme>>::Err
Source§fn approx_into(self) -> Result<Dst, <Src as ApproxInto<Dst, Scheme>>::Err>
fn approx_into(self) -> Result<Dst, <Src as ApproxInto<Dst, Scheme>>::Err>
Source§impl<T> BorrowMut<T> for Twhere
T: ?Sized,
impl<T> BorrowMut<T> for Twhere
T: ?Sized,
Source§fn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
Source§impl<T> CloneToUninit for Twhere
T: Clone,
impl<T> CloneToUninit for Twhere
T: Clone,
Source§impl<T, Dst> ConvAsUtil<Dst> for T
impl<T, Dst> ConvAsUtil<Dst> for T
Source§impl<T> ConvUtil for T
impl<T> ConvUtil for T
Source§fn approx_as<Dst>(self) -> Result<Dst, Self::Err>where
Self: Sized + ApproxInto<Dst>,
fn approx_as<Dst>(self) -> Result<Dst, Self::Err>where
Self: Sized + ApproxInto<Dst>,
Source§fn approx_as_by<Dst, Scheme>(self) -> Result<Dst, Self::Err>
fn approx_as_by<Dst, Scheme>(self) -> Result<Dst, Self::Err>
Source§impl<T> IntoEither for T
impl<T> IntoEither for T
Source§fn into_either(self, into_left: bool) -> Either<Self, Self>
fn into_either(self, into_left: bool) -> Either<Self, Self>
self
into a Left
variant of Either<Self, Self>
if into_left
is true
.
Converts self
into a Right
variant of Either<Self, Self>
otherwise. Read moreSource§fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
self
into a Left
variant of Either<Self, Self>
if into_left(&self)
returns true
.
Converts self
into a Right
variant of Either<Self, Self>
otherwise. Read moreSource§impl<T> Pointable for T
impl<T> Pointable for T
Source§impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
Source§fn to_subset(&self) -> Option<SS>
fn to_subset(&self) -> Option<SS>
self
from the equivalent element of its
superset. Read moreSource§fn is_in_subset(&self) -> bool
fn is_in_subset(&self) -> bool
self
is actually part of its subset T
(and can be converted to it).Source§fn to_subset_unchecked(&self) -> SS
fn to_subset_unchecked(&self) -> SS
self.to_subset
but without any property checks. Always succeeds.Source§fn from_subset(element: &SS) -> SP
fn from_subset(element: &SS) -> SP
self
to the equivalent element of its superset.