Skip to main content
Ctrl+K

cugraph-docs 25.06.00 documentation

  • cuGraph Introduction
  • nx-cugraph
  • Installation
  • Tutorials
  • Graph Support
    • WholeGraph
    • References
    • Developer Resources
    • API Reference
  • GitHub
  • Twitter
Home
cugraph
cucimcudf-javacudfcugraphcumlcuprojcuspatialcuvscuxfilterdask-cudadask-cudfkvikiolibcudflibcumllibcuprojlibcuspatiallibkvikiolibrmmlibucxxraftrapids-cmakermm
nightly (25.06)
nightly (25.06)stable (25.04)legacy (25.02)
  • cuGraph Introduction
  • nx-cugraph
  • Installation
  • Tutorials
  • Graph Support
  • WholeGraph
  • References
  • Developer Resources
  • API Reference
  • GitHub
  • Twitter

Section Navigation

Core Graph API Documentation

  • cugraph API Reference
    • Graph Classes
      • cugraph.Graph
      • cugraph.MultiGraph
      • cugraph.Graph.from_cudf_adjlist
      • cugraph.Graph.from_cudf_edgelist
      • cugraph.Graph.from_dask_cudf_edgelist
      • cugraph.Graph.from_pandas_adjacency
      • cugraph.Graph.from_pandas_edgelist
      • cugraph.Graph.from_numpy_array
      • cugraph.Graph.from_numpy_matrix
      • cugraph.Graph.add_internal_vertex_id
      • cugraph.Graph.add_nodes_from
      • cugraph.Graph.clear
      • cugraph.Graph.unrenumber
      • cugraph.symmetrize
      • cugraph.symmetrize_ddf
      • cugraph.symmetrize_df
      • cugraph.from_adjlist
      • cugraph.from_cudf_edgelist
      • cugraph.from_edgelist
      • cugraph.from_numpy_array
      • cugraph.from_numpy_matrix
      • cugraph.from_pandas_adjacency
      • cugraph.from_pandas_edgelist
      • cugraph.to_numpy_array
      • cugraph.to_numpy_matrix
      • cugraph.to_pandas_adjacency
      • cugraph.to_pandas_edgelist
      • cugraph.structure.NumberMap
      • cugraph.structure.NumberMap.from_internal_vertex_id
      • cugraph.structure.NumberMap.to_internal_vertex_id
      • cugraph.structure.NumberMap.add_internal_vertex_id
      • cugraph.structure.NumberMap.compute_vals
      • cugraph.structure.NumberMap.compute_vals_types
      • cugraph.structure.NumberMap.generate_unused_column_name
      • cugraph.structure.NumberMap.renumber
      • cugraph.structure.NumberMap.renumber_and_segment
      • cugraph.structure.NumberMap.set_renumbered_col_names
      • cugraph.structure.NumberMap.unrenumber
      • cugraph.structure.NumberMap.vertex_column_size
      • cugraph.hypergraph
    • Graph Implementation
      • cugraph.structure.graph_implementation.simpleGraphImpl.view_edge_list
      • cugraph.structure.graph_implementation.simpleGraphImpl.delete_edge_list
      • cugraph.structure.graph_implementation.simpleGraphImpl.view_adj_list
      • cugraph.structure.graph_implementation.simpleGraphImpl.view_transposed_adj_list
      • cugraph.structure.graph_implementation.simpleGraphImpl.delete_adj_list
      • cugraph.structure.graph_implementation.simpleGraphImpl.enable_batch
      • cugraph.structure.graph_implementation.simpleGraphImpl.get_two_hop_neighbors
      • cugraph.structure.graph_implementation.simpleGraphImpl.number_of_vertices
      • cugraph.structure.graph_implementation.simpleGraphImpl.number_of_nodes
      • cugraph.structure.graph_implementation.simpleGraphImpl.number_of_edges
      • cugraph.structure.graph_implementation.simpleGraphImpl.in_degree
      • cugraph.structure.graph_implementation.simpleGraphImpl.out_degree
      • cugraph.structure.graph_implementation.simpleGraphImpl.degree
      • cugraph.structure.graph_implementation.simpleGraphImpl.degrees
      • cugraph.structure.graph_implementation.simpleGraphImpl.has_edge
      • cugraph.structure.graph_implementation.simpleGraphImpl.has_node
      • cugraph.structure.graph_implementation.simpleGraphImpl.has_self_loop
      • cugraph.structure.graph_implementation.simpleGraphImpl.edges
      • cugraph.structure.graph_implementation.simpleGraphImpl.nodes
      • cugraph.structure.graph_implementation.simpleGraphImpl.neighbors
      • cugraph.structure.graph_implementation.simpleGraphImpl.vertex_column_size
    • Property Graph
      • cugraph.experimental.PropertySelection
      • cugraph.experimental.PropertyGraph
      • cugraph.experimental.PropertyGraph.add_edge_data
      • cugraph.experimental.PropertyGraph.add_vertex_data
      • cugraph.experimental.PropertyGraph.annotate_dataframe
      • cugraph.experimental.PropertyGraph.edge_props_to_graph
      • cugraph.experimental.PropertyGraph.extract_subgraph
      • cugraph.experimental.PropertyGraph.get_edge_data
      • cugraph.experimental.PropertyGraph.get_num_edges
      • cugraph.experimental.PropertyGraph.get_num_vertices
      • cugraph.experimental.PropertyGraph.get_vertex_data
      • cugraph.experimental.PropertyGraph.get_vertices
      • cugraph.experimental.PropertyGraph.has_duplicate_edges
      • cugraph.experimental.PropertyGraph.is_multigraph
      • cugraph.experimental.PropertyGraph.renumber_edges_by_type
      • cugraph.experimental.PropertyGraph.renumber_vertices_by_type
      • cugraph.experimental.PropertyGraph.select_edges
      • cugraph.experimental.PropertyGraph.select_vertices
    • Centrality
      • cugraph.centrality.betweenness_centrality
      • cugraph.centrality.edge_betweenness_centrality
      • cugraph.dask.centrality.betweenness_centrality
      • cugraph.centrality.katz_centrality
      • cugraph.dask.centrality.katz_centrality.katz_centrality
      • cugraph.centrality.degree_centrality
      • cugraph.centrality.eigenvector_centrality
      • cugraph.dask.centrality.eigenvector_centrality.eigenvector_centrality
    • Community
      • cugraph.batched_ego_graphs
      • cugraph.ego_graph
      • cugraph.ecg
      • cugraph.k_truss
      • cugraph.ktruss_subgraph
      • cugraph.leiden
      • cugraph.louvain
      • cugraph.dask.community.louvain.louvain
      • cugraph.analyzeClustering_edge_cut
      • cugraph.analyzeClustering_modularity
      • cugraph.analyzeClustering_ratio_cut
      • cugraph.spectralBalancedCutClustering
      • cugraph.spectralModularityMaximizationClustering
      • cugraph.subgraph
      • cugraph.triangle_count
    • Components
      • cugraph.connected_components
      • cugraph.strongly_connected_components
      • cugraph.weakly_connected_components
      • cugraph.dask.components.connectivity.weakly_connected_components
    • Cores
      • cugraph.core_number
      • cugraph.k_core
    • Layout
      • cugraph.force_atlas2
    • Linear Assignment
      • cugraph.hungarian
      • cugraph.dense_hungarian
    • Link Analysis
      • cugraph.hits
      • cugraph.dask.link_analysis.hits.hits
      • cugraph.pagerank
      • cugraph.dask.link_analysis.pagerank.pagerank
    • Link Prediction
      • cugraph.jaccard
      • cugraph.jaccard_coefficient
      • cugraph.overlap
      • cugraph.overlap_coefficient
      • cugraph.sorensen
      • cugraph.sorensen_coefficient
    • Sampling
      • cugraph.random_walks
      • cugraph.ego_graph
      • cugraph.uniform_neighbor_sample
      • cugraph.node2vec
    • Traversal
      • cugraph.bfs
      • cugraph.bfs_edges
      • cugraph.dask.traversal.bfs.bfs
      • cugraph.filter_unreachable
      • cugraph.shortest_path
      • cugraph.shortest_path_length
      • cugraph.sssp
      • cugraph.dask.traversal.sssp.sssp
    • Tree
      • cugraph.tree.minimum_spanning_tree.minimum_spanning_tree
      • cugraph.tree.minimum_spanning_tree.maximum_spanning_tree
    • Generators
      • cugraph.generators.rmat
    • DASK MG Helper functions
      • cugraph.dask.comms.comms.initialize
      • cugraph.dask.comms.comms.destroy
      • cugraph.dask.comms.comms.is_initialized
      • cugraph.dask.comms.comms.get_comms
      • cugraph.dask.comms.comms.get_workers
      • cugraph.dask.comms.comms.get_session_id
      • cugraph.dask.comms.comms.get_2D_partition
      • cugraph.dask.comms.comms.get_default_handle
      • cugraph.dask.comms.comms.get_handle
      • cugraph.dask.comms.comms.get_worker_id
      • cugraph.dask.common.read_utils.get_chunksize
    • Multi-GPU with cuGraph
  • pylibcugraph API reference
    • pylibcugraph.eigenvector_centrality
    • pylibcugraph.katz_centrality
    • pylibcugraph.strongly_connected_components
    • pylibcugraph.weakly_connected_components
    • pylibcugraph.pagerank
    • pylibcugraph.hits
    • pylibcugraph.node2vec
    • pylibcugraph.bfs
    • pylibcugraph.sssp
  • cuGraph C API documentation
    • Centrality
    • Community
    • Core
    • Components
    • Sampling
    • Similarity
    • Traversal
  • cuGraph C++ API
    • Algorithmns
      • Centrality
      • Community
      • Sampling
      • Similarity
      • Traversal
      • Linear
      • Link Analysis
      • Layout
      • Tree
      • Utility Functions
    • Graph Functions
    • Graph Generators
    • Legacy Graph Functions
    • Sampling Functions
    • Collection Wrappers
    • Low Level cuGraph C++ API

