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