xref: /petsc/src/ts/adapt/interface/tsadapt.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
184df9cb4SJed Brown 
2af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I  "petscts.h" I*/
384df9cb4SJed Brown 
41b9b13dfSLisandro Dalcin PetscClassId TSADAPT_CLASSID;
51b9b13dfSLisandro Dalcin 
6140e18c1SBarry Smith static PetscFunctionList TSAdaptList;
784df9cb4SJed Brown static PetscBool         TSAdaptPackageInitialized;
884df9cb4SJed Brown static PetscBool         TSAdaptRegisterAllCalled;
984df9cb4SJed Brown 
108cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_None(TSAdapt);
111566a47fSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt);
12cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptCreate_DSP(TSAdapt);
138cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt);
141917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt);
15d949e4a4SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptCreate_History(TSAdapt);
1684df9cb4SJed Brown 
1784df9cb4SJed Brown /*@C
181c84c290SBarry Smith    TSAdaptRegister -  adds a TSAdapt implementation
191c84c290SBarry Smith 
201c84c290SBarry Smith    Not Collective
211c84c290SBarry Smith 
221c84c290SBarry Smith    Input Parameters:
232fe279fdSBarry Smith +  sname - name of user-defined adaptivity scheme
242fe279fdSBarry Smith -  function - routine to create method context
251c84c290SBarry Smith 
26bcf0153eSBarry Smith    Level: advanced
27bcf0153eSBarry Smith 
281c84c290SBarry Smith    Notes:
29bcf0153eSBarry Smith    `TSAdaptRegister()` may be called multiple times to add several user-defined families.
301c84c290SBarry Smith 
311c84c290SBarry Smith    Sample usage:
321c84c290SBarry Smith .vb
33bdf89e91SBarry Smith    TSAdaptRegister("my_scheme", MySchemeCreate);
341c84c290SBarry Smith .ve
351c84c290SBarry Smith 
361c84c290SBarry Smith    Then, your scheme can be chosen with the procedural interface via
371c84c290SBarry Smith $     TSAdaptSetType(ts, "my_scheme")
381c84c290SBarry Smith    or at runtime via the option
391c84c290SBarry Smith $     -ts_adapt_type my_scheme
4084df9cb4SJed Brown 
41*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdaptRegisterAll()`
4284df9cb4SJed Brown @*/
43d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptRegister(const char sname[], PetscErrorCode (*function)(TSAdapt))
44d71ae5a4SJacob Faibussowitsch {
4584df9cb4SJed Brown   PetscFunctionBegin;
469566063dSJacob Faibussowitsch   PetscCall(TSAdaptInitializePackage());
479566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&TSAdaptList, sname, function));
483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4984df9cb4SJed Brown }
5084df9cb4SJed Brown 
5184df9cb4SJed Brown /*@C
522fe279fdSBarry Smith   TSAdaptRegisterAll - Registers all of the adaptivity schemes in `TSAdapt`
5384df9cb4SJed Brown 
5484df9cb4SJed Brown   Not Collective
5584df9cb4SJed Brown 
5684df9cb4SJed Brown   Level: advanced
5784df9cb4SJed Brown 
58*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdaptRegisterDestroy()`
5984df9cb4SJed Brown @*/
60d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptRegisterAll(void)
61d71ae5a4SJacob Faibussowitsch {
6284df9cb4SJed Brown   PetscFunctionBegin;
633ba16761SJacob Faibussowitsch   if (TSAdaptRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
640f51fdf8SToby Isaac   TSAdaptRegisterAllCalled = PETSC_TRUE;
659566063dSJacob Faibussowitsch   PetscCall(TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None));
669566063dSJacob Faibussowitsch   PetscCall(TSAdaptRegister(TSADAPTBASIC, TSAdaptCreate_Basic));
679566063dSJacob Faibussowitsch   PetscCall(TSAdaptRegister(TSADAPTDSP, TSAdaptCreate_DSP));
689566063dSJacob Faibussowitsch   PetscCall(TSAdaptRegister(TSADAPTCFL, TSAdaptCreate_CFL));
699566063dSJacob Faibussowitsch   PetscCall(TSAdaptRegister(TSADAPTGLEE, TSAdaptCreate_GLEE));
709566063dSJacob Faibussowitsch   PetscCall(TSAdaptRegister(TSADAPTHISTORY, TSAdaptCreate_History));
713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7284df9cb4SJed Brown }
7384df9cb4SJed Brown 
7484df9cb4SJed Brown /*@C
752fe279fdSBarry Smith   TSAdaptFinalizePackage - This function destroys everything in the `TS` package. It is
762fe279fdSBarry Smith   called from `PetscFinalize()`.
7784df9cb4SJed Brown 
7884df9cb4SJed Brown   Level: developer
7984df9cb4SJed Brown 
80*1cc06b55SBarry Smith .seealso: [](ch_ts), `PetscFinalize()`
8184df9cb4SJed Brown @*/
82d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptFinalizePackage(void)
83d71ae5a4SJacob Faibussowitsch {
8484df9cb4SJed Brown   PetscFunctionBegin;
859566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&TSAdaptList));
8684df9cb4SJed Brown   TSAdaptPackageInitialized = PETSC_FALSE;
8784df9cb4SJed Brown   TSAdaptRegisterAllCalled  = PETSC_FALSE;
883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8984df9cb4SJed Brown }
9084df9cb4SJed Brown 
9184df9cb4SJed Brown /*@C
922fe279fdSBarry Smith   TSAdaptInitializePackage - This function initializes everything in the `TSAdapt` package. It is
932fe279fdSBarry Smith   called from `TSInitializePackage()`.
9484df9cb4SJed Brown 
9584df9cb4SJed Brown   Level: developer
9684df9cb4SJed Brown 
97*1cc06b55SBarry Smith .seealso: [](ch_ts), `PetscInitialize()`
9884df9cb4SJed Brown @*/
99d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptInitializePackage(void)
100d71ae5a4SJacob Faibussowitsch {
10184df9cb4SJed Brown   PetscFunctionBegin;
1023ba16761SJacob Faibussowitsch   if (TSAdaptPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
10384df9cb4SJed Brown   TSAdaptPackageInitialized = PETSC_TRUE;
1049566063dSJacob Faibussowitsch   PetscCall(PetscClassIdRegister("TSAdapt", &TSADAPT_CLASSID));
1059566063dSJacob Faibussowitsch   PetscCall(TSAdaptRegisterAll());
1069566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSAdaptFinalizePackage));
1073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10884df9cb4SJed Brown }
10984df9cb4SJed Brown 
1107eef6055SBarry Smith /*@C
111bcf0153eSBarry Smith   TSAdaptSetType - sets the approach used for the error adapter
1127eef6055SBarry Smith 
11320f4b53cSBarry Smith   Logicially Collective
1147eef6055SBarry Smith 
115d8d19677SJose E. Roman   Input Parameters:
11658b119d1SBarry Smith + adapt - the `TS` adapter, most likely obtained with `TSGetAdapt()`
117bcf0153eSBarry Smith - type - one of the `TSAdaptType`
1187eef6055SBarry Smith 
119bcf0153eSBarry Smith   Options Database Key:
120ec18b7bcSLisandro Dalcin . -ts_adapt_type <basic or dsp or none> - to set the adapter type
1217eef6055SBarry Smith 
1227eef6055SBarry Smith   Level: intermediate
1237eef6055SBarry Smith 
124*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSGetAdapt()`, `TSAdaptDestroy()`, `TSAdaptType`, `TSAdaptGetType()`, `TSAdaptType`
1257eef6055SBarry Smith @*/
126d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetType(TSAdapt adapt, TSAdaptType type)
127d71ae5a4SJacob Faibussowitsch {
128ef749922SLisandro Dalcin   PetscBool match;
1295f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(TSAdapt);
13084df9cb4SJed Brown 
13184df9cb4SJed Brown   PetscFunctionBegin;
1324782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
133b92453a8SLisandro Dalcin   PetscValidCharPointer(type, 2);
1349566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)adapt, type, &match));
1353ba16761SJacob Faibussowitsch   if (match) PetscFunctionReturn(PETSC_SUCCESS);
1369566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(TSAdaptList, type, &r));
1373c633725SBarry Smith   PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSAdapt type \"%s\" given", type);
138dbbe0bcdSBarry Smith   PetscTryTypeMethod(adapt, destroy);
1399566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(adapt->ops, sizeof(struct _TSAdaptOps)));
1409566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)adapt, type));
1419566063dSJacob Faibussowitsch   PetscCall((*r)(adapt));
1423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14384df9cb4SJed Brown }
14484df9cb4SJed Brown 
145d0288e4fSLisandro Dalcin /*@C
146a1cb98faSBarry Smith   TSAdaptGetType - gets the `TS` adapter method type (as a string).
147d0288e4fSLisandro Dalcin 
148d0288e4fSLisandro Dalcin   Not Collective
149d0288e4fSLisandro Dalcin 
150d0288e4fSLisandro Dalcin   Input Parameter:
151a1cb98faSBarry Smith . adapt - The `TS` adapter, most likely obtained with `TSGetAdapt()`
152d0288e4fSLisandro Dalcin 
153d0288e4fSLisandro Dalcin   Output Parameter:
154a1cb98faSBarry Smith . type - The name of `TS` adapter method
155d0288e4fSLisandro Dalcin 
156d0288e4fSLisandro Dalcin   Level: intermediate
157d0288e4fSLisandro Dalcin 
158a1cb98faSBarry Smith .seealso: `TSAdapt`, `TSAdaptType`, `TSAdaptSetType()`
159d0288e4fSLisandro Dalcin @*/
160d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetType(TSAdapt adapt, TSAdaptType *type)
161d71ae5a4SJacob Faibussowitsch {
162d0288e4fSLisandro Dalcin   PetscFunctionBegin;
163d0288e4fSLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
164d0288e4fSLisandro Dalcin   PetscValidPointer(type, 2);
165d0288e4fSLisandro Dalcin   *type = ((PetscObject)adapt)->type_name;
1663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
167d0288e4fSLisandro Dalcin }
168d0288e4fSLisandro Dalcin 
169d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt, const char prefix[])
170d71ae5a4SJacob Faibussowitsch {
17184df9cb4SJed Brown   PetscFunctionBegin;
1724782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
1739566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adapt, prefix));
1743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17584df9cb4SJed Brown }
17684df9cb4SJed Brown 
177ad6bc421SBarry Smith /*@C
1782fe279fdSBarry Smith   TSAdaptLoad - Loads a TSAdapt that has been stored in binary with `TSAdaptView()`.
179ad6bc421SBarry Smith 
180c3339decSBarry Smith   Collective
181ad6bc421SBarry Smith 
182ad6bc421SBarry Smith   Input Parameters:
183bcf0153eSBarry Smith + newdm - the newly loaded `TSAdapt`, this needs to have been created with `TSAdaptCreate()` or
184bcf0153eSBarry Smith            some related function before a call to `TSAdaptLoad()`.
185bcf0153eSBarry Smith - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()` or
186bcf0153eSBarry Smith            HDF5 file viewer, obtained from `PetscViewerHDF5Open()`
187ad6bc421SBarry Smith 
188ad6bc421SBarry Smith    Level: intermediate
189ad6bc421SBarry Smith 
190bcf0153eSBarry Smith   Note:
191bcf0153eSBarry Smith    The type is determined by the data in the file, any type set into the `TSAdapt` before this call is ignored.
192ad6bc421SBarry Smith 
193*1cc06b55SBarry Smith .seealso: [](ch_ts), `PetscViewerBinaryOpen()`, `TSAdaptView()`, `MatLoad()`, `VecLoad()`, `TSAdapt`
194ad6bc421SBarry Smith @*/
195d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptLoad(TSAdapt adapt, PetscViewer viewer)
196d71ae5a4SJacob Faibussowitsch {
197ad6bc421SBarry Smith   PetscBool isbinary;
198ad6bc421SBarry Smith   char      type[256];
199ad6bc421SBarry Smith 
200ad6bc421SBarry Smith   PetscFunctionBegin;
2014782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
202ad6bc421SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
2039566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
2043c633725SBarry Smith   PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()");
205ad6bc421SBarry Smith 
2069566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR));
2079566063dSJacob Faibussowitsch   PetscCall(TSAdaptSetType(adapt, type));
208dbbe0bcdSBarry Smith   PetscTryTypeMethod(adapt, load, viewer);
2093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
210ad6bc421SBarry Smith }
211ad6bc421SBarry Smith 
212d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptView(TSAdapt adapt, PetscViewer viewer)
213d71ae5a4SJacob Faibussowitsch {
2141c167fc2SEmil Constantinescu   PetscBool iascii, isbinary, isnone, isglee;
21584df9cb4SJed Brown 
21684df9cb4SJed Brown   PetscFunctionBegin;
2174782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
2189566063dSJacob Faibussowitsch   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt), &viewer));
2194782b174SLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
2204782b174SLisandro Dalcin   PetscCheckSameComm(adapt, 1, viewer, 2);
2219566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
2229566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
22384df9cb4SJed Brown   if (iascii) {
2249566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)adapt, viewer));
2259566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTNONE, &isnone));
2269566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTGLEE, &isglee));
2271917a363SLisandro Dalcin     if (!isnone) {
2289566063dSJacob Faibussowitsch       if (adapt->always_accept) PetscCall(PetscViewerASCIIPrintf(viewer, "  always accepting steps\n"));
2299566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  safety factor %g\n", (double)adapt->safety));
2309566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  extra safety factor after step rejection %g\n", (double)adapt->reject_safety));
2319566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  clip fastest increase %g\n", (double)adapt->clip[1]));
2329566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  clip fastest decrease %g\n", (double)adapt->clip[0]));
2339566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum allowed timestep %g\n", (double)adapt->dt_max));
2349566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  minimum allowed timestep %g\n", (double)adapt->dt_min));
2359566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum solution absolute value to be ignored %g\n", (double)adapt->ignore_max));
2361c167fc2SEmil Constantinescu     }
2371c167fc2SEmil Constantinescu     if (isglee) {
2381c167fc2SEmil Constantinescu       if (adapt->glee_use_local) {
2399566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  GLEE uses local error control\n"));
2401c167fc2SEmil Constantinescu       } else {
2419566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  GLEE uses global error control\n"));
2421c167fc2SEmil Constantinescu       }
2431917a363SLisandro Dalcin     }
2449566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
245dbbe0bcdSBarry Smith     PetscTryTypeMethod(adapt, view, viewer);
2469566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
247ad6bc421SBarry Smith   } else if (isbinary) {
248ad6bc421SBarry Smith     char type[256];
249ad6bc421SBarry Smith 
250ad6bc421SBarry Smith     /* need to save FILE_CLASS_ID for adapt class */
2519566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(type, ((PetscObject)adapt)->type_name, 256));
2529566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR));
253dbbe0bcdSBarry Smith   } else PetscTryTypeMethod(adapt, view, viewer);
2543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
25584df9cb4SJed Brown }
25684df9cb4SJed Brown 
257099af0a3SLisandro Dalcin /*@
258bcf0153eSBarry Smith    TSAdaptReset - Resets a `TSAdapt` context to its defaults
259099af0a3SLisandro Dalcin 
260c3339decSBarry Smith    Collective
261099af0a3SLisandro Dalcin 
262099af0a3SLisandro Dalcin    Input Parameter:
263bcf0153eSBarry Smith .  adapt - the `TSAdapt` context obtained from `TSGetAdapt()` or `TSAdaptCreate()`
264099af0a3SLisandro Dalcin 
265099af0a3SLisandro Dalcin    Level: developer
266099af0a3SLisandro Dalcin 
267*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSGetAdapt()`, `TSAdapt`, `TSAdaptCreate()`, `TSAdaptDestroy()`
268099af0a3SLisandro Dalcin @*/
269d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptReset(TSAdapt adapt)
270d71ae5a4SJacob Faibussowitsch {
271099af0a3SLisandro Dalcin   PetscFunctionBegin;
272099af0a3SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
273dbbe0bcdSBarry Smith   PetscTryTypeMethod(adapt, reset);
2743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
275099af0a3SLisandro Dalcin }
276099af0a3SLisandro Dalcin 
277d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptDestroy(TSAdapt *adapt)
278d71ae5a4SJacob Faibussowitsch {
27984df9cb4SJed Brown   PetscFunctionBegin;
2803ba16761SJacob Faibussowitsch   if (!*adapt) PetscFunctionReturn(PETSC_SUCCESS);
28184df9cb4SJed Brown   PetscValidHeaderSpecific(*adapt, TSADAPT_CLASSID, 1);
2829371c9d4SSatish Balay   if (--((PetscObject)(*adapt))->refct > 0) {
2839371c9d4SSatish Balay     *adapt = NULL;
2843ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2859371c9d4SSatish Balay   }
286099af0a3SLisandro Dalcin 
2879566063dSJacob Faibussowitsch   PetscCall(TSAdaptReset(*adapt));
288099af0a3SLisandro Dalcin 
289dbbe0bcdSBarry Smith   PetscTryTypeMethod((*adapt), destroy);
2909566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*adapt)->monitor));
2919566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(adapt));
2923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29384df9cb4SJed Brown }
29484df9cb4SJed Brown 
2951c3436cfSJed Brown /*@
2961c3436cfSJed Brown    TSAdaptSetMonitor - Monitor the choices made by the adaptive controller
2971c3436cfSJed Brown 
298c3339decSBarry Smith    Collective
2991c3436cfSJed Brown 
3004165533cSJose E. Roman    Input Parameters:
3011c3436cfSJed Brown +  adapt - adaptive controller context
302bcf0153eSBarry Smith -  flg - `PETSC_TRUE` to active a monitor, `PETSC_FALSE` to disable
3031c3436cfSJed Brown 
304bcf0153eSBarry Smith    Options Database Key:
305ec18b7bcSLisandro Dalcin .  -ts_adapt_monitor - to turn on monitoring
306bf997491SLisandro Dalcin 
3071c3436cfSJed Brown    Level: intermediate
3081c3436cfSJed Brown 
309*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptChoose()`
3101c3436cfSJed Brown @*/
311d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt, PetscBool flg)
312d71ae5a4SJacob Faibussowitsch {
3131c3436cfSJed Brown   PetscFunctionBegin;
3144782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
3154782b174SLisandro Dalcin   PetscValidLogicalCollectiveBool(adapt, flg, 2);
3161c3436cfSJed Brown   if (flg) {
3179566063dSJacob Faibussowitsch     if (!adapt->monitor) PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt), "stdout", &adapt->monitor));
3181c3436cfSJed Brown   } else {
3199566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&adapt->monitor));
3201c3436cfSJed Brown   }
3213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3221c3436cfSJed Brown }
3231c3436cfSJed Brown 
3240873808bSJed Brown /*@C
325bf997491SLisandro Dalcin    TSAdaptSetCheckStage - Set a callback to check convergence for a stage
3260873808bSJed Brown 
3272fe279fdSBarry Smith    Logically Collective
3280873808bSJed Brown 
3294165533cSJose E. Roman    Input Parameters:
3300873808bSJed Brown +  adapt - adaptive controller context
3310873808bSJed Brown -  func - stage check function
3320873808bSJed Brown 
3332fe279fdSBarry Smith   Calling Sequence of `func`:
3340873808bSJed Brown $  PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept)
3350873808bSJed Brown +  adapt - adaptive controller context
3360873808bSJed Brown .  ts - time stepping context
3370873808bSJed Brown -  accept - pending choice of whether to accept, can be modified by this routine
3380873808bSJed Brown 
3390873808bSJed Brown    Level: advanced
3400873808bSJed Brown 
341*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptChoose()`
3420873808bSJed Brown @*/
343d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt, PetscErrorCode (*func)(TSAdapt, TS, PetscReal, Vec, PetscBool *))
344d71ae5a4SJacob Faibussowitsch {
3450873808bSJed Brown   PetscFunctionBegin;
3460873808bSJed Brown   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
34768ae941aSLisandro Dalcin   adapt->checkstage = func;
3483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3490873808bSJed Brown }
3500873808bSJed Brown 
3511c3436cfSJed Brown /*@
352bf997491SLisandro Dalcin    TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of
353bf997491SLisandro Dalcin    any error or stability condition not meeting the prescribed goal.
354bf997491SLisandro Dalcin 
3552fe279fdSBarry Smith    Logically Collective
356bf997491SLisandro Dalcin 
3574165533cSJose E. Roman    Input Parameters:
358bcf0153eSBarry Smith +  adapt - time step adaptivity context, usually gotten with `TSGetAdapt()`
359bf997491SLisandro Dalcin -  flag - whether to always accept steps
360bf997491SLisandro Dalcin 
361bcf0153eSBarry Smith    Options Database Key:
362ec18b7bcSLisandro Dalcin .  -ts_adapt_always_accept - to always accept steps
363bf997491SLisandro Dalcin 
364bf997491SLisandro Dalcin    Level: intermediate
365bf997491SLisandro Dalcin 
366*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptChoose()`
367bf997491SLisandro Dalcin @*/
368d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt, PetscBool flag)
369d71ae5a4SJacob Faibussowitsch {
370bf997491SLisandro Dalcin   PetscFunctionBegin;
371bf997491SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
372bf997491SLisandro Dalcin   PetscValidLogicalCollectiveBool(adapt, flag, 2);
373bf997491SLisandro Dalcin   adapt->always_accept = flag;
3743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
375bf997491SLisandro Dalcin }
376bf997491SLisandro Dalcin 
377bf997491SLisandro Dalcin /*@
378bcf0153eSBarry Smith    TSAdaptSetSafety - Set safety factors for time step adaptor
3791c3436cfSJed Brown 
3802fe279fdSBarry Smith    Logically Collective
3811917a363SLisandro Dalcin 
3824165533cSJose E. Roman    Input Parameters:
3831917a363SLisandro Dalcin +  adapt - adaptive controller context
3841917a363SLisandro Dalcin .  safety - safety factor relative to target error/stability goal
385ec18b7bcSLisandro Dalcin -  reject_safety - extra safety factor to apply if the last step was rejected
3861917a363SLisandro Dalcin 
3871917a363SLisandro Dalcin    Options Database Keys:
388ec18b7bcSLisandro Dalcin +  -ts_adapt_safety <safety> - to set safety factor
389ec18b7bcSLisandro Dalcin -  -ts_adapt_reject_safety <reject_safety> - to set reject safety factor
3901917a363SLisandro Dalcin 
3911917a363SLisandro Dalcin    Level: intermediate
3921917a363SLisandro Dalcin 
393*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptGetSafety()`, `TSAdaptChoose()`
3941917a363SLisandro Dalcin @*/
395d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetSafety(TSAdapt adapt, PetscReal safety, PetscReal reject_safety)
396d71ae5a4SJacob Faibussowitsch {
3971917a363SLisandro Dalcin   PetscFunctionBegin;
3981917a363SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
3991917a363SLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt, safety, 2);
4001917a363SLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt, reject_safety, 3);
40113bcc0bdSJacob Faibussowitsch   PetscCheck(safety == (PetscReal)PETSC_DEFAULT || safety >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Safety factor %g must be non negative", (double)safety);
40213bcc0bdSJacob Faibussowitsch   PetscCheck(safety == (PetscReal)PETSC_DEFAULT || safety <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Safety factor %g must be less than one", (double)safety);
40313bcc0bdSJacob Faibussowitsch   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);
40413bcc0bdSJacob Faibussowitsch   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);
40513bcc0bdSJacob Faibussowitsch   if (safety != (PetscReal)PETSC_DEFAULT) adapt->safety = safety;
40613bcc0bdSJacob Faibussowitsch   if (reject_safety != (PetscReal)PETSC_DEFAULT) adapt->reject_safety = reject_safety;
4073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4081917a363SLisandro Dalcin }
4091917a363SLisandro Dalcin 
4101917a363SLisandro Dalcin /*@
411bcf0153eSBarry Smith    TSAdaptGetSafety - Get safety factors for time step adapter
4121917a363SLisandro Dalcin 
4131917a363SLisandro Dalcin    Not Collective
4141917a363SLisandro Dalcin 
4154165533cSJose E. Roman    Input Parameter:
4161917a363SLisandro Dalcin .  adapt - adaptive controller context
4171917a363SLisandro Dalcin 
4184165533cSJose E. Roman    Output Parameters:
4191917a363SLisandro Dalcin .  safety - safety factor relative to target error/stability goal
4201917a363SLisandro Dalcin +  reject_safety - extra safety factor to apply if the last step was rejected
4211917a363SLisandro Dalcin 
4221917a363SLisandro Dalcin    Level: intermediate
4231917a363SLisandro Dalcin 
424*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptSetSafety()`, `TSAdaptChoose()`
4251917a363SLisandro Dalcin @*/
426d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetSafety(TSAdapt adapt, PetscReal *safety, PetscReal *reject_safety)
427d71ae5a4SJacob Faibussowitsch {
4281917a363SLisandro Dalcin   PetscFunctionBegin;
4291917a363SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
4301917a363SLisandro Dalcin   if (safety) PetscValidRealPointer(safety, 2);
4311917a363SLisandro Dalcin   if (reject_safety) PetscValidRealPointer(reject_safety, 3);
4321917a363SLisandro Dalcin   if (safety) *safety = adapt->safety;
4331917a363SLisandro Dalcin   if (reject_safety) *reject_safety = adapt->reject_safety;
4343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4351917a363SLisandro Dalcin }
4361917a363SLisandro Dalcin 
4371917a363SLisandro Dalcin /*@
438bcf0153eSBarry Smith    TSAdaptSetMaxIgnore - Set error estimation threshold. Solution components below this threshold value will not be considered when computing error norms
439bcf0153eSBarry Smith    for time step adaptivity (in absolute value). A negative value (default) of the threshold leads to considering all solution components.
44076cddca1SEmil Constantinescu 
4412fe279fdSBarry Smith    Logically Collective
44276cddca1SEmil Constantinescu 
4434165533cSJose E. Roman    Input Parameters:
44476cddca1SEmil Constantinescu +  adapt - adaptive controller context
44576cddca1SEmil Constantinescu -  max_ignore - threshold for solution components that are ignored during error estimation
44676cddca1SEmil Constantinescu 
447bcf0153eSBarry Smith    Options Database Key:
44876cddca1SEmil Constantinescu .  -ts_adapt_max_ignore <max_ignore> - to set the threshold
44976cddca1SEmil Constantinescu 
45076cddca1SEmil Constantinescu    Level: intermediate
45176cddca1SEmil Constantinescu 
452*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptGetMaxIgnore()`, `TSAdaptChoose()`
45376cddca1SEmil Constantinescu @*/
454d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt adapt, PetscReal max_ignore)
455d71ae5a4SJacob Faibussowitsch {
45676cddca1SEmil Constantinescu   PetscFunctionBegin;
45776cddca1SEmil Constantinescu   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
45876cddca1SEmil Constantinescu   PetscValidLogicalCollectiveReal(adapt, max_ignore, 2);
45976cddca1SEmil Constantinescu   adapt->ignore_max = max_ignore;
4603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
46176cddca1SEmil Constantinescu }
46276cddca1SEmil Constantinescu 
46376cddca1SEmil Constantinescu /*@
464bcf0153eSBarry Smith    TSAdaptGetMaxIgnore - Get error estimation threshold. Solution components below this threshold value will not be considered when computing error norms
465bcf0153eSBarry Smith    for time step adaptivity (in absolute value).
46676cddca1SEmil Constantinescu 
46776cddca1SEmil Constantinescu    Not Collective
46876cddca1SEmil Constantinescu 
4694165533cSJose E. Roman    Input Parameter:
47076cddca1SEmil Constantinescu .  adapt - adaptive controller context
47176cddca1SEmil Constantinescu 
4724165533cSJose E. Roman    Output Parameter:
47376cddca1SEmil Constantinescu .  max_ignore - threshold for solution components that are ignored during error estimation
47476cddca1SEmil Constantinescu 
47576cddca1SEmil Constantinescu    Level: intermediate
47676cddca1SEmil Constantinescu 
477*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptSetMaxIgnore()`, `TSAdaptChoose()`
47876cddca1SEmil Constantinescu @*/
479d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt adapt, PetscReal *max_ignore)
480d71ae5a4SJacob Faibussowitsch {
48176cddca1SEmil Constantinescu   PetscFunctionBegin;
48276cddca1SEmil Constantinescu   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
48376cddca1SEmil Constantinescu   PetscValidRealPointer(max_ignore, 2);
48476cddca1SEmil Constantinescu   *max_ignore = adapt->ignore_max;
4853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
48676cddca1SEmil Constantinescu }
48776cddca1SEmil Constantinescu 
48876cddca1SEmil Constantinescu /*@
489bcf0153eSBarry Smith    TSAdaptSetClip - Sets the admissible decrease/increase factor in step size in the time step adapter
4901917a363SLisandro Dalcin 
491c3339decSBarry Smith    Logically collective
4921917a363SLisandro Dalcin 
4934165533cSJose E. Roman    Input Parameters:
4941917a363SLisandro Dalcin +  adapt - adaptive controller context
4951917a363SLisandro Dalcin .  low - admissible decrease factor
496ec18b7bcSLisandro Dalcin -  high - admissible increase factor
4971917a363SLisandro Dalcin 
498bcf0153eSBarry Smith    Options Database Key:
499ec18b7bcSLisandro Dalcin .  -ts_adapt_clip <low>,<high> - to set admissible time step decrease and increase factors
5001917a363SLisandro Dalcin 
5011917a363SLisandro Dalcin    Level: intermediate
5021917a363SLisandro Dalcin 
503*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptGetClip()`, `TSAdaptSetScaleSolveFailed()`
5041917a363SLisandro Dalcin @*/
505d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetClip(TSAdapt adapt, PetscReal low, PetscReal high)
506d71ae5a4SJacob Faibussowitsch {
5071917a363SLisandro Dalcin   PetscFunctionBegin;
5081917a363SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
5091917a363SLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt, low, 2);
5101917a363SLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt, high, 3);
51113bcc0bdSJacob Faibussowitsch   PetscCheck(low == (PetscReal)PETSC_DEFAULT || low >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Decrease factor %g must be non negative", (double)low);
51213bcc0bdSJacob Faibussowitsch   PetscCheck(low == (PetscReal)PETSC_DEFAULT || low <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Decrease factor %g must be less than one", (double)low);
51313bcc0bdSJacob Faibussowitsch   PetscCheck(high == (PetscReal)PETSC_DEFAULT || high >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Increase factor %g must be greater than one", (double)high);
51413bcc0bdSJacob Faibussowitsch   if (low != (PetscReal)PETSC_DEFAULT) adapt->clip[0] = low;
51513bcc0bdSJacob Faibussowitsch   if (high != (PetscReal)PETSC_DEFAULT) adapt->clip[1] = high;
5163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5171917a363SLisandro Dalcin }
5181917a363SLisandro Dalcin 
5191917a363SLisandro Dalcin /*@
520bcf0153eSBarry Smith    TSAdaptGetClip - Gets the admissible decrease/increase factor in step size in the time step adapter
5211917a363SLisandro Dalcin 
5221917a363SLisandro Dalcin    Not Collective
5231917a363SLisandro Dalcin 
5244165533cSJose E. Roman    Input Parameter:
5251917a363SLisandro Dalcin .  adapt - adaptive controller context
5261917a363SLisandro Dalcin 
5274165533cSJose E. Roman    Output Parameters:
5281917a363SLisandro Dalcin +  low - optional, admissible decrease factor
5291917a363SLisandro Dalcin -  high - optional, admissible increase factor
5301917a363SLisandro Dalcin 
5311917a363SLisandro Dalcin    Level: intermediate
5321917a363SLisandro Dalcin 
533*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()`
5341917a363SLisandro Dalcin @*/
535d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetClip(TSAdapt adapt, PetscReal *low, PetscReal *high)
536d71ae5a4SJacob Faibussowitsch {
5371917a363SLisandro Dalcin   PetscFunctionBegin;
5381917a363SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
5391917a363SLisandro Dalcin   if (low) PetscValidRealPointer(low, 2);
5401917a363SLisandro Dalcin   if (high) PetscValidRealPointer(high, 3);
5411917a363SLisandro Dalcin   if (low) *low = adapt->clip[0];
5421917a363SLisandro Dalcin   if (high) *high = adapt->clip[1];
5433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5441917a363SLisandro Dalcin }
5451917a363SLisandro Dalcin 
5461917a363SLisandro Dalcin /*@
547bcf0153eSBarry Smith    TSAdaptSetScaleSolveFailed - Scale step size by this factor if solve fails
54862c23b28SMark 
54920f4b53cSBarry Smith    Logically Collective
55062c23b28SMark 
5514165533cSJose E. Roman    Input Parameters:
55262c23b28SMark +  adapt - adaptive controller context
553ee300463SSatish Balay -  scale - scale
55462c23b28SMark 
555bcf0153eSBarry Smith    Options Database Key:
55662c23b28SMark .  -ts_adapt_scale_solve_failed <scale> - to set scale step by this factor if solve fails
55762c23b28SMark 
55862c23b28SMark    Level: intermediate
55962c23b28SMark 
560*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptGetScaleSolveFailed()`, `TSAdaptGetClip()`
56162c23b28SMark @*/
562d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt adapt, PetscReal scale)
563d71ae5a4SJacob Faibussowitsch {
56462c23b28SMark   PetscFunctionBegin;
56562c23b28SMark   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
56662c23b28SMark   PetscValidLogicalCollectiveReal(adapt, scale, 2);
56713bcc0bdSJacob Faibussowitsch   PetscCheck(scale == (PetscReal)PETSC_DEFAULT || scale > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Scale factor %g must be positive", (double)scale);
56813bcc0bdSJacob Faibussowitsch   PetscCheck(scale == (PetscReal)PETSC_DEFAULT || scale <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Scale factor %g must be less than one", (double)scale);
56913bcc0bdSJacob Faibussowitsch   if (scale != (PetscReal)PETSC_DEFAULT) adapt->scale_solve_failed = scale;
5703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
57162c23b28SMark }
57262c23b28SMark 
57362c23b28SMark /*@
57462c23b28SMark    TSAdaptGetScaleSolveFailed - Gets the admissible decrease/increase factor in step size
57562c23b28SMark 
57662c23b28SMark    Not Collective
57762c23b28SMark 
5784165533cSJose E. Roman    Input Parameter:
57962c23b28SMark .  adapt - adaptive controller context
58062c23b28SMark 
5814165533cSJose E. Roman    Output Parameter:
582ee300463SSatish Balay .  scale - scale factor
58362c23b28SMark 
58462c23b28SMark    Level: intermediate
58562c23b28SMark 
586*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetClip()`
58762c23b28SMark @*/
588d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt adapt, PetscReal *scale)
589d71ae5a4SJacob Faibussowitsch {
59062c23b28SMark   PetscFunctionBegin;
59162c23b28SMark   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
59262c23b28SMark   if (scale) PetscValidRealPointer(scale, 2);
59362c23b28SMark   if (scale) *scale = adapt->scale_solve_failed;
5943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
59562c23b28SMark }
59662c23b28SMark 
59762c23b28SMark /*@
598bcf0153eSBarry Smith    TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the time step controller
5991917a363SLisandro Dalcin 
60020f4b53cSBarry Smith    Logically Collective
6011c3436cfSJed Brown 
6024165533cSJose E. Roman    Input Parameters:
603bcf0153eSBarry Smith +  adapt - time step adaptivity context, usually gotten with `TSGetAdapt()`
6041c3436cfSJed Brown .  hmin - minimum time step
6051c3436cfSJed Brown -  hmax - maximum time step
6061c3436cfSJed Brown 
6071c3436cfSJed Brown    Options Database Keys:
608ec18b7bcSLisandro Dalcin +  -ts_adapt_dt_min <min> - to set minimum time step
609ec18b7bcSLisandro Dalcin -  -ts_adapt_dt_max <max> - to set maximum time step
6101c3436cfSJed Brown 
6111c3436cfSJed Brown    Level: intermediate
6121c3436cfSJed Brown 
613*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptGetStepLimits()`, `TSAdaptChoose()`
6141c3436cfSJed Brown @*/
615d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt, PetscReal hmin, PetscReal hmax)
616d71ae5a4SJacob Faibussowitsch {
6171c3436cfSJed Brown   PetscFunctionBegin;
6184782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
619b1f5bccaSLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt, hmin, 2);
620b1f5bccaSLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt, hmax, 3);
62113bcc0bdSJacob Faibussowitsch   PetscCheck(hmin == (PetscReal)PETSC_DEFAULT || hmin >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum time step %g must be non negative", (double)hmin);
62213bcc0bdSJacob Faibussowitsch   PetscCheck(hmax == (PetscReal)PETSC_DEFAULT || hmax >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum time step %g must be non negative", (double)hmax);
62313bcc0bdSJacob Faibussowitsch   if (hmin != (PetscReal)PETSC_DEFAULT) adapt->dt_min = hmin;
62413bcc0bdSJacob Faibussowitsch   if (hmax != (PetscReal)PETSC_DEFAULT) adapt->dt_max = hmax;
625b1f5bccaSLisandro Dalcin   hmin = adapt->dt_min;
626b1f5bccaSLisandro Dalcin   hmax = adapt->dt_max;
6273c633725SBarry Smith   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);
6283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6291917a363SLisandro Dalcin }
6301917a363SLisandro Dalcin 
6311917a363SLisandro Dalcin /*@
632bcf0153eSBarry Smith    TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the time step controller
6331917a363SLisandro Dalcin 
6341917a363SLisandro Dalcin    Not Collective
6351917a363SLisandro Dalcin 
6364165533cSJose E. Roman    Input Parameter:
637bcf0153eSBarry Smith .  adapt - time step adaptivity context, usually gotten with `TSGetAdapt()`
6381917a363SLisandro Dalcin 
6394165533cSJose E. Roman    Output Parameters:
6401917a363SLisandro Dalcin +  hmin - minimum time step
6411917a363SLisandro Dalcin -  hmax - maximum time step
6421917a363SLisandro Dalcin 
6431917a363SLisandro Dalcin    Level: intermediate
6441917a363SLisandro Dalcin 
645*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptSetStepLimits()`, `TSAdaptChoose()`
6461917a363SLisandro Dalcin @*/
647d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt, PetscReal *hmin, PetscReal *hmax)
648d71ae5a4SJacob Faibussowitsch {
6491917a363SLisandro Dalcin   PetscFunctionBegin;
6501917a363SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
6511917a363SLisandro Dalcin   if (hmin) PetscValidRealPointer(hmin, 2);
6521917a363SLisandro Dalcin   if (hmax) PetscValidRealPointer(hmax, 3);
6531917a363SLisandro Dalcin   if (hmin) *hmin = adapt->dt_min;
6541917a363SLisandro Dalcin   if (hmax) *hmax = adapt->dt_max;
6553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6561c3436cfSJed Brown }
6571c3436cfSJed Brown 
658e55864a3SBarry Smith /*
659bcf0153eSBarry Smith    TSAdaptSetFromOptions - Sets various `TSAdapt` parameters from user options.
66084df9cb4SJed Brown 
661c3339decSBarry Smith    Collective
66284df9cb4SJed Brown 
66384df9cb4SJed Brown    Input Parameter:
6642fe279fdSBarry Smith +  adapt - the `TSAdapt` context
6652fe279fdSBarry Smith -  PetscOptionsObject - object created by `PetscOptionsBegin()`
66684df9cb4SJed Brown 
66784df9cb4SJed Brown    Options Database Keys:
6681917a363SLisandro Dalcin +  -ts_adapt_type <type> - algorithm to use for adaptivity
669bf997491SLisandro Dalcin .  -ts_adapt_always_accept - always accept steps regardless of error/stability goals
6701917a363SLisandro Dalcin .  -ts_adapt_safety <safety> - safety factor relative to target error/stability goal
6711917a363SLisandro Dalcin .  -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected
6721917a363SLisandro Dalcin .  -ts_adapt_clip <low,high> - admissible time step decrease and increase factors
67323de1d84SBarry Smith .  -ts_adapt_dt_min <min> - minimum timestep to use
67423de1d84SBarry Smith .  -ts_adapt_dt_max <max> - maximum timestep to use
67523de1d84SBarry Smith .  -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails
676de50f1caSBarry Smith .  -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates
677de50f1caSBarry Smith -  -ts_adapt_time_step_increase_delay - number of timesteps to delay increasing the time step after it has been decreased due to failed solver
67884df9cb4SJed Brown 
67984df9cb4SJed Brown    Level: advanced
68084df9cb4SJed Brown 
681bcf0153eSBarry Smith    Note:
682bcf0153eSBarry Smith    This function is automatically called by `TSSetFromOptions()`
68384df9cb4SJed Brown 
684*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptSetAlwaysAccept()`, `TSAdaptSetSafety()`,
685db781477SPatrick Sanan           `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetStepLimits()`, `TSAdaptSetMonitor()`
686e55864a3SBarry Smith */
687d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetFromOptions(TSAdapt adapt, PetscOptionItems *PetscOptionsObject)
688d71ae5a4SJacob Faibussowitsch {
68984df9cb4SJed Brown   char      type[256] = TSADAPTBASIC;
69062c23b28SMark   PetscReal safety, reject_safety, clip[2], scale, hmin, hmax;
6911c3436cfSJed Brown   PetscBool set, flg;
6921917a363SLisandro Dalcin   PetscInt  two;
69384df9cb4SJed Brown 
69484df9cb4SJed Brown   PetscFunctionBegin;
695dbbe0bcdSBarry Smith   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
69684df9cb4SJed Brown   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
6971566a47fSLisandro Dalcin    * function can only be called from inside TSSetFromOptions()  */
698d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "TS Adaptivity options");
6999566063dSJacob Faibussowitsch   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));
70048a46eb9SPierre Jolivet   if (flg || !((PetscObject)adapt)->type_name) PetscCall(TSAdaptSetType(adapt, type));
7011917a363SLisandro Dalcin 
7029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_adapt_always_accept", "Always accept the step", "TSAdaptSetAlwaysAccept", adapt->always_accept, &flg, &set));
7039566063dSJacob Faibussowitsch   if (set) PetscCall(TSAdaptSetAlwaysAccept(adapt, flg));
7041917a363SLisandro Dalcin 
7059371c9d4SSatish Balay   safety        = adapt->safety;
7069371c9d4SSatish Balay   reject_safety = adapt->reject_safety;
7079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_adapt_safety", "Safety factor relative to target error/stability goal", "TSAdaptSetSafety", safety, &safety, &set));
7089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_adapt_reject_safety", "Extra safety factor to apply if the last step was rejected", "TSAdaptSetSafety", reject_safety, &reject_safety, &flg));
7099566063dSJacob Faibussowitsch   if (set || flg) PetscCall(TSAdaptSetSafety(adapt, safety, reject_safety));
7101917a363SLisandro Dalcin 
7119371c9d4SSatish Balay   two     = 2;
7129371c9d4SSatish Balay   clip[0] = adapt->clip[0];
7139371c9d4SSatish Balay   clip[1] = adapt->clip[1];
7149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRealArray("-ts_adapt_clip", "Admissible decrease/increase factor in step size", "TSAdaptSetClip", clip, &two, &set));
7153c633725SBarry Smith   PetscCheck(!set || (two == 2), PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Must give exactly two values to -ts_adapt_clip");
7169566063dSJacob Faibussowitsch   if (set) PetscCall(TSAdaptSetClip(adapt, clip[0], clip[1]));
7171917a363SLisandro Dalcin 
7189371c9d4SSatish Balay   hmin = adapt->dt_min;
7199371c9d4SSatish Balay   hmax = adapt->dt_max;
7209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_adapt_dt_min", "Minimum time step considered", "TSAdaptSetStepLimits", hmin, &hmin, &set));
7219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_adapt_dt_max", "Maximum time step considered", "TSAdaptSetStepLimits", hmax, &hmax, &flg));
7229566063dSJacob Faibussowitsch   if (set || flg) PetscCall(TSAdaptSetStepLimits(adapt, hmin, hmax));
7231917a363SLisandro Dalcin 
7249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_adapt_max_ignore", "Adaptor ignores (absolute) solution values smaller than this value", "", adapt->ignore_max, &adapt->ignore_max, &set));
7259566063dSJacob Faibussowitsch   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));
726d580f011SEmil Constantinescu 
7279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_adapt_scale_solve_failed", "Scale step by this factor if solve fails", "TSAdaptSetScaleSolveFailed", adapt->scale_solve_failed, &scale, &set));
7289566063dSJacob Faibussowitsch   if (set) PetscCall(TSAdaptSetScaleSolveFailed(adapt, scale));
7291917a363SLisandro Dalcin 
7309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-ts_adapt_wnormtype", "Type of norm computed for error estimation", "", NormTypes, (PetscEnum)adapt->wnormtype, (PetscEnum *)&adapt->wnormtype, NULL));
7313c633725SBarry Smith   PetscCheck(adapt->wnormtype == NORM_2 || adapt->wnormtype == NORM_INFINITY, PetscObjectComm((PetscObject)adapt), PETSC_ERR_SUP, "Only 2-norm and infinite norm supported");
7321917a363SLisandro Dalcin 
7339566063dSJacob Faibussowitsch   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));
734de50f1caSBarry Smith 
7359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_adapt_monitor", "Print choices made by adaptive controller", "TSAdaptSetMonitor", adapt->monitor ? PETSC_TRUE : PETSC_FALSE, &flg, &set));
7369566063dSJacob Faibussowitsch   if (set) PetscCall(TSAdaptSetMonitor(adapt, flg));
7371917a363SLisandro Dalcin 
738dbbe0bcdSBarry Smith   PetscTryTypeMethod(adapt, setfromoptions, PetscOptionsObject);
739d0609cedSBarry Smith   PetscOptionsHeadEnd();
7403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
74184df9cb4SJed Brown }
74284df9cb4SJed Brown 
74384df9cb4SJed Brown /*@
74484df9cb4SJed Brown    TSAdaptCandidatesClear - clear any previously set candidate schemes
74584df9cb4SJed Brown 
74620f4b53cSBarry Smith    Logically Collective
74784df9cb4SJed Brown 
7484165533cSJose E. Roman    Input Parameter:
74984df9cb4SJed Brown .  adapt - adaptive controller
75084df9cb4SJed Brown 
75184df9cb4SJed Brown    Level: developer
75284df9cb4SJed Brown 
753*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptCreate()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()`
75484df9cb4SJed Brown @*/
755d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
756d71ae5a4SJacob Faibussowitsch {
75784df9cb4SJed Brown   PetscFunctionBegin;
7584782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
7599566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&adapt->candidates, sizeof(adapt->candidates)));
7603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
76184df9cb4SJed Brown }
76284df9cb4SJed Brown 
76384df9cb4SJed Brown /*@C
76484df9cb4SJed Brown    TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from
76584df9cb4SJed Brown 
76620f4b53cSBarry Smith    Logically Collective; No Fortran Support
76784df9cb4SJed Brown 
7684165533cSJose E. Roman    Input Parameters:
769bcf0153eSBarry Smith +  adapt - time step adaptivity context, obtained with `TSGetAdapt()` or `TSAdaptCreate()`
77084df9cb4SJed Brown .  name - name of the candidate scheme to add
77184df9cb4SJed Brown .  order - order of the candidate scheme
77284df9cb4SJed Brown .  stageorder - stage order of the candidate scheme
7738d59e960SJed Brown .  ccfl - stability coefficient relative to explicit Euler, used for CFL constraints
77484df9cb4SJed Brown .  cost - relative measure of the amount of work required for the candidate scheme
77584df9cb4SJed Brown -  inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme
77684df9cb4SJed Brown 
77784df9cb4SJed Brown    Level: developer
77884df9cb4SJed Brown 
779*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptChoose()`
78084df9cb4SJed Brown @*/
781d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt, const char name[], PetscInt order, PetscInt stageorder, PetscReal ccfl, PetscReal cost, PetscBool inuse)
782d71ae5a4SJacob Faibussowitsch {
78384df9cb4SJed Brown   PetscInt c;
78484df9cb4SJed Brown 
78584df9cb4SJed Brown   PetscFunctionBegin;
78684df9cb4SJed Brown   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
78763a3b9bcSJacob Faibussowitsch   PetscCheck(order >= 1, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Classical order %" PetscInt_FMT " must be a positive integer", order);
78884df9cb4SJed Brown   if (inuse) {
7893c633725SBarry Smith     PetscCheck(!adapt->candidates.inuse_set, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONGSTATE, "Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()");
79084df9cb4SJed Brown     adapt->candidates.inuse_set = PETSC_TRUE;
79184df9cb4SJed Brown   }
7921c3436cfSJed Brown   /* first slot if this is the current scheme, otherwise the next available slot */
7931c3436cfSJed Brown   c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n;
794bbd56ea5SKarl Rupp 
79584df9cb4SJed Brown   adapt->candidates.name[c]       = name;
79684df9cb4SJed Brown   adapt->candidates.order[c]      = order;
79784df9cb4SJed Brown   adapt->candidates.stageorder[c] = stageorder;
7988d59e960SJed Brown   adapt->candidates.ccfl[c]       = ccfl;
79984df9cb4SJed Brown   adapt->candidates.cost[c]       = cost;
80084df9cb4SJed Brown   adapt->candidates.n++;
8013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
80284df9cb4SJed Brown }
80384df9cb4SJed Brown 
8048d59e960SJed Brown /*@C
8058d59e960SJed Brown    TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost
8068d59e960SJed Brown 
8078d59e960SJed Brown    Not Collective
8088d59e960SJed Brown 
8094165533cSJose E. Roman    Input Parameter:
8108d59e960SJed Brown .  adapt - time step adaptivity context
8118d59e960SJed Brown 
8124165533cSJose E. Roman    Output Parameters:
8138d59e960SJed Brown +  n - number of candidate schemes, always at least 1
8148d59e960SJed Brown .  order - the order of each candidate scheme
8158d59e960SJed Brown .  stageorder - the stage order of each candidate scheme
8168d59e960SJed Brown .  ccfl - the CFL coefficient of each scheme
8178d59e960SJed Brown -  cost - the relative cost of each scheme
8188d59e960SJed Brown 
8198d59e960SJed Brown    Level: developer
8208d59e960SJed Brown 
8218d59e960SJed Brown    Note:
8228d59e960SJed Brown    The current scheme is always returned in the first slot
8238d59e960SJed Brown 
824*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()`
8258d59e960SJed Brown @*/
826d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt, PetscInt *n, const PetscInt **order, const PetscInt **stageorder, const PetscReal **ccfl, const PetscReal **cost)
827d71ae5a4SJacob Faibussowitsch {
8288d59e960SJed Brown   PetscFunctionBegin;
8298d59e960SJed Brown   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
8308d59e960SJed Brown   if (n) *n = adapt->candidates.n;
8318d59e960SJed Brown   if (order) *order = adapt->candidates.order;
8328d59e960SJed Brown   if (stageorder) *stageorder = adapt->candidates.stageorder;
8338d59e960SJed Brown   if (ccfl) *ccfl = adapt->candidates.ccfl;
8348d59e960SJed Brown   if (cost) *cost = adapt->candidates.cost;
8353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8368d59e960SJed Brown }
8378d59e960SJed Brown 
83884df9cb4SJed Brown /*@C
83984df9cb4SJed Brown    TSAdaptChoose - choose which method and step size to use for the next step
84084df9cb4SJed Brown 
841c3339decSBarry Smith    Collective
84284df9cb4SJed Brown 
8434165533cSJose E. Roman    Input Parameters:
844da81f932SPierre Jolivet +  adapt - adaptive controller
84597bb3fdcSJose E. Roman .  ts - time stepper
84684df9cb4SJed Brown -  h - current step size
84784df9cb4SJed Brown 
8484165533cSJose E. Roman    Output Parameters:
8491566a47fSLisandro Dalcin +  next_sc - optional, scheme to use for the next step
85084df9cb4SJed Brown .  next_h - step size to use for the next step
851bcf0153eSBarry Smith -  accept - `PETSC_TRUE` to accept the current step, `PETSC_FALSE` to repeat the current step with the new step size
8521c3436cfSJed Brown 
85384df9cb4SJed Brown    Level: developer
85484df9cb4SJed Brown 
855bcf0153eSBarry Smith    Note:
856bcf0153eSBarry Smith    The input value of parameter accept is retained from the last time step, so it will be `PETSC_FALSE` if the step is
857bcf0153eSBarry Smith    being retried after an initial rejection.
858bcf0153eSBarry Smith 
859*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()`
86084df9cb4SJed Brown @*/
861d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptChoose(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept)
862d71ae5a4SJacob Faibussowitsch {
8631566a47fSLisandro Dalcin   PetscInt  ncandidates = adapt->candidates.n;
8641566a47fSLisandro Dalcin   PetscInt  scheme      = 0;
8650b99f514SJed Brown   PetscReal wlte        = -1.0;
8665e4ed32fSEmil Constantinescu   PetscReal wltea       = -1.0;
8675e4ed32fSEmil Constantinescu   PetscReal wlter       = -1.0;
86884df9cb4SJed Brown 
86984df9cb4SJed Brown   PetscFunctionBegin;
87084df9cb4SJed Brown   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
87184df9cb4SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 2);
8721566a47fSLisandro Dalcin   if (next_sc) PetscValidIntPointer(next_sc, 4);
873dadcf809SJacob Faibussowitsch   PetscValidRealPointer(next_h, 5);
874064a246eSJacob Faibussowitsch   PetscValidBoolPointer(accept, 6);
8751566a47fSLisandro Dalcin   if (next_sc) *next_sc = 0;
8761566a47fSLisandro Dalcin 
8771566a47fSLisandro Dalcin   /* Do not mess with adaptivity while handling events*/
8781566a47fSLisandro Dalcin   if (ts->event && ts->event->status != TSEVENT_NONE) {
8791566a47fSLisandro Dalcin     *next_h = h;
8801566a47fSLisandro Dalcin     *accept = PETSC_TRUE;
8813ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
8821566a47fSLisandro Dalcin   }
8831566a47fSLisandro Dalcin 
884dbbe0bcdSBarry Smith   PetscUseTypeMethod(adapt, choose, ts, h, &scheme, next_h, accept, &wlte, &wltea, &wlter);
88563a3b9bcSJacob Faibussowitsch   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);
8863c633725SBarry Smith   PetscCheck(*next_h >= 0, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Computed step size %g must be positive", (double)*next_h);
8871566a47fSLisandro Dalcin   if (next_sc) *next_sc = scheme;
8881566a47fSLisandro Dalcin 
88949354f04SShri Abhyankar   if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
89036b54a69SLisandro Dalcin     /* Increase/reduce step size if end time of next step is close to or overshoots max time */
89136b54a69SLisandro Dalcin     PetscReal t = ts->ptime + ts->time_step, h = *next_h;
8924a658b32SHong Zhang     PetscReal tend = t + h, tmax, hmax;
893fe7350e0SStefano Zampini     PetscReal a    = (PetscReal)(1.0 + adapt->matchstepfac[0]);
894fe7350e0SStefano Zampini     PetscReal b    = adapt->matchstepfac[1];
8954a658b32SHong Zhang 
8964a658b32SHong Zhang     if (ts->tspan) {
897e1db57b0SHong Zhang       if (PetscIsCloseAtTol(t, ts->tspan->span_times[ts->tspan->spanctr], ts->tspan->reltol * h + ts->tspan->abstol, 0)) /* hit a span time point */
8984a658b32SHong Zhang         if (ts->tspan->spanctr + 1 < ts->tspan->num_span_times) tmax = ts->tspan->span_times[ts->tspan->spanctr + 1];
8994a658b32SHong Zhang         else tmax = ts->max_time; /* hit the last span time point */
9004a658b32SHong Zhang       else tmax = ts->tspan->span_times[ts->tspan->spanctr];
9014a658b32SHong Zhang     } else tmax = ts->max_time;
9024a658b32SHong Zhang     hmax = tmax - t;
90336b54a69SLisandro Dalcin     if (t < tmax && tend > tmax) *next_h = hmax;
904fe7350e0SStefano Zampini     if (t < tmax && tend < tmax && h * b > hmax) *next_h = hmax / 2;
90536b54a69SLisandro Dalcin     if (t < tmax && tend < tmax && h * a > hmax) *next_h = hmax;
906e1db57b0SHong Zhang     /* if step size is changed to match a span time point */
907e1db57b0SHong Zhang     if (ts->tspan && h != *next_h && !adapt->dt_span_cached) adapt->dt_span_cached = h;
908e1db57b0SHong Zhang     /* reset time step after a span time point */
909e1db57b0SHong Zhang     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)) {
910e1db57b0SHong Zhang       *next_h               = adapt->dt_span_cached;
911e1db57b0SHong Zhang       adapt->dt_span_cached = 0;
91249354f04SShri Abhyankar     }
913e1db57b0SHong Zhang   }
9141c3436cfSJed Brown   if (adapt->monitor) {
9151566a47fSLisandro Dalcin     const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : "";
9169566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
9170b99f514SJed Brown     if (wlte < 0) {
9189371c9d4SSatish Balay       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",
9199371c9d4SSatish Balay                                        (double)ts->ptime, (double)h, (double)*next_h));
9200b99f514SJed Brown     } else {
9219371c9d4SSatish Balay       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",
9229371c9d4SSatish Balay                                        (double)ts->ptime, (double)h, (double)*next_h, (double)wlte, (double)wltea, (double)wlter));
9230b99f514SJed Brown     }
9249566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
9251c3436cfSJed Brown   }
9263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
92784df9cb4SJed Brown }
92884df9cb4SJed Brown 
92997335746SJed Brown /*@
930de50f1caSBarry Smith    TSAdaptSetTimeStepIncreaseDelay - The number of timesteps to wait after a decrease in the timestep due to failed solver
931de50f1caSBarry Smith                                      before increasing the time step.
932de50f1caSBarry Smith 
933c3339decSBarry Smith    Logicially Collective
934de50f1caSBarry Smith 
9354165533cSJose E. Roman    Input Parameters:
936de50f1caSBarry Smith +  adapt - adaptive controller context
937de50f1caSBarry Smith -  cnt - the number of timesteps
938de50f1caSBarry Smith 
939de50f1caSBarry Smith    Options Database Key:
940de50f1caSBarry Smith .  -ts_adapt_time_step_increase_delay cnt - number of steps to delay the increase
941de50f1caSBarry Smith 
942de50f1caSBarry Smith    Level: advanced
943de50f1caSBarry Smith 
944bcf0153eSBarry Smith    Notes:
945bcf0153eSBarry Smith    This is to prevent an adaptor from bouncing back and forth between two nearby timesteps. The default is 0.
946bcf0153eSBarry Smith 
947bcf0153eSBarry Smith    The successful use of this option is problem dependent
948bcf0153eSBarry Smith 
949bcf0153eSBarry Smith    Developer Note:
950bcf0153eSBarry Smith    There is no theory to support this option
951bcf0153eSBarry Smith 
952*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`
953de50f1caSBarry Smith @*/
954d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt, PetscInt cnt)
955d71ae5a4SJacob Faibussowitsch {
956de50f1caSBarry Smith   PetscFunctionBegin;
957de50f1caSBarry Smith   adapt->timestepjustdecreased_delay = cnt;
9583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
959de50f1caSBarry Smith }
960de50f1caSBarry Smith 
961de50f1caSBarry Smith /*@
9626bc98fa9SBarry Smith    TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails or solution vector is infeasible)
96397335746SJed Brown 
964c3339decSBarry Smith    Collective
96597335746SJed Brown 
9664165533cSJose E. Roman    Input Parameters:
96797335746SJed Brown +  adapt - adaptive controller context
968b295832fSPierre Barbier de Reuille .  ts - time stepper
969b295832fSPierre Barbier de Reuille .  t - Current simulation time
970b295832fSPierre Barbier de Reuille -  Y - Current solution vector
97197335746SJed Brown 
9724165533cSJose E. Roman    Output Parameter:
973bcf0153eSBarry Smith .  accept - `PETSC_TRUE` to accept the stage, `PETSC_FALSE` to reject
97497335746SJed Brown 
97597335746SJed Brown    Level: developer
97697335746SJed Brown 
977*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`
97897335746SJed Brown @*/
979d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCheckStage(TSAdapt adapt, TS ts, PetscReal t, Vec Y, PetscBool *accept)
980d71ae5a4SJacob Faibussowitsch {
9811566a47fSLisandro Dalcin   SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING;
98297335746SJed Brown 
98397335746SJed Brown   PetscFunctionBegin;
9844782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
9854782b174SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 2);
986064a246eSJacob Faibussowitsch   PetscValidBoolPointer(accept, 5);
9871566a47fSLisandro Dalcin 
9889566063dSJacob Faibussowitsch   if (ts->snes) PetscCall(SNESGetConvergedReason(ts->snes, &snesreason));
98997335746SJed Brown   if (snesreason < 0) {
99097335746SJed Brown     *accept = PETSC_FALSE;
9916de24e2aSJed Brown     if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) {
99297335746SJed Brown       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
99363a3b9bcSJacob Faibussowitsch       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));
99497335746SJed Brown       if (adapt->monitor) {
9959566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
9969371c9d4SSatish Balay         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,
9979371c9d4SSatish Balay                                          (double)ts->ptime, (double)ts->time_step, ts->num_snes_failures));
9989566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
99997335746SJed Brown       }
1000cb9d8021SPierre Barbier de Reuille     }
1001cb9d8021SPierre Barbier de Reuille   } else {
10021566a47fSLisandro Dalcin     *accept = PETSC_TRUE;
10039566063dSJacob Faibussowitsch     PetscCall(TSFunctionDomainError(ts, t, Y, accept));
1004cb9d8021SPierre Barbier de Reuille     if (*accept && adapt->checkstage) {
10059566063dSJacob Faibussowitsch       PetscCall((*adapt->checkstage)(adapt, ts, t, Y, accept));
10066bc98fa9SBarry Smith       if (!*accept) {
100763a3b9bcSJacob Faibussowitsch         PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", solution rejected by user function provided by TSSetFunctionDomainError()\n", ts->steps));
10086bc98fa9SBarry Smith         if (adapt->monitor) {
10099566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
101063a3b9bcSJacob Faibussowitsch           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));
10119566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
10126bc98fa9SBarry Smith         }
10136bc98fa9SBarry Smith       }
1014cb9d8021SPierre Barbier de Reuille     }
1015cb9d8021SPierre Barbier de Reuille   }
1016cb9d8021SPierre Barbier de Reuille 
10171566a47fSLisandro Dalcin   if (!(*accept) && !ts->reason) {
10181566a47fSLisandro Dalcin     PetscReal dt, new_dt;
10199566063dSJacob Faibussowitsch     PetscCall(TSGetTimeStep(ts, &dt));
1020cb9d8021SPierre Barbier de Reuille     new_dt = dt * adapt->scale_solve_failed;
10219566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts, new_dt));
1022de50f1caSBarry Smith     adapt->timestepjustdecreased += adapt->timestepjustdecreased_delay;
102397335746SJed Brown     if (adapt->monitor) {
10249566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
102563a3b9bcSJacob Faibussowitsch       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));
10269566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
102797335746SJed Brown     }
102897335746SJed Brown   }
10293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
103097335746SJed Brown }
103197335746SJed Brown 
103284df9cb4SJed Brown /*@
103384df9cb4SJed Brown   TSAdaptCreate - create an adaptive controller context for time stepping
103484df9cb4SJed Brown 
1035d083f849SBarry Smith   Collective
103684df9cb4SJed Brown 
103784df9cb4SJed Brown   Input Parameter:
103884df9cb4SJed Brown . comm - The communicator
103984df9cb4SJed Brown 
104084df9cb4SJed Brown   Output Parameter:
1041bcf0153eSBarry Smith . adapt - new `TSAdapt` object
104284df9cb4SJed Brown 
104384df9cb4SJed Brown   Level: developer
104484df9cb4SJed Brown 
1045bcf0153eSBarry Smith   Note:
1046bcf0153eSBarry Smith   `TSAdapt` creation is handled by `TS`, so users should not need to call this function.
104784df9cb4SJed Brown 
1048*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptDestroy()`
104984df9cb4SJed Brown @*/
1050d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCreate(MPI_Comm comm, TSAdapt *inadapt)
1051d71ae5a4SJacob Faibussowitsch {
105284df9cb4SJed Brown   TSAdapt adapt;
105384df9cb4SJed Brown 
105484df9cb4SJed Brown   PetscFunctionBegin;
1055064a246eSJacob Faibussowitsch   PetscValidPointer(inadapt, 2);
10563b3bcf4cSLisandro Dalcin   *inadapt = NULL;
10579566063dSJacob Faibussowitsch   PetscCall(TSAdaptInitializePackage());
10583b3bcf4cSLisandro Dalcin 
10599566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(adapt, TSADAPT_CLASSID, "TSAdapt", "Time stepping adaptivity", "TS", comm, TSAdaptDestroy, TSAdaptView));
10601c3436cfSJed Brown 
1061bf997491SLisandro Dalcin   adapt->always_accept      = PETSC_FALSE;
10621917a363SLisandro Dalcin   adapt->safety             = 0.9;
10631917a363SLisandro Dalcin   adapt->reject_safety      = 0.5;
10641917a363SLisandro Dalcin   adapt->clip[0]            = 0.1;
10651917a363SLisandro Dalcin   adapt->clip[1]            = 10.;
10661c3436cfSJed Brown   adapt->dt_min             = 1e-20;
10671917a363SLisandro Dalcin   adapt->dt_max             = 1e+20;
10681c167fc2SEmil Constantinescu   adapt->ignore_max         = -1.0;
1069d580f011SEmil Constantinescu   adapt->glee_use_local     = PETSC_TRUE;
107097335746SJed Brown   adapt->scale_solve_failed = 0.25;
1071fe7350e0SStefano Zampini   /* these two safety factors are not public, and they are used only in the TS_EXACTFINALTIME_MATCHSTEP case
1072fe7350e0SStefano Zampini      to prevent from situations were unreasonably small time steps are taken in order to match the final time */
1073fe7350e0SStefano Zampini   adapt->matchstepfac[0]             = 0.01; /* allow 1% step size increase in the last step */
1074fe7350e0SStefano Zampini   adapt->matchstepfac[1]             = 2.0;  /* halve last step if it is greater than what remains divided this factor */
10757619abb3SShri   adapt->wnormtype                   = NORM_2;
1076de50f1caSBarry Smith   adapt->timestepjustdecreased_delay = 0;
10771917a363SLisandro Dalcin 
107884df9cb4SJed Brown   *inadapt = adapt;
10793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
108084df9cb4SJed Brown }
1081