Actual source code: tsadapt.c

petsc-3.7.7 2017-09-25
Report Typos and Errors
  2: #include <petsc/private/tsimpl.h> /*I  "petscts.h" I*/

  4: PetscClassId TSADAPT_CLASSID;

  6: static PetscFunctionList TSAdaptList;
  7: static PetscBool         TSAdaptPackageInitialized;
  8: static PetscBool         TSAdaptRegisterAllCalled;

 10: PETSC_EXTERN PetscErrorCode TSAdaptCreate_None(TSAdapt);
 11: PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt);
 12: PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt);

 16: /*@C
 17:    TSAdaptRegister -  adds a TSAdapt implementation

 19:    Not Collective

 21:    Input Parameters:
 22: +  name_scheme - name of user-defined adaptivity scheme
 23: -  routine_create - routine to create method context

 25:    Notes:
 26:    TSAdaptRegister() may be called multiple times to add several user-defined families.

 28:    Sample usage:
 29: .vb
 30:    TSAdaptRegister("my_scheme",MySchemeCreate);
 31: .ve

 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

 38:    Level: advanced

 40: .keywords: TSAdapt, register

 42: .seealso: TSAdaptRegisterAll()
 43: @*/
 44: PetscErrorCode  TSAdaptRegister(const char sname[],PetscErrorCode (*function)(TSAdapt))
 45: {

 49:   PetscFunctionListAdd(&TSAdaptList,sname,function);
 50:   return(0);
 51: }

 55: /*@C
 56:   TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt

 58:   Not Collective

 60:   Level: advanced

 62: .keywords: TSAdapt, register, all

 64: .seealso: TSAdaptRegisterDestroy()
 65: @*/
 66: PetscErrorCode  TSAdaptRegisterAll(void)
 67: {

 71:   if (TSAdaptRegisterAllCalled) return(0);
 72:   TSAdaptRegisterAllCalled = PETSC_TRUE;
 73:   TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None);
 74:   TSAdaptRegister(TSADAPTBASIC,TSAdaptCreate_Basic);
 75:   TSAdaptRegister(TSADAPTCFL,  TSAdaptCreate_CFL);
 76:   return(0);
 77: }

 81: /*@C
 82:   TSAdaptFinalizePackage - This function destroys everything in the TS package. It is
 83:   called from PetscFinalize().

 85:   Level: developer

 87: .keywords: Petsc, destroy, package
 88: .seealso: PetscFinalize()
 89: @*/
 90: PetscErrorCode  TSAdaptFinalizePackage(void)
 91: {

 95:   PetscFunctionListDestroy(&TSAdaptList);
 96:   TSAdaptPackageInitialized = PETSC_FALSE;
 97:   TSAdaptRegisterAllCalled  = PETSC_FALSE;
 98:   return(0);
 99: }

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.

108:   Level: developer

110: .keywords: TSAdapt, initialize, package
111: .seealso: PetscInitialize()
112: @*/
113: PetscErrorCode  TSAdaptInitializePackage(void)
114: {

118:   if (TSAdaptPackageInitialized) return(0);
119:   TSAdaptPackageInitialized = PETSC_TRUE;
120:   PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);
121:   TSAdaptRegisterAll();
122:   PetscRegisterFinalize(TSAdaptFinalizePackage);
123:   return(0);
124: }

128: /*@C
129:   TSAdaptSetType - sets the approach used for the error adapter, currently there is only TSADAPTBASIC and TSADAPTNONE

131:   Logicially Collective on TSAdapt

133:   Input Parameter:
134: + adapt - the TS error adapter, most likely obtained with TSGetAdapt()
135: - type - either  TSADAPTBASIC or TSADAPTNONE

137:   Options Database:
138: .  -ts_adapt_type basic or none - to setting the adapter type

140:   Level: intermediate

142: .keywords: TSAdapt, create

144: .seealso: TSGetAdapt(), TSAdaptDestroy(), TSAdaptType, TSAdaptGetType()
145: @*/
146: PetscErrorCode  TSAdaptSetType(TSAdapt adapt,TSAdaptType type)
147: {
148:   PetscBool      match;
149:   PetscErrorCode ierr,(*r)(TSAdapt);

153:   PetscObjectTypeCompare((PetscObject)adapt,type,&match);
154:   if (match) return(0);
155:   PetscFunctionListFind(TSAdaptList,type,&r);
156:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type);
157:   if (adapt->ops->destroy) {(*adapt->ops->destroy)(adapt);}
158:   PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));
159:   PetscObjectChangeTypeName((PetscObject)adapt,type);
160:   (*r)(adapt);
161:   return(0);
162: }

