Skip to content

Commit ebb87a3

Browse files
authored
Allow real adjoint cholmod solve with complex rhs (#327)
1 parent d4c36be commit ebb87a3

File tree

2 files changed

+17
-5
lines changed

2 files changed

+17
-5
lines changed

src/solvers/cholmod.jl

Lines changed: 8 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -1559,10 +1559,10 @@ end
15591559
(\)(L::Factor{T}, B::Dense{T}) where {T<:VTypes} = solve(CHOLMOD_A, L, B)
15601560
# Explicit typevars are necessary to avoid ambiguities with defs in linalg/factorizations.jl
15611561
# Likewise the two following explicit Vector and Matrix defs (rather than a single VecOrMat)
1562-
(\)(L::Factor{T}, B::Vector{Complex{T}}) where {T<:Float64} = complex.(L\real(B), L\imag(B))
1563-
(\)(L::Factor{T}, B::Matrix{Complex{T}}) where {T<:Float64} = complex.(L\real(B), L\imag(B))
1564-
(\)(L::Factor{T}, B::Adjoint{<:Any, <:Matrix{Complex{T}}}) where {T<:Float64} = complex.(L\real(B), L\imag(B))
1565-
(\)(L::Factor{T}, B::Transpose{<:Any, <:Matrix{Complex{T}}}) where {T<:Float64} = complex.(L\real(B), L\imag(B))
1562+
(\)(L::Factor{Float64}, B::Vector{Complex{Float64}}) = complex.(L\real(B), L\imag(B))
1563+
(\)(L::Factor{Float64}, B::Matrix{Complex{Float64}}) = complex.(L\real(B), L\imag(B))
1564+
(\)(L::Factor{Float64}, B::Adjoint{<:Any,Matrix{Complex{Float64}}}) = complex.(L\real(B), L\imag(B))
1565+
(\)(L::Factor{Float64}, B::Transpose{<:Any,Matrix{Complex{Float64}}}) = complex.(L\real(B), L\imag(B))
15661566

15671567
(\)(L::Factor{T}, b::StridedVector) where {T<:VTypes} = Vector(L\Dense{T}(b))
15681568
(\)(L::Factor{T}, B::StridedMatrix) where {T<:VTypes} = Matrix(L\Dense{T}(B))
@@ -1580,6 +1580,10 @@ end
15801580
\(adjL::AdjType{<:Any,<:Factor}, B::Sparse) = (L = adjL.parent; spsolve(CHOLMOD_A, L, B))
15811581
\(adjL::AdjType{<:Any,<:Factor}, B::SparseVecOrMat) = (L = adjL.parent; \(adjoint(L), Sparse(B)))
15821582

1583+
(\)(adjL::AdjType{Float64,Factor{Float64}}, B::Vector{Complex{Float64}}) = complex.(adjL\real(B), adjL\imag(B))
1584+
(\)(adjL::AdjType{Float64,Factor{Float64}}, B::Matrix{Complex{Float64}}) = complex.(adjL\real(B), adjL\imag(B))
1585+
(\)(adjL::AdjType{Float64,Factor{Float64}}, B::Adjoint{<:Any,Matrix{Complex{Float64}}}) = complex.(adjL\real(B), adjL\imag(B))
1586+
(\)(adjL::AdjType{Float64,Factor{Float64}}, B::Transpose{<:Any,Matrix{Complex{Float64}}}) = complex.(adjL\real(B), adjL\imag(B))
15831587
function \(adjL::AdjType{<:Any,<:Factor}, b::StridedVector)
15841588
L = adjL.parent
15851589
return Vector(solve(CHOLMOD_A, L, Dense(b)))

test/cholmod.jl

Lines changed: 9 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -723,8 +723,16 @@ end
723723

724724
@testset "Real factorization and complex rhs" begin
725725
A = sprandn(5, 5, 0.4) |> t -> t't + I
726-
B = complex.(randn(5, 2), randn(5, 2))
726+
B = complex.(randn(5, 5), randn(5, 5))
727+
b = B[:,1]
728+
@test cholesky(A)\b A\b
727729
@test cholesky(A)\B A\B
730+
@test cholesky(A)\B' A\B'
731+
@test cholesky(A)\transpose(B) A\transpose(B)
732+
@test cholesky(A)'\b copy(A')\b
733+
@test cholesky(A)'\B copy(A')\B
734+
@test cholesky(A)'\B' copy(A')\B'
735+
@test cholesky(A)'\transpose(B) copy(A')\transpose(B)
728736
end
729737

730738
@testset "Make sure that ldlt performs an LDLt (Issue #19032)" begin

0 commit comments

Comments
 (0)