Actual source code: tsadapt.c
petsc-3.8.0 2017-09-26
2: #include <petsc/private/tsimpl.h>
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);
13: PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(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: }
52: /*@C
53: TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt
55: Not Collective
57: Level: advanced
59: .keywords: TSAdapt, register, all
61: .seealso: TSAdaptRegisterDestroy()
62: @*/
63: PetscErrorCode TSAdaptRegisterAll(void)
64: {
68: if (TSAdaptRegisterAllCalled) return(0);
69: TSAdaptRegisterAllCalled = PETSC_TRUE;
70: TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None);
71: TSAdaptRegister(TSADAPTBASIC,TSAdaptCreate_Basic);
72: TSAdaptRegister(TSADAPTCFL, TSAdaptCreate_CFL);
73: TSAdaptRegister(TSADAPTGLEE ,TSAdaptCreate_GLEE);
74: return(0);
75: }
77: /*@C
78: TSAdaptFinalizePackage - This function destroys everything in the TS package. It is
79: called from PetscFinalize().
81: Level: developer
83: .keywords: Petsc, destroy, package
84: .seealso: PetscFinalize()
85: @*/
86: PetscErrorCode TSAdaptFinalizePackage(void)
87: {
91: PetscFunctionListDestroy(&TSAdaptList);
92: TSAdaptPackageInitialized = PETSC_FALSE;
93: TSAdaptRegisterAllCalled = PETSC_FALSE;
94: return(0);
95: }
97: /*@C
98: TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is
99: called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to
100: TSAdaptCreate() when using static libraries.
102: Level: developer
104: .keywords: TSAdapt, initialize, package
105: .seealso: PetscInitialize()
106: @*/
107: PetscErrorCode TSAdaptInitializePackage(void)
108: {
112: if (TSAdaptPackageInitialized) return(0);
113: TSAdaptPackageInitialized = PETSC_TRUE;
114: PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);
115: TSAdaptRegisterAll();
116: PetscRegisterFinalize(TSAdaptFinalizePackage);
117: return(0);
118: }
120: /*@C
121: TSAdaptSetType - sets the approach used for the error adapter, currently there is only TSADAPTBASIC and TSADAPTNONE
123: Logicially Collective on TSAdapt
125: Input Parameter:
126: + adapt - the TS adapter, most likely obtained with TSGetAdapt()
127: - type - either TSADAPTBASIC or TSADAPTNONE
129: Options Database:
130: . -ts_adapt_type basic or none - to setting the adapter type
132: Level: intermediate
134: .keywords: TSAdapt, create
136: .seealso: TSGetAdapt(), TSAdaptDestroy(), TSAdaptType, TSAdaptGetType()
137: @*/
138: PetscErrorCode TSAdaptSetType(TSAdapt adapt,TSAdaptType type)
139: {
140: PetscBool match;
141: PetscErrorCode ierr,(*r)(TSAdapt);
146: PetscObjectTypeCompare((PetscObject)adapt,type,&match);
147: if (match) return(0);
148: PetscFunctionListFind(TSAdaptList,type,&r);
149: if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type);
150: if (adapt->ops->destroy) {(*adapt->ops->destroy)(adapt);}
151: PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));
152: PetscObjectChangeTypeName((PetscObject)adapt,type);
153: (*r)(adapt);
154: return(0);
155: }
157: /*@C
158: TSAdaptGetType - gets the TS adapter method type (as a string).
160: Not Collective
162: Input Parameter:
163: . adapt - The TS adapter, most likely obtained with TSGetAdapt()
165: Output Parameter:
166: . type - The name of TS adapter method
168: Level: intermediate
170: .keywords: TSAdapt, get, type
171: .seealso TSAdaptSetType()
172: @*/
173: PetscErrorCode TSAdaptGetType(TSAdapt adapt,TSAdaptType *type)
174: {
178: *type = ((PetscObject)adapt)->type_name;
179: return(0);
180: }
182: PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[])
183: {
188: PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);
189: return(0);
190: }
192: /*@C
193: TSAdaptLoad - Loads a TSAdapt that has been stored in binary with TSAdaptView().
195: Collective on PetscViewer
197: Input Parameters:
198: + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or
199: some related function before a call to TSAdaptLoad().
200: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
201: HDF5 file viewer, obtained from PetscViewerHDF5Open()
203: Level: intermediate
205: Notes:
206: The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored.
208: Notes for advanced users:
209: Most users should not need to know the details of the binary storage
210: format, since TSAdaptLoad() and TSAdaptView() completely hide these details.
211: But for anyone who's interested, the standard binary matrix storage
212: format is
213: .vb
214: has not yet been determined
215: .ve
217: .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad()
218: @*/
219: PetscErrorCode TSAdaptLoad(TSAdapt adapt,PetscViewer viewer)
220: {
222: PetscBool isbinary;
223: char type[256];
228: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
229: if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
231: PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);
232: TSAdaptSetType(adapt,type);
233: if (adapt->ops->load) {
234: (*adapt->ops->load)(adapt,viewer);
235: }
236: return(0);
237: }
239: PetscErrorCode TSAdaptView(TSAdapt adapt,PetscViewer viewer)
240: {
242: PetscBool iascii,isbinary,isnone;
246: if (!viewer) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt),&viewer);}
249: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
250: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
251: if (iascii) {
252: PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);
253: PetscObjectTypeCompare((PetscObject)adapt,TSADAPTNONE,&isnone);
254: if (!isnone) {
255: if (adapt->always_accept) {PetscViewerASCIIPrintf(viewer," always accepting steps\n");}
256: PetscViewerASCIIPrintf(viewer," safety factor %g\n",(double)adapt->safety);
257: PetscViewerASCIIPrintf(viewer," extra safety factor after step rejection %g\n",(double)adapt->reject_safety);
258: PetscViewerASCIIPrintf(viewer," clip fastest increase %g\n",(double)adapt->clip[1]);
259: PetscViewerASCIIPrintf(viewer," clip fastest decrease %g\n",(double)adapt->clip[0]);
260: PetscViewerASCIIPrintf(viewer," maximum allowed timestep %g\n",(double)adapt->dt_max);
261: PetscViewerASCIIPrintf(viewer," minimum allowed timestep %g\n",(double)adapt->dt_min);
262: }
263: if (adapt->ops->view) {
264: PetscViewerASCIIPushTab(viewer);
265: (*adapt->ops->view)(adapt,viewer);
266: PetscViewerASCIIPopTab(viewer);
267: }
268: } else if (isbinary) {
269: char type[256];
271: /* need to save FILE_CLASS_ID for adapt class */
272: PetscStrncpy(type,((PetscObject)adapt)->type_name,256);
273: PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);
274: } else if (adapt->ops->view) {
275: (*adapt->ops->view)(adapt,viewer);
276: }
277: return(0);
278: }
280: /*@
281: TSAdaptReset - Resets a TSAdapt context.
283: Collective on TS
285: Input Parameter:
286: . adapt - the TSAdapt context obtained from TSAdaptCreate()
288: Level: developer
290: .seealso: TSAdaptCreate(), TSAdaptDestroy()
291: @*/
292: PetscErrorCode TSAdaptReset(TSAdapt adapt)
293: {
298: if (adapt->ops->reset) {(*adapt->ops->reset)(adapt);}
299: return(0);
300: }
302: PetscErrorCode TSAdaptDestroy(TSAdapt *adapt)
303: {
307: if (!*adapt) return(0);
309: if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; return(0);}
311: TSAdaptReset(*adapt);
313: if ((*adapt)->ops->destroy) {(*(*adapt)->ops->destroy)(*adapt);}
314: PetscViewerDestroy(&(*adapt)->monitor);
315: PetscHeaderDestroy(adapt);
316: return(0);
317: }
319: /*@
320: TSAdaptSetMonitor - Monitor the choices made by the adaptive controller
322: Collective on TSAdapt
324: Input Arguments:
325: + adapt - adaptive controller context
326: - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable
328: Options Database Keys:
329: . -ts_adapt_monitor
331: Level: intermediate
333: .seealso: TSAdaptChoose()
334: @*/
335: PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg)
336: {
342: if (flg) {
343: if (!adapt->monitor) {PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);}
344: } else {
345: PetscViewerDestroy(&adapt->monitor);
346: }
347: return(0);
348: }
350: /*@C
351: TSAdaptSetCheckStage - Set a callback to check convergence for a stage
353: Logically collective on TSAdapt
355: Input Arguments:
356: + adapt - adaptive controller context
357: - func - stage check function
359: Arguments of func:
360: $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept)
362: + adapt - adaptive controller context
363: . ts - time stepping context
364: - accept - pending choice of whether to accept, can be modified by this routine
366: Level: advanced
368: .seealso: TSAdaptChoose()
369: @*/
370: PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscReal,Vec,PetscBool*))
371: {
375: adapt->checkstage = func;
376: return(0);
377: }
379: /*@
380: TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of
381: any error or stability condition not meeting the prescribed goal.
383: Logically collective on TSAdapt
385: Input Arguments:
386: + adapt - time step adaptivity context, usually gotten with TSGetAdapt()
387: - flag - whether to always accept steps
389: Options Database Keys:
390: . -ts_adapt_always_accept
392: Level: intermediate
394: .seealso: TSAdapt, TSAdaptChoose()
395: @*/
396: PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt,PetscBool flag)
397: {
401: adapt->always_accept = flag;
402: return(0);
403: }
405: /*@
406: TSAdaptSetSafety - Set safety factors
408: Logically collective on TSAdapt
410: Input Arguments:
411: + adapt - adaptive controller context
412: . safety - safety factor relative to target error/stability goal
413: + reject_safety - extra safety factor to apply if the last step was rejected
415: Options Database Keys:
416: + -ts_adapt_safety
417: - -ts_adapt_reject_safety
419: Level: intermediate
421: .seealso: TSAdapt, TSAdaptGetSafety(), TSAdaptChoose()
422: @*/
423: PetscErrorCode TSAdaptSetSafety(TSAdapt adapt,PetscReal safety,PetscReal reject_safety)
424: {
429: if (safety != PETSC_DEFAULT && safety < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Safety factor %g must be non negative",(double)safety);
430: if (safety != PETSC_DEFAULT && safety > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Safety factor %g must be less than one",(double)safety);
431: if (reject_safety != PETSC_DEFAULT && reject_safety < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reject safety factor %g must be non negative",(double)reject_safety);
432: if (reject_safety != PETSC_DEFAULT && reject_safety > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Reject safety factor %g must be less than one",(double)reject_safety);
433: if (safety != PETSC_DEFAULT) adapt->safety = safety;
434: if (reject_safety != PETSC_DEFAULT) adapt->reject_safety = reject_safety;
435: return(0);
436: }
438: /*@
439: TSAdaptGetSafety - Get safety factors
441: Not Collective
443: Input Arguments:
444: . adapt - adaptive controller context
446: Ouput Arguments:
447: . safety - safety factor relative to target error/stability goal
448: + reject_safety - extra safety factor to apply if the last step was rejected
450: Level: intermediate
452: .seealso: TSAdapt, TSAdaptSetSafety(), TSAdaptChoose()
453: @*/
454: PetscErrorCode TSAdaptGetSafety(TSAdapt adapt,PetscReal *safety,PetscReal *reject_safety)
455: {
460: if (safety) *safety = adapt->safety;
461: if (reject_safety) *reject_safety = adapt->reject_safety;
462: return(0);
463: }
465: /*@
466: TSAdaptSetClip - Sets the admissible decrease/increase factor in step size
468: Logically collective on TSAdapt
470: Input Arguments:
471: + adapt - adaptive controller context
472: . low - admissible decrease factor
473: + high - admissible increase factor
475: Options Database Keys:
476: . -ts_adapt_clip
478: Level: intermediate
480: .seealso: TSAdaptChoose(), TSAdaptGetClip()
481: @*/
482: PetscErrorCode TSAdaptSetClip(TSAdapt adapt,PetscReal low,PetscReal high)
483: {
488: if (low != PETSC_DEFAULT && low < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be non negative",(double)low);
489: if (low != PETSC_DEFAULT && low > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be less than one",(double)low);
490: if (high != PETSC_DEFAULT && high < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Increase factor %g must be geather than one",(double)high);
491: if (low != PETSC_DEFAULT) adapt->clip[0] = low;
492: if (high != PETSC_DEFAULT) adapt->clip[1] = high;
493: return(0);
494: }
496: /*@
497: TSAdaptGetClip - Gets the admissible decrease/increase factor in step size
499: Not Collective
501: Input Arguments:
502: . adapt - adaptive controller context
504: Ouput Arguments:
505: + low - optional, admissible decrease factor
506: - high - optional, admissible increase factor
508: Level: intermediate
510: .seealso: TSAdaptChoose(), TSAdaptSetClip()
511: @*/
512: PetscErrorCode TSAdaptGetClip(TSAdapt adapt,PetscReal *low,PetscReal *high)
513: {
518: if (low) *low = adapt->clip[0];
519: if (high) *high = adapt->clip[1];
520: return(0);
521: }
523: /*@
524: TSAdaptSetStepLimits - Set the minimum and maximum step sizes to be considered by the controller
526: Logically collective on TSAdapt
528: Input Arguments:
529: + adapt - time step adaptivity context, usually gotten with TSGetAdapt()
530: . hmin - minimum time step
531: - hmax - maximum time step
533: Options Database Keys:
534: + -ts_adapt_dt_min - minimum time step
535: - -ts_adapt_dt_max - maximum time step
537: Level: intermediate
539: .seealso: TSAdapt, TSAdaptGetStepLimits(), TSAdaptChoose()
540: @*/
541: PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax)
542: {
548: if (hmin != PETSC_DEFAULT && hmin < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmin);
549: if (hmax != PETSC_DEFAULT && hmax < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmax);
550: if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin;
551: if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax;
552: hmin = adapt->dt_min;
553: hmax = adapt->dt_max;
554: 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);
555: return(0);
556: }
558: /*@
559: TSAdaptGetStepLimits - Get the minimum and maximum step sizes to be considered by the controller
561: Not Collective
563: Input Arguments:
564: . adapt - time step adaptivity context, usually gotten with TSGetAdapt()
566: Output Arguments:
567: + hmin - minimum time step
568: - hmax - maximum time step
570: Level: intermediate
572: .seealso: TSAdapt, TSAdaptSetStepLimits(), TSAdaptChoose()
573: @*/
574: PetscErrorCode TSAdaptGetStepLimits(TSAdapt adapt,PetscReal *hmin,PetscReal *hmax)
575: {
581: if (hmin) *hmin = adapt->dt_min;
582: if (hmax) *hmax = adapt->dt_max;
583: return(0);
584: }
586: /*
587: TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options.
589: Collective on TSAdapt
591: Input Parameter:
592: . adapt - the TSAdapt context
594: Options Database Keys:
595: + -ts_adapt_type <type> - algorithm to use for adaptivity
596: . -ts_adapt_always_accept - always accept steps regardless of error/stability goals
597: . -ts_adapt_safety <safety> - safety factor relative to target error/stability goal
598: . -ts_adapt_reject_safety <safety> - extra safety factor to apply if the last step was rejected
599: . -ts_adapt_clip <low,high> - admissible time step decrease and increase factors
600: . -ts_adapt_dt_min <min> - minimum timestep to use
601: . -ts_adapt_dt_max <max> - maximum timestep to use
602: . -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails
603: - -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates
605: Level: advanced
607: Notes:
608: This function is automatically called by TSSetFromOptions()
610: .keywords: TS, TSGetAdapt(), TSAdaptSetType(), TSAdaptSetStepLimits()
612: .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptSetAlwaysAccept(), TSAdaptSetSafety(),
613: TSAdaptSetClip(), TSAdaptSetStepLimits(), TSAdaptSetMonitor()
614: */
615: PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSAdapt adapt)
616: {
618: char type[256] = TSADAPTBASIC;
619: PetscReal safety,reject_safety,clip[2],hmin,hmax;
620: PetscBool set,flg;
621: PetscInt two;
625: /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
626: * function can only be called from inside TSSetFromOptions() */
627: PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");
628: PetscOptionsFList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList,((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);
629: if (flg || !((PetscObject)adapt)->type_name) {
630: TSAdaptSetType(adapt,type);
631: }
633: PetscOptionsBool("-ts_adapt_always_accept","Always accept the step","TSAdaptSetAlwaysAccept",adapt->always_accept,&flg,&set);
634: if (set) {TSAdaptSetAlwaysAccept(adapt,flg);}
636: safety = adapt->safety; reject_safety = adapt->reject_safety;
637: PetscOptionsReal("-ts_adapt_safety","Safety factor relative to target error/stability goal","TSAdaptSetSafety",safety,&safety,&set);
638: PetscOptionsReal("-ts_adapt_reject_safety","Extra safety factor to apply if the last step was rejected","TSAdaptSetSafety",reject_safety,&reject_safety,&flg);
639: if (set || flg) {TSAdaptSetSafety(adapt,safety,reject_safety);}
641: two = 2; clip[0] = adapt->clip[0]; clip[1] = adapt->clip[1];
642: PetscOptionsRealArray("-ts_adapt_clip","Admissible decrease/increase factor in step size","TSAdaptSetClip",clip,&two,&set);
643: if (set && (two != 2)) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_clip");
644: if (set) {TSAdaptSetClip(adapt,clip[0],clip[1]);}
646: hmin = adapt->dt_min; hmax = adapt->dt_max;
647: PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",hmin,&hmin,&set);
648: PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",hmax,&hmax,&flg);
649: if (set || flg) {TSAdaptSetStepLimits(adapt,hmin,hmax);}
651: PetscOptionsReal("-ts_adapt_scale_solve_failed","Scale step by this factor if solve fails","",adapt->scale_solve_failed,&adapt->scale_solve_failed,NULL);
653: PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);
654: if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported");
656: PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);
657: if (set) {TSAdaptSetMonitor(adapt,flg);}
659: if (adapt->ops->setfromoptions) {(*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);}
660: PetscOptionsTail();
661: return(0);
662: }
664: /*@
665: TSAdaptCandidatesClear - clear any previously set candidate schemes
667: Logically collective on TSAdapt
669: Input Argument:
670: . adapt - adaptive controller
672: Level: developer
674: .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose()
675: @*/
676: PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
677: {
682: PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));
683: return(0);
684: }
686: /*@C
687: TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from
689: Logically collective on TSAdapt
691: Input Arguments:
692: + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate()
693: . name - name of the candidate scheme to add
694: . order - order of the candidate scheme
695: . stageorder - stage order of the candidate scheme
696: . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints
697: . cost - relative measure of the amount of work required for the candidate scheme
698: - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme
700: Note:
701: This routine is not available in Fortran.
703: Level: developer
705: .seealso: TSAdaptCandidatesClear(), TSAdaptChoose()
706: @*/
707: PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse)
708: {
709: PetscInt c;
713: if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order);
714: if (inuse) {
715: if (adapt->candidates.inuse_set) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()");
716: adapt->candidates.inuse_set = PETSC_TRUE;
717: }
718: /* first slot if this is the current scheme, otherwise the next available slot */
719: c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n;
721: adapt->candidates.name[c] = name;
722: adapt->candidates.order[c] = order;
723: adapt->candidates.stageorder[c] = stageorder;
724: adapt->candidates.ccfl[c] = ccfl;
725: adapt->candidates.cost[c] = cost;
726: adapt->candidates.n++;
727: return(0);
728: }
730: /*@C
731: TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost
733: Not Collective
735: Input Arguments:
736: . adapt - time step adaptivity context
738: Output Arguments:
739: + n - number of candidate schemes, always at least 1
740: . order - the order of each candidate scheme
741: . stageorder - the stage order of each candidate scheme
742: . ccfl - the CFL coefficient of each scheme
743: - cost - the relative cost of each scheme
745: Level: developer
747: Note:
748: The current scheme is always returned in the first slot
750: .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose()
751: @*/
752: PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost)
753: {
756: if (n) *n = adapt->candidates.n;
757: if (order) *order = adapt->candidates.order;
758: if (stageorder) *stageorder = adapt->candidates.stageorder;
759: if (ccfl) *ccfl = adapt->candidates.ccfl;
760: if (cost) *cost = adapt->candidates.cost;
761: return(0);
762: }
764: /*@C
765: TSAdaptChoose - choose which method and step size to use for the next step
767: Collective on TSAdapt
769: Input Arguments:
770: + adapt - adaptive contoller
771: - h - current step size
773: Output Arguments:
774: + next_sc - optional, scheme to use for the next step
775: . next_h - step size to use for the next step
776: - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size
778: Note:
779: The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is
780: being retried after an initial rejection.
782: Level: developer
784: .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd()
785: @*/
786: PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept)
787: {
789: PetscInt ncandidates = adapt->candidates.n;
790: PetscInt scheme = 0;
791: PetscReal wlte = -1.0;
792: PetscReal wltea = -1.0;
793: PetscReal wlter = -1.0;
801: if (next_sc) *next_sc = 0;
803: /* Do not mess with adaptivity while handling events*/
804: if (ts->event && ts->event->status != TSEVENT_NONE) {
805: *next_h = h;
806: *accept = PETSC_TRUE;
807: return(0);
808: }
810: (*adapt->ops->choose)(adapt,ts,h,&scheme,next_h,accept,&wlte,&wltea,&wlter);
811: 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);
812: if (*next_h < 0) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h);
813: if (next_sc) *next_sc = scheme;
815: if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
816: /* Increase/reduce step size if end time of next step is close to or overshoots max time */
817: PetscReal t = ts->ptime + ts->time_step, h = *next_h;
818: PetscReal tend = t + h, tmax = ts->max_time, hmax = tmax - t;
819: PetscReal a = (PetscReal)1.01; /* allow 1% step size increase */
820: if (t < tmax && tend > tmax) *next_h = hmax;
821: if (t < tmax && tend < tmax && h > hmax/2) *next_h = hmax/2;
822: if (t < tmax && tend < tmax && h*a > hmax) *next_h = hmax;
823: }
825: if (adapt->monitor) {
826: const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : "";
827: PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);
828: if (wlte < 0) {
829: PetscViewerASCIIPrintf(adapt->monitor," TSAdapt %s %s %D:%s step %3D %s t=%-11g+%10.3e dt=%-10.3e\n",((PetscObject)adapt)->type_name,((PetscObject)ts)->type_name,scheme,sc_name,ts->steps,*accept ? "accepted" : "rejected",(double)ts->ptime,(double)h,(double)*next_h);
830: } else {
831: PetscViewerASCIIPrintf(adapt->monitor," TSAdapt %s %s %D:%s step %3D %s t=%-11g+%10.3e dt=%-10.3e wlte=%5.3g wltea=%5.3g wlter=%5.3g\n",((PetscObject)adapt)->type_name,((PetscObject)ts)->type_name,scheme,sc_name,ts->steps,*accept ? "accepted" : "rejected",(double)ts->ptime,(double)h,(double)*next_h,(double)wlte,(double)wltea,(double)wlter);
832: }
833: PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);
834: }
835: return(0);
836: }
838: /*@
839: TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails)
841: Collective on TSAdapt
843: Input Arguments:
844: + adapt - adaptive controller context
845: . ts - time stepper
846: . t - Current simulation time
847: - Y - Current solution vector
849: Output Arguments:
850: . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject
852: Level: developer
854: .seealso:
855: @*/
856: PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool *accept)
857: {
858: PetscErrorCode ierr;
859: SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING;
866: if (ts->snes) {SNESGetConvergedReason(ts->snes,&snesreason);}
867: if (snesreason < 0) {
868: *accept = PETSC_FALSE;
869: if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) {
870: ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
871: PetscInfo2(ts,"Step=%D, nonlinear solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);
872: if (adapt->monitor) {
873: PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);
874: 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);
875: PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);
876: }
877: }
878: } else {
879: *accept = PETSC_TRUE;
880: TSFunctionDomainError(ts,t,Y,accept);
881: if(*accept && adapt->checkstage) {
882: (*adapt->checkstage)(adapt,ts,t,Y,accept);
883: }
884: }
886: if(!(*accept) && !ts->reason) {
887: PetscReal dt,new_dt;
888: TSGetTimeStep(ts,&dt);
889: new_dt = dt * adapt->scale_solve_failed;
890: TSSetTimeStep(ts,new_dt);
891: adapt->timestepjustincreased += 4;
892: if (adapt->monitor) {
893: PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);
894: PetscViewerASCIIPrintf(adapt->monitor," TSAdapt %s step %3D stage rejected (%s) t=%-11g+%10.3e retrying with dt=%-10.3e\n",((PetscObject)adapt)->type_name,ts->steps,SNESConvergedReasons[snesreason],(double)ts->ptime,(double)dt,(double)new_dt);
895: PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);
896: }
897: }
898: return(0);
899: }
901: /*@
902: TSAdaptCreate - create an adaptive controller context for time stepping
904: Collective on MPI_Comm
906: Input Parameter:
907: . comm - The communicator
909: Output Parameter:
910: . adapt - new TSAdapt object
912: Level: developer
914: Notes:
915: TSAdapt creation is handled by TS, so users should not need to call this function.
917: .keywords: TSAdapt, create
918: .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy()
919: @*/
920: PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt)
921: {
923: TSAdapt adapt;
927: *inadapt = NULL;
928: TSAdaptInitializePackage();
930: PetscHeaderCreate(adapt,TSADAPT_CLASSID,"TSAdapt","Time stepping adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);
932: adapt->always_accept = PETSC_FALSE;
933: adapt->safety = 0.9;
934: adapt->reject_safety = 0.5;
935: adapt->clip[0] = 0.1;
936: adapt->clip[1] = 10.;
937: adapt->dt_min = 1e-20;
938: adapt->dt_max = 1e+20;
939: adapt->scale_solve_failed = 0.25;
940: adapt->wnormtype = NORM_2;
942: *inadapt = adapt;
943: return(0);
944: }