xref: /petsc/src/ts/adapt/interface/tsadapt.c (revision 0f51fdf84cc2af9c6ef13eed784e698b930e50d2)
1 
2 #include <petsc-private/tsimpl.h> /*I  "petscts.h" I*/
3 
4 static PetscFunctionList TSAdaptList;
5 static PetscBool         TSAdaptPackageInitialized;
6 static PetscBool         TSAdaptRegisterAllCalled;
7 static PetscClassId      TSADAPT_CLASSID;
8 
9 PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt);
10 PETSC_EXTERN PetscErrorCode TSAdaptCreate_None(TSAdapt);
11 PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt);
12 
13 #undef __FUNCT__
14 #define __FUNCT__ "TSAdaptRegister"
15 /*@C
16    TSAdaptRegister -  adds a TSAdapt implementation
17 
18    Not Collective
19 
20    Input Parameters:
21 +  name_scheme - name of user-defined adaptivity scheme
22 -  routine_create - routine to create method context
23 
24    Notes:
25    TSAdaptRegister() may be called multiple times to add several user-defined families.
26 
27    Sample usage:
28 .vb
29    TSAdaptRegister("my_scheme",MySchemeCreate);
30 .ve
31 
32    Then, your scheme can be chosen with the procedural interface via
33 $     TSAdaptSetType(ts,"my_scheme")
34    or at runtime via the option
35 $     -ts_adapt_type my_scheme
36 
37    Level: advanced
38 
39 .keywords: TSAdapt, register
40 
41 .seealso: TSAdaptRegisterAll()
42 @*/
43 PetscErrorCode  TSAdaptRegister(const char sname[],PetscErrorCode (*function)(TSAdapt))
44 {
45   PetscErrorCode ierr;
46 
47   PetscFunctionBegin;
48   ierr = PetscFunctionListAdd(&TSAdaptList,sname,function);CHKERRQ(ierr);
49   PetscFunctionReturn(0);
50 }
51 
52 #undef __FUNCT__
53 #define __FUNCT__ "TSAdaptRegisterAll"
54 /*@C
55   TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt
56 
57   Not Collective
58 
59   Level: advanced
60 
61 .keywords: TSAdapt, register, all
62 
63 .seealso: TSAdaptRegisterDestroy()
64 @*/
65 PetscErrorCode  TSAdaptRegisterAll(void)
66 {
67   PetscErrorCode ierr;
68 
69   PetscFunctionBegin;
70   if (TSAdaptRegisterAllCalled) PetscFunctionReturn(0);
71   TSAdaptRegisterAllCalled = PETSC_TRUE;
72   ierr = TSAdaptRegister(TSADAPTBASIC,TSAdaptCreate_Basic);CHKERRQ(ierr);
73   ierr = TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None);CHKERRQ(ierr);
74   ierr = TSAdaptRegister(TSADAPTCFL,  TSAdaptCreate_CFL);CHKERRQ(ierr);
75   PetscFunctionReturn(0);
76 }
77 
78 #undef __FUNCT__
79 #define __FUNCT__ "TSAdaptFinalizePackage"
80 /*@C
81   TSFinalizePackage - This function destroys everything in the TS package. It is
82   called from PetscFinalize().
83 
84   Level: developer
85 
86 .keywords: Petsc, destroy, package
87 .seealso: PetscFinalize()
88 @*/
89 PetscErrorCode  TSAdaptFinalizePackage(void)
90 {
91   PetscErrorCode ierr;
92 
93   PetscFunctionBegin;
94   ierr = PetscFunctionListDestroy(&TSAdaptList);CHKERRQ(ierr);
95   TSAdaptPackageInitialized = PETSC_FALSE;
96   TSAdaptRegisterAllCalled  = PETSC_FALSE;
97   PetscFunctionReturn(0);
98 }
99 
100 #undef __FUNCT__
101 #define __FUNCT__ "TSAdaptInitializePackage"
102 /*@C
103   TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is
104   called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to
105   TSCreate_GL() when using static libraries.
106 
107   Level: developer
108 
109 .keywords: TSAdapt, initialize, package
110 .seealso: PetscInitialize()
111 @*/
112 PetscErrorCode  TSAdaptInitializePackage(void)
113 {
114   PetscErrorCode ierr;
115 
116   PetscFunctionBegin;
117   if (TSAdaptPackageInitialized) PetscFunctionReturn(0);
118   TSAdaptPackageInitialized = PETSC_TRUE;
119   ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr);
120   ierr = TSAdaptRegisterAll();CHKERRQ(ierr);
121   ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr);
122   PetscFunctionReturn(0);
123 }
124 
125 #undef __FUNCT__
126 #define __FUNCT__ "TSAdaptSetType"
127 PetscErrorCode  TSAdaptSetType(TSAdapt adapt,TSAdaptType type)
128 {
129   PetscErrorCode ierr,(*r)(TSAdapt);
130 
131   PetscFunctionBegin;
132   ierr = PetscFunctionListFind(TSAdaptList,type,&r);CHKERRQ(ierr);
133   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type);
134   if (((PetscObject)adapt)->type_name) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);}
135   /* Reinitialize function pointers in TSAdaptOps structure */
136   ierr = PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));CHKERRQ(ierr);
137   ierr = (*r)(adapt);CHKERRQ(ierr);
138   ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr);
139   PetscFunctionReturn(0);
140 }
141 
142 #undef __FUNCT__
143 #define __FUNCT__ "TSAdaptSetOptionsPrefix"
144 PetscErrorCode  TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[])
145 {
146   PetscErrorCode ierr;
147 
148   PetscFunctionBegin;
149   ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr);
150   PetscFunctionReturn(0);
151 }
152 
153 #undef __FUNCT__
154 #define __FUNCT__ "TSAdaptLoad"
155 /*@C
156   TSAdaptLoad - Loads a TSAdapt that has been stored in binary  with TSAdaptView().
157 
158   Collective on PetscViewer
159 
160   Input Parameters:
161 + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or
162            some related function before a call to TSAdaptLoad().
163 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
164            HDF5 file viewer, obtained from PetscViewerHDF5Open()
165 
166    Level: intermediate
167 
168   Notes:
169    The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored.
170 
171   Notes for advanced users:
172   Most users should not need to know the details of the binary storage
173   format, since TSAdaptLoad() and TSAdaptView() completely hide these details.
174   But for anyone who's interested, the standard binary matrix storage
175   format is
176 .vb
177      has not yet been determined
178 .ve
179 
180 .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad()
181 @*/
182 PetscErrorCode  TSAdaptLoad(TSAdapt tsadapt, PetscViewer viewer)
183 {
184   PetscErrorCode ierr;
185   PetscBool      isbinary;
186   char           type[256];
187 
188   PetscFunctionBegin;
189   PetscValidHeaderSpecific(tsadapt,TSADAPT_CLASSID,1);
190   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
191   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
192   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
193 
194   ierr = PetscViewerBinaryRead(viewer,type,256,PETSC_CHAR);CHKERRQ(ierr);
195   ierr = TSAdaptSetType(tsadapt, type);CHKERRQ(ierr);
196   if (tsadapt->ops->load) {
197     ierr = (*tsadapt->ops->load)(tsadapt,viewer);CHKERRQ(ierr);
198   }
199   PetscFunctionReturn(0);
200 }
201 
202 #undef __FUNCT__
203 #define __FUNCT__ "TSAdaptView"
204 PetscErrorCode  TSAdaptView(TSAdapt adapt,PetscViewer viewer)
205 {
206   PetscErrorCode ierr;
207   PetscBool      iascii,isbinary;
208 
209   PetscFunctionBegin;
210   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
211   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
212   if (iascii) {
213     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr);
214     ierr = PetscViewerASCIIPrintf(viewer,"  number of candidates %D\n",adapt->candidates.n);CHKERRQ(ierr);
215     if (adapt->ops->view) {
216       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
217       ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr);
218       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
219     }
220   } else if (isbinary) {
221     char type[256];
222 
223     /* need to save FILE_CLASS_ID for adapt class */
224     ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr);
225     ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr);
226   } else if (adapt->ops->view) {
227     ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr);
228   }
229   PetscFunctionReturn(0);
230 }
231 
232 #undef __FUNCT__
233 #define __FUNCT__ "TSAdaptDestroy"
234 PetscErrorCode  TSAdaptDestroy(TSAdapt *adapt)
235 {
236   PetscErrorCode ierr;
237 
238   PetscFunctionBegin;
239   if (!*adapt) PetscFunctionReturn(0);
240   PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1);
241   if (--((PetscObject)(*adapt))->refct > 0) {*adapt = 0; PetscFunctionReturn(0);}
242   if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);}
243   ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr);
244   ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr);
245   PetscFunctionReturn(0);
246 }
247 
248 #undef __FUNCT__
249 #define __FUNCT__ "TSAdaptSetMonitor"
250 /*@
251    TSAdaptSetMonitor - Monitor the choices made by the adaptive controller
252 
253    Collective on TSAdapt
254 
255    Input Arguments:
256 +  adapt - adaptive controller context
257 -  flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable
258 
259    Level: intermediate
260 
261 .seealso: TSAdaptChoose()
262 @*/
263 PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg)
264 {
265   PetscErrorCode ierr;
266 
267   PetscFunctionBegin;
268   if (flg) {
269     if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);}
270   } else {
271     ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr);
272   }
273   PetscFunctionReturn(0);
274 }
275 
276 #undef __FUNCT__
277 #define __FUNCT__ "TSAdaptSetCheckStage"
278 /*@C
279    TSAdaptSetCheckStage - set a callback to check convergence for a stage
280 
281    Logically collective on TSAdapt
282 
283    Input Arguments:
284 +  adapt - adaptive controller context
285 -  func - stage check function
286 
287    Arguments of func:
288 $  PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept)
289 
290 +  adapt - adaptive controller context
291 .  ts - time stepping context
292 -  accept - pending choice of whether to accept, can be modified by this routine
293 
294    Level: advanced
295 
296 .seealso: TSAdaptChoose()
297 @*/
298 PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscBool*))
299 {
300 
301   PetscFunctionBegin;
302   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
303   adapt->ops->checkstage = func;
304   PetscFunctionReturn(0);
305 }
306 
307 #undef __FUNCT__
308 #define __FUNCT__ "TSAdaptSetStepLimits"
309 /*@
310    TSAdaptSetStepLimits - Set minimum and maximum step sizes to be considered by the controller
311 
312    Logically Collective
313 
314    Input Arguments:
315 +  adapt - time step adaptivity context, usually gotten with TSGetAdapt()
316 .  hmin - minimum time step
317 -  hmax - maximum time step
318 
319    Options Database Keys:
320 +  -ts_adapt_dt_min - minimum time step
321 -  -ts_adapt_dt_max - maximum time step
322 
323    Level: intermediate
324 
325 .seealso: TSAdapt
326 @*/
327 PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax)
328 {
329 
330   PetscFunctionBegin;
331   if (hmin != PETSC_DECIDE) adapt->dt_min = hmin;
332   if (hmax != PETSC_DECIDE) adapt->dt_max = hmax;
333   PetscFunctionReturn(0);
334 }
335 
336 #undef __FUNCT__
337 #define __FUNCT__ "TSAdaptSetFromOptions"
338 /*@
339    TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options.
340 
341    Collective on TSAdapt
342 
343    Input Parameter:
344 .  adapt - the TSAdapt context
345 
346    Options Database Keys:
347 .  -ts_adapt_type <type> - basic
348 
349    Level: advanced
350 
351    Notes:
352    This function is automatically called by TSSetFromOptions()
353 
354 .keywords: TS, TSGetAdapt(), TSAdaptSetType()
355 
356 .seealso: TSGetType()
357 @*/
358 PetscErrorCode  TSAdaptSetFromOptions(TSAdapt adapt)
359 {
360   PetscErrorCode ierr;
361   char           type[256] = TSADAPTBASIC;
362   PetscBool      set,flg;
363 
364   PetscFunctionBegin;
365   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
366   * function can only be called from inside TSSetFromOptions_GL()  */
367   ierr = PetscOptionsHead("TS Adaptivity options");CHKERRQ(ierr);
368   ierr = PetscOptionsFList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList,
369                           ((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);CHKERRQ(ierr);
370   if (flg || !((PetscObject)adapt)->type_name) {
371     ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr);
372   }
373   if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(adapt);CHKERRQ(ierr);}
374   ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",adapt->dt_min,&adapt->dt_min,NULL);CHKERRQ(ierr);
375   ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",adapt->dt_max,&adapt->dt_max,NULL);CHKERRQ(ierr);
376   ierr = PetscOptionsReal("-ts_adapt_scale_solve_failed","Scale step by this factor if solve fails","",adapt->scale_solve_failed,&adapt->scale_solve_failed,NULL);CHKERRQ(ierr);
377   ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
378   if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);}
379   ierr = PetscOptionsTail();CHKERRQ(ierr);
380   PetscFunctionReturn(0);
381 }
382 
383 #undef __FUNCT__
384 #define __FUNCT__ "TSAdaptCandidatesClear"
385 /*@
386    TSAdaptCandidatesClear - clear any previously set candidate schemes
387 
388    Logically Collective
389 
390    Input Argument:
391 .  adapt - adaptive controller
392 
393    Level: developer
394 
395 .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose()
396 @*/
397 PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
398 {
399   PetscErrorCode ierr;
400 
401   PetscFunctionBegin;
402   ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr);
403   PetscFunctionReturn(0);
404 }
405 
406 #undef __FUNCT__
407 #define __FUNCT__ "TSAdaptCandidateAdd"
408 /*@C
409    TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from
410 
411    Logically Collective
412 
413    Input Arguments:
414 +  adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate()
415 .  name - name of the candidate scheme to add
416 .  order - order of the candidate scheme
417 .  stageorder - stage order of the candidate scheme
418 .  ccfl - stability coefficient relative to explicit Euler, used for CFL constraints
419 .  cost - relative measure of the amount of work required for the candidate scheme
420 -  inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme
421 
422    Note:
423    This routine is not available in Fortran.
424 
425    Level: developer
426 
427 .seealso: TSAdaptCandidatesClear(), TSAdaptChoose()
428 @*/
429 PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse)
430 {
431   PetscInt c;
432 
433   PetscFunctionBegin;
434   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
435   if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order);
436   if (inuse) {
437     if (adapt->candidates.inuse_set) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()");
438     adapt->candidates.inuse_set = PETSC_TRUE;
439   }
440   /* first slot if this is the current scheme, otherwise the next available slot */
441   c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n;
442 
443   adapt->candidates.name[c]       = name;
444   adapt->candidates.order[c]      = order;
445   adapt->candidates.stageorder[c] = stageorder;
446   adapt->candidates.ccfl[c]       = ccfl;
447   adapt->candidates.cost[c]       = cost;
448   adapt->candidates.n++;
449   PetscFunctionReturn(0);
450 }
451 
452 #undef __FUNCT__
453 #define __FUNCT__ "TSAdaptCandidatesGet"
454 /*@C
455    TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost
456 
457    Not Collective
458 
459    Input Arguments:
460 .  adapt - time step adaptivity context
461 
462    Output Arguments:
463 +  n - number of candidate schemes, always at least 1
464 .  order - the order of each candidate scheme
465 .  stageorder - the stage order of each candidate scheme
466 .  ccfl - the CFL coefficient of each scheme
467 -  cost - the relative cost of each scheme
468 
469    Level: developer
470 
471    Note:
472    The current scheme is always returned in the first slot
473 
474 .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose()
475 @*/
476 PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost)
477 {
478   PetscFunctionBegin;
479   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
480   if (n) *n = adapt->candidates.n;
481   if (order) *order = adapt->candidates.order;
482   if (stageorder) *stageorder = adapt->candidates.stageorder;
483   if (ccfl) *ccfl = adapt->candidates.ccfl;
484   if (cost) *cost = adapt->candidates.cost;
485   PetscFunctionReturn(0);
486 }
487 
488 #undef __FUNCT__
489 #define __FUNCT__ "TSAdaptChoose"
490 /*@C
491    TSAdaptChoose - choose which method and step size to use for the next step
492 
493    Logically Collective
494 
495    Input Arguments:
496 +  adapt - adaptive contoller
497 -  h - current step size
498 
499    Output Arguments:
500 +  next_sc - scheme to use for the next step
501 .  next_h - step size to use for the next step
502 -  accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size
503 
504    Note:
505    The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is
506    being retried after an initial rejection.
507 
508    Level: developer
509 
510 .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd()
511 @*/
512 PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept)
513 {
514   PetscErrorCode ierr;
515   PetscReal      wlte = -1.0;
516 
517   PetscFunctionBegin;
518   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
519   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
520   PetscValidIntPointer(next_sc,4);
521   PetscValidPointer(next_h,5);
522   PetscValidIntPointer(accept,6);
523   if (adapt->candidates.n < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"%D candidates have been registered",adapt->candidates.n);
524   if (!adapt->candidates.inuse_set) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"The current in-use scheme is not among the %D candidates",adapt->candidates.n);
525   ierr = (*adapt->ops->choose)(adapt,ts,h,next_sc,next_h,accept,&wlte);CHKERRQ(ierr);
526   if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
527     /* Reduce time step if it overshoots max time */
528     PetscReal max_time = ts->max_time;
529     PetscReal next_dt  = 0.0;
530     if (ts->ptime + ts->time_step + *next_h >= max_time) {
531       next_dt = max_time - (ts->ptime + ts->time_step);
532       if (next_dt > PETSC_SMALL) *next_h = next_dt;
533       else ts->reason = TS_CONVERGED_TIME;
534     }
535   }
536   if (*next_sc < 0 || adapt->candidates.n <= *next_sc) SETERRQ2(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Chosen scheme %D not in valid range 0..%D",*next_sc,adapt->candidates.n-1);
537   if (!(*next_h > 0.)) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h);
538 
539   if (adapt->monitor) {
540     ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
541     if (wlte < 0) {
542       ierr = PetscViewerASCIIPrintf(adapt->monitor,"    TSAdapt '%s': step %3D %s t=%-11g+%10.3e family='%s' scheme=%D:'%s' dt=%-10g\n",((PetscObject)adapt)->type_name,ts->steps,*accept ? "accepted" : "rejected",(double)ts->ptime,(double)h,((PetscObject)ts)->type_name,*next_sc,adapt->candidates.name[*next_sc],(double)*next_h);CHKERRQ(ierr);
543     } else {
544       ierr = PetscViewerASCIIPrintf(adapt->monitor,"    TSAdapt '%s': step %3D %s t=%-11g+%10.3e wlte=%5.3g family='%s' scheme=%D:'%s' dt=%-10.3e\n",((PetscObject)adapt)->type_name,ts->steps,*accept ? "accepted" : "rejected",(double)ts->ptime,(double)h,(double)wlte,((PetscObject)ts)->type_name,*next_sc,adapt->candidates.name[*next_sc],(double)*next_h);CHKERRQ(ierr);
545     }
546     ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
547   }
548   PetscFunctionReturn(0);
549 }
550 
551 #undef __FUNCT__
552 #define __FUNCT__ "TSAdaptCheckStage"
553 /*@
554    TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails)
555 
556    Collective
557 
558    Input Arguments:
559 +  adapt - adaptive controller context
560 -  ts - time stepper
561 
562    Output Arguments:
563 .  accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject
564 
565    Level: developer
566 
567 .seealso:
568 @*/
569 PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscBool *accept)
570 {
571   PetscErrorCode      ierr;
572   SNES                snes;
573   SNESConvergedReason snesreason;
574 
575   PetscFunctionBegin;
576   *accept = PETSC_TRUE;
577   ierr    = TSGetSNES(ts,&snes);CHKERRQ(ierr);
578   ierr    = SNESGetConvergedReason(snes,&snesreason);CHKERRQ(ierr);
579   if (snesreason < 0) {
580     PetscReal dt,new_dt;
581     *accept = PETSC_FALSE;
582     ierr    = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
583     if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) {
584       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
585       ierr = PetscInfo2(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);CHKERRQ(ierr);
586       if (adapt->monitor) {
587         ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
588         ierr = PetscViewerASCIIPrintf(adapt->monitor,"    TSAdapt '%s': step %3D stage rejected t=%-11g+%10.3e, %D failures exceeds current TS allowed\n",((PetscObject)adapt)->type_name,ts->steps,(double)ts->ptime,dt,ts->num_snes_failures);CHKERRQ(ierr);
589         ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
590       }
591     } else {
592       new_dt = dt*adapt->scale_solve_failed;
593       ierr   = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr);
594       if (adapt->monitor) {
595         ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
596         ierr = PetscViewerASCIIPrintf(adapt->monitor,"    TSAdapt '%s': step %3D stage rejected t=%-11g+%10.3e retrying with dt=%-10.3e\n",((PetscObject)adapt)->type_name,ts->steps,(double)ts->ptime,(double)dt,(double)new_dt);CHKERRQ(ierr);
597         ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
598       }
599     }
600   }
601   if (adapt->ops->checkstage) {ierr = (*adapt->ops->checkstage)(adapt,ts,accept);CHKERRQ(ierr);}
602   PetscFunctionReturn(0);
603 }
604 
605 
606 
607 #undef __FUNCT__
608 #define __FUNCT__ "TSAdaptCreate"
609 /*@
610   TSAdaptCreate - create an adaptive controller context for time stepping
611 
612   Collective on MPI_Comm
613 
614   Input Parameter:
615 . comm - The communicator
616 
617   Output Parameter:
618 . adapt - new TSAdapt object
619 
620   Level: developer
621 
622   Notes:
623   TSAdapt creation is handled by TS, so users should not need to call this function.
624 
625 .keywords: TSAdapt, create
626 .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy()
627 @*/
628 PetscErrorCode  TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt)
629 {
630   PetscErrorCode ierr;
631   TSAdapt        adapt;
632 
633   PetscFunctionBegin;
634   *inadapt = 0;
635   ierr     = PetscHeaderCreate(adapt,_p_TSAdapt,struct _TSAdaptOps,TSADAPT_CLASSID,"TSAdapt","General Linear adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr);
636 
637   adapt->dt_min             = 1e-20;
638   adapt->dt_max             = 1e50;
639   adapt->scale_solve_failed = 0.25;
640 
641   *inadapt = adapt;
642   PetscFunctionReturn(0);
643 }
644