Interface to run Boost algorithms#

Wrapper for a Boost graph. The Boost graphs are Cython C++ variables, and they cannot be converted to Python objects: as a consequence, only functions defined with cdef are able to create, read, modify, and delete these graphs.

A very important feature of Boost graph library is that all object are generic: for instance, adjacency lists can be stored using different data structures, and (most of) the functions work with all implementations provided. This feature is implemented in our interface using fused types: however, Cython’s support for fused types is still experimental, and some features are missing. For instance, there cannot be nested generic function calls, and no variable can have a generic type, apart from the arguments of a generic function.

All the input functions use pointers, because otherwise we might have problems with delete().

Basic Boost Graph operations:

clustering_coeff()

Return the clustering coefficient of all vertices in the graph.

edge_connectivity()

Return the edge connectivity of the graph.

dominator_tree()

Return a dominator tree of the graph.

bandwidth_heuristics()

Use heuristics to approximate the bandwidth of the graph.

min_spanning_tree()

Compute a minimum spanning tree of a (weighted) graph.

shortest_paths()

Use Dijkstra or Bellman-Ford algorithm to compute the single-source shortest paths.

johnson_shortest_paths()

Use Johnson algorithm to compute the all-pairs shortest paths.

floyd_warshall_shortest_paths()

Use Floyd-Warshall algorithm to compute the all-pairs shortest paths.

johnson_closeness_centrality()

Use Johnson algorithm to compute the closeness centrality of all vertices.

blocks_and_cut_vertices()

Use Tarjan’s algorithm to compute the blocks and cut vertices of the graph.

min_cycle_basis()

Return a minimum weight cycle basis of the input graph.

Functions#

sage.graphs.base.boost_graph.bandwidth_heuristics(g, algorithm='cuthill_mckee')#

Use Boost heuristics to approximate the bandwidth of the input graph.

The bandwidth \(bw(M)\) of a matrix \(M\) is the smallest integer \(k\) such that all non-zero entries of \(M\) are at distance \(k\) from the diagonal. The bandwidth \(bw(g)\) of an undirected graph \(g\) is the minimum bandwidth of the adjacency matrix of \(g\), over all possible relabellings of its vertices (for more information, see the bandwidth module).

Unfortunately, exactly computing the bandwidth is NP-hard (and an exponential algorithm is implemented in Sagemath in routine bandwidth()). Here, we implement two heuristics to find good orderings: Cuthill-McKee, reverse Cuthill-McKee (also known as RCM) and King.

This function works only in undirected graphs, and its running time is \(O(md_{max}\log d_{max})\) for the Cuthill-McKee ordering, and \(O(md_{max}^2\log d_{max})\) for the King ordering, where \(m\) is the number of edges, and \(d_{max}\) is the maximum degree in the graph.

INPUT:

  • g – the input Sage graph

  • algorithm – string (default: 'cuthill_mckee'); the heuristic used to compute the ordering among 'cuthill_mckee', 'reverse_cuthill_mckee' and 'king'

OUTPUT:

A pair [bandwidth, ordering], where ordering is the ordering of vertices, bandwidth is the bandwidth of that specific ordering (which is not necessarily the bandwidth of the graph, because this is a heuristic).

EXAMPLES:

sage: from sage.graphs.base.boost_graph import bandwidth_heuristics
sage: bandwidth_heuristics(graphs.PathGraph(10))
(1, [0, 1, 2, 3, 4, 5, 6, 7, 8, 9])
sage: bandwidth_heuristics(graphs.GridGraph([3,3]))
(3, [(0, 0), (1, 0), (0, 1), (2, 0), (1, 1), (0, 2), (2, 1), (1, 2), (2, 2)])
sage: bandwidth_heuristics(graphs.GridGraph([3,3]), algorithm='reverse_cuthill_mckee')
(3, [(2, 2), (1, 2), (2, 1), (0, 2), (1, 1), (2, 0), (0, 1), (1, 0), (0, 0)])
sage: bandwidth_heuristics(graphs.GridGraph([3,3]), algorithm='king')
(3, [(0, 0), (1, 0), (0, 1), (2, 0), (1, 1), (0, 2), (2, 1), (1, 2), (2, 2)])
sage.graphs.base.boost_graph.blocks_and_cut_vertices(g)#

Compute the blocks and cut vertices of the graph.

