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" 130b99f514SJed Brown static PetscErrorCode TSAdaptChoose_Basic(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte) 1484df9cb4SJed Brown { 151c3436cfSJed Brown TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 161c3436cfSJed Brown PetscErrorCode ierr; 171c3436cfSJed Brown Vec X,Y; 187e1ba4dcSJed Brown PetscReal enorm,hfac_lte,h_lte,safety; 191c3436cfSJed Brown PetscInt order,stepno; 2084df9cb4SJed Brown 2184df9cb4SJed Brown PetscFunctionBegin; 221c3436cfSJed Brown ierr = TSGetTimeStepNumber(ts,&stepno);CHKERRQ(ierr); 231c3436cfSJed Brown ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 241c3436cfSJed Brown if (!basic->Y) {ierr = VecDuplicate(X,&basic->Y);CHKERRQ(ierr);} 251c3436cfSJed Brown Y = basic->Y; 261c3436cfSJed Brown order = adapt->candidates.order[0]; 270298fd71SBarry Smith ierr = TSEvaluateStep(ts,order-1,Y,NULL);CHKERRQ(ierr); 281c3436cfSJed Brown 297e1ba4dcSJed Brown safety = basic->safety; 30a4868fbcSLisandro Dalcin ierr = TSErrorWeightedNorm(ts,X,Y,adapt->wnormtype,&enorm);CHKERRQ(ierr); 311c3436cfSJed Brown if (enorm > 1.) { 327e1ba4dcSJed Brown if (!*accept) safety *= basic->reject_safety; /* The last attempt also failed, shorten more aggressively */ 33fd94acc0SJed Brown if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) { 3457622a8eSBarry 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); 35fd94acc0SJed Brown *accept = PETSC_TRUE; 367e1ba4dcSJed Brown } else if (basic->always_accept) { 3757622a8eSBarry 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); 387e1ba4dcSJed Brown *accept = PETSC_TRUE; 391c3436cfSJed Brown } else { 4057622a8eSBarry Smith ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, rejecting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr); 41fd94acc0SJed Brown *accept = PETSC_FALSE; 42fd94acc0SJed Brown } 43fd94acc0SJed Brown } else { 4457622a8eSBarry Smith ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr); 45fd94acc0SJed Brown *accept = PETSC_TRUE; 461c3436cfSJed Brown } 471c3436cfSJed Brown 481c3436cfSJed Brown /* The optimal new step based purely on local truncation error for this step. */ 49*bf0c4ff7SBarry Smith if (enorm == 0.0) { 50*bf0c4ff7SBarry Smith hfac_lte = safety * PETSC_INFINITY; 51*bf0c4ff7SBarry Smith } else { 52a3e8c5ccSSatish Balay hfac_lte = safety * PetscPowReal(enorm,-1./order); 53*bf0c4ff7SBarry Smith } 541c3436cfSJed Brown h_lte = h * PetscClipInterval(hfac_lte,basic->clip[0],basic->clip[1]); 551c3436cfSJed Brown 5684df9cb4SJed Brown *next_sc = 0; 571c3436cfSJed Brown *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max); 580b99f514SJed Brown *wlte = enorm; 5984df9cb4SJed Brown PetscFunctionReturn(0); 6084df9cb4SJed Brown } 6184df9cb4SJed Brown 6284df9cb4SJed Brown #undef __FUNCT__ 63099af0a3SLisandro Dalcin #define __FUNCT__ "TSAdaptReset_Basic" 64099af0a3SLisandro Dalcin static PetscErrorCode TSAdaptReset_Basic(TSAdapt adapt) 6584df9cb4SJed Brown { 661c3436cfSJed Brown TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 6784df9cb4SJed Brown PetscErrorCode ierr; 6884df9cb4SJed Brown 6984df9cb4SJed Brown PetscFunctionBegin; 701c3436cfSJed Brown ierr = VecDestroy(&basic->Y);CHKERRQ(ierr); 71099af0a3SLisandro Dalcin PetscFunctionReturn(0); 72099af0a3SLisandro Dalcin } 73099af0a3SLisandro Dalcin 74099af0a3SLisandro Dalcin #undef __FUNCT__ 75099af0a3SLisandro Dalcin #define __FUNCT__ "TSAdaptDestroy_Basic" 76099af0a3SLisandro Dalcin static PetscErrorCode TSAdaptDestroy_Basic(TSAdapt adapt) 77099af0a3SLisandro Dalcin { 78099af0a3SLisandro Dalcin PetscErrorCode ierr; 79099af0a3SLisandro Dalcin 80099af0a3SLisandro Dalcin PetscFunctionBegin; 81099af0a3SLisandro Dalcin ierr = TSAdaptReset_Basic(adapt);CHKERRQ(ierr); 8284df9cb4SJed Brown ierr = PetscFree(adapt->data);CHKERRQ(ierr); 8384df9cb4SJed Brown PetscFunctionReturn(0); 8484df9cb4SJed Brown } 8584df9cb4SJed Brown 861c3436cfSJed Brown #undef __FUNCT__ 871c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetFromOptions_Basic" 888c34d3f5SBarry Smith static PetscErrorCode TSAdaptSetFromOptions_Basic(PetscOptions *PetscOptionsObject,TSAdapt adapt) 891c3436cfSJed Brown { 901c3436cfSJed Brown TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 911c3436cfSJed Brown PetscErrorCode ierr; 921c3436cfSJed Brown PetscInt two; 931c3436cfSJed Brown PetscBool set; 941c3436cfSJed Brown 951c3436cfSJed Brown PetscFunctionBegin; 96e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Basic adaptive controller options");CHKERRQ(ierr); 971c3436cfSJed Brown two = 2; 981c3436cfSJed Brown ierr = PetscOptionsRealArray("-ts_adapt_basic_clip","Admissible decrease/increase in step size","",basic->clip,&two,&set);CHKERRQ(ierr); 99ce94432eSBarry Smith 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"); 1000298fd71SBarry Smith ierr = PetscOptionsReal("-ts_adapt_basic_safety","Safety factor relative to target error","",basic->safety,&basic->safety,NULL);CHKERRQ(ierr); 1010298fd71SBarry 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); 1020298fd71SBarry 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); 1031c3436cfSJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 1041c3436cfSJed Brown PetscFunctionReturn(0); 1051c3436cfSJed Brown } 1061c3436cfSJed Brown 10739850373SJed Brown #undef __FUNCT__ 10839850373SJed Brown #define __FUNCT__ "TSAdaptView_Basic" 10939850373SJed Brown static PetscErrorCode TSAdaptView_Basic(TSAdapt adapt,PetscViewer viewer) 11039850373SJed Brown { 11139850373SJed Brown TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 11239850373SJed Brown PetscErrorCode ierr; 11339850373SJed Brown PetscBool iascii; 11439850373SJed Brown 11539850373SJed Brown PetscFunctionBegin; 11639850373SJed Brown ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 11739850373SJed Brown if (iascii) { 11839850373SJed Brown if (basic->always_accept) {ierr = PetscViewerASCIIPrintf(viewer," Basic: always accepting steps\n");CHKERRQ(ierr);} 11957622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Basic: clip fastest decrease %g, fastest increase %g\n",(double)basic->clip[0],(double)basic->clip[1]);CHKERRQ(ierr); 12057622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Basic: safety factor %g, extra factor after step rejection %g\n",(double)basic->safety,(double)basic->reject_safety);CHKERRQ(ierr); 12139850373SJed Brown } 12239850373SJed Brown PetscFunctionReturn(0); 12339850373SJed Brown } 12439850373SJed Brown 12584df9cb4SJed Brown #undef __FUNCT__ 12684df9cb4SJed Brown #define __FUNCT__ "TSAdaptCreate_Basic" 12784df9cb4SJed Brown /*MC 12884df9cb4SJed Brown TSADAPTBASIC - Basic adaptive controller for time stepping 12984df9cb4SJed Brown 13084df9cb4SJed Brown Level: intermediate 13184df9cb4SJed Brown 13284df9cb4SJed Brown .seealso: TS, TSAdapt, TSSetAdapt() 13384df9cb4SJed Brown M*/ 1348cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt adapt) 13584df9cb4SJed Brown { 13684df9cb4SJed Brown PetscErrorCode ierr; 13784df9cb4SJed Brown TSAdapt_Basic *a; 13884df9cb4SJed Brown 13984df9cb4SJed Brown PetscFunctionBegin; 140b00a9115SJed Brown ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr); 14184df9cb4SJed Brown adapt->data = (void*)a; 1421c3436cfSJed Brown adapt->ops->choose = TSAdaptChoose_Basic; 1431c3436cfSJed Brown adapt->ops->setfromoptions = TSAdaptSetFromOptions_Basic; 14484df9cb4SJed Brown adapt->ops->destroy = TSAdaptDestroy_Basic; 14539850373SJed Brown adapt->ops->view = TSAdaptView_Basic; 1461c3436cfSJed Brown 1471c3436cfSJed Brown a->clip[0] = 0.1; 1481c3436cfSJed Brown a->clip[1] = 10.; 1491c3436cfSJed Brown a->safety = 0.9; 1507e1ba4dcSJed Brown a->reject_safety = 0.5; 1517e1ba4dcSJed Brown a->always_accept = PETSC_FALSE; 15284df9cb4SJed Brown PetscFunctionReturn(0); 15384df9cb4SJed Brown } 154