From 026783757650a617a3a55ccde314bee090011bb5 Mon Sep 17 00:00:00 2001 From: Fredrik Ekre Date: Sun, 24 Jul 2022 19:12:22 +0200 Subject: [PATCH] Add solve(!) methods for using SparseMatrixCS(C|R) directly. --- src/solvers.jl | 18 ++++++++++++++++++ test/runtests.jl | 16 ++++++++++++++++ 2 files changed, 34 insertions(+) diff --git a/src/solvers.jl b/src/solvers.jl index e134981..05989a1 100644 --- a/src/solvers.jl +++ b/src/solvers.jl @@ -51,6 +51,24 @@ function solve!(solver::HYPRESolver, x::PVector, A::PSparseMatrix, b::PVector) return x end +######################################## +# SparseMatrixCS(C|R) solver interface # +######################################## + +# TODO: This could use the HYPRE compile flag for sequential mode to avoid MPI overhead + +function solve(solver::HYPRESolver, A::Union{SparseMatrixCSC,SparseMatrixCSR}, b::Vector) + hypre_x = solve(solver, HYPREMatrix(A), HYPREVector(b)) + x = copy!(copy(b), hypre_x) + return x +end +function solve!(solver::HYPRESolver, x::Vector, A::Union{SparseMatrixCSC,SparseMatrixCSR}, b::Vector) + hypre_x = HYPREVector(x) + solve!(solver, hypre_x, HYPREMatrix(A), HYPREVector(b)) + copy!(x, hypre_x) + return x +end + ##################################### ## Concrete solver implementations ## diff --git a/test/runtests.jl b/test/runtests.jl index 6a842b1..7d64bb7 100644 --- a/test/runtests.jl +++ b/test/runtests.jl @@ -454,3 +454,19 @@ end x_p = HYPRE.solve(pcg, A_p, b_p) @test tomain(x_p) ≈ A \ b atol=tol end + +@testset "solve with SparseMatrixCS(C|R)" begin + # Setup + A = sprand(100, 100, 0.05); A = A'A + 5I + b = rand(100) + x = zeros(100) + # Solve + tol = 1e-9 + pcg = HYPRE.PCG(; Tol = tol) + ## solve! + HYPRE.solve!(pcg, x, A, b) + @test x ≈ A \ b atol=tol + ## solve + x = HYPRE.solve(pcg, A, b) + @test x ≈ A \ b atol=tol +end