This method uses the implementation of Tarjan’s algorithm available in the Boost library .

INPUT:

  • g – the input Sage graph

OUTPUT:

A 2-dimensional vector with m+1 rows (m is the number of biconnected components), where each of the first m rows correspond to vertices in a block, and the last row is the list of cut vertices.

EXAMPLES:

sage: from sage.graphs.base.boost_graph import blocks_and_cut_vertices
sage: g = graphs.KrackhardtKiteGraph()
sage: blocks_and_cut_vertices(g)
([[8, 9], [7, 8], [0, 1, 2, 3, 5, 4, 6, 7]], [8, 7])

sage: G = Graph([(0,1,{'name':'a','weight':1}), (0,2,{'name':'b','weight':3}), (1,2,{'name':'b','weight':1})])
sage: blocks_and_cut_vertices(G)
([[0, 1, 2]], [])
sage.graphs.base.boost_graph.clustering_coeff(g, vertices=None)#

Compute the clustering coefficient of the input graph, using Boost.

INPUT:

  • g – the input Sage Graph

  • vertices – list (default: None); the list of vertices to analyze (if None, compute the clustering coefficient of all vertices)

OUTPUT: a pair (average_clustering_coefficient, clust_of_v), where average_clustering_coefficient is the average clustering of the vertices in variable vertices, clust_of_v is a dictionary that associates to each vertex its clustering coefficient. If vertices is None, all vertices are considered.

EXAMPLES:

Computing the clustering coefficient of a clique:

sage: from sage.graphs.base.boost_graph import clustering_coeff
sage: g = graphs.CompleteGraph(5)
sage: clustering_coeff(g)
(1.0, {0: 1.0, 1: 1.0, 2: 1.0, 3: 1.0, 4: 1.0})
sage: clustering_coeff(g, vertices = [0,1,2])
(1.0, {0: 1.0, 1: 1.0, 2: 1.0})

Of a non-clique graph with triangles:

sage: g = graphs.IcosahedralGraph()
sage: clustering_coeff(g, vertices=[1,2,3])
(0.5, {1: 0.5, 2: 0.5, 3: 0.5})

With labels:

sage: g.relabel(list("abcdefghiklm"))
sage: clustering_coeff(g, vertices="abde")
(0.5, {'a': 0.5, 'b': 0.5, 'd': 0.5, 'e': 0.5})
sage.graphs.base.boost_graph.diameter(G, algorithm=None, source=None, weight_function=None, check_weight=True)#

Return the diameter of \(G\).

This method returns Infinity if the digraph is not strongly connected. It can also quickly return a lower bound on the diameter using the 2Dsweep scheme.

INPUT:

  • G – the input sage digraph.

  • algorithm – string (default: None); specifies the algorithm to use among:

    • '2Dsweep' – Computes lower bound on the diameter of a weighted directed graph using the weighted version of the algorithm proposed in [Broder2000]. See the code’s documentation for more details.

    • 'DiFUB' – Computes the diameter of a weighted directed graph using the weighted version of the algorithm proposed in [CGLM2012]. See the code’s documentation for more details.

  • source – (default: None) vertex from which to start the computation. If source==None, an arbitrary vertex of the graph is chosen. Raise an error if the initial vertex is not in \(G\).

  • weight_function – function (default: None); a function that associates a weight to each edge. If None (default), the weights of G are used, if G.weighted()==True, otherwise all edges have weight 1.

  • check_weight – boolean (default: True); if True, we check that the weight_function outputs a number for each edge.

EXAMPLES:

sage: from sage.graphs.base.boost_graph import diameter
sage: G = DiGraph([(0, 1, 2), (1, 0, -1)])
sage: diameter(G, algorithm='DiFUB')
1.0
sage: diameter(G, algorithm='DiFUB', weight_function=lambda e:e[2])
2.0
sage: G = DiGraph([(0, 1, -1), (1, 0, 2)])
sage: diameter(G, algorithm='DiFUB', weight_function=lambda e:e[2])
2.0
sage.graphs.base.boost_graph.diameter_DHV(g, weight_function=None, check_weight=True)#

Return the diameter of weighted graph \(g\).

This method computes the diameter of undirected graph using the algorithm proposed in [Dragan2018].

This method returns Infinity if graph is not connected.

