xref: /petsc/src/ts/adapt/interface/tsadapt.c (revision 1c167fc2a0c466620c6784b5ffb486042a891a3a)
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    TSAdaptSetClip - Sets the admissible decrease/increase factor in step size
480 
481    Logically collective on TSAdapt
482 
483    Input Arguments:
484 +  adapt - adaptive controller context
485 .  low - admissible decrease factor
486 -  high - admissible increase factor
487 
488    Options Database Keys:
489 .  -ts_adapt_clip <low>,<high> - to set admissible time step decrease and increase factors
490 
491    Level: intermediate
492 
493 .seealso: TSAdaptChoose(), TSAdaptGetClip()
494 @*/
495 PetscErrorCode TSAdaptSetClip(TSAdapt adapt,PetscReal low,PetscReal high)
496 {
497   PetscFunctionBegin;
498   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
499   PetscValidLogicalCollectiveReal(adapt,low,2);
500   PetscValidLogicalCollectiveReal(adapt,high,3);
501   if (low  != PETSC_DEFAULT && low  < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be non negative",(double)low);
502   if (low  != PETSC_DEFAULT && low  > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be less than one",(double)low);
503   if (high != PETSC_DEFAULT && high < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Increase factor %g must be geather than one",(double)high);
504   if (low  != PETSC_DEFAULT) adapt->clip[0] = low;
505   if (high != PETSC_DEFAULT) adapt->clip[1] = high;
506   PetscFunctionReturn(0);
507 }
508 
509 /*@
510    TSAdaptGetClip - Gets the admissible decrease/increase factor in step size
511 
512    Not Collective
513 
514    Input Arguments:
515 .  adapt - adaptive controller context
516 
517    Ouput Arguments:
518 +  low - optional, admissible decrease factor
519 -  high - optional, admissible increase factor
520 
521    Level: intermediate
522 
523 .seealso: TSAdaptChoose(), TSAdaptSetClip()
524 @*/
525 PetscErrorCode TSAdaptGetClip(TSAdapt adapt,PetscReal *low,PetscReal *high)
526 {
527   PetscFunctionBegin;
528   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
529   if (low)  PetscValidRealPointer(low,2);
530   if (high) PetscValidRealPointer(high,3);
531   if (low)  *low  = adapt->clip[0];
532   if (high) *high = adapt->clip[1];
533   PetscFunctionReturn(0);
534 }
535 
536 /*@
537    TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the controller
538 
539    Logically collective on TSAdapt
540 
541    Input Arguments:
542 +  adapt - time step adaptivity context, usually gotten with TSGetAdapt()
543 .  hmin - minimum time step
544 -  hmax - maximum time step
545 
546    Options Database Keys:
547 +  -ts_adapt_dt_min <min> - to set minimum time step
548 -  -ts_adapt_dt_max <max> - to set maximum time step
549 
550    Level: intermediate
551 
552 .seealso: TSAdapt, TSAdaptGetStepLimits(), TSAdaptChoose()
553 @*/
554 PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax)
555 {
556 
557   PetscFunctionBegin;
558   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
559   PetscValidLogicalCollectiveReal(adapt,hmin,2);
560   PetscValidLogicalCollectiveReal(adapt,hmax,3);
561   if (hmin != PETSC_DEFAULT && hmin < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmin);
562   if (hmax != PETSC_DEFAULT && hmax < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmax);
563   if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin;
564   if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax;
565   hmin = adapt->dt_min;
566   hmax = adapt->dt_max;
567   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);
568   PetscFunctionReturn(0);
569 }
570 
571 /*@
572    TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the controller
573 
574    Not Collective
575 
576    Input Arguments:
577 .  adapt - time step adaptivity context, usually gotten with TSGetAdapt()
578 
579    Output Arguments:
580 +  hmin - minimum time step
581 -  hmax - maximum time step
582 
583    Level: intermediate
584 
585 .seealso: TSAdapt, TSAdaptSetStepLimits(), TSAdaptChoose()
586 @*/
587 PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt,PetscReal *hmin,PetscReal *hmax)
588 {
589 
590   PetscFunctionBegin;
591   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
592   if (hmin) PetscValidRealPointer(hmin,2);
593   if (hmax) PetscValidRealPointer(hmax,3);
594   if (hmin) *hmin = adapt->dt_min;
595   if (hmax) *hmax = adapt->dt_max;
596   PetscFunctionReturn(0);
597 }
598 
599 /*
600    TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options.
601 
602    Collective on TSAdapt
603 
604    Input Parameter:
605 .  adapt - the TSAdapt context
606 
607    Options Database Keys:
608 +  -ts_adapt_type <type> - algorithm to use for adaptivity
609 .  -ts_adapt_always_accept - always accept steps regardless of error/stability goals
610 .  -ts_adapt_safety <safety> - safety factor relative to target error/stability goal
611 .  -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected
612 .  -ts_adapt_clip <low,high> - admissible time step decrease and increase factors
613 .  -ts_adapt_dt_min <min> - minimum timestep to use
614 .  -ts_adapt_dt_max <max> - maximum timestep to use
615 .  -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails
616 .  -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates
617 -  -ts_adapt_time_step_increase_delay - number of timesteps to delay increasing the time step after it has been decreased due to failed solver
618 
619    Level: advanced
620 
621    Notes:
622    This function is automatically called by TSSetFromOptions()
623 
624 .keywords: TS, TSGetAdapt(), TSAdaptSetType(), TSAdaptSetStepLimits()
625 
626 .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptSetAlwaysAccept(), TSAdaptSetSafety(),
627           TSAdaptSetClip(), TSAdaptSetStepLimits(), TSAdaptSetMonitor()
628 */
629 PetscErrorCode  TSAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSAdapt adapt)
630 {
631   PetscErrorCode ierr;
632   char           type[256] = TSADAPTBASIC;
633   PetscReal      safety,reject_safety,clip[2],hmin,hmax;
634   PetscBool      set,flg;
635   PetscInt       two;
636 
637   PetscFunctionBegin;
638   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
639   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
640    * function can only be called from inside TSSetFromOptions()  */
641   ierr = PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");CHKERRQ(ierr);
642   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);
643   if (flg || !((PetscObject)adapt)->type_name) {
644     ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr);
645   }
646 
647   ierr = PetscOptionsBool("-ts_adapt_always_accept","Always accept the step","TSAdaptSetAlwaysAccept",adapt->always_accept,&flg,&set);CHKERRQ(ierr);
648   if (set) {ierr = TSAdaptSetAlwaysAccept(adapt,flg);CHKERRQ(ierr);}
649 
650   safety = adapt->safety; reject_safety = adapt->reject_safety;
651   ierr = PetscOptionsReal("-ts_adapt_safety","Safety factor relative to target error/stability goal","TSAdaptSetSafety",safety,&safety,&set);CHKERRQ(ierr);
652   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);
653   if (set || flg) {ierr = TSAdaptSetSafety(adapt,safety,reject_safety);CHKERRQ(ierr);}
654 
655   two = 2; clip[0] = adapt->clip[0]; clip[1] = adapt->clip[1];
656   ierr = PetscOptionsRealArray("-ts_adapt_clip","Admissible decrease/increase factor in step size","TSAdaptSetClip",clip,&two,&set);CHKERRQ(ierr);
657   if (set && (two != 2)) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_clip");
658   if (set) {ierr = TSAdaptSetClip(adapt,clip[0],clip[1]);CHKERRQ(ierr);}
659 
660   hmin = adapt->dt_min; hmax = adapt->dt_max;
661   ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",hmin,&hmin,&set);CHKERRQ(ierr);
662   ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",hmax,&hmax,&flg);CHKERRQ(ierr);
663   if (set || flg) {ierr = TSAdaptSetStepLimits(adapt,hmin,hmax);CHKERRQ(ierr);}
664 
665   ierr = PetscOptionsReal("-ts_adapt_max_ignore","Adaptor ignores (absolute) solution values smaller than this value","",adapt->ignore_max,&adapt->ignore_max,&set);CHKERRQ(ierr);
666   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);
667 
668   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);
669 
670   ierr = PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);CHKERRQ(ierr);
671   if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported");
672 
673   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);
674 
675   ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
676   if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);}
677 
678   if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);}
679   ierr = PetscOptionsTail();CHKERRQ(ierr);
680   PetscFunctionReturn(0);
681 }
682 
683 /*@
684    TSAdaptCandidatesClear - clear any previously set candidate schemes
685 
686    Logically collective on TSAdapt
687 
688    Input Argument:
689 .  adapt - adaptive controller
690 
691    Level: developer
692 
693 .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose()
694 @*/
695 PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
696 {
697   PetscErrorCode ierr;
698 
699   PetscFunctionBegin;
700   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
701   ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr);
702   PetscFunctionReturn(0);
703 }
704 
705 /*@C
706    TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from
707 
708    Logically collective on TSAdapt
709 
710    Input Arguments:
711 +  adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate()
712 .  name - name of the candidate scheme to add
713 .  order - order of the candidate scheme
714 .  stageorder - stage order of the candidate scheme
715 .  ccfl - stability coefficient relative to explicit Euler, used for CFL constraints
716 .  cost - relative measure of the amount of work required for the candidate scheme
717 -  inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme
718 
719    Note:
720    This routine is not available in Fortran.
721 
722    Level: developer
723 
724 .seealso: TSAdaptCandidatesClear(), TSAdaptChoose()
725 @*/
726 PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse)
727 {
728   PetscInt c;
729 
730   PetscFunctionBegin;
731   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
732   if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order);
733   if (inuse) {
734     if (adapt->candidates.inuse_set) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()");
735     adapt->candidates.inuse_set = PETSC_TRUE;
736   }
737   /* first slot if this is the current scheme, otherwise the next available slot */
738   c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n;
739 
740   adapt->candidates.name[c]       = name;
741   adapt->candidates.order[c]      = order;
742   adapt->candidates.stageorder[c] = stageorder;
743   adapt->candidates.ccfl[c]       = ccfl;
744   adapt->candidates.cost[c]       = cost;
745   adapt->candidates.n++;
746   PetscFunctionReturn(0);
747 }
748 
749 /*@C
750    TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost
751 
752    Not Collective
753 
754    Input Arguments:
755 .  adapt - time step adaptivity context
756 
757    Output Arguments:
758 +  n - number of candidate schemes, always at least 1
759 .  order - the order of each candidate scheme
760 .  stageorder - the stage order of each candidate scheme
761 .  ccfl - the CFL coefficient of each scheme
762 -  cost - the relative cost of each scheme
763 
764    Level: developer
765 
766    Note:
767    The current scheme is always returned in the first slot
768 
769 .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose()
770 @*/
771 PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost)
772 {
773   PetscFunctionBegin;
774   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
775   if (n) *n = adapt->candidates.n;
776   if (order) *order = adapt->candidates.order;
777   if (stageorder) *stageorder = adapt->candidates.stageorder;
778   if (ccfl) *ccfl = adapt->candidates.ccfl;
779   if (cost) *cost = adapt->candidates.cost;
780   PetscFunctionReturn(0);
781 }
782 
783 /*@C
784    TSAdaptChoose - choose which method and step size to use for the next step
785 
786    Collective on TSAdapt
787 
788    Input Arguments:
789 +  adapt - adaptive contoller
790 -  h - current step size
791 
792    Output Arguments:
793 +  next_sc - optional, scheme to use for the next step
794 .  next_h - step size to use for the next step
795 -  accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size
796 
797    Note:
798    The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is
799    being retried after an initial rejection.
800 
801    Level: developer
802 
803 .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd()
804 @*/
805 PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept)
806 {
807   PetscErrorCode ierr;
808   PetscInt       ncandidates = adapt->candidates.n;
809   PetscInt       scheme = 0;
810   PetscReal      wlte = -1.0;
811   PetscReal      wltea = -1.0;
812   PetscReal      wlter = -1.0;
813 
814   PetscFunctionBegin;
815   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
816   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
817   if (next_sc) PetscValidIntPointer(next_sc,4);
818   PetscValidPointer(next_h,5);
819   PetscValidIntPointer(accept,6);
820   if (next_sc) *next_sc = 0;
821 
822   /* Do not mess with adaptivity while handling events*/
823   if (ts->event && ts->event->status != TSEVENT_NONE) {
824     *next_h = h;
825     *accept = PETSC_TRUE;
826     PetscFunctionReturn(0);
827   }
828 
829   ierr = (*adapt->ops->choose)(adapt,ts,h,&scheme,next_h,accept,&wlte,&wltea,&wlter);CHKERRQ(ierr);
830   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);
831   if (*next_h < 0) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h);
832   if (next_sc) *next_sc = scheme;
833 
834   if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
835     /* Increase/reduce step size if end time of next step is close to or overshoots max time */
836     PetscReal t = ts->ptime + ts->time_step, h = *next_h;
837     PetscReal tend = t + h, tmax = ts->max_time, hmax = tmax - t;
838     PetscReal a = (PetscReal)(1.0 + adapt->matchstepfac[0]);
839     PetscReal b = adapt->matchstepfac[1];
840     if (t < tmax && tend > tmax) *next_h = hmax;
841     if (t < tmax && tend < tmax && h*b > hmax) *next_h = hmax/2;
842     if (t < tmax && tend < tmax && h*a > hmax) *next_h = hmax;
843   }
844 
845   if (adapt->monitor) {
846     const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : "";
847     ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
848     if (wlte < 0) {
849       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);
850     } else {
851       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);
852     }
853     ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
854   }
855   PetscFunctionReturn(0);
856 }
857 
858 /*@
859    TSAdaptSetTimeStepIncreaseDelay - The number of timesteps to wait after a decrease in the timestep due to failed solver
860                                      before increasing the time step.
861 
862    Logicially Collective on TSAdapt
863 
864    Input Arguments:
865 +  adapt - adaptive controller context
866 -  cnt - the number of timesteps
867 
868    Options Database Key:
869 .  -ts_adapt_time_step_increase_delay cnt - number of steps to delay the increase
870 
871    Notes: This is to prevent an adaptor from bouncing back and forth between two nearby timesteps. The default is 0.
872           The successful use of this option is problem dependent
873 
874    Developer Note: there is no theory to support this option
875 
876    Level: advanced
877 
878 .seealso:
879 @*/
880 PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt,PetscInt cnt)
881 {
882   PetscFunctionBegin;
883   adapt->timestepjustdecreased_delay = cnt;
884   PetscFunctionReturn(0);
885 }
886 
887 
888 /*@
889    TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails)
890 
891    Collective on TSAdapt
892 
893    Input Arguments:
894 +  adapt - adaptive controller context
895 .  ts - time stepper
896 .  t - Current simulation time
897 -  Y - Current solution vector
898 
899    Output Arguments:
900 .  accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject
901 
902    Level: developer
903 
904 .seealso:
905 @*/
906 PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool *accept)
907 {
908   PetscErrorCode      ierr;
909   SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING;
910 
911   PetscFunctionBegin;
912   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
913   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
914   PetscValidIntPointer(accept,3);
915 
916   if (ts->snes) {ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);}
917   if (snesreason < 0) {
918     *accept = PETSC_FALSE;
919     if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) {
920       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
921       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);
922       if (adapt->monitor) {
923         ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
924         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);
925         ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
926       }
927     }
928   } else {
929     *accept = PETSC_TRUE;
930     ierr = TSFunctionDomainError(ts,t,Y,accept);CHKERRQ(ierr);
931     if(*accept && adapt->checkstage) {
932       ierr = (*adapt->checkstage)(adapt,ts,t,Y,accept);CHKERRQ(ierr);
933     }
934   }
935 
936   if(!(*accept) && !ts->reason) {
937     PetscReal dt,new_dt;
938     ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
939     new_dt = dt * adapt->scale_solve_failed;
940     ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr);
941     adapt->timestepjustdecreased += adapt->timestepjustdecreased_delay;
942     if (adapt->monitor) {
943       ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
944       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);
945       ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
946     }
947   }
948   PetscFunctionReturn(0);
949 }
950 
951 /*@
952   TSAdaptCreate - create an adaptive controller context for time stepping
953 
954   Collective on MPI_Comm
955 
956   Input Parameter:
957 . comm - The communicator
958 
959   Output Parameter:
960 . adapt - new TSAdapt object
961 
962   Level: developer
963 
964   Notes:
965   TSAdapt creation is handled by TS, so users should not need to call this function.
966 
967 .keywords: TSAdapt, create
968 .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy()
969 @*/
970 PetscErrorCode  TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt)
971 {
972   PetscErrorCode ierr;
973   TSAdapt        adapt;
974 
975   PetscFunctionBegin;
976   PetscValidPointer(inadapt,1);
977   *inadapt = NULL;
978   ierr = TSAdaptInitializePackage();CHKERRQ(ierr);
979 
980   ierr = PetscHeaderCreate(adapt,TSADAPT_CLASSID,"TSAdapt","Time stepping adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr);
981 
982   adapt->always_accept      = PETSC_FALSE;
983   adapt->safety             = 0.9;
984   adapt->reject_safety      = 0.5;
985   adapt->clip[0]            = 0.1;
986   adapt->clip[1]            = 10.;
987   adapt->dt_min             = 1e-20;
988   adapt->dt_max             = 1e+20;
989   adapt->ignore_max         = -1.0;
990   adapt->glee_use_local     = PETSC_TRUE;
991   adapt->scale_solve_failed = 0.25;
992   /* these two safety factors are not public, and they are used only in the TS_EXACTFINALTIME_MATCHSTEP case
993      to prevent from situations were unreasonably small time steps are taken in order to match the final time */
994   adapt->matchstepfac[0]    = 0.01; /* allow 1% step size increase in the last step */
995   adapt->matchstepfac[1]    = 2.0;  /* halve last step if it is greater than what remains divided this factor */
996   adapt->wnormtype          = NORM_2;
997   adapt->timestepjustdecreased_delay = 0;
998 
999   *inadapt = adapt;
1000   PetscFunctionReturn(0);
1001 }
1002