xref: /petsc/src/ts/impls/bdf/bdf.c (revision 21e3ffae2f3b73c0bd738cf6d0a809700fc04bb0)
1 /*
2   Code for timestepping with BDF methods
3 */
4 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
5 #include <petscdm.h>
6 
7 static PetscBool  cited      = PETSC_FALSE;
8 static const char citation[] = "@book{Brenan1995,\n"
9                                "  title     = {Numerical Solution of Initial-Value Problems in Differential-Algebraic Equations},\n"
10                                "  author    = {Brenan, K. and Campbell, S. and Petzold, L.},\n"
11                                "  publisher = {Society for Industrial and Applied Mathematics},\n"
12                                "  year      = {1995},\n"
13                                "  doi       = {10.1137/1.9781611971224},\n}\n";
14 
15 typedef struct {
16   PetscInt  k, n;
17   PetscReal time[6 + 2];
18   Vec       work[6 + 2];
19   Vec       tvwork[6 + 2];
20   PetscReal shift;
21   Vec       vec_dot; /* Xdot when !transientvar, else Cdot where C(X) is the transient variable. */
22   Vec       vec_wrk;
23   Vec       vec_lte;
24 
25   PetscBool    transientvar;
26   PetscInt     order;
27   TSStepStatus status;
28 } TS_BDF;
29 
30 /* Compute Lagrange polynomials on T[:n] evaluated at t.
31  * If one has data (T[i], Y[i]), then the interpolation/extrapolation is f(t) = \sum_i L[i]*Y[i].
32  */
33 static inline void LagrangeBasisVals(PetscInt n, PetscReal t, const PetscReal T[], PetscScalar L[])
34 {
35   PetscInt k, j;
36   for (k = 0; k < n; k++)
37     for (L[k] = 1, j = 0; j < n; j++)
38       if (j != k) L[k] *= (t - T[j]) / (T[k] - T[j]);
39 }
40 
41 static inline void LagrangeBasisDers(PetscInt n, PetscReal t, const PetscReal T[], PetscScalar dL[])
42 {
43   PetscInt k, j, i;
44   for (k = 0; k < n; k++)
45     for (dL[k] = 0, j = 0; j < n; j++)
46       if (j != k) {
47         PetscReal L = 1 / (T[k] - T[j]);
48         for (i = 0; i < n; i++)
49           if (i != j && i != k) L *= (t - T[i]) / (T[k] - T[i]);
50         dL[k] += L;
51       }
52 }
53 
54 static PetscErrorCode TSBDF_GetVecs(TS ts, DM dm, Vec *Xdot, Vec *Ydot)
55 {
56   TS_BDF *bdf = (TS_BDF *)ts->data;
57 
58   PetscFunctionBegin;
59   if (dm && dm != ts->dm) {
60     PetscCall(DMGetNamedGlobalVector(dm, "TSBDF_Vec_Xdot", Xdot));
61     PetscCall(DMGetNamedGlobalVector(dm, "TSBDF_Vec_Ydot", Ydot));
62   } else {
63     *Xdot = bdf->vec_dot;
64     *Ydot = bdf->vec_wrk;
65   }
66   PetscFunctionReturn(PETSC_SUCCESS);
67 }
68 
69 static PetscErrorCode TSBDF_RestoreVecs(TS ts, DM dm, Vec *Xdot, Vec *Ydot)
70 {
71   TS_BDF *bdf = (TS_BDF *)ts->data;
72 
73   PetscFunctionBegin;
74   if (dm && dm != ts->dm) {
75     PetscCall(DMRestoreNamedGlobalVector(dm, "TSBDF_Vec_Xdot", Xdot));
76     PetscCall(DMRestoreNamedGlobalVector(dm, "TSBDF_Vec_Ydot", Ydot));
77   } else {
78     PetscCheck(*Xdot == bdf->vec_dot, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_INCOMP, "Vec does not match the cache");
79     PetscCheck(*Ydot == bdf->vec_wrk, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_INCOMP, "Vec does not match the cache");
80     *Xdot = NULL;
81     *Ydot = NULL;
82   }
83   PetscFunctionReturn(PETSC_SUCCESS);
84 }
85 
86 static PetscErrorCode DMCoarsenHook_TSBDF(DM fine, DM coarse, void *ctx)
87 {
88   PetscFunctionBegin;
89   PetscFunctionReturn(PETSC_SUCCESS);
90 }
91 
92 static PetscErrorCode DMRestrictHook_TSBDF(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx)
93 {
94   TS  ts = (TS)ctx;
95   Vec Ydot, Ydot_c;
96   Vec Xdot, Xdot_c;
97 
98   PetscFunctionBegin;
99   PetscCall(TSBDF_GetVecs(ts, fine, &Xdot, &Ydot));
100   PetscCall(TSBDF_GetVecs(ts, coarse, &Xdot_c, &Ydot_c));
101 
102   PetscCall(MatRestrict(restrct, Ydot, Ydot_c));
103   PetscCall(VecPointwiseMult(Ydot_c, rscale, Ydot_c));
104 
105   PetscCall(TSBDF_RestoreVecs(ts, fine, &Xdot, &Ydot));
106   PetscCall(TSBDF_RestoreVecs(ts, coarse, &Xdot_c, &Ydot_c));
107   PetscFunctionReturn(PETSC_SUCCESS);
108 }
109 
110 static PetscErrorCode TSBDF_Advance(TS ts, PetscReal t, Vec X)
111 {
112   TS_BDF  *bdf = (TS_BDF *)ts->data;
113   PetscInt i, n = (PetscInt)(sizeof(bdf->work) / sizeof(Vec));
114   Vec      tail = bdf->work[n - 1], tvtail = bdf->tvwork[n - 1];
115 
116   PetscFunctionBegin;
117   for (i = n - 1; i >= 2; i--) {
118     bdf->time[i]   = bdf->time[i - 1];
119     bdf->work[i]   = bdf->work[i - 1];
120     bdf->tvwork[i] = bdf->tvwork[i - 1];
121   }
122   bdf->n         = PetscMin(bdf->n + 1, n - 1);
123   bdf->time[1]   = t;
124   bdf->work[1]   = tail;
125   bdf->tvwork[1] = tvtail;
126   PetscCall(VecCopy(X, tail));
127   PetscCall(TSComputeTransientVariable(ts, tail, tvtail));
128   PetscFunctionReturn(PETSC_SUCCESS);
129 }
130 
131 static PetscErrorCode TSBDF_VecLTE(TS ts, PetscInt order, Vec lte)
132 {
133   TS_BDF     *bdf = (TS_BDF *)ts->data;
134   PetscInt    i, n = order + 1;
135   PetscReal  *time = bdf->time;
136   Vec        *vecs = bdf->work;
137   PetscScalar a[8], b[8], alpha[8];
138 
139   PetscFunctionBegin;
140   LagrangeBasisDers(n + 0, time[0], time, a);
141   a[n] = 0;
142   LagrangeBasisDers(n + 1, time[0], time, b);
143   for (i = 0; i < n + 1; i++) alpha[i] = (a[i] - b[i]) / a[0];
144   PetscCall(VecZeroEntries(lte));
145   PetscCall(VecMAXPY(lte, n + 1, alpha, vecs));
146   PetscFunctionReturn(PETSC_SUCCESS);
147 }
148 
149 static PetscErrorCode TSBDF_Extrapolate(TS ts, PetscInt order, PetscReal t, Vec X)
150 {
151   TS_BDF     *bdf  = (TS_BDF *)ts->data;
152   PetscInt    n    = order + 1;
153   PetscReal  *time = bdf->time + 1;
154   Vec        *vecs = bdf->work + 1;
155   PetscScalar alpha[7];
156 
157   PetscFunctionBegin;
158   n = PetscMin(n, bdf->n);
159   LagrangeBasisVals(n, t, time, alpha);
160   PetscCall(VecZeroEntries(X));
161   PetscCall(VecMAXPY(X, n, alpha, vecs));
162   PetscFunctionReturn(PETSC_SUCCESS);
163 }
164 
165 static PetscErrorCode TSBDF_Interpolate(TS ts, PetscInt order, PetscReal t, Vec X)
166 {
167   TS_BDF     *bdf  = (TS_BDF *)ts->data;
168   PetscInt    n    = order + 1;
169   PetscReal  *time = bdf->time;
170   Vec        *vecs = bdf->work;
171   PetscScalar alpha[7];
172 
173   PetscFunctionBegin;
174   LagrangeBasisVals(n, t, time, alpha);
175   PetscCall(VecZeroEntries(X));
176   PetscCall(VecMAXPY(X, n, alpha, vecs));
177   PetscFunctionReturn(PETSC_SUCCESS);
178 }
179 
180 /* Compute the affine term V0 such that Xdot = shift*X + V0.
181  *
182  * When using transient variables, we're computing Cdot = shift*C(X) + V0, and thus choose a linear combination of tvwork.
183  */
184 static PetscErrorCode TSBDF_PreSolve(TS ts)
185 {
186   TS_BDF     *bdf = (TS_BDF *)ts->data;
187   PetscInt    i, n = PetscMax(bdf->k, 1) + 1;
188   Vec         V, V0;
189   Vec         vecs[7];
190   PetscScalar alpha[7];
191 
192   PetscFunctionBegin;
193   PetscCall(TSBDF_GetVecs(ts, NULL, &V, &V0));
194   LagrangeBasisDers(n, bdf->time[0], bdf->time, alpha);
195   for (i = 1; i < n; i++) vecs[i] = bdf->transientvar ? bdf->tvwork[i] : bdf->work[i];
196   PetscCall(VecZeroEntries(V0));
197   PetscCall(VecMAXPY(V0, n - 1, alpha + 1, vecs + 1));
198   bdf->shift = PetscRealPart(alpha[0]);
199   PetscCall(TSBDF_RestoreVecs(ts, NULL, &V, &V0));
200   PetscFunctionReturn(PETSC_SUCCESS);
201 }
202 
203 static PetscErrorCode TSBDF_SNESSolve(TS ts, Vec b, Vec x)
204 {
205   PetscInt nits, lits;
206 
207   PetscFunctionBegin;
208   PetscCall(TSBDF_PreSolve(ts));
209   PetscCall(SNESSolve(ts->snes, b, x));
210   PetscCall(SNESGetIterationNumber(ts->snes, &nits));
211   PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits));
212   ts->snes_its += nits;
213   ts->ksp_its += lits;
214   PetscFunctionReturn(PETSC_SUCCESS);
215 }
216 
217 static PetscErrorCode TSBDF_Restart(TS ts, PetscBool *accept)
218 {
219   TS_BDF *bdf = (TS_BDF *)ts->data;
220 
221   PetscFunctionBegin;
222   bdf->k = 1;
223   bdf->n = 0;
224   PetscCall(TSBDF_Advance(ts, ts->ptime, ts->vec_sol));
225 
226   bdf->time[0] = ts->ptime + ts->time_step / 2;
227   PetscCall(VecCopy(bdf->work[1], bdf->work[0]));
228   PetscCall(TSPreStage(ts, bdf->time[0]));
229   PetscCall(TSBDF_SNESSolve(ts, NULL, bdf->work[0]));
230   PetscCall(TSPostStage(ts, bdf->time[0], 0, &bdf->work[0]));
231   PetscCall(TSAdaptCheckStage(ts->adapt, ts, bdf->time[0], bdf->work[0], accept));
232   if (!*accept) PetscFunctionReturn(PETSC_SUCCESS);
233 
234   bdf->k = PetscMin(2, bdf->order);
235   bdf->n++;
236   PetscCall(VecCopy(bdf->work[0], bdf->work[2]));
237   bdf->time[2] = bdf->time[0];
238   PetscCall(TSComputeTransientVariable(ts, bdf->work[2], bdf->tvwork[2]));
239   PetscFunctionReturn(PETSC_SUCCESS);
240 }
241 
242 static const char *const BDF_SchemeName[] = {"", "1", "2", "3", "4", "5", "6"};
243 
244 static PetscErrorCode TSStep_BDF(TS ts)
245 {
246   TS_BDF   *bdf        = (TS_BDF *)ts->data;
247   PetscInt  rejections = 0;
248   PetscBool stageok, accept = PETSC_TRUE;
249   PetscReal next_time_step = ts->time_step;
250 
251   PetscFunctionBegin;
252   PetscCall(PetscCitationsRegister(citation, &cited));
253 
254   if (!ts->steprollback && !ts->steprestart) {
255     bdf->k = PetscMin(bdf->k + 1, bdf->order);
256     PetscCall(TSBDF_Advance(ts, ts->ptime, ts->vec_sol));
257   }
258 
259   bdf->status = TS_STEP_INCOMPLETE;
260   while (!ts->reason && bdf->status != TS_STEP_COMPLETE) {
261     if (ts->steprestart) {
262       PetscCall(TSBDF_Restart(ts, &stageok));
263       if (!stageok) goto reject_step;
264     }
265 
266     bdf->time[0] = ts->ptime + ts->time_step;
267     PetscCall(TSBDF_Extrapolate(ts, bdf->k - (accept ? 0 : 1), bdf->time[0], bdf->work[0]));
268     PetscCall(TSPreStage(ts, bdf->time[0]));
269     PetscCall(TSBDF_SNESSolve(ts, NULL, bdf->work[0]));
270     PetscCall(TSPostStage(ts, bdf->time[0], 0, &bdf->work[0]));
271     PetscCall(TSAdaptCheckStage(ts->adapt, ts, bdf->time[0], bdf->work[0], &stageok));
272     if (!stageok) goto reject_step;
273 
274     bdf->status = TS_STEP_PENDING;
275     PetscCall(TSAdaptCandidatesClear(ts->adapt));
276     PetscCall(TSAdaptCandidateAdd(ts->adapt, BDF_SchemeName[bdf->k], bdf->k, 1, 1.0, 1.0, PETSC_TRUE));
277     PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
278     bdf->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
279     if (!accept) {
280       ts->time_step = next_time_step;
281       goto reject_step;
282     }
283 
284     PetscCall(VecCopy(bdf->work[0], ts->vec_sol));
285     ts->ptime += ts->time_step;
286     ts->time_step = next_time_step;
287     break;
288 
289   reject_step:
290     ts->reject++;
291     accept = PETSC_FALSE;
292     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
293       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
294       ts->reason = TS_DIVERGED_STEP_REJECTED;
295     }
296   }
297   PetscFunctionReturn(PETSC_SUCCESS);
298 }
299 
300 static PetscErrorCode TSInterpolate_BDF(TS ts, PetscReal t, Vec X)
301 {
302   TS_BDF *bdf = (TS_BDF *)ts->data;
303 
304   PetscFunctionBegin;
305   PetscCall(TSBDF_Interpolate(ts, bdf->k, t, X));
306   PetscFunctionReturn(PETSC_SUCCESS);
307 }
308 
309 static PetscErrorCode TSEvaluateWLTE_BDF(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte)
310 {
311   TS_BDF   *bdf = (TS_BDF *)ts->data;
312   PetscInt  k   = bdf->k;
313   PetscReal wltea, wlter;
314   Vec       X = bdf->work[0], Y = bdf->vec_lte;
315 
316   PetscFunctionBegin;
317   k = PetscMin(k, bdf->n - 1);
318   PetscCall(TSBDF_VecLTE(ts, k, Y));
319   PetscCall(VecAXPY(Y, 1, X));
320   PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, wlte, &wltea, &wlter));
321   if (order) *order = k + 1;
322   PetscFunctionReturn(PETSC_SUCCESS);
323 }
324 
325 static PetscErrorCode TSRollBack_BDF(TS ts)
326 {
327   TS_BDF *bdf = (TS_BDF *)ts->data;
328 
329   PetscFunctionBegin;
330   PetscCall(VecCopy(bdf->work[1], ts->vec_sol));
331   PetscFunctionReturn(PETSC_SUCCESS);
332 }
333 
334 static PetscErrorCode SNESTSFormFunction_BDF(SNES snes, Vec X, Vec F, TS ts)
335 {
336   TS_BDF   *bdf = (TS_BDF *)ts->data;
337   DM        dm, dmsave = ts->dm;
338   PetscReal t     = bdf->time[0];
339   PetscReal shift = bdf->shift;
340   Vec       V, V0;
341 
342   PetscFunctionBegin;
343   PetscCall(SNESGetDM(snes, &dm));
344   PetscCall(TSBDF_GetVecs(ts, dm, &V, &V0));
345   if (bdf->transientvar) { /* shift*C(X) + V0 */
346     PetscCall(TSComputeTransientVariable(ts, X, V));
347     PetscCall(VecAYPX(V, shift, V0));
348   } else { /* shift*X + V0 */
349     PetscCall(VecWAXPY(V, shift, X, V0));
350   }
351 
352   /* F = Function(t,X,V) */
353   ts->dm = dm;
354   PetscCall(TSComputeIFunction(ts, t, X, V, F, PETSC_FALSE));
355   ts->dm = dmsave;
356 
357   PetscCall(TSBDF_RestoreVecs(ts, dm, &V, &V0));
358   PetscFunctionReturn(PETSC_SUCCESS);
359 }
360 
361 static PetscErrorCode SNESTSFormJacobian_BDF(SNES snes, Vec X, Mat J, Mat P, TS ts)
362 {
363   TS_BDF   *bdf = (TS_BDF *)ts->data;
364   DM        dm, dmsave = ts->dm;
365   PetscReal t     = bdf->time[0];
366   PetscReal shift = bdf->shift;
367   Vec       V, V0;
368 
369   PetscFunctionBegin;
370   PetscCall(SNESGetDM(snes, &dm));
371   PetscCall(TSBDF_GetVecs(ts, dm, &V, &V0));
372 
373   /* J,P = Jacobian(t,X,V) */
374   ts->dm = dm;
375   PetscCall(TSComputeIJacobian(ts, t, X, V, shift, J, P, PETSC_FALSE));
376   ts->dm = dmsave;
377 
378   PetscCall(TSBDF_RestoreVecs(ts, dm, &V, &V0));
379   PetscFunctionReturn(PETSC_SUCCESS);
380 }
381 
382 static PetscErrorCode TSReset_BDF(TS ts)
383 {
384   TS_BDF *bdf = (TS_BDF *)ts->data;
385   size_t  i, n = sizeof(bdf->work) / sizeof(Vec);
386 
387   PetscFunctionBegin;
388   bdf->k = bdf->n = 0;
389   for (i = 0; i < n; i++) {
390     PetscCall(VecDestroy(&bdf->work[i]));
391     PetscCall(VecDestroy(&bdf->tvwork[i]));
392   }
393   PetscCall(VecDestroy(&bdf->vec_dot));
394   PetscCall(VecDestroy(&bdf->vec_wrk));
395   PetscCall(VecDestroy(&bdf->vec_lte));
396   if (ts->dm) PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSBDF, DMRestrictHook_TSBDF, ts));
397   PetscFunctionReturn(PETSC_SUCCESS);
398 }
399 
400 static PetscErrorCode TSDestroy_BDF(TS ts)
401 {
402   PetscFunctionBegin;
403   PetscCall(TSReset_BDF(ts));
404   PetscCall(PetscFree(ts->data));
405   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBDFSetOrder_C", NULL));
406   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBDFGetOrder_C", NULL));
407   PetscFunctionReturn(PETSC_SUCCESS);
408 }
409 
410 static PetscErrorCode TSSetUp_BDF(TS ts)
411 {
412   TS_BDF   *bdf = (TS_BDF *)ts->data;
413   size_t    i, n = sizeof(bdf->work) / sizeof(Vec);
414   PetscReal low, high, two = 2;
415 
416   PetscFunctionBegin;
417   PetscCall(TSHasTransientVariable(ts, &bdf->transientvar));
418   bdf->k = bdf->n = 0;
419   for (i = 0; i < n; i++) {
420     PetscCall(VecDuplicate(ts->vec_sol, &bdf->work[i]));
421     if (i && bdf->transientvar) PetscCall(VecDuplicate(ts->vec_sol, &bdf->tvwork[i]));
422   }
423   PetscCall(VecDuplicate(ts->vec_sol, &bdf->vec_dot));
424   PetscCall(VecDuplicate(ts->vec_sol, &bdf->vec_wrk));
425   PetscCall(VecDuplicate(ts->vec_sol, &bdf->vec_lte));
426   PetscCall(TSGetDM(ts, &ts->dm));
427   PetscCall(DMCoarsenHookAdd(ts->dm, DMCoarsenHook_TSBDF, DMRestrictHook_TSBDF, ts));
428 
429   PetscCall(TSGetAdapt(ts, &ts->adapt));
430   PetscCall(TSAdaptCandidatesClear(ts->adapt));
431   PetscCall(TSAdaptGetClip(ts->adapt, &low, &high));
432   PetscCall(TSAdaptSetClip(ts->adapt, low, PetscMin(high, two)));
433 
434   PetscCall(TSGetSNES(ts, &ts->snes));
435   PetscFunctionReturn(PETSC_SUCCESS);
436 }
437 
438 static PetscErrorCode TSSetFromOptions_BDF(TS ts, PetscOptionItems *PetscOptionsObject)
439 {
440   PetscFunctionBegin;
441   PetscOptionsHeadBegin(PetscOptionsObject, "BDF ODE solver options");
442   {
443     PetscBool flg;
444     PetscInt  order;
445     PetscCall(TSBDFGetOrder(ts, &order));
446     PetscCall(PetscOptionsInt("-ts_bdf_order", "Order of the BDF method", "TSBDFSetOrder", order, &order, &flg));
447     if (flg) PetscCall(TSBDFSetOrder(ts, order));
448   }
449   PetscOptionsHeadEnd();
450   PetscFunctionReturn(PETSC_SUCCESS);
451 }
452 
453 static PetscErrorCode TSView_BDF(TS ts, PetscViewer viewer)
454 {
455   TS_BDF   *bdf = (TS_BDF *)ts->data;
456   PetscBool iascii;
457 
458   PetscFunctionBegin;
459   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
460   if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, "  Order=%" PetscInt_FMT "\n", bdf->order));
461   PetscFunctionReturn(PETSC_SUCCESS);
462 }
463 
464 /* ------------------------------------------------------------ */
465 
466 static PetscErrorCode TSBDFSetOrder_BDF(TS ts, PetscInt order)
467 {
468   TS_BDF *bdf = (TS_BDF *)ts->data;
469 
470   PetscFunctionBegin;
471   if (order == bdf->order) PetscFunctionReturn(PETSC_SUCCESS);
472   PetscCheck(order >= 1 && order <= 6, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "BDF Order %" PetscInt_FMT " not implemented", order);
473   bdf->order = order;
474   PetscFunctionReturn(PETSC_SUCCESS);
475 }
476 
477 static PetscErrorCode TSBDFGetOrder_BDF(TS ts, PetscInt *order)
478 {
479   TS_BDF *bdf = (TS_BDF *)ts->data;
480 
481   PetscFunctionBegin;
482   *order = bdf->order;
483   PetscFunctionReturn(PETSC_SUCCESS);
484 }
485 
486 /* ------------------------------------------------------------ */
487 
488 /*MC
489       TSBDF - DAE solver using BDF methods
490 
491   Level: beginner
492 
493 .seealso: [](chapter_ts), `TS`, `TSCreate()`, `TSSetType()`, `TSType`
494 M*/
495 PETSC_EXTERN PetscErrorCode TSCreate_BDF(TS ts)
496 {
497   TS_BDF *bdf;
498 
499   PetscFunctionBegin;
500   ts->ops->reset          = TSReset_BDF;
501   ts->ops->destroy        = TSDestroy_BDF;
502   ts->ops->view           = TSView_BDF;
503   ts->ops->setup          = TSSetUp_BDF;
504   ts->ops->setfromoptions = TSSetFromOptions_BDF;
505   ts->ops->step           = TSStep_BDF;
506   ts->ops->evaluatewlte   = TSEvaluateWLTE_BDF;
507   ts->ops->rollback       = TSRollBack_BDF;
508   ts->ops->interpolate    = TSInterpolate_BDF;
509   ts->ops->snesfunction   = SNESTSFormFunction_BDF;
510   ts->ops->snesjacobian   = SNESTSFormJacobian_BDF;
511   ts->default_adapt_type  = TSADAPTBASIC;
512 
513   ts->usessnes = PETSC_TRUE;
514 
515   PetscCall(PetscNew(&bdf));
516   ts->data = (void *)bdf;
517 
518   bdf->status = TS_STEP_COMPLETE;
519 
520   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBDFSetOrder_C", TSBDFSetOrder_BDF));
521   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBDFGetOrder_C", TSBDFGetOrder_BDF));
522   PetscCall(TSBDFSetOrder(ts, 2));
523   PetscFunctionReturn(PETSC_SUCCESS);
524 }
525 
526 /* ------------------------------------------------------------ */
527 
528 /*@
529   TSBDFSetOrder - Set the order of the `TSBDF` method
530 
531   Logically Collective
532 
533   Input Parameters:
534 +  ts - timestepping context
535 -  order - order of the method
536 
537   Options Database Key:
538 .  -ts_bdf_order <order> - select the order
539 
540   Level: intermediate
541 
542 .seealso: `TSBDFGetOrder()`, `TS`, `TSBDF`
543 @*/
544 PetscErrorCode TSBDFSetOrder(TS ts, PetscInt order)
545 {
546   PetscFunctionBegin;
547   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
548   PetscValidLogicalCollectiveInt(ts, order, 2);
549   PetscTryMethod(ts, "TSBDFSetOrder_C", (TS, PetscInt), (ts, order));
550   PetscFunctionReturn(PETSC_SUCCESS);
551 }
552 
553 /*@
554   TSBDFGetOrder - Get the order of the `TSBDF` method
555 
556   Not Collective
557 
558   Input Parameter:
559 .  ts - timestepping context
560 
561   Output Parameter:
562 .  order - order of the method
563 
564   Level: intermediate
565 
566 .seealso: `TSBDFSetOrder()`, `TS`, `TSBDF`
567 @*/
568 PetscErrorCode TSBDFGetOrder(TS ts, PetscInt *order)
569 {
570   PetscFunctionBegin;
571   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
572   PetscValidIntPointer(order, 2);
573   PetscUseMethod(ts, "TSBDFGetOrder_C", (TS, PetscInt *), (ts, order));
574   PetscFunctionReturn(PETSC_SUCCESS);
575 }
576