Fast compiled graphs

This implements the base class for sparse and dense graphs in Sage. It is not intended for use on its own.

Data structure

The class CGraph contains the following variables:

cdef int num_verts
cdef int num_arcs
cdef int *in_degrees
cdef int *out_degrees
cdef bitset_t active_vertices

The bitset active_vertices is a list of all available vertices for use, but only the ones which are set are considered to actually be in the graph. The variables num_verts and num_arcs are self-explanatory (note that num_verts is the number of bits set in active_vertices, not the full length of the bitset). The arrays in_degrees and out_degrees are of the same length as the bitset.

For more about active vertices, see the documentation for the realloc method.

Classes and methods

class sage.graphs.base.c_graph.CGraph

Bases: object

Compiled sparse and dense graphs.

add_arc(u, v)

This function is implemented at the level of sparse and dense graphs:

sage: from sage.graphs.base.c_graph import CGraph
sage: G = CGraph()
sage: G.add_arc(0,1)
Traceback (most recent call last):
...
NotImplementedError
add_vertex(k=-1)

Adds vertex k to the graph. If k == -1, a new vertex is added and the integer used is returned.

INPUT:

  • k – non-negative integer, or -1

EXAMPLE:

sage: from sage.graphs.base.sparse_graph import SparseGraph
sage: G = SparseGraph(3, extra_vertices=3)
sage: G.add_vertex(3)
3
sage: G.add_arc(2,5)
Traceback (most recent call last):
...
RuntimeError: Vertex (5) is not a vertex of the graph.
sage: G.add_arc(1,3)
sage: G.has_arc(1,3)
True
sage: G.has_arc(2,3)
False
sage: from sage.graphs.base.dense_graph import DenseGraph
sage: G = DenseGraph(3, extra_vertices=3)
sage: G.add_vertex(3)
3
sage: G.add_arc(2,5)
Traceback (most recent call last):
...
RuntimeError: Vertex (5) is not a vertex of the graph.
sage: G.add_arc(1,3)
sage: G.has_arc(1,3)
True
sage: G.has_arc(2,3)
False
add_vertices(verts)

Adds vertices from the iterable verts.

EXAMPLE:

sage: from sage.graphs.base.sparse_graph import SparseGraph
sage: S = SparseGraph(nverts=4, extra_vertices=4)
sage: S.verts()
[0, 1, 2, 3]
sage: S.add_vertices([3,5,7,9])
sage: S.verts()
[0, 1, 2, 3, 5, 7, 9]
sage: S.realloc(20)
sage: S.verts()
[0, 1, 2, 3, 5, 7, 9]
sage: from sage.graphs.base.dense_graph import DenseGraph
sage: D = DenseGraph(nverts=4, extra_vertices=4)
sage: D.verts()
[0, 1, 2, 3]
sage: D.add_vertices([3,5,7,9])
sage: D.verts()
[0, 1, 2, 3, 5, 7, 9]
sage: D.realloc(20)
sage: D.verts()
[0, 1, 2, 3, 5, 7, 9]
all_arcs(u, v)

This function is implemented at the level of sparse and dense graphs:

sage: from sage.graphs.base.c_graph import CGraph
sage: G = CGraph()
sage: G.all_arcs(0,1)
Traceback (most recent call last):
...
NotImplementedError 
check_vertex(n)

If n is not in self, raise an error.

EXAMPLES:

sage: from sage.graphs.base.sparse_graph import SparseGraph
sage: S = SparseGraph(nverts = 10, expected_degree = 3, extra_vertices = 10)        
sage: S.check_vertex(4)
sage: S.check_vertex(12)
Traceback (most recent call last):
...
RuntimeError: Vertex (12) is not a vertex of the graph.
sage: S.check_vertex(24)
Traceback (most recent call last):
...
RuntimeError: Vertex (24) is not a vertex of the graph.
sage: S.check_vertex(-19)
Traceback (most recent call last):
...
RuntimeError: Vertex (-19) is not a vertex of the graph.
sage: from sage.graphs.base.dense_graph import DenseGraph
sage: D = DenseGraph(nverts = 10, extra_vertices = 10)
sage: D.check_vertex(4)
sage: D.check_vertex(12)
Traceback (most recent call last):
...
RuntimeError: Vertex (12) is not a vertex of the graph.
sage: D.check_vertex(24)
Traceback (most recent call last):
...
RuntimeError: Vertex (24) is not a vertex of the graph.
sage: D.check_vertex(-19)
Traceback (most recent call last):
...
RuntimeError: Vertex (-19) is not a vertex of the graph.
current_allocation()

