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 if (adapt->timestepjustincreased){ 67 hfac_lte = PetscMin(hfac_lte,1.0); 68 adapt->timestepjustincreased--; 69 } 70 h_lte = h * PetscClipInterval(hfac_lte,basic->clip[0],basic->clip[1]); 71 72 *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max); 73 *wlte = enorm; 74 *wltea = -1; /* Weighted absolute local truncation error is not used */ 75 *wlter = -1; /* Weighted relative local truncation error is not used */ 76 77 PetscFunctionReturn(0); 78 } 79 80 static PetscErrorCode TSAdaptReset_Basic(TSAdapt adapt) 81 { 82 TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 83 PetscErrorCode ierr; 84 85 PetscFunctionBegin; 86 ierr = VecDestroy(&basic->Y);CHKERRQ(ierr); 87 PetscFunctionReturn(0); 88 } 89 90 static PetscErrorCode TSAdaptDestroy_Basic(TSAdapt adapt) 91 { 92 PetscErrorCode ierr; 93 94 PetscFunctionBegin; 95 ierr = TSAdaptReset_Basic(adapt);CHKERRQ(ierr); 96 ierr = PetscFree(adapt->data);CHKERRQ(ierr); 97 PetscFunctionReturn(0); 98 } 99 100 static PetscErrorCode TSAdaptSetFromOptions_Basic(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 101 { 102 TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 103 PetscErrorCode ierr; 104 PetscInt two; 105 PetscBool set; 106 107 PetscFunctionBegin; 108 ierr = PetscOptionsHead(PetscOptionsObject,"Basic adaptive controller options");CHKERRQ(ierr); 109 two = 2; 110 ierr = PetscOptionsRealArray("-ts_adapt_basic_clip","Admissible decrease/increase factor in step size","TSAdaptBasicSetClip",basic->clip,&two,&set);CHKERRQ(ierr); 111 if (set && (two != 2)) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_basic_clip"); 112 if (set) {ierr = TSAdaptBasicSetClip(adapt,basic->clip[0],basic->clip[1]);CHKERRQ(ierr);} 113 ierr = PetscOptionsReal("-ts_adapt_basic_safety","Safety factor relative to target error","",basic->safety,&basic->safety,NULL);CHKERRQ(ierr); 114 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); 115 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); 116 ierr = PetscOptionsTail();CHKERRQ(ierr); 117 PetscFunctionReturn(0); 118 } 119 120 static PetscErrorCode TSAdaptView_Basic(TSAdapt adapt,PetscViewer viewer) 121 { 122 TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 123 PetscErrorCode ierr; 124 PetscBool iascii; 125 126 PetscFunctionBegin; 127 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 128 if (iascii) { 129 if (basic->always_accept) {ierr = PetscViewerASCIIPrintf(viewer," Basic: always accepting steps\n");CHKERRQ(ierr);} 130 ierr = PetscViewerASCIIPrintf(viewer," Basic: clip fastest decrease %g, fastest increase %g\n",(double)basic->clip[0],(double)basic->clip[1]);CHKERRQ(ierr); 131 ierr = PetscViewerASCIIPrintf(viewer," Basic: safety factor %g, extra factor after step rejection %g\n",(double)basic->safety,(double)basic->reject_safety);CHKERRQ(ierr); 132 } 133 PetscFunctionReturn(0); 134 } 135 136 /*MC 137 TSADAPTBASIC - Basic adaptive controller for time stepping 138 139 Level: intermediate 140 141 .seealso: TS, TSAdapt, TSSetAdapt() 142 M*/ 143 PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt adapt) 144 { 145 PetscErrorCode ierr; 146 TSAdapt_Basic *a; 147 148 PetscFunctionBegin; 149 ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr); 150 adapt->data = (void*)a; 151 adapt->ops->choose = TSAdaptChoose_Basic; 152 adapt->ops->setfromoptions = TSAdaptSetFromOptions_Basic; 153 adapt->ops->destroy = TSAdaptDestroy_Basic; 154 adapt->ops->view = TSAdaptView_Basic; 155 156 a->clip[0] = 0.1; 157 a->clip[1] = 10.; 158 a->safety = 0.9; 159 a->reject_safety = 0.5; 160 a->always_accept = PETSC_FALSE; 161 PetscFunctionReturn(0); 162 } 163 164 /*@ 165 TSAdaptBasicSetClip - Sets the admissible decrease/increase factor in step size 166 167 Collective on TSAdapt 168 169 Input Arguments: 170 + adapt - adaptive controller context 171 . low - admissible decrease factor 172 + high - admissible increase factor 173 174 Level: intermediate 175 176 .seealso: TSAdaptChoose() 177 @*/ 178 PetscErrorCode TSAdaptBasicSetClip(TSAdapt adapt,PetscReal low,PetscReal high) 179 { 180 TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 181 PetscBool isbasic; 182 PetscErrorCode ierr; 183 184 PetscFunctionBegin; 185 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 186 if (low != PETSC_DEFAULT && low < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be non negative",(double)low); 187 if (low != PETSC_DEFAULT && low > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be less than one",(double)low); 188 if (high != PETSC_DEFAULT && high < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Increase factor %g must be geather than one",(double)high); 189 ierr = PetscObjectTypeCompare((PetscObject)adapt,TSADAPTBASIC,&isbasic);CHKERRQ(ierr); 190 if (isbasic) { 191 if (low != PETSC_DEFAULT) basic->clip[0] = low; 192 if (high != PETSC_DEFAULT) basic->clip[1] = high; 193 } 194 PetscFunctionReturn(0); 195 } 196 197 /*@ 198 TSAdaptBasicGetClip - Gets the admissible decrease/increase factor in step size 199 200 Collective on TSAdapt 201 202 Input Arguments: 203 . adapt - adaptive controller context 204 205 Ouput Arguments: 206 + low - optional, admissible decrease factor 207 - high - optional, admissible increase factor 208 209 Level: intermediate 210 211 .seealso: TSAdaptChoose() 212 @*/ 213 PetscErrorCode TSAdaptBasicGetClip(TSAdapt adapt,PetscReal *low,PetscReal *high) 214 { 215 TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 216 PetscBool isbasic; 217 PetscErrorCode ierr; 218 219 PetscFunctionBegin; 220 PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 221 if (low) PetscValidRealPointer(low,2); 222 if (high) PetscValidRealPointer(high,3); 223 ierr = PetscObjectTypeCompare((PetscObject)adapt,TSADAPTBASIC,&isbasic);CHKERRQ(ierr); 224 if (isbasic) { 225 if (low) *low = basic->clip[0]; 226 if (high) *high = basic->clip[1]; 227 } else { 228 if (low) *low = 0; 229 if (high) *high = PETSC_MAX_REAL; 230 } 231 PetscFunctionReturn(0); 232 } 233