diff --git a/stdlib/LinearAlgebra/src/lu.jl b/stdlib/LinearAlgebra/src/lu.jl index d2e82af5d6409..dcbc5652fff15 100644 --- a/stdlib/LinearAlgebra/src/lu.jl +++ b/stdlib/LinearAlgebra/src/lu.jl @@ -521,6 +521,14 @@ function ldiv!(adjA::AdjointFactorization{<:Any,<:LU}, B::AbstractVecOrMat) _apply_inverse_ipiv_rows!(A, B) end +# To enable rdiv! via ldiv! +ldiv!(A::TransposeFactorization{T,<:LU{T,<:StridedMatrix}}, B::Transpose{T,<:StridedVecOrMat{T}}) where {T<:Union{BlasFloat,BlasComplex}} = + LAPACK.getrs!('T', A.parent.factors, A.parent.ipiv, copy(B)) +ldiv!(A::TransposeFactorization{T,<:LU{T,<:StridedMatrix}}, B::Transpose{T,<:Hermitian{T,<:StridedMatrix}}) where {T<:BlasComplex} = + LAPACK.getrs!('T', A.parent.factors, A.parent.ipiv, Matrix(B)) +ldiv!(A::TransposeFactorization{T,<:LU{T,<:StridedMatrix}}, B::Symmetric{T,<:StridedMatrix}) where {T<:BlasFloat} = + LAPACK.getrs!('T', A.parent.factors, A.parent.ipiv, Matrix(B)) + (\)(A::AdjointFactorization{T,<:LU{T,<:StridedMatrix}}, B::Adjoint{T,<:StridedVecOrMat{T}}) where {T<:BlasComplex} = LAPACK.getrs!('C', A.parent.factors, A.parent.ipiv, copy(B)) (\)(A::TransposeFactorization{T,<:LU{T,<:StridedMatrix}}, B::Transpose{T,<:StridedVecOrMat{T}}) where {T<:BlasFloat} = diff --git a/stdlib/LinearAlgebra/test/lu.jl b/stdlib/LinearAlgebra/test/lu.jl index 56a402d70493e..d21d36d9402ef 100644 --- a/stdlib/LinearAlgebra/test/lu.jl +++ b/stdlib/LinearAlgebra/test/lu.jl @@ -499,4 +499,12 @@ end end end +@testset "lu!/rdiv! via lu/ldiv! (#55422)" begin + id(A,T) = T(Matrix{eltype(A)}(I, size(A)...)) + A = randn(1000,1000); + @test inv(A) ≈ ldiv!(lu(A), id(A,Matrix)) ≈ rdiv!(id(A,Matrix), lu(A)) ≈ rdiv!(id(A,Symmetric), lu(A)) + A = complex.(randn(1000,1000),randn(1000,1000)); + @test inv(A) ≈ ldiv!(lu(A), id(A,Matrix)) ≈ rdiv!(id(A,Matrix), lu(A)) ≈ rdiv!(id(A,Hermitian), lu(A)) +end + end # module TestLU