pub struct State<E, N: Dim = Dyn, D: DualNum<f64> + Copy = f64>where
DefaultAllocator: Allocator<N>,{
pub eos: E,
pub temperature: Temperature<D>,
pub volume: Volume<D>,
pub moles: Moles<OVector<D, N>>,
pub total_moles: Moles<D>,
pub partial_density: Density<OVector<D, N>>,
pub density: Density<D>,
pub molefracs: OVector<D, N>,
/* 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: EEquation of state
temperature: Temperature<D>Temperature $T$
volume: Volume<D>Volume $V$
moles: Moles<OVector<D, N>>Mole numbers $N_i$
total_moles: Moles<D>Total number of moles $N=\sum_iN_i$
partial_density: Density<OVector<D, N>>Partial densities $\rho_i=\frac{N_i}{V}$
density: Density<D>Total density $\rho=\frac{N}{V}=\sum_i\rho_i$
molefracs: OVector<D, N>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) -> FeosResult<bool>
pub fn is_stable(&self, options: SolverOptions) -> FeosResult<bool>
Determine if the state is stable, i.e. if a phase split should occur or not.
Sourcepub fn stability_analysis(
&self,
options: SolverOptions,
) -> FeosResult<Vec<State<E>>>
pub fn stability_analysis( &self, options: SolverOptions, ) -> FeosResult<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>>,
) -> FeosResult<PhaseEquilibrium<E, 2>>
pub fn tp_flash( &self, initial_state: Option<&PhaseEquilibrium<E, 2>>, options: SolverOptions, non_volatile_components: Option<Vec<usize>>, ) -> FeosResult<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>>, ) -> FeosResult<PhaseEquilibrium<E, 2>>
Source§impl<E: Total<N, D>, N: Gradients, D: DualNum<f64> + Copy> State<E, N, D>where
DefaultAllocator: Allocator<N>,
impl<E: Total<N, D>, N: Gradients, D: DualNum<f64> + Copy> State<E, N, D>where
DefaultAllocator: Allocator<N>,
Sourcepub fn chemical_potential(
&self,
contributions: Contributions,
) -> MolarEnergy<OVector<D, N>>
pub fn chemical_potential( &self, contributions: Contributions, ) -> MolarEnergy<OVector<D, N>>
Chemical potential: $\mu_i=\left(\frac{\partial A}{\partial N_i}\right)_{T,V,N_j}$
Sourcepub fn dmu_dt(
&self,
contributions: Contributions,
) -> MolarEntropy<OVector<D, N>>
pub fn dmu_dt( &self, contributions: Contributions, ) -> MolarEntropy<OVector<D, N>>
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<D>
pub fn molar_isochoric_heat_capacity( &self, contributions: Contributions, ) -> MolarEntropy<D>
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<D> as Div<Temperature<D>>>::Output
pub fn dc_v_dt( &self, contributions: Contributions, ) -> <MolarEntropy<D> as Div<Temperature<D>>>::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<D>
pub fn molar_isobaric_heat_capacity( &self, contributions: Contributions, ) -> MolarEntropy<D>
Molar isobaric heat capacity: $c_p=\left(\frac{\partial h}{\partial T}\right)_{p,N_i}$
Sourcepub fn entropy(&self, contributions: Contributions) -> Entropy<D>
pub fn entropy(&self, contributions: Contributions) -> Entropy<D>
Entropy: $S=-\left(\frac{\partial A}{\partial T}\right)_{V,N_i}$
Sourcepub fn molar_entropy(&self, contributions: Contributions) -> MolarEntropy<D>
pub fn molar_entropy(&self, contributions: Contributions) -> MolarEntropy<D>
Molar entropy: $s=\frac{S}{N}$
Sourcepub fn partial_molar_entropy(&self) -> MolarEntropy<OVector<D, N>>
pub fn partial_molar_entropy(&self) -> MolarEntropy<OVector<D, N>>
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<D> as Div<Temperature<D>>>::Output
pub fn ds_dt( &self, contributions: Contributions, ) -> <Entropy<D> as Div<Temperature<D>>>::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<D> as Div<Temperature<D>>>::Output as Div<Temperature<D>>>::Output
pub fn d2s_dt2( &self, contributions: Contributions, ) -> <<Entropy<D> as Div<Temperature<D>>>::Output as Div<Temperature<D>>>::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<D>
pub fn enthalpy(&self, contributions: Contributions) -> Energy<D>
Enthalpy: $H=A+TS+pV$
Sourcepub fn molar_enthalpy(&self, contributions: Contributions) -> MolarEnergy<D>
pub fn molar_enthalpy(&self, contributions: Contributions) -> MolarEnergy<D>
Molar enthalpy: $h=\frac{H}{N}$
Sourcepub fn partial_molar_enthalpy(&self) -> MolarEnergy<OVector<D, N>>
pub fn partial_molar_enthalpy(&self) -> MolarEnergy<OVector<D, N>>
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<D>
pub fn helmholtz_energy(&self, contributions: Contributions) -> Energy<D>
Helmholtz energy: $A$
Sourcepub fn molar_helmholtz_energy(
&self,
contributions: Contributions,
) -> MolarEnergy<D>
pub fn molar_helmholtz_energy( &self, contributions: Contributions, ) -> MolarEnergy<D>
Molar Helmholtz energy: $a=\frac{A}{N}$
Sourcepub fn internal_energy(&self, contributions: Contributions) -> Energy<D>
pub fn internal_energy(&self, contributions: Contributions) -> Energy<D>
Internal energy: $U=A+TS$
Sourcepub fn molar_internal_energy(
&self,
contributions: Contributions,
) -> MolarEnergy<D>
pub fn molar_internal_energy( &self, contributions: Contributions, ) -> MolarEnergy<D>
Molar internal energy: $u=\frac{U}{N}$
Sourcepub fn gibbs_energy(&self, contributions: Contributions) -> Energy<D>
pub fn gibbs_energy(&self, contributions: Contributions) -> Energy<D>
Gibbs energy: $G=A+pV$
Sourcepub fn molar_gibbs_energy(&self, contributions: Contributions) -> MolarEnergy<D>
pub fn molar_gibbs_energy(&self, contributions: Contributions) -> MolarEnergy<D>
Molar Gibbs energy: $g=\frac{G}{N}$
Sourcepub fn joule_thomson(&self) -> <Temperature<D> as Div<Pressure<D>>>::Output
pub fn joule_thomson(&self) -> <Temperature<D> as Div<Pressure<D>>>::Output
Joule Thomson coefficient: $\mu_{JT}=\left(\frac{\partial T}{\partial p}\right)_{H,N_i}$
Sourcepub fn isentropic_compressibility(
&self,
) -> Quantity<D, <_Pressure as Neg>::Output>
pub fn isentropic_compressibility( &self, ) -> Quantity<D, <_Pressure as Neg>::Output>
Isentropic compressibility: $\kappa_s=-\frac{1}{V}\left(\frac{\partial V}{\partial p}\right)_{S,N_i}$
Sourcepub fn isenthalpic_compressibility(
&self,
) -> Quantity<D, <_Pressure as Neg>::Output>
pub fn isenthalpic_compressibility( &self, ) -> Quantity<D, <_Pressure as Neg>::Output>
Isenthalpic compressibility: $\kappa_H=-\frac{1}{V}\left(\frac{\partial V}{\partial p}\right)_{H,N_i}$
Sourcepub fn thermal_expansivity(&self) -> Quantity<D, <_Temperature as Neg>::Output>
pub fn thermal_expansivity(&self) -> Quantity<D, <_Temperature as Neg>::Output>
Thermal expansivity: $\alpha_p=-\frac{1}{V}\left(\frac{\partial V}{\partial T}\right)_{p,N_i}$
Sourcepub fn grueneisen_parameter(&self) -> D
pub fn grueneisen_parameter(&self) -> D
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<(&'static str, MolarEnergy<D>)>
pub fn chemical_potential_contributions( &self, component: usize, contributions: Contributions, ) -> Vec<(&'static str, MolarEnergy<D>)>
Chemical potential $\mu_i$ evaluated for each contribution of the equation of state.
Source§impl<E: Total<N, D> + Molarweight<N, D>, N: Gradients, D: DualNum<f64> + Copy> State<E, N, D>where
DefaultAllocator: Allocator<N>,
impl<E: Total<N, D> + Molarweight<N, D>, N: Gradients, D: DualNum<f64> + Copy> State<E, N, D>where
DefaultAllocator: Allocator<N>,
Sourcepub fn specific_isochoric_heat_capacity(
&self,
contributions: Contributions,
) -> SpecificEntropy<D>
pub fn specific_isochoric_heat_capacity( &self, contributions: Contributions, ) -> SpecificEntropy<D>
Specific isochoric heat capacity: $c_v^{(m)}=\frac{C_v}{m}$
Sourcepub fn specific_isobaric_heat_capacity(
&self,
contributions: Contributions,
) -> SpecificEntropy<D>
pub fn specific_isobaric_heat_capacity( &self, contributions: Contributions, ) -> SpecificEntropy<D>
Specific isobaric heat capacity: $c_p^{(m)}=\frac{C_p}{m}$
Sourcepub fn specific_entropy(
&self,
contributions: Contributions,
) -> SpecificEntropy<D>
pub fn specific_entropy( &self, contributions: Contributions, ) -> SpecificEntropy<D>
Specific entropy: $s^{(m)}=\frac{S}{m}$
Sourcepub fn specific_enthalpy(
&self,
contributions: Contributions,
) -> SpecificEnergy<D>
pub fn specific_enthalpy( &self, contributions: Contributions, ) -> SpecificEnergy<D>
Specific enthalpy: $h^{(m)}=\frac{H}{m}$
Sourcepub fn specific_helmholtz_energy(
&self,
contributions: Contributions,
) -> SpecificEnergy<D>
pub fn specific_helmholtz_energy( &self, contributions: Contributions, ) -> SpecificEnergy<D>
Specific Helmholtz energy: $a^{(m)}=\frac{A}{m}$
Sourcepub fn specific_internal_energy(
&self,
contributions: Contributions,
) -> SpecificEnergy<D>
pub fn specific_internal_energy( &self, contributions: Contributions, ) -> SpecificEnergy<D>
Specific internal energy: $u^{(m)}=\frac{U}{m}$
Sourcepub fn specific_gibbs_energy(
&self,
contributions: Contributions,
) -> SpecificEnergy<D>
pub fn specific_gibbs_energy( &self, contributions: Contributions, ) -> SpecificEnergy<D>
Specific Gibbs energy: $g^{(m)}=\frac{G}{m}$
Source§impl<E: Total<N> + Molarweight<N>, N: Gradients> State<E, N>where
DefaultAllocator: Allocator<N>,
impl<E: Total<N> + Molarweight<N>, N: Gradients> State<E, N>where
DefaultAllocator: Allocator<N>,
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<N, D>, N: Gradients, D: DualNum<f64> + Copy> State<E, N, D>where
DefaultAllocator: Allocator<N>,
§State properties
impl<E: Residual<N, D>, N: Gradients, D: DualNum<f64> + Copy> State<E, N, D>where
DefaultAllocator: Allocator<N>,
§State properties
Sourcepub fn residual_helmholtz_energy(&self) -> Energy<D>
pub fn residual_helmholtz_energy(&self) -> Energy<D>
Residual Helmholtz energy $A^\text{res}$
Sourcepub fn residual_molar_helmholtz_energy(&self) -> MolarEnergy<D>
pub fn residual_molar_helmholtz_energy(&self) -> MolarEnergy<D>
Residual molar Helmholtz energy $a^\text{res}$
Sourcepub fn residual_entropy(&self) -> Entropy<D>
pub fn residual_entropy(&self) -> Entropy<D>
Residual entropy $S^\text{res}=\left(\frac{\partial A^\text{res}}{\partial T}\right)_{V,N_i}$
Sourcepub fn residual_molar_entropy(&self) -> MolarEntropy<D>
pub fn residual_molar_entropy(&self) -> MolarEntropy<D>
Residual entropy $s^\text{res}=\left(\frac{\partial a^\text{res}}{\partial T}\right)_{V,N_i}$
Sourcepub fn pressure(&self, contributions: Contributions) -> Pressure<D>
pub fn pressure(&self, contributions: Contributions) -> Pressure<D>
Pressure: $p=-\left(\frac{\partial A}{\partial V}\right)_{T,N_i}$
Sourcepub fn residual_chemical_potential(&self) -> MolarEnergy<OVector<D, N>>
pub fn residual_chemical_potential(&self) -> MolarEnergy<OVector<D, N>>
Residual chemical potential: $\mu_i^\text{res}=\left(\frac{\partial A^\text{res}}{\partial N_i}\right)_{T,V,N_j}$
Sourcepub fn compressibility(&self, contributions: Contributions) -> D
pub fn compressibility(&self, contributions: Contributions) -> D
Compressibility factor: $Z=\frac{pV}{NRT}$
Sourcepub fn dp_dv(
&self,
contributions: Contributions,
) -> <Pressure<D> as Div<Volume<D>>>::Output
pub fn dp_dv( &self, contributions: Contributions, ) -> <Pressure<D> as Div<Volume<D>>>::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<D> as Div<Density<D>>>::Output
pub fn dp_drho( &self, contributions: Contributions, ) -> <Pressure<D> as Div<Density<D>>>::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,
) -> Quantity<D, <_Pressure as Sub<_Temperature>>::Output>
pub fn dp_dt( &self, contributions: Contributions, ) -> Quantity<D, <_Pressure as Sub<_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,
) -> Quantity<OVector<D, N>, <_Pressure as Sub<_Moles>>::Output>
pub fn dp_dni( &self, contributions: Contributions, ) -> Quantity<OVector<D, N>, <_Pressure as Sub<_Moles>>::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<D> as Div<Volume<D>>>::Output as Div<Volume<D>>>::Output
pub fn d2p_dv2( &self, contributions: Contributions, ) -> <<Pressure<D> as Div<Volume<D>>>::Output as Div<Volume<D>>>::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<D> as Div<Density<D>>>::Output as Div<Density<D>>>::Output
pub fn d2p_drho2( &self, contributions: Contributions, ) -> <<Pressure<D> as Div<Density<D>>>::Output as Div<Density<D>>>::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) -> D
pub fn structure_factor(&self) -> D
Structure factor: $S(0)=k_BT\left(\frac{\partial\rho}{\partial p}\right)_{T,N_i}$
Sourcepub fn partial_molar_volume(&self) -> MolarVolume<OVector<D, N>>
pub fn partial_molar_volume(&self) -> MolarVolume<OVector<D, N>>
Partial molar volume: $v_i=\left(\frac{\partial V}{\partial N_i}\right)_{T,p,N_j}$
Sourcepub fn dmu_dni(
&self,
contributions: Contributions,
) -> Quantity<OMatrix<D, N, N>, <_MolarEnergy as Sub<_Moles>>::Output>where
DefaultAllocator: Allocator<N, N>,
pub fn dmu_dni(
&self,
contributions: Contributions,
) -> Quantity<OMatrix<D, N, N>, <_MolarEnergy as Sub<_Moles>>::Output>where
DefaultAllocator: Allocator<N, N>,
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,
) -> Quantity<D, <_Pressure as Neg>::Output>
pub fn isothermal_compressibility( &self, ) -> Quantity<D, <_Pressure as Neg>::Output>
Isothermal compressibility: $\kappa_T=-\frac{1}{V}\left(\frac{\partial V}{\partial p}\right)_{T,N_i}$
Sourcepub fn ds_res_dt(&self) -> <Entropy<D> as Div<Temperature<D>>>::Output
pub fn ds_res_dt(&self) -> <Entropy<D> as Div<Temperature<D>>>::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<D> as Div<Temperature<D>>>::Output as Div<Temperature<D>>>::Output
pub fn d2s_res_dt2( &self, ) -> <<Entropy<D> as Div<Temperature<D>>>::Output as Div<Temperature<D>>>::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) -> MolarEntropy<OVector<D, N>>
pub fn dmu_res_dt(&self) -> MolarEntropy<OVector<D, N>>
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) -> OVector<D, N>
pub fn ln_phi(&self) -> OVector<D, N>
Logarithm of the fugacity coefficient: $\ln\varphi_i=\beta\mu_i^\mathrm{res}\left(T,p,\lbrace N_i\rbrace\right)$
Sourcepub fn dln_phi_dt(
&self,
) -> Quantity<OVector<D, N>, <_Temperature as Neg>::Output>
pub fn dln_phi_dt( &self, ) -> Quantity<OVector<D, N>, <_Temperature as Neg>::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) -> Quantity<OVector<D, N>, <_Pressure as Neg>::Output>
pub fn dln_phi_dp(&self) -> Quantity<OVector<D, N>, <_Pressure as Neg>::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) -> Quantity<OMatrix<D, N, N>, <_Moles as Neg>::Output>where
DefaultAllocator: Allocator<N, N>,
pub fn dln_phi_dnj(&self) -> Quantity<OMatrix<D, N, N>, <_Moles as Neg>::Output>where
DefaultAllocator: Allocator<N, N>,
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}$
Source§impl<E: Residual + Subset> State<E>
impl<E: Residual + Subset> State<E>
Sourcepub fn ln_phi_pure_liquid(&self) -> FeosResult<DVector<f64>>
pub fn ln_phi_pure_liquid(&self) -> FeosResult<DVector<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) -> FeosResult<DVector<f64>>
pub fn ln_symmetric_activity_coefficient(&self) -> FeosResult<DVector<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: &E,
temperature: Temperature,
molefracs: &DVector<f64>,
) -> FeosResult<Vec<Pressure>>
pub fn henrys_law_constant( eos: &E, temperature: Temperature, molefracs: &DVector<f64>, ) -> FeosResult<Vec<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}}}$
The composition of the (possibly mixed) solvent is determined by the molefracs. All components for which the composition is 0 are treated as solutes.
For some reason the compiler is overwhelmed if returning a quantity array, therefore it is returned as list.
Sourcepub fn henrys_law_constant_binary(
eos: &E,
temperature: Temperature,
) -> FeosResult<Pressure>
pub fn henrys_law_constant_binary( eos: &E, temperature: Temperature, ) -> FeosResult<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.
Source§impl<E: Residual> State<E>
impl<E: Residual> State<E>
Sourcepub fn thermodynamic_factor(&self) -> DMatrix<f64>
pub fn thermodynamic_factor(&self) -> DMatrix<f64>
Thermodynamic factor: $\Gamma_{ij}=\delta_{ij}+x_i\left(\frac{\partial\ln\varphi_i}{\partial x_j}\right)_{T,p,\Sigma}$
Source§impl<E: Residual<N, D>, N: Gradients, D: DualNum<f64> + Copy> State<E, N, D>where
DefaultAllocator: Allocator<N>,
impl<E: Residual<N, D>, N: Gradients, D: DualNum<f64> + Copy> State<E, N, D>where
DefaultAllocator: Allocator<N>,
Sourcepub fn residual_molar_isochoric_heat_capacity(&self) -> MolarEntropy<D>
pub fn residual_molar_isochoric_heat_capacity(&self) -> MolarEntropy<D>
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<D> as Div<Temperature<D>>>::Output
pub fn dc_v_res_dt(&self) -> <MolarEntropy<D> as Div<Temperature<D>>>::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<D>
pub fn residual_molar_isobaric_heat_capacity(&self) -> MolarEntropy<D>
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<D>
pub fn residual_enthalpy(&self) -> Energy<D>
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<D>
pub fn residual_molar_enthalpy(&self) -> MolarEnergy<D>
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<D>
pub fn residual_internal_energy(&self) -> Energy<D>
Residual internal energy: $U^\text{res}(T,V,\mathbf{n})=A^\text{res}+TS^\text{res}$
Sourcepub fn residual_molar_internal_energy(&self) -> MolarEnergy<D>
pub fn residual_molar_internal_energy(&self) -> MolarEnergy<D>
Residual molar internal energy: $u^\text{res}(T,V,\mathbf{n})=a^\text{res}+Ts^\text{res}$
Sourcepub fn residual_gibbs_energy(&self) -> Energy<D>
pub fn residual_gibbs_energy(&self) -> Energy<D>
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<D>
pub fn residual_molar_gibbs_energy(&self) -> MolarEnergy<D>
Residual Gibbs energy: $g^\text{res}(T,p,\mathbf{n})=a^\text{res}+p^\text{res}v-RT \ln Z$
Sourcepub fn residual_molar_helmholtz_energy_contributions(
&self,
) -> Vec<(&'static str, MolarEnergy<D>)>
pub fn residual_molar_helmholtz_energy_contributions( &self, ) -> Vec<(&'static str, MolarEnergy<D>)>
Molar Helmholtz energy $a^\text{res}$ evaluated for each residual contribution of the equation of state.
Sourcepub fn residual_chemical_potential_contributions(
&self,
component: usize,
) -> Vec<(&'static str, MolarEnergy<D>)>
pub fn residual_chemical_potential_contributions( &self, component: usize, ) -> Vec<(&'static str, MolarEnergy<D>)>
Chemical potential $\mu_i^\text{res}$ evaluated for each residual contribution of the equation of state.
Sourcepub fn pressure_contributions(&self) -> Vec<(&'static str, Pressure<D>)>
pub fn pressure_contributions(&self) -> Vec<(&'static str, Pressure<D>)>
Pressure $p$ evaluated for each contribution of the equation of state.
Source§impl<E: Residual<N, D> + Molarweight<N, D>, N: Gradients, D: DualNum<f64> + Copy> State<E, N, D>where
DefaultAllocator: Allocator<N>,
impl<E: Residual<N, D> + Molarweight<N, D>, N: Gradients, D: DualNum<f64> + Copy> State<E, N, D>where
DefaultAllocator: Allocator<N>,
Sourcepub fn total_molar_weight(&self) -> MolarWeight<D>
pub fn total_molar_weight(&self) -> MolarWeight<D>
Total molar weight: $MW=\sum_ix_iMW_i$
Sourcepub fn total_mass(&self) -> Mass<D>
pub fn total_mass(&self) -> Mass<D>
Total mass: $m=\sum_im_i=nMW$
Sourcepub fn mass_density(&self) -> MassDensity<D>
pub fn mass_density(&self) -> MassDensity<D>
Mass density: $\rho^{(m)}=\frac{m}{V}$
Source§impl<E: Residual<N, D> + EntropyScaling<N, D>, N: Gradients, D: DualNum<f64> + Copy> State<E, N, D>where
DefaultAllocator: Allocator<N>,
§Transport properties
These properties are available for equations of state
that implement the EntropyScaling trait.
impl<E: Residual<N, D> + EntropyScaling<N, D>, N: Gradients, D: DualNum<f64> + Copy> State<E, N, D>where
DefaultAllocator: Allocator<N>,
§Transport properties
These properties are available for equations of state that implement the EntropyScaling trait.
Sourcepub fn ln_viscosity_reduced(&self) -> D
pub fn ln_viscosity_reduced(&self) -> D
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) -> Viscosity<D>
pub fn viscosity_reference(&self) -> Viscosity<D>
Return the viscosity reference as used in entropy scaling.
Sourcepub fn diffusion(&self) -> Diffusivity<D>
pub fn diffusion(&self) -> Diffusivity<D>
Return the diffusion via entropy scaling.
Sourcepub fn ln_diffusion_reduced(&self) -> D
pub fn ln_diffusion_reduced(&self) -> D
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) -> Diffusivity<D>
pub fn diffusion_reference(&self) -> Diffusivity<D>
Return the diffusion reference as used in entropy scaling.
Sourcepub fn thermal_conductivity(&self) -> ThermalConductivity<D>
pub fn thermal_conductivity(&self) -> ThermalConductivity<D>
Return the thermal conductivity via entropy scaling.
Sourcepub fn ln_thermal_conductivity_reduced(&self) -> D
pub fn ln_thermal_conductivity_reduced(&self) -> D
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) -> ThermalConductivity<D>
pub fn thermal_conductivity_reference(&self) -> ThermalConductivity<D>
Return the thermal conductivity reference as used in entropy scaling.
Source§impl<R: Residual + Subset> State<R>
§Critical points
impl<R: Residual + Subset> State<R>
§Critical points
Sourcepub fn critical_point_pure(
eos: &R,
initial_temperature: Option<Temperature>,
initial_density: Option<Density>,
options: SolverOptions,
) -> FeosResult<Vec<Self>>
pub fn critical_point_pure( eos: &R, initial_temperature: Option<Temperature>, initial_density: Option<Density>, options: SolverOptions, ) -> FeosResult<Vec<Self>>
Calculate the pure component critical point of all components.
Source§impl<E: Residual<N, D>, N: Gradients, D: DualNum<f64> + Copy> State<E, N, D>
impl<E: Residual<N, D>, N: Gradients, D: DualNum<f64> + Copy> State<E, N, D>
pub fn critical_point_binary<TP: TemperatureOrPressure<D>>( eos: &E, temperature_or_pressure: TP, initial_temperature: Option<Temperature>, initial_molefracs: Option<[f64; 2]>, initial_density: Option<Density>, options: SolverOptions, ) -> FeosResult<Self>
Sourcepub fn critical_point(
eos: &E,
molefracs: Option<&OVector<D, N>>,
initial_temperature: Option<Temperature>,
initial_density: Option<Density>,
options: SolverOptions,
) -> FeosResult<Self>
pub fn critical_point( eos: &E, molefracs: Option<&OVector<D, N>>, initial_temperature: Option<Temperature>, initial_density: Option<Density>, options: SolverOptions, ) -> FeosResult<Self>
Calculate the critical point of a system for given moles.
Source§impl<E: Residual<N>, N: Gradients> State<E, N>
impl<E: Residual<N>, N: Gradients> State<E, N>
pub fn spinodal( eos: &E, temperature: Temperature, molefracs: Option<&OVector<f64, N>>, options: SolverOptions, ) -> FeosResult<[Self; 2]>
Source§impl<E: Residual<N, D>, N: Dim, D: DualNum<f64> + Copy> State<E, N, D>where
DefaultAllocator: Allocator<N>,
impl<E: Residual<N, D>, N: Dim, D: DualNum<f64> + Copy> State<E, N, D>where
DefaultAllocator: Allocator<N>,
Sourcepub fn new_nvt(
eos: &E,
temperature: Temperature<D>,
volume: Volume<D>,
moles: &Moles<OVector<D, N>>,
) -> FeosResult<Self>
pub fn new_nvt( eos: &E, temperature: Temperature<D>, volume: Volume<D>, moles: &Moles<OVector<D, N>>, ) -> FeosResult<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_intensive(
eos: &E,
temperature: Temperature<D>,
density: Density<D>,
molefracs: &OVector<D, N>,
) -> FeosResult<Self>
pub fn new_intensive( eos: &E, temperature: Temperature<D>, density: Density<D>, molefracs: &OVector<D, N>, ) -> FeosResult<Self>
Return a new State for which the total amount of substance is unspecified.
Internally the total number of moles will be set to 1 mol.
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: &E,
temperature: Temperature<D>,
density: Density<D>,
) -> FeosResult<Self>
pub fn new_pure( eos: &E, temperature: Temperature<D>, density: Density<D>, ) -> FeosResult<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: &E,
temperature: Option<Temperature<D>>,
volume: Option<Volume<D>>,
density: Option<Density<D>>,
partial_density: Option<&Density<OVector<D, N>>>,
total_moles: Option<Moles<D>>,
moles: Option<&Moles<OVector<D, N>>>,
molefracs: Option<&OVector<D, N>>,
pressure: Option<Pressure<D>>,
density_initialization: Option<DensityInitialization>,
) -> FeosResult<Self>
pub fn new( eos: &E, temperature: Option<Temperature<D>>, volume: Option<Volume<D>>, density: Option<Density<D>>, partial_density: Option<&Density<OVector<D, N>>>, total_moles: Option<Moles<D>>, moles: Option<&Moles<OVector<D, N>>>, molefracs: Option<&OVector<D, N>>, pressure: Option<Pressure<D>>, density_initialization: Option<DensityInitialization>, ) -> FeosResult<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: &E,
temperature: Temperature<D>,
pressure: Pressure<D>,
moles: &Moles<OVector<D, N>>,
density_initialization: Option<DensityInitialization>,
) -> FeosResult<Self>
pub fn new_npt( eos: &E, temperature: Temperature<D>, pressure: Pressure<D>, moles: &Moles<OVector<D, N>>, density_initialization: Option<DensityInitialization>, ) -> FeosResult<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_xpt(
eos: &E,
temperature: Temperature<D>,
pressure: Pressure<D>,
molefracs: &OVector<D, N>,
density_initialization: Option<DensityInitialization>,
) -> FeosResult<Self>
pub fn new_xpt( eos: &E, temperature: Temperature<D>, pressure: Pressure<D>, molefracs: &OVector<D, N>, density_initialization: Option<DensityInitialization>, ) -> FeosResult<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: &E,
temperature: Temperature<D>,
pressure: Pressure<D>,
volume: Volume<D>,
molefracs: &OVector<D, N>,
density_initialization: Option<DensityInitialization>,
) -> FeosResult<Self>
pub fn new_npvx( eos: &E, temperature: Temperature<D>, pressure: Pressure<D>, volume: Volume<D>, molefracs: &OVector<D, N>, density_initialization: Option<DensityInitialization>, ) -> FeosResult<Self>
Return a new State for given pressure $p$, volume $V$, temperature $T$ and composition $x_i$.
Source§impl<E: Total<N, D>, N: Gradients, D: DualNum<f64> + Copy> State<E, N, D>where
DefaultAllocator: Allocator<N>,
impl<E: Total<N, D>, N: Gradients, D: DualNum<f64> + Copy> State<E, N, D>where
DefaultAllocator: Allocator<N>,
Sourcepub fn new_full(
eos: &E,
temperature: Option<Temperature<D>>,
volume: Option<Volume<D>>,
density: Option<Density<D>>,
partial_density: Option<&Density<OVector<D, N>>>,
total_moles: Option<Moles<D>>,
moles: Option<&Moles<OVector<D, N>>>,
molefracs: Option<&OVector<D, N>>,
pressure: Option<Pressure<D>>,
molar_enthalpy: Option<MolarEnergy<D>>,
molar_entropy: Option<MolarEntropy<D>>,
molar_internal_energy: Option<MolarEnergy<D>>,
density_initialization: Option<DensityInitialization>,
initial_temperature: Option<Temperature<D>>,
) -> FeosResult<Self>
pub fn new_full( eos: &E, temperature: Option<Temperature<D>>, volume: Option<Volume<D>>, density: Option<Density<D>>, partial_density: Option<&Density<OVector<D, N>>>, total_moles: Option<Moles<D>>, moles: Option<&Moles<OVector<D, N>>>, molefracs: Option<&OVector<D, N>>, pressure: Option<Pressure<D>>, molar_enthalpy: Option<MolarEnergy<D>>, molar_entropy: Option<MolarEntropy<D>>, molar_internal_energy: Option<MolarEnergy<D>>, density_initialization: Option<DensityInitialization>, initial_temperature: Option<Temperature<D>>, ) -> FeosResult<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: &E,
pressure: Pressure<D>,
molar_enthalpy: MolarEnergy<D>,
moles: &Moles<OVector<D, N>>,
density_initialization: Option<DensityInitialization>,
initial_temperature: Option<Temperature<D>>,
) -> FeosResult<Self>
pub fn new_nph( eos: &E, pressure: Pressure<D>, molar_enthalpy: MolarEnergy<D>, moles: &Moles<OVector<D, N>>, density_initialization: Option<DensityInitialization>, initial_temperature: Option<Temperature<D>>, ) -> FeosResult<Self>
Return a new State for given pressure $p$ and molar enthalpy $h$.
Sourcepub fn new_nth(
eos: &E,
temperature: Temperature<D>,
molar_enthalpy: MolarEnergy<D>,
moles: &Moles<OVector<D, N>>,
density_initialization: Option<DensityInitialization>,
) -> FeosResult<Self>
pub fn new_nth( eos: &E, temperature: Temperature<D>, molar_enthalpy: MolarEnergy<D>, moles: &Moles<OVector<D, N>>, density_initialization: Option<DensityInitialization>, ) -> FeosResult<Self>
Return a new State for given temperature $T$ and molar enthalpy $h$.
Sourcepub fn new_nts(
eos: &E,
temperature: Temperature<D>,
molar_entropy: MolarEntropy<D>,
moles: &Moles<OVector<D, N>>,
density_initialization: Option<DensityInitialization>,
) -> FeosResult<Self>
pub fn new_nts( eos: &E, temperature: Temperature<D>, molar_entropy: MolarEntropy<D>, moles: &Moles<OVector<D, N>>, density_initialization: Option<DensityInitialization>, ) -> FeosResult<Self>
Return a new State for given temperature $T$ and molar entropy $s$.
Sourcepub fn new_nps(
eos: &E,
pressure: Pressure<D>,
molar_entropy: MolarEntropy<D>,
moles: &Moles<OVector<D, N>>,
density_initialization: Option<DensityInitialization>,
initial_temperature: Option<Temperature<D>>,
) -> FeosResult<Self>
pub fn new_nps( eos: &E, pressure: Pressure<D>, molar_entropy: MolarEntropy<D>, moles: &Moles<OVector<D, N>>, density_initialization: Option<DensityInitialization>, initial_temperature: Option<Temperature<D>>, ) -> FeosResult<Self>
Return a new State for given pressure $p$ and molar entropy $s$.
Sourcepub fn new_nvu(
eos: &E,
volume: Volume<D>,
molar_internal_energy: MolarEnergy<D>,
moles: &Moles<OVector<D, N>>,
initial_temperature: Option<Temperature<D>>,
) -> FeosResult<Self>
pub fn new_nvu( eos: &E, volume: Volume<D>, molar_internal_energy: MolarEnergy<D>, moles: &Moles<OVector<D, N>>, initial_temperature: Option<Temperature<D>>, ) -> FeosResult<Self>
Return a new State for given volume $V$ and molar internal energy $u$.
Trait Implementations§
Source§impl<E: Clone, N: Dim, D: DualNum<f64> + Copy> Clone for State<E, N, D>where
DefaultAllocator: Allocator<N>,
impl<E: Clone, N: Dim, D: DualNum<f64> + Copy> Clone for State<E, N, D>where
DefaultAllocator: Allocator<N>,
Source§impl<E: Debug, N: Debug + Dim, D: Debug + DualNum<f64> + Copy> Debug for State<E, N, D>where
DefaultAllocator: Allocator<N>,
impl<E: Debug, N: Debug + Dim, D: Debug + DualNum<f64> + Copy> Debug for State<E, N, D>where
DefaultAllocator: Allocator<N>,
Auto Trait Implementations§
impl<E, N = Dyn, D = f64> !Freeze for State<E, N, D>
impl<E, N = Dyn, D = f64> !RefUnwindSafe for State<E, N, D>
impl<E, N = Dyn, D = f64> !Send for State<E, N, D>
impl<E, N = Dyn, D = f64> !Sync for State<E, N, D>
impl<E, N = Dyn, D = f64> !Unpin for State<E, N, D>
impl<E, N = Dyn, D = f64> !UnwindSafe for State<E, N, D>
Blanket Implementations§
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> 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<T> PropertiesAD for T
impl<T> PropertiesAD for T
fn vapor_pressure<const P: usize>( &self, temperature: Temperature, ) -> FeosResult<Pressure<DualSVec<f64, f64, P>>>
fn equilibrium_liquid_density<const P: usize>( &self, temperature: Temperature, ) -> FeosResult<(Pressure<DualSVec<f64, f64, P>>, Density<DualSVec<f64, f64, P>>)>
fn liquid_density<const P: usize>( &self, temperature: Temperature, pressure: Pressure, ) -> FeosResult<Density<DualSVec<f64, f64, P>>>
fn vapor_pressure_parallel<const P: usize>(
parameter_names: [String; P],
parameters: ArrayView2<'_, f64>,
input: ArrayView2<'_, f64>,
) -> (Array1<f64>, Array2<f64>, Array1<bool>)where
Self: ParametersAD<1>,
fn liquid_density_parallel<const P: usize>(
parameter_names: [String; P],
parameters: ArrayView2<'_, f64>,
input: ArrayView2<'_, f64>,
) -> (Array1<f64>, Array2<f64>, Array1<bool>)where
Self: ParametersAD<1>,
fn equilibrium_liquid_density_parallel<const P: usize>(
parameter_names: [String; P],
parameters: ArrayView2<'_, f64>,
input: ArrayView2<'_, f64>,
) -> (Array1<f64>, Array2<f64>, Array1<bool>)where
Self: ParametersAD<1>,
fn bubble_point_pressure<const P: usize>( &self, temperature: Temperature, pressure: Option<Pressure>, liquid_molefracs: SVector<f64, 2>, ) -> FeosResult<Pressure<DualSVec<f64, f64, P>>>
fn dew_point_pressure<const P: usize>( &self, temperature: Temperature, pressure: Option<Pressure>, vapor_molefracs: SVector<f64, 2>, ) -> FeosResult<Pressure<DualSVec<f64, f64, P>>>
fn bubble_point_pressure_parallel<const P: usize>(
parameter_names: [String; P],
parameters: ArrayView2<'_, f64>,
input: ArrayView2<'_, f64>,
) -> (Array1<f64>, Array2<f64>, Array1<bool>)where
Self: ParametersAD<2>,
fn dew_point_pressure_parallel<const P: usize>(
parameter_names: [String; P],
parameters: ArrayView2<'_, f64>,
input: ArrayView2<'_, f64>,
) -> (Array1<f64>, Array2<f64>, Array1<bool>)where
Self: ParametersAD<2>,
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.