Report the number of vertices allocated.

INPUT:

  • total - integer, the total size to make the array

Returns -1 and fails if reallocation would destroy any active vertices.

EXAMPLES:

sage: from sage.graphs.base.sparse_graph import SparseGraph
sage: S = SparseGraph(nverts=4, extra_vertices=4)
sage: S.current_allocation()
8
sage: S.add_vertex(6)
6
sage: S.current_allocation()
8
sage: S.add_vertex(10)
10
sage: S.current_allocation()
16
sage: S.add_vertex(40)
Traceback (most recent call last):
...
RuntimeError: Requested vertex is past twice the allocated range: use realloc.
sage: S.realloc(50)
sage: S.add_vertex(40)
40
sage: S.current_allocation()
50
sage: S.realloc(30)
-1
sage: S.current_allocation()
50
sage: S.del_vertex(40)
sage: S.realloc(30)
sage: S.current_allocation()
30
del_all_arcs(u, v)

This function is implemented at the level of sparse and dense graphs:

sage: from sage.graphs.base.c_graph import CGraph
sage: G = CGraph()
sage: G.del_all_arcs(0,1)
Traceback (most recent call last):
...
NotImplementedError 
del_vertex(v)

Deletes the vertex v, along with all edges incident to it. If v is not in self, fails silently.

EXAMPLES:

sage: from sage.graphs.base.sparse_graph import SparseGraph
sage: G = SparseGraph(3)
sage: G.add_arc(0,1)
sage: G.add_arc(0,2)
sage: G.add_arc(1,2)
sage: G.add_arc(2,0)
sage: G.del_vertex(2)
sage: for i in range(2):
...    for j in range(2):
...        if G.has_arc(i,j):
...            print i,j
0 1
sage: G = SparseGraph(3)
sage: G.add_arc(0,1)
sage: G.add_arc(0,2)
sage: G.add_arc(1,2)
sage: G.add_arc(2,0)
sage: G.del_vertex(1)
sage: for i in xrange(3):
...    for j in xrange(3):
...        if G.has_arc(i,j):
...            print i,j
0 2
2 0
has_arc(u, v)

This function is implemented at the level of sparse and dense graphs:

sage: from sage.graphs.base.c_graph import CGraph
sage: G = CGraph()
sage: G.has_arc(0,1)
Not Implemented!
False
has_vertex(n)

Return whether n is in self.

INPUT:
  • n - integer

EXAMPLES:

Upon initialization, a SparseGraph or DenseGraph has the first nverts vertices:

sage: from sage.graphs.base.sparse_graph import SparseGraph
sage: S = SparseGraph(nverts = 10, expected_degree = 3, extra_vertices = 10)
sage: S.has_vertex(6)
True
sage: S.has_vertex(12)
False
sage: S.has_vertex(24)
False
sage: S.has_vertex(-19)
False
sage: from sage.graphs.base.dense_graph import DenseGraph
sage: D = DenseGraph(nverts = 10, extra_vertices = 10)
sage: D.has_vertex(6)
True
sage: D.has_vertex(12)
False
sage: D.has_vertex(24)
False
sage: D.has_vertex(-19)
False
in_neighbors(v)

This function is implemented at the level of sparse and dense graphs:

sage: from sage.graphs.base.c_graph import CGraph
sage: G = CGraph()
sage: G.in_neighbors(0)
Traceback (most recent call last):
...
NotImplementedError 
out_neighbors(u)

This function is implemented at the level of sparse and dense graphs:

