TSSetCostIntegrand#

Sets the routine for evaluating the integral term in one or more cost functions

Synopsis#

#include <petscts.h>  
PetscErrorCode TSSetCostIntegrand(TS ts, PetscInt numcost, Vec costintegral, PetscErrorCode (*rf)(TS ts, PetscReal t, Vec U, Vec F, void *ctx), PetscErrorCode (*drduf)(TS ts, PetscReal t, Vec U, Vec *dRdU, void *ctx), PetscErrorCode (*drdpf)(TS ts, PetscReal t, Vec U, Vec *dRdP, void *ctx), PetscBool fwd, void *ctx)

Logically Collective

Input Parameters#

  • ts - the TS context obtained from TSCreate()

  • numcost - number of gradients to be computed, this is the number of cost functions

  • costintegral - vector that stores the integral values

  • rf - routine for evaluating the integrand function

  • drduf - function that computes the gradients of the r with respect to u

  • drdpf - function that computes the gradients of the r with respect to p, can be NULL if parametric sensitivity is not desired (mu = NULL)

  • fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run

  • ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)

Calling sequence of rf#

  • ts - the integrator

  • t - the time

  • U - the solution

  • F - the computed value of the function

  • ctx - the user context

Calling sequence of drduf#

  • ts - the integrator

  • t - the time

  • U - the solution

  • dRdU - the computed gradients of the r with respect to u

  • ctx - the user context

Calling sequence of drdpf#

  • ts - the integrator

  • t - the time

  • U - the solution

  • dRdP - the computed gradients of the r with respect to p

  • ctx - the user context

Notes#

For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions

Use TSCreateQuadratureTS() and TSForwardSetSensitivities() instead

See Also#

TS: Scalable ODE and DAE Solvers, TS, TSSetRHSJacobianP(), TSGetCostGradients(), TSSetCostGradients(), TSCreateQuadratureTS(), TSForwardSetSensitivities()

Level#

deprecated

Location#

src/ts/interface/sensitivity/tssen.c


Index of all Sensitivity routines
Table of Contents for all manual pages
Index of all manual pages