From 8ddea4c94c65deff451b8e5b7b3e956025ebd563 Mon Sep 17 00:00:00 2001 From: Fredrik Ekre Date: Wed, 27 Jul 2022 16:26:23 +0200 Subject: [PATCH] Automatically apply required default settings when using solvers as preconditioner. --- gen/solver_options.jl | 1 + src/Internals.jl | 1 + src/solver_options.jl | 3 +++ src/solvers.jl | 9 +++++++++ test/runtests.jl | 8 ++++++++ 5 files changed, 22 insertions(+) diff --git a/gen/solver_options.jl b/gen/solver_options.jl index b8a01ba..082c2f1 100644 --- a/gen/solver_options.jl +++ b/gen/solver_options.jl @@ -28,6 +28,7 @@ function generate_options(io, structname, prefixes...) print(io, " $(first ? "" : "else")if k === :$(k)") println(io) if k == "Precond" + println(io, " Internals.set_precond_defaults(v)") println(io, " Internals.set_precond(s, v)") elseif nargs == 1 println(io, " @check ", n, "(solver)") diff --git a/src/Internals.jl b/src/Internals.jl index 208e8e3..c733d5b 100644 --- a/src/Internals.jl +++ b/src/Internals.jl @@ -12,5 +12,6 @@ function set_options end function solve_func end function setup_func end function set_precond end +function set_precond_defaults end end # module Internals diff --git a/src/solver_options.jl b/src/solver_options.jl index 11248ed..9fc7287 100644 --- a/src/solver_options.jl +++ b/src/solver_options.jl @@ -18,6 +18,7 @@ function Internals.set_options(s::BiCGSTAB, kwargs) elseif k === :MinIter @check HYPRE_ParCSRBiCGSTABSetMinIter(solver, v) elseif k === :Precond + Internals.set_precond_defaults(v) Internals.set_precond(s, v) elseif k === :PrintLevel @check HYPRE_ParCSRBiCGSTABSetPrintLevel(solver, v) @@ -308,6 +309,7 @@ function Internals.set_options(s::GMRES, kwargs) elseif k === :MinIter @check HYPRE_ParCSRGMRESSetMinIter(solver, v) elseif k === :Precond + Internals.set_precond_defaults(v) Internals.set_precond(s, v) elseif k === :PrintLevel @check HYPRE_ParCSRGMRESSetPrintLevel(solver, v) @@ -341,6 +343,7 @@ function Internals.set_options(s::PCG, kwargs) elseif k === :MaxIter @check HYPRE_ParCSRPCGSetMaxIter(solver, v) elseif k === :Precond + Internals.set_precond_defaults(v) Internals.set_precond(s, v) elseif k === :PrintLevel @check HYPRE_ParCSRPCGSetPrintLevel(solver, v) diff --git a/src/solvers.jl b/src/solvers.jl index be03d97..116532b 100644 --- a/src/solvers.jl +++ b/src/solvers.jl @@ -7,6 +7,9 @@ Abstract super type of all the wrapped HYPRE solvers. """ abstract type HYPRESolver end +# Fallback for the solvers that doesn't have required defaults +Internals.set_precond_defaults(::HYPRESolver) = nothing + # Generic fallback allocating a zero vector as initial guess # TODO: This should allocate x using the owned cols instead of rows of A/b, but currently # it is assumed these are always equivalent. @@ -142,6 +145,12 @@ end Internals.solve_func(::BoomerAMG) = HYPRE_BoomerAMGSolve Internals.setup_func(::BoomerAMG) = HYPRE_BoomerAMGSetup +function Internals.set_precond_defaults(amg::BoomerAMG) + defaults = (; Tol = 0.0, MaxIter = 1) + Internals.set_options(amg, pairs(defaults)) + return nothing +end + ######### # GMRES # diff --git a/test/runtests.jl b/test/runtests.jl index b0a71e4..f07d05a 100644 --- a/test/runtests.jl +++ b/test/runtests.jl @@ -312,6 +312,14 @@ end x_h = HYPRE.solve(bicg, A_h, b_h) copy!(x, x_h) @test x ≈ A \ b atol=tol + # Tests Internals.set_precond_defaults for BoomerAMG + precond = HYPRE.BoomerAMG() + bicg = HYPRE.BiCGSTAB(; Tol = tol, Precond = precond) + x_h = HYPREVector(zeros(100)) + HYPRE.solve!(bicg, x_h, A_h, b_h) + copy!(x, x_h) + # Test result with direct solver + @test x ≈ A \ b atol=tol end @testset "BoomerAMG" begin