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