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 */ 8*657c1e31SEmil Constantinescu Vec X; 9ae2316d0SEmil Constantinescu Vec Y; 10726095e4SEmil Constantinescu } TSAdapt_GLEE; 11ae2316d0SEmil Constantinescu 12ae2316d0SEmil Constantinescu #undef __FUNCT__ 13726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptChoose_GLEE" 14726095e4SEmil Constantinescu static PetscErrorCode TSAdaptChoose_GLEE(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte) 15ae2316d0SEmil Constantinescu { 16*657c1e31SEmil Constantinescu TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 170a01e1b2SEmil Constantinescu TSType time_scheme; /* Type of time-integration scheme */ 18ae2316d0SEmil Constantinescu PetscErrorCode ierr; 19*657c1e31SEmil Constantinescu Vec X,Y,Z; 207453f775SEmil Constantinescu PetscReal enorm,enorma,enormr,hfac_lte,h_lte,safety; 21ae2316d0SEmil Constantinescu PetscInt order,stepno; 22ae2316d0SEmil Constantinescu 23ae2316d0SEmil Constantinescu PetscFunctionBegin; 24ae2316d0SEmil Constantinescu ierr = TSGetTimeStepNumber(ts,&stepno);CHKERRQ(ierr); 250a01e1b2SEmil Constantinescu 26*657c1e31SEmil Constantinescu safety = glee->safety; 270a01e1b2SEmil Constantinescu ierr = TSGetType(ts,&time_scheme);CHKERRQ(ierr); 280a01e1b2SEmil Constantinescu if (!strcmp(time_scheme,TSGLEE)){ 290a01e1b2SEmil Constantinescu /* the method is of GLEE type */ 30*657c1e31SEmil Constantinescu if (!glee->X) { 31*657c1e31SEmil Constantinescu ierr = TSGetSolution(ts,&Z);CHKERRQ(ierr); 32*657c1e31SEmil Constantinescu ierr = VecDuplicate(Z,&glee->X);CHKERRQ(ierr); 33*657c1e31SEmil Constantinescu } 34*657c1e31SEmil Constantinescu X = glee->X; 35*657c1e31SEmil Constantinescu if (!glee->Y) {ierr = VecDuplicate(X,&glee->Y);CHKERRQ(ierr);} 36*657c1e31SEmil Constantinescu Y = glee->Y; 370a01e1b2SEmil Constantinescu ierr = TSGetTimeError(ts,-1,&X);CHKERRQ(ierr); 380a01e1b2SEmil Constantinescu ierr = TSGetTimeError(ts, 0,&Y);CHKERRQ(ierr); 390a01e1b2SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr); 400a01e1b2SEmil Constantinescu } else { 410a01e1b2SEmil Constantinescu /* the method is NOT of GLEE type */ 42ae2316d0SEmil Constantinescu ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 43*657c1e31SEmil Constantinescu if (!glee->Y) {ierr = VecDuplicate(X,&glee->Y);CHKERRQ(ierr);} 44*657c1e31SEmil Constantinescu Y = glee->Y; 45ae2316d0SEmil Constantinescu order = adapt->candidates.order[0]; 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 } 49ae2316d0SEmil Constantinescu if (enorm > 1.) { 50*657c1e31SEmil Constantinescu if (!*accept) safety *= glee->reject_safety; /* The last attempt also failed, shorten more aggressively */ 51ae2316d0SEmil Constantinescu if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) { 52ae2316d0SEmil 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); 53ae2316d0SEmil Constantinescu *accept = PETSC_TRUE; 54*657c1e31SEmil Constantinescu } else if (glee->always_accept) { 55ae2316d0SEmil 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); 56ae2316d0SEmil Constantinescu *accept = PETSC_TRUE; 57ae2316d0SEmil Constantinescu } else { 58ae2316d0SEmil Constantinescu ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, rejecting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr); 59ae2316d0SEmil Constantinescu *accept = PETSC_FALSE; 60ae2316d0SEmil Constantinescu } 61ae2316d0SEmil Constantinescu } else { 62ae2316d0SEmil Constantinescu ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr); 63ae2316d0SEmil Constantinescu *accept = PETSC_TRUE; 64ae2316d0SEmil Constantinescu } 65ae2316d0SEmil Constantinescu 66ae2316d0SEmil Constantinescu /* The optimal new step based purely on local truncation error for this step. */ 67ae2316d0SEmil Constantinescu if (enorm == 0.0) { 68ae2316d0SEmil Constantinescu hfac_lte = safety * PETSC_INFINITY; 69ae2316d0SEmil Constantinescu } else { 70ae2316d0SEmil Constantinescu hfac_lte = safety * PetscPowReal(enorm,-1./order); 71ae2316d0SEmil Constantinescu } 72*657c1e31SEmil Constantinescu h_lte = h * PetscClipInterval(hfac_lte,glee->clip[0],glee->clip[1]); 73ae2316d0SEmil Constantinescu 74ae2316d0SEmil Constantinescu *next_sc = 0; 75ae2316d0SEmil Constantinescu *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max); 76ae2316d0SEmil Constantinescu *wlte = enorm; 77ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 78ae2316d0SEmil Constantinescu } 79ae2316d0SEmil Constantinescu 80ae2316d0SEmil Constantinescu #undef __FUNCT__ 81726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptReset_GLEE" 82726095e4SEmil Constantinescu static PetscErrorCode TSAdaptReset_GLEE(TSAdapt adapt) 83ae2316d0SEmil Constantinescu { 84*657c1e31SEmil Constantinescu TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 85ae2316d0SEmil Constantinescu PetscErrorCode ierr; 86ae2316d0SEmil Constantinescu 87ae2316d0SEmil Constantinescu PetscFunctionBegin; 88*657c1e31SEmil Constantinescu ierr = VecDestroy(&glee->Y);CHKERRQ(ierr); 89ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 90ae2316d0SEmil Constantinescu } 91ae2316d0SEmil Constantinescu 92ae2316d0SEmil Constantinescu #undef __FUNCT__ 93726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptDestroy_GLEE" 94726095e4SEmil Constantinescu static PetscErrorCode TSAdaptDestroy_GLEE(TSAdapt adapt) 95ae2316d0SEmil Constantinescu { 96ae2316d0SEmil Constantinescu PetscErrorCode ierr; 97ae2316d0SEmil Constantinescu 98ae2316d0SEmil Constantinescu PetscFunctionBegin; 99726095e4SEmil Constantinescu ierr = TSAdaptReset_GLEE(adapt);CHKERRQ(ierr); 100ae2316d0SEmil Constantinescu ierr = PetscFree(adapt->data);CHKERRQ(ierr); 101ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 102ae2316d0SEmil Constantinescu } 103ae2316d0SEmil Constantinescu 104ae2316d0SEmil Constantinescu #undef __FUNCT__ 105726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptSetFromOptions_GLEE" 106726095e4SEmil Constantinescu static PetscErrorCode TSAdaptSetFromOptions_GLEE(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 107ae2316d0SEmil Constantinescu { 108*657c1e31SEmil Constantinescu TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 109ae2316d0SEmil Constantinescu PetscErrorCode ierr; 110ae2316d0SEmil Constantinescu PetscInt two; 111ae2316d0SEmil Constantinescu PetscBool set; 112ae2316d0SEmil Constantinescu 113ae2316d0SEmil Constantinescu PetscFunctionBegin; 114726095e4SEmil Constantinescu ierr = PetscOptionsHead(PetscOptionsObject,"GLEE adaptive controller options");CHKERRQ(ierr); 115ae2316d0SEmil Constantinescu two = 2; 116*657c1e31SEmil Constantinescu ierr = PetscOptionsRealArray("-ts_adapt_glee_clip","Admissible decrease/increase in step size","",glee->clip,&two,&set);CHKERRQ(ierr); 117*657c1e31SEmil 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"); 118*657c1e31SEmil Constantinescu ierr = PetscOptionsReal("-ts_adapt_glee_safety","Safety factor relative to target error","",glee->safety,&glee->safety,NULL);CHKERRQ(ierr); 119*657c1e31SEmil 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); 120*657c1e31SEmil 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); 121ae2316d0SEmil Constantinescu ierr = PetscOptionsTail();CHKERRQ(ierr); 122ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 123ae2316d0SEmil Constantinescu } 124ae2316d0SEmil Constantinescu 125ae2316d0SEmil Constantinescu #undef __FUNCT__ 126726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptView_GLEE" 127726095e4SEmil Constantinescu static PetscErrorCode TSAdaptView_GLEE(TSAdapt adapt,PetscViewer viewer) 128ae2316d0SEmil Constantinescu { 129*657c1e31SEmil Constantinescu TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 130ae2316d0SEmil Constantinescu PetscErrorCode ierr; 131ae2316d0SEmil Constantinescu PetscBool iascii; 132ae2316d0SEmil Constantinescu 133ae2316d0SEmil Constantinescu PetscFunctionBegin; 134ae2316d0SEmil Constantinescu ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 135ae2316d0SEmil Constantinescu if (iascii) { 136*657c1e31SEmil Constantinescu if (glee->always_accept) {ierr = PetscViewerASCIIPrintf(viewer," GLEE: always accepting steps\n");CHKERRQ(ierr);} 137*657c1e31SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," GLEE: clip fastest decrease %g, fastest increase %g\n",(double)glee->clip[0],(double)glee->clip[1]);CHKERRQ(ierr); 138*657c1e31SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," GLEE: safety factor %g, extra factor after step rejection %g\n",(double)glee->safety,(double)glee->reject_safety);CHKERRQ(ierr); 139ae2316d0SEmil Constantinescu } 140ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 141ae2316d0SEmil Constantinescu } 142ae2316d0SEmil Constantinescu 143ae2316d0SEmil Constantinescu #undef __FUNCT__ 144726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptCreate_GLEE" 145ae2316d0SEmil Constantinescu /*MC 146*657c1e31SEmil Constantinescu TSADAPTGLEE - GLEE adaptive controller for time stepping 147ae2316d0SEmil Constantinescu 148ae2316d0SEmil Constantinescu Level: intermediate 149ae2316d0SEmil Constantinescu 150ae2316d0SEmil Constantinescu .seealso: TS, TSAdapt, TSSetAdapt() 151ae2316d0SEmil Constantinescu M*/ 152726095e4SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt adapt) 153ae2316d0SEmil Constantinescu { 154ae2316d0SEmil Constantinescu PetscErrorCode ierr; 155726095e4SEmil Constantinescu TSAdapt_GLEE *a; 156ae2316d0SEmil Constantinescu 157ae2316d0SEmil Constantinescu PetscFunctionBegin; 158ae2316d0SEmil Constantinescu ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr); 159ae2316d0SEmil Constantinescu adapt->data = (void*)a; 160726095e4SEmil Constantinescu adapt->ops->choose = TSAdaptChoose_GLEE; 161726095e4SEmil Constantinescu adapt->ops->setfromoptions = TSAdaptSetFromOptions_GLEE; 162726095e4SEmil Constantinescu adapt->ops->destroy = TSAdaptDestroy_GLEE; 163726095e4SEmil Constantinescu adapt->ops->view = TSAdaptView_GLEE; 164ae2316d0SEmil Constantinescu 165ae2316d0SEmil Constantinescu a->clip[0] = 0.1; 166ae2316d0SEmil Constantinescu a->clip[1] = 10.; 167ae2316d0SEmil Constantinescu a->safety = 0.9; 168ae2316d0SEmil Constantinescu a->reject_safety = 0.5; 169ae2316d0SEmil Constantinescu a->always_accept = PETSC_FALSE; 170ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 171ae2316d0SEmil Constantinescu } 172