sage: from sage.graphs.base.c_graph import CGraph
sage: G = CGraph()
sage: G.out_neighbors(0)
Traceback (most recent call last):
...
NotImplementedError 
realloc(total)

Reallocate the number of vertices to use, without actually adding any.

INPUT:

  • total - integer, the total size to make the array

Returns -1 and fails if reallocation would destroy any active vertices.

EXAMPLES:

First, note that realloc is implemented for SparseGraph and DenseGraph differently, and is not implemented at the CGraph level:

sage: from sage.graphs.base.c_graph import CGraph
sage: G = CGraph()
sage: G.realloc(20)
Traceback (most recent call last):
...
NotImplementedError
sage: from sage.graphs.base.sparse_graph import SparseGraph
sage: S = SparseGraph(nverts=4, extra_vertices=4)
sage: S.current_allocation()
8
sage: S.add_vertex(6)
6
sage: S.current_allocation()
8
sage: S.add_vertex(10)
10
sage: S.current_allocation()
16
sage: S.add_vertex(40)
Traceback (most recent call last):
...
RuntimeError: Requested vertex is past twice the allocated range: use realloc.
sage: S.realloc(50)
sage: S.add_vertex(40)
40
sage: S.current_allocation()
50
sage: S.realloc(30)
-1
sage: S.current_allocation()
50
sage: S.del_vertex(40)
sage: S.realloc(30)
sage: S.current_allocation()
30
sage: from sage.graphs.base.dense_graph import DenseGraph
sage: D = DenseGraph(nverts=4, extra_vertices=4)
sage: D.current_allocation()
8
sage: D.add_vertex(6)
6
sage: D.current_allocation()
8
sage: D.add_vertex(10)
10
sage: D.current_allocation()
16
sage: D.add_vertex(40)
Traceback (most recent call last):
...
RuntimeError: Requested vertex is past twice the allocated range: use realloc.
sage: D.realloc(50)
sage: D.add_vertex(40)
40
sage: D.current_allocation()
50
sage: D.realloc(30)
-1
sage: D.current_allocation()
50
sage: D.del_vertex(40)
sage: D.realloc(30)
sage: D.current_allocation()
30
verts()

Returns a list of the vertices in self.

EXAMPLE:

sage: from sage.graphs.base.sparse_graph import SparseGraph
sage: S = SparseGraph(nverts=4, extra_vertices=4)
sage: S.verts()
[0, 1, 2, 3]
sage: S.add_vertices([3,5,7,9])
sage: S.verts()
[0, 1, 2, 3, 5, 7, 9]
sage: S.realloc(20)
sage: S.verts()
[0, 1, 2, 3, 5, 7, 9]
class sage.graphs.base.c_graph.CGraphBackend

Bases: sage.graphs.base.graph_backends.GenericGraphBackend

Base class for sparse and dense graph backends.

sage: from sage.graphs.base.c_graph import CGraphBackend

This class is extended by SparseGraphBackend and DenseGraphBackend, which are fully functional backends. This class is mainly just for vertex functions, which are the same for both. A CGraphBackend will not work on its own:

sage: from sage.graphs.base.c_graph import CGraphBackend
sage: CGB = CGraphBackend()
sage: CGB.degree(0,True)
Traceback (most recent call last):
...
AttributeError: 'CGraphBackend' object has no attribute 'vertex_ints'

The appropriate way to use these backends is via Sage graphs:

sage: G = Graph(30, implementation="c_graph")
sage: G.add_edges([(0,1), (0,3), (4,5), (9, 23)])
sage: G.edges(labels=False)
[(0, 1), (0, 3), (4, 5), (9, 23)]
add_vertex(name)

Add a vertex to self.

INPUT:

  • name - the vertex to be added (must be hashable)

EXAMPLE:

sage: D = sage.graphs.base.dense_graph.DenseGraphBackend(9)
sage: D.add_vertex(10)
sage: D.add_vertex([])
Traceback (most recent call last):
...
TypeError: unhashable type: 'list'
sage: S = sage.graphs.base.sparse_graph.SparseGraphBackend(9)
sage: S.add_vertex(10)
sage: S.add_vertex([])
Traceback (most recent call last):
...
TypeError: unhashable type: 'list'
add_vertices(vertices)