INPUT:

  • g – the input Sage graph

  • weight_function – function (default: None); a function that associates a weight to each edge. If None (default), the weights of g are used, if g.weighted()==True, otherwise all edges have weight 1.

  • check_weight – boolean (default: True); if True, we check that the weight_function outputs a number for each edge

EXAMPLES:

sage: from sage.graphs.base.boost_graph import diameter_DHV
sage: G = graphs.ButterflyGraph()
sage: diameter_DHV(G)
2.0
sage.graphs.base.boost_graph.dominator_tree(g, root, return_dict=False, reverse=False)#

Use Boost to compute the dominator tree of g, rooted at root.

A node \(d\) dominates a node \(n\) if every path from the entry node root to \(n\) must go through \(d\). The immediate dominator of a node \(n\) is the unique node that strictly dominates \(n\) but does not dominate any other node that dominates \(n\). A dominator tree is a tree where each node’s children are those nodes it immediately dominates. For more information, see the Wikipedia article Dominator_(graph_theory).

If the graph is connected and undirected, the parent of a vertex \(v\) is:

  • the root if \(v\) is in the same biconnected component as the root;

  • the first cut vertex in a path from \(v\) to the root, otherwise.

If the graph is not connected, the dominator tree of the whole graph is equal to the dominator tree of the connected component of the root.

If the graph is directed, computing a dominator tree is more complicated, and it needs time \(O(m\log m)\), where \(m\) is the number of edges. The implementation provided by Boost is the most general one, so it needs time \(O(m\log m)\) even for undirected graphs.

INPUT:

  • g – the input Sage (Di)Graph

  • root – the root of the dominator tree

  • return_dict – boolean (default: False); if True, the function returns a dictionary associating to each vertex its parent in the dominator tree. If False (default), it returns the whole tree, as a Graph or a DiGraph.

  • reverse – boolean (default: False); when set to True, computes the dominator tree in the reverse graph

OUTPUT:

The dominator tree, as a graph or as a dictionary, depending on the value of return_dict. If the output is a dictionary, it will contain None in correspondence of root and of vertices that are not reachable from root. If the output is a graph, it will not contain vertices that are not reachable from root.

EXAMPLES:

An undirected grid is biconnected, and its dominator tree is a star (everyone’s parent is the root):

sage: g = graphs.GridGraph([2,2]).dominator_tree((0,0))
sage: g.to_dictionary()
{(0, 0): [(0, 1), (1, 0), (1, 1)], (0, 1): [(0, 0)], (1, 0): [(0, 0)], (1, 1): [(0, 0)]}

If the graph is made by two 3-cycles \(C_1,C_2\) connected by an edge \((v,w)\), with \(v \in C_1\), \(w \in C_2\), the cut vertices are \(v\) and \(w\), the biconnected components are \(C_1\), \(C_2\), and the edge \((v,w)\). If the root is in \(C_1\), the parent of each vertex in \(C_1\) is the root, the parent of \(w\) is \(v\), and the parent of each vertex in \(C_2\) is \(w\):

sage: G = 2 * graphs.CycleGraph(3)
sage: v = 0
sage: w = 3
sage: G.add_edge(v,w)
sage: G.dominator_tree(1, return_dict=True)
{0: 1, 1: None, 2: 1, 3: 0, 4: 3, 5: 3}

An example with a directed graph:

sage: g = digraphs.Circuit(10).dominator_tree(5)
sage: g.to_dictionary()
{0: [1], 1: [2], 2: [3], 3: [4], 4: [], 5: [6], 6: [7], 7: [8], 8: [9], 9: [0]}
sage: g = digraphs.Circuit(10).dominator_tree(5, reverse=True)
sage: g.to_dictionary()
{0: [9], 1: [0], 2: [1], 3: [2], 4: [3], 5: [4], 6: [], 7: [6], 8: [7], 9: [8]}

If the output is a dictionary:

sage: graphs.GridGraph([2,2]).dominator_tree((0,0), return_dict=True)
{(0, 0): None, (0, 1): (0, 0), (1, 0): (0, 0), (1, 1): (0, 0)}
sage.graphs.base.boost_graph.eccentricity_DHV(g, vertex_list=None, weight_function=None, check_weight=True)#

Return the vector of eccentricities using the algorithm of [Dragan2018].

