xref: /petsc/src/ts/adapt/impls/basic/adaptbasic.c (revision 1566a47f5f7ddc4d41547989035e2da0a3965616)
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;
16*1566a47fSLisandro Dalcin   PetscInt       order;
17*1566a47fSLisandro Dalcin   PetscReal      safety = basic->safety;
18*1566a47fSLisandro Dalcin   PetscReal      enorm=-1,hfac_lte,h_lte;
191c3436cfSJed Brown   PetscErrorCode ierr;
2084df9cb4SJed Brown 
2184df9cb4SJed Brown   PetscFunctionBegin;
22*1566a47fSLisandro Dalcin   *next_sc = 0; /* Reuse the same order scheme */
23*1566a47fSLisandro Dalcin 
241ebf93c6SLisandro Dalcin   if (ts->ops->evaluatewlte) {
25*1566a47fSLisandro Dalcin     order = PETSC_DECIDE;
26*1566a47fSLisandro Dalcin     ierr = TSEvaluateWLTE(ts,adapt->wnormtype,&order,&enorm);CHKERRQ(ierr);
27*1566a47fSLisandro Dalcin     if (enorm >= 0 && order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed error order %D must be positive",order);
28*1566a47fSLisandro Dalcin   } else if (ts->ops->evaluatestep) {
29*1566a47fSLisandro Dalcin     if (adapt->candidates.n < 1) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"No candidate has been registered");
30*1566a47fSLisandro 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);
31*1566a47fSLisandro Dalcin     if (!basic->Y) {ierr = VecDuplicate(ts->vec_sol,&basic->Y);CHKERRQ(ierr);}
32*1566a47fSLisandro Dalcin     order = adapt->candidates.order[0];
331ebf93c6SLisandro Dalcin     ierr = TSEvaluateStep(ts,order-1,basic->Y,NULL);CHKERRQ(ierr);
34*1566a47fSLisandro Dalcin     ierr = TSErrorWeightedNorm(ts,ts->vec_sol,basic->Y,adapt->wnormtype,&enorm);CHKERRQ(ierr);
351ebf93c6SLisandro Dalcin   }
361c3436cfSJed Brown 
37*1566a47fSLisandro Dalcin   if (enorm < 0) {
38*1566a47fSLisandro Dalcin     *accept  = PETSC_TRUE;
39*1566a47fSLisandro Dalcin     *next_h  = h;            /* Reuse the old step */
40*1566a47fSLisandro Dalcin     *wlte    = -1;           /* Weighted local truncation error was not evaluated */
41*1566a47fSLisandro Dalcin     PetscFunctionReturn(0);
42*1566a47fSLisandro Dalcin   }
43*1566a47fSLisandro Dalcin 
44*1566a47fSLisandro Dalcin   /* Determine whether the step is accepted of rejected */
45*1566a47fSLisandro Dalcin   if (enorm > 1) {
467e1ba4dcSJed Brown     if (!*accept) safety *= basic->reject_safety; /* The last attempt also failed, shorten more aggressively */
47fd94acc0SJed Brown     if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) {
4857622a8eSBarry 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);
49fd94acc0SJed Brown       *accept = PETSC_TRUE;
507e1ba4dcSJed Brown     } else if (basic->always_accept) {
5157622a8eSBarry 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);
527e1ba4dcSJed Brown       *accept = PETSC_TRUE;
531c3436cfSJed Brown     } else {
5457622a8eSBarry Smith       ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, rejecting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr);
55fd94acc0SJed Brown       *accept = PETSC_FALSE;
56fd94acc0SJed Brown     }
57fd94acc0SJed Brown   } else {
5857622a8eSBarry Smith     ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr);
59fd94acc0SJed Brown     *accept = PETSC_TRUE;
601c3436cfSJed Brown   }
611c3436cfSJed Brown 
621c3436cfSJed Brown   /* The optimal new step based purely on local truncation error for this step. */
63*1566a47fSLisandro Dalcin   if (enorm > 0)
64*1566a47fSLisandro Dalcin     hfac_lte = safety * PetscPowReal(enorm,((PetscReal)-1)/order);
65*1566a47fSLisandro Dalcin   else
66bf0c4ff7SBarry Smith     hfac_lte = safety * PETSC_INFINITY;
671c3436cfSJed Brown   h_lte = h * PetscClipInterval(hfac_lte,basic->clip[0],basic->clip[1]);
681c3436cfSJed Brown 
691c3436cfSJed Brown   *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max);
700b99f514SJed Brown   *wlte   = enorm;
7184df9cb4SJed Brown   PetscFunctionReturn(0);
7284df9cb4SJed Brown }
7384df9cb4SJed Brown 
7484df9cb4SJed Brown #undef __FUNCT__
75099af0a3SLisandro Dalcin #define __FUNCT__ "TSAdaptReset_Basic"
76099af0a3SLisandro Dalcin static PetscErrorCode TSAdaptReset_Basic(TSAdapt adapt)
7784df9cb4SJed Brown {
781c3436cfSJed Brown   TSAdapt_Basic  *basic = (TSAdapt_Basic*)adapt->data;
7984df9cb4SJed Brown   PetscErrorCode ierr;
8084df9cb4SJed Brown 
8184df9cb4SJed Brown   PetscFunctionBegin;
821c3436cfSJed Brown   ierr = VecDestroy(&basic->Y);CHKERRQ(ierr);
83099af0a3SLisandro Dalcin   PetscFunctionReturn(0);
84099af0a3SLisandro Dalcin }
85099af0a3SLisandro Dalcin 
86099af0a3SLisandro Dalcin #undef __FUNCT__
87099af0a3SLisandro Dalcin #define __FUNCT__ "TSAdaptDestroy_Basic"
88099af0a3SLisandro Dalcin static PetscErrorCode TSAdaptDestroy_Basic(TSAdapt adapt)
89099af0a3SLisandro Dalcin {
90099af0a3SLisandro Dalcin   PetscErrorCode ierr;
91099af0a3SLisandro Dalcin 
92099af0a3SLisandro Dalcin   PetscFunctionBegin;
93099af0a3SLisandro Dalcin   ierr = TSAdaptReset_Basic(adapt);CHKERRQ(ierr);
9484df9cb4SJed Brown   ierr = PetscFree(adapt->data);CHKERRQ(ierr);
9584df9cb4SJed Brown   PetscFunctionReturn(0);
9684df9cb4SJed Brown }
9784df9cb4SJed Brown 
981c3436cfSJed Brown #undef __FUNCT__
991c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetFromOptions_Basic"
1004416b707SBarry Smith static PetscErrorCode TSAdaptSetFromOptions_Basic(PetscOptionItems *PetscOptionsObject,TSAdapt adapt)
1011c3436cfSJed Brown {
1021c3436cfSJed Brown   TSAdapt_Basic  *basic = (TSAdapt_Basic*)adapt->data;
1031c3436cfSJed Brown   PetscErrorCode ierr;
1041c3436cfSJed Brown   PetscInt       two;
1051c3436cfSJed Brown   PetscBool      set;
1061c3436cfSJed Brown 
1071c3436cfSJed Brown   PetscFunctionBegin;
108e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Basic adaptive controller options");CHKERRQ(ierr);
1091c3436cfSJed Brown   two  = 2;
1101c3436cfSJed Brown   ierr = PetscOptionsRealArray("-ts_adapt_basic_clip","Admissible decrease/increase in step size","",basic->clip,&two,&set);CHKERRQ(ierr);
111ce94432eSBarry 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");
1120298fd71SBarry Smith   ierr = PetscOptionsReal("-ts_adapt_basic_safety","Safety factor relative to target error","",basic->safety,&basic->safety,NULL);CHKERRQ(ierr);
1130298fd71SBarry 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);
1140298fd71SBarry 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);
1151c3436cfSJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
1161c3436cfSJed Brown   PetscFunctionReturn(0);
1171c3436cfSJed Brown }
1181c3436cfSJed Brown 
11939850373SJed Brown #undef __FUNCT__
12039850373SJed Brown #define __FUNCT__ "TSAdaptView_Basic"
12139850373SJed Brown static PetscErrorCode TSAdaptView_Basic(TSAdapt adapt,PetscViewer viewer)
12239850373SJed Brown {
12339850373SJed Brown   TSAdapt_Basic  *basic = (TSAdapt_Basic*)adapt->data;
12439850373SJed Brown   PetscErrorCode ierr;
12539850373SJed Brown   PetscBool      iascii;
12639850373SJed Brown 
12739850373SJed Brown   PetscFunctionBegin;
12839850373SJed Brown   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
12939850373SJed Brown   if (iascii) {
13039850373SJed Brown     if (basic->always_accept) {ierr = PetscViewerASCIIPrintf(viewer,"  Basic: always accepting steps\n");CHKERRQ(ierr);}
13157622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Basic: clip fastest decrease %g, fastest increase %g\n",(double)basic->clip[0],(double)basic->clip[1]);CHKERRQ(ierr);
13257622a8eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Basic: safety factor %g, extra factor after step rejection %g\n",(double)basic->safety,(double)basic->reject_safety);CHKERRQ(ierr);
13339850373SJed Brown   }
13439850373SJed Brown   PetscFunctionReturn(0);
13539850373SJed Brown }
13639850373SJed Brown 
13784df9cb4SJed Brown #undef __FUNCT__
13884df9cb4SJed Brown #define __FUNCT__ "TSAdaptCreate_Basic"
13984df9cb4SJed Brown /*MC
14084df9cb4SJed Brown    TSADAPTBASIC - Basic adaptive controller for time stepping
14184df9cb4SJed Brown 
14284df9cb4SJed Brown    Level: intermediate
14384df9cb4SJed Brown 
14484df9cb4SJed Brown .seealso: TS, TSAdapt, TSSetAdapt()
14584df9cb4SJed Brown M*/
1468cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt adapt)
14784df9cb4SJed Brown {
14884df9cb4SJed Brown   PetscErrorCode ierr;
14984df9cb4SJed Brown   TSAdapt_Basic  *a;
15084df9cb4SJed Brown 
15184df9cb4SJed Brown   PetscFunctionBegin;
152b00a9115SJed Brown   ierr                       = PetscNewLog(adapt,&a);CHKERRQ(ierr);
15384df9cb4SJed Brown   adapt->data                = (void*)a;
1541c3436cfSJed Brown   adapt->ops->choose         = TSAdaptChoose_Basic;
1551c3436cfSJed Brown   adapt->ops->setfromoptions = TSAdaptSetFromOptions_Basic;
15684df9cb4SJed Brown   adapt->ops->destroy        = TSAdaptDestroy_Basic;
15739850373SJed Brown   adapt->ops->view           = TSAdaptView_Basic;
1581c3436cfSJed Brown 
1591c3436cfSJed Brown   a->clip[0]       = 0.1;
1601c3436cfSJed Brown   a->clip[1]       = 10.;
1611c3436cfSJed Brown   a->safety        = 0.9;
1627e1ba4dcSJed Brown   a->reject_safety = 0.5;
1637e1ba4dcSJed Brown   a->always_accept = PETSC_FALSE;
16484df9cb4SJed Brown   PetscFunctionReturn(0);
16584df9cb4SJed Brown }
166