Browse Source

Add solve function that allocates the initial guess/output vector.

This patch adds a function
solve(::HYPRESolver, ::HYPREMatrix, ::HYPREVector) which initializes an
output vector initialized to 0 and pass it to solve!.
fe/copyto
Fredrik Ekre 3 years ago
parent
commit
11ddf91fd2
  1. 24
      src/solvers.jl
  2. 13
      test/runtests.jl

24
src/solvers.jl

@ -7,6 +7,30 @@ Abstract super type of all the wrapped HYPRE solvers.
""" """
abstract type HYPRESolver end abstract type HYPRESolver end
# Generic fallback allocating a zero vector as initial guess
# TODO: This should allocate x using the owned cols instead of rows of A/b, but currently
# it is assumed these are always equivalent.
"""
solve(solver::HYPRESolver, A::HYPREMatrix, b::HYPREVector)
Solve the linear system `A x = b` using `solver` and return the solution vector.
This method allocates the initial guess/output vector `x`, initialized to 0.
See also [`solve!`](@ref).
"""
solve(solver::HYPRESolver, A::HYPREMatrix, b::HYPREVector) = solve!(solver, zero(b), A, b)
"""
solve!(solver::HYPRESolver, x::HYPREVector, A::HYPREMatrix, b::HYPREVector)
Solve the linear system `A x = b` using `solver` with `x` as the initial guess.
See also [`solve`](@ref).
"""
solve!(pcg::HYPRESolver, x::HYPREVector, A::HYPREMatrix, ::HYPREVector)
############# #############
# BoomerAMG # # BoomerAMG #
############# #############

13
test/runtests.jl

@ -11,6 +11,7 @@ using SparseMatricesCSR
using Test using Test
MPI.Init() MPI.Init()
HYPRE_Init()
@testset "HYPREMatrix" begin @testset "HYPREMatrix" begin
H = HYPREMatrix() H = HYPREMatrix()
@ -264,6 +265,10 @@ end
copy!(x, x_h) copy!(x, x_h)
# Test result with direct solver # Test result with direct solver
@test x A \ b atol=tol @test x A \ b atol=tol
# Test without passing initial guess
x_h = HYPRE.solve(amg, A_h, b_h)
copy!(x, x_h)
@test x A \ b atol=tol
end end
@testset "(ParCSR)PCG" begin @testset "(ParCSR)PCG" begin
@ -282,6 +287,10 @@ end
copy!(x, x_h) copy!(x, x_h)
# Test result with direct solver # Test result with direct solver
@test x A \ b atol=tol @test x A \ b atol=tol
# Test without passing initial guess
x_h = HYPRE.solve(pcg, A_h, b_h)
copy!(x, x_h)
@test x A \ b atol=tol
# Solve with AMG preconditioner # Solve with AMG preconditioner
precond = HYPRE.BoomerAMG(; MaxIter = 1, Tol = 0.0) precond = HYPRE.BoomerAMG(; MaxIter = 1, Tol = 0.0)
pcg = HYPRE.PCG(; Tol = tol, Precond = precond) pcg = HYPRE.PCG(; Tol = tol, Precond = precond)
@ -290,4 +299,8 @@ end
copy!(x, x_h) copy!(x, x_h)
# Test result with direct solver # Test result with direct solver
@test x A \ b atol=tol @test x A \ b atol=tol
# Test without passing initial guess
x_h = HYPRE.solve(pcg, A_h, b_h)
copy!(x, x_h)
@test x A \ b atol=tol
end end

Loading…
Cancel
Save