166: PetscErrorCode  TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[])
167: {

172:   PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);
173:   return(0);
174: }

178: /*@C
179:   TSAdaptLoad - Loads a TSAdapt that has been stored in binary  with TSAdaptView().

181:   Collective on PetscViewer

183:   Input Parameters:
184: + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or
185:            some related function before a call to TSAdaptLoad().
186: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
187:            HDF5 file viewer, obtained from PetscViewerHDF5Open()

189:    Level: intermediate

191:   Notes:
192:    The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored.

194:   Notes for advanced users:
195:   Most users should not need to know the details of the binary storage
196:   format, since TSAdaptLoad() and TSAdaptView() completely hide these details.
197:   But for anyone who's interested, the standard binary matrix storage
198:   format is
199: .vb
200:      has not yet been determined
201: .ve

203: .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad()
204: @*/
205: PetscErrorCode  TSAdaptLoad(TSAdapt adapt,PetscViewer viewer)
206: {
208:   PetscBool      isbinary;
209:   char           type[256];

214:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
215:   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");

217:   PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);
218:   TSAdaptSetType(adapt,type);
219:   if (adapt->ops->load) {
220:     (*adapt->ops->load)(adapt,viewer);
221:   }
222:   return(0);
223: }

227: PetscErrorCode  TSAdaptView(TSAdapt adapt,PetscViewer viewer)
228: {
230:   PetscBool      iascii,isbinary;

234:   if (!viewer) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt),&viewer);}
237:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
238:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
239:   if (iascii) {
240:     PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);
241:     PetscViewerASCIIPrintf(viewer,"  number of candidates %D\n",adapt->candidates.n);
242:     if (adapt->ops->view) {
243:       PetscViewerASCIIPushTab(viewer);
244:       (*adapt->ops->view)(adapt,viewer);
245:       PetscViewerASCIIPopTab(viewer);
246:     }
247:   } else if (isbinary) {
248:     char type[256];

250:     /* need to save FILE_CLASS_ID for adapt class */
251:     PetscStrncpy(type,((PetscObject)adapt)->type_name,256);
252:     PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);
253:   } else if (adapt->ops->view) {
254:     (*adapt->ops->view)(adapt,viewer);
255:   }
256:   return(0);
257: }

261: /*@
262:    TSAdaptReset - Resets a TSAdapt context.

264:    Collective on TS

266:    Input Parameter:
267: .  adapt - the TSAdapt context obtained from TSAdaptCreate()

269:    Level: developer

271: .seealso: TSAdaptCreate(), TSAdaptDestroy()
272: @*/
273: PetscErrorCode  TSAdaptReset(TSAdapt adapt)
274: {

279:   if (adapt->ops->reset) {(*adapt->ops->reset)(adapt);}
280:   return(0);
281: }

285: PetscErrorCode  TSAdaptDestroy(TSAdapt *adapt)
286: {

290:   if (!*adapt) return(0);
292:   if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; return(0);}

294:   TSAdaptReset(*adapt);

296:   if ((*adapt)->ops->destroy) {(*(*adapt)->ops->destroy)(*adapt);}
297:   PetscViewerDestroy(&(*adapt)->monitor);
298:   PetscHeaderDestroy(adapt);
299:   return(0);
300: }

304: /*@
305:    TSAdaptSetMonitor - Monitor the choices made by the adaptive controller

307:    Collective on TSAdapt

309:    Input Arguments:
310: +  adapt - adaptive controller context
311: -  flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable

313:    Level: intermediate

315: .seealso: TSAdaptChoose()
316: @*/
317: PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg)
318: {

324:   if (flg) {
325:     if (!adapt->monitor) {PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);}
326:   } else {
327:     PetscViewerDestroy(&adapt->monitor);
328:   }
329:   return(0);
330: }

334: /*@C
335:    TSAdaptSetCheckStage - set a callback to check convergence for a stage

337:    Logically collective on TSAdapt

339:    Input Arguments:
340: +  adapt - adaptive controller context
341: -  func - stage check function

343:    Arguments of func:
344: $  PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept)

346: +  adapt - adaptive controller context
347: .  ts - time stepping context
348: -  accept - pending choice of whether to accept, can be modified by this routine

350:    Level: advanced

352: .seealso: TSAdaptChoose()
353: @*/
354: PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscReal,Vec,PetscBool*))
355: {

359:   adapt->checkstage = func;
360:   return(0);
361: }