Graph Neural Networks API Documentation

  • cugraph-dgl API Reference
    • cugraph_dgl.convert.cugraph_storage_from_heterograph
    • cugraph_dgl.cugraph_storage.CuGraphStorage
  • cugraph-pyg API Reference
    • cugraph_pyg.data.graph_store.GraphStore
    • cugraph_pyg.data.feature_store.TensorDictFeatureStore
    • cugraph_pyg.data.feature_store.WholeFeatureStore
    • cugraph_pyg.loader.node_loader.NodeLoader
    • cugraph_pyg.loader.neighbor_loader.NeighborLoader
    • cugraph_pyg.sampler.sampler.BaseSampler
    • cugraph_pyg.sampler.sampler.SampleReader
    • cugraph_pyg.sampler.sampler.HomogeneousSampleReader
    • cugraph_pyg.sampler.sampler.SampleIterator

Additional Graph Packages API Documentation

  • cugraph-service API Reference
    • cugraph-service-client API Reference
    • cugraph-service-server API Reference
  • API Reference
  • cuGraph C++ API
  • Algorithmns
  • Link Analysis

Link Analysis#

template<typename vertex_t, typename edge_t, typename weight_t, typename result_t, bool multi_gpu>
void pagerank(raft::handle_t const &handle, graph_view_t<vertex_t, edge_t, true, multi_gpu> const &graph_view, std::optional<edge_property_view_t<edge_t, weight_t const*>> edge_weight_view, std::optional<weight_t const*> precomputed_vertex_out_weight_sums, std::optional<vertex_t const*> personalization_vertices, std::optional<result_t const*> personalization_values, std::optional<vertex_t> personalization_vector_size, result_t *pageranks, result_t alpha, result_t epsilon, size_t max_iterations = 500, bool has_initial_guess = false, bool do_expensive_check = false)#

