Actual source code: tsadapt.c

petsc-3.3-p2 2012-07-13
  2: #include <petsc-private/tsimpl.h> /*I  "petscts.h" I*/

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

  9: EXTERN_C_BEGIN
 10: PetscErrorCode  TSAdaptCreate_Basic(TSAdapt);
 11: PetscErrorCode  TSAdaptCreate_None(TSAdapt);
 12: PetscErrorCode  TSAdaptCreate_CFL(TSAdapt);
 13: EXTERN_C_END

 17: /*@C
 18:    TSAdaptRegister - see TSAdaptRegisterDynamic()

 20:    Level: advanced
 21: @*/
 22: PetscErrorCode  TSAdaptRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(TSAdapt))
 23: {
 25:   char           fullname[PETSC_MAX_PATH_LEN];

 28:   PetscFListConcat(path,name,fullname);
 29:   PetscFListAdd(&TSAdaptList,sname,fullname,(void(*)(void))function);
 30:   return(0);
 31: }

 35: /*@C
 36:   TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt

 38:   Not Collective

 40:   Level: advanced

 42: .keywords: TSAdapt, register, all

 44: .seealso: TSAdaptRegisterDestroy()
 45: @*/
 46: PetscErrorCode  TSAdaptRegisterAll(const char path[])
 47: {

 51:   TSAdaptRegisterDynamic(TSADAPTBASIC,path,"TSAdaptCreate_Basic",TSAdaptCreate_Basic);
 52:   TSAdaptRegisterDynamic(TSADAPTNONE, path,"TSAdaptCreate_None", TSAdaptCreate_None);
 53:   TSAdaptRegisterDynamic(TSADAPTCFL,  path,"TSAdaptCreate_CFL",  TSAdaptCreate_CFL);
 54:   return(0);
 55: }

 59: /*@C
 60:   TSFinalizePackage - This function destroys everything in the TS package. It is
 61:   called from PetscFinalize().

 63:   Level: developer

 65: .keywords: Petsc, destroy, package
 66: .seealso: PetscFinalize()
 67: @*/
 68: PetscErrorCode  TSAdaptFinalizePackage(void)
 69: {
 71:   TSAdaptPackageInitialized = PETSC_FALSE;
 72:   TSAdaptRegisterAllCalled  = PETSC_FALSE;
 73:   TSAdaptList               = PETSC_NULL;
 74:   return(0);
 75: }

 79: /*@C
 80:   TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is
 81:   called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to
 82:   TSCreate_GL() when using static libraries.

 84:   Input Parameter:
 85:   path - The dynamic library path, or PETSC_NULL

 87:   Level: developer

 89: .keywords: TSAdapt, initialize, package
 90: .seealso: PetscInitialize()
 91: @*/
 92: PetscErrorCode  TSAdaptInitializePackage(const char path[])
 93: {

 97:   if (TSAdaptPackageInitialized) return(0);
 98:   TSAdaptPackageInitialized = PETSC_TRUE;
 99:   PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);
100:   TSAdaptRegisterAll(path);
101:   PetscRegisterFinalize(TSAdaptFinalizePackage);
102:   return(0);
103: }

107: /*@C
108:    TSAdaptRegisterDestroy - Frees the list of adaptivity schemes that were registered by TSAdaptRegister()/TSAdaptRegisterDynamic().

110:    Not Collective

112:    Level: advanced

114: .keywords: TSAdapt, register, destroy
115: .seealso: TSAdaptRegister(), TSAdaptRegisterAll(), TSAdaptRegisterDynamic()
116: @*/
117: PetscErrorCode  TSAdaptRegisterDestroy(void)
118: {

122:   PetscFListDestroy(&TSAdaptList);
123:   TSAdaptRegisterAllCalled = PETSC_FALSE;
124:   return(0);
125: }


130: PetscErrorCode  TSAdaptSetType(TSAdapt adapt,const TSAdaptType type)
131: {
132:   PetscErrorCode ierr,(*r)(TSAdapt);

135:   PetscFListFind(TSAdaptList,((PetscObject)adapt)->comm,type,PETSC_TRUE,(void(**)(void))&r);
136:   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type);
137:   if (((PetscObject)adapt)->type_name) {(*adapt->ops->destroy)(adapt);}
138:   (*r)(adapt);
139:   PetscObjectChangeTypeName((PetscObject)adapt,type);
140:   return(0);
141: }

145: PetscErrorCode  TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[])
146: {

150:   PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);
151:   return(0);
152: }

156: PetscErrorCode  TSAdaptView(TSAdapt adapt,PetscViewer viewer)
157: {
159:   PetscBool      iascii;

162:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
163:   if (iascii) {
164:     PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer,"TSAdapt Object");
165:     PetscViewerASCIIPrintf(viewer,"number of candidates %D\n",adapt->candidates.n);
166:     if (adapt->ops->view) {
167:       PetscViewerASCIIPushTab(viewer);
168:       (*adapt->ops->view)(adapt,viewer);
169:       PetscViewerASCIIPopTab(viewer);
170:     }
171:   } else {
172:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported",((PetscObject)viewer)->type_name);
173:   }
174:   return(0);
175: }

