Skip to content

Kinetics (polykin.kinetics)¤

This module implements methods to create and visualize the types of kinetic coefficients most often used in polymer reactor models.

Arrhenius ¤

Arrhenius kinetic rate coefficient.

This class implements the following temperature dependence:

\[ k(T)=k_0 \exp\left[-\frac{E_a}{R}\left(\frac{1}{T}-\frac{1}{T_0}\right)\right] \]

where \(T_0\) is a reference temperature, \(E_a\) is the activation energy, and \(k_0=k(T_0)\). In the limit \(T\rightarrow+\infty\), the usual form of the Arrhenius equation with \(k_0=A\) is recovered.

PARAMETER DESCRIPTION
k0

Coefficient value at the reference temperature, \(k_0=k(T_0)\). Unit = unit.

TYPE: float | FloatArrayLike

EaR

Energy of activation, \(E_a/R\). Unit = K.

TYPE: float | FloatArrayLike

T0

Reference temperature, \(T_0\). Unit = K.

TYPE: float | FloatArrayLike DEFAULT: inf

Tmin

Lower temperature bound. Unit = K.

TYPE: float | FloatArrayLike DEFAULT: 0.0

Tmax

Upper temperature bound. Unit = K.

TYPE: float | FloatArrayLike DEFAULT: inf

unit

Unit of coefficient.

TYPE: str DEFAULT: '-'

symbol

Symbol of coefficient \(k\).

TYPE: str DEFAULT: 'k'

name

Name.

TYPE: str DEFAULT: ''

See also

Examples:

Define and evaluate the propagation rate coefficient of styrene.

>>> from polykin.kinetics import Arrhenius 
>>> kp = Arrhenius(
...     10**7.63,       # pre-exponential factor
...     32.5e3/8.314,   # Ea/R, K
...     Tmin=261., Tmax=366.,
...     symbol='k_p',
...     unit='L/mol/s',
...     name='kp of styrene')
>>> kp(25.,'C') 
86.28385101961442
Source code in src/polykin/kinetics/arrhenius.py
 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
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
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
class Arrhenius(KineticCoefficientT):
    r"""[Arrhenius](https://en.wikipedia.org/wiki/Arrhenius_equation) kinetic
    rate coefficient.

    This class implements the following temperature dependence:

    $$
    k(T)=k_0
    \exp\left[-\frac{E_a}{R}\left(\frac{1}{T}-\frac{1}{T_0}\right)\right]
    $$

    where $T_0$ is a reference temperature, $E_a$ is the activation energy,
    and $k_0=k(T_0)$. In the limit $T\rightarrow+\infty$, the usual
    form of the Arrhenius equation with $k_0=A$ is recovered.

    Parameters
    ----------
    k0 : float | FloatArrayLike
        Coefficient value at the reference temperature, $k_0=k(T_0)$.
        Unit = `unit`.
    EaR : float | FloatArrayLike
        Energy of activation, $E_a/R$.
        Unit = K.
    T0 : float | FloatArrayLike
        Reference temperature, $T_0$.
        Unit = K.
    Tmin : float | FloatArrayLike
        Lower temperature bound.
        Unit = K.
    Tmax : float | FloatArrayLike
        Upper temperature bound.
        Unit = K.
    unit : str
        Unit of coefficient.
    symbol : str
        Symbol of coefficient $k$.
    name : str
        Name.

    See also
    --------
    * [`Eyring`](Eyring.md): alternative method.

    Examples
    --------
    Define and evaluate the propagation rate coefficient of styrene.
    >>> from polykin.kinetics import Arrhenius 
    >>> kp = Arrhenius(
    ...     10**7.63,       # pre-exponential factor
    ...     32.5e3/8.314,   # Ea/R, K
    ...     Tmin=261., Tmax=366.,
    ...     symbol='k_p',
    ...     unit='L/mol/s',
    ...     name='kp of styrene')
    >>> kp(25.,'C') 
    86.28385101961442
    """

    _pinfo = {'k0': ('#', True), 'EaR': ('K', True), 'T0': ('K', False)}

    def __init__(self,
                 k0: Union[float, FloatArrayLike],
                 EaR: Union[float, FloatArrayLike],
                 T0: Union[float, FloatArrayLike] = np.inf,
                 Tmin: Union[float, FloatArrayLike] = 0.0,
                 Tmax: Union[float, FloatArrayLike] = np.inf,
                 unit: str = '-',
                 symbol: str = 'k',
                 name: str = ''
                 ) -> None:

        # Convert lists to arrays
        k0, EaR, T0, Tmin, Tmax = \
            convert_FloatOrArrayLike_to_FloatOrArray([k0, EaR, T0, Tmin, Tmax])

        # Check shapes
        self._shape = check_shapes([k0, EaR], [T0, Tmin, Tmax])

        # Check bounds
        check_bounds(k0, 0., np.inf, 'k0')
        check_bounds(EaR, -np.inf, np.inf, 'EaR')
        check_bounds(T0, 0., np.inf, 'T0')

        self.p = {'k0': k0, 'EaR': EaR, 'T0': T0}
        super().__init__((Tmin, Tmax), unit, symbol, name)

    @staticmethod
    def equation(T: Union[float, FloatArray],
                 k0: Union[float, FloatArray],
                 EaR: Union[float, FloatArray],
                 T0: Union[float, FloatArray],
                 ) -> Union[float, FloatArray]:
        r"""Arrhenius equation.

        Parameters
        ----------
        T : float | FloatArray
            Temperature.
            Unit = K.
        k0 : float | FloatArray
            Coefficient value at the reference temperature, $k_0=k(T_0)$.
            Unit = Any.
        EaR : float | FloatArray
            Energy of activation, $E_a/R$.
            Unit = K.
        T0 : float | FloatArray
            Reference temperature, $T_0$.
            Unit = K.

        Returns
        -------
        float | FloatArray
            Coefficient value. Unit = [k0].
        """
        return k0 * exp(-EaR*(1/T - 1/T0))

    @property
    def A(self) -> Union[float, FloatArray]:
        r"""Pre-exponential factor, $A=k_0 e^{E_a/(R T_0)}$."""
        return self.__call__(np.inf)

    def __mul__(self,
                other: Union[int, float, Arrhenius]
                ) -> Arrhenius:
        """Multipy Arrhenius coefficient(s).

        Create a new Arrhenius coefficient from a product of two Arrhenius
        coefficients with identical shapes or a product of an Arrhenius
        coefficient and a numerical constant.

        Parameters
        ----------
        other : int | float | Arrhenius
            Another Arrhenius coefficient or number.

        Returns
        -------
        Arrhenius
            Product coefficient.
        """
        if isinstance(other, Arrhenius):
            if self._shape == other._shape:
                return Arrhenius(k0=self.A*other.A,
                                 EaR=self.p['EaR'] + other.p['EaR'],
                                 Tmin=np.maximum(
                                     self.Trange[0], other.Trange[0]),
                                 Tmax=np.minimum(
                                     self.Trange[1], other.Trange[1]),
                                 unit=f"{self.unit}·{other.unit}",
                                 symbol=f"{self.symbol}·{other.symbol}",
                                 name=f"{self.name}·{other.name}")
            else:
                raise ShapeError(
                    "Product of array-like coefficients requires identical shapes.")  # noqa: E501
        elif isinstance(other, (int, float)):
            return Arrhenius(k0=self.p['k0']*other,
                             EaR=self.p['EaR'],
                             T0=self.p['T0'],
                             Tmin=self.Trange[0],
                             Tmax=self.Trange[1],
                             unit=self.unit,
                             symbol=f"{str(other)}·{self.symbol}",
                             name=f"{str(other)}·{self.name}")
        else:
            return NotImplemented

    def __rmul__(self,
                 other: Union[int, float, Arrhenius]
                 ) -> Arrhenius:
        return self.__mul__(other)

    def __truediv__(self,
                    other: Union[int, float, Arrhenius]
                    ) -> Arrhenius:
        """Divide Arrhenius coefficient(s).

        Create a new Arrhenius coefficient from a division of two Arrhenius
        coefficients with identical shapes or a division involving an Arrhenius
        coefficient and a numerical constant.

        Parameters
        ----------
        other : int | float | Arrhenius
            Another Arrhenius coefficient or number.

        Returns
        -------
        Arrhenius
            Quotient coefficient.
        """
        if isinstance(other, Arrhenius):
            if self._shape == other._shape:
                return Arrhenius(k0=self.A/other.A,
                                 EaR=self.p['EaR'] - other.p['EaR'],
                                 Tmin=np.maximum(
                                     self.Trange[0], other.Trange[0]),
                                 Tmax=np.minimum(
                                     self.Trange[1], other.Trange[1]),
                                 unit=f"{self.unit}/{other.unit}",
                                 symbol=f"{self.symbol}/{other.symbol}",
                                 name=f"{self.name}/{other.name}")
            else:
                raise ShapeError(
                    "Division of array-like coefficients requires identical shapes.")  # noqa: E501
        elif isinstance(other, (int, float)):
            return Arrhenius(k0=self.p['k0']/other,
                             EaR=self.p['EaR'],
                             T0=self.p['T0'],
                             Tmin=self.Trange[0],
                             Tmax=self.Trange[1],
                             unit=self.unit,
                             symbol=f"{self.symbol}/{str(other)}",
                             name=f"{self.name}/{str(other)}")
        else:
            return NotImplemented

    def __rtruediv__(self,
                     other: Union[int, float]
                     ) -> Arrhenius:
        if isinstance(other, (int, float)):
            return Arrhenius(k0=other/self.p['k0'],
                             EaR=-self.p['EaR'],
                             T0=self.p['T0'],
                             Tmin=self.Trange[0],
                             Tmax=self.Trange[1],
                             unit=f"1/{self.unit}",
                             symbol=f"{str(other)}/{self.symbol}",
                             name=f"{str(other)}/{self.name}")
        else:
            return NotImplemented

    def __pow__(self,
                other: Union[int, float]
                ) -> Arrhenius:
        """Power of an Arrhenius coefficient.

        Create a new Arrhenius coefficient from the exponentiation of an
        Arrhenius coefficient.

        Parameters
        ----------
        other : int | float
            Exponent.

        Returns
        -------
        Arrhenius
            Power coefficient.
        """
        if isinstance(other, (int, float)):
            return Arrhenius(k0=self.p['k0']**other,
                             EaR=self.p['EaR']*other,
                             T0=self.p['T0'],
                             Tmin=self.Trange[0],
                             Tmax=self.Trange[1],
                             unit=f"({self.unit})^{str(other)}",
                             symbol=f"({self.symbol})^{str(other)}",
                             name=f"{self.name}^{str(other)}"
                             )
        else:
            return NotImplemented

    def __rpow__(self, other):
        return NotImplemented