Compute PageRank scores.

Deprecated:

This API will be deprecated to replaced by the new version below that returns metadata about the algorithm.

This function computes general (if personalization_vertices is nullptr) or personalized (if personalization_vertices is not nullptr.) PageRank scores.

Throws:

cugraph::logic_error – on erroneous input arguments or if fails to converge before max_iterations.

Template Parameters:
  • vertex_t – Type of vertex identifiers. Needs to be an integral type.

  • edge_t – Type of edge identifiers. Needs to be an integral type.

  • weight_t – Type of edge weights. Needs to be a floating point type.

  • result_t – Type of PageRank scores.

  • multi_gpu – Flag indicating whether template instantiation should target single-GPU (false) or multi-GPU (true).

Parameters:
  • handle – RAFT handle object to encapsulate resources (e.g. CUDA stream, communicator, and handles to various CUDA libraries) to run graph algorithms.

  • graph_view – Graph view object.

  • edge_weight_view – Optional view object holding edge weights for graph_view. If edge_weight_view.has_value() == false, edge weights are assumed to be 1.0.

  • precomputed_vertex_out_weight_sums – Pointer to an array storing sums of out-going edge weights for the vertices (for re-use) or std::nullopt. If std::nullopt, these values are freshly computed. Computing these values outside this function reduces the number of memory allocations/deallocations and computing if a user repeatedly computes PageRank scores using the same graph with different personalization vectors.

  • personalization_vertices – Pointer to an array storing personalization vertex identifiers (compute personalized PageRank) or std::nullopt (compute general PageRank).

  • personalization_values – Pointer to an array storing personalization values for the vertices in the personalization set. Relevant only if personalization_vertices is not std::nullopt.

  • personalization_vector_size – Size of the personalization set. If @personalization_vertices is not std::nullopt, the sizes of the arrays pointed by personalization_vertices and personalization_values should be personalization_vector_size.

  • pageranks – Pointer to the output PageRank score array.

  • alpha – PageRank damping factor.

  • epsilon – Error tolerance to check convergence. Convergence is assumed if the sum of the differences in PageRank values between two consecutive iterations is less than the number of vertices in the graph multiplied by epsilon.

  • max_iterations – Maximum number of PageRank iterations.

  • has_initial_guess – If set to true, values in the PageRank output array (pointed by pageranks) is used as initial PageRank values. If false, initial PageRank values are set to 1.0 divided by the number of vertices in the graph.

  • do_expensive_check – A flag to run expensive checks for input arguments (if set to true).