365: /*@
366:    TSAdaptSetStepLimits - Set minimum and maximum step sizes to be considered by the controller

368:    Logically Collective

370:    Input Arguments:
371: +  adapt - time step adaptivity context, usually gotten with TSGetAdapt()
372: .  hmin - minimum time step
373: -  hmax - maximum time step

375:    Options Database Keys:
376: +  -ts_adapt_dt_min - minimum time step
377: -  -ts_adapt_dt_max - maximum time step

379:    Level: intermediate

381: .seealso: TSAdapt
382: @*/
383: PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax)
384: {

390:   if (hmin != PETSC_DEFAULT && hmin < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmin);
391:   if (hmax != PETSC_DEFAULT && hmax < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmax);
392:   if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin;
393:   if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax;
394: #if defined(PETSC_USE_DEBUG)
395:   hmin = adapt->dt_min;
396:   hmax = adapt->dt_max;
397:   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);
398: #endif
399:   return(0);
400: }

404: /*
405:    TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options.

407:    Collective on TSAdapt

409:    Input Parameter:
410: .  adapt - the TSAdapt context

412:    Options Database Keys:
413: .  -ts_adapt_type <type> - basic

415:    Level: advanced

417:    Notes:
418:    This function is automatically called by TSSetFromOptions()

420: .keywords: TS, TSGetAdapt(), TSAdaptSetType()

422: .seealso: TSGetType()
423: */
424: PetscErrorCode  TSAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSAdapt adapt)
425: {
427:   char           type[256] = TSADAPTBASIC;
428:   PetscBool      set,flg;

432:   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
433:    * function can only be called from inside TSSetFromOptions()  */
434:   PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");
435:   PetscOptionsFList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList,
436:                           ((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);
437:   if (flg || !((PetscObject)adapt)->type_name) {
438:     TSAdaptSetType(adapt,type);
439:   }
440:   PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",adapt->dt_min,&adapt->dt_min,NULL);
441:   PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",adapt->dt_max,&adapt->dt_max,NULL);
442:   PetscOptionsReal("-ts_adapt_scale_solve_failed","Scale step by this factor if solve fails","",adapt->scale_solve_failed,&adapt->scale_solve_failed,NULL);
443:   PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);
444:   PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);
445:   if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported");
446:   if (set) {TSAdaptSetMonitor(adapt,flg);}
447:   if (adapt->ops->setfromoptions) {(*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);}
448:   PetscOptionsTail();
449:   return(0);
450: }

454: /*@
455:    TSAdaptCandidatesClear - clear any previously set candidate schemes

457:    Logically Collective

459:    Input Argument:
460: .  adapt - adaptive controller

462:    Level: developer

464: .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose()
465: @*/
466: PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
467: {

472:   PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));
473:   return(0);
474: }

478: /*@C
479:    TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from

481:    Logically Collective

483:    Input Arguments:
484: +  adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate()
485: .  name - name of the candidate scheme to add
486: .  order - order of the candidate scheme
487: .  stageorder - stage order of the candidate scheme
488: .  ccfl - stability coefficient relative to explicit Euler, used for CFL constraints
489: .  cost - relative measure of the amount of work required for the candidate scheme
490: -  inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme

492:    Note:
493:    This routine is not available in Fortran.

495:    Level: developer

497: .seealso: TSAdaptCandidatesClear(), TSAdaptChoose()
498: @*/
499: PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse)
500: {
501:   PetscInt c;

505:   if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order);
506:   if (inuse) {
507:     if (adapt->candidates.inuse_set) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()");
508:     adapt->candidates.inuse_set = PETSC_TRUE;
509:   }
510:   /* first slot if this is the current scheme, otherwise the next available slot */
511:   c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n;

513:   adapt->candidates.name[c]       = name;
514:   adapt->candidates.order[c]      = order;
515:   adapt->candidates.stageorder[c] = stageorder;
516:   adapt->candidates.ccfl[c]       = ccfl;
517:   adapt->candidates.cost[c]       = cost;
518:   adapt->candidates.n++;
519:   return(0);
520: }

524: /*@C
525:    TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost

527:    Not Collective

529:    Input Arguments:
530: .  adapt - time step adaptivity context

532:    Output Arguments:
533: +  n - number of candidate schemes, always at least 1
534: .  order - the order of each candidate scheme
535: .  stageorder - the stage order of each candidate scheme
536: .  ccfl - the CFL coefficient of each scheme
537: -  cost - the relative cost of each scheme

539:    Level: developer

541:    Note:
542:    The current scheme is always returned in the first slot

544: .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose()
545: @*/
546: PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost)
547: {
550:   if (n) *n = adapt->candidates.n;
551:   if (order) *order = adapt->candidates.order;
552:   if (stageorder) *stageorder = adapt->candidates.stageorder;
553:   if (ccfl) *ccfl = adapt->candidates.ccfl;
554:   if (cost) *cost = adapt->candidates.cost;
555:   return(0);
556: }