A property ¤

A: Union[float, FloatArray]

Pre-exponential factor, \(A=k_0 e^{E_a/(R T_0)}\).

__call__ ¤

__call__(
    T: Union[float, FloatArrayLike],
    Tunit: Literal["C", "K"] = "K",
) -> Union[float, FloatArray]

Evaluate property equation at given temperature, including unit conversion and range check.

PARAMETER DESCRIPTION
T

Temperature. Unit = Tunit.

TYPE: float | FloatArrayLike

Tunit

Temperature unit.

TYPE: Literal['C', 'K'] DEFAULT: 'K'

RETURNS DESCRIPTION
float | FloatArray

Correlation value.

Source code in src/polykin/properties/equations/base.py
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
def __call__(self,
             T: Union[float, FloatArrayLike],
             Tunit: Literal['C', 'K'] = 'K'
             ) -> Union[float, FloatArray]:
    r"""Evaluate property equation at given temperature, including unit
    conversion and range check.

    Parameters
    ----------
    T : float | FloatArrayLike
        Temperature.
        Unit = `Tunit`.
    Tunit : Literal['C', 'K']
        Temperature unit.

    Returns
    -------
    float | FloatArray
        Correlation value.
    """
    TK = convert_check_temperature(T, Tunit, self.Trange)
    return self.equation(TK, **self.p)

__mul__ ¤

__mul__(other: Union[int, float, Arrhenius]) -> Arrhenius

Multipy Arrhenius coefficient(s).

Create a new Arrhenius coefficient from a product of two Arrhenius coefficients with identical shapes or a product of an Arrhenius coefficient and a numerical constant.

PARAMETER DESCRIPTION
other

Another Arrhenius coefficient or number.

TYPE: int | float | Arrhenius

RETURNS DESCRIPTION
Arrhenius

Product coefficient.

Source code in src/polykin/kinetics/arrhenius.py
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
def __mul__(self,
            other: Union[int, float, Arrhenius]
            ) -> Arrhenius:
    """Multipy Arrhenius coefficient(s).

    Create a new Arrhenius coefficient from a product of two Arrhenius
    coefficients with identical shapes or a product of an Arrhenius
    coefficient and a numerical constant.

    Parameters
    ----------
    other : int | float | Arrhenius
        Another Arrhenius coefficient or number.

    Returns
    -------
    Arrhenius
        Product coefficient.
    """
    if isinstance(other, Arrhenius):
        if self._shape == other._shape:
            return Arrhenius(k0=self.A*other.A,
                             EaR=self.p['EaR'] + other.p['EaR'],
                             Tmin=np.maximum(
                                 self.Trange[0], other.Trange[0]),
                             Tmax=np.minimum(
                                 self.Trange[1], other.Trange[1]),
                             unit=f"{self.unit}·{other.unit}",
                             symbol=f"{self.symbol}·{other.symbol}",
                             name=f"{self.name}·{other.name}")
        else:
            raise ShapeError(
                "Product of array-like coefficients requires identical shapes.")  # noqa: E501
    elif isinstance(other, (int, float)):
        return Arrhenius(k0=self.p['k0']*other,
                         EaR=self.p['EaR'],
                         T0=self.p['T0'],
                         Tmin=self.Trange[0],
                         Tmax=self.Trange[1],
                         unit=self.unit,
                         symbol=f"{str(other)}·{self.symbol}",
                         name=f"{str(other)}·{self.name}")
    else:
        return NotImplemented

__pow__ ¤

__pow__(other: Union[int, float]) -> Arrhenius

Power of an Arrhenius coefficient.

Create a new Arrhenius coefficient from the exponentiation of an Arrhenius coefficient.

PARAMETER DESCRIPTION
other

Exponent.

TYPE: int | float

RETURNS DESCRIPTION
Arrhenius

Power coefficient.

Source code in src/polykin/kinetics/arrhenius.py
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
def __pow__(self,
            other: Union[int, float]
            ) -> Arrhenius:
    """Power of an Arrhenius coefficient.

    Create a new Arrhenius coefficient from the exponentiation of an
    Arrhenius coefficient.

    Parameters
    ----------
    other : int | float
        Exponent.

    Returns
    -------
    Arrhenius
        Power coefficient.
    """
    if isinstance(other, (int, float)):
        return Arrhenius(k0=self.p['k0']**other,
                         EaR=self.p['EaR']*other,
                         T0=self.p['T0'],
                         Tmin=self.Trange[0],
                         Tmax=self.Trange[1],
                         unit=f"({self.unit})^{str(other)}",
                         symbol=f"({self.symbol})^{str(other)}",
                         name=f"{self.name}^{str(other)}"
                         )
    else:
        return NotImplemented

__truediv__ ¤

__truediv__(
    other: Union[int, float, Arrhenius]
) -> Arrhenius

Divide Arrhenius coefficient(s).

Create a new Arrhenius coefficient from a division of two Arrhenius coefficients with identical shapes or a division involving an Arrhenius coefficient and a numerical constant.

PARAMETER DESCRIPTION
other

Another Arrhenius coefficient or number.

TYPE: int | float | Arrhenius

RETURNS DESCRIPTION
Arrhenius

Quotient coefficient.

