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 */ 8ae2316d0SEmil Constantinescu Vec Y; 9726095e4SEmil Constantinescu } TSAdapt_GLEE; 10ae2316d0SEmil Constantinescu 11ae2316d0SEmil Constantinescu #undef __FUNCT__ 12726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptChoose_GLEE" 13726095e4SEmil Constantinescu static PetscErrorCode TSAdaptChoose_GLEE(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte) 14ae2316d0SEmil Constantinescu { 15726095e4SEmil Constantinescu TSAdapt_GLEE *basic = (TSAdapt_GLEE*)adapt->data; 16ae2316d0SEmil Constantinescu PetscErrorCode ierr; 17ae2316d0SEmil Constantinescu Vec X,Y; 18*7453f775SEmil Constantinescu PetscReal enorm,enorma,enormr,hfac_lte,h_lte,safety; 19ae2316d0SEmil Constantinescu PetscInt order,stepno; 20ae2316d0SEmil Constantinescu 21ae2316d0SEmil Constantinescu PetscFunctionBegin; 22ae2316d0SEmil Constantinescu ierr = TSGetTimeStepNumber(ts,&stepno);CHKERRQ(ierr); 23ae2316d0SEmil Constantinescu ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 24ae2316d0SEmil Constantinescu if (!basic->Y) {ierr = VecDuplicate(X,&basic->Y);CHKERRQ(ierr);} 25ae2316d0SEmil Constantinescu Y = basic->Y; 26ae2316d0SEmil Constantinescu order = adapt->candidates.order[0]; 27ae2316d0SEmil Constantinescu ierr = TSEvaluateStep(ts,order-1,Y,NULL);CHKERRQ(ierr); 28ae2316d0SEmil Constantinescu 29ae2316d0SEmil Constantinescu safety = basic->safety; 30*7453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr); 31ae2316d0SEmil Constantinescu if (enorm > 1.) { 32ae2316d0SEmil Constantinescu if (!*accept) safety *= basic->reject_safety; /* The last attempt also failed, shorten more aggressively */ 33ae2316d0SEmil Constantinescu if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) { 34ae2316d0SEmil 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); 35ae2316d0SEmil Constantinescu *accept = PETSC_TRUE; 36ae2316d0SEmil Constantinescu } else if (basic->always_accept) { 37ae2316d0SEmil 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); 38ae2316d0SEmil Constantinescu *accept = PETSC_TRUE; 39ae2316d0SEmil Constantinescu } else { 40ae2316d0SEmil Constantinescu ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, rejecting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr); 41ae2316d0SEmil Constantinescu *accept = PETSC_FALSE; 42ae2316d0SEmil Constantinescu } 43ae2316d0SEmil Constantinescu } else { 44ae2316d0SEmil Constantinescu ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr); 45ae2316d0SEmil Constantinescu *accept = PETSC_TRUE; 46ae2316d0SEmil Constantinescu } 47ae2316d0SEmil Constantinescu 48ae2316d0SEmil Constantinescu /* The optimal new step based purely on local truncation error for this step. */ 49ae2316d0SEmil Constantinescu if (enorm == 0.0) { 50ae2316d0SEmil Constantinescu hfac_lte = safety * PETSC_INFINITY; 51ae2316d0SEmil Constantinescu } else { 52ae2316d0SEmil Constantinescu hfac_lte = safety * PetscPowReal(enorm,-1./order); 53ae2316d0SEmil Constantinescu } 54ae2316d0SEmil Constantinescu h_lte = h * PetscClipInterval(hfac_lte,basic->clip[0],basic->clip[1]); 55ae2316d0SEmil Constantinescu 56ae2316d0SEmil Constantinescu *next_sc = 0; 57ae2316d0SEmil Constantinescu *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max); 58ae2316d0SEmil Constantinescu *wlte = enorm; 59ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 60ae2316d0SEmil Constantinescu } 61ae2316d0SEmil Constantinescu 62ae2316d0SEmil Constantinescu #undef __FUNCT__ 63726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptReset_GLEE" 64726095e4SEmil Constantinescu static PetscErrorCode TSAdaptReset_GLEE(TSAdapt adapt) 65ae2316d0SEmil Constantinescu { 66726095e4SEmil Constantinescu TSAdapt_GLEE *basic = (TSAdapt_GLEE*)adapt->data; 67ae2316d0SEmil Constantinescu PetscErrorCode ierr; 68ae2316d0SEmil Constantinescu 69ae2316d0SEmil Constantinescu PetscFunctionBegin; 70ae2316d0SEmil Constantinescu ierr = VecDestroy(&basic->Y);CHKERRQ(ierr); 71ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 72ae2316d0SEmil Constantinescu } 73ae2316d0SEmil Constantinescu 74ae2316d0SEmil Constantinescu #undef __FUNCT__ 75726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptDestroy_GLEE" 76726095e4SEmil Constantinescu static PetscErrorCode TSAdaptDestroy_GLEE(TSAdapt adapt) 77ae2316d0SEmil Constantinescu { 78ae2316d0SEmil Constantinescu PetscErrorCode ierr; 79ae2316d0SEmil Constantinescu 80ae2316d0SEmil Constantinescu PetscFunctionBegin; 81726095e4SEmil Constantinescu ierr = TSAdaptReset_GLEE(adapt);CHKERRQ(ierr); 82ae2316d0SEmil Constantinescu ierr = PetscFree(adapt->data);CHKERRQ(ierr); 83ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 84ae2316d0SEmil Constantinescu } 85ae2316d0SEmil Constantinescu 86ae2316d0SEmil Constantinescu #undef __FUNCT__ 87726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptSetFromOptions_GLEE" 88726095e4SEmil Constantinescu static PetscErrorCode TSAdaptSetFromOptions_GLEE(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 89ae2316d0SEmil Constantinescu { 90726095e4SEmil Constantinescu TSAdapt_GLEE *basic = (TSAdapt_GLEE*)adapt->data; 91ae2316d0SEmil Constantinescu PetscErrorCode ierr; 92ae2316d0SEmil Constantinescu PetscInt two; 93ae2316d0SEmil Constantinescu PetscBool set; 94ae2316d0SEmil Constantinescu 95ae2316d0SEmil Constantinescu PetscFunctionBegin; 96726095e4SEmil Constantinescu ierr = PetscOptionsHead(PetscOptionsObject,"GLEE adaptive controller options");CHKERRQ(ierr); 97ae2316d0SEmil Constantinescu two = 2; 98ae2316d0SEmil Constantinescu ierr = PetscOptionsRealArray("-ts_adapt_basic_clip","Admissible decrease/increase in step size","",basic->clip,&two,&set);CHKERRQ(ierr); 99ae2316d0SEmil 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"); 100ae2316d0SEmil Constantinescu ierr = PetscOptionsReal("-ts_adapt_basic_safety","Safety factor relative to target error","",basic->safety,&basic->safety,NULL);CHKERRQ(ierr); 101ae2316d0SEmil 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); 102ae2316d0SEmil 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); 103ae2316d0SEmil Constantinescu ierr = PetscOptionsTail();CHKERRQ(ierr); 104ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 105ae2316d0SEmil Constantinescu } 106ae2316d0SEmil Constantinescu 107ae2316d0SEmil Constantinescu #undef __FUNCT__ 108726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptView_GLEE" 109726095e4SEmil Constantinescu static PetscErrorCode TSAdaptView_GLEE(TSAdapt adapt,PetscViewer viewer) 110ae2316d0SEmil Constantinescu { 111726095e4SEmil Constantinescu TSAdapt_GLEE *basic = (TSAdapt_GLEE*)adapt->data; 112ae2316d0SEmil Constantinescu PetscErrorCode ierr; 113ae2316d0SEmil Constantinescu PetscBool iascii; 114ae2316d0SEmil Constantinescu 115ae2316d0SEmil Constantinescu PetscFunctionBegin; 116ae2316d0SEmil Constantinescu ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 117ae2316d0SEmil Constantinescu if (iascii) { 118726095e4SEmil Constantinescu if (basic->always_accept) {ierr = PetscViewerASCIIPrintf(viewer," GLEE: always accepting steps\n");CHKERRQ(ierr);} 119726095e4SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," GLEE: clip fastest decrease %g, fastest increase %g\n",(double)basic->clip[0],(double)basic->clip[1]);CHKERRQ(ierr); 120726095e4SEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," GLEE: safety factor %g, extra factor after step rejection %g\n",(double)basic->safety,(double)basic->reject_safety);CHKERRQ(ierr); 121ae2316d0SEmil Constantinescu } 122ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 123ae2316d0SEmil Constantinescu } 124ae2316d0SEmil Constantinescu 125ae2316d0SEmil Constantinescu #undef __FUNCT__ 126726095e4SEmil Constantinescu #define __FUNCT__ "TSAdaptCreate_GLEE" 127ae2316d0SEmil Constantinescu /*MC 128726095e4SEmil Constantinescu TSADAPTBASIC - GLEE adaptive controller for time stepping 129ae2316d0SEmil Constantinescu 130ae2316d0SEmil Constantinescu Level: intermediate 131ae2316d0SEmil Constantinescu 132ae2316d0SEmil Constantinescu .seealso: TS, TSAdapt, TSSetAdapt() 133ae2316d0SEmil Constantinescu M*/ 134726095e4SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt adapt) 135ae2316d0SEmil Constantinescu { 136ae2316d0SEmil Constantinescu PetscErrorCode ierr; 137726095e4SEmil Constantinescu TSAdapt_GLEE *a; 138ae2316d0SEmil Constantinescu 139ae2316d0SEmil Constantinescu PetscFunctionBegin; 140ae2316d0SEmil Constantinescu ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr); 141ae2316d0SEmil Constantinescu adapt->data = (void*)a; 142726095e4SEmil Constantinescu adapt->ops->choose = TSAdaptChoose_GLEE; 143726095e4SEmil Constantinescu adapt->ops->setfromoptions = TSAdaptSetFromOptions_GLEE; 144726095e4SEmil Constantinescu adapt->ops->destroy = TSAdaptDestroy_GLEE; 145726095e4SEmil Constantinescu adapt->ops->view = TSAdaptView_GLEE; 146ae2316d0SEmil Constantinescu 147ae2316d0SEmil Constantinescu a->clip[0] = 0.1; 148ae2316d0SEmil Constantinescu a->clip[1] = 10.; 149ae2316d0SEmil Constantinescu a->safety = 0.9; 150ae2316d0SEmil Constantinescu a->reject_safety = 0.5; 151ae2316d0SEmil Constantinescu a->always_accept = PETSC_FALSE; 152ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 153ae2316d0SEmil Constantinescu } 154