augur.align moduleο
Align multiple sequences from FASTA.
- augur.align.analyse_insertions(aln, ungapped, insertion_csv)ο
- augur.align.check_arguments(args)ο
- augur.align.check_duplicates(*values)ο
- augur.align.ensure_reference_strain_present(ref_name, existing_alignment, seqs)ο
- augur.align.generate_alignment_cmd(method, nthreads, existing_aln_fname, seqs_to_align_fname, aln_fname, log_fname)ο
- augur.align.make_gaps_ambiguous(aln)ο
replace all gaps by βNβ in all sequences in the alignment. TreeTime will treat them as fully ambiguous and replace then with the most likely state. This modifies the alignment in place.
- Parameters:
aln (Bio.Align.MultipleSeqAlignment) β Biopython Alignment
- augur.align.postprocess(output_file, ref_name, keep_reference, fill_gaps)ο
Postprocessing of the combined alignment file.
The modified alignment is written directly to output_file.
- Parameters:
output_file (str) β The file the new alignment was written to
ref_name (str) β If provided, the name of the reference strain used in the alignment
keep_reference (bool) β If the reference was provided, whether it should be kept in the alignment
fill_gaps (bool) β Replace all gaps in the alignment with βNβ to indicate ambiguous sites.
- augur.align.prepare(sequences, existing_aln_fname, output, ref_name, ref_seq_fname)ο
Prepare the sequences, existing alignment, and reference sequence for alignment.
- This function:
Combines all given input sequences into a single file
Checks to make sure the input sequences donβt overlap with the existing alignment, if one exists.
If given a reference name, check that sequence exists in either the existing alignment, if given, or the input sequences.
If given a reference sequence, either add it to the existing alignment or prepend it to the input seqeunces.
Write the input sequences to a single file, and write the alignment back out if we added the reference sequence to it.
- Parameters:
sequences (list of str) β List of paths to FASTA-formatted sequences to align.
existing_aln_fname (str) β Path of an existing alignment to use, or None
output (str) β Path the aligned sequences will be written out to.
ref_name (str) β The name of the reference sequence, if provided
ref_seq_fname (str) β The path to the reference sequence file. If this is provided, it overrides ref_name.
- Returns:
The existing alignment filename, the new sequences filename, and the name of the reference sequence.
- Return type:
- augur.align.prettify_alignment(aln)ο
Converts all bases to uppercase and removes auto reverse-complement prefix (_R_). This modifies the alignment in place.
- Parameters:
aln (Bio.Align.MultipleSeqAlignment) β Biopython Alignment
- augur.align.prune_seqs_matching_alignment(seqs, aln)ο
Return a set of seqs excluding those already in the alignment & print a warning message for each sequence which is exluded.
- augur.align.read_alignment(fname)ο
- augur.align.read_reference(ref_fname)ο
- augur.align.read_sequences(*fnames)ο
return list of sequences from all fnames
- augur.align.register_arguments(parser)ο
Add arguments to parser. Kept as a separate function than register_parser to continue to support unit tests that use this function to create argparser.
- augur.align.register_parser(parent_subparsers)ο
- augur.align.remove_reference_sequence(seqs, reference_name)ο
- augur.align.run(args)ο
- Parameters:
args (argparse.Namespace) β arguments passed in via the command-line from augur
- Returns:
returns 0 for success, 1 for general error
- Return type:
- augur.align.strip_non_reference(aln, reference, insertion_csv=None)ο
return sequences that have all insertions relative to the reference removed. The aligment is returned as list of sequences.
- Parameters:
aln (Bio.Align.MultipleSeqAlignment) β Biopython Alignment
reference (str) β name of reference sequence, assumed to be part of the alignment
- Returns:
list of trimmed sequences, effectively a multiple alignment
- Return type:
Examples
>>> [s.name for s in strip_non_reference(read_alignment("tests/data/align/test_aligned_sequences.fasta"), "with_gaps")] Trimmed gaps in with_gaps from the alignment ['with_gaps', 'no_gaps', 'some_other_seq', '_R_crick_strand'] >>> [s.name for s in strip_non_reference(read_alignment("tests/data/align/test_aligned_sequences.fasta"), "no_gaps")] No gaps in alignment to trim (with respect to the reference, no_gaps) ['with_gaps', 'no_gaps', 'some_other_seq', '_R_crick_strand'] >>> [s.name for s in strip_non_reference(read_alignment("tests/data/align/test_aligned_sequences.fasta"), "missing")] Traceback (most recent call last): ... augur.align.AlignmentError: ERROR: reference missing not found in alignment
- augur.align.write_seqs(seqs, fname)ο
A wrapper around SeqIO.write with error handling