xref: /petsc/src/ts/adapt/impls/dsp/adaptdsp.c (revision 95452b02e12c0ee11232c7ff2b24b568a8e07e43)
1cb7ab186SLisandro Dalcin #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
2cb7ab186SLisandro Dalcin 
3cb7ab186SLisandro Dalcin static const char *citation[] = {
4cb7ab186SLisandro Dalcin   "@article{Soderlind2003,\n"
5cb7ab186SLisandro Dalcin   " author = {S\"{o}derlind, Gustaf},\n"
6cb7ab186SLisandro Dalcin   " title = {Digital Filters in Adaptive Time-stepping},\n"
7cb7ab186SLisandro Dalcin   " journal = {ACM Transactions on Mathematical Software},\n"
8cb7ab186SLisandro Dalcin   " volume = {29},\n"
9cb7ab186SLisandro Dalcin   " number = {1},\n"
10cb7ab186SLisandro Dalcin   " pages = {1--26},\n"
11cb7ab186SLisandro Dalcin   " year = {2003},\n"
12cb7ab186SLisandro Dalcin   " issn = {0098-3500},\n"
13cb7ab186SLisandro Dalcin   " doi = {http://dx.doi.org/10.1145/641876.641877},\n"
14cb7ab186SLisandro Dalcin   "}\n",
15cb7ab186SLisandro Dalcin   "@article{Soderlind2006,\n"
16cb7ab186SLisandro Dalcin   " author = {Gustaf S\"{o}derlind and Lina Wang},\n"
17cb7ab186SLisandro Dalcin   " title = {Adaptive time-stepping and computational stability},\n"
18cb7ab186SLisandro Dalcin   " journal = {Journal of Computational and Applied Mathematics},\n"
19cb7ab186SLisandro Dalcin   " volume = {185},\n"
20cb7ab186SLisandro Dalcin   " number = {2},\n"
21cb7ab186SLisandro Dalcin   " pages = {225--243},\n"
22cb7ab186SLisandro Dalcin   " year = {2006},\n"
23cb7ab186SLisandro Dalcin   " issn = {0377-0427},\n"
24cb7ab186SLisandro Dalcin   " doi = {http://dx.doi.org/10.1016/j.cam.2005.03.008},\n"
25cb7ab186SLisandro Dalcin   "}\n",
26cb7ab186SLisandro Dalcin };
27cb7ab186SLisandro Dalcin static PetscBool cited[] = {PETSC_FALSE,PETSC_FALSE};
28cb7ab186SLisandro Dalcin 
29cb7ab186SLisandro Dalcin typedef struct {
30cb7ab186SLisandro Dalcin   Vec       Y;
31cb7ab186SLisandro Dalcin   PetscReal kBeta[3];  /* filter parameters */
32cb7ab186SLisandro Dalcin   PetscReal Alpha[2];  /* filter parameters */
33cb7ab186SLisandro Dalcin   PetscReal cerror[3]; /* control error (controller input) history */
34cb7ab186SLisandro Dalcin   PetscReal hratio[3]; /* stepsize ratio (controller output) history */
35cb7ab186SLisandro Dalcin   PetscBool rollback;
36cb7ab186SLisandro Dalcin } TSAdapt_DSP;
37cb7ab186SLisandro Dalcin 
38cb7ab186SLisandro Dalcin static PetscReal Limiter(PetscReal value,PetscReal kappa)
39cb7ab186SLisandro Dalcin {
40cb7ab186SLisandro Dalcin   return 1 + kappa*PetscAtanReal((value - 1)/kappa);
41cb7ab186SLisandro Dalcin }
42cb7ab186SLisandro Dalcin 
43cb7ab186SLisandro Dalcin static PetscErrorCode TSAdaptRestart_DSP(TSAdapt adapt)
44cb7ab186SLisandro Dalcin {
45cb7ab186SLisandro Dalcin   TSAdapt_DSP *dsp = (TSAdapt_DSP*)adapt->data;
46cb7ab186SLisandro Dalcin   PetscFunctionBegin;
47cb7ab186SLisandro Dalcin   dsp->cerror[0] = dsp->hratio[0] = 1.0;
48cb7ab186SLisandro Dalcin   dsp->cerror[1] = dsp->hratio[1] = 1.0;
49cb7ab186SLisandro Dalcin   dsp->cerror[2] = dsp->hratio[2] = 1.0;
50cb7ab186SLisandro Dalcin   PetscFunctionReturn(0);
51cb7ab186SLisandro Dalcin }
52cb7ab186SLisandro Dalcin 
53cb7ab186SLisandro Dalcin static PetscErrorCode TSAdaptRollBack_DSP(TSAdapt adapt)
54cb7ab186SLisandro Dalcin {
55cb7ab186SLisandro Dalcin   TSAdapt_DSP *dsp = (TSAdapt_DSP*)adapt->data;
56cb7ab186SLisandro Dalcin   PetscFunctionBegin;
57cb7ab186SLisandro Dalcin   dsp->cerror[0] = dsp->cerror[1];
58cb7ab186SLisandro Dalcin   dsp->cerror[1] = dsp->cerror[2];
59cb7ab186SLisandro Dalcin   dsp->cerror[2] = 1.0;
60cb7ab186SLisandro Dalcin   dsp->hratio[0] = dsp->hratio[1];
61cb7ab186SLisandro Dalcin   dsp->hratio[1] = dsp->hratio[2];
62cb7ab186SLisandro Dalcin   dsp->hratio[2] = 1.0;
63cb7ab186SLisandro Dalcin   PetscFunctionReturn(0);
64cb7ab186SLisandro Dalcin }
65cb7ab186SLisandro Dalcin 
66cb7ab186SLisandro Dalcin static PetscErrorCode TSAdaptChoose_DSP(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte,PetscReal *wltea,PetscReal *wlter)
67cb7ab186SLisandro Dalcin {
68cb7ab186SLisandro Dalcin   TSAdapt_DSP    *dsp = (TSAdapt_DSP*)adapt->data;
69cb7ab186SLisandro Dalcin   PetscInt       order = PETSC_DECIDE;
70cb7ab186SLisandro Dalcin   PetscReal      enorm = -1;
71cb7ab186SLisandro Dalcin   PetscReal      enorma,enormr;
72cb7ab186SLisandro Dalcin   PetscReal      safety = adapt->safety * (PetscReal)0.9;
73cb7ab186SLisandro Dalcin   PetscReal      hnew,hfac = PETSC_INFINITY;
74cb7ab186SLisandro Dalcin   PetscReal      hmin = adapt->dt_min*(1 + PETSC_SQRT_MACHINE_EPSILON);
75cb7ab186SLisandro Dalcin   PetscErrorCode ierr;
76cb7ab186SLisandro Dalcin 
77cb7ab186SLisandro Dalcin   PetscFunctionBegin;
78cb7ab186SLisandro Dalcin   *next_sc = 0;   /* Reuse the same order scheme */
79cb7ab186SLisandro Dalcin   *wltea   = -1;  /* Weighted absolute local truncation error is not used */
80cb7ab186SLisandro Dalcin   *wlter   = -1;  /* Weighted relative local truncation error is not used */
81cb7ab186SLisandro Dalcin 
82cb7ab186SLisandro Dalcin   if (ts->ops->evaluatewlte) {
83cb7ab186SLisandro Dalcin     ierr = TSEvaluateWLTE(ts,adapt->wnormtype,&order,&enorm);CHKERRQ(ierr);
84cb7ab186SLisandro Dalcin     if (enorm >= 0 && order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed error order %D must be positive",order);
85cb7ab186SLisandro Dalcin   } else if (ts->ops->evaluatestep) {
86cb7ab186SLisandro Dalcin     if (adapt->candidates.n < 1) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"No candidate has been registered");
87cb7ab186SLisandro 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);
88cb7ab186SLisandro Dalcin     if (!dsp->Y) {ierr = VecDuplicate(ts->vec_sol,&dsp->Y);CHKERRQ(ierr);}
89cb7ab186SLisandro Dalcin     order = adapt->candidates.order[0];
90cb7ab186SLisandro Dalcin     ierr = TSEvaluateStep(ts,order-1,dsp->Y,NULL);CHKERRQ(ierr);
91cb7ab186SLisandro Dalcin     ierr = TSErrorWeightedNorm(ts,ts->vec_sol,dsp->Y,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr);
92cb7ab186SLisandro Dalcin   }
93cb7ab186SLisandro Dalcin   if (enorm < 0) {
94cb7ab186SLisandro Dalcin     ierr = TSAdaptRestart_DSP(adapt);CHKERRQ(ierr);
95cb7ab186SLisandro Dalcin     *accept = PETSC_TRUE;  /* Accept the step */
96cb7ab186SLisandro Dalcin     *next_h = h;           /* Reuse the old step size */
97cb7ab186SLisandro Dalcin     *wlte   = -1;          /* Weighted local truncation error was not evaluated */
98cb7ab186SLisandro Dalcin     PetscFunctionReturn(0);
99cb7ab186SLisandro Dalcin   }
100cb7ab186SLisandro Dalcin 
101cb7ab186SLisandro Dalcin   ierr = PetscCitationsRegister(citation[0],&cited[0]);CHKERRQ(ierr);
102cb7ab186SLisandro Dalcin   ierr = PetscCitationsRegister(citation[1],&cited[1]);CHKERRQ(ierr);
103cb7ab186SLisandro Dalcin 
104cb7ab186SLisandro Dalcin   /* Update history after rollback */
105cb7ab186SLisandro Dalcin   if (!ts->steprollback)
106cb7ab186SLisandro Dalcin     dsp->rollback = PETSC_FALSE;
107cb7ab186SLisandro Dalcin   else if (!dsp->rollback) {
108cb7ab186SLisandro Dalcin     dsp->rollback = PETSC_TRUE;
109cb7ab186SLisandro Dalcin     ierr = TSAdaptRollBack_DSP(adapt);CHKERRQ(ierr);
110cb7ab186SLisandro Dalcin   }
111cb7ab186SLisandro Dalcin   /* Reset history after restart */
112cb7ab186SLisandro Dalcin   if (ts->steprestart) {
113cb7ab186SLisandro Dalcin     ierr = TSAdaptRestart_DSP(adapt);CHKERRQ(ierr);
114cb7ab186SLisandro Dalcin   }
115cb7ab186SLisandro Dalcin 
116cb7ab186SLisandro Dalcin   {
117cb7ab186SLisandro Dalcin     PetscReal k = (PetscReal)order;
118cb7ab186SLisandro Dalcin     PetscReal b1 = dsp->kBeta[0];
119cb7ab186SLisandro Dalcin     PetscReal b2 = dsp->kBeta[1];
120cb7ab186SLisandro Dalcin     PetscReal b3 = dsp->kBeta[2];
121cb7ab186SLisandro Dalcin     PetscReal a2 = dsp->Alpha[0];
122cb7ab186SLisandro Dalcin     PetscReal a3 = dsp->Alpha[0];
123cb7ab186SLisandro Dalcin 
124cb7ab186SLisandro Dalcin     PetscReal ctr0;
125cb7ab186SLisandro Dalcin     PetscReal ctr1 = dsp->cerror[0];
126cb7ab186SLisandro Dalcin     PetscReal ctr2 = dsp->cerror[1];
127cb7ab186SLisandro Dalcin     PetscReal rho0;
128cb7ab186SLisandro Dalcin     PetscReal rho1 = dsp->hratio[0];
129cb7ab186SLisandro Dalcin     PetscReal rho2 = dsp->hratio[1];
130cb7ab186SLisandro Dalcin 
131cb7ab186SLisandro Dalcin     /* Compute the step size ratio */
132cb7ab186SLisandro Dalcin     enorm = PetscMax(enorm,PETSC_SMALL);
133cb7ab186SLisandro Dalcin     ctr0  = PetscPowReal(1/enorm,1/k);
134cb7ab186SLisandro Dalcin     rho0  = PetscPowReal(ctr0,b1);
135cb7ab186SLisandro Dalcin     rho0 *= PetscPowReal(ctr1,b2);
136cb7ab186SLisandro Dalcin     rho0 *= PetscPowReal(ctr2,b3);
137cb7ab186SLisandro Dalcin     rho0 *= PetscPowReal(rho1,-a2);
138cb7ab186SLisandro Dalcin     rho0 *= PetscPowReal(rho2,-a3);
139cb7ab186SLisandro Dalcin     rho0  = Limiter(rho0,1);
140cb7ab186SLisandro Dalcin 
141cb7ab186SLisandro Dalcin     /* Determine whether the step is accepted or rejected */
142cb7ab186SLisandro Dalcin     if (rho0 >= safety)
143cb7ab186SLisandro Dalcin       *accept = PETSC_TRUE;
144cb7ab186SLisandro Dalcin     else if (adapt->always_accept)
145cb7ab186SLisandro Dalcin       *accept = PETSC_TRUE;
146cb7ab186SLisandro Dalcin     else if (h < hmin)
147cb7ab186SLisandro Dalcin       *accept = PETSC_TRUE;
148cb7ab186SLisandro Dalcin     else
149cb7ab186SLisandro Dalcin       *accept = PETSC_FALSE;
150cb7ab186SLisandro Dalcin 
151cb7ab186SLisandro Dalcin     /* Update history after accept */
152cb7ab186SLisandro Dalcin     if (*accept) {
153cb7ab186SLisandro Dalcin       dsp->cerror[2] = dsp->cerror[1];
154cb7ab186SLisandro Dalcin       dsp->cerror[1] = dsp->cerror[0];
155cb7ab186SLisandro Dalcin       dsp->cerror[0] = ctr0;
156cb7ab186SLisandro Dalcin       dsp->hratio[2] = dsp->hratio[1];
157cb7ab186SLisandro Dalcin       dsp->hratio[1] = dsp->hratio[0];
158cb7ab186SLisandro Dalcin       dsp->hratio[0] = rho0;
159cb7ab186SLisandro Dalcin       dsp->rollback  = PETSC_FALSE;
160cb7ab186SLisandro Dalcin     }
161cb7ab186SLisandro Dalcin 
162cb7ab186SLisandro Dalcin     hfac = rho0;
163cb7ab186SLisandro Dalcin   }
164cb7ab186SLisandro Dalcin 
165cb7ab186SLisandro Dalcin   hnew    = h * PetscClipInterval(hfac,adapt->clip[0],adapt->clip[1]);
166cb7ab186SLisandro Dalcin   *next_h = PetscClipInterval(hnew,adapt->dt_min,adapt->dt_max);
167cb7ab186SLisandro Dalcin   *wlte   = enorm;
168cb7ab186SLisandro Dalcin   PetscFunctionReturn(0);
169cb7ab186SLisandro Dalcin }
170cb7ab186SLisandro Dalcin 
171cb7ab186SLisandro Dalcin static PetscErrorCode TSAdaptReset_DSP(TSAdapt adapt)
172cb7ab186SLisandro Dalcin {
173cb7ab186SLisandro Dalcin   TSAdapt_DSP    *dsp = (TSAdapt_DSP*)adapt->data;
174cb7ab186SLisandro Dalcin   PetscErrorCode ierr;
175cb7ab186SLisandro Dalcin 
176cb7ab186SLisandro Dalcin   PetscFunctionBegin;
177cb7ab186SLisandro Dalcin   ierr = VecDestroy(&dsp->Y);CHKERRQ(ierr);
178cb7ab186SLisandro Dalcin   PetscFunctionReturn(0);
179cb7ab186SLisandro Dalcin }
180cb7ab186SLisandro Dalcin 
181cb7ab186SLisandro Dalcin static PetscErrorCode TSAdaptDestroy_DSP(TSAdapt adapt)
182cb7ab186SLisandro Dalcin {
183cb7ab186SLisandro Dalcin   PetscErrorCode ierr;
184cb7ab186SLisandro Dalcin 
185cb7ab186SLisandro Dalcin   PetscFunctionBegin;
186cb7ab186SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)adapt,"TSAdaptDSPSetFilter_C",NULL);CHKERRQ(ierr);
187cb7ab186SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)adapt,"TSAdaptDSPSetPID_C",NULL);CHKERRQ(ierr);
188cb7ab186SLisandro Dalcin   ierr = TSAdaptReset_DSP(adapt);CHKERRQ(ierr);
189cb7ab186SLisandro Dalcin   ierr = PetscFree(adapt->data);CHKERRQ(ierr);
190cb7ab186SLisandro Dalcin   PetscFunctionReturn(0);
191cb7ab186SLisandro Dalcin }
192cb7ab186SLisandro Dalcin 
193cb7ab186SLisandro Dalcin static PetscErrorCode TSAdaptView_DSP(TSAdapt adapt,PetscViewer viewer)
194cb7ab186SLisandro Dalcin {
195cb7ab186SLisandro Dalcin   TSAdapt_DSP    *dsp = (TSAdapt_DSP*)adapt->data;
196cb7ab186SLisandro Dalcin   PetscBool      iascii;
197cb7ab186SLisandro Dalcin   PetscErrorCode ierr;
198cb7ab186SLisandro Dalcin 
199cb7ab186SLisandro Dalcin   PetscFunctionBegin;
200cb7ab186SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
201cb7ab186SLisandro Dalcin   if (iascii) {
202cb7ab186SLisandro Dalcin     double a2 = (double)dsp->Alpha[0], a3 = (double)dsp->Alpha[1];
203cb7ab186SLisandro Dalcin     double b1 = (double)dsp->kBeta[0], b2 = (double)dsp->kBeta[1], b3 = (double)dsp->kBeta[2];
204cb7ab186SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"filter parameters kBeta=[%g,%g,%g] Alpha=[%g,%g]\n",b1,b2,b3,a2,a3);CHKERRQ(ierr);
205cb7ab186SLisandro Dalcin   }
206cb7ab186SLisandro Dalcin   PetscFunctionReturn(0);
207cb7ab186SLisandro Dalcin }
208cb7ab186SLisandro Dalcin 
209cb7ab186SLisandro Dalcin struct FilterTab {
210cb7ab186SLisandro Dalcin   const char *name;
211cb7ab186SLisandro Dalcin   PetscReal scale;
212cb7ab186SLisandro Dalcin   PetscReal kBeta[3];
213cb7ab186SLisandro Dalcin   PetscReal Alpha[2];
214cb7ab186SLisandro Dalcin };
215cb7ab186SLisandro Dalcin 
216cb7ab186SLisandro Dalcin static struct FilterTab filterlist[] = {
217cb7ab186SLisandro Dalcin   {"basic",    1, {  1,  0,  0 }, {   0,  0 }},
218cb7ab186SLisandro Dalcin 
219cb7ab186SLisandro Dalcin   {"PI30",     3, {  1,  0,  0 }, {   0,  0 }},
220cb7ab186SLisandro Dalcin   {"PI42",     5, {  3, -1,  0 }, {   0,  0 }},
221cb7ab186SLisandro Dalcin   {"PI33",     3, {  2, -1,  0 }, {   0,  0 }},
222cb7ab186SLisandro Dalcin   {"PI34",    10, {  7, -4,  0 }, {   0,  0 }},
223cb7ab186SLisandro Dalcin 
224cb7ab186SLisandro Dalcin   {"PC11",     1, {  2, -1,  0 }, {  -1,  0 }},
225cb7ab186SLisandro Dalcin   {"PC47",    10, { 11, -7,  0 }, { -10,  0 }},
226cb7ab186SLisandro Dalcin   {"PC36",    10, {  9, -6,  0 }, { -10,  0 }},
227cb7ab186SLisandro Dalcin 
228cb7ab186SLisandro Dalcin   {"H0211",    2, {  1,  1,  0 }, {   1,  0 }},
229cb7ab186SLisandro Dalcin   {"H211b",    4, {  1,  1,  0 }, {   1,  0 }},
230cb7ab186SLisandro Dalcin   {"H211PI",   6, {  1,  1,  0 }, {   0,  0 }},
231cb7ab186SLisandro Dalcin 
232cb7ab186SLisandro Dalcin   {"H0312",    4, {  1,  2,  1 }, {   3,  1 }},
233cb7ab186SLisandro Dalcin   {"H312b",    8, {  1,  2,  1 }, {   3,  1 }},
234cb7ab186SLisandro Dalcin   {"H312PID", 18, {  1,  2,  1 }, {   0,  0 }},
235cb7ab186SLisandro Dalcin 
236cb7ab186SLisandro Dalcin   {"H0321",    4, {  5,  2,- 3 }, {  -1, -3 }},
237cb7ab186SLisandro Dalcin   {"H321",    18, {  6,  1,- 5 }, { -15, -3 }},
238cb7ab186SLisandro Dalcin };
239cb7ab186SLisandro Dalcin 
240cb7ab186SLisandro Dalcin static PetscErrorCode TSAdaptDSPSetFilter_DSP(TSAdapt adapt,const char *name)
241cb7ab186SLisandro Dalcin {
242cb7ab186SLisandro Dalcin   TSAdapt_DSP       *dsp = (TSAdapt_DSP*)adapt->data;
243cb7ab186SLisandro Dalcin   PetscInt          i,count = (PetscInt)(sizeof(filterlist)/sizeof(filterlist[0]));
244cb7ab186SLisandro Dalcin   struct FilterTab* tab = NULL;
245cb7ab186SLisandro Dalcin   PetscBool         match;
246cb7ab186SLisandro Dalcin   PetscErrorCode    ierr;
247cb7ab186SLisandro Dalcin 
248cb7ab186SLisandro Dalcin   PetscFunctionBegin;
249cb7ab186SLisandro Dalcin   for (i=0; i<count; i++) {
250cb7ab186SLisandro Dalcin     ierr = PetscStrcasecmp(name,filterlist[i].name,&match);CHKERRQ(ierr);
251cb7ab186SLisandro Dalcin     if (match) { tab = &filterlist[i]; break; }
252cb7ab186SLisandro Dalcin   }
253cb7ab186SLisandro Dalcin   if (!tab) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_UNKNOWN_TYPE,"Filter name %s not found",name);
254cb7ab186SLisandro Dalcin   dsp->kBeta[0] = tab->kBeta[0]/tab->scale;
255cb7ab186SLisandro Dalcin   dsp->kBeta[1] = tab->kBeta[1]/tab->scale;
256cb7ab186SLisandro Dalcin   dsp->kBeta[2] = tab->kBeta[2]/tab->scale;
257cb7ab186SLisandro Dalcin   dsp->Alpha[0] = tab->Alpha[0]/tab->scale;
258cb7ab186SLisandro Dalcin   dsp->Alpha[1] = tab->Alpha[1]/tab->scale;
259cb7ab186SLisandro Dalcin   PetscFunctionReturn(0);
260cb7ab186SLisandro Dalcin }
261cb7ab186SLisandro Dalcin 
262cb7ab186SLisandro Dalcin static PetscErrorCode TSAdaptDSPSetPID_DSP(TSAdapt adapt,PetscReal kkI,PetscReal kkP,PetscReal kkD)
263cb7ab186SLisandro Dalcin {
264cb7ab186SLisandro Dalcin   TSAdapt_DSP    *dsp = (TSAdapt_DSP*)adapt->data;
265cb7ab186SLisandro Dalcin 
266cb7ab186SLisandro Dalcin   PetscFunctionBegin;
267cb7ab186SLisandro Dalcin   dsp->kBeta[0] = kkI + kkP + kkD;
268cb7ab186SLisandro Dalcin   dsp->kBeta[1] = -(kkP + 2*kkD);
269cb7ab186SLisandro Dalcin   dsp->kBeta[2] = kkD;
270cb7ab186SLisandro Dalcin   dsp->Alpha[0] = 0;
271cb7ab186SLisandro Dalcin   dsp->Alpha[1] = 0;
272cb7ab186SLisandro Dalcin   PetscFunctionReturn(0);
273cb7ab186SLisandro Dalcin }
274cb7ab186SLisandro Dalcin 
275cb7ab186SLisandro Dalcin static PetscErrorCode TSAdaptSetFromOptions_DSP(PetscOptionItems *PetscOptionsObject,TSAdapt adapt)
276cb7ab186SLisandro Dalcin {
277cb7ab186SLisandro Dalcin   TSAdapt_DSP    *dsp = (TSAdapt_DSP*)adapt->data;
278cb7ab186SLisandro Dalcin   const char     *names[sizeof(filterlist)/sizeof(filterlist[0])];
279cb7ab186SLisandro Dalcin   PetscInt       count = (PetscInt)(sizeof(filterlist)/sizeof(filterlist[0]));
280cb7ab186SLisandro Dalcin   PetscInt       index = 2; /* PI42 */
281cb7ab186SLisandro Dalcin   PetscReal      pid[3] = {1,0,0};
282cb7ab186SLisandro Dalcin   PetscInt       i,n;
283cb7ab186SLisandro Dalcin   PetscBool      set;
284cb7ab186SLisandro Dalcin   PetscErrorCode ierr;
285cb7ab186SLisandro Dalcin 
286cb7ab186SLisandro Dalcin   PetscFunctionBegin;
287cb7ab186SLisandro Dalcin   for (i=0; i<count; i++) names[i] = filterlist[i].name;
288cb7ab186SLisandro Dalcin   ierr = PetscOptionsHead(PetscOptionsObject,"DSP adaptive controller options");CHKERRQ(ierr);
289cb7ab186SLisandro Dalcin 
290cb7ab186SLisandro Dalcin   ierr = PetscOptionsEList("-ts_adapt_dsp_filter","Filter name","TSAdaptDSPSetFilter",names,count,names[index],&index,&set);CHKERRQ(ierr);
291cb7ab186SLisandro Dalcin   if (set) { ierr = TSAdaptDSPSetFilter(adapt,names[index]);CHKERRQ(ierr);}
292cb7ab186SLisandro Dalcin 
293cb7ab186SLisandro Dalcin   ierr = PetscOptionsRealArray("-ts_adapt_dsp_pid","PID parameters <kkI,kkP,kkD>","TSAdaptDSPSetPID",pid,(n=3,&n),&set);CHKERRQ(ierr);
294cb7ab186SLisandro Dalcin   if (set && !n) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONG,"Must provide at least one value for PID parameters");
295cb7ab186SLisandro Dalcin   if (set) {ierr = TSAdaptDSPSetPID(adapt,pid[0],pid[1],pid[2]);CHKERRQ(ierr);}
296cb7ab186SLisandro Dalcin 
297cb7ab186SLisandro Dalcin   ierr = PetscOptionsRealArray("-ts_adapt_dsp_kbeta","Filter parameters","",dsp->kBeta,(n=3,&n),&set);CHKERRQ(ierr);
298cb7ab186SLisandro Dalcin   if (set && !n) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONG,"Must provide at least one value for parameter kbeta");
299cb7ab186SLisandro Dalcin   if (set) for (i=n; i<3; i++) dsp->kBeta[i] = 0;
300cb7ab186SLisandro Dalcin 
301cb7ab186SLisandro Dalcin   ierr = PetscOptionsRealArray("-ts_adapt_dsp_alpha","Filter parameters","",dsp->Alpha,(n=2,&n),&set);CHKERRQ(ierr);
302cb7ab186SLisandro Dalcin   if (set && !n) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONG,"Must provide at least one value for parameter alpha");
303cb7ab186SLisandro Dalcin   if (set) for (i=n; i<2; i++) dsp->Alpha[i] = 0;
304cb7ab186SLisandro Dalcin 
305cb7ab186SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
306cb7ab186SLisandro Dalcin   PetscFunctionReturn(0);
307cb7ab186SLisandro Dalcin }
308cb7ab186SLisandro Dalcin 
309cb7ab186SLisandro Dalcin /*@C
310cb7ab186SLisandro Dalcin    TSAdaptDSPSetFilter - Sets internal parameters corresponding to the named filter
311cb7ab186SLisandro Dalcin 
312cb7ab186SLisandro Dalcin    Collective on TSAdapt
313cb7ab186SLisandro Dalcin 
314cb7ab186SLisandro Dalcin    Input Arguments:
315cb7ab186SLisandro Dalcin +  adapt - adaptive controller context
316cb7ab186SLisandro Dalcin -  name - filter name
317cb7ab186SLisandro Dalcin 
318cb7ab186SLisandro Dalcin    Level: intermediate
319cb7ab186SLisandro Dalcin 
320cb7ab186SLisandro Dalcin    References: http://dx.doi.org/10.1145/641876.641877
321cb7ab186SLisandro Dalcin 
322*95452b02SPatrick Sanan    Notes:
323*95452b02SPatrick Sanan     Valid filter names are
324cb7ab186SLisandro Dalcin +  "basic" - similar to TSADAPTBASIC but with different criteria for step rejections.
325cb7ab186SLisandro Dalcin .  "PI30", "PI42", "PI33", "PI34" - PI controlers.
326cb7ab186SLisandro Dalcin .  "PC11", "PC47", "PC36" - predictive controllers.
327cb7ab186SLisandro Dalcin .  "H0211", "H211b", "H211PI" - digital filters with orders dynamics=2, adaptivity=1, filter=1.
328cb7ab186SLisandro Dalcin .  "H0312", "H312b", "H312PID" - digital filters with orders dynamics=3, adaptivity=1, filter=2.
329cb7ab186SLisandro Dalcin -  "H0321", "H321" - digital filters with orders dynamics=3, adaptivity=2, filter=1.
330cb7ab186SLisandro Dalcin 
331cb7ab186SLisandro Dalcin    Options Database:
332cb7ab186SLisandro Dalcin .   -ts_adapt_dsp_filter <name>
333cb7ab186SLisandro Dalcin 
334cb7ab186SLisandro Dalcin .seealso: TS, TSAdapt, TSGetAdapt(), TSAdaptDSPSetFilter()
335cb7ab186SLisandro Dalcin @*/
336cb7ab186SLisandro Dalcin PetscErrorCode TSAdaptDSPSetFilter(TSAdapt adapt,const char *name)
337cb7ab186SLisandro Dalcin {
338cb7ab186SLisandro Dalcin   PetscErrorCode ierr;
339cb7ab186SLisandro Dalcin   PetscFunctionBegin;
340cb7ab186SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
341cb7ab186SLisandro Dalcin   PetscValidCharPointer(name,1);
342cb7ab186SLisandro Dalcin   ierr = PetscTryMethod(adapt,"TSAdaptDSPSetFilter_C",(TSAdapt,const char*),(adapt,name));CHKERRQ(ierr);
343cb7ab186SLisandro Dalcin   PetscFunctionReturn(0);
344cb7ab186SLisandro Dalcin }
345cb7ab186SLisandro Dalcin 
346cb7ab186SLisandro Dalcin /*@
347cb7ab186SLisandro Dalcin    TSAdaptDSPSetPID - Set the PID controller parameters
348cb7ab186SLisandro Dalcin 
349cb7ab186SLisandro Dalcin    Input Arguments:
350cb7ab186SLisandro Dalcin +  adapt - adaptive controller context
351cb7ab186SLisandro Dalcin .  kkI - Integral parameter
352cb7ab186SLisandro Dalcin .  kkP - Proportional parameter
353cb7ab186SLisandro Dalcin -  kkD - Derivative parameter
354cb7ab186SLisandro Dalcin 
355cb7ab186SLisandro Dalcin    Level: intermediate
356cb7ab186SLisandro Dalcin 
357cb7ab186SLisandro Dalcin    References: http://dx.doi.org/10.1016/j.cam.2005.03.008
358cb7ab186SLisandro Dalcin 
359cb7ab186SLisandro Dalcin    Options Database:
360cb7ab186SLisandro Dalcin .   -ts_adapt_dsp_pid <kkI,kkP,kkD>
361cb7ab186SLisandro Dalcin 
362cb7ab186SLisandro Dalcin .seealso: TS, TSAdapt, TSGetAdapt(), TSAdaptDSPSetPID()
363cb7ab186SLisandro Dalcin @*/
364cb7ab186SLisandro Dalcin PetscErrorCode TSAdaptDSPSetPID(TSAdapt adapt,PetscReal kkI,PetscReal kkP,PetscReal kkD)
365cb7ab186SLisandro Dalcin {
366cb7ab186SLisandro Dalcin   PetscErrorCode ierr;
367cb7ab186SLisandro Dalcin   PetscFunctionBegin;
368cb7ab186SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
369cb7ab186SLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt,kkI,2);
370cb7ab186SLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt,kkP,3);
371cb7ab186SLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt,kkD,4);
372cb7ab186SLisandro Dalcin   ierr = PetscTryMethod(adapt,"TSAdaptDSPSetPID_C",(TSAdapt,PetscReal,PetscReal,PetscReal),(adapt,kkI,kkP,kkD));CHKERRQ(ierr);
373cb7ab186SLisandro Dalcin   PetscFunctionReturn(0);
374cb7ab186SLisandro Dalcin }
375cb7ab186SLisandro Dalcin 
376cb7ab186SLisandro Dalcin /*MC
377cb7ab186SLisandro Dalcin    TSADAPTDSP - Adaptive controller for time-stepping based on digital signal processing (DSP)
378cb7ab186SLisandro Dalcin 
379cb7ab186SLisandro Dalcin    Level: intermediate
380cb7ab186SLisandro Dalcin 
381cb7ab186SLisandro Dalcin    References: http://dx.doi.org/10.1145/641876.641877 and http://dx.doi.org/10.1016/j.cam.2005.03.008
382cb7ab186SLisandro Dalcin 
383cb7ab186SLisandro Dalcin    Options Database:
384cb7ab186SLisandro Dalcin +   -ts_adapt_dsp_filter <name>
385cb7ab186SLisandro Dalcin .   -ts_adapt_dsp_pid <kkI,kkP,kkD>
386cb7ab186SLisandro Dalcin .   -ts_adapt_dsp_kbeta <b1,b2,b2>
387cb7ab186SLisandro Dalcin -   -ts_adapt_dsp_alpha <a2,a3>
388cb7ab186SLisandro Dalcin 
389cb7ab186SLisandro Dalcin .seealso: TS, TSAdapt, TSGetAdapt()
390cb7ab186SLisandro Dalcin M*/
391cb7ab186SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptCreate_DSP(TSAdapt adapt)
392cb7ab186SLisandro Dalcin {
393cb7ab186SLisandro Dalcin   TSAdapt_DSP    *dsp;
394cb7ab186SLisandro Dalcin   PetscErrorCode ierr;
395cb7ab186SLisandro Dalcin 
396cb7ab186SLisandro Dalcin   PetscFunctionBegin;
397cb7ab186SLisandro Dalcin   ierr = PetscNewLog(adapt,&dsp);CHKERRQ(ierr);
398cb7ab186SLisandro Dalcin   adapt->reject_safety = 1.0; /* unused */
399cb7ab186SLisandro Dalcin 
400cb7ab186SLisandro Dalcin   adapt->data                = (void*)dsp;
401cb7ab186SLisandro Dalcin   adapt->ops->choose         = TSAdaptChoose_DSP;
402cb7ab186SLisandro Dalcin   adapt->ops->setfromoptions = TSAdaptSetFromOptions_DSP;
403cb7ab186SLisandro Dalcin   adapt->ops->destroy        = TSAdaptDestroy_DSP;
404cb7ab186SLisandro Dalcin   adapt->ops->reset          = TSAdaptReset_DSP;
405cb7ab186SLisandro Dalcin   adapt->ops->view           = TSAdaptView_DSP;
406cb7ab186SLisandro Dalcin 
407cb7ab186SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)adapt,"TSAdaptDSPSetFilter_C",TSAdaptDSPSetFilter_DSP);CHKERRQ(ierr);
408cb7ab186SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)adapt,"TSAdaptDSPSetPID_C",TSAdaptDSPSetPID_DSP);CHKERRQ(ierr);
409cb7ab186SLisandro Dalcin 
410cb7ab186SLisandro Dalcin   ierr = TSAdaptDSPSetFilter_DSP(adapt,"PI42");CHKERRQ(ierr);
411cb7ab186SLisandro Dalcin   ierr = TSAdaptRestart_DSP(adapt);CHKERRQ(ierr);
412cb7ab186SLisandro Dalcin   PetscFunctionReturn(0);
413cb7ab186SLisandro Dalcin }
414