1 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2 3 typedef struct { 4 PetscReal safety; /* safety factor relative to target CFL constraint */ 5 PetscBool always_accept; 6 } TSAdapt_CFL; 7 8 static PetscErrorCode TSAdaptChoose_CFL(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte,PetscReal *wltea,PetscReal *wlter) 9 { 10 TSAdapt_CFL *cfl = (TSAdapt_CFL*)adapt->data; 11 PetscErrorCode ierr; 12 PetscReal hcfl,cfltimestep,ccfl; 13 PetscInt ncandidates; 14 const PetscReal *ccflarray; 15 16 PetscFunctionBegin; 17 ierr = TSGetCFLTime(ts,&cfltimestep);CHKERRQ(ierr); 18 ierr = TSAdaptCandidatesGet(adapt,&ncandidates,NULL,NULL,&ccflarray,NULL);CHKERRQ(ierr); 19 ccfl = (ncandidates > 0) ? ccflarray[0] : 1.0; 20 21 /* Determine whether the step is accepted of rejected */ 22 *accept = PETSC_TRUE; 23 if (h > cfltimestep * ccfl) { 24 if (cfl->always_accept) { 25 ierr = 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,(double)cfltimestep);CHKERRQ(ierr); 26 } else { 27 ierr = 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,(double)cfltimestep);CHKERRQ(ierr); 28 *accept = PETSC_FALSE; 29 } 30 } 31 32 /* The optimal new step based purely on CFL constraint for this step. */ 33 hcfl = cfl->safety * cfltimestep * ccfl; 34 if (hcfl < adapt->dt_min) { 35 ierr = PetscInfo4(adapt,"Cannot satisfy CFL constraint %g (with %g safety) at minimum time step %g with method coefficient %g, proceding anyway\n",(double)cfltimestep,(double)cfl->safety,(double)adapt->dt_min,(double)ccfl);CHKERRQ(ierr); 36 } 37 38 *next_sc = 0; 39 *next_h = PetscClipInterval(hcfl,adapt->dt_min,adapt->dt_max); 40 *wlte = -1; /* Weighted local truncation error was not evaluated */ 41 *wltea = -1; /* Weighted absolute local truncation error was not evaluated */ 42 *wlter = -1; /* Weighted relative local truncation error was not evaluated */ 43 PetscFunctionReturn(0); 44 } 45 46 static PetscErrorCode TSAdaptDestroy_CFL(TSAdapt adapt) 47 { 48 PetscErrorCode ierr; 49 50 PetscFunctionBegin; 51 ierr = PetscFree(adapt->data);CHKERRQ(ierr); 52 PetscFunctionReturn(0); 53 } 54 55 static PetscErrorCode TSAdaptSetFromOptions_CFL(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 56 { 57 TSAdapt_CFL *cfl = (TSAdapt_CFL*)adapt->data; 58 PetscErrorCode ierr; 59 60 PetscFunctionBegin; 61 ierr = PetscOptionsHead(PetscOptionsObject,"CFL adaptive controller options");CHKERRQ(ierr); 62 ierr = PetscOptionsReal("-ts_adapt_cfl_safety","Safety factor relative to target CFL constraint","",cfl->safety,&cfl->safety,NULL);CHKERRQ(ierr); 63 ierr = PetscOptionsBool("-ts_adapt_cfl_always_accept","Always accept the step regardless of whether CFL constraint meets goal","",cfl->always_accept,&cfl->always_accept,NULL);CHKERRQ(ierr); 64 /* The TS implementations do not currently communicate CFL information to the controller. There is a placeholder, but 65 * we do not believe it to provide sufficiently rich information. That means the CFL adaptor is incomplete and 66 * unusable. Do not delete the guard below unless you have finished the implementation. */ 67 if (!cfl->always_accept) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Step rejection not implemented. The CFL implementation is incomplete/unusable"); 68 ierr = PetscOptionsTail();CHKERRQ(ierr); 69 PetscFunctionReturn(0); 70 } 71 72 /*MC 73 TSADAPTCFL - CFL adaptive controller for time stepping 74 75 Level: intermediate 76 77 .seealso: TS, TSAdapt, TSSetAdapt() 78 M*/ 79 PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt adapt) 80 { 81 PetscErrorCode ierr; 82 TSAdapt_CFL *a; 83 84 PetscFunctionBegin; 85 ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr); 86 adapt->data = (void*)a; 87 88 adapt->ops->choose = TSAdaptChoose_CFL; 89 adapt->ops->setfromoptions = TSAdaptSetFromOptions_CFL; 90 adapt->ops->destroy = TSAdaptDestroy_CFL; 91 92 a->safety = 0.9; 93 a->always_accept = PETSC_FALSE; 94 PetscFunctionReturn(0); 95 } 96