Skip to content

polykin.thermo.eos¤

IdealGas ¤

Ideal gas equation of state.

This EOS is based on the following \(P(v,T)\) relationship:

\[ P = \frac{R T}{v} \]

where \(P\) is the pressure, \(T\) is the temperature, and \(v\) is the molar volume.

PARAMETER DESCRIPTION
N

Number of components.

TYPE: int

name

Name.

TYPE: str DEFAULT: ''

Source code in src/polykin/thermo/eos/idealgas.py
 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
class IdealGas(GasEoS):
    r"""Ideal gas equation of state.

    This EOS is based on the following $P(v,T)$ relationship:

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

    where $P$ is the pressure, $T$ is the temperature, and $v$ is the molar
    volume.

    Parameters
    ----------
    N : int
        Number of components.
    name : str
        Name.
    """

    def __init__(self,
                 N: int,
                 name: str = ''
                 ) -> None:

        super().__init__(N, name)

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

        $$ Z = 1 $$

        Returns
        -------
        float
            Compressibility factor.
        """
        return 1.0

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

        Parameters
        ----------
        T : float
            Temperature. Unit = K.
        v : float
            Molar volume. Unit = m³/mol.

        Returns
        -------
        float
             Pressure. Unit = Pa.
        """
        return R*T/v

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

        $$ \hat{\phi}_i = 1 $$

        Returns
        -------
        FloatVector
            Fugacity coefficients of all components.
        """
        return np.ones(self.N)

    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. Unit = K.
        V : float
            Volume. Unit = m³.
        n : FloatVector
            Mole amounts of all components. Unit = mol.
        v0 : float
            Molar volume in reference state. Unit = m³/mol.

        Returns
        -------
        float
            Helmholtz energy departure, $A - A^{\circ}$. Unit = J.
        """
        nT = n.sum()
        return -nT*R*T*log(V/(nT*v0))

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. Unit = K.

TYPE: float

V

Volume. Unit = m³.

TYPE: float

n

Mole amounts of all components. Unit = mol.

TYPE: FloatVector

v0

Molar volume in reference state. Unit = m³/mol.

TYPE: float

RETURNS DESCRIPTION
float

Helmholtz energy departure, \(A - A^{\circ}\). Unit = J.

Source code in src/polykin/thermo/eos/idealgas.py
 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
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. Unit = K.
    V : float
        Volume. Unit = m³.
    n : FloatVector
        Mole amounts of all components. Unit = mol.
    v0 : float
        Molar volume in reference state. Unit = m³/mol.

    Returns
    -------
    float
        Helmholtz energy departure, $A - A^{\circ}$. Unit = J.
    """
    nT = n.sum()
    return -nT*R*T*log(V/(nT*v0))

N property ¤

N: int

Number of components.

P ¤

P(T: float, v: float, y=None) -> float

Calculate the pressure of the fluid.

PARAMETER DESCRIPTION
T

Temperature. Unit = K.

TYPE: float

v

Molar volume. Unit = m³/mol.

TYPE: float

RETURNS DESCRIPTION
float

Pressure. Unit = Pa.

Source code in src/polykin/thermo/eos/idealgas.py
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
def P(self,
      T: float,
      v: float,
      y=None
      ) -> float:
    r"""Calculate the pressure of the fluid.

    Parameters
    ----------
    T : float
        Temperature. Unit = K.
    v : float
        Molar volume. Unit = m³/mol.

    Returns
    -------
    float
         Pressure. Unit = Pa.
    """
    return R*T/v

Z ¤

Z(T=None, P=None, y=None) -> float

Calculate the compressibility factor of the fluid.

\[ Z = 1 \]
RETURNS DESCRIPTION
float

Compressibility factor.

Source code in src/polykin/thermo/eos/idealgas.py
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
def Z(self,
      T=None,
      P=None,
      y=None
      ) -> float:
    r"""Calculate the compressibility factor of the fluid.

    $$ Z = 1 $$

    Returns
    -------
    float
        Compressibility factor.
    """
    return 1.0

f ¤

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

Calculate the fugacity of all components.

For each component, the fugacity is given by:

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

where \(\hat{\phi}_i(T,P,y)\) is the fugacity coefficient, \(P\) is the pressure, and \(y_i\) is the mole fraction.

PARAMETER DESCRIPTION
T

Temperature. Unit = K.

TYPE: float

P

Pressure. Unit = Pa.

TYPE: float

y

Mole fractions of all components. Unit = mol/mol.

TYPE: FloatVector

RETURNS DESCRIPTION
FloatVector

Fugacities of all components. Unit = Pa.

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

    For each component, the fugacity is given by:

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

    where $\hat{\phi}_i(T,P,y)$ is the fugacity coefficient, $P$ is the 
    pressure, and $y_i$ is the mole fraction.

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

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

phi ¤

phi(T=None, P=None, y=None) -> FloatVector

Calculate the fugacity coefficients of all components.

\[ \hat{\phi}_i = 1 \]
RETURNS DESCRIPTION
FloatVector

Fugacity coefficients of all components.

Source code in src/polykin/thermo/eos/idealgas.py
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
def phi(self,
        T=None,
        P=None,
        y=None
        ) -> FloatVector:
    r"""Calculate the fugacity coefficients of all components.

    $$ \hat{\phi}_i = 1 $$

    Returns
    -------
    FloatVector
        Fugacity coefficients of all components.
    """
    return np.ones(self.N)

v ¤

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

Calculate the molar volume the fluid.

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

where \(v\) is the molar volume, \(Z(T, P, y)\) is the compressibility factor, \(T\) is the temperature, \(P\) is the pressure, and \(y\) is the mole fraction vector.

PARAMETER DESCRIPTION
T

Temperature. Unit = K.

TYPE: float

P

Pressure. Unit = Pa.

TYPE: float

y

Mole fractions of all components. Unit = mol/mol.

TYPE: FloatVector

RETURNS DESCRIPTION
float

Molar volume of the fluid. Unit = m³/mol.

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

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

    where $v$ is the molar volume, $Z(T, P, y)$ is the compressibility
    factor, $T$ is the temperature, $P$ is the pressure, and $y$ is the
    mole fraction vector.

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

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

Examples¤

Estimate the molar volume of a gas at 0°C and 1 atm.

from polykin.thermo.eos import IdealGas

eos = IdealGas(1)
v = eos.v(T=273.15, P=1.01325e5, y=None)

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