#include "petscts.h"
PetscErrorCode TSSetCostIntegrand(TS ts,PetscInt numcost, PetscErrorCode (*rf)(TS,PetscReal,Vec,Vec,void*),
PetscErrorCode (*drdyf)(TS,PetscReal,Vec,Vec*,void*),
PetscErrorCode (*drdpf)(TS,PetscReal,Vec,Vec*,void*),void *ctx)
Logically Collective on TS
| ts | - the TS context obtained from TSCreate() | |
| numcost | - number of gradients to be computed, this is the number of cost functions | |
| rf | - routine for evaluating the integrand function | |
| drdyf | - function that computes the gradients of the r's with respect to y,NULL if not a function y | |
| drdpf | - function that computes the gradients of the r's with respect to p, NULL if not a function of p | |
| ctx | - [optional] user-defined context for private data for the function evaluation routine (may be NULL) |
rf(TS ts,PetscReal t,Vec y,Vec f[],void *ctx);
| t | - current timestep | |
| y | - input vector | |
| f | - function result; one vector entry for each cost function | |
| ctx | - [optional] user-defined function context |
PetscErroCode drdyf(TS ts,PetscReal t,Vec y,Vec *drdy,void *ctx);
PetscErroCode drdpf(TS ts,PetscReal t,Vec y,Vec *drdp,void *ctx);
Notes: For optimization there is generally a single cost function, numcost = 1. For sensitivities there may be multiple cost functions
Level:intermediate
Location:src/ts/interface/ts.c
Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages