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