1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
284df9cb4SJed Brown
31b9b13dfSLisandro Dalcin PetscClassId TSADAPT_CLASSID;
41b9b13dfSLisandro Dalcin
5140e18c1SBarry Smith static PetscFunctionList TSAdaptList;
684df9cb4SJed Brown static PetscBool TSAdaptPackageInitialized;
784df9cb4SJed Brown static PetscBool TSAdaptRegisterAllCalled;
884df9cb4SJed Brown
98cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_None(TSAdapt);
101566a47fSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt);
11cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptCreate_DSP(TSAdapt);
128cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt);
131917a363SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt);
14d949e4a4SStefano Zampini PETSC_EXTERN PetscErrorCode TSAdaptCreate_History(TSAdapt);
1584df9cb4SJed Brown
1684df9cb4SJed Brown /*@C
171c84c290SBarry Smith TSAdaptRegister - adds a TSAdapt implementation
181c84c290SBarry Smith
19cc4c1da9SBarry Smith Not Collective, No Fortran Support
201c84c290SBarry Smith
211c84c290SBarry Smith Input Parameters:
222fe279fdSBarry Smith + sname - name of user-defined adaptivity scheme
232fe279fdSBarry Smith - function - routine to create method context
241c84c290SBarry Smith
25bcf0153eSBarry Smith Level: advanced
26bcf0153eSBarry Smith
271c84c290SBarry Smith Notes:
28bcf0153eSBarry Smith `TSAdaptRegister()` may be called multiple times to add several user-defined families.
291c84c290SBarry Smith
30b43aa488SJacob Faibussowitsch Example Usage:
311c84c290SBarry Smith .vb
32bdf89e91SBarry Smith TSAdaptRegister("my_scheme", MySchemeCreate);
331c84c290SBarry Smith .ve
341c84c290SBarry Smith
351c84c290SBarry Smith Then, your scheme can be chosen with the procedural interface via
36b44f4de4SBarry Smith .vb
37b44f4de4SBarry Smith TSAdaptSetType(ts, "my_scheme")
38b44f4de4SBarry Smith .ve
391c84c290SBarry Smith or at runtime via the option
40b44f4de4SBarry Smith .vb
41b44f4de4SBarry Smith -ts_adapt_type my_scheme
42b44f4de4SBarry Smith .ve
4384df9cb4SJed Brown
44*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdaptRegisterAll()`
4584df9cb4SJed Brown @*/
TSAdaptRegister(const char sname[],PetscErrorCode (* function)(TSAdapt))46d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptRegister(const char sname[], PetscErrorCode (*function)(TSAdapt))
47d71ae5a4SJacob Faibussowitsch {
4884df9cb4SJed Brown PetscFunctionBegin;
499566063dSJacob Faibussowitsch PetscCall(TSAdaptInitializePackage());
509566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&TSAdaptList, sname, function));
513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5284df9cb4SJed Brown }
5384df9cb4SJed Brown
5484df9cb4SJed Brown /*@C
552fe279fdSBarry Smith TSAdaptRegisterAll - Registers all of the adaptivity schemes in `TSAdapt`
5684df9cb4SJed Brown
5784df9cb4SJed Brown Not Collective
5884df9cb4SJed Brown
5984df9cb4SJed Brown Level: advanced
6084df9cb4SJed Brown
611cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdaptRegisterDestroy()`
6284df9cb4SJed Brown @*/
TSAdaptRegisterAll(void)63d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptRegisterAll(void)
64d71ae5a4SJacob Faibussowitsch {
6584df9cb4SJed Brown PetscFunctionBegin;
663ba16761SJacob Faibussowitsch if (TSAdaptRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
670f51fdf8SToby Isaac TSAdaptRegisterAllCalled = PETSC_TRUE;
689566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None));
699566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTBASIC, TSAdaptCreate_Basic));
709566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTDSP, TSAdaptCreate_DSP));
719566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTCFL, TSAdaptCreate_CFL));
729566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTGLEE, TSAdaptCreate_GLEE));
739566063dSJacob Faibussowitsch PetscCall(TSAdaptRegister(TSADAPTHISTORY, TSAdaptCreate_History));
743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7584df9cb4SJed Brown }
7684df9cb4SJed Brown
7784df9cb4SJed Brown /*@C
782fe279fdSBarry Smith TSAdaptFinalizePackage - This function destroys everything in the `TS` package. It is
792fe279fdSBarry Smith called from `PetscFinalize()`.
8084df9cb4SJed Brown
8184df9cb4SJed Brown Level: developer
8284df9cb4SJed Brown
831cc06b55SBarry Smith .seealso: [](ch_ts), `PetscFinalize()`
8484df9cb4SJed Brown @*/
TSAdaptFinalizePackage(void)85d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptFinalizePackage(void)
86d71ae5a4SJacob Faibussowitsch {
8784df9cb4SJed Brown PetscFunctionBegin;
889566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&TSAdaptList));
8984df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_FALSE;
9084df9cb4SJed Brown TSAdaptRegisterAllCalled = PETSC_FALSE;
913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9284df9cb4SJed Brown }
9384df9cb4SJed Brown
9484df9cb4SJed Brown /*@C
952fe279fdSBarry Smith TSAdaptInitializePackage - This function initializes everything in the `TSAdapt` package. It is
962fe279fdSBarry Smith called from `TSInitializePackage()`.
9784df9cb4SJed Brown
9884df9cb4SJed Brown Level: developer
9984df9cb4SJed Brown
1001cc06b55SBarry Smith .seealso: [](ch_ts), `PetscInitialize()`
10184df9cb4SJed Brown @*/
TSAdaptInitializePackage(void)102d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptInitializePackage(void)
103d71ae5a4SJacob Faibussowitsch {
10484df9cb4SJed Brown PetscFunctionBegin;
1053ba16761SJacob Faibussowitsch if (TSAdaptPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
10684df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_TRUE;
1079566063dSJacob Faibussowitsch PetscCall(PetscClassIdRegister("TSAdapt", &TSADAPT_CLASSID));
1089566063dSJacob Faibussowitsch PetscCall(TSAdaptRegisterAll());
1099566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSAdaptFinalizePackage));
1103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11184df9cb4SJed Brown }
11284df9cb4SJed Brown
113cc4c1da9SBarry Smith /*@
114bcf0153eSBarry Smith TSAdaptSetType - sets the approach used for the error adapter
1157eef6055SBarry Smith
11620f4b53cSBarry Smith Logicially Collective
1177eef6055SBarry Smith
118d8d19677SJose E. Roman Input Parameters:
11958b119d1SBarry Smith + adapt - the `TS` adapter, most likely obtained with `TSGetAdapt()`
120bcf0153eSBarry Smith - type - one of the `TSAdaptType`
1217eef6055SBarry Smith
122bcf0153eSBarry Smith Options Database Key:
123ec18b7bcSLisandro Dalcin . -ts_adapt_type <basic or dsp or none> - to set the adapter type
1247eef6055SBarry Smith
1257eef6055SBarry Smith Level: intermediate
1267eef6055SBarry Smith
127*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSGetAdapt()`, `TSAdaptDestroy()`, `TSAdaptType`, `TSAdaptGetType()`
1287eef6055SBarry Smith @*/
TSAdaptSetType(TSAdapt adapt,TSAdaptType type)129d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetType(TSAdapt adapt, TSAdaptType type)
130d71ae5a4SJacob Faibussowitsch {
131ef749922SLisandro Dalcin PetscBool match;
1325f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(TSAdapt);
13384df9cb4SJed Brown
13484df9cb4SJed Brown PetscFunctionBegin;
1354782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
1364f572ea9SToby Isaac PetscAssertPointer(type, 2);
1379566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)adapt, type, &match));
1383ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS);
1399566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(TSAdaptList, type, &r));
1406adde796SStefano Zampini PetscCheck(r, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSAdapt type \"%s\" given", type);
141dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, destroy);
1429566063dSJacob Faibussowitsch PetscCall(PetscMemzero(adapt->ops, sizeof(struct _TSAdaptOps)));
1439566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)adapt, type));
1449566063dSJacob Faibussowitsch PetscCall((*r)(adapt));
1453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
14684df9cb4SJed Brown }
14784df9cb4SJed Brown
148cc4c1da9SBarry Smith /*@
149a1cb98faSBarry Smith TSAdaptGetType - gets the `TS` adapter method type (as a string).
150d0288e4fSLisandro Dalcin
151d0288e4fSLisandro Dalcin Not Collective
152d0288e4fSLisandro Dalcin
153d0288e4fSLisandro Dalcin Input Parameter:
154a1cb98faSBarry Smith . adapt - The `TS` adapter, most likely obtained with `TSGetAdapt()`
155d0288e4fSLisandro Dalcin
156d0288e4fSLisandro Dalcin Output Parameter:
157a1cb98faSBarry Smith . type - The name of `TS` adapter method
158d0288e4fSLisandro Dalcin
159d0288e4fSLisandro Dalcin Level: intermediate
160d0288e4fSLisandro Dalcin
161*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSAdaptType`, `TSAdaptSetType()`
162d0288e4fSLisandro Dalcin @*/
TSAdaptGetType(TSAdapt adapt,TSAdaptType * type)163d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetType(TSAdapt adapt, TSAdaptType *type)
164d71ae5a4SJacob Faibussowitsch {
165d0288e4fSLisandro Dalcin PetscFunctionBegin;
166d0288e4fSLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
1674f572ea9SToby Isaac PetscAssertPointer(type, 2);
168d0288e4fSLisandro Dalcin *type = ((PetscObject)adapt)->type_name;
1693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
170d0288e4fSLisandro Dalcin }
171d0288e4fSLisandro Dalcin
TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[])172d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt, const char prefix[])
173d71ae5a4SJacob Faibussowitsch {
17484df9cb4SJed Brown PetscFunctionBegin;
1754782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
1769566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adapt, prefix));
1773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
17884df9cb4SJed Brown }
17984df9cb4SJed Brown
180ffeef943SBarry Smith /*@
1812fe279fdSBarry Smith TSAdaptLoad - Loads a TSAdapt that has been stored in binary with `TSAdaptView()`.
182ad6bc421SBarry Smith
183c3339decSBarry Smith Collective
184ad6bc421SBarry Smith
185ad6bc421SBarry Smith Input Parameters:
186b43aa488SJacob Faibussowitsch + adapt - the newly loaded `TSAdapt`, this needs to have been created with `TSAdaptCreate()` or
187bcf0153eSBarry Smith some related function before a call to `TSAdaptLoad()`.
188bcf0153eSBarry Smith - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()` or
189bcf0153eSBarry Smith HDF5 file viewer, obtained from `PetscViewerHDF5Open()`
190ad6bc421SBarry Smith
191ad6bc421SBarry Smith Level: intermediate
192ad6bc421SBarry Smith
193bcf0153eSBarry Smith Note:
194bcf0153eSBarry Smith The type is determined by the data in the file, any type set into the `TSAdapt` before this call is ignored.
195ad6bc421SBarry Smith
1961cc06b55SBarry Smith .seealso: [](ch_ts), `PetscViewerBinaryOpen()`, `TSAdaptView()`, `MatLoad()`, `VecLoad()`, `TSAdapt`
197ad6bc421SBarry Smith @*/
TSAdaptLoad(TSAdapt adapt,PetscViewer viewer)198d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptLoad(TSAdapt adapt, PetscViewer viewer)
199d71ae5a4SJacob Faibussowitsch {
200ad6bc421SBarry Smith PetscBool isbinary;
201ad6bc421SBarry Smith char type[256];
202ad6bc421SBarry Smith
203ad6bc421SBarry Smith PetscFunctionBegin;
2044782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
205ad6bc421SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
2069566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
2073c633725SBarry Smith PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()");
208ad6bc421SBarry Smith
2099566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR));
2109566063dSJacob Faibussowitsch PetscCall(TSAdaptSetType(adapt, type));
211dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, load, viewer);
2123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
213ad6bc421SBarry Smith }
214ad6bc421SBarry Smith
215*efa39862SBarry Smith /*@
216*efa39862SBarry Smith TSAdaptView - Prints the `TSAdapt` data structure.
217*efa39862SBarry Smith
218*efa39862SBarry Smith Collective
219*efa39862SBarry Smith
220*efa39862SBarry Smith Input Parameters:
221*efa39862SBarry Smith + adapt - the `TSAdapt` context obtained from `TSGetAdapt()`
222*efa39862SBarry Smith - viewer - visualization context
223*efa39862SBarry Smith
224*efa39862SBarry Smith Options Database Key:
225*efa39862SBarry Smith . -ts_view - calls `TSView()` at end of `TSStep()`
226*efa39862SBarry Smith
227*efa39862SBarry Smith Level: advanced
228*efa39862SBarry Smith
229*efa39862SBarry Smith Notes:
230*efa39862SBarry Smith This is called by `TSView()` so rarely called directly.
231*efa39862SBarry Smith
232*efa39862SBarry Smith The available visualization contexts include
233*efa39862SBarry Smith + `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
234*efa39862SBarry Smith - `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
235*efa39862SBarry Smith output where only the first processor opens
236*efa39862SBarry Smith the file. All other processes send their
237*efa39862SBarry Smith data to the first process to print.
238*efa39862SBarry Smith
239*efa39862SBarry Smith The user can open an alternative visualization context with
240*efa39862SBarry Smith `PetscViewerASCIIOpen()` - output to a specified file.
241*efa39862SBarry Smith
242*efa39862SBarry Smith In the debugger you can do call `TSAdaptView`(adapt,0) to display the `TSAdapt`. (The same holds for any PETSc object viewer).
243*efa39862SBarry Smith
244*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSView()`, `PetscViewer`, `PetscViewerASCIIOpen()`
245*efa39862SBarry Smith @*/
TSAdaptView(TSAdapt adapt,PetscViewer viewer)246d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptView(TSAdapt adapt, PetscViewer viewer)
247d71ae5a4SJacob Faibussowitsch {
2489f196a02SMartin Diehl PetscBool isascii, isbinary, isnone, isglee;
24984df9cb4SJed Brown
25084df9cb4SJed Brown PetscFunctionBegin;
2514782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
2529566063dSJacob Faibussowitsch if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt), &viewer));
2534782b174SLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
2544782b174SLisandro Dalcin PetscCheckSameComm(adapt, 1, viewer, 2);
2559f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
2569566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
2579f196a02SMartin Diehl if (isascii) {
2589566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)adapt, viewer));
2599566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTNONE, &isnone));
2609566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTGLEE, &isglee));
2611917a363SLisandro Dalcin if (!isnone) {
2629566063dSJacob Faibussowitsch if (adapt->always_accept) PetscCall(PetscViewerASCIIPrintf(viewer, " always accepting steps\n"));
2639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " safety factor %g\n", (double)adapt->safety));
2649566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " extra safety factor after step rejection %g\n", (double)adapt->reject_safety));
2659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " clip fastest increase %g\n", (double)adapt->clip[1]));
2669566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " clip fastest decrease %g\n", (double)adapt->clip[0]));
2679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " maximum allowed timestep %g\n", (double)adapt->dt_max));
2689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " minimum allowed timestep %g\n", (double)adapt->dt_min));
2699566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " maximum solution absolute value to be ignored %g\n", (double)adapt->ignore_max));
2701c167fc2SEmil Constantinescu }
2711c167fc2SEmil Constantinescu if (isglee) {
2721c167fc2SEmil Constantinescu if (adapt->glee_use_local) {
2739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " GLEE uses local error control\n"));
2741c167fc2SEmil Constantinescu } else {
2759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " GLEE uses global error control\n"));
2761c167fc2SEmil Constantinescu }
2771917a363SLisandro Dalcin }
2789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer));
279dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, view, viewer);
2809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer));
281ad6bc421SBarry Smith } else if (isbinary) {
282ad6bc421SBarry Smith char type[256];
283ad6bc421SBarry Smith
284ad6bc421SBarry Smith /* need to save FILE_CLASS_ID for adapt class */
2859566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(type, ((PetscObject)adapt)->type_name, 256));
2869566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR));
287dbbe0bcdSBarry Smith } else PetscTryTypeMethod(adapt, view, viewer);
2883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
28984df9cb4SJed Brown }
29084df9cb4SJed Brown
291099af0a3SLisandro Dalcin /*@
292bcf0153eSBarry Smith TSAdaptReset - Resets a `TSAdapt` context to its defaults
293099af0a3SLisandro Dalcin
294c3339decSBarry Smith Collective
295099af0a3SLisandro Dalcin
296099af0a3SLisandro Dalcin Input Parameter:
297bcf0153eSBarry Smith . adapt - the `TSAdapt` context obtained from `TSGetAdapt()` or `TSAdaptCreate()`
298099af0a3SLisandro Dalcin
299099af0a3SLisandro Dalcin Level: developer
300099af0a3SLisandro Dalcin
301*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSGetAdapt()`, `TSAdapt`, `TSAdaptCreate()`, `TSAdaptDestroy()`
302099af0a3SLisandro Dalcin @*/
TSAdaptReset(TSAdapt adapt)303d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptReset(TSAdapt adapt)
304d71ae5a4SJacob Faibussowitsch {
305099af0a3SLisandro Dalcin PetscFunctionBegin;
306099af0a3SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
307dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, reset);
3083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
309099af0a3SLisandro Dalcin }
310099af0a3SLisandro Dalcin
TSAdaptDestroy(TSAdapt * adapt)311d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptDestroy(TSAdapt *adapt)
312d71ae5a4SJacob Faibussowitsch {
31384df9cb4SJed Brown PetscFunctionBegin;
3143ba16761SJacob Faibussowitsch if (!*adapt) PetscFunctionReturn(PETSC_SUCCESS);
31584df9cb4SJed Brown PetscValidHeaderSpecific(*adapt, TSADAPT_CLASSID, 1);
316f4f49eeaSPierre Jolivet if (--((PetscObject)*adapt)->refct > 0) {
3179371c9d4SSatish Balay *adapt = NULL;
3183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3199371c9d4SSatish Balay }
320099af0a3SLisandro Dalcin
3219566063dSJacob Faibussowitsch PetscCall(TSAdaptReset(*adapt));
322099af0a3SLisandro Dalcin
323f4f49eeaSPierre Jolivet PetscTryTypeMethod(*adapt, destroy);
3249566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&(*adapt)->monitor));
3259566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(adapt));
3263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
32784df9cb4SJed Brown }
32884df9cb4SJed Brown
3291c3436cfSJed Brown /*@
3301c3436cfSJed Brown TSAdaptSetMonitor - Monitor the choices made by the adaptive controller
3311c3436cfSJed Brown
332c3339decSBarry Smith Collective
3331c3436cfSJed Brown
3344165533cSJose E. Roman Input Parameters:
3351c3436cfSJed Brown + adapt - adaptive controller context
336bcf0153eSBarry Smith - flg - `PETSC_TRUE` to active a monitor, `PETSC_FALSE` to disable
3371c3436cfSJed Brown
338bcf0153eSBarry Smith Options Database Key:
339ec18b7bcSLisandro Dalcin . -ts_adapt_monitor - to turn on monitoring
340bf997491SLisandro Dalcin
3411c3436cfSJed Brown Level: intermediate
3421c3436cfSJed Brown
343*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSGetAdapt()`, `TSAdaptChoose()`
3441c3436cfSJed Brown @*/
TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg)345d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt, PetscBool flg)
346d71ae5a4SJacob Faibussowitsch {
3471c3436cfSJed Brown PetscFunctionBegin;
3484782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
3494782b174SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt, flg, 2);
3501c3436cfSJed Brown if (flg) {
3519566063dSJacob Faibussowitsch if (!adapt->monitor) PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt), "stdout", &adapt->monitor));
3521c3436cfSJed Brown } else {
3539566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&adapt->monitor));
3541c3436cfSJed Brown }
3553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3561c3436cfSJed Brown }
3571c3436cfSJed Brown
3580873808bSJed Brown /*@C
359bf997491SLisandro Dalcin TSAdaptSetCheckStage - Set a callback to check convergence for a stage
3600873808bSJed Brown
3612fe279fdSBarry Smith Logically Collective
3620873808bSJed Brown
3634165533cSJose E. Roman Input Parameters:
3640873808bSJed Brown + adapt - adaptive controller context
3650873808bSJed Brown - func - stage check function
3660873808bSJed Brown
367b43aa488SJacob Faibussowitsch Calling sequence:
3680873808bSJed Brown + adapt - adaptive controller context
3690873808bSJed Brown . ts - time stepping context
37014d0ab18SJacob Faibussowitsch . t - current time
37114d0ab18SJacob Faibussowitsch . Y - current solution vector
3720873808bSJed Brown - accept - pending choice of whether to accept, can be modified by this routine
3730873808bSJed Brown
3740873808bSJed Brown Level: advanced
3750873808bSJed Brown
376*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSGetAdapt()`, `TSAdaptChoose()`
3770873808bSJed Brown @*/
TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (* func)(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool * accept))37814d0ab18SJacob Faibussowitsch PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt, PetscErrorCode (*func)(TSAdapt adapt, TS ts, PetscReal t, Vec Y, PetscBool *accept))
379d71ae5a4SJacob Faibussowitsch {
3800873808bSJed Brown PetscFunctionBegin;
3810873808bSJed Brown PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
38268ae941aSLisandro Dalcin adapt->checkstage = func;
3833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3840873808bSJed Brown }
3850873808bSJed Brown
3861c3436cfSJed Brown /*@
387bf997491SLisandro Dalcin TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of
388bf997491SLisandro Dalcin any error or stability condition not meeting the prescribed goal.
389bf997491SLisandro Dalcin
3902fe279fdSBarry Smith Logically Collective
391bf997491SLisandro Dalcin
3924165533cSJose E. Roman Input Parameters:
393bcf0153eSBarry Smith + adapt - time step adaptivity context, usually gotten with `TSGetAdapt()`
394bf997491SLisandro Dalcin - flag - whether to always accept steps
395bf997491SLisandro Dalcin
396bcf0153eSBarry Smith Options Database Key:
397ec18b7bcSLisandro Dalcin . -ts_adapt_always_accept - to always accept steps
398bf997491SLisandro Dalcin
399bf997491SLisandro Dalcin Level: intermediate
400bf997491SLisandro Dalcin
401*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSGetAdapt()`, `TSAdaptChoose()`
402bf997491SLisandro Dalcin @*/
TSAdaptSetAlwaysAccept(TSAdapt adapt,PetscBool flag)403d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt, PetscBool flag)
404d71ae5a4SJacob Faibussowitsch {
405bf997491SLisandro Dalcin PetscFunctionBegin;
406bf997491SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
407bf997491SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt, flag, 2);
408bf997491SLisandro Dalcin adapt->always_accept = flag;
4093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
410bf997491SLisandro Dalcin }
411bf997491SLisandro Dalcin
412bf997491SLisandro Dalcin /*@
413bcf0153eSBarry Smith TSAdaptSetSafety - Set safety factors for time step adaptor
4141c3436cfSJed Brown
4152fe279fdSBarry Smith Logically Collective
4161917a363SLisandro Dalcin
4174165533cSJose E. Roman Input Parameters:
4181917a363SLisandro Dalcin + adapt - adaptive controller context
4191917a363SLisandro Dalcin . safety - safety factor relative to target error/stability goal
420ec18b7bcSLisandro Dalcin - reject_safety - extra safety factor to apply if the last step was rejected
4211917a363SLisandro Dalcin
4221917a363SLisandro Dalcin Options Database Keys:
423ec18b7bcSLisandro Dalcin + -ts_adapt_safety <safety> - to set safety factor
424ec18b7bcSLisandro Dalcin - -ts_adapt_reject_safety <reject_safety> - to set reject safety factor
4251917a363SLisandro Dalcin
4261917a363SLisandro Dalcin Level: intermediate
4271917a363SLisandro Dalcin
42809cb0f53SBarry Smith Note:
42909cb0f53SBarry Smith Use `PETSC_CURRENT` to keep the current value for either parameter
43009cb0f53SBarry Smith
43109cb0f53SBarry Smith Fortran Note:
43209cb0f53SBarry Smith Use `PETSC_CURRENT_REAL`
43309cb0f53SBarry Smith
434*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSAdaptGetSafety()`, `TSAdaptChoose()`
4351917a363SLisandro Dalcin @*/
TSAdaptSetSafety(TSAdapt adapt,PetscReal safety,PetscReal reject_safety)436d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetSafety(TSAdapt adapt, PetscReal safety, PetscReal reject_safety)
437d71ae5a4SJacob Faibussowitsch {
4381917a363SLisandro Dalcin PetscFunctionBegin;
4391917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
4401917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, safety, 2);
4411917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, reject_safety, 3);
44209cb0f53SBarry Smith PetscCheck(safety == (PetscReal)PETSC_CURRENT || safety >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Safety factor %g must be non negative", (double)safety);
44309cb0f53SBarry Smith PetscCheck(safety == (PetscReal)PETSC_CURRENT || safety <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Safety factor %g must be less than one", (double)safety);
44409cb0f53SBarry Smith PetscCheck(reject_safety == (PetscReal)PETSC_CURRENT || reject_safety >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Reject safety factor %g must be non negative", (double)reject_safety);
44509cb0f53SBarry Smith PetscCheck(reject_safety == (PetscReal)PETSC_CURRENT || reject_safety <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Reject safety factor %g must be less than one", (double)reject_safety);
44609cb0f53SBarry Smith if (safety != (PetscReal)PETSC_CURRENT) adapt->safety = safety;
44709cb0f53SBarry Smith if (reject_safety != (PetscReal)PETSC_CURRENT) adapt->reject_safety = reject_safety;
4483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4491917a363SLisandro Dalcin }
4501917a363SLisandro Dalcin
4511917a363SLisandro Dalcin /*@
452bcf0153eSBarry Smith TSAdaptGetSafety - Get safety factors for time step adapter
4531917a363SLisandro Dalcin
4541917a363SLisandro Dalcin Not Collective
4551917a363SLisandro Dalcin
4564165533cSJose E. Roman Input Parameter:
4571917a363SLisandro Dalcin . adapt - adaptive controller context
4581917a363SLisandro Dalcin
4594165533cSJose E. Roman Output Parameters:
460b43aa488SJacob Faibussowitsch + safety - safety factor relative to target error/stability goal
461b43aa488SJacob Faibussowitsch - reject_safety - extra safety factor to apply if the last step was rejected
4621917a363SLisandro Dalcin
4631917a363SLisandro Dalcin Level: intermediate
4641917a363SLisandro Dalcin
465*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSAdaptSetSafety()`, `TSAdaptChoose()`
4661917a363SLisandro Dalcin @*/
TSAdaptGetSafety(TSAdapt adapt,PetscReal * safety,PetscReal * reject_safety)467d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetSafety(TSAdapt adapt, PetscReal *safety, PetscReal *reject_safety)
468d71ae5a4SJacob Faibussowitsch {
4691917a363SLisandro Dalcin PetscFunctionBegin;
4701917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
4714f572ea9SToby Isaac if (safety) PetscAssertPointer(safety, 2);
4724f572ea9SToby Isaac if (reject_safety) PetscAssertPointer(reject_safety, 3);
4731917a363SLisandro Dalcin if (safety) *safety = adapt->safety;
4741917a363SLisandro Dalcin if (reject_safety) *reject_safety = adapt->reject_safety;
4753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4761917a363SLisandro Dalcin }
4771917a363SLisandro Dalcin
4781917a363SLisandro Dalcin /*@
479bcf0153eSBarry Smith TSAdaptSetMaxIgnore - Set error estimation threshold. Solution components below this threshold value will not be considered when computing error norms
480bcf0153eSBarry Smith for time step adaptivity (in absolute value). A negative value (default) of the threshold leads to considering all solution components.
48176cddca1SEmil Constantinescu
4822fe279fdSBarry Smith Logically Collective
48376cddca1SEmil Constantinescu
4844165533cSJose E. Roman Input Parameters:
48576cddca1SEmil Constantinescu + adapt - adaptive controller context
48676cddca1SEmil Constantinescu - max_ignore - threshold for solution components that are ignored during error estimation
48776cddca1SEmil Constantinescu
488bcf0153eSBarry Smith Options Database Key:
48976cddca1SEmil Constantinescu . -ts_adapt_max_ignore <max_ignore> - to set the threshold
49076cddca1SEmil Constantinescu
49176cddca1SEmil Constantinescu Level: intermediate
49276cddca1SEmil Constantinescu
493*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSAdaptGetMaxIgnore()`, `TSAdaptChoose()`
49476cddca1SEmil Constantinescu @*/
TSAdaptSetMaxIgnore(TSAdapt adapt,PetscReal max_ignore)495d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt adapt, PetscReal max_ignore)
496d71ae5a4SJacob Faibussowitsch {
49776cddca1SEmil Constantinescu PetscFunctionBegin;
49876cddca1SEmil Constantinescu PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
49976cddca1SEmil Constantinescu PetscValidLogicalCollectiveReal(adapt, max_ignore, 2);
50076cddca1SEmil Constantinescu adapt->ignore_max = max_ignore;
5013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
50276cddca1SEmil Constantinescu }
50376cddca1SEmil Constantinescu
50476cddca1SEmil Constantinescu /*@
505bcf0153eSBarry Smith TSAdaptGetMaxIgnore - Get error estimation threshold. Solution components below this threshold value will not be considered when computing error norms
506bcf0153eSBarry Smith for time step adaptivity (in absolute value).
50776cddca1SEmil Constantinescu
50876cddca1SEmil Constantinescu Not Collective
50976cddca1SEmil Constantinescu
5104165533cSJose E. Roman Input Parameter:
51176cddca1SEmil Constantinescu . adapt - adaptive controller context
51276cddca1SEmil Constantinescu
5134165533cSJose E. Roman Output Parameter:
51476cddca1SEmil Constantinescu . max_ignore - threshold for solution components that are ignored during error estimation
51576cddca1SEmil Constantinescu
51676cddca1SEmil Constantinescu Level: intermediate
51776cddca1SEmil Constantinescu
518*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSAdaptSetMaxIgnore()`, `TSAdaptChoose()`
51976cddca1SEmil Constantinescu @*/
TSAdaptGetMaxIgnore(TSAdapt adapt,PetscReal * max_ignore)520d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt adapt, PetscReal *max_ignore)
521d71ae5a4SJacob Faibussowitsch {
52276cddca1SEmil Constantinescu PetscFunctionBegin;
52376cddca1SEmil Constantinescu PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
5244f572ea9SToby Isaac PetscAssertPointer(max_ignore, 2);
52576cddca1SEmil Constantinescu *max_ignore = adapt->ignore_max;
5263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
52776cddca1SEmil Constantinescu }
52876cddca1SEmil Constantinescu
52976cddca1SEmil Constantinescu /*@
530bcf0153eSBarry Smith TSAdaptSetClip - Sets the admissible decrease/increase factor in step size in the time step adapter
5311917a363SLisandro Dalcin
532c3339decSBarry Smith Logically collective
5331917a363SLisandro Dalcin
5344165533cSJose E. Roman Input Parameters:
5351917a363SLisandro Dalcin + adapt - adaptive controller context
5361917a363SLisandro Dalcin . low - admissible decrease factor
537ec18b7bcSLisandro Dalcin - high - admissible increase factor
5381917a363SLisandro Dalcin
539bcf0153eSBarry Smith Options Database Key:
540ec18b7bcSLisandro Dalcin . -ts_adapt_clip <low>,<high> - to set admissible time step decrease and increase factors
5411917a363SLisandro Dalcin
5421917a363SLisandro Dalcin Level: intermediate
5431917a363SLisandro Dalcin
54409cb0f53SBarry Smith Note:
54509cb0f53SBarry Smith Use `PETSC_CURRENT` to keep the current value for either parameter
54609cb0f53SBarry Smith
54709cb0f53SBarry Smith Fortran Note:
54809cb0f53SBarry Smith Use `PETSC_CURRENT_REAL`
54909cb0f53SBarry Smith
550*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptGetClip()`, `TSAdaptSetScaleSolveFailed()`
5511917a363SLisandro Dalcin @*/
TSAdaptSetClip(TSAdapt adapt,PetscReal low,PetscReal high)552d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetClip(TSAdapt adapt, PetscReal low, PetscReal high)
553d71ae5a4SJacob Faibussowitsch {
5541917a363SLisandro Dalcin PetscFunctionBegin;
5551917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
5561917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, low, 2);
5571917a363SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, high, 3);
55809cb0f53SBarry Smith PetscCheck(low == (PetscReal)PETSC_CURRENT || low >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Decrease factor %g must be non negative", (double)low);
55909cb0f53SBarry Smith PetscCheck(low == (PetscReal)PETSC_CURRENT || low <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Decrease factor %g must be less than one", (double)low);
56009cb0f53SBarry Smith PetscCheck(high == (PetscReal)PETSC_CURRENT || high >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Increase factor %g must be greater than one", (double)high);
56109cb0f53SBarry Smith if (low != (PetscReal)PETSC_CURRENT) adapt->clip[0] = low;
56209cb0f53SBarry Smith if (high != (PetscReal)PETSC_CURRENT) adapt->clip[1] = high;
5633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5641917a363SLisandro Dalcin }
5651917a363SLisandro Dalcin
5661917a363SLisandro Dalcin /*@
567bcf0153eSBarry Smith TSAdaptGetClip - Gets the admissible decrease/increase factor in step size in the time step adapter
5681917a363SLisandro Dalcin
5691917a363SLisandro Dalcin Not Collective
5701917a363SLisandro Dalcin
5714165533cSJose E. Roman Input Parameter:
5721917a363SLisandro Dalcin . adapt - adaptive controller context
5731917a363SLisandro Dalcin
5744165533cSJose E. Roman Output Parameters:
5751917a363SLisandro Dalcin + low - optional, admissible decrease factor
5761917a363SLisandro Dalcin - high - optional, admissible increase factor
5771917a363SLisandro Dalcin
5781917a363SLisandro Dalcin Level: intermediate
5791917a363SLisandro Dalcin
580*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()`
5811917a363SLisandro Dalcin @*/
TSAdaptGetClip(TSAdapt adapt,PetscReal * low,PetscReal * high)582d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetClip(TSAdapt adapt, PetscReal *low, PetscReal *high)
583d71ae5a4SJacob Faibussowitsch {
5841917a363SLisandro Dalcin PetscFunctionBegin;
5851917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
5864f572ea9SToby Isaac if (low) PetscAssertPointer(low, 2);
5874f572ea9SToby Isaac if (high) PetscAssertPointer(high, 3);
5881917a363SLisandro Dalcin if (low) *low = adapt->clip[0];
5891917a363SLisandro Dalcin if (high) *high = adapt->clip[1];
5903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5911917a363SLisandro Dalcin }
5921917a363SLisandro Dalcin
5931917a363SLisandro Dalcin /*@
594bcf0153eSBarry Smith TSAdaptSetScaleSolveFailed - Scale step size by this factor if solve fails
59562c23b28SMark
59620f4b53cSBarry Smith Logically Collective
59762c23b28SMark
5984165533cSJose E. Roman Input Parameters:
59962c23b28SMark + adapt - adaptive controller context
600ee300463SSatish Balay - scale - scale
60162c23b28SMark
602bcf0153eSBarry Smith Options Database Key:
60362c23b28SMark . -ts_adapt_scale_solve_failed <scale> - to set scale step by this factor if solve fails
60462c23b28SMark
60562c23b28SMark Level: intermediate
60662c23b28SMark
607*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptGetScaleSolveFailed()`, `TSAdaptGetClip()`
60862c23b28SMark @*/
TSAdaptSetScaleSolveFailed(TSAdapt adapt,PetscReal scale)609d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt adapt, PetscReal scale)
610d71ae5a4SJacob Faibussowitsch {
61162c23b28SMark PetscFunctionBegin;
61262c23b28SMark PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
61362c23b28SMark PetscValidLogicalCollectiveReal(adapt, scale, 2);
61409cb0f53SBarry Smith PetscCheck(scale > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Scale factor %g must be positive", (double)scale);
61509cb0f53SBarry Smith PetscCheck(scale <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Scale factor %g must be less than one", (double)scale);
61609cb0f53SBarry Smith adapt->scale_solve_failed = scale;
6173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
61862c23b28SMark }
61962c23b28SMark
62062c23b28SMark /*@
62162c23b28SMark TSAdaptGetScaleSolveFailed - Gets the admissible decrease/increase factor in step size
62262c23b28SMark
62362c23b28SMark Not Collective
62462c23b28SMark
6254165533cSJose E. Roman Input Parameter:
62662c23b28SMark . adapt - adaptive controller context
62762c23b28SMark
6284165533cSJose E. Roman Output Parameter:
629ee300463SSatish Balay . scale - scale factor
63062c23b28SMark
63162c23b28SMark Level: intermediate
63262c23b28SMark
633*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSAdaptChoose()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetClip()`
63462c23b28SMark @*/
TSAdaptGetScaleSolveFailed(TSAdapt adapt,PetscReal * scale)635d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt adapt, PetscReal *scale)
636d71ae5a4SJacob Faibussowitsch {
63762c23b28SMark PetscFunctionBegin;
63862c23b28SMark PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
6394f572ea9SToby Isaac if (scale) PetscAssertPointer(scale, 2);
64062c23b28SMark if (scale) *scale = adapt->scale_solve_failed;
6413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
64262c23b28SMark }
64362c23b28SMark
64462c23b28SMark /*@
645bcf0153eSBarry Smith TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the time step controller
6461917a363SLisandro Dalcin
64720f4b53cSBarry Smith Logically Collective
6481c3436cfSJed Brown
6494165533cSJose E. Roman Input Parameters:
650bcf0153eSBarry Smith + adapt - time step adaptivity context, usually gotten with `TSGetAdapt()`
6511c3436cfSJed Brown . hmin - minimum time step
6521c3436cfSJed Brown - hmax - maximum time step
6531c3436cfSJed Brown
6541c3436cfSJed Brown Options Database Keys:
655ec18b7bcSLisandro Dalcin + -ts_adapt_dt_min <min> - to set minimum time step
656ec18b7bcSLisandro Dalcin - -ts_adapt_dt_max <max> - to set maximum time step
6571c3436cfSJed Brown
6581c3436cfSJed Brown Level: intermediate
6591c3436cfSJed Brown
66009cb0f53SBarry Smith Note:
66109cb0f53SBarry Smith Use `PETSC_CURRENT` to keep the current value for either parameter
66209cb0f53SBarry Smith
66309cb0f53SBarry Smith Fortran Note:
66409cb0f53SBarry Smith Use `PETSC_CURRENT_REAL`
66509cb0f53SBarry Smith
666*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSAdaptGetStepLimits()`, `TSAdaptChoose()`
6671c3436cfSJed Brown @*/
TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax)668d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt, PetscReal hmin, PetscReal hmax)
669d71ae5a4SJacob Faibussowitsch {
6701c3436cfSJed Brown PetscFunctionBegin;
6714782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
672b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, hmin, 2);
673b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, hmax, 3);
67409cb0f53SBarry Smith PetscCheck(hmin == (PetscReal)PETSC_CURRENT || hmin >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum time step %g must be non negative", (double)hmin);
67509cb0f53SBarry Smith PetscCheck(hmax == (PetscReal)PETSC_CURRENT || hmax >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum time step %g must be non negative", (double)hmax);
67609cb0f53SBarry Smith if (hmin != (PetscReal)PETSC_CURRENT) adapt->dt_min = hmin;
67709cb0f53SBarry Smith if (hmax != (PetscReal)PETSC_CURRENT) adapt->dt_max = hmax;
678b1f5bccaSLisandro Dalcin hmin = adapt->dt_min;
679b1f5bccaSLisandro Dalcin hmax = adapt->dt_max;
6803c633725SBarry 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);
6813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6821917a363SLisandro Dalcin }
6831917a363SLisandro Dalcin
6841917a363SLisandro Dalcin /*@
685bcf0153eSBarry Smith TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the time step controller
6861917a363SLisandro Dalcin
6871917a363SLisandro Dalcin Not Collective
6881917a363SLisandro Dalcin
6894165533cSJose E. Roman Input Parameter:
690bcf0153eSBarry Smith . adapt - time step adaptivity context, usually gotten with `TSGetAdapt()`
6911917a363SLisandro Dalcin
6924165533cSJose E. Roman Output Parameters:
6931917a363SLisandro Dalcin + hmin - minimum time step
6941917a363SLisandro Dalcin - hmax - maximum time step
6951917a363SLisandro Dalcin
6961917a363SLisandro Dalcin Level: intermediate
6971917a363SLisandro Dalcin
698*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSAdaptSetStepLimits()`, `TSAdaptChoose()`
6991917a363SLisandro Dalcin @*/
TSAdaptGetStepLimits(TSAdapt adapt,PetscReal * hmin,PetscReal * hmax)700d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt, PetscReal *hmin, PetscReal *hmax)
701d71ae5a4SJacob Faibussowitsch {
7021917a363SLisandro Dalcin PetscFunctionBegin;
7031917a363SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
7044f572ea9SToby Isaac if (hmin) PetscAssertPointer(hmin, 2);
7054f572ea9SToby Isaac if (hmax) PetscAssertPointer(hmax, 3);
7061917a363SLisandro Dalcin if (hmin) *hmin = adapt->dt_min;
7071917a363SLisandro Dalcin if (hmax) *hmax = adapt->dt_max;
7083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7091c3436cfSJed Brown }
7101c3436cfSJed Brown
71114d0ab18SJacob Faibussowitsch /*@C
712bcf0153eSBarry Smith TSAdaptSetFromOptions - Sets various `TSAdapt` parameters from user options.
71384df9cb4SJed Brown
714c3339decSBarry Smith Collective
71584df9cb4SJed Brown
71614d0ab18SJacob Faibussowitsch Input Parameters:
7172fe279fdSBarry Smith + adapt - the `TSAdapt` context
7182fe279fdSBarry Smith - PetscOptionsObject - object created by `PetscOptionsBegin()`
71984df9cb4SJed Brown
72084df9cb4SJed Brown Options Database Keys:
7211917a363SLisandro Dalcin + -ts_adapt_type <type> - algorithm to use for adaptivity
722bf997491SLisandro Dalcin . -ts_adapt_always_accept - always accept steps regardless of error/stability goals
7231917a363SLisandro Dalcin . -ts_adapt_safety <safety> - safety factor relative to target error/stability goal
7241917a363SLisandro Dalcin . -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected
7251917a363SLisandro Dalcin . -ts_adapt_clip <low,high> - admissible time step decrease and increase factors
72623de1d84SBarry Smith . -ts_adapt_dt_min <min> - minimum timestep to use
72723de1d84SBarry Smith . -ts_adapt_dt_max <max> - maximum timestep to use
72823de1d84SBarry Smith . -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails
729de50f1caSBarry Smith . -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates
730de50f1caSBarry 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
73184df9cb4SJed Brown
73284df9cb4SJed Brown Level: advanced
73384df9cb4SJed Brown
734bcf0153eSBarry Smith Note:
735bcf0153eSBarry Smith This function is automatically called by `TSSetFromOptions()`
73684df9cb4SJed Brown
737*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptSetAlwaysAccept()`, `TSAdaptSetSafety()`,
738db781477SPatrick Sanan `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetStepLimits()`, `TSAdaptSetMonitor()`
73914d0ab18SJacob Faibussowitsch @*/
TSAdaptSetFromOptions(TSAdapt adapt,PetscOptionItems PetscOptionsObject)740ce78bad3SBarry Smith PetscErrorCode TSAdaptSetFromOptions(TSAdapt adapt, PetscOptionItems PetscOptionsObject)
741d71ae5a4SJacob Faibussowitsch {
74284df9cb4SJed Brown char type[256] = TSADAPTBASIC;
74362c23b28SMark PetscReal safety, reject_safety, clip[2], scale, hmin, hmax;
7441c3436cfSJed Brown PetscBool set, flg;
7451917a363SLisandro Dalcin PetscInt two;
74684df9cb4SJed Brown
74784df9cb4SJed Brown PetscFunctionBegin;
748dbbe0bcdSBarry Smith PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
74984df9cb4SJed Brown /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
7501566a47fSLisandro Dalcin * function can only be called from inside TSSetFromOptions() */
751d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "TS Adaptivity options");
7529566063dSJacob 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));
75348a46eb9SPierre Jolivet if (flg || !((PetscObject)adapt)->type_name) PetscCall(TSAdaptSetType(adapt, type));
7541917a363SLisandro Dalcin
7559566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_adapt_always_accept", "Always accept the step", "TSAdaptSetAlwaysAccept", adapt->always_accept, &flg, &set));
7569566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetAlwaysAccept(adapt, flg));
7571917a363SLisandro Dalcin
7589371c9d4SSatish Balay safety = adapt->safety;
7599371c9d4SSatish Balay reject_safety = adapt->reject_safety;
7609566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_safety", "Safety factor relative to target error/stability goal", "TSAdaptSetSafety", safety, &safety, &set));
7619566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_reject_safety", "Extra safety factor to apply if the last step was rejected", "TSAdaptSetSafety", reject_safety, &reject_safety, &flg));
7629566063dSJacob Faibussowitsch if (set || flg) PetscCall(TSAdaptSetSafety(adapt, safety, reject_safety));
7631917a363SLisandro Dalcin
7649371c9d4SSatish Balay two = 2;
7659371c9d4SSatish Balay clip[0] = adapt->clip[0];
7669371c9d4SSatish Balay clip[1] = adapt->clip[1];
7679566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-ts_adapt_clip", "Admissible decrease/increase factor in step size", "TSAdaptSetClip", clip, &two, &set));
7683c633725SBarry Smith PetscCheck(!set || (two == 2), PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Must give exactly two values to -ts_adapt_clip");
7699566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetClip(adapt, clip[0], clip[1]));
7701917a363SLisandro Dalcin
7719371c9d4SSatish Balay hmin = adapt->dt_min;
7729371c9d4SSatish Balay hmax = adapt->dt_max;
7739566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_dt_min", "Minimum time step considered", "TSAdaptSetStepLimits", hmin, &hmin, &set));
7749566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_dt_max", "Maximum time step considered", "TSAdaptSetStepLimits", hmax, &hmax, &flg));
7759566063dSJacob Faibussowitsch if (set || flg) PetscCall(TSAdaptSetStepLimits(adapt, hmin, hmax));
7761917a363SLisandro Dalcin
7779566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_max_ignore", "Adaptor ignores (absolute) solution values smaller than this value", "", adapt->ignore_max, &adapt->ignore_max, &set));
7789566063dSJacob 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));
779d580f011SEmil Constantinescu
7809566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_adapt_scale_solve_failed", "Scale step by this factor if solve fails", "TSAdaptSetScaleSolveFailed", adapt->scale_solve_failed, &scale, &set));
7819566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetScaleSolveFailed(adapt, scale));
7821917a363SLisandro Dalcin
7839566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-ts_adapt_wnormtype", "Type of norm computed for error estimation", "", NormTypes, (PetscEnum)adapt->wnormtype, (PetscEnum *)&adapt->wnormtype, NULL));
7843c633725SBarry Smith PetscCheck(adapt->wnormtype == NORM_2 || adapt->wnormtype == NORM_INFINITY, PetscObjectComm((PetscObject)adapt), PETSC_ERR_SUP, "Only 2-norm and infinite norm supported");
7851917a363SLisandro Dalcin
7869566063dSJacob 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));
787de50f1caSBarry Smith
7889566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_adapt_monitor", "Print choices made by adaptive controller", "TSAdaptSetMonitor", adapt->monitor ? PETSC_TRUE : PETSC_FALSE, &flg, &set));
7899566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptSetMonitor(adapt, flg));
7901917a363SLisandro Dalcin
791dbbe0bcdSBarry Smith PetscTryTypeMethod(adapt, setfromoptions, PetscOptionsObject);
792d0609cedSBarry Smith PetscOptionsHeadEnd();
7933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
79484df9cb4SJed Brown }
79584df9cb4SJed Brown
79684df9cb4SJed Brown /*@
79784df9cb4SJed Brown TSAdaptCandidatesClear - clear any previously set candidate schemes
79884df9cb4SJed Brown
79920f4b53cSBarry Smith Logically Collective
80084df9cb4SJed Brown
8014165533cSJose E. Roman Input Parameter:
80284df9cb4SJed Brown . adapt - adaptive controller
80384df9cb4SJed Brown
80484df9cb4SJed Brown Level: developer
80584df9cb4SJed Brown
806*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSAdaptCreate()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()`
80784df9cb4SJed Brown @*/
TSAdaptCandidatesClear(TSAdapt adapt)808d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
809d71ae5a4SJacob Faibussowitsch {
81084df9cb4SJed Brown PetscFunctionBegin;
8114782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
8129566063dSJacob Faibussowitsch PetscCall(PetscMemzero(&adapt->candidates, sizeof(adapt->candidates)));
8133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
81484df9cb4SJed Brown }
81584df9cb4SJed Brown
81684df9cb4SJed Brown /*@C
81784df9cb4SJed Brown TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from
81884df9cb4SJed Brown
81920f4b53cSBarry Smith Logically Collective; No Fortran Support
82084df9cb4SJed Brown
8214165533cSJose E. Roman Input Parameters:
822bcf0153eSBarry Smith + adapt - time step adaptivity context, obtained with `TSGetAdapt()` or `TSAdaptCreate()`
82384df9cb4SJed Brown . name - name of the candidate scheme to add
82484df9cb4SJed Brown . order - order of the candidate scheme
82584df9cb4SJed Brown . stageorder - stage order of the candidate scheme
8268d59e960SJed Brown . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints
82784df9cb4SJed Brown . cost - relative measure of the amount of work required for the candidate scheme
82884df9cb4SJed Brown - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme
82984df9cb4SJed Brown
83084df9cb4SJed Brown Level: developer
83184df9cb4SJed Brown
832*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptChoose()`
83384df9cb4SJed Brown @*/
TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse)834d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt, const char name[], PetscInt order, PetscInt stageorder, PetscReal ccfl, PetscReal cost, PetscBool inuse)
835d71ae5a4SJacob Faibussowitsch {
83684df9cb4SJed Brown PetscInt c;
83784df9cb4SJed Brown
83884df9cb4SJed Brown PetscFunctionBegin;
83984df9cb4SJed Brown PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
84063a3b9bcSJacob Faibussowitsch PetscCheck(order >= 1, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Classical order %" PetscInt_FMT " must be a positive integer", order);
84184df9cb4SJed Brown if (inuse) {
8423c633725SBarry Smith PetscCheck(!adapt->candidates.inuse_set, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONGSTATE, "Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()");
84384df9cb4SJed Brown adapt->candidates.inuse_set = PETSC_TRUE;
84484df9cb4SJed Brown }
8451c3436cfSJed Brown /* first slot if this is the current scheme, otherwise the next available slot */
8461c3436cfSJed Brown c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n;
847bbd56ea5SKarl Rupp
84884df9cb4SJed Brown adapt->candidates.name[c] = name;
84984df9cb4SJed Brown adapt->candidates.order[c] = order;
85084df9cb4SJed Brown adapt->candidates.stageorder[c] = stageorder;
8518d59e960SJed Brown adapt->candidates.ccfl[c] = ccfl;
85284df9cb4SJed Brown adapt->candidates.cost[c] = cost;
85384df9cb4SJed Brown adapt->candidates.n++;
8543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
85584df9cb4SJed Brown }
85684df9cb4SJed Brown
8578d59e960SJed Brown /*@C
8588d59e960SJed Brown TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost
8598d59e960SJed Brown
8608d59e960SJed Brown Not Collective
8618d59e960SJed Brown
8624165533cSJose E. Roman Input Parameter:
8638d59e960SJed Brown . adapt - time step adaptivity context
8648d59e960SJed Brown
8654165533cSJose E. Roman Output Parameters:
8668d59e960SJed Brown + n - number of candidate schemes, always at least 1
8678d59e960SJed Brown . order - the order of each candidate scheme
8688d59e960SJed Brown . stageorder - the stage order of each candidate scheme
8698d59e960SJed Brown . ccfl - the CFL coefficient of each scheme
8708d59e960SJed Brown - cost - the relative cost of each scheme
8718d59e960SJed Brown
8728d59e960SJed Brown Level: developer
8738d59e960SJed Brown
8748d59e960SJed Brown Note:
8758d59e960SJed Brown The current scheme is always returned in the first slot
8768d59e960SJed Brown
877*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()`
8788d59e960SJed Brown @*/
TSAdaptCandidatesGet(TSAdapt adapt,PetscInt * n,const PetscInt ** order,const PetscInt ** stageorder,const PetscReal ** ccfl,const PetscReal ** cost)879d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt, PetscInt *n, const PetscInt **order, const PetscInt **stageorder, const PetscReal **ccfl, const PetscReal **cost)
880d71ae5a4SJacob Faibussowitsch {
8818d59e960SJed Brown PetscFunctionBegin;
8828d59e960SJed Brown PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
8838d59e960SJed Brown if (n) *n = adapt->candidates.n;
8848d59e960SJed Brown if (order) *order = adapt->candidates.order;
8858d59e960SJed Brown if (stageorder) *stageorder = adapt->candidates.stageorder;
8868d59e960SJed Brown if (ccfl) *ccfl = adapt->candidates.ccfl;
8878d59e960SJed Brown if (cost) *cost = adapt->candidates.cost;
8883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
8898d59e960SJed Brown }
8908d59e960SJed Brown
89184df9cb4SJed Brown /*@C
89284df9cb4SJed Brown TSAdaptChoose - choose which method and step size to use for the next step
89384df9cb4SJed Brown
894c3339decSBarry Smith Collective
89584df9cb4SJed Brown
8964165533cSJose E. Roman Input Parameters:
897da81f932SPierre Jolivet + adapt - adaptive controller
89897bb3fdcSJose E. Roman . ts - time stepper
89984df9cb4SJed Brown - h - current step size
90084df9cb4SJed Brown
9014165533cSJose E. Roman Output Parameters:
9021566a47fSLisandro Dalcin + next_sc - optional, scheme to use for the next step
90384df9cb4SJed Brown . next_h - step size to use for the next step
904bcf0153eSBarry Smith - accept - `PETSC_TRUE` to accept the current step, `PETSC_FALSE` to repeat the current step with the new step size
9051c3436cfSJed Brown
90684df9cb4SJed Brown Level: developer
90784df9cb4SJed Brown
908bcf0153eSBarry Smith Note:
909bcf0153eSBarry Smith The input value of parameter accept is retained from the last time step, so it will be `PETSC_FALSE` if the step is
910bcf0153eSBarry Smith being retried after an initial rejection.
911bcf0153eSBarry Smith
912*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()`
91384df9cb4SJed Brown @*/
TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt * next_sc,PetscReal * next_h,PetscBool * accept)914d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptChoose(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept)
915d71ae5a4SJacob Faibussowitsch {
9161566a47fSLisandro Dalcin PetscInt ncandidates = adapt->candidates.n;
9171566a47fSLisandro Dalcin PetscInt scheme = 0;
9180b99f514SJed Brown PetscReal wlte = -1.0;
9195e4ed32fSEmil Constantinescu PetscReal wltea = -1.0;
9205e4ed32fSEmil Constantinescu PetscReal wlter = -1.0;
92184df9cb4SJed Brown
92284df9cb4SJed Brown PetscFunctionBegin;
92384df9cb4SJed Brown PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
92484df9cb4SJed Brown PetscValidHeaderSpecific(ts, TS_CLASSID, 2);
9254f572ea9SToby Isaac if (next_sc) PetscAssertPointer(next_sc, 4);
9264f572ea9SToby Isaac PetscAssertPointer(next_h, 5);
9274f572ea9SToby Isaac PetscAssertPointer(accept, 6);
9281566a47fSLisandro Dalcin if (next_sc) *next_sc = 0;
9291566a47fSLisandro Dalcin
9301566a47fSLisandro Dalcin /* Do not mess with adaptivity while handling events */
931ca4445c7SIlya Fursov if (ts->event && ts->event->processing) {
9321566a47fSLisandro Dalcin *next_h = h;
9331566a47fSLisandro Dalcin *accept = PETSC_TRUE;
934fe4ad979SIlya Fursov if (adapt->monitor) {
935fe4ad979SIlya Fursov PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
936fe4ad979SIlya Fursov
937fe4ad979SIlya Fursov if (ts->event->iterctr == 0) {
938fe4ad979SIlya Fursov /*
939fe4ad979SIlya Fursov An event has been found, now finalising the event processing: performing the 1st and 2nd post-event steps.
940fe4ad979SIlya Fursov Entering this if-branch means both these steps (set to either PETSC_DECIDE or numerical value) are managed
941fe4ad979SIlya Fursov by the event handler. In this case the 1st post-event step is always accepted, without interference of TSAdapt.
942fe4ad979SIlya Fursov Note: if the 2nd post-event step is not managed by the event handler (e.g. given 1st = numerical, 2nd = PETSC_DECIDE),
943fe4ad979SIlya Fursov this if-branch is not entered, and TSAdapt may reject/adjust the proposed 1st post-event step.
944fe4ad979SIlya Fursov */
945fe4ad979SIlya Fursov PetscCall(PetscViewerASCIIPrintf(adapt->monitor, " TSAdapt does not interfere, step %3" PetscInt_FMT " accepted. Processing post-event steps: 1-st accepted just now, 2-nd yet to come\n", ts->steps));
946fe4ad979SIlya Fursov } else PetscCall(PetscViewerASCIIPrintf(adapt->monitor, " TSAdapt does not interfere, step %3" PetscInt_FMT " accepted. Event handling in progress\n", ts->steps));
947fe4ad979SIlya Fursov
948fe4ad979SIlya Fursov PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
949fe4ad979SIlya Fursov }
9503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9511566a47fSLisandro Dalcin }
9521566a47fSLisandro Dalcin
953dbbe0bcdSBarry Smith PetscUseTypeMethod(adapt, choose, ts, h, &scheme, next_h, accept, &wlte, &wltea, &wlter);
95463a3b9bcSJacob 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);
9553c633725SBarry Smith PetscCheck(*next_h >= 0, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Computed step size %g must be positive", (double)*next_h);
9561566a47fSLisandro Dalcin if (next_sc) *next_sc = scheme;
9571566a47fSLisandro Dalcin
95849354f04SShri Abhyankar if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
95936b54a69SLisandro Dalcin /* Increase/reduce step size if end time of next step is close to or overshoots max time */
960ca4445c7SIlya Fursov PetscReal t = ts->ptime + ts->time_step, tend, tmax, h1, hmax;
961fe7350e0SStefano Zampini PetscReal a = (PetscReal)(1.0 + adapt->matchstepfac[0]);
962fe7350e0SStefano Zampini PetscReal b = adapt->matchstepfac[1];
9634a658b32SHong Zhang
964fe4ad979SIlya Fursov /*
965fe4ad979SIlya Fursov Logic in using 'dt_span_cached':
966fe4ad979SIlya Fursov 1. It always overrides *next_h, except (any of):
967fe4ad979SIlya Fursov a) the current step was rejected,
968fe4ad979SIlya Fursov b) the adaptor proposed to decrease the next step,
969fe4ad979SIlya Fursov c) the adaptor proposed *next_h > dt_span_cached.
970136cf249SJames Wright 2. If *next_h was adjusted by eval_times points (or the final point):
971fe4ad979SIlya Fursov -- when dt_span_cached is filled (>0), it keeps its value,
972fe4ad979SIlya Fursov -- when dt_span_cached is clear (==0), it gets the unadjusted version of *next_h.
973fe4ad979SIlya Fursov 3. If *next_h was not adjusted as in (2), dt_span_cached is cleared.
974fe4ad979SIlya Fursov Note, if a combination (1.b || 1.c) && (3) takes place, this means that
975fe4ad979SIlya Fursov dt_span_cached remains unused at the moment of clearing.
976fe4ad979SIlya Fursov If (1.a) takes place, dt_span_cached keeps its value.
977fe4ad979SIlya Fursov Also, dt_span_cached can be updated by the event handler, see tsevent.c.
978fe4ad979SIlya Fursov */
979136cf249SJames Wright if (h <= *next_h && *next_h <= adapt->dt_eval_times_cached) *next_h = adapt->dt_eval_times_cached; /* try employing the cache */
980ca4445c7SIlya Fursov h1 = *next_h;
981ca4445c7SIlya Fursov tend = t + h1;
982ca4445c7SIlya Fursov
983136cf249SJames Wright if (ts->eval_times && ts->eval_times->time_point_idx < ts->eval_times->num_time_points) {
984136cf249SJames Wright PetscCheck(ts->eval_times->worktol == 0, PetscObjectComm((PetscObject)adapt), PETSC_ERR_PLIB, "Unexpected state (tspan->worktol != 0) in TSAdaptChoose()");
985136cf249SJames Wright ts->eval_times->worktol = ts->eval_times->reltol * h1 + ts->eval_times->abstol;
986136cf249SJames Wright if (PetscIsCloseAtTol(t, ts->eval_times->time_points[ts->eval_times->time_point_idx], ts->eval_times->worktol, 0)) /* hit a span time point */
987136cf249SJames Wright if (ts->eval_times->time_point_idx + 1 < ts->eval_times->num_time_points) tmax = ts->eval_times->time_points[ts->eval_times->time_point_idx + 1];
9884a658b32SHong Zhang else tmax = ts->max_time; /* hit the last span time point */
989136cf249SJames Wright else tmax = ts->eval_times->time_points[ts->eval_times->time_point_idx];
9904a658b32SHong Zhang } else tmax = ts->max_time;
991ca4445c7SIlya Fursov tmax = PetscMin(tmax, ts->max_time);
9924a658b32SHong Zhang hmax = tmax - t;
993ca4445c7SIlya Fursov
99436b54a69SLisandro Dalcin if (t < tmax && tend > tmax) *next_h = hmax;
995ca4445c7SIlya Fursov if (t < tmax && tend < tmax && h1 * b > hmax) *next_h = hmax / 2;
996ca4445c7SIlya Fursov if (t < tmax && tend < tmax && h1 * a > hmax) *next_h = hmax;
997136cf249SJames Wright if (ts->eval_times && h1 != *next_h && !adapt->dt_eval_times_cached) adapt->dt_eval_times_cached = h1; /* cache the step size if it is to be changed */
998136cf249SJames Wright if (ts->eval_times && h1 == *next_h && adapt->dt_eval_times_cached) adapt->dt_eval_times_cached = 0; /* clear the cache if the step size is unchanged */
999e1db57b0SHong Zhang }
10001c3436cfSJed Brown if (adapt->monitor) {
10011566a47fSLisandro Dalcin const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : "";
10029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
10030b99f514SJed Brown if (wlte < 0) {
10049371c9d4SSatish 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",
10059371c9d4SSatish Balay (double)ts->ptime, (double)h, (double)*next_h));
10060b99f514SJed Brown } else {
10079371c9d4SSatish 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",
10089371c9d4SSatish Balay (double)ts->ptime, (double)h, (double)*next_h, (double)wlte, (double)wltea, (double)wlter));
10090b99f514SJed Brown }
10109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
10111c3436cfSJed Brown }
10123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
101384df9cb4SJed Brown }
101484df9cb4SJed Brown
101597335746SJed Brown /*@
1016de50f1caSBarry Smith TSAdaptSetTimeStepIncreaseDelay - The number of timesteps to wait after a decrease in the timestep due to failed solver
1017de50f1caSBarry Smith before increasing the time step.
1018de50f1caSBarry Smith
1019c3339decSBarry Smith Logicially Collective
1020de50f1caSBarry Smith
10214165533cSJose E. Roman Input Parameters:
1022de50f1caSBarry Smith + adapt - adaptive controller context
1023de50f1caSBarry Smith - cnt - the number of timesteps
1024de50f1caSBarry Smith
1025de50f1caSBarry Smith Options Database Key:
1026de50f1caSBarry Smith . -ts_adapt_time_step_increase_delay cnt - number of steps to delay the increase
1027de50f1caSBarry Smith
1028de50f1caSBarry Smith Level: advanced
1029de50f1caSBarry Smith
1030bcf0153eSBarry Smith Notes:
1031bcf0153eSBarry Smith This is to prevent an adaptor from bouncing back and forth between two nearby timesteps. The default is 0.
1032bcf0153eSBarry Smith
1033bcf0153eSBarry Smith The successful use of this option is problem dependent
1034bcf0153eSBarry Smith
1035b43aa488SJacob Faibussowitsch Developer Notes:
1036bcf0153eSBarry Smith There is no theory to support this option
1037bcf0153eSBarry Smith
1038*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`
1039de50f1caSBarry Smith @*/
TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt,PetscInt cnt)1040d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt, PetscInt cnt)
1041d71ae5a4SJacob Faibussowitsch {
1042de50f1caSBarry Smith PetscFunctionBegin;
1043de50f1caSBarry Smith adapt->timestepjustdecreased_delay = cnt;
10443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1045de50f1caSBarry Smith }
1046de50f1caSBarry Smith
1047de50f1caSBarry Smith /*@
10486bc98fa9SBarry 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)
104997335746SJed Brown
1050c3339decSBarry Smith Collective
105197335746SJed Brown
10524165533cSJose E. Roman Input Parameters:
105397335746SJed Brown + adapt - adaptive controller context
1054b295832fSPierre Barbier de Reuille . ts - time stepper
1055b295832fSPierre Barbier de Reuille . t - Current simulation time
1056b295832fSPierre Barbier de Reuille - Y - Current solution vector
105797335746SJed Brown
10584165533cSJose E. Roman Output Parameter:
1059bcf0153eSBarry Smith . accept - `PETSC_TRUE` to accept the stage, `PETSC_FALSE` to reject
106097335746SJed Brown
106197335746SJed Brown Level: developer
106297335746SJed Brown
1063*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`
106497335746SJed Brown @*/
TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool * accept)1065d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCheckStage(TSAdapt adapt, TS ts, PetscReal t, Vec Y, PetscBool *accept)
1066d71ae5a4SJacob Faibussowitsch {
10671566a47fSLisandro Dalcin SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING;
10686c6709e3SStefano Zampini PetscBool func_accept, snes_div_func;
106997335746SJed Brown
107097335746SJed Brown PetscFunctionBegin;
10714782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
10724782b174SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 2);
10734f572ea9SToby Isaac PetscAssertPointer(accept, 5);
10741566a47fSLisandro Dalcin
10756c6709e3SStefano Zampini PetscCall(TSFunctionDomainError(ts, t, Y, &func_accept));
10769566063dSJacob Faibussowitsch if (ts->snes) PetscCall(SNESGetConvergedReason(ts->snes, &snesreason));
10776c6709e3SStefano Zampini snes_div_func = (PetscBool)(snesreason == SNES_DIVERGED_FUNCTION_DOMAIN);
10786c6709e3SStefano Zampini if (func_accept && snesreason < 0 && !snes_div_func) {
107997335746SJed Brown *accept = PETSC_FALSE;
10806c6709e3SStefano Zampini PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", nonlinear solve failure: %s\n", ts->steps, SNESConvergedReasons[snesreason]));
108109cb0f53SBarry Smith if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures != PETSC_UNLIMITED) {
108297335746SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
108363a3b9bcSJacob 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));
108497335746SJed Brown if (adapt->monitor) {
10859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
10869371c9d4SSatish 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,
10879371c9d4SSatish Balay (double)ts->ptime, (double)ts->time_step, ts->num_snes_failures));
10889566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
108997335746SJed Brown }
1090cb9d8021SPierre Barbier de Reuille }
1091cb9d8021SPierre Barbier de Reuille } else {
10926c6709e3SStefano Zampini *accept = (PetscBool)(func_accept && !snes_div_func);
10936c6709e3SStefano Zampini if (*accept && adapt->checkstage) PetscCall((*adapt->checkstage)(adapt, ts, t, Y, accept));
10946bc98fa9SBarry Smith if (!*accept) {
10956c6709e3SStefano Zampini const char *user_func = !func_accept ? "TSSetFunctionDomainError()" : "TSAdaptSetCheckStage";
10966c6709e3SStefano Zampini const char *snes_err = "SNES invalid function domain";
10976c6709e3SStefano Zampini const char *err_msg = snes_div_func && func_accept ? snes_err : user_func;
10986c6709e3SStefano Zampini PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", solution rejected by %s\n", ts->steps, err_msg));
10996bc98fa9SBarry Smith if (adapt->monitor) {
11009566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
11016c6709e3SStefano Zampini PetscCall(PetscViewerASCIIPrintf(adapt->monitor, " TSAdapt %s step %3" PetscInt_FMT " stage rejected by %s\n", ((PetscObject)adapt)->type_name, ts->steps, err_msg));
11029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
11036bc98fa9SBarry Smith }
11046bc98fa9SBarry Smith }
1105cb9d8021SPierre Barbier de Reuille }
1106cb9d8021SPierre Barbier de Reuille
110757508eceSPierre Jolivet if (!*accept && !ts->reason) {
11081566a47fSLisandro Dalcin PetscReal dt, new_dt;
11099566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt));
1110cb9d8021SPierre Barbier de Reuille new_dt = dt * adapt->scale_solve_failed;
11119566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, new_dt));
1112de50f1caSBarry Smith adapt->timestepjustdecreased += adapt->timestepjustdecreased_delay;
111397335746SJed Brown if (adapt->monitor) {
11149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
11156c6709e3SStefano 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],
11166c6709e3SStefano Zampini (double)ts->ptime, (double)dt, (double)new_dt));
11179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
111897335746SJed Brown }
111997335746SJed Brown }
11203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
112197335746SJed Brown }
112297335746SJed Brown
112384df9cb4SJed Brown /*@
112484df9cb4SJed Brown TSAdaptCreate - create an adaptive controller context for time stepping
112584df9cb4SJed Brown
1126d083f849SBarry Smith Collective
112784df9cb4SJed Brown
112884df9cb4SJed Brown Input Parameter:
112984df9cb4SJed Brown . comm - The communicator
113084df9cb4SJed Brown
113184df9cb4SJed Brown Output Parameter:
1132b43aa488SJacob Faibussowitsch . inadapt - new `TSAdapt` object
113384df9cb4SJed Brown
113484df9cb4SJed Brown Level: developer
113584df9cb4SJed Brown
1136bcf0153eSBarry Smith Note:
1137bcf0153eSBarry Smith `TSAdapt` creation is handled by `TS`, so users should not need to call this function.
113884df9cb4SJed Brown
1139*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TSAdapt`, `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptDestroy()`
114084df9cb4SJed Brown @*/
TSAdaptCreate(MPI_Comm comm,TSAdapt * inadapt)1141d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCreate(MPI_Comm comm, TSAdapt *inadapt)
1142d71ae5a4SJacob Faibussowitsch {
114384df9cb4SJed Brown TSAdapt adapt;
114484df9cb4SJed Brown
114584df9cb4SJed Brown PetscFunctionBegin;
11464f572ea9SToby Isaac PetscAssertPointer(inadapt, 2);
11479566063dSJacob Faibussowitsch PetscCall(TSAdaptInitializePackage());
11483b3bcf4cSLisandro Dalcin
11499566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(adapt, TSADAPT_CLASSID, "TSAdapt", "Time stepping adaptivity", "TS", comm, TSAdaptDestroy, TSAdaptView));
1150bf997491SLisandro Dalcin adapt->always_accept = PETSC_FALSE;
11511917a363SLisandro Dalcin adapt->safety = 0.9;
11521917a363SLisandro Dalcin adapt->reject_safety = 0.5;
11531917a363SLisandro Dalcin adapt->clip[0] = 0.1;
11541917a363SLisandro Dalcin adapt->clip[1] = 10.;
11551c3436cfSJed Brown adapt->dt_min = 1e-20;
11561917a363SLisandro Dalcin adapt->dt_max = 1e+20;
11571c167fc2SEmil Constantinescu adapt->ignore_max = -1.0;
1158d580f011SEmil Constantinescu adapt->glee_use_local = PETSC_TRUE;
115997335746SJed Brown adapt->scale_solve_failed = 0.25;
1160fe7350e0SStefano Zampini /* these two safety factors are not public, and they are used only in the TS_EXACTFINALTIME_MATCHSTEP case
1161fe7350e0SStefano Zampini to prevent from situations were unreasonably small time steps are taken in order to match the final time */
1162fe7350e0SStefano Zampini adapt->matchstepfac[0] = 0.01; /* allow 1% step size increase in the last step */
1163fe7350e0SStefano Zampini adapt->matchstepfac[1] = 2.0; /* halve last step if it is greater than what remains divided this factor */
11647619abb3SShri adapt->wnormtype = NORM_2;
1165de50f1caSBarry Smith adapt->timestepjustdecreased_delay = 0;
116684df9cb4SJed Brown *inadapt = adapt;
11673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
116884df9cb4SJed Brown }
1169