:orphan:
# KSPFGMRES
Implements the Flexible Generalized Minimal Residual method. [](sec_flexibleksp) 
## 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
- ***-ksp_fgmres_modifypcnochange -*** do not change the preconditioner between iterations
- ***-ksp_fgmres_modifypcksp -*** modify the preconditioner using `KSPFGMRESModifyPCKSP()`





## Notes
See `KSPFGMRESSetModifyPC()` for how to vary the preconditioner between iterations

Only right preconditioning is supported.

The following options -ksp_type fgmres -pc_type ksp -ksp_ksp_type bcgs -ksp_view -ksp_pc_type jacobi make the preconditioner (or inner solver)
be bi-CG-stab with a preconditioner of Jacobi.


## Developer Note
This object is subclassed off of `KSPGMRES`


## Contributed by
Allison Baker


## See Also
 [](ch_ksp), [](sec_flexibleksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `KSP`, `KSPGMRES`, `KSPLGMRES`,
`KSPGMRESSetRestart()`, `KSPGMRESSetHapTol()`, `KSPGMRESSetPreAllocateVectors()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`,
`KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`,
`KSPGMRESCGSRefinementType`, `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESMonitorKrylov()`, `KSPFGMRESSetModifyPC()`,
`KSPFGMRESModifyPCKSP()`

## Level
beginner

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

## Examples
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/dm/impls/stag/tutorials/ex2.c.html">src/dm/impls/stag/tutorials/ex2.c</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/dm/impls/stag/tutorials/ex3.c.html">src/dm/impls/stag/tutorials/ex3.c</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/dm/impls/stag/tutorials/ex4.c.html">src/dm/impls/stag/tutorials/ex4.c</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/dm/impls/stag/tutorials/ex8.c.html">src/dm/impls/stag/tutorials/ex8.c</A><BR>


---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ksp/ksp/impls/gmres/fgmres/fgmres.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)  