The array returned is of length \(n\), and by default its \(i\)-th component is the eccentricity of the \(i\)-th vertex in g.vertices(sort=True), if vertex_list is None, otherwise ecc[i] is the eccentricity of vertex vertex_list[i].

The algorithm proposed in [Dragan2018] is based on the observation that for all nodes \(v,w\in V\), we have \(\max(ecc[v]-d(v,w), d(v,w))\leq ecc[w] \leq ecc[v] + d(v,w)\). Also the algorithm iteratively improves upper and lower bounds on the eccentricity of each vertex until no further improvements can be done.

INPUT:

  • g – the input Sage graph.

  • vertex_list – list (default: None); a list of \(n\) vertices specifying a mapping from \((0, \ldots, n-1)\) to vertex labels in \(g\). When set, ecc[i] is the eccentricity of vertex vertex_list[i]. When vertex_list is None, ecc[i] is the eccentricity of vertex g.vertices(sort=True)[i].

  • weight_function – function (default: None); a function that associates a weight to each edge. If None (default), the weights of g are used, if g.weighted()==True, otherwise all edges have weight 1.

  • check_weight – boolean (default: True); if True, we check that the weight_function outputs a number for each edge

EXAMPLES:

sage: from sage.graphs.base.boost_graph import eccentricity_DHV
sage: G = graphs.BullGraph()
sage: eccentricity_DHV(G)
[2.0, 2.0, 2.0, 3.0, 3.0]
sage.graphs.base.boost_graph.edge_connectivity(g)#

Compute the edge connectivity of the input graph, using Boost.

OUTPUT: a pair (ec, edges), where ec is the edge connectivity, edges is the list of edges in a minimum cut.

EXAMPLES:

Computing the edge connectivity of a clique:

sage: from sage.graphs.base.boost_graph import edge_connectivity
sage: g = graphs.CompleteGraph(5)
sage: edge_connectivity(g)
(4, [(0, 1), (0, 2), (0, 3), (0, 4)])

Vertex-labeled graphs:

sage: from sage.graphs.base.boost_graph import edge_connectivity
sage: g = graphs.GridGraph([2,2])
sage: edge_connectivity(g)
(2, [((0, 0), (0, 1)), ((0, 0), (1, 0))])
sage.graphs.base.boost_graph.floyd_warshall_shortest_paths(g, weight_function=None, distances=True, predecessors=False)#

Use Floyd-Warshall algorithm to solve the all-pairs-shortest-paths.

This routine outputs the distance between each pair of vertices and the predecessors matrix (depending on the values of boolean distances and predecessors) using a dictionary of dictionaries. This method should be preferred only if the graph is dense. If the graph is sparse the much faster johnson_shortest_paths should be used.

The time-complexity is \(O(n^3 + nm)\), where \(n\) is the number of nodes and \(m\) the number of edges. The factor \(nm\) in the complexity is added only when predecessors is set to True.

INPUT:

  • g – the input Sage graph

  • weight_function – function (default: None); a function that associates a weight to each edge. If None (default), the weights of g are used, if g.weighted()==True, otherwise all edges have weight 1.

  • distances – boolean (default: True); whether to return the dictionary of shortest distances

  • predecessors – boolean (default: False); whether to return the predecessors matrix

OUTPUT:

Depending on the input, this function return the dictionary of predecessors, the dictionary of distances, or a pair of dictionaries (distances, predecessors) where distance[u][v] denotes the distance of a shortest path from \(u\) to \(v\) and predecessors[u][v] indicates the predecessor of \(w\) on a shortest path from \(u\) to \(v\).

EXAMPLES:

Undirected graphs:

sage: from sage.graphs.base.boost_graph import floyd_warshall_shortest_paths
sage: g = Graph([(0,1,1),(1,2,2),(1,3,4),(2,3,1)], weighted=True)
sage: floyd_warshall_shortest_paths(g)
{0: {0: 0, 1: 1, 2: 3, 3: 4},
 1: {0: 1, 1: 0, 2: 2, 3: 3},
 2: {0: 3, 1: 2, 2: 0, 3: 1},
 3: {0: 4, 1: 3, 2: 1, 3: 0}}
