1ae2316d0SEmil Constantinescu #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2ae2316d0SEmil Constantinescu 3ae2316d0SEmil Constantinescu typedef struct { 4ae2316d0SEmil Constantinescu PetscReal clip[2]; /* admissible decrease/increase factors */ 5ae2316d0SEmil Constantinescu PetscReal safety; /* safety factor relative to target error */ 6ae2316d0SEmil Constantinescu PetscReal reject_safety; /* extra safety factor if the last step was rejected */ 7657c1e31SEmil Constantinescu Vec X; 8ae2316d0SEmil Constantinescu Vec Y; 98a175baeSEmil Constantinescu Vec E; 10726095e4SEmil Constantinescu } TSAdapt_GLEE; 11ae2316d0SEmil Constantinescu 125e4ed32fSEmil 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) 13ae2316d0SEmil Constantinescu { 14657c1e31SEmil Constantinescu TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 150a01e1b2SEmil Constantinescu TSType time_scheme; /* Type of time-integration scheme */ 16ae2316d0SEmil Constantinescu PetscErrorCode ierr; 178a175baeSEmil Constantinescu Vec X,Y,E; 18df3bac28SEmil Constantinescu PetscReal enorm,enorma,enormr,hfac_lte,hfac_ltea,hfac_lter,h_lte,safety; 195e4ed32fSEmil Constantinescu PetscInt order,stepno; 20df3bac28SEmil Constantinescu PetscBool bGTEMethod=PETSC_FALSE; 21ae2316d0SEmil Constantinescu 22ae2316d0SEmil Constantinescu PetscFunctionBegin; 23df3bac28SEmil Constantinescu 245e4ed32fSEmil Constantinescu *next_sc = 0; /* Reuse the same order scheme */ 25657c1e31SEmil Constantinescu safety = glee->safety; 26df3bac28SEmil Constantinescu ierr = TSGetTimeStepNumber(ts,&stepno);CHKERRQ(ierr); 270a01e1b2SEmil Constantinescu ierr = TSGetType(ts,&time_scheme);CHKERRQ(ierr); 28df3bac28SEmil Constantinescu if (!strcmp(time_scheme,TSGLEE)) bGTEMethod=PETSC_TRUE; 29df3bac28SEmil Constantinescu order = adapt->candidates.order[0]; 30df3bac28SEmil Constantinescu 31df3bac28SEmil Constantinescu if (bGTEMethod){/* the method is of GLEE type */ 328a175baeSEmil Constantinescu ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 338a175baeSEmil Constantinescu if (!glee->E) {ierr = VecDuplicate(X,&glee->E);CHKERRQ(ierr);} 348a175baeSEmil Constantinescu E = glee->E; 358a175baeSEmil Constantinescu ierr = TSGetTimeError(ts,0,&E);CHKERRQ(ierr); 36df3bac28SEmil Constantinescu /* this should be called with Y (the solution at the beginning of the step)*/ 378a175baeSEmil Constantinescu ierr = TSErrorWeightedENorm(ts,E,X,X,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr); 380a01e1b2SEmil Constantinescu } else { 39df3bac28SEmil Constantinescu /* the method is NOT of GLEE type; use the stantard basic augmented by separate atol and rtol */ 40ae2316d0SEmil Constantinescu ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 41657c1e31SEmil Constantinescu if (!glee->Y) {ierr = VecDuplicate(X,&glee->Y);CHKERRQ(ierr);} 42657c1e31SEmil Constantinescu Y = glee->Y; 43ae2316d0SEmil Constantinescu ierr = TSEvaluateStep(ts,order-1,Y,NULL);CHKERRQ(ierr); 447453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr); 450a01e1b2SEmil Constantinescu } 465e4ed32fSEmil Constantinescu 475e4ed32fSEmil Constantinescu if (enorm < 0) { 485e4ed32fSEmil Constantinescu *accept = PETSC_TRUE; 495e4ed32fSEmil Constantinescu *next_h = h; /* Reuse the old step */ 50df3bac28SEmil Constantinescu *wlte = -1; /* Weighted error was not evaluated */ 51df3bac28SEmil Constantinescu *wltea = -1; /* Weighted absolute error was not evaluated */ 52df3bac28SEmil Constantinescu *wlter = -1; /* Weighted relative error was not evaluated */ 535e4ed32fSEmil Constantinescu PetscFunctionReturn(0); 545e4ed32fSEmil Constantinescu } 555e4ed32fSEmil Constantinescu 56df3bac28SEmil Constantinescu if (enorm > 1. || enorma > 1. || enormr > 1.) { 57657c1e31SEmil Constantinescu if (!*accept) safety *= glee->reject_safety; /* The last attempt also failed, shorten more aggressively */ 58ae2316d0SEmil Constantinescu if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) { 59df3bac28SEmil 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); 60ae2316d0SEmil Constantinescu *accept = PETSC_TRUE; 61*bf997491SLisandro Dalcin } else if (adapt->always_accept) { 62df3bac28SEmil 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); 63ae2316d0SEmil Constantinescu *accept = PETSC_TRUE; 64ae2316d0SEmil Constantinescu } else { 65df3bac28SEmil 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); 66ae2316d0SEmil Constantinescu *accept = PETSC_FALSE; 67ae2316d0SEmil Constantinescu } 68ae2316d0SEmil Constantinescu } else { 69df3bac28SEmil 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); 70ae2316d0SEmil Constantinescu *accept = PETSC_TRUE; 71ae2316d0SEmil Constantinescu } 72ae2316d0SEmil Constantinescu 73df3bac28SEmil Constantinescu if (bGTEMethod){ 74df3bac28SEmil Constantinescu /* The optimal new step based on the current global truncation error. */ 75df3bac28SEmil Constantinescu if (enorm > 0) { 76df3bac28SEmil Constantinescu /* factor based on the absolute tolerance */ 77df3bac28SEmil Constantinescu hfac_ltea = safety * PetscPowReal(1./enorma,((PetscReal)1)/order); 78df3bac28SEmil Constantinescu /* factor based on the relative tolerance */ 79df3bac28SEmil Constantinescu hfac_lter = safety * PetscPowReal(1./enormr,((PetscReal)1)/order); 80df3bac28SEmil Constantinescu /* pick the minimum time step among the relative and absolute tolerances */ 81df3bac28SEmil Constantinescu hfac_lte = PetscMin(hfac_ltea,hfac_lter); 82df3bac28SEmil Constantinescu } else { 83ae2316d0SEmil Constantinescu hfac_lte = safety * PETSC_INFINITY; 84df3bac28SEmil Constantinescu } 85657c1e31SEmil Constantinescu h_lte = h * PetscClipInterval(hfac_lte,glee->clip[0],glee->clip[1]); 86ae2316d0SEmil Constantinescu *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max); 87df3bac28SEmil Constantinescu } else { 88df3bac28SEmil Constantinescu /* The optimal new step based purely on local truncation error for this step. */ 89df3bac28SEmil Constantinescu if (enorm > 0) { 90df3bac28SEmil Constantinescu /* factor based on the absolute tolerance */ 91df3bac28SEmil Constantinescu hfac_ltea = safety * PetscPowReal(enorma,((PetscReal)-1)/order); 92df3bac28SEmil Constantinescu /* factor based on the relative tolerance */ 93df3bac28SEmil Constantinescu hfac_lter = safety * PetscPowReal(enormr,((PetscReal)-1)/order); 94df3bac28SEmil Constantinescu /* pick the minimum time step among the relative and absolute tolerances */ 95df3bac28SEmil Constantinescu hfac_lte = PetscMin(hfac_ltea,hfac_lter); 96df3bac28SEmil Constantinescu } else { 97df3bac28SEmil Constantinescu hfac_lte = safety * PETSC_INFINITY; 98df3bac28SEmil Constantinescu } 99df3bac28SEmil Constantinescu h_lte = h * PetscClipInterval(hfac_lte,glee->clip[0],glee->clip[1]); 100df3bac28SEmil Constantinescu *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max); 101df3bac28SEmil Constantinescu } 102ae2316d0SEmil Constantinescu *wlte = enorm; 1035e4ed32fSEmil Constantinescu *wltea = enorma; 1045e4ed32fSEmil Constantinescu *wlter = enormr; 105ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 106ae2316d0SEmil Constantinescu } 107ae2316d0SEmil Constantinescu 108726095e4SEmil Constantinescu static PetscErrorCode TSAdaptReset_GLEE(TSAdapt adapt) 109ae2316d0SEmil Constantinescu { 110657c1e31SEmil Constantinescu TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 111ae2316d0SEmil Constantinescu PetscErrorCode ierr; 112ae2316d0SEmil Constantinescu 113ae2316d0SEmil Constantinescu PetscFunctionBegin; 114657c1e31SEmil Constantinescu ierr = VecDestroy(&glee->Y);CHKERRQ(ierr); 115ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 116ae2316d0SEmil Constantinescu } 117ae2316d0SEmil Constantinescu 118726095e4SEmil Constantinescu static PetscErrorCode TSAdaptDestroy_GLEE(TSAdapt adapt) 119ae2316d0SEmil Constantinescu { 120ae2316d0SEmil Constantinescu PetscErrorCode ierr; 121ae2316d0SEmil Constantinescu 122ae2316d0SEmil Constantinescu PetscFunctionBegin; 123726095e4SEmil Constantinescu ierr = TSAdaptReset_GLEE(adapt);CHKERRQ(ierr); 124ae2316d0SEmil Constantinescu ierr = PetscFree(adapt->data);CHKERRQ(ierr); 125ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 126ae2316d0SEmil Constantinescu } 127ae2316d0SEmil Constantinescu 128726095e4SEmil Constantinescu static PetscErrorCode TSAdaptSetFromOptions_GLEE(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 129ae2316d0SEmil Constantinescu { 130657c1e31SEmil Constantinescu TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 131ae2316d0SEmil Constantinescu PetscErrorCode ierr; 132ae2316d0SEmil Constantinescu PetscInt two; 133ae2316d0SEmil Constantinescu PetscBool set; 134ae2316d0SEmil Constantinescu 135ae2316d0SEmil Constantinescu PetscFunctionBegin; 136726095e4SEmil Constantinescu ierr = PetscOptionsHead(PetscOptionsObject,"GLEE adaptive controller options");CHKERRQ(ierr); 137ae2316d0SEmil Constantinescu two = 2; 138657c1e31SEmil Constantinescu ierr = PetscOptionsRealArray("-ts_adapt_glee_clip","Admissible decrease/increase in step size","",glee->clip,&two,&set);CHKERRQ(ierr); 139657c1e31SEmil 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"); 140657c1e31SEmil Constantinescu ierr = PetscOptionsReal("-ts_adapt_glee_safety","Safety factor relative to target error","",glee->safety,&glee->safety,NULL);CHKERRQ(ierr); 141657c1e31SEmil 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); 142ae2316d0SEmil Constantinescu ierr = PetscOptionsTail();CHKERRQ(ierr); 143ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 144ae2316d0SEmil Constantinescu } 145ae2316d0SEmil Constantinescu 146726095e4SEmil Constantinescu static PetscErrorCode TSAdaptView_GLEE(TSAdapt adapt,PetscViewer viewer) 147ae2316d0SEmil Constantinescu { 148657c1e31SEmil Constantinescu TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 149ae2316d0SEmil Constantinescu PetscErrorCode ierr; 150ae2316d0SEmil Constantinescu PetscBool iascii; 151ae2316d0SEmil Constantinescu 152ae2316d0SEmil Constantinescu PetscFunctionBegin; 153ae2316d0SEmil Constantinescu ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 154ae2316d0SEmil Constantinescu if (iascii) { 155657c1e31SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," GLEE: clip fastest decrease %g, fastest increase %g\n",(double)glee->clip[0],(double)glee->clip[1]);CHKERRQ(ierr); 156657c1e31SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," GLEE: safety factor %g, extra factor after step rejection %g\n",(double)glee->safety,(double)glee->reject_safety);CHKERRQ(ierr); 157ae2316d0SEmil Constantinescu } 158ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 159ae2316d0SEmil Constantinescu } 160ae2316d0SEmil Constantinescu 161ae2316d0SEmil Constantinescu /*MC 162657c1e31SEmil Constantinescu TSADAPTGLEE - GLEE adaptive controller for time stepping 163ae2316d0SEmil Constantinescu 164ae2316d0SEmil Constantinescu Level: intermediate 165ae2316d0SEmil Constantinescu 166ae2316d0SEmil Constantinescu .seealso: TS, TSAdapt, TSSetAdapt() 167ae2316d0SEmil Constantinescu M*/ 168726095e4SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt adapt) 169ae2316d0SEmil Constantinescu { 170ae2316d0SEmil Constantinescu PetscErrorCode ierr; 171726095e4SEmil Constantinescu TSAdapt_GLEE *a; 172ae2316d0SEmil Constantinescu 173ae2316d0SEmil Constantinescu PetscFunctionBegin; 174ae2316d0SEmil Constantinescu ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr); 175ae2316d0SEmil Constantinescu adapt->data = (void*)a; 176726095e4SEmil Constantinescu adapt->ops->choose = TSAdaptChoose_GLEE; 177726095e4SEmil Constantinescu adapt->ops->setfromoptions = TSAdaptSetFromOptions_GLEE; 178726095e4SEmil Constantinescu adapt->ops->destroy = TSAdaptDestroy_GLEE; 179726095e4SEmil Constantinescu adapt->ops->view = TSAdaptView_GLEE; 180ae2316d0SEmil Constantinescu 181ae2316d0SEmil Constantinescu a->clip[0] = 0.1; 182ae2316d0SEmil Constantinescu a->clip[1] = 10.; 1835e4ed32fSEmil Constantinescu a->safety = 0.9; 184ae2316d0SEmil Constantinescu a->reject_safety = 0.5; 185ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 186ae2316d0SEmil Constantinescu } 187