Linear Algebra Solvers#
Eigen Decomposition#
#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)# eig decomp with divide and conquer method for the column-major symmetric matrices
- Template Parameters:
ValueType – the data-type of input and output
IntegerType – Integer used for addressing
- Parameters:
handle – raft::resources
in – [in] input raft::device_matrix_view (symmetric matrix that has real eig values and vectors)
eig_vectors – [out] eigenvectors output of type raft::device_matrix_view
eig_vals – [out] eigen values output of type raft::device_vector_view
-
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)# eig decomp to select top-n eigen values with divide and conquer method for the column-major symmetric matrices
- Template Parameters:
ValueType – the data-type of input and output
IntegerType – Integer used for addressing
- Parameters:
handle – [in] raft::resources
in – [in] input raft::device_matrix_view (symmetric matrix that has real eig values and vectors)
eig_vectors – [out] eigenvectors output of type raft::device_matrix_view
eig_vals – [out] eigen values output of type raft::device_vector_view
n_eig_vals – [in] number of eigenvectors to be generated
memUsage – [in] the memory selection for eig vector output
-
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)# overloaded function for eig decomp with Jacobi method for the column-major symmetric matrices (in parameter)
- Template Parameters:
ValueType – the data-type of input and output
IntegerType – Integer used for addressing
- Parameters:
handle – raft::resources
in – [in] input raft::device_matrix_view (symmetric matrix that has real eig values and vectors)
eig_vectors – [out] eigenvectors output of type raft::device_matrix_view
eig_vals – [out] eigen values output of type raft::device_vector_view
tol – [in] error tolerance for the jacobi method. Algorithm stops when the Frobenius norm of the absolute error is below tol
sweeps – [in] number of sweeps in the Jacobi algorithm. The more the better accuracy.
QR Decomposition#
#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)# Compute the QR decomposition of matrix M and return only the Q matrix.
- Parameters:
handle – [in] raft::resources
M – [in] Input raft::device_matrix_view
Q – [out] Output raft::device_matrix_view
-
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)# Compute the QR decomposition of matrix M and return both the Q and R matrices.
- Parameters:
handle – [in] raft::resources
M – [in] Input raft::device_matrix_view
Q – [in] Output raft::device_matrix_view
R – [out] Output raft::device_matrix_view
Randomized Singular-Value Decomposition#
#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)# randomized singular value decomposition (RSVD) on a column major rectangular matrix using QR decomposition, by specifying no. of PCs and upsamples directly
- Template Parameters:
ValueType – value type of parameters
IndexType – index type of parameters
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
- Parameters:
handle – [in] raft::resources
M – [in] input raft::device_matrix_view with layout raft::col_major of shape (M, N)
S_vec – [out] singular values raft::device_vector_view of shape (K)
p – [in] no. of upsamples
U_in – [out] std::optional left singular values of raft::device_matrix_view with layout raft::col_major
V_in – [out] std::optional right singular values of raft::device_matrix_view with layout raft::col_major
-
template<typename ...Args, typename = std::enable_if_t<sizeof...(Args) == 4>>
void rsvd_fixed_rank(Args... args)# Overload of
rsvd_fixed_rank
to help the compiler find the above overload, in case users pass instd::nullopt
for one or both of the optional arguments.Please see above for documentation of
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)# randomized singular value decomposition (RSVD) on a column major rectangular matrix using symmetric Eigen decomposition, by specifying no. of PCs and upsamples directly. The rectangular input matrix is made square and symmetric using B @ B^T
- Template Parameters:
ValueType – value type of parameters
IndexType – index type of parameters
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
- Parameters:
handle – [in] raft::resources
M – [in] input raft::device_matrix_view with layout raft::col_major of shape (M, N)
S_vec – [out] singular values raft::device_vector_view of shape (K)
p – [in] no. of upsamples
U_in – [out] std::optional left singular values of raft::device_matrix_view with layout raft::col_major
V_in – [out] std::optional right singular values of raft::device_matrix_view with layout raft::col_major
-
template<typename ...Args, typename = std::enable_if_t<sizeof...(Args) == 4>>
void rsvd_fixed_rank_symmetric(Args... args)# Overload of
rsvd_fixed_rank_symmetric
to help the compiler find the above overload, in case users pass instd::nullopt
for one or both of the optional arguments.Please see above for documentation of
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)# randomized singular value decomposition (RSVD) on a column major rectangular matrix using Jacobi method, by specifying no. of PCs and upsamples directly
- Template Parameters:
ValueType – value type of parameters
IndexType – index type of parameters
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
- Parameters:
handle – [in] raft::resources
M – [in] input raft::device_matrix_view with layout raft::col_major of shape (M, N)
S_vec – [out] singular values raft::device_vector_view of shape (K)
p – [in] no. of upsamples
tol – [in] tolerance for Jacobi-based solvers
max_sweeps – [in] maximum number of sweeps for Jacobi-based solvers
U_in – [out] std::optional left singular values of raft::device_matrix_view with layout raft::col_major
V_in – [out] std::optional right singular values of raft::device_matrix_view with layout raft::col_major
-
template<typename ...Args, typename = std::enable_if_t<sizeof...(Args) == 6>>
void rsvd_fixed_rank_jacobi(Args... args)# Overload of
rsvd_fixed_rank_jacobi
to help the compiler find the above overload, in case users pass instd::nullopt
for one or both of the optional arguments.Please see above for documentation of
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)# randomized singular value decomposition (RSVD) on a column major rectangular matrix using Jacobi method, by specifying no. of PCs and upsamples directly. The rectangular input matrix is made square and symmetric using B @ B^T
- Template Parameters:
ValueType – value type of parameters
IndexType – index type of parameters
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
- Parameters:
handle – [in] raft::resources
M – [in] input raft::device_matrix_view with layout raft::col_major of shape (M, N)
S_vec – [out] singular values raft::device_vector_view of shape (K)
p – [in] no. of upsamples
tol – [in] tolerance for Jacobi-based solvers
max_sweeps – [in] maximum number of sweeps for Jacobi-based solvers
U_in – [out] std::optional left singular values of raft::device_matrix_view with layout raft::col_major
V_in – [out] std::optional right singular values of raft::device_matrix_view with layout raft::col_major
-
template<typename ...Args, typename = std::enable_if_t<sizeof...(Args) == 6>>
void rsvd_fixed_rank_symmetric_jacobi(Args... args)# Overload of
rsvd_fixed_rank_symmetric_jacobi
to help the compiler find the above overload, in case users pass instd::nullopt
for one or both of the optional arguments.Please see above for documentation of
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)# randomized singular value decomposition (RSVD) on a column major rectangular matrix using QR decomposition, by specifying the PC and upsampling ratio
- Template Parameters:
ValueType – value type of parameters
IndexType – index type of parameters
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
- Parameters:
handle – [in] raft::resources
M – [in] input raft::device_matrix_view with layout raft::col_major of shape (M, N)
S_vec – [out] singular values raft::device_vector_view of shape (K)
PC_perc – [in] percentage of singular values to be computed
UpS_perc – [in] upsampling percentage
U_in – [out] std::optional left singular values of raft::device_matrix_view with layout raft::col_major
V_in – [out] std::optional right singular values of raft::device_matrix_view with layout raft::col_major
-
template<typename ...Args, typename = std::enable_if_t<sizeof...(Args) == 5>>
void rsvd_perc(Args... args)# Overload of
rsvd_perc
to help the compiler find the above overload, in case users pass instd::nullopt
for one or both of the optional arguments.Please see above for documentation of
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)# randomized singular value decomposition (RSVD) on a column major rectangular matrix using symmetric Eigen decomposition, by specifying the PC and upsampling ratio. The rectangular input matrix is made square and symmetric using B @ B^T
- Template Parameters:
ValueType – value type of parameters
IndexType – index type of parameters
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
- Parameters:
handle – [in] raft::resources
M – [in] input raft::device_matrix_view with layout raft::col_major of shape (M, N)
S_vec – [out] singular values raft::device_vector_view of shape (K)
PC_perc – [in] percentage of singular values to be computed
UpS_perc – [in] upsampling percentage
U_in – [out] std::optional left singular values of raft::device_matrix_view with layout raft::col_major
V_in – [out] std::optional right singular values of raft::device_matrix_view with layout raft::col_major
-
template<typename ...Args, typename = std::enable_if_t<sizeof...(Args) == 5>>
void rsvd_perc_symmetric(Args... args)# Overload of
rsvd_perc_symmetric
to help the compiler find the above overload, in case users pass instd::nullopt
for one or both of the optional arguments.Please see above for documentation of
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)# randomized singular value decomposition (RSVD) on a column major rectangular matrix using Jacobi method, by specifying the PC and upsampling ratio
- Template Parameters:
ValueType – value type of parameters
IndexType – index type of parameters
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
- Parameters:
handle – [in] raft::resources
M – [in] input raft::device_matrix_view with layout raft::col_major of shape (M, N)
S_vec – [out] singular values raft::device_vector_view of shape (K)
PC_perc – [in] percentage of singular values to be computed
UpS_perc – [in] upsampling percentage
tol – [in] tolerance for Jacobi-based solvers
max_sweeps – [in] maximum number of sweeps for Jacobi-based solvers
U_in – [out] std::optional left singular values of raft::device_matrix_view with layout raft::col_major
V_in – [out] std::optional right singular values of raft::device_matrix_view with layout raft::col_major
-
template<typename ...Args, typename = std::enable_if_t<sizeof...(Args) == 7>>
void rsvd_perc_jacobi(Args... args)# Overload of
rsvd_perc_jacobi
to help the compiler find the above overload, in case users pass instd::nullopt
for one or both of the optional arguments.Please see above for documentation of
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)# randomized singular value decomposition (RSVD) on a column major rectangular matrix using Jacobi method, by specifying the PC and upsampling ratio. The rectangular input matrix is made square and symmetric using B @ B^T
- Template Parameters:
ValueType – value type of parameters
IndexType – index type of parameters
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
- Parameters:
handle – [in] raft::resources
M – [in] input raft::device_matrix_view with layout raft::col_major of shape (M, N)
S_vec – [out] singular values raft::device_vector_view of shape (K)
PC_perc – [in] percentage of singular values to be computed
UpS_perc – [in] upsampling percentage
tol – [in] tolerance for Jacobi-based solvers
max_sweeps – [in] maximum number of sweeps for Jacobi-based solvers
U_in – [out] std::optional left singular values of raft::device_matrix_view with layout raft::col_major
V_in – [out] std::optional right singular values of raft::device_matrix_view with layout raft::col_major
-
template<typename ...Args, typename = std::enable_if_t<sizeof...(Args) == 7>>
void rsvd_perc_symmetric_jacobi(Args... args)# Overload of
rsvd_perc_symmetric_jacobi
to help the compiler find the above overload, in case users pass instd::nullopt
for one or both of the optional arguments.Please see above for documentation of
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)# randomized singular value decomposition (RSVD) using cusolver
- Template Parameters:
math_t – the data type
idx_t – index type
- Parameters:
handle – [in] raft handle
in – [in] input matrix in col-major format. Warning: the content of this matrix is modified by the cuSOLVER routines. [dim = n_rows * n_cols]
S – [out] array of singular values of input matrix. The rank k must be less than min(m,n). [dim = k]
U – [out] optional left singular values of input matrix. Use std::nullopt to not generate it. [dim = n_rows * k]
V – [out] optional right singular values of input matrix. Use std::nullopt to not generate it. [dim = k * n_cols]
p – [in] Oversampling. The size of the subspace will be (k + p). (k+p) is less than min(m,n). (Recommended to be at least 2*k)
niters – [in] Number of iteration of power method. (2 is recommended)
-
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)# Overload of
randomized_svd
to help the compiler find the above overload, in case users pass instd::nullopt
for the optional arguments.Please see above for documentation of
randomized_svd
.
Singular-Value Decomposition#
#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> 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)# singular value decomposition (SVD) on a column major matrix using QR decomposition
- Template Parameters:
ValueType – value type of parameters
IndexType – index type of parameters
- Parameters:
handle – [in] raft::resources
in – [in] input raft::device_matrix_view with layout raft::col_major of shape (M, N)
sing_vals – [out] singular values raft::device_vector_view of shape (K)
U – [out] std::optional left singular values of raft::device_matrix_view with layout raft::col_major and dimensions (m, n)
V – [out] std::optional right singular values of raft::device_matrix_view with layout raft::col_major and dimensions (n, n)
-
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)# Overload of
svd_qr
to help the compiler find the above overload, in case users pass instd::nullopt
for one or both of the optional arguments.Please see above for documentation of
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> 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)# singular value decomposition (SVD) on a column major matrix using QR decomposition. Right singular vector matrix is transposed before returning
- Template Parameters:
ValueType – value type of parameters
IndexType – index type of parameters
- Parameters:
handle – [in] raft::resources
in – [in] input raft::device_matrix_view with layout raft::col_major of shape (M, N)
sing_vals – [out] singular values raft::device_vector_view of shape (K)
U – [out] std::optional left singular values of raft::device_matrix_view with layout raft::col_major and dimensions (m, n)
V – [out] std::optional right singular values of raft::device_matrix_view with layout raft::col_major and dimensions (n, n)
-
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)# Overload of
svd_qr_transpose_right_vec
to help the compiler find the above overload, in case users pass instd::nullopt
for one or both of the optional arguments.Please see above for documentation of
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> 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)# singular value decomposition (SVD) on a column major matrix using Eigen decomposition. A square symmetric covariance matrix is constructed for the SVD
- Parameters:
handle – [in] raft::resources
in – [in] input raft::device_matrix_view with layout raft::col_major of shape (M, N)
S – [out] singular values raft::device_vector_view of shape (K)
V – [out] right singular values of raft::device_matrix_view with layout raft::col_major and dimensions (n, n)
U – [out] optional left singular values of raft::device_matrix_view with layout raft::col_major and dimensions (m, n)
-
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)# reconstruct a matrix use left and right singular vectors and singular values
- Parameters:
handle – [in] raft::resources
U – [in] left singular values of raft::device_matrix_view with layout raft::col_major and dimensions (m, k)
S – [in] square matrix with singular values on its diagonal of shape (k, k)
V – [in] right singular values of raft::device_matrix_view with layout raft::col_major and dimensions (k, n)
out – [out] output raft::device_matrix_view with layout raft::col_major of shape (m, n)
Least Squares#
#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)# Solves the linear ordinary least squares problem
Aw = b
Via SVD decomposition ofA = U S Vt
.- Template Parameters:
ValueType – the data-type of input/output
- Parameters:
handle – [in] raft::resources
A – [inout] input raft::device_matrix_view Warning: the content of this matrix is modified.
b – [inout] input target raft::device_vector_view Warning: the content of this vector is modified.
w – [out] output coefficient raft::device_vector_view
-
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> b, raft::device_vector_view<ValueType, IndexType> w)# Solves the linear ordinary least squares problem
Aw = b
Via SVD decomposition ofA = U S V^T
using Jacobi iterations.- Template Parameters:
ValueType – the data-type of input/output
- Parameters:
handle – [in] raft::resources
A – [inout] input raft::device_matrix_view Warning: the content of this matrix is modified.
b – [inout] input target raft::device_vector_view Warning: the content of this vector is modified.
w – [out] output coefficient raft::device_vector_view
-
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> b, raft::device_vector_view<ValueType, IndexType> w)# Solves the linear ordinary least squares problem
Aw = b
via eigenvalue decomposition ofA^T * A
(covariance matrix for dataset A). (w = (A^T A)^-1 A^T b
)- Template Parameters:
ValueType – the data-type of input/output
- Parameters:
handle – [in] raft::resources
A – [inout] input raft::device_matrix_view Warning: the content of this matrix is modified by the cuSOLVER routines.
b – [inout] input target raft::device_vector_view Warning: the content of this vector is modified by the cuSOLVER routines.
w – [out] output coefficient raft::device_vector_view
-
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> b, raft::device_vector_view<ValueType, IndexType> w)# Solves the linear ordinary least squares problem
Aw = b
via QR decomposition ofA = QR
. (triangular system of equationsRw = Q^T b
)- Template Parameters:
ValueType – the data-type of input/output
- Parameters:
handle – [in] raft::resources
A – [inout] input raft::device_matrix_view Warning: the content of this matrix is modified.
b – [inout] input target raft::device_vector_view Warning: the content of this vector is modified.
w – [out] output coefficient raft::device_vector_view