:orphan: # SNESGetConvergenceHistory Gets the array used to hold the convergence history. ## Synopsis ``` #include "petscsnes.h" PetscErrorCode SNESGetConvergenceHistory(SNES snes, PetscReal *a[], PetscInt *its[], PetscInt *na) ``` Not Collective ## Input Parameter - ***snes -*** iterative context obtained from `SNESCreate()` ## Output Parameters - ***a -*** array to hold history, usually was set with `SNESSetConvergenceHistory()` - ***its -*** integer array holds the number of linear iterations (or negative if not converged) for each solve. - ***na -*** size of `a` and `its` ## Note This routine is useful, e.g., when running a code for purposes of accurate performance monitoring, when no I/O should be done during the section of code that is being timed. ## Fortran Note The calling sequence for this routine in Fortran is ```none call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) ``` ## See Also [](ch_snes), `SNES`, `SNESSolve()`, `SNESSetConvergenceHistory()` ## Level intermediate ## Location src/snes/interface/snes.c --- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/snes/interface/snes.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)