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" 15*5e4ed32fSEmil 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; 217453f775SEmil Constantinescu PetscReal enorm,enorma,enormr,hfac_lte,h_lte,safety; 22*5e4ed32fSEmil Constantinescu PetscInt order,stepno; 23ae2316d0SEmil Constantinescu 24ae2316d0SEmil Constantinescu PetscFunctionBegin; 25*5e4ed32fSEmil Constantinescu *next_sc = 0; /* Reuse the same order scheme */ 26*5e4ed32fSEmil Constantinescu 27ae2316d0SEmil Constantinescu ierr = TSGetTimeStepNumber(ts,&stepno);CHKERRQ(ierr); 280a01e1b2SEmil Constantinescu 29657c1e31SEmil Constantinescu safety = glee->safety; 300a01e1b2SEmil Constantinescu ierr = TSGetType(ts,&time_scheme);CHKERRQ(ierr); 310a01e1b2SEmil Constantinescu if (!strcmp(time_scheme,TSGLEE)){ 320a01e1b2SEmil Constantinescu /* the method is of GLEE type */ 338a175baeSEmil Constantinescu ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 348a175baeSEmil Constantinescu /* ierr = TSGetPreviousSolution(ts,&Y);CHKERRQ(ierr); 358a175baeSEmil Constantinescu */ 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); 398a175baeSEmil Constantinescu /* this should be called with Y */ 408a175baeSEmil Constantinescu ierr = TSErrorWeightedENorm(ts,E,X,X,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr); 410a01e1b2SEmil Constantinescu } else { 420a01e1b2SEmil Constantinescu /* the method is NOT of GLEE type */ 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 order = adapt->candidates.order[0]; 47ae2316d0SEmil Constantinescu ierr = TSEvaluateStep(ts,order-1,Y,NULL);CHKERRQ(ierr); 487453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr); 490a01e1b2SEmil Constantinescu } 50*5e4ed32fSEmil Constantinescu 51*5e4ed32fSEmil Constantinescu if (enorm < 0) { 52*5e4ed32fSEmil Constantinescu *accept = PETSC_TRUE; 53*5e4ed32fSEmil Constantinescu *next_h = h; /* Reuse the old step */ 54*5e4ed32fSEmil Constantinescu *wlte = -1; /* Weighted local truncation error was not evaluated */ 55*5e4ed32fSEmil Constantinescu *wltea = -1; /* Weighted local truncation error was not evaluated */ 56*5e4ed32fSEmil Constantinescu *wlter = -1; /* Weighted local truncation error was not evaluated */ 57*5e4ed32fSEmil Constantinescu PetscFunctionReturn(0); 58*5e4ed32fSEmil Constantinescu } 59*5e4ed32fSEmil Constantinescu 60ae2316d0SEmil Constantinescu if (enorm > 1.) { 61657c1e31SEmil Constantinescu if (!*accept) safety *= glee->reject_safety; /* The last attempt also failed, shorten more aggressively */ 62ae2316d0SEmil Constantinescu if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) { 63ae2316d0SEmil 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); 64ae2316d0SEmil Constantinescu *accept = PETSC_TRUE; 65657c1e31SEmil Constantinescu } else if (glee->always_accept) { 66ae2316d0SEmil 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); 67ae2316d0SEmil Constantinescu *accept = PETSC_TRUE; 68ae2316d0SEmil Constantinescu } else { 69ae2316d0SEmil Constantinescu ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, rejecting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr); 70ae2316d0SEmil Constantinescu *accept = PETSC_FALSE; 71ae2316d0SEmil Constantinescu } 72ae2316d0SEmil Constantinescu } else { 73ae2316d0SEmil Constantinescu ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr); 74ae2316d0SEmil Constantinescu *accept = PETSC_TRUE; 75ae2316d0SEmil Constantinescu } 76ae2316d0SEmil Constantinescu 77ae2316d0SEmil Constantinescu /* The optimal new step based purely on local truncation error for this step. */ 78*5e4ed32fSEmil Constantinescu if (enorm > 0) 79*5e4ed32fSEmil Constantinescu hfac_lte = safety * PetscPowReal(enorm,((PetscReal)-1)/order); 80*5e4ed32fSEmil Constantinescu else 81ae2316d0SEmil Constantinescu hfac_lte = safety * PETSC_INFINITY; 82657c1e31SEmil Constantinescu h_lte = h * PetscClipInterval(hfac_lte,glee->clip[0],glee->clip[1]); 83ae2316d0SEmil Constantinescu 84ae2316d0SEmil Constantinescu *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max); 85ae2316d0SEmil Constantinescu *wlte = enorm; 86*5e4ed32fSEmil Constantinescu *wltea = enorma; 87*5e4ed32fSEmil Constantinescu *wlter = enormr; 88ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 89ae2316d0SEmil Constantinescu } 90ae2316d0SEmil Constantinescu 91ae2316d0SEmil Constantinescu #undef __FUNCT__ 92726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptReset_GLEE" 93726095e4SEmil Constantinescu static PetscErrorCode TSAdaptReset_GLEE(TSAdapt adapt) 94ae2316d0SEmil Constantinescu { 95657c1e31SEmil Constantinescu TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 96ae2316d0SEmil Constantinescu PetscErrorCode ierr; 97ae2316d0SEmil Constantinescu 98ae2316d0SEmil Constantinescu PetscFunctionBegin; 99657c1e31SEmil Constantinescu ierr = VecDestroy(&glee->Y);CHKERRQ(ierr); 100ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 101ae2316d0SEmil Constantinescu } 102ae2316d0SEmil Constantinescu 103ae2316d0SEmil Constantinescu #undef __FUNCT__ 104726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptDestroy_GLEE" 105726095e4SEmil Constantinescu static PetscErrorCode TSAdaptDestroy_GLEE(TSAdapt adapt) 106ae2316d0SEmil Constantinescu { 107ae2316d0SEmil Constantinescu PetscErrorCode ierr; 108ae2316d0SEmil Constantinescu 109ae2316d0SEmil Constantinescu PetscFunctionBegin; 110726095e4SEmil Constantinescu ierr = TSAdaptReset_GLEE(adapt);CHKERRQ(ierr); 111ae2316d0SEmil Constantinescu ierr = PetscFree(adapt->data);CHKERRQ(ierr); 112ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 113ae2316d0SEmil Constantinescu } 114ae2316d0SEmil Constantinescu 115ae2316d0SEmil Constantinescu #undef __FUNCT__ 116726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptSetFromOptions_GLEE" 117726095e4SEmil Constantinescu static PetscErrorCode TSAdaptSetFromOptions_GLEE(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 118ae2316d0SEmil Constantinescu { 119657c1e31SEmil Constantinescu TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 120ae2316d0SEmil Constantinescu PetscErrorCode ierr; 121ae2316d0SEmil Constantinescu PetscInt two; 122ae2316d0SEmil Constantinescu PetscBool set; 123ae2316d0SEmil Constantinescu 124ae2316d0SEmil Constantinescu PetscFunctionBegin; 125726095e4SEmil Constantinescu ierr = PetscOptionsHead(PetscOptionsObject,"GLEE adaptive controller options");CHKERRQ(ierr); 126ae2316d0SEmil Constantinescu two = 2; 127657c1e31SEmil Constantinescu ierr = PetscOptionsRealArray("-ts_adapt_glee_clip","Admissible decrease/increase in step size","",glee->clip,&two,&set);CHKERRQ(ierr); 128657c1e31SEmil 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"); 129657c1e31SEmil Constantinescu ierr = PetscOptionsReal("-ts_adapt_glee_safety","Safety factor relative to target error","",glee->safety,&glee->safety,NULL);CHKERRQ(ierr); 130657c1e31SEmil 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); 131657c1e31SEmil 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); 132ae2316d0SEmil Constantinescu ierr = PetscOptionsTail();CHKERRQ(ierr); 133ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 134ae2316d0SEmil Constantinescu } 135ae2316d0SEmil Constantinescu 136ae2316d0SEmil Constantinescu #undef __FUNCT__ 137726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptView_GLEE" 138726095e4SEmil Constantinescu static PetscErrorCode TSAdaptView_GLEE(TSAdapt adapt,PetscViewer viewer) 139ae2316d0SEmil Constantinescu { 140657c1e31SEmil Constantinescu TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 141ae2316d0SEmil Constantinescu PetscErrorCode ierr; 142ae2316d0SEmil Constantinescu PetscBool iascii; 143ae2316d0SEmil Constantinescu 144ae2316d0SEmil Constantinescu PetscFunctionBegin; 145ae2316d0SEmil Constantinescu ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 146ae2316d0SEmil Constantinescu if (iascii) { 147657c1e31SEmil Constantinescu if (glee->always_accept) {ierr = PetscViewerASCIIPrintf(viewer," GLEE: always accepting steps\n");CHKERRQ(ierr);} 148657c1e31SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," GLEE: clip fastest decrease %g, fastest increase %g\n",(double)glee->clip[0],(double)glee->clip[1]);CHKERRQ(ierr); 149657c1e31SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," GLEE: safety factor %g, extra factor after step rejection %g\n",(double)glee->safety,(double)glee->reject_safety);CHKERRQ(ierr); 150ae2316d0SEmil Constantinescu } 151ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 152ae2316d0SEmil Constantinescu } 153ae2316d0SEmil Constantinescu 154ae2316d0SEmil Constantinescu #undef __FUNCT__ 155726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptCreate_GLEE" 156ae2316d0SEmil Constantinescu /*MC 157657c1e31SEmil Constantinescu TSADAPTGLEE - GLEE adaptive controller for time stepping 158ae2316d0SEmil Constantinescu 159ae2316d0SEmil Constantinescu Level: intermediate 160ae2316d0SEmil Constantinescu 161ae2316d0SEmil Constantinescu .seealso: TS, TSAdapt, TSSetAdapt() 162ae2316d0SEmil Constantinescu M*/ 163726095e4SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt adapt) 164ae2316d0SEmil Constantinescu { 165ae2316d0SEmil Constantinescu PetscErrorCode ierr; 166726095e4SEmil Constantinescu TSAdapt_GLEE *a; 167ae2316d0SEmil Constantinescu 168ae2316d0SEmil Constantinescu PetscFunctionBegin; 169ae2316d0SEmil Constantinescu ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr); 170ae2316d0SEmil Constantinescu adapt->data = (void*)a; 171726095e4SEmil Constantinescu adapt->ops->choose = TSAdaptChoose_GLEE; 172726095e4SEmil Constantinescu adapt->ops->setfromoptions = TSAdaptSetFromOptions_GLEE; 173726095e4SEmil Constantinescu adapt->ops->destroy = TSAdaptDestroy_GLEE; 174726095e4SEmil Constantinescu adapt->ops->view = TSAdaptView_GLEE; 175ae2316d0SEmil Constantinescu 176ae2316d0SEmil Constantinescu a->clip[0] = 0.1; 177ae2316d0SEmil Constantinescu a->clip[1] = 10.; 178*5e4ed32fSEmil Constantinescu a->safety = 0.9; 179ae2316d0SEmil Constantinescu a->reject_safety = 0.5; 180ae2316d0SEmil Constantinescu a->always_accept = PETSC_FALSE; 181ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 182ae2316d0SEmil Constantinescu } 183