Path and Traversal

Graphs.jl provides several traversal and shortest-path algorithms, along with various utility functions. Where appropriate, edge distances may be passed in as a matrix of real number values.

Edge distances for most traversals may be passed in as a sparse or dense matrix of values, indexed by [src,dst] vertices. That is, distmx[2,4] = 2.5 assigns the distance 2.5 to the (directed) edge connecting vertex 2 and vertex 4. Note that also for undirected graphs distmx[4,2] has to be set.

Default edge distances may be passed in via the

Graphs.DefaultDistanceType
DefaultDistance

An array-like structure that provides distance values of 1 for any src, dst combination.

structure.

Any graph traversal will traverse an edge only if it is present in the graph. When a distance matrix is passed in,

  1. distance values for undefined edges will be ignored, and
  2. any unassigned values (in sparse distance matrices), for edges that are present in the graph, will be assumed to take the default value of 1.0.
  3. any zero values (in sparse/dense distance matrices), for edges that are present in the graph, will instead have an implicit edge cost of 1.0.

Graph Traversal

Graph traversal refers to a process that traverses vertices of a graph following certain order (starting from user-input sources). This package implements three traversal schemes:

  • BreadthFirst,
  • DepthFirst, and
  • MaximumAdjacency.
Graphs.bfs_treeFunction
bfs_tree(g, s[; dir=:out])

Provide a breadth-first traversal of the graph g starting with source vertex s, and return a directed acyclic graph of vertices in the order they were discovered. If dir is specified, use the corresponding edge direction (:in and :out are acceptable values).

Graphs.dfs_treeFunction
dfs_tree(g, s)

Return an ordered vector of vertices representing a directed acyclic graph based on depth-first traversal of the graph g starting with source vertex s.

Graphs.maximum_adjacency_visitFunction
maximum_adjacency_visit(g[, distmx][, log][, io])

Return the vertices in g traversed by maximum adjacency search. An optional distmx matrix may be specified; if omitted, edge distances are assumed to be 1. If log (default false) is true, visitor events will be printed to io, which defaults to STDOUT; otherwise, no event information will be displayed.

Graphs.bfs_parentsFunction
bfs_parents(g, s[; dir=:out])

Perform a breadth-first search of graph g starting from vertex s. Return a vector of parent vertices indexed by vertex. If dir is specified, use the corresponding edge direction (:in and :out are acceptable values).

Performance

This implementation is designed to perform well on large graphs. There are implementations which are marginally faster in practice for smaller graphs, but the performance improvements using this implementation on large graphs can be significant.

Graphs.has_pathFunction
has_path(g::AbstractGraph, u, v; exclude_vertices=Vector())

Return true if there is a path from u to v in g (while avoiding vertices in exclude_vertices) or u == v. Return false if there is no such path or if u or v is in excluded_vertices.

Graphs.diffusionFunction
diffusion(g, p, n)

Run diffusion simulation on g for n steps with spread probabilities based on p. Return a vector with the set of new vertices reached at each step of the simulation.

Optional Arguments

  • initial_infections=sample(vertices(g), 1): A list of vertices that

are infected at the start of the simulation.

  • watch=Vector(): While simulation is always run on the full graph,

specifying watch limits reporting to a specific set of vertices reached during the simulation. If left empty, all vertices will be watched.

  • normalize=false: if false, set the probability of spread from a vertex $i$ to

each of the outneighbors of $i$ to $p$. If true, set the probability of spread from a vertex $i$ to each of the outneighbors of $i$ to $\frac{p}{outdegreee(g, i)}$.

Graphs.diffusion_rateFunction
diffusion_rate(results)
diffusion_rate(g, p, n; ...)

Given the results of a diffusion output or the parameters to the diffusion simulation itself, (run and) return the rate of diffusion as a vector representing the cumulative number of vertices infected at each simulation step, restricted to vertices included in watch, if specified.

Graphs.mincutFunction
mincut(g, distmx=weights(g))

Return a tuple (parity, bestcut), where parity is a vector of integer values that determines the partition in g (1 or 2) and bestcut is the weight of the cut that makes this partition. An optional distmx matrix may be specified; if omitted, edge distances are assumed to be 1.

Random walks

Graphs includes uniform random walks and self avoiding walks:

Graphs.randomwalkFunction
randomwalk(g, s, niter; seed=-1)

