Core Module

yaeos Python API core module.

ArModel and GeModel abstract classes definition. Also, the implementation of the models’ thermoprops methods.

class GeModel[source]

Bases: ABC

Excess Gibbs (Ge) model abstract class.

ln_gamma(moles, temperature)[source]

Calculate activity coefficients.

Calculate \(\ln \gamma_i(n,T)\) vector.

Parameters:
  • moles (array_like) – Moles number vector [mol]

  • temperature (float) – Temperature [K]

Returns:

\(ln \gamma_i(n,T)\) vector

Return type:

np.ndarray

Example

import numpy as np

from yaeos import NRTL

a = np.array([[0, 0.3], [0.3, 0]])
b = np.array([[0, 0.4], [0.4, 0]])
c = np.array([[0, 0.5], [0.5, 0]])

nrtl = NRTL(a, b, c)

print(nrtl.ln_gamma([5.0, 5.6], 300.0))
class ArModel[source]

Bases: ABC

Residual Helmholtz (Ar) model abstract class.

lnphi_vt(moles, volume: float, temperature: float, dt: bool = False, dp: bool = False, dn: bool = False) ndarray | tuple[ndarray, dict][source]

Calculate fugacity coefficent given volume and temperature.

Calculate \(ln \phi_i(n,V,T)\) and its derivatives with respect to temperature, pressure and moles number.

Parameters:
  • moles (array_like) – Moles number vector [mol]

  • volume (float) – Volume [L]

  • temperature (float) – Temperature [K]

  • dt (bool, optional) – Calculate temperature derivative, by default False

  • dp (bool, optional) – Calculate pressure derivative, by default False

  • dn (bool, optional) – Calculate moles derivative, by default False

Returns:

\(ln \phi_i(n,V,T)\) vector or tuple with \(ln \phi_i(n,V,T)\) vector and derivatives dictionary if any derivative is asked

Return type:

Union[np.ndarray, tuple[np.ndarray, dict]]

Example

import numpy as np

from yaeos import PengRobinson76

tc = np.array([320.0, 375.0])   # critical temperatures [K]
pc = np.array([45.0, 60.0])     # critical pressures [bar]
w = np.array([0.0123, 0.045])   # acentric factors

model = PengRobinson76(tc, pc, w)

# Evaluating ln_phi only
# will print: [-1.45216274 -2.01044828]

print(model.lnphi_vt([5.0, 5.6], 1.0, 300.0))

# Asking for derivatives
# will print:
# (
# array([-1.45216274, -2.01044828]),
# {'dt': array([0.01400063, 0.01923493]), 'dp': None, 'dn': None}
# )

