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