Perform a random walk on graph g starting at vertex s and continuing for a maximum of niter steps. Return a vector of vertices visited in order.

Graphs.non_backtracking_randomwalkFunction
non_backtracking_randomwalk(g, s, niter; seed=-1)

Perform a non-backtracking random walk on directed graph g starting at vertex s and continuing for a maximum of niter steps. Return a vector of vertices visited in order.

Graphs.self_avoiding_walkFunction
self_avoiding_walk(g, s, niter; seed=-1)

Perform a self-avoiding walk on graph g starting at vertex s and continuing for a maximum of niter steps. Return a vector of vertices visited in order.

Connectivity / Bipartiteness

Graph connectivity functions are defined on both undirected and directed graphs:

Graphs.is_connectedFunction
is_connected(g)

Return true if graph g is connected. For directed graphs, return true if graph g is weakly connected.

Examples

julia> g = SimpleGraph([0 1 0; 1 0 1; 0 1 0]);

julia> is_connected(g)
true

julia> g = SimpleGraph([0 1 0 0 0; 1 0 1 0 0; 0 1 0 0 0; 0 0 0 0 1; 0 0 0 1 0]);

julia> is_connected(g)
false

julia> g = SimpleDiGraph([0 1 0; 0 0 1; 1 0 0]);

julia> is_connected(g)
true
Graphs.is_strongly_connectedFunction
is_strongly_connected(g)

Return true if directed graph g is strongly connected.

Examples

julia> g = SimpleDiGraph([0 1 0; 0 0 1; 1 0 0]);

julia> is_strongly_connected(g)
true
Graphs.is_weakly_connectedFunction
is_weakly_connected(g)

Return true if the graph g is weakly connected. If g is undirected, this function is equivalent to is_connected(g).

Examples

julia> g = SimpleDiGraph([0 1 0; 0 0 1; 1 0 0]);

julia> is_weakly_connected(g)
true

julia> g = SimpleDiGraph([0 1 0; 1 0 1; 0 0 0]);

julia> is_connected(g)
true

julia> is_strongly_connected(g)
false

julia> is_weakly_connected(g)
true
Graphs.connected_componentsFunction
connected_components(g)

Return the connected components of an undirected graph g as a vector of components, with each element a vector of vertices belonging to the component.

For directed graphs, see strongly_connected_components and weakly_connected_components.

Examples

julia> g = SimpleGraph([0 1 0; 1 0 1; 0 1 0]);

julia> connected_components(g)
1-element Array{Array{Int64,1},1}:
 [1, 2, 3]

julia> g = SimpleGraph([0 1 0 0 0; 1 0 1 0 0; 0 1 0 0 0; 0 0 0 0 1; 0 0 0 1 0]);

julia> connected_components(g)
2-element Array{Array{Int64,1},1}:
 [1, 2, 3]
 [4, 5]
Graphs.strongly_connected_componentsFunction
strongly_connected_components(g)

Compute the strongly connected components of a directed graph g.

Return an array of arrays, each of which is the entire connected component.

Implementation Notes

The order of the components is not part of the API contract.

Examples

julia> g = SimpleDiGraph([0 1 0; 1 0 1; 0 0 0]);

julia> strongly_connected_components(g)
2-element Array{Array{Int64,1},1}:
 [3]
 [1, 2]


julia> g=SimpleDiGraph(11)
{11, 0} directed simple Int64 graph

julia> edge_list=[(1,2),(2,3),(3,4),(4,1),(3,5),(5,6),(6,7),(7,5),(5,8),(8,9),(9,8),(10,11),(11,10)];

julia> g = SimpleDiGraph(Edge.(edge_list))
{11, 13} directed simple Int64 graph

julia> strongly_connected_components(g)
4-element Array{Array{Int64,1},1}:
 [8, 9]      
 [5, 6, 7]   
 [1, 2, 3, 4]
 [10, 11]    
Graphs.strongly_connected_components_kosarajuFunction
strongly_connected_components_kosaraju(g)

Compute the strongly connected components of a directed graph g using Kosaraju's Algorithm. (https://en.wikipedia.org/wiki/Kosaraju%27s_algorithm).

Return an array of arrays, each of which is the entire connected component.

Performance

Time Complexity : O(|E|+|V|) Space Complexity : O(|V|) {Excluding the memory required for storing graph}

|V| = Number of vertices |E| = Number of edges