Source code in src/polykin/kinetics/arrhenius.py
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
def __truediv__(self,
                other: Union[int, float, Arrhenius]
                ) -> Arrhenius:
    """Divide Arrhenius coefficient(s).

    Create a new Arrhenius coefficient from a division of two Arrhenius
    coefficients with identical shapes or a division involving an Arrhenius
    coefficient and a numerical constant.

    Parameters
    ----------
    other : int | float | Arrhenius
        Another Arrhenius coefficient or number.

    Returns
    -------
    Arrhenius
        Quotient coefficient.
    """
    if isinstance(other, Arrhenius):
        if self._shape == other._shape:
            return Arrhenius(k0=self.A/other.A,
                             EaR=self.p['EaR'] - other.p['EaR'],
                             Tmin=np.maximum(
                                 self.Trange[0], other.Trange[0]),
                             Tmax=np.minimum(
                                 self.Trange[1], other.Trange[1]),
                             unit=f"{self.unit}/{other.unit}",
                             symbol=f"{self.symbol}/{other.symbol}",
                             name=f"{self.name}/{other.name}")
        else:
            raise ShapeError(
                "Division of array-like coefficients requires identical shapes.")  # noqa: E501
    elif isinstance(other, (int, float)):
        return Arrhenius(k0=self.p['k0']/other,
                         EaR=self.p['EaR'],
                         T0=self.p['T0'],
                         Tmin=self.Trange[0],
                         Tmax=self.Trange[1],
                         unit=self.unit,
                         symbol=f"{self.symbol}/{str(other)}",
                         name=f"{self.name}/{str(other)}")
    else:
        return NotImplemented

equation staticmethod ¤

equation(
    T: Union[float, FloatArray],
    k0: Union[float, FloatArray],
    EaR: Union[float, FloatArray],
    T0: Union[float, FloatArray],
) -> Union[float, FloatArray]

Arrhenius equation.

PARAMETER DESCRIPTION
T

Temperature. Unit = K.

TYPE: float | FloatArray

k0

Coefficient value at the reference temperature, \(k_0=k(T_0)\). Unit = Any.

TYPE: float | FloatArray

EaR

Energy of activation, \(E_a/R\). Unit = K.

TYPE: float | FloatArray

T0

Reference temperature, \(T_0\). Unit = K.

TYPE: float | FloatArray

RETURNS DESCRIPTION
float | FloatArray

Coefficient value. Unit = [k0].

Source code in src/polykin/kinetics/arrhenius.py
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
@staticmethod
def equation(T: Union[float, FloatArray],
             k0: Union[float, FloatArray],
             EaR: Union[float, FloatArray],
             T0: Union[float, FloatArray],
             ) -> Union[float, FloatArray]:
    r"""Arrhenius equation.

    Parameters
    ----------
    T : float | FloatArray
        Temperature.
        Unit = K.
    k0 : float | FloatArray
        Coefficient value at the reference temperature, $k_0=k(T_0)$.
        Unit = Any.
    EaR : float | FloatArray
        Energy of activation, $E_a/R$.
        Unit = K.
    T0 : float | FloatArray
        Reference temperature, $T_0$.
        Unit = K.

    Returns
    -------
    float | FloatArray
        Coefficient value. Unit = [k0].
    """
    return k0 * exp(-EaR*(1/T - 1/T0))

fit ¤

fit(
    T: FloatVectorLike,
    Y: FloatVectorLike,
    sigmaY: Optional[FloatVectorLike] = None,
    fit_only: Optional[list[str]] = None,
    logY: bool = False,
    plot: bool = True,
) -> dict

Fit equation to data using non-linear regression.

PARAMETER DESCRIPTION
T

Temperature. Unit = K.

TYPE: FloatVector

Y

Property to be fitted. Unit = Any.

TYPE: FloatVector

sigmaY

Standard deviation of Y. Unit = [Y].

TYPE: FloatVector | None DEFAULT: None

fit_only

List with name of parameters to be fitted.

TYPE: list[str] | None DEFAULT: None

logY

If True, the fit will be done in terms of log(Y).

TYPE: bool DEFAULT: False

plot

If True a plot comparing data and fitted correlation will be generated.

TYPE: bool DEFAULT: True

RETURNS DESCRIPTION
dict

A dictionary of results with the following keys: 'success', 'parameters', 'covariance', and 'plot'.

Source code in src/polykin/properties/equations/base.py
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
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
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
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
def fit(self,
        T: FloatVectorLike,
        Y: FloatVectorLike,
        sigmaY: Optional[FloatVectorLike] = None,
        fit_only: Optional[list[str]] = None,
        logY: bool = False,
        plot: bool = True,
        ) -> dict:
    """Fit equation to data using non-linear regression.

    Parameters
    ----------
    T : FloatVector
        Temperature. Unit = K.
    Y : FloatVector
        Property to be fitted. Unit = Any.
    sigmaY : FloatVector | None
        Standard deviation of Y. Unit = [Y].
    fit_only : list[str] | None
        List with name of parameters to be fitted.
    logY : bool
        If `True`, the fit will be done in terms of log(Y).
    plot : bool
        If `True` a plot comparing data and fitted correlation will be
        generated.

    Returns
    -------
    dict
        A dictionary of results with the following keys: 'success',
        'parameters', 'covariance', and 'plot'.
    """

    # Current parameter values
    pdict = self.p.copy()

    # Select parameters to be fitted
    pnames_fit = [name for name, info in self._pinfo.items() if info[1]]
    if fit_only:
        pnames_fit = set(fit_only) & set(pnames_fit)
    p0 = [pdict[pname] for pname in pnames_fit]

    # Fit function
    def ffit(x, *p):
        for pname, pvalue in zip(pnames_fit, p):
            pdict[pname] = pvalue
        Yfit = self.equation(T=x, **pdict)
        if logY:
            Yfit = log(Yfit)
        return Yfit

    solution = curve_fit(ffit,
                         xdata=T,
                         ydata=log(Y) if logY else Y,
                         p0=p0,
                         sigma=sigmaY,
                         absolute_sigma=False,
                         full_output=True)
    result = {}
    result['success'] = bool(solution[4])
    if solution[4]:
        popt = solution[0]
        pcov = solution[1]
        print("Fit successful.")
        for pname, pvalue in zip(pnames_fit, popt):
            print(f"{pname}: {pvalue}")
        print("Covariance:")
        print(pcov)
        result['covariance'] = pcov

        # Update attributes
        self.Trange = (min(T), max(T))
        for pname, pvalue in zip(pnames_fit, popt):
            self.p[pname] = pvalue
        result['parameters'] = pdict

        # plot
        if plot:
            kind = 'semilogy' if logY else 'linear'
            fig, ax = self.plot(kind=kind, return_objects=True)  # ok
            ax.plot(T, Y, 'o', mfc='none')
            result['plot'] = (fig, ax)
    else:
        print("Fit error: ", solution[3])
        result['message'] = solution[3]

    return result

plot ¤

plot(
    kind: Literal[
        "linear", "semilogy", "Arrhenius"
    ] = "linear",
    Trange: Optional[tuple[float, float]] = None,
    Tunit: Literal["C", "K"] = "K",
    title: Optional[str] = None,
    axes: Optional[Axes] = None,
    return_objects: bool = False,
) -> Optional[tuple[Optional[Figure], Axes]]

Plot quantity as a function of temperature.

PARAMETER DESCRIPTION
kind

Kind of plot to be generated.

TYPE: Literal['linear', 'semilogy', 'Arrhenius'] DEFAULT: 'linear'

Trange

Temperature range for x-axis. If None, the validity range (Tmin, Tmax) will be used. If no validity range was defined, the range will default to 0-100°C.

TYPE: tuple[float, float] | None DEFAULT: None

Tunit

Temperature unit.

TYPE: Literal['C', 'K'] DEFAULT: 'K'

title

Title of plot. If None, the object name will be used.

TYPE: str | None DEFAULT: None

axes

Matplotlib Axes object.

TYPE: Axes | None DEFAULT: None

return_objects

If True, the Figure and Axes objects are returned (for saving or further manipulations).

TYPE: bool DEFAULT: False

RETURNS DESCRIPTION
tuple[Figure | None, Axes] | None

Figure and Axes objects if return_objects is True.

Source code in src/polykin/properties/equations/base.py
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
def plot(self,
         kind: Literal['linear', 'semilogy', 'Arrhenius'] = 'linear',
         Trange: Optional[tuple[float, float]] = None,
         Tunit: Literal['C', 'K'] = 'K',
         title: Optional[str] = None,
         axes: Optional[Axes] = None,
         return_objects: bool = False
         ) -> Optional[tuple[Optional[Figure], Axes]]:
    """Plot quantity as a function of temperature.

    Parameters
    ----------
    kind : Literal['linear', 'semilogy', 'Arrhenius']
        Kind of plot to be generated.
    Trange : tuple[float, float] | None
        Temperature range for x-axis. If `None`, the validity range
        (Tmin, Tmax) will be used. If no validity range was defined, the
        range will default to 0-100°C.
    Tunit : Literal['C', 'K']
        Temperature unit.
    title : str | None
        Title of plot. If `None`, the object name will be used.
    axes : Axes | None
        Matplotlib Axes object.
    return_objects : bool
        If `True`, the Figure and Axes objects are returned (for saving or
        further manipulations).

    Returns
    -------
    tuple[Figure | None, Axes] | None
        Figure and Axes objects if return_objects is `True`.
    """

    # Check inputs
    check_in_set(kind, {'linear', 'semilogy', 'Arrhenius'}, 'kind')
    check_in_set(Tunit, {'K', 'C'}, 'Tunit')
    if Trange is not None:
        Trange_min = 0.
        if Tunit == 'C':
            Trange_min = -273.15
        check_valid_range(Trange, Trange_min, np.inf, 'Trange')

    # Plot objects
    if axes is None:
        fig, ax = plt.subplots()
        if title is None:
            title = self.name
        if title:
            fig.suptitle(title)
        label = None
    else:
        fig = None
        ax = axes
        label = self.name

    # Units and xlabel
    Tunit_range = Tunit
    if kind == 'Arrhenius':
        Tunit = 'K'
    Tsymbol = Tunit
    if Tunit == 'C':
        Tsymbol = '°' + Tunit

    if kind == 'Arrhenius':
        xlabel = r"$1/T$ [" + Tsymbol + r"$^{-1}$]"
    else:
        xlabel = fr"$T$ [{Tsymbol}]"

    # ylabel
    ylabel = fr"${self.symbol}$ [{self.unit}]"
    if axes is not None:
        ylabel0 = ax.get_ylabel()
        if ylabel0 and ylabel not in ylabel0:
            ylabel = ylabel0 + ", " + ylabel

    ax.set_xlabel(xlabel)
    ax.set_ylabel(ylabel)
    ax.grid(True)

    # x-axis vector
    if Trange is not None:
        if Tunit_range == 'C':
            Trange = (Trange[0]+273.15, Trange[1]+273.15)
    else:
        Trange = (np.min(self.Trange[0]), np.max(self.Trange[1]))
        if Trange == (0.0, np.inf):
            Trange = (273.15, 373.15)

    try:
        shape = self._shape
    except AttributeError:
        shape = None
    if shape is not None:
        print("Plot method not yet implemented for array-like equations.")
    else:
        TK = np.linspace(*Trange, 100)
        y = self.__call__(TK, 'K')
        T = TK
        if Tunit == 'C':
            T -= 273.15
        if kind == 'linear':
            ax.plot(T, y, label=label)
        elif kind == 'semilogy':
            ax.semilogy(T, y, label=label)
        elif kind == 'Arrhenius':
            ax.semilogy(1/TK, y, label=label)

    if fig is None:
        ax.legend(bbox_to_anchor=(1.05, 1.0), loc="upper left")

    if return_objects:
        return (fig, ax)

shape property ¤

shape: Optional[tuple[int, ...]]

Shape of underlying parameter array.

Eyring ¤

Eyring kinetic rate coefficient.

This class implements the following temperature dependence:

\[ k(T) = \dfrac{\kappa k_B T}{h} \exp\left(\frac{\Delta S^\ddagger}{R}\right) \exp\left(-\frac{\Delta H^\ddagger}{R T}\right) \]

where \(\kappa\) is the transmission coefficient, \(\Delta S^\ddagger\) is the entropy of activation, and \(\Delta H^\ddagger\) is the enthalpy of activation. The unit of \(k\) is physically set to s\(^{-1}\).

PARAMETER DESCRIPTION
DSa

Entropy of activation, \(\Delta S^\ddagger\). Unit = J/(mol·K).

TYPE: float | FloatArrayLike

DHa

Enthalpy of activation, \(\Delta H^\ddagger\). Unit = J/mol.

TYPE: float | FloatArrayLike

kappa

Transmission coefficient.

TYPE: float | FloatArrayLike DEFAULT: 1.0

Tmin

Lower temperature bound. Unit = K.

TYPE: float | FloatArrayLike DEFAULT: 0.0

Tmax

Upper temperature bound. Unit = K.

TYPE: float | FloatArrayLike DEFAULT: inf

symbol

Symbol of coefficient \(k\).

TYPE: str DEFAULT: 'k'

name

Name.

TYPE: str DEFAULT: ''

See also

Examples:

Define and evaluate a rate coefficient from transition state properties.

>>> from polykin.kinetics import Eyring 
>>> k = Eyring(
...     DSa=20.,
...     DHa=5e4,
...     kappa=0.8,
...     Tmin=273., Tmax=373.,
...     symbol='k',
...     name='A->B')
>>> k(25.,'C')
95808.36742009166
Source code in src/polykin/kinetics/eyring.py
 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
122
123
124
125
126
127
128
129
130
131
132
133
134
class Eyring(KineticCoefficientT):
    r"""[Eyring](https://en.wikipedia.org/wiki/Eyring_equation) kinetic rate
    coefficient.

    This class implements the following temperature dependence:

    $$
    k(T) = \dfrac{\kappa k_B T}{h}
           \exp\left(\frac{\Delta S^\ddagger}{R}\right)
           \exp\left(-\frac{\Delta H^\ddagger}{R T}\right)
    $$

    where $\kappa$ is the transmission coefficient, $\Delta S^\ddagger$ is
    the entropy of activation, and $\Delta H^\ddagger$ is the enthalpy of
    activation. The unit of $k$ is physically set to s$^{-1}$.

    Parameters
    ----------
    DSa : float | FloatArrayLike
        Entropy of activation, $\Delta S^\ddagger$.
        Unit = J/(mol·K).
    DHa : float | FloatArrayLike
        Enthalpy of activation, $\Delta H^\ddagger$.
        Unit = J/mol.
    kappa : float | FloatArrayLike
        Transmission coefficient.
    Tmin : float | FloatArrayLike
        Lower temperature bound.
        Unit = K.
    Tmax : float | FloatArrayLike
        Upper temperature bound.
        Unit = K.
    symbol : str
        Symbol of coefficient $k$.
    name : str
        Name.

    See also
    --------
    * [`Arrhenius`](Arrhenius.md): alternative method.

    Examples
    --------
    Define and evaluate a rate coefficient from transition state properties.
    >>> from polykin.kinetics import Eyring 
    >>> k = Eyring(
    ...     DSa=20.,
    ...     DHa=5e4,
    ...     kappa=0.8,
    ...     Tmin=273., Tmax=373.,
    ...     symbol='k',
    ...     name='A->B')
    >>> k(25.,'C')
    95808.36742009166
    """

    _pinfo = {'DSa': ('J/(mol·K)', True),
              'DHa': ('J/mol', True), 'kappa': ('', False)}

    def __init__(self,
                 DSa: Union[float, FloatArrayLike],
                 DHa: Union[float, FloatArrayLike],
                 kappa: Union[float, FloatArrayLike] = 1.0,
                 Tmin: Union[float, FloatArrayLike] = 0.0,
                 Tmax: Union[float, FloatArrayLike] = np.inf,
                 symbol: str = 'k',
                 name: str = ''
                 ) -> None:

        # Convert lists to arrays
        DSa, DHa, kappa, Tmin, Tmax = \
            convert_FloatOrArrayLike_to_FloatOrArray(
                [DSa, DHa, kappa, Tmin, Tmax])

        # Check shapes
        self._shape = check_shapes([DSa, DHa], [kappa, Tmin, Tmax])

        # Check bounds
        check_bounds(DSa, 0., np.inf, 'DSa')
        check_bounds(DHa, 0., np.inf, 'DHa')
        check_bounds(kappa, 0., 1., 'kappa')

        self.p = {'DSa': DSa, 'DHa': DHa, 'kappa': kappa}
        super().__init__((Tmin, Tmax), '1/s', symbol, name)

    @staticmethod
    def equation(T: Union[float, FloatArray],
                 DSa: Union[float, FloatArray],
                 DHa: Union[float, FloatArray],
                 kappa: Union[float, FloatArray],
                 ) -> Union[float, FloatArray]:
        r"""Eyring equation.

        Parameters
        ----------
        T : float | FloatArray
            Temperature.
            Unit = K.
        DSa : float | FloatArray
            Entropy of activation, $\Delta S^\ddagger$.
            Unit = J/(mol·K).
        DHa : float | FloatArray
            Enthalpy of activation, $\Delta H^\ddagger$.
            Unit = J/mol.
        kappa : float | FloatArray
            Transmission coefficient.

        Returns
        -------
        float | FloatArray
            Coefficient value. Unit = 1/s.
        """
        return kappa * kB*T/h * exp((DSa - DHa/T)/R)

