Skip to content

polykin.thermo.eos¤

Virial ¤

Bases: GasEoS

Virial equation of state truncated to the second coefficient.

This EoS is based on the following \(Z(T,P)\) relationship:

\[ Z = 1 + \frac{B_m P}{R T} \]

where \(P\) is the pressure, \(T\) is the temperature, and \(B_m\) is the mixture second virial coefficient. The parameter \(B_m=B_m(T,y)\) is estimated based on the critical properties of the pure components and the mixture composition \(y\).

Important

This equation is an improvement over the ideal gas model, but it should only be used up to moderate pressures such that \(v/v_c > 2\).

References

  • RC Reid, JM Prausniz, and BE Poling. The properties of gases & liquids 4th edition, 1986, p. 37, 40, 80, 82, 145.
PARAMETER DESCRIPTION
Tc

Critical temperatures of all components [K].

TYPE: float | FloatVectorLike(N)

Pc

Critical pressures of all components [Pa].

TYPE: float | FloatVectorLike(N)

Zc

Critical compressibility factors of all components.

TYPE: float | FloatVectorLike(N)

w

Acentric factors of all components.

TYPE: float | FloatVectorLike(N)

name

Name.

TYPE: str DEFAULT: ''

Source code in src/polykin/thermo/eos/virial.py
 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
