# TSSetCostGradients
Sets the initial value of the gradients of the cost function w.r.t. initial values and w.r.t. the problem parameters for use by the `TS` adjoint routines. 
## Synopsis
```
#include "petscts.h"  
PetscErrorCode TSSetCostGradients(TS ts, PetscInt numcost, Vec *lambda, Vec *mu)
```
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
- ***lambda -*** gradients with respect to the initial condition variables, the dimension and parallel layout of these vectors is the same as the ODE solution vector
- ***mu -*** gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters





## Notes
the entries in these vectors must be correctly initialized with the values lambda_i = df/dy|finaltime  mu_i = df/dp|finaltime

After `TSAdjointSolve()` is called the lambda and the mu contain the computed sensitivities


## See Also
 `TS`, `TSAdjointSolve()`, `TSGetCostGradients()`

## Level
beginner

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/interface/sensitivity/tssen.c.html#TSSetCostGradients">src/ts/interface/sensitivity/tssen.c</A>

## Examples
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/ex20adj.c.html">src/ts/tutorials/ex20adj.c.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/ex20opt_ic.c.html">src/ts/tutorials/ex20opt_ic.c.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/ex20opt_p.c.html">src/ts/tutorials/ex20opt_p.c.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/ex20td.c.html">src/ts/tutorials/ex20td.c.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/ex23fwdadj.c.html">src/ts/tutorials/ex23fwdadj.c.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/extchem.c.html">src/ts/tutorials/extchem.c.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/ts/tutorials/extchemfield.c.html">src/ts/tutorials/extchemfield.c.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/tao/unconstrained/tutorials/burgers_spectral.c.html">src/tao/unconstrained/tutorials/burgers_spectral.c.html</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/tao/unconstrained/tutorials/spectraladjointassimilation.c.html">src/tao/unconstrained/tutorials/spectraladjointassimilation.c.html</A><BR>


---
[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](/docs/manualpages/index.md)  
[Index of all manual pages](/docs/manualpages/singleindex.md)  
