LLᵀ and LLᴴ
LinearAlgebra.cholesky
— Methodsolver = cholesky(A::CuSparseMatrixCSR{T,Cint}; view::Char='F')
Compute the LLᴴ factorization of a sparse matrix A
on an NVIDIA GPU. The type T
can be Float32
, Float64
, ComplexF32
or ComplexF64
.
Input argument
A
: a sparse Hermitian positive definite matrix stored in theCuSparseMatrixCSR
format.
Keyword argument
*view
: A character that specifies which triangle of the sparse matrix is provided. Possible options are L
for the lower triangle, U
for the upper triangle, and F
for the full matrix.
Output argument
solver
: Opaque structureCudssSolver
that stores the factors of the LLᴴ decomposition.
LinearAlgebra.cholesky!
— Methodsolver = cholesky!(solver::CudssSolver{T}, A::CuSparseMatrixCSR{T,Cint})
Compute the LLᴴ factorization of a sparse matrix A
on an NVIDIA GPU, reusing the symbolic factorization stored in solver
. The type T
can be Float32
, Float64
, ComplexF32
or ComplexF64
.
using CUDA, CUDA.CUSPARSE
+Generic API · CUDSS.jl LLᵀ and LLᴴ
LinearAlgebra.cholesky
— Methodsolver = cholesky(A::CuSparseMatrixCSR{T,Cint}; view::Char='F')
Compute the LLᴴ factorization of a sparse matrix A
on an NVIDIA GPU. The type T
can be Float32
, Float64
, ComplexF32
or ComplexF64
.
Input argument
A
: a sparse Hermitian positive definite matrix stored in the CuSparseMatrixCSR
format.
Keyword argument
*view
: A character that specifies which triangle of the sparse matrix is provided. Possible options are L
for the lower triangle, U
for the upper triangle, and F
for the full matrix.
Output argument
solver
: Opaque structure CudssSolver
that stores the factors of the LLᴴ decomposition.
sourceLinearAlgebra.cholesky!
— Methodsolver = cholesky!(solver::CudssSolver{T}, A::CuSparseMatrixCSR{T,Cint})
Compute the LLᴴ factorization of a sparse matrix A
on an NVIDIA GPU, reusing the symbolic factorization stored in solver
. The type T
can be Float32
, Float64
, ComplexF32
or ComplexF64
.
sourceusing CUDA, CUDA.CUSPARSE
using CUDSS
using LinearAlgebra
using SparseArrays
@@ -33,7 +33,7 @@
R_gpu = C_gpu - ( CuSparseMatrixCSR(A_cpu) + Diagonal(d_gpu) ) * X_gpu
norm(R_gpu)
Note If we only store one triangle of A_gpu
, we can also use the wrappers Symmetric
and Hermitian
instead of using the keyword argument view
in cholesky
. For real matrices, both wrappers are allowed but only Hermitian
can be used for complex matrices.
H_gpu = Hermitian(A_gpu, :U)
-F = cholesky(H_gpu)
LDLᵀ and LDLᴴ
LinearAlgebra.ldlt
— Methodsolver = ldlt(A::CuSparseMatrixCSR{T,Cint}; view::Char='F')
Compute the LDLᴴ factorization of a sparse matrix A
on an NVIDIA GPU. The type T
can be Float32
, Float64
, ComplexF32
or ComplexF64
.
Input argument
A
: a sparse Hermitian matrix stored in the CuSparseMatrixCSR
format.
Keyword argument
*view
: A character that specifies which triangle of the sparse matrix is provided. Possible options are L
for the lower triangle, U
for the upper triangle, and F
for the full matrix.
Output argument
solver
: Opaque structure CudssSolver
that stores the factors of the LDLᴴ decomposition.
sourceLinearAlgebra.ldlt!
— Methodsolver = ldlt!(solver::CudssSolver{T}, A::CuSparseMatrixCSR{T,Cint})
Compute the LDLᴴ factorization of a sparse matrix A
on an NVIDIA GPU, reusing the symbolic factorization stored in solver
. The type T
can be Float32
, Float64
, ComplexF32
or ComplexF64
.
sourceusing CUDA, CUDA.CUSPARSE
+F = cholesky(H_gpu)
LDLᵀ and LDLᴴ
LinearAlgebra.ldlt
— Methodsolver = ldlt(A::CuSparseMatrixCSR{T,Cint}; view::Char='F')
Compute the LDLᴴ factorization of a sparse matrix A
on an NVIDIA GPU. The type T
can be Float32
, Float64
, ComplexF32
or ComplexF64
.
Input argument
A
: a sparse Hermitian matrix stored in the CuSparseMatrixCSR
format.
Keyword argument
*view
: A character that specifies which triangle of the sparse matrix is provided. Possible options are L
for the lower triangle, U
for the upper triangle, and F
for the full matrix.
Output argument
solver
: Opaque structure CudssSolver
that stores the factors of the LDLᴴ decomposition.
sourceLinearAlgebra.ldlt!
— Methodsolver = ldlt!(solver::CudssSolver{T}, A::CuSparseMatrixCSR{T,Cint})
Compute the LDLᴴ factorization of a sparse matrix A
on an NVIDIA GPU, reusing the symbolic factorization stored in solver
. The type T
can be Float32
, Float64
, ComplexF32
or ComplexF64
.
sourceusing CUDA, CUDA.CUSPARSE
using CUDSS
using LinearAlgebra
using SparseArrays
@@ -67,7 +67,7 @@
R_gpu = C_gpu - ( CuSparseMatrixCSR(A_cpu) + Diagonal(d_gpu) ) * X_gpu
norm(R_gpu)
Note If we only store one triangle of A_gpu
, we can also use the wrappers Symmetric
and Hermitian
instead of using the keyword argument view
in ldlt
. For real matrices, both wrappers are allowed but only Hermitian
can be used for complex matrices.
S_gpu = Symmetric(A_gpu, :L)
-F = ldlt(S_gpu)
LU
LinearAlgebra.lu
— Methodsolver = lu(A::CuSparseMatrixCSR{T,Cint})
Compute the LU factorization of a sparse matrix A
on an NVIDIA GPU. The type T
can be Float32
, Float64
, ComplexF32
or ComplexF64
.
Input argument
A
: a sparse square matrix stored in the CuSparseMatrixCSR
format.
Output argument
solver
: an opaque structure CudssSolver
that stores the factors of the LU decomposition.
sourceLinearAlgebra.lu!
— Methodsolver = lu!(solver::CudssSolver{T}, A::CuSparseMatrixCSR{T,Cint})
Compute the LU factorization of a sparse matrix A
on an NVIDIA GPU, reusing the symbolic factorization stored in solver
. The type T
can be Float32
, Float64
, ComplexF32
or ComplexF64
.
sourceusing CUDA, CUDA.CUSPARSE
+F = ldlt(S_gpu)
LU
LinearAlgebra.lu
— Methodsolver = lu(A::CuSparseMatrixCSR{T,Cint})
Compute the LU factorization of a sparse matrix A
on an NVIDIA GPU. The type T
can be Float32
, Float64
, ComplexF32
or ComplexF64
.
Input argument
A
: a sparse square matrix stored in the CuSparseMatrixCSR
format.
Output argument
solver
: an opaque structure CudssSolver
that stores the factors of the LU decomposition.
sourceLinearAlgebra.lu!
— Methodsolver = lu!(solver::CudssSolver{T}, A::CuSparseMatrixCSR{T,Cint})
Compute the LU factorization of a sparse matrix A
on an NVIDIA GPU, reusing the symbolic factorization stored in solver
. The type T
can be Float32
, Float64
, ComplexF32
or ComplexF64
.
sourceusing CUDA, CUDA.CUSPARSE
using CUDSS
using LinearAlgebra
using SparseArrays
@@ -96,4 +96,4 @@
ldiv!(x_gpu, F, c_gpu)
r_gpu = c_gpu - A_gpu * x_gpu
-norm(r_gpu)
Settings
This document was generated with Documenter.jl version 1.8.0 on Thursday 12 December 2024. Using Julia version 1.11.2.
+norm(r_gpu)