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