Sparse Arrays
줄리아는 희소 벡터와 sparse matrices를 SparseArrays
표준 라이브러리 모듈에서 지원합니다. 희소 배열은 충분한 수의 0을 포함하고 있어, 이를 특별한 데이터 구조에 저장하면 밀집 배열에 비해 공간과 실행 시간에서 절약을 할 수 있습니다.
다양한 희소 저장 유형, 다차원 희소 배열 등을 구현하는 외부 패키지는 Noteworthy External Sparse Packages에서 찾을 수 있습니다.
Compressed Sparse Column (CSC) Sparse Matrix Storage
줄리아에서 희소 행렬은 Compressed Sparse Column (CSC) format에 저장됩니다. 줄리아 희소 행렬의 타입은 SparseMatrixCSC{Tv,Ti}
이며, 여기서 Tv
는 저장된 값의 타입이고, Ti
는 열 포인터와 행 인덱스를 저장하기 위한 정수 타입입니다. SparseMatrixCSC
의 내부 표현은 다음과 같습니다:
struct SparseMatrixCSC{Tv,Ti<:Integer} <: AbstractSparseMatrixCSC{Tv,Ti}
m::Int # Number of rows
n::Int # Number of columns
colptr::Vector{Ti} # Column j is in colptr[j]:(colptr[j+1]-1)
rowval::Vector{Ti} # Row indices of stored values
nzval::Vector{Tv} # Stored values, typically nonzeros
end
압축 희소 열 저장 방식은 희소 행렬의 열에 있는 요소에 쉽게 빠르게 접근할 수 있게 해주지만, 행으로 희소 행렬에 접근하는 것은 상당히 느립니다. CSC 구조에서 이전에 저장되지 않은 항목을 한 번에 하나씩 삽입하는 작업은 느린 경향이 있습니다. 이는 삽입 지점 이후의 희소 행렬의 모든 요소를 한 칸씩 이동해야 하기 때문입니다.
모든 희소 행렬에 대한 연산은 성능을 위해 CSC 데이터 구조를 활용하고 비싼 연산을 피하도록 신중하게 구현됩니다.
CSC 형식의 데이터가 다른 애플리케이션이나 라이브러리에서 온 경우, 이를 줄리아에 가져오려면 1 기반 인덱싱을 사용해야 합니다. 모든 열의 행 인덱스는 정렬되어 있어야 하며, 그렇지 않으면 행렬이 잘못 표시됩니다. SparseMatrixCSC
객체에 정렬되지 않은 행 인덱스가 포함되어 있는 경우, 이를 정렬하는 빠른 방법 중 하나는 이중 전치(transpose)를 수행하는 것입니다. 전치 연산은 지연(lazy) 방식으로 수행되므로, 각 전치를 구체화하기 위해 복사본을 만들어야 합니다.
일부 응용 프로그램에서는 SparseMatrixCSC
에 명시적인 제로 값을 저장하는 것이 편리합니다. 이러한 값은 Base
의 함수에서 허용되지만, 변형 작업에서 보존된다는 보장은 없습니다. 명시적으로 저장된 제로는 많은 루틴에서 구조적 비제로로 처리됩니다. nnz
함수는 비구조적 제로를 포함하여 희소 데이터 구조에 명시적으로 저장된 요소의 수를 반환합니다. 정확한 수치 비제로의 수를 세려면 count(!iszero, x)
를 사용하여 희소 행렬의 모든 저장된 요소를 검사하십시오. dropzeros
및 제자리에서 사용하는 dropzeros!
는 희소 행렬에서 저장된 제로를 제거하는 데 사용할 수 있습니다.
julia> A = sparse([1, 1, 2, 3], [1, 3, 2, 3], [0, 1, 2, 0])
3×3 SparseMatrixCSC{Int64, Int64} with 4 stored entries:
0 ⋅ 1
⋅ 2 ⋅
⋅ ⋅ 0
julia> dropzeros(A)
3×3 SparseMatrixCSC{Int64, Int64} with 2 stored entries:
⋅ ⋅ 1
⋅ 2 ⋅
⋅ ⋅ ⋅
Sparse Vector Storage
희소 벡터는 희소 행렬의 압축 희소 열 형식과 유사한 방식으로 저장됩니다. Julia에서 희소 벡터는 SparseVector{Tv,Ti}
유형을 가지며, 여기서 Tv
는 저장된 값의 유형이고 Ti
는 인덱스의 정수 유형입니다. 내부 표현은 다음과 같습니다:
struct SparseVector{Tv,Ti<:Integer} <: AbstractSparseVector{Tv,Ti}
n::Int # Length of the sparse vector
nzind::Vector{Ti} # Indices of stored values
nzval::Vector{Tv} # Stored values, typically nonzeros
end
SparseMatrixCSC
에 관해서는 SparseVector
유형이 명시적으로 저장된 0을 포함할 수 있습니다. (참조: Sparse Matrix Storage.)
Sparse Vector and Matrix Constructors
희소 배열을 만드는 가장 간단한 방법은 Julia가 밀집 배열을 다루기 위해 제공하는 zeros
함수와 동등한 함수를 사용하는 것입니다. 대신 희소 배열을 생성하려면 sp
접두사를 붙여 같은 이름을 사용할 수 있습니다:
julia> spzeros(3)
3-element SparseVector{Float64, Int64} with 0 stored entries
sparse
함수는 종종 희소 배열을 구성하는 유용한 방법입니다. 예를 들어, 희소 행렬을 구성하기 위해 행 인덱스의 벡터 I
, 열 인덱스의 벡터 J
, 저장된 값의 벡터 V
를 입력할 수 있습니다(이는 COO (coordinate) format로도 알려져 있습니다). 그런 다음 sparse(I,J,V)
는 S[I[k], J[k]] = V[k]
가 되도록 희소 행렬을 구성합니다. 동등한 희소 벡터 생성자는 sparsevec
로, (행) 인덱스 벡터 I
와 저장된 값이 있는 벡터 V
를 사용하여 희소 벡터 R
을 구성하며, R[I[k]] = V[k]
가 됩니다.
julia> I = [1, 4, 3, 5]; J = [4, 7, 18, 9]; V = [1, 2, -5, 3];
julia> S = sparse(I,J,V)
5×18 SparseMatrixCSC{Int64, Int64} with 4 stored entries:
⎡⠀⠈⠀⠀⠀⠀⠀⠀⢀⎤
⎣⠀⠀⠀⠂⡀⠀⠀⠀⠀⎦
julia> R = sparsevec(I,V)
5-element SparseVector{Int64, Int64} with 4 stored entries:
[1] = 1
[3] = -5
[4] = 2
[5] = 3
sparse
및 sparsevec
함수의 역은 findnz
로, 이는 희소 배열을 생성하는 데 사용된 입력값(영으로 저장된 항목 포함)을 검색합니다. findall(!iszero, x)
는 x
의 비영 항목의 카르테시안 인덱스를 반환합니다(영으로 저장된 항목은 포함하지 않음).
julia> findnz(S)
([1, 4, 5, 3], [4, 7, 9, 18], [1, 2, 3, -5])
julia> findall(!iszero, S)
4-element Vector{CartesianIndex{2}}:
CartesianIndex(1, 4)
CartesianIndex(4, 7)
CartesianIndex(5, 9)
CartesianIndex(3, 18)
julia> findnz(R)
([1, 3, 4, 5], [1, -5, 2, 3])
julia> findall(!iszero, R)
4-element Vector{Int64}:
1
3
4
5
희소 배열을 생성하는 또 다른 방법은 sparse
함수를 사용하여 밀집 배열을 희소 배열로 변환하는 것입니다:
julia> sparse(Matrix(1.0I, 5, 5))
5×5 SparseMatrixCSC{Float64, Int64} with 5 stored entries:
1.0 ⋅ ⋅ ⋅ ⋅
⋅ 1.0 ⋅ ⋅ ⋅
⋅ ⋅ 1.0 ⋅ ⋅
⋅ ⋅ ⋅ 1.0 ⋅
⋅ ⋅ ⋅ ⋅ 1.0
julia> sparse([1.0, 0.0, 1.0])
3-element SparseVector{Float64, Int64} with 2 stored entries:
[1] = 1.0
[3] = 1.0
다른 방향으로 가려면 Array
생성자를 사용할 수 있습니다. issparse
함수는 행렬이 희소한지 쿼리하는 데 사용할 수 있습니다.
julia> issparse(spzeros(5))
true
Sparse matrix operations
희소 행렬에 대한 산술 연산은 밀집 행렬에서와 동일하게 작동합니다. 희소 행렬의 인덱싱, 할당 및 연결은 밀집 행렬과 동일한 방식으로 작동합니다. 인덱싱 작업, 특히 할당은 한 번에 하나의 요소를 처리할 때 비용이 많이 듭니다. 많은 경우 희소 행렬을 findnz
형식의 (I,J,V)
로 변환하고, 밀집 벡터 (I,J,V)
에서 값이나 구조를 조작한 다음, 희소 행렬을 재구성하는 것이 더 나을 수 있습니다.
Correspondence of dense and sparse methods
다음 표는 희소 행렬의 내장 메서드와 밀집 행렬 유형의 해당 메서드 간의 대응 관계를 제공합니다. 일반적으로 희소 행렬을 생성하는 메서드는 결과 행렬이 주어진 희소 행렬 S
와 동일한 희소성 패턴을 따르거나, 결과 희소 행렬이 밀도 d
를 가지며, 즉 각 행렬 요소가 0이 아닐 확률이 d
인 점에서 밀집 행렬의 대응 메서드와 다릅니다.
자세한 내용은 표준 라이브러리 참조의 Sparse Vectors and Matrices 섹션에서 확인할 수 있습니다.
Sparse | Dense | Description |
---|---|---|
spzeros(m,n) | zeros(m,n) | Creates a m-by-n matrix of zeros. (spzeros(m,n) is empty.) |
sparse(I,n,n) | Matrix(I,n,n) | Creates a n-by-n identity matrix. |
sparse(A) | Array(S) | Interconverts between dense and sparse formats. |
sprand(m,n,d) | rand(m,n) | Creates a m-by-n random matrix (of density d) with iid non-zero elements distributed uniformly on the half-open interval $[0, 1)$. |
sprandn(m,n,d) | randn(m,n) | Creates a m-by-n random matrix (of density d) with iid non-zero elements distributed according to the standard normal (Gaussian) distribution. |
sprandn(rng,m,n,d) | randn(rng,m,n) | Creates a m-by-n random matrix (of density d) with iid non-zero elements generated with the rng random number generator |
Sparse Linear Algebra
희소 행렬 솔버는 SuiteSparse에서 함수를 호출합니다. 다음과 같은 분해가 가능합니다:
Type | Description |
---|---|
CHOLMOD.Factor | Cholesky and LDLt factorizations |
UMFPACK.UmfpackLU | LU factorization |
SPQR.QRSparse | QR factorization |
이러한 인수분해는 Sparse Linear Algebra API section에서 더 자세히 설명되어 있습니다:
SparseArrays API
SparseArrays.AbstractSparseArray
— TypeAbstractSparseArray{Tv,Ti,N}
Tv
유형의 요소와 Ti
유형의 인덱스를 가진 N
차원 희소 배열(또는 배열 유사 유형)의 슈퍼타입입니다. SparseMatrixCSC
, SparseVector
및 SuiteSparse.CHOLMOD.Sparse
는 이의 하위 유형입니다.
SparseArrays.AbstractSparseVector
— TypeAbstractSparseVector{Tv,Ti}
타입 Tv
와 인덱스 타입 Ti
의 요소를 가진 1차원 희소 배열(또는 배열 유사 타입)의 슈퍼타입. AbstractSparseArray{Tv,Ti,1}
의 별칭입니다.
SparseArrays.AbstractSparseMatrix
— TypeAbstractSparseMatrix{Tv,Ti}
Tv
유형의 요소와 Ti
인덱스 유형을 가진 2차원 희소 배열(또는 배열 유사 유형)의 슈퍼타입입니다. AbstractSparseArray{Tv,Ti,2}
의 별칭입니다.
SparseArrays.SparseVector
— TypeSparseVector{Tv,Ti<:Integer} <: AbstractSparseVector{Tv,Ti}
희소 벡터를 저장하기 위한 벡터 유형입니다. 벡터의 길이, 정렬된 비영(非零) 인덱스 벡터, 비영(非零) 값 벡터를 전달하여 생성할 수 있습니다.
예를 들어, 벡터 [5, 6, 0, 7]
는 다음과 같이 표현할 수 있습니다.
SparseVector(4, [1, 2, 4], [5, 6, 7])
이는 인덱스 1의 요소가 5, 인덱스 2의 요소가 6, 인덱스 3의 요소가 zero(Int)
, 인덱스 4의 요소가 7임을 나타냅니다.
dense
벡터에서 직접 희소 벡터를 생성하는 것이 더 편리할 수 있으며, sparse
를 사용하여 다음과 같이 할 수 있습니다.
sparse([5, 6, 0, 7])
이는 동일한 희소 벡터를 생성합니다.
SparseArrays.SparseMatrixCSC
— TypeSparseMatrixCSC{Tv,Ti<:Integer} <: AbstractSparseMatrixCSC{Tv,Ti}
희소 행렬을 압축 희소 열 형식으로 저장하기 위한 행렬 유형입니다. SparseMatrixCSC를 구성하는 표준 방법은 sparse
함수를 사용하는 것입니다. 또한 spzeros
, spdiagm
및 sprand
를 참조하십시오.
SparseArrays.sparse
— Functionsparse(A)
AbstractMatrix A
를 희소 행렬로 변환합니다.
예제
julia> A = Matrix(1.0I, 3, 3)
3×3 Matrix{Float64}:
1.0 0.0 0.0
0.0 1.0 0.0
0.0 0.0 1.0
julia> sparse(A)
3×3 SparseMatrixCSC{Float64, Int64} with 3 stored entries:
1.0 ⋅ ⋅
⋅ 1.0 ⋅
⋅ ⋅ 1.0
sparse(I, J, V,[ m, n, combine])
S
라는 희소 행렬을 m x n
크기로 생성하여 S[I[k], J[k]] = V[k]
가 되도록 합니다. combine
함수는 중복을 결합하는 데 사용됩니다. m
과 n
이 지정되지 않으면 각각 maximum(I)
와 maximum(J)
로 설정됩니다. combine
함수가 제공되지 않으면, combine
은 기본적으로 +
로 설정되지만, V
의 요소가 Boolean인 경우 combine
은 기본적으로 |
로 설정됩니다. I
의 모든 요소는 1 <= I[k] <= m
을 만족해야 하며, J
의 모든 요소는 1 <= J[k] <= n
을 만족해야 합니다. (I
, J
, V
)의 수치적 0은 구조적 비제로로 유지됩니다; 수치적 0을 제거하려면 dropzeros!
를 사용하십시오.
추가 문서 및 전문가 드라이버는 SparseArrays.sparse!
를 참조하십시오.
예제
julia> Is = [1; 2; 3];
julia> Js = [1; 2; 3];
julia> Vs = [1; 2; 3];
julia> sparse(Is, Js, Vs)
3×3 SparseMatrixCSC{Int64, Int64} with 3 stored entries:
1 ⋅ ⋅
⋅ 2 ⋅
⋅ ⋅ 3
SparseArrays.sparse!
— Functionsparse!(I::AbstractVector{Ti}, J::AbstractVector{Ti}, V::AbstractVector{Tv},
m::Integer, n::Integer, combine, klasttouch::Vector{Ti},
csrrowptr::Vector{Ti}, csrcolval::Vector{Ti}, csrnzval::Vector{Tv},
[csccolptr::Vector{Ti}], [cscrowval::Vector{Ti}, cscnzval::Vector{Tv}] ) where {Tv,Ti<:Integer}
sparse
의 부모 및 전문가 드라이버; 기본 사용법은 sparse
를 참조하십시오. 이 메서드는 사용자가 아래에 설명된 대로 sparse
의 중간 객체와 결과를 위한 미리 할당된 저장소를 제공할 수 있도록 합니다. 이 기능은 좌표 표현에서 SparseMatrixCSC
로의 보다 효율적인 연속 구성을 가능하게 하며, 추가 비용 없이 결과의 전치의 정렬되지 않은 열 표현을 추출할 수 있게 합니다.
이 메서드는 세 가지 주요 단계로 구성됩니다: (1) 제공된 좌표 표현을 반복 항목을 포함한 정렬되지 않은 행 CSR 형식으로 카운팅 정렬합니다. (2) CSR 형식을 통해 스윕하면서 원하는 CSC 형식의 열 포인터 배열을 동시에 계산하고, 반복 항목을 감지하며, 반복 항목이 결합된 CSR 형식을 다시 포장합니다. 이 단계는 반복 항목이 없는 정렬되지 않은 행 CSR 형식을 생성합니다. (3) 이전 CSR 형식을 반복 항목이 없는 완전히 정렬된 CSC 형식으로 카운팅 정렬합니다.
입력 배열 csrrowptr
, csrcolval
, 및 csrnzval
은 중간 CSR 형식의 저장소를 구성하며 length(csrrowptr) >= m + 1
, length(csrcolval) >= length(I)
, 및 length(csrnzval >= length(I))
를 요구합니다. 입력 배열 klasttouch
, 두 번째 단계의 작업 공간,은 length(klasttouch) >= n
을 요구합니다. 선택적 입력 배열 csccolptr
, cscrowval
, 및 cscnzval
은 반환된 CSC 형식 S
의 저장소를 구성합니다. 필요할 경우, 이들은 length(csccolptr) = n + 1
, length(cscrowval) = nnz(S)
및 length(cscnzval) = nnz(S)
를 만족하도록 자동으로 크기가 조정됩니다. 따라서 nnz(S)
가 처음에 알려지지 않은 경우, 적절한 유형의 빈 벡터(Vector{Ti}()
및 Vector{Tv}()
)를 전달하거나 cscrowval
및 cscnzval
을 무시하고 sparse!
메서드를 호출하면 충분합니다.
반환 시, csrrowptr
, csrcolval
, 및 csrnzval
은 결과의 전치의 정렬되지 않은 열 표현을 포함합니다.
입력 배열의 저장소(I
, J
, V
)를 출력 배열(csccolptr
, cscrowval
, cscnzval
)에 재사용할 수 있습니다. 예를 들어, sparse!(I, J, V, csrrowptr, csrcolval, csrnzval, I, J, V)
를 호출할 수 있습니다. 이들은 위의 조건을 만족하도록 크기가 조정됩니다.
효율성을 위해, 이 메서드는 1 <= I[k] <= m
및 1 <= J[k] <= n
을 넘어서는 인수 검사를 수행하지 않습니다. 주의해서 사용하십시오. --check-bounds=yes
로 테스트하는 것이 좋습니다.
이 메서드는 O(m, n, length(I))
시간에 실행됩니다. F. Gustavson의 "Two fast algorithms for sparse matrices: multiplication and permuted transposition," ACM TOMS 4(3), 250-269 (1978)에서 설명된 HALFPERM 알고리즘이 이 메서드의 카운팅 정렬 쌍 사용에 영감을 주었습니다.
SparseArrays.sparse!(I, J, V, [m, n, combine]) -> SparseMatrixCSC
입력 벡터(I
, J
, V
)를 최종 행렬 저장소에 재사용하는 sparse!
의 변형입니다. 생성 후 입력 벡터는 행렬 버퍼를 참조하게 되며, S.colptr === I
, S.rowval === J
, S.nzval === V
가 성립하며, 필요에 따라 resize!
됩니다.
일부 작업 버퍼는 여전히 할당됩니다. 구체적으로, 이 메서드는 sparse!(I, J, V, m, n, combine, klasttouch, csrrowptr, csrcolval, csrnzval, csccolptr, cscrowval, cscnzval)
에 대한 편의 래퍼로, 이 메서드는 적절한 크기의 klasttouch
, csrrowptr
, csrcolval
, 및 csrnzval
을 할당하지만, csccolptr
, cscrowval
, 및 cscnzval
에 대해 I
, J
, 및 V
를 재사용합니다.
인수 m
, n
, 및 combine
의 기본값은 각각 maximum(I)
, maximum(J)
, 및 +
입니다.
이 메서드는 Julia 버전 1.10 이상이 필요합니다.
SparseArrays.sparsevec
— Functionsparsevec(I, V, [m, combine])
길이 m
인 희소 벡터 S
를 생성하여 S[I[k]] = V[k]
가 되도록 합니다. 중복 항목은 combine
함수를 사용하여 결합되며, combine
인수가 제공되지 않으면 기본값은 +
입니다. 단, V
의 요소가 부울인 경우 combine
의 기본값은 |
입니다.
예제
julia> II = [1, 3, 3, 5]; V = [0.1, 0.2, 0.3, 0.2];
julia> sparsevec(II, V)
5-element SparseVector{Float64, Int64} with 3 stored entries:
[1] = 0.1
[3] = 0.5
[5] = 0.2
julia> sparsevec(II, V, 8, -)
8-element SparseVector{Float64, Int64} with 3 stored entries:
[1] = 0.1
[3] = -0.1
[5] = 0.2
julia> sparsevec([1, 3, 1, 2, 2], [true, true, false, false, false])
3-element SparseVector{Bool, Int64} with 3 stored entries:
[1] = 1
[2] = 0
[3] = 1
sparsevec(d::Dict, [m])
길이 m
의 희소 벡터를 생성하며, 비영(非零) 인덱스는 사전의 키이고, 비영 값은 사전의 값입니다.
예제
julia> sparsevec(Dict(1 => 3, 2 => 2))
2-element SparseVector{Int64, Int64} with 2 stored entries:
[1] = 3
[2] = 2
sparsevec(A)
벡터 A
를 길이 m
의 희소 벡터로 변환합니다.
예제
julia> sparsevec([1.0, 2.0, 0.0, 0.0, 3.0, 0.0])
6-element SparseVector{Float64, Int64} with 3 stored entries:
[1] = 1.0
[2] = 2.0
[5] = 3.0
Base.similar
— Method유사한(A::AbstractSparseMatrixCSC{Tv,Ti}, [::Type{TvNew}, ::Type{TiNew}, m::Integer, n::Integer]) where {Tv,Ti}
주어진 요소 유형, 인덱스 유형 및 크기를 기반으로 주어진 소스 SparseMatrixCSC
로 초기화되지 않은 가변 배열을 생성합니다. 새로운 희소 행렬은 원래 희소 행렬의 구조를 유지하지만 출력 행렬의 차원이 다를 경우를 제외합니다.
출력 행렬은 입력과 동일한 위치에 0이 있지만 비영(非零) 위치에는 초기화되지 않은 값이 있습니다.
SparseArrays.issparse
— Functionissparse(S)
S
가 희소 행렬이면 true
를 반환하고, 그렇지 않으면 false
를 반환합니다.
예제
julia> sv = sparsevec([1, 4], [2.3, 2.2], 10)
10-element SparseVector{Float64, Int64} with 2 stored entries:
[1] = 2.3
[4] = 2.2
julia> issparse(sv)
true
julia> issparse(Array(sv))
false
SparseArrays.nnz
— Functionnnz(A)
희소 배열에 저장된(채워진) 요소의 수를 반환합니다.
예제
julia> A = sparse(2I, 3, 3)
3×3 SparseMatrixCSC{Int64, Int64} with 3 stored entries:
2 ⋅ ⋅
⋅ 2 ⋅
⋅ ⋅ 2
julia> nnz(A)
3
SparseArrays.findnz
— Functionfindnz(A::SparseMatrixCSC)
희소 행렬 A
에 저장된 ("구조적으로 비영(zero)") 값의 행 및 열 인덱스인 튜플 (I, J, V)
를 반환합니다. 여기서 I
와 J
는 인덱스이고, V
는 값의 벡터입니다.
예제
julia> A = sparse([1 2 0; 0 0 3; 0 4 0])
3×3 SparseMatrixCSC{Int64, Int64} with 4 stored entries:
1 2 ⋅
⋅ ⋅ 3
⋅ 4 ⋅
julia> findnz(A)
([1, 1, 3, 2], [1, 2, 2, 3], [1, 2, 4, 3])
SparseArrays.spzeros
— Functionspzeros([type,]m[,n])
길이 m
의 희소 벡터 또는 크기 m x n
의 희소 행렬을 생성합니다. 이 희소 배열은 비영 값이 포함되지 않습니다. 생성 중에 비영 값에 대한 저장소는 할당되지 않습니다. 타입이 지정되지 않으면 기본값은 Float64
입니다.
예제
julia> spzeros(3, 3)
3×3 SparseMatrixCSC{Float64, Int64} with 0 stored entries:
⋅ ⋅ ⋅
⋅ ⋅ ⋅
⋅ ⋅ ⋅
julia> spzeros(Float32, 4)
4-element SparseVector{Float32, Int64} with 0 stored entries
spzeros([type], I::AbstractVector, J::AbstractVector, [m, n])
구조적 제로가 있는 차원 m x n
의 희소 행렬 S
를 생성합니다. S[I[k], J[k]]
에 구조적 제로가 있습니다.
이 방법은 행렬의 희소성 패턴을 구성하는 데 사용될 수 있으며, 예를 들어 sparse(I, J, zeros(length(I)))
를 사용하는 것보다 더 효율적입니다.
추가 문서 및 전문가 드라이버에 대한 내용은 SparseArrays.spzeros!
를 참조하십시오.
이 방법은 Julia 버전 1.10 이상이 필요합니다.
SparseArrays.spzeros!
— Functionspzeros!(::Type{Tv}, I::AbstractVector{Ti}, J::AbstractVector{Ti}, m::Integer, n::Integer,
klasttouch::Vector{Ti}, csrrowptr::Vector{Ti}, csrcolval::Vector{Ti},
[csccolptr::Vector{Ti}], [cscrowval::Vector{Ti}, cscnzval::Vector{Tv}]) where {Tv,Ti<:Integer}
spzeros(I, J)
의 부모이자 전문가 드라이버로, 사용자가 중간 객체에 대한 미리 할당된 저장소를 제공할 수 있도록 합니다. 이 메서드는 spzeros
에 대해 SparseArrays.sparse!
가 sparse
에 해당하는 것과 같습니다. 필요한 버퍼 길이에 대한 자세한 내용은 SparseArrays.sparse!
문서를 참조하십시오.
이 메서드는 Julia 버전 1.10 이상이 필요합니다.
SparseArrays.spzeros!(::Type{Tv}, I, J, [m, n]) -> SparseMatrixCSC{Tv}
입력 벡터 I
와 J
를 최종 행렬 저장소에 재사용하는 spzeros!
의 변형입니다. 생성 후 입력 벡터는 행렬 버퍼를 별칭으로 사용하게 되며, S.colptr === I
및 S.rowval === J
가 성립하고 필요에 따라 resize!
됩니다.
일부 작업 버퍼는 여전히 할당될 것임에 유의하십시오. 구체적으로, 이 메서드는 spzeros!(Tv, I, J, m, n, klasttouch, csrrowptr, csrcolval, csccolptr, cscrowval)
에 대한 편의 래퍼이며, 이 메서드는 적절한 크기의 klasttouch
, csrrowptr
, 및 csrcolval
을 할당하지만 csccolptr
및 cscrowval
에 대해 I
와 J
를 재사용합니다.
인수 m
과 n
은 각각 maximum(I)
및 maximum(J)
로 기본값이 설정됩니다.
이 메서드는 Julia 버전 1.10 이상이 필요합니다.
SparseArrays.spdiagm
— Functionspdiagm(kv::Pair{<:Integer,<:AbstractVector}...)
spdiagm(m::Integer, n::Integer, kv::Pair{<:Integer,<:AbstractVector}...)
Pair
의 벡터와 대각선으로부터 희소 대각 행렬을 구성합니다. 각 벡터 kv.second
는 kv.first
대각선에 배치됩니다. 기본적으로 행렬은 정사각형이며 그 크기는 kv
로부터 유추되지만, 필요에 따라 0으로 패딩된 비정사각형 크기 m
×n
을 첫 번째 인수로 m,n
을 전달하여 지정할 수 있습니다.
예제
julia> spdiagm(-1 => [1,2,3,4], 1 => [4,3,2,1])
5×5 SparseMatrixCSC{Int64, Int64} with 8 stored entries:
⋅ 4 ⋅ ⋅ ⋅
1 ⋅ 3 ⋅ ⋅
⋅ 2 ⋅ 2 ⋅
⋅ ⋅ 3 ⋅ 1
⋅ ⋅ ⋅ 4 ⋅
spdiagm(v::AbstractVector)
spdiagm(m::Integer, n::Integer, v::AbstractVector)
벡터의 요소를 대각선 요소로 갖는 희소 행렬을 생성합니다. 기본적으로(m
과 n
이 주어지지 않은 경우) 행렬은 정사각형이며 크기는 length(v)
로 주어지지만, m
과 n
을 첫 번째 인수로 전달하여 비정사각형 크기 m
×n
을 지정할 수 있습니다.
이 함수는 최소한 Julia 1.6이 필요합니다.
예제
julia> spdiagm([1,2,3])
3×3 SparseMatrixCSC{Int64, Int64} with 3 stored entries:
1 ⋅ ⋅
⋅ 2 ⋅
⋅ ⋅ 3
julia> spdiagm(sparse([1,0,3]))
3×3 SparseMatrixCSC{Int64, Int64} with 2 stored entries:
1 ⋅ ⋅
⋅ ⋅ ⋅
⋅ ⋅ 3
SparseArrays.sparse_hcat
— Functionsparse_hcat(A...)
차원 2를 따라 연결합니다. SparseMatrixCSC 객체를 반환합니다.
이 메서드는 Julia 1.8에 추가되었습니다. LinearAlgebra.jl의 특수화된 "희소" 행렬 유형과의 연결이 SparseArray 인수가 없더라도 자동으로 희소 출력을 생성하는 이전 연결 동작을 모방합니다.
SparseArrays.sparse_vcat
— Functionsparse_vcat(A...)
차원 1을 따라 연결합니다. SparseMatrixCSC 객체를 반환합니다.
이 메서드는 Julia 1.8에 추가되었습니다. LinearAlgebra.jl의 특수화된 "희소" 행렬 유형과의 연결이 SparseArray 인수가 없더라도 자동으로 희소 출력을 생성하는 이전 연결 동작을 모방합니다.
SparseArrays.sparse_hvcat
— Functionsparse_hvcat(rows::Tuple{Vararg{Int}}, values...)
한 번의 호출로 희소 수평 및 수직 연결. 이 함수는 블록 행렬 구문에 대해 호출됩니다. 첫 번째 인수는 각 블록 행에서 연결할 인수의 수를 지정합니다.
이 메서드는 Julia 1.8에 추가되었습니다. 이는 LinearAlgebra.jl의 특수 "희소" 행렬 유형과의 연결이 SparseArray 인수가 없더라도 자동으로 희소 출력을 생성하는 이전 연결 동작을 모방합니다.
SparseArrays.blockdiag
— Functionblockdiag(A...)
행렬을 블록 대각선으로 연결합니다. 현재는 희소 행렬에 대해서만 구현되어 있습니다.
예제
julia> blockdiag(sparse(2I, 3, 3), sparse(4I, 2, 2))
5×5 SparseMatrixCSC{Int64, Int64} with 5 stored entries:
2 ⋅ ⋅ ⋅ ⋅
⋅ 2 ⋅ ⋅ ⋅
⋅ ⋅ 2 ⋅ ⋅
⋅ ⋅ ⋅ 4 ⋅
⋅ ⋅ ⋅ ⋅ 4
SparseArrays.sprand
— Functionsprand([rng],[T::Type],m,[n],p::AbstractFloat)
sprand([rng],m,[n],p::AbstractFloat,[rfn=rand])
길이가 m
인 랜덤 희소 벡터 또는 m
x n
희소 행렬을 생성합니다. 이 행렬에서 어떤 요소가 0이 아닐 확률은 독립적으로 p
로 주어지며 (따라서 비영 값의 평균 밀도도 정확히 p
입니다). 선택적 rng
인자는 난수 생성기를 지정하며, Random Numbers를 참조하십시오. 선택적 T
인자는 요소 유형을 지정하며, 기본값은 Float64
입니다.
기본적으로 비영 값은 rand
함수를 사용하여 균일 분포에서 샘플링됩니다. 즉, rand(T)
또는 rng
가 제공된 경우 rand(rng, T)
로 샘플링됩니다. 기본 T=Float64
의 경우, 이는 [0,1)
에서 균일하게 샘플링된 비영 값에 해당합니다.
다른 분포에서 비영 값을 샘플링하려면 rand
대신 사용자 정의 rfn
함수를 전달하면 됩니다. 이 함수는 원하는 분포에서 샘플링된 k
개의 랜덤 숫자의 배열을 반환하는 rfn(k)
이어야 하며, rng
가 제공된 경우에는 대신 rfn(rng, k)
가 되어야 합니다.
예제
julia> sprand(Bool, 2, 2, 0.5)
2×2 SparseMatrixCSC{Bool, Int64} with 2 stored entries:
1 1
⋅ ⋅
julia> sprand(Float64, 3, 0.75)
3-element SparseVector{Float64, Int64} with 2 stored entries:
[1] = 0.795547
[2] = 0.49425
SparseArrays.sprandn
— Functionsprandn([rng][,Type],m[,n],p::AbstractFloat)
길이가 m
인 임의의 희소 벡터 또는 크기가 m
x n
인 희소 행렬을 생성하며, 각 항목이 0이 아닐 확률 p
를 지정합니다. 0이 아닌 값은 정규 분포에서 샘플링됩니다. 선택적 rng
인자는 난수 생성기를 지정하며, 난수를 참조하십시오.
출력 요소 유형 Type
을 지정하려면 최소한 Julia 1.1이 필요합니다.
예제
julia> sprandn(2, 2, 0.75)
2×2 SparseMatrixCSC{Float64, Int64} with 3 stored entries:
-1.20577 ⋅
0.311817 -0.234641
SparseArrays.nonzeros
— Functionnonzeros(A)
희소 배열 A
의 구조적 비영(非零) 값을 벡터로 반환합니다. 여기에는 희소 배열에 명시적으로 저장된 영(零) 값도 포함됩니다. 반환된 벡터는 A
의 내부 비영 저장소를 직접 가리키며, 반환된 벡터에 대한 수정은 A
도 변형시킵니다. rowvals
및 nzrange
를 참조하세요.
예제
julia> A = sparse(2I, 3, 3)
3×3 SparseMatrixCSC{Int64, Int64} with 3 stored entries:
2 ⋅ ⋅
⋅ 2 ⋅
⋅ ⋅ 2
julia> nonzeros(A)
3-element Vector{Int64}:
2
2
2
SparseArrays.rowvals
— Functionrowvals(A::AbstractSparseMatrixCSC)
A
의 행 인덱스 벡터를 반환합니다. 반환된 벡터에 대한 수정은 A
도 변형시킵니다. 행 인덱스가 내부적으로 어떻게 저장되는지에 대한 접근을 제공하는 것은 구조적 비영 값에 대해 반복하는 것과 함께 유용할 수 있습니다. nonzeros
및 nzrange
도 참조하십시오.
예제
julia> A = sparse(2I, 3, 3)
3×3 SparseMatrixCSC{Int64, Int64} with 3 stored entries:
2 ⋅ ⋅
⋅ 2 ⋅
⋅ ⋅ 2
julia> rowvals(A)
3-element Vector{Int64}:
1
2
3
SparseArrays.nzrange
— Functionnzrange(A::AbstractSparseMatrixCSC, col::Integer)
희소 행렬 열의 구조적 비영(非零) 값에 대한 인덱스 범위를 반환합니다. nonzeros
및 rowvals
와 함께 사용하면 희소 행렬을 편리하게 반복(iterate)할 수 있습니다:
A = sparse(I,J,V)
rows = rowvals(A)
vals = nonzeros(A)
m, n = size(A)
for j = 1:n
for i in nzrange(A, j)
row = rows[i]
val = vals[i]
# 희소 마법 수행...
end
end
!!! 경고 행렬에 비영(非零) 요소를 추가하거나 제거하면 nzrange
가 무효화될 수 있으므로, 반복(iterating)하는 동안 행렬을 변경해서는 안 됩니다.
nzrange(x::SparseVectorUnion, col)
희소 벡터의 구조적 비영 값에 대한 인덱스 범위를 제공합니다. 열 인덱스 col
은 무시됩니다(기본값은 1
으로 가정됨).
SparseArrays.droptol!
— Functiondroptol!(A::AbstractSparseMatrixCSC, tol)
A
에서 절대값이 tol
보다 작거나 같은 저장된 값을 제거합니다.
droptol!(x::AbstractCompressedVector, tol)
x
에서 절대값이 tol
보다 작거나 같은 저장된 값을 제거합니다.
SparseArrays.dropzeros!
— Functiondropzeros!(x::AbstractCompressedVector)
저장된 숫자 0을 x
에서 제거합니다.
제자리에서 제거하지 않는 버전은 dropzeros
를 참조하세요. 알고리즘 정보는 fkeep!
를 참조하세요.
SparseArrays.dropzeros
— Functiondropzeros(A::AbstractSparseMatrixCSC;)
A
의 복사본을 생성하고 해당 복사본에서 저장된 수치 0을 제거합니다.
제자리 버전 및 알고리즘 정보는 dropzeros!
를 참조하십시오.
예제
julia> A = sparse([1, 2, 3], [1, 2, 3], [1.0, 0.0, 1.0])
3×3 SparseMatrixCSC{Float64, Int64} with 3 stored entries:
1.0 ⋅ ⋅
⋅ 0.0 ⋅
⋅ ⋅ 1.0
julia> dropzeros(A)
3×3 SparseMatrixCSC{Float64, Int64} with 2 stored entries:
1.0 ⋅ ⋅
⋅ ⋅ ⋅
⋅ ⋅ 1.0
dropzeros(x::AbstractCompressedVector)
x
의 복사본을 생성하고 해당 복사본에서 숫자 0을 제거합니다.
제자리 버전 및 알고리즘 정보는 dropzeros!
를 참조하십시오.
예제
julia> A = sparsevec([1, 2, 3], [1.0, 0.0, 1.0])
3-element SparseVector{Float64, Int64} with 3 stored entries:
[1] = 1.0
[2] = 0.0
[3] = 1.0
julia> dropzeros(A)
3-element SparseVector{Float64, Int64} with 2 stored entries:
[1] = 1.0
[3] = 1.0
SparseArrays.permute
— Functionpermute(A::AbstractSparseMatrixCSC{Tv,Ti}, p::AbstractVector{<:Integer},
q::AbstractVector{<:Integer}) where {Tv,Ti}
행렬 A
를 양방향으로 순열하여 PAQ
(A[p,q]
)를 반환합니다. 열 순열 q
의 길이는 A
의 열 수와 일치해야 합니다 (length(q) == size(A, 2)
). 행 순열 p
의 길이는 A
의 행 수와 일치해야 합니다 (length(p) == size(A, 1)
).
전문가 드라이버 및 추가 정보는 permute!
를 참조하십시오.
예제
julia> A = spdiagm(0 => [1, 2, 3, 4], 1 => [5, 6, 7])
4×4 SparseMatrixCSC{Int64, Int64} with 7 stored entries:
1 5 ⋅ ⋅
⋅ 2 6 ⋅
⋅ ⋅ 3 7
⋅ ⋅ ⋅ 4
julia> permute(A, [4, 3, 2, 1], [1, 2, 3, 4])
4×4 SparseMatrixCSC{Int64, Int64} with 7 stored entries:
⋅ ⋅ ⋅ 4
⋅ ⋅ 3 7
⋅ 2 6 ⋅
1 5 ⋅ ⋅
julia> permute(A, [1, 2, 3, 4], [4, 3, 2, 1])
4×4 SparseMatrixCSC{Int64, Int64} with 7 stored entries:
⋅ ⋅ 5 1
⋅ 6 2 ⋅
7 3 ⋅ ⋅
4 ⋅ ⋅ ⋅
Base.permute!
— Methodpermute!(X::AbstractSparseMatrixCSC{Tv,Ti}, A::AbstractSparseMatrixCSC{Tv,Ti},
p::AbstractVector{<:Integer}, q::AbstractVector{<:Integer},
[C::AbstractSparseMatrixCSC{Tv,Ti}]) where {Tv,Ti}
A
를 양방향으로 치환하여 결과 PAQ
(A[p,q]
)를 X
에 저장합니다. 선택적 인수 C
가 있는 경우 중간 결과 (AQ)^T
(transpose(A[:,q])
)를 저장합니다. X
, A
, 그리고, 만약 존재한다면 C
가 서로 별개이어야 합니다; 결과 PAQ
를 A
에 다시 저장하려면 X
없이 다음 방법을 사용하십시오:
permute!(A::AbstractSparseMatrixCSC{Tv,Ti}, p::AbstractVector{<:Integer},
q::AbstractVector{<:Integer}[, C::AbstractSparseMatrixCSC{Tv,Ti},
[workcolptr::Vector{Ti}]]) where {Tv,Ti}
X
의 차원은 A
의 차원과 일치해야 합니다 (size(X, 1) == size(A, 1)
및 size(X, 2) == size(A, 2)
), 그리고 X
는 A
의 모든 할당된 항목을 수용할 수 있는 충분한 저장 공간을 가져야 합니다 (length(rowvals(X)) >= nnz(A)
및 length(nonzeros(X)) >= nnz(A)
). 열 치환 q
의 길이는 A
의 열 수와 일치해야 합니다 (length(q) == size(A, 2)
). 행 치환 p
의 길이는 A
의 행 수와 일치해야 합니다 (length(p) == size(A, 1)
).
C
의 차원은 transpose(A)
의 차원과 일치해야 합니다 (size(C, 1) == size(A, 2)
및 size(C, 2) == size(A, 1)
), 그리고 C
는 A
의 모든 할당된 항목을 수용할 수 있는 충분한 저장 공간을 가져야 합니다 (length(rowvals(C)) >= nnz(A)
및 length(nonzeros(C)) >= nnz(A)
).
추가적인 (알고리즘적) 정보와 인수 검사를 생략한 이러한 메서드의 버전은 (내보내지 않은) 부모 메서드 unchecked_noalias_permute!
및 unchecked_aliasing_permute!
를 참조하십시오.
또한 permute
를 참조하십시오.
SparseArrays.halfperm!
— Functionhalfperm!(X::AbstractSparseMatrixCSC{Tv,Ti}, A::AbstractSparseMatrixCSC{TvA,Ti},
q::AbstractVector{<:Integer}, f::Function = identity) where {Tv,TvA,Ti}
열을 재배열하고 A
를 전치하며, 동시에 f
를 A
의 각 항목에 적용하고 결과 (f(A)Q)^T
(map(f, transpose(A[:,q]))
)를 X
에 저장합니다.
X
의 요소 유형 Tv
는 f(::TvA)
와 일치해야 하며, 여기서 TvA
는 A
의 요소 유형입니다. X
의 차원은 transpose(A)
의 차원과 일치해야 하며 (size(X, 1) == size(A, 2)
및 size(X, 2) == size(A, 1)
), X
는 A
의 모든 할당된 항목을 수용할 수 있는 충분한 저장 공간을 가져야 합니다 (length(rowvals(X)) >= nnz(A)
및 length(nonzeros(X)) >= nnz(A)
). 열 재배열 q
의 길이는 A
의 열 수와 일치해야 합니다 (length(q) == size(A, 2)
).
이 메서드는 SparseMatrixCSC
에서 전치 및 재배열 작업을 수행하는 여러 메서드의 부모입니다. 이 메서드는 인수 검사를 수행하지 않으므로, 직접 사용하기보다는 더 안전한 자식 메서드 ([c]transpose[!]
, permute[!]
)를 사용하는 것이 좋습니다.
이 메서드는 F. Gustavson의 "희소 행렬을 위한 두 가지 빠른 알고리즘: 곱셈 및 재배열된 전치," ACM TOMS 4(3), 250-269 (1978)에서 설명된 HALFPERM
알고리즘을 구현합니다. 이 알고리즘은 O(size(A, 1), size(A, 2), nnz(A))
시간에 실행되며, 전달된 것 외에 추가 공간이 필요하지 않습니다.
SparseArrays.ftranspose!
— Functionftranspose!(X::AbstractSparseMatrixCSC{Tv,Ti}, A::AbstractSparseMatrixCSC{Tv,Ti}, f::Function) where {Tv,Ti}
행렬 A
를 전치하고 비영(非零) 요소에 함수 f
를 적용하여 결과를 X
에 저장합니다. f
에 의해 생성된 영(零) 요소는 제거되지 않습니다. size(X)
는 size(transpose(A))
와 같아야 합니다. 필요할 경우 X
의 rowval 및 nzval 크기를 조정하는 것 외에는 추가 메모리가 할당되지 않습니다.
halfperm!
를 참조하세요.
Sparse Linear Algebra API
LinearAlgebra.cholesky
— Functioncholesky(A, NoPivot(); check = true) -> Cholesky
밀집 대칭 양의 정부호 행렬 A
의 Cholesky 분해를 계산하고 Cholesky
분해를 반환합니다. 행렬 A
는 Symmetric
또는 Hermitian
AbstractMatrix
일 수 있으며, 완벽하게 대칭 또는 에르미트 AbstractMatrix
일 수 있습니다.
삼각형 Cholesky 인자는 분해 F
에서 F.L
및 F.U
를 통해 얻을 수 있으며, 여기서 A ≈ F.U' * F.U ≈ F.L * F.L'
입니다.
Cholesky
객체에 대해 사용할 수 있는 함수는 size
, \
, inv
, det
, logdet
및 isposdef
입니다.
구성 과정에서 반올림 오류로 인해 약간 비에르미트인 행렬 A
가 있는 경우, cholesky
에 전달하기 전에 Hermitian(A)
로 감싸서 완벽한 에르미트로 처리하십시오.
check = true
일 때 분해가 실패하면 오류가 발생합니다. check = false
일 때 분해의 유효성을 확인하는 책임은 사용자에게 있습니다( issuccess
사용).
예제
julia> A = [4. 12. -16.; 12. 37. -43.; -16. -43. 98.]
3×3 Matrix{Float64}:
4.0 12.0 -16.0
12.0 37.0 -43.0
-16.0 -43.0 98.0
julia> C = cholesky(A)
Cholesky{Float64, Matrix{Float64}}
U 인자:
3×3 UpperTriangular{Float64, Matrix{Float64}}:
2.0 6.0 -8.0
⋅ 1.0 5.0
⋅ ⋅ 3.0
julia> C.U
3×3 UpperTriangular{Float64, Matrix{Float64}}:
2.0 6.0 -8.0
⋅ 1.0 5.0
⋅ ⋅ 3.0
julia> C.L
3×3 LowerTriangular{Float64, Matrix{Float64}}:
2.0 ⋅ ⋅
6.0 1.0 ⋅
-8.0 5.0 3.0
julia> C.L * C.U == A
true
cholesky(A, RowMaximum(); tol = 0.0, check = true) -> CholeskyPivoted
밀집 대칭 양의 준정칙 행렬 A
의 피벗된 Cholesky 분해를 계산하고 CholeskyPivoted
분해를 반환합니다. 행렬 A
는 Symmetric
또는 Hermitian
AbstractMatrix
일 수 있으며, 완벽하게 대칭 또는 허미티안 AbstractMatrix
일 수 있습니다.
삼각형 Cholesky 인자는 분해 F
에서 F.L
및 F.U
를 통해 얻을 수 있으며, 순열은 F.p
를 통해 얻을 수 있습니다. 여기서 A[F.p, F.p] ≈ Ur' * Ur ≈ Lr * Lr'
이며, Ur = F.U[1:F.rank, :]
및 Lr = F.L[:, 1:F.rank]
입니다. 또는 대안으로 A ≈ Up' * Up ≈ Lp * Lp'
이며, Up = F.U[1:F.rank, invperm(F.p)]
및 Lp = F.L[invperm(F.p), 1:F.rank]
입니다.
다음 함수는 CholeskyPivoted
객체에 대해 사용할 수 있습니다: size
, \
, inv
, det
, 및 rank
.
인자 tol
은 랭크를 결정하기 위한 허용 오차를 결정합니다. 음수 값의 경우, 허용 오차는 기계 정밀도입니다.
구성 과정에서 반올림 오류로 인해 약간 비허미티안인 행렬 A
가 있는 경우, cholesky
에 전달하기 전에 Hermitian(A)
로 감싸서 완벽한 허미티안으로 처리하십시오.
check = true
일 때, 분해가 실패하면 오류가 발생합니다. check = false
일 때, 분해의 유효성을 확인하는 책임( issuccess
사용)은 사용자에게 있습니다.
예제
julia> X = [1.0, 2.0, 3.0, 4.0];
julia> A = X * X';
julia> C = cholesky(A, RowMaximum(), check = false)
CholeskyPivoted{Float64, Matrix{Float64}, Vector{Int64}}
U 인자는 랭크 1:
4×4 UpperTriangular{Float64, Matrix{Float64}}:
4.0 2.0 3.0 1.0
⋅ 0.0 6.0 2.0
⋅ ⋅ 9.0 3.0
⋅ ⋅ ⋅ 1.0
순열:
4-element Vector{Int64}:
4
2
3
1
julia> C.U[1:C.rank, :]' * C.U[1:C.rank, :] ≈ A[C.p, C.p]
true
julia> l, u = C; # 반복을 통한 구조 분해
julia> l == C.L && u == C.U
true
cholesky(A::SparseMatrixCSC; shift = 0.0, check = true, perm = nothing) -> CHOLMOD.Factor
희소 양의 정부호 행렬 A
의 Cholesky 분해를 계산합니다. A
는 SparseMatrixCSC
또는 SparseMatrixCSC
의 Symmetric
/Hermitian
뷰여야 합니다. A
가 타입 태그가 없더라도 여전히 대칭 또는 에르미트여야 합니다. perm
이 주어지지 않으면, 채우기 감소 순열이 사용됩니다. F = cholesky(A)
는 F\b
로 방정식 시스템을 해결하는 데 가장 자주 사용되지만, diag
, det
, 및 logdet
와 같은 메서드도 F
에 대해 정의되어 있습니다. F
에서 개별 요소를 추출할 수도 있으며, F.L
을 사용할 수 있습니다. 그러나 피벗이 기본적으로 활성화되어 있기 때문에, 분해는 내부적으로 A == P'*L*L'*P
로 표현되며, 여기서 P
는 순열 행렬입니다. P
를 고려하지 않고 단순히 L
을 사용하면 잘못된 결과가 나옵니다. 순열의 효과를 포함하려면, 일반적으로 PtL = F.PtL
(즉, P'*L
) 및 LtP = F.UP
(즉, L'*P
)와 같은 "결합된" 요소를 추출하는 것이 바람직합니다.
check = true
일 때, 분해가 실패하면 오류가 발생합니다. check = false
일 때, 분해의 유효성을 확인하는 책임은 사용자에게 있습니다(즉, issuccess
사용).
선택적 shift
키워드 인수를 설정하면 A
대신 A+shift*I
의 분해를 계산합니다. perm
인수가 제공되면, 이는 사용해야 할 순서를 제공하는 1:size(A,1)
의 순열이어야 합니다(즉, CHOLMOD의 기본 AMD 순서 대신).
예제
다음 예제에서 사용된 채우기 감소 순열은 [3, 2, 1]
입니다. perm
이 1:3
으로 설정되어 순열이 없도록 강제하면, 요소의 수는 6입니다.
julia> A = [2 1 1; 1 2 0; 1 0 2]
3×3 Matrix{Int64}:
2 1 1
1 2 0
1 0 2
julia> C = cholesky(sparse(A))
SparseArrays.CHOLMOD.Factor{Float64, Int64}
type: LLt
method: simplicial
maxnnz: 5
nnz: 5
success: true
julia> C.p
3-element Vector{Int64}:
3
2
1
julia> L = sparse(C.L);
julia> Matrix(L)
3×3 Matrix{Float64}:
1.41421 0.0 0.0
0.0 1.41421 0.0
0.707107 0.707107 1.0
julia> L * L' ≈ A[C.p, C.p]
true
julia> P = sparse(1:3, C.p, ones(3))
3×3 SparseMatrixCSC{Float64, Int64} with 3 stored entries:
⋅ ⋅ 1.0
⋅ 1.0 ⋅
1.0 ⋅ ⋅
julia> P' * L * L' * P ≈ A
true
julia> C = cholesky(sparse(A), perm=1:3)
SparseArrays.CHOLMOD.Factor{Float64, Int64}
type: LLt
method: simplicial
maxnnz: 6
nnz: 6
success: true
julia> L = sparse(C.L);
julia> Matrix(L)
3×3 Matrix{Float64}:
1.41421 0.0 0.0
0.707107 1.22474 0.0
0.707107 -0.408248 1.1547
julia> L * L' ≈ A
true
이 방법은 SuiteSparse에서 CHOLMOD[ACM887][DavisHager2009] 라이브러리를 사용합니다. CHOLMOD는 단일 또는 이중 정밀도의 실수 또는 복소수 유형만 지원합니다. 이러한 요소 유형이 아닌 입력 행렬은 적절하게 이러한 유형으로 변환됩니다.
CHOLMOD의 많은 다른 함수가 래핑되었지만 Base.SparseArrays.CHOLMOD
모듈에서 내보내지지는 않았습니다.
LinearAlgebra.cholesky!
— Functioncholesky!(A::AbstractMatrix, NoPivot(); check = true) -> Cholesky
cholesky
와 동일하지만, 입력 A
를 복사하는 대신 덮어써서 공간을 절약합니다. 인수 분해가 A
의 요소 유형으로 표현할 수 없는 숫자를 생성하는 경우, 예를 들어 정수 유형의 경우, InexactError
예외가 발생합니다.
예제
julia> A = [1 2; 2 50]
2×2 Matrix{Int64}:
1 2
2 50
julia> cholesky!(A)
ERROR: InexactError: Int64(6.782329983125268)
Stacktrace:
[...]
cholesky!(A::AbstractMatrix, RowMaximum(); tol = 0.0, check = true) -> CholeskyPivoted
cholesky
와 동일하지만, 복사본을 생성하는 대신 입력 A
를 덮어써서 공간을 절약합니다. 인수 분해가 A
의 요소 유형으로 표현할 수 없는 숫자를 생성하는 경우, 예를 들어 정수 유형의 경우 InexactError
예외가 발생합니다.
cholesky!(F::CHOLMOD.Factor, A::SparseMatrixCSC; shift = 0.0, check = true) -> CHOLMOD.Factor
A
의 Cholesky ($LL'$) 분해를 계산하고, 기호 분해 F
를 재사용합니다. A
는 SparseMatrixCSC
또는 SparseMatrixCSC
의 Symmetric
/Hermitian
뷰여야 합니다. A
가 타입 태그가 없더라도 여전히 대칭 또는 에르미트여야 합니다.
자세한 내용은 cholesky
를 참조하세요.
이 방법은 SuiteSparse의 CHOLMOD 라이브러리를 사용하며, 이는 단일 또는 이중 정밀도의 실수 또는 복소수 유형만 지원합니다. 해당 요소 유형이 아닌 입력 행렬은 적절하게 이러한 유형으로 변환됩니다.
LinearAlgebra.lowrankupdate
— Functionlowrankupdate(C::Cholesky, v::AbstractVector) -> CC::Cholesky
벡터 v
로 Cholesky 분해 C
를 업데이트합니다. 만약 A = C.U'C.U
라면 CC = cholesky(C.U'C.U + v*v')
이지만, CC
의 계산은 O(n^2)
연산만 사용합니다.
lowrankupdate(F::CHOLMOD.Factor, C::AbstractArray) -> FF::CHOLMOD.Factor
주어진 A
의 LDLt
또는 LLt
분해 F
에 대해 A + C*C'
의 LDLt
분해를 가져옵니다.
반환된 인자는 항상 LDLt
분해입니다.
또한 lowrankupdate!
, lowrankdowndate
, lowrankdowndate!
를 참조하십시오.
LinearAlgebra.lowrankupdate!
— Functionlowrankupdate!(C::Cholesky, v::AbstractVector) -> CC::Cholesky
벡터 v
로 Cholesky 분해 C
를 업데이트합니다. 만약 A = C.U'C.U
라면 CC = cholesky(C.U'C.U + v*v')
이지만, CC
의 계산은 O(n^2)
연산만 사용합니다. 입력 분해 C
는 제자리에서 업데이트되어 종료 시 C == CC
가 됩니다. 벡터 v
는 계산 중에 파괴됩니다.
lowrankupdate!(F::CHOLMOD.Factor, C::AbstractArray)
A
의 LDLt
또는 LLt
분해 F
를 A + C*C'
의 분해로 업데이트합니다.
LLt
분해는 LDLt
로 변환됩니다.
자세한 내용은 lowrankupdate
, lowrankdowndate
, lowrankdowndate!
를 참조하세요.
LinearAlgebra.lowrankdowndate
— Functionlowrankdowndate(C::Cholesky, v::AbstractVector) -> CC::Cholesky
Cholesky 분해 C
를 벡터 v
로 다운데이트합니다. 만약 A = C.U'C.U
라면 CC = cholesky(C.U'C.U - v*v')
이지만, CC
의 계산은 O(n^2)
연산만 사용합니다.
lowrankdowndate(F::CHOLMOD.Factor, C::AbstractArray) -> FF::CHOLMOD.Factor
주어진 A
의 LDLt
또는 LLt
분해 F
에 대해 A + C*C'
의 LDLt
분해를 가져옵니다.
반환된 분해는 항상 LDLt
분해입니다.
또한 lowrankdowndate!
, lowrankupdate
, lowrankupdate!
를 참조하십시오.
LinearAlgebra.lowrankdowndate!
— Functionlowrankdowndate!(C::Cholesky, v::AbstractVector) -> CC::Cholesky
Cholesky 분해 C
를 벡터 v
로 다운데이트합니다. 만약 A = C.U'C.U
라면 CC = cholesky(C.U'C.U - v*v')
이지만, CC
의 계산은 O(n^2)
연산만 사용합니다. 입력 분해 C
는 제자리에서 업데이트되어 종료 시 C == CC
가 됩니다. 벡터 v
는 계산 중에 파괴됩니다.
lowrankdowndate!(F::CHOLMOD.Factor, C::AbstractArray)
A
의 LDLt
또는 LLt
분해 F
를 A - C*C'
의 분해로 업데이트합니다.
LLt
분해는 LDLt
로 변환됩니다.
자세한 내용은 lowrankdowndate
, lowrankupdate
, lowrankupdate!
를 참조하세요.
SparseArrays.CHOLMOD.lowrankupdowndate!
— Functionlowrankupdowndate!(F::CHOLMOD.Factor, C::Sparse, update::Cint)
A
의 LDLt
또는 LLt
분해 F
를 A ± C*C'
의 분해로 업데이트합니다.
희소성을 유지하는 분해가 사용되는 경우, 즉 L*L' == P*A*P'
이면 새로운 분해는 L*L' == P*A*P' + C'*C
가 됩니다.
update
: Cint(1)
은 A + CC'
에 대한 것이고, Cint(0)
은 A - CC'
에 대한 것입니다.
LinearAlgebra.ldlt
— Functionldlt(S::SymTridiagonal) -> LDLt
실수 대칭 삼중 대각 행렬 S
의 LDLt
(즉, $LDL^T$) 분해를 계산합니다. 여기서 S = L*Diagonal(d)*L'
이며, L
은 단위 하삼각 행렬이고 d
는 벡터입니다. LDLt
분해 F = ldlt(S)
의 주요 용도는 선형 방정식 시스템 Sx = b
를 F\b
로 푸는 것입니다.
유사하지만 피벗된 임의의 대칭 또는 에르미트 행렬의 분해에 대한 bunchkaufman
도 참조하십시오.
예제
julia> S = SymTridiagonal([3., 4., 5.], [1., 2.])
3×3 SymTridiagonal{Float64, Vector{Float64}}:
3.0 1.0 ⋅
1.0 4.0 2.0
⋅ 2.0 5.0
julia> ldltS = ldlt(S);
julia> b = [6., 7., 8.];
julia> ldltS \ b
3-element Vector{Float64}:
1.7906976744186047
0.627906976744186
1.3488372093023255
julia> S \ b
3-element Vector{Float64}:
1.7906976744186047
0.627906976744186
1.3488372093023255
ldlt(A::SparseMatrixCSC; shift = 0.0, check = true, perm=nothing) -> CHOLMOD.Factor
희소 행렬 A
의 $LDL'$ 분해를 계산합니다. A
는 SparseMatrixCSC
또는 SparseMatrixCSC
의 Symmetric
/Hermitian
뷰여야 합니다. A
가 타입 태그가 없더라도 여전히 대칭 또는 에르미트여야 합니다. 채우기 감소 순열이 사용됩니다. F = ldlt(A)
는 시스템의 방정식 A*x = b
를 F\b
로 해결하는 데 가장 자주 사용됩니다. 반환된 분해 객체 F
는 또한 diag
, det
, logdet
, 및 inv
메서드를 지원합니다. F
에서 개별 요소를 추출하려면 F.L
을 사용하면 됩니다. 그러나 피벗이 기본적으로 활성화되어 있기 때문에, 분해는 내부적으로 A == P'*L*D*L'*P
로 표현되며, 여기서 P
는 순열 행렬입니다; P
를 고려하지 않고 단순히 L
을 사용하면 잘못된 결과가 나옵니다. 순열의 효과를 포함하려면, 일반적으로 PtL = F.PtL
(즉, P'*L
에 해당) 및 LtP = F.UP
(즉, L'*P
에 해당)와 같은 "결합된" 요소를 추출하는 것이 바람직합니다. 지원되는 요소의 전체 목록은 :L, :PtL, :D, :UP, :U, :LD, :DU, :PtLD, :DUP
입니다.
check = true
일 때, 분해가 실패하면 오류가 발생합니다. check = false
일 때, 분해의 유효성을 확인하는 책임은 사용자에게 있습니다 (via issuccess
).
선택적 shift
키워드 인수를 설정하면 A
대신 A+shift*I
의 분해를 계산합니다. perm
인수가 제공되면, 이는 사용해야 할 순서를 제공하는 1:size(A,1)
의 순열이어야 합니다 (CHOLMOD의 기본 AMD 순서 대신).
이 방법은 SuiteSparse에서 제공하는 CHOLMOD[ACM887][DavisHager2009] 라이브러리를 사용합니다. CHOLMOD는 단일 또는 이중 정밀도의 실수 또는 복소수 타입만 지원합니다. 해당 요소 타입이 아닌 입력 행렬은 적절하게 이러한 타입으로 변환됩니다.
CHOLMOD의 많은 다른 함수가 래핑되었지만 Base.SparseArrays.CHOLMOD
모듈에서 내보내지지는 않았습니다.
LinearAlgebra.lu
— Functionlu(A::AbstractSparseMatrixCSC; check = true, q = nothing, control = get_umfpack_control()) -> F::UmfpackLU
희소 행렬 A
의 LU 분해를 계산합니다.
실수 또는 복소수 요소 유형을 가진 희소 A
의 경우, F
의 반환 유형은 UmfpackLU{Tv, Ti}
이며, 여기서 Tv
는 Float64
또는 ComplexF64
이고, Ti
는 정수 유형 (Int32
또는 Int64
)입니다.
check = true
일 때, 분해가 실패하면 오류가 발생합니다. check = false
일 때, 분해의 유효성을 확인하는 책임은 사용자에게 있습니다 (via issuccess
).
치환 q
는 치환 벡터이거나 nothing
일 수 있습니다. 치환 벡터가 제공되지 않거나 q
가 nothing
인 경우, UMFPACK의 기본값이 사용됩니다. 치환이 0 기반이 아닐 경우, 0 기반 복사가 이루어집니다.
control
벡터는 UMFPACK에 대한 Julia SparseArrays 패키지의 기본 구성으로 기본값이 설정됩니다 (참고: 이는 반복 정제를 비활성화하기 위해 UMFPACK 기본값에서 수정됨), 그러나 UMFPACK_CONTROL
길이의 벡터를 전달하여 변경할 수 있습니다. 가능한 구성에 대해서는 UMFPACK 매뉴얼을 참조하십시오. 예를 들어 반복 정제를 다시 활성화하려면:
umfpack_control = SparseArrays.UMFPACK.get_umfpack_control(Float64, Int64) # Float64 희소 행렬에 대한 Julia 기본 구성 읽기
SparseArrays.UMFPACK.show_umf_ctrl(umfpack_control) # 선택 사항 - 값 표시
umfpack_control[SparseArrays.UMFPACK.JL_UMFPACK_IRSTEP] = 2.0 # 반복 정제 다시 활성화 (2는 UMFPACK 기본 최대 반복 정제 단계)
Alu = lu(A; control = umfpack_control)
x = Alu \ b # Ax = b를 해결하며, UMFPACK 반복 정제를 포함
분해 F
의 개별 구성 요소는 인덱싱하여 접근할 수 있습니다:
구성 요소 | 설명 |
---|---|
L | LU 의 L (하삼각) 부분 |
U | LU 의 U (상삼각) 부분 |
p | 오른쪽 치환 Vector |
q | 왼쪽 치환 Vector |
Rs | 스케일링 계수의 Vector |
: | (L,U,p,q,Rs) 구성 요소 |
F
와 A
간의 관계는 다음과 같습니다.
F.L*F.U == (F.Rs .* A)[F.p, F.q]
F
는 다음 함수를 추가로 지원합니다:
또한 lu!
를 참조하십시오.
lu(A::AbstractSparseMatrixCSC)
는 SuiteSparse의 일부인 UMFPACK[ACM832] 라이브러리를 사용합니다. 이 라이브러리는 Float64
또는 ComplexF64
요소를 가진 희소 행렬만 지원하므로, lu
는 A
를 적절하게 SparseMatrixCSC{Float64}
또는 SparseMatrixCSC{ComplexF64}
유형의 복사본으로 변환합니다.
lu(A, pivot = RowMaximum(); check = true, allowsingular = false) -> F::LU
행렬 A
의 LU 분해를 계산합니다.
check = true
일 때, 분해가 실패하면 오류가 발생합니다. check = false
일 때, 분해의 유효성을 확인하는 책임은 사용자에게 있습니다 (via issuccess
).
기본적으로, check = true
일 때, 분해가 유효한 요소를 생성하더라도 상삼각 행렬 U
가 랭크 결핍인 경우에도 오류가 발생합니다. 이는 allowsingular = true
를 전달하여 변경할 수 있습니다.
대부분의 경우, A
가 AbstractMatrix{T}
의 하위 유형 S
이고 요소 유형 T
가 +
, -
, *
및 /
를 지원하는 경우, 반환 유형은 LU{T,S{T}}
입니다.
일반적으로 LU 분해는 행렬의 행을 순열하는 것을 포함하며 (아래에 설명된 F.p
출력에 해당), 이를 "피벗팅"이라고 합니다 (이는 "피벗"이 포함된 행, 즉 F.U
의 대각선 항목을 선택하는 것과 관련이 있습니다). 다음의 피벗팅 전략 중 하나를 선택할 수 있습니다:
RowMaximum()
(기본값): 표준 피벗팅 전략; 피벗은 남아 있는, 분해될 행 중 최대 절대값을 가지는 요소에 해당합니다. 이 피벗팅 전략은 요소 유형이abs
및<
도 지원해야 합니다. (이는 일반적으로 부동 소수점 행렬에 대해 유일하게 수치적으로 안정적인 옵션입니다.)RowNonZero()
: 피벗은 남아 있는, 분해될 행 중 첫 번째 비영(非零) 요소에 해당합니다. (이는 수작업 계산에서 일반적인 선택에 해당하며,iszero
를 지원하지만abs
또는<
를 지원하지 않는 보다 일반적인 대수적 수 유형에도 유용합니다.)NoPivot()
: 피벗팅이 꺼집니다 (피벗 위치에서 0 항목이 발견되면 실패하며,allowsingular = true
일 때도 마찬가지입니다).
분해 F
의 개별 구성 요소는 getproperty
를 통해 접근할 수 있습니다:
구성 요소 | 설명 |
---|---|
F.L | LU 의 L (하삼각) 부분 |
F.U | LU 의 U (상삼각) 부분 |
F.p | (오른쪽) 순열 Vector |
F.P | (오른쪽) 순열 Matrix |
분해를 반복하면 구성 요소 F.L
, F.U
, 및 F.p
가 생성됩니다.
F
와 A
의 관계는 다음과 같습니다.
F.L*F.U == A[F.p, :]
F
는 다음 함수를 추가로 지원합니다:
지원 함수 | LU | LU{T,Tridiagonal{T}} |
---|---|---|
/ | ✓ | |
\ | ✓ | ✓ |
inv | ✓ | ✓ |
det | ✓ | ✓ |
logdet | ✓ | ✓ |
logabsdet | ✓ | ✓ |
size | ✓ | ✓ |
allowsingular
키워드 인자는 Julia 1.11에서 추가되었습니다.
예제
julia> A = [4 3; 6 3]
2×2 Matrix{Int64}:
4 3
6 3
julia> F = lu(A)
LU{Float64, Matrix{Float64}, Vector{Int64}}
L factor:
2×2 Matrix{Float64}:
1.0 0.0
0.666667 1.0
U factor:
2×2 Matrix{Float64}:
6.0 3.0
0.0 1.0
julia> F.L * F.U == A[F.p, :]
true
julia> l, u, p = lu(A); # 반복을 통한 구조 분해
julia> l == F.L && u == F.U && p == F.p
true
julia> lu([1 2; 1 2], allowsingular = true)
LU{Float64, Matrix{Float64}, Vector{Int64}}
L factor:
2×2 Matrix{Float64}:
1.0 0.0
1.0 1.0
U factor (rank-deficient):
2×2 Matrix{Float64}:
1.0 2.0
0.0 0.0
LinearAlgebra.qr
— Functionqr(A::SparseMatrixCSC; tol=_default_tol(A), ordering=ORDERING_DEFAULT) -> QRSparse
희소 행렬 A
의 QR
분해를 계산합니다. F.R = F.Q'*A[F.prow,F.pcol]
가 되도록 행과 열의 채우기 감소 순열이 사용됩니다. 이 유형의 주요 응용 프로그램은 \
를 사용하여 최소 제곱 또는 미결정 문제를 해결하는 것입니다. 이 함수는 C 라이브러리 SPQR[ACM933]를 호출합니다.
qr(A::SparseMatrixCSC)
는 SuiteSparse의 일부인 SPQR 라이브러리를 사용합니다. 이 라이브러리는 Float64
또는 ComplexF64
요소를 가진 희소 행렬만 지원하므로, Julia v1.4부터 qr
는 A
를 SparseMatrixCSC{Float64}
또는 SparseMatrixCSC{ComplexF64}
유형의 복사본으로 변환합니다.
예제
julia> A = sparse([1,2,3,4], [1,1,2,2], [1.0,1.0,1.0,1.0])
4×2 SparseMatrixCSC{Float64, Int64} with 4 stored entries:
1.0 ⋅
1.0 ⋅
⋅ 1.0
⋅ 1.0
julia> qr(A)
SparseArrays.SPQR.QRSparse{Float64, Int64}
Q factor:
4×4 SparseArrays.SPQR.QRSparseQ{Float64, Int64}
R factor:
2×2 SparseMatrixCSC{Float64, Int64} with 2 stored entries:
-1.41421 ⋅
⋅ -1.41421
Row permutation:
4-element Vector{Int64}:
1
3
4
2
Column permutation:
2-element Vector{Int64}:
1
2
qr(A, pivot = NoPivot(); blocksize) -> F
행렬 A
의 QR 분해를 계산합니다: 직교(또는 A
가 복소수 값일 경우 유니타리) 행렬 Q
와 상삼각 행렬 R
이 있어 다음과 같은 관계를 만족합니다.
\[A = Q R\]
반환된 객체 F
는 패킹된 형식으로 분해를 저장합니다:
pivot == ColumnNorm()
인 경우F
는QRPivoted
객체입니다.- 그렇지 않고
A
의 요소 유형이 BLAS 유형(Float32
,Float64
,ComplexF32
또는ComplexF64
)인 경우,F
는QRCompactWY
객체입니다. - 그렇지 않으면
F
는QR
객체입니다.
분해의 개별 구성 요소 F
는 속성 접근자를 통해 검색할 수 있습니다:
F.R
를 통해 상삼각 요소에 대한 각 참조는 새로운 배열을 할당합니다. 따라서 R = F.R
와 같이 해당 배열을 캐시하고 R
로 계속 작업하는 것이 좋습니다.
분해를 반복하면 구성 요소 Q
, R
, 그리고 존재하는 경우 p
를 생성합니다.
QR
객체에 대해 사용할 수 있는 함수는 다음과 같습니다: inv
, size
, 및 \
. A
가 직사각형인 경우, \
는 최소 제곱 해를 반환하며, 해가 유일하지 않은 경우 가장 작은 노름을 가진 해가 반환됩니다. A
가 풀 랭크가 아닐 경우, 최소 노름 해를 얻기 위해 (열) 피벗을 사용한 분해가 필요합니다.
전체/정사각형 또는 비전체/비정사각형 Q
에 대한 곱셈이 허용됩니다. 즉, F.Q*F.R
와 F.Q*A
모두 지원됩니다. Q
행렬은 Matrix
로 일반 행렬로 변환할 수 있습니다. 이 작업은 "얇은" Q 요소를 반환합니다. 즉, A
가 m
×n
이고 m>=n
인 경우, Matrix(F.Q)
는 직교 정규 열을 가진 m
×n
행렬을 생성합니다. "전체" Q 요소를 검색하려면 F.Q*I
또는 collect(F.Q)
를 사용하십시오. m<=n
인 경우, Matrix(F.Q)
는 m
×m
직교 행렬을 생성합니다.
QR 분해를 위한 블록 크기는 pivot == NoPivot()
이고 A isa StridedMatrix{<:BlasFloat}
일 때 키워드 인수 blocksize :: Integer
로 지정할 수 있습니다. blocksize > minimum(size(A))
인 경우 무시됩니다. QRCompactWY
를 참조하십시오.
blocksize
키워드 인수는 Julia 1.4 이상이 필요합니다.
예제
julia> A = [3.0 -6.0; 4.0 -8.0; 0.0 1.0]
3×2 Matrix{Float64}:
3.0 -6.0
4.0 -8.0
0.0 1.0
julia> F = qr(A)
LinearAlgebra.QRCompactWY{Float64, Matrix{Float64}, Matrix{Float64}}
Q 요소: 3×3 LinearAlgebra.QRCompactWYQ{Float64, Matrix{Float64}, Matrix{Float64}}
R 요소:
2×2 Matrix{Float64}:
-5.0 10.0
0.0 -1.0
julia> F.Q * F.R == A
true
qr
는 LAPACK가 Householder 기본 반사체의 곱의 메모리 저장 요구 사항을 최소화하는 여러 표현을 사용하기 때문에 여러 유형을 반환합니다. 따라서 Q
및 R
행렬을 두 개의 개별 밀집 행렬 대신 압축하여 저장할 수 있습니다.
Noteworthy External Sparse Packages
여러 다른 Julia 패키지가 언급해야 할 희소 행렬 구현을 제공합니다:
SuiteSparseGraphBLAS.jl는 빠르고 다중 스레드 지원을 하는 SuiteSparse:GraphBLAS C 라이브러리에 대한 래퍼입니다. CPU에서 이는 일반적으로 가장 빠른 옵션으로, 종종 MKLSparse보다 상당히 뛰어난 성능을 발휘합니다.
SparseMatricesCSR.jl는 압축 희소 행렬(CSR) 형식의 Julia 네이티브 구현을 제공합니다.
MKLSparse.jl는 Intel의 MKL 라이브러리를 사용하여 SparseArrays의 희소-밀집 행렬 연산을 가속화합니다.
SparseArrayKit.jl 다차원 희소 배열에 사용할 수 있습니다.
LuxurySparse.jl는 정적 희소 배열 형식과 좌표 형식을 제공합니다.
ExtendableSparse.jl는 새로운 저장된 인덱스에 대한 지연 접근 방식을 사용하여 희소 행렬에 빠른 삽입을 가능하게 합니다.
Finch.jl는 미니 텐서 언어와 컴파일러를 통해 광범위한 다차원 희소 배열 형식 및 연산을 지원하며, 모두 네이티브 줄리아로 작성되었습니다. COO, CSF, CSR, CSC 등과 같은 형식과 브로드캐스트, 리듀스 등과 같은 연산 및 사용자 정의 연산을 지원합니다.
외부 패키지에서 제공하는 희소 직접 솔버:
외부 패키지들은 고유 시스템과 특이값 분해의 반복적 해결을 위한 솔버를 제공합니다:
그래프 작업을 위한 외부 패키지:
- ACM887Chen, Y., Davis, T. A., Hager, W. W., & Rajamanickam, S. (2008). Algorithm 887: CHOLMOD, Supernodal Sparse Cholesky Factorization and Update/Downdate. ACM Trans. Math. Softw., 35(3). doi:10.1145/1391989.1391995
- DavisHager2009Davis, Timothy A., & Hager, W. W. (2009). Dynamic Supernodes in Sparse Cholesky Update/Downdate and Triangular Solves. ACM Trans. Math. Softw., 35(4). doi:10.1145/1462173.1462176
- ACM832Davis, Timothy A. (2004b). Algorithm 832: UMFPACK V4.3–-an Unsymmetric-Pattern Multifrontal Method. ACM Trans. Math. Softw., 30(2), 196–199. doi:10.1145/992200.992206
- ACM933Foster, L. V., & Davis, T. A. (2013). Algorithm 933: Reliable Calculation of Numerical Rank, Null Space Bases, Pseudoinverse Solutions, and Basic Solutions Using SuitesparseQR. ACM Trans. Math. Softw., 40(1). doi:10.1145/2513109.2513116