Examples


julia> g=SimpleDiGraph(3)
{3, 0} directed simple Int64 graph

julia> g = SimpleDiGraph([0 1 0 ; 0 0 1; 0 0 0])
{3, 2} directed simple Int64 graph

julia> strongly_connected_components_kosaraju(g)
3-element Array{Array{Int64,1},1}:
 [1]
 [2]
 [3]


julia> g=SimpleDiGraph(11)
{11, 0} directed simple Int64 graph

julia> edge_list=[(1,2),(2,3),(3,4),(4,1),(3,5),(5,6),(6,7),(7,5),(5,8),(8,9),(9,8),(10,11),(11,10)]
13-element Array{Tuple{Int64,Int64},1}:
 (1, 2)  
 (2, 3)  
 (3, 4)  
 (4, 1)  
 (3, 5)  
 (5, 6)  
 (6, 7)  
 (7, 5)  
 (5, 8)  
 (8, 9)  
 (9, 8)  
 (10, 11)
 (11, 10)

julia> g = SimpleDiGraph(Edge.(edge_list))
{11, 13} directed simple Int64 graph

julia> strongly_connected_components_kosaraju(g)
4-element Array{Array{Int64,1},1}:
 [11, 10]    
 [2, 3, 4, 1]
 [6, 7, 5]   
 [9, 8]      
Graphs.weakly_connected_componentsFunction
weakly_connected_components(g)

Return the weakly connected components of the graph g. This is equivalent to the connected components of the undirected equivalent of g. For undirected graphs this is equivalent to the connected_components of g.

Examples

julia> g = SimpleDiGraph([0 1 0; 1 0 1; 0 0 0]);

julia> weakly_connected_components(g)
1-element Array{Array{Int64,1},1}:
 [1, 2, 3]
Graphs.has_self_loopsFunction
has_self_loops(g)

Return true if g has any self loops.

Examples

julia> using Graphs

julia> g = SimpleGraph(2);

julia> add_edge!(g, 1, 2);

julia> has_self_loops(g)
false

julia> add_edge!(g, 1, 1);

julia> has_self_loops(g)
true
Graphs.attracting_componentsFunction
attracting_components(g)

Return a vector of vectors of integers representing lists of attracting components in the directed graph g.

The attracting components are a subset of the strongly connected components in which the components do not have any leaving edges.

Examples

julia> g = SimpleDiGraph([0 1 0 0 0; 0 0 1 0 0; 1 0 0 1 0; 0 0 0 0 1; 0 0 0 1 0])
{5, 6} directed simple Int64 graph

julia> strongly_connected_components(g)
2-element Array{Array{Int64,1},1}:
 [4, 5]
 [1, 2, 3]

julia> attracting_components(g)
1-element Array{Array{Int64,1},1}:
 [4, 5]
Graphs.is_bipartiteFunction
is_bipartite(g)

Return true if graph g is bipartite.

Examples

julia> using Graphs

julia> g = SimpleGraph(3);

julia> add_edge!(g, 1, 2);

julia> add_edge!(g, 2, 3);

julia> is_bipartite(g)
true

julia> add_edge!(g, 1, 3);

julia> is_bipartite(g)
false
Graphs.bipartite_mapFunction
bipartite_map(g) -> Vector{UInt8}

For a bipartite graph g, return a vector c of size $|V|$ containing the assignment of each vertex to one of the two sets ($c_i == 1$ or $c_i == 2$). If g is not bipartite, return an empty vector.

Implementation Notes

Note that an empty vector does not necessarily indicate non-bipartiteness. An empty graph will return an empty vector but is bipartite.

Examples

julia> using Graphs

julia> g = SimpleGraph(3);

julia> bipartite_map(g)
3-element Array{UInt8,1}:
 0x01
 0x01
 0x01

julia> add_vertices!(g, 3);

julia> add_edge!(g, 1, 2);

julia> add_edge!(g, 2, 3);

julia> bipartite_map(g)
3-element Array{UInt8,1}:
 0x01
 0x02
 0x01
Graphs.biconnected_componentsFunction
biconnected_components(g) -> Vector{Vector{Edge{eltype(g)}}}

Compute the biconnected components of an undirected graph gand return a vector of vectors containing each biconnected component.

Performance: Time complexity is $\mathcal{O}(|V|)$.

Examples

julia> using Graphs

