Skip to content

Provide interface to the precs API of LinearSolve #116

New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Merged
merged 1 commit into from
Mar 18, 2025
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
20 changes: 20 additions & 0 deletions README.md
Original file line number Diff line number Diff line change
Expand Up @@ -56,6 +56,26 @@ p = aspreconditioner(ml)
c = cg(A, A*ones(1000), Pl = p)
```


### As a preconditioner with LinearSolve.jl

`RugeStubenPreconBuilder` and `SmoothedAggregationPreconBuilder` work with the
[`precs` API](https://docs.sciml.ai/LinearSolve/stable/basics/Preconditioners/#Specifying-Preconditioners)
of LinearSolve.jl

```julia
A = poisson( (100,100) )
u0= rand(size(A,1))
b=A*u0

prob = LinearProblem(A, b)
strategy = KrylovJL_CG(precs = RugeStubenPreconBuilder())
sol = solve(prob, strategy, atol=1.0e-14)

strategy = KrylovJL_CG(precs = SmoothedAggregationPreconBuilder())
sol = solve(prob, strategy, atol=1.0e-14)
```

## Features and Roadmap

This package currently supports:
Expand Down
7 changes: 6 additions & 1 deletion src/AlgebraicMultigrid.jl
Original file line number Diff line number Diff line change
Expand Up @@ -3,7 +3,9 @@ module AlgebraicMultigrid
using Reexport
using LinearAlgebra
using LinearSolve
using SparseArrays, Printf
using SparseArrays
using SparseArrays: AbstractSparseMatrixCSC
using Printf
@reexport import CommonSolve: solve, solve!, init
using Reexport

Expand Down Expand Up @@ -40,4 +42,7 @@ export fit_candidates, smoothed_aggregation
include("preconditioner.jl")
export aspreconditioner

include("precs.jl")
export SmoothedAggregationPreconBuilder, RugeStubenPreconBuilder

end # module
38 changes: 38 additions & 0 deletions src/precs.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,38 @@
"""
SmoothedAggregationPreconBuilder(;blocksize=1, kwargs...)

Return callable object constructing a left smoothed aggregation algebraic multigrid preconditioner
to be used with the `precs` API of LinearSolve.
"""
struct SmoothedAggregationPreconBuilder{Tk}
blocksize::Int
kwargs::Tk
end

function SmoothedAggregationPreconBuilder(; blocksize = 1, kwargs...)
return SmoothedAggregationPreconBuilder(blocksize, kwargs)
end

function (b::SmoothedAggregationPreconBuilder)(A::AbstractSparseMatrixCSC, p)
return (aspreconditioner(smoothed_aggregation(SparseMatrixCSC(A), Val{b.blocksize}; b.kwargs...)), I)
end


"""
RugeStubenPreconBuilder(;blocksize=1, kwargs...)

Return callable object constructing a left algebraic multigrid preconditioner after Ruge & Stüben
to be used with the `precs` API of LinearSolve.
"""
struct RugeStubenPreconBuilder{Tk}
blocksize::Int
kwargs::Tk
end

function RugeStubenPreconBuilder(; blocksize = 1, kwargs...)
return RugeStubenPreconBuilder(blocksize, kwargs)
end

function (b::RugeStubenPreconBuilder)(A::AbstractSparseMatrixCSC, p)
return (aspreconditioner(ruge_stuben(SparseMatrixCSC(A), Val{b.blocksize}; b.kwargs...)), I)
end
25 changes: 24 additions & 1 deletion test/runtests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -330,4 +330,27 @@ end
X = poisson(27_000)+24.0*I
ml = ruge_stuben(X)
b = rand(27_000)
@test AlgebraicMultigrid._solve(ml, b, reltol = 1e-10) ≈ X \ b rtol = 1e-10
@test AlgebraicMultigrid._solve(ml, b, reltol = 1e-10) ≈ X \ b rtol = 1e-10




# LinearSolve precs interface
@testset "LinearSolvePrecs" begin

for sz in [ (10,10), (20,20), (50,50)]
A = poisson(sz)
u0= ones(size(A,1))
b=A*u0

prob = LinearProblem(A, b)
strategy = KrylovJL_CG(precs = RugeStubenPreconBuilder())
@test solve(prob, strategy, atol=1.0e-14) ≈ u0 rtol = 1.0e-8


strategy = KrylovJL_CG(precs = SmoothedAggregationPreconBuilder())
@test solve(prob, strategy, atol=1.0e-14) ≈ u0 rtol = 1.0e-8

end

end
Loading