:orphan: # KSPPGMRES Implements the Pipelined Generalized Minimal Residual method. [](sec_pipelineksp) ## 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 -*** 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 [](doc_faq_pipelined) ## Reference Ghysels, Ashby, Meerbergen, Vanroose, Hiding global communication latencies in the GMRES algorithm on massively parallel machines, 2012. ## Developer Note This object is subclassed off of `KSPGMRES` ## See Also [](ch_ksp), [](sec_pipelineksp), [](doc_faq_pipelined), `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 --- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ksp/ksp/impls/gmres/pgmres/pgmres.c) [Index of all KSP routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)