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