sage: expected = {0: {0: None, 1: 0, 2: 1, 3: 2},
....:             1: {0: 1, 1: None, 2: 1, 3: 2},
....:             2: {0: 1, 1: 2, 2: None, 3: 2},
....:             3: {0: 1, 1: 2, 2: 3, 3: None}}
sage: floyd_warshall_shortest_paths(g, distances=False, predecessors=True) == expected
True

Directed graphs:

sage: g = DiGraph([(0,1,1),(1,2,-2),(1,3,4),(2,3,1)], weighted=True)
sage: floyd_warshall_shortest_paths(g)
{0: {0: 0, 1: 1, 2: -1, 3: 0},
 1: {1: 0, 2: -2, 3: -1},
 2: {2: 0, 3: 1},
 3: {3: 0}}
sage: g = DiGraph([(1,2,3),(2,3,2),(1,4,1),(4,2,1)], weighted=True)
sage: floyd_warshall_shortest_paths(g, distances=False, predecessors=True)
{1: {1: None, 2: 4, 3: 2, 4: 1},
 2: {2: None, 3: 2},
 3: {3: None},
 4: {2: 4, 3: 2, 4: None}}
sage.graphs.base.boost_graph.johnson_closeness_centrality(g, weight_function=None)#

Use Johnson algorithm to compute the closeness centrality of all vertices.

This routine is preferable to johnson_shortest_paths() because it does not create a doubly indexed dictionary of distances, saving memory.

The time-complexity is \(O(mn\log n)\), where \(n\) is the number of nodes and \(m\) is the number of edges.

INPUT:

  • g – the input Sage graph

  • weight_function – function (default: None); a function that associates a weight to each edge. If None (default), the weights of g are used, if g.weighted()==True, otherwise all edges have weight 1.

OUTPUT:

A dictionary associating each vertex v to its closeness centrality.

EXAMPLES:

Undirected graphs:

sage: from sage.graphs.base.boost_graph import johnson_closeness_centrality
sage: g = Graph([(0,1,1),(1,2,2),(1,3,4),(2,3,1)], weighted=True)
sage: johnson_closeness_centrality(g)
{0: 0.375, 1: 0.5, 2: 0.5, 3: 0.375}

Directed graphs:

sage: from sage.graphs.base.boost_graph import johnson_closeness_centrality
sage: g = DiGraph([(0,1,1),(1,2,-2),(1,3,4),(2,3,1)], weighted=True)
sage: johnson_closeness_centrality(g)
{0: inf, 1: -0.4444444444444444, 2: 0.3333333333333333}
sage.graphs.base.boost_graph.johnson_shortest_paths(g, weight_function=None, distances=True, predecessors=False)#

Use Johnson algorithm to solve the all-pairs-shortest-paths.

This routine outputs the distance between each pair of vertices and the predecessors matrix (depending on the values of boolean distances and predecessors) using a dictionary of dictionaries. It works on all kinds of graphs, but it is designed specifically for graphs with negative weights (otherwise there are more efficient algorithms, like Dijkstra).

The time-complexity is \(O(mn\log n)\), where \(n\) is the number of nodes and \(m\) is the number of edges.

INPUT:

  • g – the input Sage graph

  • weight_function – function (default: None); a function that associates a weight to each edge. If None (default), the weights of g are used, if g.weighted()==True, otherwise all edges have weight 1.

  • distances – boolean (default: True); whether to return the dictionary of shortest distances

  • predecessors – boolean (default: False); whether to return the predecessors matrix

OUTPUT:

Depending on the input, this function return the dictionary of predecessors, the dictionary of distances, or a pair of dictionaries (distances, predecessors) where distance[u][v] denotes the distance of a shortest path from \(u\) to \(v\) and predecessors[u][v] indicates the predecessor of \(w\) on a shortest path from \(u\) to \(v\).

EXAMPLES:

Undirected graphs:

sage: from sage.graphs.base.boost_graph import johnson_shortest_paths
sage: g = Graph([(0,1,1),(1,2,2),(1,3,4),(2,3,1)], weighted=True)
sage: johnson_shortest_paths(g)
{0: {0: 0, 1: 1, 2: 3, 3: 4},
 1: {0: 1, 1: 0, 2: 2, 3: 3},
 2: {0: 3, 1: 2, 2: 0, 3: 1},
 3: {0: 4, 1: 3, 2: 1, 3: 0}}
