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