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; 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; 503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51cb7ab186SLisandro Dalcin } 52cb7ab186SLisandro Dalcin 53d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptRollBack_DSP(TSAdapt adapt) 54d71ae5a4SJacob Faibussowitsch { 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; 633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 64cb7ab186SLisandro Dalcin } 65cb7ab186SLisandro Dalcin 66d71ae5a4SJacob 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) 67d71ae5a4SJacob Faibussowitsch { 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 76cb7ab186SLisandro Dalcin PetscFunctionBegin; 77cb7ab186SLisandro Dalcin *next_sc = 0; /* Reuse the same order scheme */ 78cb7ab186SLisandro Dalcin *wltea = -1; /* Weighted absolute local truncation error is not used */ 79cb7ab186SLisandro Dalcin *wlter = -1; /* Weighted relative local truncation error is not used */ 80cb7ab186SLisandro Dalcin 81cb7ab186SLisandro Dalcin if (ts->ops->evaluatewlte) { 829566063dSJacob Faibussowitsch PetscCall(TSEvaluateWLTE(ts, adapt->wnormtype, &order, &enorm)); 83cad9d221SBarry Smith PetscCheck(enorm < 0 || order >= 1, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Computed error order %" PetscInt_FMT " must be positive", order); 84cb7ab186SLisandro Dalcin } else if (ts->ops->evaluatestep) { 8586e171c2SStefano Zampini DM dm; 8686e171c2SStefano Zampini Vec Y; 8786e171c2SStefano Zampini 8808401ef6SPierre Jolivet PetscCheck(adapt->candidates.n >= 1, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONGSTATE, "No candidate has been registered"); 8963a3b9bcSJacob 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); 90cb7ab186SLisandro Dalcin order = adapt->candidates.order[0]; 919566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 929566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &Y)); 939566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep(ts, order - 1, Y, NULL)); 949566063dSJacob Faibussowitsch PetscCall(TSErrorWeightedNorm(ts, ts->vec_sol, Y, adapt->wnormtype, &enorm, &enorma, &enormr)); 959566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &Y)); 96cb7ab186SLisandro Dalcin } 97cb7ab186SLisandro Dalcin if (enorm < 0) { 989566063dSJacob Faibussowitsch PetscCall(TSAdaptRestart_DSP(adapt)); 99cb7ab186SLisandro Dalcin *accept = PETSC_TRUE; /* Accept the step */ 100cb7ab186SLisandro Dalcin *next_h = h; /* Reuse the old step size */ 101cb7ab186SLisandro Dalcin *wlte = -1; /* Weighted local truncation error was not evaluated */ 1023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 103cb7ab186SLisandro Dalcin } 104cb7ab186SLisandro Dalcin 1059566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(citation[0], &cited[0])); 1069566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(citation[1], &cited[1])); 107cb7ab186SLisandro Dalcin 108cb7ab186SLisandro Dalcin /* Update history after rollback */ 1099371c9d4SSatish Balay if (!ts->steprollback) dsp->rollback = PETSC_FALSE; 110cb7ab186SLisandro Dalcin else if (!dsp->rollback) { 111cb7ab186SLisandro Dalcin dsp->rollback = PETSC_TRUE; 1129566063dSJacob Faibussowitsch PetscCall(TSAdaptRollBack_DSP(adapt)); 113cb7ab186SLisandro Dalcin } 114cb7ab186SLisandro Dalcin /* Reset history after restart */ 1151baa6e33SBarry Smith if (ts->steprestart) PetscCall(TSAdaptRestart_DSP(adapt)); 116cb7ab186SLisandro Dalcin 117cb7ab186SLisandro Dalcin { 118cb7ab186SLisandro Dalcin PetscReal k = (PetscReal)order; 119cb7ab186SLisandro Dalcin PetscReal b1 = dsp->kBeta[0]; 120cb7ab186SLisandro Dalcin PetscReal b2 = dsp->kBeta[1]; 121cb7ab186SLisandro Dalcin PetscReal b3 = dsp->kBeta[2]; 122cb7ab186SLisandro Dalcin PetscReal a2 = dsp->Alpha[0]; 1230ad35642SHendrik Ranocha PetscReal a3 = dsp->Alpha[1]; 124cb7ab186SLisandro Dalcin 125cb7ab186SLisandro Dalcin PetscReal ctr0; 126cb7ab186SLisandro Dalcin PetscReal ctr1 = dsp->cerror[0]; 127cb7ab186SLisandro Dalcin PetscReal ctr2 = dsp->cerror[1]; 128cb7ab186SLisandro Dalcin PetscReal rho0; 129cb7ab186SLisandro Dalcin PetscReal rho1 = dsp->hratio[0]; 130cb7ab186SLisandro Dalcin PetscReal rho2 = dsp->hratio[1]; 131cb7ab186SLisandro Dalcin 132cb7ab186SLisandro Dalcin /* Compute the step size ratio */ 133cb7ab186SLisandro Dalcin enorm = PetscMax(enorm, PETSC_SMALL); 134cb7ab186SLisandro Dalcin ctr0 = PetscPowReal(1 / enorm, 1 / k); 135cb7ab186SLisandro Dalcin rho0 = PetscPowReal(ctr0, b1); 136cb7ab186SLisandro Dalcin rho0 *= PetscPowReal(ctr1, b2); 137cb7ab186SLisandro Dalcin rho0 *= PetscPowReal(ctr2, b3); 138cb7ab186SLisandro Dalcin rho0 *= PetscPowReal(rho1, -a2); 139cb7ab186SLisandro Dalcin rho0 *= PetscPowReal(rho2, -a3); 140cb7ab186SLisandro Dalcin rho0 = Limiter(rho0, 1); 141cb7ab186SLisandro Dalcin 142cb7ab186SLisandro Dalcin /* Determine whether the step is accepted or rejected */ 1439371c9d4SSatish Balay if (rho0 >= safety) *accept = PETSC_TRUE; 1449371c9d4SSatish Balay else if (adapt->always_accept) *accept = PETSC_TRUE; 1459371c9d4SSatish Balay else if (h < hmin) *accept = PETSC_TRUE; 1469371c9d4SSatish Balay else *accept = PETSC_FALSE; 147cb7ab186SLisandro Dalcin 148cb7ab186SLisandro Dalcin /* Update history after accept */ 149cb7ab186SLisandro Dalcin if (*accept) { 150cb7ab186SLisandro Dalcin dsp->cerror[2] = dsp->cerror[1]; 151cb7ab186SLisandro Dalcin dsp->cerror[1] = dsp->cerror[0]; 152cb7ab186SLisandro Dalcin dsp->cerror[0] = ctr0; 153cb7ab186SLisandro Dalcin dsp->hratio[2] = dsp->hratio[1]; 154cb7ab186SLisandro Dalcin dsp->hratio[1] = dsp->hratio[0]; 155cb7ab186SLisandro Dalcin dsp->hratio[0] = rho0; 156cb7ab186SLisandro Dalcin dsp->rollback = PETSC_FALSE; 157cb7ab186SLisandro Dalcin } 158cb7ab186SLisandro Dalcin 159cb7ab186SLisandro Dalcin hfac = rho0; 160cb7ab186SLisandro Dalcin } 161cb7ab186SLisandro Dalcin 162cb7ab186SLisandro Dalcin hnew = h * PetscClipInterval(hfac, adapt->clip[0], adapt->clip[1]); 163cb7ab186SLisandro Dalcin *next_h = PetscClipInterval(hnew, adapt->dt_min, adapt->dt_max); 164cb7ab186SLisandro Dalcin *wlte = enorm; 1653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 166cb7ab186SLisandro Dalcin } 167cb7ab186SLisandro Dalcin 168d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptDestroy_DSP(TSAdapt adapt) 169d71ae5a4SJacob Faibussowitsch { 170cb7ab186SLisandro Dalcin PetscFunctionBegin; 1719566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)adapt, "TSAdaptDSPSetFilter_C", NULL)); 1729566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)adapt, "TSAdaptDSPSetPID_C", NULL)); 1739566063dSJacob Faibussowitsch PetscCall(PetscFree(adapt->data)); 1743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 175cb7ab186SLisandro Dalcin } 176cb7ab186SLisandro Dalcin 177d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptView_DSP(TSAdapt adapt, PetscViewer viewer) 178d71ae5a4SJacob Faibussowitsch { 179cb7ab186SLisandro Dalcin TSAdapt_DSP *dsp = (TSAdapt_DSP *)adapt->data; 180cb7ab186SLisandro Dalcin PetscBool iascii; 181cb7ab186SLisandro Dalcin 182cb7ab186SLisandro Dalcin PetscFunctionBegin; 1839566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 184cb7ab186SLisandro Dalcin if (iascii) { 185cb7ab186SLisandro Dalcin double a2 = (double)dsp->Alpha[0], a3 = (double)dsp->Alpha[1]; 186cb7ab186SLisandro Dalcin double b1 = (double)dsp->kBeta[0], b2 = (double)dsp->kBeta[1], b3 = (double)dsp->kBeta[2]; 1879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "filter parameters kBeta=[%g,%g,%g] Alpha=[%g,%g]\n", b1, b2, b3, a2, a3)); 188cb7ab186SLisandro Dalcin } 1893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 190cb7ab186SLisandro Dalcin } 191cb7ab186SLisandro Dalcin 192cb7ab186SLisandro Dalcin struct FilterTab { 193cb7ab186SLisandro Dalcin const char *name; 194cb7ab186SLisandro Dalcin PetscReal scale; 195cb7ab186SLisandro Dalcin PetscReal kBeta[3]; 196cb7ab186SLisandro Dalcin PetscReal Alpha[2]; 197cb7ab186SLisandro Dalcin }; 198cb7ab186SLisandro Dalcin 199cb7ab186SLisandro Dalcin static struct FilterTab filterlist[] = { 200cb7ab186SLisandro Dalcin {"basic", 1, {1, 0, 0}, {0, 0} }, 201cb7ab186SLisandro Dalcin 202cb7ab186SLisandro Dalcin {"PI30", 3, {1, 0, 0}, {0, 0} }, 203cb7ab186SLisandro Dalcin {"PI42", 5, {3, -1, 0}, {0, 0} }, 204cb7ab186SLisandro Dalcin {"PI33", 3, {2, -1, 0}, {0, 0} }, 205cb7ab186SLisandro Dalcin {"PI34", 10, {7, -4, 0}, {0, 0} }, 206cb7ab186SLisandro Dalcin 207cb7ab186SLisandro Dalcin {"PC11", 1, {2, -1, 0}, {-1, 0} }, 208cb7ab186SLisandro Dalcin {"PC47", 10, {11, -7, 0}, {-10, 0} }, 209cb7ab186SLisandro Dalcin {"PC36", 10, {9, -6, 0}, {-10, 0} }, 210cb7ab186SLisandro Dalcin 211cb7ab186SLisandro Dalcin {"H0211", 2, {1, 1, 0}, {1, 0} }, 212cb7ab186SLisandro Dalcin {"H211b", 4, {1, 1, 0}, {1, 0} }, 213cb7ab186SLisandro Dalcin {"H211PI", 6, {1, 1, 0}, {0, 0} }, 214cb7ab186SLisandro Dalcin 215cb7ab186SLisandro Dalcin {"H0312", 4, {1, 2, 1}, {3, 1} }, 216cb7ab186SLisandro Dalcin {"H312b", 8, {1, 2, 1}, {3, 1} }, 217cb7ab186SLisandro Dalcin {"H312PID", 18, {1, 2, 1}, {0, 0} }, 218cb7ab186SLisandro Dalcin 219cb7ab186SLisandro Dalcin {"H0321", 4, {5, 2, -3}, {-1, -3} }, 220cb7ab186SLisandro Dalcin {"H321", 18, {6, 1, -5}, {-15, -3}}, 221cb7ab186SLisandro Dalcin }; 222cb7ab186SLisandro Dalcin 223d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptDSPSetFilter_DSP(TSAdapt adapt, const char *name) 224d71ae5a4SJacob Faibussowitsch { 225cb7ab186SLisandro Dalcin TSAdapt_DSP *dsp = (TSAdapt_DSP *)adapt->data; 226400f6f02SBarry Smith PetscInt i, count = PETSC_STATIC_ARRAY_LENGTH(filterlist); 227cb7ab186SLisandro Dalcin struct FilterTab *tab = NULL; 228cb7ab186SLisandro Dalcin PetscBool match; 229cb7ab186SLisandro Dalcin 230cb7ab186SLisandro Dalcin PetscFunctionBegin; 231cb7ab186SLisandro Dalcin for (i = 0; i < count; i++) { 2329566063dSJacob Faibussowitsch PetscCall(PetscStrcasecmp(name, filterlist[i].name, &match)); 2339371c9d4SSatish Balay if (match) { 2349371c9d4SSatish Balay tab = &filterlist[i]; 2359371c9d4SSatish Balay break; 2369371c9d4SSatish Balay } 237cb7ab186SLisandro Dalcin } 2383c633725SBarry Smith PetscCheck(tab, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_UNKNOWN_TYPE, "Filter name %s not found", name); 239cb7ab186SLisandro Dalcin dsp->kBeta[0] = tab->kBeta[0] / tab->scale; 240cb7ab186SLisandro Dalcin dsp->kBeta[1] = tab->kBeta[1] / tab->scale; 241cb7ab186SLisandro Dalcin dsp->kBeta[2] = tab->kBeta[2] / tab->scale; 242cb7ab186SLisandro Dalcin dsp->Alpha[0] = tab->Alpha[0] / tab->scale; 243cb7ab186SLisandro Dalcin dsp->Alpha[1] = tab->Alpha[1] / tab->scale; 2443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 245cb7ab186SLisandro Dalcin } 246cb7ab186SLisandro Dalcin 247d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptDSPSetPID_DSP(TSAdapt adapt, PetscReal kkI, PetscReal kkP, PetscReal kkD) 248d71ae5a4SJacob Faibussowitsch { 249cb7ab186SLisandro Dalcin TSAdapt_DSP *dsp = (TSAdapt_DSP *)adapt->data; 250cb7ab186SLisandro Dalcin 251cb7ab186SLisandro Dalcin PetscFunctionBegin; 252cb7ab186SLisandro Dalcin dsp->kBeta[0] = kkI + kkP + kkD; 253cb7ab186SLisandro Dalcin dsp->kBeta[1] = -(kkP + 2 * kkD); 254cb7ab186SLisandro Dalcin dsp->kBeta[2] = kkD; 255cb7ab186SLisandro Dalcin dsp->Alpha[0] = 0; 256cb7ab186SLisandro Dalcin dsp->Alpha[1] = 0; 2573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 258cb7ab186SLisandro Dalcin } 259cb7ab186SLisandro Dalcin 260d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptSetFromOptions_DSP(TSAdapt adapt, PetscOptionItems *PetscOptionsObject) 261d71ae5a4SJacob Faibussowitsch { 262cb7ab186SLisandro Dalcin TSAdapt_DSP *dsp = (TSAdapt_DSP *)adapt->data; 263400f6f02SBarry Smith const char *names[PETSC_STATIC_ARRAY_LENGTH(filterlist)]; 264400f6f02SBarry Smith PetscInt count = PETSC_STATIC_ARRAY_LENGTH(filterlist); 265cb7ab186SLisandro Dalcin PetscInt index = 2; /* PI42 */ 266cb7ab186SLisandro Dalcin PetscReal pid[3] = {1, 0, 0}; 267cb7ab186SLisandro Dalcin PetscInt i, n; 268cb7ab186SLisandro Dalcin PetscBool set; 269cb7ab186SLisandro Dalcin 270cb7ab186SLisandro Dalcin PetscFunctionBegin; 271cb7ab186SLisandro Dalcin for (i = 0; i < count; i++) names[i] = filterlist[i].name; 272d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "DSP adaptive controller options"); 273cb7ab186SLisandro Dalcin 2749566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_adapt_dsp_filter", "Filter name", "TSAdaptDSPSetFilter", names, count, names[index], &index, &set)); 2759566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptDSPSetFilter(adapt, names[index])); 276cb7ab186SLisandro Dalcin 2779566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-ts_adapt_dsp_pid", "PID parameters <kkI,kkP,kkD>", "TSAdaptDSPSetPID", pid, (n = 3, &n), &set)); 27808401ef6SPierre Jolivet PetscCheck(!set || n, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONG, "Must provide at least one value for PID parameters"); 2799566063dSJacob Faibussowitsch if (set) PetscCall(TSAdaptDSPSetPID(adapt, pid[0], pid[1], pid[2])); 280cb7ab186SLisandro Dalcin 2819566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-ts_adapt_dsp_kbeta", "Filter parameters", "", dsp->kBeta, (n = 3, &n), &set)); 28208401ef6SPierre Jolivet PetscCheck(!set || n, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONG, "Must provide at least one value for parameter kbeta"); 2839371c9d4SSatish Balay if (set) 2849371c9d4SSatish Balay for (i = n; i < 3; i++) dsp->kBeta[i] = 0; 285cb7ab186SLisandro Dalcin 2869566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-ts_adapt_dsp_alpha", "Filter parameters", "", dsp->Alpha, (n = 2, &n), &set)); 28708401ef6SPierre Jolivet PetscCheck(!set || n, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONG, "Must provide at least one value for parameter alpha"); 2889371c9d4SSatish Balay if (set) 2899371c9d4SSatish Balay for (i = n; i < 2; i++) dsp->Alpha[i] = 0; 290cb7ab186SLisandro Dalcin 291d0609cedSBarry Smith PetscOptionsHeadEnd(); 2923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 293cb7ab186SLisandro Dalcin } 294cb7ab186SLisandro Dalcin 295cb7ab186SLisandro Dalcin /*@C 296*1d27aa22SBarry Smith TSAdaptDSPSetFilter - Sets internal parameters corresponding to the named filter {cite}`soderlind2006adaptive` {cite}`soderlind2003digital` 297cb7ab186SLisandro Dalcin 298c3339decSBarry Smith Collective 299cb7ab186SLisandro Dalcin 3004165533cSJose E. Roman Input Parameters: 301cb7ab186SLisandro Dalcin + adapt - adaptive controller context 302cb7ab186SLisandro Dalcin - name - filter name 303cb7ab186SLisandro Dalcin 304bcf0153eSBarry Smith Options Database Key: 30514d0ab18SJacob Faibussowitsch + -ts_adapt_dsp_filter <name> - Sets predefined controller by name; use -help for a list of available controllers 306bcf0153eSBarry Smith 30714d0ab18SJacob Faibussowitsch Filter names\: 30814d0ab18SJacob Faibussowitsch . basic - similar to `TSADAPTBASIC` but with different criteria for step rejections. 309bcf0153eSBarry Smith . PI30, PI42, PI33, PI34 - PI controllers. 310bcf0153eSBarry Smith . PC11, PC47, PC36 - predictive controllers. 311bcf0153eSBarry Smith . H0211, H211b, H211PI - digital filters with orders dynamics=2, adaptivity=1, filter=1. 312bcf0153eSBarry Smith . H0312, H312b, H312PID - digital filters with orders dynamics=3, adaptivity=1, filter=2. 313bcf0153eSBarry Smith - H0321, H321 - digital filters with orders dynamics=3, adaptivity=2, filter=1. 314bcf0153eSBarry Smith 315cb7ab186SLisandro Dalcin Level: intermediate 316cb7ab186SLisandro Dalcin 31714d0ab18SJacob Faibussowitsch .seealso: [](ch_ts), `TSADAPTDSP`, `TS`, `TSAdapt`, `TSGetAdapt()`, `TSAdaptDSPSetPID()` 318cb7ab186SLisandro Dalcin @*/ 319d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptDSPSetFilter(TSAdapt adapt, const char *name) 320d71ae5a4SJacob Faibussowitsch { 321cb7ab186SLisandro Dalcin PetscFunctionBegin; 322cb7ab186SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 3234f572ea9SToby Isaac PetscAssertPointer(name, 2); 324cac4c232SBarry Smith PetscTryMethod(adapt, "TSAdaptDSPSetFilter_C", (TSAdapt, const char *), (adapt, name)); 3253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 326cb7ab186SLisandro Dalcin } 327cb7ab186SLisandro Dalcin 328cb7ab186SLisandro Dalcin /*@ 329*1d27aa22SBarry Smith TSAdaptDSPSetPID - Set the PID controller parameters {cite}`soderlind2006adaptive` {cite}`soderlind2003digital` 330cb7ab186SLisandro Dalcin 3314165533cSJose E. Roman Input Parameters: 332cb7ab186SLisandro Dalcin + adapt - adaptive controller context 333cb7ab186SLisandro Dalcin . kkI - Integral parameter 334cb7ab186SLisandro Dalcin . kkP - Proportional parameter 335cb7ab186SLisandro Dalcin - kkD - Derivative parameter 336cb7ab186SLisandro Dalcin 337bcf0153eSBarry Smith Options Database Key: 338bcf0153eSBarry Smith . -ts_adapt_dsp_pid <kkI,kkP,kkD> - Sets PID controller parameters 339bcf0153eSBarry Smith 340cb7ab186SLisandro Dalcin Level: intermediate 341cb7ab186SLisandro Dalcin 3421cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdapt`, `TSGetAdapt()`, `TSAdaptDSPSetFilter()` 343cb7ab186SLisandro Dalcin @*/ 344d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdaptDSPSetPID(TSAdapt adapt, PetscReal kkI, PetscReal kkP, PetscReal kkD) 345d71ae5a4SJacob Faibussowitsch { 346cb7ab186SLisandro Dalcin PetscFunctionBegin; 347cb7ab186SLisandro Dalcin PetscValidHeaderSpecific(adapt, TSADAPT_CLASSID, 1); 348cb7ab186SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, kkI, 2); 349cb7ab186SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, kkP, 3); 350cb7ab186SLisandro Dalcin PetscValidLogicalCollectiveReal(adapt, kkD, 4); 351cac4c232SBarry Smith PetscTryMethod(adapt, "TSAdaptDSPSetPID_C", (TSAdapt, PetscReal, PetscReal, PetscReal), (adapt, kkI, kkP, kkD)); 3523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 353cb7ab186SLisandro Dalcin } 354cb7ab186SLisandro Dalcin 355cb7ab186SLisandro Dalcin /*MC 356cb7ab186SLisandro Dalcin TSADAPTDSP - Adaptive controller for time-stepping based on digital signal processing (DSP) 357*1d27aa22SBarry Smith {cite}`soderlind2006adaptive` {cite}`soderlind2003digital` 358cb7ab186SLisandro Dalcin 359*1d27aa22SBarry Smith Options Database Keys: 360bcf0153eSBarry Smith + -ts_adapt_dsp_filter <name> - Sets predefined controller by name; use -help for a list of available controllers 361bcf0153eSBarry Smith . -ts_adapt_dsp_pid <kkI,kkP,kkD> - Sets PID controller parameters 362bcf0153eSBarry Smith . -ts_adapt_dsp_kbeta <b1,b2,b2> - Sets general filter parameters 363bcf0153eSBarry Smith - -ts_adapt_dsp_alpha <a2,a3> - Sets general filter parameters 364bcf0153eSBarry Smith 365cb7ab186SLisandro Dalcin Level: intermediate 366cb7ab186SLisandro Dalcin 3671cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdapt`, `TSGetAdapt()`, `TSAdaptDSPSetPID()`, `TSAdaptDSPSetFilter()`, `TSAdaptType` 368cb7ab186SLisandro Dalcin M*/ 369d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSAdaptCreate_DSP(TSAdapt adapt) 370d71ae5a4SJacob Faibussowitsch { 371cb7ab186SLisandro Dalcin TSAdapt_DSP *dsp; 372cb7ab186SLisandro Dalcin 373cb7ab186SLisandro Dalcin PetscFunctionBegin; 3744dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&dsp)); 375cb7ab186SLisandro Dalcin adapt->reject_safety = 1.0; /* unused */ 376cb7ab186SLisandro Dalcin 377cb7ab186SLisandro Dalcin adapt->data = (void *)dsp; 378cb7ab186SLisandro Dalcin adapt->ops->choose = TSAdaptChoose_DSP; 379cb7ab186SLisandro Dalcin adapt->ops->setfromoptions = TSAdaptSetFromOptions_DSP; 380cb7ab186SLisandro Dalcin adapt->ops->destroy = TSAdaptDestroy_DSP; 381cb7ab186SLisandro Dalcin adapt->ops->view = TSAdaptView_DSP; 382cb7ab186SLisandro Dalcin 3839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)adapt, "TSAdaptDSPSetFilter_C", TSAdaptDSPSetFilter_DSP)); 3849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)adapt, "TSAdaptDSPSetPID_C", TSAdaptDSPSetPID_DSP)); 385cb7ab186SLisandro Dalcin 3869566063dSJacob Faibussowitsch PetscCall(TSAdaptDSPSetFilter_DSP(adapt, "PI42")); 3879566063dSJacob Faibussowitsch PetscCall(TSAdaptRestart_DSP(adapt)); 3883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 389cb7ab186SLisandro Dalcin } 390