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 #undef __FUNCT__ 9 #define __FUNCT__ "TSAdaptChoose_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,cfltimestep,ccfl; 15 PetscInt ncandidates; 16 const PetscReal *ccflarray; 17 18 PetscFunctionBegin; 19 ierr = TSGetCFLTime(ts,&cfltimestep);CHKERRQ(ierr); 20 ierr = TSAdaptCandidatesGet(adapt,&ncandidates,NULL,NULL,&ccflarray,NULL);CHKERRQ(ierr); 21 ccfl = (ncandidates > 0) ? ccflarray[0] : 1.0; 22 23 /* Determine whether the step is accepted of rejected */ 24 *accept = PETSC_TRUE; 25 if (h > cfltimestep * ccfl) { 26 if (cfl->always_accept) { 27 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); 28 } else { 29 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); 30 *accept = PETSC_FALSE; 31 } 32 } 33 34 /* The optimal new step based purely on CFL constraint for this step. */ 35 hcfl = cfl->safety * cfltimestep * ccfl; 36 if (hcfl < adapt->dt_min) { 37 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); 38 } 39 40 *next_sc = 0; 41 *next_h = PetscClipInterval(hcfl,adapt->dt_min,adapt->dt_max); 42 *wlte = -1; /* Weighted local truncation error was not evaluated */ 43 PetscFunctionReturn(0); 44 } 45 46 #undef __FUNCT__ 47 #define __FUNCT__ "TSAdaptDestroy_CFL" 48 static PetscErrorCode TSAdaptDestroy_CFL(TSAdapt adapt) 49 { 50 PetscErrorCode ierr; 51 52 PetscFunctionBegin; 53 ierr = PetscFree(adapt->data);CHKERRQ(ierr); 54 PetscFunctionReturn(0); 55 } 56 57 #undef __FUNCT__ 58 #define __FUNCT__ "TSAdaptSetFromOptions_CFL" 59 static PetscErrorCode TSAdaptSetFromOptions_CFL(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 60 { 61 TSAdapt_CFL *cfl = (TSAdapt_CFL*)adapt->data; 62 PetscErrorCode ierr; 63 64 PetscFunctionBegin; 65 ierr = PetscOptionsHead(PetscOptionsObject,"CFL adaptive controller options");CHKERRQ(ierr); 66 ierr = PetscOptionsReal("-ts_adapt_cfl_safety","Safety factor relative to target CFL constraint","",cfl->safety,&cfl->safety,NULL);CHKERRQ(ierr); 67 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); 68 ierr = PetscOptionsTail();CHKERRQ(ierr); 69 PetscFunctionReturn(0); 70 } 71 72 #undef __FUNCT__ 73 #define __FUNCT__ "TSAdaptCreate_CFL" 74 /*MC 75 TSADAPTCFL - CFL adaptive controller for time stepping 76 77 Level: intermediate 78 79 .seealso: TS, TSAdapt, TSSetAdapt() 80 M*/ 81 PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt adapt) 82 { 83 PetscErrorCode ierr; 84 TSAdapt_CFL *a; 85 86 PetscFunctionBegin; 87 ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr); 88 adapt->data = (void*)a; 89 90 adapt->ops->choose = TSAdaptChoose_CFL; 91 adapt->ops->setfromoptions = TSAdaptSetFromOptions_CFL; 92 adapt->ops->destroy = TSAdaptDestroy_CFL; 93 94 a->safety = 0.9; 95 a->always_accept = PETSC_FALSE; 96 PetscFunctionReturn(0); 97 } 98