:orphan:
# TSSetIFunction
Set the function to compute F(t,U,U_t) where F() = 0 is the DAE to be solved.
## Synopsis
```
#include "petscts.h"
PetscErrorCode TSSetIFunction(TS ts, Vec r, TSIFunction f, void *ctx)
```
Logically Collective
## Input Parameters
- ***ts -*** the `TS` context obtained from `TSCreate()`
- ***r -*** vector to hold the residual (or `NULL` to have it created internally)
- ***f -*** the function evaluation routine
- ***ctx -*** user-defined context for private data for the function evaluation routine (may be `NULL`)
## Calling sequence of `f`
```none
PetscErrorCode f(TS ts, PetscReal t, Vec u, Vec u_t, Vec F, void *ctx)
```
- ***ts -*** the `TS` context obtained from `TSCreate()`
- ***t -*** time at step/stage being solved
- ***u -*** state vector
- ***u_t -*** time derivative of state vector
- ***F -*** function vector
- ***ctx -*** [optional] user-defined context for matrix evaluation routine
## Note
The user MUST call either this routine or `TSSetRHSFunction()` to define the ODE. When solving DAEs you must use this function.
## See Also
[](ch_ts), `TS`, `TSSetRHSJacobian()`, `TSSetRHSFunction()`, `TSSetIJacobian()`
## Level
beginner
## Location
src/ts/interface/ts.c
## Examples
src/ts/tutorials/ex10.c
src/ts/tutorials/ex14.c
src/ts/tutorials/ex15.c
src/ts/tutorials/ex16.c
src/ts/tutorials/ex17.c
src/ts/tutorials/ex19.c
src/ts/tutorials/ex20.c
src/ts/tutorials/ex20adj.c
src/ts/tutorials/ex20fwd.c
src/ts/tutorials/ex20opt_ic.c
src/ts/tutorials/ex20opt_p.c
---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ts/interface/ts.c)
[Index of all TS routines](index.md)
[Table of Contents for all manual pages](/manualpages/index.md)
[Index of all manual pages](/manualpages/singleindex.md)