Skip to main content
Ctrl+K
Logo image

graph-tool documentation (2.59)

  • graph-tool
  • Welcome to graph-tool’s documentation!
  • Quick start guide
  • Cookbook guides
    • Inferring modular network structure
    • Animations with graph-tool
    • Integration with matplotlib
    • Writing extensions in C++
  • Submodules
    • graph_tool
      • graph_tool.Graph
      • graph_tool.GraphView
      • graph_tool.Vertex
      • graph_tool.Edge
      • graph_tool.PropertyMap
      • graph_tool.VertexPropertyMap
      • graph_tool.EdgePropertyMap
      • graph_tool.GraphPropertyMap
      • graph_tool.PropertyArray
      • graph_tool.group_vector_property
      • graph_tool.ungroup_vector_property
      • graph_tool.map_property_values
      • graph_tool.infect_vertex_property
      • graph_tool.edge_endpoint_property
      • graph_tool.incident_edges_op
      • graph_tool.perfect_prop_hash
      • graph_tool.value_types
      • graph_tool.load_graph
      • graph_tool.load_graph_from_csv
      • graph_tool.openmp_enabled
      • graph_tool.openmp_get_num_threads
      • graph_tool.openmp_set_num_threads
      • graph_tool.openmp_get_schedule
      • graph_tool.openmp_set_schedule
      • graph_tool.show_config
    • graph_tool.centrality
      • graph_tool.centrality.pagerank
      • graph_tool.centrality.betweenness
      • graph_tool.centrality.central_point_dominance
      • graph_tool.centrality.closeness
      • graph_tool.centrality.eigenvector
      • graph_tool.centrality.katz
      • graph_tool.centrality.hits
      • graph_tool.centrality.eigentrust
      • graph_tool.centrality.trust_transitivity
    • graph_tool.clustering
      • graph_tool.clustering.local_clustering
      • graph_tool.clustering.global_clustering
      • graph_tool.clustering.extended_clustering
      • graph_tool.clustering.motifs
      • graph_tool.clustering.motif_significance
    • graph_tool.collection
      • graph_tool.collection.LCF_graph
      • graph_tool.collection.bull_graph
      • graph_tool.collection.chvatal_graph
      • graph_tool.collection.cubical_graph
      • graph_tool.collection.desargues_graph
      • graph_tool.collection.diamond_graph
      • graph_tool.collection.dodecahedral_graph
      • graph_tool.collection.frucht_graph
      • graph_tool.collection.heawood_graph
      • graph_tool.collection.hoffman_singleton_graph
      • graph_tool.collection.house_graph
      • graph_tool.collection.icosahedral_graph
      • graph_tool.collection.krackhardt_kite_graph
      • graph_tool.collection.moebius_kantor_graph
      • graph_tool.collection.octahedral_graph
      • graph_tool.collection.pappus_graph
      • graph_tool.collection.petersen_graph
      • graph_tool.collection.sedgewick_maze_graph
      • graph_tool.collection.tetrahedral_graph
      • graph_tool.collection.truncated_cube_graph
      • graph_tool.collection.truncated_tetrahedron_graph
      • graph_tool.collection.tutte_graph
    • graph_tool.correlations
      • graph_tool.correlations.assortativity
      • graph_tool.correlations.scalar_assortativity
      • graph_tool.correlations.corr_hist
      • graph_tool.correlations.combined_corr_hist
      • graph_tool.correlations.avg_neighbor_corr
      • graph_tool.correlations.avg_combined_corr
    • graph_tool.dynamics
      • graph_tool.dynamics.DiscreteStateBase
      • graph_tool.dynamics.EpidemicStateBase
      • graph_tool.dynamics.SIState
      • graph_tool.dynamics.SISState
      • graph_tool.dynamics.SIRState
      • graph_tool.dynamics.SIRSState
      • graph_tool.dynamics.VoterState
      • graph_tool.dynamics.MajorityVoterState
      • graph_tool.dynamics.BinaryThresholdState
      • graph_tool.dynamics.IsingGlauberState
      • graph_tool.dynamics.CIsingGlauberState
      • graph_tool.dynamics.IsingMetropolisState
      • graph_tool.dynamics.PottsGlauberState
      • graph_tool.dynamics.PottsMetropolisState
      • graph_tool.dynamics.AxelrodState
      • graph_tool.dynamics.BooleanState
      • graph_tool.dynamics.KirmanState
      • graph_tool.dynamics.ContinuousStateBase
      • graph_tool.dynamics.KuramotoState
    • graph_tool.draw
      • graph_tool.draw.sfdp_layout
      • graph_tool.draw.fruchterman_reingold_layout
      • graph_tool.draw.arf_layout
      • graph_tool.draw.radial_tree_layout
      • graph_tool.draw.planar_layout
      • graph_tool.draw.random_layout
      • graph_tool.draw.graph_draw
      • graph_tool.draw.draw_hierarchy
      • graph_tool.draw.graphviz_draw
      • graph_tool.draw.prop_to_size
      • graph_tool.draw.get_hierarchy_control_points
      • graph_tool.draw.cairo_draw
      • graph_tool.draw.interactive_window
      • graph_tool.draw.GraphWidget
      • graph_tool.draw.GraphWindow
    • graph_tool.flow
      • graph_tool.flow.edmonds_karp_max_flow
      • graph_tool.flow.push_relabel_max_flow
      • graph_tool.flow.boykov_kolmogorov_max_flow
      • graph_tool.flow.min_st_cut
      • graph_tool.flow.min_cut
    • graph_tool.generation
      • graph_tool.generation.random_graph
      • graph_tool.generation.random_rewire
      • graph_tool.generation.add_random_edges
      • graph_tool.generation.remove_random_edges
      • graph_tool.generation.generate_triadic_closure
      • graph_tool.generation.price_network
      • graph_tool.generation.generate_sbm
      • graph_tool.generation.generate_maxent_sbm
      • graph_tool.generation.solve_sbm_fugacities
      • graph_tool.generation.generate_knn
      • graph_tool.generation.geometric_graph
      • graph_tool.generation.triangulation
      • graph_tool.generation.predecessor_tree
      • graph_tool.generation.line_graph
      • graph_tool.generation.condensation_graph
      • graph_tool.generation.contract_parallel_edges
      • graph_tool.generation.remove_parallel_edges
      • graph_tool.generation.expand_parallel_edges
      • graph_tool.generation.label_parallel_edges
      • graph_tool.generation.remove_self_loops
      • graph_tool.generation.label_self_loops
      • graph_tool.generation.graph_union
      • graph_tool.generation.lattice
      • graph_tool.generation.complete_graph
      • graph_tool.generation.circular_graph
    • graph_tool.inference
      • graph_tool.inference.minimize_blockmodel_dl
      • graph_tool.inference.minimize_nested_blockmodel_dl
      • graph_tool.inference.BlockState
      • graph_tool.inference.OverlapBlockState
      • graph_tool.inference.LayeredBlockState
      • graph_tool.inference.NestedBlockState
      • graph_tool.inference.PPBlockState
      • graph_tool.inference.RankedBlockState
      • graph_tool.inference.ModularityState
      • graph_tool.inference.NormCutState
      • graph_tool.inference.TemperingState
      • graph_tool.inference.CliqueState
      • graph_tool.inference.MCMCState
      • graph_tool.inference.MultiflipMCMCState
      • graph_tool.inference.MultilevelMCMCState
      • graph_tool.inference.GibbsMCMCState
      • graph_tool.inference.MulticanonicalMCMCState
      • graph_tool.inference.ExhaustiveSweepState
      • graph_tool.inference.DrawBlockState
      • graph_tool.inference.mcmc_equilibrate
      • graph_tool.inference.mcmc_anneal
      • graph_tool.inference.multicanonical_equilibrate
      • graph_tool.inference.MulticanonicalState
      • graph_tool.inference.PartitionModeState
      • graph_tool.inference.ModeClusterState
      • graph_tool.inference.PartitionCentroidState
      • graph_tool.inference.partition_overlap
      • graph_tool.inference.nested_partition_overlap
      • graph_tool.inference.variation_information
      • graph_tool.inference.mutual_information
      • graph_tool.inference.reduced_mutual_information
      • graph_tool.inference.contingency_graph
      • graph_tool.inference.shuffle_partition_labels
      • graph_tool.inference.order_partition_labels
      • graph_tool.inference.order_nested_partition_labels
      • graph_tool.inference.align_partition_labels
      • graph_tool.inference.align_nested_partition_labels
      • graph_tool.inference.partition_overlap_center
      • graph_tool.inference.nested_partition_overlap_center
      • graph_tool.inference.nested_partition_clear_null
      • graph_tool.inference.contiguous_map
      • graph_tool.inference.nested_contiguous_map
      • graph_tool.inference.mf_entropy
      • graph_tool.inference.bethe_entropy
      • graph_tool.inference.microstate_entropy
      • graph_tool.inference.marginal_multigraph_entropy
      • graph_tool.inference.half_edge_graph
      • graph_tool.inference.get_block_edge_gradient
      • graph_tool.inference.PartitionHist
      • graph_tool.inference.BlockPairHist
      • graph_tool.inference.LatentLayerBaseState
      • graph_tool.inference.LatentMultigraphBlockState
      • graph_tool.inference.LatentClosureBlockState
      • graph_tool.inference.MeasuredBlockState
      • graph_tool.inference.MeasuredClosureBlockState
      • graph_tool.inference.MixedMeasuredBlockState
      • graph_tool.inference.UncertainBlockState
      • graph_tool.inference.UncertainBaseState
      • graph_tool.inference.DynamicsBlockStateBase
      • graph_tool.inference.EpidemicsBlockState
      • graph_tool.inference.IsingBaseBlockState
      • graph_tool.inference.IsingGlauberBlockState
      • graph_tool.inference.CIsingGlauberBlockState
      • graph_tool.inference.PseudoIsingBlockState
      • graph_tool.inference.PseudoCIsingBlockState
      • graph_tool.inference.latent_multigraph
      • graph_tool.inference.EMBlockState
      • graph_tool.inference.em_infer
      • graph_tool.inference.modularity
    • graph_tool.search
      • graph_tool.search.bfs_search
      • graph_tool.search.bfs_iterator
      • graph_tool.search.dfs_search
      • graph_tool.search.dfs_iterator
      • graph_tool.search.dijkstra_search
      • graph_tool.search.dijkstra_iterator
      • graph_tool.search.astar_search
      • graph_tool.search.astar_iterator
      • graph_tool.search.bellman_ford_search
      • graph_tool.search.BFSVisitor
      • graph_tool.search.DFSVisitor
      • graph_tool.search.DijkstraVisitor
      • graph_tool.search.BellmanFordVisitor
      • graph_tool.search.AStarVisitor
      • graph_tool.search.StopSearch
    • graph_tool.spectral
      • graph_tool.spectral.adjacency
      • graph_tool.spectral.laplacian
      • graph_tool.spectral.incidence
      • graph_tool.spectral.transition
      • graph_tool.spectral.modularity_matrix
      • graph_tool.spectral.hashimoto
      • graph_tool.spectral.AdjacencyOperator
      • graph_tool.spectral.LaplacianOperator
      • graph_tool.spectral.IncidenceOperator
      • graph_tool.spectral.TransitionOperator
      • graph_tool.spectral.HashimotoOperator
      • graph_tool.spectral.CompactHashimotoOperator
    • graph_tool.stats
      • graph_tool.stats.vertex_hist
      • graph_tool.stats.edge_hist
      • graph_tool.stats.vertex_average
      • graph_tool.stats.edge_average
      • graph_tool.stats.distance_histogram
    • graph_tool.topology
      • graph_tool.topology.shortest_distance
      • graph_tool.topology.shortest_path
      • graph_tool.topology.random_shortest_path
      • graph_tool.topology.count_shortest_paths
      • graph_tool.topology.all_shortest_paths
      • graph_tool.topology.all_predecessors
      • graph_tool.topology.all_paths
      • graph_tool.topology.all_circuits
      • graph_tool.topology.pseudo_diameter
      • graph_tool.topology.similarity
      • graph_tool.topology.vertex_similarity
      • graph_tool.topology.isomorphism
      • graph_tool.topology.subgraph_isomorphism
      • graph_tool.topology.mark_subgraph
      • graph_tool.topology.max_cliques
      • graph_tool.topology.max_cardinality_matching
      • graph_tool.topology.max_independent_vertex_set
      • graph_tool.topology.min_spanning_tree
      • graph_tool.topology.random_spanning_tree
      • graph_tool.topology.dominator_tree
      • graph_tool.topology.topological_sort
      • graph_tool.topology.transitive_closure
      • graph_tool.topology.label_components
      • graph_tool.topology.label_biconnected_components
      • graph_tool.topology.label_largest_component
      • graph_tool.topology.extract_largest_component
      • graph_tool.topology.label_out_component
      • graph_tool.topology.vertex_percolation
      • graph_tool.topology.edge_percolation
      • graph_tool.topology.kcore_decomposition
      • graph_tool.topology.is_bipartite
      • graph_tool.topology.is_DAG
      • graph_tool.topology.is_planar
      • graph_tool.topology.make_maximal_planar
      • graph_tool.topology.edge_reciprocity
      • graph_tool.topology.tsp_tour
      • graph_tool.topology.sequential_vertex_coloring
    • graph_tool.util
      • graph_tool.util.find_vertex
      • graph_tool.util.find_vertex_range
      • graph_tool.util.find_edge
      • graph_tool.util.find_edge_range
  • The gt file format
  • FAQ

