Skip to content

Commit 6c66f3e

Browse files
authored
Adjust to AdjointFactorization (#313)
1 parent 992ab32 commit 6c66f3e

File tree

5 files changed

+27
-13
lines changed

5 files changed

+27
-13
lines changed

Project.toml

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,6 @@
11
name = "BandedMatrices"
22
uuid = "aae01518-5342-5314-be14-df237901396f"
3-
version = "0.17.11"
3+
version = "0.17.12"
44

55
[deps]
66
ArrayLayouts = "4c555306-a7a7-4459-81d9-ec55ddd5c99a"

src/BandedMatrices.jl

Lines changed: 6 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -40,6 +40,12 @@ import FillArrays: AbstractFill, getindex_value, _broadcasted_zeros, unique_valu
4040

4141
const libblas = Base.libblas_name
4242
const liblapack = Base.liblapack_name
43+
const AdjointFact = isdefined(LinearAlgebra, :AdjointFactorization) ?
44+
LinearAlgebra.AdjointFactorization :
45+
Adjoint
46+
const TransposeFact = isdefined(LinearAlgebra, :TransposeFactorization) ?
47+
LinearAlgebra.TransposeFactorization :
48+
Transpose
4349

4450
export BandedMatrix,
4551
bandrange,

src/banded/BandedLU.jl

Lines changed: 4 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -32,8 +32,10 @@ Base.iterate(S::BandedLU, ::Val{:U}) = (S.U, Val(:p))
3232
Base.iterate(S::BandedLU, ::Val{:p}) = (S.p, Val(:done))
3333
Base.iterate(S::BandedLU, ::Val{:done}) = nothing
3434

35-
adjoint(F::BandedLU) = Adjoint(F)
36-
transpose(F::BandedLU) = Transpose(F)
35+
if !(isdefined(LinearAlgebra, :AdjointFactorization)) # VERSION < v"1.10-"
36+
adjoint(F::BandedLU) = Adjoint(F)
37+
end
38+
transpose(F::BandedLU) = TransposeFact(F)
3739

3840
lu(S::BandedLU) = S
3941

src/banded/linalg.jl

Lines changed: 7 additions & 7 deletions
Original file line numberDiff line numberDiff line change
@@ -38,38 +38,38 @@ function ldiv!(A::BandedLU{T}, B::AbstractVecOrMat{Complex{T}}) where T<:Real
3838
B .= a .+ im.*b
3939
end
4040

41-
function ldiv!(transA::Transpose{T,<:BandedLU{T,<:BandedMatrix}}, B::StridedVecOrMat{T}) where {T<:BlasFloat}
41+
function ldiv!(transA::TransposeFact{T,<:BandedLU{T,<:BandedMatrix}}, B::StridedVecOrMat{T}) where {T<:BlasFloat}
4242
A = transA.parent
4343
m = size(A.factors,1)
4444
l,u = bandwidths(A.factors)
4545
data = bandeddata(A.factors)
4646
LAPACK.gbtrs!('T', l, u-l, m, data, A.ipiv, B)
4747
end
4848

49-
function ldiv!(transA::Transpose{<:Any,<:BandedLU{<:Any,<:BandedMatrix}}, B::AbstractVecOrMat)
49+
function ldiv!(transA::TransposeFact{<:Any,<:BandedLU{<:Any,<:BandedMatrix}}, B::AbstractVecOrMat)
5050
A = transA.parent
5151
ldiv!(transpose(UnitLowerTriangular(A.factors)), ldiv!(transpose(UpperTriangular(A.factors)), B))
5252
_apply_inverse_ipiv_rows!(A, B)
5353
end
5454

55-
ldiv!(adjF::Adjoint{T,<:BandedLU{T,<:BandedMatrix}}, B::AbstractVecOrMat{T}) where {T<:Real} =
55+
ldiv!(adjF::AdjointFact{T,<:BandedLU{T,<:BandedMatrix}}, B::AbstractVecOrMat{T}) where {T<:Real} =
5656
(F = adjF.parent; ldiv!(transpose(F), B))
57-
function ldiv!(adjA::Adjoint{T,<:BandedLU{T,<:BandedMatrix}}, B::StridedVecOrMat{T}) where {T<:BlasComplex}
57+
function ldiv!(adjA::AdjointFact{T,<:BandedLU{T,<:BandedMatrix}}, B::StridedVecOrMat{T}) where {T<:BlasComplex}
5858
A = adjA.parent
5959
m = size(A.factors,1)
6060
l,u = bandwidths(A.factors)
6161
data = bandeddata(A.factors)
6262
LAPACK.gbtrs!('C', l, u-l, m, data, A.ipiv, B)
6363
end
6464

65-
function ldiv!(adjA::Adjoint{<:Any,<:BandedLU{<:Any,<:BandedMatrix}}, B::AbstractVecOrMat)
65+
function ldiv!(adjA::AdjointFact{<:Any,<:BandedLU{<:Any,<:BandedMatrix}}, B::AbstractVecOrMat)
6666
error("Implement")
6767
A = adjA.parent
6868
ldiv!(adjoint(UnitLowerTriangular(A.factors)), ldiv!(adjoint(UpperTriangular(A.factors)), B))
6969
_apply_inverse_ipiv!(A, B)
7070
end
7171

72-
\(A::Adjoint{<:Any,<:BandedLU}, B::Adjoint{<:Any,<:AbstractVecOrMat}) = A \ copy(B)
73-
\(A::Transpose{<:Any,<:BandedLU}, B::Transpose{<:Any,<:AbstractVecOrMat}) = A \ copy(B)
72+
\(A::AdjointFact{<:Any,<:BandedLU}, B::Adjoint{<:Any,<:AbstractVecOrMat}) = A \ copy(B)
73+
\(A::TransposeFact{<:Any,<:BandedLU}, B::Transpose{<:Any,<:AbstractVecOrMat}) = A \ copy(B)
7474

7575
_factorize(::AbstractBandedLayout, _, A) = size(A,1) == size(A,2) ? lu(A) : qr(A)

src/symbanded/SplitCholesky.jl

Lines changed: 9 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -45,7 +45,13 @@ function splitcholesky!(::SymmetricLayout{<:BandedColumnMajor},
4545
SplitCholesky(A, A.uplo)
4646
end
4747

48-
adjoint(S::SplitCholesky) = Adjoint(S)
48+
if !(isdefined(LinearAlgebra, :AdjointFactorization)) # VERSION < v"1.10-"
49+
adjoint(S::SplitCholesky) = Adjoint(S)
50+
else
51+
transpose(S::SplitCholesky{<:Real}) = S'
52+
transpose(::SplitCholesky) =
53+
throw(ArgumentError("transpose of SplitCholesky decomposition is not supported, consider using adjoint"))
54+
end
4955

5056
function lmul!(S::SplitCholesky{T,Symmetric{T,M}}, B::AbstractVecOrMat{T}) where {T,M<:BandedMatrix{T}}
5157
require_one_based_indexing(B)
@@ -103,7 +109,7 @@ function ldiv!(S::SplitCholesky{T,Symmetric{T,M}}, B::AbstractVecOrMat{T}) where
103109
B
104110
end
105111

106-
function lmul!(S::Adjoint{T,SplitCholesky{T,Symmetric{T,M}}}, B::AbstractVecOrMat{T}) where {T,M<:BandedMatrix{T}}
112+
function lmul!(S::AdjointFact{T,SplitCholesky{T,Symmetric{T,M}}}, B::AbstractVecOrMat{T}) where {T,M<:BandedMatrix{T}}
107113
require_one_based_indexing(B)
108114
n, nrhs = size(B, 1), size(B, 2)
109115
if size(S, 1) != n
@@ -138,7 +144,7 @@ function lmul!(S::Adjoint{T,SplitCholesky{T,Symmetric{T,M}}}, B::AbstractVecOrMa
138144
B
139145
end
140146

141-
function ldiv!(S::Adjoint{T,SplitCholesky{T,Symmetric{T,M}}}, B::AbstractVecOrMat{T}) where {T,M<:BandedMatrix{T}}
147+
function ldiv!(S::AdjointFact{T,SplitCholesky{T,Symmetric{T,M}}}, B::AbstractVecOrMat{T}) where {T,M<:BandedMatrix{T}}
142148
require_one_based_indexing(B)
143149
n, nrhs = size(B, 1), size(B, 2)
144150
if size(S, 1) != n

0 commit comments

Comments
 (0)