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