feos-pets 0.1.0

Implementation of PeTS equation of state and corresponding Helmholtz energy functional.
Documentation
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
use crate::parameters::PetsParameters;
use feos_core::joback::Joback;
use feos_core::parameter::Parameter;
use feos_core::{
    Contributions, EntropyScaling, EosError, EosResult, EosUnit, EquationOfState, HelmholtzEnergy,
    IdealGasContribution, MolarWeight, State,
};
use ndarray::Array1;
use quantity::si::*;
use std::f64::consts::{FRAC_PI_6, PI};
use std::rc::Rc;

pub(crate) mod dispersion;
pub(crate) mod hard_sphere;
mod qspr;
use dispersion::Dispersion;
use hard_sphere::HardSphere;
use qspr::QSPR;

#[allow(clippy::upper_case_acronyms)]
enum IdealGasContributions {
    QSPR(QSPR),
    Joback(Joback),
}

#[derive(Copy, Clone)]
pub struct PetsOptions {
    pub max_eta: f64,
}

impl Default for PetsOptions {
    fn default() -> Self {
        Self { max_eta: 0.5 }
    }
}

pub struct Pets {
    parameters: Rc<PetsParameters>,
    options: PetsOptions,
    contributions: Vec<Box<dyn HelmholtzEnergy>>,
    ideal_gas: IdealGasContributions,
}

impl Pets {
    pub fn new(parameters: Rc<PetsParameters>) -> Self {
        Self::with_options(parameters, PetsOptions::default())
    }

    pub fn with_options(parameters: Rc<PetsParameters>, options: PetsOptions) -> Self {
        let mut contributions: Vec<Box<dyn HelmholtzEnergy>> = Vec::with_capacity(2);
        contributions.push(Box::new(HardSphere {
            parameters: parameters.clone(),
        }));
        contributions.push(Box::new(Dispersion {
            parameters: parameters.clone(),
        }));

        let joback_records = parameters.joback_records.clone();

        Self {
            parameters: parameters.clone(),
            options,
            contributions,
            ideal_gas: joback_records.map_or(
                IdealGasContributions::QSPR(QSPR { parameters }),
                |joback_records| IdealGasContributions::Joback(Joback::new(joback_records)),
            ),
        }
    }
}

impl EquationOfState for Pets {
    fn components(&self) -> usize {
        self.parameters.pure_records.len()
    }

    fn subset(&self, component_list: &[usize]) -> Self {
        Self::with_options(
            Rc::new(self.parameters.subset(component_list)),
            self.options,
        )
    }

    fn compute_max_density(&self, moles: &Array1<f64>) -> f64 {
        self.options.max_eta * moles.sum()
            / (FRAC_PI_6 * self.parameters.sigma.mapv(|v| v.powi(3)) * moles).sum()
    }

    fn residual(&self) -> &[Box<dyn HelmholtzEnergy>] {
        &self.contributions
    }

    fn ideal_gas(&self) -> &dyn IdealGasContribution {
        match &self.ideal_gas {
            IdealGasContributions::QSPR(qspr) => qspr,
            IdealGasContributions::Joback(joback) => joback,
        }
    }
}

impl MolarWeight<SIUnit> for Pets {
    fn molar_weight(&self) -> SIArray1 {
        self.parameters.molarweight.clone() * GRAM / MOL
    }
}

fn omega11(t: f64) -> f64 {
    1.06036 * t.powf(-0.15610)
        + 0.19300 * (-0.47635 * t).exp()
        + 1.03587 * (-1.52996 * t).exp()
        + 1.76474 * (-3.89411 * t).exp()
}

fn omega22(t: f64) -> f64 {
    1.16145 * t.powf(-0.14874) + 0.52487 * (-0.77320 * t).exp() + 2.16178 * (-2.43787 * t).exp()
        - 6.435e-4 * t.powf(0.14874) * (18.0323 * t.powf(-0.76830) - 7.27371).sin()
}

