Sorted Glossary of Methods¶
Mutator Methods¶
The following lists sorts methods in :class:odgi.graph by what types of objects they modify.
odgi.edge
¶
-
odgi.graph.
create_edge
(*args, **kwargs)¶ Overloaded function.
create_edge(self: odgi.graph, arg0: odgi.handle, arg1: odgi.handle) -> None
Create an edge connecting the given handles in the given order and orientations.
create_edge(self: odgi.graph, arg0: Tuple[odgi.handle, odgi.handle]) -> None
Create an edge connecting the given handles in the given order and orientations.
-
odgi.graph.
destroy_edge
(*args, **kwargs)¶ Overloaded function.
destroy_edge(self: odgi.graph, arg0: Tuple[odgi.handle, odgi.handle]) -> None
Remove the edge connecting the given handles in the given order and orientations.
destroy_edge(self: odgi.graph, arg0: odgi.handle, arg1: odgi.handle) -> None
Remove the edge connecting the given handles in the given order and orientations.
odgi.graph
¶
-
odgi.graph.
apply_ordering
(self: odgi.graph, order: List[odgi.handle], compact_ids: bool = False) → None¶ Reorder the graph’s internal structure to match that given. Optionally compact the id space of the graph to match the ordering, from 1->|ordering|.
-
odgi.graph.
apply_path_ordering
(self: odgi.graph, arg0: List[odgi.path_handle]) → None¶ Reorder the graph’s paths as given.
-
odgi.graph.
clear
(self: odgi.graph) → None¶ Remove all nodes and edges. Does not update any stored paths.
-
odgi.graph.
clear_paths
(self: odgi.graph) → None¶ Remove all stored paths.
-
odgi.graph.
load
(self: odgi.graph, arg0: str) → None¶ Load the graph from the given file.
-
odgi.graph.
optimize
(self: odgi.graph, allow_id_reassignment: bool = False) → None¶ Organize the graph for better performance and memory use.
-
odgi.graph.
serialize
(self: odgi.graph, arg0: str) → None¶ Save the graph to the given file, returning the number of bytes written.
odgi.handle
¶
-
odgi.graph.
apply_orientation
(self: odgi.graph, arg0: odgi.handle) → odgi.handle¶ Alter the node that the given handle corresponds to so the orientation indicated by the handle becomes the node’s local forward orientation. Updates all links and path steps to match the new orientation.
-
odgi.graph.
combine_handles
(self: odgi.graph, arg0: List[odgi.handle]) → odgi.handle¶ Join handles into a new node, returning the handle of the new node.
-
odgi.graph.
create_handle
(*args, **kwargs)¶ Overloaded function.
create_handle(self: odgi.graph, arg0: str) -> odgi.handle
Create a new node with the given sequence and return the handle.
create_handle(self: odgi.graph, arg0: str, arg1: int) -> odgi.handle
Create a new node with the given sequence and return the handle.
-
odgi.graph.
destroy_handle
(self: odgi.graph, arg0: odgi.handle) → None¶ Remove the node belonging to the given handle and all of its edges. Does not update any stored paths. Invalidates the destroyed handle.
-
odgi.graph.
divide_handle
(*args, **kwargs)¶ Overloaded function.
divide_handle(self: odgi.graph, arg0: odgi.handle, arg1: List[int]) -> List[odgi.handle]
Split a handle’s underlying node at the given offsets in the handle’s orientation. Returns the handles to the new parts.
divide_handle(self: odgi.graph, arg0: odgi.handle, arg1: int) -> Tuple[odgi.handle, odgi.handle]
Split a handle’s underlying node at the given offset in the handle’s orientation. Returns the handles to the new parts.
-
odgi.graph.
flip
(self: odgi.graph, handle: odgi.handle) → odgi.handle¶ Flip the handle to the opposite orientation.
odgi.path_handle
¶
-
odgi.graph.
append_step
(self: odgi.graph, arg0: odgi.path_handle, arg1: odgi.handle) → odgi.step_handle¶ Append a visit to a node to the given path. Returns a handle to the new final step on the path which is appended.
-
odgi.graph.
create_path_handle
(self: odgi.graph, name: str, is_circular: bool = False) → odgi.path_handle¶ Create a path with the given name. The caller must ensure that no path with the given name already exists.
-
odgi.graph.
destroy_path
(self: odgi.graph, arg0: odgi.path_handle) → None¶ Destroy the given path. Invalidates handles to the path and its node steps.
-
odgi.graph.
prepend_step
(self: odgi.graph, arg0: odgi.path_handle, arg1: odgi.handle) → odgi.step_handle¶ Append a visit to a node to the given path. Returns a handle to the new final step on the path which is appended.
-
odgi.graph.
set_circularity
(self: odgi.graph, arg0: odgi.path_handle, arg1: bool) → None¶ Set if the path is circular or not.
odgi.step_handle
¶
-
odgi.graph.
rewrite_segment
(self: odgi.graph, arg0: odgi.step_handle, arg1: odgi.step_handle, arg2: List[odgi.handle]) → Tuple[odgi.step_handle, odgi.step_handle]¶ Replace the path range with the new segment, returning the new start and end step handles for the segment.
-
odgi.graph.
set_step
(self: odgi.graph, arg0: odgi.step_handle, arg1: odgi.handle) → odgi.step_handle¶ Set the step to the given handle, possibly re-linking and cleaning up if needed.
Accessor Methods¶
The following list sorts methods in :class:odgi.graph by what object they return information about.
odgi.edge
¶
-
odgi.graph.
edge_handle
(self: odgi.graph, arg0: odgi.handle, arg1: odgi.handle) → Tuple[odgi.handle, odgi.handle]¶ Return the edge handle for the given pair of handles.
odgi.graph
¶
-
odgi.graph.
get_node_count
(self: odgi.graph) → int¶ Return the number of nodes in the graph.
-
odgi.graph.
get_path_count
(self: odgi.graph) → int¶ Return the path count of the graph
-
odgi.graph.
has_node
(self: odgi.graph, node_id: int) → bool¶ Return true if the given node is in the graph.
-
odgi.graph.
has_path
(self: odgi.graph, arg0: str) → bool¶ Return if a path with the givenv name exists in the graph.
-
odgi.graph.
max_node_id
(self: odgi.graph) → int¶ Return the maximum node id in the graph.
-
odgi.graph.
min_node_id
(self: odgi.graph) → int¶ Return the minimum node id in the graph.
-
odgi.graph.
to_gfa
(self: odgi.graph) → None¶ Display as GFA
odgi.handle
¶
-
odgi.graph.
forward
(self: odgi.graph, arg0: odgi.handle) → odgi.handle¶ Return the forward version of the handle.
-
odgi.graph.
get_degree
(self: odgi.graph, arg0: odgi.handle, arg1: bool) → int¶ Return the degree of the given node.
-
odgi.graph.
get_handle
(self: odgi.graph, node_id: int, is_reverse: bool = False) → odgi.handle¶ Return the handle for the given node id.
-
odgi.graph.
get_id
(self: odgi.graph, handle: odgi.handle) → int¶ Return the id of the given handle.
-
odgi.graph.
get_is_reverse
(self: odgi.graph, handle: odgi.handle) → bool¶ Return true if the handle refers to the node reverse complement.
-
odgi.graph.
get_length
(self: odgi.graph, handle: odgi.handle) → int¶ Return the length of the node referred to by the handle.
-
odgi.graph.
get_sequence
(self: odgi.graph, handle: odgi.handle) → str¶
-
odgi.graph.
get_step_count
(*args, **kwargs)¶ Overloaded function.
get_step_count(self: odgi.graph, arg0: odgi.path_handle) -> int
Return the step count of a given path.
get_step_count(self: odgi.graph, arg0: odgi.handle) -> int
Return the number of steps on the given handle.
-
odgi.graph.
has_edge
(self: odgi.graph, arg0: odgi.handle, arg1: odgi.handle) → bool¶ Returns true if the given edge exists
-
odgi.graph.
steps_of_handle
(self: odgi.graph, arg0: odgi.handle, arg1: bool) → List[odgi.step_handle]¶ Obtain the steps on a given handle.
odgi.path_handle
¶
-
odgi.graph.
get_handle_of_step
(self: odgi.graph, arg0: odgi.step_handle) → odgi.handle¶ Return the handle that a given step occurs on.
-
odgi.graph.
get_is_circular
(self: odgi.graph, arg0: odgi.path_handle) → bool¶ Returns true if the path is circular.
-
odgi.graph.
get_path_handle
(self: odgi.graph, arg0: str) → odgi.path_handle¶ Return the path handle for the named path.
-
odgi.graph.
get_path_name
(self: odgi.graph, arg0: odgi.path_handle) → str¶ Return the path name for a given path handle.
-
odgi.graph.
get_step_count
(*args, **kwargs)¶ Overloaded function.
get_step_count(self: odgi.graph, arg0: odgi.path_handle) -> int
Return the step count of a given path.
get_step_count(self: odgi.graph, arg0: odgi.handle) -> int
Return the number of steps on the given handle.
-
odgi.graph.
is_empty
(self: odgi.graph, arg0: odgi.path_handle) → bool¶ Returns true if the given path is empty, and false otherwise.
-
odgi.graph.
path_back
(self: odgi.graph, arg0: odgi.path_handle) → odgi.step_handle¶ Return a step handle to the last step, which is arbitrary in the case of a circular path.
-
odgi.graph.
path_begin
(self: odgi.graph, arg0: odgi.path_handle) → odgi.step_handle¶ Return the step handle for the first step in the given path.
-
odgi.graph.
path_end
(self: odgi.graph, arg0: odgi.path_handle) → odgi.step_handle¶ Return a step handle to a fictitious handle one past the end of the path.
-
odgi.graph.
path_front_end
(self: odgi.graph, arg0: odgi.path_handle) → odgi.step_handle¶ Return a step handle to a fictitious handle one past the start of the path.
odgi.step_handle
¶
-
odgi.graph.
get_next_step
(self: odgi.graph, arg0: odgi.step_handle) → odgi.step_handle¶ Returns a handle to the next step on the path. Calling on an end marker step returns the same end marker.
-
odgi.graph.
get_path
(self: odgi.graph, arg0: odgi.step_handle) → odgi.path_handle¶ Return the path of a given step handle.
-
odgi.graph.
get_path_handle_of_step
(self: odgi.graph, arg0: odgi.step_handle) → odgi.path_handle¶ Returns a handle to the path that an step is on.
-
odgi.graph.
get_previous_step
(self: odgi.graph, arg0: odgi.step_handle) → odgi.step_handle¶ Returns a handle to the previous step on the path. Calling on a front end marker step returns the same end marker.
-
odgi.graph.
has_next_step
(self: odgi.graph, arg0: odgi.step_handle) → bool¶ Returns true if the step is not the last step on the path, else false.
-
odgi.graph.
has_previous_step
(self: odgi.graph, arg0: odgi.step_handle) → bool¶ Returns true if the step is not the first step on the path, else false.
-
odgi.graph.
is_path_end
(self: odgi.graph, arg0: odgi.step_handle) → bool¶ Returns true if the step handle is an end magic handle.
-
odgi.graph.
is_path_front_end
(self: odgi.graph, arg0: odgi.step_handle) → bool¶ Returns true if the step handle is a front end magic handle.
Iteratator Methods¶
The following list sorts methods in :class:odgi.graph by what kind of iteratee they operate on.
odgi.edge
¶
-
odgi.graph.
follow_edges
(self: odgi.graph, arg0: odgi.handle, arg1: bool, arg2: Callable[[odgi.handle], bool]) → bool¶ Follow edges starting at a given node.
odgi.handle
¶
-
odgi.graph.
for_each_handle
(self: odgi.graph, iteratee: Callable[[odgi.handle], bool], parallel: bool = False) → bool¶ Iterate over all the nodes in the graph.
odgi.path_handle
¶
-
odgi.graph.
for_each_path_handle
(self: odgi.graph, arg0: Callable[[odgi.path_handle], bool]) → bool¶ Invoke the callback for each path in the graph.
odgi.step_handle
¶
-
odgi.graph.
for_each_step_in_path
(self: odgi.graph, arg0: odgi.path_handle, arg1: Callable[[odgi.step_handle], None]) → None¶ Invoke the callback for each step in a given path.
-
odgi.graph.
for_each_step_on_handle
(self: odgi.graph, arg0: odgi.handle, arg1: Callable[[odgi.step_handle], bool]) → bool¶ Invoke the callback for each of the steps on a given handle.