Expand description
Safe Rust wrappers for NVIDIA cuSOLVER.
Covers the dense API (Dn) for all four BLAS scalar types:
- LU factorization:
getrf+getrs - QR factorization:
geqrf - Cholesky:
potrf+potrs - SVD:
gesvd - Symmetric / Hermitian eigendecomposition:
syevd/heevd
The generic 64-bit X… API (xgetrf, xgeqrf, xpotrf) gives
type-erased data pointers and is exposed under xapi. The sparse API
(cusolverSp*) is under sparse. The refactor API (cusolverRf*) is
under refactor.
Modules§
- mg
- Multi-GPU dense solvers via
libcusolverMg. Shares dimensions with the single-GPU API but takes arrays of device pointers (one per physical GPU afterHandle::device_select). - refactor
cusolverRf*— fast re-factorization given a sparsity pattern, for solving many systems that differ only in numeric values.- sparse
cusolverSp*— solve sparse linear systems via Cholesky or QR.- xapi
- The generic 64-bit cuSOLVER API (
cusolverDnX*). Matrix dimensions arei64; element types are passed at call-time ascudaDataType. Workspace sizes are split between on-device and on-host buffers.
Structs§
- DnHandle
- Dense cuSOLVER handle.
- Gesvdj
Info - Jacobi-SVD tuning handle.
- Syevj
Info - Jacobi-eigen tuning handle (tolerance + max sweeps).
Enums§
Traits§
- Solver
Scalar - Scalars supported by cuSOLVER’s Dn S/D/C/Z API.
Functions§
- gels
- Solve
A * X = Bin the least-squares sense (iterative-refinement).Aism × n,Bism × nrhs,Xisn × nrhs.AandBmay be overwritten. Returnsiter: number of refinement iterations used (-1 = fallback to full precision). - geqrf
- QR factorization:
A = Q * R. Overwritesa(upper triangle = R, lower = Householder reflectors);taureceives reflector scalars. - gesvd
- Full SVD:
A = U * Σ * Vᵀ.jobu/jobvtare LAPACK-style single-byte selectors (b’A’ = all, b’S’ = economy, b’N’ = none, b’O’ = overwrite A). - gesvdj
- Jacobi SVD:
A = U * diag(s) * Vᴴ.econselects thin-SVD when set. - gesvdj_
batched - Batched Jacobi SVD: batch of
m × nmatrices with stridem×n. - getrf
- In-place LU factorization of a column-major matrix. Overwrites
a. - getrs
- Solve
op(A) * X = Busing the LU factorization fromgetrf. - orgqr
- Generate the orthogonal matrix
Qfrom the factorization produced bygeqrf. After this,aholds the firstncolumns ofQ. - ormqr
- Apply
op(Q)toC:C = op(Q) * C(Left) orC = C * op(Q)(Right), whereQis packed ina+taufromgeqrf. - potrf
- Cholesky factorization:
A = L * Lᵀ(orUᵀ * U). Overwritesa. - potri
- Compute
A = (Lᵀ * L)⁻¹orA = (U * Uᵀ)⁻¹given the Cholesky factor already stored in the triangle selected byuplo.amust hold the output ofpotrfin-place. - potrs
- Solve
A * X = Busing the Cholesky factorization frompotrf. - sgetrf
- Shortcut for
getrfonf32. - sgetrs
- Shortcut for
getrsonf32. - syevd
- Symmetric / Hermitian eigenvalue decomposition:
A = Q * diag(w) * Qᵀ. - syevj
- Jacobi symmetric/Hermitian eigendecomposition (smaller matrices than
syevd, faster convergence on well-conditioned problems). - syevj_
batched - Batched Jacobi symmetric/Hermitian eigendecomposition. Every matrix in
the batch is
n × nand striden × n.wholdsn * batch_sizeeigenvalues, strided byn.
Type Aliases§
- Error
- Error type for cuSOLVER operations.
- Gesvdj
Info Raw - Result
- Result alias.
- Syevj
Info Raw