xref: /petsc/src/ts/adapt/interface/tsadapt.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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:
231c84c290SBarry Smith +  name_scheme - name of user-defined adaptivity scheme
241c84c290SBarry Smith -  routine_create - routine to create method context
251c84c290SBarry Smith 
261c84c290SBarry Smith    Notes:
271c84c290SBarry Smith    TSAdaptRegister() may be called multiple times to add several user-defined families.
281c84c290SBarry Smith 
291c84c290SBarry Smith    Sample usage:
301c84c290SBarry Smith .vb
31bdf89e91SBarry Smith    TSAdaptRegister("my_scheme",MySchemeCreate);
321c84c290SBarry Smith .ve
331c84c290SBarry Smith 
341c84c290SBarry Smith    Then, your scheme can be chosen with the procedural interface via
351c84c290SBarry Smith $     TSAdaptSetType(ts,"my_scheme")
361c84c290SBarry Smith    or at runtime via the option
371c84c290SBarry Smith $     -ts_adapt_type my_scheme
3884df9cb4SJed Brown 
3984df9cb4SJed Brown    Level: advanced
401c84c290SBarry Smith 
41db781477SPatrick Sanan .seealso: `TSAdaptRegisterAll()`
4284df9cb4SJed Brown @*/
43*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptRegister(const char sname[], PetscErrorCode (*function)(TSAdapt))
44*d71ae5a4SJacob Faibussowitsch {
4584df9cb4SJed Brown   PetscFunctionBegin;
469566063dSJacob Faibussowitsch   PetscCall(TSAdaptInitializePackage());
479566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&TSAdaptList, sname, function));
4884df9cb4SJed Brown   PetscFunctionReturn(0);
4984df9cb4SJed Brown }
5084df9cb4SJed Brown 
5184df9cb4SJed Brown /*@C
5284df9cb4SJed Brown   TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt
5384df9cb4SJed Brown 
5484df9cb4SJed Brown   Not Collective
5584df9cb4SJed Brown 
5684df9cb4SJed Brown   Level: advanced
5784df9cb4SJed Brown 
58db781477SPatrick Sanan .seealso: `TSAdaptRegisterDestroy()`
5984df9cb4SJed Brown @*/
60*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptRegisterAll(void)
61*d71ae5a4SJacob Faibussowitsch {
6284df9cb4SJed Brown   PetscFunctionBegin;
630f51fdf8SToby Isaac   if (TSAdaptRegisterAllCalled) PetscFunctionReturn(0);
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));
7184df9cb4SJed Brown   PetscFunctionReturn(0);
7284df9cb4SJed Brown }
7384df9cb4SJed Brown 
7484df9cb4SJed Brown /*@C
75560360afSLisandro Dalcin   TSAdaptFinalizePackage - This function destroys everything in the TS package. It is
7684df9cb4SJed Brown   called from PetscFinalize().
7784df9cb4SJed Brown 
7884df9cb4SJed Brown   Level: developer
7984df9cb4SJed Brown 
80db781477SPatrick Sanan .seealso: `PetscFinalize()`
8184df9cb4SJed Brown @*/
82*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptFinalizePackage(void)
83*d71ae5a4SJacob Faibussowitsch {
8484df9cb4SJed Brown   PetscFunctionBegin;
859566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&TSAdaptList));
8684df9cb4SJed Brown   TSAdaptPackageInitialized = PETSC_FALSE;
8784df9cb4SJed Brown   TSAdaptRegisterAllCalled  = PETSC_FALSE;
8884df9cb4SJed Brown   PetscFunctionReturn(0);
8984df9cb4SJed Brown }
9084df9cb4SJed Brown 
9184df9cb4SJed Brown /*@C
9284df9cb4SJed Brown   TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is
938a690491SBarry Smith   called from TSInitializePackage().
9484df9cb4SJed Brown 
9584df9cb4SJed Brown   Level: developer
9684df9cb4SJed Brown 
97db781477SPatrick Sanan .seealso: `PetscInitialize()`
9884df9cb4SJed Brown @*/
99*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptInitializePackage(void)
100*d71ae5a4SJacob Faibussowitsch {
10184df9cb4SJed Brown   PetscFunctionBegin;
10284df9cb4SJed Brown   if (TSAdaptPackageInitialized) PetscFunctionReturn(0);
10384df9cb4SJed Brown   TSAdaptPackageInitialized = PETSC_TRUE;
1049566063dSJacob Faibussowitsch   PetscCall(PetscClassIdRegister("TSAdapt", &TSADAPT_CLASSID));
1059566063dSJacob Faibussowitsch   PetscCall(TSAdaptRegisterAll());
1069566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(TSAdaptFinalizePackage));
10784df9cb4SJed Brown   PetscFunctionReturn(0);
10884df9cb4SJed Brown }
10984df9cb4SJed Brown 
1107eef6055SBarry Smith /*@C
1117eef6055SBarry Smith   TSAdaptSetType - sets the approach used for the error adapter, currently there is only TSADAPTBASIC and TSADAPTNONE
1127eef6055SBarry Smith 
1137eef6055SBarry Smith   Logicially Collective on TSAdapt
1147eef6055SBarry Smith 
115d8d19677SJose E. Roman   Input Parameters:
116d0288e4fSLisandro Dalcin + adapt - the TS adapter, most likely obtained with TSGetAdapt()
1177eef6055SBarry Smith - type - either  TSADAPTBASIC or TSADAPTNONE
1187eef6055SBarry Smith 
1197eef6055SBarry Smith   Options Database:
120ec18b7bcSLisandro Dalcin . -ts_adapt_type <basic or dsp or none> - to set the adapter type
1217eef6055SBarry Smith 
1227eef6055SBarry Smith   Level: intermediate
1237eef6055SBarry Smith 
124db781477SPatrick Sanan .seealso: `TSGetAdapt()`, `TSAdaptDestroy()`, `TSAdaptType`, `TSAdaptGetType()`
1257eef6055SBarry Smith @*/
126*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetType(TSAdapt adapt, TSAdaptType type)
127*d71ae5a4SJacob Faibussowitsch {
128ef749922SLisandro Dalcin   PetscBool match;
1295f80ce2aSJacob Faibussowitsch   PetscErrorCode (*r)(TSAdapt);
13084df9cb4SJed Brown 
13184df9cb4SJed Brown   PetscFunctionBegin;
1324782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
133b92453a8SLisandro Dalcin   PetscValidCharPointer(type, 2);
1349566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)adapt, type, &match));
135ef749922SLisandro Dalcin   if (match) PetscFunctionReturn(0);
1369566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(TSAdaptList, type, &r));
1373c633725SBarry Smith   PetscCheck(r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSAdapt type \"%s\" given", type);
138dbbe0bcdSBarry Smith   PetscTryTypeMethod(adapt, destroy);
1399566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(adapt->ops, sizeof(struct _TSAdaptOps)));
1409566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)adapt, type));
1419566063dSJacob Faibussowitsch   PetscCall((*r)(adapt));
14284df9cb4SJed Brown   PetscFunctionReturn(0);
14384df9cb4SJed Brown }
14484df9cb4SJed Brown 
145d0288e4fSLisandro Dalcin /*@C
146d0288e4fSLisandro Dalcin   TSAdaptGetType - gets the TS adapter method type (as a string).
147d0288e4fSLisandro Dalcin 
148d0288e4fSLisandro Dalcin   Not Collective
149d0288e4fSLisandro Dalcin 
150d0288e4fSLisandro Dalcin   Input Parameter:
151d0288e4fSLisandro Dalcin . adapt - The TS adapter, most likely obtained with TSGetAdapt()
152d0288e4fSLisandro Dalcin 
153d0288e4fSLisandro Dalcin   Output Parameter:
154d0288e4fSLisandro Dalcin . type - The name of TS adapter method
155d0288e4fSLisandro Dalcin 
156d0288e4fSLisandro Dalcin   Level: intermediate
157d0288e4fSLisandro Dalcin 
158db781477SPatrick Sanan .seealso `TSAdaptSetType()`
159d0288e4fSLisandro Dalcin @*/
160*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetType(TSAdapt adapt, TSAdaptType *type)
161*d71ae5a4SJacob Faibussowitsch {
162d0288e4fSLisandro Dalcin   PetscFunctionBegin;
163d0288e4fSLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
164d0288e4fSLisandro Dalcin   PetscValidPointer(type, 2);
165d0288e4fSLisandro Dalcin   *type = ((PetscObject)adapt)->type_name;
166d0288e4fSLisandro Dalcin   PetscFunctionReturn(0);
167d0288e4fSLisandro Dalcin }
168d0288e4fSLisandro Dalcin 
169*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt, const char prefix[])
170*d71ae5a4SJacob Faibussowitsch {
17184df9cb4SJed Brown   PetscFunctionBegin;
1724782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
1739566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)adapt, prefix));
17484df9cb4SJed Brown   PetscFunctionReturn(0);
17584df9cb4SJed Brown }
17684df9cb4SJed Brown 
177ad6bc421SBarry Smith /*@C
178ad6bc421SBarry Smith   TSAdaptLoad - Loads a TSAdapt that has been stored in binary  with TSAdaptView().
179ad6bc421SBarry Smith 
180ad6bc421SBarry Smith   Collective on PetscViewer
181ad6bc421SBarry Smith 
182ad6bc421SBarry Smith   Input Parameters:
183ad6bc421SBarry Smith + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or
184ad6bc421SBarry Smith            some related function before a call to TSAdaptLoad().
185ad6bc421SBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
186ad6bc421SBarry Smith            HDF5 file viewer, obtained from PetscViewerHDF5Open()
187ad6bc421SBarry Smith 
188ad6bc421SBarry Smith    Level: intermediate
189ad6bc421SBarry Smith 
190ad6bc421SBarry Smith   Notes:
191ad6bc421SBarry Smith    The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored.
192ad6bc421SBarry Smith 
193ad6bc421SBarry Smith   Notes for advanced users:
194ad6bc421SBarry Smith   Most users should not need to know the details of the binary storage
195ad6bc421SBarry Smith   format, since TSAdaptLoad() and TSAdaptView() completely hide these details.
196ad6bc421SBarry Smith   But for anyone who's interested, the standard binary matrix storage
197ad6bc421SBarry Smith   format is
198ad6bc421SBarry Smith .vb
199ad6bc421SBarry Smith      has not yet been determined
200ad6bc421SBarry Smith .ve
201ad6bc421SBarry Smith 
202db781477SPatrick Sanan .seealso: `PetscViewerBinaryOpen()`, `TSAdaptView()`, `MatLoad()`, `VecLoad()`
203ad6bc421SBarry Smith @*/
204*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptLoad(TSAdapt adapt, PetscViewer viewer)
205*d71ae5a4SJacob Faibussowitsch {
206ad6bc421SBarry Smith   PetscBool isbinary;
207ad6bc421SBarry Smith   char      type[256];
208ad6bc421SBarry Smith 
209ad6bc421SBarry Smith   PetscFunctionBegin;
2104782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
211ad6bc421SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
2129566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
2133c633725SBarry Smith   PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()");
214ad6bc421SBarry Smith 
2159566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR));
2169566063dSJacob Faibussowitsch   PetscCall(TSAdaptSetType(adapt, type));
217dbbe0bcdSBarry Smith   PetscTryTypeMethod(adapt, load, viewer);
218ad6bc421SBarry Smith   PetscFunctionReturn(0);
219ad6bc421SBarry Smith }
220ad6bc421SBarry Smith 
221*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptView(TSAdapt adapt, PetscViewer viewer)
222*d71ae5a4SJacob Faibussowitsch {
2231c167fc2SEmil Constantinescu   PetscBool iascii, isbinary, isnone, isglee;
22484df9cb4SJed Brown 
22584df9cb4SJed Brown   PetscFunctionBegin;
2264782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
2279566063dSJacob Faibussowitsch   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt), &viewer));
2284782b174SLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
2294782b174SLisandro Dalcin   PetscCheckSameComm(adapt, 1, viewer, 2);
2309566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
2319566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
23284df9cb4SJed Brown   if (iascii) {
2339566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)adapt, viewer));
2349566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTNONE, &isnone));
2359566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)adapt, TSADAPTGLEE, &isglee));
2361917a363SLisandro Dalcin     if (!isnone) {
2379566063dSJacob Faibussowitsch       if (adapt->always_accept) PetscCall(PetscViewerASCIIPrintf(viewer, "  always accepting steps\n"));
2389566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  safety factor %g\n", (double)adapt->safety));
2399566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  extra safety factor after step rejection %g\n", (double)adapt->reject_safety));
2409566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  clip fastest increase %g\n", (double)adapt->clip[1]));
2419566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  clip fastest decrease %g\n", (double)adapt->clip[0]));
2429566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum allowed timestep %g\n", (double)adapt->dt_max));
2439566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  minimum allowed timestep %g\n", (double)adapt->dt_min));
2449566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  maximum solution absolute value to be ignored %g\n", (double)adapt->ignore_max));
2451c167fc2SEmil Constantinescu     }
2461c167fc2SEmil Constantinescu     if (isglee) {
2471c167fc2SEmil Constantinescu       if (adapt->glee_use_local) {
2489566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  GLEE uses local error control\n"));
2491c167fc2SEmil Constantinescu       } else {
2509566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  GLEE uses global error control\n"));
2511c167fc2SEmil Constantinescu       }
2521917a363SLisandro Dalcin     }
2539566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
254dbbe0bcdSBarry Smith     PetscTryTypeMethod(adapt, view, viewer);
2559566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
256ad6bc421SBarry Smith   } else if (isbinary) {
257ad6bc421SBarry Smith     char type[256];
258ad6bc421SBarry Smith 
259ad6bc421SBarry Smith     /* need to save FILE_CLASS_ID for adapt class */
2609566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(type, ((PetscObject)adapt)->type_name, 256));
2619566063dSJacob Faibussowitsch     PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR));
262dbbe0bcdSBarry Smith   } else PetscTryTypeMethod(adapt, view, viewer);
26384df9cb4SJed Brown   PetscFunctionReturn(0);
26484df9cb4SJed Brown }
26584df9cb4SJed Brown 
266099af0a3SLisandro Dalcin /*@
267099af0a3SLisandro Dalcin    TSAdaptReset - Resets a TSAdapt context.
268099af0a3SLisandro Dalcin 
269099af0a3SLisandro Dalcin    Collective on TS
270099af0a3SLisandro Dalcin 
271099af0a3SLisandro Dalcin    Input Parameter:
272099af0a3SLisandro Dalcin .  adapt - the TSAdapt context obtained from TSAdaptCreate()
273099af0a3SLisandro Dalcin 
274099af0a3SLisandro Dalcin    Level: developer
275099af0a3SLisandro Dalcin 
276db781477SPatrick Sanan .seealso: `TSAdaptCreate()`, `TSAdaptDestroy()`
277099af0a3SLisandro Dalcin @*/
278*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptReset(TSAdapt adapt)
279*d71ae5a4SJacob Faibussowitsch {
280099af0a3SLisandro Dalcin   PetscFunctionBegin;
281099af0a3SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
282dbbe0bcdSBarry Smith   PetscTryTypeMethod(adapt, reset);
283099af0a3SLisandro Dalcin   PetscFunctionReturn(0);
284099af0a3SLisandro Dalcin }
285099af0a3SLisandro Dalcin 
286*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptDestroy(TSAdapt *adapt)
287*d71ae5a4SJacob Faibussowitsch {
28884df9cb4SJed Brown   PetscFunctionBegin;
28984df9cb4SJed Brown   if (!*adapt) PetscFunctionReturn(0);
29084df9cb4SJed Brown   PetscValidHeaderSpecific(*adapt, TSADAPT_CLASSID, 1);
2919371c9d4SSatish Balay   if (--((PetscObject)(*adapt))->refct > 0) {
2929371c9d4SSatish Balay     *adapt = NULL;
2939371c9d4SSatish Balay     PetscFunctionReturn(0);
2949371c9d4SSatish Balay   }
295099af0a3SLisandro Dalcin 
2969566063dSJacob Faibussowitsch   PetscCall(TSAdaptReset(*adapt));
297099af0a3SLisandro Dalcin 
298dbbe0bcdSBarry Smith   PetscTryTypeMethod((*adapt), destroy);
2999566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*adapt)->monitor));
3009566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(adapt));
30184df9cb4SJed Brown   PetscFunctionReturn(0);
30284df9cb4SJed Brown }
30384df9cb4SJed Brown 
3041c3436cfSJed Brown /*@
3051c3436cfSJed Brown    TSAdaptSetMonitor - Monitor the choices made by the adaptive controller
3061c3436cfSJed Brown 
3071c3436cfSJed Brown    Collective on TSAdapt
3081c3436cfSJed Brown 
3094165533cSJose E. Roman    Input Parameters:
3101c3436cfSJed Brown +  adapt - adaptive controller context
3111c3436cfSJed Brown -  flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable
3121c3436cfSJed Brown 
313bf997491SLisandro Dalcin    Options Database Keys:
314ec18b7bcSLisandro Dalcin .  -ts_adapt_monitor - to turn on monitoring
315bf997491SLisandro Dalcin 
3161c3436cfSJed Brown    Level: intermediate
3171c3436cfSJed Brown 
318db781477SPatrick Sanan .seealso: `TSAdaptChoose()`
3191c3436cfSJed Brown @*/
320*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt, PetscBool flg)
321*d71ae5a4SJacob Faibussowitsch {
3221c3436cfSJed Brown   PetscFunctionBegin;
3234782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
3244782b174SLisandro Dalcin   PetscValidLogicalCollectiveBool(adapt, flg, 2);
3251c3436cfSJed Brown   if (flg) {
3269566063dSJacob Faibussowitsch     if (!adapt->monitor) PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt), "stdout", &adapt->monitor));
3271c3436cfSJed Brown   } else {
3289566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&adapt->monitor));
3291c3436cfSJed Brown   }
3301c3436cfSJed Brown   PetscFunctionReturn(0);
3311c3436cfSJed Brown }
3321c3436cfSJed Brown 
3330873808bSJed Brown /*@C
334bf997491SLisandro Dalcin    TSAdaptSetCheckStage - Set a callback to check convergence for a stage
3350873808bSJed Brown 
3360873808bSJed Brown    Logically collective on TSAdapt
3370873808bSJed Brown 
3384165533cSJose E. Roman    Input Parameters:
3390873808bSJed Brown +  adapt - adaptive controller context
3400873808bSJed Brown -  func - stage check function
3410873808bSJed Brown 
3420873808bSJed Brown    Arguments of func:
3430873808bSJed Brown $  PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept)
3440873808bSJed Brown 
3450873808bSJed Brown +  adapt - adaptive controller context
3460873808bSJed Brown .  ts - time stepping context
3470873808bSJed Brown -  accept - pending choice of whether to accept, can be modified by this routine
3480873808bSJed Brown 
3490873808bSJed Brown    Level: advanced
3500873808bSJed Brown 
351db781477SPatrick Sanan .seealso: `TSAdaptChoose()`
3520873808bSJed Brown @*/
353*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt, PetscErrorCode (*func)(TSAdapt, TS, PetscReal, Vec, PetscBool *))
354*d71ae5a4SJacob Faibussowitsch {
3550873808bSJed Brown   PetscFunctionBegin;
3560873808bSJed Brown   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
35768ae941aSLisandro Dalcin   adapt->checkstage = func;
3580873808bSJed Brown   PetscFunctionReturn(0);
3590873808bSJed Brown }
3600873808bSJed Brown 
3611c3436cfSJed Brown /*@
362bf997491SLisandro Dalcin    TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of
363bf997491SLisandro Dalcin    any error or stability condition not meeting the prescribed goal.
364bf997491SLisandro Dalcin 
3651917a363SLisandro Dalcin    Logically collective on TSAdapt
366bf997491SLisandro Dalcin 
3674165533cSJose E. Roman    Input Parameters:
368bf997491SLisandro Dalcin +  adapt - time step adaptivity context, usually gotten with TSGetAdapt()
369bf997491SLisandro Dalcin -  flag - whether to always accept steps
370bf997491SLisandro Dalcin 
371bf997491SLisandro Dalcin    Options Database Keys:
372ec18b7bcSLisandro Dalcin .  -ts_adapt_always_accept - to always accept steps
373bf997491SLisandro Dalcin 
374bf997491SLisandro Dalcin    Level: intermediate
375bf997491SLisandro Dalcin 
376db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptChoose()`
377bf997491SLisandro Dalcin @*/
378*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt, PetscBool flag)
379*d71ae5a4SJacob Faibussowitsch {
380bf997491SLisandro Dalcin   PetscFunctionBegin;
381bf997491SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
382bf997491SLisandro Dalcin   PetscValidLogicalCollectiveBool(adapt, flag, 2);
383bf997491SLisandro Dalcin   adapt->always_accept = flag;
384bf997491SLisandro Dalcin   PetscFunctionReturn(0);
385bf997491SLisandro Dalcin }
386bf997491SLisandro Dalcin 
387bf997491SLisandro Dalcin /*@
3881917a363SLisandro Dalcin    TSAdaptSetSafety - Set safety factors
3891c3436cfSJed Brown 
3901917a363SLisandro Dalcin    Logically collective on TSAdapt
3911917a363SLisandro Dalcin 
3924165533cSJose E. Roman    Input Parameters:
3931917a363SLisandro Dalcin +  adapt - adaptive controller context
3941917a363SLisandro Dalcin .  safety - safety factor relative to target error/stability goal
395ec18b7bcSLisandro Dalcin -  reject_safety - extra safety factor to apply if the last step was rejected
3961917a363SLisandro Dalcin 
3971917a363SLisandro Dalcin    Options Database Keys:
398ec18b7bcSLisandro Dalcin +  -ts_adapt_safety <safety> - to set safety factor
399ec18b7bcSLisandro Dalcin -  -ts_adapt_reject_safety <reject_safety> - to set reject safety factor
4001917a363SLisandro Dalcin 
4011917a363SLisandro Dalcin    Level: intermediate
4021917a363SLisandro Dalcin 
403db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptGetSafety()`, `TSAdaptChoose()`
4041917a363SLisandro Dalcin @*/
405*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetSafety(TSAdapt adapt, PetscReal safety, PetscReal reject_safety)
406*d71ae5a4SJacob Faibussowitsch {
4071917a363SLisandro Dalcin   PetscFunctionBegin;
4081917a363SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
4091917a363SLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt, safety, 2);
4101917a363SLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt, reject_safety, 3);
4113c633725SBarry Smith   PetscCheck(safety == PETSC_DEFAULT || safety >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Safety factor %g must be non negative", (double)safety);
4123c633725SBarry Smith   PetscCheck(safety == PETSC_DEFAULT || safety <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Safety factor %g must be less than one", (double)safety);
4133c633725SBarry Smith   PetscCheck(reject_safety == PETSC_DEFAULT || reject_safety >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Reject safety factor %g must be non negative", (double)reject_safety);
4143c633725SBarry Smith   PetscCheck(reject_safety == PETSC_DEFAULT || reject_safety <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Reject safety factor %g must be less than one", (double)reject_safety);
4151917a363SLisandro Dalcin   if (safety != PETSC_DEFAULT) adapt->safety = safety;
4161917a363SLisandro Dalcin   if (reject_safety != PETSC_DEFAULT) adapt->reject_safety = reject_safety;
4171917a363SLisandro Dalcin   PetscFunctionReturn(0);
4181917a363SLisandro Dalcin }
4191917a363SLisandro Dalcin 
4201917a363SLisandro Dalcin /*@
4211917a363SLisandro Dalcin    TSAdaptGetSafety - Get safety factors
4221917a363SLisandro Dalcin 
4231917a363SLisandro Dalcin    Not Collective
4241917a363SLisandro Dalcin 
4254165533cSJose E. Roman    Input Parameter:
4261917a363SLisandro Dalcin .  adapt - adaptive controller context
4271917a363SLisandro Dalcin 
4284165533cSJose E. Roman    Output Parameters:
4291917a363SLisandro Dalcin .  safety - safety factor relative to target error/stability goal
4301917a363SLisandro Dalcin +  reject_safety - extra safety factor to apply if the last step was rejected
4311917a363SLisandro Dalcin 
4321917a363SLisandro Dalcin    Level: intermediate
4331917a363SLisandro Dalcin 
434db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptSetSafety()`, `TSAdaptChoose()`
4351917a363SLisandro Dalcin @*/
436*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetSafety(TSAdapt adapt, PetscReal *safety, PetscReal *reject_safety)
437*d71ae5a4SJacob Faibussowitsch {
4381917a363SLisandro Dalcin   PetscFunctionBegin;
4391917a363SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
4401917a363SLisandro Dalcin   if (safety) PetscValidRealPointer(safety, 2);
4411917a363SLisandro Dalcin   if (reject_safety) PetscValidRealPointer(reject_safety, 3);
4421917a363SLisandro Dalcin   if (safety) *safety = adapt->safety;
4431917a363SLisandro Dalcin   if (reject_safety) *reject_safety = adapt->reject_safety;
4441917a363SLisandro Dalcin   PetscFunctionReturn(0);
4451917a363SLisandro Dalcin }
4461917a363SLisandro Dalcin 
4471917a363SLisandro Dalcin /*@
44876cddca1SEmil Constantinescu    TSAdaptSetMaxIgnore - Set error estimation threshold. Solution components below this threshold value will not be considered when computing error norms for time step adaptivity (in absolute value). A negative value (default) of the threshold leads to considering all solution components.
44976cddca1SEmil Constantinescu 
45076cddca1SEmil Constantinescu    Logically collective on TSAdapt
45176cddca1SEmil Constantinescu 
4524165533cSJose E. Roman    Input Parameters:
45376cddca1SEmil Constantinescu +  adapt - adaptive controller context
45476cddca1SEmil Constantinescu -  max_ignore - threshold for solution components that are ignored during error estimation
45576cddca1SEmil Constantinescu 
45676cddca1SEmil Constantinescu    Options Database Keys:
45776cddca1SEmil Constantinescu .  -ts_adapt_max_ignore <max_ignore> - to set the threshold
45876cddca1SEmil Constantinescu 
45976cddca1SEmil Constantinescu    Level: intermediate
46076cddca1SEmil Constantinescu 
461db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptGetMaxIgnore()`, `TSAdaptChoose()`
46276cddca1SEmil Constantinescu @*/
463*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetMaxIgnore(TSAdapt adapt, PetscReal max_ignore)
464*d71ae5a4SJacob Faibussowitsch {
46576cddca1SEmil Constantinescu   PetscFunctionBegin;
46676cddca1SEmil Constantinescu   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
46776cddca1SEmil Constantinescu   PetscValidLogicalCollectiveReal(adapt, max_ignore, 2);
46876cddca1SEmil Constantinescu   adapt->ignore_max = max_ignore;
46976cddca1SEmil Constantinescu   PetscFunctionReturn(0);
47076cddca1SEmil Constantinescu }
47176cddca1SEmil Constantinescu 
47276cddca1SEmil Constantinescu /*@
47376cddca1SEmil Constantinescu    TSAdaptGetMaxIgnore - Get error estimation threshold. Solution components below this threshold value will not be considered when computing error norms for time step adaptivity (in absolute value).
47476cddca1SEmil Constantinescu 
47576cddca1SEmil Constantinescu    Not Collective
47676cddca1SEmil Constantinescu 
4774165533cSJose E. Roman    Input Parameter:
47876cddca1SEmil Constantinescu .  adapt - adaptive controller context
47976cddca1SEmil Constantinescu 
4804165533cSJose E. Roman    Output Parameter:
48176cddca1SEmil Constantinescu .  max_ignore - threshold for solution components that are ignored during error estimation
48276cddca1SEmil Constantinescu 
48376cddca1SEmil Constantinescu    Level: intermediate
48476cddca1SEmil Constantinescu 
485db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptSetMaxIgnore()`, `TSAdaptChoose()`
48676cddca1SEmil Constantinescu @*/
487*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetMaxIgnore(TSAdapt adapt, PetscReal *max_ignore)
488*d71ae5a4SJacob Faibussowitsch {
48976cddca1SEmil Constantinescu   PetscFunctionBegin;
49076cddca1SEmil Constantinescu   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
49176cddca1SEmil Constantinescu   PetscValidRealPointer(max_ignore, 2);
49276cddca1SEmil Constantinescu   *max_ignore = adapt->ignore_max;
49376cddca1SEmil Constantinescu   PetscFunctionReturn(0);
49476cddca1SEmil Constantinescu }
49576cddca1SEmil Constantinescu 
49676cddca1SEmil Constantinescu /*@
4971917a363SLisandro Dalcin    TSAdaptSetClip - Sets the admissible decrease/increase factor in step size
4981917a363SLisandro Dalcin 
4991917a363SLisandro Dalcin    Logically collective on TSAdapt
5001917a363SLisandro Dalcin 
5014165533cSJose E. Roman    Input Parameters:
5021917a363SLisandro Dalcin +  adapt - adaptive controller context
5031917a363SLisandro Dalcin .  low - admissible decrease factor
504ec18b7bcSLisandro Dalcin -  high - admissible increase factor
5051917a363SLisandro Dalcin 
5061917a363SLisandro Dalcin    Options Database Keys:
507ec18b7bcSLisandro Dalcin .  -ts_adapt_clip <low>,<high> - to set admissible time step decrease and increase factors
5081917a363SLisandro Dalcin 
5091917a363SLisandro Dalcin    Level: intermediate
5101917a363SLisandro Dalcin 
511db781477SPatrick Sanan .seealso: `TSAdaptChoose()`, `TSAdaptGetClip()`, `TSAdaptSetScaleSolveFailed()`
5121917a363SLisandro Dalcin @*/
513*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetClip(TSAdapt adapt, PetscReal low, PetscReal high)
514*d71ae5a4SJacob Faibussowitsch {
5151917a363SLisandro Dalcin   PetscFunctionBegin;
5161917a363SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
5171917a363SLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt, low, 2);
5181917a363SLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt, high, 3);
5193c633725SBarry Smith   PetscCheck(low == PETSC_DEFAULT || low >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Decrease factor %g must be non negative", (double)low);
5203c633725SBarry Smith   PetscCheck(low == PETSC_DEFAULT || low <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Decrease factor %g must be less than one", (double)low);
5213c633725SBarry Smith   PetscCheck(high == PETSC_DEFAULT || high >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Increase factor %g must be greater than one", (double)high);
5221917a363SLisandro Dalcin   if (low != PETSC_DEFAULT) adapt->clip[0] = low;
5231917a363SLisandro Dalcin   if (high != PETSC_DEFAULT) adapt->clip[1] = high;
5241917a363SLisandro Dalcin   PetscFunctionReturn(0);
5251917a363SLisandro Dalcin }
5261917a363SLisandro Dalcin 
5271917a363SLisandro Dalcin /*@
5281917a363SLisandro Dalcin    TSAdaptGetClip - Gets the admissible decrease/increase factor in step size
5291917a363SLisandro Dalcin 
5301917a363SLisandro Dalcin    Not Collective
5311917a363SLisandro Dalcin 
5324165533cSJose E. Roman    Input Parameter:
5331917a363SLisandro Dalcin .  adapt - adaptive controller context
5341917a363SLisandro Dalcin 
5354165533cSJose E. Roman    Output Parameters:
5361917a363SLisandro Dalcin +  low - optional, admissible decrease factor
5371917a363SLisandro Dalcin -  high - optional, admissible increase factor
5381917a363SLisandro Dalcin 
5391917a363SLisandro Dalcin    Level: intermediate
5401917a363SLisandro Dalcin 
541db781477SPatrick Sanan .seealso: `TSAdaptChoose()`, `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()`
5421917a363SLisandro Dalcin @*/
543*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetClip(TSAdapt adapt, PetscReal *low, PetscReal *high)
544*d71ae5a4SJacob Faibussowitsch {
5451917a363SLisandro Dalcin   PetscFunctionBegin;
5461917a363SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
5471917a363SLisandro Dalcin   if (low) PetscValidRealPointer(low, 2);
5481917a363SLisandro Dalcin   if (high) PetscValidRealPointer(high, 3);
5491917a363SLisandro Dalcin   if (low) *low = adapt->clip[0];
5501917a363SLisandro Dalcin   if (high) *high = adapt->clip[1];
5511917a363SLisandro Dalcin   PetscFunctionReturn(0);
5521917a363SLisandro Dalcin }
5531917a363SLisandro Dalcin 
5541917a363SLisandro Dalcin /*@
55562c23b28SMark    TSAdaptSetScaleSolveFailed - Scale step by this factor if solve fails
55662c23b28SMark 
55762c23b28SMark    Logically collective on TSAdapt
55862c23b28SMark 
5594165533cSJose E. Roman    Input Parameters:
56062c23b28SMark +  adapt - adaptive controller context
561ee300463SSatish Balay -  scale - scale
56262c23b28SMark 
56362c23b28SMark    Options Database Keys:
56462c23b28SMark .  -ts_adapt_scale_solve_failed <scale> - to set scale step by this factor if solve fails
56562c23b28SMark 
56662c23b28SMark    Level: intermediate
56762c23b28SMark 
568db781477SPatrick Sanan .seealso: `TSAdaptChoose()`, `TSAdaptGetScaleSolveFailed()`, `TSAdaptGetClip()`
56962c23b28SMark @*/
570*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetScaleSolveFailed(TSAdapt adapt, PetscReal scale)
571*d71ae5a4SJacob Faibussowitsch {
57262c23b28SMark   PetscFunctionBegin;
57362c23b28SMark   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
57462c23b28SMark   PetscValidLogicalCollectiveReal(adapt, scale, 2);
5753c633725SBarry Smith   PetscCheck(scale == PETSC_DEFAULT || scale > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Scale factor %g must be positive", (double)scale);
5763c633725SBarry Smith   PetscCheck(scale == PETSC_DEFAULT || scale <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Scale factor %g must be less than one", (double)scale);
57762c23b28SMark   if (scale != PETSC_DEFAULT) adapt->scale_solve_failed = scale;
57862c23b28SMark   PetscFunctionReturn(0);
57962c23b28SMark }
58062c23b28SMark 
58162c23b28SMark /*@
58262c23b28SMark    TSAdaptGetScaleSolveFailed - Gets the admissible decrease/increase factor in step size
58362c23b28SMark 
58462c23b28SMark    Not Collective
58562c23b28SMark 
5864165533cSJose E. Roman    Input Parameter:
58762c23b28SMark .  adapt - adaptive controller context
58862c23b28SMark 
5894165533cSJose E. Roman    Output Parameter:
590ee300463SSatish Balay .  scale - scale factor
59162c23b28SMark 
59262c23b28SMark    Level: intermediate
59362c23b28SMark 
594db781477SPatrick Sanan .seealso: `TSAdaptChoose()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetClip()`
59562c23b28SMark @*/
596*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetScaleSolveFailed(TSAdapt adapt, PetscReal *scale)
597*d71ae5a4SJacob Faibussowitsch {
59862c23b28SMark   PetscFunctionBegin;
59962c23b28SMark   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
60062c23b28SMark   if (scale) PetscValidRealPointer(scale, 2);
60162c23b28SMark   if (scale) *scale = adapt->scale_solve_failed;
60262c23b28SMark   PetscFunctionReturn(0);
60362c23b28SMark }
60462c23b28SMark 
60562c23b28SMark /*@
6061917a363SLisandro Dalcin    TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the controller
6071917a363SLisandro Dalcin 
6081917a363SLisandro Dalcin    Logically collective on TSAdapt
6091c3436cfSJed Brown 
6104165533cSJose E. Roman    Input Parameters:
611552698daSJed Brown +  adapt - time step adaptivity context, usually gotten with TSGetAdapt()
6121c3436cfSJed Brown .  hmin - minimum time step
6131c3436cfSJed Brown -  hmax - maximum time step
6141c3436cfSJed Brown 
6151c3436cfSJed Brown    Options Database Keys:
616ec18b7bcSLisandro Dalcin +  -ts_adapt_dt_min <min> - to set minimum time step
617ec18b7bcSLisandro Dalcin -  -ts_adapt_dt_max <max> - to set maximum time step
6181c3436cfSJed Brown 
6191c3436cfSJed Brown    Level: intermediate
6201c3436cfSJed Brown 
621db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptGetStepLimits()`, `TSAdaptChoose()`
6221c3436cfSJed Brown @*/
623*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt, PetscReal hmin, PetscReal hmax)
624*d71ae5a4SJacob Faibussowitsch {
6251c3436cfSJed Brown   PetscFunctionBegin;
6264782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
627b1f5bccaSLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt, hmin, 2);
628b1f5bccaSLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt, hmax, 3);
6293c633725SBarry Smith   PetscCheck(hmin == PETSC_DEFAULT || hmin >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum time step %g must be non negative", (double)hmin);
6303c633725SBarry Smith   PetscCheck(hmax == PETSC_DEFAULT || hmax >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Minimum time step %g must be non negative", (double)hmax);
631b1f5bccaSLisandro Dalcin   if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin;
632b1f5bccaSLisandro Dalcin   if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax;
633b1f5bccaSLisandro Dalcin   hmin = adapt->dt_min;
634b1f5bccaSLisandro Dalcin   hmax = adapt->dt_max;
6353c633725SBarry 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);
6361917a363SLisandro Dalcin   PetscFunctionReturn(0);
6371917a363SLisandro Dalcin }
6381917a363SLisandro Dalcin 
6391917a363SLisandro Dalcin /*@
6401917a363SLisandro Dalcin    TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the controller
6411917a363SLisandro Dalcin 
6421917a363SLisandro Dalcin    Not Collective
6431917a363SLisandro Dalcin 
6444165533cSJose E. Roman    Input Parameter:
6451917a363SLisandro Dalcin .  adapt - time step adaptivity context, usually gotten with TSGetAdapt()
6461917a363SLisandro Dalcin 
6474165533cSJose E. Roman    Output Parameters:
6481917a363SLisandro Dalcin +  hmin - minimum time step
6491917a363SLisandro Dalcin -  hmax - maximum time step
6501917a363SLisandro Dalcin 
6511917a363SLisandro Dalcin    Level: intermediate
6521917a363SLisandro Dalcin 
653db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptSetStepLimits()`, `TSAdaptChoose()`
6541917a363SLisandro Dalcin @*/
655*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt, PetscReal *hmin, PetscReal *hmax)
656*d71ae5a4SJacob Faibussowitsch {
6571917a363SLisandro Dalcin   PetscFunctionBegin;
6581917a363SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
6591917a363SLisandro Dalcin   if (hmin) PetscValidRealPointer(hmin, 2);
6601917a363SLisandro Dalcin   if (hmax) PetscValidRealPointer(hmax, 3);
6611917a363SLisandro Dalcin   if (hmin) *hmin = adapt->dt_min;
6621917a363SLisandro Dalcin   if (hmax) *hmax = adapt->dt_max;
6631c3436cfSJed Brown   PetscFunctionReturn(0);
6641c3436cfSJed Brown }
6651c3436cfSJed Brown 
666e55864a3SBarry Smith /*
66784df9cb4SJed Brown    TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options.
66884df9cb4SJed Brown 
66984df9cb4SJed Brown    Collective on TSAdapt
67084df9cb4SJed Brown 
67184df9cb4SJed Brown    Input Parameter:
67284df9cb4SJed Brown .  adapt - the TSAdapt context
67384df9cb4SJed Brown 
67484df9cb4SJed Brown    Options Database Keys:
6751917a363SLisandro Dalcin +  -ts_adapt_type <type> - algorithm to use for adaptivity
676bf997491SLisandro Dalcin .  -ts_adapt_always_accept - always accept steps regardless of error/stability goals
6771917a363SLisandro Dalcin .  -ts_adapt_safety <safety> - safety factor relative to target error/stability goal
6781917a363SLisandro Dalcin .  -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected
6791917a363SLisandro Dalcin .  -ts_adapt_clip <low,high> - admissible time step decrease and increase factors
68023de1d84SBarry Smith .  -ts_adapt_dt_min <min> - minimum timestep to use
68123de1d84SBarry Smith .  -ts_adapt_dt_max <max> - maximum timestep to use
68223de1d84SBarry Smith .  -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails
683de50f1caSBarry Smith .  -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates
684de50f1caSBarry 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
68584df9cb4SJed Brown 
68684df9cb4SJed Brown    Level: advanced
68784df9cb4SJed Brown 
68884df9cb4SJed Brown    Notes:
68984df9cb4SJed Brown    This function is automatically called by TSSetFromOptions()
69084df9cb4SJed Brown 
691db781477SPatrick Sanan .seealso: `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptSetAlwaysAccept()`, `TSAdaptSetSafety()`,
692db781477SPatrick Sanan           `TSAdaptSetClip()`, `TSAdaptSetScaleSolveFailed()`, `TSAdaptSetStepLimits()`, `TSAdaptSetMonitor()`
693e55864a3SBarry Smith */
694*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetFromOptions(TSAdapt adapt, PetscOptionItems *PetscOptionsObject)
695*d71ae5a4SJacob Faibussowitsch {
69684df9cb4SJed Brown   char      type[256] = TSADAPTBASIC;
69762c23b28SMark   PetscReal safety, reject_safety, clip[2], scale, hmin, hmax;
6981c3436cfSJed Brown   PetscBool set, flg;
6991917a363SLisandro Dalcin   PetscInt  two;
70084df9cb4SJed Brown 
70184df9cb4SJed Brown   PetscFunctionBegin;
702dbbe0bcdSBarry Smith   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
70384df9cb4SJed Brown   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
7041566a47fSLisandro Dalcin    * function can only be called from inside TSSetFromOptions()  */
705d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "TS Adaptivity options");
7069566063dSJacob 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));
70748a46eb9SPierre Jolivet   if (flg || !((PetscObject)adapt)->type_name) PetscCall(TSAdaptSetType(adapt, type));
7081917a363SLisandro Dalcin 
7099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_adapt_always_accept", "Always accept the step", "TSAdaptSetAlwaysAccept", adapt->always_accept, &flg, &set));
7109566063dSJacob Faibussowitsch   if (set) PetscCall(TSAdaptSetAlwaysAccept(adapt, flg));
7111917a363SLisandro Dalcin 
7129371c9d4SSatish Balay   safety        = adapt->safety;
7139371c9d4SSatish Balay   reject_safety = adapt->reject_safety;
7149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_adapt_safety", "Safety factor relative to target error/stability goal", "TSAdaptSetSafety", safety, &safety, &set));
7159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_adapt_reject_safety", "Extra safety factor to apply if the last step was rejected", "TSAdaptSetSafety", reject_safety, &reject_safety, &flg));
7169566063dSJacob Faibussowitsch   if (set || flg) PetscCall(TSAdaptSetSafety(adapt, safety, reject_safety));
7171917a363SLisandro Dalcin 
7189371c9d4SSatish Balay   two     = 2;
7199371c9d4SSatish Balay   clip[0] = adapt->clip[0];
7209371c9d4SSatish Balay   clip[1] = adapt->clip[1];
7219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRealArray("-ts_adapt_clip", "Admissible decrease/increase factor in step size", "TSAdaptSetClip", clip, &two, &set));
7223c633725SBarry Smith   PetscCheck(!set || (two == 2), PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Must give exactly two values to -ts_adapt_clip");
7239566063dSJacob Faibussowitsch   if (set) PetscCall(TSAdaptSetClip(adapt, clip[0], clip[1]));
7241917a363SLisandro Dalcin 
7259371c9d4SSatish Balay   hmin = adapt->dt_min;
7269371c9d4SSatish Balay   hmax = adapt->dt_max;
7279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_adapt_dt_min", "Minimum time step considered", "TSAdaptSetStepLimits", hmin, &hmin, &set));
7289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_adapt_dt_max", "Maximum time step considered", "TSAdaptSetStepLimits", hmax, &hmax, &flg));
7299566063dSJacob Faibussowitsch   if (set || flg) PetscCall(TSAdaptSetStepLimits(adapt, hmin, hmax));
7301917a363SLisandro Dalcin 
7319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_adapt_max_ignore", "Adaptor ignores (absolute) solution values smaller than this value", "", adapt->ignore_max, &adapt->ignore_max, &set));
7329566063dSJacob 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));
733d580f011SEmil Constantinescu 
7349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ts_adapt_scale_solve_failed", "Scale step by this factor if solve fails", "TSAdaptSetScaleSolveFailed", adapt->scale_solve_failed, &scale, &set));
7359566063dSJacob Faibussowitsch   if (set) PetscCall(TSAdaptSetScaleSolveFailed(adapt, scale));
7361917a363SLisandro Dalcin 
7379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-ts_adapt_wnormtype", "Type of norm computed for error estimation", "", NormTypes, (PetscEnum)adapt->wnormtype, (PetscEnum *)&adapt->wnormtype, NULL));
7383c633725SBarry Smith   PetscCheck(adapt->wnormtype == NORM_2 || adapt->wnormtype == NORM_INFINITY, PetscObjectComm((PetscObject)adapt), PETSC_ERR_SUP, "Only 2-norm and infinite norm supported");
7391917a363SLisandro Dalcin 
7409566063dSJacob 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));
741de50f1caSBarry Smith 
7429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ts_adapt_monitor", "Print choices made by adaptive controller", "TSAdaptSetMonitor", adapt->monitor ? PETSC_TRUE : PETSC_FALSE, &flg, &set));
7439566063dSJacob Faibussowitsch   if (set) PetscCall(TSAdaptSetMonitor(adapt, flg));
7441917a363SLisandro Dalcin 
745dbbe0bcdSBarry Smith   PetscTryTypeMethod(adapt, setfromoptions, PetscOptionsObject);
746d0609cedSBarry Smith   PetscOptionsHeadEnd();
74784df9cb4SJed Brown   PetscFunctionReturn(0);
74884df9cb4SJed Brown }
74984df9cb4SJed Brown 
75084df9cb4SJed Brown /*@
75184df9cb4SJed Brown    TSAdaptCandidatesClear - clear any previously set candidate schemes
75284df9cb4SJed Brown 
7531917a363SLisandro Dalcin    Logically collective on TSAdapt
75484df9cb4SJed Brown 
7554165533cSJose E. Roman    Input Parameter:
75684df9cb4SJed Brown .  adapt - adaptive controller
75784df9cb4SJed Brown 
75884df9cb4SJed Brown    Level: developer
75984df9cb4SJed Brown 
760db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptCreate()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()`
76184df9cb4SJed Brown @*/
762*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
763*d71ae5a4SJacob Faibussowitsch {
76484df9cb4SJed Brown   PetscFunctionBegin;
7654782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
7669566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(&adapt->candidates, sizeof(adapt->candidates)));
76784df9cb4SJed Brown   PetscFunctionReturn(0);
76884df9cb4SJed Brown }
76984df9cb4SJed Brown 
77084df9cb4SJed Brown /*@C
77184df9cb4SJed Brown    TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from
77284df9cb4SJed Brown 
7731917a363SLisandro Dalcin    Logically collective on TSAdapt
77484df9cb4SJed Brown 
7754165533cSJose E. Roman    Input Parameters:
776552698daSJed Brown +  adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate()
77784df9cb4SJed Brown .  name - name of the candidate scheme to add
77884df9cb4SJed Brown .  order - order of the candidate scheme
77984df9cb4SJed Brown .  stageorder - stage order of the candidate scheme
7808d59e960SJed Brown .  ccfl - stability coefficient relative to explicit Euler, used for CFL constraints
78184df9cb4SJed Brown .  cost - relative measure of the amount of work required for the candidate scheme
78284df9cb4SJed Brown -  inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme
78384df9cb4SJed Brown 
78484df9cb4SJed Brown    Note:
78584df9cb4SJed Brown    This routine is not available in Fortran.
78684df9cb4SJed Brown 
78784df9cb4SJed Brown    Level: developer
78884df9cb4SJed Brown 
789db781477SPatrick Sanan .seealso: `TSAdaptCandidatesClear()`, `TSAdaptChoose()`
79084df9cb4SJed Brown @*/
791*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt, const char name[], PetscInt order, PetscInt stageorder, PetscReal ccfl, PetscReal cost, PetscBool inuse)
792*d71ae5a4SJacob Faibussowitsch {
79384df9cb4SJed Brown   PetscInt c;
79484df9cb4SJed Brown 
79584df9cb4SJed Brown   PetscFunctionBegin;
79684df9cb4SJed Brown   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
79763a3b9bcSJacob Faibussowitsch   PetscCheck(order >= 1, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Classical order %" PetscInt_FMT " must be a positive integer", order);
79884df9cb4SJed Brown   if (inuse) {
7993c633725SBarry Smith     PetscCheck(!adapt->candidates.inuse_set, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONGSTATE, "Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()");
80084df9cb4SJed Brown     adapt->candidates.inuse_set = PETSC_TRUE;
80184df9cb4SJed Brown   }
8021c3436cfSJed Brown   /* first slot if this is the current scheme, otherwise the next available slot */
8031c3436cfSJed Brown   c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n;
804bbd56ea5SKarl Rupp 
80584df9cb4SJed Brown   adapt->candidates.name[c]       = name;
80684df9cb4SJed Brown   adapt->candidates.order[c]      = order;
80784df9cb4SJed Brown   adapt->candidates.stageorder[c] = stageorder;
8088d59e960SJed Brown   adapt->candidates.ccfl[c]       = ccfl;
80984df9cb4SJed Brown   adapt->candidates.cost[c]       = cost;
81084df9cb4SJed Brown   adapt->candidates.n++;
81184df9cb4SJed Brown   PetscFunctionReturn(0);
81284df9cb4SJed Brown }
81384df9cb4SJed Brown 
8148d59e960SJed Brown /*@C
8158d59e960SJed Brown    TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost
8168d59e960SJed Brown 
8178d59e960SJed Brown    Not Collective
8188d59e960SJed Brown 
8194165533cSJose E. Roman    Input Parameter:
8208d59e960SJed Brown .  adapt - time step adaptivity context
8218d59e960SJed Brown 
8224165533cSJose E. Roman    Output Parameters:
8238d59e960SJed Brown +  n - number of candidate schemes, always at least 1
8248d59e960SJed Brown .  order - the order of each candidate scheme
8258d59e960SJed Brown .  stageorder - the stage order of each candidate scheme
8268d59e960SJed Brown .  ccfl - the CFL coefficient of each scheme
8278d59e960SJed Brown -  cost - the relative cost of each scheme
8288d59e960SJed Brown 
8298d59e960SJed Brown    Level: developer
8308d59e960SJed Brown 
8318d59e960SJed Brown    Note:
8328d59e960SJed Brown    The current scheme is always returned in the first slot
8338d59e960SJed Brown 
834db781477SPatrick Sanan .seealso: `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()`, `TSAdaptChoose()`
8358d59e960SJed Brown @*/
836*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt, PetscInt *n, const PetscInt **order, const PetscInt **stageorder, const PetscReal **ccfl, const PetscReal **cost)
837*d71ae5a4SJacob Faibussowitsch {
8388d59e960SJed Brown   PetscFunctionBegin;
8398d59e960SJed Brown   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
8408d59e960SJed Brown   if (n) *n = adapt->candidates.n;
8418d59e960SJed Brown   if (order) *order = adapt->candidates.order;
8428d59e960SJed Brown   if (stageorder) *stageorder = adapt->candidates.stageorder;
8438d59e960SJed Brown   if (ccfl) *ccfl = adapt->candidates.ccfl;
8448d59e960SJed Brown   if (cost) *cost = adapt->candidates.cost;
8458d59e960SJed Brown   PetscFunctionReturn(0);
8468d59e960SJed Brown }
8478d59e960SJed Brown 
84884df9cb4SJed Brown /*@C
84984df9cb4SJed Brown    TSAdaptChoose - choose which method and step size to use for the next step
85084df9cb4SJed Brown 
8511917a363SLisandro Dalcin    Collective on TSAdapt
85284df9cb4SJed Brown 
8534165533cSJose E. Roman    Input Parameters:
85484df9cb4SJed Brown +  adapt - adaptive contoller
85597bb3fdcSJose E. Roman .  ts - time stepper
85684df9cb4SJed Brown -  h - current step size
85784df9cb4SJed Brown 
8584165533cSJose E. Roman    Output Parameters:
8591566a47fSLisandro Dalcin +  next_sc - optional, scheme to use for the next step
86084df9cb4SJed Brown .  next_h - step size to use for the next step
86184df9cb4SJed Brown -  accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size
86284df9cb4SJed Brown 
8631c3436cfSJed Brown    Note:
8641c3436cfSJed Brown    The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is
8651c3436cfSJed Brown    being retried after an initial rejection.
8661c3436cfSJed Brown 
86784df9cb4SJed Brown    Level: developer
86884df9cb4SJed Brown 
869db781477SPatrick Sanan .seealso: `TSAdapt`, `TSAdaptCandidatesClear()`, `TSAdaptCandidateAdd()`
87084df9cb4SJed Brown @*/
871*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptChoose(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept)
872*d71ae5a4SJacob Faibussowitsch {
8731566a47fSLisandro Dalcin   PetscInt  ncandidates = adapt->candidates.n;
8741566a47fSLisandro Dalcin   PetscInt  scheme      = 0;
8750b99f514SJed Brown   PetscReal wlte        = -1.0;
8765e4ed32fSEmil Constantinescu   PetscReal wltea       = -1.0;
8775e4ed32fSEmil Constantinescu   PetscReal wlter       = -1.0;
87884df9cb4SJed Brown 
87984df9cb4SJed Brown   PetscFunctionBegin;
88084df9cb4SJed Brown   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
88184df9cb4SJed Brown   PetscValidHeaderSpecific(ts, TS_CLASSID, 2);
8821566a47fSLisandro Dalcin   if (next_sc) PetscValidIntPointer(next_sc, 4);
883dadcf809SJacob Faibussowitsch   PetscValidRealPointer(next_h, 5);
884064a246eSJacob Faibussowitsch   PetscValidBoolPointer(accept, 6);
8851566a47fSLisandro Dalcin   if (next_sc) *next_sc = 0;
8861566a47fSLisandro Dalcin 
8871566a47fSLisandro Dalcin   /* Do not mess with adaptivity while handling events*/
8881566a47fSLisandro Dalcin   if (ts->event && ts->event->status != TSEVENT_NONE) {
8891566a47fSLisandro Dalcin     *next_h = h;
8901566a47fSLisandro Dalcin     *accept = PETSC_TRUE;
8911566a47fSLisandro Dalcin     PetscFunctionReturn(0);
8921566a47fSLisandro Dalcin   }
8931566a47fSLisandro Dalcin 
894dbbe0bcdSBarry Smith   PetscUseTypeMethod(adapt, choose, ts, h, &scheme, next_h, accept, &wlte, &wltea, &wlter);
89563a3b9bcSJacob 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);
8963c633725SBarry Smith   PetscCheck(*next_h >= 0, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Computed step size %g must be positive", (double)*next_h);
8971566a47fSLisandro Dalcin   if (next_sc) *next_sc = scheme;
8981566a47fSLisandro Dalcin 
89949354f04SShri Abhyankar   if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
90036b54a69SLisandro Dalcin     /* Increase/reduce step size if end time of next step is close to or overshoots max time */
90136b54a69SLisandro Dalcin     PetscReal t = ts->ptime + ts->time_step, h = *next_h;
9024a658b32SHong Zhang     PetscReal tend = t + h, tmax, hmax;
903fe7350e0SStefano Zampini     PetscReal a    = (PetscReal)(1.0 + adapt->matchstepfac[0]);
904fe7350e0SStefano Zampini     PetscReal b    = adapt->matchstepfac[1];
9054a658b32SHong Zhang 
9064a658b32SHong Zhang     if (ts->tspan) {
907e1db57b0SHong Zhang       if (PetscIsCloseAtTol(t, ts->tspan->span_times[ts->tspan->spanctr], ts->tspan->reltol * h + ts->tspan->abstol, 0)) /* hit a span time point */
9084a658b32SHong Zhang         if (ts->tspan->spanctr + 1 < ts->tspan->num_span_times) tmax = ts->tspan->span_times[ts->tspan->spanctr + 1];
9094a658b32SHong Zhang         else tmax = ts->max_time; /* hit the last span time point */
9104a658b32SHong Zhang       else tmax = ts->tspan->span_times[ts->tspan->spanctr];
9114a658b32SHong Zhang     } else tmax = ts->max_time;
9124a658b32SHong Zhang     hmax = tmax - t;
91336b54a69SLisandro Dalcin     if (t < tmax && tend > tmax) *next_h = hmax;
914fe7350e0SStefano Zampini     if (t < tmax && tend < tmax && h * b > hmax) *next_h = hmax / 2;
91536b54a69SLisandro Dalcin     if (t < tmax && tend < tmax && h * a > hmax) *next_h = hmax;
916e1db57b0SHong Zhang     /* if step size is changed to match a span time point */
917e1db57b0SHong Zhang     if (ts->tspan && h != *next_h && !adapt->dt_span_cached) adapt->dt_span_cached = h;
918e1db57b0SHong Zhang     /* reset time step after a span time point */
919e1db57b0SHong 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)) {
920e1db57b0SHong Zhang       *next_h               = adapt->dt_span_cached;
921e1db57b0SHong Zhang       adapt->dt_span_cached = 0;
92249354f04SShri Abhyankar     }
923e1db57b0SHong Zhang   }
9241c3436cfSJed Brown   if (adapt->monitor) {
9251566a47fSLisandro Dalcin     const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : "";
9269566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
9270b99f514SJed Brown     if (wlte < 0) {
9289371c9d4SSatish 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",
9299371c9d4SSatish Balay                                        (double)ts->ptime, (double)h, (double)*next_h));
9300b99f514SJed Brown     } else {
9319371c9d4SSatish 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",
9329371c9d4SSatish Balay                                        (double)ts->ptime, (double)h, (double)*next_h, (double)wlte, (double)wltea, (double)wlter));
9330b99f514SJed Brown     }
9349566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
9351c3436cfSJed Brown   }
93684df9cb4SJed Brown   PetscFunctionReturn(0);
93784df9cb4SJed Brown }
93884df9cb4SJed Brown 
93997335746SJed Brown /*@
940de50f1caSBarry Smith    TSAdaptSetTimeStepIncreaseDelay - The number of timesteps to wait after a decrease in the timestep due to failed solver
941de50f1caSBarry Smith                                      before increasing the time step.
942de50f1caSBarry Smith 
943de50f1caSBarry Smith    Logicially Collective on TSAdapt
944de50f1caSBarry Smith 
9454165533cSJose E. Roman    Input Parameters:
946de50f1caSBarry Smith +  adapt - adaptive controller context
947de50f1caSBarry Smith -  cnt - the number of timesteps
948de50f1caSBarry Smith 
949de50f1caSBarry Smith    Options Database Key:
950de50f1caSBarry Smith .  -ts_adapt_time_step_increase_delay cnt - number of steps to delay the increase
951de50f1caSBarry Smith 
952de50f1caSBarry Smith    Notes: This is to prevent an adaptor from bouncing back and forth between two nearby timesteps. The default is 0.
953de50f1caSBarry Smith           The successful use of this option is problem dependent
954de50f1caSBarry Smith 
955de50f1caSBarry Smith    Developer Note: there is no theory to support this option
956de50f1caSBarry Smith 
957de50f1caSBarry Smith    Level: advanced
958de50f1caSBarry Smith 
959de50f1caSBarry Smith .seealso:
960de50f1caSBarry Smith @*/
961*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptSetTimeStepIncreaseDelay(TSAdapt adapt, PetscInt cnt)
962*d71ae5a4SJacob Faibussowitsch {
963de50f1caSBarry Smith   PetscFunctionBegin;
964de50f1caSBarry Smith   adapt->timestepjustdecreased_delay = cnt;
965de50f1caSBarry Smith   PetscFunctionReturn(0);
966de50f1caSBarry Smith }
967de50f1caSBarry Smith 
968de50f1caSBarry Smith /*@
9696bc98fa9SBarry 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)
97097335746SJed Brown 
9711917a363SLisandro Dalcin    Collective on TSAdapt
97297335746SJed Brown 
9734165533cSJose E. Roman    Input Parameters:
97497335746SJed Brown +  adapt - adaptive controller context
975b295832fSPierre Barbier de Reuille .  ts - time stepper
976b295832fSPierre Barbier de Reuille .  t - Current simulation time
977b295832fSPierre Barbier de Reuille -  Y - Current solution vector
97897335746SJed Brown 
9794165533cSJose E. Roman    Output Parameter:
98097335746SJed Brown .  accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject
98197335746SJed Brown 
98297335746SJed Brown    Level: developer
98397335746SJed Brown 
98497335746SJed Brown .seealso:
98597335746SJed Brown @*/
986*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCheckStage(TSAdapt adapt, TS ts, PetscReal t, Vec Y, PetscBool *accept)
987*d71ae5a4SJacob Faibussowitsch {
9881566a47fSLisandro Dalcin   SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING;
98997335746SJed Brown 
99097335746SJed Brown   PetscFunctionBegin;
9914782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1);
9924782b174SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 2);
993064a246eSJacob Faibussowitsch   PetscValidBoolPointer(accept, 5);
9941566a47fSLisandro Dalcin 
9959566063dSJacob Faibussowitsch   if (ts->snes) PetscCall(SNESGetConvergedReason(ts->snes, &snesreason));
99697335746SJed Brown   if (snesreason < 0) {
99797335746SJed Brown     *accept = PETSC_FALSE;
9986de24e2aSJed Brown     if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) {
99997335746SJed Brown       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
100063a3b9bcSJacob 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));
100197335746SJed Brown       if (adapt->monitor) {
10029566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
10039371c9d4SSatish 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,
10049371c9d4SSatish Balay                                          (double)ts->ptime, (double)ts->time_step, ts->num_snes_failures));
10059566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
100697335746SJed Brown       }
1007cb9d8021SPierre Barbier de Reuille     }
1008cb9d8021SPierre Barbier de Reuille   } else {
10091566a47fSLisandro Dalcin     *accept = PETSC_TRUE;
10109566063dSJacob Faibussowitsch     PetscCall(TSFunctionDomainError(ts, t, Y, accept));
1011cb9d8021SPierre Barbier de Reuille     if (*accept && adapt->checkstage) {
10129566063dSJacob Faibussowitsch       PetscCall((*adapt->checkstage)(adapt, ts, t, Y, accept));
10136bc98fa9SBarry Smith       if (!*accept) {
101463a3b9bcSJacob Faibussowitsch         PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", solution rejected by user function provided by TSSetFunctionDomainError()\n", ts->steps));
10156bc98fa9SBarry Smith         if (adapt->monitor) {
10169566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
101763a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(adapt->monitor, "    TSAdapt %s step %3" PetscInt_FMT " stage rejected by user function provided by TSSetFunctionDomainError()\n", ((PetscObject)adapt)->type_name, ts->steps));
10189566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
10196bc98fa9SBarry Smith         }
10206bc98fa9SBarry Smith       }
1021cb9d8021SPierre Barbier de Reuille     }
1022cb9d8021SPierre Barbier de Reuille   }
1023cb9d8021SPierre Barbier de Reuille 
10241566a47fSLisandro Dalcin   if (!(*accept) && !ts->reason) {
10251566a47fSLisandro Dalcin     PetscReal dt, new_dt;
10269566063dSJacob Faibussowitsch     PetscCall(TSGetTimeStep(ts, &dt));
1027cb9d8021SPierre Barbier de Reuille     new_dt = dt * adapt->scale_solve_failed;
10289566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts, new_dt));
1029de50f1caSBarry Smith     adapt->timestepjustdecreased += adapt->timestepjustdecreased_delay;
103097335746SJed Brown     if (adapt->monitor) {
10319566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIAddTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
103263a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(adapt->monitor, "    TSAdapt %s step %3" PetscInt_FMT " stage rejected (%s) t=%-11g+%10.3e retrying with dt=%-10.3e\n", ((PetscObject)adapt)->type_name, ts->steps, SNESConvergedReasons[snesreason], (double)ts->ptime, (double)dt, (double)new_dt));
10339566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISubtractTab(adapt->monitor, ((PetscObject)adapt)->tablevel));
103497335746SJed Brown     }
103597335746SJed Brown   }
103697335746SJed Brown   PetscFunctionReturn(0);
103797335746SJed Brown }
103897335746SJed Brown 
103984df9cb4SJed Brown /*@
104084df9cb4SJed Brown   TSAdaptCreate - create an adaptive controller context for time stepping
104184df9cb4SJed Brown 
1042d083f849SBarry Smith   Collective
104384df9cb4SJed Brown 
104484df9cb4SJed Brown   Input Parameter:
104584df9cb4SJed Brown . comm - The communicator
104684df9cb4SJed Brown 
104784df9cb4SJed Brown   Output Parameter:
104884df9cb4SJed Brown . adapt - new TSAdapt object
104984df9cb4SJed Brown 
105084df9cb4SJed Brown   Level: developer
105184df9cb4SJed Brown 
105284df9cb4SJed Brown   Notes:
105384df9cb4SJed Brown   TSAdapt creation is handled by TS, so users should not need to call this function.
105484df9cb4SJed Brown 
1055db781477SPatrick Sanan .seealso: `TSGetAdapt()`, `TSAdaptSetType()`, `TSAdaptDestroy()`
105684df9cb4SJed Brown @*/
1057*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptCreate(MPI_Comm comm, TSAdapt *inadapt)
1058*d71ae5a4SJacob Faibussowitsch {
105984df9cb4SJed Brown   TSAdapt adapt;
106084df9cb4SJed Brown 
106184df9cb4SJed Brown   PetscFunctionBegin;
1062064a246eSJacob Faibussowitsch   PetscValidPointer(inadapt, 2);
10633b3bcf4cSLisandro Dalcin   *inadapt = NULL;
10649566063dSJacob Faibussowitsch   PetscCall(TSAdaptInitializePackage());
10653b3bcf4cSLisandro Dalcin 
10669566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(adapt, TSADAPT_CLASSID, "TSAdapt", "Time stepping adaptivity", "TS", comm, TSAdaptDestroy, TSAdaptView));
10671c3436cfSJed Brown 
1068bf997491SLisandro Dalcin   adapt->always_accept      = PETSC_FALSE;
10691917a363SLisandro Dalcin   adapt->safety             = 0.9;
10701917a363SLisandro Dalcin   adapt->reject_safety      = 0.5;
10711917a363SLisandro Dalcin   adapt->clip[0]            = 0.1;
10721917a363SLisandro Dalcin   adapt->clip[1]            = 10.;
10731c3436cfSJed Brown   adapt->dt_min             = 1e-20;
10741917a363SLisandro Dalcin   adapt->dt_max             = 1e+20;
10751c167fc2SEmil Constantinescu   adapt->ignore_max         = -1.0;
1076d580f011SEmil Constantinescu   adapt->glee_use_local     = PETSC_TRUE;
107797335746SJed Brown   adapt->scale_solve_failed = 0.25;
1078fe7350e0SStefano Zampini   /* these two safety factors are not public, and they are used only in the TS_EXACTFINALTIME_MATCHSTEP case
1079fe7350e0SStefano Zampini      to prevent from situations were unreasonably small time steps are taken in order to match the final time */
1080fe7350e0SStefano Zampini   adapt->matchstepfac[0]             = 0.01; /* allow 1% step size increase in the last step */
1081fe7350e0SStefano Zampini   adapt->matchstepfac[1]             = 2.0;  /* halve last step if it is greater than what remains divided this factor */
10827619abb3SShri   adapt->wnormtype                   = NORM_2;
1083de50f1caSBarry Smith   adapt->timestepjustdecreased_delay = 0;
10841917a363SLisandro Dalcin 
108584df9cb4SJed Brown   *inadapt = adapt;
108684df9cb4SJed Brown   PetscFunctionReturn(0);
108784df9cb4SJed Brown }
1088