1b07a2398SLisandro Dalcin /* 2b07a2398SLisandro Dalcin Code for timestepping with implicit generalized-\alpha method 3b07a2398SLisandro Dalcin for first order systems. 4b07a2398SLisandro Dalcin */ 5b07a2398SLisandro Dalcin #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 6b07a2398SLisandro Dalcin 7b07a2398SLisandro Dalcin static PetscBool cited = PETSC_FALSE; 89371c9d4SSatish Balay static const char citation[] = "@article{Jansen2000,\n" 9b07a2398SLisandro Dalcin " title = {A generalized-$\\alpha$ method for integrating the filtered {N}avier--{S}tokes equations with a stabilized finite element method},\n" 10b07a2398SLisandro Dalcin " author = {Kenneth E. Jansen and Christian H. Whiting and Gregory M. Hulbert},\n" 11b07a2398SLisandro Dalcin " journal = {Computer Methods in Applied Mechanics and Engineering},\n" 12b07a2398SLisandro Dalcin " volume = {190},\n" 13b07a2398SLisandro Dalcin " number = {3--4},\n" 14b07a2398SLisandro Dalcin " pages = {305--319},\n" 15b07a2398SLisandro Dalcin " year = {2000},\n" 16b07a2398SLisandro Dalcin " issn = {0045-7825},\n" 17b07a2398SLisandro Dalcin " doi = {http://dx.doi.org/10.1016/S0045-7825(00)00203-6}\n}\n"; 18b07a2398SLisandro Dalcin 19b07a2398SLisandro Dalcin typedef struct { 20b07a2398SLisandro Dalcin PetscReal stage_time; 21b07a2398SLisandro Dalcin PetscReal shift_V; 22b07a2398SLisandro Dalcin PetscReal scale_F; 23b07a2398SLisandro Dalcin Vec X0, Xa, X1; 24b07a2398SLisandro Dalcin Vec V0, Va, V1; 25b07a2398SLisandro Dalcin 26b07a2398SLisandro Dalcin PetscReal Alpha_m; 27b07a2398SLisandro Dalcin PetscReal Alpha_f; 28b07a2398SLisandro Dalcin PetscReal Gamma; 29b07a2398SLisandro Dalcin PetscInt order; 30b07a2398SLisandro Dalcin 31b07a2398SLisandro Dalcin Vec vec_sol_prev; 321566a47fSLisandro Dalcin Vec vec_lte_work; 33b07a2398SLisandro Dalcin 34b07a2398SLisandro Dalcin TSStepStatus status; 35b07a2398SLisandro Dalcin } TS_Alpha; 36b07a2398SLisandro Dalcin 378ec9177eSStefano Zampini static PetscErrorCode TSResizeRegister_Alpha(TS ts, PetscBool reg) 388ec9177eSStefano Zampini { 398ec9177eSStefano Zampini TS_Alpha *th = (TS_Alpha *)ts->data; 408ec9177eSStefano Zampini 418ec9177eSStefano Zampini PetscFunctionBegin; 42ecc87898SStefano Zampini if (reg) { 43ecc87898SStefano Zampini PetscCall(TSResizeRegisterVec(ts, "ts:theta:sol_prev", th->vec_sol_prev)); 44ecc87898SStefano Zampini PetscCall(TSResizeRegisterVec(ts, "ts:theta:X0", th->X0)); 45ecc87898SStefano Zampini } else { 46ecc87898SStefano Zampini PetscCall(TSResizeRetrieveVec(ts, "ts:theta:sol_prev", &th->vec_sol_prev)); 47ecc87898SStefano Zampini PetscCall(PetscObjectReference((PetscObject)th->vec_sol_prev)); 48ecc87898SStefano Zampini PetscCall(TSResizeRetrieveVec(ts, "ts:theta:X0", &th->X0)); 498ec9177eSStefano Zampini PetscCall(PetscObjectReference((PetscObject)th->X0)); 508ec9177eSStefano Zampini } 518ec9177eSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 528ec9177eSStefano Zampini } 538ec9177eSStefano Zampini 54d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_StageTime(TS ts) 55d71ae5a4SJacob Faibussowitsch { 56b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 57b07a2398SLisandro Dalcin PetscReal t = ts->ptime; 58b07a2398SLisandro Dalcin PetscReal dt = ts->time_step; 59b07a2398SLisandro Dalcin PetscReal Alpha_m = th->Alpha_m; 60b07a2398SLisandro Dalcin PetscReal Alpha_f = th->Alpha_f; 61b07a2398SLisandro Dalcin PetscReal Gamma = th->Gamma; 62b07a2398SLisandro Dalcin 63b07a2398SLisandro Dalcin PetscFunctionBegin; 64b07a2398SLisandro Dalcin th->stage_time = t + Alpha_f * dt; 65b07a2398SLisandro Dalcin th->shift_V = Alpha_m / (Alpha_f * Gamma * dt); 66b07a2398SLisandro Dalcin th->scale_F = 1 / Alpha_f; 673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 68b07a2398SLisandro Dalcin } 69b07a2398SLisandro Dalcin 70d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_StageVecs(TS ts, Vec X) 71d71ae5a4SJacob Faibussowitsch { 72b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 73b07a2398SLisandro Dalcin Vec X1 = X, V1 = th->V1; 74b07a2398SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va; 75b07a2398SLisandro Dalcin Vec X0 = th->X0, V0 = th->V0; 76b07a2398SLisandro Dalcin PetscReal dt = ts->time_step; 77b07a2398SLisandro Dalcin PetscReal Alpha_m = th->Alpha_m; 78b07a2398SLisandro Dalcin PetscReal Alpha_f = th->Alpha_f; 79b07a2398SLisandro Dalcin PetscReal Gamma = th->Gamma; 80b07a2398SLisandro Dalcin 81b07a2398SLisandro Dalcin PetscFunctionBegin; 82b07a2398SLisandro Dalcin /* V1 = 1/(Gamma*dT)*(X1-X0) + (1-1/Gamma)*V0 */ 839566063dSJacob Faibussowitsch PetscCall(VecWAXPY(V1, -1.0, X0, X1)); 849566063dSJacob Faibussowitsch PetscCall(VecAXPBY(V1, 1 - 1 / Gamma, 1 / (Gamma * dt), V0)); 85b07a2398SLisandro Dalcin /* Xa = X0 + Alpha_f*(X1-X0) */ 869566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Xa, -1.0, X0, X1)); 879566063dSJacob Faibussowitsch PetscCall(VecAYPX(Xa, Alpha_f, X0)); 88b07a2398SLisandro Dalcin /* Va = V0 + Alpha_m*(V1-V0) */ 899566063dSJacob Faibussowitsch PetscCall(VecWAXPY(Va, -1.0, V0, V1)); 909566063dSJacob Faibussowitsch PetscCall(VecAYPX(Va, Alpha_m, V0)); 913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 92b07a2398SLisandro Dalcin } 93b07a2398SLisandro Dalcin 94d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_SNESSolve(TS ts, Vec b, Vec x) 95d71ae5a4SJacob Faibussowitsch { 96b07a2398SLisandro Dalcin PetscInt nits, lits; 97b07a2398SLisandro Dalcin 98b07a2398SLisandro Dalcin PetscFunctionBegin; 999566063dSJacob Faibussowitsch PetscCall(SNESSolve(ts->snes, b, x)); 1009566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(ts->snes, &nits)); 1019566063dSJacob Faibussowitsch PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits)); 1029371c9d4SSatish Balay ts->snes_its += nits; 1039371c9d4SSatish Balay ts->ksp_its += lits; 1043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 105b07a2398SLisandro Dalcin } 106b07a2398SLisandro Dalcin 107b07a2398SLisandro Dalcin /* 108b07a2398SLisandro Dalcin Compute a consistent initial state for the generalized-alpha method. 109b07a2398SLisandro Dalcin - Solve two successive backward Euler steps with halved time step. 110b07a2398SLisandro Dalcin - Compute the initial time derivative using backward differences. 111b07a2398SLisandro Dalcin - If using adaptivity, estimate the LTE of the initial step. 112b07a2398SLisandro Dalcin */ 113d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_Restart(TS ts, PetscBool *initok) 114d71ae5a4SJacob Faibussowitsch { 115b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 116b07a2398SLisandro Dalcin PetscReal time_step; 117b07a2398SLisandro Dalcin PetscReal alpha_m, alpha_f, gamma; 118b07a2398SLisandro Dalcin Vec X0 = ts->vec_sol, X1, X2 = th->X1; 119b07a2398SLisandro Dalcin PetscBool stageok; 120b07a2398SLisandro Dalcin 121b07a2398SLisandro Dalcin PetscFunctionBegin; 1229566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X0, &X1)); 123b07a2398SLisandro Dalcin 124b07a2398SLisandro Dalcin /* Setup backward Euler with halved time step */ 1259566063dSJacob Faibussowitsch PetscCall(TSAlphaGetParams(ts, &alpha_m, &alpha_f, &gamma)); 1269566063dSJacob Faibussowitsch PetscCall(TSAlphaSetParams(ts, 1, 1, 1)); 1279566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &time_step)); 128b07a2398SLisandro Dalcin ts->time_step = time_step / 2; 1299566063dSJacob Faibussowitsch PetscCall(TSAlpha_StageTime(ts)); 130b07a2398SLisandro Dalcin th->stage_time = ts->ptime; 1319566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->V0)); 132b07a2398SLisandro Dalcin 133b07a2398SLisandro Dalcin /* First BE step, (t0,X0) -> (t1,X1) */ 134b07a2398SLisandro Dalcin th->stage_time += ts->time_step; 1359566063dSJacob Faibussowitsch PetscCall(VecCopy(X0, th->X0)); 1369566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, th->stage_time)); 1379566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0, X1)); 1389566063dSJacob Faibussowitsch PetscCall(TSAlpha_SNESSolve(ts, NULL, X1)); 1399566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, th->stage_time, 0, &X1)); 1409566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, X1, &stageok)); 141b07a2398SLisandro Dalcin if (!stageok) goto finally; 142b07a2398SLisandro Dalcin 143b07a2398SLisandro Dalcin /* Second BE step, (t1,X1) -> (t2,X2) */ 144b07a2398SLisandro Dalcin th->stage_time += ts->time_step; 1459566063dSJacob Faibussowitsch PetscCall(VecCopy(X1, th->X0)); 1469566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, th->stage_time)); 1479566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0, X2)); 1489566063dSJacob Faibussowitsch PetscCall(TSAlpha_SNESSolve(ts, NULL, X2)); 1499566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, th->stage_time, 0, &X2)); 1509566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, X2, &stageok)); 151b07a2398SLisandro Dalcin if (!stageok) goto finally; 152b07a2398SLisandro Dalcin 153b07a2398SLisandro Dalcin /* Compute V0 ~ dX/dt at t0 with backward differences */ 1549566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->V0)); 1559566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->V0, -3 / ts->time_step, X0)); 1569566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->V0, +4 / ts->time_step, X1)); 1579566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->V0, -1 / ts->time_step, X2)); 158b07a2398SLisandro Dalcin 159b07a2398SLisandro Dalcin /* Rough, lower-order estimate LTE of the initial step */ 1602ffb9264SLisandro Dalcin if (th->vec_lte_work) { 1619566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(th->vec_lte_work)); 1629566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->vec_lte_work, +2, X2)); 1639566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->vec_lte_work, -4, X1)); 1649566063dSJacob Faibussowitsch PetscCall(VecAXPY(th->vec_lte_work, +2, X0)); 165b07a2398SLisandro Dalcin } 166b07a2398SLisandro Dalcin 167b07a2398SLisandro Dalcin finally: 168b07a2398SLisandro Dalcin /* Revert TSAlpha to the initial state (t0,X0) */ 169b07a2398SLisandro Dalcin if (initok) *initok = stageok; 1709566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, time_step)); 1719566063dSJacob Faibussowitsch PetscCall(TSAlphaSetParams(ts, alpha_m, alpha_f, gamma)); 1729566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, th->X0)); 173b07a2398SLisandro Dalcin 1749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X1)); 1753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 176b07a2398SLisandro Dalcin } 177b07a2398SLisandro Dalcin 178d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Alpha(TS ts) 179d71ae5a4SJacob Faibussowitsch { 180b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 181b07a2398SLisandro Dalcin PetscInt rejections = 0; 182b07a2398SLisandro Dalcin PetscBool stageok, accept = PETSC_TRUE; 183b07a2398SLisandro Dalcin PetscReal next_time_step = ts->time_step; 184b07a2398SLisandro Dalcin 185b07a2398SLisandro Dalcin PetscFunctionBegin; 1869566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(citation, &cited)); 187b07a2398SLisandro Dalcin 188b07a2398SLisandro Dalcin if (!ts->steprollback) { 1899566063dSJacob Faibussowitsch if (th->vec_sol_prev) PetscCall(VecCopy(th->X0, th->vec_sol_prev)); 1909566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, th->X0)); 1919566063dSJacob Faibussowitsch PetscCall(VecCopy(th->V1, th->V0)); 192b07a2398SLisandro Dalcin } 193b07a2398SLisandro Dalcin 1941566a47fSLisandro Dalcin th->status = TS_STEP_INCOMPLETE; 195b07a2398SLisandro Dalcin while (!ts->reason && th->status != TS_STEP_COMPLETE) { 196fecfb714SLisandro Dalcin if (ts->steprestart) { 1979566063dSJacob Faibussowitsch PetscCall(TSAlpha_Restart(ts, &stageok)); 198fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 199b07a2398SLisandro Dalcin } 200b07a2398SLisandro Dalcin 2019566063dSJacob Faibussowitsch PetscCall(TSAlpha_StageTime(ts)); 2029566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0, th->X1)); 2039566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, th->stage_time)); 2049566063dSJacob Faibussowitsch PetscCall(TSAlpha_SNESSolve(ts, NULL, th->X1)); 2059566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, th->stage_time, 0, &th->Xa)); 2069566063dSJacob Faibussowitsch PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->Xa, &stageok)); 207fecfb714SLisandro Dalcin if (!stageok) goto reject_step; 208b07a2398SLisandro Dalcin 2091566a47fSLisandro Dalcin th->status = TS_STEP_PENDING; 2109566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X1, ts->vec_sol)); 2119566063dSJacob Faibussowitsch PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept)); 2121566a47fSLisandro Dalcin th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 213b07a2398SLisandro Dalcin if (!accept) { 2149566063dSJacob Faibussowitsch PetscCall(VecCopy(th->X0, ts->vec_sol)); 215be5899b3SLisandro Dalcin ts->time_step = next_time_step; 216be5899b3SLisandro Dalcin goto reject_step; 217b07a2398SLisandro Dalcin } 218b07a2398SLisandro Dalcin 219b07a2398SLisandro Dalcin ts->ptime += ts->time_step; 220b07a2398SLisandro Dalcin ts->time_step = next_time_step; 221b07a2398SLisandro Dalcin break; 222b07a2398SLisandro Dalcin 223b07a2398SLisandro Dalcin reject_step: 2249371c9d4SSatish Balay ts->reject++; 2259371c9d4SSatish Balay accept = PETSC_FALSE; 226b07a2398SLisandro Dalcin if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 227b07a2398SLisandro Dalcin ts->reason = TS_DIVERGED_STEP_REJECTED; 22863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections)); 229b07a2398SLisandro Dalcin } 230b07a2398SLisandro Dalcin } 2313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 232b07a2398SLisandro Dalcin } 233b07a2398SLisandro Dalcin 234d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateWLTE_Alpha(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte) 235d71ae5a4SJacob Faibussowitsch { 236b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 2379808bdc1SLisandro Dalcin Vec X = th->X1; /* X = solution */ 2381566a47fSLisandro Dalcin Vec Y = th->vec_lte_work; /* Y = X + LTE */ 2397453f775SEmil Constantinescu PetscReal wltea, wlter; 240b07a2398SLisandro Dalcin 241b07a2398SLisandro Dalcin PetscFunctionBegin; 2429371c9d4SSatish Balay if (!th->vec_sol_prev) { 2439371c9d4SSatish Balay *wlte = -1; 2443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2459371c9d4SSatish Balay } 2469371c9d4SSatish Balay if (!th->vec_lte_work) { 2479371c9d4SSatish Balay *wlte = -1; 2483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2499371c9d4SSatish Balay } 250fecfb714SLisandro Dalcin if (ts->steprestart) { 251fecfb714SLisandro Dalcin /* th->vec_lte_work is set to the LTE in TSAlpha_Restart() */ 2529566063dSJacob Faibussowitsch PetscCall(VecAXPY(Y, 1, X)); 253b07a2398SLisandro Dalcin } else { 254b07a2398SLisandro Dalcin /* Compute LTE using backward differences with non-constant time step */ 255be5899b3SLisandro Dalcin PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 256be5899b3SLisandro Dalcin PetscReal a = 1 + h_prev / h; 2579371c9d4SSatish Balay PetscScalar scal[3]; 2589371c9d4SSatish Balay Vec vecs[3]; 2599371c9d4SSatish Balay scal[0] = +1 / a; 2609371c9d4SSatish Balay scal[1] = -1 / (a - 1); 2619371c9d4SSatish Balay scal[2] = +1 / (a * (a - 1)); 2629371c9d4SSatish Balay vecs[0] = th->X1; 2639371c9d4SSatish Balay vecs[1] = th->X0; 2649371c9d4SSatish Balay vecs[2] = th->vec_sol_prev; 2659566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Y)); 2669566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Y, 3, scal, vecs)); 267b07a2398SLisandro Dalcin } 2689566063dSJacob Faibussowitsch PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, wlte, &wltea, &wlter)); 2699808bdc1SLisandro Dalcin if (order) *order = 2; 2703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 271b07a2398SLisandro Dalcin } 272b07a2398SLisandro Dalcin 273d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_Alpha(TS ts, PetscReal t, Vec X) 274d71ae5a4SJacob Faibussowitsch { 275b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 276b07a2398SLisandro Dalcin PetscReal dt = t - ts->ptime; 277*51adc54cSHong Zhang PetscReal Gamma = th->Gamma; 278b07a2398SLisandro Dalcin 279b07a2398SLisandro Dalcin PetscFunctionBegin; 280*51adc54cSHong Zhang PetscCall(VecWAXPY(th->V1, -1.0, th->X0, ts->vec_sol)); 281*51adc54cSHong Zhang PetscCall(VecAXPBY(th->V1, 1 - 1 / Gamma, 1 / (Gamma * ts->time_step), th->V0)); 2829566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 283*51adc54cSHong Zhang /* X = X + Gamma*dT*V1 */ 2849566063dSJacob Faibussowitsch PetscCall(VecAXPY(X, th->Gamma * dt, th->V1)); 285*51adc54cSHong Zhang /* X = X + (1-Gamma)*dT*V0 */ 2869566063dSJacob Faibussowitsch PetscCall(VecAXPY(X, (1 - th->Gamma) * dt, th->V0)); 2873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 288b07a2398SLisandro Dalcin } 289b07a2398SLisandro Dalcin 290d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes, Vec X, Vec F, TS ts) 291d71ae5a4SJacob Faibussowitsch { 292b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 293b07a2398SLisandro Dalcin PetscReal ta = th->stage_time; 294b07a2398SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va; 295b07a2398SLisandro Dalcin 296b07a2398SLisandro Dalcin PetscFunctionBegin; 2979566063dSJacob Faibussowitsch PetscCall(TSAlpha_StageVecs(ts, X)); 298b07a2398SLisandro Dalcin /* F = Function(ta,Xa,Va) */ 2999566063dSJacob Faibussowitsch PetscCall(TSComputeIFunction(ts, ta, Xa, Va, F, PETSC_FALSE)); 3009566063dSJacob Faibussowitsch PetscCall(VecScale(F, th->scale_F)); 3013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 302b07a2398SLisandro Dalcin } 303b07a2398SLisandro Dalcin 304d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes, PETSC_UNUSED Vec X, Mat J, Mat P, TS ts) 305d71ae5a4SJacob Faibussowitsch { 306b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 307b07a2398SLisandro Dalcin PetscReal ta = th->stage_time; 308b07a2398SLisandro Dalcin Vec Xa = th->Xa, Va = th->Va; 309b07a2398SLisandro Dalcin PetscReal dVdX = th->shift_V; 310b07a2398SLisandro Dalcin 311b07a2398SLisandro Dalcin PetscFunctionBegin; 312b07a2398SLisandro Dalcin /* J,P = Jacobian(ta,Xa,Va) */ 3139566063dSJacob Faibussowitsch PetscCall(TSComputeIJacobian(ts, ta, Xa, Va, dVdX, J, P, PETSC_FALSE)); 3143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 315b07a2398SLisandro Dalcin } 316b07a2398SLisandro Dalcin 317d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Alpha(TS ts) 318d71ae5a4SJacob Faibussowitsch { 319b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 320b07a2398SLisandro Dalcin 321b07a2398SLisandro Dalcin PetscFunctionBegin; 3229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->X0)); 3239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->Xa)); 3249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->X1)); 3259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->V0)); 3269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->Va)); 3279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->V1)); 3289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_sol_prev)); 3299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&th->vec_lte_work)); 3303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 331b07a2398SLisandro Dalcin } 332b07a2398SLisandro Dalcin 333d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Alpha(TS ts) 334d71ae5a4SJacob Faibussowitsch { 335b07a2398SLisandro Dalcin PetscFunctionBegin; 3369566063dSJacob Faibussowitsch PetscCall(TSReset_Alpha(ts)); 3379566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 338b07a2398SLisandro Dalcin 3399566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetRadius_C", NULL)); 3409566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetParams_C", NULL)); 3419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaGetParams_C", NULL)); 3423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 343b07a2398SLisandro Dalcin } 344b07a2398SLisandro Dalcin 345d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Alpha(TS ts) 346d71ae5a4SJacob Faibussowitsch { 347b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 3482ffb9264SLisandro Dalcin PetscBool match; 349b07a2398SLisandro Dalcin 350b07a2398SLisandro Dalcin PetscFunctionBegin; 3518ec9177eSStefano Zampini if (!th->X0) PetscCall(VecDuplicate(ts->vec_sol, &th->X0)); 3529566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->Xa)); 3539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->X1)); 3549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->V0)); 3559566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->Va)); 3569566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &th->V1)); 3571566a47fSLisandro Dalcin 3589566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &ts->adapt)); 3599566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesClear(ts->adapt)); 3609566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match)); 3612ffb9264SLisandro Dalcin if (!match) { 362ecc87898SStefano Zampini if (!th->vec_sol_prev) PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev)); 363ecc87898SStefano Zampini if (!th->vec_lte_work) PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work)); 364b07a2398SLisandro Dalcin } 3651566a47fSLisandro Dalcin 3669566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &ts->snes)); 3673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 368b07a2398SLisandro Dalcin } 369b07a2398SLisandro Dalcin 370d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_Alpha(TS ts, PetscOptionItems *PetscOptionsObject) 371d71ae5a4SJacob Faibussowitsch { 372b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 373b07a2398SLisandro Dalcin 374b07a2398SLisandro Dalcin PetscFunctionBegin; 375d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Generalized-Alpha ODE solver options"); 376b07a2398SLisandro Dalcin { 377b07a2398SLisandro Dalcin PetscBool flg; 378b07a2398SLisandro Dalcin PetscReal radius = 1; 3799566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_radius", "Spectral radius (high-frequency dissipation)", "TSAlphaSetRadius", radius, &radius, &flg)); 3809566063dSJacob Faibussowitsch if (flg) PetscCall(TSAlphaSetRadius(ts, radius)); 3819566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_alpha_m", "Algorithmic parameter alpha_m", "TSAlphaSetParams", th->Alpha_m, &th->Alpha_m, NULL)); 3829566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_alpha_f", "Algorithmic parameter alpha_f", "TSAlphaSetParams", th->Alpha_f, &th->Alpha_f, NULL)); 3839566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ts_alpha_gamma", "Algorithmic parameter gamma", "TSAlphaSetParams", th->Gamma, &th->Gamma, NULL)); 3849566063dSJacob Faibussowitsch PetscCall(TSAlphaSetParams(ts, th->Alpha_m, th->Alpha_f, th->Gamma)); 385b07a2398SLisandro Dalcin } 386d0609cedSBarry Smith PetscOptionsHeadEnd(); 3873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 388b07a2398SLisandro Dalcin } 389b07a2398SLisandro Dalcin 390d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Alpha(TS ts, PetscViewer viewer) 391d71ae5a4SJacob Faibussowitsch { 392b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 3939c334d8fSLisandro Dalcin PetscBool iascii; 394b07a2398SLisandro Dalcin 395b07a2398SLisandro Dalcin PetscFunctionBegin; 3969566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 39748a46eb9SPierre Jolivet if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, " Alpha_m=%g, Alpha_f=%g, Gamma=%g\n", (double)th->Alpha_m, (double)th->Alpha_f, (double)th->Gamma)); 3983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 399b07a2398SLisandro Dalcin } 400b07a2398SLisandro Dalcin 401d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlphaSetRadius_Alpha(TS ts, PetscReal radius) 402d71ae5a4SJacob Faibussowitsch { 403b07a2398SLisandro Dalcin PetscReal alpha_m, alpha_f, gamma; 404b07a2398SLisandro Dalcin 405b07a2398SLisandro Dalcin PetscFunctionBegin; 406cad9d221SBarry Smith PetscCheck(radius >= 0 && radius <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius); 407b07a2398SLisandro Dalcin alpha_m = (PetscReal)0.5 * (3 - radius) / (1 + radius); 408b07a2398SLisandro Dalcin alpha_f = 1 / (1 + radius); 409b07a2398SLisandro Dalcin gamma = (PetscReal)0.5 + alpha_m - alpha_f; 4109566063dSJacob Faibussowitsch PetscCall(TSAlphaSetParams(ts, alpha_m, alpha_f, gamma)); 4113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 412b07a2398SLisandro Dalcin } 413b07a2398SLisandro Dalcin 414d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlphaSetParams_Alpha(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma) 415d71ae5a4SJacob Faibussowitsch { 416b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 417b07a2398SLisandro Dalcin PetscReal tol = 100 * PETSC_MACHINE_EPSILON; 418b07a2398SLisandro Dalcin PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma; 419b07a2398SLisandro Dalcin 420b07a2398SLisandro Dalcin PetscFunctionBegin; 421b07a2398SLisandro Dalcin th->Alpha_m = alpha_m; 422b07a2398SLisandro Dalcin th->Alpha_f = alpha_f; 423b07a2398SLisandro Dalcin th->Gamma = gamma; 424b07a2398SLisandro Dalcin th->order = (PetscAbsReal(res) < tol) ? 2 : 1; 4253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 426b07a2398SLisandro Dalcin } 427b07a2398SLisandro Dalcin 428d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlphaGetParams_Alpha(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma) 429d71ae5a4SJacob Faibussowitsch { 430b07a2398SLisandro Dalcin TS_Alpha *th = (TS_Alpha *)ts->data; 431b07a2398SLisandro Dalcin 432b07a2398SLisandro Dalcin PetscFunctionBegin; 433b07a2398SLisandro Dalcin if (alpha_m) *alpha_m = th->Alpha_m; 434b07a2398SLisandro Dalcin if (alpha_f) *alpha_f = th->Alpha_f; 435b07a2398SLisandro Dalcin if (gamma) *gamma = th->Gamma; 4363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 437b07a2398SLisandro Dalcin } 438b07a2398SLisandro Dalcin 439b07a2398SLisandro Dalcin /*MC 4401d27aa22SBarry Smith TSALPHA - ODE/DAE solver using the implicit Generalized-Alpha method {cite}`jansen_2000` {cite}`chung1993` for first-order systems 441b07a2398SLisandro Dalcin 442b07a2398SLisandro Dalcin Level: beginner 443b07a2398SLisandro Dalcin 4441cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSCreate()`, `TSSetType()`, `TSAlphaSetRadius()`, `TSAlphaSetParams()` 445b07a2398SLisandro Dalcin M*/ 446d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Alpha(TS ts) 447d71ae5a4SJacob Faibussowitsch { 448b07a2398SLisandro Dalcin TS_Alpha *th; 449b07a2398SLisandro Dalcin 450b07a2398SLisandro Dalcin PetscFunctionBegin; 451b07a2398SLisandro Dalcin ts->ops->reset = TSReset_Alpha; 452b07a2398SLisandro Dalcin ts->ops->destroy = TSDestroy_Alpha; 453b07a2398SLisandro Dalcin ts->ops->view = TSView_Alpha; 454b07a2398SLisandro Dalcin ts->ops->setup = TSSetUp_Alpha; 455b07a2398SLisandro Dalcin ts->ops->setfromoptions = TSSetFromOptions_Alpha; 456b07a2398SLisandro Dalcin ts->ops->step = TSStep_Alpha; 4579808bdc1SLisandro Dalcin ts->ops->evaluatewlte = TSEvaluateWLTE_Alpha; 458b07a2398SLisandro Dalcin ts->ops->interpolate = TSInterpolate_Alpha; 4598ec9177eSStefano Zampini ts->ops->resizeregister = TSResizeRegister_Alpha; 460b07a2398SLisandro Dalcin ts->ops->snesfunction = SNESTSFormFunction_Alpha; 461b07a2398SLisandro Dalcin ts->ops->snesjacobian = SNESTSFormJacobian_Alpha; 4622ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 463b07a2398SLisandro Dalcin 464efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 465efd4aadfSBarry Smith 4664dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&th)); 467b07a2398SLisandro Dalcin ts->data = (void *)th; 468b07a2398SLisandro Dalcin 469b07a2398SLisandro Dalcin th->Alpha_m = 0.5; 470b07a2398SLisandro Dalcin th->Alpha_f = 0.5; 471b07a2398SLisandro Dalcin th->Gamma = 0.5; 472b07a2398SLisandro Dalcin th->order = 2; 473b07a2398SLisandro Dalcin 4749566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetRadius_C", TSAlphaSetRadius_Alpha)); 4759566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetParams_C", TSAlphaSetParams_Alpha)); 4769566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaGetParams_C", TSAlphaGetParams_Alpha)); 4773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 478b07a2398SLisandro Dalcin } 479b07a2398SLisandro Dalcin 480b07a2398SLisandro Dalcin /*@ 481bcf0153eSBarry Smith TSAlphaSetRadius - sets the desired spectral radius of the method for `TSALPHA` 482b07a2398SLisandro Dalcin (i.e. high-frequency numerical damping) 483b07a2398SLisandro Dalcin 484c3339decSBarry Smith Logically Collective 485b07a2398SLisandro Dalcin 486d8d19677SJose E. Roman Input Parameters: 487b07a2398SLisandro Dalcin + ts - timestepping context 488b07a2398SLisandro Dalcin - radius - the desired spectral radius 489b07a2398SLisandro Dalcin 490bcf0153eSBarry Smith Options Database Key: 49167b8a455SSatish Balay . -ts_alpha_radius <radius> - set alpha radius 492b07a2398SLisandro Dalcin 493b07a2398SLisandro Dalcin Level: intermediate 494b07a2398SLisandro Dalcin 49514d0ab18SJacob Faibussowitsch Notes: 49614d0ab18SJacob Faibussowitsch The algorithmic parameters $\alpha_m$ and $\alpha_f$ of the generalized-$\alpha$ method can 49714d0ab18SJacob Faibussowitsch be computed in terms of a specified spectral radius $\rho$ in [0, 1] for infinite time step 49814d0ab18SJacob Faibussowitsch in order to control high-frequency numerical damping\: 499562efe2eSBarry Smith 50014d0ab18SJacob Faibussowitsch $$ 501562efe2eSBarry Smith \begin{align*} 502562efe2eSBarry Smith \alpha_m = 0.5*(3-\rho)/(1+\rho) \\ 50314d0ab18SJacob Faibussowitsch \alpha_f = 1/(1+\rho) 504562efe2eSBarry Smith \end{align*} 50514d0ab18SJacob Faibussowitsch $$ 50614d0ab18SJacob Faibussowitsch 5071cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA`, `TSAlphaSetParams()`, `TSAlphaGetParams()` 508b07a2398SLisandro Dalcin @*/ 509d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlphaSetRadius(TS ts, PetscReal radius) 510d71ae5a4SJacob Faibussowitsch { 511b07a2398SLisandro Dalcin PetscFunctionBegin; 512b07a2398SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 513b07a2398SLisandro Dalcin PetscValidLogicalCollectiveReal(ts, radius, 2); 514cad9d221SBarry Smith PetscCheck(radius >= 0 && radius <= 1, ((PetscObject)ts)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius); 515cac4c232SBarry Smith PetscTryMethod(ts, "TSAlphaSetRadius_C", (TS, PetscReal), (ts, radius)); 5163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 517b07a2398SLisandro Dalcin } 518b07a2398SLisandro Dalcin 519b07a2398SLisandro Dalcin /*@ 520bcf0153eSBarry Smith TSAlphaSetParams - sets the algorithmic parameters for `TSALPHA` 521b07a2398SLisandro Dalcin 522c3339decSBarry Smith Logically Collective 523b07a2398SLisandro Dalcin 524d8d19677SJose E. Roman Input Parameters: 525b07a2398SLisandro Dalcin + ts - timestepping context 5262fe279fdSBarry Smith . alpha_m - algorithmic parameter 5272fe279fdSBarry Smith . alpha_f - algorithmic parameter 5282fe279fdSBarry Smith - gamma - algorithmic parameter 529b07a2398SLisandro Dalcin 530bcf0153eSBarry Smith Options Database Keys: 53167b8a455SSatish Balay + -ts_alpha_alpha_m <alpha_m> - set alpha_m 53267b8a455SSatish Balay . -ts_alpha_alpha_f <alpha_f> - set alpha_f 53367b8a455SSatish Balay - -ts_alpha_gamma <gamma> - set gamma 534b07a2398SLisandro Dalcin 535bcf0153eSBarry Smith Level: advanced 536bcf0153eSBarry Smith 537b07a2398SLisandro Dalcin Note: 538562efe2eSBarry Smith Second-order accuracy can be obtained so long as\: $\gamma = 0.5 + \alpha_m - \alpha_f$ 53914d0ab18SJacob Faibussowitsch 540562efe2eSBarry Smith Unconditional stability requires\: $\alpha_m >= \alpha_f >= 0.5$ 54114d0ab18SJacob Faibussowitsch 542562efe2eSBarry Smith Backward Euler method is recovered with\: $\alpha_m = \alpha_f = \gamma = 1$ 54314d0ab18SJacob Faibussowitsch 54414d0ab18SJacob Faibussowitsch Use of this function is normally only required to hack `TSALPHA` to use a modified 54514d0ab18SJacob Faibussowitsch integration scheme. Users should call `TSAlphaSetRadius()` to set the desired spectral radius 54614d0ab18SJacob Faibussowitsch of the methods (i.e. high-frequency damping) in order so select optimal values for these 54714d0ab18SJacob Faibussowitsch parameters. 548b07a2398SLisandro Dalcin 5491cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA`, `TSAlphaSetRadius()`, `TSAlphaGetParams()` 550b07a2398SLisandro Dalcin @*/ 551d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlphaSetParams(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma) 552d71ae5a4SJacob Faibussowitsch { 553b07a2398SLisandro Dalcin PetscFunctionBegin; 554b07a2398SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 555b07a2398SLisandro Dalcin PetscValidLogicalCollectiveReal(ts, alpha_m, 2); 556b07a2398SLisandro Dalcin PetscValidLogicalCollectiveReal(ts, alpha_f, 3); 557b07a2398SLisandro Dalcin PetscValidLogicalCollectiveReal(ts, gamma, 4); 558cac4c232SBarry Smith PetscTryMethod(ts, "TSAlphaSetParams_C", (TS, PetscReal, PetscReal, PetscReal), (ts, alpha_m, alpha_f, gamma)); 5593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 560b07a2398SLisandro Dalcin } 561b07a2398SLisandro Dalcin 562b07a2398SLisandro Dalcin /*@ 563bcf0153eSBarry Smith TSAlphaGetParams - gets the algorithmic parameters for `TSALPHA` 564b07a2398SLisandro Dalcin 565b07a2398SLisandro Dalcin Not Collective 566b07a2398SLisandro Dalcin 567b07a2398SLisandro Dalcin Input Parameter: 568b07a2398SLisandro Dalcin . ts - timestepping context 569b07a2398SLisandro Dalcin 570b07a2398SLisandro Dalcin Output Parameters: 5712fe279fdSBarry Smith + alpha_m - algorithmic parameter 5722fe279fdSBarry Smith . alpha_f - algorithmic parameter 5732fe279fdSBarry Smith - gamma - algorithmic parameter 574b07a2398SLisandro Dalcin 575bcf0153eSBarry Smith Level: advanced 576bcf0153eSBarry Smith 577b07a2398SLisandro Dalcin Note: 57814d0ab18SJacob Faibussowitsch Use of this function is normally only required to hack `TSALPHA` to use a modified 57914d0ab18SJacob Faibussowitsch integration scheme. Users should call `TSAlphaSetRadius()` to set the high-frequency damping 58014d0ab18SJacob Faibussowitsch (i.e. spectral radius of the method) in order so select optimal values for these parameters. 581b07a2398SLisandro Dalcin 5821cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA`, `TSAlphaSetRadius()`, `TSAlphaSetParams()` 583b07a2398SLisandro Dalcin @*/ 584d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlphaGetParams(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma) 585d71ae5a4SJacob Faibussowitsch { 586b07a2398SLisandro Dalcin PetscFunctionBegin; 587b07a2398SLisandro Dalcin PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 5884f572ea9SToby Isaac if (alpha_m) PetscAssertPointer(alpha_m, 2); 5894f572ea9SToby Isaac if (alpha_f) PetscAssertPointer(alpha_f, 3); 5904f572ea9SToby Isaac if (gamma) PetscAssertPointer(gamma, 4); 591cac4c232SBarry Smith PetscUseMethod(ts, "TSAlphaGetParams_C", (TS, PetscReal *, PetscReal *, PetscReal *), (ts, alpha_m, alpha_f, gamma)); 5923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 593b07a2398SLisandro Dalcin } 594