1use crate::errors::{FeosError, FeosResult};
2use crate::phase_equilibria::PhaseEquilibrium;
3use crate::state::{
4 Contributions,
5 DensityInitialization::{InitialDensity, Liquid, Vapor},
6};
7use crate::{ReferenceSystem, Residual, SolverOptions, State, Verbosity};
8use nalgebra::allocator::Allocator;
9use nalgebra::{DMatrix, DVector, DefaultAllocator, Dim, Dyn, OVector, U1};
10#[cfg(feature = "ndarray")]
11use ndarray::Array1;
12use num_dual::linalg::LU;
13use num_dual::{DualNum, DualStruct, Gradients};
14use quantity::{Density, Dimensionless, Moles, Pressure, Quantity, RGAS, SIUnit, Temperature};
15use typenum::{N1, N2, P1, Z0};
16
17const MAX_ITER_INNER: usize = 5;
18const TOL_INNER: f64 = 1e-9;
19const MAX_ITER_OUTER: usize = 400;
20const TOL_OUTER: f64 = 1e-10;
21
22const MAX_TSTEP: f64 = 20.0;
23const MAX_LNPSTEP: f64 = 0.1;
24const NEWTON_TOL: f64 = 1e-3;
25
26pub trait TemperatureOrPressure<D: DualNum<f64> + Copy = f64>: Copy {
28 type Other: Copy;
29
30 const IDENTIFIER: &'static str;
31
32 fn temperature(&self) -> Option<Temperature<D>>;
33 fn pressure(&self) -> Option<Pressure<D>>;
34
35 fn temperature_pressure(
36 &self,
37 tp_init: Option<Self::Other>,
38 ) -> (Option<Temperature<D>>, Option<Pressure<D>>, bool);
39
40 fn from_state<E: Residual<N, D>, N: Gradients>(state: &State<E, N, D>) -> Self::Other
41 where
42 DefaultAllocator: Allocator<N>;
43
44 #[cfg(feature = "ndarray")]
45 fn linspace(
46 &self,
47 start: Self::Other,
48 end: Self::Other,
49 n: usize,
50 ) -> (Temperature<Array1<f64>>, Pressure<Array1<f64>>);
51}
52
53impl<D: DualNum<f64> + Copy> TemperatureOrPressure<D> for Temperature<D> {
54 type Other = Pressure<D>;
55 const IDENTIFIER: &'static str = "temperature";
56
57 fn temperature(&self) -> Option<Temperature<D>> {
58 Some(*self)
59 }
60
61 fn pressure(&self) -> Option<Pressure<D>> {
62 None
63 }
64
65 fn temperature_pressure(
66 &self,
67 tp_init: Option<Self::Other>,
68 ) -> (Option<Temperature<D>>, Option<Pressure<D>>, bool) {
69 (Some(*self), tp_init, true)
70 }
71
72 fn from_state<E: Residual<N, D>, N: Gradients>(state: &State<E, N, D>) -> Self::Other
73 where
74 DefaultAllocator: Allocator<N>,
75 {
76 state.pressure(Contributions::Total)
77 }
78
79 #[cfg(feature = "ndarray")]
80 fn linspace(
81 &self,
82 start: Pressure<D>,
83 end: Pressure<D>,
84 n: usize,
85 ) -> (Temperature<Array1<f64>>, Pressure<Array1<f64>>) {
86 (
87 Temperature::linspace(self.re(), self.re(), n),
88 Pressure::linspace(start.re(), end.re(), n),
89 )
90 }
91}
92
93impl<D: DualNum<f64> + Copy> TemperatureOrPressure<D>
97 for Quantity<D, SIUnit<N2, N1, P1, Z0, Z0, Z0, Z0>>
98{
99 type Other = Temperature<D>;
100 const IDENTIFIER: &'static str = "pressure";
101
102 fn temperature(&self) -> Option<Temperature<D>> {
103 None
104 }
105
106 fn pressure(&self) -> Option<Pressure<D>> {
107 Some(*self)
108 }
109
110 fn temperature_pressure(
111 &self,
112 tp_init: Option<Self::Other>,
113 ) -> (Option<Temperature<D>>, Option<Pressure<D>>, bool) {
114 (tp_init, Some(*self), false)
115 }
116
117 fn from_state<E: Residual<N, D>, N: Dim>(state: &State<E, N, D>) -> Self::Other
118 where
119 DefaultAllocator: Allocator<N>,
120 {
121 state.temperature
122 }
123
124 #[cfg(feature = "ndarray")]
125 fn linspace(
126 &self,
127 start: Temperature<D>,
128 end: Temperature<D>,
129 n: usize,
130 ) -> (Temperature<Array1<f64>>, Pressure<Array1<f64>>) {
131 (
132 Temperature::linspace(start.re(), end.re(), n),
133 Pressure::linspace(self.re(), self.re(), n),
134 )
135 }
136}
137
138impl<E: Residual<N, D>, N: Gradients, D: DualNum<f64> + Copy> PhaseEquilibrium<E, 2, N, D>
140where
141 DefaultAllocator: Allocator<N> + Allocator<N, N> + Allocator<U1, N>,
142{
143 pub fn bubble_point<TP: TemperatureOrPressure<D>>(
146 eos: &E,
147 temperature_or_pressure: TP,
148 liquid_molefracs: &OVector<D, N>,
149 tp_init: Option<TP::Other>,
150 vapor_molefracs: Option<&OVector<f64, N>>,
151 options: (SolverOptions, SolverOptions),
152 ) -> FeosResult<Self> {
153 Self::bubble_dew_point(
154 eos,
155 temperature_or_pressure,
156 liquid_molefracs,
157 tp_init,
158 vapor_molefracs,
159 true,
160 options,
161 )
162 }
163
164 pub fn dew_point<TP: TemperatureOrPressure<D>>(
167 eos: &E,
168 temperature_or_pressure: TP,
169 vapor_molefracs: &OVector<D, N>,
170 tp_init: Option<TP::Other>,
171 liquid_molefracs: Option<&OVector<f64, N>>,
172 options: (SolverOptions, SolverOptions),
173 ) -> FeosResult<Self> {
174 Self::bubble_dew_point(
175 eos,
176 temperature_or_pressure,
177 vapor_molefracs,
178 tp_init,
179 liquid_molefracs,
180 false,
181 options,
182 )
183 }
184
185 pub(super) fn bubble_dew_point<TP: TemperatureOrPressure<D>>(
186 eos: &E,
187 temperature_or_pressure: TP,
188 vapor_molefracs: &OVector<D, N>,
189 tp_init: Option<TP::Other>,
190 liquid_molefracs: Option<&OVector<f64, N>>,
191 bubble: bool,
192 options: (SolverOptions, SolverOptions),
193 ) -> FeosResult<Self> {
194 let (temperature, pressure, iterate_p) =
195 temperature_or_pressure.temperature_pressure(tp_init);
196 Self::bubble_dew_point_tp(
197 eos,
198 temperature,
199 pressure,
200 vapor_molefracs,
201 liquid_molefracs,
202 bubble,
203 iterate_p,
204 options,
205 )
206 }
207
208 #[expect(clippy::too_many_arguments)]
209 fn bubble_dew_point_tp(
210 eos: &E,
211 temperature: Option<Temperature<D>>,
212 pressure: Option<Pressure<D>>,
213 molefracs_spec: &OVector<D, N>,
214 molefracs_init: Option<&OVector<f64, N>>,
215 bubble: bool,
216 iterate_p: bool,
217 options: (SolverOptions, SolverOptions),
218 ) -> FeosResult<Self> {
219 let eos_re = eos.re();
220 let mut temperature_re = temperature.map(|t| t.re());
221 let mut pressure_re = pressure.map(|p| p.re());
222 let molefracs_spec_re = molefracs_spec.map(|x| x.re());
223 let (v1, rho2) = if iterate_p {
224 let temperature_re = temperature_re.as_mut().unwrap();
226
227 if let Some(p) = pressure_re.as_mut() {
229 PhaseEquilibrium::iterate_bubble_dew(
230 &eos_re,
231 temperature_re,
232 p,
233 &molefracs_spec_re,
234 molefracs_init,
235 bubble,
236 iterate_p,
237 options,
238 )?
239 } else {
240 let x2 = PhaseEquilibrium::starting_pressure_ideal_gas(
242 &eos_re,
243 *temperature_re,
244 &molefracs_spec_re,
245 bubble,
246 )
247 .and_then(|(p, x)| {
248 pressure_re = Some(p);
249 PhaseEquilibrium::iterate_bubble_dew(
250 &eos_re,
251 temperature_re,
252 pressure_re.as_mut().unwrap(),
253 &molefracs_spec_re,
254 molefracs_init.or(Some(&x)),
255 bubble,
256 iterate_p,
257 options,
258 )
259 });
260
261 x2.or_else(|_| {
263 PhaseEquilibrium::starting_pressure_spinodal(
264 &eos_re,
265 *temperature_re,
266 &molefracs_spec_re,
267 )
268 .and_then(|p| {
269 pressure_re = Some(p);
270 PhaseEquilibrium::iterate_bubble_dew(
271 &eos_re,
272 temperature_re,
273 pressure_re.as_mut().unwrap(),
274 &molefracs_spec_re,
275 molefracs_init,
276 bubble,
277 iterate_p,
278 options,
279 )
280 })
281 })?
282 }
283 } else {
284 let pressure_re = pressure_re.as_mut().unwrap();
286
287 let temperature_re = temperature_re.as_mut().expect("An initial temperature is required for the calculation of bubble/dew points at given pressure!");
288 PhaseEquilibrium::iterate_bubble_dew(
289 &eos.re(),
290 temperature_re,
291 pressure_re,
292 &molefracs_spec_re,
293 molefracs_init,
294 bubble,
295 iterate_p,
296 options,
297 )?
298 };
299
300 let (mut t, mut p) = if iterate_p {
302 (
303 temperature.unwrap().into_reduced(),
304 D::from(pressure_re.unwrap().into_reduced()),
305 )
306 } else {
307 (
308 D::from(temperature_re.unwrap().into_reduced()),
309 pressure.unwrap().into_reduced(),
310 )
311 };
312 let mut molar_volume = D::from(v1);
313 let mut rho2 = rho2.map(D::from);
314 for _ in 0..D::NDERIV {
315 if iterate_p {
316 Self::newton_step_t(
317 eos,
318 t,
319 molefracs_spec,
320 &mut p,
321 &mut molar_volume,
322 &mut rho2,
323 Verbosity::None,
324 )
325 } else {
326 Self::newton_step_p(
327 eos,
328 &mut t,
329 molefracs_spec,
330 p,
331 &mut molar_volume,
332 &mut rho2,
333 Verbosity::None,
334 )
335 };
336 }
337 let state1 = State::new_intensive(
338 eos,
339 Temperature::from_reduced(t),
340 Density::from_reduced(molar_volume.recip()),
341 molefracs_spec,
342 )?;
343 let rho2_total = rho2.sum();
344 let x2 = rho2 / rho2_total;
345 let state2 = State::new_intensive(
346 eos,
347 Temperature::from_reduced(t),
348 Density::from_reduced(rho2_total),
349 &x2,
350 )?;
351
352 Ok(PhaseEquilibrium(if bubble {
353 [state2, state1]
354 } else {
355 [state1, state2]
356 }))
357 }
358
359 fn newton_step_t(
360 eos: &E,
361 temperature: D,
362 molefracs: &OVector<D, N>,
363 pressure: &mut D,
364 molar_volume: &mut D,
365 partial_density_other_phase: &mut OVector<D, N>,
366 verbosity: Verbosity,
367 ) -> f64 {
368 let (p_1, mu_res_1, dp_1, dmu_1) = eos.dmu_drho(temperature, partial_density_other_phase);
370 let (p_2, mu_res_2, dp_2, dmu_2) = eos.dmu_dv(temperature, *molar_volume, molefracs);
371
372 let n = molefracs.len();
374 let f = DVector::from_fn(n + 2, |i, _| {
375 if i == n {
376 p_1 - *pressure
377 } else if i == n + 1 {
378 p_2 - *pressure
379 } else {
380 mu_res_1[i] - mu_res_2[i]
381 + (partial_density_other_phase[i] * *molar_volume / molefracs[i]).ln()
382 * temperature
383 }
384 });
385
386 let jac = DMatrix::from_fn(n + 2, n + 2, |i, j| {
388 if i < n && j < n {
389 dmu_1[(i, j)]
390 } else if i < n && j == n {
391 -dmu_2[i]
392 } else if i == n && j < n {
393 dp_1[j]
394 } else if i == n + 1 && j == n {
395 dp_2
396 } else if i >= n && j == n + 1 {
397 -D::one()
398 } else {
399 D::zero()
400 }
401 });
402
403 let dx = LU::<_, _, Dyn>::new(jac).unwrap().solve(&f);
405
406 for i in 0..n {
408 partial_density_other_phase[i] -= dx[i];
409 }
410 *molar_volume -= dx[n];
411 *pressure -= dx[n + 1];
412
413 let error = f.map(|r| r.re()).norm();
414
415 let x = partial_density_other_phase.map(|r| r.re());
416 let x = &x / x.sum();
417 log_iteration(
418 verbosity,
419 Some(error),
420 Temperature::from_reduced(temperature.re()),
421 Pressure::from_reduced(pressure.re()),
422 x.as_slice(),
423 true,
424 );
425 error
426 }
427
428 fn newton_step_p(
429 eos: &E,
430 temperature: &mut D,
431 molefracs: &OVector<D, N>,
432 pressure: D,
433 molar_volume: &mut D,
434 partial_density_other_phase: &mut OVector<D, N>,
435 verbosity: Verbosity,
436 ) -> f64 {
437 let (p_1, mu_res_1, dp_1, dmu_1) = eos.dmu_drho(*temperature, partial_density_other_phase);
439 let (p_2, mu_res_2, dp_2, dmu_2) = eos.dmu_dv(*temperature, *molar_volume, molefracs);
440 let (dp_dt_1, dmu_res_dt_1) = eos.dmu_dt(*temperature, partial_density_other_phase);
441 let (dp_dt_2, dmu_res_dt_2) = eos.dmu_dt(*temperature, &(molefracs / *molar_volume));
442
443 let n = molefracs.len();
445 let f = DVector::from_fn(n + 2, |i, _| {
446 if i == n {
447 p_1 - pressure
448 } else if i == n + 1 {
449 p_2 - pressure
450 } else {
451 mu_res_1[i] - mu_res_2[i]
452 + (partial_density_other_phase[i] * *molar_volume / molefracs[i]).ln()
453 * *temperature
454 }
455 });
456
457 let jac = DMatrix::from_fn(n + 2, n + 2, |i, j| {
459 if i < n && j < n {
460 dmu_1[(i, j)]
461 } else if i < n && j == n {
462 -dmu_2[i]
463 } else if i < n && j == n + 1 {
464 dmu_res_dt_1[i] - dmu_res_dt_2[i]
465 + (partial_density_other_phase[i] * *molar_volume / molefracs[i]).ln()
466 } else if i == n && j < n {
467 dp_1[j]
468 } else if i == n && j == n + 1 {
469 dp_dt_1
470 } else if i == n + 1 && j == n {
471 dp_2
472 } else if i == n + 1 && j == n + 1 {
473 dp_dt_2
474 } else {
475 D::zero()
476 }
477 });
478
479 let dx = LU::<_, _, Dyn>::new(jac).unwrap().solve(&f);
481
482 for i in 0..n {
484 partial_density_other_phase[i] -= dx[i];
485 }
486 *molar_volume -= dx[n];
487 *temperature -= dx[n + 1];
488
489 let error = f.map(|r| r.re()).norm();
490
491 let x = partial_density_other_phase.map(|r| r.re());
492 let x = &x / x.sum();
493 log_iteration(
494 verbosity,
495 Some(error),
496 Temperature::from_reduced(temperature.re()),
497 Pressure::from_reduced(pressure.re()),
498 x.as_slice(),
499 true,
500 );
501 error
502 }
503}
504
505impl<E: Residual<N>, N: Gradients> PhaseEquilibrium<E, 2, N>
507where
508 DefaultAllocator: Allocator<N> + Allocator<N, N> + Allocator<U1, N>,
509{
510 #[expect(clippy::too_many_arguments)]
511 fn iterate_bubble_dew(
512 eos: &E,
513 temperature: &mut Temperature,
514 pressure: &mut Pressure,
515 molefracs_spec: &OVector<f64, N>,
516 molefracs_init: Option<&OVector<f64, N>>,
517 bubble: bool,
518 iterate_p: bool,
519 options: (SolverOptions, SolverOptions),
520 ) -> FeosResult<(f64, OVector<f64, N>)> {
521 let [mut state1, mut state2] = if bubble {
522 Self::starting_x2_bubble(eos, *temperature, *pressure, molefracs_spec, molefracs_init)
523 } else {
524 Self::starting_x2_dew(eos, *temperature, *pressure, molefracs_spec, molefracs_init)
525 }?;
526 let (options_inner, options_outer) = options;
527
528 let mut err_out = 1.0;
530 let mut k_out = 0;
531
532 if PhaseEquilibrium::is_trivial_solution(&state1, &state2) {
533 log_iter!(options_outer.verbosity, "Trivial solution encountered!");
534 return Err(FeosError::TrivialSolution);
535 }
536
537 log_iter!(
538 options_outer.verbosity,
539 "res outer loop | res inner loop | temperature | pressure | molefracs second phase",
540 );
541 log_iter!(options_outer.verbosity, "{:-<104}", "");
542 log_iteration(
543 options_outer.verbosity,
544 None,
545 *temperature,
546 *pressure,
547 state2.molefracs.as_slice(),
548 false,
549 );
550
551 for ko in 0..options_outer.max_iter.unwrap_or(MAX_ITER_OUTER) {
553 err_out = if err_out > NEWTON_TOL {
555 for _ in 0..options_inner.max_iter.unwrap_or(MAX_ITER_INNER) {
557 let res = if iterate_p {
558 Self::adjust_p(
559 *temperature,
560 pressure,
561 &mut state1,
562 &mut state2,
563 options_inner.verbosity,
564 )?
565 } else {
566 Self::adjust_t(
567 temperature,
568 *pressure,
569 &mut state1,
570 &mut state2,
571 options_inner.verbosity,
572 )?
573 };
574 if res < options_inner.tol.unwrap_or(TOL_INNER) {
575 break;
576 }
577 }
578 Self::adjust_x2(&state1, &mut state2, options_outer.verbosity)
579 } else {
580 let mut t = temperature.into_reduced();
581 let mut p = pressure.into_reduced();
582 let mut molar_volume = state1.density.into_reduced().recip();
583 let mut rho2 = state2.partial_density.to_reduced();
584 let err = if iterate_p {
585 Self::newton_step_t(
586 &state1.eos,
587 t,
588 &state1.molefracs,
589 &mut p,
590 &mut molar_volume,
591 &mut rho2,
592 options_outer.verbosity,
593 )
594 } else {
595 Self::newton_step_p(
596 &state1.eos,
597 &mut t,
598 &state1.molefracs,
599 p,
600 &mut molar_volume,
601 &mut rho2,
602 options_outer.verbosity,
603 )
604 };
605 *temperature = Temperature::from_reduced(t);
606 *pressure = Pressure::from_reduced(p);
607 state1.density = Density::from_reduced(molar_volume.recip());
608 state2.partial_density = Density::from_reduced(rho2);
609 Ok(err)
610 }?;
611
612 if Self::is_trivial_solution(&state1, &state2) {
613 log_iter!(options_outer.verbosity, "Trivial solution encountered!");
614 return Err(FeosError::TrivialSolution);
615 }
616
617 if err_out < options_outer.tol.unwrap_or(TOL_OUTER) {
618 k_out = ko + 1;
619 break;
620 }
621 }
622
623 if err_out < options_outer.tol.unwrap_or(TOL_OUTER) {
624 log_result!(
625 options_outer.verbosity,
626 "Bubble/dew point: calculation converged in {} step(s)\n",
627 k_out
628 );
629 Ok((
630 state1.density.into_reduced().recip(),
631 state2.partial_density.to_reduced(),
632 ))
633 } else {
634 Err(FeosError::NotConverged(String::from(
636 "bubble-dew-iteration",
637 )))
638 }
639 }
640
641 fn adjust_p(
642 temperature: Temperature,
643 pressure: &mut Pressure,
644 state1: &mut State<E, N>,
645 state2: &mut State<E, N>,
646 verbosity: Verbosity,
647 ) -> FeosResult<f64> {
648 let ln_phi_1 = state1.ln_phi();
650 let ln_phi_2 = state2.ln_phi();
651 let k = (&ln_phi_1 - &ln_phi_2).map(f64::exp);
652
653 let xk = state1.molefracs.component_mul(&k);
655 let f = xk.sum() - 1.0;
656
657 let ln_phi_1_dp = state1.dln_phi_dp();
659 let ln_phi_2_dp = state2.dln_phi_dp();
660 let df = ((ln_phi_1_dp - ln_phi_2_dp) * *pressure)
661 .into_value()
662 .component_mul(&xk)
663 .sum();
664 let mut lnpstep = -f / df;
665
666 lnpstep = lnpstep.clamp(-MAX_LNPSTEP, MAX_LNPSTEP);
668
669 *pressure *= lnpstep.exp();
671
672 Self::adjust_states(temperature, *pressure, state1, state2, None)?;
674
675 log_iteration(
677 verbosity,
678 Some(f),
679 temperature,
680 *pressure,
681 state2.molefracs.as_slice(),
682 false,
683 );
684
685 Ok(f.abs())
686 }
687
688 fn adjust_t(
689 temperature: &mut Temperature,
690 pressure: Pressure,
691 state1: &mut State<E, N>,
692 state2: &mut State<E, N>,
693 verbosity: Verbosity,
694 ) -> FeosResult<f64> {
695 let ln_phi_1 = state1.ln_phi();
697 let ln_phi_2 = state2.ln_phi();
698 let k = (&ln_phi_1 - &ln_phi_2).map(f64::exp);
699
700 let f = state1.molefracs.dot(&k) - 1.0;
702
703 let ln_phi_1_dt = state1.dln_phi_dt();
705 let ln_phi_2_dt = state2.dln_phi_dt();
706 let df = ((ln_phi_1_dt - ln_phi_2_dt)
707 .component_mul(&Dimensionless::new(state1.molefracs.component_mul(&k))))
708 .sum();
709 let mut tstep = -f / df;
710
711 if tstep < -Temperature::from_reduced(MAX_TSTEP) {
713 tstep = -Temperature::from_reduced(MAX_TSTEP);
714 } else if tstep > Temperature::from_reduced(MAX_TSTEP) {
715 tstep = Temperature::from_reduced(MAX_TSTEP);
716 }
717
718 *temperature += tstep;
720
721 Self::adjust_states(*temperature, pressure, state1, state2, None)?;
723
724 log_iteration(
726 verbosity,
727 Some(f),
728 *temperature,
729 pressure,
730 state2.molefracs.as_slice(),
731 false,
732 );
733
734 Ok(f.abs())
735 }
736
737 fn starting_pressure_ideal_gas(
738 eos: &E,
739 temperature: Temperature,
740 molefracs_spec: &OVector<f64, N>,
741 bubble: bool,
742 ) -> FeosResult<(Pressure, OVector<f64, N>)> {
743 if bubble {
744 Self::starting_pressure_ideal_gas_bubble(eos, temperature, molefracs_spec)
745 } else {
746 Self::starting_pressure_ideal_gas_dew(eos, temperature, molefracs_spec)
747 }
748 }
749
750 pub(super) fn starting_pressure_ideal_gas_bubble(
751 eos: &E,
752 temperature: Temperature,
753 liquid_molefracs: &OVector<f64, N>,
754 ) -> FeosResult<(Pressure, OVector<f64, N>)> {
755 let density = 0.75 * Density::from_reduced(eos.compute_max_density(liquid_molefracs));
756 let liquid = State::new_intensive(eos, temperature, density, liquid_molefracs)?;
757 let v_l = liquid.partial_molar_volume();
758 let p_l = liquid.pressure(Contributions::Total);
759 let mu_l = liquid.residual_chemical_potential();
760 let k_i = liquid_molefracs.component_mul(
761 &((mu_l - v_l * p_l) / (RGAS * temperature))
762 .into_value()
763 .map(f64::exp),
764 );
765 let p = k_i.sum() * RGAS * temperature * density;
766 let y = &k_i / k_i.sum();
767 Ok((p, y))
768 }
769
770 fn starting_pressure_ideal_gas_dew(
771 eos: &E,
772 temperature: Temperature,
773 vapor_molefracs: &OVector<f64, N>,
774 ) -> FeosResult<(Pressure, OVector<f64, N>)> {
775 let mut p: Option<Pressure> = None;
776
777 let mut x = vapor_molefracs.clone();
778 for _ in 0..5 {
779 let density = Density::from_reduced(0.75 * eos.compute_max_density(&x));
780 let liquid = State::new_intensive(eos, temperature, density, &x)?;
781 let v_l = liquid.partial_molar_volume();
782 let p_l = liquid.pressure(Contributions::Total);
783 let mu_l = liquid.residual_chemical_potential();
784 let k = vapor_molefracs.clone().component_div(
785 &((mu_l - v_l * p_l) / (RGAS * temperature))
786 .into_value()
787 .map(f64::exp),
788 );
789 let k_sum = k.sum();
790 let p_new = RGAS * temperature * density / k_sum;
791 x = k / k_sum;
792 if let Some(p_old) = p
793 && ((p_new - p_old) / p_old).into_value().abs() < 1e-5
794 {
795 p = Some(p_new);
796 break;
797 }
798 p = Some(p_new);
799 }
800 Ok((p.unwrap(), x))
801 }
802
803 pub(super) fn starting_pressure_spinodal(
804 eos: &E,
805 temperature: Temperature,
806 molefracs: &OVector<f64, N>,
807 ) -> FeosResult<Pressure> {
808 let [sp_v, sp_l] = State::spinodal(eos, temperature, Some(molefracs), Default::default())?;
809 let pv = sp_v.pressure(Contributions::Total);
810 let pl = sp_l.pressure(Contributions::Total);
811 Ok(0.5 * (Pressure::from_reduced(0.0).max(pl) + pv))
812 }
813
814 fn starting_x2_bubble(
815 eos: &E,
816 temperature: Temperature,
817 pressure: Pressure,
818 liquid_molefracs: &OVector<f64, N>,
819 vapor_molefracs: Option<&OVector<f64, N>>,
820 ) -> FeosResult<[State<E, N>; 2]> {
821 let liquid_state =
822 State::new_xpt(eos, temperature, pressure, liquid_molefracs, Some(Liquid))?;
823 let xv = match vapor_molefracs {
824 Some(xv) => xv.clone(),
825 None => liquid_state
826 .ln_phi()
827 .map(f64::exp)
828 .component_mul(liquid_molefracs),
829 };
830 let vapor_state = State::new_xpt(eos, temperature, pressure, &xv, Some(Vapor))?;
831 Ok([liquid_state, vapor_state])
832 }
833
834 fn starting_x2_dew(
835 eos: &E,
836 temperature: Temperature,
837 pressure: Pressure,
838 vapor_molefracs: &OVector<f64, N>,
839 liquid_molefracs: Option<&OVector<f64, N>>,
840 ) -> FeosResult<[State<E, N>; 2]> {
841 let vapor_state = State::new_npt(
842 eos,
843 temperature,
844 pressure,
845 &Moles::from_reduced(vapor_molefracs.clone()),
846 Some(Vapor),
847 )?;
848 let xl = match liquid_molefracs {
849 Some(xl) => xl.clone(),
850 None => {
851 let xl = vapor_state
852 .ln_phi()
853 .map(f64::exp)
854 .component_mul(vapor_molefracs);
855 let liquid_state = State::new_xpt(eos, temperature, pressure, &xl, Some(Liquid))?;
856 (vapor_state.ln_phi() - liquid_state.ln_phi())
857 .map(f64::exp)
858 .component_mul(vapor_molefracs)
859 }
860 };
861 let liquid_state = State::new_xpt(eos, temperature, pressure, &xl, Some(Liquid))?;
862 Ok([vapor_state, liquid_state])
863 }
864
865 fn adjust_states(
866 temperature: Temperature,
867 pressure: Pressure,
868 state1: &mut State<E, N>,
869 state2: &mut State<E, N>,
870 moles_state2: Option<&Moles<OVector<f64, N>>>,
871 ) -> FeosResult<()> {
872 *state1 = State::new_npt(
873 &state1.eos,
874 temperature,
875 pressure,
876 &state1.moles,
877 Some(InitialDensity(state1.density)),
878 )?;
879 *state2 = State::new_npt(
880 &state2.eos,
881 temperature,
882 pressure,
883 moles_state2.unwrap_or(&state2.moles),
884 Some(InitialDensity(state2.density)),
885 )?;
886 Ok(())
887 }
888
889 fn adjust_x2(
890 state1: &State<E, N>,
891 state2: &mut State<E, N>,
892 verbosity: Verbosity,
893 ) -> FeosResult<f64> {
894 let x1 = &state1.molefracs;
895 let ln_phi_1 = state1.ln_phi();
896 let ln_phi_2 = state2.ln_phi();
897 let k = (ln_phi_1 - ln_phi_2).map(f64::exp);
898 let kx1 = k.component_mul(x1);
899 let err_out = kx1
900 .component_div(&state2.molefracs)
901 .map(|e| (e - 1.0).abs())
902 .sum();
903 let x2 = &kx1 / kx1.sum();
904 log_iter!(
905 verbosity,
906 "{:<14.8e} | {:14} | {:14} | {:16} |",
907 err_out,
908 "",
909 "",
910 ""
911 );
912 *state2 = State::new_xpt(
913 &state2.eos,
914 state2.temperature,
915 state2.pressure(Contributions::Total),
916 &x2,
917 Some(InitialDensity(state2.density)),
918 )?;
919 Ok(err_out)
920 }
921}
922
923fn log_iteration(
924 verbosity: Verbosity,
925 error: Option<f64>,
926 temperature: Temperature,
927 pressure: Pressure,
928 x2: &[f64],
929 newton: bool,
930) {
931 let error = error.map_or_else(|| format!("{:14}", ""), |e| format!("{:<14.8e}", e.abs()));
932 log_iter!(
933 verbosity,
934 "{:14} | {} | {:12.8} | {:12.8} | {:.8?} {}",
935 "",
936 error,
937 temperature,
938 pressure,
939 x2,
940 if newton { "NEWTON" } else { "" }
941 );
942}