petsc-3.9.4 2018-09-11
Report Typos and Errors

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,PetscReal,Vec,Vec,void*),
                                                          PetscErrorCode (*drdyf)(TS,PetscReal,Vec,Vec*,void*),
                                                          PetscErrorCode (*drdpf)(TS,PetscReal,Vec,Vec*,void*),
                                                          PetscBool fwd,void *ctx)
Logically Collective on TS

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
drdyf - function that computes the gradients of the r's with respect to y
drdpf - function that computes the gradients of the r's 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

  PetscErrorCode rf(TS ts,PetscReal t,Vec y,Vec f,void *ctx);

Calling sequence of drdyf

  PetscErroCode drdyf(TS ts,PetscReal t,Vec y,Vec *drdy,void *ctx);

Calling sequence of drdpf

  PetscErroCode drdpf(TS ts,PetscReal t,Vec y,Vec *drdp,void *ctx);

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

Keywords

TS, sensitivity analysis, timestep, set, quadrature, function

See Also

TSSetRHSJacobianP(), TSGetCostGradients(), TSSetCostGradients()

Level

intermediate

Location

src/ts/interface/tssen.c
Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages