xref: /petsc/src/ts/adapt/interface/tsadapt.c (revision 5dc64a9727a0cc1147601b21df2b9dcd374fa7c9)
1 
2 #include <petsc/private/tsimpl.h> /*I  "petscts.h" I*/
3 
4 PetscClassId TSADAPT_CLASSID;
5 
6 static PetscFunctionList TSAdaptList;
7 static PetscBool         TSAdaptPackageInitialized;
8 static PetscBool         TSAdaptRegisterAllCalled;
9 
10 PETSC_EXTERN PetscErrorCode TSAdaptCreate_None(TSAdapt);
11 PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt);
12 PETSC_EXTERN PetscErrorCode TSAdaptCreate_DSP(TSAdapt);
13 PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt);
14 PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt);
15 PETSC_EXTERN PetscErrorCode TSAdaptCreate_History(TSAdapt);
16 
17 /*@C
18    TSAdaptRegister -  adds a TSAdapt implementation
19 
20    Not Collective
21 
22    Input Parameters:
23 +  name_scheme - name of user-defined adaptivity scheme
24 -  routine_create - routine to create method context
25 
26    Notes:
27    TSAdaptRegister() may be called multiple times to add several user-defined families.
28 
29    Sample usage:
30 .vb
31    TSAdaptRegister("my_scheme",MySchemeCreate);
32 .ve
33 
34    Then, your scheme can be chosen with the procedural interface via
35 $     TSAdaptSetType(ts,"my_scheme")
36    or at runtime via the option
37 $     -ts_adapt_type my_scheme
38 
39    Level: advanced
40 
41 .keywords: TSAdapt, register
42 
43 .seealso: TSAdaptRegisterAll()
44 @*/
45 PetscErrorCode  TSAdaptRegister(const char sname[],PetscErrorCode (*function)(TSAdapt))
46 {
47   PetscErrorCode ierr;
48 
49   PetscFunctionBegin;
50   ierr = TSAdaptInitializePackage();CHKERRQ(ierr);
51   ierr = PetscFunctionListAdd(&TSAdaptList,sname,function);CHKERRQ(ierr);
52   PetscFunctionReturn(0);
53 }
54 
55 /*@C
56   TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt
57 
58   Not Collective
59 
60   Level: advanced
61 
62 .keywords: TSAdapt, register, all
63 
64 .seealso: TSAdaptRegisterDestroy()
65 @*/
66 PetscErrorCode  TSAdaptRegisterAll(void)
67 {
68   PetscErrorCode ierr;
69 
70   PetscFunctionBegin;
71   if (TSAdaptRegisterAllCalled) PetscFunctionReturn(0);
72   TSAdaptRegisterAllCalled = PETSC_TRUE;
73   ierr = TSAdaptRegister(TSADAPTNONE,   TSAdaptCreate_None);CHKERRQ(ierr);
74   ierr = TSAdaptRegister(TSADAPTBASIC,  TSAdaptCreate_Basic);CHKERRQ(ierr);
75   ierr = TSAdaptRegister(TSADAPTDSP,    TSAdaptCreate_DSP);CHKERRQ(ierr);
76   ierr = TSAdaptRegister(TSADAPTCFL,    TSAdaptCreate_CFL);CHKERRQ(ierr);
77   ierr = TSAdaptRegister(TSADAPTGLEE,   TSAdaptCreate_GLEE);CHKERRQ(ierr);
78   ierr = TSAdaptRegister(TSADAPTHISTORY,TSAdaptCreate_History);CHKERRQ(ierr);
79   PetscFunctionReturn(0);
80 }
81 
82 /*@C
83   TSAdaptFinalizePackage - This function destroys everything in the TS package. It is
84   called from PetscFinalize().
85 
86   Level: developer
87 
88 .keywords: Petsc, destroy, package
89 .seealso: PetscFinalize()
90 @*/
91 PetscErrorCode  TSAdaptFinalizePackage(void)
92 {
93   PetscErrorCode ierr;
94 
95   PetscFunctionBegin;
96   ierr = PetscFunctionListDestroy(&TSAdaptList);CHKERRQ(ierr);
97   TSAdaptPackageInitialized = PETSC_FALSE;
98   TSAdaptRegisterAllCalled  = PETSC_FALSE;
99   PetscFunctionReturn(0);
100 }
101 
102 /*@C
103   TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is
104   called from TSInitializePackage().
105 
106   Level: developer
107 
108 .keywords: TSAdapt, initialize, package
109 .seealso: PetscInitialize()
110 @*/
111 PetscErrorCode  TSAdaptInitializePackage(void)
112 {
113   PetscErrorCode ierr;
114 
115   PetscFunctionBegin;
116   if (TSAdaptPackageInitialized) PetscFunctionReturn(0);
117   TSAdaptPackageInitialized = PETSC_TRUE;
118   ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr);
119   ierr = TSAdaptRegisterAll();CHKERRQ(ierr);
120   ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr);
121   PetscFunctionReturn(0);
122 }
123 
124 /*@C
125   TSAdaptSetType - sets the approach used for the error adapter, currently there is only TSADAPTBASIC and TSADAPTNONE
126 
127   Logicially Collective on TSAdapt
128 
129   Input Parameter:
130 + adapt - the TS adapter, most likely obtained with TSGetAdapt()
131 - type - either  TSADAPTBASIC or TSADAPTNONE
132 
133   Options Database:
134 . -ts_adapt_type <basic or dsp or none> - to set the adapter type
135 
136   Level: intermediate
137 
138 .keywords: TSAdapt, create
139 
140 .seealso: TSGetAdapt(), TSAdaptDestroy(), TSAdaptType, TSAdaptGetType()
141 @*/
142 PetscErrorCode  TSAdaptSetType(TSAdapt adapt,TSAdaptType type)
143 {
144   PetscBool      match;
145   PetscErrorCode ierr,(*r)(TSAdapt);
146 
147   PetscFunctionBegin;
148   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
149   PetscValidCharPointer(type,2);
150   ierr = PetscObjectTypeCompare((PetscObject)adapt,type,&match);CHKERRQ(ierr);
151   if (match) PetscFunctionReturn(0);
152   ierr = PetscFunctionListFind(TSAdaptList,type,&r);CHKERRQ(ierr);
153   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type);
154   if (adapt->ops->destroy) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);}
155   ierr = PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));CHKERRQ(ierr);
156   ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr);
157   ierr = (*r)(adapt);CHKERRQ(ierr);
158   PetscFunctionReturn(0);
159 }
160 
161 /*@C
162   TSAdaptGetType - gets the TS adapter method type (as a string).
163 
164   Not Collective
165 
166   Input Parameter:
167 . adapt - The TS adapter, most likely obtained with TSGetAdapt()
168 
169   Output Parameter:
170 . type - The name of TS adapter method
171 
172   Level: intermediate
173 
174 .keywords: TSAdapt, get, type
175 .seealso TSAdaptSetType()
176 @*/
177 PetscErrorCode TSAdaptGetType(TSAdapt adapt,TSAdaptType *type)
178 {
179   PetscFunctionBegin;
180   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
181   PetscValidPointer(type,2);
182   *type = ((PetscObject)adapt)->type_name;
183   PetscFunctionReturn(0);
184 }
185 
186 PetscErrorCode  TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[])
187 {
188   PetscErrorCode ierr;
189 
190   PetscFunctionBegin;
191   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
192   ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr);
193   PetscFunctionReturn(0);
194 }
195 
196 /*@C
197   TSAdaptLoad - Loads a TSAdapt that has been stored in binary  with TSAdaptView().
198 
199   Collective on PetscViewer
200 
201   Input Parameters:
202 + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or
203            some related function before a call to TSAdaptLoad().
204 - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
205            HDF5 file viewer, obtained from PetscViewerHDF5Open()
206 
207    Level: intermediate
208 
209   Notes:
210    The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored.
211 
212   Notes for advanced users:
213   Most users should not need to know the details of the binary storage
214   format, since TSAdaptLoad() and TSAdaptView() completely hide these details.
215   But for anyone who's interested, the standard binary matrix storage
216   format is
217 .vb
218      has not yet been determined
219 .ve
220 
221 .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad()
222 @*/
223 PetscErrorCode  TSAdaptLoad(TSAdapt adapt,PetscViewer viewer)
224 {
225   PetscErrorCode ierr;
226   PetscBool      isbinary;
227   char           type[256];
228 
229   PetscFunctionBegin;
230   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
231   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
232   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
233   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
234 
235   ierr = PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);CHKERRQ(ierr);
236   ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr);
237   if (adapt->ops->load) {
238     ierr = (*adapt->ops->load)(adapt,viewer);CHKERRQ(ierr);
239   }
240   PetscFunctionReturn(0);
241 }
242 
243 PetscErrorCode  TSAdaptView(TSAdapt adapt,PetscViewer viewer)
244 {
245   PetscErrorCode ierr;
246   PetscBool      iascii,isbinary,isnone,isglee;
247 
248   PetscFunctionBegin;
249   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
250   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt),&viewer);CHKERRQ(ierr);}
251   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
252   PetscCheckSameComm(adapt,1,viewer,2);
253   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
254   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
255   if (iascii) {
256     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr);
257     ierr = PetscObjectTypeCompare((PetscObject)adapt,TSADAPTNONE,&isnone);CHKERRQ(ierr);
258     ierr = PetscObjectTypeCompare((PetscObject)adapt,TSADAPTGLEE,&isglee);CHKERRQ(ierr);
259     if (!isnone) {
260       if (adapt->always_accept) {ierr = PetscViewerASCIIPrintf(viewer,"  always accepting steps\n");CHKERRQ(ierr);}
261       ierr = PetscViewerASCIIPrintf(viewer,"  safety factor %g\n",(double)adapt->safety);CHKERRQ(ierr);
262       ierr = PetscViewerASCIIPrintf(viewer,"  extra safety factor after step rejection %g\n",(double)adapt->reject_safety);CHKERRQ(ierr);
263       ierr = PetscViewerASCIIPrintf(viewer,"  clip fastest increase %g\n",(double)adapt->clip[1]);CHKERRQ(ierr);
264       ierr = PetscViewerASCIIPrintf(viewer,"  clip fastest decrease %g\n",(double)adapt->clip[0]);CHKERRQ(ierr);
265       ierr = PetscViewerASCIIPrintf(viewer,"  maximum allowed timestep %g\n",(double)adapt->dt_max);CHKERRQ(ierr);
266       ierr = PetscViewerASCIIPrintf(viewer,"  minimum allowed timestep %g\n",(double)adapt->dt_min);CHKERRQ(ierr);
267       ierr = PetscViewerASCIIPrintf(viewer,"  maximum solution absolute value to be ignored %g\n",(double)adapt->ignore_max);CHKERRQ(ierr);
268     }
269     if (isglee) {
270       if (adapt->glee_use_local) {
271         ierr = PetscViewerASCIIPrintf(viewer,"  GLEE uses local error control\n");CHKERRQ(ierr);
272       } else {
273         ierr = PetscViewerASCIIPrintf(viewer,"  GLEE uses global error control\n");CHKERRQ(ierr);
274       }
275     }
276     if (adapt->ops->view) {
277       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
278       ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr);
279       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
280     }
281   } else if (isbinary) {
282     char type[256];
283 
284     /* need to save FILE_CLASS_ID for adapt class */
285     ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr);
286     ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr);
287   } else if (adapt->ops->view) {
288     ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr);
289   }
290   PetscFunctionReturn(0);
291 }
292 
293 /*@
294    TSAdaptReset - Resets a TSAdapt context.
295 
296    Collective on TS
297 
298    Input Parameter:
299 .  adapt - the TSAdapt context obtained from TSAdaptCreate()
300 
301    Level: developer
302 
303 .seealso: TSAdaptCreate(), TSAdaptDestroy()
304 @*/
305 PetscErrorCode  TSAdaptReset(TSAdapt adapt)
306 {
307   PetscErrorCode ierr;
308 
309   PetscFunctionBegin;
310   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
311   if (adapt->ops->reset) {ierr = (*adapt->ops->reset)(adapt);CHKERRQ(ierr);}
312   PetscFunctionReturn(0);
313 }
314 
315 PetscErrorCode  TSAdaptDestroy(TSAdapt *adapt)
316 {
317   PetscErrorCode ierr;
318 
319   PetscFunctionBegin;
320   if (!*adapt) PetscFunctionReturn(0);
321   PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1);
322   if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);}
323 
324   ierr = TSAdaptReset(*adapt);CHKERRQ(ierr);
325 
326   if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);}
327   ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr);
328   ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr);
329   PetscFunctionReturn(0);
330 }
331 
332 /*@
333    TSAdaptSetMonitor - Monitor the choices made by the adaptive controller
334 
335    Collective on TSAdapt
336 
337    Input Arguments:
338 +  adapt - adaptive controller context
339 -  flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable
340 
341    Options Database Keys:
342 .  -ts_adapt_monitor - to turn on monitoring
343 
344    Level: intermediate
345 
346 .seealso: TSAdaptChoose()
347 @*/
348 PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg)
349 {
350   PetscErrorCode ierr;
351 
352   PetscFunctionBegin;
353   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
354   PetscValidLogicalCollectiveBool(adapt,flg,2);
355   if (flg) {
356     if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);}
357   } else {
358     ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr);
359   }
360   PetscFunctionReturn(0);
361 }
362 
363 /*@C
364    TSAdaptSetCheckStage - Set a callback to check convergence for a stage
365 
366    Logically collective on TSAdapt
367 
368    Input Arguments:
369 +  adapt - adaptive controller context
370 -  func - stage check function
371 
372    Arguments of func:
373 $  PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept)
374 
375 +  adapt - adaptive controller context
376 .  ts - time stepping context
377 -  accept - pending choice of whether to accept, can be modified by this routine
378 
379    Level: advanced
380 
381 .seealso: TSAdaptChoose()
382 @*/
383 PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscReal,Vec,PetscBool*))
384 {
385 
386   PetscFunctionBegin;
387   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
388   adapt->checkstage = func;
389   PetscFunctionReturn(0);
390 }
391 
392 /*@
393    TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of
394    any error or stability condition not meeting the prescribed goal.
395 
396    Logically collective on TSAdapt
397 
398    Input Arguments:
399 +  adapt - time step adaptivity context, usually gotten with TSGetAdapt()
400 -  flag - whether to always accept steps
401 
402    Options Database Keys:
403 .  -ts_adapt_always_accept - to always accept steps
404 
405    Level: intermediate
406 
407 .seealso: TSAdapt, TSAdaptChoose()
408 @*/
409 PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt,PetscBool flag)
410 {
411   PetscFunctionBegin;
412   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
413   PetscValidLogicalCollectiveBool(adapt,flag,2);
414   adapt->always_accept = flag;
415   PetscFunctionReturn(0);
416 }
417 
418 /*@
419    TSAdaptSetSafety - Set safety factors
420 
421    Logically collective on TSAdapt
422 
423    Input Arguments:
424 +  adapt - adaptive controller context
425 .  safety - safety factor relative to target error/stability goal
426 -  reject_safety - extra safety factor to apply if the last step was rejected
427 
428    Options Database Keys:
429 +  -ts_adapt_safety <safety> - to set safety factor
430 -  -ts_adapt_reject_safety <reject_safety> - to set reject safety factor
431 
432    Level: intermediate
433 
434 .seealso: TSAdapt, TSAdaptGetSafety(), TSAdaptChoose()
435 @*/
436 PetscErrorCode TSAdaptSetSafety(TSAdapt adapt,PetscReal safety,PetscReal reject_safety)
437 {
438   PetscFunctionBegin;
439   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
440   PetscValidLogicalCollectiveReal(adapt,safety,2);
441   PetscValidLogicalCollectiveReal(adapt,reject_safety,3);
442   if (safety != PETSC_DEFAULT && safety < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Safety factor %g must be non negative",(double)safety);
443   if (safety != PETSC_DEFAULT && safety > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Safety factor %g must be less than one",(double)safety);
444   if (reject_safety != PETSC_DEFAULT && reject_safety < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reject safety factor %g must be non negative",(double)reject_safety);
445   if (reject_safety != PETSC_DEFAULT && reject_safety > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reject safety factor %g must be less than one",(double)reject_safety);
446   if (safety != PETSC_DEFAULT) adapt->safety = safety;
447   if (reject_safety != PETSC_DEFAULT) adapt->reject_safety = reject_safety;
448   PetscFunctionReturn(0);
449 }
450 
451 /*@
452    TSAdaptGetSafety - Get safety factors
453 
454    Not Collective
455 
456    Input Arguments:
457 .  adapt - adaptive controller context
458 
459    Ouput Arguments:
460 .  safety - safety factor relative to target error/stability goal
461 +  reject_safety - extra safety factor to apply if the last step was rejected
462 
463    Level: intermediate
464 
465 .seealso: TSAdapt, TSAdaptSetSafety(), TSAdaptChoose()
466 @*/
467 PetscErrorCode TSAdaptGetSafety(TSAdapt adapt,PetscReal *safety,PetscReal *reject_safety)
468 {
469   PetscFunctionBegin;
470   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
471   if (safety)        PetscValidRealPointer(safety,2);
472   if (reject_safety) PetscValidRealPointer(reject_safety,3);
473   if (safety)        *safety        = adapt->safety;
474   if (reject_safety) *reject_safety = adapt->reject_safety;
475   PetscFunctionReturn(0);
476 }
477 
478 /*@
479    TSAdaptSetMaxIgnore - Set error estimation threshold. Solution components below this threshold value will not be considered when computing error norms for time step adaptivity (in absolute value). A negative value (default) of the threshold leads to considering all solution components.
480 
481    Logically collective on TSAdapt
482 
483    Input Arguments:
484 +  adapt - adaptive controller context
485 -  max_ignore - threshold for solution components that are ignored during error estimation
486 
487    Options Database Keys:
488 .  -ts_adapt_max_ignore <max_ignore> - to set the threshold
489 
490    Level: intermediate
491 
492 .seealso: TSAdapt, TSAdaptGetMaxIgnore(), TSAdaptChoose()
493 @*/
494 PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt adapt,PetscReal max_ignore)
495 {
496   PetscFunctionBegin;
497   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
498   PetscValidLogicalCollectiveReal(adapt,max_ignore,2);
499   adapt->ignore_max = max_ignore;
500   PetscFunctionReturn(0);
501 }
502 
503 /*@
504    TSAdaptGetMaxIgnore - Get error estimation threshold. Solution components below this threshold value will not be considered when computing error norms for time step adaptivity (in absolute value).
505 
506    Not Collective
507 
508    Input Arguments:
509 .  adapt - adaptive controller context
510 
511    Ouput Arguments:
512 .  max_ignore - threshold for solution components that are ignored during error estimation
513 
514    Level: intermediate
515 
516 .seealso: TSAdapt, TSAdaptSetMaxIgnore(), TSAdaptChoose()
517 @*/
518 PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt adapt,PetscReal *max_ignore)
519 {
520   PetscFunctionBegin;
521   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
522   PetscValidRealPointer(max_ignore,2);
523   *max_ignore = adapt->ignore_max;
524   PetscFunctionReturn(0);
525 }
526 
527 
528 /*@
529    TSAdaptSetClip - Sets the admissible decrease/increase factor in step size
530 
531    Logically collective on TSAdapt
532 
533    Input Arguments:
534 +  adapt - adaptive controller context
535 .  low - admissible decrease factor
536 -  high - admissible increase factor
537 
538    Options Database Keys:
539 .  -ts_adapt_clip <low>,<high> - to set admissible time step decrease and increase factors
540 
541    Level: intermediate
542 
543 .seealso: TSAdaptChoose(), TSAdaptGetClip()
544 @*/
545 PetscErrorCode TSAdaptSetClip(TSAdapt adapt,PetscReal low,PetscReal high)
546 {
547   PetscFunctionBegin;
548   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
549   PetscValidLogicalCollectiveReal(adapt,low,2);
550   PetscValidLogicalCollectiveReal(adapt,high,3);
551   if (low  != PETSC_DEFAULT && low  < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be non negative",(double)low);
552   if (low  != PETSC_DEFAULT && low  > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be less than one",(double)low);
553   if (high != PETSC_DEFAULT && high < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Increase factor %g must be geather than one",(double)high);
554   if (low  != PETSC_DEFAULT) adapt->clip[0] = low;
555   if (high != PETSC_DEFAULT) adapt->clip[1] = high;
556   PetscFunctionReturn(0);
557 }
558 
559 /*@
560    TSAdaptGetClip - Gets the admissible decrease/increase factor in step size
561 
562    Not Collective
563 
564    Input Arguments:
565 .  adapt - adaptive controller context
566 
567    Ouput Arguments:
568 +  low - optional, admissible decrease factor
569 -  high - optional, admissible increase factor
570 
571    Level: intermediate
572 
573 .seealso: TSAdaptChoose(), TSAdaptSetClip()
574 @*/
575 PetscErrorCode TSAdaptGetClip(TSAdapt adapt,PetscReal *low,PetscReal *high)
576 {
577   PetscFunctionBegin;
578   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
579   if (low)  PetscValidRealPointer(low,2);
580   if (high) PetscValidRealPointer(high,3);
581   if (low)  *low  = adapt->clip[0];
582   if (high) *high = adapt->clip[1];
583   PetscFunctionReturn(0);
584 }
585 
586 /*@
587    TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the controller
588 
589    Logically collective on TSAdapt
590 
591    Input Arguments:
592 +  adapt - time step adaptivity context, usually gotten with TSGetAdapt()
593 .  hmin - minimum time step
594 -  hmax - maximum time step
595 
596    Options Database Keys:
597 +  -ts_adapt_dt_min <min> - to set minimum time step
598 -  -ts_adapt_dt_max <max> - to set maximum time step
599 
600    Level: intermediate
601 
602 .seealso: TSAdapt, TSAdaptGetStepLimits(), TSAdaptChoose()
603 @*/
604 PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax)
605 {
606 
607   PetscFunctionBegin;
608   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
609   PetscValidLogicalCollectiveReal(adapt,hmin,2);
610   PetscValidLogicalCollectiveReal(adapt,hmax,3);
611   if (hmin != PETSC_DEFAULT && hmin < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmin);
612   if (hmax != PETSC_DEFAULT && hmax < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmax);
613   if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin;
614   if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax;
615   hmin = adapt->dt_min;
616   hmax = adapt->dt_max;
617   if (hmax <= hmin) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Maximum time step %g must geather than minimum time step %g",(double)hmax,(double)hmin);
618   PetscFunctionReturn(0);
619 }
620 
621 /*@
622    TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the controller
623 
624    Not Collective
625 
626    Input Arguments:
627 .  adapt - time step adaptivity context, usually gotten with TSGetAdapt()
628 
629    Output Arguments:
630 +  hmin - minimum time step
631 -  hmax - maximum time step
632 
633    Level: intermediate
634 
635 .seealso: TSAdapt, TSAdaptSetStepLimits(), TSAdaptChoose()
636 @*/
637 PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt,PetscReal *hmin,PetscReal *hmax)
638 {
639 
640   PetscFunctionBegin;
641   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
642   if (hmin) PetscValidRealPointer(hmin,2);
643   if (hmax) PetscValidRealPointer(hmax,3);
644   if (hmin) *hmin = adapt->dt_min;
645   if (hmax) *hmax = adapt->dt_max;
646   PetscFunctionReturn(0);
647 }
648 
649 /*
650    TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options.
651 
652    Collective on TSAdapt
653 
654    Input Parameter:
655 .  adapt - the TSAdapt context
656 
657    Options Database Keys:
658 +  -ts_adapt_type <type> - algorithm to use for adaptivity
659 .  -ts_adapt_always_accept - always accept steps regardless of error/stability goals
660 .  -ts_adapt_safety <safety> - safety factor relative to target error/stability goal
661 .  -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected
662 .  -ts_adapt_clip <low,high> - admissible time step decrease and increase factors
663 .  -ts_adapt_dt_min <min> - minimum timestep to use
664 .  -ts_adapt_dt_max <max> - maximum timestep to use
665 .  -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails
666 .  -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates
667 -  -ts_adapt_time_step_increase_delay - number of timesteps to delay increasing the time step after it has been decreased due to failed solver
668 
669    Level: advanced
670 
671    Notes:
672    This function is automatically called by TSSetFromOptions()
673 
674 .keywords: TS, TSGetAdapt(), TSAdaptSetType(), TSAdaptSetStepLimits()
675 
676 .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptSetAlwaysAccept(), TSAdaptSetSafety(),
677           TSAdaptSetClip(), TSAdaptSetStepLimits(), TSAdaptSetMonitor()
678 */
679 PetscErrorCode  TSAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSAdapt adapt)
680 {
681   PetscErrorCode ierr;
682   char           type[256] = TSADAPTBASIC;
683   PetscReal      safety,reject_safety,clip[2],hmin,hmax;
684   PetscBool      set,flg;
685   PetscInt       two;
686 
687   PetscFunctionBegin;
688   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
689   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
690    * function can only be called from inside TSSetFromOptions()  */
691   ierr = PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");CHKERRQ(ierr);
692   ierr = PetscOptionsFList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList,((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);CHKERRQ(ierr);
693   if (flg || !((PetscObject)adapt)->type_name) {
694     ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr);
695   }
696 
697   ierr = PetscOptionsBool("-ts_adapt_always_accept","Always accept the step","TSAdaptSetAlwaysAccept",adapt->always_accept,&flg,&set);CHKERRQ(ierr);
698   if (set) {ierr = TSAdaptSetAlwaysAccept(adapt,flg);CHKERRQ(ierr);}
699 
700   safety = adapt->safety; reject_safety = adapt->reject_safety;
701   ierr = PetscOptionsReal("-ts_adapt_safety","Safety factor relative to target error/stability goal","TSAdaptSetSafety",safety,&safety,&set);CHKERRQ(ierr);
702   ierr = PetscOptionsReal("-ts_adapt_reject_safety","Extra safety factor to apply if the last step was rejected","TSAdaptSetSafety",reject_safety,&reject_safety,&flg);CHKERRQ(ierr);
703   if (set || flg) {ierr = TSAdaptSetSafety(adapt,safety,reject_safety);CHKERRQ(ierr);}
704 
705   two = 2; clip[0] = adapt->clip[0]; clip[1] = adapt->clip[1];
706   ierr = PetscOptionsRealArray("-ts_adapt_clip","Admissible decrease/increase factor in step size","TSAdaptSetClip",clip,&two,&set);CHKERRQ(ierr);
707   if (set && (two != 2)) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_clip");
708   if (set) {ierr = TSAdaptSetClip(adapt,clip[0],clip[1]);CHKERRQ(ierr);}
709 
710   hmin = adapt->dt_min; hmax = adapt->dt_max;
711   ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",hmin,&hmin,&set);CHKERRQ(ierr);
712   ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",hmax,&hmax,&flg);CHKERRQ(ierr);
713   if (set || flg) {ierr = TSAdaptSetStepLimits(adapt,hmin,hmax);CHKERRQ(ierr);}
714 
715   ierr = PetscOptionsReal("-ts_adapt_max_ignore","Adaptor ignores (absolute) solution values smaller than this value","",adapt->ignore_max,&adapt->ignore_max,&set);CHKERRQ(ierr);
716   ierr = PetscOptionsBool("-ts_adapt_glee_use_local","GLEE adaptor uses local error estimation for step control","",adapt->glee_use_local,&adapt->glee_use_local,&set);CHKERRQ(ierr);
717 
718   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);
719 
720   ierr = PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);CHKERRQ(ierr);
721   if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported");
722 
723   ierr = PetscOptionsInt("-ts_adapt_time_step_increase_delay","Number of timesteps to delay increasing the time step after it has been decreased due to failed solver","TSAdaptSetTimeStepIncreaseDelay",adapt->timestepjustdecreased_delay,&adapt->timestepjustdecreased_delay,NULL);CHKERRQ(ierr);
724 
725   ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
726   if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);}
727 
728   if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);}
729   ierr = PetscOptionsTail();CHKERRQ(ierr);
730   PetscFunctionReturn(0);
731 }
732 
733 /*@
734    TSAdaptCandidatesClear - clear any previously set candidate schemes
735 
736    Logically collective on TSAdapt
737 
738    Input Argument:
739 .  adapt - adaptive controller
740 
741    Level: developer
742 
743 .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose()
744 @*/
745 PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
746 {
747   PetscErrorCode ierr;
748 
749   PetscFunctionBegin;
750   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
751   ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr);
752   PetscFunctionReturn(0);
753 }
754 
755 /*@C
756    TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from
757 
758    Logically collective on TSAdapt
759 
760    Input Arguments:
761 +  adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate()
762 .  name - name of the candidate scheme to add
763 .  order - order of the candidate scheme
764 .  stageorder - stage order of the candidate scheme
765 .  ccfl - stability coefficient relative to explicit Euler, used for CFL constraints
766 .  cost - relative measure of the amount of work required for the candidate scheme
767 -  inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme
768 
769    Note:
770    This routine is not available in Fortran.
771 
772    Level: developer
773 
774 .seealso: TSAdaptCandidatesClear(), TSAdaptChoose()
775 @*/
776 PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse)
777 {
778   PetscInt c;
779 
780   PetscFunctionBegin;
781   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
782   if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order);
783   if (inuse) {
784     if (adapt->candidates.inuse_set) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()");
785     adapt->candidates.inuse_set = PETSC_TRUE;
786   }
787   /* first slot if this is the current scheme, otherwise the next available slot */
788   c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n;
789 
790   adapt->candidates.name[c]       = name;
791   adapt->candidates.order[c]      = order;
792   adapt->candidates.stageorder[c] = stageorder;
793   adapt->candidates.ccfl[c]       = ccfl;
794   adapt->candidates.cost[c]       = cost;
795   adapt->candidates.n++;
796   PetscFunctionReturn(0);
797 }
798 
799 /*@C
800    TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost
801 
802    Not Collective
803 
804    Input Arguments:
805 .  adapt - time step adaptivity context
806 
807    Output Arguments:
808 +  n - number of candidate schemes, always at least 1
809 .  order - the order of each candidate scheme
810 .  stageorder - the stage order of each candidate scheme
811 .  ccfl - the CFL coefficient of each scheme
812 -  cost - the relative cost of each scheme
813 
814    Level: developer
815 
816    Note:
817    The current scheme is always returned in the first slot
818 
819 .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose()
820 @*/
821 PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost)
822 {
823   PetscFunctionBegin;
824   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
825   if (n) *n = adapt->candidates.n;
826   if (order) *order = adapt->candidates.order;
827   if (stageorder) *stageorder = adapt->candidates.stageorder;
828   if (ccfl) *ccfl = adapt->candidates.ccfl;
829   if (cost) *cost = adapt->candidates.cost;
830   PetscFunctionReturn(0);
831 }
832 
833 /*@C
834    TSAdaptChoose - choose which method and step size to use for the next step
835 
836    Collective on TSAdapt
837 
838    Input Arguments:
839 +  adapt - adaptive contoller
840 -  h - current step size
841 
842    Output Arguments:
843 +  next_sc - optional, scheme to use for the next step
844 .  next_h - step size to use for the next step
845 -  accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size
846 
847    Note:
848    The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is
849    being retried after an initial rejection.
850 
851    Level: developer
852 
853 .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd()
854 @*/
855 PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept)
856 {
857   PetscErrorCode ierr;
858   PetscInt       ncandidates = adapt->candidates.n;
859   PetscInt       scheme = 0;
860   PetscReal      wlte = -1.0;
861   PetscReal      wltea = -1.0;
862   PetscReal      wlter = -1.0;
863 
864   PetscFunctionBegin;
865   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
866   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
867   if (next_sc) PetscValidIntPointer(next_sc,4);
868   PetscValidPointer(next_h,5);
869   PetscValidIntPointer(accept,6);
870   if (next_sc) *next_sc = 0;
871 
872   /* Do not mess with adaptivity while handling events*/
873   if (ts->event && ts->event->status != TSEVENT_NONE) {
874     *next_h = h;
875     *accept = PETSC_TRUE;
876     PetscFunctionReturn(0);
877   }
878 
879   ierr = (*adapt->ops->choose)(adapt,ts,h,&scheme,next_h,accept,&wlte,&wltea,&wlter);CHKERRQ(ierr);
880   if (scheme < 0 || (ncandidates > 0 && scheme >= ncandidates)) SETERRQ2(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Chosen scheme %D not in valid range 0..%D",scheme,ncandidates-1);
881   if (*next_h < 0) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h);
882   if (next_sc) *next_sc = scheme;
883 
884   if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
885     /* Increase/reduce step size if end time of next step is close to or overshoots max time */
886     PetscReal t = ts->ptime + ts->time_step, h = *next_h;
887     PetscReal tend = t + h, tmax = ts->max_time, hmax = tmax - t;
888     PetscReal a = (PetscReal)(1.0 + adapt->matchstepfac[0]);
889     PetscReal b = adapt->matchstepfac[1];
890     if (t < tmax && tend > tmax) *next_h = hmax;
891     if (t < tmax && tend < tmax && h*b > hmax) *next_h = hmax/2;
892     if (t < tmax && tend < tmax && h*a > hmax) *next_h = hmax;
893   }
894 
895   if (adapt->monitor) {
896     const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : "";
897     ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
898     if (wlte < 0) {
899       ierr = PetscViewerASCIIPrintf(adapt->monitor,"    TSAdapt %s %s %D:%s step %3D %s t=%-11g+%10.3e dt=%-10.3e\n",((PetscObject)adapt)->type_name,((PetscObject)ts)->type_name,scheme,sc_name,ts->steps,*accept ? "accepted" : "rejected",(double)ts->ptime,(double)h,(double)*next_h);CHKERRQ(ierr);
900     } else {
901       ierr = PetscViewerASCIIPrintf(adapt->monitor,"    TSAdapt %s %s %D:%s step %3D %s t=%-11g+%10.3e dt=%-10.3e wlte=%5.3g  wltea=%5.3g wlter=%5.3g\n",((PetscObject)adapt)->type_name,((PetscObject)ts)->type_name,scheme,sc_name,ts->steps,*accept ? "accepted" : "rejected",(double)ts->ptime,(double)h,(double)*next_h,(double)wlte,(double)wltea,(double)wlter);CHKERRQ(ierr);
902     }
903     ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
904   }
905   PetscFunctionReturn(0);
906 }
907 
908 /*@
909    TSAdaptSetTimeStepIncreaseDelay - The number of timesteps to wait after a decrease in the timestep due to failed solver
910                                      before increasing the time step.
911 
912    Logicially Collective on TSAdapt
913 
914    Input Arguments:
915 +  adapt - adaptive controller context
916 -  cnt - the number of timesteps
917 
918    Options Database Key:
919 .  -ts_adapt_time_step_increase_delay cnt - number of steps to delay the increase
920 
921    Notes: This is to prevent an adaptor from bouncing back and forth between two nearby timesteps. The default is 0.
922           The successful use of this option is problem dependent
923 
924    Developer Note: there is no theory to support this option
925 
926    Level: advanced
927 
928 .seealso:
929 @*/
930 PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt,PetscInt cnt)
931 {
932   PetscFunctionBegin;
933   adapt->timestepjustdecreased_delay = cnt;
934   PetscFunctionReturn(0);
935 }
936 
937 
938 /*@
939    TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails)
940 
941    Collective on TSAdapt
942 
943    Input Arguments:
944 +  adapt - adaptive controller context
945 .  ts - time stepper
946 .  t - Current simulation time
947 -  Y - Current solution vector
948 
949    Output Arguments:
950 .  accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject
951 
952    Level: developer
953 
954 .seealso:
955 @*/
956 PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool *accept)
957 {
958   PetscErrorCode      ierr;
959   SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING;
960 
961   PetscFunctionBegin;
962   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
963   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
964   PetscValidIntPointer(accept,3);
965 
966   if (ts->snes) {ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);}
967   if (snesreason < 0) {
968     *accept = PETSC_FALSE;
969     if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) {
970       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
971       ierr = PetscInfo2(ts,"Step=%D, nonlinear solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);CHKERRQ(ierr);
972       if (adapt->monitor) {
973         ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
974         ierr = PetscViewerASCIIPrintf(adapt->monitor,"    TSAdapt %s step %3D stage rejected t=%-11g+%10.3e, nonlinear solve failures %D greater than current TS allowed\n",((PetscObject)adapt)->type_name,ts->steps,(double)ts->ptime,(double)ts->time_step,ts->num_snes_failures);CHKERRQ(ierr);
975         ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
976       }
977     }
978   } else {
979     *accept = PETSC_TRUE;
980     ierr = TSFunctionDomainError(ts,t,Y,accept);CHKERRQ(ierr);
981     if(*accept && adapt->checkstage) {
982       ierr = (*adapt->checkstage)(adapt,ts,t,Y,accept);CHKERRQ(ierr);
983     }
984   }
985 
986   if(!(*accept) && !ts->reason) {
987     PetscReal dt,new_dt;
988     ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
989     new_dt = dt * adapt->scale_solve_failed;
990     ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr);
991     adapt->timestepjustdecreased += adapt->timestepjustdecreased_delay;
992     if (adapt->monitor) {
993       ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
994       ierr = PetscViewerASCIIPrintf(adapt->monitor,"    TSAdapt %s step %3D stage rejected (%s) t=%-11g+%10.3e retrying with dt=%-10.3e\n",((PetscObject)adapt)->type_name,ts->steps,SNESConvergedReasons[snesreason],(double)ts->ptime,(double)dt,(double)new_dt);CHKERRQ(ierr);
995       ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
996     }
997   }
998   PetscFunctionReturn(0);
999 }
1000 
1001 /*@
1002   TSAdaptCreate - create an adaptive controller context for time stepping
1003 
1004   Collective on MPI_Comm
1005 
1006   Input Parameter:
1007 . comm - The communicator
1008 
1009   Output Parameter:
1010 . adapt - new TSAdapt object
1011 
1012   Level: developer
1013 
1014   Notes:
1015   TSAdapt creation is handled by TS, so users should not need to call this function.
1016 
1017 .keywords: TSAdapt, create
1018 .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy()
1019 @*/
1020 PetscErrorCode  TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt)
1021 {
1022   PetscErrorCode ierr;
1023   TSAdapt        adapt;
1024 
1025   PetscFunctionBegin;
1026   PetscValidPointer(inadapt,1);
1027   *inadapt = NULL;
1028   ierr = TSAdaptInitializePackage();CHKERRQ(ierr);
1029 
1030   ierr = PetscHeaderCreate(adapt,TSADAPT_CLASSID,"TSAdapt","Time stepping adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr);
1031 
1032   adapt->always_accept      = PETSC_FALSE;
1033   adapt->safety             = 0.9;
1034   adapt->reject_safety      = 0.5;
1035   adapt->clip[0]            = 0.1;
1036   adapt->clip[1]            = 10.;
1037   adapt->dt_min             = 1e-20;
1038   adapt->dt_max             = 1e+20;
1039   adapt->ignore_max         = -1.0;
1040   adapt->glee_use_local     = PETSC_TRUE;
1041   adapt->scale_solve_failed = 0.25;
1042   /* these two safety factors are not public, and they are used only in the TS_EXACTFINALTIME_MATCHSTEP case
1043      to prevent from situations were unreasonably small time steps are taken in order to match the final time */
1044   adapt->matchstepfac[0]    = 0.01; /* allow 1% step size increase in the last step */
1045   adapt->matchstepfac[1]    = 2.0;  /* halve last step if it is greater than what remains divided this factor */
1046   adapt->wnormtype          = NORM_2;
1047   adapt->timestepjustdecreased_delay = 0;
1048 
1049   *inadapt = adapt;
1050   PetscFunctionReturn(0);
1051 }
1052