注意

RAFT 中的向量搜索和聚类算法正在迁移到专门用于向量搜索的新库 cuVS。在迁移过程中,我们将继续支持 RAFT 中的向量搜索算法,但在 RAPIDS 24.06 (June) 版本后将不再更新它们。我们计划在 RAPIDS 24.10 (October) 版本之前完成迁移,并在 24.12 (December) 版本中将它们完全从 RAFT 中移除。

线性代数求解器#

特征分解#

#include <raft/linalg/eig.cuh>

namespace raft::linalg

template<typename ValueType, typename IndexType>
void eig_dc(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> in, raft::device_matrix_view<ValueType, IndexType, raft::col_major> eig_vectors, raft::device_vector_view<ValueType, IndexType> eig_vals)#

使用分治法对列主序对称矩阵进行特征分解

模板参数:
  • ValueType – 输入和输出的数据类型

  • IntegerType – 用于寻址的整数

参数:
template<typename ValueType, typename IndexType>
void eig_dc_selective(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> in, raft::device_matrix_view<ValueType, IndexType, raft::col_major> eig_vectors, raft::device_vector_view<ValueType, IndexType> eig_vals, std::size_t n_eig_vals, EigVecMemUsage memUsage)#

使用分治法对列主序对称矩阵进行特征分解,以选择前 n 个特征值

模板参数:
  • ValueType – 输入和输出的数据类型

  • IntegerType – 用于寻址的整数

参数:
  • handle[in] raft::resources

  • in[in] 输入 raft::device_matrix_view (具有实特征值和特征向量的对称矩阵)

  • eig_vectors[out] 类型为 raft::device_matrix_view 的特征向量输出

  • eig_vals[out] 类型为 raft::device_vector_view 的特征值输出

  • n_eig_vals[in] 要生成的特征向量数量

  • memUsage[in] 特征向量输出的内存选择

template<typename ValueType, typename IndexType>
void eig_jacobi(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> in, raft::device_matrix_view<ValueType, IndexType, raft::col_major> eig_vectors, raft::device_vector_view<ValueType, IndexType> eig_vals, ValueType tol = 1.e-7, int sweeps = 15)#

使用 Jacobi 方法对列主序对称矩阵(作为输入参数)进行特征分解的重载函数

模板参数:
  • ValueType – 输入和输出的数据类型

  • IntegerType – 用于寻址的整数

参数:
  • handle – raft::resources

  • in[in] 输入 raft::device_matrix_view (具有实特征值和特征向量的对称矩阵)

  • eig_vectors[out] 类型为 raft::device_matrix_view 的特征向量输出

  • eig_vals[out] 类型为 raft::device_vector_view 的特征值输出

  • tol[in] Jacobi 方法的误差容限。当绝对误差的 Frobenius 范数低于 tol 时,算法停止。

  • sweeps[in] Jacobi 算法中的扫描次数。次数越多精度越高。

QR 分解#

#include <raft/linalg/qr.cuh>

namespace raft::linalg

template<typename ElementType, typename IndexType>
void qr_get_q(raft::resources const &handle, raft::device_matrix_view<const ElementType, IndexType, raft::col_major> M, raft::device_matrix_view<ElementType, IndexType, raft::col_major> Q)#

计算矩阵 M 的 QR 分解并仅返回矩阵 Q。

参数:
template<typename ElementType, typename IndexType>
void qr_get_qr(raft::resources const &handle, raft::device_matrix_view<const ElementType, IndexType, raft::col_major> M, raft::device_matrix_view<ElementType, IndexType, raft::col_major> Q, raft::device_matrix_view<ElementType, IndexType, raft::col_major> R)#

计算矩阵 M 的 QR 分解并返回矩阵 Q 和 R。

参数:

随机奇异值分解#

#include <raft/linalg/rsvd.cuh>

namespace raft::linalg

template<typename ValueType, typename IndexType, typename UType, typename VType>
void rsvd_fixed_rank(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> M, raft::device_vector_view<ValueType, IndexType> S_vec, IndexType p, UType &&U_in, VType &&V_in)#

使用 QR 分解对列主序矩形矩阵进行随机奇异值分解 (RSVD),通过直接指定主成分 (PC) 和上采样数量

模板参数:
  • ValueType – 参数的值类型

  • IndexType – 参数的索引类型

  • UType – std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> U_in

  • VType – std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> V_in

参数:
template<typename ...Args, typename = std::enable_if_t<sizeof...(Args) == 4>>
void rsvd_fixed_rank(Args... args)#

rsvd_fixed_rank 的重载,用于帮助编译器找到上述重载,以防用户为可选参数之一或全部传入 std::nullopt

