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 Y; 9 } TSAdapt_Basic; 10 11 static PetscErrorCode TSAdaptChoose_Basic(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte,PetscReal *wltea,PetscReal *wlter) 12 { 13 TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 14 PetscInt order = PETSC_DECIDE; 15 PetscReal enorm = -1; 16 PetscReal enorma,enormr; 17 PetscReal safety = basic->safety; 18 PetscReal hfac_lte,h_lte; 19 PetscErrorCode ierr; 20 21 PetscFunctionBegin; 22 *next_sc = 0; /* Reuse the same order scheme */ 23 24 if (ts->ops->evaluatewlte) { 25 ierr = TSEvaluateWLTE(ts,adapt->wnormtype,&order,&enorm);CHKERRQ(ierr); 26 if (enorm >= 0 && order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed error order %D must be positive",order); 27 } else if (ts->ops->evaluatestep) { 28 if (adapt->candidates.n < 1) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"No candidate has been registered"); 29 if (!adapt->candidates.inuse_set) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"The current in-use scheme is not among the %D candidates",adapt->candidates.n); 30 if (!basic->Y) {ierr = VecDuplicate(ts->vec_sol,&basic->Y);CHKERRQ(ierr);} 31 order = adapt->candidates.order[0]; 32 ierr = TSEvaluateStep(ts,order-1,basic->Y,NULL);CHKERRQ(ierr); 33 ierr = TSErrorWeightedNorm(ts,ts->vec_sol,basic->Y,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr); 34 } 35 36 if (enorm < 0) { 37 *accept = PETSC_TRUE; 38 *next_h = h; /* Reuse the old step */ 39 *wlte = -1; /* Weighted local truncation error was not evaluated */ 40 PetscFunctionReturn(0); 41 } 42 43 /* Determine whether the step is accepted of rejected */ 44 if (enorm > 1) { 45 if (!*accept) safety *= basic->reject_safety; /* The last attempt also failed, shorten more aggressively */ 46 if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) { 47 ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting because step size %g is at minimum\n",(double)enorm,(double)h);CHKERRQ(ierr); 48 *accept = PETSC_TRUE; 49 } else if (basic->always_accept) { 50 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); 51 *accept = PETSC_TRUE; 52 } else { 53 ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, rejecting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr); 54 *accept = PETSC_FALSE; 55 } 56 } else { 57 ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr); 58 *accept = PETSC_TRUE; 59 } 60 61 /* The optimal new step based purely on local truncation error for this step. */ 62 if (enorm > 0) 63 hfac_lte = safety * PetscPowReal(enorm,((PetscReal)-1)/order); 64 else 65 hfac_lte = safety * PETSC_INFINITY; 66 h_lte = h * PetscClipInterval(hfac_lte,basic->clip[0],basic->clip[1]); 67 68 *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max); 69 *wlte = enorm; 70 *wltea = -1; /* Weighted absolute local truncation error is not used */ 71 *wlter = -1; /* Weighted relative local truncation error is not used */ 72 73 PetscFunctionReturn(0); 74 } 75 76 static PetscErrorCode TSAdaptReset_Basic(TSAdapt adapt) 77 { 78 TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 79 PetscErrorCode ierr; 80 81 PetscFunctionBegin; 82 ierr = VecDestroy(&basic->Y);CHKERRQ(ierr); 83 PetscFunctionReturn(0); 84 } 85 86 static PetscErrorCode TSAdaptDestroy_Basic(TSAdapt adapt) 87 { 88 PetscErrorCode ierr; 89 90 PetscFunctionBegin; 91 ierr = TSAdaptReset_Basic(adapt);CHKERRQ(ierr); 92 ierr = PetscFree(adapt->data);CHKERRQ(ierr); 93 PetscFunctionReturn(0); 94 } 95 96 static PetscErrorCode TSAdaptSetFromOptions_Basic(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 97 { 98 TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 99 PetscErrorCode ierr; 100 PetscInt two; 101 PetscBool set; 102 103 PetscFunctionBegin; 104 ierr = PetscOptionsHead(PetscOptionsObject,"Basic adaptive controller options");CHKERRQ(ierr); 105 two = 2; 106 ierr = PetscOptionsRealArray("-ts_adapt_basic_clip","Admissible decrease/increase factor in step size","TSAdaptBasicSetClip",basic->clip,&two,&set);CHKERRQ(ierr); 107 if (set && (two != 2)) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_basic_clip"); 108 if (set) {ierr = TSAdaptBasicSetClip(adapt,basic->clip[0],basic->clip[1]);CHKERRQ(ierr);} 109 ierr = PetscOptionsReal("-ts_adapt_basic_safety","Safety factor relative to target error","",basic->safety,&basic->safety,NULL);CHKERRQ(ierr); 110 ierr = PetscOptionsReal("-ts_adapt_basic_reject_safety","Extra safety factor to apply if the last step was rejected","",basic->reject_safety,&basic->reject_safety,NULL);CHKERRQ(ierr); 111 ierr = PetscOptionsBool("-ts_adapt_basic_always_accept","Always accept the step regardless of whether local truncation error meets goal","",basic->always_accept,&basic->always_accept,NULL);CHKERRQ(ierr); 112 ierr = PetscOptionsTail();CHKERRQ(ierr); 113 PetscFunctionReturn(0); 114 } 115 116 static PetscErrorCode TSAdaptView_Basic(TSAdapt adapt,PetscViewer viewer) 117 { 118 TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 119 PetscErrorCode ierr; 120 PetscBool iascii; 121 122 PetscFunctionBegin; 123 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 124 if (iascii) { 125 if (basic->always_accept) {ierr = PetscViewerASCIIPrintf(viewer," Basic: always accepting steps\n");CHKERRQ(ierr);} 126 ierr = PetscViewerASCIIPrintf(viewer," Basic: clip fastest decrease %g, fastest increase %g\n",(double)basic->clip[0],(double)basic->clip[1]);CHKERRQ(ierr); 127 ierr = PetscViewerASCIIPrintf(viewer," Basic: safety factor %g, extra factor after step rejection %g\n",(double)basic->safety,(double)basic->reject_safety);CHKERRQ(ierr); 128 } 129 PetscFunctionReturn(0); 130 } 131 132 /*MC 133 TSADAPTBASIC - Basic adaptive controller for time stepping 134 135 Level: intermediate 136 137 .seealso: TS, TSAdapt, TSSetAdapt() 138 M*/ 139 PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt adapt) 140 { 141 PetscErrorCode ierr; 142 TSAdapt_Basic *a; 143 144 PetscFunctionBegin; 145 ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr); 146 adapt->data = (void*)a; 147 adapt->ops->choose = TSAdaptChoose_Basic; 148 adapt->ops->setfromoptions = TSAdaptSetFromOptions_Basic; 149 adapt->ops->destroy = TSAdaptDestroy_Basic; 150 adapt->ops->view = TSAdaptView_Basic; 151 152 a->clip[0] = 0.1; 153 a->clip[1] = 10.; 154 a->safety = 0.9; 155 a->reject_safety = 0.5; 156 a->always_accept = PETSC_FALSE; 157 PetscFunctionReturn(0); 158 } 159 160 /*@ 161 TSAdaptBasicSetClip - Sets the admissible decrease/increase factor in step size 162 163 Collective on TSAdapt 164 165 Input Arguments: 166 + adapt - adaptive controller context 167 . low - admissible decrease factor 168 + high - admissible increase factor 169 170 Level: intermediate 171 172 .seealso: TSAdaptChoose() 173 @*/ 174 PetscErrorCode TSAdaptBasicSetClip(TSAdapt adapt,PetscReal low,PetscReal high) 175 { 176 TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 177 PetscBool isbasic; 178 PetscErrorCode ierr; 179 180 PetscFunctionBegin; 181 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 182 if (low != PETSC_DEFAULT && low < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be non negative",(double)low); 183 if (low != PETSC_DEFAULT && low > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be less than one",(double)low); 184 if (high != PETSC_DEFAULT && high < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Increase factor %g must be geather than one",(double)high); 185 ierr = PetscObjectTypeCompare((PetscObject)adapt,TSADAPTBASIC,&isbasic);CHKERRQ(ierr); 186 if (isbasic) { 187 if (low != PETSC_DEFAULT) basic->clip[0] = low; 188 if (high != PETSC_DEFAULT) basic->clip[1] = high; 189 } 190 PetscFunctionReturn(0); 191 } 192 193 /*@ 194 TSAdaptBasicGetClip - Gets the admissible decrease/increase factor in step size 195 196 Collective on TSAdapt 197 198 Input Arguments: 199 . adapt - adaptive controller context 200 201 Ouput Arguments: 202 + low - optional, admissible decrease factor 203 - high - optional, admissible increase factor 204 205 Level: intermediate 206 207 .seealso: TSAdaptChoose() 208 @*/ 209 PetscErrorCode TSAdaptBasicGetClip(TSAdapt adapt,PetscReal *low,PetscReal *high) 210 { 211 TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 212 PetscBool isbasic; 213 PetscErrorCode ierr; 214 215 PetscFunctionBegin; 216 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 217 if (low) PetscValidRealPointer(low,2); 218 if (high) PetscValidRealPointer(high,3); 219 ierr = PetscObjectTypeCompare((PetscObject)adapt,TSADAPTBASIC,&isbasic);CHKERRQ(ierr); 220 if (isbasic) { 221 if (low) *low = basic->clip[0]; 222 if (high) *high = basic->clip[1]; 223 } else { 224 if (low) *low = 0; 225 if (high) *high = PETSC_MAX_REAL; 226 } 227 PetscFunctionReturn(0); 228 } 229