179: PetscErrorCode  TSAdaptDestroy(TSAdapt *adapt)
180: {

184:   if (!*adapt) return(0);
186:   if (--((PetscObject)(*adapt))->refct > 0) {*adapt = 0; return(0);}
187:   if ((*adapt)->ops->destroy) {(*(*adapt)->ops->destroy)(*adapt);}
188:   PetscViewerDestroy(&(*adapt)->monitor);
189:   PetscHeaderDestroy(adapt);
190:   return(0);
191: }

195: /*@
196:    TSAdaptSetMonitor - Monitor the choices made by the adaptive controller

198:    Collective on TSAdapt

200:    Input Arguments:
201: +  adapt - adaptive controller context
202: -  flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable

204:    Level: intermediate

206: .seealso: TSAdaptChoose()
207: @*/
208: PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg)
209: {

213:   if (flg) {
214:     if (!adapt->monitor) {PetscViewerASCIIOpen(((PetscObject)adapt)->comm,"stdout",&adapt->monitor);}
215:   } else {
216:     PetscViewerDestroy(&adapt->monitor);
217:   }
218:   return(0);
219: }

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

226:    Logically Collective

228:    Input Arguments:
229: +  adapt - time step adaptivity context, usually gotten with TSGetAdapt()
230: .  hmin - minimum time step
231: -  hmax - maximum time step

233:    Options Database Keys:
234: +  -ts_adapt_dt_min - minimum time step
235: -  -ts_adapt_dt_max - maximum time step

237:    Level: intermediate

239: .seealso: TSAdapt
240: @*/
241: PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax)
242: {

245:   if (hmin != PETSC_DECIDE) adapt->dt_min = hmin;
246:   if (hmax != PETSC_DECIDE) adapt->dt_max = hmax;
247:   return(0);
248: }

252: /*@
253:    TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options.

255:    Collective on TSAdapt

257:    Input Parameter:
258: .  adapt - the TSAdapt context

260:    Options Database Keys:
261: .  -ts_adapt_type <type> - basic

263:    Level: advanced

265:    Notes:
266:    This function is automatically called by TSSetFromOptions()

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

270: .seealso: TSGetType()
271: @*/
272: PetscErrorCode  TSAdaptSetFromOptions(TSAdapt adapt)
273: {
275:   char           type[256] = TSADAPTBASIC;
276:   PetscBool      set,flg;

279:   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
280:   * function can only be called from inside TSSetFromOptions_GL()  */
281:   PetscOptionsHead("TS Adaptivity options");
282:   PetscOptionsList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList,
283:                           ((PetscObject)adapt)->type_name?((PetscObject)adapt)->type_name:type,type,sizeof type,&flg);
284:   if (flg || !((PetscObject)adapt)->type_name) {
285:     TSAdaptSetType(adapt,type);
286:   }
287:   if (adapt->ops->setfromoptions) {(*adapt->ops->setfromoptions)(adapt);}
288:   PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",adapt->dt_min,&adapt->dt_min,PETSC_NULL);
289:   PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",adapt->dt_max,&adapt->dt_max,PETSC_NULL);
290:   PetscOptionsReal("-ts_adapt_scale_solve_failed","Scale step by this factor if solve fails","",adapt->scale_solve_failed,&adapt->scale_solve_failed,PETSC_NULL);
291:   PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);
292:   if (set) {TSAdaptSetMonitor(adapt,flg);}
293:   PetscOptionsTail();
294:   return(0);
295: }

299: /*@
300:    TSAdaptCandidatesClear - clear any previously set candidate schemes

302:    Logically Collective

304:    Input Argument:
305: .  adapt - adaptive controller

307:    Level: developer

309: .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose()
310: @*/
311: PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
312: {

316:   PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));
317:   return(0);
318: }

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

325:    Logically Collective

327:    Input Arguments:
328: +  adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate()
329: .  name - name of the candidate scheme to add
330: .  order - order of the candidate scheme
331: .  stageorder - stage order of the candidate scheme
332: .  ccfl - stability coefficient relative to explicit Euler, used for CFL constraints
333: .  cost - relative measure of the amount of work required for the candidate scheme
334: -  inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme

336:    Note:
337:    This routine is not available in Fortran.

339:    Level: developer

341: .seealso: TSAdaptCandidatesClear(), TSAdaptChoose()
342: @*/
343: PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse)
344: {
345:   PetscInt c;

349:   if (order < 1) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order);
350:   if (inuse) {
351:     if (adapt->candidates.inuse_set) SETERRQ(((PetscObject)adapt)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()");
352:     adapt->candidates.inuse_set = PETSC_TRUE;
353:   }
354:   /* first slot if this is the current scheme, otherwise the next available slot */
355:   c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n;
356:   adapt->candidates.name[c]         = name;
357:   adapt->candidates.order[c]        = order;
358:   adapt->candidates.stageorder[c]   = stageorder;
359:   adapt->candidates.ccfl[c]         = ccfl;
360:   adapt->candidates.cost[c]         = cost;
361:   adapt->candidates.n++;
362:   return(0);
363: }

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

