API reference

Random access of Cortex graphs

This module contains classes for inspecting Cortex graphs with random access to their kmers.

class cortexpy.graph.parser.random_access.RandomAccess(graph_handle, kmer_cache_size=None)[source]

Provide fast k-mer access to Cortex graph in log(n) time (n = number of kmers in graph)

__getitem__(lexlo_string)[source]

Return kmer associated with kmer string

No check is performed to make sure that the input string is a lexicographically-lowest kmer string. Use get_kmer_for_string() in order to convert a kmer string to its lexlo form before retrieving it from the cortex object.

__iter__()[source]

Iterate over kmer strings in graph in order stored in graph

get_kmer_for_string(string)[source]

Will compute the revcomp of kmer string before getting a kmer

items()[source]

Iterate over kmer strings and kmers in graph in order stored in graph

values()[source]

Iterate over kmers in cortex graph

Cortex graph headers

This module contains classes for parsing and representing a Cortex file header

class cortexpy.graph.parser.header.Header(version=6, kmer_size=1, kmer_container_size=None, num_colors=None, mean_read_lengths=None, total_sequences=None, sample_names=None, error_rates=None, color_info_blocks=NOTHING)[source]

Cortex header object

This object allows access to header information contained in a cortex file

classmethod from_stream(stream)[source]

Extract a cortex header from a file handle

Cortex kmers

This module provides classes and functions for working with Cortex kmers.

class cortexpy.graph.parser.kmer.Kmer(kmer_data, num_colors, kmer_size, revcomp=None)[source]

Represents a Cortex kmer

This class wraps a kmer data object with attributes and methods for inspecting and manipulating the underlying kmer data object.

increment_color_coverage(color)[source]

Increment the coverage of a color by one

class cortexpy.graph.parser.kmer.StringKmerConverter(kmer_size)[source]

Converts kmer strings to various binary representations

to_uints(kmer_string)[source]

Converts kmer_string to big-endian uint64 array

cortexpy.graph.parser.kmer.connect_kmers(first, second, color, identical_kmer_check=True)[source]

Connect two kmers

cortexpy.graph.parser.kmer.disconnect_kmers(first, second, colors)[source]

Disconnect two kmers

cortexpy.graph.parser.kmer.find_all_neighbors(first, second)[source]

Return kmers and letters to get from first kmer to second

Utility functions

This module contains utility functions that are used inside cortexpy. These functions may also be useful outside of cortexpy.

cortexpy.utils.kmerize_contig(contig, kmer_size)[source]

Return generator of kmers in contig

The returned kmers are not lexicographically lowest.

>>> list(kmerize_contig('ATTT', 3))
['ATT', 'TTT']
cortexpy.utils.kmerize_fasta(fasta, kmer_size)[source]

Return generator to all kmers in fasta

cortexpy.utils.lexlo[source]

Return lexicographically lowest version of a kmer string and its reverse complement

The reverse complement of a kmer string is generated and the lexicographically-lowest kmer string is returned.

>>> lexlo('AAA')
'AAA'
>>> lexlo('TTT')
'AAA'