Skip to content

descriptors

aliphatic_index(peptide)

Calculate the aliphatic index of a peptide.

The aliphatic index is a measure of the thermal stability of a peptide. It is defined as the volume of a protein that is occupied by aliphatic side chains of amino acids, such as alanine, valine, isoleucine, and leucine (Ikai, 1980).

Aliphatic index is calculated as:

\(\textit{aliphatic_index} = 100 * (f_A + 2.9 f_V + 3.9 f_I + 3.9 f_L)\),

where \(f_A\), \(f_V\), \(f_I\), and \(f_L\) are the frequencies of alanine, valine, isoleucine, and leucine, respectively.

Parameters:

Name Type Description Default
peptide str

Amino acid sequence of the peptide.

required

Returns:

Type Description
float

The aliphatic index of the peptide.

Raises:

Type Description
ValueError

If the peptide sequence contains unknown amino acids or a syntax error.

Examples:

>>> aliphatic_index("AVIL")
292.5
>>> # doubling the length of the peptide (with no aliphatic amino acids) halves the aliphatic index
>>> aliphatic_index('AVILMNPS_p')
146.25
>>> aliphatic_index("AKLVT")
156.0
>>> aliphatic_index('ACD')
33.3...
>>> # Equals to 0, if the peptide contains no aliphatic amino acids
>>> aliphatic_index("WYGHP")
0.0
>>> aliphatic_index('DEFR_mSC')
0.0
Source code in peptidy/descriptors.py
 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
def aliphatic_index(peptide: str) -> float:
    """
    Calculate the aliphatic index of a peptide.

    The aliphatic index is a measure of the thermal stability of a peptide. It is defined as the volume
    of a protein that is occupied by aliphatic side chains of amino acids, such as alanine, valine,
    isoleucine, and leucine ([Ikai, 1980](https://academic.oup.com/jb/article-abstract/88/6/1895/773432?redirectedFrom=fulltext)).

    Aliphatic index is calculated as:

    $\\textit{aliphatic_index} = 100 * (f_A + 2.9 f_V + 3.9 f_I + 3.9 f_L)$,

    where $f_A$, $f_V$, $f_I$, and $f_L$ are the frequencies of alanine, valine, isoleucine, and leucine, respectively.

    Parameters
    ----------
    peptide : str
        Amino acid sequence of the peptide.

    Returns
    -------
    float
        The aliphatic index of the peptide.

    Raises
    ------
    ValueError
        If the peptide sequence contains unknown amino acids or a syntax error.

    Examples
    --------
    >>> aliphatic_index("AVIL")
    292.5
    >>> # doubling the length of the peptide (with no aliphatic amino acids) halves the aliphatic index
    >>> aliphatic_index('AVILMNPS_p')
    146.25
    >>> aliphatic_index("AKLVT")
    156.0
    >>> aliphatic_index('ACD')
    33.3...
    >>> # Equals to 0, if the peptide contains no aliphatic amino acids
    >>> aliphatic_index("WYGHP")
    0.0
    >>> aliphatic_index('DEFR_mSC')
    0.0
    """
    tokenized_peptide = tokenizer.tokenize_peptide(peptide)
    aa_counts = dict(Counter(tokenized_peptide))
    aa_counts = {aa: count / length(peptide) * 100 for aa, count in aa_counts.items()}
    return (
        aa_counts.get("A", 0)
        + 2.9 * aa_counts.get("V", 0)
        + 3.9 * (aa_counts.get("I", 0) + aa_counts.get("L", 0))
    )

aminoacid_frequencies(peptide)

Calculate the frequency (count / peptide length) of all amino acids in the input sequence.

Parameters:

Name Type Description Default
peptide str

Amino acid sequence of the peptide.

required

Returns:

Type Description
Dict[str, float]

A dictionary that contains the frequencies of all amino acids and post-translations, not only the ones present in the sequence. The keys have the format freq_<AA> (AA = letter code of the amino acid).

Raises:

Type Description
ValueError

If the peptide sequence contains unknown amino acids or a syntax error.

Examples:

>>> freqs = aminoacid_frequencies("AVIL")
>>> freqs["freq_A"]
0.25
>>> freqs["freq_V"]
0.25
>>> freqs["freq_C_m"]
0.0
>>> freqs["freq_R"]
0.0
>>> freqs = aminoacid_frequencies('AC_mD')
>>> freqs["freq_C_m"]
0.3...
>>> freqs["freq_D"]
0.3...
>>> freqs["freq_C"]
0.0
>>> freqs["freq_W"]
0.0
>>> aminoacid_frequencies('AXR')
Traceback (most recent call last):
    ...
