KSPPGMRES#
Implements the Pipelined Generalized Minimal Residual method [GAMV13]. Pipelined Krylov Methods
Options Database Keys#
-ksp_gmres_restart
- the number of Krylov directions to orthogonalize against-ksp_gmres_haptol
- sets the tolerance for “happy ending” (exact convergence)-ksp_gmres_preallocate - preallocate all the Krylov search directions initially (otherwise groups of vectors are allocated as needed)
-ksp_gmres_classicalgramschmidt - use classical (unmodified) Gram-Schmidt to orthogonalize against the Krylov space (fast) (the default)
-ksp_gmres_modifiedgramschmidt - use modified Gram-Schmidt in the orthogonalization (more stable, but slower)
-ksp_gmres_cgs_refinement_type <refine_never,refine_ifneeded,refine_always> - determine if iterative refinement is used to increase the stability of the classical Gram-Schmidt orthogonalization.
-ksp_gmres_krylov_monitor - plot the Krylov space generated
Note#
MPI configuration may be necessary for reductions to make asynchronous progress, which is important for performance of pipelined methods. See What steps are necessary to make the pipelined solvers execute efficiently?
Developer Note#
This object is subclassed off of KSPGMRES
, see the source code in src/ksp/ksp/impls/gmres for comments on the structure of the code
References#
P. Ghysels, T.J. Ashby, K. Meerbergen, and W. Vanroose. Hiding global communication latency in the GMRES algorithm on massively parallel machines. SIAM Journal on Scientific Computing, 35(1):C48–C71, 2013.
See Also#
KSP: Linear System Solvers, Pipelined Krylov Methods, What steps are necessary to make the pipelined solvers execute efficiently?, KSPCreate()
, KSPSetType()
, KSPType
, KSP
, KSPGMRES
, KSPLGMRES
, KSPPIPECG
, KSPPIPECR
,
KSPGMRESSetRestart()
, KSPGMRESSetHapTol()
, KSPGMRESSetPreAllocateVectors()
, KSPGMRESSetOrthogonalization()
, KSPGMRESGetOrthogonalization()
,
KSPGMRESClassicalGramSchmidtOrthogonalization()
, KSPGMRESModifiedGramSchmidtOrthogonalization()
,
KSPGMRESCGSRefinementType
, KSPGMRESSetCGSRefinementType()
, KSPGMRESGetCGSRefinementType()
, KSPGMRESMonitorKrylov()
Level#
beginner
Location#
src/ksp/ksp/impls/gmres/pgmres/pgmres.c
Index of all KSP routines
Table of Contents for all manual pages
Index of all manual pages