__call__ ¤

__call__(
    T: Union[float, FloatArrayLike],
    Tunit: Literal["C", "K"] = "K",
) -> Union[float, FloatArray]

Evaluate property equation at given temperature, including unit conversion and range check.

PARAMETER DESCRIPTION
T

Temperature. Unit = Tunit.

TYPE: float | FloatArrayLike

Tunit

Temperature unit.

TYPE: Literal['C', 'K'] DEFAULT: 'K'

RETURNS DESCRIPTION
float | FloatArray

Correlation value.

Source code in src/polykin/properties/equations/base.py
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
def __call__(self,
             T: Union[float, FloatArrayLike],
             Tunit: Literal['C', 'K'] = 'K'
             ) -> Union[float, FloatArray]:
    r"""Evaluate property equation at given temperature, including unit
    conversion and range check.

    Parameters
    ----------
    T : float | FloatArrayLike
        Temperature.
        Unit = `Tunit`.
    Tunit : Literal['C', 'K']
        Temperature unit.

    Returns
    -------
    float | FloatArray
        Correlation value.
    """
    TK = convert_check_temperature(T, Tunit, self.Trange)
    return self.equation(TK, **self.p)

equation staticmethod ¤

equation(
    T: Union[float, FloatArray],
    DSa: Union[float, FloatArray],
    DHa: Union[float, FloatArray],
    kappa: Union[float, FloatArray],
) -> Union[float, FloatArray]

Eyring equation.

PARAMETER DESCRIPTION
T

Temperature. Unit = K.

TYPE: float | FloatArray

DSa

Entropy of activation, \(\Delta S^\ddagger\). Unit = J/(mol·K).

TYPE: float | FloatArray

DHa

Enthalpy of activation, \(\Delta H^\ddagger\). Unit = J/mol.

TYPE: float | FloatArray

kappa

Transmission coefficient.

TYPE: float | FloatArray

RETURNS DESCRIPTION
float | FloatArray

Coefficient value. Unit = 1/s.

Source code in src/polykin/kinetics/eyring.py
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
@staticmethod
def equation(T: Union[float, FloatArray],
             DSa: Union[float, FloatArray],
             DHa: Union[float, FloatArray],
             kappa: Union[float, FloatArray],
             ) -> Union[float, FloatArray]:
    r"""Eyring equation.

    Parameters
    ----------
    T : float | FloatArray
        Temperature.
        Unit = K.
    DSa : float | FloatArray
        Entropy of activation, $\Delta S^\ddagger$.
        Unit = J/(mol·K).
    DHa : float | FloatArray
        Enthalpy of activation, $\Delta H^\ddagger$.
        Unit = J/mol.
    kappa : float | FloatArray
        Transmission coefficient.

    Returns
    -------
    float | FloatArray
        Coefficient value. Unit = 1/s.
    """
    return kappa * kB*T/h * exp((DSa - DHa/T)/R)

fit ¤

fit(
    T: FloatVectorLike,
    Y: FloatVectorLike,
    sigmaY: Optional[FloatVectorLike] = None,
    fit_only: Optional[list[str]] = None,
    logY: bool = False,
    plot: bool = True,
) -> dict

Fit equation to data using non-linear regression.

PARAMETER DESCRIPTION
T

Temperature. Unit = K.

TYPE: FloatVector

Y

Property to be fitted. Unit = Any.

TYPE: FloatVector

sigmaY

Standard deviation of Y. Unit = [Y].

TYPE: FloatVector | None DEFAULT: None

fit_only

List with name of parameters to be fitted.

TYPE: list[str] | None DEFAULT: None

logY

If True, the fit will be done in terms of log(Y).

TYPE: bool DEFAULT: False

plot

If True a plot comparing data and fitted correlation will be generated.

TYPE: bool DEFAULT: True

RETURNS DESCRIPTION
dict

A dictionary of results with the following keys: 'success', 'parameters', 'covariance', and 'plot'.

Source code in src/polykin/properties/equations/base.py
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
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
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
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
def fit(self,
        T: FloatVectorLike,
        Y: FloatVectorLike,
        sigmaY: Optional[FloatVectorLike] = None,
        fit_only: Optional[list[str]] = None,
        logY: bool = False,
        plot: bool = True,
        ) -> dict:
    """Fit equation to data using non-linear regression.

    Parameters
    ----------
    T : FloatVector
        Temperature. Unit = K.
    Y : FloatVector
        Property to be fitted. Unit = Any.
    sigmaY : FloatVector | None
        Standard deviation of Y. Unit = [Y].
    fit_only : list[str] | None
        List with name of parameters to be fitted.
    logY : bool
        If `True`, the fit will be done in terms of log(Y).
    plot : bool
        If `True` a plot comparing data and fitted correlation will be
        generated.

    Returns
    -------
    dict
        A dictionary of results with the following keys: 'success',
        'parameters', 'covariance', and 'plot'.
    """

    # Current parameter values
    pdict = self.p.copy()

    # Select parameters to be fitted
    pnames_fit = [name for name, info in self._pinfo.items() if info[1]]
    if fit_only:
        pnames_fit = set(fit_only) & set(pnames_fit)
    p0 = [pdict[pname] for pname in pnames_fit]

    # Fit function
    def ffit(x, *p):
        for pname, pvalue in zip(pnames_fit, p):
            pdict[pname] = pvalue
        Yfit = self.equation(T=x, **pdict)
        if logY:
            Yfit = log(Yfit)
        return Yfit

    solution = curve_fit(ffit,
                         xdata=T,
                         ydata=log(Y) if logY else Y,
                         p0=p0,
                         sigma=sigmaY,
                         absolute_sigma=False,
                         full_output=True)
    result = {}
    result['success'] = bool(solution[4])
    if solution[4]:
        popt = solution[0]
        pcov = solution[1]
        print("Fit successful.")
        for pname, pvalue in zip(pnames_fit, popt):
            print(f"{pname}: {pvalue}")
        print("Covariance:")
        print(pcov)
        result['covariance'] = pcov

        # Update attributes
        self.Trange = (min(T), max(T))
        for pname, pvalue in zip(pnames_fit, popt):
            self.p[pname] = pvalue
        result['parameters'] = pdict

        # plot
        if plot:
            kind = 'semilogy' if logY else 'linear'
            fig, ax = self.plot(kind=kind, return_objects=True)  # ok
            ax.plot(T, Y, 'o', mfc='none')
            result['plot'] = (fig, ax)
    else:
        print("Fit error: ", solution[3])
        result['message'] = solution[3]

    return result

plot ¤

plot(
    kind: Literal[
        "linear", "semilogy", "Arrhenius"
    ] = "linear",
    Trange: Optional[tuple[float, float]] = None,
    Tunit: Literal["C", "K"] = "K",
    title: Optional[str] = None,
    axes: Optional[Axes] = None,
    return_objects: bool = False,
) -> Optional[tuple[Optional[Figure], Axes]]

Plot quantity as a function of temperature.

PARAMETER DESCRIPTION
kind

Kind of plot to be generated.

TYPE: Literal['linear', 'semilogy', 'Arrhenius'] DEFAULT: 'linear'

Trange

Temperature range for x-axis. If None, the validity range (Tmin, Tmax) will be used. If no validity range was defined, the range will default to 0-100°C.

TYPE: tuple[float, float] | None DEFAULT: None

Tunit

Temperature unit.

TYPE: Literal['C', 'K'] DEFAULT: 'K'

title

Title of plot. If None, the object name will be used.

TYPE: str | None DEFAULT: None

axes

Matplotlib Axes object.

TYPE: Axes | None DEFAULT: None

return_objects

If True, the Figure and Axes objects are returned (for saving or further manipulations).

TYPE: bool DEFAULT: False

RETURNS DESCRIPTION
tuple[Figure | None, Axes] | None

Figure and Axes objects if return_objects is True.

