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