请参阅上方关于 rsvd_fixed_rank 的文档。

template<typename ValueType, typename IndexType, typename UType, typename VType>
void rsvd_fixed_rank_symmetric(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> M, raft::device_vector_view<ValueType, IndexType> S_vec, IndexType p, UType &&U_in, VType &&V_in)#

使用对称特征分解对列主序矩形矩阵进行随机奇异值分解 (RSVD),通过直接指定主成分 (PC) 和上采样数量。通过 B @ B^T 使矩形输入矩阵变为方阵和对称矩阵

模板参数:
  • ValueType – 参数的值类型

  • IndexType – 参数的索引类型

  • UType – std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> U_in

  • VType – std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> V_in

参数:
template<typename ...Args, typename = std::enable_if_t<sizeof...(Args) == 4>>
void rsvd_fixed_rank_symmetric(Args... args)#

rsvd_fixed_rank_symmetric 的重载,用于帮助编译器找到上述重载,以防用户为可选参数之一或全部传入 std::nullopt

请参阅上方关于 rsvd_fixed_rank_symmetric 的文档。

template<typename ValueType, typename IndexType, typename UType, typename VType>
void rsvd_fixed_rank_jacobi(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> M, raft::device_vector_view<ValueType, IndexType> S_vec, IndexType p, ValueType tol, int max_sweeps, UType &&U_in, VType &&V_in)#

使用 Jacobi 方法对列主序矩形矩阵进行随机奇异值分解 (RSVD),通过直接指定主成分 (PC) 和上采样数量

模板参数:
  • ValueType – 参数的值类型

  • IndexType – 参数的索引类型

  • UType – std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> U_in

  • VType – std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> V_in

参数:
template<typename ...Args, typename = std::enable_if_t<sizeof...(Args) == 6>>
void rsvd_fixed_rank_jacobi(Args... args)#

rsvd_fixed_rank_jacobi 的重载,用于帮助编译器找到上述重载,以防用户为可选参数之一或全部传入 std::nullopt

请参阅上方关于 rsvd_fixed_rank_jacobi 的文档。

template<typename ValueType, typename IndexType, typename UType, typename VType>
void rsvd_fixed_rank_symmetric_jacobi(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> M, raft::device_vector_view<ValueType, IndexType> S_vec, IndexType p, ValueType tol, int max_sweeps, UType &&U_in, VType &&V_in)#

使用 Jacobi 方法对列主序矩形矩阵进行随机奇异值分解 (RSVD),通过直接指定主成分 (PC) 和上采样数量。通过 B @ B^T 使矩形输入矩阵变为方阵和对称矩阵

模板参数:
  • ValueType – 参数的值类型

  • IndexType – 参数的索引类型

  • UType – std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> U_in

  • VType – std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> V_in

参数:
template<typename ...Args, typename = std::enable_if_t<sizeof...(Args) == 6>>
void rsvd_fixed_rank_symmetric_jacobi(Args... args)#

rsvd_fixed_rank_symmetric_jacobi 的重载,用于帮助编译器找到上述重载,以防用户为可选参数之一或全部传入 std::nullopt

请参阅上方关于 rsvd_fixed_rank_symmetric_jacobi 的文档。

template<typename ValueType, typename IndexType, typename UType, typename VType>
void rsvd_perc(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> M, raft::device_vector_view<ValueType, IndexType> S_vec, ValueType PC_perc, ValueType UpS_perc, UType &&U_in, VType &&V_in)#

使用 QR 分解对列主序矩形矩阵进行随机奇异值分解 (RSVD),通过指定主成分 (PC) 和上采样比例

模板参数:
  • ValueType – 参数的值类型

  • IndexType – 参数的索引类型

  • UType – std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> U_in

  • VType – std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> V_in

参数:
template<typename ...Args, typename = std::enable_if_t<sizeof...(Args) == 5>>
void rsvd_perc(Args... args)#

rsvd_perc 的重载,用于帮助编译器找到上述重载,以防用户为可选参数之一或全部传入 std::nullopt

请参阅上方关于 rsvd_perc 的文档。

template<typename ValueType, typename IndexType, typename UType, typename VType>
void rsvd_perc_symmetric(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> M, raft::device_vector_view<ValueType, IndexType> S_vec, ValueType PC_perc, ValueType UpS_perc, UType &&U_in, VType &&V_in)#

通过指定 PC 百分比和过采样率,使用对称特征分解对列主序矩形矩阵进行随机奇异值分解 (RSVD)。矩形输入矩阵通过 B @ B^T 变换为方阵和对称矩阵

模板参数:
  • ValueType – 参数的值类型

  • IndexType – 参数的索引类型

  • UType – std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> U_in

  • VType – std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> V_in