print(model.lnphi_vt([5.0, 5.6], 1.0, 300.0, dt=True)
lnphi_pt(moles, pressure: float, temperature: float, root: str = 'stable', dt: bool = False, dp: bool = False, dn: bool = False) ndarray | tuple[ndarray, dict][source]

Calculate fugacity coefficent given pressure and temperature.

Calculate \(ln \phi_i(n,P,T)\) and its derivatives with respect to temperature, pressure and moles number.

Parameters:
  • moles (array_like) – Moles number vector [mol]

  • pressure (float) – Pressure [bar]

  • temperature (float) – Temperature [K]

  • root (str, optional) – Volume root, use: “liquid”, “vapor” or “stable”, by default “stable”

  • dt (bool, optional) – Calculate temperature derivative, by default False

  • dp (bool, optional) – Calculate pressure derivative, by default False

  • dn (bool, optional) – Calculate moles derivative, by default False

Returns:

\(ln \phi_i(n,P,T)\) vector or tuple with \(ln \phi_i(n,P,T)\) vector and derivatives dictionary if any derivative is asked

Return type:

Union[np.ndarray, tuple[np.ndarray, dict]]

Example

import numpy as np

from yaeos import PengRobinson76

tc = np.array([320.0, 375.0])   # critical temperatures [K]
pc = np.array([45.0, 60.0])     # critical pressures [bar]
w = np.array([0.0123, 0.045])   # acentric factors

model = PengRobinson76(tc, pc, w)

# Evaluating ln_phi only
# will print: [-0.10288733 -0.11909807]

print(model.lnphi_pt([5.0, 5.6], 10.0, 300.0))

# Asking for derivatives
# will print:
# (
# array([-0.10288733, -0.11909807]),
# {'dt': array([0.00094892, 0.00108809]), 'dp': None, 'dn': None}
# )

print(model.lnphi_pt([5.0, 5.6], 10.0, 300.0, dt=True)
pressure(moles, volume: float, temperature: float, dv: bool = False, dt: bool = False, dn: bool = False) float | tuple[float, dict][source]

Calculate pressure given volume and temperature [bar].

Calculate \(P(n,V,T)\) and its derivatives with respect to volume, temperature and moles number.

Parameters:
  • moles (array_like) – Moles number vector [mol]

  • volume (float) – Volume [L]

  • temperature (float) – Temperature [K]

  • dv (bool, optional) – Calculate volume derivative, by default False

  • dt (bool, optional) – Calculate temperature derivative, by default False

  • dn (bool, optional) – Calculate moles derivative, by default False

Returns:

Pressure or tuple with Presure and derivatives dictionary if any derivative is asked [bar]

Return type:

Union[float, tuple[float, dict]]

Example

import numpy as np

from yaeos import PengRobinson76

tc = np.array([320.0, 375.0])   # critical temperatures [K]
pc = np.array([45.0, 60.0])     # critical pressures [bar]
w = np.array([0.0123, 0.045])   # acentric factors

model = PengRobinson76(tc, pc, w)

# Evaluating pressure only
# will print: 16.011985733846956

print(model.pressure(np.array([5.0, 5.6]), 2.0, 300.0))

# Asking for derivatives
# will print:
# (
# 16.011985733846956,
# {'dv': None, 'dt': np.float64(0.7664672352866752), 'dn': None}
# )

print(model.pressure(np.array([5.0, 5.6]), 2.0, 300.0, dt=True))
volume(moles, pressure: float, temperature: float, root: str = 'stable') float[source]

Calculate volume given pressure and temperature [L].

Parameters:
  • moles (array_like) – Moles number vector [mol]

  • pressure (float) – Pressure [bar]

  • temperature (float) – Temperature [K]

  • root (str, optional) – Volume root, use: “liquid”, “vapor” or “stable”, by default “stable”

Returns:

Volume [L]

Return type:

float

Example

import numpy as np

from yaeos import PengRobinson76

tc = np.array([320.0, 375.0])   # critical temperatures [K]
pc = np.array([45.0, 60.0])     # critical pressures [bar]
w = np.array([0.0123, 0.045])   # acentric factors

model = PengRobinson76(tc, pc, w)

# Evaluating stable root volume
# will print: 23.373902973572587

print(model.volume(np.array([5.0, 5.6]), 10.0, 300.0))

# Liquid root volume (not stable)
# will print: 0.8156388756398074

print(model.volume(np.array([5.0, 5.6]), 10.0, 300.0, "liquid"))
enthalpy_residual_vt(moles, volume: float, temperature: float, dt: bool = False, dv: bool = False, dn: bool = False) float | tuple[float, dict][source]

Calculate residual enthalpy given volume and temperature [bar L].

Parameters:
  • moles (array_like) – Moles number vector [mol]

  • volume (float) – Volume [L]

  • temperature (float) – Temperature [K]

Returns:

Residual enthalpy or tuple with Residual enthalpy and derivatives dictionary if any derivative is asked [bar L]

Return type:

Union[float, tuple[float, dict]]

Example

import numpy as np

from yaeos import PengRobinson76

tc = np.array([320.0, 375.0])   # critical temperatures [K]
pc = np.array([45.0, 60.0])     # critical pressures [bar]
w = np.array([0.0123, 0.045])   # acentric factors

model = PengRobinson76(tc, pc, w)

# Evaluating residual enthalpy only
# will print: -182.50424367123696

print(
    model.enthalpy_residual_vt(np.array([5.0, 5.6]), 10.0, 300.0)
)

# Asking for derivatives
# will print:
# (
# -182.50424367123696,
# {'dt': 0.21542452742588686, 'dv': None, 'dn': None}
# )

print(
    model.enthalpy_residual_vt(
        np.array([5.0, 5.6]),
        10.0,
        300.0,
        dt=True)
    )
)
gibbs_residual_vt(moles, volume: float, temperature: float, dt: bool = False, dv: bool = False, dn: bool = False) float | tuple[float, dict][source]

Calculate residual Gibbs energy at volume and temperature [bar L].

Parameters:
  • moles (array_like) – Moles number vector [mol]

  • volume (float) – Volume [L]

  • temperature (float) – Temperature [K]

Returns:

Residual Gibbs energy or tuple with Residual Gibbs energy and derivatives dictionary if any derivative is asked [bar L]

Return type:

Union[float, tuple[float, dict]]

Example

import numpy as np

from yaeos import PengRobinson76

tc = np.array([320.0, 375.0])   # critical temperatures [K]
pc = np.array([45.0, 60.0])     # critical pressures [bar]
w = np.array([0.0123, 0.045])   # acentric factors

model = PengRobinson76(tc, pc, w)

# Evaluating residual gibbs energy only
# will print: -138.60374582274

print(model.gibbs_residual_vt(np.array([5.0, 5.6]), 10.0, 300.0))

# Asking for derivatives
# will print:
# (
# -138.60374582274,
# {'dt': 0.289312908265414, 'dv': None, 'dn': None}
# )

print(
    model.gibbs_residual_vt(
        np.array([5.0, 5.6]),
        10.0,
        300.0,
        dt=True
    )
)
entropy_residual_vt(moles, volume: float, temperature: float, dt: bool = False, dv: bool = False, dn: bool = False) float | tuple[float, dict][source]

Calculate residual entropy given volume and temperature [bar L / K].

Parameters:
  • moles (array_like) – Moles number vector [mol]

  • volume (float) – Volume [L]

  • temperature (float) – Temperature [K]

Returns:

Residual entropy or tuple with Residual entropy and derivatives dictionary if any derivative is asked [bar L]

Return type:

Union[float, tuple[float, dict]]

Example

import numpy as np

from yaeos import PengRobinson76

tc = np.array([320.0, 375.0])   # critical temperatures [K]
pc = np.array([45.0, 60.0])     # critical pressures [bar]
w = np.array([0.0123, 0.045])   # acentric factors

model = PengRobinson76(tc, pc, w)

# Evaluating residual entropy only
# will print: -0.1463349928283233

print(model.entropy_residual_vt(np.array([5.0, 5.6]), 10.0, 300.0))

# Asking for derivatives
# will print:
# (
# (-0.1463349928283233,
# {'dt': 0.00024148870662932045, 'dv': None, 'dn': None})
# )

print(
    model.entropy_residual_vt(
        np.array([5.0, 5.6]),
        10.0,
        300.0,
        dt=True
    )
)
cv_residual_vt(moles, volume: float, temperature: float) float[source]

Residual isochoric heat capacity given V and T [bar L / K].

Parameters:
  • moles (array_like) – Moles number vector [mol]

  • volume (float) – Volume [L]

  • temperature (float) – Temperature [K]

Returns:

Residual isochoric heat capacity [bar L / K]

Return type:

float

Example


import numpy as np

from yaeos import PengRobinson76

tc = np.array([320.0, 375.0]) # critical temperatures [K] pc = np.array([45.0, 60.0]) # critical pressures [bar] w = np.array([0.0123, 0.045]) # acentric factors

model = PengRobinson76(tc, pc, w)

# Evaluating residual isochoric heat capacity only # will print: 0.07244661198879614

print(model.cv_residual_vt(np.array([5.0, 5.6]), 10.0, 300.0))

cp_residual_vt(moles, volume: float, temperature: float) float[source]

Calculate residual isobaric heat capacity given V and T [bar L / K].

Parameters:
  • moles (array_like) – Moles number vector [mol]

  • volume (float) – Volume [L]

  • temperature (float) – Temperature [K]

Returns:

Residual isochoric heat capacity [bar L / K]

Return type:

float

Example

import numpy as np

from yaeos import PengRobinson76

tc = np.array([320.0, 375.0])   # critical temperatures [K]
pc = np.array([45.0, 60.0])     # critical pressures [bar]
w = np.array([0.0123, 0.045])   # acentric factors

model = PengRobinson76(tc, pc, w)

# Evaluating residual isobaric heat capacity only
# will print: 1.4964025088916886

print(model.cp_residual_vt(np.array([5.0, 5.6]), 10.0, 300.0))
flash_pt(z, pressure: float, temperature: float, k0=None) dict[source]

Two-phase split with specification of temperature and pressure.

Parameters:
  • z (array_like) – Global mole fractions

  • pressure (float) – Pressure [bar]

  • temperature (float) – Temperature [K]

  • k0 (array_like, optional) – Initial guess for the split, by default None (will use k_wilson)

Returns:

Flash result dictionary with keys:
  • x: heavy phase mole fractions

  • y: light phase mole fractions

  • Vx: heavy phase volume [L]

  • Vy: light phase volume [L]

  • P: pressure [bar]

  • T: temperature [K]

  • beta: light phase fraction

Return type:

dict

Example

import numpy as np

from yaeos import PengRobinson76

tc = np.array([369.83, 507.6])       # critical temperatures [K]
pc = np.array([42.48, 30.25])        # critical pressures [bar]
w = np.array([0.152291, 0.301261])   # acentric factors

model = PengRobinson76(tc, pc, w)

# Flash calculation
# will print:
# {
#   'x': array([0.3008742, 0.6991258]),
#   'y': array([0.85437317, 0.14562683]),
#   'Vx': 0.12742569165483714,
#   'Vy': 3.218831515959867,
#   'P': 8.0,
#   'T': 350.0,
#   'beta': 0.35975821044266726
# }

print(model.flash_pt([0.5, 0.5], 8.0, 350.0))
flash_pt_grid(z, pressures, temperatures, parallel=False) dict[source]

Two-phase split with specification of temperature and pressure grid.

Parameters:
  • z (array_like) – Global mole fractions

  • pressures (array_like) – Pressures grid [bar]

  • temperatures (array_like) – Temperatures grid [K]

  • parallel (bool, optional) – Use parallel processing, by default False

Returns:

Flash grid result dictionary with keys:
  • x: heavy phase mole fractions

  • y: light phase mole fractions

  • Vx: heavy phase volume [L]

  • Vy: light phase volume [L]

  • P: pressure [bar]

  • T: temperature [K]

  • beta: light phase fraction

Return type:

dict

Example

import numpy as np

from yaeos import PengRobinson76

tc = np.array([369.83, 507.6])       # critical temperatures [K]
pc = np.array([42.48, 30.25])        # critical pressures [bar]
w = np.array([0.152291, 0.301261])

temperatures = [350.0, 360.0, 400.0]
pressures = [10, 20, 30]
saturation_pressure(z, temperature: float, kind: str = 'bubble', p0: float = 0) dict[source]

Saturation pressure at specified temperature.

Parameters:
  • z (array_like) – Global molar fractions

  • temperature (float) – Temperature [K]

  • kind (str, optional) –

    Kind of saturation point, defaults to “bubble”. Options are
    • ”bubble”

    • ”dew”

    • ”liquid-liquid”

Returns:

Saturation pressure calculation result dictionary with keys:
  • x: heavy phase mole fractions

  • y: light phase mole fractions

  • Vx: heavy phase volume [L]

  • Vy: light phase volume [L]

  • P: pressure [bar]

  • T: temperature [K]

  • beta: light phase fraction

Return type:

dict

Example

import numpy as np

from yaeos import PengRobinson76

tc = np.array([369.83, 507.6])       # critical temperatures [K]
pc = np.array([42.48, 30.25])        # critical pressures [bar]
w = np.array([0.152291, 0.301261])   # acentric factors

model = PengRobinson76(tc, pc, w)

# Saturation pressure calculation
# will print:
# {
# 'x': array([0.5, 0.5]),
# 'y': array([0.9210035 , 0.07899651]),
# 'Vx': 0.11974125553488875,
# 'Vy': 1.849650524323853,
# 'T': 350.0,
# 'P': 12.990142036059941,
# 'beta': 0.0
# }

print(model.saturation_pressure(np.array([0.5, 0.5]), 350.0))
saturation_temperature(z, pressure: float, kind: str = 'bubble', t0: float = 0) dict[source]

Saturation pressure at specified temperature.

Parameters:
  • z (array_like) – Global molar fractions

  • pressure (float) – Pressure [bar]

  • kind (str, optional) –

    Kind of saturation point, defaults to “bubble”. Options are
    • ”bubble”

    • ”dew”

    • ”liquid-liquid”

Returns:

Saturation pressure calculation result dictionary with keys:
  • x: heavy phase mole fractions

  • y: light phase mole fractions

  • Vx: heavy phase volume [L]

  • Vy: light phase volume [L]

  • P: pressure [bar]

  • T: temperature [K]

  • beta: light phase fraction

Return type:

dict

Example

import numpy as np

from yaeos import PengRobinson76

tc = np.array([369.83, 507.6])       # critical temperatures [K]
pc = np.array([42.48, 30.25])        # critical pressures [bar]
w = np.array([0.152291, 0.301261])   # acentric factors

model = PengRobinson76(tc, pc, w)

# Saturation pressure calculation
# will print:
# {
# 'x': array([0.5, 0.5]),
# 'y': array([0.9210035 , 0.07899651]),
# 'Vx': 0.11974125553488875,
# 'Vy': 1.849650524323853,
# 'T': 350.0,
# 'P': 12.99,
# 'beta': 0.0
# }

print(model.saturation_temperature(np.array([0.5, 0.5]), 12.99))
phase_envelope_pt(z, kind: str = 'bubble', max_points: int = 300, t0: float = 150.0, p0: float = 1.0) dict[source]

Two phase envelope calculation (PT).

Parameters:
  • z (array_like) – Global mole fractions

  • kind (str, optional) – Kind of saturation point to start the envelope calculation, defaults to “bubble”. Options are - “bubble” - “dew”

  • max_points (int, optional) – Envelope’s maximum points to calculate (T, P), by default 300

  • t0 (float, optional) – Initial guess for temperature [K] for the saturation point of kind: kind, by default 150.0

  • p0 (float, optional) – Initial guess for pressure [bar] for the saturation point of kind: kind, by default 1.0

Returns:

Envelope calculation result dictionary with keys:
  • Ts: temperatures [K]

  • Ps: pressures [bar]

  • Tcs: critical temperatures [K]

  • Pcs: critical pressures [bar]

Return type:

dict

Example

import numpy as np

import matplotlib.pyplot as plt

from yaeos import PengRobinson76

tc = np.array([369.83, 507.6])       # critical temperatures [K]
pc = np.array([42.48, 30.25])        # critical pressures [bar]
w = np.array([0.152291, 0.301261])   # acentric factors

model = PengRobinson76(tc, pc, w)

# Two phase envelope calculation and plot
env = model.phase_envelope_pt(
    np.array([0.5, 0.5]),
    t0=150.0,
    p0=1.0
)

plt.plot(env["Ts"], env["Ps"])
plt.scatter(env["Tcs"], env["Pcs"])