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