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