Add vertices to self.

INPUT:

  • vertices - iterator of vertex labels

EXAMPLE:

sage: D = sage.graphs.base.sparse_graph.SparseGraphBackend(1)
sage: D.add_vertices([1,2,3])
sage: G = sage.graphs.base.sparse_graph.SparseGraphBackend(0)
sage: G.add_vertices([0,1])
sage: list(G.iterator_verts(None))
[0, 1]
sage: list(G.iterator_edges([0,1], True))
[]
sage: import sage.graphs.base.dense_graph
sage: D = sage.graphs.base.dense_graph.DenseGraphBackend(9)
sage: D.add_vertices([10,11,12])
bidirectional_dijkstra(x, y)

Returns the shortest path between x and y using a bidirectional version of Dijkstra

EXAMPLE:

sage: G = Graph(graphs.PetersenGraph(), implementation='c_graph')
sage: for (u,v) in G.edges(labels=None):
...      G.set_edge_label(u,v,1)
sage: G.shortest_path(0,1,by_weight=True)
[0, 1]

TEST:

Bugfix from #7673

sage: G = Graph(implementation="networkx")
sage: G.add_edges([(0,1,9),(0,2,8),(1,2,7)])
sage: Gc = G.copy(implementation='c_graph')
sage: sp = G.shortest_path_length(0,1,by_weight=True)
sage: spc = Gc.shortest_path_length(0,1,by_weight=True)
sage: sp == spc
True

Returns a breadth-first search from vertex v.

INPUT:

  • v – a vertex from which to start the breadth-first search.
  • reverse – boolean (default: False). This is only relevant to digraphs. If this is a digraph, consider the reversed graph in which the out-neighbors become the in-neighbors and vice versa.
  • ignore_direction – boolean (default: False). This is only relevant to digraphs. If this is a digraph, ignore all orientations and consider the graph as undirected.

ALGORITHM:

Below is a general template for breadth-first search.

  • Input: A directed or undirected graph G = (V, E) of order n > 0. A vertex s from which to start the search. The vertices are numbered from 1 to n = |V|, i.e. V = \{1, 2, \dots, n\}.
  • Output: A list D of distances of all vertices from s. A tree T rooted at s.
  1. Q \leftarrow [s] # a queue of nodes to visit
  2. D \leftarrow [\infty, \infty, \dots, \infty] # n copies of \infty
  3. D[s] \leftarrow 0
  4. T \leftarrow [\,]
  5. while \text{length}(Q) > 0 do
    1. v \leftarrow \text{dequeue}(Q)
    2. for each w \in \text{adj}(v) do # for digraphs, use out-neighbor set \text{oadj}(v)
      1. if D[w] = \infty then
        1. D[w] \leftarrow D[v] + 1
        2. \text{enqueue}(Q, w)
        3. \text{append}(T, vw)
  6. return (D, T)

See also

EXAMPLES:

Breadth-first search of the Petersen graph starting at vertex 0:

sage: G = Graph(graphs.PetersenGraph(), implementation="c_graph")
sage: list(G.breadth_first_search(0))
[0, 1, 4, 5, 2, 6, 3, 9, 7, 8]

Visiting German cities using breadth-first search:

sage: G = Graph({"Mannheim": ["Frankfurt","Karlsruhe"],
...   "Frankfurt": ["Mannheim","Wurzburg","Kassel"],
...   "Kassel": ["Frankfurt","Munchen"],
...   "Munchen": ["Kassel","Nurnberg","Augsburg"],
...   "Augsburg": ["Munchen","Karlsruhe"],
...   "Karlsruhe": ["Mannheim","Augsburg"],
...   "Wurzburg": ["Frankfurt","Erfurt","Nurnberg"],
...   "Nurnberg": ["Wurzburg","Stuttgart","Munchen"],
...   "Stuttgart": ["Nurnberg"],
...   "Erfurt": ["Wurzburg"]}, implementation="c_graph")
sage: list(G.breadth_first_search("Frankfurt"))
['Frankfurt', 'Mannheim', 'Kassel', 'Wurzburg', 'Karlsruhe', 'Munchen', 'Erfurt', 'Nurnberg', 'Augsburg', 'Stuttgart']
degree(v, directed)

