:orphan: # TSAdjointMonitorSet Sets an ADDITIONAL function that is to be used at every timestep to display the iteration's progress. ## Synopsis ``` #include "petscts.h" PetscErrorCode TSAdjointMonitorSet(TS ts, PetscErrorCode (*adjointmonitor)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *), void *adjointmctx, PetscErrorCode (*adjointmdestroy)(void **)) ``` Logically Collective ## Input Parameters - ***ts -*** the `TS` context obtained from `TSCreate()` - ***adjointmonitor -*** monitoring routine - ***adjointmctx -*** [optional] user-defined context for private data for the monitor routine (use `NULL` if no context is desired) - ***adjointmonitordestroy -*** [optional] routine that frees monitor context (may be `NULL`) ## Calling sequence of `adjointmonitor` ```none PetscErrorCode adjointmonitor(TS ts, PetscInt steps, PetscReal time, Vec u, PetscInt numcost, Vec *lambda, Vec *mu, void *adjointmctx) ``` - ***ts -*** the `TS` context - ***steps -*** iteration number (after the final time step the monitor routine is called with a step of -1, this is at the final time which may have been interpolated to) - ***time -*** current time - ***u -*** current iterate - ***numcost -*** number of cost functionos - ***lambda -*** sensitivities to initial conditions - ***mu -*** sensitivities to parameters - ***adjointmctx -*** [optional] adjoint monitoring context ## Note This routine adds an additional monitor to the list of monitors that already has been loaded. ## Fortran Note Only a single monitor function can be set for each `TS` object ## See Also [](ch_ts), `TS`, `TSAdjointSolve()`, `TSAdjointMonitorCancel()` ## Level intermediate ## Location src/ts/interface/sensitivity/tssen.c ## Examples src/ts/tutorials/ex20td.c
--- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ts/interface/sensitivity/tssen.c) [Index of all Sensitivity routines](index.md) [Table of Contents for all manual pages](/manualpages/index.md) [Index of all manual pages](/manualpages/singleindex.md)