Community#
- template<typename vertex_t, typename edge_t, typename weight_t, bool multi_gpu>
std::pair<size_t, weight_t> louvain(raft::handle_t const &handle, std::optional<std::reference_wrapper<raft::random::RngState>> rng_state, graph_view_t<vertex_t, edge_t, false, multi_gpu> const &graph_view, std::optional<edge_property_view_t<edge_t, weight_t const*>> edge_weight_view, vertex_t *clustering, size_t max_level = 100, weight_t threshold = weight_t{1e-7}, weight_t resolution = weight_t{1})#Louvain implementation.
Compute a clustering of the graph by maximizing modularity
Computed using the Louvain method described in:
VD Blondel, J-L Guillaume, R Lambiotte and E Lefebvre: Fast unfolding of community hierarchies in large networks, J Stat Mech P10008 (2008), http://arxiv.org/abs/0803.0476
- Throws:
cugraph::logic_error – when an error occurs.
- 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.
multi_gpu – Flag indicating whether template instantiation should target single-GPU (false)
- Parameters:
handle – [in] Library handle (RAFT). If a communicator is set in the handle,
rng_state – [in] The RngState instance holding pseudo-random number generator state.
graph_view – [in] Input graph view object.
edge_weight_view – [in] Optional view object holding edge weights for
graph_view
. If @pedge_weight_view.has_value() == false, edge weights are assumed to be 1.0.clustering – [out] Pointer to device array where the clustering should be stored
max_level – [in] (optional) maximum number of levels to run (default 100)
threshold – [in] (optional) threshold for convergence at each level (default 1e-7)
resolution – [in] (optional) The value of the resolution parameter to use. Called gamma in the modularity formula, this changes the size of the communities. Higher resolutions lead to more smaller communities, lower resolutions lead to fewer larger communities. (default 1)
- Returns:
a pair containing: 1) number of levels of the returned clustering 2) modularity of the returned clustering
- template<typename vertex_t, typename edge_t, typename weight_t, bool multi_gpu>
std::pair<std::unique_ptr<Dendrogram<vertex_t>>, weight_t> louvain(raft::handle_t const &handle, std::optional<std::reference_wrapper<raft::random::RngState>> rng_state, graph_view_t<vertex_t, edge_t, false, multi_gpu> const &graph_view, std::optional<edge_property_view_t<edge_t, weight_t const*>> edge_weight_view, size_t max_level = 100, weight_t threshold = weight_t{1e-7}, weight_t resolution = weight_t{1})#Louvain implementation, returning dendrogram.
Compute a clustering of the graph by maximizing modularity
Computed using the Louvain method described in:
VD Blondel, J-L Guillaume, R Lambiotte and E Lefebvre: Fast unfolding of community hierarchies in large networks, J Stat Mech P10008 (2008), http://arxiv.org/abs/0803.0476
- Throws:
cugraph::logic_error – when an error occurs.
- 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.
multi_gpu – Flag indicating whether template instantiation should target single-GPU (false)
- Parameters:
handle – [in] Library handle (RAFT). If a communicator is set in the handle,
rng_state – [in] The RngState instance holding pseudo-random number generator state.
graph_view – [in] Input graph view object.
edge_weight_view – [in] Optional view object holding edge weights for
graph_view
. If @pedge_weight_view.has_value() == false, edge weights are assumed to be 1.0.max_level – [in] (optional) maximum number of levels to run (default 100)
threshold – [in] (optional) threshold for convergence at each level (default 1e-7)
resolution – [in] (optional) The value of the resolution parameter to use. Called gamma in the modularity formula, this changes the size of the communities. Higher resolutions lead to more smaller communities, lower resolutions lead to fewer larger communities. (default 1)
- Returns:
a pair containing: 1) unique pointer to dendrogram 2) modularity of the returned clustering
- template<typename graph_view_t>
void flatten_dendrogram(raft::handle_t const &handle, graph_view_t const &graph_view, Dendrogram<typename graph_view_t::vertex_type> const &dendrogram, typename graph_view_t::vertex_type *clustering)#Flatten a Dendrogram at a particular level.
A Dendrogram represents a hierarchical clustering/partitioning of a graph. This function will flatten the hierarchical clustering into a label for each vertex representing the final cluster/partition to which it is assigned
- Throws:
cugraph::logic_error – when an error occurs.
- Template Parameters:
graph_view_t – Type of graph
- Parameters:
handle – [in] Library handle (RAFT). If a communicator is set in the handle,
graph – [in] input graph object
dendrogram – [in] input dendrogram object
clustering – [out] Pointer to device array where the clustering should be stored
- template<typename vertex_t, typename edge_t, typename weight_t, bool multi_gpu>
std::pair<std::unique_ptr<Dendrogram<vertex_t>>, weight_t> leiden(raft::handle_t const &handle, raft::random::RngState &rng_state, graph_view_t<vertex_t, edge_t, false, multi_gpu> const &graph_view, std::optional<edge_property_view_t<edge_t, weight_t const*>> edge_weight_view, size_t max_level = 100, weight_t resolution = weight_t{1}, weight_t theta = weight_t{1})#Leiden implementation.
Compute a clustering of the graph by maximizing modularity using the Leiden improvements to the Louvain method.
Computed using the Leiden method described in:
Traag, V. A., Waltman, L., & van Eck, N. J. (2019). From Louvain to Leiden: guaranteeing well-connected communities. Scientific reports, 9(1), 5233. doi: 10.1038/s41598-019-41695-z
- Throws:
cugraph::logic_error – when an error occurs.
- Template Parameters:
vertex_t – Type of vertex identifiers. Supported value : int (signed, 32-bit)
edge_t – Type of edge identifiers. Supported value : int (signed, 32-bit)
weight_t – Type of edge weights. Supported values : float or double.
- Parameters:
handle – RAFT handle object to encapsulate resources (e.g. CUDA stream, communicator, and handles to various CUDA libraries) to run graph algorithms.
rng_state – The RngState instance holding pseudo-random number generator state.
graph_view – Graph view object.
edge_weight_view – Optional view object holding edge weights for
graph_view
. Ifedge_weight_view.has_value()
== false, edge weights are assumed to be 1.0.max_level – [in] (optional) maximum number of levels to run (default 100)
resolution – [in] (optional) The value of the resolution parameter to use. Called gamma in the modularity formula, this changes the size of the communities. Higher resolutions lead to more smaller communities, lower resolutions lead to fewer larger communities. (default 1)
theta – [in] (optional) The value of the parameter to scale modularity gain in Leiden refinement phase. It is used to compute the probability of joining a random leiden community. Called theta in the Leiden algorithm.
- Returns:
a pair containing: 1) unique pointer to dendrogram 2) modularity of the returned clustering
- template<typename vertex_t, typename edge_t, typename weight_t, bool multi_gpu>
std::pair<size_t, weight_t> leiden(raft::handle_t const &handle, raft::random::RngState &rng_state, graph_view_t<vertex_t, edge_t, false, multi_gpu> const &graph_view, std::optional<edge_property_view_t<edge_t, weight_t const*>> edge_weight_view, vertex_t *clustering, size_t max_level = 100, weight_t resolution = weight_t{1}, weight_t theta = weight_t{1})#Leiden implementation.
.*
Compute a clustering of the graph by maximizing modularity using the Leiden improvements to the Louvain method.
Computed using the Leiden method described in:
Traag, V. A., Waltman, L., & van Eck, N. J. (2019). From Louvain to Leiden: guaranteeing well-connected communities. Scientific reports, 9(1), 5233. doi: 10.1038/s41598-019-41695-z
- Throws:
cugraph::logic_error – when an error occurs.
- Template Parameters:
vertex_t – Type of vertex identifiers. Supported value : int (signed, 32-bit)
edge_t – Type of edge identifiers. Supported value : int (signed, 32-bit)
weight_t – Type of edge weights. Supported values : float or double.
- Parameters:
handle – RAFT handle object to encapsulate resources (e.g. CUDA stream, communicator, and handles to various CUDA libraries) to run graph algorithms.
rng_state – The RngState instance holding pseudo-random number generator state.
graph_view – Graph view object.
edge_weight_view – Optional view object holding edge weights for
graph_view
. Ifedge_weight_view.has_value()
== false, edge weights are assumed to be 1.0.max_level – [in] (optional) maximum number of levels to run (default 100)
resolution – [in] (optional) The value of the resolution parameter to use. Called gamma in the modularity formula, this changes the size of the communities. Higher resolutions lead to more smaller communities, lower resolutions lead to fewer larger communities. (default 1)
theta – [in] (optional) The value of the parameter to scale modularity gain in Leiden refinement phase. It is used to compute the probability of joining a random leiden community. Called theta in the Leiden algorithm. communities. (default 1)
- Returns:
a pair containing: 1) number of levels of the returned clustering 2) modularity of the returned clustering
- template<typename vertex_t, typename edge_t, typename weight_t, bool multi_gpu>
std::tuple<rmm::device_uvector<vertex_t>, size_t, weight_t> ecg(raft::handle_t const &handle, raft::random::RngState &rng_state, graph_view_t<vertex_t, edge_t, false, multi_gpu> const &graph_view, std::optional<edge_property_view_t<edge_t, weight_t const*>> edge_weight_view, weight_t min_weight, size_t ensemble_size, size_t max_level = 100, weight_t threshold = weight_t{1e-7}, weight_t resolution = weight_t{1})#Computes the ecg clustering of the given graph.
.*
ECG runs truncated Louvain on an ensemble of permutations of the input graph, then uses the ensemble partitions to determine weights for the input graph. The final result is found by running full Louvain on the input graph using the determined weights. See https://arxiv.org/abs/1809.05578 for further information.
- Throws:
cugraph::logic_error – when an error occurs.
- 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.
multi_gpu – Flag indicating whether template instantiation should target single-GPU (false)
- Parameters:
handle – [in] Library handle (RAFT). If a communicator is set in the handle,
rng_state – [in] The RngState instance holding pseudo-random number generator state.
graph_view – [in] Input graph view object
edge_weight_view – [in] View object holding edge weights for
graph_view
.min_weight – [in] Minimum edge weight to use in the final call of the clustering algorithm if an edge does not appear in any of the ensemble runs.
ensemble_size – [in] The ensemble size parameter
max_level – [in] (optional) maximum number of levels to run (default 100)
threshold – [in] (optional) threshold for convergence at each level (default 1e-7)
resolution – [in] (optional) The value of the resolution parameter to use. Called gamma in the modularity formula, this changes the size of the communities. Higher resolutions lead to more smaller communities, lower resolutions lead to fewer larger communities. (default 1)
- Returns:
a tuple containing: 1) Device vector containing clustering result 2) number of levels of the returned clustering 3) modularity of the returned clustering
- template<typename VT, typename ET, typename WT>
void spectralModularityMaximization(legacy::GraphCSRView<VT, ET, WT> const &graph, VT n_clusters, VT n_eig_vects, WT evs_tolerance, int evs_max_iter, WT kmean_tolerance, int kmean_max_iter, VT *clustering)#Wrapper function for Nvgraph spectral modularity maximization algorithm.
- Throws:
cugraph::logic_error – when an error occurs.
- Template Parameters:
VT – Type of vertex identifiers. Supported value : int (signed, 32-bit)
ET – Type of edge identifiers. Supported value : int (signed, 32-bit)
WT – Type of edge weights. Supported values : float or double.
- Parameters:
graph – [in] input graph object (CSR)
num_clusters – [in] The desired number of clusters
num_eigen_vects – [in] The number of eigenvectors to use
evs_tolerance – [in] The tolerance to use for the eigenvalue solver
evs_max_iter – [in] The maximum number of iterations of the eigenvalue solver
kmean_tolerance – [in] The tolerance to use for the kmeans solver
kmean_max_iter – [in] The maximum number of iteration of the k-means solver
clustering – [out] Pointer to device memory where the resulting clustering will be stored
- template<typename VT, typename ET, typename WT>
void analyzeClustering_modularity(legacy::GraphCSRView<VT, ET, WT> const &graph, int n_clusters, VT const *clustering, WT *score)#Wrapper function for Nvgraph clustering modularity metric.
- Throws:
cugraph::logic_error – when an error occurs.
- Template Parameters:
VT – Type of vertex identifiers. Supported value : int (signed, 32-bit)
ET – Type of edge identifiers. Supported value : int (signed, 32-bit)
WT – Type of edge weights. Supported values : float or double.
- Parameters:
graph – [in] input graph object (CSR)
n_clusters – [in] Number of clusters in the clustering
clustering – [in] Pointer to device array containing the clustering to analyze
score – [out] Pointer to a float in which the result will be written
- template<typename VT, typename ET, typename WT>
void analyzeClustering_edge_cut(legacy::GraphCSRView<VT, ET, WT> const &graph, int n_clusters, VT const *clustering, WT *score)#Wrapper function for Nvgraph clustering edge cut metric.
- Throws:
cugraph::logic_error – when an error occurs.
- Template Parameters:
VT – Type of vertex identifiers. Supported value : int (signed, 32-bit)
ET – Type of edge identifiers. Supported value : int (signed, 32-bit)
WT – Type of edge weights. Supported values : float or double.
- Parameters:
graph – [in] input graph object (CSR)
n_clusters – [in] Number of clusters in the clustering
clustering – [in] Pointer to device array containing the clustering to analyze
score – [out] Pointer to a float in which the result will be written
- template<typename VT, typename ET, typename WT>
void analyzeClustering_ratio_cut(legacy::GraphCSRView<VT, ET, WT> const &graph, int n_clusters, VT const *clustering, WT *score)#Wrapper function for Nvgraph clustering ratio cut metric.
- Throws:
cugraph::logic_error – when an error occurs.
- Template Parameters:
VT – Type of vertex identifiers. Supported value : int (signed, 32-bit)
ET – Type of edge identifiers. Supported value : int (signed, 32-bit)
WT – Type of edge weights. Supported values : float or double.
- Parameters:
graph – [in] input graph object (CSR)
n_clusters – [in] Number of clusters in the clustering
clustering – [in] Pointer to device array containing the clustering to analyze
score – [out] Pointer to a float in which the result will be written
- template<typename vertex_t, typename edge_t, typename weight_t, bool multi_gpu>
std::tuple<rmm::device_uvector<vertex_t>, rmm::device_uvector<vertex_t>, std::optional<rmm::device_uvector<weight_t>>, rmm::device_uvector<size_t>> extract_ego(raft::handle_t const &handle, graph_view_t<vertex_t, edge_t, false, multi_gpu> const &graph_view, std::optional<edge_property_view_t<edge_t, weight_t const*>> edge_weight_view, vertex_t *source_vertex, vertex_t n_subgraphs, vertex_t radius)#returns induced EgoNet subgraph(s) of neighbors centered at nodes in source_vertex within a given radius.
.*
- Deprecated:
This algorithm will be deprecated to replaced by the new version that uses the raft::device_span.
- 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.
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. Must have at least one worker stream.
graph_view – Graph view object of, we extract induced egonet subgraphs from
graph_view
.edge_weight_view – Optional view object holding edge weights for
graph_view
.source_vertex – Pointer to egonet center vertices (size ==
n_subgraphs
).n_subgraphs – Number of induced EgoNet subgraphs to extract (ie. number of elements in
source_vertex
).radius – Include all neighbors of distance <= radius from
source_vertex
.- Returns:
Quadraplet of edge source vertices, edge destination vertices, edge weights (if
edge_weight_view.has_value()
== true), and edge offsets for each induced EgoNet subgraph.
- template<typename vertex_t, typename edge_t, typename weight_t, bool multi_gpu>
std::tuple<rmm::device_uvector<vertex_t>, rmm::device_uvector<vertex_t>, std::optional<rmm::device_uvector<weight_t>>, rmm::device_uvector<size_t>> extract_ego(raft::handle_t const &handle, graph_view_t<vertex_t, edge_t, false, multi_gpu> const &graph_view, std::optional<edge_property_view_t<edge_t, weight_t const*>> edge_weight_view, raft::device_span<vertex_t const> source_vertices, vertex_t radius, bool do_expensive_check = false)#returns induced EgoNet subgraph(s) of neighbors centered at nodes in source_vertex within a given radius.
.*
- 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.
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. Must have at least one worker stream.
graph_view – Graph view object of, we extract induced egonet subgraphs from
graph_view
.edge_weight_view – Optional view object holding edge weights for
graph_view
.source_vertex – Pointer to egonet center vertices (size ==
n_subgraphs
).n_subgraphs – Number of induced EgoNet subgraphs to extract (ie. number of elements in
source_vertex
).radius – Include all neighbors of distance <= radius from
source_vertex
.- Returns:
std::tuple<rmm::device_uvector<vertex_t>, rmm::device_uvector<vertex_t>, rmm::device_uvector<weight_t>, rmm::device_uvector<size_t>> Quadraplet of edge source vertices, edge destination vertices, edge weights, and edge offsets for each induced EgoNet subgraph.
- template<typename vertex_t, typename edge_t, bool multi_gpu>
void triangle_count(raft::handle_t const &handle, graph_view_t<vertex_t, edge_t, false, multi_gpu> const &graph_view, std::optional<raft::device_span<vertex_t const>> vertices, raft::device_span<edge_t> counts, bool do_expensive_check = false)#Compute triangle counts.
Compute triangle counts for the entire set of vertices (if
vertices
is std::nullopt) or the given vertices (vertices.has_value()
is true).
- 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)
- 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.
vertices – Vertices to compute triangle counts. If
vertices.has_value()
is false, compute triangle counts for the entire set of vertices.counts – Output triangle count array. The size of the array should be the local vertex partition range size (if
vertices
is std::nullopt) or the size ofvertices
(ifvertices.has_value()
is true).do_expensive_check – A flag to run expensive checks for input arguments (if set to
true
).
- template<typename vertex_t, typename edge_t, bool multi_gpu>
edge_property_t<graph_view_t<vertex_t, edge_t, false, multi_gpu>, edge_t> edge_triangle_count(raft::handle_t const &handle, graph_view_t<vertex_t, edge_t, false, multi_gpu> const &graph_view, bool do_expensive_check = false)#Compute edge triangle counts.
.*
Compute edge triangle counts for the entire set of edges.
- 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)
- 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.
do_expensive_check – A flag to run expensive checks for input arguments (if set to
true
).- Returns:
edge_property_t containing the edge triangle count
- template<typename vertex_t, typename edge_t, typename weight_t, bool multi_gpu>
std::tuple<rmm::device_uvector<vertex_t>, rmm::device_uvector<vertex_t>, std::optional<rmm::device_uvector<weight_t>>> k_truss(raft::handle_t const &handle, graph_view_t<vertex_t, edge_t, false, multi_gpu> const &graph_view, std::optional<edge_property_view_t<edge_t, weight_t const*>> edge_weight_view, edge_t k, bool do_expensive_check = false)#Compute K-Truss.
.*
Extract the K-Truss subgraph of a graph
- 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)
- 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
.k – The desired k to be used for extracting the K-Truss subgraph
do_expensive_check – A flag to run expensive checks for input arguments (if set to
true
).- Returns:
edge list of the K-Truss subgraph