ArrayAlignment

class ArrayAlignment(*args, **kwargs)

Holds a dense array representing a multiple sequence alignment.

An Alignment is _often_, but not necessarily, an array of chars. You might want to use some other data type for the alignment if you have a large number of symbols. For example, codons on an ungapped DNA alphabet has 4*4*4=64 entries so can fit in a standard char data type, but tripeptides on the 20-letter ungapped protein alphabet has 20*20*20=8000 entries so can _not_ fit in a char and values will wrap around (i.e. you will get an unpredictable, wrong value for any item whose index is greater than the max value, e.g. 255 for uint8), so in this case you would need to use UInt16, which can hold 65536 values. DO NOT USE SIGNED DATA TYPES FOR YOUR ALIGNMENT ARRAY UNLESS YOU LOVE MISERY AND HARD-TO-DEBUG PROBLEMS.

Implementation: aln[i] returns position i in the alignment.

aln.positions[i] returns the same as aln[i] – usually, users think of this as a ‘column’, because alignment editors such as Clustal typically display each sequence as a row so a position that cuts across sequences is a column.

aln.seqs[i] returns a sequence, or ‘row’ of the alignment in standard terminology.

WARNING: aln.seqs and aln.positions are different views of the same array, so if you change one you will change the other. This will no longer be true if you assign to seqs or positions directly, so don’t do it. If you want to change the data in the whole array, always assign to a slice so that both views update: aln.seqs[:] = x instead of aln.seqs = x. If you get the two views out of sync, you will get all sorts of exceptions. No validation is performed on aln.seqs and aln.positions for performance reasons, so this can really get you into trouble.

Alignments are immutable, though this is not enforced. If you change the data after the alignment is created, all sorts of bad things might happen.

Class properties: alphabet: should be an Alphabet object. Must provide mapping between items (possibly, but not necessarily, characters) in the alignment and indices of those characters in the resulting Alignment object.

SequenceType: Constructor to use when building sequences. Default: Sequence

_input_handlers: dict of {input_type:input_handler} where input_handler is from the _input_handlers above and input_type is a result of the method self._guess_input_type (should always be a string).

Creating a new array will always result in a new object unless you use the force_same_object=True parameter.

WARNING: Rebinding the names attribute in a ArrayAlignment is not recommended because not all methods will use the updated name order. This is because the original sequence and name order are used to produce data structures that are cached for efficiency, and are not updated if you change the names attribute.

WARNING: ArrayAlignment strips off info objects from sequences that have them, primarily for efficiency.

add_from_ref_aln(ref_aln, before_name=None, after_name=None)

Insert sequence(s) to self based on their alignment to a reference sequence. Assumes the first sequence in ref_aln.names[0] is the reference.

By default the sequence is appended to the end of the alignment, this can be changed by using either before_name or after_name arguments.

Returns Alignment object of the same class.

