pub struct State<U, E> {
pub eos: Rc<E>,
pub temperature: QuantityScalar<U>,
pub volume: QuantityScalar<U>,
pub moles: QuantityArray1<U>,
pub total_moles: QuantityScalar<U>,
pub partial_density: QuantityArray1<U>,
pub density: QuantityScalar<U>,
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: Rc<E>
Equation of state
temperature: QuantityScalar<U>
Temperature $T$
volume: QuantityScalar<U>
Volume $V$
moles: QuantityArray1<U>
Mole numbers $N_i$
total_moles: QuantityScalar<U>
Total number of moles $N=\sum_iN_i$
partial_density: QuantityArray1<U>
Partial densities $\rho_i=\frac{N_i}{V}$
density: QuantityScalar<U>
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
sourceimpl<U: EosUnit, E: EquationOfState> State<U, E>
impl<U: EosUnit, E: EquationOfState> State<U, E>
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<U, E>>>
pub fn stability_analysis(
&self,
options: SolverOptions
) -> EosResult<Vec<State<U, 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.
sourceimpl<U: EosUnit, E: EquationOfState> State<U, E>
impl<U: EosUnit, E: EquationOfState> State<U, E>
sourcepub fn tp_flash(
&self,
initial_state: Option<&PhaseEquilibrium<U, E, 2>>,
options: SolverOptions,
non_volatile_components: Option<Vec<usize>>
) -> EosResult<PhaseEquilibrium<U, E, 2>>
pub fn tp_flash(
&self,
initial_state: Option<&PhaseEquilibrium<U, E, 2>>,
options: SolverOptions,
non_volatile_components: Option<Vec<usize>>
) -> EosResult<PhaseEquilibrium<U, 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).
sourceimpl<U: EosUnit, E: EquationOfState> State<U, E>
impl<U: EosUnit, E: EquationOfState> State<U, E>
sourcepub fn pressure(&self, contributions: Contributions) -> QuantityScalar<U>
pub fn pressure(&self, contributions: Contributions) -> QuantityScalar<U>
Pressure: $p=-\left(\frac{\partial A}{\partial V}\right)_{T,N_i}$
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) -> QuantityScalar<U>
pub fn dp_dv(&self, contributions: Contributions) -> QuantityScalar<U>
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) -> QuantityScalar<U>
pub fn dp_drho(&self, contributions: Contributions) -> QuantityScalar<U>
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) -> QuantityScalar<U>
pub fn dp_dt(&self, contributions: Contributions) -> QuantityScalar<U>
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) -> QuantityArray1<U>
pub fn dp_dni(&self, contributions: Contributions) -> QuantityArray1<U>
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) -> QuantityScalar<U>
pub fn d2p_dv2(&self, contributions: Contributions) -> QuantityScalar<U>
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) -> QuantityScalar<U>
pub fn d2p_drho2(&self, contributions: Contributions) -> QuantityScalar<U>
Second partial derivative of pressure w.r.t. density: $\left(\frac{\partial^2 p}{\partial \rho^2}\right)_{T,N_j}$
sourcepub fn molar_volume(&self, contributions: Contributions) -> QuantityArray1<U>
pub fn molar_volume(&self, contributions: Contributions) -> QuantityArray1<U>
Partial molar volume: $v_i=\left(\frac{\partial V}{\partial N_i}\right)_{T,p,N_j}$
sourcepub fn chemical_potential(
&self,
contributions: Contributions
) -> QuantityArray1<U>
pub fn chemical_potential(
&self,
contributions: Contributions
) -> QuantityArray1<U>
Chemical potential: $\mu_i=\left(\frac{\partial A}{\partial N_i}\right)_{T,V,N_j}$
sourcepub fn dmu_dt(&self, contributions: Contributions) -> QuantityArray1<U>
pub fn dmu_dt(&self, contributions: Contributions) -> QuantityArray1<U>
Partial derivative of chemical potential w.r.t. temperature: $\left(\frac{\partial\mu_i}{\partial T}\right)_{V,N_i}$
sourcepub fn dmu_dni(&self, contributions: Contributions) -> QuantityArray2<U>
pub fn dmu_dni(&self, contributions: Contributions) -> QuantityArray2<U>
Partial derivative of chemical potential w.r.t. moles: $\left(\frac{\partial\mu_i}{\partial \mu_j}\right)_{T,V,N_k}$
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(&self) -> EosResult<Array1<f64>>
pub fn ln_phi_pure(&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(T, p)$
sourcepub fn dln_phi_dt(&self) -> QuantityArray1<U>
pub fn dln_phi_dt(&self) -> QuantityArray1<U>
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) -> QuantityArray1<U>
pub fn dln_phi_dp(&self) -> QuantityArray1<U>
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) -> QuantityArray2<U>
pub fn dln_phi_dnj(&self) -> QuantityArray2<U>
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 c_v(&self, contributions: Contributions) -> QuantityScalar<U>
pub fn c_v(&self, contributions: Contributions) -> QuantityScalar<U>
Molar isochoric heat capacity: $c_v=\left(\frac{\partial u}{\partial T}\right)_{V,N_i}$
sourcepub fn dc_v_dt(&self, contributions: Contributions) -> QuantityScalar<U>
pub fn dc_v_dt(&self, contributions: Contributions) -> QuantityScalar<U>
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 c_p(&self, contributions: Contributions) -> QuantityScalar<U>
pub fn c_p(&self, contributions: Contributions) -> QuantityScalar<U>
Molar isobaric heat capacity: $c_p=\left(\frac{\partial h}{\partial T}\right)_{p,N_i}$
sourcepub fn entropy(&self, contributions: Contributions) -> QuantityScalar<U>
pub fn entropy(&self, contributions: Contributions) -> QuantityScalar<U>
Entropy: $S=-\left(\frac{\partial A}{\partial T}\right)_{V,N_i}$
sourcepub fn ds_dt(&self, contributions: Contributions) -> QuantityScalar<U>
pub fn ds_dt(&self, contributions: Contributions) -> QuantityScalar<U>
Partial derivative of the entropy w.r.t. temperature: $\left(\frac{\partial S}{\partial T}\right)_{V,N_i}$
sourcepub fn molar_entropy(&self, contributions: Contributions) -> QuantityScalar<U>
pub fn molar_entropy(&self, contributions: Contributions) -> QuantityScalar<U>
molar entropy: $s=\frac{S}{N}$
sourcepub fn enthalpy(&self, contributions: Contributions) -> QuantityScalar<U>
pub fn enthalpy(&self, contributions: Contributions) -> QuantityScalar<U>
Enthalpy: $H=A+TS+pV$
sourcepub fn molar_enthalpy(&self, contributions: Contributions) -> QuantityScalar<U>
pub fn molar_enthalpy(&self, contributions: Contributions) -> QuantityScalar<U>
molar enthalpy: $h=\frac{H}{N}$
sourcepub fn helmholtz_energy(
&self,
contributions: Contributions
) -> QuantityScalar<U>
pub fn helmholtz_energy(
&self,
contributions: Contributions
) -> QuantityScalar<U>
Helmholtz energy: $A$
sourcepub fn molar_helmholtz_energy(
&self,
contributions: Contributions
) -> QuantityScalar<U>
pub fn molar_helmholtz_energy(
&self,
contributions: Contributions
) -> QuantityScalar<U>
molar Helmholtz energy: $a=\frac{A}{N}$
sourcepub fn internal_energy(&self, contributions: Contributions) -> QuantityScalar<U>
pub fn internal_energy(&self, contributions: Contributions) -> QuantityScalar<U>
Internal energy: $U=A+TS$
sourcepub fn molar_internal_energy(
&self,
contributions: Contributions
) -> QuantityScalar<U>
pub fn molar_internal_energy(
&self,
contributions: Contributions
) -> QuantityScalar<U>
Molar internal energy: $u=\frac{U}{N}$
sourcepub fn gibbs_energy(&self, contributions: Contributions) -> QuantityScalar<U>
pub fn gibbs_energy(&self, contributions: Contributions) -> QuantityScalar<U>
Gibbs energy: $G=A+pV$
sourcepub fn molar_gibbs_energy(
&self,
contributions: Contributions
) -> QuantityScalar<U>
pub fn molar_gibbs_energy(
&self,
contributions: Contributions
) -> QuantityScalar<U>
Molar Gibbs energy: $g=\frac{G}{N}$
sourcepub fn partial_molar_entropy(
&self,
contributions: Contributions
) -> QuantityArray1<U>
pub fn partial_molar_entropy(
&self,
contributions: Contributions
) -> QuantityArray1<U>
Partial molar entropy: $s_i=\left(\frac{\partial S}{\partial N_i}\right)_{T,p,N_j}$
sourcepub fn partial_molar_enthalpy(
&self,
contributions: Contributions
) -> QuantityArray1<U>
pub fn partial_molar_enthalpy(
&self,
contributions: Contributions
) -> QuantityArray1<U>
Partial molar enthalpy: $h_i=\left(\frac{\partial H}{\partial N_i}\right)_{T,p,N_j}$
sourcepub fn joule_thomson(&self) -> QuantityScalar<U>
pub fn joule_thomson(&self) -> QuantityScalar<U>
Joule Thomson coefficient: $\mu_{JT}=\left(\frac{\partial T}{\partial p}\right)_{H,N_i}$
sourcepub fn isentropic_compressibility(&self) -> QuantityScalar<U>
pub fn isentropic_compressibility(&self) -> QuantityScalar<U>
Isentropic compressibility: $\kappa_s=-\frac{1}{V}\left(\frac{\partial V}{\partial p}\right)_{S,N_i}$
sourcepub fn isothermal_compressibility(&self) -> QuantityScalar<U>
pub fn isothermal_compressibility(&self) -> QuantityScalar<U>
Isothermal compressibility: $\kappa_T=-\frac{1}{V}\left(\frac{\partial V}{\partial p}\right)_{T,N_i}$
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 helmholtz_energy_contributions(&self) -> Vec<(String, QuantityScalar<U>)>ⓘNotable traits for Vec<u8, A>impl<A> Write for Vec<u8, A> where
A: Allocator,
pub fn helmholtz_energy_contributions(&self) -> Vec<(String, QuantityScalar<U>)>ⓘNotable traits for Vec<u8, A>impl<A> Write for Vec<u8, A> where
A: Allocator,
A: Allocator,
Helmholtz energy $A$ evaluated for each contribution of the equation of state.
sourceimpl<U: EosUnit, E: EquationOfState + MolarWeight<U>> State<U, E>
impl<U: EosUnit, E: EquationOfState + MolarWeight<U>> State<U, E>
Mass specific state properties
These properties are available for equations of state that implement the MolarWeight trait.
sourcepub fn total_molar_weight(&self) -> QuantityScalar<U>
pub fn total_molar_weight(&self) -> QuantityScalar<U>
Total molar weight: $MW=\sum_ix_iMW_i$
sourcepub fn mass(&self) -> QuantityArray1<U>
pub fn mass(&self) -> QuantityArray1<U>
Mass of each component: $m_i=n_iMW_i$
sourcepub fn total_mass(&self) -> QuantityScalar<U>
pub fn total_mass(&self) -> QuantityScalar<U>
Total mass: $m=\sum_im_i=nMW$
sourcepub fn mass_density(&self) -> QuantityScalar<U>
pub fn mass_density(&self) -> QuantityScalar<U>
Mass density: $\rho^{(m)}=\frac{m}{V}$
sourcepub fn specific_entropy(
&self,
contributions: Contributions
) -> QuantityScalar<U>
pub fn specific_entropy(
&self,
contributions: Contributions
) -> QuantityScalar<U>
Specific entropy: $s^{(m)}=\frac{S}{m}$
sourcepub fn specific_enthalpy(
&self,
contributions: Contributions
) -> QuantityScalar<U>
pub fn specific_enthalpy(
&self,
contributions: Contributions
) -> QuantityScalar<U>
Specific enthalpy: $h^{(m)}=\frac{H}{m}$
sourcepub fn specific_helmholtz_energy(
&self,
contributions: Contributions
) -> QuantityScalar<U>
pub fn specific_helmholtz_energy(
&self,
contributions: Contributions
) -> QuantityScalar<U>
Specific Helmholtz energy: $a^{(m)}=\frac{A}{m}$
sourcepub fn specific_internal_energy(
&self,
contributions: Contributions
) -> QuantityScalar<U>
pub fn specific_internal_energy(
&self,
contributions: Contributions
) -> QuantityScalar<U>
Specific internal energy: $u^{(m)}=\frac{U}{m}$
sourcepub fn specific_gibbs_energy(
&self,
contributions: Contributions
) -> QuantityScalar<U>
pub fn specific_gibbs_energy(
&self,
contributions: Contributions
) -> QuantityScalar<U>
Specific Gibbs energy: $g^{(m)}=\frac{G}{m}$
sourcepub fn speed_of_sound(&self) -> QuantityScalar<U>
pub fn speed_of_sound(&self) -> QuantityScalar<U>
Speed of sound: $c=\sqrt{\left(\frac{\partial p}{\partial\rho}\right)_{S,N_i}}$
sourceimpl<U: EosUnit, E: EquationOfState + EntropyScaling<U>> State<U, E>
impl<U: EosUnit, E: EquationOfState + EntropyScaling<U>> State<U, E>
Transport properties
These properties are available for equations of state that implement the EntropyScaling trait.
sourcepub fn viscosity(&self) -> EosResult<QuantityScalar<U>>
pub fn viscosity(&self) -> EosResult<QuantityScalar<U>>
Return the viscosity via entropy scaling.
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<QuantityScalar<U>>
pub fn viscosity_reference(&self) -> EosResult<QuantityScalar<U>>
Return the viscosity reference as used in entropy scaling.
sourcepub fn diffusion(&self) -> EosResult<QuantityScalar<U>>
pub fn diffusion(&self) -> EosResult<QuantityScalar<U>>
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<QuantityScalar<U>>
pub fn diffusion_reference(&self) -> EosResult<QuantityScalar<U>>
Return the diffusion reference as used in entropy scaling.
sourcepub fn thermal_conductivity(&self) -> EosResult<QuantityScalar<U>>
pub fn thermal_conductivity(&self) -> EosResult<QuantityScalar<U>>
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<QuantityScalar<U>>
pub fn thermal_conductivity_reference(&self) -> EosResult<QuantityScalar<U>>
Return the thermal conductivity reference as used in entropy scaling.
sourceimpl<U: EosUnit, E: EquationOfState> State<U, E>
impl<U: EosUnit, E: EquationOfState> State<U, E>
sourcepub fn critical_point_pure(
eos: &Rc<E>,
initial_temperature: Option<QuantityScalar<U>>,
options: SolverOptions
) -> EosResult<Vec<Self>> where
QuantityScalar<U>: Display,
pub fn critical_point_pure(
eos: &Rc<E>,
initial_temperature: Option<QuantityScalar<U>>,
options: SolverOptions
) -> EosResult<Vec<Self>> where
QuantityScalar<U>: Display,
Calculate the pure component critical point of all components.
pub fn critical_point_binary(
eos: &Rc<E>,
temperature_or_pressure: QuantityScalar<U>,
initial_temperature: Option<QuantityScalar<U>>,
initial_molefracs: Option<[f64; 2]>,
options: SolverOptions
) -> EosResult<Self> where
QuantityScalar<U>: Display,
sourcepub fn critical_point(
eos: &Rc<E>,
moles: Option<&QuantityArray1<U>>,
initial_temperature: Option<QuantityScalar<U>>,
options: SolverOptions
) -> EosResult<Self> where
QuantityScalar<U>: Display,
pub fn critical_point(
eos: &Rc<E>,
moles: Option<&QuantityArray1<U>>,
initial_temperature: Option<QuantityScalar<U>>,
options: SolverOptions
) -> EosResult<Self> where
QuantityScalar<U>: Display,
Calculate the critical point of a system for given moles.
sourceimpl<U: EosUnit, E: EquationOfState> State<U, E>
impl<U: EosUnit, E: EquationOfState> State<U, E>
sourcepub fn new_nvt(
eos: &Rc<E>,
temperature: QuantityScalar<U>,
volume: QuantityScalar<U>,
moles: &QuantityArray1<U>
) -> EosResult<Self>
pub fn new_nvt(
eos: &Rc<E>,
temperature: QuantityScalar<U>,
volume: QuantityScalar<U>,
moles: &QuantityArray1<U>
) -> 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: &Rc<E>,
temperature: QuantityScalar<U>,
density: QuantityScalar<U>
) -> EosResult<Self>
pub fn new_pure(
eos: &Rc<E>,
temperature: QuantityScalar<U>,
density: QuantityScalar<U>
) -> 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: &Rc<E>,
temperature: Option<QuantityScalar<U>>,
volume: Option<QuantityScalar<U>>,
density: Option<QuantityScalar<U>>,
partial_density: Option<&QuantityArray1<U>>,
total_moles: Option<QuantityScalar<U>>,
moles: Option<&QuantityArray1<U>>,
molefracs: Option<&Array1<f64>>,
pressure: Option<QuantityScalar<U>>,
molar_enthalpy: Option<QuantityScalar<U>>,
molar_entropy: Option<QuantityScalar<U>>,
molar_internal_energy: Option<QuantityScalar<U>>,
density_initialization: DensityInitialization<U>,
initial_temperature: Option<QuantityScalar<U>>
) -> EosResult<Self>
pub fn new(
eos: &Rc<E>,
temperature: Option<QuantityScalar<U>>,
volume: Option<QuantityScalar<U>>,
density: Option<QuantityScalar<U>>,
partial_density: Option<&QuantityArray1<U>>,
total_moles: Option<QuantityScalar<U>>,
moles: Option<&QuantityArray1<U>>,
molefracs: Option<&Array1<f64>>,
pressure: Option<QuantityScalar<U>>,
molar_enthalpy: Option<QuantityScalar<U>>,
molar_entropy: Option<QuantityScalar<U>>,
molar_internal_energy: Option<QuantityScalar<U>>,
density_initialization: DensityInitialization<U>,
initial_temperature: Option<QuantityScalar<U>>
) -> 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_npt(
eos: &Rc<E>,
temperature: QuantityScalar<U>,
pressure: QuantityScalar<U>,
moles: &QuantityArray1<U>,
density_initialization: DensityInitialization<U>
) -> EosResult<Self>
pub fn new_npt(
eos: &Rc<E>,
temperature: QuantityScalar<U>,
pressure: QuantityScalar<U>,
moles: &QuantityArray1<U>,
density_initialization: DensityInitialization<U>
) -> 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: &Rc<E>,
temperature: QuantityScalar<U>,
pressure: QuantityScalar<U>,
volume: QuantityScalar<U>,
molefracs: &Array1<f64>,
density_initialization: DensityInitialization<U>
) -> EosResult<Self>
pub fn new_npvx(
eos: &Rc<E>,
temperature: QuantityScalar<U>,
pressure: QuantityScalar<U>,
volume: QuantityScalar<U>,
molefracs: &Array1<f64>,
density_initialization: DensityInitialization<U>
) -> EosResult<Self>
Return a new State
for given pressure $p$, volume $V$, temperature $T$ and composition $x_i$.
sourcepub fn new_nph(
eos: &Rc<E>,
pressure: QuantityScalar<U>,
molar_enthalpy: QuantityScalar<U>,
moles: &QuantityArray1<U>,
density_initialization: DensityInitialization<U>,
initial_temperature: Option<QuantityScalar<U>>
) -> EosResult<Self>
pub fn new_nph(
eos: &Rc<E>,
pressure: QuantityScalar<U>,
molar_enthalpy: QuantityScalar<U>,
moles: &QuantityArray1<U>,
density_initialization: DensityInitialization<U>,
initial_temperature: Option<QuantityScalar<U>>
) -> EosResult<Self>
Return a new State
for given pressure $p$ and molar enthalpy $h$.
sourcepub fn new_nth(
eos: &Rc<E>,
temperature: QuantityScalar<U>,
molar_enthalpy: QuantityScalar<U>,
moles: &QuantityArray1<U>,
density_initialization: DensityInitialization<U>
) -> EosResult<Self>
pub fn new_nth(
eos: &Rc<E>,
temperature: QuantityScalar<U>,
molar_enthalpy: QuantityScalar<U>,
moles: &QuantityArray1<U>,
density_initialization: DensityInitialization<U>
) -> EosResult<Self>
Return a new State
for given temperature $T$ and molar enthalpy $h$.
sourcepub fn new_nts(
eos: &Rc<E>,
temperature: QuantityScalar<U>,
molar_entropy: QuantityScalar<U>,
moles: &QuantityArray1<U>,
density_initialization: DensityInitialization<U>
) -> EosResult<Self>
pub fn new_nts(
eos: &Rc<E>,
temperature: QuantityScalar<U>,
molar_entropy: QuantityScalar<U>,
moles: &QuantityArray1<U>,
density_initialization: DensityInitialization<U>
) -> EosResult<Self>
Return a new State
for given temperature $T$ and molar entropy $s$.
sourcepub fn new_nps(
eos: &Rc<E>,
pressure: QuantityScalar<U>,
molar_entropy: QuantityScalar<U>,
moles: &QuantityArray1<U>,
density_initialization: DensityInitialization<U>,
initial_temperature: Option<QuantityScalar<U>>
) -> EosResult<Self>
pub fn new_nps(
eos: &Rc<E>,
pressure: QuantityScalar<U>,
molar_entropy: QuantityScalar<U>,
moles: &QuantityArray1<U>,
density_initialization: DensityInitialization<U>,
initial_temperature: Option<QuantityScalar<U>>
) -> EosResult<Self>
Return a new State
for given pressure $p$ and molar entropy $s$.
sourcepub fn new_nvu(
eos: &Rc<E>,
volume: QuantityScalar<U>,
molar_internal_energy: QuantityScalar<U>,
moles: &QuantityArray1<U>,
initial_temperature: Option<QuantityScalar<U>>
) -> EosResult<Self>
pub fn new_nvu(
eos: &Rc<E>,
volume: QuantityScalar<U>,
molar_internal_energy: QuantityScalar<U>,
moles: &QuantityArray1<U>,
initial_temperature: Option<QuantityScalar<U>>
) -> EosResult<Self>
Return a new State
for given volume $V$ and molar internal energy $u$.
sourcepub fn update_temperature(
&self,
temperature: QuantityScalar<U>
) -> EosResult<Self>
pub fn update_temperature(
&self,
temperature: QuantityScalar<U>
) -> EosResult<Self>
Update the state with the given temperature
sourcepub fn update_chemical_potential(
&mut self,
chemical_potential: &QuantityArray1<U>
) -> EosResult<()>
pub fn update_chemical_potential(
&mut self,
chemical_potential: &QuantityArray1<U>
) -> EosResult<()>
Update the state with the given chemical potential.
sourcepub fn update_gibbs_energy(
self,
molar_gibbs_energy: QuantityScalar<U>
) -> EosResult<Self>
pub fn update_gibbs_energy(
self,
molar_gibbs_energy: QuantityScalar<U>
) -> EosResult<Self>
Update the state with the given molar Gibbs energy.
Trait Implementations
sourceimpl<U, E> Display for State<U, E> where
QuantityScalar<U>: Display,
QuantityArray1<U>: Display,
E: EquationOfState,
impl<U, E> Display for State<U, E> where
QuantityScalar<U>: Display,
QuantityArray1<U>: Display,
E: EquationOfState,
sourceimpl<'a, U, E> FromIterator<&'a State<U, E>> for StateVec<'a, U, E>
impl<'a, U, E> FromIterator<&'a State<U, E>> for StateVec<'a, U, E>
sourcefn from_iter<I: IntoIterator<Item = &'a State<U, E>>>(iter: I) -> Self
fn from_iter<I: IntoIterator<Item = &'a State<U, E>>>(iter: I) -> Self
Creates a value from an iterator. Read more
Auto Trait Implementations
impl<U, E> !RefUnwindSafe for State<U, E>
impl<U, E> !Send for State<U, E>
impl<U, E> !Sync for State<U, E>
impl<U, E> Unpin for State<U, E> where
U: Unpin,
impl<U, E> UnwindSafe for State<U, E> where
E: RefUnwindSafe,
U: UnwindSafe,
Blanket Implementations
sourceimpl<T> BorrowMut<T> for T where
T: ?Sized,
impl<T> BorrowMut<T> for T where
T: ?Sized,
const: unstable · sourcefn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
Mutably borrows from an owned value. Read more
sourceimpl<T> ToOwned for T where
T: Clone,
impl<T> ToOwned for T where
T: Clone,
type Owned = T
type Owned = T
The resulting type after obtaining ownership.
sourcefn clone_into(&self, target: &mut T)
fn clone_into(&self, target: &mut T)
toowned_clone_into
)Uses borrowed data to replace owned data, usually by cloning. Read more