# KSPGMRESClassicalGramSchmidtOrthogonalization
This is the basic orthogonalization routine using classical Gram-Schmidt with possible iterative refinement to improve the stability 
## Synopsis
```
PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP ksp, PetscInt it)
```
Collective


## Input Parameters

- ***ksp -*** KSP object, must be associated with `KSPGMRES`, `KSPFGMRES`, or `KSPLGMRES` Krylov method
- ***its -*** one less then the current GMRES restart iteration, i.e. the size of the Krylov space



## Options Database Keys

- ***-ksp_gmres_classicalgramschmidt -*** Activates `KSPGMRESClassicalGramSchmidtOrthogonalization()`
- ***-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.





## Notes
Use `KSPGMRESSetCGSRefinementType()` to determine if iterative refinement is to be used.
This is much faster than `KSPGMRESModifiedGramSchmidtOrthogonalization()` but has the small possibility of stability issues
that can usually be handled by using a a single step of iterative refinement with `KSPGMRESSetCGSRefinementType()`


## See Also
 [](chapter_ksp), `KSPGMRESSetOrthogonalization()`, `KSPGMRESSetCGSRefinementType()`,
`KSPGMRESGetCGSRefinementType()`, `KSPGMRESGetOrthogonalization()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()`

## Level
intermediate

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

## Examples
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex5f.F90.html">src/ksp/ksp/tutorials/ex5f.F90.html</A><BR>


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