ValueError: Unknown amino acid(s) in peptide: {'X'}
Source code in peptidy/descriptors.py
 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
def aminoacid_frequencies(peptide: str) -> Dict[str, float]:
    """
    Calculate the frequency (count / peptide length) of all amino acids in the input sequence.

    Parameters
    ----------
    peptide : str
        Amino acid sequence of the peptide.

    Returns
    -------
    Dict[str, float]
        A dictionary that contains the frequencies of all amino acids and post-translations,
        not only the ones present in the sequence. The keys have the format `freq_<AA>` (`AA` = letter code of the
        amino acid).

    Raises
    ------
    ValueError
        If the peptide sequence contains unknown amino acids or a syntax error.

    Examples
    --------
    >>> freqs = aminoacid_frequencies("AVIL")
    >>> freqs["freq_A"]
    0.25
    >>> freqs["freq_V"]
    0.25
    >>> freqs["freq_C_m"]
    0.0
    >>> freqs["freq_R"]
    0.0
    >>> freqs = aminoacid_frequencies('AC_mD')
    >>> freqs["freq_C_m"]
    0.3...
    >>> freqs["freq_D"]
    0.3...
    >>> freqs["freq_C"]
    0.0
    >>> freqs["freq_W"]
    0.0
    >>> aminoacid_frequencies('AXR')
    Traceback (most recent call last):
        ...
    ValueError: Unknown amino acid(s) in peptide: {'X'}
    """
    tokenized_peptide = tokenizer.tokenize_peptide(peptide)
    aa_counts = dict(Counter(tokenized_peptide))
    peptide_len = length(tokenized_peptide)
    return {
        f"freq_{aa}": aa_counts.get(aa, 0) / peptide_len for aa in biology.aminoacids
    }

aromaticity(peptide)

Calculate the sum of the frequencies of aromatic amino-acids ("F", "W", "Y", and "Y_p") as a measure of aromaticity of a peptide (Lobry,1994).

Parameters:

Name Type Description Default
peptide str

Amino acid sequence of the peptide.

required

Returns:

Type Description
float

Frequency of aromatic residues in the peptide.

Raises:

Type Description
ValueError

If the peptide sequence contains unknown amino acids or a syntax error.

Examples:

>>> aromaticity("AVIL")
0.0
>>> aromaticity("WYGHP")
0.4
>>> aromaticity('ACDEF')
0.2
>>> aromaticity('DAY_pCDFWY')
0.5
Source code in peptidy/descriptors.py
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
def aromaticity(peptide: str) -> float:
    """
    Calculate the sum of the frequencies of aromatic amino-acids
    ("F", "W", "Y", and "Y_p") as a measure of aromaticity of a peptide
    ([Lobry,1994](https://academic.oup.com/nar/article-abstract/22/15/3174/1087817)).

    Parameters
    ----------
    peptide : str
        Amino acid sequence of the peptide.

    Returns
    -------
    float
        Frequency of aromatic residues in the peptide.

    Raises
    ------
    ValueError
        If the peptide sequence contains unknown amino acids or a syntax error.

    Examples
    --------
    >>> aromaticity("AVIL")
    0.0
    >>> aromaticity("WYGHP")
    0.4
    >>> aromaticity('ACDEF')
    0.2
    >>> aromaticity('DAY_pCDFWY')
    0.5
    """
    peptide = tokenizer.tokenize_peptide(peptide)
    counts = dict(Counter(peptide))
    return sum([counts.get(aa, 0) for aa in biology.aromatic_aas]) / length(peptide)

average_n_rotatable_bonds(peptide)

Calculate the number of total rotatable bonds divided by the number of amino acids in the peptide.

