API reference

class minotaor.MinotaorTranslator

Bases: object

Please install dna_features_viewer to use this class.

minotaor.add_aa_content(seqrecord, aa, window_size, cutoff, name=None)

Compute and annotate global or local content of selected amino acids.

Parameters:
  • seqrecord (SeqRecord) – The amino acid SeqRecord to annotate.

  • aa (list) – List of amino acids to search for.

  • window_size (int) – The search window size.

  • cutoff (float) – Annotate section with at least this proportion (between 0 and 1).

  • name (str, optional) – Annotation label. Default: >=#% X/Y/Z.

minotaor.add_elm(seqrecord, elm)

Annotate SeqRecord with Eukaryotic Linear Motif (ELM) database search results.

Parameters:
  • seqrecord (SeqRecord) – SeqRecord to annotate.

  • elm (pandas.DataFrame) – Dataframe of results.

minotaor.add_elm_tsv(seqrecord, elm_tsv)

Annotate SeqRecord with Eukaryotic Linear Motif (ELM) database search results.

Parameters:
  • seqrecord (SeqRecord) – SeqRecord to annotate.

  • elm_tsv (str) – Path to TSV file of results.

minotaor.add_interpro(seqrecord, interpro, hit_types=None, include_description=True)

Annotate SeqRecord with InterPro results.

Parameters:
  • seqrecord (SeqRecord) – SeqRecord to annotate.

  • interpro (QueryResult) – QueryResult object output of Bio.SearchIO.read(handle, “interproscan-xml”).

  • hit_types (list, optional) – The InterProScan hit types to filter for. Default includes all.

  • include_description (bool) – If True, includes description in the label, otherwise only in the note qualifier of the SeqRecord.

minotaor.add_scanprosite_results(seqrecord, scanprosite_record)
minotaor.annotate_record(seqrecord, seq_dataset=None)

Annotate a record with entries of a reference sequence dataset.

Note that the search is case sensitive.

Parameters:
  • seqrecord (SeqRecord) – SeqRecord to annotate.

  • seq_dataset (pandas.DataFrame) – A minotaor sequence dataset. Default uses the built-in data. The sequence and name columns are used for search and naming of the motifs. If there is a class column, then it is used for the SeqFeature’s mino_class qualifier, which is used to determine the color during plotting. If the dataframe has a description column, then its entry is added as a note qualifier.

minotaor.convert_dna_to_aa_pattern(dna)

Convert a DNA string to a list of patterns representing its translations.

Parameters:

dna (str) – DNA (ATCG)

minotaor.convert_prosite_to_regex(prosite_string)

Convert a PROSITE motif string to a regex string.

Parameters:

prosite_string (str) – The PROSITE string.

minotaor.convert_regex_to_prosite(regex)

Convert a compatible regex string to a PROSITE motif.

Parameters:

regex (str) – The regex string.

minotaor.convert_tokens_to_prosite(tokens)
minotaor.create_and_annotate_record(sequence, seq_dataset=None)

Create a SeqRecord from an amino acid sequence string.

Parameters:

sequence (str) – The amino acid sequence.

minotaor.create_postfix_regex(postfix)
minotaor.create_prefix_regex(prefix)
minotaor.evaluate_content(sequence, aa, window_size, cutoff)

Compute global or local content of selected amino acids.

minotaor.generate_postfix_codons(postfix)
minotaor.generate_prefix_codons(prefix)
minotaor.get_content(sequence, aa, window_size)

Compute proportion of selected amino acids in string.

minotaor.make_regex_from_dna(dna, prefix, postfix)

Convert three DNA strings into a regex.

The first DNA string (dna) must be divisible by 3, the length of the second (prefix) and third (postfix) must be 1 or 2.

Parameters:
  • dna (str) – DNA (ATCG).

  • prefix (str) – DNA (ATCG).

  • postfix (str) – DNA (ATCG).

minotaor.tokenize_simple_regex(regex)

Lex regex into list of tokens.