# KSPGMRES
Implements the Generalized Minimal Residual method [1] with restart 
## Options Database Keys

- ***-ksp_gmres_restart <restart> -*** the number of Krylov directions to orthogonalize against
- ***-ksp_gmres_haptol <tol> -*** 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
Left and right preconditioning are supported, but not symmetric preconditioning.


## Reference

- ***[1] -*** YOUCEF SAAD AND MARTIN H. SCHULTZ, GMRES: A GENERALIZED MINIMAL RESIDUAL ALGORITHM FOR SOLVING NONSYMMETRIC LINEAR SYSTEMS.
SIAM J. ScI. STAT. COMPUT. Vo|. 7, No. 3, July 1986.



## See Also
 [](chapter_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPFGMRES`, `KSPLGMRES`,
`KSPGMRESSetRestart()`, `KSPGMRESSetHapTol()`, `KSPGMRESSetPreAllocateVectors()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
`KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`,
`KSPGMRESCGSRefinementType`, `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESMonitorKrylov()`, `KSPSetPCSide()`

## Level
beginner

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/impls/gmres/gmres.c.html#KSPGMRES">src/ksp/ksp/impls/gmres/gmres.c</A>

## Examples
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex62.c.html">src/ksp/ksp/tutorials/ex62.c.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex7.c.html">src/ksp/ksp/tutorials/ex7.c.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex78.c.html">src/ksp/ksp/tutorials/ex78.c.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex7f.F90.html">src/ksp/ksp/tutorials/ex7f.F90.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex8.c.html">src/ksp/ksp/tutorials/ex8.c.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/tao/pde_constrained/tutorials/hyperbolic.c.html">src/tao/pde_constrained/tutorials/hyperbolic.c.html</A><BR>


---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ksp/ksp/impls/gmres/gmres.c)


[Index of all KSP routines](index.md)  
[Table of Contents for all manual pages](/docs/manualpages/index.md)  
[Index of all manual pages](/docs/manualpages/singleindex.md)  
