TSAdjointMonitor#
Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet()
Synopsis#
#include "petscts.h"
PetscErrorCode TSAdjointMonitor(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscInt numcost, Vec *lambda, Vec *mu)
Collective
Input Parameters#
ts - time stepping context obtained from
TSCreate()
step - step number that has just completed
ptime - model time of the state
u - state at the current model time
numcost - number of cost functions (dimension of lambda or mu)
lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
mu - vectors containing the gradients of the cost functions with respect to the problem parameters
Note#
TSAdjointMonitor()
is typically used automatically within the time stepping implementations.
Users would almost never call this routine directly.
See Also#
Level#
developer
Location#
src/ts/interface/sensitivity/tssen.c
Index of all Sensitivity routines
Table of Contents for all manual pages
Index of all manual pages