Skip to content

Commit af48e04

Browse files
Add user-defined near null space to aggregation-based AMG (#118)
* add near null space (B) as kwarg to solve * add test coverage for near null space via beam equation discretization
1 parent cf58288 commit af48e04

File tree

9 files changed

+304
-60
lines changed

9 files changed

+304
-60
lines changed

.gitignore

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -2,3 +2,4 @@
22
*.jl.*.cov
33
*.jl.mem
44
/Manifest.toml
5+
.vscode

src/AlgebraicMultigrid.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -9,7 +9,7 @@ using Printf
99
@reexport import CommonSolve: solve, solve!, init
1010
using Reexport
1111

12-
using LinearAlgebra: rmul!
12+
using LinearAlgebra: rmul!, qr
1313

1414
include("utils.jl")
1515
export approximate_spectral_radius

src/aggregation.jl

Lines changed: 36 additions & 56 deletions
Original file line numberDiff line numberDiff line change
@@ -1,5 +1,6 @@
1-
function smoothed_aggregation(A::TA,
1+
function smoothed_aggregation(A::TA,
22
::Type{Val{bs}}=Val{1};
3+
B = nothing,
34
symmetry = HermitianSymmetry(),
45
strength = SymmetricStrength(),
56
aggregate = StandardAggregation(),
@@ -12,10 +13,9 @@ function smoothed_aggregation(A::TA,
1213
diagonal_dominance = false,
1314
keep = false,
1415
coarse_solver = Pinv, kwargs...) where {T,V,bs,TA<:SparseMatrixCSC{T,V}}
15-
1616
n = size(A, 1)
17-
# B = kron(ones(n, 1), eye(1))
18-
B = ones(T,n)
17+
B = isnothing(B) ? ones(T,n) : copy(B)
18+
@assert size(A, 1) == size(B, 1)
1919

2020
#=max_levels, max_coarse, strength =
2121
levelize_strength_or_aggregation(max_levels, max_coarse, strength)
@@ -70,7 +70,7 @@ function extend_hierarchy!(levels, strength, aggregate, smooth,
7070
# b = zeros(eltype(A), size(A, 1))
7171

7272
# Improve candidates
73-
b = zeros(size(A,1))
73+
b = zeros(size(A,1),size(B,2))
7474
improve_candidates(A, B, b)
7575
T, B = fit_candidates(AggOp, B)
7676

@@ -88,59 +88,39 @@ function extend_hierarchy!(levels, strength, aggregate, smooth,
8888
end
8989
construct_R(::HermitianSymmetry, P) = P'
9090

91-
function fit_candidates(AggOp, B, tol = 1e-10)
92-
91+
function fit_candidates(AggOp, B; tol=1e-10)
9392
A = adjoint(AggOp)
94-
n_fine, n_coarse = size(A)
95-
n_col = n_coarse
96-
97-
R = zeros(eltype(B), n_coarse)
98-
Qx = zeros(eltype(B), nnz(A))
99-
# copy!(Qx, B)
100-
for i = 1:size(Qx, 1)
101-
Qx[i] = B[i]
102-
end
103-
# copy!(A.nzval, B)
104-
for i = 1:n_col
105-
for j in nzrange(A,i)
106-
row = A.rowval[j]
107-
A.nzval[j] = B[row]
108-
end
109-
end
110-
k = 1
111-
for i = 1:n_col
112-
norm_i = norm_col(A, Qx, i)
113-
threshold_i = tol * norm_i
114-
if norm_i > threshold_i
115-
scale = 1 / norm_i
116-
R[i] = norm_i
117-
else
118-
scale = 0
119-
R[i] = 0
93+
n_fine, m = ndims(B) == 1 ? (length(B), 1) : size(B)
94+
n_fine2, n_agg = size(A)
95+
@assert n_fine2 == n_fine
96+
n_coarse = m * n_agg
97+
T = eltype(B)
98+
Qs = spzeros(T, n_fine, n_coarse)
99+
R = zeros(T, n_coarse, m)
100+
101+
for agg in 1:n_agg
102+
rows = A.rowval[A.colptr[agg]:A.colptr[agg+1]-1]
103+
M = @view B[rows, :] # size(rows) × m
104+
105+
106+
F = qr(M)
107+
r = min(length(rows), m)
108+
Qfull = Matrix(F.Q)
109+
Qj = Qfull[:, 1:r]
110+
Rj = F.R
111+
112+
offset = (agg - 1) * m
113+
114+
for local_i in 1:length(rows), local_j in 1:r
115+
val = Qj[local_i, local_j]
116+
if abs(val) >= tol
117+
Qs[rows[local_i], offset+local_j] = val
118+
end
120119
end
121-
for j in nzrange(A, i)
122-
row = A.rowval[j]
123-
# Qx[row] *= scale
124-
#@show k
125-
# Qx[k] *= scale
126-
# k += 1
127-
A.nzval[j] *= scale
128-
end
129-
end
120+
dropzeros!(Qs)
130121

131-
# SparseMatrixCSC(size(A)..., A.colptr, A.rowval, Qx), R
132-
A, R
133-
end
134-
function norm_col(A, Qx, i)
135-
s = zero(eltype(A))
136-
for j in nzrange(A, i)
137-
if A.rowval[j] > length(Qx)
138-
val = 1
139-
else
140-
val = Qx[A.rowval[j]]
141-
end
142-
# val = A.nzval[A.rowval[j]]
143-
s += val*val
122+
R[offset+1:offset+r, :] .= Rj[1:r, :]
144123
end
145-
sqrt(s)
124+
125+
return Qs, R
146126
end

src/classical.jl

Lines changed: 4 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -9,6 +9,10 @@ function ruge_stuben(_A::Union{TA, Symmetric{Ti, TA}, Hermitian{Ti, TA}},
99
max_coarse = 10,
1010
coarse_solver = Pinv, kwargs...) where {Ti,Tv,bs,TA<:SparseMatrixCSC{Ti,Tv}}
1111

12+
13+
# fails if near null space `B` is provided
14+
haskey(kwargs, :B) && kwargs[:B] !== nothing && error("near null space `B` is only supported for smoothed aggregation AMG, not Ruge-Stüben AMG.")
15+
1216
if _A isa Symmetric && Ti <: Real || _A isa Hermitian
1317
A = _A.data
1418
symmetric = true

src/multilevel.jl

Lines changed: 8 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -187,7 +187,13 @@ Keyword Arguments
187187
* maxiter::Int64 - maximum number of iterations to execute
188188
* verbose::Bool - display residual at each iteration
189189
* log::Bool - return vector of residuals along with solution
190+
* B::AbstractArray - the **near null space** in SA-AMG, which represents the low energy that cannot be attenuated by relaxtion, and thus needs to be perserved across the coarse grid.
190191
192+
!!! note
193+
`B` can be:
194+
- a `Vector` (e.g., for scalar PDEs),
195+
- a `Matrix` (e.g., for vector PDEs or systems with multiple equations),
196+
If `B` is not provided, it defaults to a vector of ones.
191197
"""
192198
function _solve(ml::MultiLevel, b::AbstractArray, args...; kwargs...)
193199
n = length(ml) == 1 ? size(ml.final_A, 1) : size(ml.levels[1].A, 1)
@@ -296,9 +302,9 @@ end
296302
function init(::RugeStubenAMG, A, b, args...; kwargs...)
297303
AMGSolver(ruge_stuben(A; kwargs...), b)
298304
end
299-
function init(::SmoothedAggregationAMG, A, b; kwargs...)
305+
function init(sa::SmoothedAggregationAMG, A, b; kwargs...)
300306
AMGSolver(smoothed_aggregation(A; kwargs...), b)
301307
end
302308
function solve!(solt::AMGSolver, args...; kwargs...)
303309
_solve(solt.ml, solt.b, args...; kwargs...)
304-
end
310+
end

src/smoother.jl

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -29,6 +29,7 @@ end
2929
function gs!(A, b, x, start, step, stop)
3030
n = size(A, 1)
3131
z = zero(eltype(A))
32+
@assert size(x,2) == size(b, 2) "x and b must have the same number of columns"
3233
@inbounds for col in 1:size(x, 2)
3334
for i in start:step:stop
3435
rsum = z

test/lin_elastic_2d.jld2

55.2 KB
Binary file not shown.

0 commit comments

Comments
 (0)