1 #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 2 3 typedef struct { 4 PetscBool always_accept; 5 PetscReal safety; /* safety factor relative to target error */ 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,cfltime; 15 PetscInt stepno,ncandidates; 16 const PetscInt *order; 17 const PetscReal *ccfl; 18 19 PetscFunctionBegin; 20 ierr = TSGetTimeStepNumber(ts,&stepno);CHKERRQ(ierr); 21 ierr = TSGetCFLTime(ts,&cfltime);CHKERRQ(ierr); 22 ierr = TSAdaptCandidatesGet(adapt,&ncandidates,&order,PETSC_NULL,&ccfl,PETSC_NULL);CHKERRQ(ierr); 23 24 hcfl = cfl->safety * cfltime * ccfl[0]; 25 if (hcfl < adapt->dt_min) { 26 ierr = PetscInfo4(adapt,"Cannot satisfy CFL constraint %G (with %G safety) at minimum time step %G with method coefficient %G, proceding anyway\n",cfltime,cfl->safety,adapt->dt_min,ccfl[0]);CHKERRQ(ierr); 27 } 28 29 if (h > cfltime * ccfl[0]) { 30 if (cfl->always_accept) { 31 ierr = PetscInfo3(adapt,"Step length %G with scheme of CFL coefficient %G did not satisfy user-provided CFL constraint %G, proceeding anyway\n",h,ccfl[0],cfltime);CHKERRQ(ierr); 32 } else { 33 ierr = PetscInfo3(adapt,"Step length %G with scheme of CFL coefficient %G did not satisfy user-provided CFL constraint %G, step REJECTED\n",h,ccfl[0],cfltime);CHKERRQ(ierr); 34 *next_sc = 0; 35 *next_h = PetscClipInterval(hcfl,adapt->dt_min,adapt->dt_max); 36 *accept = PETSC_FALSE; 37 } 38 } 39 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 PetscFunctionReturn(0); 45 } 46 47 #undef __FUNCT__ 48 #define __FUNCT__ "TSAdaptDestroy_CFL" 49 static PetscErrorCode TSAdaptDestroy_CFL(TSAdapt adapt) 50 { 51 PetscErrorCode ierr; 52 53 PetscFunctionBegin; 54 ierr = PetscFree(adapt->data);CHKERRQ(ierr); 55 PetscFunctionReturn(0); 56 } 57 58 #undef __FUNCT__ 59 #define __FUNCT__ "TSAdaptSetFromOptions_CFL" 60 static PetscErrorCode TSAdaptSetFromOptions_CFL(TSAdapt adapt) 61 { 62 TSAdapt_CFL *cfl = (TSAdapt_CFL*)adapt->data; 63 PetscErrorCode ierr; 64 65 PetscFunctionBegin; 66 ierr = PetscOptionsHead("CFL adaptive controller options");CHKERRQ(ierr); 67 ierr = PetscOptionsReal("-ts_adapt_cfl_safety","Safety factor relative to target error","",cfl->safety,&cfl->safety,PETSC_NULL);CHKERRQ(ierr); 68 ierr = PetscOptionsBool("-ts_adapt_cfl_always_accept","Always accept the step regardless of whether local truncation error meets goal","",cfl->always_accept,&cfl->always_accept,PETSC_NULL);CHKERRQ(ierr); 69 if (!cfl->always_accept) SETERRQ(((PetscObject)adapt)->comm,PETSC_ERR_SUP,"step rejection not implemented yet"); 70 ierr = PetscOptionsTail();CHKERRQ(ierr); 71 PetscFunctionReturn(0); 72 } 73 74 EXTERN_C_BEGIN 75 #undef __FUNCT__ 76 #define __FUNCT__ "TSAdaptCreate_CFL" 77 /*MC 78 TSADAPTCFL - CFL adaptive controller for time stepping 79 80 Level: intermediate 81 82 .seealso: TS, TSAdapt, TSSetAdapt() 83 M*/ 84 PetscErrorCode TSAdaptCreate_CFL(TSAdapt adapt) 85 { 86 PetscErrorCode ierr; 87 TSAdapt_CFL *a; 88 89 PetscFunctionBegin; 90 ierr = PetscNewLog(adapt,TSAdapt_CFL,&a);CHKERRQ(ierr); 91 adapt->data = (void*)a; 92 adapt->ops->choose = TSAdaptChoose_CFL; 93 adapt->ops->setfromoptions = TSAdaptSetFromOptions_CFL; 94 adapt->ops->destroy = TSAdaptDestroy_CFL; 95 96 a->safety = 0.9; 97 a->always_accept = PETSC_FALSE; 98 PetscFunctionReturn(0); 99 } 100 EXTERN_C_END 101