Attention
The vector search and clustering algorithms in RAFT are being migrated to a new library dedicated to vector search called cuVS. We will continue to support the vector search algorithms in RAFT during this move, but will no longer update them after the RAPIDS 24.06 (June) release. We plan to complete the migration by RAPIDS 24.10 (October) release and they will be removed from RAFT altogether in the 24.12 (December) release.
Sparse Solvers#
-
enum raft::sparse::solver::LANCZOS_WHICH#
Enumeration specifying which eigenvalues to compute in the Lanczos algorithm.
Values:
-
enumerator LA#
LA: Largest (algebraic) eigenvalues.
-
enumerator LM#
LM: Largest (in magnitude) eigenvalues.
-
enumerator SA#
SA: Smallest (algebraic) eigenvalues.
-
enumerator SM#
SM: Smallest (in magnitude) eigenvalues.
-
enumerator LA#
-
template<typename IndexTypeT, typename ValueTypeT>
int raft::sparse::solver::lanczos_compute_smallest_eigenvectors(raft::resources const &handle, lanczos_solver_config<ValueTypeT> const &config, raft::device_csr_matrix_view<ValueTypeT, IndexTypeT, IndexTypeT, IndexTypeT> A, std::optional<raft::device_vector_view<ValueTypeT, uint32_t, raft::row_major>> v0, raft::device_vector_view<ValueTypeT, uint32_t, raft::col_major> eigenvalues, raft::device_matrix_view<ValueTypeT, uint32_t, raft::col_major> eigenvectors)# Find the eigenpairs using lanczos solver.
- Template Parameters:
IndexTypeT – the type of data used for indexing.
ValueTypeT – the type of data used for weights, distances.
- Parameters:
handle – the raft handle.
config – lanczos config used to set hyperparameters
A – Sparse matrix in CSR format.
v0 – Optional Initial lanczos vector
eigenvalues – output eigenvalues
eigenvectors – output eigenvectors
- Returns:
Zero if successful. Otherwise non-zero.
-
template<typename IndexTypeT, typename ValueTypeT>
int raft::sparse::solver::lanczos_compute_smallest_eigenvectors(raft::resources const &handle, lanczos_solver_config<ValueTypeT> const &config, raft::device_coo_matrix_view<ValueTypeT, IndexTypeT, IndexTypeT, IndexTypeT> A, std::optional<raft::device_vector_view<ValueTypeT, uint32_t, raft::row_major>> v0, raft::device_vector_view<ValueTypeT, uint32_t, raft::col_major> eigenvalues, raft::device_matrix_view<ValueTypeT, uint32_t, raft::col_major> eigenvectors)# Find the eigenpairs using lanczos solver.
- Template Parameters:
IndexTypeT – the type of data used for indexing.
ValueTypeT – the type of data used for weights, distances.
- Parameters:
handle – the raft handle.
config – lanczos config used to set hyperparameters
A – Sparse matrix in COO format.
v0 – Optional Initial lanczos vector
eigenvalues – output eigenvalues
eigenvectors – output eigenvectors
- Returns:
Zero if successful. Otherwise non-zero.
-
template<typename IndexTypeT, typename ValueTypeT>
int raft::sparse::solver::lanczos_compute_smallest_eigenvectors(raft::resources const &handle, lanczos_solver_config<ValueTypeT> const &config, raft::device_vector_view<IndexTypeT, uint32_t, raft::row_major> rows, raft::device_vector_view<IndexTypeT, uint32_t, raft::row_major> cols, raft::device_vector_view<ValueTypeT, uint32_t, raft::row_major> vals, std::optional<raft::device_vector_view<ValueTypeT, uint32_t, raft::row_major>> v0, raft::device_vector_view<ValueTypeT, uint32_t, raft::col_major> eigenvalues, raft::device_matrix_view<ValueTypeT, uint32_t, raft::col_major> eigenvectors)# Find the eigenpairs using lanczos solver.
- Template Parameters:
index_type_t – the type of data used for indexing.
value_type_t – the type of data used for weights, distances.
- Parameters:
handle – the raft handle.
config – lanczos config used to set hyperparameters
rows – Vector view of the rows of the sparse CSR matrix.
cols – Vector view of the cols of the sparse CSR matrix.
vals – Vector view of the vals of the sparse CSR matrix.
v0 – Optional Initial lanczos vector
eigenvalues – output eigenvalues
eigenvectors – output eigenvectors
- Returns:
Zero if successful. Otherwise non-zero.
-
template<typename vertex_t, typename edge_t, typename weight_t, typename alteration_t = weight_t>
Graph_COO<vertex_t, edge_t, weight_t> raft::sparse::solver::mst(raft::resources const &handle, edge_t const *offsets, vertex_t const *indices, weight_t const *weights, vertex_t const v, edge_t const e, vertex_t *color, cudaStream_t stream, bool symmetrize_output = true, bool initialize_colors = true, int iterations = 0) Compute the minimum spanning tree (MST) or minimum spanning forest (MSF) depending on the connected components of the given graph.
- Template Parameters:
vertex_t – integral type for precision of vertex indexing
edge_t – integral type for precision of edge indexing
weight_t – type of weights array
alteration_t – type to use for random alteration
- Parameters:
handle –
offsets – csr inptr array of row offsets (size v+1)
indices – csr array of column indices (size e)
weights – csr array of weights (size e)
v – number of vertices in graph
e – number of edges in graph
color – array to store resulting colors for MSF
stream – cuda stream for ordering operations
symmetrize_output – should the resulting output edge list should be symmetrized?
initialize_colors – should the colors array be initialized inside the MST?
iterations – maximum number of iterations to perform
- Returns:
a list of edges containing the mst (or a subset of the edges guaranteed to be in the mst when an msf is encountered)
-
template<typename vertex_t, typename edge_t, typename weight_t>
struct Graph_COO#
-
template<typename ValueTypeT>
struct lanczos_solver_config# - #include <lanczos_types.hpp>
Configuration parameters for the Lanczos eigensolver.
This structure encapsulates all configuration parameters needed to run the Lanczos algorithm for computing eigenvalues and eigenvectors of large sparse matrices.
- Template Parameters:
ValueTypeT – Data type for values (float or double)
Public Members
-
int n_components#
The number of eigenvalues and eigenvectors to compute.
Note
Must be 1 <= n_components < n, where n is the matrix dimension
-
int max_iterations#
Maximum number of iterations allowed for the algorithm to converge.
-
int ncv#
The number of Lanczos vectors to generate.
Note
Must satisfy n_components + 1 < ncv < n, where n is the matrix dimension
-
ValueTypeT tolerance#
Convergence tolerance for residuals.
Note
Used to determine when to stop iteration based on ||Ax - wx|| < tolerance
-
LANCZOS_WHICH which#
Specifies which eigenvalues to compute in the Lanczos algorithm.
See also
LANCZOS_WHICH for possible values (SA, LA, SM, LM)
-
std::optional<uint64_t> seed = std::nullopt#
Random seed for initialization of the algorithm.
Note
Controls reproducibility of results
-
template<typename vertex_t, typename edge_t, typename weight_t, typename alteration_t>
class MST_solver#