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