Mixing Rules

Cubic EoS mixing rules implementations module.

class CubicMixRule[source]

Bases: ABC

Cubic mix rule abstract class.

abstractmethod set_mixrule(ar_model_id: int) None[source]

Set mix rule abstract method.

Changes the default mixing rule of the given cubic EoS model.

Parameters:

ar_model_id (int) – ID of the cubic EoS model

Raises:

NotImplementedError – Abstract error, this method must be implemented in the subclass

class QMR(kij, lij)[source]

Bases: CubicMixRule

Quadratic mixing rule.

Parameters:
  • kij (array_like) – kij binary interaction parameters matrix

  • lij (array_like) – lij binary interaction parameters matrix

kij

kij binary interaction parameters matrix

Type:

array_like

lij

lij binary interaction parameters matrix

Type:

array_like

Example

from yaeos import QMR, SoaveRedlichKwong

kij = [[0.0, 0.1], [0.1, 0.0]]
lij = [[0.0, 0.02], [0.02, 0.0]]

mixrule = QMR(kij, lij)     # Quadratic mixing rule instance

tc = [305.32, 469.7]        # critical temperature [K]
pc = [48.72, 33.7]          # critical pressure [bar]
w = [0.0995, 0.152]         # acentric factor

model = SoaveRedlichKwong(tc, pc, w, mixrule)
name = 'QMR'
set_mixrule(ar_model_id: int) None[source]

Set quadratic mix rule method.

Parameters:

ar_model_id (int) – ID of the cubic EoS model

class QMRTD(kij_0, kij_inf, t_ref, lij)[source]

Bases: CubicMixRule

Quadratic mixing rule, with temperature dependence.

..math::

k_{ij}(T) = k_{ij}^{infty} + k_{ij}^0 exp{left(-T/T^{ref}right)}

Parameters:
  • kij_0 (matrix_like) – kij_0 binary interaction parameters matrix

  • kij_inf (matrix_like) – kij_inf binary interaction parameters matrix

  • t_ref (matrix_like) – Reference temperature

  • lij (matrix_like) – lij binary interaction parameters matrix

kij_0

kij_0 binary interaction parameters matrix

Type:

matrix_like

kij_inf

kij_inf binary interaction parameters matrix

Type:

matrix_like

t_ref

Reference temperature

Type:

matrix_like

lij

lij binary interaction parameters matrix

Type:

matrix_like

Example

from yaeos import QMRTD, SoaveRedlichKwong

kij_0 = [[0.0, 0.1], [0.1, 0.0]]
kij_inf = [[0.0, 0.1], [0.1, 0.0]]
Tref = [[0.0, 390], [390, 0.0]]
lij = [[0.0, 0.02], [0.02, 0.0]]

# Quadratic mixing rule instance
mixrule = QMRTD(kij_0, kij_inf, Tref, lij)

tc = [305.32, 469.7]        # critical temperature [K]
pc = [48.72, 33.7]          # critical pressure [bar]
w = [0.0995, 0.152]         # acentric factor

model = SoaveRedlichKwong(tc, pc, w, mixrule)
name = 'QMRTD'
set_mixrule(ar_model_id: int) None[source]

Set mix rule method.

class MHV(ge: GeModel, q: float, lij=None)[source]

Bases: CubicMixRule

Modified Huron-Vidal mixing rule.

Parameters:
  • ge (GeModel) – Excess Gibbs energy model

  • q (float) –

    q parameter. Use:

    q = -0.594 for Soave-Redlich-Kwong q = -0.53 for Peng-Robinson q = -0.85 for Van der Waals

  • lij (array_like, optional) – lij binary interaction parameters matrix, by default None

ge

Excess Gibbs energy model

Type:

GeModel

q

q parameter

Type:

float

lij

lij binary interaction parameters matrix

Type:

array_like

Example

from yaeos import MHV, SoaveRedlichKwong, NRTL

tc = [647.14, 513.92]               # critical temperature [K]
pc = [220.64, 61.48]                # critical pressure [bar]
w =  [0.344, 0.649]                 # acentric factor

a = [[0, 3.458], [-0.801, 0]]       # NRTL aij parameters
b = [[0, -586.1], [246.2, 0]]       # NRTL bij parameters
c = [[0, 0.3], [0.3, 0]]            # NRTL cij parameters

ge_model = NRTL(a, b, c)
mixrule = MHV(ge_model, q=-0.53)

model_mhv = PengRobinson76(tc, pc, w, mixrule)
name = 'MHV'
set_mixrule(ar_model_id: int) None[source]

Set modified Huron-Vidal mix rule method.

Parameters:

ar_model_id (int) – ID of the cubic EoS model

class HV(ge: GeModel)[source]

Bases: CubicMixRule

Huron-Vidal mixing rule.

Parameters:

ge (GeModel) – Excess Gibbs energy model

ge

Excess Gibbs energy model

Type:

GeModel

Example

from yaeos import HV, SoaveRedlichKwong, NRTL