class Virial(GasEoS):
    r"""Virial equation of state truncated to the second coefficient.

    This EoS is based on the following $Z(T,P)$ relationship:

    $$ Z = 1 + \frac{B_m P}{R T} $$

    where $P$ is the pressure, $T$ is the temperature, and $B_m$ is the mixture
    second virial coefficient. The parameter $B_m=B_m(T,y)$ is estimated based
    on the critical properties of the pure components and the mixture composition
    $y$.

    !!! important

        This equation is an improvement over the ideal gas model, but it
        should only be used up to moderate pressures such that $v/v_c > 2$.

    **References**

    *   RC Reid, JM Prausniz, and BE Poling. The properties of gases &
        liquids 4th edition, 1986, p. 37, 40, 80, 82, 145.

    Parameters
    ----------
    Tc : float | FloatVectorLike (N)
        Critical temperatures of all components [K].
    Pc : float | FloatVectorLike (N)
        Critical pressures of all components [Pa].
    Zc : float | FloatVectorLike (N)
        Critical compressibility factors of all components.
    w : float | FloatVectorLike (N)
        Acentric factors of all components.
    name : str
        Name.
    """

    Tc: FloatVector
    Pc: FloatVector
    Zc: FloatVector
    w: FloatVector

    def __init__(
        self,
        Tc: float | FloatVectorLike,
        Pc: float | FloatVectorLike,
        Zc: float | FloatVectorLike,
        w: float | FloatVectorLike,
        name: str = "",
    ) -> None:

        Tc, Pc, Zc, w = convert_FloatOrVectorLike_to_FloatVector([Tc, Pc, Zc, w])

        N = len(Tc)
        super().__init__(N, name)
        self.Tc = Tc
        self.Pc = Pc
        self.Zc = Zc
        self.w = w

    @functools.cache
    def Bij(self, T: float) -> FloatSquareMatrix:
        r"""Calculate the matrix of interaction virial coefficients.

        The calculation is handled by [`B_mixture`](../../properties/pvt/B_mixture.md).

        Parameters
        ----------
        T : float
            Temperature [K].

        Returns
        -------
        FloatSquareMatrix (N,N)
            Matrix of interaction virial coefficients, $B_{ij}$ [m³/mol].
        """
        return B_mixture(T, self.Tc, self.Pc, self.Zc, self.w)

    def Bm(self, T: float, y: FloatVector) -> float:
        r"""Calculate the second virial coefficient of the mixture.

        $$ B_m = \sum_i \sum_j y_i y_j B_{ij} $$

        **References**

        *   RC Reid, JM Prausniz, and BE Poling. The properties of gases &
            liquids 4th edition, 1986, p. 79.

        Parameters
        ----------
        T : float
            Temperature [K].
        y : FloatVector (N)
            Mole fractions of all components [mol/mol].

        Returns
        -------
        float
            Mixture second virial coefficient, $B_m$ [m³/mol].
        """
        return quadratic_mixing(y, self.Bij(T))

    def Z(self, T: float, P: float, y: FloatVector) -> float:
        r"""Calculate the compressibility factor of the fluid.

        $$ Z = 1 + \frac{B_m P}{R T} $$

        Parameters
        ----------
        T : float
            Temperature [K].
        P : float
            Pressure [Pa].
        y : FloatVector (N)
            Mole fractions of all components [mol/mol].

        Returns
        -------
        float
            Compressibility factor.
        """
        Bm = self.Bm(T, y)
        return 1.0 + Bm * P / (R * T)

    @override
    def gR(self, T: float, P: float, y: FloatVector) -> float:
        """Calculate the molar residual Gibbs energy of the fluid.

        $$ g^R = B_m P $$

        Parameters
        ----------
        T : float
            Temperature [K].
        P : float
            Pressure [Pa].
        y : FloatVector
            Mole fractions of all components [mol/mol].

        Returns
        -------
        float
            Molar residual Gibbs energy [J/mol].
        """
        return self.Bm(T, y) * P

    @override
    def P(self, T: float, v: float, y: FloatVector) -> float:
        r"""Calculate the pressure of the fluid.

        $$ P = \frac{R T}{v - B_m} $$

        Parameters
        ----------
        T : float
            Temperature [K].
        v : float
            Molar volume [m³/mol].
        y : FloatVector (N)
            Mole fractions of all components [mol/mol].

        Returns
        -------
        float
            Pressure [Pa].
        """
        Bm = self.Bm(T, y)
        return R * T / (v - Bm)

    @override
    def phi(self, T: float, P: float, y: FloatVector) -> FloatVector:
        r"""Calculate the fugacity coefficients of all components.

        $$ \ln \hat{\phi}_i
           = \left(2\sum_j {y_jB_{ij}} -B_m \right)\frac{P}{RT} $$

        **References**

        *   RC Reid, JM Prausniz, and BE Poling. The properties of gases &
            liquids 4th edition, 1986, p. 145.

        Parameters
        ----------
        T : float
            Temperature [K].
        P : float
            Pressure [Pa].
        y : FloatVector (N)
            Mole fractions of all components [mol/mol].

        Returns
        -------
        FloatVector (N)
            Fugacity coefficients of all components.
        """
        B = self.Bij(T)
        Bm = self.Bm(T, y)
        return exp((2 * dot(B, y) - Bm) * P / (R * T))

    @override
    def kappa(self, T: float, P: float, y: FloatVector) -> float:
        r"""Calculate the isothermal compressibility coefficient.

        $$ \kappa \equiv
           - \frac{1}{v} \left( \frac{\partial v}{\partial P} \right)_T
           = \frac{1}{P Z} $$

        Parameters
        ----------
        T : float
            Temperature [K].
        P : float
            Pressure [Pa].
        y : FloatVector (N)
            Mole fractions of all components [mol/mol].

        Returns
        -------
        float
            Isothermal compressibility coefficient, $\kappa$ [Pa⁻¹].
        """
        Z = self.Z(T, P, y)
        return 1 / (P * Z)

Bij cached ¤

Bij(T: float) -> FloatSquareMatrix

Calculate the matrix of interaction virial coefficients.

The calculation is handled by B_mixture.

PARAMETER DESCRIPTION
T

Temperature [K].

TYPE: float

RETURNS DESCRIPTION
FloatSquareMatrix(N, N)

Matrix of interaction virial coefficients, \(B_{ij}\) [m³/mol].

Source code in src/polykin/thermo/eos/virial.py
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
@functools.cache
def Bij(self, T: float) -> FloatSquareMatrix:
    r"""Calculate the matrix of interaction virial coefficients.

    The calculation is handled by [`B_mixture`](../../properties/pvt/B_mixture.md).

    Parameters
    ----------
    T : float
        Temperature [K].

    Returns
    -------
    FloatSquareMatrix (N,N)
        Matrix of interaction virial coefficients, $B_{ij}$ [m³/mol].
    """
    return B_mixture(T, self.Tc, self.Pc, self.Zc, self.w)

Bm ¤

Bm(T: float, y: FloatVector) -> float

Calculate the second virial coefficient of the mixture.