Indices

  • General Index
  • Python Module Index
  • .rst

graph_tool.search.dijkstra_search

Contents

  • dijkstra_search()

graph_tool.search.dijkstra_search#

graph_tool.search.dijkstra_search(g, weight, source=None, visitor=<graph_tool.search.DijkstraVisitor object>, dist_map=None, pred_map=None, combine=<function <lambda>>, compare=<function <lambda>>, zero=0, infinity=inf)[source]#

Dijkstra traversal of a directed or undirected graph, with non-negative weights.

Parameters:
gGraph

Graph to be used.

weightEdgePropertyMap

Edge property map with weight values.

sourceVertex (optional, default: None)

Source vertex. If unspecified, all vertices will be traversed, by iterating over starting vertices according to their index in increasing order.

visitorDijkstraVisitor (optional, default: DijkstraVisitor())

A visitor object that is invoked at the event points inside the algorithm. This should be a subclass of DijkstraVisitor.

dist_mapVertexPropertyMap (optional, default: None)

A vertex property map where the distances from the source will be stored.

pred_mapVertexPropertyMap (optional, default: None)

A vertex property map where the predecessor map will be stored (must have value type “int64_t”).

combinebinary function (optional, default: lambda a, b: a + b)

This function is used to combine distances to compute the distance of a path.

comparebinary function (optional, default: lambda a, b: a < b)

This function is use to compare distances to determine which vertex is closer to the source vertex.

zeroint or float (optional, default: 0)

Value assumed to correspond to a distance of zero by the combine and compare functions.

infinityint or float (optional, default: numpy.inf)

Value assumed to correspond to a distance of infinity by the combine and compare functions.

Returns:
dist_mapVertexPropertyMap

A vertex property map with the computed distances from the source.

pred_mapVertexPropertyMap

A vertex property map with the predecessor tree.

See also

bfs_search

Breadth-first search

dfs_search

Depth-first search

astar_search

\(A^*\) heuristic search algorithm

Notes

Dijkstra’s algorithm finds all the shortest paths from the source vertex to every other vertex by iteratively “growing” the set of vertices S to which it knows the shortest path. At each step of the algorithm, the next vertex added to S is determined by a priority queue. The queue contains the vertices in V - S prioritized by their distance label, which is the length of the shortest path seen so far for each vertex. The vertex u at the top of the priority queue is then added to S, and each of its out-edges is relaxed: if the distance to u plus the weight of the out-edge (u,v) is less than the distance label for v then the estimated distance for vertex v is reduced. The algorithm then loops back, processing the next vertex at the top of the priority queue. The algorithm finishes when the priority queue is empty.

