petsc-3.14.6 2021-03-30
Report Typos and Errors

TSDiscGradSetFormulation

Set the construction method for S, F, and grad F from the formulation u_t = S(u) grad F(u)

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

Calling sequence of Sfunc

PetscErrorCode func(TS ts, PetscReal time, Vec u, Mat S, void *)

Calling sequence of Ffunc

PetscErrorCode func(TS ts, PetscReal time, Vec u, PetscScalar *F, void *)

Calling sequence of Gfunc

PetscErrorCode func(TS ts, PetscReal time, Vec u, Vec G, void *)

See Also

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

Index of all TS routines
Table of Contents for all manual pages
Index of all manual pages