参数:
template<typename ...Args, typename = std::enable_if_t<sizeof...(Args) == 5>>
void rsvd_perc_symmetric(Args... args)#

该函数是 rsvd_perc_symmetric 的重载,用于帮助编译器找到上方的重载定义,以防用户为可选参数中的一个或两个传入 std::nullopt

请参阅上方 rsvd_perc_symmetric 的文档。

template<typename ValueType, typename IndexType, typename UType, typename VType>
void rsvd_perc_jacobi(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> M, raft::device_vector_view<ValueType, IndexType> S_vec, ValueType PC_perc, ValueType UpS_perc, ValueType tol, int max_sweeps, UType &&U_in, VType &&V_in)#

通过指定 PC 百分比和过采样率,使用 Jacobi 方法对列主序矩形矩阵进行随机奇异值分解 (RSVD)

模板参数:
  • ValueType – 参数的值类型

  • IndexType – 参数的索引类型

  • UType – std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> U_in

  • VType – std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> V_in

参数:
template<typename ...Args, typename = std::enable_if_t<sizeof...(Args) == 7>>
void rsvd_perc_jacobi(Args... args)#

该函数是 rsvd_perc_jacobi 的重载,用于帮助编译器找到上方的重载定义,以防用户为可选参数中的一个或两个传入 std::nullopt

请参阅上方 rsvd_perc_jacobi 的文档。

template<typename ValueType, typename IndexType, typename UType, typename VType>
void rsvd_perc_symmetric_jacobi(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> M, raft::device_vector_view<ValueType, IndexType> S_vec, ValueType PC_perc, ValueType UpS_perc, ValueType tol, int max_sweeps, UType &&U_in, VType &&V_in)#

通过指定 PC 百分比和过采样率,使用 Jacobi 方法对列主序矩形矩阵进行随机奇异值分解 (RSVD)。矩形输入矩阵通过 B @ B^T 变换为方阵和对称矩阵

模板参数:
  • ValueType – 参数的值类型

  • IndexType – 参数的索引类型

  • UType – std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> U_in

  • VType – std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> V_in

参数:
template<typename ...Args, typename = std::enable_if_t<sizeof...(Args) == 7>>
void rsvd_perc_symmetric_jacobi(Args... args)#

该函数是 rsvd_perc_symmetric_jacobi 的重载,用于帮助编译器找到上方的重载定义,以防用户为可选参数中的一个或两个传入 std::nullopt

请参阅上方 rsvd_perc_symmetric_jacobi 的文档。

template<typename math_t, typename idx_t>
void randomized_svd(const raft::resources &handle, raft::device_matrix_view<const math_t, idx_t, raft::col_major> in, raft::device_vector_view<math_t, idx_t> S, std::optional<raft::device_matrix_view<math_t, idx_t, raft::col_major>> U, std::optional<raft::device_matrix_view<math_t, idx_t, raft::col_major>> V, std::size_t p, std::size_t niters)#

使用 cusolver 进行随机奇异值分解 (RSVD)

模板参数:
  • math_t – 数据类型

  • idx_t – 索引类型

参数:
  • handle[输入] raft 句柄

  • in[输入] 列主序格式的输入矩阵。警告:此矩阵的内容会被 cuSOLVER 例程修改。[维度 = 行数 * 列数]

  • S[输出] 输入矩阵的奇异值数组。秩 k 必须小于 min(m,n)。[维度 = k]

  • U[输出] 可选的输入矩阵左奇异值。使用 std::nullopt 表示不生成它。[维度 = 行数 * k]

  • V[输出] 可选的输入矩阵右奇异值。使用 std::nullopt 表示不生成它。[维度 = k * 列数]

  • p[输入] 过采样。子空间的尺寸将是 (k + p)。(k+p) 小于 min(m,n)。(建议至少为 2*k)

  • niters[输入] 幂方法的迭代次数。(建议为 2)

template<typename math_t, typename idx_t, typename opt_u_vec_t, typename opt_v_vec_t>
void randomized_svd(const raft::resources &handle, raft::device_matrix_view<const math_t, idx_t, raft::col_major> in, raft::device_vector_view<math_t, idx_t> S, opt_u_vec_t &&U, opt_v_vec_t &&V, std::size_t p, std::size_t niters)#

该函数是 randomized_svd 的重载,用于帮助编译器找到上方的重载定义,以防用户为可选参数传入 std::nullopt

请参阅上方 randomized_svd 的文档。

奇异值分解#

#include <raft/linalg/svd.cuh>

namespace raft::linalg

template<typename ValueType, typename IndexType>
void svd_qr(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> in, raft::device_vector_view<ValueType, IndexType, raft::col_major>> sing_vals, std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> U = std::nullopt, std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> V = std::nullopt)#