The time complexity is \(O(E + V \log V)\).

The pseudo-code for Dijkstra’s algorithm is listed below, with the annotated event points, for which the given visitor object will be called with the appropriate method.

DIJKSTRA(G, source, weight)
  for each vertex u in V                      initialize vertex u
    d[u] := infinity
    p[u] := u
  end for
  d[source] := 0
  INSERT(Q, source)                           discover vertex s
  while (Q != Ø)
    u := EXTRACT-MIN(Q)                       examine vertex u
    for each vertex v in Adj[u]               examine edge (u,v)
      if (weight[(u,v)] + d[u] < d[v])        edge (u,v) relaxed
        d[v] := weight[(u,v)] + d[u]
        p[v] := u
        DECREASE-KEY(Q, v)
      else                                    edge (u,v) not relaxed
        ...
      if (d[v] was originally infinity)
        INSERT(Q, v)                          discover vertex v
    end for                                   finish vertex u
  end while
  return d

References

[dijkstra]

E. Dijkstra, “A note on two problems in connexion with graphs”, Numerische Mathematik, 1:269-271, 1959.

[dijkstra-bgl]

http://www.boost.org/doc/libs/release/libs/graph/doc/dijkstra_shortest_paths_no_color_map.html

[dijkstra-wikipedia]

