xref: /petsc/src/ts/impls/bdf/bdf.c (revision 609bdbee21ea3be08735c64dbe00a9ab27759925)
1 /*
2   Code for timestepping with BDF methods
3 */
4 #include <petsc/private/tsimpl.h>  /*I "petscts.h" I*/
5 
6 static PetscBool  cited = PETSC_FALSE;
7 static const char citation[] =
8   "@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   PetscReal shift;
20   Vec       vec_dot;
21   Vec       vec_lte;
22 
23   PetscInt     order;
24   PetscBool    adapt;
25   TSStepStatus status;
26 } TS_BDF;
27 
28 
29 PETSC_STATIC_INLINE void LagrangeBasisVals(PetscInt n,PetscReal t,const PetscReal T[],PetscScalar L[])
30 {
31   PetscInt k,j;
32   for (k=0; k<n; k++)
33     for (L[k]=1, j=0; j<n; j++)
34       if (j != k)
35         L[k] *= (t - T[j])/(T[k] - T[j]);
36 }
37 
38 PETSC_STATIC_INLINE void LagrangeBasisDers(PetscInt n,PetscReal t,const PetscReal T[],PetscScalar dL[])
39 {
40   PetscInt  k,j,i;
41   for (k=0; k<n; k++)
42     for (dL[k]=0, j=0; j<n; j++)
43       if (j != k) {
44         PetscReal L = 1/(T[k] - T[j]);
45         for (i=0; i<n; i++)
46           if (i != j && i != k)
47             L *= (t - T[i])/(T[k] - T[i]);
48         dL[k] += L;
49       }
50 }
51 
52 static PetscErrorCode TSBDF_Advance(TS ts,PetscReal t,Vec X)
53 {
54   TS_BDF         *bdf = (TS_BDF*)ts->data;
55   PetscInt       i,n = (PetscInt)(sizeof(bdf->work)/sizeof(Vec));
56   Vec            tail = bdf->work[n-1];
57   PetscErrorCode ierr;
58 
59   PetscFunctionBegin;
60   for (i=n-1; i>=2; i--) {
61     bdf->time[i] = bdf->time[i-1];
62     bdf->work[i] = bdf->work[i-1];
63   }
64   bdf->n       = PetscMin(bdf->n+1,n-1);
65   bdf->time[1] = t;
66   bdf->work[1] = tail;
67   ierr = VecCopy(X,tail);CHKERRQ(ierr);
68   PetscFunctionReturn(0);
69 }
70 
71 static PetscErrorCode TSBDF_VecDot(TS ts,PetscInt order,PetscReal t,Vec X,Vec Xdot,PetscReal *shift)
72 {
73   TS_BDF         *bdf = (TS_BDF*)ts->data;
74   PetscInt       i,n = order+1;
75   PetscReal      time[7];
76   Vec            vecs[7];
77   PetscScalar    alpha[7];
78   PetscErrorCode ierr;
79 
80   PetscFunctionBegin;
81   n = PetscMax(n,2);
82   time[0] = t; for (i=1; i<n; i++) time[i] = bdf->time[i];
83   vecs[0] = X; for (i=1; i<n; i++) vecs[i] = bdf->work[i];
84   LagrangeBasisDers(n,t,time,alpha);
85   ierr = VecZeroEntries(Xdot);CHKERRQ(ierr);
86   ierr = VecMAXPY(Xdot,n,alpha,vecs);CHKERRQ(ierr);
87   if (shift) *shift = PetscRealPart(alpha[0]);
88   PetscFunctionReturn(0);
89 }
90 
91 static PetscErrorCode TSBDF_VecLTE(TS ts,PetscInt order,Vec lte)
92 {
93   TS_BDF         *bdf = (TS_BDF*)ts->data;
94   PetscInt       i,n = order+1;
95   PetscReal      *time = bdf->time;
96   Vec            *vecs = bdf->work;
97   PetscScalar    a[8],b[8],alpha[8];
98   PetscErrorCode ierr;
99 
100   PetscFunctionBegin;
101   LagrangeBasisDers(n+0,time[0],time,a); a[n] =0;
102   LagrangeBasisDers(n+1,time[0],time,b);
103   for (i=0; i<n+1; i++) alpha[i] = (a[i]-b[i])/a[0];
104   ierr = VecZeroEntries(lte);CHKERRQ(ierr);
105   ierr = VecMAXPY(lte,n+1,alpha,vecs);CHKERRQ(ierr);
106   PetscFunctionReturn(0);
107 }
108 
109 static PetscErrorCode TSBDF_Extrapolate(TS ts,PetscInt order,PetscReal t,Vec X)
110 {
111   TS_BDF         *bdf = (TS_BDF*)ts->data;
112   PetscInt       n = order+1;
113   PetscReal      *time = bdf->time+1;
114   Vec            *vecs = bdf->work+1;
115   PetscScalar    alpha[7];
116   PetscErrorCode ierr;
117 
118   PetscFunctionBegin;
119   n = PetscMin(n,bdf->n);
120   LagrangeBasisVals(n,t,time,alpha);
121   ierr = VecZeroEntries(X);CHKERRQ(ierr);
122   ierr = VecMAXPY(X,n,alpha,vecs);CHKERRQ(ierr);
123   PetscFunctionReturn(0);
124 }
125 
126 static PetscErrorCode TSBDF_Interpolate(TS ts,PetscInt order,PetscReal t,Vec X)
127 {
128   TS_BDF         *bdf = (TS_BDF*)ts->data;
129   PetscInt       n = order+1;
130   PetscReal      *time = bdf->time;
131   Vec            *vecs = bdf->work;
132   PetscScalar    alpha[7];
133   PetscErrorCode ierr;
134 
135   PetscFunctionBegin;
136   LagrangeBasisVals(n,t,time,alpha);
137   ierr = VecZeroEntries(X);CHKERRQ(ierr);
138   ierr = VecMAXPY(X,n,alpha,vecs);CHKERRQ(ierr);
139   PetscFunctionReturn(0);
140 }
141 
142 static PetscErrorCode TS_SNESSolve(TS ts,Vec b,Vec x)
143 {
144   PetscInt       nits,lits;
145   PetscErrorCode ierr;
146 
147   PetscFunctionBegin;
148   ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr);
149   ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr);
150   ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
151   ts->snes_its += nits; ts->ksp_its += lits;
152   PetscFunctionReturn(0);
153 }
154 
155 static PetscErrorCode TSBDF_Restart(TS ts,PetscBool *accept)
156 {
157   TS_BDF         *bdf = (TS_BDF*)ts->data;
158   PetscErrorCode ierr;
159 
160   PetscFunctionBegin;
161   bdf->k = 1; bdf->n = 0;
162   ierr = TSBDF_Advance(ts,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
163 
164   bdf->time[0] = ts->ptime + ts->time_step/2;
165   ierr = VecCopy(bdf->work[1],bdf->work[0]);CHKERRQ(ierr);
166   ierr = TSPreStage(ts,bdf->time[0]);CHKERRQ(ierr);
167   ierr = TS_SNESSolve(ts,NULL,bdf->work[0]);CHKERRQ(ierr);
168   ierr = TSPostStage(ts,bdf->time[0],0,&bdf->work[0]);CHKERRQ(ierr);
169   ierr = TSAdaptCheckStage(ts->adapt,ts,bdf->time[0],bdf->work[0],accept);CHKERRQ(ierr);
170   if (!*accept) PetscFunctionReturn(0);
171 
172   bdf->k = PetscMin(2,bdf->order); bdf->n++;
173   ierr = VecCopy(bdf->work[0],bdf->work[2]);CHKERRQ(ierr);
174   bdf->time[2] = bdf->time[0];
175   PetscFunctionReturn(0);
176 }
177 
178 static const char *const BDF_SchemeName[] = {"", "1", "2", "3", "4", "5", "6"};
179 
180 static PetscErrorCode TSStep_BDF(TS ts)
181 {
182   TS_BDF         *bdf = (TS_BDF*)ts->data;
183   PetscInt       rejections = 0;
184   PetscBool      stageok,accept = PETSC_TRUE;
185   PetscReal      next_time_step = ts->time_step;
186   PetscErrorCode ierr;
187 
188   PetscFunctionBegin;
189   ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr);
190 
191   if (!ts->steprollback && !ts->steprestart) {
192     bdf->k = PetscMin(bdf->k+1,bdf->order);
193     ierr = TSBDF_Advance(ts,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
194   }
195 
196   bdf->status = TS_STEP_INCOMPLETE;
197   while (!ts->reason && bdf->status != TS_STEP_COMPLETE) {
198 
199     if (ts->steprestart) {
200       ierr = TSBDF_Restart(ts,&stageok);CHKERRQ(ierr);
201       if (!stageok) goto reject_step;
202     }
203 
204     bdf->time[0] = ts->ptime + ts->time_step;
205     ierr = TSBDF_Extrapolate(ts,bdf->k-(accept?0:1),bdf->time[0],bdf->work[0]);CHKERRQ(ierr);
206     ierr = TSPreStage(ts,bdf->time[0]);CHKERRQ(ierr);
207     ierr = TS_SNESSolve(ts,NULL,bdf->work[0]);CHKERRQ(ierr);
208     ierr = TSPostStage(ts,bdf->time[0],0,&bdf->work[0]);CHKERRQ(ierr);
209     ierr = TSAdaptCheckStage(ts->adapt,ts,bdf->time[0],bdf->work[0],&stageok);CHKERRQ(ierr);
210     if (!stageok) goto reject_step;
211 
212     bdf->status = TS_STEP_PENDING;
213     ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
214     ierr = TSAdaptCandidateAdd(ts->adapt,BDF_SchemeName[bdf->k],bdf->k,1,1.0,1.0,PETSC_TRUE);CHKERRQ(ierr);
215     ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
216     bdf->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
217     if (!accept) { ts->time_step = next_time_step; goto reject_step; }
218 
219     ierr = VecCopy(bdf->work[0],ts->vec_sol);CHKERRQ(ierr);
220     ts->ptime += ts->time_step;
221     ts->time_step = next_time_step;
222     break;
223 
224   reject_step:
225     ts->reject++; accept = PETSC_FALSE;
226     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
227       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
228       ts->reason = TS_DIVERGED_STEP_REJECTED;
229     }
230   }
231   PetscFunctionReturn(0);
232 }
233 
234 static PetscErrorCode TSInterpolate_BDF(TS ts,PetscReal t,Vec X)
235 {
236   TS_BDF         *bdf = (TS_BDF*)ts->data;
237   PetscErrorCode ierr;
238 
239   PetscFunctionBegin;
240   ierr = TSBDF_Interpolate(ts,bdf->k,t,X);CHKERRQ(ierr);
241   PetscFunctionReturn(0);
242 }
243 
244 static PetscErrorCode TSEvaluateWLTE_BDF(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
245 {
246   TS_BDF         *bdf = (TS_BDF*)ts->data;
247   PetscInt       k = bdf->k;
248   PetscReal      wltea,wlter;
249   Vec            X = bdf->work[0], Y = bdf->vec_lte;
250   PetscErrorCode ierr;
251 
252   PetscFunctionBegin;
253   k = PetscMin(k,bdf->n-1);
254   ierr = TSBDF_VecLTE(ts,k,Y);CHKERRQ(ierr);
255   ierr = VecAXPY(Y,1,X);CHKERRQ(ierr);
256   ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr);
257   if (order) *order = k + 1;
258   PetscFunctionReturn(0);
259 }
260 
261 static PetscErrorCode TSRollBack_BDF(TS ts)
262 {
263   TS_BDF         *bdf = (TS_BDF*)ts->data;
264   PetscErrorCode ierr;
265 
266   PetscFunctionBegin;
267   ierr = VecCopy(bdf->work[1],ts->vec_sol);CHKERRQ(ierr);
268   PetscFunctionReturn(0);
269 }
270 
271 static PetscErrorCode SNESTSFormFunction_BDF(PETSC_UNUSED SNES snes,Vec X,Vec F,TS ts)
272 {
273   TS_BDF         *bdf = (TS_BDF*)ts->data;
274   PetscInt       k = bdf->k;
275   PetscReal      t = bdf->time[0];
276   Vec            V = bdf->vec_dot;
277   PetscErrorCode ierr;
278 
279   PetscFunctionBegin;
280   ierr = TSBDF_VecDot(ts,k,t,X,V,&bdf->shift);CHKERRQ(ierr);
281   /* F = Function(t,X,V) */
282   ierr = TSComputeIFunction(ts,t,X,V,F,PETSC_FALSE);CHKERRQ(ierr);
283   PetscFunctionReturn(0);
284 }
285 
286 static PetscErrorCode SNESTSFormJacobian_BDF(PETSC_UNUSED SNES snes,
287                                              PETSC_UNUSED Vec X,
288                                              Mat J,Mat P,
289                                              TS ts)
290 {
291   TS_BDF         *bdf = (TS_BDF*)ts->data;
292   PetscReal      t = bdf->time[0];
293   Vec            V = bdf->vec_dot;
294   PetscReal      dVdX = bdf->shift;
295   PetscErrorCode ierr;
296 
297   PetscFunctionBegin;
298   /* J,P = Jacobian(t,X,V) */
299   ierr = TSComputeIJacobian(ts,t,X,V,dVdX,J,P,PETSC_FALSE);CHKERRQ(ierr);
300   PetscFunctionReturn(0);
301 }
302 
303 static PetscErrorCode TSReset_BDF(TS ts)
304 {
305   TS_BDF         *bdf = (TS_BDF*)ts->data;
306   size_t         i,n = sizeof(bdf->work)/sizeof(Vec);
307   PetscErrorCode ierr;
308 
309   PetscFunctionBegin;
310   for (i=0; i<n; i++) {ierr = VecDestroy(&bdf->work[i]);CHKERRQ(ierr);}
311   ierr = VecDestroy(&bdf->vec_dot);CHKERRQ(ierr);
312   ierr = VecDestroy(&bdf->vec_lte);CHKERRQ(ierr);
313   PetscFunctionReturn(0);
314 }
315 
316 static PetscErrorCode TSDestroy_BDF(TS ts)
317 {
318   PetscErrorCode ierr;
319 
320   PetscFunctionBegin;
321   ierr = TSReset_BDF(ts);CHKERRQ(ierr);
322   ierr = PetscFree(ts->data);CHKERRQ(ierr);
323   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBDFSetOrder_C",NULL);CHKERRQ(ierr);
324   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBDFGetOrder_C",NULL);CHKERRQ(ierr);
325   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBDFUseAdapt_C",NULL);CHKERRQ(ierr);
326   PetscFunctionReturn(0);
327 }
328 
329 static PetscErrorCode TSSetUp_BDF(TS ts)
330 {
331   TS_BDF         *bdf = (TS_BDF*)ts->data;
332   size_t         i,n = sizeof(bdf->work)/sizeof(Vec);
333   PetscErrorCode ierr;
334 
335   PetscFunctionBegin;
336   bdf->k = bdf->n = 0;
337   for (i=0; i<n; i++) {ierr = VecDuplicate(ts->vec_sol,&bdf->work[i]);CHKERRQ(ierr);}
338   ierr = VecDuplicate(ts->vec_sol,&bdf->vec_dot);CHKERRQ(ierr);
339   ierr = VecDuplicate(ts->vec_sol,&bdf->vec_lte);CHKERRQ(ierr);
340 
341   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
342   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
343   if (!bdf->adapt) {
344     ierr = TSAdaptSetType(ts->adapt,TSADAPTNONE);CHKERRQ(ierr);
345   } else {
346     PetscReal low,high;
347     ierr = TSAdaptBasicGetClip(ts->adapt,&low,&high);CHKERRQ(ierr);
348     high = PetscMin(high,2.0);
349     ierr = TSAdaptBasicSetClip(ts->adapt,low,high);CHKERRQ(ierr);
350     if (ts->exact_final_time == TS_EXACTFINALTIME_UNSPECIFIED)
351       ts->exact_final_time = TS_EXACTFINALTIME_MATCHSTEP;
352   }
353 
354   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
355   PetscFunctionReturn(0);
356 }
357 
358 static PetscErrorCode TSSetFromOptions_BDF(PetscOptionItems *PetscOptionsObject,TS ts)
359 {
360   TS_BDF         *bdf = (TS_BDF*)ts->data;
361   PetscErrorCode ierr;
362 
363   PetscFunctionBegin;
364   ierr = PetscOptionsHead(PetscOptionsObject,"BDF ODE solver options");CHKERRQ(ierr);
365   {
366     PetscBool flg;
367     PetscInt  order = bdf->order;
368     PetscBool adapt = bdf->adapt;
369     ierr = PetscOptionsInt("-ts_bdf_order","Order of the BDF method","TSBDFSetOrder",order,&order,&flg);CHKERRQ(ierr);
370     if (flg) {ierr = TSBDFSetOrder(ts,order);CHKERRQ(ierr);}
371     ierr = PetscOptionsBool("-ts_bdf_adapt","Use time-step adaptivity with the BDF method","TSBDFUseAdapt",adapt,&adapt,&flg);CHKERRQ(ierr);
372     if (flg) {ierr = TSBDFUseAdapt(ts,adapt);CHKERRQ(ierr);}
373   }
374   ierr = PetscOptionsTail();CHKERRQ(ierr);
375   PetscFunctionReturn(0);
376 }
377 
378 static PetscErrorCode TSView_BDF(TS ts,PetscViewer viewer)
379 {
380   TS_BDF         *bdf = (TS_BDF*)ts->data;
381   PetscBool      iascii;
382   PetscErrorCode ierr;
383 
384   PetscFunctionBegin;
385   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
386   if (iascii)    {
387     ierr = PetscViewerASCIIPrintf(viewer,"  Order=%D\n",bdf->order);CHKERRQ(ierr);
388   }
389   if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);}
390   if (ts->snes)  {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
391   PetscFunctionReturn(0);
392 }
393 
394 /* ------------------------------------------------------------ */
395 
396 static PetscErrorCode TSBDFSetOrder_BDF(TS ts,PetscInt order)
397 {
398   TS_BDF *bdf = (TS_BDF*)ts->data;
399 
400   PetscFunctionBegin;
401   if (order == bdf->order) PetscFunctionReturn(0);
402   if (order < 1 || order > 6) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"BDF Order %D not implemented",order);
403   bdf->order = order;
404   PetscFunctionReturn(0);
405 }
406 
407 static PetscErrorCode TSBDFGetOrder_BDF(TS ts,PetscInt *order)
408 {
409   TS_BDF *bdf = (TS_BDF*)ts->data;
410 
411   PetscFunctionBegin;
412   *order = bdf->order;
413   PetscFunctionReturn(0);
414 }
415 
416 static PetscErrorCode TSBDFUseAdapt_BDF(TS ts,PetscBool use)
417 {
418   TS_BDF *bdf = (TS_BDF*)ts->data;
419 
420   PetscFunctionBegin;
421   if (use == bdf->adapt) PetscFunctionReturn(0);
422   if (ts->setupcalled) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ORDER,"Cannot change adaptivity after TSSetUp()");
423   bdf->adapt = use;
424   PetscFunctionReturn(0);
425 }
426 
427 /* ------------------------------------------------------------ */
428 
429 /*MC
430       TSBDF - DAE solver using BDF methods
431 
432   Level: beginner
433 
434 .seealso:  TS, TSCreate(), TSSetType()
435 M*/
436 PETSC_EXTERN PetscErrorCode TSCreate_BDF(TS ts)
437 {
438   TS_BDF         *bdf;
439   PetscErrorCode ierr;
440 
441   PetscFunctionBegin;
442   ts->ops->reset          = TSReset_BDF;
443   ts->ops->destroy        = TSDestroy_BDF;
444   ts->ops->view           = TSView_BDF;
445   ts->ops->setup          = TSSetUp_BDF;
446   ts->ops->setfromoptions = TSSetFromOptions_BDF;
447   ts->ops->step           = TSStep_BDF;
448   ts->ops->evaluatewlte   = TSEvaluateWLTE_BDF;
449   ts->ops->rollback       = TSRollBack_BDF;
450   ts->ops->interpolate    = TSInterpolate_BDF;
451   ts->ops->snesfunction   = SNESTSFormFunction_BDF;
452   ts->ops->snesjacobian   = SNESTSFormJacobian_BDF;
453 
454   ierr = PetscNewLog(ts,&bdf);CHKERRQ(ierr);
455   ts->data = (void*)bdf;
456 
457   bdf->order  = 2;
458   bdf->adapt  = PETSC_FALSE;
459   bdf->status = TS_STEP_COMPLETE;
460 
461   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBDFSetOrder_C",TSBDFSetOrder_BDF);CHKERRQ(ierr);
462   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBDFGetOrder_C",TSBDFGetOrder_BDF);CHKERRQ(ierr);
463   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBDFUseAdapt_C",TSBDFUseAdapt_BDF);CHKERRQ(ierr);
464   PetscFunctionReturn(0);
465 }
466 
467 /* ------------------------------------------------------------ */
468 
469 /*@
470   TSBDFSetOrder - Set the order of the BDF method
471 
472   Logically Collective on TS
473 
474   Input Parameter:
475 +  ts - timestepping context
476 -  order - order of the method
477 
478   Options Database:
479 .  -ts_bdf_order <order>
480 
481   Level: intermediate
482 
483 @*/
484 PetscErrorCode TSBDFSetOrder(TS ts,PetscInt order)
485 {
486   PetscErrorCode ierr;
487 
488   PetscFunctionBegin;
489   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
490   PetscValidLogicalCollectiveInt(ts,order,2);
491   ierr = PetscTryMethod(ts,"TSBDFSetOrder_C",(TS,PetscInt),(ts,order));CHKERRQ(ierr);
492   PetscFunctionReturn(0);
493 }
494 
495 /*@
496   TSBDFGetOrder - Get the order of the BDF method
497 
498   Not Collective
499 
500   Input Parameter:
501 .  ts - timestepping context
502 
503   Output Parameter:
504 .  order - order of the method
505 
506   Level: intermediate
507 
508 @*/
509 PetscErrorCode TSBDFGetOrder(TS ts,PetscInt *order)
510 {
511   PetscErrorCode ierr;
512 
513   PetscFunctionBegin;
514   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
515   PetscValidIntPointer(order,2);
516   ierr = PetscUseMethod(ts,"TSBDFGetOrder_C",(TS,PetscInt*),(ts,order));CHKERRQ(ierr);
517   PetscFunctionReturn(0);
518 }
519 
520 /*@
521   TSBDFUseAdapt - Use time-step adaptivity with the BDF method
522 
523   Logically Collective on TS
524 
525   Input Parameter:
526 +  ts - timestepping context
527 -  use - flag to use adaptivity
528 
529   Options Database:
530 .  -ts_bdf_adapt
531 
532   Level: intermediate
533 
534 .seealso: TSAdapt
535 @*/
536 PetscErrorCode TSBDFUseAdapt(TS ts,PetscBool use)
537 {
538   PetscErrorCode ierr;
539 
540   PetscFunctionBegin;
541   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
542   PetscValidLogicalCollectiveBool(ts,use,2);
543   ierr = PetscTryMethod(ts,"TSBDFUseAdapt_C",(TS,PetscBool),(ts,use));CHKERRQ(ierr);
544   PetscFunctionReturn(0);
545 }
546