xref: /petsc/src/ts/adapt/interface/tsadapt.c (revision 27f49a208b01d2e827ab9db411a2d16003fe9262)
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    Level: advanced
27 
28    Notes:
29    `TSAdaptRegister()` may be called multiple times to add several user-defined families.
30 
31    Sample 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: [](chapter_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: [](chapter_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: [](chapter_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: [](chapter_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: [](chapter_ts), `TSGetAdapt()`, `TSAdaptDestroy()`, `TSAdaptType`, `TSAdaptGetType()`, `TSAdaptType`
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   PetscValidCharPointer(type, 2);
134   PetscCall(PetscObjectTypeCompare((PetscObject)adapt, type, &match));
135   if (match) PetscFunctionReturn(PETSC_SUCCESS);
136   PetscCall(PetscFunctionListFind(TSAdaptList, type, &r));
137   PetscCheck(r, PETSC_COMM_SELF, 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   PetscValidPointer(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 + newdm - 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: [](chapter_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: [](chapter_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: [](chapter_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    Arguments of func:
334 $  PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept)
335 
336 +  adapt - adaptive controller context
337 .  ts - time stepping context
338 -  accept - pending choice of whether to accept, can be modified by this routine
339 
340    Level: advanced
341 
342 .seealso: [](chapter_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptChoose()`
343 @*/
344 PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt, PetscErrorCode (*func)(TSAdapt, TS, PetscReal, Vec, PetscBool *))
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: [](chapter_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: [](chapter_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: [](chapter_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) PetscValidRealPointer(safety, 2);
432   if (reject_safety) PetscValidRealPointer(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: [](chapter_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: [](chapter_ts), `TSAdapt`, `TSAdaptSetMaxIgnore()`, `TSAdaptChoose()`
479 @*/
480 PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt adapt, PetscReal *max_ignore)
481 {
482   PetscFunctionBegin;
483   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
484   PetscValidRealPointer(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: [](chapter_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: [](chapter_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) PetscValidRealPointer(low, 2);
541   if (high) PetscValidRealPointer(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: [](chapter_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: [](chapter_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) PetscValidRealPointer(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: [](chapter_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: [](chapter_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) PetscValidRealPointer(hmin, 2);
653   if (hmax) PetscValidRealPointer(hmax, 3);
654   if (hmin) *hmin = adapt->dt_min;
655   if (hmax) *hmax = adapt->dt_max;
656   PetscFunctionReturn(PETSC_SUCCESS);
657 }
658 
659 /*
660    TSAdaptSetFromOptions - Sets various `TSAdapt` parameters from user options.
661 
662    Collective
663 
664    Input Parameter:
665 .  adapt - the `TSAdapt` context
666 
667    Options Database Keys:
668 +  -ts_adapt_type <type> - algorithm to use for adaptivity
669 .  -ts_adapt_always_accept - always accept steps regardless of error/stability goals
670 .  -ts_adapt_safety <safety> - safety factor relative to target error/stability goal
671 .  -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected
672 .  -ts_adapt_clip <low,high> - admissible time step decrease and increase factors
673 .  -ts_adapt_dt_min <min> - minimum timestep to use
674 .  -ts_adapt_dt_max <max> - maximum timestep to use
675 .  -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails
676 .  -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates
677 -  -ts_adapt_time_step_increase_delay - number of timesteps to delay increasing the time step after it has been decreased due to failed solver
678 
679    Level: advanced
680 
681    Note:
682    This function is automatically called by `TSSetFromOptions()`
683 
684 .seealso: [](chapter_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptSetAlwaysAccept()`, `TSAdaptSetSafety()`,
685           `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetStepLimits()`, `TSAdaptSetMonitor()`
686 */
687 PetscErrorCode TSAdaptSetFromOptions(TSAdapt adapt, PetscOptionItems *PetscOptionsObject)
688 {
689   char      type[256] = TSADAPTBASIC;
690   PetscReal safety, reject_safety, clip[2], scale, hmin, hmax;
691   PetscBool set, flg;
692   PetscInt  two;
693 
694   PetscFunctionBegin;
695   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
696   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
697    * function can only be called from inside TSSetFromOptions()  */
698   PetscOptionsHeadBegin(PetscOptionsObject, "TS Adaptivity options");
699   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));
700   if (flg || !((PetscObject)adapt)->type_name) PetscCall(TSAdaptSetType(adapt, type));
701 
702   PetscCall(PetscOptionsBool("-ts_adapt_always_accept", "Always accept the step", "TSAdaptSetAlwaysAccept", adapt->always_accept, &flg, &set));
703   if (set) PetscCall(TSAdaptSetAlwaysAccept(adapt, flg));
704 
705   safety        = adapt->safety;
706   reject_safety = adapt->reject_safety;
707   PetscCall(PetscOptionsReal("-ts_adapt_safety", "Safety factor relative to target error/stability goal", "TSAdaptSetSafety", safety, &safety, &set));
708   PetscCall(PetscOptionsReal("-ts_adapt_reject_safety", "Extra safety factor to apply if the last step was rejected", "TSAdaptSetSafety", reject_safety, &reject_safety, &flg));
709   if (set || flg) PetscCall(TSAdaptSetSafety(adapt, safety, reject_safety));
710 
711   two     = 2;
712   clip[0] = adapt->clip[0];
713   clip[1] = adapt->clip[1];
714   PetscCall(PetscOptionsRealArray("-ts_adapt_clip", "Admissible decrease/increase factor in step size", "TSAdaptSetClip", clip, &two, &set));
715   PetscCheck(!set || (two == 2), PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Must give exactly two values to -ts_adapt_clip");
716   if (set) PetscCall(TSAdaptSetClip(adapt, clip[0], clip[1]));
717 
718   hmin = adapt->dt_min;
719   hmax = adapt->dt_max;
720   PetscCall(PetscOptionsReal("-ts_adapt_dt_min", "Minimum time step considered", "TSAdaptSetStepLimits", hmin, &hmin, &set));
721   PetscCall(PetscOptionsReal("-ts_adapt_dt_max", "Maximum time step considered", "TSAdaptSetStepLimits", hmax, &hmax, &flg));
722   if (set || flg) PetscCall(TSAdaptSetStepLimits(adapt, hmin, hmax));
723 
724   PetscCall(PetscOptionsReal("-ts_adapt_max_ignore", "Adaptor ignores (absolute) solution values smaller than this value", "", adapt->ignore_max, &adapt->ignore_max, &set));
725   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));
726 
727   PetscCall(PetscOptionsReal("-ts_adapt_scale_solve_failed", "Scale step by this factor if solve fails", "TSAdaptSetScaleSolveFailed", adapt->scale_solve_failed, &scale, &set));
728   if (set) PetscCall(TSAdaptSetScaleSolveFailed(adapt, scale));
729 
730   PetscCall(PetscOptionsEnum("-ts_adapt_wnormtype", "Type of norm computed for error estimation", "", NormTypes, (PetscEnum)adapt->wnormtype, (PetscEnum *)&adapt->wnormtype, NULL));
731   PetscCheck(adapt->wnormtype == NORM_2 || adapt->wnormtype == NORM_INFINITY, PetscObjectComm((PetscObject)adapt), PETSC_ERR_SUP, "Only 2-norm and infinite norm supported");
732 
733   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));
734 
735   PetscCall(PetscOptionsBool("-ts_adapt_monitor", "Print choices made by adaptive controller", "TSAdaptSetMonitor", adapt->monitor ? PETSC_TRUE : PETSC_FALSE, &flg, &set));
736   if (set) PetscCall(TSAdaptSetMonitor(adapt, flg));
737 
738   PetscTryTypeMethod(adapt, setfromoptions, PetscOptionsObject);
739   PetscOptionsHeadEnd();
740   PetscFunctionReturn(PETSC_SUCCESS);
741 }
742 
743 /*@
744    TSAdaptCandidatesClear - clear any previously set candidate schemes
745 
746    Logically Collective
747 
748    Input Parameter:
749 .  adapt - adaptive controller
750 
751    Level: developer
752 
753 .seealso: [](chapter_ts), `TSAdapt`, `TSAdaptCreate()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()`
754 @*/
755 PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
756 {
757   PetscFunctionBegin;
758   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
759   PetscCall(PetscMemzero(&adapt->candidates, sizeof(adapt->candidates)));
760   PetscFunctionReturn(PETSC_SUCCESS);
761 }
762 
763 /*@C
764    TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from
765 
766    Logically Collective; No Fortran Support
767 
768    Input Parameters:
769 +  adapt - time step adaptivity context, obtained with `TSGetAdapt()` or `TSAdaptCreate()`
770 .  name - name of the candidate scheme to add
771 .  order - order of the candidate scheme
772 .  stageorder - stage order of the candidate scheme
773 .  ccfl - stability coefficient relative to explicit Euler, used for CFL constraints
774 .  cost - relative measure of the amount of work required for the candidate scheme
775 -  inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme
776 
777    Level: developer
778 
779 .seealso: [](chapter_ts), `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptChoose()`
780 @*/
781 PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt, const char name[], PetscInt order, PetscInt stageorder, PetscReal ccfl, PetscReal cost, PetscBool inuse)
782 {
783   PetscInt c;
784 
785   PetscFunctionBegin;
786   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
787   PetscCheck(order >= 1, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Classical order %" PetscInt_FMT " must be a positive integer", order);
788   if (inuse) {
789     PetscCheck(!adapt->candidates.inuse_set, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONGSTATE, "Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()");
790     adapt->candidates.inuse_set = PETSC_TRUE;
791   }
792   /* first slot if this is the current scheme, otherwise the next available slot */
793   c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n;
794 
795   adapt->candidates.name[c]       = name;
796   adapt->candidates.order[c]      = order;
797   adapt->candidates.stageorder[c] = stageorder;
798   adapt->candidates.ccfl[c]       = ccfl;
799   adapt->candidates.cost[c]       = cost;
800   adapt->candidates.n++;
801   PetscFunctionReturn(PETSC_SUCCESS);
802 }
803 
804 /*@C
805    TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost
806 
807    Not Collective
808 
809    Input Parameter:
810 .  adapt - time step adaptivity context
811 
812    Output Parameters:
813 +  n - number of candidate schemes, always at least 1
814 .  order - the order of each candidate scheme
815 .  stageorder - the stage order of each candidate scheme
816 .  ccfl - the CFL coefficient of each scheme
817 -  cost - the relative cost of each scheme
818 
819    Level: developer
820 
821    Note:
822    The current scheme is always returned in the first slot
823 
824 .seealso: [](chapter_ts), `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()`
825 @*/
826 PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt, PetscInt *n, const PetscInt **order, const PetscInt **stageorder, const PetscReal **ccfl, const PetscReal **cost)
827 {
828   PetscFunctionBegin;
829   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
830   if (n) *n = adapt->candidates.n;
831   if (order) *order = adapt->candidates.order;
832   if (stageorder) *stageorder = adapt->candidates.stageorder;
833   if (ccfl) *ccfl = adapt->candidates.ccfl;
834   if (cost) *cost = adapt->candidates.cost;
835   PetscFunctionReturn(PETSC_SUCCESS);
836 }
837 
838 /*@C
839    TSAdaptChoose - choose which method and step size to use for the next step
840 
841    Collective
842 
843    Input Parameters:
844 +  adapt - adaptive controller
845 .  ts - time stepper
846 -  h - current step size
847 
848    Output Parameters:
849 +  next_sc - optional, scheme to use for the next step
850 .  next_h - step size to use for the next step
851 -  accept - `PETSC_TRUE` to accept the current step, `PETSC_FALSE` to repeat the current step with the new step size
852 
853    Level: developer
854 
855    Note:
856    The input value of parameter accept is retained from the last time step, so it will be `PETSC_FALSE` if the step is
857    being retried after an initial rejection.
858 
859 .seealso: [](chapter_ts), `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()`
860 @*/
861 PetscErrorCode TSAdaptChoose(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept)
862 {
863   PetscInt  ncandidates = adapt->candidates.n;
864   PetscInt  scheme      = 0;
865   PetscReal wlte        = -1.0;
866   PetscReal wltea       = -1.0;
867   PetscReal wlter       = -1.0;
868 
869   PetscFunctionBegin;
870   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
871   PetscValidHeaderSpecific(ts, TS_CLASSID, 2);
872   if (next_sc) PetscValidIntPointer(next_sc, 4);
873   PetscValidRealPointer(next_h, 5);
874   PetscValidBoolPointer(accept, 6);
875   if (next_sc) *next_sc = 0;
876 
877   /* Do not mess with adaptivity while handling events*/
878   if (ts->event && ts->event->status != TSEVENT_NONE) {
879     *next_h = h;
880     *accept = PETSC_TRUE;
881     PetscFunctionReturn(PETSC_SUCCESS);
882   }
883 
884   PetscUseTypeMethod(adapt, choose, ts, h, &scheme, next_h, accept, &wlte, &wltea, &wlter);
885   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);
886   PetscCheck(*next_h >= 0, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Computed step size %g must be positive", (double)*next_h);
887   if (next_sc) *next_sc = scheme;
888 
889   if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
890     /* Increase/reduce step size if end time of next step is close to or overshoots max time */
891     PetscReal t = ts->ptime + ts->time_step, h = *next_h;
892     PetscReal tend = t + h, tmax, hmax;
893     PetscReal a    = (PetscReal)(1.0 + adapt->matchstepfac[0]);
894     PetscReal b    = adapt->matchstepfac[1];
895 
896     if (ts->tspan) {
897       if (PetscIsCloseAtTol(t, ts->tspan->span_times[ts->tspan->spanctr], ts->tspan->reltol * h + ts->tspan->abstol, 0)) /* hit a span time point */
898         if (ts->tspan->spanctr + 1 < ts->tspan->num_span_times) tmax = ts->tspan->span_times[ts->tspan->spanctr + 1];
899         else tmax = ts->max_time; /* hit the last span time point */
900       else tmax = ts->tspan->span_times[ts->tspan->spanctr];
901     } else tmax = ts->max_time;
902     hmax = tmax - t;
903     if (t < tmax && tend > tmax) *next_h = hmax;
904     if (t < tmax && tend < tmax && h * b > hmax) *next_h = hmax / 2;
905     if (t < tmax && tend < tmax && h * a > hmax) *next_h = hmax;
906     /* if step size is changed to match a span time point */
907     if (ts->tspan && h != *next_h && !adapt->dt_span_cached) adapt->dt_span_cached = h;
908     /* reset time step after a span time point */
909     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)) {
910       *next_h               = adapt->dt_span_cached;
911       adapt->dt_span_cached = 0;
912     }
913   }
914   if (adapt->monitor) {
915     const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : "";
916     PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
917     if (wlte < 0) {
918       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",
919                                        (double)ts->ptime, (double)h, (double)*next_h));
920     } else {
921       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",
922                                        (double)ts->ptime, (double)h, (double)*next_h, (double)wlte, (double)wltea, (double)wlter));
923     }
924     PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
925   }
926   PetscFunctionReturn(PETSC_SUCCESS);
927 }
928 
929 /*@
930    TSAdaptSetTimeStepIncreaseDelay - The number of timesteps to wait after a decrease in the timestep due to failed solver
931                                      before increasing the time step.
932 
933    Logicially Collective
934 
935    Input Parameters:
936 +  adapt - adaptive controller context
937 -  cnt - the number of timesteps
938 
939    Options Database Key:
940 .  -ts_adapt_time_step_increase_delay cnt - number of steps to delay the increase
941 
942    Level: advanced
943 
944    Notes:
945    This is to prevent an adaptor from bouncing back and forth between two nearby timesteps. The default is 0.
946 
947    The successful use of this option is problem dependent
948 
949    Developer Note:
950    There is no theory to support this option
951 
952 .seealso: [](chapter_ts), `TSAdapt`
953 @*/
954 PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt, PetscInt cnt)
955 {
956   PetscFunctionBegin;
957   adapt->timestepjustdecreased_delay = cnt;
958   PetscFunctionReturn(PETSC_SUCCESS);
959 }
960 
961 /*@
962    TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails or solution vector is infeasible)
963 
964    Collective
965 
966    Input Parameters:
967 +  adapt - adaptive controller context
968 .  ts - time stepper
969 .  t - Current simulation time
970 -  Y - Current solution vector
971 
972    Output Parameter:
973 .  accept - `PETSC_TRUE` to accept the stage, `PETSC_FALSE` to reject
974 
975    Level: developer
976 
977 .seealso: [](chapter_ts), `TSAdapt`
978 @*/
979 PetscErrorCode TSAdaptCheckStage(TSAdapt adapt, TS ts, PetscReal t, Vec Y, PetscBool *accept)
980 {
981   SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING;
982 
983   PetscFunctionBegin;
984   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
985   PetscValidHeaderSpecific(ts, TS_CLASSID, 2);
986   PetscValidBoolPointer(accept, 5);
987 
988   if (ts->snes) PetscCall(SNESGetConvergedReason(ts->snes, &snesreason));
989   if (snesreason < 0) {
990     *accept = PETSC_FALSE;
991     if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) {
992       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
993       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));
994       if (adapt->monitor) {
995         PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
996         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,
997                                          (double)ts->ptime, (double)ts->time_step, ts->num_snes_failures));
998         PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
999       }
1000     }
1001   } else {
1002     *accept = PETSC_TRUE;
1003     PetscCall(TSFunctionDomainError(ts, t, Y, accept));
1004     if (*accept && adapt->checkstage) {
1005       PetscCall((*adapt->checkstage)(adapt, ts, t, Y, accept));
1006       if (!*accept) {
1007         PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", solution rejected by user function provided by TSSetFunctionDomainError()\n", ts->steps));
1008         if (adapt->monitor) {
1009           PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
1010           PetscCall(PetscViewerASCIIPrintf(adapt->monitor, "    TSAdapt %s step %3" PetscInt_FMT " stage rejected by user function provided by TSSetFunctionDomainError()\n", ((PetscObject)adapt)->type_name, ts->steps));
1011           PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
1012         }
1013       }
1014     }
1015   }
1016 
1017   if (!(*accept) && !ts->reason) {
1018     PetscReal dt, new_dt;
1019     PetscCall(TSGetTimeStep(ts, &dt));
1020     new_dt = dt * adapt->scale_solve_failed;
1021     PetscCall(TSSetTimeStep(ts, new_dt));
1022     adapt->timestepjustdecreased += adapt->timestepjustdecreased_delay;
1023     if (adapt->monitor) {
1024       PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
1025       PetscCall(PetscViewerASCIIPrintf(adapt->monitor, "    TSAdapt %s step %3" PetscInt_FMT " 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));
1026       PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
1027     }
1028   }
1029   PetscFunctionReturn(PETSC_SUCCESS);
1030 }
1031 
1032 /*@
1033   TSAdaptCreate - create an adaptive controller context for time stepping
1034 
1035   Collective
1036 
1037   Input Parameter:
1038 . comm - The communicator
1039 
1040   Output Parameter:
1041 . adapt - new `TSAdapt` object
1042 
1043   Level: developer
1044 
1045   Note:
1046   `TSAdapt` creation is handled by `TS`, so users should not need to call this function.
1047 
1048 .seealso: [](chapter_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptDestroy()`
1049 @*/
1050 PetscErrorCode TSAdaptCreate(MPI_Comm comm, TSAdapt *inadapt)
1051 {
1052   TSAdapt adapt;
1053 
1054   PetscFunctionBegin;
1055   PetscValidPointer(inadapt, 2);
1056   *inadapt = NULL;
1057   PetscCall(TSAdaptInitializePackage());
1058 
1059   PetscCall(PetscHeaderCreate(adapt, TSADAPT_CLASSID, "TSAdapt", "Time stepping adaptivity", "TS", comm, TSAdaptDestroy, TSAdaptView));
1060 
1061   adapt->always_accept      = PETSC_FALSE;
1062   adapt->safety             = 0.9;
1063   adapt->reject_safety      = 0.5;
1064   adapt->clip[0]            = 0.1;
1065   adapt->clip[1]            = 10.;
1066   adapt->dt_min             = 1e-20;
1067   adapt->dt_max             = 1e+20;
1068   adapt->ignore_max         = -1.0;
1069   adapt->glee_use_local     = PETSC_TRUE;
1070   adapt->scale_solve_failed = 0.25;
1071   /* these two safety factors are not public, and they are used only in the TS_EXACTFINALTIME_MATCHSTEP case
1072      to prevent from situations were unreasonably small time steps are taken in order to match the final time */
1073   adapt->matchstepfac[0]             = 0.01; /* allow 1% step size increase in the last step */
1074   adapt->matchstepfac[1]             = 2.0;  /* halve last step if it is greater than what remains divided this factor */
1075   adapt->wnormtype                   = NORM_2;
1076   adapt->timestepjustdecreased_delay = 0;
1077 
1078   *inadapt = adapt;
1079   PetscFunctionReturn(PETSC_SUCCESS);
1080 }
1081