julia> biconnected_components(star_graph(5))
4-element Array{Array{Graphs.SimpleGraphs.SimpleEdge,1},1}:
 [Edge 1 => 3]
 [Edge 1 => 4]
 [Edge 1 => 5]
 [Edge 1 => 2]

julia> biconnected_components(cycle_graph(5))
1-element Array{Array{Graphs.SimpleGraphs.SimpleEdge,1},1}:
 [Edge 1 => 5, Edge 4 => 5, Edge 3 => 4, Edge 2 => 3, Edge 1 => 2]
Graphs.condensationFunction
condensation(g[, scc])

Return the condensation graph of the strongly connected components scc in the directed graph g. If scc is missing, generate the strongly connected components first.

Examples

julia> g = SimpleDiGraph([0 1 0 0 0; 0 0 1 0 0; 1 0 0 1 0; 0 0 0 0 1; 0 0 0 1 0])
{5, 6} directed simple Int64 graph

julia> strongly_connected_components(g)
2-element Array{Array{Int64,1},1}:
 [4, 5]
 [1, 2, 3]

julia> foreach(println, edges(condensation(g)))
Edge 2 => 1
Graphs.neighborhoodFunction
neighborhood(g, v, d, distmx=weights(g))

Return a vector of each vertex in g at a geodesic distance less than or equal to d, where distances may be specified by distmx.

Optional Arguments

  • dir=:out: If g is directed, this argument specifies the edge direction

with respect to v of the edges to be considered. Possible values: :in or :out.

Examples

julia> g = SimpleDiGraph([0 1 0 0 0; 0 0 1 0 0; 1 0 0 1 0; 0 0 0 0 1; 0 0 0 1 0]);

julia> neighborhood(g, 1, 2)
3-element Array{Int64,1}:
 1
 2
 3

julia> neighborhood(g, 1, 3)
4-element Array{Int64,1}:
 1
 2
 3
 4

julia> neighborhood(g, 1, 3, [0 1 0 0 0; 0 0 1 0 0; 1 0 0 0.25 0; 0 0 0 0 0.25; 0 0 0 0.25 0])
5-element Array{Int64,1}:
 1
 2
 3
 4
 5
Graphs.neighborhood_distsFunction
neighborhood_dists(g, v, d, distmx=weights(g))

Return a a vector of tuples representing each vertex which is at a geodesic distance less than or equal to d, along with its distance from v. Non-negative distances may be specified by distmx.

Optional Arguments

  • dir=:out: If g is directed, this argument specifies the edge direction

with respect to v of the edges to be considered. Possible values: :in or :out.

Examples

julia> g = SimpleDiGraph([0 1 0 0 0; 0 0 1 0 0; 1 0 0 1 0; 0 0 0 0 1; 0 0 0 1 0]);

julia> neighborhood_dists(g, 1, 3)
4-element Array{Tuple{Int64,Int64},1}:
 (1, 0)
 (2, 1)
 (3, 2)
 (4, 3)

julia> neighborhood_dists(g, 1, 3, [0 1 0 0 0; 0 0 1 0 0; 1 0 0 0.25 0; 0 0 0 0 0.25; 0 0 0 0.25 0])
5-element Array{Tuple{Int64,Float64},1}:
 (1, 0.0)
 (2, 1.0)
 (3, 2.0)
 (4, 2.25)
 (5, 2.5)

julia> neighborhood_dists(g, 4, 3)
2-element Array{Tuple{Int64,Int64},1}:
 (4, 0)
 (5, 1)

julia> neighborhood_dists(g, 4, 3, dir=:in)
5-element Array{Tuple{Int64,Int64},1}:
 (4, 0)
 (3, 1)
 (5, 1)
 (2, 2)
 (1, 3)
Graphs.articulationFunction
articulation(g)

Compute the articulation points of a connected graph g and return an array containing all cut vertices.

Examples

julia> using Graphs

julia> articulation(star_graph(5))
1-element Array{Int64,1}:
 1

julia> articulation(path_graph(5))
3-element Array{Int64,1}:
 2
 3
 4
Graphs.bridgesFunction
bridges(g)

Compute the bridges of a connected graph g and return an array containing all bridges, i.e edges whose deletion increases the number of connected components of the graph.

Examples

julia> using Graphs