Source code in src/polykin/properties/equations/base.py
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
def plot(self,
         kind: Literal['linear', 'semilogy', 'Arrhenius'] = 'linear',
         Trange: Optional[tuple[float, float]] = None,
         Tunit: Literal['C', 'K'] = 'K',
         title: Optional[str] = None,
         axes: Optional[Axes] = None,
         return_objects: bool = False
         ) -> Optional[tuple[Optional[Figure], Axes]]:
    """Plot quantity as a function of temperature.

    Parameters
    ----------
    kind : Literal['linear', 'semilogy', 'Arrhenius']
        Kind of plot to be generated.
    Trange : tuple[float, float] | None
        Temperature range for x-axis. If `None`, the validity range
        (Tmin, Tmax) will be used. If no validity range was defined, the
        range will default to 0-100°C.
    Tunit : Literal['C', 'K']
        Temperature unit.
    title : str | None
        Title of plot. If `None`, the object name will be used.
    axes : Axes | None
        Matplotlib Axes object.
    return_objects : bool
        If `True`, the Figure and Axes objects are returned (for saving or
        further manipulations).

    Returns
    -------
    tuple[Figure | None, Axes] | None
        Figure and Axes objects if return_objects is `True`.
    """

    # Check inputs
    check_in_set(kind, {'linear', 'semilogy', 'Arrhenius'}, 'kind')
    check_in_set(Tunit, {'K', 'C'}, 'Tunit')
    if Trange is not None:
        Trange_min = 0.
        if Tunit == 'C':
            Trange_min = -273.15
        check_valid_range(Trange, Trange_min, np.inf, 'Trange')

    # Plot objects
    if axes is None:
        fig, ax = plt.subplots()
        if title is None:
            title = self.name
        if title:
            fig.suptitle(title)
        label = None
    else:
        fig = None
        ax = axes
        label = self.name

    # Units and xlabel
    Tunit_range = Tunit
    if kind == 'Arrhenius':
        Tunit = 'K'
    Tsymbol = Tunit
    if Tunit == 'C':
        Tsymbol = '°' + Tunit

    if kind == 'Arrhenius':
        xlabel = r"$1/T$ [" + Tsymbol + r"$^{-1}$]"
    else:
        xlabel = fr"$T$ [{Tsymbol}]"

    # ylabel
    ylabel = fr"${self.symbol}$ [{self.unit}]"
    if axes is not None:
        ylabel0 = ax.get_ylabel()
        if ylabel0 and ylabel not in ylabel0:
            ylabel = ylabel0 + ", " + ylabel

    ax.set_xlabel(xlabel)
    ax.set_ylabel(ylabel)
    ax.grid(True)

    # x-axis vector
    if Trange is not None:
        if Tunit_range == 'C':
            Trange = (Trange[0]+273.15, Trange[1]+273.15)
    else:
        Trange = (np.min(self.Trange[0]), np.max(self.Trange[1]))
        if Trange == (0.0, np.inf):
            Trange = (273.15, 373.15)

    try:
        shape = self._shape
    except AttributeError:
        shape = None
    if shape is not None:
        print("Plot method not yet implemented for array-like equations.")
    else:
        TK = np.linspace(*Trange, 100)
        y = self.__call__(TK, 'K')
        T = TK
        if Tunit == 'C':
            T -= 273.15
        if kind == 'linear':
            ax.plot(T, y, label=label)
        elif kind == 'semilogy':
            ax.semilogy(T, y, label=label)
        elif kind == 'Arrhenius':
            ax.semilogy(1/TK, y, label=label)

    if fig is None:
        ax.legend(bbox_to_anchor=(1.05, 1.0), loc="upper left")

    if return_objects:
        return (fig, ax)

shape property ¤

shape: Optional[tuple[int, ...]]

Shape of underlying parameter array.

PropagationHalfLength ¤

Half-length model for the decay of the propagation rate coefficient with chain length.

This model implements the chain-length dependence:

\[ k_p(i) = k_p \left[ 1+ (C - 1)\exp{\left (-\frac{\ln 2}{i_{1/2}} (i-1) \right )} \right] \]

where \(k_p=k_p(\infty)\) is the long-chain value of the propagation rate coefficient, \(C\ge 1\) is the ratio \(k_p(1)/k_p\) and \((i_{1/2}+1)\) is the hypothetical chain-length at which the difference \(k_p(1) - k_p\) is halved.

References

  • Smith, Gregory B., et al. "The effects of chain length dependent propagation and termination on the kinetics of free-radical polymerization at low chain lengths." European polymer journal 41.2 (2005): 225-230.
PARAMETER DESCRIPTION
kp

Long-chain value of the propagation rate coefficient, \(k_p\).

TYPE: Arrhenius | Eyring

C

Ratio of the propagation coefficients of a monomeric radical and a long-chain radical, \(C\).

TYPE: float DEFAULT: 10.0

ihalf

Half-length, \(i_{1/2}\).

TYPE: float DEFAULT: 1.0

name

Name.

TYPE: str DEFAULT: ''

Examples:

>>> from polykin.kinetics import PropagationHalfLength, Arrhenius
>>> kp = Arrhenius(
...    10**7.63, 32.5e3/8.314, Tmin=261., Tmax=366.,
...    symbol='k_p', unit='L/mol/s', name='kp of styrene')
>>> kpi = PropagationHalfLength(kp, C=10, ihalf=0.5,
...    name='kp(T,i) of styrene')
>>> kpi
name:    kp(T,i) of styrene
C:       10
ihalf:   0.5
kp:
  name:            kp of styrene
  symbol:          k_p
  unit:            L/mol/s
  Trange [K]:      (261.0, 366.0)
  k0 [L/mol/s]:    42657951.88015926
  EaR [K]:         3909.0690401732018
  T0 [K]:          inf
>>> kpi(T=50., i=3, Tunit='C')
371.75986615653215
Source code in src/polykin/kinetics/cldpropagation.py
 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
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
class PropagationHalfLength(KineticCoefficientCLD):
    r"""Half-length model for the decay of the propagation rate coefficient
    with chain length.

    This model implements the chain-length dependence:

    $$ k_p(i) = k_p \left[ 1+ (C - 1)\exp{\left (-\frac{\ln 2}{i_{1/2}} (i-1)
                \right )} \right] $$

    where $k_p=k_p(\infty)$ is the long-chain value of the propagation rate
    coefficient, $C\ge 1$ is the ratio $k_p(1)/k_p$ and $(i_{1/2}+1)$ is the
    hypothetical chain-length at which the difference $k_p(1) - k_p$ is halved.

    **References**

    *   Smith, Gregory B., et al. "The effects of chain length dependent
        propagation and termination on the kinetics of free-radical
        polymerization at low chain lengths." European polymer journal
        41.2 (2005): 225-230.

    Parameters
    ----------
    kp : Arrhenius | Eyring
        Long-chain value of the propagation rate coefficient, $k_p$.
    C : float
        Ratio of the propagation coefficients of a monomeric radical and a
        long-chain radical, $C$.
    ihalf : float
        Half-length, $i_{1/2}$.
    name : str
        Name.

    Examples
    --------
    >>> from polykin.kinetics import PropagationHalfLength, Arrhenius
    >>> kp = Arrhenius(
    ...    10**7.63, 32.5e3/8.314, Tmin=261., Tmax=366.,
    ...    symbol='k_p', unit='L/mol/s', name='kp of styrene')

    >>> kpi = PropagationHalfLength(kp, C=10, ihalf=0.5,
    ...    name='kp(T,i) of styrene')
    >>> kpi
    name:    kp(T,i) of styrene
    C:       10
    ihalf:   0.5
    kp:
      name:            kp of styrene
      symbol:          k_p
      unit:            L/mol/s
      Trange [K]:      (261.0, 366.0)
      k0 [L/mol/s]:    42657951.88015926
      EaR [K]:         3909.0690401732018
      T0 [K]:          inf

    >>> kpi(T=50., i=3, Tunit='C')
    371.75986615653215
    """

    kp: Union[Arrhenius, Eyring]
    C: float
    ihalf: float
    symbol: str = 'k_p(i)'

    def __init__(self,
                 kp: Union[Arrhenius, Eyring],
                 C: float = 10.0,
                 ihalf: float = 1.0,
                 name: str = ''
                 ) -> None:

        check_type(kp, (Arrhenius, Eyring), 'kp')
        check_bounds(C, 1., 100., 'C')
        check_bounds(ihalf, 0.1, 10, 'ihalf')

        self.kp = kp
        self.C = C
        self.ihalf = ihalf
        self.name = name

    def __repr__(self) -> str:
        return custom_repr(self, ('name', 'C', 'ihalf', 'kp'))

    @staticmethod
    def equation(i: Union[int, IntArray],
                 kp: Union[float, FloatArray],
                 C: float,
                 ihalf: float
                 ) -> Union[float, FloatArray]:
        r"""Half-length model chain-length dependence equation.

        Parameters
        ----------
        i : int | IntArray
            Chain length of radical.
        kp : float | FloatArray
            Long-chain value of the propagation rate coefficient, $k_p$.
        C : float
            Ratio of the propagation coefficients of a monomeric radical and a
            long-chain radical, $C$.
        ihalf : float
            Half-length, $i_{i/2}$.

        Returns
        -------
        float | FloatArray
            Coefficient value.
        """

        return kp*(1 + (C - 1)*exp(-log(2)*(i - 1)/ihalf))

    def __call__(self,
                 T: Union[float, FloatArrayLike],
                 i: Union[int, IntArrayLike],
                 Tunit: Literal['C', 'K'] = 'K'
                 ) -> Union[float, FloatArray]:
        r"""Evaluate kinetic coefficient at given conditions, including unit
        conversion and range check.

        Parameters
        ----------
        T : float | FloatArrayLike
            Temperature.
            Unit = `Tunit`.
        i : int | IntArrayLike
            Chain length of radical.
        Tunit : Literal['C', 'K']
            Temperature unit.

        Returns
        -------
        float | FloatArray
            Coefficient value.
        """
        TK = convert_check_temperature(T, Tunit, self.kp.Trange)
        if isinstance(i, (list, tuple)):
            i = np.array(i, dtype=np.int32)
        return self.equation(i, self.kp(TK), self.C, self.ihalf)

