1cb7ab186SLisandro Dalcin #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 286e171c2SStefano Zampini #include <petscdm.h> 3cb7ab186SLisandro Dalcin 4cb7ab186SLisandro Dalcin static const char *citation[] = { 5cb7ab186SLisandro Dalcin "@article{Soderlind2003,\n" 6cb7ab186SLisandro Dalcin " author = {S\"{o}derlind, Gustaf},\n" 7cb7ab186SLisandro Dalcin " title = {Digital Filters in Adaptive Time-stepping},\n" 8cb7ab186SLisandro Dalcin " journal = {ACM Transactions on Mathematical Software},\n" 9cb7ab186SLisandro Dalcin " volume = {29},\n" 10cb7ab186SLisandro Dalcin " number = {1},\n" 11cb7ab186SLisandro Dalcin " pages = {1--26},\n" 12cb7ab186SLisandro Dalcin " year = {2003},\n" 13cb7ab186SLisandro Dalcin " issn = {0098-3500},\n" 14cb7ab186SLisandro Dalcin " doi = {http://dx.doi.org/10.1145/641876.641877},\n" 15cb7ab186SLisandro Dalcin "}\n", 16cb7ab186SLisandro Dalcin "@article{Soderlind2006,\n" 17cb7ab186SLisandro Dalcin " author = {Gustaf S\"{o}derlind and Lina Wang},\n" 18cb7ab186SLisandro Dalcin " title = {Adaptive time-stepping and computational stability},\n" 19cb7ab186SLisandro Dalcin " journal = {Journal of Computational and Applied Mathematics},\n" 20cb7ab186SLisandro Dalcin " volume = {185},\n" 21cb7ab186SLisandro Dalcin " number = {2},\n" 22cb7ab186SLisandro Dalcin " pages = {225--243},\n" 23cb7ab186SLisandro Dalcin " year = {2006},\n" 24cb7ab186SLisandro Dalcin " issn = {0377-0427},\n" 25cb7ab186SLisandro Dalcin " doi = {http://dx.doi.org/10.1016/j.cam.2005.03.008},\n" 26cb7ab186SLisandro Dalcin "}\n", 27cb7ab186SLisandro Dalcin }; 28cb7ab186SLisandro Dalcin static PetscBool cited[] = {PETSC_FALSE, PETSC_FALSE}; 29cb7ab186SLisandro Dalcin 30cb7ab186SLisandro Dalcin typedef struct { 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 38d71ae5a4SJacob Faibussowitsch static PetscReal Limiter(PetscReal value, PetscReal kappa) 39d71ae5a4SJacob Faibussowitsch { 40cb7ab186SLisandro Dalcin return 1 + kappa * PetscAtanReal((value - 1) / kappa); 41cb7ab186SLisandro Dalcin } 42cb7ab186SLisandro Dalcin 43d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptRestart_DSP(TSAdapt adapt) 44d71ae5a4SJacob Faibussowitsch { 45cb7ab186SLisandro Dalcin TSAdapt_DSP *dsp = (TSAdapt_DSP *)adapt->data; 46*4d86920dSPierre Jolivet 47cb7ab186SLisandro Dalcin PetscFunctionBegin; 48cb7ab186SLisandro Dalcin dsp->cerror[0] = dsp->hratio[0] = 1.0; 49cb7ab186SLisandro Dalcin dsp->cerror[1] = dsp->hratio[1] = 1.0; 50cb7ab186SLisandro Dalcin dsp->cerror[2] = dsp->hratio[2] = 1.0; 513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 52cb7ab186SLisandro Dalcin } 53cb7ab186SLisandro Dalcin 54d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptRollBack_DSP(TSAdapt adapt) 55d71ae5a4SJacob Faibussowitsch { 56cb7ab186SLisandro Dalcin TSAdapt_DSP *dsp = (TSAdapt_DSP *)adapt->data; 57*4d86920dSPierre Jolivet 58cb7ab186SLisandro Dalcin PetscFunctionBegin; 59cb7ab186SLisandro Dalcin dsp->cerror[0] = dsp->cerror[1]; 60cb7ab186SLisandro Dalcin dsp->cerror[1] = dsp->cerror[2]; 61cb7ab186SLisandro Dalcin dsp->cerror[2] = 1.0; 62cb7ab186SLisandro Dalcin dsp->hratio[0] = dsp->hratio[1]; 63cb7ab186SLisandro Dalcin dsp->hratio[1] = dsp->hratio[2]; 64cb7ab186SLisandro Dalcin dsp->hratio[2] = 1.0; 653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 66cb7ab186SLisandro Dalcin } 67cb7ab186SLisandro Dalcin 68d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptChoose_DSP(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept, PetscReal *wlte, PetscReal *wltea, PetscReal *wlter) 69d71ae5a4SJacob Faibussowitsch { 70cb7ab186SLisandro Dalcin TSAdapt_DSP *dsp = (TSAdapt_DSP *)adapt->data; 71cb7ab186SLisandro Dalcin PetscInt order = PETSC_DECIDE; 72cb7ab186SLisandro Dalcin PetscReal enorm = -1; 73cb7ab186SLisandro Dalcin PetscReal enorma, enormr; 74cb7ab186SLisandro Dalcin PetscReal safety = adapt->safety * (PetscReal)0.9; 75cb7ab186SLisandro Dalcin PetscReal hnew, hfac = PETSC_INFINITY; 76cb7ab186SLisandro Dalcin PetscReal hmin = adapt->dt_min * (1 + PETSC_SQRT_MACHINE_EPSILON); 77cb7ab186SLisandro Dalcin 78cb7ab186SLisandro Dalcin PetscFunctionBegin; 79cb7ab186SLisandro Dalcin *next_sc = 0; /* Reuse the same order scheme */ 80cb7ab186SLisandro Dalcin *wltea = -1; /* Weighted absolute local truncation error is not used */ 81cb7ab186SLisandro Dalcin *wlter = -1; /* Weighted relative local truncation error is not used */ 82cb7ab186SLisandro Dalcin 83cb7ab186SLisandro Dalcin if (ts->ops->evaluatewlte) { 849566063dSJacob Faibussowitsch PetscCall(TSEvaluateWLTE(ts, adapt->wnormtype, &order, &enorm)); 85cad9d221SBarry Smith PetscCheck(enorm < 0 || order >= 1, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Computed error order %" PetscInt_FMT " must be positive", order); 86cb7ab186SLisandro Dalcin } else if (ts->ops->evaluatestep) { 8786e171c2SStefano Zampini DM dm; 8886e171c2SStefano Zampini Vec Y; 8986e171c2SStefano Zampini 9008401ef6SPierre Jolivet PetscCheck(adapt->candidates.n >= 1, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONGSTATE, "No candidate has been registered"); 9163a3b9bcSJacob Faibussowitsch PetscCheck(adapt->candidates.inuse_set, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONGSTATE, "The current in-use scheme is not among the %" PetscInt_FMT " candidates", adapt->candidates.n); 92cb7ab186SLisandro Dalcin order = adapt->candidates.order[0]; 939566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 949566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &Y)); 959566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep(ts, order - 1, Y, NULL)); 969566063dSJacob Faibussowitsch PetscCall(TSErrorWeightedNorm(ts, ts->vec_sol, Y, adapt->wnormtype, &enorm, &enorma, &enormr)); 979566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &Y)); 98cb7ab186SLisandro Dalcin } 99cb7ab186SLisandro Dalcin if (enorm < 0) { 1009566063dSJacob Faibussowitsch PetscCall(TSAdaptRestart_DSP(adapt)); 101cb7ab186SLisandro Dalcin *accept = PETSC_TRUE; /* Accept the step */ 102cb7ab186SLisandro Dalcin *next_h = h; /* Reuse the old step size */ 103cb7ab186SLisandro Dalcin *wlte = -1; /* Weighted local truncation error was not evaluated */ 1043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 105cb7ab186SLisandro Dalcin } 106cb7ab186SLisandro Dalcin 1079566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(citation[0], &cited[0])); 1089566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(citation[1], &cited[1])); 109cb7ab186SLisandro Dalcin 110cb7ab186SLisandro Dalcin /* Update history after rollback */ 1119371c9d4SSatish Balay if (!ts->steprollback) dsp->rollback = PETSC_FALSE; 112cb7ab186SLisandro Dalcin else if (!dsp->rollback) { 113cb7ab186SLisandro Dalcin dsp->rollback = PETSC_TRUE; 1149566063dSJacob Faibussowitsch PetscCall(TSAdaptRollBack_DSP(adapt)); 115cb7ab186SLisandro Dalcin } 116cb7ab186SLisandro Dalcin /* Reset history after restart */ 1171baa6e33SBarry Smith if (ts->steprestart) PetscCall(TSAdaptRestart_DSP(adapt)); 118cb7ab186SLisandro Dalcin 119cb7ab186SLisandro Dalcin { 120cb7ab186SLisandro Dalcin PetscReal k = (PetscReal)order; 121cb7ab186SLisandro Dalcin PetscReal b1 = dsp->kBeta[0]; 122cb7ab186SLisandro Dalcin PetscReal b2 = dsp->kBeta[1]; 123cb7ab186SLisandro Dalcin PetscReal b3 = dsp->kBeta[2]; 124cb7ab186SLisandro Dalcin PetscReal a2 = dsp->Alpha[0]; 1250ad35642SHendrik Ranocha PetscReal a3 = dsp->Alpha[1]; 126cb7ab186SLisandro Dalcin 127cb7ab186SLisandro Dalcin PetscReal ctr0; 128cb7ab186SLisandro Dalcin PetscReal ctr1 = dsp->cerror[0]; 129cb7ab186SLisandro Dalcin PetscReal ctr2 = dsp->cerror[1]; 130cb7ab186SLisandro Dalcin PetscReal rho0; 131cb7ab186SLisandro Dalcin PetscReal rho1 = dsp->hratio[0]; 132cb7ab186SLisandro Dalcin PetscReal rho2 = dsp->hratio[1]; 133cb7ab186SLisandro Dalcin 134cb7ab186SLisandro Dalcin /* Compute the step size ratio */ 135cb7ab186SLisandro Dalcin enorm = PetscMax(enorm, PETSC_SMALL); 136cb7ab186SLisandro Dalcin ctr0 = PetscPowReal(1 / enorm, 1 / k); 137cb7ab186SLisandro Dalcin rho0 = PetscPowReal(ctr0, b1); 138cb7ab186SLisandro Dalcin rho0 *= PetscPowReal(ctr1, b2); 139cb7ab186SLisandro Dalcin rho0 *= PetscPowReal(ctr2, b3); 140cb7ab186SLisandro Dalcin rho0 *= PetscPowReal(rho1, -a2); 141cb7ab186SLisandro Dalcin rho0 *= PetscPowReal(rho2, -a3); 142cb7ab186SLisandro Dalcin rho0 = Limiter(rho0, 1); 143cb7ab186SLisandro Dalcin 144cb7ab186SLisandro Dalcin /* Determine whether the step is accepted or rejected */ 1459371c9d4SSatish Balay if (rho0 >= safety) *accept = PETSC_TRUE; 1469371c9d4SSatish Balay else if (adapt->always_accept) *accept = PETSC_TRUE; 1479371c9d4SSatish Balay else if (h < hmin) *accept = PETSC_TRUE; 1489371c9d4SSatish Balay else *accept = PETSC_FALSE; 149cb7ab186SLisandro Dalcin 150cb7ab186SLisandro Dalcin /* Update history after accept */ 151cb7ab186SLisandro Dalcin if (*accept) { 152cb7ab186SLisandro Dalcin dsp->cerror[2] = dsp->cerror[1]; 153cb7ab186SLisandro Dalcin dsp->cerror[1] = dsp->cerror[0]; 154cb7ab186SLisandro Dalcin dsp->cerror[0] = ctr0; 155cb7ab186SLisandro Dalcin dsp->hratio[2] = dsp->hratio[1]; 156cb7ab186SLisandro Dalcin dsp->hratio[1] = dsp->hratio[0]; 157cb7ab186SLisandro Dalcin dsp->hratio[0] = rho0; 158cb7ab186SLisandro Dalcin dsp->rollback = PETSC_FALSE; 159cb7ab186SLisandro Dalcin } 160cb7ab186SLisandro Dalcin 161cb7ab186SLisandro Dalcin hfac = rho0; 162cb7ab186SLisandro Dalcin } 163cb7ab186SLisandro Dalcin 164cb7ab186SLisandro Dalcin hnew = h * PetscClipInterval(hfac, adapt->clip[0], adapt->clip[1]); 165cb7ab186SLisandro Dalcin *next_h = PetscClipInterval(hnew, adapt->dt_min, adapt->dt_max); 166cb7ab186SLisandro Dalcin *wlte = enorm; 1673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 168cb7ab186SLisandro Dalcin } 169cb7ab186SLisandro Dalcin 170d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptDestroy_DSP(TSAdapt adapt) 171d71ae5a4SJacob Faibussowitsch { 172cb7ab186SLisandro Dalcin PetscFunctionBegin; 1739566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)adapt, "TSAdaptDSPSetFilter_C", NULL)); 1749566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)adapt, "TSAdaptDSPSetPID_C", NULL)); 1759566063dSJacob Faibussowitsch PetscCall(PetscFree(adapt->data)); 1763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 177cb7ab186SLisandro Dalcin } 178cb7ab186SLisandro Dalcin 179d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptView_DSP(TSAdapt adapt, PetscViewer viewer) 180d71ae5a4SJacob Faibussowitsch { 181cb7ab186SLisandro Dalcin TSAdapt_DSP *dsp = (TSAdapt_DSP *)adapt->data; 182cb7ab186SLisandro Dalcin PetscBool iascii; 183cb7ab186SLisandro Dalcin 184cb7ab186SLisandro Dalcin PetscFunctionBegin; 1859566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 186cb7ab186SLisandro Dalcin if (iascii) { 187cb7ab186SLisandro Dalcin double a2 = (double)dsp->Alpha[0], a3 = (double)dsp->Alpha[1]; 188cb7ab186SLisandro Dalcin double b1 = (double)dsp->kBeta[0], b2 = (double)dsp->kBeta[1], b3 = (double)dsp->kBeta[2]; 1899566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "filter parameters kBeta=[%g,%g,%g] Alpha=[%g,%g]\n", b1, b2, b3, a2, a3)); 190cb7ab186SLisandro Dalcin } 1913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 192cb7ab186SLisandro Dalcin } 193cb7ab186SLisandro Dalcin 194cb7ab186SLisandro Dalcin struct FilterTab { 195cb7ab186SLisandro Dalcin const char *name; 196cb7ab186SLisandro Dalcin PetscReal scale; 197cb7ab186SLisandro Dalcin PetscReal kBeta[3]; 198cb7ab186SLisandro Dalcin PetscReal Alpha[2]; 199cb7ab186SLisandro Dalcin }; 200cb7ab186SLisandro Dalcin 201cb7ab186SLisandro Dalcin static struct FilterTab filterlist[] = { 202cb7ab186SLisandro Dalcin {"basic", 1, {1, 0, 0}, {0, 0} }, 203cb7ab186SLisandro Dalcin 204cb7ab186SLisandro Dalcin {"PI30", 3, {1, 0, 0}, {0, 0} }, 205cb7ab186SLisandro Dalcin {"PI42", 5, {3, -1, 0}, {0, 0} }, 206cb7ab186SLisandro Dalcin {"PI33", 3, {2, -1, 0}, {0, 0} }, 207cb7ab186SLisandro Dalcin {"PI34", 10, {7, -4, 0}, {0, 0} }, 208cb7ab186SLisandro Dalcin 209cb7ab186SLisandro Dalcin {"PC11", 1, {2, -1, 0}, {-1, 0} }, 210cb7ab186SLisandro Dalcin {"PC47", 10, {11, -7, 0}, {-10, 0} }, 211cb7ab186SLisandro Dalcin {"PC36", 10, {9, -6, 0}, {-10, 0} }, 212cb7ab186SLisandro Dalcin 213cb7ab186SLisandro Dalcin {"H0211", 2, {1, 1, 0}, {1, 0} }, 214cb7ab186SLisandro Dalcin {"H211b", 4, {1, 1, 0}, {1, 0} }, 215cb7ab186SLisandro Dalcin {"H211PI", 6, {1, 1, 0}, {0, 0} }, 216cb7ab186SLisandro Dalcin 217cb7ab186SLisandro Dalcin {"H0312", 4, {1, 2, 1}, {3, 1} }, 218cb7ab186SLisandro Dalcin {"H312b", 8, {1, 2, 1}, {3, 1} }, 219cb7ab186SLisandro Dalcin {"H312PID", 18, {1, 2, 1}, {0, 0} }, 220cb7ab186SLisandro Dalcin 221cb7ab186SLisandro Dalcin {"H0321", 4, {5, 2, -3}, {-1, -3} }, 222cb7ab186SLisandro Dalcin {"H321", 18, {6, 1, -5}, {-15, -3}}, 223cb7ab186SLisandro Dalcin }; 224cb7ab186SLisandro Dalcin 225d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptDSPSetFilter_DSP(TSAdapt adapt, const char *name) 226d71ae5a4SJacob Faibussowitsch { 227cb7ab186SLisandro Dalcin TSAdapt_DSP *dsp = (TSAdapt_DSP *)adapt->data; 228400f6f02SBarry Smith PetscInt i, count = PETSC_STATIC_ARRAY_LENGTH(filterlist); 229cb7ab186SLisandro Dalcin struct FilterTab *tab = NULL; 230cb7ab186SLisandro Dalcin PetscBool match; 231cb7ab186SLisandro Dalcin 232cb7ab186SLisandro Dalcin PetscFunctionBegin; 233cb7ab186SLisandro Dalcin for (i = 0; i < count; i++) { 2349566063dSJacob Faibussowitsch PetscCall(PetscStrcasecmp(name, filterlist[i].name, &match)); 2359371c9d4SSatish Balay if (match) { 2369371c9d4SSatish Balay tab = &filterlist[i]; 2379371c9d4SSatish Balay break; 2389371c9d4SSatish Balay } 239cb7ab186SLisandro Dalcin } 2403c633725SBarry Smith PetscCheck(tab, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_UNKNOWN_TYPE, "Filter name %s not found", name); 241cb7ab186SLisandro Dalcin dsp->kBeta[0] = tab->kBeta[0] / tab->scale; 242cb7ab186SLisandro Dalcin dsp->kBeta[1] = tab->kBeta[1] / tab->scale; 243cb7ab186SLisandro Dalcin dsp->kBeta[2] = tab->kBeta[2] / tab->scale; 244cb7ab186SLisandro Dalcin dsp->Alpha[0] = tab->Alpha[0] / tab->scale; 245cb7ab186SLisandro Dalcin dsp->Alpha[1] = tab->Alpha[1] / tab->scale; 2463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 247cb7ab186SLisandro Dalcin } 248cb7ab186SLisandro Dalcin 249d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptDSPSetPID_DSP(TSAdapt adapt, PetscReal kkI, PetscReal kkP, PetscReal kkD) 250d71ae5a4SJacob Faibussowitsch { 251cb7ab186SLisandro Dalcin TSAdapt_DSP *dsp = (TSAdapt_DSP *)adapt->data; 252cb7ab186SLisandro Dalcin 253cb7ab186SLisandro Dalcin PetscFunctionBegin; 254cb7ab186SLisandro Dalcin dsp->kBeta[0] = kkI + kkP + kkD; 255cb7ab186SLisandro Dalcin dsp->kBeta[1] = -(kkP + 2 * kkD); 256cb7ab186SLisandro Dalcin dsp->kBeta[2] = kkD; 257cb7ab186SLisandro Dalcin dsp->Alpha[0] = 0; 258cb7ab186SLisandro Dalcin dsp->Alpha[1] = 0; 2593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 260cb7ab186SLisandro Dalcin } 261cb7ab186SLisandro Dalcin 262d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptSetFromOptions_DSP(TSAdapt adapt, PetscOptionItems *PetscOptionsObject) 263d71ae5a4SJacob Faibussowitsch { 264cb7ab186SLisandro Dalcin TSAdapt_DSP *dsp = (TSAdapt_DSP *)adapt->data; 265400f6f02SBarry Smith const char *names[PETSC_STATIC_ARRAY_LENGTH(filterlist)]; 266400f6f02SBarry Smith PetscInt count = PETSC_STATIC_ARRAY_LENGTH(filterlist); 267cb7ab186SLisandro Dalcin PetscInt index = 2; /* PI42 */ 268cb7ab186SLisandro Dalcin PetscReal pid[3] = {1, 0, 0}; 269cb7ab186SLisandro Dalcin PetscInt i, n; 270cb7ab186SLisandro Dalcin PetscBool set; 271cb7ab186SLisandro Dalcin 272cb7ab186SLisandro Dalcin PetscFunctionBegin; 273cb7ab186SLisandro Dalcin for (i = 0; i < count; i++) names[i] = filterlist[i].name; 274d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "DSP adaptive controller options"); 275cb7ab186SLisandro Dalcin 2769566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_adapt_dsp_filter", "Filter name", "TSAdaptDSPSetFilter", names, count, names[index], &index, &set)); 2779566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptDSPSetFilter(adapt, names[index])); 278cb7ab186SLisandro Dalcin 2799566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-ts_adapt_dsp_pid", "PID parameters <kkI,kkP,kkD>", "TSAdaptDSPSetPID", pid, (n = 3, &n), &set)); 28008401ef6SPierre Jolivet PetscCheck(!set || n, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONG, "Must provide at least one value for PID parameters"); 2819566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptDSPSetPID(adapt, pid[0], pid[1], pid[2])); 282cb7ab186SLisandro Dalcin 2839566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-ts_adapt_dsp_kbeta", "Filter parameters", "", dsp->kBeta, (n = 3, &n), &set)); 28408401ef6SPierre Jolivet PetscCheck(!set || n, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONG, "Must provide at least one value for parameter kbeta"); 2859371c9d4SSatish Balay if (set) 2869371c9d4SSatish Balay for (i = n; i < 3; i++) dsp->kBeta[i] = 0; 287cb7ab186SLisandro Dalcin 2889566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-ts_adapt_dsp_alpha", "Filter parameters", "", dsp->Alpha, (n = 2, &n), &set)); 28908401ef6SPierre Jolivet PetscCheck(!set || n, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONG, "Must provide at least one value for parameter alpha"); 2909371c9d4SSatish Balay if (set) 2919371c9d4SSatish Balay for (i = n; i < 2; i++) dsp->Alpha[i] = 0; 292cb7ab186SLisandro Dalcin 293d0609cedSBarry Smith PetscOptionsHeadEnd(); 2943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 295cb7ab186SLisandro Dalcin } 296cb7ab186SLisandro Dalcin 297cb7ab186SLisandro Dalcin /*@C 2981d27aa22SBarry Smith TSAdaptDSPSetFilter - Sets internal parameters corresponding to the named filter {cite}`soderlind2006adaptive` {cite}`soderlind2003digital` 299cb7ab186SLisandro Dalcin 300c3339decSBarry Smith Collective 301cb7ab186SLisandro Dalcin 3024165533cSJose E. Roman Input Parameters: 303cb7ab186SLisandro Dalcin + adapt - adaptive controller context 304cb7ab186SLisandro Dalcin - name - filter name 305cb7ab186SLisandro Dalcin 306bcf0153eSBarry Smith Options Database Key: 30714d0ab18SJacob Faibussowitsch + -ts_adapt_dsp_filter <name> - Sets predefined controller by name; use -help for a list of available controllers 308bcf0153eSBarry Smith 30914d0ab18SJacob Faibussowitsch Filter names\: 31014d0ab18SJacob Faibussowitsch . basic - similar to `TSADAPTBASIC` but with different criteria for step rejections. 311bcf0153eSBarry Smith . PI30, PI42, PI33, PI34 - PI controllers. 312bcf0153eSBarry Smith . PC11, PC47, PC36 - predictive controllers. 313bcf0153eSBarry Smith . H0211, H211b, H211PI - digital filters with orders dynamics=2, adaptivity=1, filter=1. 314bcf0153eSBarry Smith . H0312, H312b, H312PID - digital filters with orders dynamics=3, adaptivity=1, filter=2. 315bcf0153eSBarry Smith - H0321, H321 - digital filters with orders dynamics=3, adaptivity=2, filter=1. 316bcf0153eSBarry Smith 317cb7ab186SLisandro Dalcin Level: intermediate 318cb7ab186SLisandro Dalcin 31914d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TSADAPTDSP`, `TS`, `TSAdapt`, `TSGetAdapt()`, `TSAdaptDSPSetPID()` 320cb7ab186SLisandro Dalcin @*/ 321d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptDSPSetFilter(TSAdapt adapt, const char *name) 322d71ae5a4SJacob Faibussowitsch { 323cb7ab186SLisandro Dalcin PetscFunctionBegin; 324cb7ab186SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 3254f572ea9SToby Isaac PetscAssertPointer(name, 2); 326cac4c232SBarry Smith PetscTryMethod(adapt, "TSAdaptDSPSetFilter_C", (TSAdapt, const char *), (adapt, name)); 3273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 328cb7ab186SLisandro Dalcin } 329cb7ab186SLisandro Dalcin 330cb7ab186SLisandro Dalcin /*@ 3311d27aa22SBarry Smith TSAdaptDSPSetPID - Set the PID controller parameters {cite}`soderlind2006adaptive` {cite}`soderlind2003digital` 332cb7ab186SLisandro Dalcin 3334165533cSJose E. Roman Input Parameters: 334cb7ab186SLisandro Dalcin + adapt - adaptive controller context 335cb7ab186SLisandro Dalcin . kkI - Integral parameter 336cb7ab186SLisandro Dalcin . kkP - Proportional parameter 337cb7ab186SLisandro Dalcin - kkD - Derivative parameter 338cb7ab186SLisandro Dalcin 339bcf0153eSBarry Smith Options Database Key: 340bcf0153eSBarry Smith . -ts_adapt_dsp_pid <kkI,kkP,kkD> - Sets PID controller parameters 341bcf0153eSBarry Smith 342cb7ab186SLisandro Dalcin Level: intermediate 343cb7ab186SLisandro Dalcin 3441cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdapt`, `TSGetAdapt()`, `TSAdaptDSPSetFilter()` 345cb7ab186SLisandro Dalcin @*/ 346d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptDSPSetPID(TSAdapt adapt, PetscReal kkI, PetscReal kkP, PetscReal kkD) 347d71ae5a4SJacob Faibussowitsch { 348cb7ab186SLisandro Dalcin PetscFunctionBegin; 349cb7ab186SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 350cb7ab186SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, kkI, 2); 351cb7ab186SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, kkP, 3); 352cb7ab186SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, kkD, 4); 353cac4c232SBarry Smith PetscTryMethod(adapt, "TSAdaptDSPSetPID_C", (TSAdapt, PetscReal, PetscReal, PetscReal), (adapt, kkI, kkP, kkD)); 3543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 355cb7ab186SLisandro Dalcin } 356cb7ab186SLisandro Dalcin 357cb7ab186SLisandro Dalcin /*MC 358cb7ab186SLisandro Dalcin TSADAPTDSP - Adaptive controller for time-stepping based on digital signal processing (DSP) 3591d27aa22SBarry Smith {cite}`soderlind2006adaptive` {cite}`soderlind2003digital` 360cb7ab186SLisandro Dalcin 3611d27aa22SBarry Smith Options Database Keys: 362bcf0153eSBarry Smith + -ts_adapt_dsp_filter <name> - Sets predefined controller by name; use -help for a list of available controllers 363bcf0153eSBarry Smith . -ts_adapt_dsp_pid <kkI,kkP,kkD> - Sets PID controller parameters 364bcf0153eSBarry Smith . -ts_adapt_dsp_kbeta <b1,b2,b2> - Sets general filter parameters 365bcf0153eSBarry Smith - -ts_adapt_dsp_alpha <a2,a3> - Sets general filter parameters 366bcf0153eSBarry Smith 367cb7ab186SLisandro Dalcin Level: intermediate 368cb7ab186SLisandro Dalcin 3691cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdapt`, `TSGetAdapt()`, `TSAdaptDSPSetPID()`, `TSAdaptDSPSetFilter()`, `TSAdaptType` 370cb7ab186SLisandro Dalcin M*/ 371d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSAdaptCreate_DSP(TSAdapt adapt) 372d71ae5a4SJacob Faibussowitsch { 373cb7ab186SLisandro Dalcin TSAdapt_DSP *dsp; 374cb7ab186SLisandro Dalcin 375cb7ab186SLisandro Dalcin PetscFunctionBegin; 3764dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&dsp)); 377cb7ab186SLisandro Dalcin adapt->reject_safety = 1.0; /* unused */ 378cb7ab186SLisandro Dalcin 379cb7ab186SLisandro Dalcin adapt->data = (void *)dsp; 380cb7ab186SLisandro Dalcin adapt->ops->choose = TSAdaptChoose_DSP; 381cb7ab186SLisandro Dalcin adapt->ops->setfromoptions = TSAdaptSetFromOptions_DSP; 382cb7ab186SLisandro Dalcin adapt->ops->destroy = TSAdaptDestroy_DSP; 383cb7ab186SLisandro Dalcin adapt->ops->view = TSAdaptView_DSP; 384cb7ab186SLisandro Dalcin 3859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)adapt, "TSAdaptDSPSetFilter_C", TSAdaptDSPSetFilter_DSP)); 3869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)adapt, "TSAdaptDSPSetPID_C", TSAdaptDSPSetPID_DSP)); 387cb7ab186SLisandro Dalcin 3889566063dSJacob Faibussowitsch PetscCall(TSAdaptDSPSetFilter_DSP(adapt, "PI42")); 3899566063dSJacob Faibussowitsch PetscCall(TSAdaptRestart_DSP(adapt)); 3903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 391cb7ab186SLisandro Dalcin } 392