Return the degree of the vertex v.

INPUT:

  • v - a vertex of the graph

EXAMPLE:

sage: from sage.graphs.base.sparse_graph import SparseGraphBackend
sage: B = SparseGraphBackend(7)
sage: B.degree(3, False)
0
del_vertex(v)

Delete a vertex in self, failing silently if the vertex is not in the graph.

INPUT:

  • v - vertex to be deleted

EXAMPLE:

sage: D = sage.graphs.base.dense_graph.DenseGraphBackend(9)
sage: D.del_vertex(0)
sage: D.has_vertex(0)
False
sage: S = sage.graphs.base.sparse_graph.SparseGraphBackend(9)
sage: S.del_vertex(0)
sage: S.has_vertex(0)
False
del_vertices(vertices)

Delete vertices from an iterable container.

INPUT:

  • vertices - iterator of vertex labels

EXAMPLE:

sage: import sage.graphs.base.dense_graph
sage: D = sage.graphs.base.dense_graph.DenseGraphBackend(9)
sage: D.del_vertices([7,8])
sage: D.has_vertex(7)
False
sage: D.has_vertex(6)
True
sage: D = sage.graphs.base.sparse_graph.SparseGraphBackend(9)
sage: D.del_vertices([1,2,3])
sage: D.has_vertex(1)
False
sage: D.has_vertex(0)
True

Returns a depth-first search from vertex v.

INPUT:

  • v – a vertex from which to start the depth-first search.
  • reverse – boolean (default: False). This is only relevant to digraphs. If this is a digraph, consider the reversed graph in which the out-neighbors become the in-neighbors and vice versa.
  • ignore_direction – boolean (default: False). This is only relevant to digraphs. If this is a digraph, ignore all orientations and consider the graph as undirected.

ALGORITHM:

Below is a general template for depth-first search.

  • Input: A directed or undirected graph G = (V, E) of order n > 0. A vertex s from which to start the search. The vertices are numbered from 1 to n = |V|, i.e. V = \{1, 2, \dots, n\}.
  • Output: A list D of distances of all vertices from s. A tree T rooted at s.
  1. S \leftarrow [s] # a stack of nodes to visit
  2. D \leftarrow [\infty, \infty, \dots, \infty] # n copies of \infty
  3. D[s] \leftarrow 0
  4. T \leftarrow [\,]
  5. while \text{length}(S) > 0 do
    1. v \leftarrow \text{pop}(S)
    2. for each w \in \text{adj}(v) do # for digraphs, use out-neighbor set \text{oadj}(v)
      1. if D[w] = \infty then
        1. D[w] \leftarrow D[v] + 1
        2. \text{push}(S, w)
        3. \text{append}(T, vw)
  6. return (D, T)

See also

EXAMPLES:

Traversing the Petersen graph using depth-first search:

sage: G = Graph(graphs.PetersenGraph(), implementation="c_graph")
sage: list(G.depth_first_search(0))
[0, 5, 8, 6, 9, 7, 2, 3, 4, 1]

Visiting German cities using depth-first search:

sage: G = Graph({"Mannheim": ["Frankfurt","Karlsruhe"],
...   "Frankfurt": ["Mannheim","Wurzburg","Kassel"],
...   "Kassel": ["Frankfurt","Munchen"],
...   "Munchen": ["Kassel","Nurnberg","Augsburg"],
...   "Augsburg": ["Munchen","Karlsruhe"],
...   "Karlsruhe": ["Mannheim","Augsburg"],
...   "Wurzburg": ["Frankfurt","Erfurt","Nurnberg"],
...   "Nurnberg": ["Wurzburg","Stuttgart","Munchen"],
...   "Stuttgart": ["Nurnberg"],
...   "Erfurt": ["Wurzburg"]}, implementation="c_graph")
sage: list(G.depth_first_search("Frankfurt"))
['Frankfurt', 'Wurzburg', 'Nurnberg', 'Munchen', 'Kassel', 'Augsburg', 'Karlsruhe', 'Mannheim', 'Stuttgart', 'Erfurt']
has_vertex(v)

