diff --git a/src/algebra/LinearSolver.cpp b/src/algebra/LinearSolver.cpp index a275cd5816d6ce04715fd1d160008b4dd8508051..6a6beaa05bd3ac46b1a81f6e8645721003fe5577 100644 --- a/src/algebra/LinearSolver.cpp +++ b/src/algebra/LinearSolver.cpp @@ -88,6 +88,7 @@ struct LinearSolver::Internals { switch (precond) { case LSPrecond::none: + case LSPrecond::amg: case LSPrecond::diagonal: case LSPrecond::incomplete_choleski: case LSPrecond::incomplete_LU: { @@ -244,6 +245,10 @@ struct LinearSolver::Internals if (not direct_solver) { switch (options.precond()) { + case LSPrecond::amg: { + PCSetType(pc, PCGAMG); + break; + } case LSPrecond::diagonal: { PCSetType(pc, PCJACOBI); break; diff --git a/src/algebra/LinearSolverOptions.hpp b/src/algebra/LinearSolverOptions.hpp index bd7538b53df8740ee6c151e0d2e72a106fdc3c8c..014ee86d64bee4d2a128e820ab8876450a195369 100644 --- a/src/algebra/LinearSolverOptions.hpp +++ b/src/algebra/LinearSolverOptions.hpp @@ -37,6 +37,7 @@ enum class LSPrecond : int8_t diagonal, incomplete_choleski, incomplete_LU, + amg, // LS__end }; @@ -101,6 +102,9 @@ name(const LSPrecond precond) case LSPrecond::incomplete_LU: { return "ILU"; } + case LSPrecond::amg: { + return "AMG"; + } case LSPrecond::LS__end: { } }