xref: /petsc/src/ts/adapt/impls/basic/adaptbasic.c (revision fe998a80077c9ee0917a39496df43fc256e1b478)
1 #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/
2 
3 typedef struct {
4   PetscBool always_accept;
5   PetscReal clip[2];            /* admissible decrease/increase factors */
6   PetscReal safety;             /* safety factor relative to target error */
7   PetscReal reject_safety;      /* extra safety factor if the last step was rejected */
8   Vec       Y;
9 } TSAdapt_Basic;
10 
11 #undef __FUNCT__
12 #define __FUNCT__ "TSAdaptChoose_Basic"
13 static PetscErrorCode TSAdaptChoose_Basic(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte)
14 {
15   TSAdapt_Basic  *basic = (TSAdapt_Basic*)adapt->data;
16   PetscErrorCode ierr;
17   Vec            X,Y;
18   PetscReal      enorm,hfac_lte,h_lte,safety;
19   PetscInt       order,stepno;
20 
21   PetscFunctionBegin;
22   ierr = TSGetTimeStepNumber(ts,&stepno);CHKERRQ(ierr);
23   ierr = TSGetSolution(ts,&X);CHKERRQ(ierr);
24   if (!basic->Y) {ierr = VecDuplicate(X,&basic->Y);CHKERRQ(ierr);}
25   Y     = basic->Y;
26   order = adapt->candidates.order[0];
27   ierr  = TSEvaluateStep(ts,order-1,Y,NULL);CHKERRQ(ierr);
28 
29   safety = basic->safety;
30   ierr   = TSErrorNormWRMS(ts,Y,&enorm);CHKERRQ(ierr);
31   if (enorm > 1.) {
32     if (!*accept) safety *= basic->reject_safety; /* The last attempt also failed, shorten more aggressively */
33     if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) {
34       ierr    = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting because step size %g is at minimum\n",(double)enorm,(double)h);CHKERRQ(ierr);
35       *accept = PETSC_TRUE;
36     } else if (basic->always_accept) {
37       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);
38       *accept = PETSC_TRUE;
39     } else {
40       ierr    = PetscInfo2(adapt,"Estimated scaled local truncation error %g, rejecting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr);
41       *accept = PETSC_FALSE;
42     }
43   } else {
44     ierr    = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr);
45     *accept = PETSC_TRUE;
46   }
47 
48   /* The optimal new step based purely on local truncation error for this step. */
49   hfac_lte = safety * PetscRealPart(PetscPowScalar((PetscScalar)enorm,(PetscReal)(-1./order)));
50   h_lte    = h * PetscClipInterval(hfac_lte,basic->clip[0],basic->clip[1]);
51 
52   *next_sc = 0;
53   *next_h  = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max);
54   *wlte    = enorm;
55   PetscFunctionReturn(0);
56 }
57 
58 #undef __FUNCT__
59 #define __FUNCT__ "TSAdaptDestroy_Basic"
60 static PetscErrorCode TSAdaptDestroy_Basic(TSAdapt adapt)
61 {
62   TSAdapt_Basic  *basic = (TSAdapt_Basic*)adapt->data;
63   PetscErrorCode ierr;
64 
65   PetscFunctionBegin;
66   ierr = VecDestroy(&basic->Y);CHKERRQ(ierr);
67   ierr = PetscFree(adapt->data);CHKERRQ(ierr);
68   PetscFunctionReturn(0);
69 }
70 
71 #undef __FUNCT__
72 #define __FUNCT__ "TSAdaptSetFromOptions_Basic"
73 static PetscErrorCode TSAdaptSetFromOptions_Basic(TSAdapt adapt)
74 {
75   TSAdapt_Basic  *basic = (TSAdapt_Basic*)adapt->data;
76   PetscErrorCode ierr;
77   PetscInt       two;
78   PetscBool      set;
79 
80   PetscFunctionBegin;
81   ierr = PetscOptionsHead("Basic adaptive controller options");CHKERRQ(ierr);
82   two  = 2;
83   ierr = PetscOptionsRealArray("-ts_adapt_basic_clip","Admissible decrease/increase in step size","",basic->clip,&two,&set);CHKERRQ(ierr);
84   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");
85   ierr = PetscOptionsReal("-ts_adapt_basic_safety","Safety factor relative to target error","",basic->safety,&basic->safety,NULL);CHKERRQ(ierr);
86   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);
87   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);
88   ierr = PetscOptionsTail();CHKERRQ(ierr);
89   PetscFunctionReturn(0);
90 }
91 
92 #undef __FUNCT__
93 #define __FUNCT__ "TSAdaptView_Basic"
94 static PetscErrorCode TSAdaptView_Basic(TSAdapt adapt,PetscViewer viewer)
95 {
96   TSAdapt_Basic  *basic = (TSAdapt_Basic*)adapt->data;
97   PetscErrorCode ierr;
98   PetscBool      iascii;
99 
100   PetscFunctionBegin;
101   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
102   if (iascii) {
103     if (basic->always_accept) {ierr = PetscViewerASCIIPrintf(viewer,"  Basic: always accepting steps\n");CHKERRQ(ierr);}
104     ierr = PetscViewerASCIIPrintf(viewer,"  Basic: clip fastest decrease %g, fastest increase %g\n",(double)basic->clip[0],(double)basic->clip[1]);CHKERRQ(ierr);
105     ierr = PetscViewerASCIIPrintf(viewer,"  Basic: safety factor %g, extra factor after step rejection %g\n",(double)basic->safety,(double)basic->reject_safety);CHKERRQ(ierr);
106   }
107   PetscFunctionReturn(0);
108 }
109 
110 #undef __FUNCT__
111 #define __FUNCT__ "TSAdaptCreate_Basic"
112 /*MC
113    TSADAPTBASIC - Basic adaptive controller for time stepping
114 
115    Level: intermediate
116 
117 .seealso: TS, TSAdapt, TSSetAdapt()
118 M*/
119 PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt adapt)
120 {
121   PetscErrorCode ierr;
122   TSAdapt_Basic  *a;
123 
124   PetscFunctionBegin;
125   ierr                       = PetscNewLog(adapt,&a);CHKERRQ(ierr);
126   adapt->data                = (void*)a;
127   adapt->ops->choose         = TSAdaptChoose_Basic;
128   adapt->ops->setfromoptions = TSAdaptSetFromOptions_Basic;
129   adapt->ops->destroy        = TSAdaptDestroy_Basic;
130   adapt->ops->view           = TSAdaptView_Basic;
131 
132   a->clip[0]       = 0.1;
133   a->clip[1]       = 10.;
134   a->safety        = 0.9;
135   a->reject_safety = 0.5;
136   a->always_accept = PETSC_FALSE;
137   PetscFunctionReturn(0);
138 }
139