\[ B_m = \sum_i \sum_j y_i y_j B_{ij} \]

References

  • RC Reid, JM Prausniz, and BE Poling. The properties of gases & liquids 4th edition, 1986, p. 79.
PARAMETER DESCRIPTION
T

Temperature [K].

TYPE: float

y

Mole fractions of all components [mol/mol].

TYPE: FloatVector(N)

RETURNS DESCRIPTION
float

Mixture second virial coefficient, \(B_m\) [m³/mol].

Source code in src/polykin/thermo/eos/virial.py
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
def Bm(self, T: float, y: FloatVector) -> float:
    r"""Calculate the second virial coefficient of the mixture.

    $$ B_m = \sum_i \sum_j y_i y_j B_{ij} $$

    **References**

    *   RC Reid, JM Prausniz, and BE Poling. The properties of gases &
        liquids 4th edition, 1986, p. 79.

    Parameters
    ----------
    T : float
        Temperature [K].
    y : FloatVector (N)
        Mole fractions of all components [mol/mol].

    Returns
    -------
    float
        Mixture second virial coefficient, $B_m$ [m³/mol].
    """
    return quadratic_mixing(y, self.Bij(T))

DA ¤

DA(T: float, V: float, n: FloatVector, v0: float) -> float

Calculate the departure of Helmholtz energy.

\[ A(T,V,n) - A^{\circ}(T,V,n)\]
PARAMETER DESCRIPTION
T

Temperature [K].

TYPE: float

V

Volume [m³].

TYPE: float

n

Mole amounts of all components [mol].

TYPE: FloatVector(N)

v0

Molar volume in reference state [m³/mol].

TYPE: float

RETURNS DESCRIPTION
float

Helmholtz energy departure, \(A - A^{\circ}\) [J].

Source code in src/polykin/thermo/eos/base.py
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
def DA(
    self,
    T: float,
    V: float,
    n: FloatVector,
    v0: float,
) -> float:
    r"""Calculate the departure of Helmholtz energy.

    $$ A(T,V,n) - A^{\circ}(T,V,n)$$

    Parameters
    ----------
    T : float
        Temperature [K].
    V : float
        Volume [m³].
    n : FloatVector (N)
        Mole amounts of all components [mol].
    v0 : float
        Molar volume in reference state [m³/mol].

    Returns
    -------
    float
        Helmholtz energy departure, $A - A^{\circ}$ [J].
    """

N property ¤

N: int

Number of components.

P ¤

P(T: float, v: float, y: FloatVector) -> float

Calculate the pressure of the fluid.

\[ P = \frac{R T}{v - B_m} \]
PARAMETER DESCRIPTION
T

Temperature [K].

TYPE: float

v

Molar volume [m³/mol].

TYPE: float

y

Mole fractions of all components [mol/mol].

TYPE: FloatVector(N)

RETURNS DESCRIPTION
float

Pressure [Pa].

Source code in src/polykin/thermo/eos/virial.py
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
@override
def P(self, T: float, v: float, y: FloatVector) -> float:
    r"""Calculate the pressure of the fluid.

    $$ P = \frac{R T}{v - B_m} $$

    Parameters
    ----------
    T : float
        Temperature [K].
    v : float
        Molar volume [m³/mol].
    y : FloatVector (N)
        Mole fractions of all components [mol/mol].

    Returns
    -------
    float
        Pressure [Pa].
    """
    Bm = self.Bm(T, y)
    return R * T / (v - Bm)

Z ¤

Z(T: float, P: float, y: FloatVector) -> float

Calculate the compressibility factor of the fluid.

\[ Z = 1 + \frac{B_m P}{R T} \]
PARAMETER DESCRIPTION
T

Temperature [K].

TYPE: float

P

Pressure [Pa].

TYPE: float

y

Mole fractions of all components [mol/mol].

TYPE: FloatVector(N)

RETURNS DESCRIPTION
float

Compressibility factor.

Source code in src/polykin/thermo/eos/virial.py
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
def Z(self, T: float, P: float, y: FloatVector) -> float:
    r"""Calculate the compressibility factor of the fluid.

    $$ Z = 1 + \frac{B_m P}{R T} $$

    Parameters
    ----------
    T : float
        Temperature [K].
    P : float
        Pressure [Pa].
    y : FloatVector (N)
        Mole fractions of all components [mol/mol].

    Returns
    -------
    float
        Compressibility factor.
    """
    Bm = self.Bm(T, y)
    return 1.0 + Bm * P / (R * T)