370:    Not Collective

372:    Input Arguments:
373: .  adapt - time step adaptivity context

375:    Output Arguments:
376: +  n - number of candidate schemes, always at least 1
377: .  order - the order of each candidate scheme
378: .  stageorder - the stage order of each candidate scheme
379: .  ccfl - the CFL coefficient of each scheme
380: -  cost - the relative cost of each scheme

382:    Level: developer

384:    Note:
385:    The current scheme is always returned in the first slot

387: .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose()
388: @*/
389: PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost)
390: {
393:   if (n) *n = adapt->candidates.n;
394:   if (order) *order = adapt->candidates.order;
395:   if (stageorder) *stageorder = adapt->candidates.stageorder;
396:   if (ccfl) *ccfl = adapt->candidates.ccfl;
397:   if (cost) *cost = adapt->candidates.cost;
398:   return(0);
399: }

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

406:    Logically Collective

408:    Input Arguments:
409: +  adapt - adaptive contoller
410: -  h - current step size

412:    Output Arguments:
413: +  next_sc - scheme to use for the next step
414: .  next_h - step size to use for the next step
415: -  accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size

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

421:    Level: developer

423: .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd()
424: @*/
425: PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept)
426: {
428:   PetscReal wlte = -1.0;

436:   if (adapt->candidates.n < 1) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_WRONGSTATE,"%D candidates have been registered",adapt->candidates.n);
437:   if (!adapt->candidates.inuse_set) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_WRONGSTATE,"The current in-use scheme is not among the %D candidates",adapt->candidates.n);
438:   (*adapt->ops->choose)(adapt,ts,h,next_sc,next_h,accept,&wlte);
439:   if (*next_sc < 0 || adapt->candidates.n <= *next_sc) SETERRQ2(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Chosen scheme %D not in valid range 0..%D",*next_sc,adapt->candidates.n-1);
440:   if (!(*next_h > 0.)) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %G must be positive",*next_h);

442:   if (adapt->monitor) {
443:     PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);
444:     if (wlte < 0) {
445:       PetscViewerASCIIPrintf(adapt->monitor,"    TSAdapt '%s': step %3D %s t=%-11g+%10.3e family='%s' scheme=%D:'%s' dt=%-10g\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);
446:     } else {
447:       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);
448:     }
449:     PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);
450:   }
451:   return(0);
452: }

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

459:    Collective

461:    Input Arguments:
462: +  adapt - adaptive controller context
463: -  ts - time stepper

465:    Output Arguments:
466: .  accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject

468:    Level: developer

470: .seealso:
471: @*/
472: PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscBool *accept)
473: {
474:   PetscErrorCode      ierr;
475:   SNES                snes;
476:   SNESConvergedReason snesreason;

479:   *accept = PETSC_TRUE;
480:   TSGetSNES(ts,&snes);
481:   SNESGetConvergedReason(snes,&snesreason);
482:   if (snesreason < 0) {
483:     PetscReal dt,new_dt;
484:     *accept = PETSC_FALSE;
485:     TSGetTimeStep(ts,&dt);
486:     if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) {
487:       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
488:       PetscInfo2(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);
489:       if (adapt->monitor) {
490:         PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);
491:         PetscViewerASCIIPrintf(adapt->monitor,"    TSAdapt '%s': step %3D stage rejected t=%-11g+%10.3e, %D failures exceeds current TS allowed\n",((PetscObject)adapt)->type_name,ts->steps,(double)ts->ptime,dt,ts->num_snes_failures);
492:         PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);
493:       }
494:     } else {
495:       new_dt = dt*adapt->scale_solve_failed;
496:       TSSetTimeStep(ts,new_dt);
497:       if (adapt->monitor) {
498:         PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);
499:         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);
500:         PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);
501:       }
502:     }
503:   }
504:   return(0);
505: }

509: /*@
510:   TSAdaptCreate - create an adaptive controller context for time stepping

512:   Collective on MPI_Comm

514:   Input Parameter:
515: . comm - The communicator

517:   Output Parameter:
518: . adapt - new TSAdapt object

520:   Level: developer

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

525: .keywords: TSAdapt, create
526: .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy()
527: @*/
528: PetscErrorCode  TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt)
529: {
531:   TSAdapt adapt;

534:   *inadapt = 0;
535:   PetscHeaderCreate(adapt,_p_TSAdapt,struct _TSAdaptOps,TSADAPT_CLASSID,0,"TSAdapt","General Linear adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);

537:   adapt->dt_min             = 1e-20;
538:   adapt->dt_max             = 1e50;
539:   adapt->scale_solve_failed = 0.25;

541:   *inadapt = adapt;
542:   return(0);
543: }