Skip to content

polykin.copolymerization¤

inst_copolymer_multi ¤

inst_copolymer_multi(
    f: Optional[FloatVectorLike],
    r: Optional[FloatSquareMatrix],
    P: Optional[FloatSquareMatrix] = None,
) -> FloatVector

Calculate the instantaneous copolymer composition for a multicomponent system.

In a multicomponent system, the instantaneous copolymer composition \(F_i\) can be determined by solving the following set of linear algebraic equations:

\[ \begin{bmatrix} P_{11}-1 & P_{21} & \cdots & P_{N1} \\ P_{12} & P_{22}-1 & \cdots & P_{N2} \\ \vdots & \vdots & \vdots & \vdots \\ 1 & 1 & \cdots & 1 \end{bmatrix} \begin{bmatrix} F_1 \\ F_2 \\ \vdots \\ F_N \end{bmatrix} = \begin{bmatrix} 0 \\ 0 \\ \vdots \\ 1 \end{bmatrix} \]

where \(P_{ij}\) are the transition probabilitites, which can be computed from the instantaneous monomer composition and the reactivity matrix.

References

  • H. K. Frensdorff, R. Pariser; Copolymerization as a Markov Chain. J. Chem. Phys. 1 November 1963; 39 (9): 2303-2309.
PARAMETER DESCRIPTION
f

Vector of instantaneous monomer compositions, \(f_i\).

TYPE: FloatVectorLike(N) | None

r

Matrix of reactivity ratios, \(r_{ij}=k_{ii}/k_{ij}\).

TYPE: FloatSquareMatrix(N, N) | None

P

Matrix of transition probabilities, \(P_{ij}\). If None, it will be computed internally. When calculating other quantities (e.g., sequence lengths, tuples) that also depend on \(P\), it is more efficient to precompute \(P\) once and use it in all cases.

TYPE: FloatSquareMatrix(N, N) | None DEFAULT: None

RETURNS DESCRIPTION
FloatVector(N)

Vector of instantaneous copolymer compositions, \(F_i\).

See also

Examples:

>>> from polykin.copolymerization import inst_copolymer_multi
>>> import numpy as np

Define the reactivity ratio matrix.

>>> r = np.ones((3, 3))
>>> r[0, 1] = 0.2
>>> r[1, 0] = 2.3
>>> r[0, 2] = 3.0
>>> r[2, 0] = 0.9
>>> r[1, 2] = 0.4
>>> r[2, 1] = 1.5

Evaluate the instantaneous copolymer composition at f1=0.5, f2=0.3, f3=0.2.

>>> f = [0.5, 0.3, 0.2]
>>> F = inst_copolymer_multi(f, r)
>>> F
array([0.32138111, 0.41041608, 0.26820282])
Source code in src/polykin/copolymerization/multicomponent.py
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
def inst_copolymer_multi(f: Optional[FloatVectorLike],
                         r: Optional[FloatSquareMatrix],
                         P: Optional[FloatSquareMatrix] = None
                         ) -> FloatVector:
    r"""Calculate the instantaneous copolymer composition for a multicomponent
    system.

    In a multicomponent system, the instantaneous copolymer composition $F_i$
    can be determined by solving the following set of linear algebraic
    equations:

    $$ \begin{bmatrix}
    P_{11}-1  & P_{21}   & \cdots  & P_{N1} \\
    P_{12}    & P_{22}-1 & \cdots  & P_{N2} \\
    \vdots    & \vdots   & \vdots  & \vdots \\
    1         & 1        & \cdots  & 1
    \end{bmatrix}
    \begin{bmatrix}
    F_1    \\
    F_2    \\
    \vdots \\
    F_N
    \end{bmatrix} =
    \begin{bmatrix}
    0      \\
    0      \\
    \vdots \\
    1
    \end{bmatrix} $$

    where $P_{ij}$ are the transition probabilitites, which can be computed
    from the instantaneous monomer composition and the reactivity matrix.

    **References**

    *   H. K. Frensdorff, R. Pariser; Copolymerization as a Markov Chain.
        J. Chem. Phys. 1 November 1963; 39 (9): 2303-2309.

    Parameters
    ----------
    f : FloatVectorLike (N) | None
        Vector of instantaneous monomer compositions, $f_i$.
    r : FloatSquareMatrix (N, N) | None
        Matrix of reactivity ratios, $r_{ij}=k_{ii}/k_{ij}$.
    P : FloatSquareMatrix (N, N) | None
        Matrix of transition probabilities, $P_{ij}$. If `None`, it will
        be computed internally. When calculating other quantities (e.g.,
        sequence lengths, tuples) that also depend on $P$, it is more efficient
        to precompute $P$ once and use it in all cases.

    Returns
    -------
    FloatVector (N)
        Vector of instantaneous copolymer compositions, $F_i$.

    See also
    --------
    * [`inst_copolymer_binary`](inst_copolymer_binary.md):
      specific method for binary systems.
    * [`inst_copolymer_ternary`](inst_copolymer_ternary.md):
      specific method for terpolymer systems.
    * [`monomer_drift_multi`](monomer_drift_multi.md):
      monomer composition drift.
    * [`transitions_multi`](transitions_multi.md):
      instantaneous transition probabilities.

    Examples
    --------
    >>> from polykin.copolymerization import inst_copolymer_multi
    >>> import numpy as np

    Define the reactivity ratio matrix.
    >>> r = np.ones((3, 3))
    >>> r[0, 1] = 0.2
    >>> r[1, 0] = 2.3
    >>> r[0, 2] = 3.0
    >>> r[2, 0] = 0.9
    >>> r[1, 2] = 0.4
    >>> r[2, 1] = 1.5

    Evaluate the instantaneous copolymer composition at f1=0.5, f2=0.3, f3=0.2.
    >>> f = [0.5, 0.3, 0.2]
    >>> F = inst_copolymer_multi(f, r)
    >>> F
    array([0.32138111, 0.41041608, 0.26820282])
    """

    if P is not None and (f is None and r is None):
        pass
    elif P is None and (f is not None and r is not None):
        P = transitions_multi(f, r)
    else:
        raise ValueError("Invalid combination of inputs.")

    N = P.shape[0]
    A = P.T - np.eye(N)
    A[-1, :] = 1.
    b = np.zeros(N)
    b[-1] = 1.
    F = np.linalg.solve(A, b)

    return F