petsc-3.6.4 2016-04-12
Report Typos and Errors

TaoGetConvergenceHistory

Gets the arrays used to hold the convergence history.

Synopsis

#include "petsctao.h" 
PetscErrorCode TaoGetConvergenceHistory(Tao tao, PetscReal **obj, PetscReal **resid, PetscReal **cnorm, PetscInt **lits, PetscInt *nhist)
Collective on Tao

Input Parameter

tao -the Tao context

Output Parameters

obj - array used to hold objective value history
resid - array used to hold residual history
cnorm - array used to hold constraint violation history
lits - integer array used to hold linear solver iteration count
nhist - size of obj, resid, cnorm, and lits (will be less than or equal to na given in TaoSetHistory)

Notes

This routine must be preceded by calls to TaoSetConvergenceHistory() and TaoSolve(), otherwise it returns useless information.

The calling sequence for this routine in Fortran is

  call TaoGetConvergenceHistory(Tao tao, PetscInt nhist, PetscErrorCode ierr)

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

TaoSetConvergenceHistory()

Level:advanced
Location:
src/tao/interface/taosolver.c
Index of all Tao routines
Table of Contents for all manual pages
Index of all manual pages

Examples

src/tao/leastsquares/examples/tutorials/chwirut1.c.html
src/tao/leastsquares/examples/tutorials/chwirut1f.F.html