http://en.wikipedia.org/wiki/Dijkstra’s_algorithm

Examples

We must define what should be done during the search by subclassing DijkstraVisitor, and specializing the appropriate methods. In the following we will keep track of the discover time, and the predecessor tree.

class VisitorExample(gt.DijkstraVisitor):

    def __init__(self, name, time):
        self.name = name
        self.time = time
        self.last_time = 0

    def discover_vertex(self, u):
        print("-->", self.name[u], "has been discovered!")
        self.time[u] = self.last_time
        self.last_time += 1

    def examine_edge(self, e):
        print("edge (%s, %s) has been examined..." % \
            (self.name[e.source()], self.name[e.target()]))

    def edge_relaxed(self, e):
        print("edge (%s, %s) has been relaxed..." % \
            (self.name[e.source()], self.name[e.target()]))

With the above class defined, we can perform the Dijkstra search as follows.

>>> g = gt.load_graph("search_example.xml")
>>> name = g.vp["name"]
>>> weight = g.ep["weight"]
>>> time = g.new_vertex_property("int")
>>> dist, pred = gt.dijkstra_search(g, weight, g.vertex(0), VisitorExample(name, time))
--> Bob has been discovered!
edge (Bob, Eve) has been examined...
edge (Bob, Eve) has been relaxed...
--> Eve has been discovered!
edge (Bob, Chuck) has been examined...
edge (Bob, Chuck) has been relaxed...
--> Chuck has been discovered!
edge (Bob, Carlos) has been examined...
edge (Bob, Carlos) has been relaxed...
--> Carlos has been discovered!
edge (Bob, Isaac) has been examined...
edge (Bob, Isaac) has been relaxed...
--> Isaac has been discovered!
edge (Eve, Isaac) has been examined...
edge (Eve, Imothep) has been examined...
edge (Eve, Imothep) has been relaxed...
--> Imothep has been discovered!
edge (Eve, Carlos) has been examined...
edge (Eve, Chuck) has been examined...
edge (Eve, Bob) has been examined...
edge (Eve, Carol) has been examined...
edge (Eve, Carol) has been relaxed...
--> Carol has been discovered!
edge (Isaac, Bob) has been examined...
edge (Isaac, Chuck) has been examined...
edge (Isaac, Eve) has been examined...
edge (Chuck, Eve) has been examined...
edge (Chuck, Isaac) has been examined...
edge (Chuck, Imothep) has been examined...
edge (Chuck, Bob) has been examined...
edge (Carlos, Eve) has been examined...
edge (Carlos, Imothep) has been examined...
edge (Carlos, Bob) has been examined...
edge (Carlos, Alice) has been examined...
edge (Carlos, Alice) has been relaxed...
--> Alice has been discovered!
edge (Imothep, Carol) has been examined...
edge (Imothep, Carlos) has been examined...
edge (Imothep, Chuck) has been examined...
edge (Imothep, Eve) has been examined...
edge (Alice, Oscar) has been examined...
edge (Alice, Oscar) has been relaxed...
--> Oscar has been discovered!
edge (Alice, Dave) has been examined...
edge (Alice, Dave) has been relaxed...
--> Dave has been discovered!
edge (Alice, Carlos) has been examined...
edge (Carol, Eve) has been examined...
edge (Carol, Imothep) has been examined...
edge (Oscar, Alice) has been examined...
edge (Oscar, Dave) has been examined...
edge (Dave, Oscar) has been examined...
edge (Dave, Alice) has been examined...
>>> print(time.a)
[0 7 6 3 2 9 1 4 8 5]
>>> print(pred.a)
[0 3 6 0 0 1 0 0 1 6]
>>> print(dist.a)
[ 0.          8.91915887  9.27141329  4.29277116  4.02118246 12.23513866
  3.23790211  3.45487436 11.04391549  7.74858396]

previous

graph_tool.search.dfs_iterator

next

graph_tool.search.dijkstra_iterator

Contents
  • dijkstra_search()

© Copyright 2006-2023, Tiago de Paula Peixoto <tiago@skewed.de>.

Last updated on Jan 07, 2024.