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