使用 QR 分解对列主序矩阵进行奇异值分解 (SVD)

模板参数:
  • ValueType – 参数的值类型

  • IndexType – 参数的索引类型

参数:
template<typename ValueType, typename IndexType, typename UType, typename VType>
void svd_qr(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> in, raft::device_vector_view<ValueType, IndexType> sing_vals, UType &&U_in = std::nullopt, VType &&V_in = std::nullopt)#

该函数是 svd_qr 的重载,用于帮助编译器找到上方的重载定义,以防用户为可选参数中的一个或两个传入 std::nullopt

请参阅上方 svd_qr 的文档。

template<typename ValueType, typename IndexType>
void svd_qr_transpose_right_vec(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> in, raft::device_vector_view<ValueType, IndexType, raft::col_major> sing_vals, std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> U = std::nullopt, std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> V = std::nullopt)#

使用 QR 分解对列主序矩阵进行奇异值分解 (SVD)。返回前右奇异向量矩阵会被转置

模板参数:
  • ValueType – 参数的值类型

  • IndexType – 参数的索引类型

参数:
template<typename ValueType, typename IndexType, typename UType, typename VType>
void svd_qr_transpose_right_vec(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> in, raft::device_vector_view<ValueType, IndexType> sing_vals, UType &&U_in = std::nullopt, VType &&V_in = std::nullopt)#

该函数是 svd_qr_transpose_right_vec 的重载,用于帮助编译器找到上方的重载定义,以防用户为可选参数中的一个或两个传入 std::nullopt

请参阅上方 svd_qr_transpose_right_vec 的文档。

template<typename ValueType, typename IndexType>
void svd_eig(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> in, raft::device_vector_view<ValueType, IndexType, raft::col_major>> S, raft::device_matrix_view<ValueType, IndexType, raft::col_major> V, std::optional<raft::device_matrix_view<ValueType, IndexType, raft::col_major>> U = std::nullopt)#

使用特征分解对列主序矩阵进行奇异值分解 (SVD)。为进行 SVD 构建一个方阵对称协方差矩阵

参数:
template<typename ValueType, typename IndexType, typename UType>
void svd_eig(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> in, raft::device_vector_view<ValueType, IndexType> S, raft::device_matrix_view<ValueType, IndexType, raft::col_major> V, UType &&U = std::nullopt)#
template<typename ValueType, typename IndexType>
void svd_reconstruction(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> U, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> S, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> V, raft::device_matrix_view<ValueType, IndexType, raft::col_major> out)#

使用左奇异向量、右奇异向量和奇异值重建矩阵

参数:

最小二乘#

#include <raft/linalg/lstsq.cuh>

namespace raft::linalg

template<typename ValueType, typename IndexType>
void lstsq_svd_qr(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> A, raft::device_vector_view<const ValueType, IndexType> b, raft::device_vector_view<ValueType, IndexType> w)#

通过 A = U S Vt 的 SVD 分解求解线性普通最小二乘问题 Aw = b

模板参数:

ValueType – 输入/输出的数据类型

参数:
template<typename ValueType, typename IndexType>
void lstsq_svd_jacobi(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> A, raft::device_vector_view<const ValueType, IndexType, raft::col_major> A, raft::device_vector_view<const ValueType, IndexType> b, raft::device_vector_view<ValueType, IndexType> w)#

通过使用 Jacobi 迭代的 A = U S V^T 的 SVD 分解求解线性普通最小二乘问题 Aw = b

模板参数:

ValueType – 输入/输出的数据类型

参数:
template<typename ValueType, typename IndexType>
void lstsq_eig(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> A, raft::device_vector_view<const ValueType, IndexType, raft::col_major> A, raft::device_vector_view<const ValueType, IndexType> b, raft::device_vector_view<ValueType, IndexType> w)#

通过 A^T * A 的特征值分解(数据集 A 的协方差矩阵)求解线性普通最小二乘问题 Aw = b。(w = (A^T A)^-1 A^T b

模板参数:

ValueType – 输入/输出的数据类型

参数:
template<typename ValueType, typename IndexType>
void lstsq_qr(raft::resources const &handle, raft::device_matrix_view<const ValueType, IndexType, raft::col_major> A, raft::device_vector_view<const ValueType, IndexType, raft::col_major> A, raft::device_vector_view<const ValueType, IndexType> b, raft::device_vector_view<ValueType, IndexType> w)#

通过 A = QR 的 QR 分解求解线性普通最小二乘问题 Aw = b。(三角方程组 Rw = Q^T b

模板参数:

ValueType – 输入/输出的数据类型

参数: