Python API¶
Handle Graph API¶
The bdsg.handlegraph
module defines the handle graph interface.
Bindings for ::handlegraph namespace
Handles¶
The module contains definitions for different types of handles. THese are references to graph elements. A basic bdsg.handelgraph.handle_t
is a reference to a strand or orientation of a node in the graph.
- class bdsg.handlegraph.handle_t¶
Bases:
pybind11_object
Represents a traversal of a node in a graph in a particular direction
- assign(self: bdsg.handlegraph.handle_t, : bdsg.handlegraph.handle_t) bdsg.handlegraph.handle_t ¶
C++: handlegraph::handle_t::operator=(const struct handlegraph::handle_t &) –> struct handlegraph::handle_t &
- class bdsg.handlegraph.path_handle_t¶
Bases:
pybind11_object
Represents the internal id of a path entity
- assign(self: bdsg.handlegraph.path_handle_t, : bdsg.handlegraph.path_handle_t) bdsg.handlegraph.path_handle_t ¶
C++: handlegraph::path_handle_t::operator=(const struct handlegraph::path_handle_t &) –> struct handlegraph::path_handle_t &
- class bdsg.handlegraph.step_handle_t¶
Bases:
pybind11_object
A step handle is an opaque reference to a single step of an oriented node on a path in a graph
- assign(self: bdsg.handlegraph.step_handle_t, : bdsg.handlegraph.step_handle_t) bdsg.handlegraph.step_handle_t ¶
C++: handlegraph::step_handle_t::operator=(const struct handlegraph::step_handle_t &) –> struct handlegraph::step_handle_t &
- class bdsg.handlegraph.net_handle_t¶
Bases:
pybind11_object
A net handle is an opaque reference to a category of traversals of a single node, a chain, or the interior of a snarl, in the snarl decomposition of a graph.
Snarls and chains are bounded by two particular points, but the traversal may not visit both or any of them (as is the case for traversals between internal tips).
The handle refers to the snarl or chain itself and also a particular category of traversals of it. Each of the start and end of the traversal can be the start of the snarl/chain, the end of the snarl/chain, or some internal tip, for 6 distinct combinations.
For single nodes, we only have forward and reverse.
Graph Interfaces¶
The bdsg.handlegraph
module also defines a hierarchy of interfaces for graph implementations that provide different levels of features.
HandleGraph¶
The most basic is the bdsg.handlegraph.HandleGraph
, a completely immutable, unannotated graph.
- class bdsg.handlegraph.HandleGraph¶
Bases:
pybind11_object
This is the interface that a graph that uses handles needs to support. It is also the interface that users should code against.
- assign(self: bdsg.handlegraph.HandleGraph, : bdsg.handlegraph.HandleGraph) bdsg.handlegraph.HandleGraph ¶
C++: handlegraph::HandleGraph::operator=(const class handlegraph::HandleGraph &) –> class handlegraph::HandleGraph &
- edge_handle(self: bdsg.handlegraph.HandleGraph, left: bdsg.handlegraph.handle_t, right: bdsg.handlegraph.handle_t) Tuple[bdsg.handlegraph.handle_t, bdsg.handlegraph.handle_t] ¶
- A pair of handles can be used as an edge. When so used, the handles have a
canonical order and orientation.
C++: handlegraph::HandleGraph::edge_handle(const struct handlegraph::handle_t &, const struct handlegraph::handle_t &) const –> struct std::pair<struct handlegraph::handle_t, struct handlegraph::handle_t>
- flip(self: bdsg.handlegraph.HandleGraph, handle: bdsg.handlegraph.handle_t) bdsg.handlegraph.handle_t ¶
Invert the orientation of a handle (potentially without getting its ID)
C++: handlegraph::HandleGraph::flip(const struct handlegraph::handle_t &) const –> struct handlegraph::handle_t
- follow_edges(self: bdsg.handlegraph.HandleGraph, handle: bdsg.handlegraph.handle_t, go_left: bool, iteratee: Callable[[bdsg.handlegraph.handle_t], bool]) bool ¶
C++: handlegraph::HandleGraph::follow_edges(const struct handlegraph::handle_t &, bool, const class std::function<bool (const struct handlegraph::handle_t &)> &) const –> bool
- for_each_edge(*args, **kwargs)¶
Overloaded function.
for_each_edge(self: bdsg.handlegraph.HandleGraph, iteratee: Callable[[Tuple[bdsg.handlegraph.handle_t, bdsg.handlegraph.handle_t]], bool]) -> bool
for_each_edge(self: bdsg.handlegraph.HandleGraph, iteratee: Callable[[Tuple[bdsg.handlegraph.handle_t, bdsg.handlegraph.handle_t]], bool], parallel: bool) -> bool
C++: handlegraph::HandleGraph::for_each_edge(const class std::function<bool (const struct std::pair<struct handlegraph::handle_t, struct handlegraph::handle_t> &)> &, bool) const –> bool
- for_each_handle(*args, **kwargs)¶
Overloaded function.
for_each_handle(self: bdsg.handlegraph.HandleGraph, iteratee: Callable[[bdsg.handlegraph.handle_t], bool]) -> bool
for_each_handle(self: bdsg.handlegraph.HandleGraph, iteratee: Callable[[bdsg.handlegraph.handle_t], bool], parallel: bool) -> bool
C++: handlegraph::HandleGraph::for_each_handle(const class std::function<bool (const struct handlegraph::handle_t &)> &, bool) const –> bool
- forward(self: bdsg.handlegraph.HandleGraph, handle: bdsg.handlegraph.handle_t) bdsg.handlegraph.handle_t ¶
Get the locally forward version of a handle
C++: handlegraph::HandleGraph::forward(const struct handlegraph::handle_t &) const –> struct handlegraph::handle_t
- get_base(self: bdsg.handlegraph.HandleGraph, handle: bdsg.handlegraph.handle_t, index: int) str ¶
- Returns one base of a handle’s sequence, in the orientation of the
handle.
C++: handlegraph::HandleGraph::get_base(const struct handlegraph::handle_t &, unsigned long) const –> char
- get_degree(self: bdsg.handlegraph.HandleGraph, handle: bdsg.handlegraph.handle_t, go_left: bool) int ¶
- Get the number of edges on the right (go_left = false) or left (go_left
= true) side of the given handle. The default implementation is O(n) in the number of edges returned, but graph implementations that track this information more efficiently can override this method.
C++: handlegraph::HandleGraph::get_degree(const struct handlegraph::handle_t &, bool) const –> unsigned long
- get_edge_count(self: bdsg.handlegraph.HandleGraph) int ¶
- Return the total number of edges in the graph. If not overridden,
counts them all in linear time.
C++: handlegraph::HandleGraph::get_edge_count() const –> unsigned long
- get_handle(*args, **kwargs)¶
Overloaded function.
get_handle(self: bdsg.handlegraph.HandleGraph, node_id: int) -> bdsg.handlegraph.handle_t
get_handle(self: bdsg.handlegraph.HandleGraph, node_id: int, is_reverse: bool) -> bdsg.handlegraph.handle_t
Look up the handle for the node with the given ID in the given orientation
C++: handlegraph::HandleGraph::get_handle(const long long &, bool) const –> struct handlegraph::handle_t
- get_id(self: bdsg.handlegraph.HandleGraph, handle: bdsg.handlegraph.handle_t) int ¶
Get the ID from a handle
C++: handlegraph::HandleGraph::get_id(const struct handlegraph::handle_t &) const –> long long
- get_is_reverse(self: bdsg.handlegraph.HandleGraph, handle: bdsg.handlegraph.handle_t) bool ¶
Get the orientation of a handle
C++: handlegraph::HandleGraph::get_is_reverse(const struct handlegraph::handle_t &) const –> bool
- get_length(self: bdsg.handlegraph.HandleGraph, handle: bdsg.handlegraph.handle_t) int ¶
Get the length of a node
C++: handlegraph::HandleGraph::get_length(const struct handlegraph::handle_t &) const –> unsigned long
- get_node_count(self: bdsg.handlegraph.HandleGraph) int ¶
Return the number of nodes in the graph
C++: handlegraph::HandleGraph::get_node_count() const –> unsigned long
- get_sequence(self: bdsg.handlegraph.HandleGraph, handle: bdsg.handlegraph.handle_t) str ¶
- Get the sequence of a node, presented in the handle’s local forward
orientation.
C++: handlegraph::HandleGraph::get_sequence(const struct handlegraph::handle_t &) const –> std::string
- get_subsequence(self: bdsg.handlegraph.HandleGraph, handle: bdsg.handlegraph.handle_t, index: int, size: int) str ¶
- Returns a substring of a handle’s sequence, in the orientation of the
handle. If the indicated substring would extend beyond the end of the handle’s sequence, the return value is truncated to the sequence’s end. By default O(n) in the size of the handle’s sequence, but can be overriden.
C++: handlegraph::HandleGraph::get_subsequence(const struct handlegraph::handle_t &, unsigned long, unsigned long) const –> std::string
- get_total_length(self: bdsg.handlegraph.HandleGraph) int ¶
- Return the total length of all nodes in the graph, in bp. If not
overridden, loops over all nodes in linear time.
C++: handlegraph::HandleGraph::get_total_length() const –> unsigned long
- has_edge(*args, **kwargs)¶
Overloaded function.
has_edge(self: bdsg.handlegraph.HandleGraph, left: bdsg.handlegraph.handle_t, right: bdsg.handlegraph.handle_t) -> bool
- Returns true if there is an edge that allows traversal from the left
handle to the right handle. By default O(n) in the number of edges on left, but can be overridden with more efficient implementations.
C++: handlegraph::HandleGraph::has_edge(const struct handlegraph::handle_t &, const struct handlegraph::handle_t &) const –> bool
has_edge(self: bdsg.handlegraph.HandleGraph, edge: Tuple[bdsg.handlegraph.handle_t, bdsg.handlegraph.handle_t]) -> bool
Convenient wrapper of has_edge for edge_t argument.
C++: handlegraph::HandleGraph::has_edge(const struct std::pair<struct handlegraph::handle_t, struct handlegraph::handle_t> &) const –> bool
- has_node(self: bdsg.handlegraph.HandleGraph, node_id: int) bool ¶
Method to check if a node exists by ID
C++: handlegraph::HandleGraph::has_node(long long) const –> bool
- max_node_id(self: bdsg.handlegraph.HandleGraph) int ¶
- Return the largest ID in the graph, or some larger number if the
largest ID is unavailable. Return value is unspecified if the graph is empty.
C++: handlegraph::HandleGraph::max_node_id() const –> long long
- min_node_id(self: bdsg.handlegraph.HandleGraph) int ¶
- Return the smallest ID in the graph, or some smaller number if the
smallest ID is unavailable. Return value is unspecified if the graph is empty.
C++: handlegraph::HandleGraph::min_node_id() const –> long long
- traverse_edge_handle(self: bdsg.handlegraph.HandleGraph, edge: Tuple[bdsg.handlegraph.handle_t, bdsg.handlegraph.handle_t], left: bdsg.handlegraph.handle_t) bdsg.handlegraph.handle_t ¶
- Such a pair can be viewed from either inward end handle and produce the
outward handle you would arrive at.
C++: handlegraph::HandleGraph::traverse_edge_handle(const struct std::pair<struct handlegraph::handle_t, struct handlegraph::handle_t> &, const struct handlegraph::handle_t &) const –> struct handlegraph::handle_t
PathHandleGraph¶
On top of this, there is the bdsg.handlegraph.PathHandleGraph
, which allows for embedded, named paths in the graph.
- class bdsg.handlegraph.PathHandleGraph¶
Bases:
HandleGraph
,PathMetadata
This is the interface for a handle graph that stores embedded paths.
- assign(self: bdsg.handlegraph.PathHandleGraph, : bdsg.handlegraph.PathHandleGraph) bdsg.handlegraph.PathHandleGraph ¶
C++: handlegraph::PathHandleGraph::operator=(const class handlegraph::PathHandleGraph &) –> class handlegraph::PathHandleGraph &
- for_each_path_handle(self: bdsg.handlegraph.PathHandleGraph, iteratee: Callable[[bdsg.handlegraph.path_handle_t], bool]) bool ¶
C++: handlegraph::PathHandleGraph::for_each_path_handle(const class std::function<bool (const struct handlegraph::path_handle_t &)> &) const –> bool
- for_each_step_in_path(self: bdsg.handlegraph.PathHandleGraph, path: bdsg.handlegraph.path_handle_t, iteratee: Callable[[bdsg.handlegraph.step_handle_t], bool]) bool ¶
C++: handlegraph::PathHandleGraph::for_each_step_in_path(const struct handlegraph::path_handle_t &, const class std::function<bool (const struct handlegraph::step_handle_t &)> &) const –> bool
- for_each_step_on_handle(self: bdsg.handlegraph.PathHandleGraph, handle: bdsg.handlegraph.handle_t, iteratee: Callable[[bdsg.handlegraph.step_handle_t], bool]) bool ¶
C++: handlegraph::PathHandleGraph::for_each_step_on_handle(const struct handlegraph::handle_t &, const class std::function<bool (const struct handlegraph::step_handle_t &)> &) const –> bool
- get_handle_of_step(self: bdsg.handlegraph.PathHandleGraph, step_handle: bdsg.handlegraph.step_handle_t) bdsg.handlegraph.handle_t ¶
Get a node handle (node ID and orientation) from a handle to an step on a path
C++: handlegraph::PathHandleGraph::get_handle_of_step(const struct handlegraph::step_handle_t &) const –> struct handlegraph::handle_t
- get_is_circular(self: bdsg.handlegraph.PathHandleGraph, path_handle: bdsg.handlegraph.path_handle_t) bool ¶
Look up whether a path is circular
C++: handlegraph::PathHandleGraph::get_is_circular(const struct handlegraph::path_handle_t &) const –> bool
- get_next_step(self: bdsg.handlegraph.PathHandleGraph, step_handle: bdsg.handlegraph.step_handle_t) bdsg.handlegraph.step_handle_t ¶
- Returns a handle to the next step on the path. If the given step is the final step
of a non-circular path, this method has undefined behavior. In a circular path, the “last” step will loop around to the “first” step.
C++: handlegraph::PathHandleGraph::get_next_step(const struct handlegraph::step_handle_t &) const –> struct handlegraph::step_handle_t
- get_path_count(self: bdsg.handlegraph.PathHandleGraph) int ¶
Returns the number of paths stored in the graph
C++: handlegraph::PathHandleGraph::get_path_count() const –> unsigned long
- get_path_handle(self: bdsg.handlegraph.PathHandleGraph, path_name: str) bdsg.handlegraph.path_handle_t ¶
- Look up the path handle for the given path name.
The path with that name must exist.
C++: handlegraph::PathHandleGraph::get_path_handle(const std::string &) const –> struct handlegraph::path_handle_t
- get_path_handle_of_step(self: bdsg.handlegraph.PathHandleGraph, step_handle: bdsg.handlegraph.step_handle_t) bdsg.handlegraph.path_handle_t ¶
Returns a handle to the path that an step is on
C++: handlegraph::PathHandleGraph::get_path_handle_of_step(const struct handlegraph::step_handle_t &) const –> struct handlegraph::path_handle_t
- get_path_name(self: bdsg.handlegraph.PathHandleGraph, path_handle: bdsg.handlegraph.path_handle_t) str ¶
Look up the name of a path from a handle to it
C++: handlegraph::PathHandleGraph::get_path_name(const struct handlegraph::path_handle_t &) const –> std::string
- get_previous_step(self: bdsg.handlegraph.PathHandleGraph, step_handle: bdsg.handlegraph.step_handle_t) bdsg.handlegraph.step_handle_t ¶
- Returns a handle to the previous step on the path. If the given step is the first
step of a non-circular path, this method has undefined behavior. In a circular path, it will loop around from the “first” step (i.e. the one returned by path_begin) to the “last” step.
C++: handlegraph::PathHandleGraph::get_previous_step(const struct handlegraph::step_handle_t &) const –> struct handlegraph::step_handle_t
- get_step_count(*args, **kwargs)¶
Overloaded function.
get_step_count(self: bdsg.handlegraph.PathHandleGraph, path_handle: bdsg.handlegraph.path_handle_t) -> int
Returns the number of node steps in the path
C++: handlegraph::PathHandleGraph::get_step_count(const struct handlegraph::path_handle_t &) const –> unsigned long
get_step_count(self: bdsg.handlegraph.PathHandleGraph, handle: bdsg.handlegraph.handle_t) -> int
Returns the number of node steps on a handle
C++: handlegraph::PathHandleGraph::get_step_count(const struct handlegraph::handle_t &) const –> unsigned long
- has_next_step(self: bdsg.handlegraph.PathHandleGraph, step_handle: bdsg.handlegraph.step_handle_t) bool ¶
Returns true if the step is not the last step in a non-circular path.
C++: handlegraph::PathHandleGraph::has_next_step(const struct handlegraph::step_handle_t &) const –> bool
- has_path(self: bdsg.handlegraph.PathHandleGraph, path_name: str) bool ¶
Determine if a path name exists and is legal to get a path handle for.
C++: handlegraph::PathHandleGraph::has_path(const std::string &) const –> bool
- has_previous_step(self: bdsg.handlegraph.PathHandleGraph, step_handle: bdsg.handlegraph.step_handle_t) bool ¶
Returns true if the step is not the first step in a non-circular path.
C++: handlegraph::PathHandleGraph::has_previous_step(const struct handlegraph::step_handle_t &) const –> bool
- is_empty(self: bdsg.handlegraph.PathHandleGraph, path_handle: bdsg.handlegraph.path_handle_t) bool ¶
Returns true if the given path is empty, and false otherwise
C++: handlegraph::PathHandleGraph::is_empty(const struct handlegraph::path_handle_t &) const –> bool
- path_back(self: bdsg.handlegraph.PathHandleGraph, path_handle: bdsg.handlegraph.path_handle_t) bdsg.handlegraph.step_handle_t ¶
- Get a handle to the last step, which will be an arbitrary step in a circular path that
we consider “last” based on our construction of the path. If the path is empty then the implementation must return the same value as path_front_end().
C++: handlegraph::PathHandleGraph::path_back(const struct handlegraph::path_handle_t &) const –> struct handlegraph::step_handle_t
- path_begin(self: bdsg.handlegraph.PathHandleGraph, path_handle: bdsg.handlegraph.path_handle_t) bdsg.handlegraph.step_handle_t ¶
- Get a handle to the first step, which will be an arbitrary step in a circular path
that we consider “first” based on our construction of the path. If the path is empty, then the implementation must return the same value as path_end().
C++: handlegraph::PathHandleGraph::path_begin(const struct handlegraph::path_handle_t &) const –> struct handlegraph::step_handle_t
- path_end(self: bdsg.handlegraph.PathHandleGraph, path_handle: bdsg.handlegraph.path_handle_t) bdsg.handlegraph.step_handle_t ¶
- Get a handle to a fictitious position past the end of a path. This position is
returned by get_next_step for the final step in a path in a non-circular path. Note: get_next_step will NEVER return this value for a circular path.
C++: handlegraph::PathHandleGraph::path_end(const struct handlegraph::path_handle_t &) const –> struct handlegraph::step_handle_t
- path_front_end(self: bdsg.handlegraph.PathHandleGraph, path_handle: bdsg.handlegraph.path_handle_t) bdsg.handlegraph.step_handle_t ¶
- Get a handle to a fictitious position before the beginning of a path. This position is
return by get_previous_step for the first step in a path in a non-circular path. Note: get_previous_step will NEVER return this value for a circular path.
C++: handlegraph::PathHandleGraph::path_front_end(const struct handlegraph::path_handle_t &) const –> struct handlegraph::step_handle_t
- scan_path(self: bdsg.handlegraph.PathHandleGraph, path: bdsg.handlegraph.path_handle_t) handlegraph::PathForEachSocket ¶
- Returns a class with an STL-style iterator interface that can be used directly
in a for each loop like: for (handle_t handle : graph->scan_path(path)) { }
C++: handlegraph::PathHandleGraph::scan_path(const struct handlegraph::path_handle_t &) const –> class handlegraph::PathForEachSocket
- steps_of_handle(*args, **kwargs)¶
Overloaded function.
steps_of_handle(self: bdsg.handlegraph.PathHandleGraph, handle: bdsg.handlegraph.handle_t) -> std::vector<handlegraph::step_handle_t, std::allocator<handlegraph::step_handle_t> >
steps_of_handle(self: bdsg.handlegraph.PathHandleGraph, handle: bdsg.handlegraph.handle_t, match_orientation: bool) -> std::vector<handlegraph::step_handle_t, std::allocator<handlegraph::step_handle_t> >
- Returns a vector of all steps of a node on paths. Optionally restricts to
steps that match the handle in orientation.
C++: handlegraph::PathHandleGraph::steps_of_handle(const struct handlegraph::handle_t &, bool) const –> class std::vector<handlegraph::step_handle_t>
Mutable and Deletable Interfaces¶
Then for each there are versions where the underlying graph is “mutable” (meaning that material can be added to it and nodes can be split) and “deletable” (meaning that nodes and edges can actually be removed from the graph), and for bdsg.handlegraph.PathHandleGraph
there are versions where the paths can be altered.
- class bdsg.handlegraph.MutableHandleGraph¶
Bases:
HandleGraph
This is the interface for a handle graph that supports addition of new graph material.
- apply_ordering(*args, **kwargs)¶
Overloaded function.
apply_ordering(self: bdsg.handlegraph.MutableHandleGraph, order: bdsg.std.vector_handlegraph_handle_t) -> bool
apply_ordering(self: bdsg.handlegraph.MutableHandleGraph, order: bdsg.std.vector_handlegraph_handle_t, compact_ids: bool) -> bool
- Reorder the graph’s internal structure to match that given.
This sets the order that is used for iteration in functions like for_each_handle. If compact_ids is true, may (but will not necessarily) compact the id space of the graph to match the ordering, from 1->|ordering|. In other cases, node IDs will be preserved. This may be a no-op in the case of graph implementations that do not have any mechanism to maintain an ordering. This may invalidate outstanding handles. Returns true if node IDs actually were adjusted to match the given order, and false if they remain unchanged.
C++: handlegraph::MutableHandleGraph::apply_ordering(const class std::vector<handlegraph::handle_t> &, bool) –> bool
- apply_orientation(self: bdsg.handlegraph.MutableHandleGraph, handle: bdsg.handlegraph.handle_t) bdsg.handlegraph.handle_t ¶
- Alter the node that the given handle corresponds to so the orientation
indicated by the handle becomes the node’s local forward orientation. Rewrites all edges pointing to the node and the node’s sequence to reflect this. Invalidates all handles to the node (including the one passed). Returns a new, valid handle to the node in its new forward orientation. Note that it is possible for the node’s ID to change. Does not update any stored paths. May change the ordering of the underlying graph.
C++: handlegraph::MutableHandleGraph::apply_orientation(const struct handlegraph::handle_t &) –> struct handlegraph::handle_t
- assign(self: bdsg.handlegraph.MutableHandleGraph, : bdsg.handlegraph.MutableHandleGraph) bdsg.handlegraph.MutableHandleGraph ¶
C++: handlegraph::MutableHandleGraph::operator=(const class handlegraph::MutableHandleGraph &) –> class handlegraph::MutableHandleGraph &
- create_edge(*args, **kwargs)¶
Overloaded function.
create_edge(self: bdsg.handlegraph.MutableHandleGraph, left: bdsg.handlegraph.handle_t, right: bdsg.handlegraph.handle_t) -> None
- Create an edge connecting the given handles in the given order and orientations.
Ignores existing edges.
C++: handlegraph::MutableHandleGraph::create_edge(const struct handlegraph::handle_t &, const struct handlegraph::handle_t &) –> void
create_edge(self: bdsg.handlegraph.MutableHandleGraph, edge: Tuple[bdsg.handlegraph.handle_t, bdsg.handlegraph.handle_t]) -> None
Convenient wrapper for create_edge.
C++: handlegraph::MutableHandleGraph::create_edge(const struct std::pair<struct handlegraph::handle_t, struct handlegraph::handle_t> &) –> void
- create_handle(*args, **kwargs)¶
Overloaded function.
create_handle(self: bdsg.handlegraph.MutableHandleGraph, sequence: str) -> bdsg.handlegraph.handle_t
- Create a new node with the given sequence and return the handle.
The sequence may not be empty.
C++: handlegraph::MutableHandleGraph::create_handle(const std::string &) –> struct handlegraph::handle_t
create_handle(self: bdsg.handlegraph.MutableHandleGraph, sequence: str, id: int) -> bdsg.handlegraph.handle_t
- Create a new node with the given id and sequence, then return the handle.
The sequence may not be empty. The ID must be strictly greater than 0.
C++: handlegraph::MutableHandleGraph::create_handle(const std::string &, const long long &) –> struct handlegraph::handle_t
- divide_handle(*args, **kwargs)¶
Overloaded function.
divide_handle(self: bdsg.handlegraph.MutableHandleGraph, handle: bdsg.handlegraph.handle_t, offsets: bdsg.std.vector_unsigned_long) -> bdsg.std.vector_handlegraph_handle_t
- Split a handle’s underlying node at the given offsets in the handle’s
orientation. Returns all of the handles to the parts. Other handles to the node being split may be invalidated. The split pieces stay in the same local forward orientation as the original node, but the returned handles come in the order and orientation appropriate for the handle passed in. Updates stored paths. This invalidates the passed handles, but not other handles.
C++: handlegraph::MutableHandleGraph::divide_handle(const struct handlegraph::handle_t &, const class std::vector<unsigned long> &) –> class std::vector<handlegraph::handle_t>
divide_handle(self: bdsg.handlegraph.MutableHandleGraph, handle: bdsg.handlegraph.handle_t, offset: int) -> Tuple[bdsg.handlegraph.handle_t, bdsg.handlegraph.handle_t]
Specialization of divide_handle for a single division point
C++: handlegraph::MutableHandleGraph::divide_handle(const struct handlegraph::handle_t &, unsigned long) –> struct std::pair<struct handlegraph::handle_t, struct handlegraph::handle_t>
- increment_node_ids(*args, **kwargs)¶
Overloaded function.
increment_node_ids(self: bdsg.handlegraph.MutableHandleGraph, increment: int) -> None
- Add the given value to all node IDs.
Has a default implementation in terms of reassign_node_ids, but can be implemented more efficiently in some graphs.
C++: handlegraph::MutableHandleGraph::increment_node_ids(long long) –> void
increment_node_ids(self: bdsg.handlegraph.MutableHandleGraph, increment: int) -> None
This specialization for long appears to be needed to avoid confusion about nid_t
C++: handlegraph::MutableHandleGraph::increment_node_ids(long) –> void
- optimize(*args, **kwargs)¶
Overloaded function.
optimize(self: bdsg.handlegraph.MutableHandleGraph) -> None
optimize(self: bdsg.handlegraph.MutableHandleGraph, allow_id_reassignment: bool) -> None
- Adjust the representation of the graph in memory to improve performance.
Optionally, allow the node IDs to be reassigned to further improve performance. Note: Ideally, this method is called one time once there is expected to be few graph modifications in the future. This may invalidate outstanding handles.
C++: handlegraph::MutableHandleGraph::optimize(bool) –> void
- reassign_node_ids(self: bdsg.handlegraph.MutableHandleGraph, get_new_id: Callable[[int], int]) None ¶
- Renumber all node IDs using the given function, which, given an old ID, returns the new ID.
Modifies the graph in place. Invalidates all outstanding handles. If the graph supports paths, they also must be updated. The mapping function may return 0. In this case, the input ID will remain unchanged. The mapping function should not return any ID for which it would return 0.
C++: handlegraph::MutableHandleGraph::reassign_node_ids(const class std::function<long long (const long long &)> &) –> void
- set_id_increment(self: bdsg.handlegraph.MutableHandleGraph, min_id: int) None ¶
- Set a minimum id to increment the id space by, used as a hint during construction.
May have no effect on a backing implementation.
C++: handlegraph::MutableHandleGraph::set_id_increment(const long long &) –> void
- class bdsg.handlegraph.DeletableHandleGraph¶
Bases:
MutableHandleGraph
This is the interface for a handle graph that supports both addition of new nodes and edges as well as deletion of nodes and edges.
- assign(self: bdsg.handlegraph.DeletableHandleGraph, : bdsg.handlegraph.DeletableHandleGraph) bdsg.handlegraph.DeletableHandleGraph ¶
C++: handlegraph::DeletableHandleGraph::operator=(const class handlegraph::DeletableHandleGraph &) –> class handlegraph::DeletableHandleGraph &
- clear(self: bdsg.handlegraph.DeletableHandleGraph) None ¶
Remove all nodes and edges. May also remove all paths, if applicable.
C++: handlegraph::DeletableHandleGraph::clear() –> void
- destroy_edge(*args, **kwargs)¶
Overloaded function.
destroy_edge(self: bdsg.handlegraph.DeletableHandleGraph, left: bdsg.handlegraph.handle_t, right: bdsg.handlegraph.handle_t) -> None
- Remove the edge connecting the given handles in the given order and orientations.
Ignores nonexistent edges. Does not update any stored paths.
C++: handlegraph::DeletableHandleGraph::destroy_edge(const struct handlegraph::handle_t &, const struct handlegraph::handle_t &) –> void
destroy_edge(self: bdsg.handlegraph.DeletableHandleGraph, edge: Tuple[bdsg.handlegraph.handle_t, bdsg.handlegraph.handle_t]) -> None
Convenient wrapper for destroy_edge.
C++: handlegraph::DeletableHandleGraph::destroy_edge(const struct std::pair<struct handlegraph::handle_t, struct handlegraph::handle_t> &) –> void
- destroy_handle(self: bdsg.handlegraph.DeletableHandleGraph, handle: bdsg.handlegraph.handle_t) None ¶
- Remove the node belonging to the given handle and all of its edges.
Either destroys any paths in which the node participates, or leaves a “hidden”, un-iterateable handle in the path to represent the sequence of the removed node. Invalidates the destroyed handle. May be called during serial for_each_handle iteration ONLY on the node being iterated. May NOT be called during parallel for_each_handle iteration. May NOT be called on the node from which edges are being followed during follow_edges. May NOT be called during iteration over paths, if it could destroy a path. May NOT be called during iteration along a path, if it could destroy that path.
C++: handlegraph::DeletableHandleGraph::destroy_handle(const struct handlegraph::handle_t &) –> void
- truncate_handle(self: bdsg.handlegraph.DeletableHandleGraph, handle: bdsg.handlegraph.handle_t, trunc_left: bool, offset: int) bdsg.handlegraph.handle_t ¶
- Shorten a node by truncating either the left or right side of the node, relative to the orientation
of the handle, starting from a given offset along the nodes sequence. Any edges on the truncated end of the node are deleted. Returns a (possibly altered) handle to the truncated node. May invalid stored paths.
C++: handlegraph::DeletableHandleGraph::truncate_handle(const struct handlegraph::handle_t &, bool, unsigned long) –> struct handlegraph::handle_t
- class bdsg.handlegraph.MutablePathHandleGraph¶
Bases:
PathHandleGraph
,MutablePathMetadata
This is the interface for a handle graph with embedded paths where the paths can be modified. Note that if the graph can also be modified, the implementation will also need to inherit from MutableHandleGraph, via the combination MutablePathMutableHandleGraph interface. TODO: This is a very limited interface at the moment. It will probably need to be extended.
- append_step(self: bdsg.handlegraph.MutablePathHandleGraph, path: bdsg.handlegraph.path_handle_t, to_append: bdsg.handlegraph.handle_t) bdsg.handlegraph.step_handle_t ¶
- Append a visit to a node to the given path. Returns a handle to the new
final step on the path which is appended. If the path is cirular, the new step is placed between the steps considered “last” and “first” by the method path_begin. Handles to prior steps on the path, and to other paths, must remain valid.
C++: handlegraph::MutablePathHandleGraph::append_step(const struct handlegraph::path_handle_t &, const struct handlegraph::handle_t &) –> struct handlegraph::step_handle_t
- assign(self: bdsg.handlegraph.MutablePathHandleGraph, : bdsg.handlegraph.MutablePathHandleGraph) bdsg.handlegraph.MutablePathHandleGraph ¶
C++: handlegraph::MutablePathHandleGraph::operator=(const class handlegraph::MutablePathHandleGraph &) –> class handlegraph::MutablePathHandleGraph &
- create_path_handle(*args, **kwargs)¶
Overloaded function.
create_path_handle(self: bdsg.handlegraph.MutablePathHandleGraph, name: str) -> bdsg.handlegraph.path_handle_t
create_path_handle(self: bdsg.handlegraph.MutablePathHandleGraph, name: str, is_circular: bool) -> bdsg.handlegraph.path_handle_t
- Create a path with the given name. The caller must ensure that no path
with the given name exists already, or the behavior is undefined. Returns a handle to the created empty path. Handles to other paths must remain valid.
C++: handlegraph::MutablePathHandleGraph::create_path_handle(const std::string &, bool) –> struct handlegraph::path_handle_t
- destroy_path(self: bdsg.handlegraph.MutablePathHandleGraph, path_handle: bdsg.handlegraph.path_handle_t) None ¶
Destroy the given path. Invalidates handles to the path and its steps.
C++: handlegraph::MutablePathHandleGraph::destroy_path(const struct handlegraph::path_handle_t &) –> void
- pop_back_step(self: bdsg.handlegraph.MutablePathHandleGraph, path_handle: bdsg.handlegraph.path_handle_t) None ¶
Remove the last step in a path. Undefined behavior if path is empty.
C++: handlegraph::MutablePathHandleGraph::pop_back_step(const struct handlegraph::path_handle_t &) –> void
- pop_front_step(self: bdsg.handlegraph.MutablePathHandleGraph, path_handle: bdsg.handlegraph.path_handle_t) None ¶
Remove the first step in a path. Undefined behavior if path is empty.
C++: handlegraph::MutablePathHandleGraph::pop_front_step(const struct handlegraph::path_handle_t &) –> void
- prepend_step(self: bdsg.handlegraph.MutablePathHandleGraph, path: bdsg.handlegraph.path_handle_t, to_prepend: bdsg.handlegraph.handle_t) bdsg.handlegraph.step_handle_t ¶
- Prepend a visit to a node to the given path. Returns a handle to the new
first step on the path which is appended. If the path is cirular, the new step is placed between the steps considered “last” and “first” by the method path_begin. Handles to later steps on the path, and to other paths, must remain valid.
C++: handlegraph::MutablePathHandleGraph::prepend_step(const struct handlegraph::path_handle_t &, const struct handlegraph::handle_t &) –> struct handlegraph::step_handle_t
- rename_path(self: bdsg.handlegraph.MutablePathHandleGraph, path_handle: bdsg.handlegraph.path_handle_t, new_name: str) bdsg.handlegraph.path_handle_t ¶
Renames a path. Existing path_handle_t’s may become invalidated..
C++: handlegraph::MutablePathHandleGraph::rename_path(const struct handlegraph::path_handle_t &, const std::string &) –> struct handlegraph::path_handle_t
- rewrite_segment(self: bdsg.handlegraph.MutablePathHandleGraph, segment_begin: bdsg.handlegraph.step_handle_t, segment_end: bdsg.handlegraph.step_handle_t, new_segment: bdsg.std.vector_handlegraph_handle_t) Tuple[bdsg.handlegraph.step_handle_t, bdsg.handlegraph.step_handle_t] ¶
- Delete a segment of a path and rewrite it as some other sequence of
steps. Returns a pair of step_handle_t’s that indicate the range of the new segment in the path. The segment to delete should be designated by the first (begin) and past-last (end) step handles. If the step that is returned by path_begin is deleted, path_begin will now return the first step from the new segment or, in the case that the new segment is empty, the step used as segment_end. Empty ranges consist of two copies of the same step handle. Empty ranges in empty paths consist of two copies of the end sentinel handle for the path. Rewriting an empty range inserts before the provided end handle.
C++: handlegraph::MutablePathHandleGraph::rewrite_segment(const struct handlegraph::step_handle_t &, const struct handlegraph::step_handle_t &, const class std::vector<handlegraph::handle_t> &) –> struct std::pair<struct handlegraph::step_handle_t, struct handlegraph::step_handle_t>
- set_circularity(self: bdsg.handlegraph.MutablePathHandleGraph, path: bdsg.handlegraph.path_handle_t, circular: bool) None ¶
- Make a path circular or non-circular. If the path is becoming circular, the
last step is joined to the first step. If the path is becoming linear, the step considered “last” is unjoined from the step considered “first” according to the method path_begin.
C++: handlegraph::MutablePathHandleGraph::set_circularity(const struct handlegraph::path_handle_t &, bool) –> void
- class bdsg.handlegraph.MutablePathMutableHandleGraph¶
Bases:
MutablePathHandleGraph
,MutableHandleGraph
This is the interface for a graph which is mutable and which has paths which are also mutable.
- assign(self: bdsg.handlegraph.MutablePathMutableHandleGraph, : bdsg.handlegraph.MutablePathMutableHandleGraph) bdsg.handlegraph.MutablePathMutableHandleGraph ¶
C++: handlegraph::MutablePathMutableHandleGraph::operator=(const class handlegraph::MutablePathMutableHandleGraph &) –> class handlegraph::MutablePathMutableHandleGraph &
- class bdsg.handlegraph.MutablePathDeletableHandleGraph¶
Bases:
MutablePathMutableHandleGraph
,DeletableHandleGraph
This is the interface for a graph which is deletable and which has paths which are also mutable.
- assign(self: bdsg.handlegraph.MutablePathDeletableHandleGraph, : bdsg.handlegraph.MutablePathDeletableHandleGraph) bdsg.handlegraph.MutablePathDeletableHandleGraph ¶
C++: handlegraph::MutablePathDeletableHandleGraph::operator=(const class handlegraph::MutablePathDeletableHandleGraph &) –> class handlegraph::MutablePathDeletableHandleGraph &
Note that there is no bdsg.handlegraph.PathMutableHandleGraph
or bdsg.handlegraph.PathDeletableHandleGraph
; it does not make sense for the paths to be static while the graph can be modified.
Position and Ordering Interfaces¶
For paths, there is also the bdsg.handlegraph.PathPositionHandleGraph
which provides efficient random access by or lookup of base offset along each embedded path. Additionally, there is bdsg.handlegraph.VectorizableHandleGraph
which provides the same operations for a linearization of all of the graph’s bases. There is also a bdsg.handlegraph.RankedHandleGraph
interface, which provides an ordering, though not necessarily a base-level linearization, of nodes and edges.
- class bdsg.handlegraph.PathPositionHandleGraph¶
Bases:
PathHandleGraph
This is the interface for a path handle graph with path positions
- assign(self: bdsg.handlegraph.PathPositionHandleGraph, : bdsg.handlegraph.PathPositionHandleGraph) bdsg.handlegraph.PathPositionHandleGraph ¶
C++: handlegraph::PathPositionHandleGraph::operator=(const class handlegraph::PathPositionHandleGraph &) –> class handlegraph::PathPositionHandleGraph &
- for_each_step_position_on_handle(self: bdsg.handlegraph.PathPositionHandleGraph, handle: bdsg.handlegraph.handle_t, iteratee: Callable[[bdsg.handlegraph.step_handle_t, bool, int], bool]) bool ¶
- Execute an iteratee on each step on a path, along with its orientation relative to
the path (true if it is reverse the orientation of the handle on the path), and its position measured in bases of sequence along the path. Positions are always measured on the forward strand.
Iteration will stop early if the iteratee returns false. This method returns false if iteration was stopped early, else true
C++: handlegraph::PathPositionHandleGraph::for_each_step_position_on_handle(const struct handlegraph::handle_t &, const class std::function<bool (const struct handlegraph::step_handle_t &, const bool &, const unsigned long &)> &) const –> bool
- get_path_length(self: bdsg.handlegraph.PathPositionHandleGraph, path_handle: bdsg.handlegraph.path_handle_t) int ¶
Returns the length of a path measured in bases of sequence.
C++: handlegraph::PathPositionHandleGraph::get_path_length(const struct handlegraph::path_handle_t &) const –> unsigned long
- get_position_of_step(self: bdsg.handlegraph.PathPositionHandleGraph, step: bdsg.handlegraph.step_handle_t) int ¶
- Returns the position along the path of the beginning of this step measured in
bases of sequence. In a circular path, positions start at the step returned by path_begin().
C++: handlegraph::PathPositionHandleGraph::get_position_of_step(const struct handlegraph::step_handle_t &) const –> unsigned long
- get_step_at_position(self: bdsg.handlegraph.PathPositionHandleGraph, path: bdsg.handlegraph.path_handle_t, position: int) bdsg.handlegraph.step_handle_t ¶
- Returns the step at this position, measured in bases of sequence starting at
the step returned by path_begin(). If the position is past the end of the path, returns path_end().
C++: handlegraph::PathPositionHandleGraph::get_step_at_position(const struct handlegraph::path_handle_t &, const unsigned long &) const –> struct handlegraph::step_handle_t
- class bdsg.handlegraph.VectorizableHandleGraph¶
Bases:
RankedHandleGraph
Defines an interface providing a vectorization of the graph nodes and edges, which can be co-inherited alongside HandleGraph.
- assign(self: bdsg.handlegraph.VectorizableHandleGraph, : bdsg.handlegraph.VectorizableHandleGraph) bdsg.handlegraph.VectorizableHandleGraph ¶
C++: handlegraph::VectorizableHandleGraph::operator=(const class handlegraph::VectorizableHandleGraph &) –> class handlegraph::VectorizableHandleGraph &
- edge_index(self: bdsg.handlegraph.VectorizableHandleGraph, edge: Tuple[bdsg.handlegraph.handle_t, bdsg.handlegraph.handle_t]) int ¶
Return a unique index among edges in the graph
C++: handlegraph::VectorizableHandleGraph::edge_index(const struct std::pair<struct handlegraph::handle_t, struct handlegraph::handle_t> &) const –> unsigned long
- node_at_vector_offset(self: bdsg.handlegraph.VectorizableHandleGraph, offset: int) int ¶
Return the node overlapping the given offset in the implicit node vector
C++: handlegraph::VectorizableHandleGraph::node_at_vector_offset(const unsigned long &) const –> long long
- node_vector_offset(self: bdsg.handlegraph.VectorizableHandleGraph, node_id: int) int ¶
- Return the start position of the node in a (possibly implict) sorted array
constructed from the concatenation of the node sequences
C++: handlegraph::VectorizableHandleGraph::node_vector_offset(const long long &) const –> unsigned long
- class bdsg.handlegraph.RankedHandleGraph¶
Bases:
HandleGraph
Defines an interface for a handle graph that can rank its nodes and handles.
Doesn’t actually require an efficient node lookup by sequence position as in VectorizableHandleGraph.
- assign(self: bdsg.handlegraph.RankedHandleGraph, : bdsg.handlegraph.RankedHandleGraph) bdsg.handlegraph.RankedHandleGraph ¶
C++: handlegraph::RankedHandleGraph::operator=(const class handlegraph::RankedHandleGraph &) –> class handlegraph::RankedHandleGraph &
- handle_to_rank(self: bdsg.handlegraph.RankedHandleGraph, handle: bdsg.handlegraph.handle_t) int ¶
- Return the rank of a handle (ranks start at 1 and are dense, and each
orientation has its own rank). Handle ranks may not have anything to do with node ranks.
C++: handlegraph::RankedHandleGraph::handle_to_rank(const struct handlegraph::handle_t &) const –> unsigned long
- id_to_rank(self: bdsg.handlegraph.RankedHandleGraph, node_id: int) int ¶
Return the rank of a node (ranks start at 1 and are dense).
C++: handlegraph::RankedHandleGraph::id_to_rank(const long long &) const –> unsigned long
- rank_to_handle(self: bdsg.handlegraph.RankedHandleGraph, rank: int) bdsg.handlegraph.handle_t ¶
Return the handle with a given rank.
C++: handlegraph::RankedHandleGraph::rank_to_handle(const unsigned long &) const –> struct handlegraph::handle_t
- rank_to_id(self: bdsg.handlegraph.RankedHandleGraph, rank: int) int ¶
Return the node with a given rank.
C++: handlegraph::RankedHandleGraph::rank_to_id(const unsigned long &) const –> long long
Algorithm implementers are encouraged to take the least capable graph type necessary for their algorithm to function.
SerializableHandleGraph¶
Orthogonal to the mutability and paths hierarchy, there is a bdsg.handlegraph.SerializableHandleGraph
interface that is implemented by graphs that can be saved to and loaded from disk. The C++ API supports saving to and loading from C++ streams, but the Python API provides only the ability to save to or load from filenames.
- class bdsg.handlegraph.SerializableHandleGraph¶
Bases:
Serializable
- assign(self: bdsg.handlegraph.SerializableHandleGraph, : bdsg.handlegraph.SerializableHandleGraph) bdsg.handlegraph.SerializableHandleGraph ¶
C++: handlegraph::SerializableHandleGraph::operator=(const class handlegraph::SerializableHandleGraph &) –> class handlegraph::SerializableHandleGraph &
- deserialize(self: bdsg.handlegraph.Serializable, filename: str) None ¶
- Sets the contents of this object to the contents of a serialized object
from a file. The serialized object must be from the same implementation of this interface as is calling deserialize(). Can only be called on an empty object.
C++: handlegraph::Serializable::deserialize(const std::string &) –> void
- get_magic_number(self: bdsg.handlegraph.Serializable) int ¶
- Returns a number that is specific to the serialized implementation for type
checking. Does not depend on the contents of any particular instantiation (i.e. behaves as if static, but cannot be static and virtual).
C++: handlegraph::Serializable::get_magic_number() const –> unsigned int
- serialize(self: bdsg.handlegraph.Serializable, filename: str) None ¶
- Write the contents of this object to a named file. Makes sure to include
a leading magic number.
C++: handlegraph::Serializable::serialize(const std::string &) –> void
SnarlDecomposition¶
A “snarl decomposition” describes the decomposition of the graph into nested substructures known as snarls and chains. The bdsg.handlegraph.SnarlDecomposition
interface defines methods for traversing the snarl decomposition of a graph using bdsg.handelgraph.net_handle_t
.
- class bdsg.handlegraph.SnarlDecomposition¶
Bases:
pybind11_object
- assign(self: bdsg.handlegraph.SnarlDecomposition, : bdsg.handlegraph.SnarlDecomposition) bdsg.handlegraph.SnarlDecomposition ¶
C++: handlegraph::SnarlDecomposition::operator=(const class handlegraph::SnarlDecomposition &) –> class handlegraph::SnarlDecomposition &
- canonical(self: bdsg.handlegraph.SnarlDecomposition, net: bdsg.handlegraph.net_handle_t) bdsg.handlegraph.net_handle_t ¶
- Get a canonical traversal handle from any net handle. All handles to the
same net graph element have the same canonical traversal. That canonical traversal must be realizable, and might not always be start-to-end or even consistently be the same kind of traversal for different snarls, chains, or nodes. Mostly useful to normalize for equality comparisons.
C++: handlegraph::SnarlDecomposition::canonical(const struct handlegraph::net_handle_t &) const –> struct handlegraph::net_handle_t
- class endpoint_t¶
Bases:
pybind11_object
- Represents a place that a traversal can start or end. Traversals can start
or end at the start, end, or an internal tip of the thing they traverse.
Members:
START
END
TIP
- property name¶
- ends_at(self: bdsg.handlegraph.SnarlDecomposition, traversal: bdsg.handlegraph.net_handle_t) bdsg.handlegraph.SnarlDecomposition.endpoint_t ¶
Return the kind of location at which the given traversal ends.
C++: handlegraph::SnarlDecomposition::ends_at(const struct handlegraph::net_handle_t &) const –> enum handlegraph::SnarlDecomposition::endpoint_t
- ends_at_end(self: bdsg.handlegraph.SnarlDecomposition, net: bdsg.handlegraph.net_handle_t) bool ¶
Return true if the given net handle describes a category of traversal that ends at the local end of the snarl/chain/node.
C++: handlegraph::SnarlDecomposition::ends_at_end(const struct handlegraph::net_handle_t &) const –> bool
- ends_at_start(self: bdsg.handlegraph.SnarlDecomposition, net: bdsg.handlegraph.net_handle_t) bool ¶
Return true if the given net handle describes a category of traversal that ends at the local start of the snarl/chain/node.
C++: handlegraph::SnarlDecomposition::ends_at_start(const struct handlegraph::net_handle_t &) const –> bool
- ends_at_tip(self: bdsg.handlegraph.SnarlDecomposition, net: bdsg.handlegraph.net_handle_t) bool ¶
Return true if the given net handle describes a category of traversal that ends at the an internal tip in the snarl/chain. Never true for nodes.
C++: handlegraph::SnarlDecomposition::ends_at_tip(const struct handlegraph::net_handle_t &) const –> bool
- flip(self: bdsg.handlegraph.SnarlDecomposition, net: bdsg.handlegraph.net_handle_t) bdsg.handlegraph.net_handle_t ¶
- Return a net handle to the same snarl/chain/node in the opposite orientation.
No effect on tip-to-tip, start-to-start, or end-to-end net handles. Flips all the others.
C++: handlegraph::SnarlDecomposition::flip(const struct handlegraph::net_handle_t &) const –> struct handlegraph::net_handle_t
- follow_net_edges(self: bdsg.handlegraph.SnarlDecomposition, here: bdsg.handlegraph.net_handle_t, graph: bdsg.handlegraph.HandleGraph, go_left: bool, iteratee: Callable[[bdsg.handlegraph.net_handle_t], bool]) bool ¶
C++: handlegraph::SnarlDecomposition::follow_net_edges(const struct handlegraph::net_handle_t &, const class handlegraph::HandleGraph *, bool, const class std::function<bool (const struct handlegraph::net_handle_t &)> &) const –> bool
- for_each_child(self: bdsg.handlegraph.SnarlDecomposition, parent: bdsg.handlegraph.net_handle_t, iteratee: Callable[[bdsg.handlegraph.net_handle_t], bool]) bool ¶
C++: handlegraph::SnarlDecomposition::for_each_child(const struct handlegraph::net_handle_t &, const class std::function<bool (const struct handlegraph::net_handle_t &)> &) const –> bool
- for_each_tippy_child(self: bdsg.handlegraph.SnarlDecomposition, parent: bdsg.handlegraph.net_handle_t, iteratee: Callable[[bdsg.handlegraph.net_handle_t], bool]) bool ¶
C++: handlegraph::SnarlDecomposition::for_each_tippy_child(const struct handlegraph::net_handle_t &, const class std::function<bool (const struct handlegraph::net_handle_t &)> &) const –> bool
- for_each_traversal(self: bdsg.handlegraph.SnarlDecomposition, item: bdsg.handlegraph.net_handle_t, iteratee: Callable[[bdsg.handlegraph.net_handle_t], bool]) bool ¶
C++: handlegraph::SnarlDecomposition::for_each_traversal(const struct handlegraph::net_handle_t &, const class std::function<bool (const struct handlegraph::net_handle_t &)> &) const –> bool
- for_each_traversal_end(self: bdsg.handlegraph.SnarlDecomposition, traversal: bdsg.handlegraph.net_handle_t, iteratee: Callable[[bdsg.handlegraph.net_handle_t], bool]) bool ¶
C++: handlegraph::SnarlDecomposition::for_each_traversal_end(const struct handlegraph::net_handle_t &, const class std::function<bool (const struct handlegraph::net_handle_t &)> &) const –> bool
- for_each_traversal_start(self: bdsg.handlegraph.SnarlDecomposition, traversal: bdsg.handlegraph.net_handle_t, iteratee: Callable[[bdsg.handlegraph.net_handle_t], bool]) bool ¶
C++: handlegraph::SnarlDecomposition::for_each_traversal_start(const struct handlegraph::net_handle_t &, const class std::function<bool (const struct handlegraph::net_handle_t &)> &) const –> bool
- get_bound(self: bdsg.handlegraph.SnarlDecomposition, snarl: bdsg.handlegraph.net_handle_t, get_end: bool, face_in: bool) bdsg.handlegraph.net_handle_t ¶
- Get the bounding handle for the snarl or chain referenced by the given
net handle, getting the start or end facing in or out as appropriate.
For snarls, returns the bounding sentinel net handles. For chains, returns net handles for traversals of the bounding nodes of the chain.
Ignores traversal type.
May not be called on traversals of individual nodes.
C++: handlegraph::SnarlDecomposition::get_bound(const struct handlegraph::net_handle_t &, bool, bool) const –> struct handlegraph::net_handle_t
- get_end_bound(self: bdsg.handlegraph.SnarlDecomposition, parent: bdsg.handlegraph.net_handle_t) bdsg.handlegraph.net_handle_t ¶
- Get a handle to the outward-facing traversal of the last node in a chain
or the end boundary in a snarl.
This isn’t necessarily where the traversal specified by the given handle actually ends (it may be start to start or tip to tip, for example).
C++: handlegraph::SnarlDecomposition::get_end_bound(const struct handlegraph::net_handle_t &) const –> struct handlegraph::net_handle_t
- get_handle(self: bdsg.handlegraph.SnarlDecomposition, net: bdsg.handlegraph.net_handle_t, graph: bdsg.handlegraph.HandleGraph) bdsg.handlegraph.handle_t ¶
- For a net handle to a traversal of a single node, get the handle for that node in the orientation it is traversed.
May not be called for other net handles.
C++: handlegraph::SnarlDecomposition::get_handle(const struct handlegraph::net_handle_t &, const class handlegraph::HandleGraph *) const –> struct handlegraph::handle_t
- get_net(self: bdsg.handlegraph.SnarlDecomposition, handle: bdsg.handlegraph.handle_t, graph: bdsg.handlegraph.HandleGraph) bdsg.handlegraph.net_handle_t ¶
Turn a handle to an oriented node into a net handle for a start-to-end or end-to-start traversal of the node, as appropriate.
C++: handlegraph::SnarlDecomposition::get_net(const struct handlegraph::handle_t &, const class handlegraph::HandleGraph *) const –> struct handlegraph::net_handle_t
- get_parent(self: bdsg.handlegraph.SnarlDecomposition, child: bdsg.handlegraph.net_handle_t) bdsg.handlegraph.net_handle_t ¶
- Get the parent snarl of a chain, or the parent chain of a snarl or node.
If the child is start-to-end or end-to-start, and the parent is a chain, the chain comes out facing the same way, accounting for the relative orientation of the child snarl or node in the chain. Otherwise, everything is produced as start-to-end, even if that is not actually a realizable traversal of a snarl or chain. May not be called on the root snarl.
Also works on snarl boundary sentinels.
C++: handlegraph::SnarlDecomposition::get_parent(const struct handlegraph::net_handle_t &) const –> struct handlegraph::net_handle_t
- get_parent_traversal(self: bdsg.handlegraph.SnarlDecomposition, traversal_start: bdsg.handlegraph.net_handle_t, traversal_end: bdsg.handlegraph.net_handle_t) bdsg.handlegraph.net_handle_t ¶
- Get a net handle for traversals of a the snarl or chain that contains
the given oriented bounding node traversals or sentinels. Given two sentinels for a snarl, produces a net handle to a start-to-end, end-to-end, end-to-start, or start-to-start traversal of that snarl. Given handles to traversals of the bounding nodes of a chain, similarly produces a net handle to a traversal of the chain.
For a chain, either or both handles can also be a snarl containing tips, for a tip-to-start, tip-to-end, start-to-tip, end-to-tip, or tip-to-tip traversal. Similarly, for a snarl, either or both handles can be a chain in the snarl that contains internal tips, or that has no edges on the appropriate end.
May only be called if a path actually exists between the given start and end.
C++: handlegraph::SnarlDecomposition::get_parent_traversal(const struct handlegraph::net_handle_t &, const struct handlegraph::net_handle_t &) const –> struct handlegraph::net_handle_t
- get_root(self: bdsg.handlegraph.SnarlDecomposition) bdsg.handlegraph.net_handle_t ¶
- Get a net handle referring to a tip-to-tip traversal of the contents of the root snarl.
TODO: Special handling for circular things in the root snarl? Circular traversal type?
C++: handlegraph::SnarlDecomposition::get_root() const –> struct handlegraph::net_handle_t
- get_start_bound(self: bdsg.handlegraph.SnarlDecomposition, parent: bdsg.handlegraph.net_handle_t) bdsg.handlegraph.net_handle_t ¶
- Get a handle to the inward-facing traversal of the first node in a chain
or the start boundary in a snarl.
This isn’t necessarily where the traversal specified by the given handle actually starts (it may be end to end or tip to tip, for example).
C++: handlegraph::SnarlDecomposition::get_start_bound(const struct handlegraph::net_handle_t &) const –> struct handlegraph::net_handle_t
- is_chain(self: bdsg.handlegraph.SnarlDecomposition, net: bdsg.handlegraph.net_handle_t) bool ¶
Returns true if the given net handle refers to (a traversal of) a chain.
C++: handlegraph::SnarlDecomposition::is_chain(const struct handlegraph::net_handle_t &) const –> bool
- is_node(self: bdsg.handlegraph.SnarlDecomposition, net: bdsg.handlegraph.net_handle_t) bool ¶
Returns true if the given net handle refers to (a traversal of) a single node, and thus has a corresponding handle_t.
C++: handlegraph::SnarlDecomposition::is_node(const struct handlegraph::net_handle_t &) const –> bool
- is_root(self: bdsg.handlegraph.SnarlDecomposition, net: bdsg.handlegraph.net_handle_t) bool ¶
- Return true if the given handle refers to (a traversal of) the root
snarl, and false otherwise.
C++: handlegraph::SnarlDecomposition::is_root(const struct handlegraph::net_handle_t &) const –> bool
- is_sentinel(self: bdsg.handlegraph.SnarlDecomposition, net: bdsg.handlegraph.net_handle_t) bool ¶
- Return true if the given net handle is a snarl bound sentinel (in either
inward or outward orientation), and false otherwise.
C++: handlegraph::SnarlDecomposition::is_sentinel(const struct handlegraph::net_handle_t &) const –> bool
- is_snarl(self: bdsg.handlegraph.SnarlDecomposition, net: bdsg.handlegraph.net_handle_t) bool ¶
Returns true if the given net handle refers to (a traversal of) a snarl.
C++: handlegraph::SnarlDecomposition::is_snarl(const struct handlegraph::net_handle_t &) const –> bool
- starts_at(self: bdsg.handlegraph.SnarlDecomposition, traversal: bdsg.handlegraph.net_handle_t) bdsg.handlegraph.SnarlDecomposition.endpoint_t ¶
Return the kind of location at which the given traversal starts.
C++: handlegraph::SnarlDecomposition::starts_at(const struct handlegraph::net_handle_t &) const –> enum handlegraph::SnarlDecomposition::endpoint_t
- starts_at_end(self: bdsg.handlegraph.SnarlDecomposition, net: bdsg.handlegraph.net_handle_t) bool ¶
Return true if the given net handle describes a category of traversal that starts at the local end of the snarl/chain/node.
C++: handlegraph::SnarlDecomposition::starts_at_end(const struct handlegraph::net_handle_t &) const –> bool
- starts_at_start(self: bdsg.handlegraph.SnarlDecomposition, net: bdsg.handlegraph.net_handle_t) bool ¶
Return true if the given net handle describes a category of traversal that starts at the local start of the snarl/chain/node.
C++: handlegraph::SnarlDecomposition::starts_at_start(const struct handlegraph::net_handle_t &) const –> bool
- starts_at_tip(self: bdsg.handlegraph.SnarlDecomposition, net: bdsg.handlegraph.net_handle_t) bool ¶
Return true if the given net handle describes a category of traversal that starts at an internal tip in the snarl/chain. Never true for nodes.
C++: handlegraph::SnarlDecomposition::starts_at_tip(const struct handlegraph::net_handle_t &) const –> bool
libbdsg Handle Graph Implementations¶
The bdsg.bdsg
module provides useful implementations of the Handle Graph API.
Bindings for ::bdsg namespace
Full Graph Implementations¶
There are two full graph implementations in the module: bdsg.bdsg.PackedGraph
and bdsg.bdsg.HashGraph
. Previously, a third implementation, ODGI, was provided, but that implementation is now part of its own odgi project.
PackedGraph¶
- class bdsg.bdsg.PackedGraph¶
Bases:
GraphProxy_bdsg_BasePackedGraph_t
- assign(self: bdsg.bdsg.PackedGraph, : bdsg.bdsg.PackedGraph) bdsg.bdsg.PackedGraph ¶
C++: bdsg::PackedGraph::operator=(const class bdsg::PackedGraph &) –> class bdsg::PackedGraph &
HashGraph¶
- class bdsg.bdsg.HashGraph¶
Bases:
MutablePathDeletableHandleGraph
,SerializableHandleGraph
HashGraph is a HandleGraph implementation designed for simplicity. Nodes are plain C++ objects stored in a hash map, which contain C++ vectors representing their adjacencies.
HashGraph is a good choice when fast access to or modification of a graph is required, but can use more memory than other graph implementations.
- append_step(self: bdsg.bdsg.HashGraph, path: bdsg.handlegraph.path_handle_t, to_append: bdsg.handlegraph.handle_t) bdsg.handlegraph.step_handle_t ¶
- Append a visit to a node to the given path. Returns a handle to the new
final step on the path which is appended. If the path is cirular, the new step is placed between the steps considered “last” and “first” by the method path_begin. Handles to prior steps on the path, and to other paths, must remain valid.
C++: bdsg::HashGraph::append_step(const struct handlegraph::path_handle_t &, const struct handlegraph::handle_t &) –> struct handlegraph::step_handle_t
- apply_ordering(*args, **kwargs)¶
Overloaded function.
apply_ordering(self: bdsg.bdsg.HashGraph, order: bdsg.std.vector_handlegraph_handle_t) -> bool
apply_ordering(self: bdsg.bdsg.HashGraph, order: bdsg.std.vector_handlegraph_handle_t, compact_ids: bool) -> bool
- Reorder the graph’s internal structure to match that given.
This sets the order that is used for iteration in functions like for_each_handle. If compact_ids is true, may (but will not necessarily) compact the id space of the graph to match the ordering, from 1->|ordering|. In other cases, node IDs will be preserved. This may be a no-op in the case of graph implementations that do not have any mechanism to maintain an ordering. This may invalidate outstanding handles. Returns true if node IDs actually were adjusted to match the given order, and false if they remain unchanged.
C++: bdsg::HashGraph::apply_ordering(const class std::vector<handlegraph::handle_t> &, bool) –> bool
- apply_orientation(self: bdsg.bdsg.HashGraph, handle: bdsg.handlegraph.handle_t) bdsg.handlegraph.handle_t ¶
- Alter the node that the given handle corresponds to so the orientation
indicated by the handle becomes the node’s local forward orientation. Rewrites all edges pointing to the node and the node’s sequence to reflect this. Invalidates all handles to the node (including the one passed). Returns a new, valid handle to the node in its new forward orientation. Note that it is possible for the node’s ID to change. Does not update any stored paths. May change the ordering of the underlying graph.
C++: bdsg::HashGraph::apply_orientation(const struct handlegraph::handle_t &) –> struct handlegraph::handle_t
- assign(self: bdsg.bdsg.HashGraph, other: bdsg.bdsg.HashGraph) bdsg.bdsg.HashGraph ¶
C++: bdsg::HashGraph::operator=(const class bdsg::HashGraph &) –> class bdsg::HashGraph &
- clear(self: bdsg.bdsg.HashGraph) None ¶
Remove all nodes and edges. Does not update any stored paths.
C++: bdsg::HashGraph::clear() –> void
- create_edge(self: bdsg.bdsg.HashGraph, left: bdsg.handlegraph.handle_t, right: bdsg.handlegraph.handle_t) None ¶
- Create an edge connecting the given handles in the given order and orientations.
Ignores existing edges.
C++: bdsg::HashGraph::create_edge(const struct handlegraph::handle_t &, const struct handlegraph::handle_t &) –> void
- create_handle(*args, **kwargs)¶
Overloaded function.
create_handle(self: bdsg.bdsg.HashGraph, sequence: str) -> bdsg.handlegraph.handle_t
- Create a new node with the given sequence and return the handle.
The sequence may not be empty.
C++: bdsg::HashGraph::create_handle(const std::string &) –> struct handlegraph::handle_t
create_handle(self: bdsg.bdsg.HashGraph, sequence: str, id: int) -> bdsg.handlegraph.handle_t
- Create a new node with the given id and sequence, then return the handle.
The sequence may not be empty. The ID must be strictly greater than 0.
C++: bdsg::HashGraph::create_handle(const std::string &, const long long &) –> struct handlegraph::handle_t
- create_path_handle(*args, **kwargs)¶
Overloaded function.
create_path_handle(self: bdsg.bdsg.HashGraph, name: str) -> bdsg.handlegraph.path_handle_t
create_path_handle(self: bdsg.bdsg.HashGraph, name: str, is_circular: bool) -> bdsg.handlegraph.path_handle_t
- Create a path with the given name. The caller must ensure that no path
with the given name exists already, or the behavior is undefined. Returns a handle to the created empty path. Handles to other paths must remain valid.
C++: bdsg::HashGraph::create_path_handle(const std::string &, bool) –> struct handlegraph::path_handle_t
- destroy_edge(self: bdsg.bdsg.HashGraph, left: bdsg.handlegraph.handle_t, right: bdsg.handlegraph.handle_t) None ¶
- Remove the edge connecting the given handles in the given order and orientations.
Ignores nonexistent edges. Does not update any stored paths.
C++: bdsg::HashGraph::destroy_edge(const struct handlegraph::handle_t &, const struct handlegraph::handle_t &) –> void
- destroy_handle(self: bdsg.bdsg.HashGraph, handle: bdsg.handlegraph.handle_t) None ¶
- Remove the node belonging to the given handle and all of its edges.
Destroys any paths in which the node participates. Invalidates the destroyed handle. May be called during serial for_each_handle iteration ONLY on the node being iterated. May NOT be called during parallel for_each_handle iteration. May NOT be called on the node from which edges are being followed during follow_edges. May NOT be called during iteration over paths, if it would destroy a path. May NOT be called during iteration along a path, if it would destroy that path.
C++: bdsg::HashGraph::destroy_handle(const struct handlegraph::handle_t &) –> void
- destroy_path(self: bdsg.bdsg.HashGraph, path: bdsg.handlegraph.path_handle_t) None ¶
Destroy the given path. Invalidates handles to the path and its node steps.
C++: bdsg::HashGraph::destroy_path(const struct handlegraph::path_handle_t &) –> void
- divide_handle(self: bdsg.bdsg.HashGraph, handle: bdsg.handlegraph.handle_t, offsets: bdsg.std.vector_unsigned_long) bdsg.std.vector_handlegraph_handle_t ¶
- Split a handle’s underlying node at the given offsets in the handle’s
orientation. Returns all of the handles to the parts. Other handles to the node being split may be invalidated. The split pieces stay in the same local forward orientation as the original node, but the returned handles come in the order and orientation appropriate for the handle passed in. Updates stored paths.
C++: bdsg::HashGraph::divide_handle(const struct handlegraph::handle_t &, const class std::vector<unsigned long> &) –> class std::vector<handlegraph::handle_t>
- flip(self: bdsg.bdsg.HashGraph, handle: bdsg.handlegraph.handle_t) bdsg.handlegraph.handle_t ¶
Invert the orientation of a handle (potentially without getting its ID)
C++: bdsg::HashGraph::flip(const struct handlegraph::handle_t &) const –> struct handlegraph::handle_t
- follow_edges_impl(self: bdsg.bdsg.HashGraph, handle: bdsg.handlegraph.handle_t, go_left: bool, iteratee: Callable[[bdsg.handlegraph.handle_t], bool]) bool ¶
- Loop over all the handles to next/previous (right/left) nodes. Passes
them to a callback which returns false to stop iterating and true to continue. Returns true if we finished and false if we stopped early.
C++: bdsg::HashGraph::follow_edges_impl(const struct handlegraph::handle_t &, bool, const class std::function<bool (const struct handlegraph::handle_t &)> &) const –> bool
- for_each_handle_impl(*args, **kwargs)¶
Overloaded function.
for_each_handle_impl(self: bdsg.bdsg.HashGraph, iteratee: Callable[[bdsg.handlegraph.handle_t], bool]) -> bool
for_each_handle_impl(self: bdsg.bdsg.HashGraph, iteratee: Callable[[bdsg.handlegraph.handle_t], bool], parallel: bool) -> bool
- Loop over all the nodes in the graph in their local forward
orientations, in their internal stored order. Stop if the iteratee returns false. Can be told to run in parallel, in which case stopping after a false return value is on a best-effort basis and iteration order is not defined.
C++: bdsg::HashGraph::for_each_handle_impl(const class std::function<bool (const struct handlegraph::handle_t &)> &, bool) const –> bool
- for_each_path_handle_impl(self: bdsg.bdsg.HashGraph, iteratee: Callable[[bdsg.handlegraph.path_handle_t], bool]) bool ¶
Execute a function on each path in the graph
C++: bdsg::HashGraph::for_each_path_handle_impl(const class std::function<bool (const struct handlegraph::path_handle_t &)> &) const –> bool
- for_each_step_on_handle_impl(self: bdsg.bdsg.HashGraph, handle: bdsg.handlegraph.handle_t, iteratee: Callable[[bdsg.handlegraph.step_handle_t], bool]) bool ¶
Calls a function with all steps of a node on paths.
C++: bdsg::HashGraph::for_each_step_on_handle_impl(const struct handlegraph::handle_t &, const class std::function<bool (const struct handlegraph::step_handle_t &)> &) const –> bool
- get_base(self: bdsg.bdsg.HashGraph, handle: bdsg.handlegraph.handle_t, index: int) str ¶
- Returns one base of a handle’s sequence, in the orientation of the
handle.
C++: bdsg::HashGraph::get_base(const struct handlegraph::handle_t &, unsigned long) const –> char
- get_degree(self: bdsg.bdsg.HashGraph, handle: bdsg.handlegraph.handle_t, go_left: bool) int ¶
Efficiently get the number of edges attached to one side of a handle.
C++: bdsg::HashGraph::get_degree(const struct handlegraph::handle_t &, bool) const –> unsigned long
- get_handle(*args, **kwargs)¶
Overloaded function.
get_handle(self: bdsg.bdsg.HashGraph, node_id: int) -> bdsg.handlegraph.handle_t
get_handle(self: bdsg.bdsg.HashGraph, node_id: int, is_reverse: bool) -> bdsg.handlegraph.handle_t
Look up the handle for the node with the given ID in the given orientation
C++: bdsg::HashGraph::get_handle(const long long &, bool) const –> struct handlegraph::handle_t
- get_handle_of_step(self: bdsg.bdsg.HashGraph, step_handle: bdsg.handlegraph.step_handle_t) bdsg.handlegraph.handle_t ¶
Get a node handle (node ID and orientation) from a handle to a step on a path
C++: bdsg::HashGraph::get_handle_of_step(const struct handlegraph::step_handle_t &) const –> struct handlegraph::handle_t
- get_id(self: bdsg.bdsg.HashGraph, handle: bdsg.handlegraph.handle_t) int ¶
Get the ID from a handle
C++: bdsg::HashGraph::get_id(const struct handlegraph::handle_t &) const –> long long
- get_is_circular(self: bdsg.bdsg.HashGraph, path_handle: bdsg.handlegraph.path_handle_t) bool ¶
Look up whether a path is circular
C++: bdsg::HashGraph::get_is_circular(const struct handlegraph::path_handle_t &) const –> bool
- get_is_reverse(self: bdsg.bdsg.HashGraph, handle: bdsg.handlegraph.handle_t) bool ¶
Get the orientation of a handle
C++: bdsg::HashGraph::get_is_reverse(const struct handlegraph::handle_t &) const –> bool
- get_length(self: bdsg.bdsg.HashGraph, handle: bdsg.handlegraph.handle_t) int ¶
Get the length of a node
C++: bdsg::HashGraph::get_length(const struct handlegraph::handle_t &) const –> unsigned long
- get_magic_number(self: bdsg.bdsg.HashGraph) int ¶
Returns a static high-entropy number to indicate the class
C++: bdsg::HashGraph::get_magic_number() const –> unsigned int
- get_next_step(self: bdsg.bdsg.HashGraph, step_handle: bdsg.handlegraph.step_handle_t) bdsg.handlegraph.step_handle_t ¶
- Returns a handle to the next step on the path. If the given step is the final step
of a non-circular path, returns the past-the-last step that is also returned by path_end. In a circular path, the “last” step will loop around to the “first” (i.e. the one returned by path_begin). Note: to iterate over each step one time, even in a circular path, consider for_each_step_in_path.
C++: bdsg::HashGraph::get_next_step(const struct handlegraph::step_handle_t &) const –> struct handlegraph::step_handle_t
- get_node_count(self: bdsg.bdsg.HashGraph) int ¶
- Return the number of nodes in the graph
TODO: can’t be node_count because XG has a field named node_count.
C++: bdsg::HashGraph::get_node_count() const –> unsigned long
- get_path_count(self: bdsg.bdsg.HashGraph) int ¶
Returns the number of paths stored in the graph
C++: bdsg::HashGraph::get_path_count() const –> unsigned long
- get_path_handle(self: bdsg.bdsg.HashGraph, path_name: str) bdsg.handlegraph.path_handle_t ¶
- Look up the path handle for the given path name.
The path with that name must exist.
C++: bdsg::HashGraph::get_path_handle(const std::string &) const –> struct handlegraph::path_handle_t
- get_path_handle_of_step(self: bdsg.bdsg.HashGraph, step_handle: bdsg.handlegraph.step_handle_t) bdsg.handlegraph.path_handle_t ¶
Returns a handle to the path that a step is on
C++: bdsg::HashGraph::get_path_handle_of_step(const struct handlegraph::step_handle_t &) const –> struct handlegraph::path_handle_t
- get_path_name(self: bdsg.bdsg.HashGraph, path_handle: bdsg.handlegraph.path_handle_t) str ¶
Look up the name of a path from a handle to it
C++: bdsg::HashGraph::get_path_name(const struct handlegraph::path_handle_t &) const –> std::string
- get_previous_step(self: bdsg.bdsg.HashGraph, step_handle: bdsg.handlegraph.step_handle_t) bdsg.handlegraph.step_handle_t ¶
- Returns a handle to the previous step on the path. If the given step is the first
step of a non-circular path, this method has undefined behavior. In a circular path, it will loop around from the “first” step (i.e. the one returned by path_begin) to the “last” step. Note: to iterate over each step one time, even in a circular path, consider for_each_step_in_path.
C++: bdsg::HashGraph::get_previous_step(const struct handlegraph::step_handle_t &) const –> struct handlegraph::step_handle_t
- get_sequence(self: bdsg.bdsg.HashGraph, handle: bdsg.handlegraph.handle_t) str ¶
Get the sequence of a node, presented in the handle’s local forward orientation.
C++: bdsg::HashGraph::get_sequence(const struct handlegraph::handle_t &) const –> std::string
- get_step_count(self: bdsg.bdsg.HashGraph, path_handle: bdsg.handlegraph.path_handle_t) int ¶
Returns the number of node steps in the path
C++: bdsg::HashGraph::get_step_count(const struct handlegraph::path_handle_t &) const –> unsigned long
- get_subsequence(self: bdsg.bdsg.HashGraph, handle: bdsg.handlegraph.handle_t, index: int, size: int) str ¶
- Returns a substring of a handle’s sequence, in the orientation of the
handle. If the indicated substring would extend beyond the end of the handle’s sequence, the return value is truncated to the sequence’s end.
C++: bdsg::HashGraph::get_subsequence(const struct handlegraph::handle_t &, unsigned long, unsigned long) const –> std::string
- has_next_step(self: bdsg.bdsg.HashGraph, step_handle: bdsg.handlegraph.step_handle_t) bool ¶
Returns true if the step is not the last step in a non-circular path.
C++: bdsg::HashGraph::has_next_step(const struct handlegraph::step_handle_t &) const –> bool
- has_node(self: bdsg.bdsg.HashGraph, node_id: int) bool ¶
Method to check if a node exists by ID
C++: bdsg::HashGraph::has_node(long long) const –> bool
- has_path(self: bdsg.bdsg.HashGraph, path_name: str) bool ¶
Determine if a path name exists and is legal to get a path handle for.
C++: bdsg::HashGraph::has_path(const std::string &) const –> bool
- has_previous_step(self: bdsg.bdsg.HashGraph, step_handle: bdsg.handlegraph.step_handle_t) bool ¶
Returns true if the step is not the first step in a non-circular path.
C++: bdsg::HashGraph::has_previous_step(const struct handlegraph::step_handle_t &) const –> bool
- increment_node_ids(self: bdsg.bdsg.HashGraph, increment: int) None ¶
Add the given value to all node IDs
C++: bdsg::HashGraph::increment_node_ids(long long) –> void
- max_node_id(self: bdsg.bdsg.HashGraph) int ¶
- Return the largest ID in the graph, or some larger number if the
largest ID is unavailable. Return value is unspecified if the graph is empty.
C++: bdsg::HashGraph::max_node_id() const –> long long
- min_node_id(self: bdsg.bdsg.HashGraph) int ¶
- Return the smallest ID in the graph, or some smaller number if the
smallest ID is unavailable. Return value is unspecified if the graph is empty.
C++: bdsg::HashGraph::min_node_id() const –> long long
- optimize(*args, **kwargs)¶
Overloaded function.
optimize(self: bdsg.bdsg.HashGraph) -> None
optimize(self: bdsg.bdsg.HashGraph, allow_id_reassignment: bool) -> None
- Adjust the representation of the graph in memory to improve performance.
Optionally, allow the node IDs to be reassigned to further improve performance. Note: Ideally, this method is called one time once there is expected to be few graph modifications in the future.
C++: bdsg::HashGraph::optimize(bool) –> void
- path_back(self: bdsg.bdsg.HashGraph, path_handle: bdsg.handlegraph.path_handle_t) bdsg.handlegraph.step_handle_t ¶
- Get a handle to the last step, which will be an arbitrary step in a circular path that
we consider “last” based on our construction of the path. If the path is empty then the implementation must return the same value as path_front_end().
C++: bdsg::HashGraph::path_back(const struct handlegraph::path_handle_t &) const –> struct handlegraph::step_handle_t
- path_begin(self: bdsg.bdsg.HashGraph, path_handle: bdsg.handlegraph.path_handle_t) bdsg.handlegraph.step_handle_t ¶
- Get a handle to the first step, or in a circular path to an arbitrary step
considered “first”. If the path is empty, returns the past-the-last step returned by path_end.
C++: bdsg::HashGraph::path_begin(const struct handlegraph::path_handle_t &) const –> struct handlegraph::step_handle_t
- path_end(self: bdsg.bdsg.HashGraph, path_handle: bdsg.handlegraph.path_handle_t) bdsg.handlegraph.step_handle_t ¶
- Get a handle to a fictitious position past the end of a path. This position is
return by get_next_step for the final step in a path in a non-circular path. Note that get_next_step will NEVER return this value for a circular path.
C++: bdsg::HashGraph::path_end(const struct handlegraph::path_handle_t &) const –> struct handlegraph::step_handle_t
- path_front_end(self: bdsg.bdsg.HashGraph, path_handle: bdsg.handlegraph.path_handle_t) bdsg.handlegraph.step_handle_t ¶
- Get a handle to a fictitious position before the beginning of a path. This position is
return by get_previous_step for the first step in a path in a non-circular path. Note: get_previous_step will NEVER return this value for a circular path.
C++: bdsg::HashGraph::path_front_end(const struct handlegraph::path_handle_t &) const –> struct handlegraph::step_handle_t
- prepend_step(self: bdsg.bdsg.HashGraph, path: bdsg.handlegraph.path_handle_t, to_prepend: bdsg.handlegraph.handle_t) bdsg.handlegraph.step_handle_t ¶
- Prepend a visit to a node to the given path. Returns a handle to the new
first step on the path which is appended. If the path is cirular, the new step is placed between the steps considered “last” and “first” by the method path_begin. Handles to later steps on the path, and to other paths, must remain valid.
C++: bdsg::HashGraph::prepend_step(const struct handlegraph::path_handle_t &, const struct handlegraph::handle_t &) –> struct handlegraph::step_handle_t
- reassign_node_ids(self: bdsg.bdsg.HashGraph, get_new_id: Callable[[int], int]) None ¶
Reassign all node IDs as specified by the old->new mapping function.
C++: bdsg::HashGraph::reassign_node_ids(const class std::function<long long (const long long &)> &) –> void
- rewrite_segment(self: bdsg.bdsg.HashGraph, segment_begin: bdsg.handlegraph.step_handle_t, segment_end: bdsg.handlegraph.step_handle_t, new_segment: bdsg.std.vector_handlegraph_handle_t) Tuple[bdsg.handlegraph.step_handle_t, bdsg.handlegraph.step_handle_t] ¶
- Delete a segment of a path and rewrite it as some other sequence of
steps. Returns a pair of step_handle_t’s that indicate the range of the new segment in the path. The segment to delete should be designated by the first (begin) and past-last (end) step handles. If the step that is returned by path_begin is deleted, path_begin will now return the first step from the new segment or, in the case that the new segment is empty, the step used as segment_end. Empty ranges consist of two copies of the same step handle. Empty ranges in empty paths consist of two copies of the end sentinel handle for the path. Rewriting an empty range inserts before the provided end handle.
C++: bdsg::HashGraph::rewrite_segment(const struct handlegraph::step_handle_t &, const struct handlegraph::step_handle_t &, const class std::vector<handlegraph::handle_t> &) –> struct std::pair<struct handlegraph::step_handle_t, struct handlegraph::step_handle_t>
- set_circularity(self: bdsg.bdsg.HashGraph, path: bdsg.handlegraph.path_handle_t, circular: bool) None ¶
- Make a path circular or non-circular. If the path is becoming circular, the
last step is joined to the first step. If the path is becoming linear, the step considered “last” is unjoined from the step considered “first” according to the method path_begin.
C++: bdsg::HashGraph::set_circularity(const struct handlegraph::path_handle_t &, bool) –> void
- set_id_increment(self: bdsg.bdsg.HashGraph, min_id: int) None ¶
- Set a minimum id to increment the id space by, used as a hint during construction.
May have no effect on a backing implementation.
C++: bdsg::HashGraph::set_id_increment(const long long &) –> void
- truncate_handle(self: bdsg.bdsg.HashGraph, handle: bdsg.handlegraph.handle_t, trunc_left: bool, offset: int) bdsg.handlegraph.handle_t ¶
- Shorten a node by truncating either the left or right side of the node, relative to the orientation
of the handle, starting from a given offset along the nodes sequence. Any edges on the truncated end of the node are deleted. Returns a (possibly altered) handle to the truncated node. May invalid stored paths.
C++: bdsg::HashGraph::truncate_handle(const struct handlegraph::handle_t &, bool, unsigned long) –> struct handlegraph::handle_t
SnarlDecomposition Implementations¶
There is one implementation for a snarl decomposition
SnarlDistanceIndex¶
- class bdsg.bdsg.SnarlDistanceIndex¶
Bases:
SnarlDecomposition
,TriviallySerializable
The distance index, which also acts as a snarl decomposition.
The distance index provides an interface to traverse the snarl tree and to find minimum distances between two sibling nodes in the snarl tree (eg between two chains that are children of the same snarl).
It also provides a method for quickly calculating the minimum distance between two positions on the graph.
The implementation here is tightly coupled with the filling-in code in vg (see vg::fill_in_distance_index()). To make a SnarlDistanceIndex that actually works, you have to construct the object, and then call get_snarl_tree_records() with zero or more TemporaryDistanceIndex objects for connected components, and a graph.
The TemporaryDistanceIndex needs to have a variety of TemporaryRecord implementation classes (TemporaryChainRecord, TemporarySnarlRecord, TemporaryNodeRecord) set up and added to it; this all has to be done “by hand”, as it were, because no code is in this library to help you do it.
- static bit_width(value: int) int ¶
C++: bdsg::SnarlDistanceIndex::bit_width(unsigned long) –> unsigned long
- canonical(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) bdsg.handlegraph.net_handle_t ¶
Get a canonical traversal handle from any net handle. All handles to the same net graph element have the same canonical traversal. That canonical traversal must be realizable, and might not always be start-to-end or even consistently be the same kind of traversal for different snarls, chains, or nodes. Mostly useful to normalize for equality comparisons.
Any root snarl will become just the root Anything without connectivity will get START_END
C++: bdsg::SnarlDistanceIndex::canonical(const struct handlegraph::net_handle_t &) const –> struct handlegraph::net_handle_t
- chain_minimum_length(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) int ¶
The length of a chain. If it is a multicomponent chain, then the length of the last component, which is used for calculating distance, instead of inf
C++: bdsg::SnarlDistanceIndex::chain_minimum_length(const struct handlegraph::net_handle_t &) const –> unsigned long
- connected_component_count(self: bdsg.bdsg.SnarlDistanceIndex) int ¶
How many connected components are in this graph? This returns the number of topological connected components, not necessarily the number of nodes in the top-level snarl
C++: bdsg::SnarlDistanceIndex::connected_component_count() const –> unsigned long
- class connectivity_t¶
Bases:
pybind11_object
The connectivity of a net_handle- this defines the direction that the net_handle is traversed
Members:
START_START
START_END
START_TIP
END_START
END_END
END_TIP
TIP_START
TIP_END
TIP_TIP
- property name¶
- static connectivity_to_endpoints(connectivity: bdsg.bdsg.SnarlDistanceIndex.connectivity_t) Tuple[bdsg.handlegraph.SnarlDecomposition.endpoint_t, bdsg.handlegraph.SnarlDecomposition.endpoint_t] ¶
C++: bdsg::SnarlDistanceIndex::connectivity_to_endpoints(const enum bdsg::SnarlDistanceIndex::connectivity_t &) –> const struct std::pair<enum handlegraph::SnarlDecomposition::endpoint_t, enum handlegraph::SnarlDecomposition::endpoint_t>
- deserialize(*args, **kwargs)¶
Overloaded function.
deserialize(self: bdsg.bdsg.SnarlDistanceIndex, filename: str) -> None
deserialize(self: bdsg.bdsg.SnarlDistanceIndex, fd: int) -> None
C++: bdsg::SnarlDistanceIndex::deserialize(int) –> void
- dissociate(self: bdsg.bdsg.SnarlDistanceIndex) None ¶
C++: bdsg::SnarlDistanceIndex::dissociate() –> void
- distance_in_parent(*args, **kwargs)¶
Overloaded function.
distance_in_parent(self: bdsg.bdsg.SnarlDistanceIndex, parent: bdsg.handlegraph.net_handle_t, child1: bdsg.handlegraph.net_handle_t, child2: bdsg.handlegraph.net_handle_t) -> int
distance_in_parent(self: bdsg.bdsg.SnarlDistanceIndex, parent: bdsg.handlegraph.net_handle_t, child1: bdsg.handlegraph.net_handle_t, child2: bdsg.handlegraph.net_handle_t, graph: bdsg.handlegraph.HandleGraph) -> int
distance_in_parent(self: bdsg.bdsg.SnarlDistanceIndex, parent: bdsg.handlegraph.net_handle_t, child1: bdsg.handlegraph.net_handle_t, child2: bdsg.handlegraph.net_handle_t, graph: bdsg.handlegraph.HandleGraph, distance_limit: int) -> int
C++: bdsg::SnarlDistanceIndex::distance_in_parent(const struct handlegraph::net_handle_t &, const struct handlegraph::net_handle_t &, const struct handlegraph::net_handle_t &, const class handlegraph::HandleGraph *, unsigned long) const –> unsigned long
- distance_to_parent_bound(*args, **kwargs)¶
Overloaded function.
distance_to_parent_bound(self: bdsg.bdsg.SnarlDistanceIndex, parent: bdsg.handlegraph.net_handle_t, to_start: bool, child: bdsg.handlegraph.net_handle_t) -> int
distance_to_parent_bound(self: bdsg.bdsg.SnarlDistanceIndex, parent: bdsg.handlegraph.net_handle_t, to_start: bool, child: bdsg.handlegraph.net_handle_t, parent_and_child_types: Tuple[bdsg.bdsg.SnarlDistanceIndex.net_handle_record_t, bdsg.bdsg.SnarlDistanceIndex.net_handle_record_t, bdsg.bdsg.SnarlDistanceIndex.net_handle_record_t, bdsg.bdsg.SnarlDistanceIndex.net_handle_record_t]) -> int
Get the distance from the child to the start or end bound of the parent. parent_and_child_types are hints to figure out the type of snarl/chain records the parent and child are. tuple of parent record type, parent handle type, child record type, child handle type. This is really just used to see if the parent and child are trivial chains, so it might not be exactly what the actual record is.
C++: bdsg::SnarlDistanceIndex::distance_to_parent_bound(const struct handlegraph::net_handle_t &, bool, struct handlegraph::net_handle_t, class std::tuple<enum bdsg::SnarlDistanceIndex::net_handle_record_t, enum bdsg::SnarlDistanceIndex::net_handle_record_t, enum bdsg::SnarlDistanceIndex::net_handle_record_t, enum bdsg::SnarlDistanceIndex::net_handle_record_t>) const –> unsigned long
- static endpoints_to_connectivity(start: bdsg.handlegraph.SnarlDecomposition.endpoint_t, end: bdsg.handlegraph.SnarlDecomposition.endpoint_t) bdsg.bdsg.SnarlDistanceIndex.connectivity_t ¶
C++: bdsg::SnarlDistanceIndex::endpoints_to_connectivity(enum handlegraph::SnarlDecomposition::endpoint_t, enum handlegraph::SnarlDecomposition::endpoint_t) –> const enum bdsg::SnarlDistanceIndex::connectivity_t
- ends_at(self: bdsg.bdsg.SnarlDistanceIndex, traversal: bdsg.handlegraph.net_handle_t) bdsg.handlegraph.SnarlDecomposition.endpoint_t ¶
Return the kind of location at which the given traversal ends.
C++: bdsg::SnarlDistanceIndex::ends_at(const struct handlegraph::net_handle_t &) const –> enum handlegraph::SnarlDecomposition::endpoint_t
- flip(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) bdsg.handlegraph.net_handle_t ¶
Return a net handle to the same snarl/chain/node in the opposite orientation. No effect on tip-to-tip, start-to-start, or end-to-end net handles. Flips all the others.
C++: bdsg::SnarlDistanceIndex::flip(const struct handlegraph::net_handle_t &) const –> struct handlegraph::net_handle_t
- get_bound(self: bdsg.bdsg.SnarlDistanceIndex, snarl: bdsg.handlegraph.net_handle_t, get_end: bool, face_in: bool) bdsg.handlegraph.net_handle_t ¶
Get the bounding handle for the snarl or chain referenced by the given net handle, getting the start or end facing in or out as appropriate.
For snarls, returns the bounding sentinel net handles. For chains, returns net handles for traversals of the bounding nodes of the chain. If the chain is a looping chain, then the start and end of the chain are the same, so the connectivity of the bound indicates which we’re looking at; the connectivity will be start-start if it is going backwards in the node, and end-end if it is going forwards.
Ignores traversal type.
May not be called on traversals of individual nodes.
C++: bdsg::SnarlDistanceIndex::get_bound(const struct handlegraph::net_handle_t &, bool, bool) const –> struct handlegraph::net_handle_t
- get_chain_component(*args, **kwargs)¶
Overloaded function.
get_chain_component(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) -> int
get_chain_component(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t, get_end: bool) -> int
C++: bdsg::SnarlDistanceIndex::get_chain_component(const struct handlegraph::net_handle_t, bool) const –> unsigned long
- get_connected_component_number(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) int ¶
C++: bdsg::SnarlDistanceIndex::get_connected_component_number(const struct handlegraph::net_handle_t &) const –> unsigned long
- static get_connectivity(net_handle: bdsg.handlegraph.net_handle_t) bdsg.bdsg.SnarlDistanceIndex.connectivity_t ¶
C++: bdsg::SnarlDistanceIndex::get_connectivity(const struct handlegraph::net_handle_t &) –> const enum bdsg::SnarlDistanceIndex::connectivity_t
- get_depth(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) int ¶
What is the depth of this net handle? Nodes and snarls get the depth of their parent. The depth of the root is 0, the depth of its child chains is 1, the depth of the nodes and snarls that are children of those chains is also 1, and the chains that are children of those snarls have depth 2
C++: bdsg::SnarlDistanceIndex::get_depth(const struct handlegraph::net_handle_t &) const –> unsigned long
- static get_end_endpoint(*args, **kwargs)¶
Overloaded function.
get_end_endpoint(connectivity: bdsg.bdsg.SnarlDistanceIndex.connectivity_t) -> bdsg.handlegraph.SnarlDecomposition.endpoint_t
C++: bdsg::SnarlDistanceIndex::get_end_endpoint(enum bdsg::SnarlDistanceIndex::connectivity_t) –> const enum handlegraph::SnarlDecomposition::endpoint_t
get_end_endpoint(net: bdsg.handlegraph.net_handle_t) -> bdsg.handlegraph.SnarlDecomposition.endpoint_t
C++: bdsg::SnarlDistanceIndex::get_end_endpoint(const struct handlegraph::net_handle_t &) –> const enum handlegraph::SnarlDecomposition::endpoint_t
- get_forward_loop_value(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) int ¶
Get the prefix sum value for a node in a chain. Fails if the parent of net is not a chain
C++: bdsg::SnarlDistanceIndex::get_forward_loop_value(const struct handlegraph::net_handle_t) const –> unsigned long
- get_handle(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t, graph: bdsg.handlegraph.HandleGraph) bdsg.handlegraph.handle_t ¶
For a net handle to a traversal of a single node, get the handle for that node in the orientation it is traversed. May not be called for other net handles.
C++: bdsg::SnarlDistanceIndex::get_handle(const struct handlegraph::net_handle_t &, const class handlegraph::HandleGraph *) const –> struct handlegraph::handle_t
- get_handle_from_connected_component(self: bdsg.bdsg.SnarlDistanceIndex, num: int) bdsg.handlegraph.net_handle_t ¶
Given the connected component number (from get_connected_component_number), get the root-level handle pointing to it. If the connected component is a root-level snarl, then this may return a “root” handle, but it will actually point to the snarl
C++: bdsg::SnarlDistanceIndex::get_handle_from_connected_component(unsigned long) const –> struct handlegraph::net_handle_t
- static get_handle_type(net_handle: bdsg.handlegraph.net_handle_t) bdsg.bdsg.SnarlDistanceIndex.net_handle_record_t ¶
C++: bdsg::SnarlDistanceIndex::get_handle_type(const struct handlegraph::net_handle_t &) –> const enum bdsg::SnarlDistanceIndex::net_handle_record_t
- get_magic_number(self: bdsg.bdsg.SnarlDistanceIndex) int ¶
C++: bdsg::SnarlDistanceIndex::get_magic_number() const –> unsigned int
- get_max_tree_depth(self: bdsg.bdsg.SnarlDistanceIndex) int ¶
How deep is the snarl tree? The root is 0, top-level chain is 1, etc Only counts chains
C++: bdsg::SnarlDistanceIndex::get_max_tree_depth() const –> unsigned long
- get_net(self: bdsg.bdsg.SnarlDistanceIndex, handle: bdsg.handlegraph.handle_t, graph: bdsg.handlegraph.HandleGraph) bdsg.handlegraph.net_handle_t ¶
Turn a handle to an oriented node into a net handle for a start-to-end or end-to-start traversal of the node, as appropriate.
C++: bdsg::SnarlDistanceIndex::get_net(const struct handlegraph::handle_t &, const class handlegraph::HandleGraph *) const –> struct handlegraph::net_handle_t
- get_net_handle(*args, **kwargs)¶
Overloaded function.
get_net_handle(self: bdsg.bdsg.SnarlDistanceIndex, pointer: int, connectivity: bdsg.bdsg.SnarlDistanceIndex.connectivity_t) -> bdsg.handlegraph.net_handle_t
C++: bdsg::SnarlDistanceIndex::get_net_handle(unsigned long, enum bdsg::SnarlDistanceIndex::connectivity_t) const –> struct handlegraph::net_handle_t
get_net_handle(self: bdsg.bdsg.SnarlDistanceIndex, pointer: int) -> bdsg.handlegraph.net_handle_t
C++: bdsg::SnarlDistanceIndex::get_net_handle(unsigned long) const –> struct handlegraph::net_handle_t
- static get_net_handle_from_values(*args, **kwargs)¶
Overloaded function.
get_net_handle_from_values(pointer: int, connectivity: bdsg.bdsg.SnarlDistanceIndex.connectivity_t, type: bdsg.bdsg.SnarlDistanceIndex.net_handle_record_t) -> bdsg.handlegraph.net_handle_t
get_net_handle_from_values(pointer: int, connectivity: bdsg.bdsg.SnarlDistanceIndex.connectivity_t, type: bdsg.bdsg.SnarlDistanceIndex.net_handle_record_t, node_offset: int) -> bdsg.handlegraph.net_handle_t
C++: bdsg::SnarlDistanceIndex::get_net_handle_from_values(unsigned long, enum bdsg::SnarlDistanceIndex::connectivity_t, enum bdsg::SnarlDistanceIndex::net_handle_record_t, unsigned long) –> const struct handlegraph::net_handle_t
- get_node_from_sentinel(self: bdsg.bdsg.SnarlDistanceIndex, sentinel: bdsg.handlegraph.net_handle_t) bdsg.handlegraph.net_handle_t ¶
Given the sentinel of a snarl, return a handle to the node representing it
C++: bdsg::SnarlDistanceIndex::get_node_from_sentinel(const struct handlegraph::net_handle_t &) const –> struct handlegraph::net_handle_t
- get_node_net_handle(*args, **kwargs)¶
Overloaded function.
get_node_net_handle(self: bdsg.bdsg.SnarlDistanceIndex, id: int) -> bdsg.handlegraph.net_handle_t
get_node_net_handle(self: bdsg.bdsg.SnarlDistanceIndex, id: int, rev: bool) -> bdsg.handlegraph.net_handle_t
Get a net handle from a node and, optionally, an orientation
C++: bdsg::SnarlDistanceIndex::get_node_net_handle(const long long, bool) const –> struct handlegraph::net_handle_t
- static get_node_pointer_offset(id: int, min_node_id: int, component_count: int) int ¶
Get the offset into snarl_tree_records for the pointer to a node record.
C++: bdsg::SnarlDistanceIndex::get_node_pointer_offset(const long long &, const long long &, unsigned long) –> const unsigned long
- static get_node_record_offset(net_handle: bdsg.handlegraph.net_handle_t) int ¶
The offset of a node in a trivial snarl (0 if it isn’t a node in a trivial snarl)
C++: bdsg::SnarlDistanceIndex::get_node_record_offset(const struct handlegraph::net_handle_t &) –> const unsigned long
- get_parent(self: bdsg.bdsg.SnarlDistanceIndex, child: bdsg.handlegraph.net_handle_t) bdsg.handlegraph.net_handle_t ¶
Get the parent snarl of a chain, or the parent chain of a snarl or node. If the child is start-to-end or end-to-start, and the parent is a chain, the chain comes out facing the same way, accounting for the relative orientation of the child snarl or node in the chain. Otherwise, everything is produced as start-to-end, even if that is not actually a realizable traversal of a snarl or chain. May not be called on the root snarl.
Also works on snarl boundary sentinels.
C++: bdsg::SnarlDistanceIndex::get_parent(const struct handlegraph::net_handle_t &) const –> struct handlegraph::net_handle_t
- get_parent_traversal(self: bdsg.bdsg.SnarlDistanceIndex, traversal_start: bdsg.handlegraph.net_handle_t, traversal_end: bdsg.handlegraph.net_handle_t) bdsg.handlegraph.net_handle_t ¶
Get a net handle for traversals of a snarl or chain that contains the given oriented bounding node traversals or sentinels. Given two sentinels for a snarl, produces a net handle to a start-to-end, end-to-end, end-to-start, or start-to-start traversal of that snarl. Given handles to traversals of the bounding nodes of a chain, similarly produces a net handle to a traversal of the chain.
For a chain, either or both handles can also be a snarl containing tips, for a tip-to-start, tip-to-end, start-to-tip, end-to-tip, or tip-to-tip traversal. Similarly, for a snarl, either or both handles can be a chain in the snarl that contains internal tips, or that has no edges on the appropriate end.
May only be called if a path actually exists between the given start and end.
C++: bdsg::SnarlDistanceIndex::get_parent_traversal(const struct handlegraph::net_handle_t &, const struct handlegraph::net_handle_t &) const –> struct handlegraph::net_handle_t
- get_prefix(self: bdsg.bdsg.SnarlDistanceIndex) str ¶
C++: bdsg::SnarlDistanceIndex::get_prefix() const –> std::string
- get_prefix_sum_value(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) int ¶
Get the prefix sum value for a node in a chain. Fails if the parent of net is not a chain
C++: bdsg::SnarlDistanceIndex::get_prefix_sum_value(const struct handlegraph::net_handle_t) const –> unsigned long
- get_rank_in_parent(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) int ¶
For a child of a snarl, the rank is used to calculate the distance
C++: bdsg::SnarlDistanceIndex::get_rank_in_parent(const struct handlegraph::net_handle_t &) const –> unsigned long
- static get_record_handle_type(type: bdsg.bdsg.SnarlDistanceIndex.record_t) bdsg.bdsg.SnarlDistanceIndex.net_handle_record_t ¶
Given the type of the record, return the handle type. Some record types can represent multiple things, for example a simple snarl record is used to represent a snarl, and the nodes/trivial chains in it. This will return whatever is higher on the snarl tree. A simple snarl will be considered a snarl, a root snarl will be considered a root, etc
C++: bdsg::SnarlDistanceIndex::get_record_handle_type(enum bdsg::SnarlDistanceIndex::record_t) –> const enum bdsg::SnarlDistanceIndex::net_handle_record_t
- static get_record_offset(net_handle: bdsg.handlegraph.net_handle_t) int ¶
The offset into records that this handle points to
C++: bdsg::SnarlDistanceIndex::get_record_offset(const struct handlegraph::net_handle_t &) –> const unsigned long
- get_reverse_loop_value(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) int ¶
Get the prefix sum value for a node in a chain. Fails if the parent of net is not a chain
C++: bdsg::SnarlDistanceIndex::get_reverse_loop_value(const struct handlegraph::net_handle_t) const –> unsigned long
- get_root(self: bdsg.bdsg.SnarlDistanceIndex) bdsg.handlegraph.net_handle_t ¶
Get a net handle referring to a tip-to-tip traversal of the contents of the root snarl.
C++: bdsg::SnarlDistanceIndex::get_root() const –> struct handlegraph::net_handle_t
- static get_start_endpoint(*args, **kwargs)¶
Overloaded function.
get_start_endpoint(connectivity: bdsg.bdsg.SnarlDistanceIndex.connectivity_t) -> bdsg.handlegraph.SnarlDecomposition.endpoint_t
C++: bdsg::SnarlDistanceIndex::get_start_endpoint(enum bdsg::SnarlDistanceIndex::connectivity_t) –> const enum handlegraph::SnarlDecomposition::endpoint_t
get_start_endpoint(net: bdsg.handlegraph.net_handle_t) -> bdsg.handlegraph.SnarlDecomposition.endpoint_t
C++: bdsg::SnarlDistanceIndex::get_start_endpoint(struct handlegraph::net_handle_t) –> const enum handlegraph::SnarlDecomposition::endpoint_t
- get_usage(self: bdsg.bdsg.SnarlDistanceIndex) Tuple[int, int, int] ¶
C++: bdsg::SnarlDistanceIndex::get_usage() –> class std::tuple<unsigned long, unsigned long, unsigned long>
- has_connectivity(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t, start: bdsg.handlegraph.SnarlDecomposition.endpoint_t, end: bdsg.handlegraph.SnarlDecomposition.endpoint_t) bool ¶
Is there a path between the start and end endpoints within the net handle?
C++: bdsg::SnarlDistanceIndex::has_connectivity(const struct handlegraph::net_handle_t &, enum handlegraph::SnarlDecomposition::endpoint_t, enum handlegraph::SnarlDecomposition::endpoint_t) const –> bool
- has_external_connectivity(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t, start: bdsg.handlegraph.SnarlDecomposition.endpoint_t, end: bdsg.handlegraph.SnarlDecomposition.endpoint_t) bool ¶
Is there a path between the start and end endpoints outside the net handle? This is used for children of the root
C++: bdsg::SnarlDistanceIndex::has_external_connectivity(const struct handlegraph::net_handle_t &, enum handlegraph::SnarlDecomposition::endpoint_t, enum handlegraph::SnarlDecomposition::endpoint_t) const –> bool
- has_node(self: bdsg.bdsg.SnarlDistanceIndex, id: int) bool ¶
Does the graph have this node?
C++: bdsg::SnarlDistanceIndex::has_node(const long long) const –> bool
- into_which_snarl(self: bdsg.bdsg.SnarlDistanceIndex, id: int, reverse: bool) Tuple[int, bool, bool] ¶
If this node id and orientation is pointing into a snarl, then return the start. node id and orientation pointing into the snarl, and if the snarl is trivial. Returns <0, false, false> if it doesn’t point into a snarl.
C++: bdsg::SnarlDistanceIndex::into_which_snarl(const long long &, const bool &) const –> class std::tuple<long long, bool, bool>
- is_chain(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) bool ¶
Returns true if the given net handle refers to (a traversal of) a chain.
C++: bdsg::SnarlDistanceIndex::is_chain(const struct handlegraph::net_handle_t &) const –> bool
- is_externally_end_end_connected(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) bool ¶
C++: bdsg::SnarlDistanceIndex::is_externally_end_end_connected(const struct handlegraph::net_handle_t) const –> bool
- is_externally_start_end_connected(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) bool ¶
C++: bdsg::SnarlDistanceIndex::is_externally_start_end_connected(const struct handlegraph::net_handle_t) const –> bool
- is_externally_start_start_connected(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) bool ¶
C++: bdsg::SnarlDistanceIndex::is_externally_start_start_connected(const struct handlegraph::net_handle_t) const –> bool
- is_looping_chain(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) bool ¶
Returns true if the given net handle refers to (a traversal of) a chain that loops (a chain where the first and last node are the same).
C++: bdsg::SnarlDistanceIndex::is_looping_chain(const struct handlegraph::net_handle_t &) const –> bool
- is_multicomponent_chain(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) bool ¶
Returns true if the given net handle refers to (a traversal of) a chain that is not start-end connected
C++: bdsg::SnarlDistanceIndex::is_multicomponent_chain(const struct handlegraph::net_handle_t &) const –> bool
- is_node(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) bool ¶
Returns true if the given net handle refers to (a traversal of) a single node, and thus has a corresponding handle_t.
C++: bdsg::SnarlDistanceIndex::is_node(const struct handlegraph::net_handle_t &) const –> bool
- is_ordered_in_chain(self: bdsg.bdsg.SnarlDistanceIndex, child1: bdsg.handlegraph.net_handle_t, child2: bdsg.handlegraph.net_handle_t) bool ¶
Return true if child1 comes before child2 in the chain.
C++: bdsg::SnarlDistanceIndex::is_ordered_in_chain(const struct handlegraph::net_handle_t &, const struct handlegraph::net_handle_t &) const –> bool
- is_reversed_in_parent(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) bool ¶
Only really relevant for nodes in chains, is the node traversed backwards relative to the orientation of the chain
C++: bdsg::SnarlDistanceIndex::is_reversed_in_parent(const struct handlegraph::net_handle_t &) const –> bool
- is_root(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) bool ¶
Return true if the given handle refers to (a traversal of) the root snarl, and false otherwise.
C++: bdsg::SnarlDistanceIndex::is_root(const struct handlegraph::net_handle_t &) const –> bool
- is_root_snarl(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) bool ¶
Return true if the given handle refers to (a traversal of) a snarl of the root, which is considered to be the root but actually refers to a subset of the children of the root that are connected
C++: bdsg::SnarlDistanceIndex::is_root_snarl(const struct handlegraph::net_handle_t &) const –> bool
- is_sentinel(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) bool ¶
Return true if the given net handle is a snarl bound sentinel (in either inward or outward orientation), and false otherwise.
C++: bdsg::SnarlDistanceIndex::is_sentinel(const struct handlegraph::net_handle_t &) const –> bool
- is_simple_snarl(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) bool ¶
Returns true if the given net handle refers to (a traversal of) a simple snarl A simple snarl is a bubble where each child node can only reach the boundary nodes, and each side of a node reaches a different boundary node
C++: bdsg::SnarlDistanceIndex::is_simple_snarl(const struct handlegraph::net_handle_t &) const –> bool
- is_snarl(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) bool ¶
Returns true if the given net handle refers to (a traversal of) a snarl.
C++: bdsg::SnarlDistanceIndex::is_snarl(const struct handlegraph::net_handle_t &) const –> bool
- is_trivial_chain(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) bool ¶
Returns true if the given net handle refers to (a traversal of) a trivial chain that represents a single node.
C++: bdsg::SnarlDistanceIndex::is_trivial_chain(const struct handlegraph::net_handle_t &) const –> bool
- lowest_common_ancestor(self: bdsg.bdsg.SnarlDistanceIndex, net1: bdsg.handlegraph.net_handle_t, net2: bdsg.handlegraph.net_handle_t) Tuple[bdsg.handlegraph.net_handle_t, bool] ¶
For two net handles, get a net handle lowest common ancestor. If the lowest common ancestor is the root, then the two handles may be in different connected components. In this case, return false.
C++: bdsg::SnarlDistanceIndex::lowest_common_ancestor(const struct handlegraph::net_handle_t &, const struct handlegraph::net_handle_t &) const –> struct std::pair<struct handlegraph::net_handle_t, bool>
- max_distance_in_parent(*args, **kwargs)¶
Overloaded function.
max_distance_in_parent(self: bdsg.bdsg.SnarlDistanceIndex, parent: bdsg.handlegraph.net_handle_t, child1: bdsg.handlegraph.net_handle_t, child2: bdsg.handlegraph.net_handle_t) -> int
max_distance_in_parent(self: bdsg.bdsg.SnarlDistanceIndex, parent: bdsg.handlegraph.net_handle_t, child1: bdsg.handlegraph.net_handle_t, child2: bdsg.handlegraph.net_handle_t, graph: bdsg.handlegraph.HandleGraph) -> int
max_distance_in_parent(self: bdsg.bdsg.SnarlDistanceIndex, parent: bdsg.handlegraph.net_handle_t, child1: bdsg.handlegraph.net_handle_t, child2: bdsg.handlegraph.net_handle_t, graph: bdsg.handlegraph.HandleGraph, distance_limit: int) -> int
Find the maximum distance between two children in the parent. This is the same as distance_in_parent for everything except children of chains
C++: bdsg::SnarlDistanceIndex::max_distance_in_parent(const struct handlegraph::net_handle_t &, const struct handlegraph::net_handle_t &, const struct handlegraph::net_handle_t &, const class handlegraph::HandleGraph *, unsigned long) const –> unsigned long
- static maximum(x: int, y: int) int ¶
C++: bdsg::SnarlDistanceIndex::maximum(unsigned long, unsigned long) –> unsigned long
- maximum_distance(*args, **kwargs)¶
Overloaded function.
maximum_distance(self: bdsg.bdsg.SnarlDistanceIndex, id1: int, rev1: bool, offset1: int, id2: int, rev2: bool, offset2: int) -> int
maximum_distance(self: bdsg.bdsg.SnarlDistanceIndex, id1: int, rev1: bool, offset1: int, id2: int, rev2: bool, offset2: int, unoriented_distance: bool) -> int
maximum_distance(self: bdsg.bdsg.SnarlDistanceIndex, id1: int, rev1: bool, offset1: int, id2: int, rev2: bool, offset2: int, unoriented_distance: bool, graph: bdsg.handlegraph.HandleGraph) -> int
Find an approximation of the maximum distance between two positions. This isn’t a true maximum- the only guarantee is that it’s greater than or equal to the minimum distance.
C++: bdsg::SnarlDistanceIndex::maximum_distance(const long long, const bool, const unsigned long, const long long, const bool, const unsigned long, bool, const class handlegraph::HandleGraph *) const –> unsigned long
- maximum_length(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) int ¶
C++: bdsg::SnarlDistanceIndex::maximum_length(const struct handlegraph::net_handle_t &) const –> unsigned long
- minimum_length(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) int ¶
This is also the length of a net, but it can also be a snarl or chain. The length of a chain includes the boundary nodes, a snarl does not. A looping chain only includes the start/end node once
C++: bdsg::SnarlDistanceIndex::minimum_length(const struct handlegraph::net_handle_t &) const –> unsigned long
- static minus(x: int, y: int) int ¶
C++: bdsg::SnarlDistanceIndex::minus(unsigned long, unsigned long) –> unsigned long
- net_handle_as_string(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) str ¶
C++: bdsg::SnarlDistanceIndex::net_handle_as_string(const struct handlegraph::net_handle_t &) const –> std::string
- class net_handle_record_t¶
Bases:
pybind11_object
Type of a net_handle_t, which may not be the type of the record This is to allow a node record to be seen as a chain from the perspective of a handle
Members:
ROOT_HANDLE
NODE_HANDLE
SNARL_HANDLE
CHAIN_HANDLE
SENTINEL_HANDLE
- property name¶
- node_id(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) int ¶
What is the node id of the node represented by this net handle. net must be a node or a sentinel
C++: bdsg::SnarlDistanceIndex::node_id(const struct handlegraph::net_handle_t &) const –> long long
- node_length(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) int ¶
Return the length of the net, which must represent a node (or sentinel of a snarl)
C++: bdsg::SnarlDistanceIndex::node_length(const struct handlegraph::net_handle_t &) const –> unsigned long
- preload(*args, **kwargs)¶
Overloaded function.
preload(self: bdsg.bdsg.SnarlDistanceIndex) -> None
preload(self: bdsg.bdsg.SnarlDistanceIndex, blocking: bool) -> None
- Allow for preloading the index for more accurate timing of algorithms
that use it, if it fits in memory. If blocking is true, waits for the index to be paged in. Otherwise, just tells the OS that we will want to use it.
C++: bdsg::SnarlDistanceIndex::preload(bool) const –> void
- print_descendants_of(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) None ¶
C++: bdsg::SnarlDistanceIndex::print_descendants_of(const struct handlegraph::net_handle_t) const –> void
- print_self(self: bdsg.bdsg.SnarlDistanceIndex) None ¶
C++: bdsg::SnarlDistanceIndex::print_self() const –> void
- print_snarl_stats(self: bdsg.bdsg.SnarlDistanceIndex) None ¶
C++: bdsg::SnarlDistanceIndex::print_snarl_stats() const –> void
- class record_t¶
Bases:
pybind11_object
A record_t is the type of structure that a record can be.
NODE, SNARL, and CHAIN indicate that they don’t store distances. SIMPLE_SNARL is a snarl with all children connecting only to the boundary nodes in one direction. OVERSIZED_SNARL only stores distances to the boundaries. ROOT_SNARL represents a connected component of the root. It has no start or end node so
its children technically belong to the root.
- MULTICOMPONENT_CHAIN can represent a chain with snarls that are not start-end connected.
The chain is split up into components between these snarls, each node is tagged with which component it belongs to.
Members:
ROOT
NODE
DISTANCED_NODE
TRIVIAL_SNARL
DISTANCED_TRIVIAL_SNARL
SIMPLE_SNARL
DISTANCED_SIMPLE_SNARL
SNARL
DISTANCED_SNARL
OVERSIZED_SNARL
ROOT_SNARL
DISTANCED_ROOT_SNARL
CHAIN
DISTANCED_CHAIN
MULTICOMPONENT_CHAIN
CHILDREN
- property name¶
- serialize(*args, **kwargs)¶
Overloaded function.
serialize(self: bdsg.bdsg.SnarlDistanceIndex, filename: str) -> None
serialize(self: bdsg.bdsg.SnarlDistanceIndex, iteratee: Callable[[capsule, int], None]) -> None
C++: bdsg::SnarlDistanceIndex::serialize(const class std::function<void (const void *, unsigned long)> &) const –> void
serialize(self: bdsg.bdsg.SnarlDistanceIndex, fd: int) -> None
C++: bdsg::SnarlDistanceIndex::serialize(int) –> void
- set_snarl_size_limit(self: bdsg.bdsg.SnarlDistanceIndex, size: int) None ¶
C++: bdsg::SnarlDistanceIndex::set_snarl_size_limit(unsigned long) –> void
- start_end_traversal_of(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) bdsg.handlegraph.net_handle_t ¶
Makes a start-end traversal of the net. Faster than canonical because it doesn’t check the index for anything
C++: bdsg::SnarlDistanceIndex::start_end_traversal_of(const struct handlegraph::net_handle_t &) const –> struct handlegraph::net_handle_t
- starts_at(self: bdsg.bdsg.SnarlDistanceIndex, traversal: bdsg.handlegraph.net_handle_t) bdsg.handlegraph.SnarlDecomposition.endpoint_t ¶
Return the kind of location at which the given traversal starts.
C++: bdsg::SnarlDistanceIndex::starts_at(const struct handlegraph::net_handle_t &) const –> enum handlegraph::SnarlDecomposition::endpoint_t
- static sum(val1: int, val2: int) int ¶
Add integers, returning max() if any of them are max()
C++: bdsg::SnarlDistanceIndex::sum(const unsigned long &, const unsigned long &) –> unsigned long
- class temp_record_t¶
Bases:
pybind11_object
Members:
TEMP_CHAIN
TEMP_SNARL
TEMP_NODE
TEMP_ROOT
- property name¶
- time_accesses(self: bdsg.bdsg.SnarlDistanceIndex) None ¶
C++: bdsg::SnarlDistanceIndex::time_accesses() –> void
- traverse_decomposition(self: bdsg.bdsg.SnarlDistanceIndex, snarl_iteratee: Callable[[bdsg.handlegraph.net_handle_t], bool], chain_iteratee: Callable[[bdsg.handlegraph.net_handle_t], bool], node_iteratee: Callable[[bdsg.handlegraph.net_handle_t], bool]) bool ¶
C++: bdsg::SnarlDistanceIndex::traverse_decomposition(const class std::function<bool (const struct handlegraph::net_handle_t &)> &, const class std::function<bool (const struct handlegraph::net_handle_t &)> &, const class std::function<bool (const struct handlegraph::net_handle_t &)> &) const –> bool
- traverse_decomposition_helper(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t, snarl_iteratee: Callable[[bdsg.handlegraph.net_handle_t], bool], chain_iteratee: Callable[[bdsg.handlegraph.net_handle_t], bool], node_iteratee: Callable[[bdsg.handlegraph.net_handle_t], bool]) bool ¶
C++: bdsg::SnarlDistanceIndex::traverse_decomposition_helper(const struct handlegraph::net_handle_t &, const class std::function<bool (const struct handlegraph::net_handle_t &)> &, const class std::function<bool (const struct handlegraph::net_handle_t &)> &, const class std::function<bool (const struct handlegraph::net_handle_t &)> &) const –> bool
- validate_ancestors_of(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) None ¶
C++: bdsg::SnarlDistanceIndex::validate_ancestors_of(const struct handlegraph::net_handle_t) const –> void
- validate_descendants_of(self: bdsg.bdsg.SnarlDistanceIndex, net: bdsg.handlegraph.net_handle_t) None ¶
C++: bdsg::SnarlDistanceIndex::validate_descendants_of(const struct handlegraph::net_handle_t) const –> void
- validate_index(self: bdsg.bdsg.SnarlDistanceIndex) None ¶
C++: bdsg::SnarlDistanceIndex::validate_index() const –> void
- write_snarls_to_json(self: bdsg.bdsg.SnarlDistanceIndex) None ¶
C++: bdsg::SnarlDistanceIndex::write_snarls_to_json() const –> void
Graph Overlays¶
In addition to these basic implementations, there are several “overlays”. These overlays are graphs that wrap other graphs, providing features not avialable in the backing graph, or otherwise transforming it.
- class bdsg.bdsg.PositionOverlay¶
- class bdsg.bdsg.PackedPositionOverlay¶
- class bdsg.bdsg.MutablePositionOverlay¶
- class bdsg.bdsg.VectorizableOverlay¶
- class bdsg.bdsg.PathVectorizableOverlay¶
Bases:
VectorizableOverlay
,PathHandleGraph
- class bdsg.bdsg.PathPositionVectorizableOverlay¶
Many of these are based on the bdsg.handlegraph.ExpandingOverlayGraph
interface, which guarantees that the overlay does not remove any graph material, and allows handles form the backing graph and the overlay graph to be interconverted.
- class bdsg.handlegraph.ExpandingOverlayGraph¶
Bases:
HandleGraph
This is the interface for a graph that represents a transformation of some underlying HandleGraph where every node in the overlay corresponds to a node in the underlying graph, but where more than one node in the overlay can map to the same underlying node.
- assign(self: bdsg.handlegraph.ExpandingOverlayGraph, : bdsg.handlegraph.ExpandingOverlayGraph) bdsg.handlegraph.ExpandingOverlayGraph ¶
C++: handlegraph::ExpandingOverlayGraph::operator=(const class handlegraph::ExpandingOverlayGraph &) –> class handlegraph::ExpandingOverlayGraph &
- get_underlying_handle(self: bdsg.handlegraph.ExpandingOverlayGraph, handle: bdsg.handlegraph.handle_t) bdsg.handlegraph.handle_t ¶
- Returns the handle in the underlying graph that corresponds to a handle in the
overlay
C++: handlegraph::ExpandingOverlayGraph::get_underlying_handle(const struct handlegraph::handle_t &) const –> struct handlegraph::handle_t
Typed Collections¶
Some methods, such as bdsg.handlegraph.MutableHandleGraph.divide_handle()
, produce or consume collections of typed objects: C++ STL vectors with Python bindings. The typed collection classes are available in bdsg.std
. They are convertible from and to Python lists via their constructors and the list constructor, respectively.
Here is an example of how to use these typed collections:
import bdsg
g = bdsg.bdsg.HashGraph()
h = g.create_handle("GATTACA")
v = bdsg.std.vector_unsigned_long([1, 3])
parts = g.divide_handle(h, v)
# parts is a bdsg.std.vector_handlegraph_handle_t
print(list(parts))
Bindings for ::std namespace
- class bdsg.std.vector_handlegraph_handle_t¶
- append(self: bdsg.std.vector_handlegraph_handle_t, x: bdsg.handlegraph.handle_t) None ¶
Add an item to the end of the list
- capacity(self: bdsg.std.vector_handlegraph_handle_t) int ¶
returns the number of elements that can be held in currently allocated storage
- clear(*args, **kwargs)¶
Overloaded function.
clear(self: bdsg.std.vector_handlegraph_handle_t) -> None
Clear the contents
clear(self: bdsg.std.vector_handlegraph_handle_t) -> None
clears the contents
- count(self: bdsg.std.vector_handlegraph_handle_t, x: bdsg.handlegraph.handle_t) int ¶
Return the number of times
x
appears in the list
- empty(self: bdsg.std.vector_handlegraph_handle_t) bool ¶
checks whether the container is empty
- extend(*args, **kwargs)¶
Overloaded function.
extend(self: bdsg.std.vector_handlegraph_handle_t, L: bdsg.std.vector_handlegraph_handle_t) -> None
Extend the list by appending all the items in the given list
extend(self: bdsg.std.vector_handlegraph_handle_t, L: Iterable) -> None
Extend the list by appending all the items in the given list
- insert(self: bdsg.std.vector_handlegraph_handle_t, i: int, x: bdsg.handlegraph.handle_t) None ¶
Insert an item at a given position.
- max_size(self: bdsg.std.vector_handlegraph_handle_t) int ¶
returns the maximum possible number of elements
- pop(*args, **kwargs)¶
Overloaded function.
pop(self: bdsg.std.vector_handlegraph_handle_t) -> bdsg.handlegraph.handle_t
Remove and return the last item
pop(self: bdsg.std.vector_handlegraph_handle_t, i: int) -> bdsg.handlegraph.handle_t
Remove and return the item at index
i
- remove(self: bdsg.std.vector_handlegraph_handle_t, x: bdsg.handlegraph.handle_t) None ¶
Remove the first item from the list whose value is x. It is an error if there is no such item.
- reserve(self: bdsg.std.vector_handlegraph_handle_t, arg0: int) None ¶
reserves storage
- shrink_to_fit(self: bdsg.std.vector_handlegraph_handle_t) None ¶
reduces memory usage by freeing unused memory
- swap(self: bdsg.std.vector_handlegraph_handle_t, arg0: bdsg.std.vector_handlegraph_handle_t) None ¶
swaps the contents
- class bdsg.std.vector_handlegraph_step_handle_t¶
- append(self: bdsg.std.vector_handlegraph_step_handle_t, x: bdsg.handlegraph.step_handle_t) None ¶
Add an item to the end of the list
- capacity(self: bdsg.std.vector_handlegraph_step_handle_t) int ¶
returns the number of elements that can be held in currently allocated storage
- clear(*args, **kwargs)¶
Overloaded function.
clear(self: bdsg.std.vector_handlegraph_step_handle_t) -> None
Clear the contents
clear(self: bdsg.std.vector_handlegraph_step_handle_t) -> None
clears the contents
- count(self: bdsg.std.vector_handlegraph_step_handle_t, x: bdsg.handlegraph.step_handle_t) int ¶
Return the number of times
x
appears in the list
- empty(self: bdsg.std.vector_handlegraph_step_handle_t) bool ¶
checks whether the container is empty
- extend(*args, **kwargs)¶
Overloaded function.
extend(self: bdsg.std.vector_handlegraph_step_handle_t, L: bdsg.std.vector_handlegraph_step_handle_t) -> None
Extend the list by appending all the items in the given list
extend(self: bdsg.std.vector_handlegraph_step_handle_t, L: Iterable) -> None
Extend the list by appending all the items in the given list
- insert(self: bdsg.std.vector_handlegraph_step_handle_t, i: int, x: bdsg.handlegraph.step_handle_t) None ¶
Insert an item at a given position.
- max_size(self: bdsg.std.vector_handlegraph_step_handle_t) int ¶
returns the maximum possible number of elements
- pop(*args, **kwargs)¶
Overloaded function.
pop(self: bdsg.std.vector_handlegraph_step_handle_t) -> bdsg.handlegraph.step_handle_t
Remove and return the last item
pop(self: bdsg.std.vector_handlegraph_step_handle_t, i: int) -> bdsg.handlegraph.step_handle_t
Remove and return the item at index
i
- remove(self: bdsg.std.vector_handlegraph_step_handle_t, x: bdsg.handlegraph.step_handle_t) None ¶
Remove the first item from the list whose value is x. It is an error if there is no such item.
- reserve(self: bdsg.std.vector_handlegraph_step_handle_t, arg0: int) None ¶
reserves storage
- shrink_to_fit(self: bdsg.std.vector_handlegraph_step_handle_t) None ¶
reduces memory usage by freeing unused memory
- swap(self: bdsg.std.vector_handlegraph_step_handle_t, arg0: bdsg.std.vector_handlegraph_step_handle_t) None ¶
swaps the contents
- class bdsg.std.vector_unsigned_long¶
- append(self: bdsg.std.vector_unsigned_long, x: int) None ¶
Add an item to the end of the list
- capacity(self: bdsg.std.vector_unsigned_long) int ¶
returns the number of elements that can be held in currently allocated storage
- clear(*args, **kwargs)¶
Overloaded function.
clear(self: bdsg.std.vector_unsigned_long) -> None
Clear the contents
clear(self: bdsg.std.vector_unsigned_long) -> None
clears the contents
- count(self: bdsg.std.vector_unsigned_long, x: int) int ¶
Return the number of times
x
appears in the list
- empty(self: bdsg.std.vector_unsigned_long) bool ¶
checks whether the container is empty
- extend(*args, **kwargs)¶
Overloaded function.
extend(self: bdsg.std.vector_unsigned_long, L: bdsg.std.vector_unsigned_long) -> None
Extend the list by appending all the items in the given list
extend(self: bdsg.std.vector_unsigned_long, L: Iterable) -> None
Extend the list by appending all the items in the given list
- insert(self: bdsg.std.vector_unsigned_long, i: int, x: int) None ¶
Insert an item at a given position.
- max_size(self: bdsg.std.vector_unsigned_long) int ¶
returns the maximum possible number of elements
- pop(*args, **kwargs)¶
Overloaded function.
pop(self: bdsg.std.vector_unsigned_long) -> int
Remove and return the last item
pop(self: bdsg.std.vector_unsigned_long, i: int) -> int
Remove and return the item at index
i
- remove(self: bdsg.std.vector_unsigned_long, x: int) None ¶
Remove the first item from the list whose value is x. It is an error if there is no such item.
- reserve(self: bdsg.std.vector_unsigned_long, arg0: int) None ¶
reserves storage
- shrink_to_fit(self: bdsg.std.vector_unsigned_long) None ¶
reduces memory usage by freeing unused memory
- swap(self: bdsg.std.vector_unsigned_long, arg0: bdsg.std.vector_unsigned_long) None ¶
swaps the contents