:orphan:
# KSPConvergedReason
reason a Krylov method was determined to have converged or diverged 
## Synopsis
```
typedef enum { /* converged */
  KSP_CONVERGED_RTOL_NORMAL               = 1,
  KSP_CONVERGED_ATOL_NORMAL               = 9,
  KSP_CONVERGED_RTOL                      = 2,
  KSP_CONVERGED_ATOL                      = 3,
  KSP_CONVERGED_ITS                       = 4,
  KSP_CONVERGED_NEG_CURVE                 = 5,
  KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED   = 5,
  KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED = 6,
  KSP_CONVERGED_STEP_LENGTH               = 6,
  KSP_CONVERGED_HAPPY_BREAKDOWN           = 7,
  /* diverged */
  KSP_DIVERGED_NULL                      = -2,
  KSP_DIVERGED_ITS                       = -3,
  KSP_DIVERGED_DTOL                      = -4,
  KSP_DIVERGED_BREAKDOWN                 = -5,
  KSP_DIVERGED_BREAKDOWN_BICG            = -6,
  KSP_DIVERGED_NONSYMMETRIC              = -7,
  KSP_DIVERGED_INDEFINITE_PC             = -8,
  KSP_DIVERGED_NANORINF                  = -9,
  KSP_DIVERGED_INDEFINITE_MAT            = -10,
  KSP_DIVERGED_PC_FAILED                 = -11,
  KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11,

  KSP_CONVERGED_ITERATING = 0
} KSPConvergedReason;
```

## Values

- ***`KSP_CONVERGED_RTOL_NORMAL` -*** requested decrease in the residual for the normal equations
- ***`KSP_CONVERGED_ATOL_NORMAL` -*** requested absolute value in the residual for the normal equations
- ***`KSP_CONVERGED_RTOL` -*** requested decrease in the residual
- ***`KSP_CONVERGED_ATOL` -*** requested absolute value in the residual
- ***`KSP_CONVERGED_ITS` -*** requested number of iterations
- ***`KSP_CONVERGED_NEG_CURVE` -*** see note below
- ***`KSP_CONVERGED_STEP_LENGTH` -*** see note below
- ***`KSP_CONVERGED_HAPPY_BREAKDOWN` -*** happy breakdown (meaning early convergence of the `KSPType` occurred.
- ***`KSP_DIVERGED_NULL` -*** breakdown when solving the Hessenberg system within GMRES
- ***`KSP_DIVERGED_ITS` -*** requested number of iterations
- ***`KSP_DIVERGED_DTOL` -*** large increase in the residual norm
- ***`KSP_DIVERGED_BREAKDOWN` -*** breakdown in the Krylov method
- ***`KSP_DIVERGED_BREAKDOWN_BICG` -*** breakdown in the `KSPBGCS` Krylov method
- ***`KSP_DIVERGED_NONSYMMETRIC` -*** the operator or preonditioner was not symmetric for a `KSPType` that requires symmetry
- ***`KSP_DIVERGED_INDEFINITE_PC` -*** the preconditioner was indefinite for a `KSPType` that requires it be definite
- ***`KSP_DIVERGED_NANORINF` -*** a not a number of infinity was detected in a vector during the computation
- ***`KSP_DIVERGED_INDEFINITE_MAT` -*** the operator was indefinite for a `KSPType` that requires it be definite
- ***`KSP_DIVERGED_PC_FAILED` -*** the action of the preconditioner failed for some reason





## Note
The values  `KSP_CONVERGED_NEG_CURVE`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by the special `KSPNASH`,
`KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver.


## Developer Notes
This must match the values in petsc/finclude/petscksp.h

The string versions of these are `KSPConvergedReasons`; if you change
any of the values here also change them that array of names.


## See Also
 [](ch_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()`

## Level
beginner

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/include/petscksp.h.html#KSPConvergedReason">include/petscksp.h</A>

## Examples
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex72.c.html">src/ksp/ksp/tutorials/ex72.c</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ksp/ksp/tutorials/ex77.c.html">src/ksp/ksp/tutorials/ex77.c</A><BR>


---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/include/petscksp.h)


[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)  