Returns whether v is a vertex of self.

INPUT:

  • v - any object

EXAMPLE:

sage: from sage.graphs.base.sparse_graph import SparseGraphBackend
sage: B = SparseGraphBackend(7)
sage: B.has_vertex(6)
True
sage: B.has_vertex(7)
False
is_connected()

Returns whether the graph is connected.

EXAMPLE:

Petersen’s graph is connected

sage: DiGraph(graphs.PetersenGraph(),implementation="c_graph").is_connected()
True

While the disjoint union of two of them is not:

sage: DiGraph(2*graphs.PetersenGraph(),implementation="c_graph").is_connected()
False
A graph with non-integer vertex labels::
sage: Graph(graphs.CubeGraph(3), implementation=’c_graph’).is_connected() True
is_strongly_connected()

Returns whether the graph is strongly connected.

EXAMPLE:

The circuit on 3 vertices is obviously strongly connected

sage: g = DiGraph({ 0 : [1], 1 : [2], 2 : [0]},implementation="c_graph")
sage: g.is_strongly_connected()
True

But a transitive triangle is not:

sage: g = DiGraph({ 0 : [1,2], 1 : [2]},implementation="c_graph")
sage: g.is_strongly_connected()
False
iterator_in_nbrs(v)

Returns an iterator over the incoming neighbors of v.

INPUT:

  • v - a vertex

EXAMPLE:

sage: P = DiGraph(graphs.PetersenGraph().to_directed(), implementation='c_graph')
sage: list(P._backend.iterator_in_nbrs(0))
[1, 4, 5]
iterator_nbrs(v)

Returns an iterator over the neighbors of v.

INPUT:

  • v - a vertex

EXAMPLE:

sage: P = Graph(graphs.PetersenGraph(), implementation='c_graph')
sage: list(P._backend.iterator_nbrs(0))
[1, 4, 5]
iterator_out_nbrs(v)

Returns an iterator over the outgoing neighbors of v.

INPUT:

  • v - a vertex

EXAMPLE:

sage: P = DiGraph(graphs.PetersenGraph().to_directed(), implementation='c_graph')
sage: list(P._backend.iterator_out_nbrs(0))
[1, 4, 5]
iterator_verts(verts)

Returns an iterator over the vertices of self intersected with verts.

INPUT:
  • verts - an iterable container of objects

EXAMPLE:

sage: P = Graph(graphs.PetersenGraph(), implementation='c_graph')
sage: list(P._backend.iterator_verts(P))
[0, 1, 2, 3, 4, 5, 6, 7, 8, 9]
loops(new)

Returns whether loops are allowed in this graph.

INPUT:

  • new - boolean (to set) or None (to get)

EXAMPLE:

sage: G = Graph(implementation='c_graph')
sage: G._backend.loops(None)
False
sage: G._backend.loops(True)
sage: G._backend.loops(None)
True
name(new)

Returns the name of this graph.

INPUT:

  • new - boolean (to set) or None (to get)

EXAMPLE:

sage: G = Graph(graphs.PetersenGraph(), implementation='c_graph')
sage: G._backend.name(None)
'Petersen graph'
num_edges(directed)

Returns the number of edges in self.

INPUT:

  • directed - whether to count (u,v) and (v,u) as one or two edges

EXAMPLE:

sage: G = Graph(graphs.PetersenGraph(), implementation='c_graph')
sage: G._backend.num_edges(False)
15
num_verts()

Returns the number of vertices in self.

EXAMPLE:

sage: G = Graph(graphs.PetersenGraph(), implementation='c_graph')
sage: G._backend.num_verts()
10
relabel(perm, directed)

Relabels the graph according to perm.

