:orphan:
# PetscConvEstGetConvRate
Returns an estimate of the convergence rate for the discretization 
## Synopsis
```
#include "petscconvest.h" 
PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[])
```
Not Collective


## Input Parameter

- ***ce   -*** The `PetscConvEst` object



## Output Parameter

- ***alpha -*** The convergence rate for each field



## Options Database Keys

- ***-snes_convergence_estimate -*** Execute convergence estimation inside `SNESSolve()` and print out the rate
- ***-ts_convergence_estimate -*** Execute convergence estimation inside `TSSolve()` and print out the rate





## Notes
The convergence rate alpha is defined by
```none
|| u_\Delta - u_exact || < C \Delta^alpha
```
where u_\Delta is the discrete solution, and Delta is a measure of the discretization size. We usually use h for the
spatial resolution and \Delta t for the temporal resolution.

We solve a series of problems using increasing resolution (refined meshes or decreased timesteps), calculate an error
based upon the exact solution in the `PetscDS`, and then fit the result to our model above using linear regression.


## See Also
 `PetscConvEstSetSolver()`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`, `SNESSolve()`, `TSSolve()`

## Level
intermediate

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/utils/convest.c.html#PetscConvEstGetConvRate">src/snes/utils/convest.c</A>


---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/snes/utils/convest.c)


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