Parameters:
  • ref_aln – reference alignment (Alignment object/series) of reference sequence and sequences to add. New sequences in ref_aln (ref_aln.names[1:] are sequences to add. If series is used as ref_aln, it must have the structure [[‘ref_name’, SEQ], [‘name’, SEQ]]

  • before_name – name of the sequence before which sequence is added

  • after_name – name of the sequence after which sequence is added If both before_name and after_name are specified seqs will be inserted using before_name.

Examples

Aln1: -AC-DEFGHI (name: seq1) XXXXXX–XX (name: seq2) YYYY-YYYYY (name: seq3)

Aln2: ACDEFGHI (name: seq1) KL–MNPR (name: seqX) KLACMNPR (name: seqY) KL–MNPR (name: seqZ)

Out: -AC-DEFGHI (name: seq1) XXXXXX–XX (name: seq2) YYYY-YYYYY (name: seq3) -KL—MNPR (name: seqX) -KL-ACMNPR (name: seqY) -KL—MNPR (name: seqZ)

add_seqs(other, before_name=None, after_name=None)

Returns new object of class self with sequences from other added.

Parameters:
  • other – same class as self or coerceable to that class

  • before_name (str) – which sequence is added

  • after_name (str) – which sequence is added

Notes

If both before_name and after_name are specified, the seqs will be inserted using before_name.

By default the sequence is appended to the end of the alignment, this can be changed by using either before_name or after_name arguments.

alignment_quality(equifreq_mprobs=True)

Computes the alignment quality for an alignment based on eq. (2) in noted reference.

Parameters:

equifreq_mprobs (bool) – If true, specifies equally frequent motif probabilities.

Notes

    1. Hertz, G. D. Stormo - Published 1999, Bioinformatics, vol. 15 pg. 563-577.

The alignment quality statistic is a log-likelihood ratio (computed using log2) of the observed alignment column freqs versus the expected.

alphabet = ('\x00', '\x01', '\x02', '\x03', '\x04', '\x05', '\x06', '\x07', '\x08', '\t', '\n', '\x0b', '\x0c', '\r', '\x0e', '\x0f', '\x10', '\x11', '\x12', '\x13', '\x14', '\x15', '\x16', '\x17', '\x18', '\x19', '\x1a', '\x1b', '\x1c', '\x1d', '\x1e', '\x1f', ' ', '!', '"', '#', '$', '%', '&', "'", '(', ')', '*', '+', ',', '-', '.', '/', '0', '1', '2', '3', '4', '5', '6', '7', '8', '9', ':', ';', '<', '=', '>', '?', '@', 'A', 'B', 'C', 'D', 'E', 'F', 'G', 'H', 'I', 'J', 'K', 'L', 'M', 'N', 'O', 'P', 'Q', 'R', 'S', 'T', 'U', 'V', 'W', 'X', 'Y', 'Z', '[', '\\', ']', '^', '_', '`', 'a', 'b', 'c', 'd', 'e', 'f', 'g', 'h', 'i', 'j', 'k', 'l', 'm', 'n', 'o', 'p', 'q', 'r', 's', 't', 'u', 'v', 'w', 'x', 'y', 'z', '{', '|', '}', '~', '\x7f', '\x80', '\x81', '\x82', '\x83', '\x84', '\x85', '\x86', '\x87', '\x88', '\x89', '\x8a', '\x8b', '\x8c', '\x8d', '\x8e', '\x8f', '\x90', '\x91', '\x92', '\x93', '\x94', '\x95', '\x96', '\x97', '\x98', '\x99', '\x9a', '\x9b', '\x9c', '\x9d', '\x9e', '\x9f', '\xa0', '¡', '¢', '£', '¤', '¥', '¦', '§', '¨', '©', 'ª', '«', '¬', '\xad', '®', '¯', '°', '±', '²', '³', '´', 'µ', '¶', '·', '¸', '¹', 'º', '»', '¼', '½', '¾', '¿', 'À', 'Á', 'Â', 'Ã', 'Ä', 'Å', 'Æ', 'Ç', 'È', 'É', 'Ê', 'Ë', 'Ì', 'Í', 'Î', 'Ï', 'Ð', 'Ñ', 'Ò', 'Ó', 'Ô', 'Õ', 'Ö', '×', 'Ø', 'Ù', 'Ú', 'Û', 'Ü', 'Ý', 'Þ', 'ß', 'à', 'á', 'â', 'ã', 'ä', 'å', 'æ', 'ç', 'è', 'é', 'ê', 'ë', 'ì', 'í', 'î', 'ï', 'ð', 'ñ', 'ò', 'ó', 'ô', 'õ', 'ö', '÷', 'ø', 'ù', 'ú', 'û', 'ü', 'ý', 'þ', 'ÿ')
apply_pssm(pssm=None, path=None, background=None, pseudocount=0, names=None, ui=None)

scores sequences using the specified pssm

Parameters:
  • pssm (profile.PSSM) – if not provided, will be loaded from path

  • path – path to either a jaspar or cisbp matrix (path must end have a suffix matching the format).

  • pseudocount – adjustment for zero in matrix

  • names – returns only scores for these sequences and in the name order

Return type:

numpy array of log2 based scores at every position

coevolution(method='nmi', segments=None, drawable=None, show_progress=False, ui=None)

performs pairwise coevolution measurement

Parameters:
  • method (str) – coevolution metric, defaults to ‘nmi’ (Normalized Mutual Information). Valid choices are ‘rmi’ (Resampled Mutual Information) and ‘mi’, mutual information.

  • segments (coordinate series) – coordinates of the form [(start, end), …] where all possible pairs of alignment positions within and between segments are examined.

  • drawable (None or str) – Result object is capable of plotting data specified type. str value must be one of plot type ‘box’, ‘heatmap’, ‘violin’.

  • show_progress (bool) – shows a progress bar

Returns:

  • DictArray of results with lower-triangular values. Upper triangular

  • elements and estimates that could not be computed for numerical reasons

  • are set as nan

copy()

Returns deep copy of self.

count_gaps_per_pos(include_ambiguity=True)

return counts of gaps per position as a DictArray

Parameters:

include_ambiguity (bool) – if True, ambiguity characters that include the gap state are included

count_gaps_per_seq(induced_by=False, unique=False, include_ambiguity=True, drawable=False)

return counts of gaps per sequence as a DictArray

Parameters:
  • induced_by (bool) – a gapped column is considered to be induced by a seq if the seq has a non-gap character in that column.

  • unique (bool) – count is limited to gaps uniquely induced by each sequence

  • include_ambiguity (bool) – if True, ambiguity characters that include the gap state are included

  • drawable (bool or str) – if True, resulting object is capable of plotting data via specified plot type ‘bar’, ‘box’ or ‘violin’

counts(motif_length=1, include_ambiguity=False, allow_gap=False, exclude_unobserved=False)

counts of motifs

Parameters:
  • motif_length – number of elements per character.

  • include_ambiguity – if True, motifs containing ambiguous characters from the seq moltype are included. No expansion of those is attempted.

  • allow_gaps – if True, motifs containing a gap character are included.

  • exclude_unobserved – if True, unobserved motif combinations are excluded.

Notes

only non-overlapping motifs are counted

counts_per_pos(motif_length=1, include_ambiguity=False, allow_gap=False, alert=False)

return DictArray of counts per position

Parameters:

alert – warns if motif_length > 1 and alignment trimmed to produce motif columns

counts_per_seq(motif_length=1, include_ambiguity=False, allow_gap=False, exclude_unobserved=False, alert=False)

counts of non-overlapping motifs per sequence

Parameters:
  • motif_length – number of elements per character.

  • include_ambiguity – if True, motifs containing ambiguous characters from the seq moltype are included. No expansion of those is attempted.

  • allow_gaps – if True, motifs containing a gap character are included.

  • exclude_unobserved – if False, all canonical states included

  • alert – warns if motif_length > 1 and alignment trimmed to produce motif columns

Return type:

MotifCountsArray

deepcopy(sliced=True)

Returns deep copy of self.

default_gap = '-'
degap(**kwargs)

Returns copy in which sequences have no gaps.

distance_matrix(calc='percent', show_progress=False, drop_invalid=False)

Returns pairwise distances between sequences.

Parameters:
  • calc (str) – a pairwise distance calculator or name of one. For options see cogent3.evolve.fast_distance.available_distances

  • show_progress (bool) – controls progress display for distance calculation

  • drop_invalid (bool) – If True, sequences for which a pairwise distance could not be calculated are excluded. If False, an ArithmeticError is raised if a distance could not be computed on observed data.

dotplot(name1=None, name2=None, window=20, threshold=None, k=None, min_gap=0, width=500, title=None, rc=False, show_progress=False)

make a dotplot between specified sequences. Random sequences chosen if names not provided.

Parameters:
  • name1 (str or None) – names of sequences. If one is not provided, a random choice is made

  • name2 (str or None) – names of sequences. If one is not provided, a random choice is made

  • window (int) – k-mer size for comparison between sequences

  • threshold (int) – windows where the sequences are identical >= threshold are a match

  • k (int) – size of k-mer to break sequences into. Larger values increase speed but reduce resolution. If not specified, is computed as the maximum of (window-threshold), (window % k) * k <= threshold.

  • min_gap (int) – permitted gap for joining adjacent line segments, default is no gap joining

  • width (int) – figure width. Figure height is computed based on the ratio of len(seq1) / len(seq2)

  • title – title for the plot

  • rc (bool or None) – include dotplot of reverse compliment also. Only applies to Nucleic acids moltypes

Return type:

a Drawable or AnnotatedDrawable

entropy_per_pos(motif_length=1, include_ambiguity=False, allow_gap=False, alert=False)

returns shannon entropy per position

entropy_per_seq(motif_length=1, include_ambiguity=False, allow_gap=False, exclude_unobserved=True, alert=False)

returns the Shannon entropy per sequence

Parameters:
  • motif_length – number of characters per tuple.

  • include_ambiguity – if True, motifs containing ambiguous characters from the seq moltype are included. No expansion of those is attempted.

  • allow_gap – if True, motifs containing a gap character are included.

  • exclude_unobserved

    if True, unobserved motif combinations are excluded.

    For motif_length > 1, it’s advisable to specify exclude_unobserved=True, this avoids unnecessary calculations.

filtered(predicate, motif_length=1, drop_remainder=True, **kwargs)

The alignment positions where predicate(column) is true.

Parameters:
  • predicate (callable) – a callback function that takes an tuple of motifs and returns True/False

  • motif_length (int) – length of the motifs the sequences should be split into, eg. 3 for filtering aligned codons.

  • drop_remainder (bool) – If length is not modulo motif_length, allow dropping the terminal remaining columns

gap_chars = {'-': None, '?': None}
get_ambiguous_positions()

Returns dict of seq:{position:char} for ambiguous chars.

Used in likelihood calculations.

get_degapped_relative_to(name)

Remove all columns with gaps in sequence with given name.

Returns Alignment object of the same class. Note that the seqs in the new Alignment are always new objects.

Parameters:

name – sequence name

get_gap_array(include_ambiguity=True)

returns bool array with gap state True, False otherwise

Parameters:

include_ambiguity (bool) – if True, ambiguity characters that include the gap state are included

get_gapped_seq(seq_name, recode_gaps=False, moltype=None)

Return a gapped Sequence object for the specified seqname.

Note: return type may depend on what data was loaded into the SequenceCollection or Alignment.

get_identical_sets(mask_degen=False)

returns sets of names for sequences that are identical

Parameters:

mask_degen – if True, degenerate characters are ignored

get_lengths(include_ambiguity=False, allow_gap=False)

returns {name: seq length, …}

Parameters:
  • include_ambiguity – if True, motifs containing ambiguous characters from the seq moltype are included. No expansion of those is attempted.

  • allow_gaps – if True, motifs containing a gap character are included.

get_motif_probs(alphabet=None, include_ambiguity=False, exclude_unobserved=False, allow_gap=False, pseudocount=0)

Return a dictionary of motif probs, calculated as the averaged frequency across sequences.

Parameters:
  • include_ambiguity – if True resolved ambiguous codes are included in estimation of frequencies, default is False.

  • exclude_unobserved – if True, motifs that are not present in the alignment are excluded from the returned dictionary, default is False.

  • allow_gap – allow gap motif

Notes

only non-overlapping motifs are counted

get_position_indices(f, native=False, negate=False)

Returns list of column indices for which f(col) is True.

fcallable

function that returns true/false given an alignment position

nativeboolean

if True, and ArrayAlignment, f is provided with slice of array otherwise the string is used

negateboolean

if True, not f() is used

get_seq(seqname)

Return a sequence object for the specified seqname.

get_seq_indices(f, negate=False)

Returns list of keys of seqs where f(row) is True.

List will be in the same order as self.names, if present.

get_similar(target, min_similarity=0.0, max_similarity=1.0, metric=<cogent3.util.transform.for_seq object>, transform=None)

Returns new Alignment containing sequences similar to target.

Parameters:
  • target – sequence object to compare to. Can be in the alignment.

  • min_similarity – minimum similarity that will be kept. Default 0.0.

  • max_similarity – maximum similarity that will be kept. Default 1.0. (Note that both min_similarity and max_similarity are inclusive.) metric similarity function to use. Must be f(first_seq, second_seq).

  • similarity (The default metric is fraction) –

  • (0% (ranging from 0.0) –

  • lots (identical) to 1.0 (100% identical). The Sequence classes have) –

  • the (of methods that can be passed in as unbound methods to act as) –

  • metric

  • frac_same_gaps. (e.g.) –

  • transform – transformation function to use on the sequences before the metric is calculated. If None, uses the whole sequences in each case. A frequent transformation is a function that returns a specified range of a sequence, e.g. eliminating the ends. Note that the transform applies to both the real sequence and the target sequence.

  • WARNING (if the transformation changes the type of the sequence (e.g.) –

  • object) (extracting a string from an RnaSequence) –

  • that (distance metrics) –

  • fail. (depend on instance data of the original class may) –

get_sub_alignment(seqs=None, pos=None, invert_seqs=False, invert_pos=False)

Returns subalignment of specified sequences and positions.

seqs and pos can be passed in as lists of sequence indices to keep or positions to keep.

invert_seqs: if True (default False), gets everything _except_ the specified sequences.

invert_pos: if True (default False), gets everything _except_ the specified positions.

Unlike most of the other code that gets things out of an alignment, this method returns a new alignment that does NOT share data with the original alignment.

get_translation(gc=None, incomplete_ok=False, **kwargs)

translate from nucleic acid to protein

Parameters:
  • gc – genetic code, either the number or name (use cogent3.core.genetic_code.available_codes)

  • incomplete_ok (bool) – codons that are mixes of nucleotide and gaps converted to ‘?’. raises a ValueError if False

  • kwargs – related to construction of the resulting object

Return type:

A new instance of self translated into protein

has_terminal_stops(gc=None, allow_partial=False)

Returns True if any sequence has a terminal stop codon.

Parameters:
  • gc – genetic code object

  • allow_partial – if True and the sequence length is not divisible by 3, ignores the 3’ terminal incomplete codon

information_plot(width=None, height=None, window=None, stat='median', include_gap=True)

plot information per position

Parameters:
  • width (int) – figure width in pixels

  • height (int) – figure height in pixels

  • window (int or None) – used for smoothing line, defaults to sqrt(length)

  • stat (str) – ‘mean’ or ‘median, used as the summary statistic for each window

  • include_gap – whether to include gap counts, shown on right y-axis

is_array = {'array', 'array_seqs'}
is_ragged()

Returns True if alignment has sequences of different lengths.

iter_positions(pos_order=None)

Iterates over positions in the alignment, in order.

pos_order refers to a list of indices (ints) specifying the column order. This lets you rearrange positions if you want to (e.g. to pull out individual codon positions).

Note that self.iter_positions() always returns new objects, by default lists of elements. Use map(f, self.iter_positions) to apply the constructor or function f to the resulting lists (f must take a single list as a parameter).

Will raise IndexError if one of the indices in order exceeds the sequence length. This will always happen on ragged alignments: assign to self.seq_len to set all sequences to the same length.

iter_selected(seq_order=None, pos_order=None)

Iterates over elements in the alignment.

seq_order (names) can be used to select a subset of seqs. pos_order (positions) can be used to select a subset of positions.

Always iterates along a seq first, then down a position (transposes normal order of a[i][j]; possibly, this should change)..

WARNING: Alignment.iter_selected() is not the same as alignment.iteritems() (which is the built-in dict iteritems that iterates over key-value pairs).

iter_seqs(seq_order=None)

Iterates over values (sequences) in the alignment, in order.

seq_order: list of keys giving the order in which seqs will be returned. Defaults to self.Names. Note that only these sequences will be returned, and that KeyError will be raised if there are sequences in order that have been deleted from the Alignment. If self.Names is None, returns the sequences in the same order as self.named_seqs.values().

Use map(f, self.seqs()) to apply the constructor f to each seq. f must accept a single list as an argument.

Always returns references to the same objects that are values of the alignment.

iupac_consensus(alphabet=None, allow_gaps=True)

Returns string containing IUPAC consensus sequence of the alignment.

majority_consensus()

Returns list containing most frequent item at each position.

Optional parameter transform gives constructor for type to which result will be converted (useful when consensus should be same type as originals).

matching_ref(ref_name, gap_fraction, gap_run)

Returns new alignment with seqs well aligned with a reference.

gap_fraction = fraction of positions that either have a gap in the

template but not in the seq or in the seq but not in the template

gap_run = number of consecutive gaps tolerated in query relative to

sequence or sequence relative to query

moltype = MolType(('\x00', '\x01', '\x02', '\x03', '\x04', '\x05', '\x06', '\x07', '\x08', '\t', '\n', '\x0b', '\x0c', '\r', '\x0e', '\x0f', '\x10', '\x11', '\x12', '\x13', '\x14', '\x15', '\x16', '\x17', '\x18', '\x19', '\x1a', '\x1b', '\x1c', '\x1d', '\x1e', '\x1f', ' ', '!', '"', '#', '$', '%', '&', "'", '(', ')', '*', '+', ',', '-', '.', '/', '0', '1', '2', '3', '4', '5', '6', '7', '8', '9', ':', ';', '<', '=', '>', '?', '@', 'A', 'B', 'C', 'D', 'E', 'F', 'G', 'H', 'I', 'J', 'K', 'L', 'M', 'N', 'O', 'P', 'Q', 'R', 'S', 'T', 'U', 'V', 'W', 'X', 'Y', 'Z', '[', '\\', ']', '^', '_', '`', 'a', 'b', 'c', 'd', 'e', 'f', 'g', 'h', 'i', 'j', 'k', 'l', 'm', 'n', 'o', 'p', 'q', 'r', 's', 't', 'u', 'v', 'w', 'x', 'y', 'z', '{', '|', '}', '~', '\x7f', '\x80', '\x81', '\x82', '\x83', '\x84', '\x85', '\x86', '\x87', '\x88', '\x89', '\x8a', '\x8b', '\x8c', '\x8d', '\x8e', '\x8f', '\x90', '\x91', '\x92', '\x93', '\x94', '\x95', '\x96', '\x97', '\x98', '\x99', '\x9a', '\x9b', '\x9c', '\x9d', '\x9e', '\x9f', '\xa0', '¡', '¢', '£', '¤', '¥', '¦', '§', '¨', '©', 'ª', '«', '¬', '\xad', '®', '¯', '°', '±', '²', '³', '´', 'µ', '¶', '·', '¸', '¹', 'º', '»', '¼', '½', '¾', '¿', 'À', 'Á', 'Â', 'Ã', 'Ä', 'Å', 'Æ', 'Ç', 'È', 'É', 'Ê', 'Ë', 'Ì', 'Í', 'Î', 'Ï', 'Ð', 'Ñ', 'Ò', 'Ó', 'Ô', 'Õ', 'Ö', '×', 'Ø', 'Ù', 'Ú', 'Û', 'Ü', 'Ý', 'Þ', 'ß', 'à', 'á', 'â', 'ã', 'ä', 'å', 'æ', 'ç', 'è', 'é', 'ê', 'ë', 'ì', 'í', 'î', 'ï', 'ð', 'ñ', 'ò', 'ó', 'ô', 'õ', 'ö', '÷', 'ø', 'ù', 'ú', 'û', 'ü', 'ý', 'þ', 'ÿ'))
property named_seqs
no_degenerates(motif_length=1, allow_gap=False)

returns new alignment without degenerate characters

Parameters:
  • motif_length – sequences are segmented into units of this size

  • allow_gaps – whether gaps are to be treated as a degenerate character (default, most evolutionary modelling treats gaps as N) or not.

property num_seqs

Returns the number of sequences in the alignment.

omit_bad_seqs(quantile=None)

Returns new alignment without sequences with a number of uniquely introduced gaps exceeding quantile

Uses count_gaps_per_seq(unique=True) to obtain the counts of gaps uniquely introduced by a sequence. The cutoff is the quantile of this distribution.

Parameters:

quantile (float or None) – sequences whose unique gap count is in a quantile larger than this cutoff are excluded. The default quantile is (num_seqs - 1) / num_seqs

omit_gap_pos(allowed_gap_frac=0.999999, motif_length=1)

Returns new alignment where all cols (motifs) have <= allowed_gap_frac gaps.

Parameters:
  • allowed_gap_frac – specifies proportion of gaps is allowed in each column (default is just < 1, i.e. only cols with at least one gap character are preserved). Set to 1 - e-6 to exclude strictly gapped columns.

  • motif_length – set’s the “column” width, e.g. setting to 3 corresponds to codons. A motif that includes a gap at any position is included in the counting. Default is 1.

omit_gap_runs(allowed_run=1)

Returns new alignment where all seqs have runs of gaps <=allowed_run.

Note that seqs with exactly allowed_run gaps are not deleted. Default is for allowed_run to be 1 (i.e. no consecutive gaps allowed).

Because the test for whether the current gap run exceeds the maximum allowed gap run is only triggered when there is at least one gap, even negative values for allowed_run will still let sequences with no gaps through.

omit_gap_seqs(allowed_gap_frac=0)

Returns new alignment with seqs that have <= allowed_gap_frac.

allowed_gap_frac should be a fraction between 0 and 1 inclusive. Default is 0.

pad_seqs(pad_length=None, **kwargs)

Returns copy in which sequences are padded to same length.

Parameters:

pad_length – Length all sequences are to be padded to. Will pad to max sequence length if pad_length is None or less than max length.

property positions

Override superclass positions to return positions as symbols.

probs_per_pos(motif_length=1, include_ambiguity=False, allow_gap=False, alert=False)

returns MotifFreqsArray per position

probs_per_seq(motif_length=1, include_ambiguity=False, allow_gap=False, exclude_unobserved=False, alert=False)

return MotifFreqsArray per sequence

Parameters:
  • motif_length – number of characters per tuple.

  • include_ambiguity – if True, motifs containing ambiguous characters from the seq moltype are included. No expansion of those is attempted.

  • allow_gap – if True, motifs containing a gap character are included.

  • exclude_unobserved – if True, unobserved motif combinations are excluded.

quick_tree(calc='percent', bootstrap=None, drop_invalid=False, show_progress=False, ui=None)

Returns pairwise distances between sequences.

Parameters:
  • calc (str) – a pairwise distance calculator or name of one. For options see cogent3.evolve.fast_distance.available_distances

  • show_progress (bool) – controls progress display for distance calculation

  • drop_invalid (bool) – If True, sequences for which a pairwise distance could not be calculated are excluded. If False, an ArithmeticError is raised if a distance could not be computed on observed data.

  • bootstrap (int or None) – Number of non-parametric bootstrap replicates. Resamples alignment columns with replacement and builds a phylogeny for each such resampling.

  • drop_invalid – If True, sequences for which a pairwise distance could not be calculated are excluded. If False, an ArithmeticError is raised if a distance could not be computed on observed data.

Returns:

  • a phylogenetic tree. If bootstrap specified, returns the weighted

  • majority consensus. Support for each node is stored as

  • edge.params[‘params’].

Notes

Sequences in the observed alignment for which distances could not be computed are omitted. Bootstrap replicates are required to have distances for all seqs present in the observed data distance matrix.

rc()

Returns the reverse complement alignment

rename_seqs(renamer)

returns new instance with sequences renamed

Parameters:

renamer (callable) – function that will take current sequences and return the new one

replace_seqs(seqs, aa_to_codon=True)

Returns new alignment with same shape but with data taken from seqs.

Parameters:
  • aa_to_codon – If True (default) aligns codons from protein alignment, or, more generally, substituting in codons from a set of protein sequences (not necessarily aligned). For this reason, it takes characters from seqs three at a time rather than one at a time (i.e. 3 characters in seqs are put in place of 1 character in self). If False, seqs must be the same lengths.

  • alignment (If seqs is an) –

  • ignored. (any gaps in it will be) –

reverse_complement()

Returns the reverse complement alignment. A synonymn for rc.

sample(n=None, with_replacement=False, motif_length=1, randint=<built-in method randint of numpy.random.mtrand.RandomState object>, permutation=<built-in method permutation of numpy.random.mtrand.RandomState object>)

Returns random sample of positions from self, e.g. to bootstrap.

Parameters:
  • n – the number of positions to sample from the alignment. Default is alignment length

  • with_replacement – boolean flag for determining if sampled positions

  • permutation (randint and) – functions for random integer in a specified range, and permutation, respectively.

Notes

By default (resampling all positions without replacement), generates a permutation of the positions of the alignment.

Setting with_replacement to True and otherwise leaving parameters as defaults generates a standard bootstrap resampling of the alignment.

returns Drawable sequence logo using mutual information

Parameters:
  • width (float) – plot dimensions in pixels

  • height (float) – plot dimensions in pixels

  • wrap (int) – number of alignment columns per row

  • vspace (float) – vertical separation between rows, as a proportion of total plot

  • colours (dict) – mapping of characters to colours. If note provided, defaults to custom for everything ecept protein, which uses protein moltype colours.

Notes

Computes MI based on log2 and includes the gap state, so the maximum possible value is -log2(1/num_states)

property seqs
set_repr_policy(num_seqs=None, num_pos=None, ref_name=None, wrap=None)

specify policy for repr(self)

Parameters:
  • num_seqs (int or None) – number of sequences to include in represented display.

  • num_pos (int or None) – length of sequences to include in represented display.

  • ref_name (str or None) – name of sequence to be placed first, or “longest” (default). If latter, indicates longest sequence will be chosen.

  • wrap (int or None) – number of printed bases per row

sliding_windows(window, step, start=None, end=None)

Generator yielding new alignments of given length and interval.

Parameters:
  • window – The length of each returned alignment.

  • step – The interval between the start of the successive alignment objects returned.

  • start – first window start position

  • end – last window start position

strand_symmetry(motif_length=1)

returns dict of strand symmetry test results per seq

take_positions(cols, negate=False)

Returns new Alignment containing only specified positions.

By default, the seqs will be lists, but an alternative constructor can be specified.

Note that take_positions will fail on ragged positions.

take_positions_if(f, negate=False)

Returns new Alignment containing cols where f(col) is True.

take_seqs(seqs, negate=False, **kwargs)

Returns new Alignment containing only specified seqs.

Note that the seqs in the new alignment will be references to the same objects as the seqs in the old alignment.

take_seqs_if(f, negate=False, **kwargs)

Returns new Alignment containing seqs where f(row) is True.

Note that the seqs in the new Alignment are the same objects as the seqs in the old Alignment, not copies.

to_dict()

Returns the alignment as dict of names -> strings.

Note: returns strings, NOT Sequence objects.

to_dna()

returns copy of self as an alignment of DNA moltype seqs

to_fasta()

Return alignment in Fasta format

Parameters:

make_seqlabel – callback function that takes the seq object and returns a label str

to_html(name_order=None, wrap=60, limit=None, ref_name='longest', colors=None, font_size=12, font_family='Lucida Console')

returns html with embedded styles for sequence colouring

Parameters:
  • name_order – order of names for display.

  • wrap – number of alignment columns per row

  • limit – truncate alignment to this length

  • ref_name – Name of an existing sequence or ‘longest’. If the latter, the longest sequence (excluding gaps and ambiguities) is selected as the reference.

  • colors – {character moltype.

  • font_size – in points. Affects labels and sequence and line spacing (proportional to value)

  • font_family – string denoting font family

  • notebook (To display in jupyter) –

    >>> from IPython.core.display import HTML
    >>> HTML(aln.to_html())
    

to_json()

returns json formatted string

to_moltype(moltype)

returns copy of self with moltype seqs

to_nexus(seq_type, wrap=50)

Return alignment in NEXUS format and mapping to sequence ids

NOTE Not that every sequence in the alignment MUST come from

a different species!! (You can concatenate multiple sequences from same species together before building tree)

seq_type: dna, rna, or protein

Raises exception if invalid alignment

to_phylip()

Return alignment in PHYLIP format and mapping to sequence ids

raises exception if invalid alignment

to_pretty(name_order=None, wrap=None)

returns a string representation of the alignment in pretty print format

Parameters:
  • name_order – order of names for display.

  • wrap – maximum number of printed bases

to_protein()

returns copy of self as an alignment of PROTEIN moltype seqs

to_rich_dict()

returns detailed content including info and moltype attributes

to_rna()

returns copy of self as an alignment of RNA moltype seqs

to_type(array_align=False, moltype=None, alphabet=None)

returns alignment of type indicated by array_align

Parameters:
  • array_align (bool) – if True, returns as ArrayAlignment. Otherwise as “standard” Alignment class. Conversion to ArrayAlignment loses annotations.

  • moltype (MolType instance) – overrides self.moltype

  • alphabet (Alphabet instance) – overrides self.alphabet

  • self) (If array_align would result in no change (class is same as) –

:param : :param returns self:

trim_stop_codons(gc=1, allow_partial=False, **kwargs)

Removes any terminal stop codons from the sequences

Parameters:
  • gc – genetic code object

  • allow_partial – if True and the sequence length is not divisible by 3, ignores the 3’ terminal incomplete codon

variable_positions(include_gap_motif=True)

Return a list of variable position indexes.

Parameters:

include_gap_motif – if False, sequences with a gap motif in a column are ignored.

with_modified_termini()

Changes the termini to include termini char instead of gapmotif.

Useful to correct the standard gap char output by most alignment programs when aligned sequences have different ends.

write(filename=None, format=None, **kwargs)

Write the alignment to a file, preserving order of sequences.

Parameters:
  • filename – name of the sequence file

  • format – format of the sequence file

Notes

If format is None, will attempt to infer format from the filename suffix.