Actual source code: adaptcfl.c
petsc-3.5.4 2015-05-23
1: #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/
3: typedef struct {
4: PetscBool always_accept;
5: PetscReal safety; /* safety factor relative to target error */
6: } TSAdapt_CFL;
10: static PetscErrorCode TSAdaptChoose_CFL(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte)
11: {
12: TSAdapt_CFL *cfl = (TSAdapt_CFL*)adapt->data;
13: PetscErrorCode ierr;
14: PetscReal hcfl,cfltime;
15: PetscInt stepno,ncandidates;
16: const PetscInt *order;
17: const PetscReal *ccfl;
20: TSGetTimeStepNumber(ts,&stepno);
21: TSGetCFLTime(ts,&cfltime);
22: TSAdaptCandidatesGet(adapt,&ncandidates,&order,NULL,&ccfl,NULL);
24: hcfl = cfl->safety * cfltime * ccfl[0];
25: if (hcfl < adapt->dt_min) {
26: PetscInfo4(adapt,"Cannot satisfy CFL constraint %g (with %g safety) at minimum time step %g with method coefficient %g, proceding anyway\n",(double)cfltime,(double)cfl->safety,(double)adapt->dt_min,(double)ccfl[0]);
27: }
29: if (h > cfltime * ccfl[0]) {
30: if (cfl->always_accept) {
31: PetscInfo3(adapt,"Step length %g with scheme of CFL coefficient %g did not satisfy user-provided CFL constraint %g, proceeding anyway\n",(double)h,(double)ccfl[0],(double)cfltime);
32: } else {
33: PetscInfo3(adapt,"Step length %g with scheme of CFL coefficient %g did not satisfy user-provided CFL constraint %g, step REJECTED\n",(double)h,(double)ccfl[0],(double)cfltime);
34: *next_sc = 0;
35: *next_h = PetscClipInterval(hcfl,adapt->dt_min,adapt->dt_max);
36: *accept = PETSC_FALSE;
37: }
38: }
40: *next_sc = 0;
41: *next_h = PetscClipInterval(hcfl,adapt->dt_min,adapt->dt_max);
42: *accept = PETSC_TRUE;
43: *wlte = -1; /* Weighted local truncation error was not evaluated */
44: return(0);
45: }
49: static PetscErrorCode TSAdaptDestroy_CFL(TSAdapt adapt)
50: {
54: PetscFree(adapt->data);
55: return(0);
56: }
60: static PetscErrorCode TSAdaptSetFromOptions_CFL(TSAdapt adapt)
61: {
62: TSAdapt_CFL *cfl = (TSAdapt_CFL*)adapt->data;
66: PetscOptionsHead("CFL adaptive controller options");
67: PetscOptionsReal("-ts_adapt_cfl_safety","Safety factor relative to target error","",cfl->safety,&cfl->safety,NULL);
68: PetscOptionsBool("-ts_adapt_cfl_always_accept","Always accept the step regardless of whether local truncation error meets goal","",cfl->always_accept,&cfl->always_accept,NULL);
69: if (!cfl->always_accept) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"step rejection not implemented yet");
70: PetscOptionsTail();
71: return(0);
72: }
76: /*MC
77: TSADAPTCFL - CFL adaptive controller for time stepping
79: Level: intermediate
81: .seealso: TS, TSAdapt, TSSetAdapt()
82: M*/
83: PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt adapt)
84: {
86: TSAdapt_CFL *a;
89: PetscNewLog(adapt,&a);
90: adapt->data = (void*)a;
91: adapt->ops->choose = TSAdaptChoose_CFL;
92: adapt->ops->setfromoptions = TSAdaptSetFromOptions_CFL;
93: adapt->ops->destroy = TSAdaptDestroy_CFL;
95: a->safety = 0.9;
96: a->always_accept = PETSC_FALSE;
97: return(0);
98: }