julia> bridges(star_graph(5))
8-element Array{Graphs.SimpleGraphs.SimpleEdge{Int64},1}:
 Edge 1 => 2
 Edge 1 => 3
 Edge 1 => 4
 Edge 1 => 5

julia> bridges(path_graph(5))
8-element Array{Graphs.SimpleGraphs.SimpleEdge{Int64},1}:
 Edge 4 => 5
 Edge 3 => 4
 Edge 2 => 3
 Edge 1 => 2
Graphs.periodFunction
period(g)

Return the (common) period for all vertices in a strongly connected directed graph. Will throw an error if the graph is not strongly connected.

Examples

julia> g = SimpleDiGraph([0 1 0; 0 0 1; 1 0 0]);

julia> period(g)
3
Graphs.isgraphicalFunction
isgraphical(degs)

Return true if the degree sequence degs is graphical. A sequence of integers is called graphical, if there exists a graph where the degrees of its vertices form that same sequence.

Performance

Time complexity: $\mathcal{O}(|degs|*\log(|degs|))$.

Implementation Notes

According to Erdös-Gallai theorem, a degree sequence $\{d_1, ...,d_n\}$ (sorted in descending order) is graphic iff the sum of vertex degrees is even and the sequence obeys the property -

\[\sum_{i=1}^{r} d_i \leq r(r-1) + \sum_{i=r+1}^n min(r,d_i)\]

for each integer r <= n-1

Cycle Detection

In graph theory, a cycle is defined to be a path that starts from some vertex v and ends up at v.

Graphs.is_cyclicFunction
is_cyclic(g)

Return true if graph g contains a cycle.

Implementation Notes

Uses DFS.

Graphs.maxsimplecyclesFunction
maxsimplecycles(dg::::IsDirected, byscc::Bool = true)

Compute the theoretical maximum number of cycles in the directed graph dg.

The computation can be performed assuming the graph is complete or taking into account the decomposition in strongly connected components (byscc parameter).

Performance

A more efficient version is possible.

References

Graphs.simplecyclesFunction
simplecycles(dg::::IsDirected)

Compute and return all cycles of the given directed graph using Johnson's algorithm.

Performance

The number of cycles grows more than exponentially with the number of vertices, you might want to use the algorithm with a ceiling – simplecycles_iter – on large directed graphs (slightly slower). If you want to have an idea of the possible number of cycles, look at function maxsimplecycles(dg::DiGraph, byscc::Bool = true). If you only need short cycles of a limited length, simplecycles_limited_length can be more efficient.

References

Examples

julia> simplecycles(complete_digraph(3))
5-element Array{Array{Int64,1},1}:
 [1, 2]
 [1, 2, 3]
 [1, 3]
 [1, 3, 2]
 [2, 3]
Graphs.simplecycles_iterFunction
simplecycles_iter(dg::DiGraph, ceiling = 10^6)

Search all cycles of the given directed graph, using Johnson's algorithm, up to the ceiling (to avoid memory overload).

Implementation Notes

If the graph is small, the ceiling will not be reached and simplecycles(dg::DiGraph) is more efficient. It avoids the overhead of the counting and testing if the ceiling is reached. It returns all the cycles of the directed graph if the ceiling is not reached, a subset of them otherwise.

