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