sage: expected = {0: {0: None, 1: 0, 2: 1, 3: 2},
....:             1: {0: 1, 1: None, 2: 1, 3: 2},
....:             2: {0: 1, 1: 2, 2: None, 3: 2},
....:             3: {0: 1, 1: 2, 2: 3, 3: None}}
sage: johnson_shortest_paths(g, distances=False, predecessors=True) == expected
True

Directed graphs:

sage: g = DiGraph([(0,1,1),(1,2,-2),(1,3,4),(2,3,1)], weighted=True)
sage: johnson_shortest_paths(g)
{0: {0: 0, 1: 1, 2: -1, 3: 0},
 1: {1: 0, 2: -2, 3: -1},
 2: {2: 0, 3: 1},
 3: {3: 0}}
sage: g = DiGraph([(1,2,3),(2,3,2),(1,4,1),(4,2,1)], weighted=True)
sage: johnson_shortest_paths(g, distances=False, predecessors=True)
{1: {1: None, 2: 4, 3: 2, 4: 1},
 2: {2: None, 3: 2},
 3: {3: None},
 4: {2: 4, 3: 2, 4: None}}
sage.graphs.base.boost_graph.min_cycle_basis(g_sage, weight_function=None, by_weight=False)#

Return a minimum weight cycle basis of the input graph g_sage.

A cycle basis is a list of cycles (list of vertices forming a cycle) of g_sage. Note that the vertices are not necessarily returned in the order in which they appear in the cycle.

A minimum weight cycle basis is a cycle basis that minimizes the sum of the weights (length for unweighted graphs) of its cycles.

Not implemented for directed graphs and multigraphs.

INPUT:

  • g_sage – a Sage Graph

  • weight_function – function (default: None); a function that takes as input an edge (u, v, l) and outputs its weight. If not None, by_weight is automatically set to True. If None and by_weight is True, the weights of g_sage are used, if g_sage.weighted()==True, otherwise all edges have weight 1.

  • by_weight – boolean (default: False); if True, the edges in the graph are weighted, otherwise all edges have weight 1

EXAMPLES:

sage: g = Graph([(1, 2, 3), (2, 3, 5), (3, 4, 8), (4, 1, 13), (1, 3, 250), (5, 6, 9), (6, 7, 17), (7, 5, 20)])
sage: sorted(g.minimum_cycle_basis(by_weight=True))
[[1, 2, 3], [1, 2, 3, 4], [5, 6, 7]]
sage: sorted(g.minimum_cycle_basis())
[[1, 2, 3], [1, 3, 4], [5, 6, 7]]
sage.graphs.base.boost_graph.min_spanning_tree(g, weight_function=None, algorithm='Kruskal')#

Use Boost to compute the minimum spanning tree of the input graph.

INPUT:

  • g – the input Sage graph

  • weight_function – function (default: None); a function that inputs an edge e and outputs its weight. An edge has the form (u,v,l), where u and v are vertices, l is a label (that can be of any kind). The weight_function can be used to transform the label into a weight (see the example below). In particular:

    • if weight_function is not None, the weight of an edge e is weight_function(e);

    • if weight_function is None (default) and g is weighted (that is, g.weighted()==True), for each edge e=(u,v,l), we set weight l;

    • if weight_function is None and g is not weighted, we set all weights to 1 (hence, the output can be any spanning tree).

    Note that, if the weight is not convertible to a number with function float(), an error is raised (see tests below).

  • algorithm – string (default: 'Kruskal'); the algorithm to use among 'Kruskal' and 'Prim'

OUTPUT:

The edges of a minimum spanning tree of g, if one exists, otherwise the empty list.

EXAMPLES:

sage: from sage.graphs.base.boost_graph import min_spanning_tree
sage: min_spanning_tree(graphs.PathGraph(4))
[(0, 1, None), (1, 2, None), (2, 3, None)]

sage: G = Graph([(0,1,{'name':'a','weight':1}), (0,2,{'name':'b','weight':3}), (1,2,{'name':'b','weight':1})])
sage: min_spanning_tree(G, weight_function=lambda e: e[2]['weight'])
[(0, 1, {'name': 'a', 'weight': 1}), (1, 2, {'name': 'b', 'weight': 1})]
sage.graphs.base.boost_graph.radius_DHV(g, weight_function=None, check_weight=True)#

Return the radius of weighted graph \(g\).

This method computes the radius of undirected graph using the algorithm given in [Dragan2018].

