Actual source code: tsadapt.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  2: #include <petsc/private/tsimpl.h> /*I  "petscts.h" I*/

  4: static PetscFunctionList TSAdaptList;
  5: static PetscBool         TSAdaptPackageInitialized;
  6: static PetscBool         TSAdaptRegisterAllCalled;
  7: static PetscClassId      TSADAPT_CLASSID;

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

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

 18:    Not Collective

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

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

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

 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

 37:    Level: advanced

 39: .keywords: TSAdapt, register

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

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

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

 57:   Not Collective

 59:   Level: advanced

 61: .keywords: TSAdapt, register, all

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

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

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

 84:   Level: developer

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

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

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.

107:   Level: developer

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

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

127: PetscErrorCode  TSAdaptSetType(TSAdapt adapt,TSAdaptType type)
128: {
129:   PetscBool      match;
130:   PetscErrorCode ierr,(*r)(TSAdapt);

134:   PetscObjectTypeCompare((PetscObject)adapt,type,&match);
135:   if (match) return(0);
136:   PetscFunctionListFind(TSAdaptList,type,&r);
137:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type);
138:   if (adapt->ops->destroy) {(*adapt->ops->destroy)(adapt);}
139:   PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));
140:   PetscObjectChangeTypeName((PetscObject)adapt,type);
141:   (*r)(adapt);
142:   return(0);
143: }

147: PetscErrorCode  TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[])
148: {

153:   PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);
154:   return(0);
155: }

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

162:   Collective on PetscViewer

164:   Input Parameters:
165: + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or
166:            some related function before a call to TSAdaptLoad().
167: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
168:            HDF5 file viewer, obtained from PetscViewerHDF5Open()

170:    Level: intermediate

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

175:   Notes for advanced users:
176:   Most users should not need to know the details of the binary storage
177:   format, since TSAdaptLoad() and TSAdaptView() completely hide these details.
178:   But for anyone who's interested, the standard binary matrix storage
179:   format is
180: .vb
181:      has not yet been determined
182: .ve

184: .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad()
185: @*/
186: PetscErrorCode  TSAdaptLoad(TSAdapt adapt,PetscViewer viewer)
187: {
189:   PetscBool      isbinary;
190:   char           type[256];

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

198:   PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);
199:   TSAdaptSetType(adapt, type);
200:   if (adapt->ops->load) {
201:     (*adapt->ops->load)(adapt,viewer);
202:   }
203:   return(0);
204: }

208: PetscErrorCode  TSAdaptView(TSAdapt adapt,PetscViewer viewer)
209: {
211:   PetscBool      iascii,isbinary;

215:   if (!viewer) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt),&viewer);}
218:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
219:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
220:   if (iascii) {
221:     PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);
222:     PetscViewerASCIIPrintf(viewer,"  number of candidates %D\n",adapt->candidates.n);
223:     if (adapt->ops->view) {
224:       PetscViewerASCIIPushTab(viewer);
225:       (*adapt->ops->view)(adapt,viewer);
226:       PetscViewerASCIIPopTab(viewer);
227:     }
228:   } else if (isbinary) {
229:     char type[256];

231:     /* need to save FILE_CLASS_ID for adapt class */
232:     PetscStrncpy(type,((PetscObject)adapt)->type_name,256);
233:     PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);
234:   } else if (adapt->ops->view) {
235:     (*adapt->ops->view)(adapt,viewer);
236:   }
237:   return(0);
238: }

242: /*@
243:    TSAdaptReset - Resets a TSAdapt context.

245:    Collective on TS

247:    Input Parameter:
248: .  adapt - the TSAdapt context obtained from TSAdaptCreate()

250:    Level: developer

252: .seealso: TSAdaptCreate(), TSAdaptDestroy()
253: @*/
254: PetscErrorCode  TSAdaptReset(TSAdapt adapt)
255: {

260:   if (adapt->ops->reset) {(*adapt->ops->reset)(adapt);}
261:   return(0);
262: }

266: PetscErrorCode  TSAdaptDestroy(TSAdapt *adapt)
267: {

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

275:   TSAdaptReset(*adapt);

277:   if ((*adapt)->ops->destroy) {(*(*adapt)->ops->destroy)(*adapt);}
278:   PetscViewerDestroy(&(*adapt)->monitor);
279:   PetscHeaderDestroy(adapt);
280:   return(0);
281: }