impl EntropyScaling<SIUnit> for Pets {
    fn viscosity_reference(
        &self,
        temperature: SINumber,
        _: SINumber,
        moles: &SIArray1,
    ) -> EosResult<SINumber> {
        let x = moles.to_reduced(moles.sum())?;
        let p = &self.parameters;
        let mw = &p.molarweight;
        let ce: Array1<SINumber> = (0..self.components())
            .map(|i| {
                let tr = (temperature / p.epsilon_k[i] / KELVIN)
                    .into_value()
                    .unwrap();
                5.0 / 16.0
                    * (mw[i] * GRAM / MOL * KB / NAV * temperature / PI)
                        .sqrt()
                        .unwrap()
                    / omega22(tr)
                    / (p.sigma[i] * ANGSTROM).powi(2)
            })
            .collect();
        let mut ce_mix = 0.0 * MILLI * PASCAL * SECOND;
        for i in 0..self.components() {
            let denom: f64 = (0..self.components())
                .map(|j| {
                    x[j] * (1.0
                        + (ce[i] / ce[j]).into_value().unwrap().sqrt()
                            * (mw[j] / mw[i]).powf(1.0 / 4.0))
                    .powi(2)
                        / (8.0 * (1.0 + mw[i] / mw[j])).sqrt()
                })
                .sum();
            ce_mix += ce[i] * x[i] / denom
        }
        Ok(ce_mix)
    }

    fn viscosity_correlation(&self, s_res: f64, x: &Array1<f64>) -> EosResult<f64> {
        let coefficients = self
            .parameters
            .viscosity
            .as_ref()
            .expect("Missing viscosity coefficients.");
        let a: f64 = (&coefficients.row(0) * x).sum();
        let b: f64 = (&coefficients.row(1) * x).sum();
        let c: f64 = (&coefficients.row(2) * x).sum();
        let d: f64 = (&coefficients.row(3) * x).sum();
        Ok(a + b * s_res + c * s_res.powi(2) + d * s_res.powi(3))
    }

    fn diffusion_reference(
        &self,
        temperature: SINumber,
        volume: SINumber,
        moles: &SIArray1,
    ) -> EosResult<SINumber> {
        if self.components() != 1 {
            return Err(EosError::IncompatibleComponents(self.components(), 1));
        }
        let p = &self.parameters;
        let density = moles.sum() / volume;
        let res: Array1<SINumber> = (0..self.components())
            .map(|i| {
                let tr = (temperature / p.epsilon_k[i] / KELVIN)
                    .into_value()
                    .unwrap();
                3.0 / 8.0 / (p.sigma[i] * ANGSTROM).powi(2) / omega11(tr) / (density * NAV)
                    * (temperature * RGAS / PI / (p.molarweight[i] * GRAM / MOL))
                        .sqrt()
                        .unwrap()
            })
            .collect();
        Ok(res[0])
    }

    fn diffusion_correlation(&self, s_res: f64, x: &Array1<f64>) -> EosResult<f64> {
        if self.components() != 1 {
            return Err(EosError::IncompatibleComponents(self.components(), 1));
        }
        let coefficients = self
            .parameters
            .diffusion
            .as_ref()
            .expect("Missing diffusion coefficients.");
        let a: f64 = (&coefficients.row(0) * x).sum();
        let b: f64 = (&coefficients.row(1) * x).sum();
        let c: f64 = (&coefficients.row(2) * x).sum();
        let d: f64 = (&coefficients.row(3) * x).sum();
        let e: f64 = (&coefficients.row(4) * x).sum();
        Ok(a + b * s_res
            - c * (1.0 - s_res.exp()) * s_res.powi(2)
            - d * s_res.powi(4)
            - e * s_res.powi(8))
    }

