Skip to content
Open
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
4 changes: 4 additions & 0 deletions CHANGELOG.md
Original file line number Diff line number Diff line change
@@ -1,5 +1,9 @@
# Changelog

## [2.1.0]

- Add `SchurComplementPreconditioner` and `SchurComplementPreconBuilder`

## [2.0.0] - 2026-01-06

- Remove solver + precon API which is not based on precs or directly overloading `\`.
Expand Down
2 changes: 1 addition & 1 deletion Project.toml
Original file line number Diff line number Diff line change
@@ -1,6 +1,6 @@
name = "ExtendableSparse"
uuid = "95c220a8-a1cf-11e9-0c77-dbfce5f500b3"
version = "2.0.1"
version = "2.1.0"
authors = ["Juergen Fuhrmann <juergen.fuhrmann@wias-berlin.de>", "Daniel Runge"]

[deps]
Expand Down
4 changes: 2 additions & 2 deletions src/ExtendableSparse.jl
Original file line number Diff line number Diff line change
Expand Up @@ -7,7 +7,7 @@ module ExtendableSparse

using DocStringExtensions: DocStringExtensions, SIGNATURES, TYPEDEF, TYPEDFIELDS, TYPEDSIGNATURES
using ILUZero: ILUZero
using LinearAlgebra: LinearAlgebra, Diagonal, Hermitian, Symmetric, Tridiagonal, convert, mul!, ldiv!
using LinearAlgebra: LinearAlgebra, Diagonal, Hermitian, Symmetric, Tridiagonal, convert, mul!, ldiv!, lu
using SparseArrays: SparseArrays, AbstractSparseMatrix, AbstractSparseMatrixCSC, SparseMatrixCSC
using SparseArrays: dropzeros!, findnz, nzrange, sparse, spzeros, rowvals, getcolptr, nonzeros, nnz
using Sparspak: sparspaklu
Expand Down Expand Up @@ -42,7 +42,7 @@ export eliminate_dirichlet, eliminate_dirichlet!, mark_dirichlet


include("preconbuilders.jl")
export LinearSolvePreconBuilder, BlockPreconBuilder, JacobiPreconBuilder
export LinearSolvePreconBuilder, BlockPreconBuilder, JacobiPreconBuilder, SchurComplementPreconBuilder
@public ILUZeroPreconBuilder, ILUTPreconBuilder


Expand Down
92 changes: 92 additions & 0 deletions src/preconbuilders.jl
Original file line number Diff line number Diff line change
Expand Up @@ -294,3 +294,95 @@ LinearAlgebra.ldiv!(fact::JacobiPreconditioner, v) = ldiv!(fact.factorization, v

allow_views(::JacobiPreconditioner) = true
allow_views(::Type{JacobiPreconditioner}) = true


"""
SchurComplementPreconditioner{AT, ST}

Preconditioner for saddle-point problems of the form:
```
[ A B ]
[ Bᵀ 0 ]
```

# Fields
- `partitions`: Vector of two ranges defining matrix partitioning
- `A_fac::AT`: Factorization of the A-block (top-left)
- `S_fac::ST`: Factorization of the Schur complement (approximation of -BᵀA⁻¹B)
"""
struct SchurComplementPreconditioner{AT, ST}
partitions
A_fac::AT
S_fac::ST
end


"""
SchurComplementPreconBuilder(dofs_first_block, A_factorization, S_factorization = lu)

Factory function creating preconditioners for saddle-point problems.

# Arguments
- `dofs_first_block`: Number of degrees of freedom in the first block (size of A matrix)
- `A_factorization`: Factorization method for the A-block (e.g., `ilu0, Diagonal, lu`)
- `S_factorization`: Factorization method for the Schur complement (defaults to `lu`)

Returns a function `prec(M, p)` that creates a `SchurComplementPreconditioner`.
"""
function SchurComplementPreconBuilder(dofs_first_block, A_factorization, S_factorization = lu; verbosity = 0)

# this is the resulting preconditioner
function prec(M)

# We have
# M = [ A B ] → n1 dofs
# [ Bᵀ 0 ] → n2 dofs

n1 = dofs_first_block
n2 = size(M, 1) - n1

# I see no other way than creating this expensive copy
A = M[1:n1, 1:n1]
B = M[1:n1, (n1 + 1):end]

# first factorization
A_fac = A_factorization(A)

verbosity > 0 && @info "SchurComplementPreconBuilder: A ($n1×$n1) is factorized"

# compute the Schur Matrix
# S ≈ -BᵀA⁻¹B
# we use the diagonal of A: this is _very_ performant and creates a sparse result
S = - B' * (Diagonal(A) \ B)

verbosity > 0 && @info "SchurComplementPreconBuilder: S ($n2×$n2) is computed"

# factorize S
S_fac = S_factorization(S)

verbosity > 0 && @info "SchurComplementPreconBuilder: S ($n2×$n2) is factorized"

return SchurComplementPreconditioner([1:n1, (n1 + 1):(n1 + n2)], A_fac, S_fac)
end

return prec
end


function LinearAlgebra.ldiv!(u, p::SchurComplementPreconditioner, v)

(part1, part2) = p.partitions
@views ldiv!(u[part1], p.A_fac, v[part1])
@views ldiv!(u[part2], p.S_fac, v[part2])

return u
end

function LinearAlgebra.ldiv!(p::SchurComplementPreconditioner, v)

(part1, part2) = p.partitions
@views ldiv!(p.A_fac, v[part1])
@views ldiv!(p.S_fac, v[part2])

return v
end
20 changes: 17 additions & 3 deletions test/test_block.jl
Original file line number Diff line number Diff line change
@@ -1,12 +1,13 @@
module test_block
using Test
using AMGCLWrap
using ExtendableSparse
using ExtendableSparse: BlockPreconditioner, jacobi
using ExtendableSparse: BlockPreconditioner, jacobi, SchurComplementPreconBuilder
using ILUZero, AlgebraicMultigrid
using IterativeSolvers
using LinearAlgebra
using SparseArrays
using Sparspak
using AMGCLWrap
using Test

ExtendableSparse.allow_views(::typeof(ilu0)) = true

Expand All @@ -29,6 +30,19 @@ function main(; n = 100)
sol = cg(A, b, Pl = BlockPreconditioner(A; partitioning, factorization = sparspaklu))
@test sol ≈ sol0

# Schur complement: create a saddle point system
let
m = n ÷ 10
B = I[1:(n^2), 1:(m^2)]
M = [ A B; B' spzeros(m^2, m^2)]

sol1 = ones(n^2 + m^2)
c = M * sol1

sol = cg(M, c, Pl = SchurComplementPreconBuilder(n^2, lu)(M))
@test sol ≈ sol1
end

return
end

Expand Down
Loading