aR ¤

aR(T: float, P: float, y: FloatVector) -> float

Calculate the molar residual Helmholtz energy of the fluid.

\[ a^R = g^R - R T (Z - 1) \]
PARAMETER DESCRIPTION
T

Temperature [K].

TYPE: float

P

Pressure [Pa].

TYPE: float

y

Mole fractions of all components [mol/mol].

TYPE: FloatVector

RETURNS DESCRIPTION
float

Molar residual Helmholtz energy [J/mol].

Source code in src/polykin/thermo/eos/base.py
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
def aR(self, T: float, P: float, y: FloatVector) -> float:
    r"""Calculate the molar residual Helmholtz energy of the fluid.

    $$ a^R = g^R - R T (Z - 1) $$

    Parameters
    ----------
    T : float
        Temperature [K].
    P : float
        Pressure [Pa].
    y : FloatVector
        Mole fractions of all components [mol/mol].

    Returns
    -------
    float
        Molar residual Helmholtz energy [J/mol].
    """
    return self.gR(T, P, y) - R * T * (self.Z(T, P, y) - 1.0)

beta ¤

beta(T: float, P: float, y: FloatVector) -> float

Calculate the thermal expansion coefficient.

\[ \beta \equiv \frac{1}{v} \left( \frac{\partial v}{\partial T} \right)_{P,y_i} = \frac{1}{T} + \frac{1}{Z} \left( \frac{\partial Z}{\partial T} \right)_{P,y_i} \]
PARAMETER DESCRIPTION
T

Temperature [K].

TYPE: float

P

Pressure [Pa].

TYPE: float

y

Mole fractions of all components [mol/mol].

TYPE: FloatVector(N)

RETURNS DESCRIPTION
float

Thermal expansion coefficient, \(\beta\) [K⁻¹].

Source code in src/polykin/thermo/eos/base.py
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
def beta(self, T: float, P: float, y: FloatVector) -> float:
    r"""Calculate the thermal expansion coefficient.

    $$ \beta \equiv
     \frac{1}{v} \left( \frac{\partial v}{\partial T} \right)_{P,y_i}
     = \frac{1}{T}
     + \frac{1}{Z} \left( \frac{\partial Z}{\partial T} \right)_{P,y_i} $$

    Parameters
    ----------
    T : float
        Temperature [K].
    P : float
        Pressure [Pa].
    y : FloatVector (N)
        Mole fractions of all components [mol/mol].

    Returns
    -------
    float
        Thermal expansion coefficient, $\beta$ [K⁻¹].
    """
    dZdT, Z = derivative_centered(lambda T_: self.Z(T_, P, y), T)

    return 1 / T + dZdT / Z

f ¤

f(T: float, P: float, y: FloatVector) -> FloatVector

Calculate the fugacity of all components.

\[ \hat{f}_i = \hat{\phi}_i y_i P \]
PARAMETER DESCRIPTION
T

Temperature [K].

TYPE: float

P

Pressure [Pa].

TYPE: float

y

Mole fractions of all components [mol/mol].

TYPE: FloatVector(N)

RETURNS DESCRIPTION
FloatVector(N)

Fugacities of all components [Pa].

Source code in src/polykin/thermo/eos/base.py
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
def f(self, T: float, P: float, y: FloatVector) -> FloatVector:
    r"""Calculate the fugacity of all components.

    $$ \hat{f}_i = \hat{\phi}_i y_i P $$

    Parameters
    ----------
    T : float
        Temperature [K].
    P : float
        Pressure [Pa].
    y : FloatVector (N)
        Mole fractions of all components [mol/mol].

    Returns
    -------
    FloatVector (N)
        Fugacities of all components [Pa].
    """
    return self.phi(T, P, y) * y * P

gR ¤

gR(T: float, P: float, y: FloatVector) -> float

Calculate the molar residual Gibbs energy of the fluid.

\[ g^R = B_m P \]
PARAMETER DESCRIPTION
T

Temperature [K].

TYPE: float

P

Pressure [Pa].