This method returns Infinity if graph is not connected.

INPUT:

  • g – the input Sage graph

  • weight_function – function (default: None); a function that associates a weight to each edge. If None (default), the weights of g are used, if g.weighted()==True, otherwise all edges have weight 1.

  • check_weight – boolean (default: True); if True, we check that the weight_function outputs a number for each edge.

EXAMPLES:

sage: from sage.graphs.base.boost_graph import radius_DHV
sage: G = Graph([(0,1,1), (1,2,1), (0,2,3)])
sage: radius_DHV(G)
1.0
sage: G = graphs.PathGraph(7)
sage: radius_DHV(G) == G.radius(algorithm='Dijkstra_Boost')
True
sage.graphs.base.boost_graph.shortest_paths(g, start, weight_function=None, algorithm=None)#

Compute the shortest paths from start to all other vertices.

This routine outputs all shortest paths from node start to any other node in the graph. The input graph can be weighted: if the algorithm is Dijkstra, no negative weights are allowed, while if the algorithm is Bellman-Ford, negative weights are allowed, but there must be no negative cycle (otherwise, the shortest paths might not exist).

However, Dijkstra algorithm is more efficient: for this reason, we suggest to use Bellman-Ford only if necessary (which is also the default option). Note that, if the graph is undirected, a negative edge automatically creates a negative cycle: for this reason, in this case, Dijkstra algorithm is always better.

The running-time is \(O(n \log n+m)\) for Dijkstra algorithm and \(O(mn)\) for Bellman-Ford algorithm, where \(n\) is the number of nodes and \(m\) is the number of edges.

INPUT:

  • g – the input Sage graph

  • start – the starting vertex to compute shortest paths

  • weight_function – function (default: None); a function that associates a weight to each edge. If None (default), the weights of g are used, if g.weighted()==True, otherwise all edges have weight 1.

  • algorithm – string (default: None); one of the following algorithms:

    • 'Dijkstra', 'Dijkstra_Boost': the Dijkstra algorithm implemented in Boost (works only with positive weights)

    • 'Bellman-Ford', 'Bellman-Ford_Boost': the Bellman-Ford algorithm implemented in Boost (works also with negative weights, if there is no negative cycle)

OUTPUT:

A pair of dictionaries (distances, predecessors) such that, for each vertex v, distances[v] is the distance from start to v, predecessors[v] is the last vertex in a shortest path from start to v.

EXAMPLES:

Undirected graphs:

sage: from sage.graphs.base.boost_graph import shortest_paths
sage: g = Graph([(0,1,1),(1,2,2),(1,3,4),(2,3,1)], weighted=True)
sage: shortest_paths(g, 1)
({0: 1, 1: 0, 2: 2, 3: 3}, {0: 1, 1: None, 2: 1, 3: 2})
sage: g = graphs.GridGraph([2,2])
sage: shortest_paths(g,(0,0),weight_function=lambda e:2)
({(0, 0): 0, (0, 1): 2, (1, 0): 2, (1, 1): 4},
 {(0, 0): None, (0, 1): (0, 0), (1, 0): (0, 0), (1, 1): (0, 1)})

Directed graphs:

sage: g = DiGraph([(0,1,1),(1,2,2),(1,3,4),(2,3,1)], weighted=True)
sage: shortest_paths(g, 1)
({1: 0, 2: 2, 3: 3}, {1: None, 2: 1, 3: 2})
sage.graphs.base.boost_graph.shortest_paths_from_vertices(g, vertex_list=None, order=None, weight_function=None, algorithm=None)#

Compute the shortest paths to all vertices from each vertex in vertex_list.

The input graph can be weighted: if the algorithm is Dijkstra, no negative weights are allowed, while if the algorithm is Bellman-Ford, negative weights are allowed, but there must be no negative cycle (otherwise, the shortest paths might not exist).

However, Dijkstra algorithm is more efficient: for this reason, we suggest to use Bellman-Ford only if necessary (which is also the default option).

The running-time for each vertex is \(O(n \log n+m)\) for Dijkstra algorithm and \(O(mn)\) for Bellman-Ford algorithm, where \(n\) is the number of nodes and \(m\) is the number of edges.

