:orphan:
# PCApplyRichardson
Applies several steps of Richardson iteration with the particular preconditioner. This routine is usually used by the Krylov solvers and not the application code directly.
## Synopsis
```
#include "petscksp.h"
PetscErrorCode PCApplyRichardson(PC pc, Vec b, Vec y, Vec w, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt its, PetscBool guesszero, PetscInt *outits, PCRichardsonConvergedReason *reason)
```
Collective
## Input Parameters
- ***pc -*** the preconditioner context
- ***b -*** the right hand side
- ***w -*** one work vector
- ***rtol -*** relative decrease in residual norm convergence criteria
- ***abstol -*** absolute residual norm convergence criteria
- ***dtol -*** divergence residual norm increase criteria
- ***its -*** the number of iterations to apply.
- ***guesszero -*** if the input x contains nonzero initial guess
## Output Parameters
- ***outits -*** number of iterations actually used (for SOR this always equals its)
- ***reason -*** the reason the apply terminated
- ***y -*** the solution (also contains initial guess if guesszero is `PETSC_FALSE`
## Notes
Most preconditioners do not support this function. Use the command
`PCApplyRichardsonExists()` to determine if one does.
Except for the `PCMG` this routine ignores the convergence tolerances
and always runs for the number of iterations
## See Also
`PC`, `PCApplyRichardsonExists()`
## Level
developer
## Location
src/ksp/pc/interface/precon.c
## Implementations
PCApplyRichardson_PFMG in src/ksp/pc/impls/hypre/hypre.c
PCApplyRichardson_SysPFMG in src/ksp/pc/impls/hypre/hypre.c
PCApplyRichardson_SMG in src/ksp/pc/impls/hypre/hypre.c
PCApplyRichardson_MG in src/ksp/pc/impls/mg/mg.c
PCApplyRichardson_Shell in src/ksp/pc/impls/shell/shellpc.c
PCApplyRichardson_SOR in src/ksp/pc/impls/sor/sor.c
PCApplyRichardson_Telescope in src/ksp/pc/impls/telescope/telescope.c
---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ksp/pc/interface/precon.c)
[Index of all PC routines](index.md)
[Table of Contents for all manual pages](/manualpages/index.md)
[Index of all manual pages](/manualpages/singleindex.md)