    // fn thermal_conductivity_reference(
    //     &self,
    //     state: &State<SIUnit, E>,
    // ) -> EosResult<SINumber> {
    //     if self.components() != 1 {
    //         return Err(EosError::IncompatibleComponents(self.components(), 1));
    //     }
    //     let p = &self.parameters;
    //     let res: Array1<SINumber> = (0..self.components())
    //         .map(|i| {
    //             let tr = (state.temperature / p.epsilon_k[i] / KELVIN)
    //                 .into_value()
    //                 .unwrap();
    //             let cp = State::critical_point_pure(&state.eos, Some(state.temperature)).unwrap();
    //             let s_res_cp_reduced = cp
    //                 .entropy(Contributions::Residual)
    //                 .to_reduced(SIUnit::reference_entropy())
    //                 .unwrap();
    //             let s_res_reduced = cp
    //                 .entropy(Contributions::Residual)
    //                 .to_reduced(SIUnit::reference_entropy())
    //                 .unwrap();
    //             let ref_ce = 0.083235
    //                 * ((state.temperature / KELVIN).into_value().unwrap()
    //                     / (p.molarweight[0]))
    //                     .sqrt()
    //                 / p.sigma[0]
    //                 / p.sigma[0]
    //                 / omega22(tr);
    //             let alpha_visc = (-s_res_reduced / s_res_cp_reduced).exp();
    //             let ref_ts = (-0.0167141 * tr + 0.0470581 * (tr).powi(2))
    //                 * (p.sigma[i].powi(3) * p.epsilon_k[0])
    //                 / 100000.0;
    //             (ref_ce + ref_ts * alpha_visc) * WATT / METER / KELVIN
    //         })
    //         .collect();
    //     Ok(res[0])
    // }

    // Equation 11 of DOI: 10.1021/acs.iecr.9b03998
    fn thermal_conductivity_reference(
        &self,
        temperature: SINumber,
        volume: SINumber,
        moles: &SIArray1,
    ) -> EosResult<SINumber> {
        if self.components() != 1 {
            return Err(EosError::IncompatibleComponents(self.components(), 1));
        }
        let p = &self.parameters;
        let state = State::new_nvt(
            &Rc::new(Self::new(self.parameters.clone())),
            temperature,
            volume,
            moles,
        )?;
        let res: Array1<SINumber> = (0..self.components())
            .map(|i| {
                let tr = (temperature / p.epsilon_k[i] / KELVIN)
                    .into_value()
                    .unwrap();
                let ce = 83.235
                    * f64::powf(10.0, -1.5)
                    * ((temperature / KELVIN).into_value().unwrap() / p.molarweight[0]).sqrt()
                    / (p.sigma[0] * p.sigma[0])
                    / omega22(tr);
                ce * WATT / METER / KELVIN
                    + state.density
                        * self
                            .diffusion_reference(temperature, volume, moles)
                            .unwrap()
                        * self
                            .diffusion_correlation(
                                state
                                    .molar_entropy(Contributions::ResidualNvt)
                                    .to_reduced(SIUnit::reference_molar_entropy())
                                    .unwrap(),
                                &state.molefracs,
                            )
                            .unwrap()
                        * (state.c_v(Contributions::Total) - 1.5 * RGAS)
            })
            .collect();
        Ok(res[0])
    }

    fn thermal_conductivity_correlation(&self, s_res: f64, x: &Array1<f64>) -> EosResult<f64> {
        if self.components() != 1 {
            return Err(EosError::IncompatibleComponents(self.components(), 1));
        }
        let coefficients = self
            .parameters
            .thermal_conductivity
            .as_ref()
            .expect("Missing thermal conductivity coefficients");
        let a: f64 = (&coefficients.row(0) * x).sum();
        let b: f64 = (&coefficients.row(1) * x).sum();
        let c: f64 = (&coefficients.row(2) * x).sum();
        let d: f64 = (&coefficients.row(3) * x).sum();
        Ok(a + b * s_res + c * (1.0 - s_res.exp()) + d * s_res.powi(2))
    }
}

