NRTL Gibbs Excess Model (thermo.nrtl)

This module contains a class NRTL for performing activity coefficient calculations with the NRTL model. An older, functional calculation for activity coefficients only is also present, NRTL_gammas.

For reporting bugs, adding feature requests, or submitting pull requests, please use the GitHub issue tracker.

NRTL Class

class thermo.nrtl.NRTL(T, xs, tau_coeffs=None, alpha_coeffs=None, ABEFGHCD=None, tau_as=None, tau_bs=None, tau_es=None, tau_fs=None, tau_gs=None, tau_hs=None, alpha_cs=None, alpha_ds=None)[source]

Bases: thermo.activity.GibbsExcess

Class for representing an a liquid with excess gibbs energy represented by the NRTL equation. This model is capable of representing VL and LL behavior. [1] and [2] are good references on this model.

gE=RTixijτjiGjixjjGjixjg^E = RT\sum_i x_i \frac{\sum_j \tau_{ji} G_{ji} x_j} {\sum_j G_{ji}x_j}
Gij=exp(αijτij)G_{ij}=\exp(-\alpha_{ij}\tau_{ij})
αij=cij+dijT\alpha_{ij}=c_{ij} + d_{ij}T
τij=Aij+BijT+EijlnT+FijT+GijT2+HijT2\tau_{ij}=A_{ij}+\frac{B_{ij}}{T}+E_{ij}\ln T + F_{ij}T + \frac{G_{ij}}{T^2} + H_{ij}{T^2}
Parameters
Tfloat

Temperature, [K]

xslist[float]

Mole fractions, [-]

tau_coeffslist[list[list[float]]], optional

NRTL parameters, indexed by [i][j] and then each value is a 6 element list with parameters [a, b, e, f, g, h]; either (tau_coeffs and alpha_coeffs) or ABEFGHCD are required, [-]

alpha_coeffslist[list[float]], optional

NRTL alpha parameters, []

ABEFGHCDtuple[list[list[float]], 8], optional

Contains the following. One of (tau_coeffs and alpha_coeffs) or ABEFGHCD or some of the tau or alpha parameters are required, [-]

tau_aslist[list[float]], optional

a parameters used in calculating NRTL.taus, [-]

tau_bslist[list[float]], optional

b parameters used in calculating NRTL.taus, [K]

tau_eslist[list[float]], optional

e parameters used in calculating NRTL.taus, [-]

tau_fslist[list[float]], optional

f paraemeters used in calculating NRTL.taus, [1/K]

tau_gslist[list[float]], optional

e parameters used in calculating NRTL.taus, [K^2]

tau_hslist[list[float]], optional

f parameters used in calculating NRTL.taus, [1/K^2]

alpha_cslist[list[float]], optional

c parameters used in calculating NRTL.alphas, [-]

alpha_dslist[list[float]], optional

d paraemeters used in calculating NRTL.alphas, [1/K]

Notes

In addition to the methods presented here, the methods of its base class thermo.activity.GibbsExcess are available as well.

References

1

