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