1ae2316d0SEmil Constantinescu #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2ae2316d0SEmil Constantinescu 3ae2316d0SEmil Constantinescu typedef struct { 4ae2316d0SEmil Constantinescu PetscBool always_accept; 5ae2316d0SEmil Constantinescu PetscReal clip[2]; /* admissible decrease/increase factors */ 6ae2316d0SEmil Constantinescu PetscReal safety; /* safety factor relative to target error */ 7ae2316d0SEmil Constantinescu PetscReal reject_safety; /* extra safety factor if the last step was rejected */ 8657c1e31SEmil Constantinescu Vec X; 9ae2316d0SEmil Constantinescu Vec Y; 108a175baeSEmil Constantinescu Vec E; 11726095e4SEmil Constantinescu } TSAdapt_GLEE; 12ae2316d0SEmil Constantinescu 13ae2316d0SEmil Constantinescu #undef __FUNCT__ 14726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptChoose_GLEE" 155e4ed32fSEmil Constantinescu static PetscErrorCode TSAdaptChoose_GLEE(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte,PetscReal *wltea,PetscReal *wlter) 16ae2316d0SEmil Constantinescu { 17657c1e31SEmil Constantinescu TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 180a01e1b2SEmil Constantinescu TSType time_scheme; /* Type of time-integration scheme */ 19ae2316d0SEmil Constantinescu PetscErrorCode ierr; 208a175baeSEmil Constantinescu Vec X,Y,E; 21*df3bac28SEmil Constantinescu PetscReal enorm,enorma,enormr,hfac_lte,hfac_ltea,hfac_lter,h_lte,safety; 225e4ed32fSEmil Constantinescu PetscInt order,stepno; 23*df3bac28SEmil Constantinescu PetscBool bGTEMethod=PETSC_FALSE; 24ae2316d0SEmil Constantinescu 25ae2316d0SEmil Constantinescu PetscFunctionBegin; 26*df3bac28SEmil Constantinescu 275e4ed32fSEmil Constantinescu *next_sc = 0; /* Reuse the same order scheme */ 28657c1e31SEmil Constantinescu safety = glee->safety; 29*df3bac28SEmil Constantinescu ierr = TSGetTimeStepNumber(ts,&stepno);CHKERRQ(ierr); 300a01e1b2SEmil Constantinescu ierr = TSGetType(ts,&time_scheme);CHKERRQ(ierr); 31*df3bac28SEmil Constantinescu if (!strcmp(time_scheme,TSGLEE)) bGTEMethod=PETSC_TRUE; 32*df3bac28SEmil Constantinescu order = adapt->candidates.order[0]; 33*df3bac28SEmil Constantinescu 34*df3bac28SEmil Constantinescu if (bGTEMethod){/* the method is of GLEE type */ 358a175baeSEmil Constantinescu ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 368a175baeSEmil Constantinescu if (!glee->E) {ierr = VecDuplicate(X,&glee->E);CHKERRQ(ierr);} 378a175baeSEmil Constantinescu E = glee->E; 388a175baeSEmil Constantinescu ierr = TSGetTimeError(ts,0,&E);CHKERRQ(ierr); 39*df3bac28SEmil Constantinescu /* this should be called with Y (the solution at the beginning of the step)*/ 408a175baeSEmil Constantinescu ierr = TSErrorWeightedENorm(ts,E,X,X,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr); 410a01e1b2SEmil Constantinescu } else { 42*df3bac28SEmil Constantinescu /* the method is NOT of GLEE type; use the stantard basic augmented by separate atol and rtol */ 43ae2316d0SEmil Constantinescu ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 44657c1e31SEmil Constantinescu if (!glee->Y) {ierr = VecDuplicate(X,&glee->Y);CHKERRQ(ierr);} 45657c1e31SEmil Constantinescu Y = glee->Y; 46ae2316d0SEmil Constantinescu ierr = TSEvaluateStep(ts,order-1,Y,NULL);CHKERRQ(ierr); 477453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr); 480a01e1b2SEmil Constantinescu } 495e4ed32fSEmil Constantinescu 505e4ed32fSEmil Constantinescu if (enorm < 0) { 515e4ed32fSEmil Constantinescu *accept = PETSC_TRUE; 525e4ed32fSEmil Constantinescu *next_h = h; /* Reuse the old step */ 53*df3bac28SEmil Constantinescu *wlte = -1; /* Weighted error was not evaluated */ 54*df3bac28SEmil Constantinescu *wltea = -1; /* Weighted absolute error was not evaluated */ 55*df3bac28SEmil Constantinescu *wlter = -1; /* Weighted relative error was not evaluated */ 565e4ed32fSEmil Constantinescu PetscFunctionReturn(0); 575e4ed32fSEmil Constantinescu } 585e4ed32fSEmil Constantinescu 59*df3bac28SEmil Constantinescu if (enorm > 1. || enorma > 1. || enormr > 1.) { 60657c1e31SEmil Constantinescu if (!*accept) safety *= glee->reject_safety; /* The last attempt also failed, shorten more aggressively */ 61ae2316d0SEmil Constantinescu if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) { 62*df3bac28SEmil Constantinescu ierr = PetscInfo4(adapt,"Estimated scaled truncation error [combined, absolute, relative]] [%g, %g, %g], accepting because step size %g is at minimum\n",(double)enorm,(double)enorma,(double)enormr,(double)h);CHKERRQ(ierr); 63ae2316d0SEmil Constantinescu *accept = PETSC_TRUE; 64657c1e31SEmil Constantinescu } else if (glee->always_accept) { 65*df3bac28SEmil Constantinescu ierr = PetscInfo4(adapt,"Estimated scaled truncation error [combined, absolute, relative]] [%g, %g, %g], accepting step of size %g because always_accept is set\n",(double)enorm,(double)enorma,(double)enormr,(double)h);CHKERRQ(ierr); 66ae2316d0SEmil Constantinescu *accept = PETSC_TRUE; 67ae2316d0SEmil Constantinescu } else { 68*df3bac28SEmil Constantinescu ierr = PetscInfo4(adapt,"Estimated scaled truncation error [combined, absolute, relative]] [%g, %g, %g], rejecting step of size %g\n",(double)enorm,(double)enorma,(double)enormr,(double)h);CHKERRQ(ierr); 69ae2316d0SEmil Constantinescu *accept = PETSC_FALSE; 70ae2316d0SEmil Constantinescu } 71ae2316d0SEmil Constantinescu } else { 72*df3bac28SEmil Constantinescu ierr = PetscInfo4(adapt,"Estimated scaled truncation error [combined, absolute, relative] [%g, %g, %g], accepting step of size %g\n",(double)enorm,(double)enorma,(double)enormr,(double)h);CHKERRQ(ierr); 73ae2316d0SEmil Constantinescu *accept = PETSC_TRUE; 74ae2316d0SEmil Constantinescu } 75ae2316d0SEmil Constantinescu 76*df3bac28SEmil Constantinescu if (bGTEMethod){ 77*df3bac28SEmil Constantinescu /* The optimal new step based on the current global truncation error. */ 78*df3bac28SEmil Constantinescu if (enorm > 0) { 79*df3bac28SEmil Constantinescu /* factor based on the absolute tolerance */ 80*df3bac28SEmil Constantinescu hfac_ltea = safety * PetscPowReal(1./enorma,((PetscReal)1)/order); 81*df3bac28SEmil Constantinescu /* factor based on the relative tolerance */ 82*df3bac28SEmil Constantinescu hfac_lter = safety * PetscPowReal(1./enormr,((PetscReal)1)/order); 83*df3bac28SEmil Constantinescu /* pick the minimum time step among the relative and absolute tolerances */ 84*df3bac28SEmil Constantinescu hfac_lte = PetscMin(hfac_ltea,hfac_lter); 85*df3bac28SEmil Constantinescu } else { 86ae2316d0SEmil Constantinescu hfac_lte = safety * PETSC_INFINITY; 87*df3bac28SEmil Constantinescu } 88657c1e31SEmil Constantinescu h_lte = h * PetscClipInterval(hfac_lte,glee->clip[0],glee->clip[1]); 89ae2316d0SEmil Constantinescu *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max); 90*df3bac28SEmil Constantinescu } else { 91*df3bac28SEmil Constantinescu /* The optimal new step based purely on local truncation error for this step. */ 92*df3bac28SEmil Constantinescu if (enorm > 0) { 93*df3bac28SEmil Constantinescu /* factor based on the absolute tolerance */ 94*df3bac28SEmil Constantinescu hfac_ltea = safety * PetscPowReal(enorma,((PetscReal)-1)/order); 95*df3bac28SEmil Constantinescu /* factor based on the relative tolerance */ 96*df3bac28SEmil Constantinescu hfac_lter = safety * PetscPowReal(enormr,((PetscReal)-1)/order); 97*df3bac28SEmil Constantinescu /* pick the minimum time step among the relative and absolute tolerances */ 98*df3bac28SEmil Constantinescu hfac_lte = PetscMin(hfac_ltea,hfac_lter); 99*df3bac28SEmil Constantinescu } else { 100*df3bac28SEmil Constantinescu hfac_lte = safety * PETSC_INFINITY; 101*df3bac28SEmil Constantinescu } 102*df3bac28SEmil Constantinescu h_lte = h * PetscClipInterval(hfac_lte,glee->clip[0],glee->clip[1]); 103*df3bac28SEmil Constantinescu *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max); 104*df3bac28SEmil Constantinescu } 105ae2316d0SEmil Constantinescu *wlte = enorm; 1065e4ed32fSEmil Constantinescu *wltea = enorma; 1075e4ed32fSEmil Constantinescu *wlter = enormr; 108ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 109ae2316d0SEmil Constantinescu } 110ae2316d0SEmil Constantinescu 111ae2316d0SEmil Constantinescu #undef __FUNCT__ 112726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptReset_GLEE" 113726095e4SEmil Constantinescu static PetscErrorCode TSAdaptReset_GLEE(TSAdapt adapt) 114ae2316d0SEmil Constantinescu { 115657c1e31SEmil Constantinescu TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 116ae2316d0SEmil Constantinescu PetscErrorCode ierr; 117ae2316d0SEmil Constantinescu 118ae2316d0SEmil Constantinescu PetscFunctionBegin; 119657c1e31SEmil Constantinescu ierr = VecDestroy(&glee->Y);CHKERRQ(ierr); 120ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 121ae2316d0SEmil Constantinescu } 122ae2316d0SEmil Constantinescu 123ae2316d0SEmil Constantinescu #undef __FUNCT__ 124726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptDestroy_GLEE" 125726095e4SEmil Constantinescu static PetscErrorCode TSAdaptDestroy_GLEE(TSAdapt adapt) 126ae2316d0SEmil Constantinescu { 127ae2316d0SEmil Constantinescu PetscErrorCode ierr; 128ae2316d0SEmil Constantinescu 129ae2316d0SEmil Constantinescu PetscFunctionBegin; 130726095e4SEmil Constantinescu ierr = TSAdaptReset_GLEE(adapt);CHKERRQ(ierr); 131ae2316d0SEmil Constantinescu ierr = PetscFree(adapt->data);CHKERRQ(ierr); 132ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 133ae2316d0SEmil Constantinescu } 134ae2316d0SEmil Constantinescu 135ae2316d0SEmil Constantinescu #undef __FUNCT__ 136726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptSetFromOptions_GLEE" 137726095e4SEmil Constantinescu static PetscErrorCode TSAdaptSetFromOptions_GLEE(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 138ae2316d0SEmil Constantinescu { 139657c1e31SEmil Constantinescu TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 140ae2316d0SEmil Constantinescu PetscErrorCode ierr; 141ae2316d0SEmil Constantinescu PetscInt two; 142ae2316d0SEmil Constantinescu PetscBool set; 143ae2316d0SEmil Constantinescu 144ae2316d0SEmil Constantinescu PetscFunctionBegin; 145726095e4SEmil Constantinescu ierr = PetscOptionsHead(PetscOptionsObject,"GLEE adaptive controller options");CHKERRQ(ierr); 146ae2316d0SEmil Constantinescu two = 2; 147657c1e31SEmil Constantinescu ierr = PetscOptionsRealArray("-ts_adapt_glee_clip","Admissible decrease/increase in step size","",glee->clip,&two,&set);CHKERRQ(ierr); 148657c1e31SEmil Constantinescu if (set && (two != 2 || glee->clip[0] > glee->clip[1])) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_glee_clip"); 149657c1e31SEmil Constantinescu ierr = PetscOptionsReal("-ts_adapt_glee_safety","Safety factor relative to target error","",glee->safety,&glee->safety,NULL);CHKERRQ(ierr); 150657c1e31SEmil Constantinescu ierr = PetscOptionsReal("-ts_adapt_glee_reject_safety","Extra safety factor to apply if the last step was rejected","",glee->reject_safety,&glee->reject_safety,NULL);CHKERRQ(ierr); 151657c1e31SEmil Constantinescu ierr = PetscOptionsBool("-ts_adapt_glee_always_accept","Always accept the step regardless of whether local truncation error meets goal","",glee->always_accept,&glee->always_accept,NULL);CHKERRQ(ierr); 152ae2316d0SEmil Constantinescu ierr = PetscOptionsTail();CHKERRQ(ierr); 153ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 154ae2316d0SEmil Constantinescu } 155ae2316d0SEmil Constantinescu 156ae2316d0SEmil Constantinescu #undef __FUNCT__ 157726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptView_GLEE" 158726095e4SEmil Constantinescu static PetscErrorCode TSAdaptView_GLEE(TSAdapt adapt,PetscViewer viewer) 159ae2316d0SEmil Constantinescu { 160657c1e31SEmil Constantinescu TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 161ae2316d0SEmil Constantinescu PetscErrorCode ierr; 162ae2316d0SEmil Constantinescu PetscBool iascii; 163ae2316d0SEmil Constantinescu 164ae2316d0SEmil Constantinescu PetscFunctionBegin; 165ae2316d0SEmil Constantinescu ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 166ae2316d0SEmil Constantinescu if (iascii) { 167657c1e31SEmil Constantinescu if (glee->always_accept) {ierr = PetscViewerASCIIPrintf(viewer," GLEE: always accepting steps\n");CHKERRQ(ierr);} 168657c1e31SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," GLEE: clip fastest decrease %g, fastest increase %g\n",(double)glee->clip[0],(double)glee->clip[1]);CHKERRQ(ierr); 169657c1e31SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," GLEE: safety factor %g, extra factor after step rejection %g\n",(double)glee->safety,(double)glee->reject_safety);CHKERRQ(ierr); 170ae2316d0SEmil Constantinescu } 171ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 172ae2316d0SEmil Constantinescu } 173ae2316d0SEmil Constantinescu 174ae2316d0SEmil Constantinescu #undef __FUNCT__ 175726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptCreate_GLEE" 176ae2316d0SEmil Constantinescu /*MC 177657c1e31SEmil Constantinescu TSADAPTGLEE - GLEE adaptive controller for time stepping 178ae2316d0SEmil Constantinescu 179ae2316d0SEmil Constantinescu Level: intermediate 180ae2316d0SEmil Constantinescu 181ae2316d0SEmil Constantinescu .seealso: TS, TSAdapt, TSSetAdapt() 182ae2316d0SEmil Constantinescu M*/ 183726095e4SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt adapt) 184ae2316d0SEmil Constantinescu { 185ae2316d0SEmil Constantinescu PetscErrorCode ierr; 186726095e4SEmil Constantinescu TSAdapt_GLEE *a; 187ae2316d0SEmil Constantinescu 188ae2316d0SEmil Constantinescu PetscFunctionBegin; 189ae2316d0SEmil Constantinescu ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr); 190ae2316d0SEmil Constantinescu adapt->data = (void*)a; 191726095e4SEmil Constantinescu adapt->ops->choose = TSAdaptChoose_GLEE; 192726095e4SEmil Constantinescu adapt->ops->setfromoptions = TSAdaptSetFromOptions_GLEE; 193726095e4SEmil Constantinescu adapt->ops->destroy = TSAdaptDestroy_GLEE; 194726095e4SEmil Constantinescu adapt->ops->view = TSAdaptView_GLEE; 195ae2316d0SEmil Constantinescu 196ae2316d0SEmil Constantinescu a->clip[0] = 0.1; 197ae2316d0SEmil Constantinescu a->clip[1] = 10.; 1985e4ed32fSEmil Constantinescu a->safety = 0.9; 199ae2316d0SEmil Constantinescu a->reject_safety = 0.5; 200ae2316d0SEmil Constantinescu a->always_accept = PETSC_FALSE; 201ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 202ae2316d0SEmil Constantinescu } 203