285: /*@
286:    TSAdaptSetMonitor - Monitor the choices made by the adaptive controller

288:    Collective on TSAdapt

290:    Input Arguments:
291: +  adapt - adaptive controller context
292: -  flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable

294:    Level: intermediate

296: .seealso: TSAdaptChoose()
297: @*/
298: PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg)
299: {

305:   if (flg) {
306:     if (!adapt->monitor) {PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);}
307:   } else {
308:     PetscViewerDestroy(&adapt->monitor);
309:   }
310:   return(0);
311: }

315: /*@C
316:    TSAdaptSetCheckStage - set a callback to check convergence for a stage

318:    Logically collective on TSAdapt

320:    Input Arguments:
321: +  adapt - adaptive controller context
322: -  func - stage check function

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

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

331:    Level: advanced

333: .seealso: TSAdaptChoose()
334: @*/
335: PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscBool*))
336: {

340:   adapt->checkstage = func;
341:   return(0);
342: }

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

349:    Logically Collective

351:    Input Arguments:
352: +  adapt - time step adaptivity context, usually gotten with TSGetAdapt()
353: .  hmin - minimum time step
354: -  hmax - maximum time step

356:    Options Database Keys:
357: +  -ts_adapt_dt_min - minimum time step
358: -  -ts_adapt_dt_max - maximum time step

360:    Level: intermediate

362: .seealso: TSAdapt
363: @*/
364: PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax)
365: {

369:   if (hmin != PETSC_DECIDE) adapt->dt_min = hmin;
370:   if (hmax != PETSC_DECIDE) adapt->dt_max = hmax;
371:   return(0);
372: }

376: /*
377:    TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options.

379:    Collective on TSAdapt

381:    Input Parameter:
382: .  adapt - the TSAdapt context

384:    Options Database Keys:
385: .  -ts_adapt_type <type> - basic

387:    Level: advanced

389:    Notes:
390:    This function is automatically called by TSSetFromOptions()

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

394: .seealso: TSGetType()
395: */
396: PetscErrorCode  TSAdaptSetFromOptions(PetscOptions *PetscOptionsObject,TSAdapt adapt)
397: {
399:   char           type[256] = TSADAPTBASIC;
400:   PetscBool      set,flg;

404:   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
405:   * function can only be called from inside TSSetFromOptions_GL()  */
406:   PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");
407:   PetscOptionsFList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList,
408:                           ((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);
409:   if (flg || !((PetscObject)adapt)->type_name) {
410:     TSAdaptSetType(adapt,type);
411:   }
412:   PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",adapt->dt_min,&adapt->dt_min,NULL);
413:   PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",adapt->dt_max,&adapt->dt_max,NULL);
414:   PetscOptionsReal("-ts_adapt_scale_solve_failed","Scale step by this factor if solve fails","",adapt->scale_solve_failed,&adapt->scale_solve_failed,NULL);
415:   PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);
416:   PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);
417:   if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported");
418:   if (set) {TSAdaptSetMonitor(adapt,flg);}
419:   if (adapt->ops->setfromoptions) {(*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);}
420:   PetscOptionsTail();
421:   return(0);
422: }

426: /*@
427:    TSAdaptCandidatesClear - clear any previously set candidate schemes

429:    Logically Collective

431:    Input Argument:
432: .  adapt - adaptive controller

434:    Level: developer

436: .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose()
437: @*/
438: PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
439: {

444:   PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));
445:   return(0);
446: }

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

453:    Logically Collective

455:    Input Arguments:
456: +  adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate()
457: .  name - name of the candidate scheme to add
458: .  order - order of the candidate scheme
459: .  stageorder - stage order of the candidate scheme
460: .  ccfl - stability coefficient relative to explicit Euler, used for CFL constraints
461: .  cost - relative measure of the amount of work required for the candidate scheme
462: -  inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme

464:    Note:
465:    This routine is not available in Fortran.

467:    Level: developer

469: .seealso: TSAdaptCandidatesClear(), TSAdaptChoose()
470: @*/
471: PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse)
472: {
473:   PetscInt c;

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

485:   adapt->candidates.name[c]       = name;
486:   adapt->candidates.order[c]      = order;
487:   adapt->candidates.stageorder[c] = stageorder;
488:   adapt->candidates.ccfl[c]       = ccfl;
489:   adapt->candidates.cost[c]       = cost;
490:   adapt->candidates.n++;
491:   return(0);
492: }

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