tc = [647.14, 513.92]               # critical temperature [K]
pc = [220.64, 61.48]                # critical pressure [bar]
w =  [0.344, 0.649]                 # acentric factor

a = [[0, 3.458], [-0.801, 0]]       # NRTL aij parameters
b = [[0, -586.1], [246.2, 0]]       # NRTL bij parameters
c = [[0, 0.3], [0.3, 0]]            # NRTL cij parameters

ge_model = NRTL(a, b, c)
mixrule = HV(ge_model)

model_hv = PengRobinson76(tc, pc, w, mixrule)
name = 'HV'
set_mixrule(ar_model_id: int) None[source]

Set modified Huron-Vidal mix rule method.

Parameters:

ar_model_id (int) – ID of the cubic EoS model

class HVNRTL(alpha, gji, gjiT=None, use_kij=None, kij=None)[source]

Bases: CubicMixRule

Huron-Vidal mixing rule using an NRTL model.

Huron Vidal mixing rule coupled with the excess Gibbs energy model defined by Huron-Vidal. This model can use the classic VdW1f mixing rules kij parameters when desired.

Parameters:
  • alpha (matrix_like) – NRTL alpha parameters matrix

  • gji (matrix_like) – NRTL gij parameters matrix

  • gjiT (matrix_like) – NRTL gij temperature coefficient parameters matrix

  • use_kij (matrix_like) – Boolean matrix indicating whether to use kij parameters

  • kij (matrix_like) – kij binary interaction parameters matrix

alpha

NRTL alpha parameters matrix

Type:

matrix_like

gji

NRTL gij parameters matrix

Type:

matrix_like

gjiT

NRTL gij temperature coefficient parameters matrix

Type:

matrix_like

use_kij

Boolean matrix indicating whether to use kij parameters

Type:

matrix_like

kij

kij binary interaction parameters matrix

Type:

matrix_like

name = 'HVNRTL'
set_mixrule(ar_model_id: int) None[source]

Set modified Huron-Vidal mix rule method.

Parameters:

ar_model_id (int) – ID of the cubic EoS model

classmethod setup_from_kij(model, temperatures, kij)[source]

Generalizes the calculation of gji parameters and returns linear correlation matrices gji0 and gjiT for any number of components.

Warning

Internally this function depends on the \(\delta_1\) of the model. For the case of the RKPR EoS where there are different \(\delta_1\) parameters for each component the method will only use a mean value. Which will make the predictions further away of the predictions being made with \(k_{ij}\) compared to using classic Cubic EoS. This method is intended only as a way to ease the initialization of parameters on optimization methods.

Parameters:
  • model (yaeos.CubicEoS) – Thermodynamic model intended to be used as reference model.

  • temperatures (array_like) – Array of temperature values on which make the correlation

  • kij (float) – kij value to replicate.

Returns:

  • yaeos.HVNRTL (Instance of a HVNRTL mixing rule with parameters)

  • set to be close the parameters that would replicate the classic

  • QMR mixing rule.

class sDDLC(qs, kij_0=None, kij_inf=None, t_ref=None, lij=None)[source]

Bases: CubicMixRule

segmented Density Dependent Local Composition mixing rule.

..math::

k_{ij}(T) = k_{ij}^{infty} + k_{ij}^0 exp{left(-T/T^{ref}right)}

Parameters:
  • q (array_like) – Segment size for each component

  • kij_0 (matrix_like) – kij_0 binary interaction parameters matrix

  • kij_inf (matrix_like) – kij_inf binary interaction parameters matrix

  • t_ref (matrix_like) – Reference temperature

  • lij (matrix_like) – lij binary interaction parameters matrix

q

Segment size for each component

Type:

array_like

kij_0, optional

kij_0 binary interaction parameters matrix. Defaults to zero.

Type:

matrix_like

kij_inf, optional

kij_inf binary interaction parameters matrix. Defaults to zero.

Type:

matrix_like

t_ref, optional

Reference temperature. Defaults to zero.

Type:

matrix_like

lij, optional

lij binary interaction parameters matrix. Defaults to zero.

Type:

matrix_like

Example

from yaeos import sDDLC, SoaveRedlichKwong

qs = [1.16, 6.0]
kij_0 = [[0.0, 0.1], [0.1, 0.0]]
kij_inf = [[0.0, 0.1], [0.1, 0.0]]
Tref = [[0.0, 390], [390, 0.0]]
lij = [[0.0, 0.02], [0.02, 0.0]]

# sDDLC
mixrule = sDDLC(qs, kij_0, kij_inf, Tref, lij)

tc = [305.32, 469.7]        # critical temperature [K]
pc = [48.72, 33.7]          # critical pressure [bar]
w = [0.0995, 0.152]         # acentric factor

model = SoaveRedlichKwong(tc, pc, w, mixrule)
name = 'sDDLC'
set_mixrule(ar_model_id: int) None[source]

Set mix rule method.