__call__ ¤

__call__(
    T: Union[float, FloatArrayLike],
    i: Union[int, IntArrayLike],
    Tunit: Literal["C", "K"] = "K",
) -> Union[float, FloatArray]

Evaluate kinetic coefficient at given conditions, including unit conversion and range check.

PARAMETER DESCRIPTION
T

Temperature. Unit = Tunit.

TYPE: float | FloatArrayLike

i

Chain length of radical.

TYPE: int | IntArrayLike

Tunit

Temperature unit.

TYPE: Literal['C', 'K'] DEFAULT: 'K'

RETURNS DESCRIPTION
float | FloatArray

Coefficient value.

Source code in src/polykin/kinetics/cldpropagation.py
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
def __call__(self,
             T: Union[float, FloatArrayLike],
             i: Union[int, IntArrayLike],
             Tunit: Literal['C', 'K'] = 'K'
             ) -> Union[float, FloatArray]:
    r"""Evaluate kinetic coefficient at given conditions, including unit
    conversion and range check.

    Parameters
    ----------
    T : float | FloatArrayLike
        Temperature.
        Unit = `Tunit`.
    i : int | IntArrayLike
        Chain length of radical.
    Tunit : Literal['C', 'K']
        Temperature unit.

    Returns
    -------
    float | FloatArray
        Coefficient value.
    """
    TK = convert_check_temperature(T, Tunit, self.kp.Trange)
    if isinstance(i, (list, tuple)):
        i = np.array(i, dtype=np.int32)
    return self.equation(i, self.kp(TK), self.C, self.ihalf)

equation staticmethod ¤

equation(
    i: Union[int, IntArray],
    kp: Union[float, FloatArray],
    C: float,
    ihalf: float,
) -> Union[float, FloatArray]

Half-length model chain-length dependence equation.

PARAMETER DESCRIPTION
i

Chain length of radical.

TYPE: int | IntArray

kp

Long-chain value of the propagation rate coefficient, \(k_p\).

TYPE: float | FloatArray

C

Ratio of the propagation coefficients of a monomeric radical and a long-chain radical, \(C\).

TYPE: float

ihalf

Half-length, \(i_{i/2}\).

TYPE: float

RETURNS DESCRIPTION
float | FloatArray

Coefficient value.

Source code in src/polykin/kinetics/cldpropagation.py
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
@staticmethod
def equation(i: Union[int, IntArray],
             kp: Union[float, FloatArray],
             C: float,
             ihalf: float
             ) -> Union[float, FloatArray]:
    r"""Half-length model chain-length dependence equation.

    Parameters
    ----------
    i : int | IntArray
        Chain length of radical.
    kp : float | FloatArray
        Long-chain value of the propagation rate coefficient, $k_p$.
    C : float
        Ratio of the propagation coefficients of a monomeric radical and a
        long-chain radical, $C$.
    ihalf : float
        Half-length, $i_{i/2}$.

    Returns
    -------
    float | FloatArray
        Coefficient value.
    """

    return kp*(1 + (C - 1)*exp(-log(2)*(i - 1)/ihalf))

TerminationCompositeModel ¤

Composite model for the termination rate coefficient between two radicals.

This model implements the chain-length dependence:

\[ k_t(i,j)=\sqrt{k_t(i,i) k_t(j,j)} \]

with:

\[ k_t(i,i)=\begin{cases} k_t(1,1)i^{-\alpha_S}& \text{if } i \leq i_{crit} \\ k_t(1,1)i_{crit}^{-(\alpha_S-\alpha_L)}i^{-\alpha_L} & \text{if }i>i_{crit} \end{cases} \]

where \(k_t(1,1)\) is the temperature-dependent termination rate coefficient between two monomeric radicals, \(i_{crit}\) is the critical chain length, \(\alpha_S\) is the short-chain exponent, and \(\alpha_L\) is the long-chain exponent.

References

  • Smith, Gregory B., Gregory T. Russell, and Johan PA Heuts. "Termination in dilute-solution free-radical polymerization: a composite model." Macromolecular theory and simulations 12.5 (2003): 299-314.
PARAMETER DESCRIPTION
kt11

Temperature-dependent termination rate coefficient between two monomeric radicals, \(k_t(1,1)\).

TYPE: Arrhenius | Eyring

icrit

Critical chain length, \(i_{crit}\).

TYPE: int

aS

Short-chain exponent, \(\alpha_S\).

TYPE: float DEFAULT: 0.5

aL

Long-chain exponent, \(\alpha_L\).

TYPE: float DEFAULT: 0.2

name

Name.

TYPE: str DEFAULT: ''

Examples:

>>> from polykin.kinetics import TerminationCompositeModel, Arrhenius
>>> kt11 = Arrhenius(1e9, 2e3, T0=298.,
...        symbol='k_t(T,1,1)', unit='L/mol/s', name='kt11 of Y')
>>> ktij = TerminationCompositeModel(kt11, icrit=30, name='ktij of Y')
>>> ktij
name:    ktij of Y
icrit:   30
aS:      0.5
aL:      0.2
kt11:
  name:            kt11 of Y
  symbol:          k_t(T,1,1)
  unit:            L/mol/s
  Trange [K]:      (0.0, inf)
  k0 [L/mol/s]:    1000000000.0
  EaR [K]:         2000.0
  T0 [K]:          298.0