INPUT:
  • perm - anything which represents a permutation as v --> perm[v], for example a dict or a list
  • directed - ignored (this is here for compatibility with other backends)

EXAMPLE:

sage: G = Graph(graphs.PetersenGraph(), implementation='c_graph')
sage: G._backend.relabel(range(9,-1,-1), False)
sage: G.edges()
[(0, 2, None),
 (0, 3, None),
 (0, 5, None),
 (1, 3, None),
 (1, 4, None),
 (1, 6, None),
 (2, 4, None),
 (2, 7, None),
 (3, 8, None),
 (4, 9, None),
 (5, 6, None),
 (5, 9, None),
 (6, 7, None),
 (7, 8, None),
 (8, 9, None)]
shortest_path(x, y)

Returns the shortest path between x and y

EXAMPLE:

sage: G = Graph(graphs.PetersenGraph(), implementation='c_graph')
sage: G.shortest_path(0,1)
[0, 1]
shortest_path_all_vertices(v, cutoff=None)

Returns for each vertex u a shortest v-u path.

INPUT:

  • v – a vertex
  • cutoff – maximal distance. Longer paths will not be returned

OUTPUT:

A list which associates to each vertex u the shortest path between u and v if there is one.

NOTE:

  • The weight of edges is not taken into account.

ALGORITHM:

This is just a breadth-first search.

EXAMPLES:

On the Petersen Graph:

sage: g = graphs.PetersenGraph()
sage: paths = g._backend.shortest_path_all_vertices(0)
sage: all([ len(paths[v]) == 0 or len(paths[v])-1 == g.distance(0,v) for v in g])
True

On a disconnected graph

sage: g = 2*graphs.RandomGNP(20,.3)
sage: paths = g._backend.shortest_path_all_vertices(0)
sage: all([ (not paths.has_key(v) and g.distance(0,v) == +Infinity) or len(paths[v])-1 == g.distance(0,v) for v in g])
True
strongly_connected_component_containing_vertex(v)

Returns the strongly connected component containing the given vertex

INPUT:

  • v – a vertex

EXAMPLE:

The digraph obtained from the PetersenGraph has an unique strongly connected component

sage: g = DiGraph(graphs.PetersenGraph())
sage: g.strongly_connected_component_containing_vertex(0)
[0, 1, 2, 3, 4, 5, 6, 7, 8, 9]

In the Butterfly DiGraph, each vertex is a strongly connected component:

sage: g = digraphs.ButterflyGraph(3)
sage: all([[v] == g.strongly_connected_component_containing_vertex(v) for v in g])
True
class sage.graphs.base.c_graph.Search_iterator

Bases: object

An iterator for traversing a (di)graph.

This class is commonly used to perform a depth-first or breadth-first search. The class does not build all at once in memory the whole list of visited vertices. The class maintains the following variables:

  • graph – a graph whose vertices are to be iterated over.
  • direction – integer; this determines the position at which vertices to be visited are removed from the list stack. For breadth-first search (BFS), element removal occurs at the start of the list, as signified by the value direction=0. This is because in implementations of BFS, the list of vertices to visit are usually maintained by a queue, so element insertion and removal follow a first-in first-out (FIFO) protocol. For depth-first search (DFS), element removal occurs at the end of the list, as signified by the value direction=-1. The reason is that DFS is usually implemented using a stack to maintain the list of vertices to visit. Hence, element insertion and removal follow a last-in first-out (LIFO) protocol.
  • stack – a list of vertices to visit.
  • seen – a list of vertices that are already visited.
  • test_out – boolean; whether we want to consider the out-neighbors of the graph to be traversed. For undirected graphs, we consider both the in- and out-neighbors. However, for digraphs we only traverse along out-neighbors.
  • test_in – boolean; whether we want to consider the in-neighbors of the graph to be traversed. For undirected graphs, we consider both the in- and out-neighbors.
next()
x.next() -> the next value, or raise StopIteration

Table Of Contents

Previous topic

A module for dealing with lists of graphs.

Next topic

Fast sparse graphs

This Page