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