xref: /petsc/src/ts/adapt/interface/tsadapt.c (revision fa59e80594ae4f5b96fb1a0d88f82abc4bbaec69)
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,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 = PetscObjectTypeCompare((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     if (adapt->ops->view) {
173       ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr);
174     }
175   }
176   PetscFunctionReturn(0);
177 }
178 
179 #undef __FUNCT__
180 #define __FUNCT__ "TSAdaptDestroy"
181 PetscErrorCode  TSAdaptDestroy(TSAdapt *adapt)
182 {
183   PetscErrorCode ierr;
184 
185   PetscFunctionBegin;
186   if (!*adapt) PetscFunctionReturn(0);
187   PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1);
188   if (--((PetscObject)(*adapt))->refct > 0) {*adapt = 0; PetscFunctionReturn(0);}
189   if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);}
190   ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr);
191   ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr);
192   PetscFunctionReturn(0);
193 }
194 
195 #undef __FUNCT__
196 #define __FUNCT__ "TSAdaptSetMonitor"
197 /*@
198    TSAdaptSetMonitor - Monitor the choices made by the adaptive controller
199 
200    Collective on TSAdapt
201 
202    Input Arguments:
203 +  adapt - adaptive controller context
204 -  flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable
205 
206    Level: intermediate
207 
208 .seealso: TSAdaptChoose()
209 @*/
210 PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg)
211 {
212   PetscErrorCode ierr;
213 
214   PetscFunctionBegin;
215   if (flg) {
216     if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(((PetscObject)adapt)->comm,"stdout",&adapt->monitor);CHKERRQ(ierr);}
217   } else {
218     ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr);
219   }
220   PetscFunctionReturn(0);
221 }
222 
223 #undef __FUNCT__
224 #define __FUNCT__ "TSAdaptSetCheckStage"
225 /*@C
226    TSAdaptSetCheckStage - set a callback to check convergence for a stage
227 
228    Logically collective on TSAdapt
229 
230    Input Arguments:
231 +  adapt - adaptive controller context
232 -  func - stage check function
233 
234    Arguments of func:
235 $  PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept)
236 
237 +  adapt - adaptive controller context
238 .  ts - time stepping context
239 -  accept - pending choice of whether to accept, can be modified by this routine
240 
241    Level: advanced
242 
243 .seealso: TSAdaptChoose()
244 @*/
245 PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscBool*))
246 {
247 
248   PetscFunctionBegin;
249   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
250   adapt->ops->checkstage = func;
251   PetscFunctionReturn(0);
252 }
253 
254 #undef __FUNCT__
255 #define __FUNCT__ "TSAdaptSetStepLimits"
256 /*@
257    TSAdaptSetStepLimits - Set minimum and maximum step sizes to be considered by the controller
258 
259    Logically Collective
260 
261    Input Arguments:
262 +  adapt - time step adaptivity context, usually gotten with TSGetAdapt()
263 .  hmin - minimum time step
264 -  hmax - maximum time step
265 
266    Options Database Keys:
267 +  -ts_adapt_dt_min - minimum time step
268 -  -ts_adapt_dt_max - maximum time step
269 
270    Level: intermediate
271 
272 .seealso: TSAdapt
273 @*/
274 PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax)
275 {
276 
277   PetscFunctionBegin;
278   if (hmin != PETSC_DECIDE) adapt->dt_min = hmin;
279   if (hmax != PETSC_DECIDE) adapt->dt_max = hmax;
280   PetscFunctionReturn(0);
281 }
282 
283 #undef __FUNCT__
284 #define __FUNCT__ "TSAdaptSetFromOptions"
285 /*@
286    TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options.
287 
288    Collective on TSAdapt
289 
290    Input Parameter:
291 .  adapt - the TSAdapt context
292 
293    Options Database Keys:
294 .  -ts_adapt_type <type> - basic
295 
296    Level: advanced
297 
298    Notes:
299    This function is automatically called by TSSetFromOptions()
300 
301 .keywords: TS, TSGetAdapt(), TSAdaptSetType()
302 
303 .seealso: TSGetType()
304 @*/
305 PetscErrorCode  TSAdaptSetFromOptions(TSAdapt adapt)
306 {
307   PetscErrorCode ierr;
308   char           type[256] = TSADAPTBASIC;
309   PetscBool      set,flg;
310 
311   PetscFunctionBegin;
312   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
313   * function can only be called from inside TSSetFromOptions_GL()  */
314   ierr = PetscOptionsHead("TS Adaptivity options");CHKERRQ(ierr);
315   ierr = PetscOptionsList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList,
316                           ((PetscObject)adapt)->type_name?((PetscObject)adapt)->type_name:type,type,sizeof(type),&flg);CHKERRQ(ierr);
317   if (flg || !((PetscObject)adapt)->type_name) {
318     ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr);
319   }
320   if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(adapt);CHKERRQ(ierr);}
321   ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",adapt->dt_min,&adapt->dt_min,PETSC_NULL);CHKERRQ(ierr);
322   ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",adapt->dt_max,&adapt->dt_max,PETSC_NULL);CHKERRQ(ierr);
323   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);
324   ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
325   if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);}
326   ierr = PetscOptionsTail();CHKERRQ(ierr);
327   PetscFunctionReturn(0);
328 }
329 
330 #undef __FUNCT__
331 #define __FUNCT__ "TSAdaptCandidatesClear"
332 /*@
333    TSAdaptCandidatesClear - clear any previously set candidate schemes
334 
335    Logically Collective
336 
337    Input Argument:
338 .  adapt - adaptive controller
339 
340    Level: developer
341 
342 .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose()
343 @*/
344 PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
345 {
346   PetscErrorCode ierr;
347 
348   PetscFunctionBegin;
349   ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr);
350   PetscFunctionReturn(0);
351 }
352 
353 #undef __FUNCT__
354 #define __FUNCT__ "TSAdaptCandidateAdd"
355 /*@C
356    TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from
357 
358    Logically Collective
359 
360    Input Arguments:
361 +  adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate()
362 .  name - name of the candidate scheme to add
363 .  order - order of the candidate scheme
364 .  stageorder - stage order of the candidate scheme
365 .  ccfl - stability coefficient relative to explicit Euler, used for CFL constraints
366 .  cost - relative measure of the amount of work required for the candidate scheme
367 -  inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme
368 
369    Note:
370    This routine is not available in Fortran.
371 
372    Level: developer
373 
374 .seealso: TSAdaptCandidatesClear(), TSAdaptChoose()
375 @*/
376 PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse)
377 {
378   PetscInt c;
379 
380   PetscFunctionBegin;
381   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
382   if (order < 1) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order);
383   if (inuse) {
384     if (adapt->candidates.inuse_set) SETERRQ(((PetscObject)adapt)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()");
385     adapt->candidates.inuse_set = PETSC_TRUE;
386   }
387   /* first slot if this is the current scheme, otherwise the next available slot */
388   c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n;
389   adapt->candidates.name[c]         = name;
390   adapt->candidates.order[c]        = order;
391   adapt->candidates.stageorder[c]   = stageorder;
392   adapt->candidates.ccfl[c]         = ccfl;
393   adapt->candidates.cost[c]         = cost;
394   adapt->candidates.n++;
395   PetscFunctionReturn(0);
396 }
397 
398 #undef __FUNCT__
399 #define __FUNCT__ "TSAdaptCandidatesGet"
400 /*@C
401    TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost
402 
403    Not Collective
404 
405    Input Arguments:
406 .  adapt - time step adaptivity context
407 
408    Output Arguments:
409 +  n - number of candidate schemes, always at least 1
410 .  order - the order of each candidate scheme
411 .  stageorder - the stage order of each candidate scheme
412 .  ccfl - the CFL coefficient of each scheme
413 -  cost - the relative cost of each scheme
414 
415    Level: developer
416 
417    Note:
418    The current scheme is always returned in the first slot
419 
420 .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose()
421 @*/
422 PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost)
423 {
424   PetscFunctionBegin;
425   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
426   if (n) *n = adapt->candidates.n;
427   if (order) *order = adapt->candidates.order;
428   if (stageorder) *stageorder = adapt->candidates.stageorder;
429   if (ccfl) *ccfl = adapt->candidates.ccfl;
430   if (cost) *cost = adapt->candidates.cost;
431   PetscFunctionReturn(0);
432 }
433 
434 #undef __FUNCT__
435 #define __FUNCT__ "TSAdaptChoose"
436 /*@C
437    TSAdaptChoose - choose which method and step size to use for the next step
438 
439    Logically Collective
440 
441    Input Arguments:
442 +  adapt - adaptive contoller
443 -  h - current step size
444 
445    Output Arguments:
446 +  next_sc - scheme to use for the next step
447 .  next_h - step size to use for the next step
448 -  accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size
449 
450    Note:
451    The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is
452    being retried after an initial rejection.
453 
454    Level: developer
455 
456 .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd()
457 @*/
458 PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept)
459 {
460   PetscErrorCode ierr;
461   PetscReal wlte = -1.0;
462 
463   PetscFunctionBegin;
464   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
465   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
466   PetscValidIntPointer(next_sc,4);
467   PetscValidPointer(next_h,5);
468   PetscValidIntPointer(accept,6);
469   if (adapt->candidates.n < 1) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_WRONGSTATE,"%D candidates have been registered",adapt->candidates.n);
470   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);
471   ierr = (*adapt->ops->choose)(adapt,ts,h,next_sc,next_h,accept,&wlte);CHKERRQ(ierr);
472   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);
473   if (!(*next_h > 0.)) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %G must be positive",*next_h);
474 
475   if (adapt->monitor) {
476     ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
477     if (wlte < 0) {
478       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);
479     } else {
480       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);
481     }
482     ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
483   }
484   PetscFunctionReturn(0);
485 }
486 
487 #undef __FUNCT__
488 #define __FUNCT__ "TSAdaptCheckStage"
489 /*@
490    TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails)
491 
492    Collective
493 
494    Input Arguments:
495 +  adapt - adaptive controller context
496 -  ts - time stepper
497 
498    Output Arguments:
499 .  accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject
500 
501    Level: developer
502 
503 .seealso:
504 @*/
505 PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscBool *accept)
506 {
507   PetscErrorCode      ierr;
508   SNES                snes;
509   SNESConvergedReason snesreason;
510 
511   PetscFunctionBegin;
512   *accept = PETSC_TRUE;
513   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
514   ierr = SNESGetConvergedReason(snes,&snesreason);CHKERRQ(ierr);
515   if (snesreason < 0) {
516     PetscReal dt,new_dt;
517     *accept = PETSC_FALSE;
518     ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
519     if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) {
520       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
521       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);
522       if (adapt->monitor) {
523         ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
524         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);
525         ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
526       }
527     } else {
528       new_dt = dt*adapt->scale_solve_failed;
529       ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr);
530       if (adapt->monitor) {
531         ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
532         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);
533         ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
534       }
535     }
536   }
537   if (adapt->ops->checkstage) {ierr = (*adapt->ops->checkstage)(adapt,ts,accept);CHKERRQ(ierr);}
538   PetscFunctionReturn(0);
539 }
540 
541 
542 
543 #undef __FUNCT__
544 #define __FUNCT__ "TSAdaptCreate"
545 /*@
546   TSAdaptCreate - create an adaptive controller context for time stepping
547 
548   Collective on MPI_Comm
549 
550   Input Parameter:
551 . comm - The communicator
552 
553   Output Parameter:
554 . adapt - new TSAdapt object
555 
556   Level: developer
557 
558   Notes:
559   TSAdapt creation is handled by TS, so users should not need to call this function.
560 
561 .keywords: TSAdapt, create
562 .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy()
563 @*/
564 PetscErrorCode  TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt)
565 {
566   PetscErrorCode ierr;
567   TSAdapt adapt;
568 
569   PetscFunctionBegin;
570   *inadapt = 0;
571   ierr = PetscHeaderCreate(adapt,_p_TSAdapt,struct _TSAdaptOps,TSADAPT_CLASSID,0,"TSAdapt","General Linear adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr);
572 
573   adapt->dt_min             = 1e-20;
574   adapt->dt_max             = 1e50;
575   adapt->scale_solve_failed = 0.25;
576 
577   *inadapt = adapt;
578   PetscFunctionReturn(0);
579 }
580