xref: /petsc/src/ts/impls/implicit/discgrad/tsdiscgrad.c (revision e0f5bfbec699682fa3e8b8532b1176849ea4e12a)
1 /*
2   Code for timestepping with discrete gradient integrators
3 */
4 #include <petsc/private/tsimpl.h> /*I   "petscts.h"   I*/
5 #include <petscdm.h>
6 
7 PetscBool  DGCite       = PETSC_FALSE;
8 const char DGCitation[] = "@article{Gonzalez1996,\n"
9                           "  title   = {Time integration and discrete Hamiltonian systems},\n"
10                           "  author  = {Oscar Gonzalez},\n"
11                           "  journal = {Journal of Nonlinear Science},\n"
12                           "  volume  = {6},\n"
13                           "  pages   = {449--467},\n"
14                           "  doi     = {10.1007/978-1-4612-1246-1_10},\n"
15                           "  year    = {1996}\n}\n";
16 
17 typedef struct {
18   PetscReal stage_time;
19   Vec       X0, X, Xdot;
20   void     *funcCtx;
21   PetscBool gonzalez;
22   PetscErrorCode (*Sfunc)(TS, PetscReal, Vec, Mat, void *);
23   PetscErrorCode (*Ffunc)(TS, PetscReal, Vec, PetscScalar *, void *);
24   PetscErrorCode (*Gfunc)(TS, PetscReal, Vec, Vec, void *);
25 } TS_DiscGrad;
26 
27 static PetscErrorCode TSDiscGradGetX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) {
28   TS_DiscGrad *dg = (TS_DiscGrad *)ts->data;
29 
30   PetscFunctionBegin;
31   if (X0) {
32     if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSDiscGrad_X0", X0));
33     else *X0 = ts->vec_sol;
34   }
35   if (Xdot) {
36     if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSDiscGrad_Xdot", Xdot));
37     else *Xdot = dg->Xdot;
38   }
39   PetscFunctionReturn(0);
40 }
41 
42 static PetscErrorCode TSDiscGradRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot) {
43   PetscFunctionBegin;
44   if (X0) {
45     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSDiscGrad_X0", X0));
46   }
47   if (Xdot) {
48     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSDiscGrad_Xdot", Xdot));
49   }
50   PetscFunctionReturn(0);
51 }
52 
53 static PetscErrorCode DMCoarsenHook_TSDiscGrad(DM fine, DM coarse, void *ctx) {
54   PetscFunctionBegin;
55   PetscFunctionReturn(0);
56 }
57 
58 static PetscErrorCode DMRestrictHook_TSDiscGrad(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) {
59   TS  ts = (TS)ctx;
60   Vec X0, Xdot, X0_c, Xdot_c;
61 
62   PetscFunctionBegin;
63   PetscCall(TSDiscGradGetX0AndXdot(ts, fine, &X0, &Xdot));
64   PetscCall(TSDiscGradGetX0AndXdot(ts, coarse, &X0_c, &Xdot_c));
65   PetscCall(MatRestrict(restrct, X0, X0_c));
66   PetscCall(MatRestrict(restrct, Xdot, Xdot_c));
67   PetscCall(VecPointwiseMult(X0_c, rscale, X0_c));
68   PetscCall(VecPointwiseMult(Xdot_c, rscale, Xdot_c));
69   PetscCall(TSDiscGradRestoreX0AndXdot(ts, fine, &X0, &Xdot));
70   PetscCall(TSDiscGradRestoreX0AndXdot(ts, coarse, &X0_c, &Xdot_c));
71   PetscFunctionReturn(0);
72 }
73 
74 static PetscErrorCode DMSubDomainHook_TSDiscGrad(DM dm, DM subdm, void *ctx) {
75   PetscFunctionBegin;
76   PetscFunctionReturn(0);
77 }
78 
79 static PetscErrorCode DMSubDomainRestrictHook_TSDiscGrad(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) {
80   TS  ts = (TS)ctx;
81   Vec X0, Xdot, X0_sub, Xdot_sub;
82 
83   PetscFunctionBegin;
84   PetscCall(TSDiscGradGetX0AndXdot(ts, dm, &X0, &Xdot));
85   PetscCall(TSDiscGradGetX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub));
86 
87   PetscCall(VecScatterBegin(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD));
88   PetscCall(VecScatterEnd(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD));
89 
90   PetscCall(VecScatterBegin(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD));
91   PetscCall(VecScatterEnd(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD));
92 
93   PetscCall(TSDiscGradRestoreX0AndXdot(ts, dm, &X0, &Xdot));
94   PetscCall(TSDiscGradRestoreX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub));
95   PetscFunctionReturn(0);
96 }
97 
98 static PetscErrorCode TSSetUp_DiscGrad(TS ts) {
99   TS_DiscGrad *dg = (TS_DiscGrad *)ts->data;
100   DM           dm;
101 
102   PetscFunctionBegin;
103   if (!dg->X) PetscCall(VecDuplicate(ts->vec_sol, &dg->X));
104   if (!dg->X0) PetscCall(VecDuplicate(ts->vec_sol, &dg->X0));
105   if (!dg->Xdot) PetscCall(VecDuplicate(ts->vec_sol, &dg->Xdot));
106 
107   PetscCall(TSGetDM(ts, &dm));
108   PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSDiscGrad, DMRestrictHook_TSDiscGrad, ts));
109   PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSDiscGrad, DMSubDomainRestrictHook_TSDiscGrad, ts));
110   PetscFunctionReturn(0);
111 }
112 
113 static PetscErrorCode TSSetFromOptions_DiscGrad(TS ts, PetscOptionItems *PetscOptionsObject) {
114   TS_DiscGrad *dg = (TS_DiscGrad *)ts->data;
115 
116   PetscFunctionBegin;
117   PetscOptionsHeadBegin(PetscOptionsObject, "Discrete Gradients ODE solver options");
118   { PetscCall(PetscOptionsBool("-ts_discgrad_gonzalez", "Use Gonzalez term in discrete gradients formulation", "TSDiscGradUseGonzalez", dg->gonzalez, &dg->gonzalez, NULL)); }
119   PetscOptionsHeadEnd();
120   PetscFunctionReturn(0);
121 }
122 
123 static PetscErrorCode TSView_DiscGrad(TS ts, PetscViewer viewer) {
124   PetscBool iascii;
125 
126   PetscFunctionBegin;
127   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
128   if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, "  Discrete Gradients\n"));
129   PetscFunctionReturn(0);
130 }
131 
132 static PetscErrorCode TSDiscGradIsGonzalez_DiscGrad(TS ts, PetscBool *gonzalez) {
133   TS_DiscGrad *dg = (TS_DiscGrad *)ts->data;
134 
135   PetscFunctionBegin;
136   *gonzalez = dg->gonzalez;
137   PetscFunctionReturn(0);
138 }
139 
140 static PetscErrorCode TSDiscGradUseGonzalez_DiscGrad(TS ts, PetscBool flg) {
141   TS_DiscGrad *dg = (TS_DiscGrad *)ts->data;
142 
143   PetscFunctionBegin;
144   dg->gonzalez = flg;
145   PetscFunctionReturn(0);
146 }
147 
148 static PetscErrorCode TSReset_DiscGrad(TS ts) {
149   TS_DiscGrad *dg = (TS_DiscGrad *)ts->data;
150 
151   PetscFunctionBegin;
152   PetscCall(VecDestroy(&dg->X));
153   PetscCall(VecDestroy(&dg->X0));
154   PetscCall(VecDestroy(&dg->Xdot));
155   PetscFunctionReturn(0);
156 }
157 
158 static PetscErrorCode TSDestroy_DiscGrad(TS ts) {
159   DM dm;
160 
161   PetscFunctionBegin;
162   PetscCall(TSReset_DiscGrad(ts));
163   PetscCall(TSGetDM(ts, &dm));
164   if (dm) {
165     PetscCall(DMCoarsenHookRemove(dm, DMCoarsenHook_TSDiscGrad, DMRestrictHook_TSDiscGrad, ts));
166     PetscCall(DMSubDomainHookRemove(dm, DMSubDomainHook_TSDiscGrad, DMSubDomainRestrictHook_TSDiscGrad, ts));
167   }
168   PetscCall(PetscFree(ts->data));
169   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradGetFormulation_C", NULL));
170   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradSetFormulation_C", NULL));
171   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradIsGonzalez_C", NULL));
172   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradUseGonzalez_C", NULL));
173   PetscFunctionReturn(0);
174 }
175 
176 static PetscErrorCode TSInterpolate_DiscGrad(TS ts, PetscReal t, Vec X) {
177   TS_DiscGrad *dg = (TS_DiscGrad *)ts->data;
178   PetscReal    dt = t - ts->ptime;
179 
180   PetscFunctionBegin;
181   PetscCall(VecCopy(ts->vec_sol, dg->X));
182   PetscCall(VecWAXPY(X, dt, dg->Xdot, dg->X));
183   PetscFunctionReturn(0);
184 }
185 
186 static PetscErrorCode TSDiscGrad_SNESSolve(TS ts, Vec b, Vec x) {
187   SNES     snes;
188   PetscInt nits, lits;
189 
190   PetscFunctionBegin;
191   PetscCall(TSGetSNES(ts, &snes));
192   PetscCall(SNESSolve(snes, b, x));
193   PetscCall(SNESGetIterationNumber(snes, &nits));
194   PetscCall(SNESGetLinearSolveIterations(snes, &lits));
195   ts->snes_its += nits;
196   ts->ksp_its += lits;
197   PetscFunctionReturn(0);
198 }
199 
200 static PetscErrorCode TSStep_DiscGrad(TS ts) {
201   TS_DiscGrad *dg = (TS_DiscGrad *)ts->data;
202   TSAdapt      adapt;
203   TSStepStatus status     = TS_STEP_INCOMPLETE;
204   PetscInt     rejections = 0;
205   PetscBool    stageok, accept = PETSC_TRUE;
206   PetscReal    next_time_step = ts->time_step;
207 
208   PetscFunctionBegin;
209   PetscCall(TSGetAdapt(ts, &adapt));
210   if (!ts->steprollback) PetscCall(VecCopy(ts->vec_sol, dg->X0));
211 
212   while (!ts->reason && status != TS_STEP_COMPLETE) {
213     PetscReal shift = 1 / (0.5 * ts->time_step);
214 
215     dg->stage_time = ts->ptime + 0.5 * ts->time_step;
216 
217     PetscCall(VecCopy(dg->X0, dg->X));
218     PetscCall(TSPreStage(ts, dg->stage_time));
219     PetscCall(TSDiscGrad_SNESSolve(ts, NULL, dg->X));
220     PetscCall(TSPostStage(ts, dg->stage_time, 0, &dg->X));
221     PetscCall(TSAdaptCheckStage(adapt, ts, dg->stage_time, dg->X, &stageok));
222     if (!stageok) goto reject_step;
223 
224     status = TS_STEP_PENDING;
225     PetscCall(VecAXPBYPCZ(dg->Xdot, -shift, shift, 0, dg->X0, dg->X));
226     PetscCall(VecAXPY(ts->vec_sol, ts->time_step, dg->Xdot));
227     PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
228     status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
229     if (!accept) {
230       PetscCall(VecCopy(dg->X0, ts->vec_sol));
231       ts->time_step = next_time_step;
232       goto reject_step;
233     }
234     ts->ptime += ts->time_step;
235     ts->time_step = next_time_step;
236     break;
237 
238   reject_step:
239     ts->reject++;
240     accept = PETSC_FALSE;
241     if (!ts->reason && ts->max_reject >= 0 && ++rejections > ts->max_reject) {
242       ts->reason = TS_DIVERGED_STEP_REJECTED;
243       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
244     }
245   }
246   PetscFunctionReturn(0);
247 }
248 
249 static PetscErrorCode TSGetStages_DiscGrad(TS ts, PetscInt *ns, Vec **Y) {
250   TS_DiscGrad *dg = (TS_DiscGrad *)ts->data;
251 
252   PetscFunctionBegin;
253   if (ns) *ns = 1;
254   if (Y) *Y = &(dg->X);
255   PetscFunctionReturn(0);
256 }
257 
258 /*
259   This defines the nonlinear equation that is to be solved with SNES
260     G(U) = F[t0 + 0.5*dt, U, (U-U0)/dt] = 0
261 */
262 
263 /* x = (x+x')/2 */
264 /* NEED TO CALCULATE x_{n+1} from x and x_{n}*/
265 static PetscErrorCode SNESTSFormFunction_DiscGrad(SNES snes, Vec x, Vec y, TS ts) {
266   TS_DiscGrad *dg = (TS_DiscGrad *)ts->data;
267   PetscReal    norm, shift = 1 / (0.5 * ts->time_step);
268   PetscInt     n;
269   Vec          X0, Xdot, Xp, Xdiff;
270   Mat          S;
271   PetscScalar  F = 0, F0 = 0, Gp;
272   Vec          G, SgF;
273   DM           dm, dmsave;
274 
275   PetscFunctionBegin;
276   PetscCall(SNESGetDM(snes, &dm));
277 
278   PetscCall(VecDuplicate(y, &Xp));
279   PetscCall(VecDuplicate(y, &Xdiff));
280   PetscCall(VecDuplicate(y, &SgF));
281   PetscCall(VecDuplicate(y, &G));
282 
283   PetscCall(VecGetLocalSize(y, &n));
284   PetscCall(MatCreate(PETSC_COMM_WORLD, &S));
285   PetscCall(MatSetSizes(S, PETSC_DECIDE, PETSC_DECIDE, n, n));
286   PetscCall(MatSetFromOptions(S));
287   PetscCall(MatSetUp(S));
288   PetscCall(MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY));
289   PetscCall(MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY));
290 
291   PetscCall(TSDiscGradGetX0AndXdot(ts, dm, &X0, &Xdot));
292   PetscCall(VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x)); /* Xdot = shift (x - X0) */
293 
294   PetscCall(VecAXPBYPCZ(Xp, -1, 2, 0, X0, x));     /* Xp = 2*x - X0 + (0)*Xmid */
295   PetscCall(VecAXPBYPCZ(Xdiff, -1, 1, 0, X0, Xp)); /* Xdiff = xp - X0 + (0)*Xdiff */
296 
297   if (dg->gonzalez) {
298     PetscCall((*dg->Sfunc)(ts, dg->stage_time, x, S, dg->funcCtx));
299     PetscCall((*dg->Ffunc)(ts, dg->stage_time, Xp, &F, dg->funcCtx));
300     PetscCall((*dg->Ffunc)(ts, dg->stage_time, X0, &F0, dg->funcCtx));
301     PetscCall((*dg->Gfunc)(ts, dg->stage_time, x, G, dg->funcCtx));
302 
303     /* Adding Extra Gonzalez Term */
304     PetscCall(VecDot(Xdiff, G, &Gp));
305     PetscCall(VecNorm(Xdiff, NORM_2, &norm));
306     if (norm < PETSC_SQRT_MACHINE_EPSILON) {
307       Gp = 0;
308     } else {
309       /* Gp = (1/|xn+1 - xn|^2) * (F(xn+1) - F(xn) - Gp) */
310       Gp = (F - F0 - Gp) / PetscSqr(norm);
311     }
312     PetscCall(VecAXPY(G, Gp, Xdiff));
313     PetscCall(MatMult(S, G, SgF)); /* S*gradF */
314 
315   } else {
316     PetscCall((*dg->Sfunc)(ts, dg->stage_time, x, S, dg->funcCtx));
317     PetscCall((*dg->Gfunc)(ts, dg->stage_time, x, G, dg->funcCtx));
318 
319     PetscCall(MatMult(S, G, SgF)); /* Xdot = S*gradF */
320   }
321   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
322   dmsave = ts->dm;
323   ts->dm = dm;
324   PetscCall(VecAXPBYPCZ(y, 1, -1, 0, Xdot, SgF));
325   ts->dm = dmsave;
326   PetscCall(TSDiscGradRestoreX0AndXdot(ts, dm, &X0, &Xdot));
327 
328   PetscCall(VecDestroy(&Xp));
329   PetscCall(VecDestroy(&Xdiff));
330   PetscCall(VecDestroy(&SgF));
331   PetscCall(VecDestroy(&G));
332   PetscCall(MatDestroy(&S));
333 
334   PetscFunctionReturn(0);
335 }
336 
337 static PetscErrorCode SNESTSFormJacobian_DiscGrad(SNES snes, Vec x, Mat A, Mat B, TS ts) {
338   TS_DiscGrad *dg    = (TS_DiscGrad *)ts->data;
339   PetscReal    shift = 1 / (0.5 * ts->time_step);
340   Vec          Xdot;
341   DM           dm, dmsave;
342 
343   PetscFunctionBegin;
344   PetscCall(SNESGetDM(snes, &dm));
345   /* Xdot has already been computed in SNESTSFormFunction_DiscGrad (SNES guarantees this) */
346   PetscCall(TSDiscGradGetX0AndXdot(ts, dm, NULL, &Xdot));
347 
348   dmsave = ts->dm;
349   ts->dm = dm;
350   PetscCall(TSComputeIJacobian(ts, dg->stage_time, x, Xdot, shift, A, B, PETSC_FALSE));
351   ts->dm = dmsave;
352   PetscCall(TSDiscGradRestoreX0AndXdot(ts, dm, NULL, &Xdot));
353   PetscFunctionReturn(0);
354 }
355 
356 static PetscErrorCode TSDiscGradGetFormulation_DiscGrad(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) {
357   TS_DiscGrad *dg = (TS_DiscGrad *)ts->data;
358 
359   PetscFunctionBegin;
360   *Sfunc = dg->Sfunc;
361   *Ffunc = dg->Ffunc;
362   *Gfunc = dg->Gfunc;
363   PetscFunctionReturn(0);
364 }
365 
366 static PetscErrorCode TSDiscGradSetFormulation_DiscGrad(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) {
367   TS_DiscGrad *dg = (TS_DiscGrad *)ts->data;
368 
369   PetscFunctionBegin;
370   dg->Sfunc   = Sfunc;
371   dg->Ffunc   = Ffunc;
372   dg->Gfunc   = Gfunc;
373   dg->funcCtx = ctx;
374   PetscFunctionReturn(0);
375 }
376 
377 /*MC
378   TSDISCGRAD - ODE solver using the discrete gradients version of the implicit midpoint method
379 
380   Notes:
381   This is the implicit midpoint rule, with an optional term that guarantees the discrete gradient property. This
382   timestepper applies to systems of the form
383 $ u_t = S(u) grad F(u)
384   where S(u) is a linear operator, and F is a functional of u.
385 
386   Level: intermediate
387 
388 .seealso: `TSCreate()`, `TSSetType()`, `TS`, `TSDISCGRAD`, `TSDiscGradSetFormulation()`
389 M*/
390 PETSC_EXTERN PetscErrorCode TSCreate_DiscGrad(TS ts) {
391   TS_DiscGrad *th;
392 
393   PetscFunctionBegin;
394   PetscCall(PetscCitationsRegister(DGCitation, &DGCite));
395   ts->ops->reset          = TSReset_DiscGrad;
396   ts->ops->destroy        = TSDestroy_DiscGrad;
397   ts->ops->view           = TSView_DiscGrad;
398   ts->ops->setfromoptions = TSSetFromOptions_DiscGrad;
399   ts->ops->setup          = TSSetUp_DiscGrad;
400   ts->ops->step           = TSStep_DiscGrad;
401   ts->ops->interpolate    = TSInterpolate_DiscGrad;
402   ts->ops->getstages      = TSGetStages_DiscGrad;
403   ts->ops->snesfunction   = SNESTSFormFunction_DiscGrad;
404   ts->ops->snesjacobian   = SNESTSFormJacobian_DiscGrad;
405   ts->default_adapt_type  = TSADAPTNONE;
406 
407   ts->usessnes = PETSC_TRUE;
408 
409   PetscCall(PetscNew(&th));
410   ts->data = (void *)th;
411 
412   th->gonzalez = PETSC_FALSE;
413 
414   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradGetFormulation_C", TSDiscGradGetFormulation_DiscGrad));
415   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradSetFormulation_C", TSDiscGradSetFormulation_DiscGrad));
416   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradIsGonzalez_C", TSDiscGradIsGonzalez_DiscGrad));
417   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSDiscGradUseGonzalez_C", TSDiscGradUseGonzalez_DiscGrad));
418   PetscFunctionReturn(0);
419 }
420 
421 /*@C
422   TSDiscGradGetFormulation - Get the construction method for S, F, and grad F from the formulation u_t = S grad F
423 
424   Not Collective
425 
426   Input Parameter:
427 . ts - timestepping context
428 
429   Output Parameters:
430 + Sfunc - constructor for the S matrix from the formulation
431 . Ffunc - functional F from the formulation
432 . Gfunc - constructor for the gradient of F from the formulation
433 - ctx   - the user context
434 
435   Calling sequence of Sfunc:
436 $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Mat S, void *)
437 
438   Calling sequence of Ffunc:
439 $ PetscErrorCode func(TS ts, PetscReal time, Vec u, PetscScalar *F, void *)
440 
441   Calling sequence of Gfunc:
442 $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Vec G, void *)
443 
444   Level: intermediate
445 
446 .seealso: `TSDiscGradSetFormulation()`
447 @*/
448 PetscErrorCode TSDiscGradGetFormulation(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) {
449   PetscFunctionBegin;
450   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
451   PetscValidPointer(Sfunc, 2);
452   PetscValidPointer(Ffunc, 3);
453   PetscValidPointer(Gfunc, 4);
454   PetscUseMethod(ts, "TSDiscGradGetFormulation_C", (TS, PetscErrorCode(**Sfunc)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode(**Ffunc)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode(**Gfunc)(TS, PetscReal, Vec, Vec, void *), void *), (ts, Sfunc, Ffunc, Gfunc, ctx));
455   PetscFunctionReturn(0);
456 }
457 
458 /*@C
459   TSDiscGradSetFormulation - Set the construction method for S, F, and grad F from the formulation u_t = S(u) grad F(u)
460 
461   Not Collective
462 
463   Input Parameters:
464 + ts    - timestepping context
465 . Sfunc - constructor for the S matrix from the formulation
466 . Ffunc - functional F from the formulation
467 - Gfunc - constructor for the gradient of F from the formulation
468   Calling sequence of Sfunc:
469 $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Mat S, void *)
470 
471   Calling sequence of Ffunc:
472 $ PetscErrorCode func(TS ts, PetscReal time, Vec u, PetscScalar *F, void *)
473 
474   Calling sequence of Gfunc:
475 $ PetscErrorCode func(TS ts, PetscReal time, Vec u, Vec G, void *)
476 
477   Level: Intermediate
478 
479 .seealso: `TSDiscGradGetFormulation()`
480 @*/
481 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) {
482   PetscFunctionBegin;
483   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
484   PetscValidFunction(Sfunc, 2);
485   PetscValidFunction(Ffunc, 3);
486   PetscValidFunction(Gfunc, 4);
487   PetscTryMethod(ts, "TSDiscGradSetFormulation_C", (TS, PetscErrorCode(*Sfunc)(TS, PetscReal, Vec, Mat, void *), PetscErrorCode(*Ffunc)(TS, PetscReal, Vec, PetscScalar *, void *), PetscErrorCode(*Gfunc)(TS, PetscReal, Vec, Vec, void *), void *), (ts, Sfunc, Ffunc, Gfunc, ctx));
488   PetscFunctionReturn(0);
489 }
490 
491 /*@
492   TSDiscGradIsGonzalez - Checks flag for whether to use additional conservative terms in discrete gradient formulation.
493 
494   Not Collective
495 
496   Input Parameter:
497 .  ts - timestepping context
498 
499   Output Parameter:
500 .  gonzalez - PETSC_TRUE when using the Gonzalez term
501 
502   Level: Advanced
503 
504 .seealso: `TSDiscGradUseGonzalez()`, `TSDISCGRAD`
505 @*/
506 PetscErrorCode TSDiscGradIsGonzalez(TS ts, PetscBool *gonzalez) {
507   PetscFunctionBegin;
508   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
509   PetscValidBoolPointer(gonzalez, 2);
510   PetscUseMethod(ts, "TSDiscGradIsGonzalez_C", (TS, PetscBool *), (ts, gonzalez));
511   PetscFunctionReturn(0);
512 }
513 
514 /*@
515   TSDiscGradUseGonzalez - Sets discrete gradient formulation with or without additional conservative terms.  Without flag, the discrete gradients timestepper is just backwards euler
516 
517   Not Collective
518 
519   Input Parameters:
520 + ts - timestepping context
521 - flg - PETSC_TRUE to use the Gonzalez term
522 
523   Options Database:
524 . -ts_discgrad_gonzalez <flg> - use the Gonzalez term for the discrete gradient formulation
525 
526   Level: Intermediate
527 
528 .seealso: `TSDISCGRAD`
529 @*/
530 PetscErrorCode TSDiscGradUseGonzalez(TS ts, PetscBool flg) {
531   PetscFunctionBegin;
532   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
533   PetscTryMethod(ts, "TSDiscGradUseGonzalez_C", (TS, PetscBool), (ts, flg));
534   PetscFunctionReturn(0);
535 }
536