:orphan: # TaoSetConvergenceHistory Sets the array used to hold the convergence history. ## Synopsis ``` #include "petsctao.h" PetscErrorCode TaoSetConvergenceHistory(Tao tao, PetscReal obj[], PetscReal resid[], PetscReal cnorm[], PetscInt lits[], PetscInt na, PetscBool reset) ``` Logically Collective ## Input Parameters - ***tao -*** the `Tao` solver context - ***obj -*** array to hold objective value history - ***resid -*** array to hold residual history - ***cnorm -*** array to hold constraint violation history - ***lits -*** integer array holds the number of linear iterations for each Tao iteration - ***na -*** size of `obj`, `resid`, and `cnorm` - ***reset -*** `PETSC_TRUE` indicates each new minimization resets the history counter to zero, else it continues storing new values for new minimizations after the old ones ## Notes If set, `Tao` will fill the given arrays with the indicated information at each iteration. If 'obj','resid','cnorm','lits' are *all* `NULL` then space (using size `na`, or 1000 if na is `PETSC_DECIDE` or `PETSC_DEFAULT`) is allocated for the history. If not all are `NULL`, then only the non-`NULL` information categories will be stored, the others will be ignored. Any convergence information after iteration number 'na' will not be stored. 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. ## See Also [](ch_tao), `TaoGetConvergenceHistory()` ## Level intermediate ## Location src/tao/interface/taosolver.c ## Examples src/tao/leastsquares/tutorials/chwirut1.c
src/tao/leastsquares/tutorials/chwirut1f.F90
src/tao/leastsquares/tutorials/cs1.c
src/tao/leastsquares/tutorials/tomography.c
--- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/tao/interface/taosolver.c) [Index of all Tao routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)