Skip to content

Commit

Permalink
fix AMG with matrix input (#41)
Browse files Browse the repository at this point in the history
* fix AMG with matrix input

* bump version
  • Loading branch information
mohamed82008 authored Sep 1, 2023
1 parent 9c7deca commit 3e61526
Show file tree
Hide file tree
Showing 3 changed files with 25 additions and 5 deletions.
2 changes: 1 addition & 1 deletion Project.toml
Original file line number Diff line number Diff line change
@@ -1,6 +1,6 @@
name = "Preconditioners"
uuid = "af69fa37-3177-5a40-98ee-561f696e4fcd"
version = "0.6"
version = "0.6.1"

[deps]
AlgebraicMultigrid = "2169fc97-5a83-5252-b627-83903c6c433c"
Expand Down
10 changes: 8 additions & 2 deletions src/amg.jl
Original file line number Diff line number Diff line change
Expand Up @@ -25,13 +25,19 @@ AMGPreconditioner(A::AbstractMatrix; kwargs...) = AMGPreconditioner(RugeStuben,
AMGPreconditioner{T}(A::AbstractMatrix; kwargs...) where T = AMGPreconditioner(T, A; kwargs...)

@inline function \(p::AMGPreconditioner, b)
x = copy(b);
x = copy(b)
return ldiv!(x, AMG.Preconditioner(p.ml, p.cycle), b)
end
@inline *(p::AMGPreconditioner, b) = AMG.Preconditioner(p.ml, p.cycle) * b
@inline ldiv!(p::AMGPreconditioner, b) = b .= p \ b
@inline function ldiv!(x, p::AMGPreconditioner, b)
@inline function ldiv!(x::AbstractVector, p::AMGPreconditioner, b::AbstractVector)
x .= b
return ldiv!(x, AMG.Preconditioner(p.ml, p.cycle), b)
end
@inline function ldiv!(x::AbstractMatrix, p::AMGPreconditioner, b::AbstractMatrix)
foreach(zip(eachcol(x), eachcol(b))) do (_x, _b)
ldiv!(_x, p, _b)
end
return x
end
@inline mul!(b, p::AMGPreconditioner, x) = mul!(b, AMG.Preconditioner(p.ml, p.cycle), x)
18 changes: 16 additions & 2 deletions test/runtests.jl
Original file line number Diff line number Diff line change
@@ -1,9 +1,9 @@
using LinearAlgebra: I, diag, ldiv!, norm, Symmetric, Hermitian
using LinearAlgebra: I, diag, ldiv!, norm, Symmetric, Hermitian, eigen
using SparseArrays: sprand
using Preconditioners: CholeskyPreconditioner, DiagonalPreconditioner
using Preconditioners: RugeStuben, SmoothedAggregation
using Preconditioners: UpdatePreconditioner!, AMGPreconditioner
using IterativeSolvers: cg
using IterativeSolvers: cg, lobpcg
using Random: seed!
using Test: @test, @testset, @test_throws
using OffsetArrays: OffsetMatrix, OffsetVector
Expand Down Expand Up @@ -68,3 +68,17 @@ end
test_matrix(Symmetric(A), F, atol)
test_matrix(Hermitian(A), F, atol)
end

@testset "AMG + LOBPCG" begin
A = sprand(1000, 1000, 0.01)
A = A + A' + 30I
X0 = randn(1000,1)
tol = 1e-3
F = lobpcg(A, false, X0, tol=tol, maxiter=200, P=AMGPreconditioner{SmoothedAggregation}(A))
@test eigen(Matrix(A)).values[1] F.λ[1] atol = 1e-3

X0 = randn(1000,2)
tol = 1e-3
F = lobpcg(A, false, X0, tol=tol, maxiter=200, P=AMGPreconditioner{SmoothedAggregation}(A))
@test eigen(Matrix(A)).values[1:2] F.λ atol = 1e-3
end

2 comments on commit 3e61526

@mohamed82008
Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@JuliaRegistrator
Copy link

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Registration pull request created: JuliaRegistries/General/90656

After the above pull request is merged, it is recommended that a tag is created on this repository for the registered package version.

This will be done automatically if the Julia TagBot GitHub Action is installed, or can be done manually through the github interface, or via:

git tag -a v0.6.1 -m "<description of version>" 3e61526211524542584022aa9712851e28c657d0
git push origin v0.6.1

Please sign in to comment.