499:    Not Collective

501:    Input Arguments:
502: .  adapt - time step adaptivity context

504:    Output Arguments:
505: +  n - number of candidate schemes, always at least 1
506: .  order - the order of each candidate scheme
507: .  stageorder - the stage order of each candidate scheme
508: .  ccfl - the CFL coefficient of each scheme
509: -  cost - the relative cost of each scheme

511:    Level: developer

513:    Note:
514:    The current scheme is always returned in the first slot

516: .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose()
517: @*/
518: PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost)
519: {
522:   if (n) *n = adapt->candidates.n;
523:   if (order) *order = adapt->candidates.order;
524:   if (stageorder) *stageorder = adapt->candidates.stageorder;
525:   if (ccfl) *ccfl = adapt->candidates.ccfl;
526:   if (cost) *cost = adapt->candidates.cost;
527:   return(0);
528: }

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

535:    Logically Collective

537:    Input Arguments:
538: +  adapt - adaptive contoller
539: -  h - current step size

541:    Output Arguments:
542: +  next_sc - scheme to use for the next step
543: .  next_h - step size to use for the next step
544: -  accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size

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

550:    Level: developer

552: .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd()
553: @*/
554: PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept)
555: {
557:   PetscReal      wlte = -1.0;

565:   if (adapt->candidates.n < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"%D candidates have been registered",adapt->candidates.n);
566:   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);
567:   (*adapt->ops->choose)(adapt,ts,h,next_sc,next_h,accept,&wlte);
568:   if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
569:     /* Reduce time step if it overshoots max time */
570:     PetscReal max_time = ts->max_time;
571:     PetscReal next_dt  = 0.0;
572:     if (ts->ptime + ts->time_step + *next_h >= max_time) {
573:       next_dt = max_time - (ts->ptime + ts->time_step);
574:       if (next_dt > PETSC_SMALL) *next_h = next_dt;
575:       else ts->reason = TS_CONVERGED_TIME;
576:     }
577:   }
578:   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);
579:   if (!(*next_h > 0.)) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h);

581:   if (adapt->monitor) {
582:     PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);
583:     if (wlte < 0) {
584:       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);
585:     } else {
586:       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);
587:     }
588:     PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);
589:   }
590:   return(0);
591: }

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

598:    Collective

600:    Input Arguments:
601: +  adapt - adaptive controller context
602: -  ts - time stepper

604:    Output Arguments:
605: .  accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject

607:    Level: developer

609: .seealso:
610: @*/
611: PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscBool *accept)
612: {
613:   PetscErrorCode      ierr;
614:   SNES                snes;
615:   SNESConvergedReason snesreason;

621:   *accept = PETSC_TRUE;
622:   TSGetSNES(ts,&snes);
623:   SNESGetConvergedReason(snes,&snesreason);
624:   if (snesreason < 0) {
625:     PetscReal dt,new_dt;
626:     *accept = PETSC_FALSE;
627:     TSGetTimeStep(ts,&dt);
628:     if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) {
629:       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
630:       PetscInfo2(ts,"Step=%D, nonlinear solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);
631:       if (adapt->monitor) {
632:         PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);
633:         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);
634:         PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);
635:       }
636:     } else {
637:       new_dt = dt*adapt->scale_solve_failed;
638:       TSSetTimeStep(ts,new_dt);
639:       if (adapt->monitor) {
640:         PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);
641:         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);
642:         PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);
643:       }
644:     }
645:   }
646:   if (adapt->checkstage) {(*adapt->checkstage)(adapt,ts,accept);}
647:   return(0);
648: }



654: /*@
655:   TSAdaptCreate - create an adaptive controller context for time stepping

657:   Collective on MPI_Comm

659:   Input Parameter:
660: . comm - The communicator

662:   Output Parameter:
663: . adapt - new TSAdapt object

665:   Level: developer

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

670: .keywords: TSAdapt, create
671: .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy()
672: @*/
673: PetscErrorCode  TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt)
674: {
676:   TSAdapt        adapt;

680:   *inadapt = NULL;
681:   TSAdaptInitializePackage();

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

685:   adapt->dt_min             = 1e-20;
686:   adapt->dt_max             = 1e50;
687:   adapt->scale_solve_failed = 0.25;
688:   adapt->wnormtype          = NORM_2;

690:   *inadapt = adapt;
691:   return(0);
692: }