xref: /petsc/src/ts/adapt/interface/tsadapt.c (revision 14d0ab18b70ee075d422f67a0c1395817de67fab)
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 
31b43aa488SJacob Faibussowitsch   Example 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 
411cc06b55SBarry 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 
581cc06b55SBarry 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 
801cc06b55SBarry 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 
971cc06b55SBarry 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 
124b43aa488SJacob Faibussowitsch .seealso: [](ch_ts), `TSGetAdapt()`, `TSAdaptDestroy()`, `TSAdaptType`, `TSAdaptGetType()`
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);
1334f572ea9SToby Isaac   PetscAssertPointer(type, 2);
1349566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)adapt, type, &match));
1353ba16761SJacob Faibussowitsch   if (match) PetscFunctionReturn(PETSC_SUCCESS);
1369566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(TSAdaptList, type, &r));
1376adde796SStefano Zampini   PetscCheck(r, PetscObjectComm((PetscObject)adapt), 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);
1644f572ea9SToby Isaac   PetscAssertPointer(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:
183b43aa488SJacob Faibussowitsch + adapt  - 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 
1931cc06b55SBarry 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 
2671cc06b55SBarry 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 
3091cc06b55SBarry 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 
333b43aa488SJacob Faibussowitsch   Calling sequence:
3340873808bSJed Brown + adapt  - adaptive controller context
3350873808bSJed Brown . ts     - time stepping context
336*14d0ab18SJacob Faibussowitsch . t      - current time
337*14d0ab18SJacob Faibussowitsch . Y      - current solution vector
3380873808bSJed Brown - accept - pending choice of whether to accept, can be modified by this routine
3390873808bSJed Brown 
3400873808bSJed Brown   Level: advanced
3410873808bSJed Brown 
3421cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptChoose()`
3430873808bSJed Brown @*/
344*14d0ab18SJacob Faibussowitsch PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt, PetscErrorCode (*func)(TSAdapt adapt, TS ts, PetscReal t, Vec Y, PetscBool *accept))
345d71ae5a4SJacob Faibussowitsch {
3460873808bSJed Brown   PetscFunctionBegin;
3470873808bSJed Brown   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
34868ae941aSLisandro Dalcin   adapt->checkstage = func;
3493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3500873808bSJed Brown }
3510873808bSJed Brown 
3521c3436cfSJed Brown /*@
353bf997491SLisandro Dalcin   TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of
354bf997491SLisandro Dalcin   any error or stability condition not meeting the prescribed goal.
355bf997491SLisandro Dalcin 
3562fe279fdSBarry Smith   Logically Collective
357bf997491SLisandro Dalcin 
3584165533cSJose E. Roman   Input Parameters:
359bcf0153eSBarry Smith + adapt - time step adaptivity context, usually gotten with `TSGetAdapt()`
360bf997491SLisandro Dalcin - flag  - whether to always accept steps
361bf997491SLisandro Dalcin 
362bcf0153eSBarry Smith   Options Database Key:
363ec18b7bcSLisandro Dalcin . -ts_adapt_always_accept - to always accept steps
364bf997491SLisandro Dalcin 
365bf997491SLisandro Dalcin   Level: intermediate
366bf997491SLisandro Dalcin 
3671cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptChoose()`
368bf997491SLisandro Dalcin @*/
369d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt, PetscBool flag)
370d71ae5a4SJacob Faibussowitsch {
371bf997491SLisandro Dalcin   PetscFunctionBegin;
372bf997491SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
373bf997491SLisandro Dalcin   PetscValidLogicalCollectiveBool(adapt, flag, 2);
374bf997491SLisandro Dalcin   adapt->always_accept = flag;
3753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
376bf997491SLisandro Dalcin }
377bf997491SLisandro Dalcin 
378bf997491SLisandro Dalcin /*@
379bcf0153eSBarry Smith   TSAdaptSetSafety - Set safety factors for time step adaptor
3801c3436cfSJed Brown 
3812fe279fdSBarry Smith   Logically Collective
3821917a363SLisandro Dalcin 
3834165533cSJose E. Roman   Input Parameters:
3841917a363SLisandro Dalcin + adapt         - adaptive controller context
3851917a363SLisandro Dalcin . safety        - safety factor relative to target error/stability goal
386ec18b7bcSLisandro Dalcin - reject_safety - extra safety factor to apply if the last step was rejected
3871917a363SLisandro Dalcin 
3881917a363SLisandro Dalcin   Options Database Keys:
389ec18b7bcSLisandro Dalcin + -ts_adapt_safety <safety>               - to set safety factor
390ec18b7bcSLisandro Dalcin - -ts_adapt_reject_safety <reject_safety> - to set reject safety factor
3911917a363SLisandro Dalcin 
3921917a363SLisandro Dalcin   Level: intermediate
3931917a363SLisandro Dalcin 
3941cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptGetSafety()`, `TSAdaptChoose()`
3951917a363SLisandro Dalcin @*/
396d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetSafety(TSAdapt adapt, PetscReal safety, PetscReal reject_safety)
397d71ae5a4SJacob Faibussowitsch {
3981917a363SLisandro Dalcin   PetscFunctionBegin;
3991917a363SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
4001917a363SLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt, safety, 2);
4011917a363SLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt, reject_safety, 3);
40213bcc0bdSJacob Faibussowitsch   PetscCheck(safety == (PetscReal)PETSC_DEFAULT || safety >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Safety factor %g must be non negative", (double)safety);
40313bcc0bdSJacob 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);
40413bcc0bdSJacob 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);
40513bcc0bdSJacob 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);
40613bcc0bdSJacob Faibussowitsch   if (safety != (PetscReal)PETSC_DEFAULT) adapt->safety = safety;
40713bcc0bdSJacob Faibussowitsch   if (reject_safety != (PetscReal)PETSC_DEFAULT) adapt->reject_safety = reject_safety;
4083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4091917a363SLisandro Dalcin }
4101917a363SLisandro Dalcin 
4111917a363SLisandro Dalcin /*@
412bcf0153eSBarry Smith   TSAdaptGetSafety - Get safety factors for time step adapter
4131917a363SLisandro Dalcin 
4141917a363SLisandro Dalcin   Not Collective
4151917a363SLisandro Dalcin 
4164165533cSJose E. Roman   Input Parameter:
4171917a363SLisandro Dalcin . adapt - adaptive controller context
4181917a363SLisandro Dalcin 
4194165533cSJose E. Roman   Output Parameters:
420b43aa488SJacob Faibussowitsch + safety        - safety factor relative to target error/stability goal
421b43aa488SJacob Faibussowitsch - reject_safety - extra safety factor to apply if the last step was rejected
4221917a363SLisandro Dalcin 
4231917a363SLisandro Dalcin   Level: intermediate
4241917a363SLisandro Dalcin 
4251cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptSetSafety()`, `TSAdaptChoose()`
4261917a363SLisandro Dalcin @*/
427d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetSafety(TSAdapt adapt, PetscReal *safety, PetscReal *reject_safety)
428d71ae5a4SJacob Faibussowitsch {
4291917a363SLisandro Dalcin   PetscFunctionBegin;
4301917a363SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
4314f572ea9SToby Isaac   if (safety) PetscAssertPointer(safety, 2);
4324f572ea9SToby Isaac   if (reject_safety) PetscAssertPointer(reject_safety, 3);
4331917a363SLisandro Dalcin   if (safety) *safety = adapt->safety;
4341917a363SLisandro Dalcin   if (reject_safety) *reject_safety = adapt->reject_safety;
4353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4361917a363SLisandro Dalcin }
4371917a363SLisandro Dalcin 
4381917a363SLisandro Dalcin /*@
439bcf0153eSBarry Smith   TSAdaptSetMaxIgnore - Set error estimation threshold. Solution components below this threshold value will not be considered when computing error norms
440bcf0153eSBarry Smith   for time step adaptivity (in absolute value). A negative value (default) of the threshold leads to considering all solution components.
44176cddca1SEmil Constantinescu 
4422fe279fdSBarry Smith   Logically Collective
44376cddca1SEmil Constantinescu 
4444165533cSJose E. Roman   Input Parameters:
44576cddca1SEmil Constantinescu + adapt      - adaptive controller context
44676cddca1SEmil Constantinescu - max_ignore - threshold for solution components that are ignored during error estimation
44776cddca1SEmil Constantinescu 
448bcf0153eSBarry Smith   Options Database Key:
44976cddca1SEmil Constantinescu . -ts_adapt_max_ignore <max_ignore> - to set the threshold
45076cddca1SEmil Constantinescu 
45176cddca1SEmil Constantinescu   Level: intermediate
45276cddca1SEmil Constantinescu 
4531cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptGetMaxIgnore()`, `TSAdaptChoose()`
45476cddca1SEmil Constantinescu @*/
455d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt adapt, PetscReal max_ignore)
456d71ae5a4SJacob Faibussowitsch {
45776cddca1SEmil Constantinescu   PetscFunctionBegin;
45876cddca1SEmil Constantinescu   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
45976cddca1SEmil Constantinescu   PetscValidLogicalCollectiveReal(adapt, max_ignore, 2);
46076cddca1SEmil Constantinescu   adapt->ignore_max = max_ignore;
4613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
46276cddca1SEmil Constantinescu }
46376cddca1SEmil Constantinescu 
46476cddca1SEmil Constantinescu /*@
465bcf0153eSBarry Smith   TSAdaptGetMaxIgnore - Get error estimation threshold. Solution components below this threshold value will not be considered when computing error norms
466bcf0153eSBarry Smith   for time step adaptivity (in absolute value).
46776cddca1SEmil Constantinescu 
46876cddca1SEmil Constantinescu   Not Collective
46976cddca1SEmil Constantinescu 
4704165533cSJose E. Roman   Input Parameter:
47176cddca1SEmil Constantinescu . adapt - adaptive controller context
47276cddca1SEmil Constantinescu 
4734165533cSJose E. Roman   Output Parameter:
47476cddca1SEmil Constantinescu . max_ignore - threshold for solution components that are ignored during error estimation
47576cddca1SEmil Constantinescu 
47676cddca1SEmil Constantinescu   Level: intermediate
47776cddca1SEmil Constantinescu 
4781cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptSetMaxIgnore()`, `TSAdaptChoose()`
47976cddca1SEmil Constantinescu @*/
480d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt adapt, PetscReal *max_ignore)
481d71ae5a4SJacob Faibussowitsch {
48276cddca1SEmil Constantinescu   PetscFunctionBegin;
48376cddca1SEmil Constantinescu   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
4844f572ea9SToby Isaac   PetscAssertPointer(max_ignore, 2);
48576cddca1SEmil Constantinescu   *max_ignore = adapt->ignore_max;
4863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
48776cddca1SEmil Constantinescu }
48876cddca1SEmil Constantinescu 
48976cddca1SEmil Constantinescu /*@
490bcf0153eSBarry Smith   TSAdaptSetClip - Sets the admissible decrease/increase factor in step size in the time step adapter
4911917a363SLisandro Dalcin 
492c3339decSBarry Smith   Logically collective
4931917a363SLisandro Dalcin 
4944165533cSJose E. Roman   Input Parameters:
4951917a363SLisandro Dalcin + adapt - adaptive controller context
4961917a363SLisandro Dalcin . low   - admissible decrease factor
497ec18b7bcSLisandro Dalcin - high  - admissible increase factor
4981917a363SLisandro Dalcin 
499bcf0153eSBarry Smith   Options Database Key:
500ec18b7bcSLisandro Dalcin . -ts_adapt_clip <low>,<high> - to set admissible time step decrease and increase factors
5011917a363SLisandro Dalcin 
5021917a363SLisandro Dalcin   Level: intermediate
5031917a363SLisandro Dalcin 
5041cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptGetClip()`, `TSAdaptSetScaleSolveFailed()`
5051917a363SLisandro Dalcin @*/
506d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetClip(TSAdapt adapt, PetscReal low, PetscReal high)
507d71ae5a4SJacob Faibussowitsch {
5081917a363SLisandro Dalcin   PetscFunctionBegin;
5091917a363SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
5101917a363SLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt, low, 2);
5111917a363SLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt, high, 3);
51213bcc0bdSJacob Faibussowitsch   PetscCheck(low == (PetscReal)PETSC_DEFAULT || low >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Decrease factor %g must be non negative", (double)low);
51313bcc0bdSJacob 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);
51413bcc0bdSJacob 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);
51513bcc0bdSJacob Faibussowitsch   if (low != (PetscReal)PETSC_DEFAULT) adapt->clip[0] = low;
51613bcc0bdSJacob Faibussowitsch   if (high != (PetscReal)PETSC_DEFAULT) adapt->clip[1] = high;
5173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5181917a363SLisandro Dalcin }
5191917a363SLisandro Dalcin 
5201917a363SLisandro Dalcin /*@
521bcf0153eSBarry Smith   TSAdaptGetClip - Gets the admissible decrease/increase factor in step size in the time step adapter
5221917a363SLisandro Dalcin 
5231917a363SLisandro Dalcin   Not Collective
5241917a363SLisandro Dalcin 
5254165533cSJose E. Roman   Input Parameter:
5261917a363SLisandro Dalcin . adapt - adaptive controller context
5271917a363SLisandro Dalcin 
5284165533cSJose E. Roman   Output Parameters:
5291917a363SLisandro Dalcin + low  - optional, admissible decrease factor
5301917a363SLisandro Dalcin - high - optional, admissible increase factor
5311917a363SLisandro Dalcin 
5321917a363SLisandro Dalcin   Level: intermediate
5331917a363SLisandro Dalcin 
5341cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()`
5351917a363SLisandro Dalcin @*/
536d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetClip(TSAdapt adapt, PetscReal *low, PetscReal *high)
537d71ae5a4SJacob Faibussowitsch {
5381917a363SLisandro Dalcin   PetscFunctionBegin;
5391917a363SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
5404f572ea9SToby Isaac   if (low) PetscAssertPointer(low, 2);
5414f572ea9SToby Isaac   if (high) PetscAssertPointer(high, 3);
5421917a363SLisandro Dalcin   if (low) *low = adapt->clip[0];
5431917a363SLisandro Dalcin   if (high) *high = adapt->clip[1];
5443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5451917a363SLisandro Dalcin }
5461917a363SLisandro Dalcin 
5471917a363SLisandro Dalcin /*@
548bcf0153eSBarry Smith   TSAdaptSetScaleSolveFailed - Scale step size by this factor if solve fails
54962c23b28SMark 
55020f4b53cSBarry Smith   Logically Collective
55162c23b28SMark 
5524165533cSJose E. Roman   Input Parameters:
55362c23b28SMark + adapt - adaptive controller context
554ee300463SSatish Balay - scale - scale
55562c23b28SMark 
556bcf0153eSBarry Smith   Options Database Key:
55762c23b28SMark . -ts_adapt_scale_solve_failed <scale> - to set scale step by this factor if solve fails
55862c23b28SMark 
55962c23b28SMark   Level: intermediate
56062c23b28SMark 
5611cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptGetScaleSolveFailed()`, `TSAdaptGetClip()`
56262c23b28SMark @*/
563d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt adapt, PetscReal scale)
564d71ae5a4SJacob Faibussowitsch {
56562c23b28SMark   PetscFunctionBegin;
56662c23b28SMark   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
56762c23b28SMark   PetscValidLogicalCollectiveReal(adapt, scale, 2);
56813bcc0bdSJacob Faibussowitsch   PetscCheck(scale == (PetscReal)PETSC_DEFAULT || scale > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Scale factor %g must be positive", (double)scale);
56913bcc0bdSJacob 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);
57013bcc0bdSJacob Faibussowitsch   if (scale != (PetscReal)PETSC_DEFAULT) adapt->scale_solve_failed = scale;
5713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
57262c23b28SMark }
57362c23b28SMark 
57462c23b28SMark /*@
57562c23b28SMark   TSAdaptGetScaleSolveFailed - Gets the admissible decrease/increase factor in step size
57662c23b28SMark 
57762c23b28SMark   Not Collective
57862c23b28SMark 
5794165533cSJose E. Roman   Input Parameter:
58062c23b28SMark . adapt - adaptive controller context
58162c23b28SMark 
5824165533cSJose E. Roman   Output Parameter:
583ee300463SSatish Balay . scale - scale factor
58462c23b28SMark 
58562c23b28SMark   Level: intermediate
58662c23b28SMark 
5871cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetClip()`
58862c23b28SMark @*/
589d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt adapt, PetscReal *scale)
590d71ae5a4SJacob Faibussowitsch {
59162c23b28SMark   PetscFunctionBegin;
59262c23b28SMark   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
5934f572ea9SToby Isaac   if (scale) PetscAssertPointer(scale, 2);
59462c23b28SMark   if (scale) *scale = adapt->scale_solve_failed;
5953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
59662c23b28SMark }
59762c23b28SMark 
59862c23b28SMark /*@
599bcf0153eSBarry Smith   TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the time step controller
6001917a363SLisandro Dalcin 
60120f4b53cSBarry Smith   Logically Collective
6021c3436cfSJed Brown 
6034165533cSJose E. Roman   Input Parameters:
604bcf0153eSBarry Smith + adapt - time step adaptivity context, usually gotten with `TSGetAdapt()`
6051c3436cfSJed Brown . hmin  - minimum time step
6061c3436cfSJed Brown - hmax  - maximum time step
6071c3436cfSJed Brown 
6081c3436cfSJed Brown   Options Database Keys:
609ec18b7bcSLisandro Dalcin + -ts_adapt_dt_min <min> - to set minimum time step
610ec18b7bcSLisandro Dalcin - -ts_adapt_dt_max <max> - to set maximum time step
6111c3436cfSJed Brown 
6121c3436cfSJed Brown   Level: intermediate
6131c3436cfSJed Brown 
6141cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptGetStepLimits()`, `TSAdaptChoose()`
6151c3436cfSJed Brown @*/
616d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt, PetscReal hmin, PetscReal hmax)
617d71ae5a4SJacob Faibussowitsch {
6181c3436cfSJed Brown   PetscFunctionBegin;
6194782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
620b1f5bccaSLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt, hmin, 2);
621b1f5bccaSLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt, hmax, 3);
62213bcc0bdSJacob 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);
62313bcc0bdSJacob 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);
62413bcc0bdSJacob Faibussowitsch   if (hmin != (PetscReal)PETSC_DEFAULT) adapt->dt_min = hmin;
62513bcc0bdSJacob Faibussowitsch   if (hmax != (PetscReal)PETSC_DEFAULT) adapt->dt_max = hmax;
626b1f5bccaSLisandro Dalcin   hmin = adapt->dt_min;
627b1f5bccaSLisandro Dalcin   hmax = adapt->dt_max;
6283c633725SBarry 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);
6293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6301917a363SLisandro Dalcin }
6311917a363SLisandro Dalcin 
6321917a363SLisandro Dalcin /*@
633bcf0153eSBarry Smith   TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the time step controller
6341917a363SLisandro Dalcin 
6351917a363SLisandro Dalcin   Not Collective
6361917a363SLisandro Dalcin 
6374165533cSJose E. Roman   Input Parameter:
638bcf0153eSBarry Smith . adapt - time step adaptivity context, usually gotten with `TSGetAdapt()`
6391917a363SLisandro Dalcin 
6404165533cSJose E. Roman   Output Parameters:
6411917a363SLisandro Dalcin + hmin - minimum time step
6421917a363SLisandro Dalcin - hmax - maximum time step
6431917a363SLisandro Dalcin 
6441917a363SLisandro Dalcin   Level: intermediate
6451917a363SLisandro Dalcin 
6461cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptSetStepLimits()`, `TSAdaptChoose()`
6471917a363SLisandro Dalcin @*/
648d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt, PetscReal *hmin, PetscReal *hmax)
649d71ae5a4SJacob Faibussowitsch {
6501917a363SLisandro Dalcin   PetscFunctionBegin;
6511917a363SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
6524f572ea9SToby Isaac   if (hmin) PetscAssertPointer(hmin, 2);
6534f572ea9SToby Isaac   if (hmax) PetscAssertPointer(hmax, 3);
6541917a363SLisandro Dalcin   if (hmin) *hmin = adapt->dt_min;
6551917a363SLisandro Dalcin   if (hmax) *hmax = adapt->dt_max;
6563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6571c3436cfSJed Brown }
6581c3436cfSJed Brown 
659*14d0ab18SJacob Faibussowitsch /*@C
660bcf0153eSBarry Smith   TSAdaptSetFromOptions - Sets various `TSAdapt` parameters from user options.
66184df9cb4SJed Brown 
662c3339decSBarry Smith   Collective
66384df9cb4SJed Brown 
664*14d0ab18SJacob Faibussowitsch   Input Parameters:
6652fe279fdSBarry Smith + adapt              - the `TSAdapt` context
6662fe279fdSBarry Smith - PetscOptionsObject - object created by `PetscOptionsBegin()`
66784df9cb4SJed Brown 
66884df9cb4SJed Brown   Options Database Keys:
6691917a363SLisandro Dalcin + -ts_adapt_type <type>                - algorithm to use for adaptivity
670bf997491SLisandro Dalcin . -ts_adapt_always_accept              - always accept steps regardless of error/stability goals
6711917a363SLisandro Dalcin . -ts_adapt_safety <safety>            - safety factor relative to target error/stability goal
6721917a363SLisandro Dalcin . -ts_adapt_reject_safety <safety>     - extra safety factor to apply if the last step was rejected
6731917a363SLisandro Dalcin . -ts_adapt_clip <low,high>            - admissible time step decrease and increase factors
67423de1d84SBarry Smith . -ts_adapt_dt_min <min>               - minimum timestep to use
67523de1d84SBarry Smith . -ts_adapt_dt_max <max>               - maximum timestep to use
67623de1d84SBarry Smith . -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails
677de50f1caSBarry Smith . -ts_adapt_wnormtype <2 or infinity>  - type of norm for computing error estimates
678de50f1caSBarry 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
67984df9cb4SJed Brown 
68084df9cb4SJed Brown   Level: advanced
68184df9cb4SJed Brown 
682bcf0153eSBarry Smith   Note:
683bcf0153eSBarry Smith   This function is automatically called by `TSSetFromOptions()`
68484df9cb4SJed Brown 
6851cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptSetAlwaysAccept()`, `TSAdaptSetSafety()`,
686db781477SPatrick Sanan           `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetStepLimits()`, `TSAdaptSetMonitor()`
687*14d0ab18SJacob Faibussowitsch @*/
688d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetFromOptions(TSAdapt adapt, PetscOptionItems *PetscOptionsObject)
689d71ae5a4SJacob Faibussowitsch {
69084df9cb4SJed Brown   char      type[256] = TSADAPTBASIC;
69162c23b28SMark   PetscReal safety, reject_safety, clip[2], scale, hmin, hmax;
6921c3436cfSJed Brown   PetscBool set, flg;
6931917a363SLisandro Dalcin   PetscInt  two;
69484df9cb4SJed Brown 
69584df9cb4SJed Brown   PetscFunctionBegin;
696dbbe0bcdSBarry Smith   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
69784df9cb4SJed Brown   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
6981566a47fSLisandro Dalcin    * function can only be called from inside TSSetFromOptions()  */
699d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "TS Adaptivity options");
7009566063dSJacob 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));
70148a46eb9SPierre Jolivet   if (flg || !((PetscObject)adapt)->type_name) PetscCall(TSAdaptSetType(adapt, type));
7021917a363SLisandro Dalcin 
7039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_adapt_always_accept", "Always accept the step", "TSAdaptSetAlwaysAccept", adapt->always_accept, &flg, &set));
7049566063dSJacob Faibussowitsch   if (set) PetscCall(TSAdaptSetAlwaysAccept(adapt, flg));
7051917a363SLisandro Dalcin 
7069371c9d4SSatish Balay   safety        = adapt->safety;
7079371c9d4SSatish Balay   reject_safety = adapt->reject_safety;
7089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_adapt_safety", "Safety factor relative to target error/stability goal", "TSAdaptSetSafety", safety, &safety, &set));
7099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_adapt_reject_safety", "Extra safety factor to apply if the last step was rejected", "TSAdaptSetSafety", reject_safety, &reject_safety, &flg));
7109566063dSJacob Faibussowitsch   if (set || flg) PetscCall(TSAdaptSetSafety(adapt, safety, reject_safety));
7111917a363SLisandro Dalcin 
7129371c9d4SSatish Balay   two     = 2;
7139371c9d4SSatish Balay   clip[0] = adapt->clip[0];
7149371c9d4SSatish Balay   clip[1] = adapt->clip[1];
7159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRealArray("-ts_adapt_clip", "Admissible decrease/increase factor in step size", "TSAdaptSetClip", clip, &two, &set));
7163c633725SBarry Smith   PetscCheck(!set || (two == 2), PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Must give exactly two values to -ts_adapt_clip");
7179566063dSJacob Faibussowitsch   if (set) PetscCall(TSAdaptSetClip(adapt, clip[0], clip[1]));
7181917a363SLisandro Dalcin 
7199371c9d4SSatish Balay   hmin = adapt->dt_min;
7209371c9d4SSatish Balay   hmax = adapt->dt_max;
7219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_adapt_dt_min", "Minimum time step considered", "TSAdaptSetStepLimits", hmin, &hmin, &set));
7229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_adapt_dt_max", "Maximum time step considered", "TSAdaptSetStepLimits", hmax, &hmax, &flg));
7239566063dSJacob Faibussowitsch   if (set || flg) PetscCall(TSAdaptSetStepLimits(adapt, hmin, hmax));
7241917a363SLisandro Dalcin 
7259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_adapt_max_ignore", "Adaptor ignores (absolute) solution values smaller than this value", "", adapt->ignore_max, &adapt->ignore_max, &set));
7269566063dSJacob 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));
727d580f011SEmil Constantinescu 
7289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_adapt_scale_solve_failed", "Scale step by this factor if solve fails", "TSAdaptSetScaleSolveFailed", adapt->scale_solve_failed, &scale, &set));
7299566063dSJacob Faibussowitsch   if (set) PetscCall(TSAdaptSetScaleSolveFailed(adapt, scale));
7301917a363SLisandro Dalcin 
7319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-ts_adapt_wnormtype", "Type of norm computed for error estimation", "", NormTypes, (PetscEnum)adapt->wnormtype, (PetscEnum *)&adapt->wnormtype, NULL));
7323c633725SBarry Smith   PetscCheck(adapt->wnormtype == NORM_2 || adapt->wnormtype == NORM_INFINITY, PetscObjectComm((PetscObject)adapt), PETSC_ERR_SUP, "Only 2-norm and infinite norm supported");
7331917a363SLisandro Dalcin 
7349566063dSJacob 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));
735de50f1caSBarry Smith 
7369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_adapt_monitor", "Print choices made by adaptive controller", "TSAdaptSetMonitor", adapt->monitor ? PETSC_TRUE : PETSC_FALSE, &flg, &set));
7379566063dSJacob Faibussowitsch   if (set) PetscCall(TSAdaptSetMonitor(adapt, flg));
7381917a363SLisandro Dalcin 
739dbbe0bcdSBarry Smith   PetscTryTypeMethod(adapt, setfromoptions, PetscOptionsObject);
740d0609cedSBarry Smith   PetscOptionsHeadEnd();
7413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
74284df9cb4SJed Brown }
74384df9cb4SJed Brown 
74484df9cb4SJed Brown /*@
74584df9cb4SJed Brown   TSAdaptCandidatesClear - clear any previously set candidate schemes
74684df9cb4SJed Brown 
74720f4b53cSBarry Smith   Logically Collective
74884df9cb4SJed Brown 
7494165533cSJose E. Roman   Input Parameter:
75084df9cb4SJed Brown . adapt - adaptive controller
75184df9cb4SJed Brown 
75284df9cb4SJed Brown   Level: developer
75384df9cb4SJed Brown 
7541cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptCreate()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()`
75584df9cb4SJed Brown @*/
756d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
757d71ae5a4SJacob Faibussowitsch {
75884df9cb4SJed Brown   PetscFunctionBegin;
7594782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
7609566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&adapt->candidates, sizeof(adapt->candidates)));
7613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
76284df9cb4SJed Brown }
76384df9cb4SJed Brown 
76484df9cb4SJed Brown /*@C
76584df9cb4SJed Brown   TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from
76684df9cb4SJed Brown 
76720f4b53cSBarry Smith   Logically Collective; No Fortran Support
76884df9cb4SJed Brown 
7694165533cSJose E. Roman   Input Parameters:
770bcf0153eSBarry Smith + adapt      - time step adaptivity context, obtained with `TSGetAdapt()` or `TSAdaptCreate()`
77184df9cb4SJed Brown . name       - name of the candidate scheme to add
77284df9cb4SJed Brown . order      - order of the candidate scheme
77384df9cb4SJed Brown . stageorder - stage order of the candidate scheme
7748d59e960SJed Brown . ccfl       - stability coefficient relative to explicit Euler, used for CFL constraints
77584df9cb4SJed Brown . cost       - relative measure of the amount of work required for the candidate scheme
77684df9cb4SJed Brown - inuse      - indicates that this scheme is the one currently in use, this flag can only be set for one scheme
77784df9cb4SJed Brown 
77884df9cb4SJed Brown   Level: developer
77984df9cb4SJed Brown 
7801cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptChoose()`
78184df9cb4SJed Brown @*/
782d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt, const char name[], PetscInt order, PetscInt stageorder, PetscReal ccfl, PetscReal cost, PetscBool inuse)
783d71ae5a4SJacob Faibussowitsch {
78484df9cb4SJed Brown   PetscInt c;
78584df9cb4SJed Brown 
78684df9cb4SJed Brown   PetscFunctionBegin;
78784df9cb4SJed Brown   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
78863a3b9bcSJacob Faibussowitsch   PetscCheck(order >= 1, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Classical order %" PetscInt_FMT " must be a positive integer", order);
78984df9cb4SJed Brown   if (inuse) {
7903c633725SBarry Smith     PetscCheck(!adapt->candidates.inuse_set, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONGSTATE, "Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()");
79184df9cb4SJed Brown     adapt->candidates.inuse_set = PETSC_TRUE;
79284df9cb4SJed Brown   }
7931c3436cfSJed Brown   /* first slot if this is the current scheme, otherwise the next available slot */
7941c3436cfSJed Brown   c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n;
795bbd56ea5SKarl Rupp 
79684df9cb4SJed Brown   adapt->candidates.name[c]       = name;
79784df9cb4SJed Brown   adapt->candidates.order[c]      = order;
79884df9cb4SJed Brown   adapt->candidates.stageorder[c] = stageorder;
7998d59e960SJed Brown   adapt->candidates.ccfl[c]       = ccfl;
80084df9cb4SJed Brown   adapt->candidates.cost[c]       = cost;
80184df9cb4SJed Brown   adapt->candidates.n++;
8023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
80384df9cb4SJed Brown }
80484df9cb4SJed Brown 
8058d59e960SJed Brown /*@C
8068d59e960SJed Brown   TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost
8078d59e960SJed Brown 
8088d59e960SJed Brown   Not Collective
8098d59e960SJed Brown 
8104165533cSJose E. Roman   Input Parameter:
8118d59e960SJed Brown . adapt - time step adaptivity context
8128d59e960SJed Brown 
8134165533cSJose E. Roman   Output Parameters:
8148d59e960SJed Brown + n          - number of candidate schemes, always at least 1
8158d59e960SJed Brown . order      - the order of each candidate scheme
8168d59e960SJed Brown . stageorder - the stage order of each candidate scheme
8178d59e960SJed Brown . ccfl       - the CFL coefficient of each scheme
8188d59e960SJed Brown - cost       - the relative cost of each scheme
8198d59e960SJed Brown 
8208d59e960SJed Brown   Level: developer
8218d59e960SJed Brown 
8228d59e960SJed Brown   Note:
8238d59e960SJed Brown   The current scheme is always returned in the first slot
8248d59e960SJed Brown 
8251cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()`
8268d59e960SJed Brown @*/
827d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt, PetscInt *n, const PetscInt **order, const PetscInt **stageorder, const PetscReal **ccfl, const PetscReal **cost)
828d71ae5a4SJacob Faibussowitsch {
8298d59e960SJed Brown   PetscFunctionBegin;
8308d59e960SJed Brown   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
8318d59e960SJed Brown   if (n) *n = adapt->candidates.n;
8328d59e960SJed Brown   if (order) *order = adapt->candidates.order;
8338d59e960SJed Brown   if (stageorder) *stageorder = adapt->candidates.stageorder;
8348d59e960SJed Brown   if (ccfl) *ccfl = adapt->candidates.ccfl;
8358d59e960SJed Brown   if (cost) *cost = adapt->candidates.cost;
8363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8378d59e960SJed Brown }
8388d59e960SJed Brown 
83984df9cb4SJed Brown /*@C
84084df9cb4SJed Brown   TSAdaptChoose - choose which method and step size to use for the next step
84184df9cb4SJed Brown 
842c3339decSBarry Smith   Collective
84384df9cb4SJed Brown 
8444165533cSJose E. Roman   Input Parameters:
845da81f932SPierre Jolivet + adapt - adaptive controller
84697bb3fdcSJose E. Roman . ts    - time stepper
84784df9cb4SJed Brown - h     - current step size
84884df9cb4SJed Brown 
8494165533cSJose E. Roman   Output Parameters:
8501566a47fSLisandro Dalcin + next_sc - optional, scheme to use for the next step
85184df9cb4SJed Brown . next_h  - step size to use for the next step
852bcf0153eSBarry Smith - accept  - `PETSC_TRUE` to accept the current step, `PETSC_FALSE` to repeat the current step with the new step size
8531c3436cfSJed Brown 
85484df9cb4SJed Brown   Level: developer
85584df9cb4SJed Brown 
856bcf0153eSBarry Smith   Note:
857bcf0153eSBarry Smith   The input value of parameter accept is retained from the last time step, so it will be `PETSC_FALSE` if the step is
858bcf0153eSBarry Smith   being retried after an initial rejection.
859bcf0153eSBarry Smith 
8601cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()`
86184df9cb4SJed Brown @*/
862d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptChoose(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept)
863d71ae5a4SJacob Faibussowitsch {
8641566a47fSLisandro Dalcin   PetscInt  ncandidates = adapt->candidates.n;
8651566a47fSLisandro Dalcin   PetscInt  scheme      = 0;
8660b99f514SJed Brown   PetscReal wlte        = -1.0;
8675e4ed32fSEmil Constantinescu   PetscReal wltea       = -1.0;
8685e4ed32fSEmil Constantinescu   PetscReal wlter       = -1.0;
86984df9cb4SJed Brown 
87084df9cb4SJed Brown   PetscFunctionBegin;
87184df9cb4SJed Brown   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
87284df9cb4SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 2);
8734f572ea9SToby Isaac   if (next_sc) PetscAssertPointer(next_sc, 4);
8744f572ea9SToby Isaac   PetscAssertPointer(next_h, 5);
8754f572ea9SToby Isaac   PetscAssertPointer(accept, 6);
8761566a47fSLisandro Dalcin   if (next_sc) *next_sc = 0;
8771566a47fSLisandro Dalcin 
8781566a47fSLisandro Dalcin   /* Do not mess with adaptivity while handling events*/
8791566a47fSLisandro Dalcin   if (ts->event && ts->event->status != TSEVENT_NONE) {
8801566a47fSLisandro Dalcin     *next_h = h;
8811566a47fSLisandro Dalcin     *accept = PETSC_TRUE;
8823ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
8831566a47fSLisandro Dalcin   }
8841566a47fSLisandro Dalcin 
885dbbe0bcdSBarry Smith   PetscUseTypeMethod(adapt, choose, ts, h, &scheme, next_h, accept, &wlte, &wltea, &wlter);
88663a3b9bcSJacob 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);
8873c633725SBarry Smith   PetscCheck(*next_h >= 0, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Computed step size %g must be positive", (double)*next_h);
8881566a47fSLisandro Dalcin   if (next_sc) *next_sc = scheme;
8891566a47fSLisandro Dalcin 
89049354f04SShri Abhyankar   if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
89136b54a69SLisandro Dalcin     /* Increase/reduce step size if end time of next step is close to or overshoots max time */
89236b54a69SLisandro Dalcin     PetscReal t = ts->ptime + ts->time_step, h = *next_h;
8934a658b32SHong Zhang     PetscReal tend = t + h, tmax, hmax;
894fe7350e0SStefano Zampini     PetscReal a    = (PetscReal)(1.0 + adapt->matchstepfac[0]);
895fe7350e0SStefano Zampini     PetscReal b    = adapt->matchstepfac[1];
8964a658b32SHong Zhang 
8974a658b32SHong Zhang     if (ts->tspan) {
898e1db57b0SHong Zhang       if (PetscIsCloseAtTol(t, ts->tspan->span_times[ts->tspan->spanctr], ts->tspan->reltol * h + ts->tspan->abstol, 0)) /* hit a span time point */
8994a658b32SHong Zhang         if (ts->tspan->spanctr + 1 < ts->tspan->num_span_times) tmax = ts->tspan->span_times[ts->tspan->spanctr + 1];
9004a658b32SHong Zhang         else tmax = ts->max_time; /* hit the last span time point */
9014a658b32SHong Zhang       else tmax = ts->tspan->span_times[ts->tspan->spanctr];
9024a658b32SHong Zhang     } else tmax = ts->max_time;
9034a658b32SHong Zhang     hmax = tmax - t;
90436b54a69SLisandro Dalcin     if (t < tmax && tend > tmax) *next_h = hmax;
905fe7350e0SStefano Zampini     if (t < tmax && tend < tmax && h * b > hmax) *next_h = hmax / 2;
90636b54a69SLisandro Dalcin     if (t < tmax && tend < tmax && h * a > hmax) *next_h = hmax;
907e1db57b0SHong Zhang     /* if step size is changed to match a span time point */
908e1db57b0SHong Zhang     if (ts->tspan && h != *next_h && !adapt->dt_span_cached) adapt->dt_span_cached = h;
909e1db57b0SHong Zhang     /* reset time step after a span time point */
910e1db57b0SHong 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)) {
911e1db57b0SHong Zhang       *next_h               = adapt->dt_span_cached;
912e1db57b0SHong Zhang       adapt->dt_span_cached = 0;
91349354f04SShri Abhyankar     }
914e1db57b0SHong Zhang   }
9151c3436cfSJed Brown   if (adapt->monitor) {
9161566a47fSLisandro Dalcin     const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : "";
9179566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
9180b99f514SJed Brown     if (wlte < 0) {
9199371c9d4SSatish 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",
9209371c9d4SSatish Balay                                        (double)ts->ptime, (double)h, (double)*next_h));
9210b99f514SJed Brown     } else {
9229371c9d4SSatish 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",
9239371c9d4SSatish Balay                                        (double)ts->ptime, (double)h, (double)*next_h, (double)wlte, (double)wltea, (double)wlter));
9240b99f514SJed Brown     }
9259566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
9261c3436cfSJed Brown   }
9273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
92884df9cb4SJed Brown }
92984df9cb4SJed Brown 
93097335746SJed Brown /*@
931de50f1caSBarry Smith   TSAdaptSetTimeStepIncreaseDelay - The number of timesteps to wait after a decrease in the timestep due to failed solver
932de50f1caSBarry Smith   before increasing the time step.
933de50f1caSBarry Smith 
934c3339decSBarry Smith   Logicially Collective
935de50f1caSBarry Smith 
9364165533cSJose E. Roman   Input Parameters:
937de50f1caSBarry Smith + adapt - adaptive controller context
938de50f1caSBarry Smith - cnt   - the number of timesteps
939de50f1caSBarry Smith 
940de50f1caSBarry Smith   Options Database Key:
941de50f1caSBarry Smith . -ts_adapt_time_step_increase_delay cnt - number of steps to delay the increase
942de50f1caSBarry Smith 
943de50f1caSBarry Smith   Level: advanced
944de50f1caSBarry Smith 
945bcf0153eSBarry Smith   Notes:
946bcf0153eSBarry Smith   This is to prevent an adaptor from bouncing back and forth between two nearby timesteps. The default is 0.
947bcf0153eSBarry Smith 
948bcf0153eSBarry Smith   The successful use of this option is problem dependent
949bcf0153eSBarry Smith 
950b43aa488SJacob Faibussowitsch   Developer Notes:
951bcf0153eSBarry Smith   There is no theory to support this option
952bcf0153eSBarry Smith 
9531cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`
954de50f1caSBarry Smith @*/
955d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt, PetscInt cnt)
956d71ae5a4SJacob Faibussowitsch {
957de50f1caSBarry Smith   PetscFunctionBegin;
958de50f1caSBarry Smith   adapt->timestepjustdecreased_delay = cnt;
9593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
960de50f1caSBarry Smith }
961de50f1caSBarry Smith 
962de50f1caSBarry Smith /*@
9636bc98fa9SBarry 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)
96497335746SJed Brown 
965c3339decSBarry Smith   Collective
96697335746SJed Brown 
9674165533cSJose E. Roman   Input Parameters:
96897335746SJed Brown + adapt - adaptive controller context
969b295832fSPierre Barbier de Reuille . ts    - time stepper
970b295832fSPierre Barbier de Reuille . t     - Current simulation time
971b295832fSPierre Barbier de Reuille - Y     - Current solution vector
97297335746SJed Brown 
9734165533cSJose E. Roman   Output Parameter:
974bcf0153eSBarry Smith . accept - `PETSC_TRUE` to accept the stage, `PETSC_FALSE` to reject
97597335746SJed Brown 
97697335746SJed Brown   Level: developer
97797335746SJed Brown 
9781cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`
97997335746SJed Brown @*/
980d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCheckStage(TSAdapt adapt, TS ts, PetscReal t, Vec Y, PetscBool *accept)
981d71ae5a4SJacob Faibussowitsch {
9821566a47fSLisandro Dalcin   SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING;
9836c6709e3SStefano Zampini   PetscBool           func_accept, snes_div_func;
98497335746SJed Brown 
98597335746SJed Brown   PetscFunctionBegin;
9864782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
9874782b174SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 2);
9884f572ea9SToby Isaac   PetscAssertPointer(accept, 5);
9891566a47fSLisandro Dalcin 
9906c6709e3SStefano Zampini   PetscCall(TSFunctionDomainError(ts, t, Y, &func_accept));
9919566063dSJacob Faibussowitsch   if (ts->snes) PetscCall(SNESGetConvergedReason(ts->snes, &snesreason));
9926c6709e3SStefano Zampini   snes_div_func = (PetscBool)(snesreason == SNES_DIVERGED_FUNCTION_DOMAIN);
9936c6709e3SStefano Zampini   if (func_accept && snesreason < 0 && !snes_div_func) {
99497335746SJed Brown     *accept = PETSC_FALSE;
9956c6709e3SStefano Zampini     PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", nonlinear solve failure: %s\n", ts->steps, SNESConvergedReasons[snesreason]));
9966de24e2aSJed Brown     if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) {
99797335746SJed Brown       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
99863a3b9bcSJacob 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));
99997335746SJed Brown       if (adapt->monitor) {
10009566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
10019371c9d4SSatish 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,
10029371c9d4SSatish Balay                                          (double)ts->ptime, (double)ts->time_step, ts->num_snes_failures));
10039566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
100497335746SJed Brown       }
1005cb9d8021SPierre Barbier de Reuille     }
1006cb9d8021SPierre Barbier de Reuille   } else {
10076c6709e3SStefano Zampini     *accept = (PetscBool)(func_accept && !snes_div_func);
10086c6709e3SStefano Zampini     if (*accept && adapt->checkstage) PetscCall((*adapt->checkstage)(adapt, ts, t, Y, accept));
10096bc98fa9SBarry Smith     if (!*accept) {
10106c6709e3SStefano Zampini       const char *user_func = !func_accept ? "TSSetFunctionDomainError()" : "TSAdaptSetCheckStage";
10116c6709e3SStefano Zampini       const char *snes_err  = "SNES invalid function domain";
10126c6709e3SStefano Zampini       const char *err_msg   = snes_div_func && func_accept ? snes_err : user_func;
10136c6709e3SStefano Zampini       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", solution rejected by %s\n", ts->steps, err_msg));
10146bc98fa9SBarry Smith       if (adapt->monitor) {
10159566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
10166c6709e3SStefano Zampini         PetscCall(PetscViewerASCIIPrintf(adapt->monitor, "    TSAdapt %s step %3" PetscInt_FMT " stage rejected by %s\n", ((PetscObject)adapt)->type_name, ts->steps, err_msg));
10179566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
10186bc98fa9SBarry Smith       }
10196bc98fa9SBarry Smith     }
1020cb9d8021SPierre Barbier de Reuille   }
1021cb9d8021SPierre Barbier de Reuille 
10221566a47fSLisandro Dalcin   if (!(*accept) && !ts->reason) {
10231566a47fSLisandro Dalcin     PetscReal dt, new_dt;
10249566063dSJacob Faibussowitsch     PetscCall(TSGetTimeStep(ts, &dt));
1025cb9d8021SPierre Barbier de Reuille     new_dt = dt * adapt->scale_solve_failed;
10269566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts, new_dt));
1027de50f1caSBarry Smith     adapt->timestepjustdecreased += adapt->timestepjustdecreased_delay;
102897335746SJed Brown     if (adapt->monitor) {
10299566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
10306c6709e3SStefano Zampini       PetscCall(PetscViewerASCIIPrintf(adapt->monitor, "    TSAdapt %s step %3" PetscInt_FMT " stage rejected (SNES reason %s) t=%-11g+%10.3e retrying with dt=%-10.3e\n", ((PetscObject)adapt)->type_name, ts->steps, SNESConvergedReasons[snesreason],
10316c6709e3SStefano Zampini                                        (double)ts->ptime, (double)dt, (double)new_dt));
10329566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
103397335746SJed Brown     }
103497335746SJed Brown   }
10353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
103697335746SJed Brown }
103797335746SJed Brown 
103884df9cb4SJed Brown /*@
103984df9cb4SJed Brown   TSAdaptCreate - create an adaptive controller context for time stepping
104084df9cb4SJed Brown 
1041d083f849SBarry Smith   Collective
104284df9cb4SJed Brown 
104384df9cb4SJed Brown   Input Parameter:
104484df9cb4SJed Brown . comm - The communicator
104584df9cb4SJed Brown 
104684df9cb4SJed Brown   Output Parameter:
1047b43aa488SJacob Faibussowitsch . inadapt - new `TSAdapt` object
104884df9cb4SJed Brown 
104984df9cb4SJed Brown   Level: developer
105084df9cb4SJed Brown 
1051bcf0153eSBarry Smith   Note:
1052bcf0153eSBarry Smith   `TSAdapt` creation is handled by `TS`, so users should not need to call this function.
105384df9cb4SJed Brown 
10541cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdapt`, `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptDestroy()`
105584df9cb4SJed Brown @*/
1056d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCreate(MPI_Comm comm, TSAdapt *inadapt)
1057d71ae5a4SJacob Faibussowitsch {
105884df9cb4SJed Brown   TSAdapt adapt;
105984df9cb4SJed Brown 
106084df9cb4SJed Brown   PetscFunctionBegin;
10614f572ea9SToby Isaac   PetscAssertPointer(inadapt, 2);
10623b3bcf4cSLisandro Dalcin   *inadapt = NULL;
10639566063dSJacob Faibussowitsch   PetscCall(TSAdaptInitializePackage());
10643b3bcf4cSLisandro Dalcin 
10659566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(adapt, TSADAPT_CLASSID, "TSAdapt", "Time stepping adaptivity", "TS", comm, TSAdaptDestroy, TSAdaptView));
10661c3436cfSJed Brown 
1067bf997491SLisandro Dalcin   adapt->always_accept      = PETSC_FALSE;
10681917a363SLisandro Dalcin   adapt->safety             = 0.9;
10691917a363SLisandro Dalcin   adapt->reject_safety      = 0.5;
10701917a363SLisandro Dalcin   adapt->clip[0]            = 0.1;
10711917a363SLisandro Dalcin   adapt->clip[1]            = 10.;
10721c3436cfSJed Brown   adapt->dt_min             = 1e-20;
10731917a363SLisandro Dalcin   adapt->dt_max             = 1e+20;
10741c167fc2SEmil Constantinescu   adapt->ignore_max         = -1.0;
1075d580f011SEmil Constantinescu   adapt->glee_use_local     = PETSC_TRUE;
107697335746SJed Brown   adapt->scale_solve_failed = 0.25;
1077fe7350e0SStefano Zampini   /* these two safety factors are not public, and they are used only in the TS_EXACTFINALTIME_MATCHSTEP case
1078fe7350e0SStefano Zampini      to prevent from situations were unreasonably small time steps are taken in order to match the final time */
1079fe7350e0SStefano Zampini   adapt->matchstepfac[0]             = 0.01; /* allow 1% step size increase in the last step */
1080fe7350e0SStefano Zampini   adapt->matchstepfac[1]             = 2.0;  /* halve last step if it is greater than what remains divided this factor */
10817619abb3SShri   adapt->wnormtype                   = NORM_2;
1082de50f1caSBarry Smith   adapt->timestepjustdecreased_delay = 0;
10831917a363SLisandro Dalcin 
108484df9cb4SJed Brown   *inadapt = adapt;
10853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
108684df9cb4SJed Brown }
1087