KSPAGMRES#

Newton basis GMRES implementation with adaptive augmented eigenvectors The techniques used are best described in [WE11]. The contribution of this work is that it combines many of the previous work to reduce the amount of MPI messages and improve the robustness of the global approach by using deflation techniques. Please see [WE11] for numerical experiments and [WP13] for a description of these problems. There are many ongoing work that aim at avoiding (or minimizing) the communication in Krylov subspace methods. This code can be used as an experimental framework to combine several techniques in the particular case of GMRES. For instance, the computation of the shifts can be improved with techniques described in [PR12]. The orthogonalization technique can be replaced by TSQR [DGHL12]. The generation of the basis can be done using s-steps approaches[MHDY09]. See also [Sid97] and [BHR94].

Options Database Keys#

  • -ksp_gmres_restart - the number of Krylov directions

  • -ksp_gmres_krylov_monitor - plot the Krylov space generated

  • -ksp_agmres_eigen - Number of eigenvalues to deflate (Number of vectors to augment)

  • -ksp_agmres_maxeigen <max_neig> - Maximum number of eigenvalues to deflate

  • -ksp_agmres_MinRatio <1> - Relaxation parameter in the adaptive strategy; smallest multiple of the remaining number of steps allowed

  • -ksp_agmres_MaxRatio <1> - Relaxation parameter in the adaptive strategy; Largest multiple of the remaining number of steps allowed

  • -ksp_agmres_DeflPrecond - Apply deflation as a preconditioner, this is similar to KSPDGMRES but it rather builds a Newton basis.

  • -ksp_dgmres_force <0, 1> - Force the deflation at each restart.

Note#

Left and right preconditioning are supported, but not symmetric preconditioning. Complex arithmetic is not supported

Developer Note#

This object is subclassed off of KSPDGMRES, see the source code in src/ksp/ksp/impls/gmres for comments on the structure of the code

Contributed by#

Desire NUENTSA WAKAM, INRIA desire.nuentsa_wakam@inria.fr with inputs from Guy Atenekeng atenekeng@yahoo.com and R.B. Sidje roger.b.sidje@ua.edu

References#

[BHR94]

Zhaojun Bai, Dan Hu, and Lothar Reichel. A Newton basis GMRES implementation. IMA Journal of Numerical Analysis, 14(4):563–581, 1994.

[DGHL12]

James Demmel, Laura Grigori, Mark Hoemmen, and Julien Langou. Communication-optimal parallel and sequential QR and LU factorizations. SIAM Journal on Scientific Computing, 34(1):A206–A239, 2012.

[MHDY09]

Marghoob Mohiyuddin, Mark Hoemmen, James Demmel, and Katherine Yelick. Minimizing communication in sparse matrix solvers. In Proceedings of the Conference on High Performance Computing Networking, Storage and Analysis, 1–12. 2009.

[PR12]

Bernard Philippe and Lothar Reichel. On the generation of Krylov subspace bases. Applied Numerical Mathematics, 62(9):1171–1186, 2012.

[Sid97]

Roger B Sidje. Alternatives for parallel krylov subspace basis computation. Numerical linear algebra with applications, 4(4):305–331, 1997.

[WE11] (1,2)

Désiré Nuentsa Wakam and Jocelyne Erhel. Parallelism and robustness in GMRES with the Newton basis and the deflated restarting. 2011. URL: https://hal.inria.fr/inria-00638247/en.

[WP13]

Désiré Nuentsa Wakam and François Pacull. Memory efficient hybrid algebraic solvers for linear systems arising from compressible flows. Computers & Fluids, 80:158–167, 2013.

See Also#

KSP: Linear System Solvers, KSPCreate(), KSPSetType(), KSPType, KSP, KSPDGMRES, KSPPGMRES, KSPGMRESSetRestart(), KSPGMRESSetHapTol(), KSPGMRESSetPreAllocateVectors(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESModifiedGramSchmidtOrthogonalization(), KSPGMRESCGSRefinementType, KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESMonitorKrylov(), KSPSetPCSide()

Level#

intermediate

Location#

src/ksp/ksp/impls/gmres/agmres/agmres.c


Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages