Skip to content

Commit 210c5b5

Browse files
Fix rdiv of complex lhs by real factorizations (#50671)
Co-authored-by: Martin Holters <martin@holters.name>
1 parent c664a57 commit 210c5b5

File tree

5 files changed

+42
-30
lines changed

5 files changed

+42
-30
lines changed

stdlib/LinearAlgebra/src/factorization.jl

Lines changed: 30 additions & 26 deletions
Original file line numberDiff line numberDiff line change
@@ -137,6 +137,12 @@ function Base.show(io::IO, ::MIME"text/plain", x::TransposeFactorization)
137137
show(io, MIME"text/plain"(), parent(x))
138138
end
139139

140+
function (\)(F::Factorization, B::AbstractVecOrMat)
141+
require_one_based_indexing(B)
142+
TFB = typeof(oneunit(eltype(F)) \ oneunit(eltype(B)))
143+
ldiv!(F, copy_similar(B, TFB))
144+
end
145+
(\)(F::TransposeFactorization, B::AbstractVecOrMat) = conj!(adjoint(F.parent) \ conj.(B))
140146
# With a real lhs and complex rhs with the same precision, we can reinterpret
141147
# the complex rhs as a real rhs with twice the number of columns or rows
142148
function (\)(F::Factorization{T}, B::VecOrMat{Complex{T}}) where {T<:BlasReal}
@@ -151,32 +157,6 @@ end
151157
(\)(F::AdjointFactorization{T}, B::VecOrMat{Complex{T}}) where {T<:BlasReal} =
152158
@invoke \(F::typeof(F), B::VecOrMat)
153159

154-
function (/)(B::VecOrMat{Complex{T}}, F::Factorization{T}) where {T<:BlasReal}
155-
require_one_based_indexing(B)
156-
x = rdiv!(copy(reinterpret(T, B)), F)
157-
return copy(reinterpret(Complex{T}, x))
158-
end
159-
# don't do the reinterpretation for [Adjoint/Transpose]Factorization
160-
(/)(B::VecOrMat{Complex{T}}, F::TransposeFactorization{T}) where {T<:BlasReal} =
161-
conj!(adjoint(parent(F)) \ conj.(B))
162-
(/)(B::VecOrMat{Complex{T}}, F::AdjointFactorization{T}) where {T<:BlasReal} =
163-
@invoke /(B::VecOrMat{Complex{T}}, F::Factorization{T})
164-
165-
function (\)(F::Factorization, B::AbstractVecOrMat)
166-
require_one_based_indexing(B)
167-
TFB = typeof(oneunit(eltype(F)) \ oneunit(eltype(B)))
168-
ldiv!(F, copy_similar(B, TFB))
169-
end
170-
(\)(F::TransposeFactorization, B::AbstractVecOrMat) = conj!(adjoint(F.parent) \ conj.(B))
171-
172-
function (/)(B::AbstractMatrix, F::Factorization)
173-
require_one_based_indexing(B)
174-
TFB = typeof(oneunit(eltype(B)) / oneunit(eltype(F)))
175-
rdiv!(copy_similar(B, TFB), F)
176-
end
177-
(/)(A::AbstractMatrix, F::AdjointFactorization) = adjoint(adjoint(F) \ adjoint(A))
178-
(/)(A::AbstractMatrix, F::TransposeFactorization) = transpose(transpose(F) \ transpose(A))
179-
180160
function ldiv!(Y::AbstractVector, A::Factorization, B::AbstractVector)
181161
require_one_based_indexing(Y, B)
182162
m, n = size(A)
@@ -200,3 +180,27 @@ function ldiv!(Y::AbstractMatrix, A::Factorization, B::AbstractMatrix)
200180
return ldiv!(A, Y)
201181
end
202182
end
183+
184+
function (/)(B::AbstractMatrix, F::Factorization)
185+
require_one_based_indexing(B)
186+
TFB = typeof(oneunit(eltype(B)) / oneunit(eltype(F)))
187+
rdiv!(copy_similar(B, TFB), F)
188+
end
189+
# reinterpretation trick for complex lhs and real factorization
190+
function (/)(B::Union{Matrix{Complex{T}},AdjOrTrans{Complex{T},Vector{Complex{T}}}}, F::Factorization{T}) where {T<:BlasReal}
191+
require_one_based_indexing(B)
192+
x = rdiv!(copy(reinterpret(T, B)), F)
193+
return copy(reinterpret(Complex{T}, x))
194+
end
195+
# don't do the reinterpretation for [Adjoint/Transpose]Factorization
196+
(/)(B::Union{Matrix{Complex{T}},AdjOrTrans{Complex{T},Vector{Complex{T}}}}, F::TransposeFactorization{T}) where {T<:BlasReal} =
197+
@invoke /(B::AbstractMatrix, F::Factorization)
198+
(/)(B::Matrix{Complex{T}}, F::AdjointFactorization{T}) where {T<:BlasReal} =
199+
@invoke /(B::AbstractMatrix, F::Factorization)
200+
(/)(B::Adjoint{Complex{T},Vector{Complex{T}}}, F::AdjointFactorization{T}) where {T<:BlasReal} =
201+
(F' \ B')'
202+
(/)(B::Transpose{Complex{T},Vector{Complex{T}}}, F::TransposeFactorization{T}) where {T<:BlasReal} =
203+
transpose(transpose(F) \ transpose(B))
204+
205+
rdiv!(B::AbstractMatrix, A::TransposeFactorization) = transpose(ldiv!(A.parent, transpose(B)))
206+
rdiv!(B::AbstractMatrix, A::AdjointFactorization) = adjoint(ldiv!(A.parent, adjoint(B)))

stdlib/LinearAlgebra/src/hessenberg.jl

Lines changed: 0 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -600,7 +600,6 @@ function rdiv!(B::AbstractVecOrMat{<:Complex}, F::Hessenberg{<:Complex,<:Any,<:A
600600
end
601601

602602
ldiv!(F::AdjointFactorization{<:Any,<:Hessenberg}, B::AbstractVecOrMat) = rdiv!(B', F')'
603-
rdiv!(B::AbstractMatrix, F::AdjointFactorization{<:Any,<:Hessenberg}) = ldiv!(F', B')'
604603

605604
det(F::Hessenberg) = det(F.H; shift=F.μ)
606605
logabsdet(F::Hessenberg) = logabsdet(F.H; shift=F.μ)

stdlib/LinearAlgebra/src/lu.jl

Lines changed: 0 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -709,8 +709,6 @@ function ldiv!(adjA::AdjointFactorization{<:Any,<:LU{T,Tridiagonal{T,V}}}, B::Ab
709709
end
710710

711711
rdiv!(B::AbstractMatrix, A::LU) = transpose(ldiv!(transpose(A), transpose(B)))
712-
rdiv!(B::AbstractMatrix, A::TransposeFactorization{<:Any,<:LU}) = transpose(ldiv!(A.parent, transpose(B)))
713-
rdiv!(B::AbstractMatrix, A::AdjointFactorization{<:Any,<:LU}) = adjoint(ldiv!(A.parent, adjoint(B)))
714712

715713
# Conversions
716714
AbstractMatrix(F::LU) = (F.L * F.U)[invperm(F.p),:]

stdlib/LinearAlgebra/test/hessenberg.jl

Lines changed: 3 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -178,8 +178,10 @@ let n = 10
178178
@test H \ B A \ B H \ complex(B)
179179
@test (H - I) \ B (A - I) \ B
180180
@test (H - (3+4im)I) \ B (A - (3+4im)I) \ B
181-
@test b' / H b' / A complex.(b') / H
181+
@test b' / H b' / A complex(b') / H
182+
@test transpose(b) / H transpose(b) / A transpose(complex(b)) / H
182183
@test B' / H B' / A complex(B') / H
184+
@test b' / H' complex(b)' / H'
183185
@test B' / (H - I) B' / (A - I)
184186
@test B' / (H - (3+4im)I) B' / (A - (3+4im)I)
185187
@test (H - (3+4im)I)' \ B (A - (3+4im)I)' \ B

stdlib/LinearAlgebra/test/lu.jl

Lines changed: 9 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -391,6 +391,15 @@ end
391391
B = randn(elty, 5, 5)
392392
@test rdiv!(transform(A), transform(lu(B))) transform(C) / transform(B)
393393
end
394+
for elty in (Float32, Float64, ComplexF64), transF in (identity, transpose),
395+
transB in (transpose, adjoint), transT in (identity, complex)
396+
A = randn(elty, 5, 5)
397+
F = lu(A)
398+
b = randn(transT(elty), 5)
399+
@test rdiv!(transB(copy(b)), transF(F)) transB(b) / transF(F) transB(b) / transF(A)
400+
B = randn(transT(elty), 5, 5)
401+
@test rdiv!(copy(B), transF(F)) B / transF(F) B / transF(A)
402+
end
394403
end
395404

396405
@testset "transpose(A) / lu(B)' should not overwrite A (#36657)" begin

0 commit comments

Comments
 (0)