From: Timo Heister Date: Thu, 27 Feb 2020 20:17:04 +0000 (-0500) Subject: amend FGMRES documentation X-Git-Tag: v9.2.0-rc1~489^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=fc9f6b7c7c97bd584252e0ddcff8f2c9515c3da7;p=dealii.git amend FGMRES documentation --- diff --git a/doc/doxygen/references.bib b/doc/doxygen/references.bib index 8817762b27..b3ea97f3d6 100644 --- a/doc/doxygen/references.bib +++ b/doc/doxygen/references.bib @@ -483,3 +483,13 @@ MRREVIEWER = {Jose Luis Gracia}, pages={3596--3604}, year={2018} } + +@TechReport{Saad1991, + Title = {{A {F}lexible {I}nner-{O}uter {P}reconditioned {GMRES} {A}lgorithm}}, + Author = {Y. Saad}, + Institution = {Minnesota Supercomputer Institute}, + Year = {1991}, + Address = {University of Minnesota}, + Number = {91-279}, + Url = {citeseer.ist.psu.edu/saad93flexible.html} +} diff --git a/include/deal.II/lac/solver_gmres.h b/include/deal.II/lac/solver_gmres.h index 09eec56065..42c4f6f2a7 100644 --- a/include/deal.II/lac/solver_gmres.h +++ b/include/deal.II/lac/solver_gmres.h @@ -440,19 +440,23 @@ protected: /** * Implementation of the Generalized minimal residual method with flexible - * preconditioning method. + * preconditioning (flexible GMRES or FGMRES). * - * This version of the GMRES method allows for the use of a different + * This flexible version of the GMRES method allows for the use of a different * preconditioner in each iteration step. Therefore, it is also more robust * with respect to inaccurate evaluation of the preconditioner. An important - * application is also the use of a Krylov space method inside the + * application is the use of a Krylov space method inside the * preconditioner. As opposed to SolverGMRES which allows one to choose * between left and right preconditioning, this solver always applies the * preconditioner from the right. * * FGMRES needs two vectors in each iteration steps yielding a total of * 2*SolverFGMRES::AdditionalData::max_basis_size+1 auxiliary - * vectors. + * vectors. Otherwise, FGMRES requires roughly the same number of operations + * per iteration compared to GMRES, except one application of the + * preconditioner less at each restart and at the end of solve(). + * + * For more details see @cite Saad1991. * * @author Guido Kanschat, 2003 */