1 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2 3 typedef struct { 4 PetscBool always_accept; 5 PetscReal clip[2]; /* admissible decrease/increase factors */ 6 PetscReal safety; /* safety factor relative to target error */ 7 PetscReal reject_safety; /* extra safety factor if the last step was rejected */ 8 Vec X; 9 Vec Y; 10 Vec E; 11 } TSAdapt_GLEE; 12 13 #undef __FUNCT__ 14 #define __FUNCT__ "TSAdaptChoose_GLEE" 15 static PetscErrorCode TSAdaptChoose_GLEE(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte) 16 { 17 TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 18 TSType time_scheme; /* Type of time-integration scheme */ 19 PetscErrorCode ierr; 20 Vec X,Y,E; 21 PetscReal enorm,enorma,enormr,hfac_lte,h_lte,safety; 22 PetscInt order,order_add,stepno; 23 24 PetscFunctionBegin; 25 ierr = TSGetTimeStepNumber(ts,&stepno);CHKERRQ(ierr); 26 27 safety = glee->safety; 28 ierr = TSGetType(ts,&time_scheme);CHKERRQ(ierr); 29 if (!strcmp(time_scheme,TSGLEE)){ 30 /* the method is of GLEE type */ 31 order_add=0; /* typically same order estimates */ 32 ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 33 /* ierr = TSGetPreviousSolution(ts,&Y);CHKERRQ(ierr); 34 */ 35 if (!glee->E) {ierr = VecDuplicate(X,&glee->E);CHKERRQ(ierr);} 36 E = glee->E; 37 ierr = TSGetTimeError(ts,0,&E);CHKERRQ(ierr); 38 /* this should be called with Y */ 39 ierr = TSErrorWeightedENorm(ts,E,X,X,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr); 40 } else { 41 /* the method is NOT of GLEE type */ 42 order_add=0; /* typically lower order estimates */ 43 ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 44 if (!glee->Y) {ierr = VecDuplicate(X,&glee->Y);CHKERRQ(ierr);} 45 Y = glee->Y; 46 order = adapt->candidates.order[0]; 47 ierr = TSEvaluateStep(ts,order-1,Y,NULL);CHKERRQ(ierr); 48 ierr = TSErrorWeightedNorm(ts,X,Y,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr); 49 } 50 if (enorm > 1.) { 51 if (!*accept) safety *= glee->reject_safety; /* The last attempt also failed, shorten more aggressively */ 52 if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) { 53 ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting because step size %g is at minimum\n",(double)enorm,(double)h);CHKERRQ(ierr); 54 *accept = PETSC_TRUE; 55 } else if (glee->always_accept) { 56 ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g because always_accept is set\n",(double)enorm,(double)h);CHKERRQ(ierr); 57 *accept = PETSC_TRUE; 58 } else { 59 ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, rejecting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr); 60 *accept = PETSC_FALSE; 61 } 62 } else { 63 ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr); 64 *accept = PETSC_TRUE; 65 } 66 67 /* The optimal new step based purely on local truncation error for this step. */ 68 if (enorm == 0.0) { 69 hfac_lte = safety * PETSC_INFINITY; 70 } else { 71 hfac_lte = safety * PetscPowReal(enorm,-1./(order+order_add)); 72 } 73 h_lte = h * PetscClipInterval(hfac_lte,glee->clip[0],glee->clip[1]); 74 75 *next_sc = 0; 76 *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max); 77 *wlte = enorm; 78 PetscFunctionReturn(0); 79 } 80 81 #undef __FUNCT__ 82 #define __FUNCT__ "TSAdaptReset_GLEE" 83 static PetscErrorCode TSAdaptReset_GLEE(TSAdapt adapt) 84 { 85 TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 86 PetscErrorCode ierr; 87 88 PetscFunctionBegin; 89 ierr = VecDestroy(&glee->Y);CHKERRQ(ierr); 90 PetscFunctionReturn(0); 91 } 92 93 #undef __FUNCT__ 94 #define __FUNCT__ "TSAdaptDestroy_GLEE" 95 static PetscErrorCode TSAdaptDestroy_GLEE(TSAdapt adapt) 96 { 97 PetscErrorCode ierr; 98 99 PetscFunctionBegin; 100 ierr = TSAdaptReset_GLEE(adapt);CHKERRQ(ierr); 101 ierr = PetscFree(adapt->data);CHKERRQ(ierr); 102 PetscFunctionReturn(0); 103 } 104 105 #undef __FUNCT__ 106 #define __FUNCT__ "TSAdaptSetFromOptions_GLEE" 107 static PetscErrorCode TSAdaptSetFromOptions_GLEE(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 108 { 109 TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 110 PetscErrorCode ierr; 111 PetscInt two; 112 PetscBool set; 113 114 PetscFunctionBegin; 115 ierr = PetscOptionsHead(PetscOptionsObject,"GLEE adaptive controller options");CHKERRQ(ierr); 116 two = 2; 117 ierr = PetscOptionsRealArray("-ts_adapt_glee_clip","Admissible decrease/increase in step size","",glee->clip,&two,&set);CHKERRQ(ierr); 118 if (set && (two != 2 || glee->clip[0] > glee->clip[1])) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_glee_clip"); 119 ierr = PetscOptionsReal("-ts_adapt_glee_safety","Safety factor relative to target error","",glee->safety,&glee->safety,NULL);CHKERRQ(ierr); 120 ierr = PetscOptionsReal("-ts_adapt_glee_reject_safety","Extra safety factor to apply if the last step was rejected","",glee->reject_safety,&glee->reject_safety,NULL);CHKERRQ(ierr); 121 ierr = PetscOptionsBool("-ts_adapt_glee_always_accept","Always accept the step regardless of whether local truncation error meets goal","",glee->always_accept,&glee->always_accept,NULL);CHKERRQ(ierr); 122 ierr = PetscOptionsTail();CHKERRQ(ierr); 123 PetscFunctionReturn(0); 124 } 125 126 #undef __FUNCT__ 127 #define __FUNCT__ "TSAdaptView_GLEE" 128 static PetscErrorCode TSAdaptView_GLEE(TSAdapt adapt,PetscViewer viewer) 129 { 130 TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 131 PetscErrorCode ierr; 132 PetscBool iascii; 133 134 PetscFunctionBegin; 135 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 136 if (iascii) { 137 if (glee->always_accept) {ierr = PetscViewerASCIIPrintf(viewer," GLEE: always accepting steps\n");CHKERRQ(ierr);} 138 ierr = PetscViewerASCIIPrintf(viewer," GLEE: clip fastest decrease %g, fastest increase %g\n",(double)glee->clip[0],(double)glee->clip[1]);CHKERRQ(ierr); 139 ierr = PetscViewerASCIIPrintf(viewer," GLEE: safety factor %g, extra factor after step rejection %g\n",(double)glee->safety,(double)glee->reject_safety);CHKERRQ(ierr); 140 } 141 PetscFunctionReturn(0); 142 } 143 144 #undef __FUNCT__ 145 #define __FUNCT__ "TSAdaptCreate_GLEE" 146 /*MC 147 TSADAPTGLEE - GLEE adaptive controller for time stepping 148 149 Level: intermediate 150 151 .seealso: TS, TSAdapt, TSSetAdapt() 152 M*/ 153 PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt adapt) 154 { 155 PetscErrorCode ierr; 156 TSAdapt_GLEE *a; 157 158 PetscFunctionBegin; 159 ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr); 160 adapt->data = (void*)a; 161 adapt->ops->choose = TSAdaptChoose_GLEE; 162 adapt->ops->setfromoptions = TSAdaptSetFromOptions_GLEE; 163 adapt->ops->destroy = TSAdaptDestroy_GLEE; 164 adapt->ops->view = TSAdaptView_GLEE; 165 166 a->clip[0] = 0.1; 167 a->clip[1] = 10.; 168 a->safety = 0.99; 169 a->reject_safety = 0.5; 170 a->always_accept = PETSC_FALSE; 171 PetscFunctionReturn(0); 172 } 173