Chain flexibility is known to play a role in binding [Francesca Peccati & Gonzalo Jiménez-Osés, 2021]{https://pubs.acs.org/doi/10.1021/acsomega.1c00485}. The number of rotatable bonds per amino acid was retrieved from PubChem Kim et al., 2023.

Parameters:

Name Type Description Default
peptide str

Amino acid sequence of the peptide.

required

Returns:

Type Description
float

Average number of rotatable bonds in the peptide.

Raises:

Type Description
ValueError

If the peptide sequence contains unknown amino acids or a syntax error.

Examples:

>>> average_n_rotatable_bonds("AVIL")
2.25
>>> average_n_rotatable_bonds("WYGHP")
2.2
>>> average_n_rotatable_bonds('ACD')
2.0
Source code in peptidy/descriptors.py
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
def average_n_rotatable_bonds(
    peptide: str,
) -> float:
    """
    Calculate the number of total rotatable bonds divided by the number of amino acids in the peptide.

    Chain flexibility is known to play a role in binding [Francesca Peccati & Gonzalo Jiménez-Osés, 2021]{https://pubs.acs.org/doi/10.1021/acsomega.1c00485}.
    The number of rotatable bonds per amino acid was retrieved from PubChem [Kim et al., 2023](https://academic.oup.com/nar/article/51/D1/D1373/6777787?login=true).

    Parameters
    ----------
    peptide : str
        Amino acid sequence of the peptide.

    Returns
    -------
    float
        Average number of rotatable bonds in the peptide.

    Raises
    ------
    ValueError
        If the peptide sequence contains unknown amino acids or a syntax error.

    Examples
    --------
    >>> average_n_rotatable_bonds("AVIL")
    2.25
    >>> average_n_rotatable_bonds("WYGHP")
    2.2
    >>> average_n_rotatable_bonds('ACD')
    2.0
    """
    tokenized_peptide = tokenizer.tokenize_peptide(peptide)
    return sum([biology.n_rotatable_bonds[aa] for aa in tokenized_peptide]) / length(
        peptide
    )

charge(peptide, pH=7)

Calculate the total charge of the sequence.

The method used is first described by Bjellqvist [Bjellqvist et al., 1993][https://analyticalsciencejournals.onlinelibrary.wiley.com/doi/10.1002/elps.11501401163]. The total charge is then calculated based on the Henderson-Hasselbach equation Aronson, 1983. Pka of phosphoserine and phosphothreonine were retrieved from Xie,Jiang & Ben-Amotz, 2005. Pka of phosphotyrosine was taken from Wojciechowski M et al., 2003 Pka of the posttranslational of arginine was kept equal to the pka of arginine, based on Evich M et al.,2015

Parameters:

Name Type Description Default
peptide str

Amino acid sequence of the peptide.

required
pH float, optional

pH at which to calculate charge, by default 7.

7

Returns:

Type Description
float

The total charge of the sequence.

Examples:

>>> charge('ACD', pH=13)
-2.999...
>>> charge('NNNNRKTNGDDSLF')
-0.238...
Source code in peptidy/descriptors.py
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
def charge(peptide: str, pH: float = 7) -> float:
    """
    Calculate the total charge of the sequence.

    The method used is first described by Bjellqvist [Bjellqvist et al., 1993][https://analyticalsciencejournals.onlinelibrary.wiley.com/doi/10.1002/elps.11501401163].
    The total charge is then calculated based on the Henderson-Hasselbach equation [Aronson, 1983](https://www.sciencedirect.com/science/article/pii/0307441283900468?via%3Dihub).
    Pka of phosphoserine and phosphothreonine were retrieved from [Xie,Jiang & Ben-Amotz, 2005](https://www.sciencedirect.com/science/article/pii/S0003269705004124?via%3Dihub).
    Pka of phosphotyrosine was taken from [Wojciechowski M et al., 2003](https://www.ncbi.nlm.nih.gov/pmc/articles/PMC1302655/#:~:text=The%20pKa%20value%20for,et%20al.%2C%201994)
    Pka of the posttranslational of arginine was kept equal to the pka of arginine, based on [Evich M et al.,2015](https://onlinelibrary.wiley.com/doi/full/10.1002/pro.2838)
    Parameters
    ----------
    peptide : str
        Amino acid sequence of the peptide.
    pH : float, optional
        pH at which to calculate charge, by default 7.

    Returns
    -------
    float
        The total charge of the sequence.

    Examples
    --------
    >>> charge('ACD', pH=13)
    -2.999...
    >>> charge('NNNNRKTNGDDSLF')
    -0.238...
    """
    tokenized_peptide = tokenizer.tokenize_peptide(peptide)
    aa_counts = dict(Counter(tokenized_peptide))
    aa_counts["Nterm"] = 1.0
    aa_counts["Cterm"] = 1.0

    pos_crs = {aa: 10 ** (pk - pH) for aa, pk in biology.pos_pks.items()}
    pos_partial_charges = {aa: cr / (cr + 1) for aa, cr in pos_crs.items()}
    pos_charge = sum(
        [aa_counts.get(aa, 0) * pc for aa, pc in pos_partial_charges.items()]
    )

    neg_crs = {aa: 10 ** (pH - pk) for aa, pk in biology.neg_pks.items()}
    neg_partial_charges = {aa: cr / (cr + 1) for aa, cr in neg_crs.items()}
    neg_charge = sum(
        [aa_counts.get(aa, 0) * pc for aa, pc in neg_partial_charges.items()]
    )

    return pos_charge - neg_charge

charge_density(peptide, pH=7)

Calculate the charge of the peptide normalized by weight, i.e., charge / molecular weight.

Parameters:

Name Type Description Default
peptide str

Amino acid sequence of the peptide.

required
pH float, optional

pH at which to calculate charge, by default 7.

7

Returns:

Type Description
float

Charge density.

Examples:

>>> charge_density('KTTENGD')
-0.00161...
>>> charge_density('FPAL', pH=13)
-0.00223...
Source code in peptidy/descriptors.py
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
def charge_density(peptide: str, pH: float = 7) -> float:
    """
    Calculate the charge of the peptide normalized by weight, *i.e.,* charge / molecular weight.

    Parameters
    ----------
    peptide : str
        Amino acid sequence of the peptide.
    pH : float, optional
        pH at which to calculate charge, by default 7.

    Returns
    -------
    float
        Charge density.

    Examples
    --------
    >>> charge_density('KTTENGD')
    -0.00161...
    >>> charge_density('FPAL', pH=13)
    -0.00223...
    """
    return charge(peptide, pH) / molecular_weight(peptide)

compute_descriptors(peptide, descriptor_names=None, pH=7)

Computes multiple descriptors of the peptide.

Parameters:

Name Type Description Default
peptide str

Amino acid sequence of the peptide.

required
descriptor_names List[str], optional

A List of descriptor names. If any of the descriptor names is invalid, ValueError is raised. Set to None by default, and computes all descriptors in this case.

None
pH float, optional

pH value to compute charge and charge density (if requested), by default 7

7

Returns:

Type Description
Dict[str, Union[float, int]]

A dictionary that maps descriptor names to their values.

Raises:

Type Description
ValueError

ValueError is raised if any of the descriptor names is invalid.

Examples:

>>> compute_descriptors('ACD', ['charge', 'charge_density'], 13)
{'charge': -2.99..., 'charge_density': -0.00976...}
>>> compute_descriptors('STY', ['molecular_formula'])
{'n_C': 16, 'n_H': 23, 'n_N': 3, 'n_O': 7, 'n_S': 0, 'n_P': 0}
Source code in peptidy/descriptors.py
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
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
def compute_descriptors(
    peptide: str,
    descriptor_names: List[str] = None,
    pH: float = 7,
) -> Dict[str, Union[float, int]]:
    """Computes multiple descriptors of the peptide.

    Parameters
    ----------
    peptide : str
        Amino acid sequence of the peptide.
    descriptor_names : List[str], optional
        A List of descriptor names. If any of the descriptor names is invalid, ValueError is raised.
        Set to None by default, and computes all descriptors in this case.
    pH : float, optional
        pH value to compute charge and charge density (if requested), by default 7

    Returns
    ----------
    Dict[str, Union[float, int]]
        A dictionary that maps descriptor names to their values.

    Raises
    ----------
    ValueError
        ValueError is raised if any of the descriptor names is invalid.

    Examples
    ----------
    >>> compute_descriptors('ACD', ['charge', 'charge_density'], 13)
    {'charge': -2.99..., 'charge_density': -0.00976...}
    >>> compute_descriptors('STY', ['molecular_formula'])
    {'n_C': 16, 'n_H': 23, 'n_N': 3, 'n_O': 7, 'n_S': 0, 'n_P': 0}
    """
    if descriptor_names is None:
        descriptor_names = list(__DESCRIPTOR_FNS.keys())

    diff = set(descriptor_names) - set(__DESCRIPTOR_FNS.keys())
    if len(diff) > 0:
        raise ValueError(
            f"Invalid descriptor names: {diff}. Possible names are: {list(__DESCRIPTOR_FNS.keys())}"
        )

    name_to_descriptor = dict()
    for name in descriptor_names:
        if name in ["charge", "charge_density"]:
            name_to_descriptor[name] = __DESCRIPTOR_FNS[name](peptide, pH=pH)
        elif name in ["aminoacid_frequencies", "molecular_formula"]:
            name_to_descriptor = {
                **name_to_descriptor,
                **__DESCRIPTOR_FNS[name](peptide),
            }
        else:
            name_to_descriptor[name] = __DESCRIPTOR_FNS[name](peptide)
    return name_to_descriptor

hydrophobic_aa_ratio(peptide)

Calculate the total ratio of hydrophobic amino-acids (A, C, C_m, F, I, L, M, and V) in a peptide. (Nelson & Cox, 2004)

Parameters:

Name Type Description Default
peptide str

Amino acid sequence of the peptide.

required

Returns:

Type Description
float

Total ratio of hydrophobic amino-acids of the peptide.

Examples:

>>> hydrophobic_aa_ratio('FC_mPR_mS_pA')
0.5
>>> hydrophobic_aa_ratio('FPR_mXS_pA')
Traceback (most recent call last):
    ...
ValueError: Unknown amino acid(s) in peptide: {'X'}
Source code in peptidy/descriptors.py
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
def hydrophobic_aa_ratio(peptide: str) -> float:
    """
    Calculate the total ratio of hydrophobic amino-acids (A, C, C_m, F, I, L, M, and V) in a peptide.
    ([Nelson & Cox, 2004](https://mis.kp.ac.rw/admin/admin_panel/kp_lms/files/digital/Core%20Books/Core%20Books%20In%20Nursing%20%20And%20%20Midwifery/H106_%20Biochemistry_Lehninger%20Principles%20of%20Biochemistry,%20Fourth%20Edition%20-%20David%20L.%20Nelson,%20Michael%20M.%20Cox.pdf))

    Parameters
    ----------
    peptide : str
        Amino acid sequence of the peptide.

    Returns
    -------
    float
        Total ratio of hydrophobic amino-acids of the peptide.

    Examples
    --------
    >>> hydrophobic_aa_ratio('FC_mPR_mS_pA')
    0.5
    >>> hydrophobic_aa_ratio('FPR_mXS_pA')
    Traceback (most recent call last):
        ...
    ValueError: Unknown amino acid(s) in peptide: {'X'}
    """
    tokenized_peptide = tokenizer.tokenize_peptide(peptide)
    aa_counts = dict(Counter(tokenized_peptide))
    return sum([aa_counts.get(aa, 0) for aa in biology.hydrophobic_aas]) / length(
        peptide
    )

instability_index(peptide)

Calculate the instability index of the peptide.

The instability index is based on amino acid compositions and computed by summing the instability coefficient of all dipeptide combinations in the peptide. It is based on the frequency of the dipeptide occurring in stable versus unstable proteins Guruprasad, Reddy & Pandit, 1990. A value of 1 is used for amino acid pairs whose instability coefficient is unavailable.

Parameters:

Name Type Description Default
peptide str

Amino acid sequence of the peptide.

required

Returns:

Type Description
float

Examples:

>>> # from https://rdrr.io/cran/Peptides/src/R/boman.R
>>> instability_index('ACFEGM')
81.566...
>>> instability_index('FPP_mS_pA')
Traceback (most recent call last):
    ...
ValueError: Unknown amino acid(s) in peptide: {'P_m'}
Source code in peptidy/descriptors.py
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
def instability_index(peptide: str) -> float:
    """
    Calculate the instability index of the peptide.

    The instability index is based on amino acid compositions and computed by summing the instability coefficient
    of all dipeptide combinations in the peptide. It is based on the frequency of the dipeptide
    occurring in stable versus unstable proteins [Guruprasad, Reddy & Pandit, 1990](https://academic.oup.com/peds/article-abstract/4/2/155/1491271).
    A value of 1 is used for amino acid pairs whose instability coefficient is unavailable.

    Parameters
    ----------
    peptide : str
        Amino acid sequence of the peptide.

    Returns
    -------
    float

    Examples
    --------
    >>> # from https://rdrr.io/cran/Peptides/src/R/boman.R
    >>> instability_index('ACFEGM')
    81.566...
    >>> instability_index('FPP_mS_pA')
    Traceback (most recent call last):
        ...
    ValueError: Unknown amino acid(s) in peptide: {'P_m'}
    """
    tokenized_peptide = tokenizer.tokenize_peptide(peptide)
    aa_pairs = zip(tokenized_peptide, tokenized_peptide[1:])
    instabilities = [biology.instabilities[aa1][aa2] for aa1, aa2 in aa_pairs]
    return sum(instabilities) * 10 / length(peptide)

isoelectric_point(peptide)

Calculate the isoelectric point (pH that the peptide carries no net charge) of the peptide.

The isoelectric point is calculated using the peptide charge at different pH values. The method used is based on the Henderson-Hasselbach equation Aronson, 1983.

Parameters:

Name Type Description Default
peptide str

Amino acid sequence of the peptide.

required

Returns:

Type Description
float

Isoelectric point of the peptide.

Examples:

>>> isoelectric_point('ADEFGHI')
4.35...
>>> isoelectric_point('AMSTV')
5.5234375
Source code in peptidy/descriptors.py
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
def isoelectric_point(peptide: str) -> float:
    """
    Calculate the isoelectric point (pH that the peptide carries no net charge) of the peptide.

    The isoelectric point is calculated using the peptide charge at different pH values.
    The method used is based on the Henderson-Hasselbach equation [Aronson, 1983](https://www.sciencedirect.com/science/article/pii/0307441283900468?via%3Dihub).

    Parameters
    ----------
    peptide : str
        Amino acid sequence of the peptide.

    Returns
    -------
    float
        Isoelectric point of the peptide.

    Examples
    --------
    >>> isoelectric_point('ADEFGHI')
    4.35...
    >>> isoelectric_point('AMSTV')
    5.5234375
    """
    test_ph = 7
    peptide_charge = charge(peptide, test_ph)
    if peptide_charge < 0:
        lower_limit, upper_limit = 0, 7
    else:
        lower_limit, upper_limit = 7, 14

    precision = 10**-4
    while (upper_limit - lower_limit) > precision and abs(peptide_charge) > precision:
        test_ph = (upper_limit + lower_limit) / 2
        peptide_charge = charge(peptide, test_ph)
        if peptide_charge < 0:
            upper_limit = test_ph
        else:
            lower_limit = test_ph

    return test_ph

length(peptide)

Calculate the number of amino acids in the peptide.

Parameters:

Name Type Description Default
peptide str

Amino acid sequence of the peptide.

required

Returns:

Type Description
int

Number of amino acids in the peptide.

Examples:

>>> length('ACD')
3
>>> length('C_mS_pADDWY')
7
>>> length('FP_mS_pA')
Traceback (most recent call last):
    ...
ValueError: Unknown amino acid(s) in peptide: {'P_m'}
Source code in peptidy/descriptors.py
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
def length(peptide: str) -> int:
    """
    Calculate the number of amino acids in the peptide.

    Parameters
    ----------
    peptide : str
        Amino acid sequence of the peptide.

    Returns
    -------
    int
        Number of amino acids in the peptide.

    Examples
    --------
    >>> length('ACD')
    3
    >>> length('C_mS_pADDWY')
    7
    >>> length('FP_mS_pA')
    Traceback (most recent call last):
        ...
    ValueError: Unknown amino acid(s) in peptide: {'P_m'}
    """
    return len(tokenizer.tokenize_peptide(peptide))

molecular_formula(peptide)

Determine the closed molecular formula of the amino acid sequence of the peptide.

The peptide bonds between amino acids are included in the formula. Molecular formulas were retrieved from PubChem.

Parameters:

Name Type Description Default
peptide str

Amino acid sequence of the peptide.

required

Returns:

Type Description
Dict[str, int]

Count of each element in the peptide where the element symbols are keys and counts are values. The keys have the format n_<element> (element = element symbol).

Examples:

>>> formula = molecular_formula('ADEF')
>>> formula['n_H']
28
>>> formula['n_O']
9
>>> formula['n_S']
0
>>> formula = molecular_formula('PTHRAAPDES')
>>> formula['n_H']
69
>>> formula['n_O']
17
>>> formula['n_S']
0
Source code in peptidy/descriptors.py
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
def molecular_formula(
    peptide: str,
) -> Dict[str, int]:
    """
    Determine the closed molecular formula of the amino acid sequence of the peptide.

    The peptide bonds between amino acids are included in the formula. Molecular formulas
    were retrieved from [PubChem](https://pubchem.ncbi.nlm.nih.gov/).

    Parameters
    ----------
    peptide : str
        Amino acid sequence of the peptide.

    Returns
    -------
    Dict[str, int]
        Count of each element in the peptide where the element symbols are keys and counts are values.
        The keys have the format `n_<element>` (`element` = element symbol).

    Examples
    --------
    >>> formula = molecular_formula('ADEF')
    >>> formula['n_H']
    28
    >>> formula['n_O']
    9
    >>> formula['n_S']
    0
    >>> formula = molecular_formula('PTHRAAPDES')
    >>> formula['n_H']
    69
    >>> formula['n_O']
    17
    >>> formula['n_S']
    0
    """
    tokenized_peptide = tokenizer.tokenize_peptide(peptide)
    aa_formulas = [Counter(biology.formulas[aa]) for aa in tokenized_peptide]
    elements = ["C", "H", "N", "O", "S", "P"]
    peptide_formula = dict()
    for element in elements:
        peptide_formula[f"n_{element}"] = sum(
            [formula[element] for formula in aa_formulas]
        )

    peptide_formula["n_H"] = peptide_formula["n_H"] - 2 * length(peptide) + 2
    peptide_formula["n_O"] = peptide_formula["n_O"] - length(peptide) + 1

    return peptide_formula

molecular_weight(peptide)

Calculate the weight (g/mol) of the peptide without peptide bonds.

The molecular weight of the peptide is calculated by summing the weights of the amino acids and subtracting the weight of water for each peptide bond. Molecular weight were retrieved from PubChem.

Parameters:

Name Type Description Default
peptide str

Amino acid sequence of the peptide.

required

Returns:

Type Description
float

Molecular weight of the peptide.

Examples:

>>> molecular_weight('RMK_aS_pCD')
860.885
>>> molecular_weight('DEGHI')
569.56
Source code in peptidy/descriptors.py
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
def molecular_weight(peptide: str) -> float:
    """
    Calculate the weight (g/mol) of the peptide without peptide bonds.

    The molecular weight of the peptide is calculated by summing the weights of the amino acids and
    subtracting the weight of water for each peptide bond. Molecular weight were retrieved from [PubChem](https://pubchem.ncbi.nlm.nih.gov/).


    Parameters
    ----------
    peptide : str
        Amino acid sequence of the peptide.

    Returns
    -------
    float
        Molecular weight of the peptide.

    Examples
    --------
    >>> molecular_weight('RMK_aS_pCD')
    860.885
    >>> molecular_weight('DEGHI')
    569.56
    """

    tokenized_peptide = tokenizer.tokenize_peptide(peptide)
    total_weight = sum([biology.weights[aa] for aa in tokenized_peptide])
    return total_weight - 18.015 * (length(peptide) - 1)

n_h_acceptors(peptide)

Calculate the total number of hydrogen bond acceptors in the peptide.

The number of hydrogen bond acceptors in the peptide is calculated by summing the number of hydrogen bond acceptors. Hydrogen bonds are important in protein-protein interactions Hubbard & Haider, 2010. The number of hydrogen bond acceptors for each amino acid were retrieved from PubChem.

Parameters:

Name Type Description Default
peptide str

Amino acid sequence of the peptide.

required

Returns:

Type Description
float

Total number of hydrogen bond acceptors in the peptide.

Examples:

>>> n_h_acceptors('FSCA')
14
>>> n_h_acceptors('FXS_pGNM')
Traceback (most recent call last):
    ...
ValueError: Unknown amino acid(s) in peptide: {'X'}
Source code in peptidy/descriptors.py
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
def n_h_acceptors(peptide: str) -> float:
    """
    Calculate the total number of hydrogen bond acceptors in the peptide.

    The number of hydrogen bond acceptors in the peptide is calculated by summing the number of hydrogen bond acceptors.
    Hydrogen bonds are important in protein-protein interactions [Hubbard & Haider, 2010](https://doi.org/10.1002/9780470015902.a0003011.pub2).
    The number of hydrogen bond acceptors for each amino acid were retrieved from [PubChem](https://pubchem.ncbi.nlm.nih.gov/).

    Parameters
    ----------
    peptide : str
        Amino acid sequence of the peptide.

    Returns
    -------
    float
        Total number of hydrogen bond acceptors in the peptide.

    Examples
    --------
    >>> n_h_acceptors('FSCA')
    14
    >>> n_h_acceptors('FXS_pGNM')
    Traceback (most recent call last):
        ...
    ValueError: Unknown amino acid(s) in peptide: {'X'}
    """
    tokenized_peptide = tokenizer.tokenize_peptide(peptide)
    return sum([biology.n_h_acceptors[aa] for aa in tokenized_peptide])

n_h_donors(peptide)

Calculate the total number of hydrogen bond donors in the peptide.

Hydrogen bonds are important in protein-protein interactions Hubbard & Haider, 2010. The number of hydrogen bond donors for each amino acid were retrieved from PubChem.

Parameters:

Name Type Description Default
peptide str

Amino acid sequence of the peptide.

required

Returns:

Type Description
float

Total number of hydrogen bond donors in the peptide.

Examples:

>>> n_h_donors('VYP')
7
>>> n_h_donors('PGU')
Traceback (most recent call last):
    ...
ValueError: Unknown amino acid(s) in peptide: {'U'}
Source code in peptidy/descriptors.py
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
def n_h_donors(peptide: str) -> float:
    """
    Calculate the total number of hydrogen bond donors in the peptide.

    Hydrogen bonds are important in protein-protein interactions [Hubbard & Haider, 2010](https://doi.org/10.1002/9780470015902.a0003011.pub2).
    The number of hydrogen bond donors for each amino acid were retrieved from [PubChem](https://pubchem.ncbi.nlm.nih.gov/).

    Parameters
    ----------
    peptide : str
        Amino acid sequence of the peptide.

    Returns
    -------
    float
        Total number of hydrogen bond donors in the peptide.

    Examples
    --------
    >>> n_h_donors('VYP')
    7
    >>> n_h_donors('PGU')
    Traceback (most recent call last):
        ...
    ValueError: Unknown amino acid(s) in peptide: {'U'}
    """
    tokenized_peptide = tokenizer.tokenize_peptide(peptide)
    return sum([biology.n_h_donors[aa] for aa in tokenized_peptide])

topological_polar_surface_area(peptide)

Calculate the total topological polar surface area of the peptide.

The topological polar surface area relates to the Van der Waals forces Adhav & Saikrishnan, 2023.

Parameters:

Name Type Description Default
peptide str

Amino acid sequence of the peptide.

required

Returns:

Type Description
float

Total topological polar surface area of the peptide.

Examples:

>>> topological_polar_surface_area('R_dPSRMNPAWE')
853.19...
>>> topological_polar_surface_area('AYZ')
Traceback (most recent call last):
    ...
ValueError: Unknown amino acid(s) in peptide: {'Z'}
Source code in peptidy/descriptors.py
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
def topological_polar_surface_area(
    peptide: str,
) -> float:
    """
    Calculate the total topological polar surface area of the peptide.

    The topological polar surface area relates to the Van der Waals forces [Adhav & Saikrishnan, 2023](https://pubs.acs.org/doi/10.1021/acsomega.3c00205).

    Parameters
    ----------
    peptide : str
        Amino acid sequence of the peptide.

    Returns
    ----------
    float
        Total topological polar surface area of the peptide.

    Examples
    ----------
    >>> topological_polar_surface_area('R_dPSRMNPAWE')
    853.19...
    >>> topological_polar_surface_area('AYZ')
    Traceback (most recent call last):
        ...
    ValueError: Unknown amino acid(s) in peptide: {'Z'}
    """
    tokenized_peptide = tokenizer.tokenize_peptide(peptide)
    total_surface = sum([biology.tpsas[aa] for aa in tokenized_peptide])
    return total_surface

x_logp_energy(peptide)

Calculate the sum of xlogP index of the peptide divided by the length of the peptide.

The xlogP index is a measure of the hydrophobicity of the peptide Chen et al, 2007 and the indices per amino acid was retrieved from PubChem.

Parameters:

Name Type Description Default
peptide str

Amino acid sequence of the peptide.

required

Returns:

Type Description
float

xlogP index of the peptide.

Examples:

>>> x_logp_energy('R_dPSRMNPAWE')
2.9
>>> x_logp_energy('BCAF')
Traceback (most recent call last):
    ...
ValueError: Unknown amino acid(s) in peptide: {'B'}
Source code in peptidy/descriptors.py
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
def x_logp_energy(peptide: str) -> float:
    """
    Calculate the sum of xlogP index of the peptide divided by the length of the peptide.

    The xlogP index is a measure of the hydrophobicity of the peptide [Chen et al, 2007](https://pubs.acs.org/doi/10.1021/ci700257y)
    and the indices per amino acid was retrieved from [PubChem](https://pubchem.ncbi.nlm.nih.gov/).

    Parameters
    ----------
    peptide : str
        Amino acid sequence of the peptide.

    Returns
    ----------
    float
        xlogP index of the peptide.

    Examples
    ----------
    >>> x_logp_energy('R_dPSRMNPAWE')
    2.9
    >>> x_logp_energy('BCAF')
    Traceback (most recent call last):
        ...
    ValueError: Unknown amino acid(s) in peptide: {'B'}
    """
    tokenized_peptide = tokenizer.tokenize_peptide(peptide)
    energy = sum([-biology.x_logps[aa] for aa in tokenized_peptide])
    return energy / length(peptide)