template<typename vertex_t, typename edge_t, typename weight_t, typename result_t, bool multi_gpu>
std::tuple<rmm::device_uvector<result_t>, centrality_algorithm_metadata_t> pagerank(raft::handle_t const &handle, graph_view_t<vertex_t, edge_t, true, multi_gpu> const &graph_view, std::optional<edge_property_view_t<edge_t, weight_t const*>> edge_weight_view, std::optional<raft::device_span<weight_t const>> precomputed_vertex_out_weight_sums, std::optional<std::tuple<raft::device_span<vertex_t const>, raft::device_span<result_t const>>> personalization, std::optional<raft::device_span<result_t const>> initial_pageranks, result_t alpha, result_t epsilon, size_t max_iterations = 500, bool do_expensive_check = false)#

Compute PageRank scores.

.*

This function computes general (if personalization_vertices is nullptr) or personalized (if personalization_vertices is not nullptr.) PageRank scores.

Throws:

cugraph::logic_error – on erroneous input arguments or if fails to converge before max_iterations.

Template Parameters:
  • vertex_t – Type of vertex identifiers. Needs to be an integral type.

  • edge_t – Type of edge identifiers. Needs to be an integral type.

  • weight_t – Type of edge weights. Needs to be a floating point type.

  • result_t – Type of PageRank scores.

  • multi_gpu – Flag indicating whether template instantiation should target single-GPU (false) or multi-GPU (true).

