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 #undef __FUNCT__ 12 #define __FUNCT__ "TSAdaptChoose_Basic" 13 static PetscErrorCode TSAdaptChoose_Basic(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte) 14 { 15 TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 16 PetscErrorCode ierr; 17 Vec X,Y; 18 PetscReal enorm,hfac_lte,h_lte,safety; 19 PetscInt order,stepno; 20 21 PetscFunctionBegin; 22 ierr = TSGetTimeStepNumber(ts,&stepno);CHKERRQ(ierr); 23 ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 24 if (!basic->Y) {ierr = VecDuplicate(X,&basic->Y);CHKERRQ(ierr);} 25 Y = basic->Y; 26 order = adapt->candidates.order[0]; 27 ierr = TSEvaluateStep(ts,order-1,Y,NULL);CHKERRQ(ierr); 28 29 safety = basic->safety; 30 ierr = TSErrorWeightedNorm(ts,X,Y,adapt->wnormtype,&enorm);CHKERRQ(ierr); 31 if (enorm > 1.) { 32 if (!*accept) safety *= basic->reject_safety; /* The last attempt also failed, shorten more aggressively */ 33 if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) { 34 ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting because step size %g is at minimum\n",(double)enorm,(double)h);CHKERRQ(ierr); 35 *accept = PETSC_TRUE; 36 } else if (basic->always_accept) { 37 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); 38 *accept = PETSC_TRUE; 39 } else { 40 ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, rejecting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr); 41 *accept = PETSC_FALSE; 42 } 43 } else { 44 ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr); 45 *accept = PETSC_TRUE; 46 } 47 48 /* The optimal new step based purely on local truncation error for this step. */ 49 if (enorm == 0.0) { 50 hfac_lte = safety * PETSC_INFINITY; 51 } else { 52 hfac_lte = safety * PetscPowReal(enorm,-1./order); 53 } 54 h_lte = h * PetscClipInterval(hfac_lte,basic->clip[0],basic->clip[1]); 55 56 *next_sc = 0; 57 *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max); 58 *wlte = enorm; 59 PetscFunctionReturn(0); 60 } 61 62 #undef __FUNCT__ 63 #define __FUNCT__ "TSAdaptReset_Basic" 64 static PetscErrorCode TSAdaptReset_Basic(TSAdapt adapt) 65 { 66 TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 67 PetscErrorCode ierr; 68 69 PetscFunctionBegin; 70 ierr = VecDestroy(&basic->Y);CHKERRQ(ierr); 71 PetscFunctionReturn(0); 72 } 73 74 #undef __FUNCT__ 75 #define __FUNCT__ "TSAdaptDestroy_Basic" 76 static PetscErrorCode TSAdaptDestroy_Basic(TSAdapt adapt) 77 { 78 PetscErrorCode ierr; 79 80 PetscFunctionBegin; 81 ierr = TSAdaptReset_Basic(adapt);CHKERRQ(ierr); 82 ierr = PetscFree(adapt->data);CHKERRQ(ierr); 83 PetscFunctionReturn(0); 84 } 85 86 #undef __FUNCT__ 87 #define __FUNCT__ "TSAdaptSetFromOptions_Basic" 88 static PetscErrorCode TSAdaptSetFromOptions_Basic(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 89 { 90 TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 91 PetscErrorCode ierr; 92 PetscInt two; 93 PetscBool set; 94 95 PetscFunctionBegin; 96 ierr = PetscOptionsHead(PetscOptionsObject,"Basic adaptive controller options");CHKERRQ(ierr); 97 two = 2; 98 ierr = PetscOptionsRealArray("-ts_adapt_basic_clip","Admissible decrease/increase in step size","",basic->clip,&two,&set);CHKERRQ(ierr); 99 if (set && (two != 2 || basic->clip[0] > basic->clip[1])) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_basic_clip"); 100 ierr = PetscOptionsReal("-ts_adapt_basic_safety","Safety factor relative to target error","",basic->safety,&basic->safety,NULL);CHKERRQ(ierr); 101 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); 102 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); 103 ierr = PetscOptionsTail();CHKERRQ(ierr); 104 PetscFunctionReturn(0); 105 } 106 107 #undef __FUNCT__ 108 #define __FUNCT__ "TSAdaptView_Basic" 109 static PetscErrorCode TSAdaptView_Basic(TSAdapt adapt,PetscViewer viewer) 110 { 111 TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 112 PetscErrorCode ierr; 113 PetscBool iascii; 114 115 PetscFunctionBegin; 116 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 117 if (iascii) { 118 if (basic->always_accept) {ierr = PetscViewerASCIIPrintf(viewer," Basic: always accepting steps\n");CHKERRQ(ierr);} 119 ierr = PetscViewerASCIIPrintf(viewer," Basic: clip fastest decrease %g, fastest increase %g\n",(double)basic->clip[0],(double)basic->clip[1]);CHKERRQ(ierr); 120 ierr = PetscViewerASCIIPrintf(viewer," Basic: safety factor %g, extra factor after step rejection %g\n",(double)basic->safety,(double)basic->reject_safety);CHKERRQ(ierr); 121 } 122 PetscFunctionReturn(0); 123 } 124 125 #undef __FUNCT__ 126 #define __FUNCT__ "TSAdaptCreate_Basic" 127 /*MC 128 TSADAPTBASIC - Basic adaptive controller for time stepping 129 130 Level: intermediate 131 132 .seealso: TS, TSAdapt, TSSetAdapt() 133 M*/ 134 PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt adapt) 135 { 136 PetscErrorCode ierr; 137 TSAdapt_Basic *a; 138 139 PetscFunctionBegin; 140 ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr); 141 adapt->data = (void*)a; 142 adapt->ops->choose = TSAdaptChoose_Basic; 143 adapt->ops->setfromoptions = TSAdaptSetFromOptions_Basic; 144 adapt->ops->destroy = TSAdaptDestroy_Basic; 145 adapt->ops->view = TSAdaptView_Basic; 146 147 a->clip[0] = 0.1; 148 a->clip[1] = 10.; 149 a->safety = 0.9; 150 a->reject_safety = 0.5; 151 a->always_accept = PETSC_FALSE; 152 PetscFunctionReturn(0); 153 } 154