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