To get an idea of the possible number of cycles, use function `maxsimplecycles(dg::DiGraph, byscc::Bool = true) on the directed graph.

References

Graphs.simplecycles_hawick_jamesFunction
simplecycles_hawick_james(g)

Find circuits (including self-loops) in g using the algorithm of Hawick & James.

References

  • Hawick & James, "Enumerating Circuits and Loops in Graphs with Self-Arcs and Multiple-Arcs", 2008
Graphs.simplecyclescountFunction
simplecyclescount(dg::DiGraph, ceiling = 10^6)

Count the number of cycles in a directed graph, using Johnson's algorithm. Return the minimum of the ceiling and the number of cycles.

Implementation Notes

The ceiling is here to avoid memory overload if there are a lot of cycles in the graph. Default value is 10^6, but it can be higher or lower. You can use the function maxsimplecycles(dg::DiGraph, byscc::Bool = true) to get an idea of the theoretical maximum number or cycles.

References

Examples

julia> simplecyclescount(complete_digraph(6))
409
Graphs.simplecycleslengthFunction
simplecycleslength(dg::DiGraph, ceiling = 10^6)

Search all cycles of the given directed graph, using Johnson's algorithm, and return a tuple representing the cycle length and the number of cycles.

Implementation Notes

To get an idea of the possible number of cycles, using function maxsimplecycles(dg::DiGraph, byscc::Bool = true) on the directed graph.

If the ceiling is reached (ncycles = ceiling), the output is only a subset of the cycles lengths.

References

Examples

julia> simplecycleslength(complete_digraph(16))
([0, 1, 1, 1, 1, 1, 2, 10, 73, 511, 3066, 15329, 61313, 183939, 367876, 367876], 1000000)

julia> simplecycleslength(wheel_digraph(16))
([0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0], 1)
Graphs.simplecycles_limited_lengthFunction
simplecycles_limited_length(g, n, ceiling=10^6)

Compute and return at most ceiling cycles of length at most n of the given graph. Both directed and undirected graphs are supported.

Performance

The number of cycles grows very fast with the number of vertices and the allowed length of the cycles. This function is intended for finding short cycles. If you want to find cycles of any length in a directed graph, simplecycles or simplecycles_iter may be more efficient.

Graphs.karp_minimum_cycle_meanFunction
karp_minimum_cycle_mean(g[, distmx])

Return minimum cycle mean of the directed graph g with optional edge weights contained in distmx.

References

Minimum Spanning Trees (MST) Algorithms

A Minimum Spanning Tree (MST) is a subset of the edges of a connected, edge-weighted (un)directed graph that connects all the vertices together, without any cycles and with the minimum possible total edge weight.

Graphs.kruskal_mstFunction
kruskal_mst(g, distmx=weights(g); minimize=true)

Return a vector of edges representing the minimum (by default) spanning tree of a connected, undirected graph g with optional distance matrix distmx using Kruskal's algorithm.

Optional Arguments

  • minimize=true: if set to false, calculate the maximum spanning tree.
Graphs.prim_mstFunction
prim_mst(g, distmx=weights(g))

Return a vector of edges representing the minimum spanning tree of a connected, undirected graph g with optional distance matrix distmx using Prim's algorithm. Return a vector of edges.

Shortest-Path Algorithms

General properties of shortest path algorithms

  • The distance from a vertex to itself is always 0.
  • The distance between two vertices with no connecting edge is always Inf.
Graphs.a_starFunction
a_star(g, s, t[, distmx][, heuristic])

Return a vector of edges comprising the shortest path between vertices s and t using the A* search algorithm. An optional heuristic function and edge distance matrix may be supplied. If missing, the distance matrix is set to Graphs.DefaultDistance and the heuristic is set to n -> 0.

Graphs.dijkstra_shortest_pathsFunction
dijkstra_shortest_paths(g, srcs, distmx=weights(g));

Perform Dijkstra's algorithm on a graph, computing shortest distances between srcs and all other vertices. Return a Graphs.DijkstraState that contains various traversal information.

Optional Arguments

predecessors of a given vertex.

Performance

If using a sparse matrix for distmx, you may achieve better performance by passing in a transpose of its sparse transpose. That is, assuming D is the sparse distance matrix:

D = transpose(sparse(transpose(D)))

Be aware that realizing the sparse transpose of D incurs a heavy one-time penalty, so this strategy should only be used when multiple calls to dijkstra_shortest_paths with the distance matrix are planned.

Examples

julia> using Graphs

julia> ds = dijkstra_shortest_paths(cycle_graph(5), 2);

julia> ds.dists
5-element Array{Int64,1}:
 1
 0
 1
 2
 2

julia> ds = dijkstra_shortest_paths(path_graph(5), 2);

julia> ds.dists
5-element Array{Int64,1}:
 1
 0
 1
 2
 3
Graphs.desopo_pape_shortest_pathsFunction
desopo_pape_shortest_paths(g, src, distmx=weights(g))

Compute shortest paths between a source src and all other nodes in graph g using the D'Esopo-Pape algorithm. Return a Graphs.DEsopoPapeState with relevant traversal information.

Examples

julia> using Graphs

julia> ds = desopo_pape_shortest_paths(cycle_graph(5), 2);

julia> ds.dists
5-element Array{Int64,1}:
 1
 0
 1
 2
 2

julia> ds = desopo_pape_shortest_paths(path_graph(5), 2);

julia> ds.dists
5-element Array{Int64,1}:
 1
 0
 1
 2
 3
Graphs.yen_k_shortest_pathsFunction
yen_k_shortest_paths(g, source, target, distmx=weights(g), K=1; maxdist=Inf);

Perform Yen's algorithm on a graph, computing k-shortest distances between source and target other vertices. Return a YenState that contains distances and paths.

Graphs.spfa_shortest_pathsFunction
spfa_shortest_paths(g, s, distmx=weights(g))

Compute shortest paths between a source s and all other nodes in graph g using the Shortest Path Faster Algorithm.

Examples

julia> g = complete_graph(3);

julia> d = [1 -3 1; -3 1 1; 1 1 1];

julia> spfa_shortest_paths(g, 1, d)

ERROR: Graphs.NegativeCycleError()

julia> g = complete_graph(4);

julia> d = [1 1 -1 1; 1 1 -1 1; 1 1 1 1; 1 1 1 1];

julia> spfa_shortest_paths(gx, 1, d)

4-element Array{Int64,1}:
  0
  0
 -1
  0

Path discovery / enumeration

Graphs.gdistancesFunction
gdistances(g, source; sort_alg=QuickSort)

Return a vector filled with the geodesic distances of vertices in g from source. If source is a collection of vertices each element should be unique. For vertices in disconnected components the default distance is typemax(T).

An optional sorting algorithm may be specified (see Performance section).

Performance

gdistances uses QuickSort internally for its default sorting algorithm, since it performs the best of the algorithms built into Julia Base. However, passing a RadixSort (available via SortingAlgorithms.jl) will provide significant performance improvements on larger graphs.

Graphs.gdistances!Function
gdistances!(g, source, dists; sort_alg=QuickSort)

Fill dists with the geodesic distances of vertices in g from source vertex (or collection of vertices) source. dists should be a vector of length nv(g) filled with typemax(T). Return dists.

For vertices in disconnected components the default distance is typemax(T).

An optional sorting algorithm may be specified (see Performance section).

Performance

gdistances uses QuickSort internally for its default sorting algorithm, since it performs the best of the algorithms built into Julia Base. However, passing a RadixSort (available via SortingAlgorithms.jl) will provide significant performance improvements on larger graphs.

Graphs.enumerate_pathsFunction
enumerate_paths(state[, vs])

Given a path state state of type AbstractPathState, return a vector (indexed by vertex) of the paths between the source vertex used to compute the path state and a single destination vertex, a list of destination vertices, or the entire graph. For multiple destination vertices, each path is represented by a vector of vertices on the path between the source and the destination. Nonexistent paths will be indicated by an empty vector. For single destinations, the path is represented by a single vector of vertices, and will be length 0 if the path does not exist.

Implementation Notes

For Floyd-Warshall path states, please note that the output is a bit different, since this algorithm calculates all shortest paths for all pairs of vertices: enumerate_paths(state) will return a vector (indexed by source vertex) of vectors (indexed by destination vertex) of paths. enumerate_paths(state, v) will return a vector (indexed by destination vertex) of paths from source v to all other vertices. In addition, enumerate_paths(state, v, d) will return a vector representing the path from vertex v to vertex d.

Path States

All path states derive from

Graphs.AbstractPathStateType
AbstractPathState

An abstract type that provides information from shortest paths calculations.

The dijkstra_shortest_paths, floyd_warshall_shortest_paths, bellman_ford_shortest_paths, and yen_shortest_paths functions return states that contain various information about the graph learned during traversal.

Graphs.BellmanFordStateType
BellmanFordState{T, U}

An AbstractPathState designed for Bellman-Ford shortest-paths calculations.

Graphs.YenStateType
struct YenState{T, U}

Designed for yen k-shortest-paths calculations.

The above state types (with the exception of YenState) have the following common information, accessible via the type:

.dists Holds a vector of distances computed, indexed by source vertex.

.parents Holds a vector of parents of each vertex on the paths. The parent of a source vertex is always 0.

(YenState substitutes .paths for .parents.)

In addition, the following information may be populated with the appropriate arguments to dijkstra_shortest_paths:

.predecessors Holds a vector, indexed by vertex, of all the predecessors discovered during shortest-path calculations. This keeps track of all parents when there are multiple shortest paths available from the source.

.pathcounts Holds a vector, indexed by vertex, of the path counts discovered during traversal. This equals the length of each subvector in the .predecessors output above.