:orphan: # TSDiscGradSetFormulation Set the construction method for S, F, and grad F from the formulation u_t = S(u) grad F(u) for `TSDISCGRAD` ## Synopsis ``` #include "petscts.h" PetscErrorCode TSDiscGradSetFormulation(TS ts, PetscErrorCode (*Sfunc)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (*Ffunc)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (*Gfunc)(TS, PetscReal, Vec, Vec, void *), void *ctx) ``` Not Collective ## Input Parameters - ***ts -*** timestepping context - ***Sfunc -*** constructor for the S matrix from the formulation - ***Ffunc -*** functional F from the formulation - ***Gfunc -*** constructor for the gradient of F from the formulation - ***ctx -*** optional context for the functions ## Calling sequence of `Sfunc` ```none PetscErrorCode Sfunc(TS ts, PetscReal time, Vec u, Mat S, void *ctx) ``` ## Calling sequence of `Ffunc` ```none PetscErrorCode Ffunc(TS ts, PetscReal time, Vec u, PetscScalar *F, void *ctx) ``` ## Calling sequence of `Gfunc` ```none PetscErrorCode Gfunc(TS ts, PetscReal time, Vec u, Vec G, void *ctx) ``` ## See Also [](ch_ts), `TSDISCGRAD`, `TSDiscGradGetFormulation()` ## Level Intermediate ## Location src/ts/impls/implicit/discgrad/tsdiscgrad.c ## Implementations TSDiscGradSetFormulation_DiscGrad(TS ts, PetscErrorCode (*Sfunc)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode (*Ffunc)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode (*Gfunc) in src/ts/impls/implicit/discgrad/tsdiscgrad.c
--- [Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/ts/impls/implicit/discgrad/tsdiscgrad.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)