Skip to content

polykin.copolymerization¤

tuples_multi ¤

tuples_multi(
    P: FloatSquareMatrix,
    n: int,
    F: Optional[FloatVectorLike] = None,
) -> dict[tuple[int, ...], float]

Calculate the instantaneous n-tuple fractions.

For a multicomponent system, the probability of finding a specific sequence \(ijk \cdots rs\) of repeating units is:

\[ A_{ijk \cdots rs} = F_i P_{ij} P_{jk} \cdots P_{rs} \]

where \(F_i\) is the instantaneous copolymer composition, and \(P_{ij}\) is the transition probability \(i \rightarrow j\). Since the direction of chain growth does not play a role, symmetric sequences are combined under the sequence with lower index (e.g., \(A_{112} \leftarrow A_{112} + A_{211}\)).

References

  • NA Dotson, R Galván, RL Laurence, and M Tirrel. Polymerization process modeling, Wiley, 1996, p. 179.
PARAMETER DESCRIPTION
P

Matrix of transition probabilities, \(P_{ij}\).

TYPE: FloatSquareMatrix(N, N)

n

Tuple length, e.g. monads (1), diads (2), triads (3), etc.

TYPE: int

F

Vector of instantaneous copolymer composition, \(F_i\). If None, the value will be computed internally. When calculating tuples of various lengths, it is more efficient to precompute \(F\) and use it in all tuple cases.

TYPE: FloatVectorLike(N) | None DEFAULT: None

RETURNS DESCRIPTION
dict[tuple[int, ...], float]

Tuple of molar fractions.

See also

Examples:

>>> from polykin.copolymerization import tuples_multi
>>> from polykin.copolymerization import transitions_multi
>>> import numpy as np

Define 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 transition probabilities.

>>> f = [0.5, 0.3, 0.2]
>>> P = transitions_multi(f, r)

Evaluate triad fractions.

>>> A = tuples_multi(P, 3)
>>> A[(0, 0, 0)]
0.018811329044450834
>>> A[(1, 0, 1)]
0.06365013630778116
Source code in src/polykin/copolymerization/multicomponent.py
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
def tuples_multi(P: FloatSquareMatrix,
                 n: int,
                 F: Optional[FloatVectorLike] = None
                 ) -> dict[tuple[int, ...], float]:
    r"""Calculate the instantaneous n-tuple fractions.

    For a multicomponent system, the probability of finding a specific sequence
    $ijk \cdots rs$ of repeating units is:

    $$ A_{ijk \cdots rs} = F_i P_{ij} P_{jk} \cdots P_{rs} $$

    where $F_i$ is the instantaneous copolymer composition, and $P_{ij}$ is
    the transition probability $i \rightarrow j$. Since the direction of chain
    growth does not play a role, symmetric sequences are combined under the
    sequence with lower index (e.g., $A_{112} \leftarrow A_{112} + A_{211}$).

    **References**

    *   NA Dotson, R Galván, RL Laurence, and M Tirrel. Polymerization
        process modeling, Wiley, 1996, p. 179.

    Parameters
    ----------
    P : FloatSquareMatrix (N, N)
        Matrix of transition probabilities, $P_{ij}$.
    n : int
        Tuple length, e.g. monads (1), diads (2), triads (3), etc.
    F : FloatVectorLike (N) | None
        Vector of instantaneous copolymer composition, $F_i$. If `None`,
        the value will be computed internally. When calculating tuples of
        various lengths, it is more efficient to precompute $F$ and use it in
        all tuple cases.

    Returns
    -------
    dict[tuple[int, ...], float]
        Tuple of molar fractions.

    See also
    --------
    * [`sequence_multi`](sequence_multi.md):
      instantaneous sequence lengths.
    * [`transitions_multi`](transitions_multi.md):
      instantaneous transition probabilities.

    Examples
    --------
    >>> from polykin.copolymerization import tuples_multi
    >>> from polykin.copolymerization import transitions_multi
    >>> import numpy as np

    Define 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 transition probabilities.
    >>> f = [0.5, 0.3, 0.2]
    >>> P = transitions_multi(f, r)

    Evaluate triad fractions.
    >>> A = tuples_multi(P, 3)
    >>> A[(0, 0, 0)]
    0.018811329044450834
    >>> A[(1, 0, 1)]
    0.06365013630778116

    """

    # Compute F if not given
    if F is None:
        F = inst_copolymer_multi(f=None, r=None, P=P)

    # Generate all tuple combinations
    N = P.shape[0]
    indexes = list(itertools.product(range(N), repeat=n))

    result = {}
    for idx in indexes:
        # Compute tuple probability
        P_product = 1.
        for j in range(n-1):
            P_product *= P[idx[j], idx[j+1]]
        A = F[idx[0]]*P_product
        # Add probability to dict, but combine symmetric tuples: 12 <- 12 + 21
        reversed_idx = idx[::-1]
        if reversed_idx in result.keys():
            result[reversed_idx] += A
        else:
            result[idx] = A

    return result