INPUT:

  • g – the input Sage graph

  • vertex_list – list (default: None); list of vertices to compute shortest paths from. By default (None), compute shortest paths from all vertices.

  • order – list (default: None); order of vertices of \(g\)

  • weight_function – function (default: None); a function that associates a weight to each edge. If None (default), the weights of g are used, if g.weighted()==True, otherwise all edges have weight 1.

  • algorithm – string (default: None); one of the following algorithms:

    • 'Dijkstra', 'Dijkstra_Boost' - the Dijkstra algorithm implemented in Boost (works only with positive weights)

    • 'Bellman-Ford', 'Bellman-Ford_Boost' - the Bellman-Ford algorithm implemented in Boost (works also with negative weights, if there is no negative cycle)

OUTPUT:

The type of output depends on the input. More precisely -

  • A pair of dictionaries of list (distances, predecessors), when order is not None, such that for each vertex v in vertex_list, distances[v][i] store the shortest distance between v and order[i] and predecessors[v][i] store the last vertex in the shortest path from v to order[i].

  • A pair of dictionaries of dictionaries (distances, predecessors) such that for each vertex v in vertex_list, distances[v] store the shortest distances of all the other vertices from v, predecessors[v] store the last vertices in the shortest path from v to all the other vertices.

EXAMPLES:

Undirected graphs:

sage: from sage.graphs.base.boost_graph import shortest_paths_from_vertices
sage: g = Graph([(0,1,1),(1,2,2),(1,3,4),(2,3,1)], weighted=True)
sage: shortest_paths_from_vertices(g,[1,2])
({1: {0: 1.0, 1: 0.0, 2: 2.0, 3: 3.0}, 2: {0: 3.0, 1: 2.0, 2: 0.0, 3: 1.0}},
 {1: {0: 1, 1: None, 2: 1, 3: 2}, 2: {0: 1, 1: 2, 2: None, 3: 2}})

Directed graphs:

sage: g = DiGraph([(0,1,1),(1,2,-1),(2,0,2),(2,3,1)], weighted=True)
sage: shortest_paths_from_vertices(g,1)
({1: {0: 1.0, 1: 0.0, 2: -1.0, 3: 0.0}}, {1: {0: 2, 1: None, 2: 1, 3: 2}})
sage: shortest_paths_from_vertices(g, 1, [0,1,2,3])
({1: [1.0, 0.0, -1.0, 0.0]}, {1: [2, None, 1, 2]})
sage.graphs.base.boost_graph.wiener_index(g, algorithm=None, weight_function=None, check_weight=True)#

Return the Wiener index of the graph.

The Wiener index of an undirected graph \(G\) is defined as \(W(G) = \frac{1}{2} \sum_{u,v\in G} d(u,v)\) where \(d(u,v)\) denotes the distance between vertices \(u\) and \(v\) (see [KRG1996]).

The Wiener index of a directed graph \(G\) is defined as the sum of the distances between each pairs of vertices, \(W(G) = \sum_{u,v\in G} d(u,v)\).

INPUT:

  • g – the input Sage graph

  • algorithm – string (default: None); one of the following algorithms:

    • 'Dijkstra', 'Dijkstra_Boost': the Dijkstra algorithm implemented in Boost (works only with positive weights)

    • 'Bellman-Ford', 'Bellman-Ford_Boost': the Bellman-Ford algorithm implemented in Boost (works also with negative weights, if there is no negative cycle)

  • weight_function – function (default: None); a function that associates a weight to each edge. If None (default), the weights of g are used, if g.weighted()==True, otherwise all edges have weight 1.

  • check_weight – boolean (default: True); if True, we check that the weight_function outputs a number for each edge.

EXAMPLES:

sage: from sage.graphs.base.boost_graph import wiener_index sage: g = Graph([(0,1,9), (1,2,7), (2,3,4), (3,0,3)]) sage: wiener_index(g) 8.0 sage: g.weighted(True) sage: wiener_index(g) 41.0

Wiener index of circuit digraphs:

sage: n = 10
sage: g = digraphs.Circuit(n)
sage: w = lambda x: (x*x*(x-1))/2
sage: wiener_index(g) == w(n)
True

Wiener index of a graph of order 1:

sage: wiener_index(Graph(1))
0

The Wiener index is not defined on the empty graph:

sage: wiener_index(Graph())
Traceback (most recent call last):
...
ValueError: Wiener index is not defined for the empty graph