184df9cb4SJed Brown #include <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 */ 71c3436cfSJed Brown Vec Y; 884df9cb4SJed Brown } TSAdapt_Basic; 984df9cb4SJed Brown 1084df9cb4SJed Brown #undef __FUNCT__ 111c3436cfSJed Brown #define __FUNCT__ "TSAdaptChoose_Basic" 12*0b99f514SJed Brown static PetscErrorCode TSAdaptChoose_Basic(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte) 1384df9cb4SJed Brown { 141c3436cfSJed Brown TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 151c3436cfSJed Brown PetscErrorCode ierr; 161c3436cfSJed Brown Vec X,Y; 171c3436cfSJed Brown PetscReal enorm,hfac_lte,h_lte; 181c3436cfSJed Brown PetscInt order,stepno; 1984df9cb4SJed Brown 2084df9cb4SJed Brown PetscFunctionBegin; 211c3436cfSJed Brown ierr = TSGetTimeStepNumber(ts,&stepno);CHKERRQ(ierr); 221c3436cfSJed Brown ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 231c3436cfSJed Brown if (!basic->Y) {ierr = VecDuplicate(X,&basic->Y);CHKERRQ(ierr);} 241c3436cfSJed Brown Y = basic->Y; 251c3436cfSJed Brown order = adapt->candidates.order[0]; 261c3436cfSJed Brown ierr = TSEvaluateStep(ts,order-1,Y,PETSC_NULL);CHKERRQ(ierr); 271c3436cfSJed Brown 281c3436cfSJed Brown ierr = TSErrorNormWRMS(ts,Y,&enorm);CHKERRQ(ierr); 291c3436cfSJed Brown if (enorm > 1.) { 30fd94acc0SJed Brown if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) { 31fd94acc0SJed Brown ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %G, accepting because step size %G is at minimum\n",enorm,h);CHKERRQ(ierr); 32fd94acc0SJed Brown *accept = PETSC_TRUE; 331c3436cfSJed Brown } else { 34fd94acc0SJed Brown ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %G, rejecting step of size %G\n",enorm,h);CHKERRQ(ierr); 35fd94acc0SJed Brown *accept = PETSC_FALSE; 36fd94acc0SJed Brown } 37fd94acc0SJed Brown } else { 38fd94acc0SJed Brown ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %G, accepting step of size %G\n",enorm,h);CHKERRQ(ierr); 39fd94acc0SJed Brown *accept = PETSC_TRUE; 401c3436cfSJed Brown } 411c3436cfSJed Brown 421c3436cfSJed Brown /* The optimal new step based purely on local truncation error for this step. */ 433116ef58SSatish Balay hfac_lte = basic->safety * PetscRealPart(PetscPowScalar((PetscScalar)enorm,(PetscReal)(-1./order))); 441c3436cfSJed Brown h_lte = h * PetscClipInterval(hfac_lte,basic->clip[0],basic->clip[1]); 451c3436cfSJed Brown 4684df9cb4SJed Brown *next_sc = 0; 471c3436cfSJed Brown *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max); 48*0b99f514SJed Brown *wlte = enorm; 4984df9cb4SJed Brown PetscFunctionReturn(0); 5084df9cb4SJed Brown } 5184df9cb4SJed Brown 5284df9cb4SJed Brown #undef __FUNCT__ 5384df9cb4SJed Brown #define __FUNCT__ "TSAdaptDestroy_Basic" 5484df9cb4SJed Brown static PetscErrorCode TSAdaptDestroy_Basic(TSAdapt adapt) 5584df9cb4SJed Brown { 561c3436cfSJed Brown TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 5784df9cb4SJed Brown PetscErrorCode ierr; 5884df9cb4SJed Brown 5984df9cb4SJed Brown PetscFunctionBegin; 601c3436cfSJed Brown ierr = VecDestroy(&basic->Y);CHKERRQ(ierr); 6184df9cb4SJed Brown ierr = PetscFree(adapt->data);CHKERRQ(ierr); 6284df9cb4SJed Brown PetscFunctionReturn(0); 6384df9cb4SJed Brown } 6484df9cb4SJed Brown 651c3436cfSJed Brown #undef __FUNCT__ 661c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetFromOptions_Basic" 671c3436cfSJed Brown static PetscErrorCode TSAdaptSetFromOptions_Basic(TSAdapt adapt) 681c3436cfSJed Brown { 691c3436cfSJed Brown TSAdapt_Basic *basic = (TSAdapt_Basic*)adapt->data; 701c3436cfSJed Brown PetscErrorCode ierr; 711c3436cfSJed Brown PetscInt two; 721c3436cfSJed Brown PetscBool set; 731c3436cfSJed Brown 741c3436cfSJed Brown PetscFunctionBegin; 751c3436cfSJed Brown ierr = PetscOptionsHead("Basic adaptive controller options");CHKERRQ(ierr); 761c3436cfSJed Brown two = 2; 771c3436cfSJed Brown ierr = PetscOptionsRealArray("-ts_adapt_basic_clip","Admissible decrease/increase in step size","",basic->clip,&two,&set);CHKERRQ(ierr); 781c3436cfSJed Brown if (set && (two != 2 || basic->clip[0] > basic->clip[1])) 791c3436cfSJed Brown SETERRQ(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_basic_clip"); 801c3436cfSJed Brown ierr = PetscOptionsReal("-ts_adapt_basic_safety","Safety factor relative to target error","",basic->safety,&basic->safety,PETSC_NULL);CHKERRQ(ierr); 811c3436cfSJed Brown 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,PETSC_NULL);CHKERRQ(ierr); 821c3436cfSJed Brown if (!basic->always_accept) SETERRQ(((PetscObject)adapt)->comm,PETSC_ERR_SUP,"step rejection not implemented yet"); 831c3436cfSJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 841c3436cfSJed Brown PetscFunctionReturn(0); 851c3436cfSJed Brown } 861c3436cfSJed Brown 8784df9cb4SJed Brown EXTERN_C_BEGIN 8884df9cb4SJed Brown #undef __FUNCT__ 8984df9cb4SJed Brown #define __FUNCT__ "TSAdaptCreate_Basic" 9084df9cb4SJed Brown /*MC 9184df9cb4SJed Brown TSADAPTBASIC - Basic adaptive controller for time stepping 9284df9cb4SJed Brown 9384df9cb4SJed Brown Level: intermediate 9484df9cb4SJed Brown 9584df9cb4SJed Brown .seealso: TS, TSAdapt, TSSetAdapt() 9684df9cb4SJed Brown M*/ 9784df9cb4SJed Brown PetscErrorCode TSAdaptCreate_Basic(TSAdapt adapt) 9884df9cb4SJed Brown { 9984df9cb4SJed Brown PetscErrorCode ierr; 10084df9cb4SJed Brown TSAdapt_Basic *a; 10184df9cb4SJed Brown 10284df9cb4SJed Brown PetscFunctionBegin; 10384df9cb4SJed Brown ierr = PetscNewLog(adapt,TSAdapt_Basic,&a);CHKERRQ(ierr); 10484df9cb4SJed Brown adapt->data = (void*)a; 1051c3436cfSJed Brown adapt->ops->choose = TSAdaptChoose_Basic; 1061c3436cfSJed Brown adapt->ops->setfromoptions = TSAdaptSetFromOptions_Basic; 10784df9cb4SJed Brown adapt->ops->destroy = TSAdaptDestroy_Basic; 1081c3436cfSJed Brown 1091c3436cfSJed Brown a->clip[0] = 0.1; 1101c3436cfSJed Brown a->clip[1] = 10.; 1111c3436cfSJed Brown a->safety = 0.9; 1121c3436cfSJed Brown a->always_accept = PETSC_TRUE; /* fix this */ 11384df9cb4SJed Brown PetscFunctionReturn(0); 11484df9cb4SJed Brown } 11584df9cb4SJed Brown EXTERN_C_END 116