TYPE: float

y

Mole fractions of all components [mol/mol].

TYPE: FloatVector

RETURNS DESCRIPTION
float

Molar residual Gibbs energy [J/mol].

Source code in src/polykin/thermo/eos/virial.py
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
@override
def gR(self, T: float, P: float, y: FloatVector) -> float:
    """Calculate the molar residual Gibbs energy of the fluid.

    $$ g^R = B_m P $$

    Parameters
    ----------
    T : float
        Temperature [K].
    P : float
        Pressure [Pa].
    y : FloatVector
        Mole fractions of all components [mol/mol].

    Returns
    -------
    float
        Molar residual Gibbs energy [J/mol].
    """
    return self.Bm(T, y) * P

hR ¤

hR(T: float, P: float, y: FloatVector) -> float

Calculate the molar residual enthalpy of the fluid.

\[ h^R = g^R + T s^R \]
PARAMETER DESCRIPTION
T

Temperature [K].

TYPE: float

P

Pressure [Pa].

TYPE: float

y

Mole fractions of all components [mol/mol].

TYPE: FloatVector

RETURNS DESCRIPTION
float

Molar residual enthalpy [J/mol].

Source code in src/polykin/thermo/eos/base.py
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
def hR(self, T: float, P: float, y: FloatVector) -> float:
    r"""Calculate the molar residual enthalpy of the fluid.

    $$ h^R = g^R + T s^R $$

    Parameters
    ----------
    T : float
        Temperature [K].
    P : float
        Pressure [Pa].
    y : FloatVector
        Mole fractions of all components [mol/mol].

    Returns
    -------
    float
        Molar residual enthalpy [J/mol].
    """
    return self.gR(T, P, y) + T * self.sR(T, P, y)

kappa ¤

kappa(T: float, P: float, y: FloatVector) -> float

Calculate the isothermal compressibility coefficient.

\[ \kappa \equiv - \frac{1}{v} \left( \frac{\partial v}{\partial P} \right)_T = \frac{1}{P Z} \]
PARAMETER DESCRIPTION
T

Temperature [K].

TYPE: float

P

Pressure [Pa].

TYPE: float

y

Mole fractions of all components [mol/mol].

TYPE: FloatVector(N)

RETURNS DESCRIPTION
float

Isothermal compressibility coefficient, \(\kappa\) [Pa⁻¹].

Source code in src/polykin/thermo/eos/virial.py
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
@override
def kappa(self, T: float, P: float, y: FloatVector) -> float:
    r"""Calculate the isothermal compressibility coefficient.

    $$ \kappa \equiv
       - \frac{1}{v} \left( \frac{\partial v}{\partial P} \right)_T
       = \frac{1}{P Z} $$

    Parameters
    ----------
    T : float
        Temperature [K].
    P : float
        Pressure [Pa].
    y : FloatVector (N)
        Mole fractions of all components [mol/mol].

    Returns
    -------
    float
        Isothermal compressibility coefficient, $\kappa$ [Pa⁻¹].
    """
    Z = self.Z(T, P, y)
    return 1 / (P * Z)

phi ¤

phi(T: float, P: float, y: FloatVector) -> FloatVector

Calculate the fugacity coefficients of all components.

\[ \ln \hat{\phi}_i = \left(2\sum_j {y_jB_{ij}} -B_m \right)\frac{P}{RT} \]

References

  • RC Reid, JM Prausniz, and BE Poling. The properties of gases & liquids 4th edition, 1986, p. 145.
PARAMETER DESCRIPTION
T

Temperature [K].

TYPE: float

P

Pressure [Pa].

TYPE: float

y

Mole fractions of all components [mol/mol].

TYPE: FloatVector(N)

RETURNS DESCRIPTION
FloatVector(N)

Fugacity coefficients of all components.

Source code in src/polykin/thermo/eos/virial.py
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
@override
def phi(self, T: float, P: float, y: FloatVector) -> FloatVector:
    r"""Calculate the fugacity coefficients of all components.

    $$ \ln \hat{\phi}_i
       = \left(2\sum_j {y_jB_{ij}} -B_m \right)\frac{P}{RT} $$

    **References**

    *   RC Reid, JM Prausniz, and BE Poling. The properties of gases &
        liquids 4th edition, 1986, p. 145.

    Parameters
    ----------
    T : float
        Temperature [K].
    P : float
        Pressure [Pa].
    y : FloatVector (N)
        Mole fractions of all components [mol/mol].

    Returns
    -------
    FloatVector (N)
        Fugacity coefficients of all components.
    """
    B = self.Bij(T)
    Bm = self.Bm(T, y)
    return exp((2 * dot(B, y) - Bm) * P / (R * T))