Poling, Bruce E., John M. Prausnitz, and John P. O`Connell. The Properties of Gases and Liquids. 5th edition. New York: McGraw-Hill Professional, 2000.

2(1,2)

Gmehling, Jürgen, Michael Kleiber, Bärbel Kolbe, and Jürgen Rarey. Chemical Thermodynamics for Process Simulation. John Wiley & Sons, 2019.

Examples

The DDBST has published numerous problems showing this model a simple binary system, Example P05.01b in [2], shows how to use parameters from the DDBST which are in units of calorie and need the gas constant as a multiplier:

>>> from scipy.constants import calorie, R
>>> N = 2
>>> T = 70.0 + 273.15
>>> xs = [0.252, 0.748]
>>> tausA = tausE = tausF = tausG = tausH = alphaD = [[0.0]*N for i in range(N)]
>>> tausB = [[0, -121.2691/R*calorie], [1337.8574/R*calorie, 0]]
>>> alphaC =  [[0, 0.2974],[.2974, 0]]
>>> ABEFGHCD = (tausA, tausB, tausE, tausF, tausG, tausH, alphaC, alphaD)
>>> GE = NRTL(T=T, xs=xs, ABEFGHCD=ABEFGHCD)
>>> GE.gammas()
[1.93605165145, 1.15366304520]
>>> GE
NRTL(T=343.15, xs=[0.252, 0.748], tau_bs=[[0, -61.0249799309399], [673.2359767282798, 0]], alpha_cs=[[0, 0.2974], [0.2974, 0]])
>>> GE.GE(), GE.dGE_dT(), GE.d2GE_dT2()
(780.053057219, 0.5743500022, -0.003584843605528)
>>> GE.HE(), GE.SE(), GE.dHE_dT(), GE.dSE_dT()
(582.964853938, -0.57435000227, 1.230139083237, 0.0035848436055)

The solution given by the DDBST has the same values [1.936, 1.154], and can be found here: http://chemthermo.ddbst.com/Problems_Solutions/Mathcad_Files/P05.01b%20VLE%20Behavior%20of%20Ethanol%20-%20Water%20Using%20NRTL.xps

Attributes
Tfloat

Temperature, [K]

xslist[float]

Mole fractions, [-]

Methods

GE()

Calculate and return the excess Gibbs energy of a liquid phase represented by the NRTL model.

Gs()

Calculates and return the G terms in the NRTL model for a specified temperature.

alphas()

Calculates and return the alpha terms in the NRTL model for a specified temperature.

d2GE_dT2()

Calculate and return the second tempreature derivative of excess Gibbs energy of a liquid phase represented by the NRTL model.

d2GE_dTdxs()

Calculate and return the temperature derivative of mole fraction derivatives of excess Gibbs energy of a liquid represented by the NRTL model.

d2GE_dxixjs()

Calculate and return the second mole fraction derivatives of excess Gibbs energy of a liquid represented by the NRTL model.

d2Gs_dT2()

Calculates and return the second temperature derivative of G terms in the NRTL model for a specified temperature.

d2taus_dT2()

Calculate and return the second temperature derivative of the tau terms for the NRTL model for a specified temperature.

d3Gs_dT3()

Calculates and return the third temperature derivative of G terms in the NRTL model for a specified temperature.

d3taus_dT3()

Calculate and return the third temperature derivative of the tau terms for the NRTL model for a specified temperature.

dGE_dT()

Calculate and return the first tempreature derivative of excess Gibbs energy of a liquid phase represented by the NRTL model.

dGE_dxs()

Calculate and return the mole fraction derivatives of excess Gibbs energy of a liquid represented by the NRTL model.

dGs_dT()

Calculates and return the first temperature derivative of G terms in the NRTL model for a specified temperature.

dtaus_dT()

Calculate and return the temperature derivative of the tau terms for the NRTL model for a specified temperature.

taus()

Calculate and return the tau terms for the NRTL model for a specified temperature.

to_T_xs(T, xs)

Method to construct a new NRTL instance at temperature T, and mole fractions xs with the same parameters as the existing object.

to_T_xs(T, xs)[source]

Method to construct a new NRTL instance at temperature T, and mole fractions xs with the same parameters as the existing object.

Parameters
Tfloat

Temperature, [K]

xslist[float]

Mole fractions of each component, [-]

Returns
objNRTL

New NRTL object at the specified conditions [-]

Notes

If the new temperature is the same temperature as the existing temperature, if the tau, Gs, or alphas terms or their derivatives have been calculated, they will be set to the new object as well.

taus()[source]

Calculate and return the tau terms for the NRTL model for a specified temperature.

τij=Aij+BijT+EijlnT+FijT+GijT2+HijT2\tau_{ij}=A_{ij}+\frac{B_{ij}}{T}+E_{ij}\ln T + F_{ij}T + \frac{G_{ij}}{T^2} + H_{ij}{T^2}
Returns
tauslist[list[float]]

tau terms, asymmetric matrix [-]

Notes

These tau ij values (and the coefficients) are NOT symmetric.

dtaus_dT()[source]

Calculate and return the temperature derivative of the tau terms for the NRTL model for a specified temperature.

τijTP,xi=BijT2+EijT+Fij2GijT3+2HijT\frac{\partial \tau_{ij}} {\partial T}_{P, x_i} = - \frac{B_{ij}}{T^{2}} + \frac{E_{ij}}{T} + F_{ij} - \frac{2 G_{ij}}{T^{3}} + 2 H_{ij} T
Returns
dtaus_dTlist[list[float]]

First temperature derivative of tau terms, asymmetric matrix [1/K]

d2taus_dT2()[source]

Calculate and return the second temperature derivative of the tau terms for the NRTL model for a specified temperature.

2τijT2P,xi=2BijT3EijT2+6GijT4+2Hij\frac{\partial^2 \tau_{ij}} {\partial T^2}_{P, x_i} = \frac{2 B_{ij}}{T^{3}} - \frac{E_{ij}}{T^{2}} + \frac{6 G_{ij}} {T^{4}} + 2 H_{ij}
Returns
d2taus_dT2list[list[float]]

Second temperature derivative of tau terms, asymmetric matrix [1/K^2]

d3taus_dT3()[source]

Calculate and return the third temperature derivative of the tau terms for the NRTL model for a specified temperature.

3τijT3P,xi=6BijT4+2EijT324GijT5\frac{\partial^3 \tau_{ij}} {\partial T^3}_{P, x_i} = - \frac{6 B_{ij}}{T^{4}} + \frac{2 E_{ij}}{T^{3}} - \frac{24 G_{ij}}{T^{5}}
Returns
d3taus_dT3list[list[float]]

Third temperature derivative of tau terms, asymmetric matrix [1/K^3]

alphas()[source]

Calculates and return the alpha terms in the NRTL model for a specified temperature.

αij=cij+dijT\alpha_{ij}=c_{ij} + d_{ij}T
Returns
alphaslist[list[float]]

alpha terms, possibly asymmetric matrix [-]

Notes

alpha values (and therefore cij and dij are normally symmetrical; but this is not strictly required.

Some sources suggest the c term should be fit to a given system; but the d term should be fit for an entire chemical family to avoid overfitting.

Recommended values for cij according to one source are:

0.30 Nonpolar substances with nonpolar substances; low deviation from ideality. 0.20 Hydrocarbons that are saturated interacting with polar liquids that do not associate, or systems that for multiple liquid phases which are immiscible 0.47 Strongly self associative systems, interacting with non-polar substances

alpha_coeffs should be a list[list[cij, dij]] so a 3d array

Gs()[source]

Calculates and return the G terms in the NRTL model for a specified temperature.

Gij=exp(αijτij)G_{ij}=\exp(-\alpha_{ij}\tau_{ij})
Returns
Gslist[list[float]]

G terms, asymmetric matrix [-]

dGs_dT()[source]

Calculates and return the first temperature derivative of G terms in the NRTL model for a specified temperature.

GijT=(α(T)ddTτ(T)τ(T)ddTα(T))eα(T)τ(T)\frac{\partial G_{ij}}{\partial T} = \left(- \alpha{\left(T \right)} \frac{d}{d T} \tau{\left(T \right)} - \tau{\left(T \right)} \frac{d}{d T} \alpha{\left(T \right)}\right) e^{- \alpha{\left(T \right)} \tau{\left(T \right)}}
Returns
dGs_dTlist[list[float]]

Temperature derivative of G terms, asymmetric matrix [1/K]

Notes

Derived with SymPy:

>>> from sympy import * 
>>> T = symbols('T') 
>>> alpha, tau = symbols('alpha, tau', cls=Function) 
>>> diff(exp(-alpha(T)*tau(T)), T) 
d2Gs_dT2()[source]

Calculates and return the second temperature derivative of G terms in the NRTL model for a specified temperature.

2GijT2=((α(T)ddTτ(T)+τ(T)ddTα(T))2α(T)d2dT2τ(T)2ddTα(T)ddTτ(T))eα(T)τ(T)\frac{\partial^2 G_{ij}}{\partial T^2} = \left(\left(\alpha{\left(T \right)} \frac{d}{d T} \tau{\left(T \right)} + \tau{\left(T \right)} \frac{d}{d T} \alpha{\left(T \right)}\right)^{2} - \alpha{\left(T \right)} \frac{d^{2}}{d T^{2}} \tau{\left(T \right)} - 2 \frac{d}{d T} \alpha{\left(T \right)} \frac{d}{d T} \tau{\left(T \right)}\right) e^{- \alpha{\left(T \right)} \tau{\left(T \right)}}
Returns
d2Gs_dT2list[list[float]]

Second temperature derivative of G terms, asymmetric matrix [1/K^2]

Notes

Derived with SymPy:

>>> from sympy import * 
>>> T = symbols('T') 
>>> alpha, tau = symbols('alpha, tau', cls=Function) 
>>> diff(exp(-alpha(T)*tau(T)), T, 2) 
d3Gs_dT3()[source]

Calculates and return the third temperature derivative of G terms in the NRTL model for a specified temperature.

3GijT3=(α(T)ddTτ(T)+τ(T)ddTα(T))3+(3α(T)ddTτ(T)+3τ(T)ddTα(T))(α(T)d2dT2τ(T)+2ddTα(T)ddTτ(T))α(T)d3dT3τ(T)3ddTα(T)d2dT2τ(T)\frac{\partial^3 G_{ij}}{\partial T^3} = \left(\alpha{\left(T \right)} \frac{d}{d T} \tau{\left(T \right)} + \tau{\left(T \right)} \frac{d}{d T} \alpha{\left(T \right)}\right)^{3} + \left(3 \alpha{\left(T \right)} \frac{d}{d T} \tau{\left(T \right)} + 3 \tau{\left(T \right)} \frac{d}{d T} \alpha{\left(T \right)} \right) \left(\alpha{\left(T \right)} \frac{d^{2}}{d T^{2}} \tau{\left(T \right)} + 2 \frac{d}{d T} \alpha{\left(T \right)} \frac{d}{d T} \tau{\left(T \right)}\right) - \alpha{\left(T \right)} \frac{d^{3}}{d T^{3}} \tau{\left(T \right)} - 3 \frac{d}{d T} \alpha{\left(T \right)} \frac{d^{2}}{d T^{2}} \tau{\left(T \right)}
Returns
d3Gs_dT3list[list[float]]

Third temperature derivative of G terms, asymmetric matrix [1/K^3]

Notes

Derived with SymPy:

>>> from sympy import * 
>>> T = symbols('T') 
>>> alpha, tau = symbols('alpha, tau', cls=Function) 
>>> diff(exp(-alpha(T)*tau(T)), T, 3) 
GE()[source]

Calculate and return the excess Gibbs energy of a liquid phase represented by the NRTL model.

gE=RTixijτjiGjixjjGjixjg^E = RT\sum_i x_i \frac{\sum_j \tau_{ji} G_{ji} x_j} {\sum_j G_{ji}x_j}
Returns
GEfloat

Excess Gibbs energy, [J/mol]

dGE_dT()[source]

Calculate and return the first tempreature derivative of excess Gibbs energy of a liquid phase represented by the NRTL model.

Returns
dGE_dTfloat

First temperature derivative of excess Gibbs energy, [J/(mol*K)]

d2GE_dT2()[source]

Calculate and return the second tempreature derivative of excess Gibbs energy of a liquid phase represented by the NRTL model.

Returns
d2GE_dT2float

Second temperature derivative of excess Gibbs energy, [J/(mol*K^2)]

dGE_dxs()[source]

Calculate and return the mole fraction derivatives of excess Gibbs energy of a liquid represented by the NRTL model.

gExi\frac{\partial g^E}{\partial x_i}
Returns
dGE_dxslist[float]

Mole fraction derivatives of excess Gibbs energy, [J/mol]

d2GE_dxixjs()[source]

Calculate and return the second mole fraction derivatives of excess Gibbs energy of a liquid represented by the NRTL model.

2gExixj=RT[+GijτijmxmGmj+GjiτjiijmxmGmi(mxmGmjτmj)Gij(mxmGmj)2(mxmGmiτmi)Gji(mxmGmi)2k(2xk(mxmτmkGmk)GikGjk(mxmGmk)3xkGikGjk(τjk+τik)(mxmGmk)2)]\frac{\partial^2 g^E}{\partial x_i \partial x_j} = RT\left[ + \frac{G_{ij}\tau_{ij}}{\sum_m x_m G_{mj}} + \frac{G_{ji}\tau_{jiij}}{\sum_m x_m G_{mi}} -\frac{(\sum_m x_m G_{mj}\tau_{mj})G_{ij}}{(\sum_m x_m G_{mj})^2} -\frac{(\sum_m x_m G_{mi}\tau_{mi})G_{ji}}{(\sum_m x_m G_{mi})^2} \sum_k \left(\frac{2x_k(\sum_m x_m \tau_{mk}G_{mk})G_{ik}G_{jk}}{(\sum_m x_m G_{mk})^3} - \frac{x_k G_{ik}G_{jk}(\tau_{jk} + \tau_{ik})}{(\sum_m x_m G_{mk})^2} \right) \right]
Returns
d2GE_dxixjslist[list[float]]

Second mole fraction derivatives of excess Gibbs energy, [J/mol]

d2GE_dTdxs()[source]

Calculate and return the temperature derivative of mole fraction derivatives of excess Gibbs energy of a liquid represented by the NRTL model.

2gExiT=R[T(j(xj(GijτijT+τijGijT)kxkGkj+xjGijτij(kxkGkjT)(kxkGkj)2+xjGijT(kxkGkjτkj)(kxkGkj)2+xjGij(kxk(GkjτkjT+τkjGkjT))(kxkGkj)22xjGij(kxkGkjT)(kxkGkjτkj)(kxkGkj)3)k(xkGkiτkiT+xkτkiGkiTkxkGki+(kxkGkiT)(kxkGkiτki)(kxkGki)2)+jxjGjiτjijxjGji+j(xjGij(kxkGkjτkj)(kxkGkj)2+xjGijτijkxkGkj)]\frac{\partial^2 g^E}{\partial x_i \partial T} = R\left[-T\left( \sum_j \left( -\frac{x_j(G_{ij}\frac{\partial \tau_{ij}}{\partial T} + \tau_{ij}\frac{\partial G_{ij}}{\partial T})}{\sum_k x_k G_{kj}} + \frac{x_j G_{ij}\tau_{ij}(\sum_k x_k \frac{\partial G_{kj}}{\partial T})}{(\sum_k x_k G_{kj})^2} +\frac{x_j \frac{\partial G_{ij}}{\partial T}(\sum_k x_k G_{kj}\tau_{kj})}{(\sum_k x_k G_{kj})^2} + \frac{x_jG_{ij}(\sum_k x_k (G_{kj} \frac{\partial \tau_{kj}}{\partial T} + \tau_{kj} \frac{\partial G_{kj}}{\partial T} ))}{(\sum_k x_k G_{kj})^2} -2\frac{x_j G_{ij} (\sum_k x_k \frac{\partial G_{kj}}{\partial T})(\sum_k x_k G_{kj}\tau_{kj})}{(\sum_k x_k G_{kj})^3} \right) - \frac{\sum_k (x_k G_{ki} \frac{\partial \tau_{ki}}{\partial T} + x_k \tau_{ki} \frac{\partial G_{ki}}{\partial T}} {\sum_k x_k G_{ki}} + \frac{(\sum_k x_k \frac{\partial G_{ki}}{\partial T})(\sum_k x_k G_{ki} \tau_{ki})}{(\sum_k x_k G_{ki})^2} \right) + \frac{\sum_j x_j G_{ji}\tau_{ji}}{\sum_j x_j G_{ji}} + \sum_j \left( \frac{x_j G_{ij}(\sum_k x_k G_{kj}\tau_{kj})}{(\sum_k x_k G_{kj})^2} + \frac{x_j G_{ij}\tau_{ij}}{\sum_k x_k G_{kj}} \right) \right]
Returns
d2GE_dTdxslist[float]

Temperature derivative of mole fraction derivatives of excess Gibbs energy, [J/(mol*K)]

NRTL Functional Calculations

thermo.nrtl.NRTL_gammas(xs, taus, alphas)[source]

Calculates the activity coefficients of each species in a mixture using the Non-Random Two-Liquid (NRTL) method, given their mole fractions, dimensionless interaction parameters, and nonrandomness constants. Those are normally correlated with temperature in some form, and need to be calculated separately.

ln(γi)=j=1nxjτjiGjik=1nxkGki+j=1nxjGijk=1nxkGkj(τijm=1nxmτmjGmjk=1nxkGkj)\ln(\gamma_i)=\frac{\displaystyle\sum_{j=1}^{n}{x_{j}\tau_{ji}G_{ji}}} {\displaystyle\sum_{k=1}^{n}{x_{k}G_{ki}}}+\sum_{j=1}^{n} {\frac{x_{j}G_{ij}}{\displaystyle\sum_{k=1}^{n}{x_{k}G_{kj}}}} {\left ({\tau_{ij}-\frac{\displaystyle\sum_{m=1}^{n}{x_{m}\tau_{mj} G_{mj}}}{\displaystyle\sum_{k=1}^{n}{x_{k}G_{kj}}}}\right )}
Gij=exp(αijτij)G_{ij}=\text{exp}\left ({-\alpha_{ij}\tau_{ij}}\right )
Parameters
xslist[float]

Liquid mole fractions of each species, [-]

tauslist[list[float]]

Dimensionless interaction parameters of each compound with each other, [-]

alphaslist[list[float]]

Nonrandomness constants of each compound interacting with each other, [-]

Returns
gammaslist[float]

Activity coefficient for each species in the liquid mixture, [-]

Notes

This model needs N^2 parameters.

One common temperature dependence of the nonrandomness constants is:

αij=cij+dijT\alpha_{ij}=c_{ij}+d_{ij}T

Most correlations for the interaction parameters include some of the terms shown in the following form:

τij=Aij+BijT+CijT2+Dijln(T)+EijTFij\tau_{ij}=A_{ij}+\frac{B_{ij}}{T}+\frac{C_{ij}}{T^{2}}+D_{ij} \ln{\left ({T}\right )}+E_{ij}T^{F_{ij}}

The original form of this model used the temperature dependence of taus in the form (values can be found in the literature, often with units of calories/mol):

τij=bijRT\tau_{ij}=\frac{b_{ij}}{RT}

For this model to produce ideal acitivty coefficients (gammas = 1), all interaction parameters should be 0; the value of alpha does not impact the calculation when that is the case.

References

1

Renon, Henri, and J. M. Prausnitz. “Local Compositions in Thermodynamic Excess Functions for Liquid Mixtures.” AIChE Journal 14, no. 1 (1968): 135-144. doi:10.1002/aic.690140124.

2

Gmehling, Jurgen, Barbel Kolbe, Michael Kleiber, and Jurgen Rarey. Chemical Thermodynamics for Process Simulation. 1st edition. Weinheim: Wiley-VCH, 2012.

Examples

Ethanol-water example, at 343.15 K and 1 MPa:

>>> NRTL_gammas(xs=[0.252, 0.748], taus=[[0, -0.178], [1.963, 0]],
... alphas=[[0, 0.2974],[.2974, 0]])
[1.9363183763514304, 1.1537609663170014]

NRTL Regression Calculations

thermo.nrtl.NRTL_gammas_binaries(xs, tau12, tau21, alpha12, alpha21, gammas=None)[source]

Calculates activity coefficients at fixed tau and alpha values for a binary system at a series of mole fractions. This is used for regression of tau and alpha parameters. This function is highly optimized, and operates on multiple points at a time.

lnγ1=x22[τ21(G21x1+x2G21)2+τ12G12(x2+x1G12)2]\ln \gamma_1 = x_2^2\left[\tau_{21}\left(\frac{G_{21}}{x_1 + x_2 G_{21}} \right)^2 + \frac{\tau_{12}G_{12}}{(x_2 + x_1G_{12})^2} \right]
lnγ2=x12[τ12(G12x2+x1G12)2+τ21G21(x1+x2G21)2]\ln \gamma_2 = x_1^2\left[\tau_{12}\left(\frac{G_{12}}{x_2 + x_1 G_{12}} \right)^2 + \frac{\tau_{21}G_{21}}{(x_1 + x_2 G_{21})^2} \right]
Gij=exp(αijτij)G_{ij}=\exp(-\alpha_{ij}\tau_{ij})
Parameters
xslist[float]

Liquid mole fractions of each species in the format x0_0, x1_0, (component 1 point1, component 2 point 1), x0_1, x1_1, (component 1 point2, component 2 point 2), … [-]

tau12float

tau parameter for 12, [-]

tau21float

tau parameter for 21, [-]

alpha12float

alpha parameter for 12, [-]

alpha21float

alpha parameter for 21, [-]

gammaslist[float], optional

Array to store the activity coefficient for each species in the liquid mixture, indexed the same as xs; can be omitted or provided for slightly better performance [-]

Returns
gammaslist[float]

Activity coefficient for each species in the liquid mixture, indexed the same as xs, [-]

Examples

>>> NRTL_gammas_binaries([.1, .9, 0.3, 0.7, .85, .15], 0.1759, 0.7991, .2, .3)
[2.121421, 1.011342, 1.52177, 1.09773, 1.016062, 1.841391]