560: /*@C
561:    TSAdaptChoose - choose which method and step size to use for the next step

563:    Logically Collective

565:    Input Arguments:
566: +  adapt - adaptive contoller
567: -  h - current step size

569:    Output Arguments:
570: +  next_sc - optional, scheme to use for the next step
571: .  next_h - step size to use for the next step
572: -  accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size

574:    Note:
575:    The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is
576:    being retried after an initial rejection.

578:    Level: developer

580: .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd()
581: @*/
582: PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept)
583: {
585:   PetscInt       ncandidates = adapt->candidates.n;
586:   PetscInt       scheme = 0;
587:   PetscReal      wlte = -1.0;

595:   if (next_sc) *next_sc = 0;

597:   /* Do not mess with adaptivity while handling events*/
598:   if (ts->event && ts->event->status != TSEVENT_NONE) {
599:     *next_h = h;
600:     *accept = PETSC_TRUE;
601:     return(0);
602:   }

604:   (*adapt->ops->choose)(adapt,ts,h,&scheme,next_h,accept,&wlte);
605:   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);
606:   if (*next_h < 0) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h);
607:   if (next_sc) *next_sc = scheme;

609:   if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
610:     /* Reduce time step if it overshoots max time */
611:     if (ts->ptime + ts->time_step + *next_h >= ts->max_time) {
612:       PetscReal next_dt = ts->max_time - (ts->ptime + ts->time_step);
613:       if (next_dt > PETSC_SMALL) *next_h = next_dt;
614:       else ts->reason = TS_CONVERGED_TIME;
615:     }
616:   }

618:   if (adapt->monitor) {
619:     const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : "";
620:     PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);
621:     if (wlte < 0) {
622:       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);
623:     } else {
624:       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);
625:     }
626:     PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);
627:   }
628:   return(0);
629: }

633: /*@
634:    TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails)

636:    Collective

638:    Input Arguments:
639: +  adapt - adaptive controller context
640: .  ts - time stepper
641: .  t - Current simulation time
642: -  Y - Current solution vector

644:    Output Arguments:
645: .  accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject

647:    Level: developer

649: .seealso:
650: @*/
651: PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool *accept)
652: {
653:   PetscErrorCode      ierr;
654:   SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING;


661:   if (ts->snes) {SNESGetConvergedReason(ts->snes,&snesreason);}
662:   if (snesreason < 0) {
663:     *accept = PETSC_FALSE;
664:     if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) {
665:       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
666:       PetscInfo2(ts,"Step=%D, nonlinear solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);
667:       if (adapt->monitor) {
668:         PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);
669:         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);
670:         PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);
671:       }
672:     }
673:   } else {
674:     *accept = PETSC_TRUE;
675:     TSFunctionDomainError(ts,t,Y,accept);
676:     if(*accept && adapt->checkstage) {
677:       (*adapt->checkstage)(adapt,ts,t,Y,accept);
678:     }
679:   }

681:   if(!(*accept) && !ts->reason) {
682:     PetscReal dt,new_dt;
683:     TSGetTimeStep(ts,&dt);
684:     new_dt = dt * adapt->scale_solve_failed;
685:     TSSetTimeStep(ts,new_dt);
686:     if (adapt->monitor) {
687:       PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);
688:       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);
689:       PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);
690:     }
691:   }
692:   return(0);
693: }

697: /*@
698:   TSAdaptCreate - create an adaptive controller context for time stepping

700:   Collective on MPI_Comm

702:   Input Parameter:
703: . comm - The communicator

705:   Output Parameter:
706: . adapt - new TSAdapt object

708:   Level: developer

710:   Notes:
711:   TSAdapt creation is handled by TS, so users should not need to call this function.

713: .keywords: TSAdapt, create
714: .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy()
715: @*/
716: PetscErrorCode  TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt)
717: {
719:   TSAdapt        adapt;

723:   *inadapt = NULL;
724:   TSAdaptInitializePackage();

726:   PetscHeaderCreate(adapt,TSADAPT_CLASSID,"TSAdapt","Time stepping adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);

728:   adapt->dt_min             = 1e-20;
729:   adapt->dt_max             = 1e50;
730:   adapt->scale_solve_failed = 0.25;
731:   adapt->wnormtype          = NORM_2;
732:   TSAdaptSetType(adapt,TSADAPTBASIC);

734:   *inadapt = adapt;
735:   return(0);
736: }