sR ¤

sR(T: float, P: float, y: FloatVector) -> float

Calculate the molar residual entropy of the fluid.

\[ s^R = -\left( \frac{\partial g^R}{\partial T} \right)_{P,y_i} \]
PARAMETER DESCRIPTION
T

Temperature [K].

TYPE: float

P

Pressure [Pa].

TYPE: float

y

Mole fractions of all components [mol/mol].

TYPE: FloatVector

RETURNS DESCRIPTION
float

Molar residual entropy [J/(mol·K)].

Source code in src/polykin/thermo/eos/base.py
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
def sR(self, T: float, P: float, y: FloatVector) -> float:
    r"""Calculate the molar residual entropy of the fluid.

    $$ s^R = -\left( \frac{\partial g^R}{\partial T} \right)_{P,y_i} $$

    Parameters
    ----------
    T : float
        Temperature [K].
    P : float
        Pressure [Pa].
    y : FloatVector
        Mole fractions of all components [mol/mol].

    Returns
    -------
    float
        Molar residual entropy [J/(mol·K)].
    """
    return -1 * derivative_centered(lambda T_: self.gR(T_, P, y), T)[0]

v ¤

v(T: float, P: float, y: FloatVector) -> float

Calculate the molar volume the fluid.

\[ v = Z \frac{R T}{P} \]
PARAMETER DESCRIPTION
T

Temperature [K].

TYPE: float

P

Pressure [Pa].

TYPE: float

y

Mole fractions of all components [mol/mol].

TYPE: FloatVector(N)

RETURNS DESCRIPTION
float

Molar volume of the fluid [m³/mol].

Source code in src/polykin/thermo/eos/base.py
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
def v(self, T: float, P: float, y: FloatVector) -> float:
    r"""Calculate the molar volume the fluid.

    $$ v = Z \frac{R T}{P} $$

    Parameters
    ----------
    T : float
        Temperature [K].
    P : float
        Pressure [Pa].
    y : FloatVector (N)
        Mole fractions of all components [mol/mol].

    Returns
    -------
    float
        Molar volume of the fluid [m³/mol].
    """
    return self.Z(T, P, y) * R * T / P

vR ¤

vR(T: float, P: float, y: FloatVector) -> float

Calculate the molar residual volume of the fluid.

\[ v^R = \left(Z - 1 \right) \frac{R T}{P} \]
PARAMETER DESCRIPTION
T

Temperature [K].

TYPE: float

P

Pressure [Pa].

TYPE: float

y

Mole fractions of all components [mol/mol].

TYPE: FloatVector

RETURNS DESCRIPTION
float

Molar residual volume [m³/mol].

Source code in src/polykin/thermo/eos/base.py
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
def vR(self, T: float, P: float, y: FloatVector) -> float:
    r"""Calculate the molar residual volume of the fluid.

    $$ v^R = \left(Z - 1 \right) \frac{R T}{P} $$

    Parameters
    ----------
    T : float
        Temperature [K].
    P : float
        Pressure [Pa].
    y : FloatVector
        Mole fractions of all components [mol/mol].

    Returns
    -------
    float
        Molar residual volume [m³/mol].
    """
    return R * T / P * (self.Z(T, P, y) - 1.0)

Examples¤

Estimate the molar volume of a 50 mol% ethylene/nitrogen at 350 K and 10 bar.

from polykin.thermo.eos import Virial
import numpy as np

Tc = [282.4, 126.2]    # K
Pc = [50.4e5, 33.9e5]  # Pa
Zc = [0.280, 0.290]
w = [0.089, 0.039]

eos = Virial(Tc, Pc, Zc, w)
v = eos.v(T=350., P=10e5, y=np.array([0.5, 0.5]))

print(f"{v:.2e} m³/mol")
2.87e-03 m³/mol