xref: /petsc/src/ts/adapt/impls/basic/adaptbasic.c (revision af0996ce37bc06907c37d8d91773840993d61e62)
1*af0996ceSBarry 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. */
497e1ba4dcSJed Brown   hfac_lte = safety * PetscRealPart(PetscPowScalar((PetscScalar)enorm,(PetscReal)(-1./order)));
501c3436cfSJed Brown   h_lte    = h * PetscClipInterval(hfac_lte,basic->clip[0],basic->clip[1]);
511c3436cfSJed Brown 
5284df9cb4SJed Brown   *next_sc = 0;
531c3436cfSJed Brown   *next_h  = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max);
540b99f514SJed Brown   *wlte    = enorm;
5584df9cb4SJed Brown   PetscFunctionReturn(0);
5684df9cb4SJed Brown }
5784df9cb4SJed Brown 
5884df9cb4SJed Brown #undef __FUNCT__
59099af0a3SLisandro Dalcin #define __FUNCT__ "TSAdaptReset_Basic"
60099af0a3SLisandro Dalcin static PetscErrorCode TSAdaptReset_Basic(TSAdapt adapt)
6184df9cb4SJed Brown {
621c3436cfSJed Brown   TSAdapt_Basic  *basic = (TSAdapt_Basic*)adapt->data;
6384df9cb4SJed Brown   PetscErrorCode ierr;
6484df9cb4SJed Brown 
6584df9cb4SJed Brown   PetscFunctionBegin;
661c3436cfSJed Brown   ierr = VecDestroy(&basic->Y);CHKERRQ(ierr);
67099af0a3SLisandro Dalcin   PetscFunctionReturn(0);
68099af0a3SLisandro Dalcin }
69099af0a3SLisandro Dalcin 
70099af0a3SLisandro Dalcin #undef __FUNCT__
71099af0a3SLisandro Dalcin #define __FUNCT__ "TSAdaptDestroy_Basic"
72099af0a3SLisandro Dalcin static PetscErrorCode TSAdaptDestroy_Basic(TSAdapt adapt)
73099af0a3SLisandro Dalcin {
74099af0a3SLisandro Dalcin   PetscErrorCode ierr;
75099af0a3SLisandro Dalcin 
76099af0a3SLisandro Dalcin   PetscFunctionBegin;
77099af0a3SLisandro Dalcin   ierr = TSAdaptReset_Basic(adapt);CHKERRQ(ierr);
7884df9cb4SJed Brown   ierr = PetscFree(adapt->data);CHKERRQ(ierr);
7984df9cb4SJed Brown   PetscFunctionReturn(0);
8084df9cb4SJed Brown }
8184df9cb4SJed Brown 
821c3436cfSJed Brown #undef __FUNCT__
831c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetFromOptions_Basic"
848c34d3f5SBarry Smith static PetscErrorCode TSAdaptSetFromOptions_Basic(PetscOptions *PetscOptionsObject,TSAdapt adapt)
851c3436cfSJed Brown {
861c3436cfSJed Brown   TSAdapt_Basic  *basic = (TSAdapt_Basic*)adapt->data;
871c3436cfSJed Brown   PetscErrorCode ierr;
881c3436cfSJed Brown   PetscInt       two;
891c3436cfSJed Brown   PetscBool      set;
901c3436cfSJed Brown 
911c3436cfSJed Brown   PetscFunctionBegin;
92e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Basic adaptive controller options");CHKERRQ(ierr);
931c3436cfSJed Brown   two  = 2;
941c3436cfSJed Brown   ierr = PetscOptionsRealArray("-ts_adapt_basic_clip","Admissible decrease/increase in step size","",basic->clip,&two,&set);CHKERRQ(ierr);
95ce94432eSBarry 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");
960298fd71SBarry Smith   ierr = PetscOptionsReal("-ts_adapt_basic_safety","Safety factor relative to target error","",basic->safety,&basic->safety,NULL);CHKERRQ(ierr);
970298fd71SBarry 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);
980298fd71SBarry 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);
991c3436cfSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
1001c3436cfSJed Brown   PetscFunctionReturn(0);
1011c3436cfSJed Brown }
1021c3436cfSJed Brown 
10339850373SJed Brown #undef __FUNCT__
10439850373SJed Brown #define __FUNCT__ "TSAdaptView_Basic"
10539850373SJed Brown static PetscErrorCode TSAdaptView_Basic(TSAdapt adapt,PetscViewer viewer)
10639850373SJed Brown {
10739850373SJed Brown   TSAdapt_Basic  *basic = (TSAdapt_Basic*)adapt->data;
10839850373SJed Brown   PetscErrorCode ierr;
10939850373SJed Brown   PetscBool      iascii;
11039850373SJed Brown 
11139850373SJed Brown   PetscFunctionBegin;
11239850373SJed Brown   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
11339850373SJed Brown   if (iascii) {
11439850373SJed Brown     if (basic->always_accept) {ierr = PetscViewerASCIIPrintf(viewer,"  Basic: always accepting steps\n");CHKERRQ(ierr);}
11557622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Basic: clip fastest decrease %g, fastest increase %g\n",(double)basic->clip[0],(double)basic->clip[1]);CHKERRQ(ierr);
11657622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Basic: safety factor %g, extra factor after step rejection %g\n",(double)basic->safety,(double)basic->reject_safety);CHKERRQ(ierr);
11739850373SJed Brown   }
11839850373SJed Brown   PetscFunctionReturn(0);
11939850373SJed Brown }
12039850373SJed Brown 
12184df9cb4SJed Brown #undef __FUNCT__
12284df9cb4SJed Brown #define __FUNCT__ "TSAdaptCreate_Basic"
12384df9cb4SJed Brown /*MC
12484df9cb4SJed Brown    TSADAPTBASIC - Basic adaptive controller for time stepping
12584df9cb4SJed Brown 
12684df9cb4SJed Brown    Level: intermediate
12784df9cb4SJed Brown 
12884df9cb4SJed Brown .seealso: TS, TSAdapt, TSSetAdapt()
12984df9cb4SJed Brown M*/
1308cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt adapt)
13184df9cb4SJed Brown {
13284df9cb4SJed Brown   PetscErrorCode ierr;
13384df9cb4SJed Brown   TSAdapt_Basic  *a;
13484df9cb4SJed Brown 
13584df9cb4SJed Brown   PetscFunctionBegin;
136b00a9115SJed Brown   ierr                       = PetscNewLog(adapt,&a);CHKERRQ(ierr);
13784df9cb4SJed Brown   adapt->data                = (void*)a;
1381c3436cfSJed Brown   adapt->ops->choose         = TSAdaptChoose_Basic;
1391c3436cfSJed Brown   adapt->ops->setfromoptions = TSAdaptSetFromOptions_Basic;
14084df9cb4SJed Brown   adapt->ops->destroy        = TSAdaptDestroy_Basic;
14139850373SJed Brown   adapt->ops->view           = TSAdaptView_Basic;
1421c3436cfSJed Brown 
1431c3436cfSJed Brown   a->clip[0]       = 0.1;
1441c3436cfSJed Brown   a->clip[1]       = 10.;
1451c3436cfSJed Brown   a->safety        = 0.9;
1467e1ba4dcSJed Brown   a->reject_safety = 0.5;
1477e1ba4dcSJed Brown   a->always_accept = PETSC_FALSE;
14884df9cb4SJed Brown   PetscFunctionReturn(0);
14984df9cb4SJed Brown }
150