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