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