File: _protein.py

package info (click to toggle)
python-skbio 0.6.2-4
  • links: PTS, VCS
  • area: main
  • in suites: trixie
  • size: 9,312 kB
  • sloc: python: 60,482; ansic: 672; makefile: 224
file content (188 lines) | stat: -rw-r--r-- 5,330 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
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
# ----------------------------------------------------------------------------
# Copyright (c) 2013--, scikit-bio development team.
#
# Distributed under the terms of the Modified BSD License.
#
# The full license is in the file LICENSE.txt, distributed with this software.
# ----------------------------------------------------------------------------

import numpy as np
import pandas as pd

from skbio.sequence import Protein
from skbio.stats.ordination import OrdinationResults
from skbio.embedding._embedding import SequenceEmbedding, SequenceVector, _repr_helper


def _validate_protein(sequence):
    if isinstance(sequence, Protein):
        sequence = str(sequence)
    elif isinstance(sequence, str):
        if " " in sequence:
            sequence = sequence.replace(" ", "")

        # perform a check to make sure the sequence is a valid protein sequence
        _ = Protein(sequence)
    return sequence


class ProteinEmbedding(SequenceEmbedding):
    r"""Embedding of a protein sequence.

    Parameters
    ----------
    embedding : array_like
        The embedding of the protein sequence. Row vectors correspond to
        the latent residues coordinates.
    sequence : str, Protein, or 1D ndarray
        Characters representing the protein sequence itself.
    clip_head : bool, optional
        If ``True``, then the first row of the embedding will be removed.
        Some language models specify start tokens, and this parameter can
        be used to account for this.
    clip_tail : bool, optional
        If ``True``, then the last row of the embedding will be removed.
        Some language models specify end tokens, and this parameter can
        be used to account for this.

    See Also
    --------
    skbio.sequence.Protein

    Examples
    --------
    >>> from skbio.embedding import ProteinEmbedding
    >>> import numpy as np
    >>> embedding = np.random.rand(10, 3)
    >>> sequence = "ACDEFGHIKL"
    >>> ProteinEmbedding(embedding, sequence)
    ProteinEmbedding
    --------------------------
    Stats:
        length: 10
        embedding dimension: 3
        has gaps: False
        has degenerates: False
        has definites: True
        has stops: False
    --------------------------
    0 ACDEFGHIKL

    """

    default_write_format = "embed"

    def __init__(self, embedding, sequence, clip_head=False, clip_tail=False, **kwargs):
        embedding = np.asarray(embedding)
        if clip_head:
            embedding = embedding[1:]
        if clip_tail:
            embedding = embedding[:-1]

        sequence = _validate_protein(sequence)
        super(ProteinEmbedding, self).__init__(
            embedding=embedding, sequence=sequence, **kwargs
        )

    @property
    def residues(self):
        r"""Array containing underlying residue characters."""
        return self._ids.view("|S1")

    def __repr__(self):
        r"""Return a string representation of the ProteinEmbedding object.

        Returns
        -------
        str
            String representation of the ProteinEmbedding object.

        See Also
        --------
        skbio.sequence.Protein

        """
        seq = Protein(str(self))
        rstr = _repr_helper(
            repr(seq),
            "Protein",
            "ProteinEmbedding",
            "embedding",
            regex_match="has gaps",
            shape=self.embedding.shape[1],
        )
        return rstr


example_protein_embedding = ProteinEmbedding(
    np.random.default_rng(0).normal(size=(62, 1024)),
    "IGKEEIQQRLAQFVDHWKELKQLAAARGQRLEESLEYQQFVANVEEEEAWINEKMTLVASED",
)


class ProteinVector(SequenceVector):
    r"""Vector representation for a protein sequence.

    Parameters
    ----------
    vector : 1D or 2D array_like
        The vector representation of the protein sequence. Typically a 1D array. Can
        also be a 2D array with only one row.
    sequence : str, Sequence, or 1D ndarray
        Characters representing the protein sequence itself.

    See Also
    --------
    SequenceVector
    skbio.sequence.Protein

    Examples
    --------
    >>> from skbio.embedding import ProteinVector
    >>> import numpy as np
    >>> vector = np.random.rand(10)
    >>> sequence = "ACDEFGHIKL"
    >>> ProteinVector(vector, sequence)
    ProteinVector
    --------------------------
    Stats:
        length: 10
        vector dimension: 10
        has gaps: False
        has degenerates: False
        has definites: True
        has stops: False
    --------------------------
    0 ACDEFGHIKL

    """

    default_write_format = "embed"

    def __init__(self, vector, sequence: str, **kwargs):
        sequence = _validate_protein(sequence)
        super(ProteinVector, self).__init__(vector, sequence=sequence, **kwargs)

    def __repr__(self):
        r"""Return a string representation of the ProteinVector object.

        Returns
        -------
        str
            String representation of the ProteinEmbedding object.

        See Also
        --------
        skbio.sequence.Protein

        """
        seq = Protein(str(self))
        rstr = _repr_helper(
            repr(seq),
            "Protein",
            "ProteinVector",
            "vector",
            regex_match="has gaps",
            shape=self.embedding.shape[1],
        )
        return rstr