Parameters:
  • handle – RAFT handle object to encapsulate resources (e.g. CUDA stream, communicator, and handles to various CUDA libraries) to run graph algorithms.

  • graph_view – Graph view object.

  • edge_weight_view – Optional view object holding edge weights for graph_view. If edge_weight_view.has_value() == false, edge weights are assumed to be 1.0.

  • precomputed_vertex_out_weight_sums – Pointer to an array storing sums of out-going edge weights for the vertices (for re-use) or std::nullopt. If std::nullopt, these values are freshly computed. Computing these values outside this function reduces the number of memory allocations/deallocations and computing if a user repeatedly computes PageRank scores using the same graph with different personalization vectors.

  • personalization – Optional tuple containing device spans of vertex identifiers and personalization values for the vertices (compute personalized PageRank) or std::nullopt (compute general PageRank).

  • initial_pageranks – Optional device span containing initial PageRank values. If specified this array will be used as the initial values and the PageRank values will be updated in place. If not specified then the initial values will be set to 1.0 divided by the number of vertices in the graph and the return value will contain an rmm::device_uvector containing the resulting PageRank values.

  • alpha – PageRank damping factor.

  • epsilon – Error tolerance to check convergence. Convergence is assumed if the sum of the differences in PageRank values between two consecutive iterations is less than the number of vertices in the graph multiplied by epsilon.

  • max_iterations – Maximum number of PageRank iterations.

  • do_expensive_check – A flag to run expensive checks for input arguments (if set to true).

Returns:

tuple containing the optional pagerank results (populated if initial_pageranks is set to std::nullopt) and a metadata structure with metadata indicating how many iterations were run and whether the algorithm converged or not.

template<typename vertex_t, typename edge_t, typename result_t, bool multi_gpu>
std::tuple<result_t, size_t> hits(raft::handle_t const &handle, graph_view_t<vertex_t, edge_t, true, multi_gpu> const &graph_view, result_t *hubs, result_t *authorities, result_t epsilon, size_t max_iterations, bool has_initial_hubs_guess, bool normalize, bool do_expensive_check)#

Compute HITS scores.

.*

This function computes HITS scores for the vertices of a graph

Throws:

cugraph::logic_error – on erroneous input arguments

Template Parameters:
  • vertex_t – Type of vertex identifiers. Needs to be an integral type.

  • edge_t – Type of edge identifiers. Needs to be an integral type.

  • multi_gpu – Flag indicating whether template instantiation should target single-GPU (false) or multi-GPU (true).

Parameters:
  • handle – RAFT handle object to encapsulate resources (e.g. CUDA stream, communicator, and handles to various CUDA libraries) to run graph algorithms.

  • graph_view – Graph view object.

  • hubs – Pointer to the input/output hub score array.

  • authorities – Pointer to the output authorities score array.

  • epsilon – Error tolerance to check convergence. Convergence is assumed if the sum of the differences in hub values between two consecutive iterations is less than epsilon

  • max_iterations – Maximum number of HITS iterations.

  • has_initial_guess – If set to true, values in the hubs output array (pointed by hubs) is used as initial hub values. If false, initial hub values are set to 1.0 divided by the number of vertices in the graph.

  • normalize – If set to true, final hub and authority scores are normalized (the L1-norm of the returned hub and authority score arrays is 1.0) before returning.

  • do_expensive_check – A flag to run expensive checks for input arguments (if set to true).

Returns:

std::tuple<result_t, size_t> A tuple of sum of the differences of hub scores of the last two iterations and the total number of iterations taken to reach the final result

previous

Linear

next

Layout

On this page
  • pagerank()
  • pagerank()
  • hits()

This Page

  • Show Source

© Copyright 2024, NVIDIA Corporation.

Created using Sphinx 8.1.3.

Built with the PyData Sphinx Theme 0.16.1.