1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 284df9cb4SJed Brown 384df9cb4SJed Brown typedef struct { 41c3436cfSJed Brown PetscBool always_accept; 51c3436cfSJed Brown PetscReal clip[2]; /* admissible decrease/increase factors */ 61c3436cfSJed Brown PetscReal safety; /* safety factor relative to target error */ 77e1ba4dcSJed Brown PetscReal reject_safety; /* extra safety factor if the last step was rejected */ 81c3436cfSJed Brown Vec Y; 984df9cb4SJed Brown } TSAdapt_Basic; 1084df9cb4SJed Brown 1184df9cb4SJed Brown #undef __FUNCT__ 121c3436cfSJed Brown #define __FUNCT__ "TSAdaptChoose_Basic" 13*5e4ed32fSEmil Constantinescu static PetscErrorCode TSAdaptChoose_Basic(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte,PetscReal *wltea,PetscReal *wlter) 1484df9cb4SJed Brown { 151c3436cfSJed Brown TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 16b1f5bccaSLisandro Dalcin PetscInt order = PETSC_DECIDE; 17b1f5bccaSLisandro Dalcin PetscReal enorm = -1; 187453f775SEmil Constantinescu PetscReal enorma,enormr; 191566a47fSLisandro Dalcin PetscReal safety = basic->safety; 20b1f5bccaSLisandro Dalcin PetscReal hfac_lte,h_lte; 211c3436cfSJed Brown PetscErrorCode ierr; 2284df9cb4SJed Brown 2384df9cb4SJed Brown PetscFunctionBegin; 241566a47fSLisandro Dalcin *next_sc = 0; /* Reuse the same order scheme */ 251566a47fSLisandro Dalcin 261ebf93c6SLisandro Dalcin if (ts->ops->evaluatewlte) { 271566a47fSLisandro Dalcin ierr = TSEvaluateWLTE(ts,adapt->wnormtype,&order,&enorm);CHKERRQ(ierr); 281566a47fSLisandro Dalcin if (enorm >= 0 && order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed error order %D must be positive",order); 291566a47fSLisandro Dalcin } else if (ts->ops->evaluatestep) { 301566a47fSLisandro Dalcin if (adapt->candidates.n < 1) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"No candidate has been registered"); 311566a47fSLisandro Dalcin if (!adapt->candidates.inuse_set) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"The current in-use scheme is not among the %D candidates",adapt->candidates.n); 321566a47fSLisandro Dalcin if (!basic->Y) {ierr = VecDuplicate(ts->vec_sol,&basic->Y);CHKERRQ(ierr);} 331566a47fSLisandro Dalcin order = adapt->candidates.order[0]; 341ebf93c6SLisandro Dalcin ierr = TSEvaluateStep(ts,order-1,basic->Y,NULL);CHKERRQ(ierr); 357453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,ts->vec_sol,basic->Y,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr); 361ebf93c6SLisandro Dalcin } 371c3436cfSJed Brown 381566a47fSLisandro Dalcin if (enorm < 0) { 391566a47fSLisandro Dalcin *accept = PETSC_TRUE; 401566a47fSLisandro Dalcin *next_h = h; /* Reuse the old step */ 411566a47fSLisandro Dalcin *wlte = -1; /* Weighted local truncation error was not evaluated */ 421566a47fSLisandro Dalcin PetscFunctionReturn(0); 431566a47fSLisandro Dalcin } 441566a47fSLisandro Dalcin 451566a47fSLisandro Dalcin /* Determine whether the step is accepted of rejected */ 461566a47fSLisandro Dalcin if (enorm > 1) { 477e1ba4dcSJed Brown if (!*accept) safety *= basic->reject_safety; /* The last attempt also failed, shorten more aggressively */ 48fd94acc0SJed Brown if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) { 4957622a8eSBarry Smith ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting because step size %g is at minimum\n",(double)enorm,(double)h);CHKERRQ(ierr); 50fd94acc0SJed Brown *accept = PETSC_TRUE; 517e1ba4dcSJed Brown } else if (basic->always_accept) { 5257622a8eSBarry Smith 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); 537e1ba4dcSJed Brown *accept = PETSC_TRUE; 541c3436cfSJed Brown } else { 5557622a8eSBarry Smith ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, rejecting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr); 56fd94acc0SJed Brown *accept = PETSC_FALSE; 57fd94acc0SJed Brown } 58fd94acc0SJed Brown } else { 5957622a8eSBarry Smith ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr); 60fd94acc0SJed Brown *accept = PETSC_TRUE; 611c3436cfSJed Brown } 621c3436cfSJed Brown 631c3436cfSJed Brown /* The optimal new step based purely on local truncation error for this step. */ 641566a47fSLisandro Dalcin if (enorm > 0) 651566a47fSLisandro Dalcin hfac_lte = safety * PetscPowReal(enorm,((PetscReal)-1)/order); 661566a47fSLisandro Dalcin else 67bf0c4ff7SBarry Smith hfac_lte = safety * PETSC_INFINITY; 681c3436cfSJed Brown h_lte = h * PetscClipInterval(hfac_lte,basic->clip[0],basic->clip[1]); 691c3436cfSJed Brown 701c3436cfSJed Brown *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max); 710b99f514SJed Brown *wlte = enorm; 72*5e4ed32fSEmil Constantinescu *wltea = -1; /* Weighted absolute local truncation error is not used */ 73*5e4ed32fSEmil Constantinescu *wlter = -1; /* Weighted relative local truncation error is not used */ 74*5e4ed32fSEmil Constantinescu 7584df9cb4SJed Brown PetscFunctionReturn(0); 7684df9cb4SJed Brown } 7784df9cb4SJed Brown 7884df9cb4SJed Brown #undef __FUNCT__ 79099af0a3SLisandro Dalcin #define __FUNCT__ "TSAdaptReset_Basic" 80099af0a3SLisandro Dalcin static PetscErrorCode TSAdaptReset_Basic(TSAdapt adapt) 8184df9cb4SJed Brown { 821c3436cfSJed Brown TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 8384df9cb4SJed Brown PetscErrorCode ierr; 8484df9cb4SJed Brown 8584df9cb4SJed Brown PetscFunctionBegin; 861c3436cfSJed Brown ierr = VecDestroy(&basic->Y);CHKERRQ(ierr); 87099af0a3SLisandro Dalcin PetscFunctionReturn(0); 88099af0a3SLisandro Dalcin } 89099af0a3SLisandro Dalcin 90099af0a3SLisandro Dalcin #undef __FUNCT__ 91099af0a3SLisandro Dalcin #define __FUNCT__ "TSAdaptDestroy_Basic" 92099af0a3SLisandro Dalcin static PetscErrorCode TSAdaptDestroy_Basic(TSAdapt adapt) 93099af0a3SLisandro Dalcin { 94099af0a3SLisandro Dalcin PetscErrorCode ierr; 95099af0a3SLisandro Dalcin 96099af0a3SLisandro Dalcin PetscFunctionBegin; 97099af0a3SLisandro Dalcin ierr = TSAdaptReset_Basic(adapt);CHKERRQ(ierr); 9884df9cb4SJed Brown ierr = PetscFree(adapt->data);CHKERRQ(ierr); 9984df9cb4SJed Brown PetscFunctionReturn(0); 10084df9cb4SJed Brown } 10184df9cb4SJed Brown 1021c3436cfSJed Brown #undef __FUNCT__ 1031c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetFromOptions_Basic" 1044416b707SBarry Smith static PetscErrorCode TSAdaptSetFromOptions_Basic(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 1051c3436cfSJed Brown { 1061c3436cfSJed Brown TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 1071c3436cfSJed Brown PetscErrorCode ierr; 1081c3436cfSJed Brown PetscInt two; 1091c3436cfSJed Brown PetscBool set; 1101c3436cfSJed Brown 1111c3436cfSJed Brown PetscFunctionBegin; 112e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Basic adaptive controller options");CHKERRQ(ierr); 1131c3436cfSJed Brown two = 2; 1141b9b13dfSLisandro Dalcin ierr = PetscOptionsRealArray("-ts_adapt_basic_clip","Admissible decrease/increase factor in step size","TSAdaptBasicSetClip",basic->clip,&two,&set);CHKERRQ(ierr); 1151b9b13dfSLisandro Dalcin if (set && (two != 2)) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_basic_clip"); 1161b9b13dfSLisandro Dalcin if (set) {ierr = TSAdaptBasicSetClip(adapt,basic->clip[0],basic->clip[1]);CHKERRQ(ierr);} 1170298fd71SBarry Smith ierr = PetscOptionsReal("-ts_adapt_basic_safety","Safety factor relative to target error","",basic->safety,&basic->safety,NULL);CHKERRQ(ierr); 1180298fd71SBarry Smith 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); 1190298fd71SBarry Smith 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); 1201c3436cfSJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 1211c3436cfSJed Brown PetscFunctionReturn(0); 1221c3436cfSJed Brown } 1231c3436cfSJed Brown 12439850373SJed Brown #undef __FUNCT__ 12539850373SJed Brown #define __FUNCT__ "TSAdaptView_Basic" 12639850373SJed Brown static PetscErrorCode TSAdaptView_Basic(TSAdapt adapt,PetscViewer viewer) 12739850373SJed Brown { 12839850373SJed Brown TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 12939850373SJed Brown PetscErrorCode ierr; 13039850373SJed Brown PetscBool iascii; 13139850373SJed Brown 13239850373SJed Brown PetscFunctionBegin; 13339850373SJed Brown ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 13439850373SJed Brown if (iascii) { 13539850373SJed Brown if (basic->always_accept) {ierr = PetscViewerASCIIPrintf(viewer," Basic: always accepting steps\n");CHKERRQ(ierr);} 13657622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Basic: clip fastest decrease %g, fastest increase %g\n",(double)basic->clip[0],(double)basic->clip[1]);CHKERRQ(ierr); 13757622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Basic: safety factor %g, extra factor after step rejection %g\n",(double)basic->safety,(double)basic->reject_safety);CHKERRQ(ierr); 13839850373SJed Brown } 13939850373SJed Brown PetscFunctionReturn(0); 14039850373SJed Brown } 14139850373SJed Brown 14284df9cb4SJed Brown #undef __FUNCT__ 14384df9cb4SJed Brown #define __FUNCT__ "TSAdaptCreate_Basic" 14484df9cb4SJed Brown /*MC 14584df9cb4SJed Brown TSADAPTBASIC - Basic adaptive controller for time stepping 14684df9cb4SJed Brown 14784df9cb4SJed Brown Level: intermediate 14884df9cb4SJed Brown 14984df9cb4SJed Brown .seealso: TS, TSAdapt, TSSetAdapt() 15084df9cb4SJed Brown M*/ 1518cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt adapt) 15284df9cb4SJed Brown { 15384df9cb4SJed Brown PetscErrorCode ierr; 15484df9cb4SJed Brown TSAdapt_Basic *a; 15584df9cb4SJed Brown 15684df9cb4SJed Brown PetscFunctionBegin; 157b00a9115SJed Brown ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr); 15884df9cb4SJed Brown adapt->data = (void*)a; 1591c3436cfSJed Brown adapt->ops->choose = TSAdaptChoose_Basic; 1601c3436cfSJed Brown adapt->ops->setfromoptions = TSAdaptSetFromOptions_Basic; 16184df9cb4SJed Brown adapt->ops->destroy = TSAdaptDestroy_Basic; 16239850373SJed Brown adapt->ops->view = TSAdaptView_Basic; 1631c3436cfSJed Brown 1641c3436cfSJed Brown a->clip[0] = 0.1; 1651c3436cfSJed Brown a->clip[1] = 10.; 1661c3436cfSJed Brown a->safety = 0.9; 1677e1ba4dcSJed Brown a->reject_safety = 0.5; 1687e1ba4dcSJed Brown a->always_accept = PETSC_FALSE; 16984df9cb4SJed Brown PetscFunctionReturn(0); 17084df9cb4SJed Brown } 1711b9b13dfSLisandro Dalcin 1721b9b13dfSLisandro Dalcin #undef __FUNCT__ 1731b9b13dfSLisandro Dalcin #define __FUNCT__ "TSAdaptBasicSetClip" 1741b9b13dfSLisandro Dalcin /*@ 1751b9b13dfSLisandro Dalcin TSAdaptBasicSetClip - Sets the admissible decrease/increase factor in step size 1761b9b13dfSLisandro Dalcin 1771b9b13dfSLisandro Dalcin Collective on TSAdapt 1781b9b13dfSLisandro Dalcin 1791b9b13dfSLisandro Dalcin Input Arguments: 1801b9b13dfSLisandro Dalcin + adapt - adaptive controller context 1811b9b13dfSLisandro Dalcin . low - admissible decrease factor 1821b9b13dfSLisandro Dalcin + high - admissible increase factor 1831b9b13dfSLisandro Dalcin 1841b9b13dfSLisandro Dalcin Level: intermediate 1851b9b13dfSLisandro Dalcin 1861b9b13dfSLisandro Dalcin .seealso: TSAdaptChoose() 1871b9b13dfSLisandro Dalcin @*/ 1881b9b13dfSLisandro Dalcin PetscErrorCode TSAdaptBasicSetClip(TSAdapt adapt,PetscReal low,PetscReal high) 1891b9b13dfSLisandro Dalcin { 1901b9b13dfSLisandro Dalcin TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 1911b9b13dfSLisandro Dalcin PetscBool isbasic; 1921b9b13dfSLisandro Dalcin PetscErrorCode ierr; 1931b9b13dfSLisandro Dalcin 1941b9b13dfSLisandro Dalcin PetscFunctionBegin; 1951b9b13dfSLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 196b1f5bccaSLisandro Dalcin if (low != PETSC_DEFAULT && low < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be non negative",(double)low); 197b1f5bccaSLisandro Dalcin if (low != PETSC_DEFAULT && low > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Decrease factor %g must be less than one",(double)low); 198b1f5bccaSLisandro Dalcin if (high != PETSC_DEFAULT && high < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Increase factor %g must be geather than one",(double)high); 1991b9b13dfSLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)adapt,TSADAPTBASIC,&isbasic);CHKERRQ(ierr); 2001b9b13dfSLisandro Dalcin if (isbasic) { 201b1f5bccaSLisandro Dalcin if (low != PETSC_DEFAULT) basic->clip[0] = low; 202b1f5bccaSLisandro Dalcin if (high != PETSC_DEFAULT) basic->clip[1] = high; 2031b9b13dfSLisandro Dalcin } 2041b9b13dfSLisandro Dalcin PetscFunctionReturn(0); 2051b9b13dfSLisandro Dalcin } 2061b9b13dfSLisandro Dalcin 2071b9b13dfSLisandro Dalcin #undef __FUNCT__ 2081b9b13dfSLisandro Dalcin #define __FUNCT__ "TSAdaptBasicGetClip" 2091b9b13dfSLisandro Dalcin /*@ 2101b9b13dfSLisandro Dalcin TSAdaptBasicGetClip - Gets the admissible decrease/increase factor in step size 2111b9b13dfSLisandro Dalcin 2121b9b13dfSLisandro Dalcin Collective on TSAdapt 2131b9b13dfSLisandro Dalcin 2141b9b13dfSLisandro Dalcin Input Arguments: 2151b9b13dfSLisandro Dalcin . adapt - adaptive controller context 2161b9b13dfSLisandro Dalcin 2171b9b13dfSLisandro Dalcin Ouput Arguments: 2181b9b13dfSLisandro Dalcin + low - optional, admissible decrease factor 2191b9b13dfSLisandro Dalcin - high - optional, admissible increase factor 2201b9b13dfSLisandro Dalcin 2211b9b13dfSLisandro Dalcin Level: intermediate 2221b9b13dfSLisandro Dalcin 2231b9b13dfSLisandro Dalcin .seealso: TSAdaptChoose() 2241b9b13dfSLisandro Dalcin @*/ 2251b9b13dfSLisandro Dalcin PetscErrorCode TSAdaptBasicGetClip(TSAdapt adapt,PetscReal *low,PetscReal *high) 2261b9b13dfSLisandro Dalcin { 2271b9b13dfSLisandro Dalcin TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 2281b9b13dfSLisandro Dalcin PetscBool isbasic; 2291b9b13dfSLisandro Dalcin PetscErrorCode ierr; 2301b9b13dfSLisandro Dalcin 2311b9b13dfSLisandro Dalcin PetscFunctionBegin; 2321b9b13dfSLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 2331b9b13dfSLisandro Dalcin if (low) PetscValidRealPointer(low,2); 2341b9b13dfSLisandro Dalcin if (high) PetscValidRealPointer(high,3); 2351b9b13dfSLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)adapt,TSADAPTBASIC,&isbasic);CHKERRQ(ierr); 2361b9b13dfSLisandro Dalcin if (isbasic) { 2371b9b13dfSLisandro Dalcin if (low) *low = basic->clip[0]; 2381b9b13dfSLisandro Dalcin if (high) *high = basic->clip[1]; 2391b9b13dfSLisandro Dalcin } else { 2401b9b13dfSLisandro Dalcin if (low) *low = 0; 2411b9b13dfSLisandro Dalcin if (high) *high = PETSC_MAX_REAL; 2421b9b13dfSLisandro Dalcin } 2431b9b13dfSLisandro Dalcin PetscFunctionReturn(0); 2441b9b13dfSLisandro Dalcin } 245