:orphan:
# TSSetRHSFunction
Sets the routine for evaluating the function, where U_t = G(t,u).
## Synopsis
```
#include "petscts.h"
PetscErrorCode TSSetRHSFunction(TS ts, Vec r, PetscErrorCode (*f)(TS, PetscReal, Vec, Vec, void *), void *ctx)
```
Logically Collective
## Input Parameters
- ***ts -*** the `TS` context obtained from `TSCreate()`
- ***r -*** vector to put the computed right hand side (or `NULL` to have it created)
- ***f -*** routine for evaluating the right-hand-side function
- ***ctx -*** [optional] 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 F, void *ctx)
```
- ***ts -*** timestep context
- ***t -*** current timestep
- ***u -*** input vector
- ***F -*** function vector
- ***ctx -*** [optional] user-defined function context
## Note
You must call this function or `TSSetIFunction()` to define your ODE. You cannot use this function when solving a DAE.
## See Also
[](ch_ts), `TS`, `TSSetRHSJacobian()`, `TSSetIJacobian()`, `TSSetIFunction()`
## Level
beginner
## Location
src/ts/interface/ts.c
## Examples
src/tao/unconstrained/tutorials/burgers_spectral.c
src/tao/unconstrained/tutorials/spectraladjointassimilation.c
src/ts/tutorials/ex1.c
src/ts/tutorials/ex12.c
src/ts/tutorials/ex13.c
src/ts/tutorials/ex16.c
src/ts/tutorials/ex16fwd.c
src/ts/tutorials/ex1f.F90
src/ts/tutorials/ex2.c
src/ts/tutorials/ex20.c
src/ts/tutorials/ex20adj.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)