#[cfg(test)]
mod tests {
    use super::*;
    use crate::parameters::utils::{
        argon_krypton_parameters, argon_parameters, krypton_parameters,
    };
    use approx::assert_relative_eq;
    use feos_core::{Contributions, DensityInitialization, PhaseEquilibrium, State};
    use ndarray::arr1;
    use quantity::si::{BAR, KELVIN, METER, PASCAL, RGAS};

    #[test]
    fn ideal_gas_pressure() {
        let e = Rc::new(Pets::new(argon_parameters()));
        let t = 200.0 * KELVIN;
        let v = 1e-3 * METER.powi(3);
        let n = arr1(&[1.0]) * MOL;
        let s = State::new_nvt(&e, t, v, &n).unwrap();
        let p_ig = s.total_moles * RGAS * t / v;
        assert_relative_eq!(s.pressure(Contributions::IdealGas), p_ig, epsilon = 1e-10);
        assert_relative_eq!(
            s.pressure(Contributions::IdealGas) + s.pressure(Contributions::ResidualNvt),
            s.pressure(Contributions::Total),
            epsilon = 1e-10
        );
    }

    #[test]
    fn ideal_gas_heat_capacity_joback() {
        let e = Rc::new(Pets::new(argon_parameters()));
        let t = 200.0 * KELVIN;
        let v = 1e-3 * METER.powi(3);
        let n = arr1(&[1.0]) * MOL;
        let s = State::new_nvt(&e, t, v, &n).unwrap();
        let p_ig = s.total_moles * RGAS * t / v;
        assert_relative_eq!(s.pressure(Contributions::IdealGas), p_ig, epsilon = 1e-10);
        assert_relative_eq!(
            s.pressure(Contributions::IdealGas) + s.pressure(Contributions::ResidualNvt),
            s.pressure(Contributions::Total),
            epsilon = 1e-10
        );
    }

    #[test]
    fn new_tpn() {
        let e = Rc::new(Pets::new(argon_parameters()));
        let t = 300.0 * KELVIN;
        let p = BAR;
        let m = arr1(&[1.0]) * MOL;
        let s = State::new_npt(&e, t, p, &m, DensityInitialization::None);
        let p_calc = if let Ok(state) = s {
            state.pressure(Contributions::Total)
        } else {
            0.0 * PASCAL
        };
        assert_relative_eq!(p, p_calc, epsilon = 1e-6);
    }

    #[test]
    fn vle_pure_t() {
        let e = Rc::new(Pets::new(argon_parameters()));
        let t = 300.0 * KELVIN;
        let vle = PhaseEquilibrium::pure(&e, t, None, Default::default());
        if let Ok(v) = vle {
            assert_relative_eq!(
                v.vapor().pressure(Contributions::Total),
                v.liquid().pressure(Contributions::Total),
                epsilon = 1e-6
            )
        }
    }

    #[test]
    fn mix_single() {
        let e1 = Rc::new(Pets::new(argon_parameters()));
        let e2 = Rc::new(Pets::new(krypton_parameters()));
        let e12 = Rc::new(Pets::new(argon_krypton_parameters()));
        let t = 300.0 * KELVIN;
        let v = 0.02456883872966545 * METER.powi(3);
        let m1 = arr1(&[2.0]) * MOL;
        let m1m = arr1(&[2.0, 0.0]) * MOL;
        let m2m = arr1(&[0.0, 2.0]) * MOL;
        let s1 = State::new_nvt(&e1, t, v, &m1).unwrap();
        let s2 = State::new_nvt(&e2, t, v, &m1).unwrap();
        let s1m = State::new_nvt(&e12, t, v, &m1m).unwrap();
        let s2m = State::new_nvt(&e12, t, v, &m2m).unwrap();
        assert_relative_eq!(
            s1.pressure(Contributions::Total),
            s1m.pressure(Contributions::Total),
            epsilon = 1e-12
        );
        assert_relative_eq!(
            s2.pressure(Contributions::Total),
            s2m.pressure(Contributions::Total),
            epsilon = 1e-12
        )
    }
}