>>> ktij(T=25., i=150, j=200, Tunit='C')
129008375.03821689
Source code in src/polykin/kinetics/cldtermination.py
 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
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
class TerminationCompositeModel(KineticCoefficientCLD):
    r"""Composite model for the termination rate coefficient between two
    radicals.

    This model implements the chain-length dependence:

    $$ k_t(i,j)=\sqrt{k_t(i,i) k_t(j,j)} $$

    with:

    $$ k_t(i,i)=\begin{cases}
    k_t(1,1)i^{-\alpha_S}& \text{if } i \leq i_{crit} \\
    k_t(1,1)i_{crit}^{-(\alpha_S-\alpha_L)}i^{-\alpha_L} & \text{if }i>i_{crit}
    \end{cases} $$

    where $k_t(1,1)$ is the temperature-dependent termination rate coefficient
    between two monomeric radicals, $i_{crit}$ is the critical chain length,
    $\alpha_S$ is the short-chain exponent, and $\alpha_L$ is the long-chain
    exponent.

    **References**

    *   Smith, Gregory B., Gregory T. Russell, and Johan PA Heuts. "Termination
        in dilute-solution free-radical polymerization: a composite model."
        Macromolecular theory and simulations 12.5 (2003): 299-314.

    Parameters
    ----------
    kt11 : Arrhenius | Eyring
        Temperature-dependent termination rate coefficient between two
        monomeric radicals, $k_t(1,1)$.
    icrit : int
        Critical chain length, $i_{crit}$.
    aS : float
        Short-chain exponent, $\alpha_S$.
    aL : float
        Long-chain exponent, $\alpha_L$.
    name : str
        Name.

    Examples
    --------
    >>> from polykin.kinetics import TerminationCompositeModel, Arrhenius
    >>> kt11 = Arrhenius(1e9, 2e3, T0=298.,
    ...        symbol='k_t(T,1,1)', unit='L/mol/s', name='kt11 of Y')

    >>> ktij = TerminationCompositeModel(kt11, icrit=30, name='ktij of Y')
    >>> ktij
    name:    ktij of Y
    icrit:   30
    aS:      0.5
    aL:      0.2
    kt11:
      name:            kt11 of Y
      symbol:          k_t(T,1,1)
      unit:            L/mol/s
      Trange [K]:      (0.0, inf)
      k0 [L/mol/s]:    1000000000.0
      EaR [K]:         2000.0
      T0 [K]:          298.0

    >>> ktij(T=25., i=150, j=200, Tunit='C')
    129008375.03821689
    """

    kt11: Union[Arrhenius, Eyring]
    icrit: int
    aS: float
    aL: float
    symbol: str = 'k_t (i,j)'

    def __init__(self,
                 kt11: Union[Arrhenius, Eyring],
                 icrit: int,
                 aS: float = 0.5,
                 aL: float = 0.2,
                 name: str = ''
                 ) -> None:

        check_type(kt11, (Arrhenius, Eyring), 'k11')
        check_bounds(icrit, 1, 200, 'icrit')
        check_bounds(aS, 0, 1, 'alpha_short')
        check_bounds(aL, 0, 0.5, 'alpha_long')

        self.kt11 = kt11
        self.icrit = icrit
        self.aS = aS
        self.aL = aL
        self.name = name

    def __repr__(self) -> str:
        return custom_repr(self, ('name', 'icrit', 'aS', 'aL', 'kt11'))

    @staticmethod
    def equation(i: Union[int, IntArray],
                 j: Union[int, IntArray],
                 kt11: Union[float, FloatArray],
                 icrit: int,
                 aS: float,
                 aL: float,
                 ) -> Union[float, FloatArray]:
        r"""Composite model chain-length dependence equation.

        Parameters
        ----------
        i : int | IntArray
            Chain length of 1st radical.
        j : int | IntArray
            Chain length of 2nd radical.
        kt11 : float | FloatArray
            Temperature-dependent termination rate coefficient between two
            monomeric radicals, $k_t(1,1)$.
        icrit : int
            Critical chain length, $i_{crit}$.
        aS : float
            Short-chain exponent, $\alpha_S$.
        aL : float
            Long-chain exponent, $\alpha_L$.

        Returns
        -------
        float | FloatArray
            Coefficient value.
        """

        def ktii(i):
            return np.where(i <= icrit,
                            kt11*i**(-aS),
                            kt11*icrit**(-aS+aL)*i**(-aL))

        return sqrt(ktii(i)*ktii(j))

    def __call__(self,
                 T: Union[float, FloatArrayLike],
                 i: Union[int, IntArrayLike],
                 j: Union[int, IntArrayLike],
                 Tunit: Literal['C', 'K'] = 'K'
                 ) -> Union[float, FloatArray]:
        r"""Evaluate kinetic coefficient at given conditions, including unit
        conversion and range check.

        Parameters
        ----------
        T : float | FloatArrayLike
            Temperature.
            Unit = `Tunit`.
        i : int | IntArrayLike
            Chain length of 1st radical.
        j : int | IntArrayLike
            Chain length of 2nd radical.
        Tunit : Literal['C', 'K']
            Temperature unit.

        Returns
        -------
        float | FloatArray
            Coefficient value.
        """
        TK = convert_check_temperature(T, Tunit, self.kt11.Trange)
        if isinstance(i, (list, tuple)):
            i = np.array(i, dtype=np.int32)
        if isinstance(j, (list, tuple)):
            j = np.array(j, dtype=np.int32)
        return self.equation(i, j, self.kt11(TK), self.icrit, self.aS, self.aL)

__call__ ¤

__call__(
    T: Union[float, FloatArrayLike],
    i: Union[int, IntArrayLike],
    j: Union[int, IntArrayLike],
    Tunit: Literal["C", "K"] = "K",
) -> Union[float, FloatArray]

Evaluate kinetic coefficient at given conditions, including unit conversion and range check.

PARAMETER DESCRIPTION
T

Temperature. Unit = Tunit.

TYPE: float | FloatArrayLike

i

Chain length of 1st radical.

TYPE: int | IntArrayLike

j

Chain length of 2nd radical.

TYPE: int | IntArrayLike

Tunit

Temperature unit.

TYPE: Literal['C', 'K'] DEFAULT: 'K'

RETURNS DESCRIPTION
float | FloatArray

Coefficient value.

Source code in src/polykin/kinetics/cldtermination.py
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
def __call__(self,
             T: Union[float, FloatArrayLike],
             i: Union[int, IntArrayLike],
             j: Union[int, IntArrayLike],
             Tunit: Literal['C', 'K'] = 'K'
             ) -> Union[float, FloatArray]:
    r"""Evaluate kinetic coefficient at given conditions, including unit
    conversion and range check.

    Parameters
    ----------
    T : float | FloatArrayLike
        Temperature.
        Unit = `Tunit`.
    i : int | IntArrayLike
        Chain length of 1st radical.
    j : int | IntArrayLike
        Chain length of 2nd radical.
    Tunit : Literal['C', 'K']
        Temperature unit.

    Returns
    -------
    float | FloatArray
        Coefficient value.
    """
    TK = convert_check_temperature(T, Tunit, self.kt11.Trange)
    if isinstance(i, (list, tuple)):
        i = np.array(i, dtype=np.int32)
    if isinstance(j, (list, tuple)):
        j = np.array(j, dtype=np.int32)
    return self.equation(i, j, self.kt11(TK), self.icrit, self.aS, self.aL)

equation staticmethod ¤

equation(
    i: Union[int, IntArray],
    j: Union[int, IntArray],
    kt11: Union[float, FloatArray],
    icrit: int,
    aS: float,
    aL: float,
) -> Union[float, FloatArray]

Composite model chain-length dependence equation.

PARAMETER DESCRIPTION
i

Chain length of 1st radical.

TYPE: int | IntArray

j

Chain length of 2nd radical.

TYPE: int | IntArray

kt11

Temperature-dependent termination rate coefficient between two monomeric radicals, \(k_t(1,1)\).

TYPE: float | FloatArray

icrit

Critical chain length, \(i_{crit}\).

TYPE: int

aS

Short-chain exponent, \(\alpha_S\).

TYPE: float

aL

Long-chain exponent, \(\alpha_L\).

TYPE: float

RETURNS DESCRIPTION
float | FloatArray

Coefficient value.

Source code in src/polykin/kinetics/cldtermination.py
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
@staticmethod
def equation(i: Union[int, IntArray],
             j: Union[int, IntArray],
             kt11: Union[float, FloatArray],
             icrit: int,
             aS: float,
             aL: float,
             ) -> Union[float, FloatArray]:
    r"""Composite model chain-length dependence equation.

    Parameters
    ----------
    i : int | IntArray
        Chain length of 1st radical.
    j : int | IntArray
        Chain length of 2nd radical.
    kt11 : float | FloatArray
        Temperature-dependent termination rate coefficient between two
        monomeric radicals, $k_t(1,1)$.
    icrit : int
        Critical chain length, $i_{crit}$.
    aS : float
        Short-chain exponent, $\alpha_S$.
    aL : float
        Long-chain exponent, $\alpha_L$.

    Returns
    -------
    float | FloatArray
        Coefficient value.
    """

    def ktii(i):
        return np.where(i <= icrit,
                        kt11*i**(-aS),
                        kt11*icrit**(-aS+aL)*i**(-aL))

    return sqrt(ktii(i)*ktii(j))

Tutorial