SciML / LinearSolve.jl

LinearSolve.jl: High-Performance Unified Interface for Linear Solvers in Julia. Easily switch between factorization and Krylov methods, add preconditioners, and all in one interface.
https://docs.sciml.ai/LinearSolve/stable/
Other
235 stars 51 forks source link

make preconditioners part of the solver rather than a random extra #514

Open oscardssmith opened 2 weeks ago

oscardssmith commented 2 weeks ago

This still needs tests, but I'm putting it up as a PR first to get feedback. The intended goal here is to make the dolinsolve functionality in OrdinaryDiffEq of being able to pass in a precs function that takes in A and possibly some parameters and remakes your preconditioner for the new A value.

Here's an example of the new functionality in action:

using LinearSolve, Krylov, KrylovPreconditioners, LinearAlgebra, SparseArrays
A, b = sprand(10,10, .5), rand(10);
prob = LinearProblem(A,b)
solver = KrylovJL_CG(precs=(A,p=nothing) -> (BlockJacobiPreconditioner(A, 2), I), ldiv=false)
cache = init(prob, solver)
solve!(cache)
reinit!(cache, A = sprand(10,10, .2), b = rand(10))
solve!(cache)
j-fu commented 2 weeks ago

... This is the one functionality I missed with LinearSolve. I had a slightly different view on this: see https://github.com/SciML/LinearSolve.jl/issues/308 . Unfortunately I had no time to work on this, so thank you for this PR!

As for the functionality: Sometimes, e.g. in Newton solvers, it makes sense to update the matrix and keep the preconditioner. The "old" preconditioner may still be good enough for the "new" matrix in the sense that few additional iterations are cheaper than the re-creation of the preconditioner (thinking e.g. of Algebraic multigrid etc.)

May be a kwarg keep_precs (with default value false) in reinit! can help to achieve this.

oscardssmith commented 2 weeks ago

@j-fu yeah, the only reason I thought about this API is that OrdinaryDiffEq has a (slightly broken very hacky) version of this API already (dolinsolve) which someone pointed out to me a few weeks ago that it's slightly broken.

j-fu commented 2 weeks ago

A couple of thoughts for this API extension:

May one could achieve both ways via a dispatch between Function and something like AbstractPreconditioningAlgorithm. Also this would give the possibility to have a couple of "curated" preconditioners predefined (e.g. in package extensions).

UPDATE: Well, this already works with this PR :)

struct ILUZero_ILU0 end
(::ILUZero_ILU0)(A, p=nothing)=(ilu0(A),I)

p=LinearProblem(A,b)
solve(p,KrylovJL_GMRES(precs=ILUZero_ILU0(),ldiv=true) )
oscardssmith commented 2 weeks ago

Yeah, we likely would want to add aliases to create iterative algorithms with builtin preconditioners. My reasoning behind precs rather than Pr/Pl is that if you want both preconditioners, it's likely that they will share computational work. We could make it so you are allowed to pass up to one of Pr, Pl, or precs (where passing Pr or Pl would set the other one to I)

j-fu commented 1 week ago

Let's have a chat on this at JuliaCon!