xref: /petsc/src/ts/impls/implicit/alpha/alpha1.c (revision 562efe2ef922487c6beae96ba39e19afd4eefbe6)
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 /* We need to transfer X0 which will be copied into sol_prev */
388ec9177eSStefano Zampini static PetscErrorCode TSResizeRegister_Alpha(TS ts, PetscBool reg)
398ec9177eSStefano Zampini {
408ec9177eSStefano Zampini   TS_Alpha  *th     = (TS_Alpha *)ts->data;
418ec9177eSStefano Zampini   const char name[] = "ts:alpha:X0";
428ec9177eSStefano Zampini 
438ec9177eSStefano Zampini   PetscFunctionBegin;
448ec9177eSStefano Zampini   if (reg && th->vec_sol_prev) {
458ec9177eSStefano Zampini     PetscCall(TSResizeRegisterVec(ts, name, th->X0));
468ec9177eSStefano Zampini   } else if (!reg) {
478ec9177eSStefano Zampini     PetscCall(TSResizeRetrieveVec(ts, name, &th->X0));
488ec9177eSStefano Zampini     PetscCall(PetscObjectReference((PetscObject)th->X0));
498ec9177eSStefano Zampini   }
508ec9177eSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
518ec9177eSStefano Zampini }
528ec9177eSStefano Zampini 
53d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_StageTime(TS ts)
54d71ae5a4SJacob Faibussowitsch {
55b07a2398SLisandro Dalcin   TS_Alpha *th      = (TS_Alpha *)ts->data;
56b07a2398SLisandro Dalcin   PetscReal t       = ts->ptime;
57b07a2398SLisandro Dalcin   PetscReal dt      = ts->time_step;
58b07a2398SLisandro Dalcin   PetscReal Alpha_m = th->Alpha_m;
59b07a2398SLisandro Dalcin   PetscReal Alpha_f = th->Alpha_f;
60b07a2398SLisandro Dalcin   PetscReal Gamma   = th->Gamma;
61b07a2398SLisandro Dalcin 
62b07a2398SLisandro Dalcin   PetscFunctionBegin;
63b07a2398SLisandro Dalcin   th->stage_time = t + Alpha_f * dt;
64b07a2398SLisandro Dalcin   th->shift_V    = Alpha_m / (Alpha_f * Gamma * dt);
65b07a2398SLisandro Dalcin   th->scale_F    = 1 / Alpha_f;
663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
67b07a2398SLisandro Dalcin }
68b07a2398SLisandro Dalcin 
69d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_StageVecs(TS ts, Vec X)
70d71ae5a4SJacob Faibussowitsch {
71b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
72b07a2398SLisandro Dalcin   Vec       X1 = X, V1 = th->V1;
73b07a2398SLisandro Dalcin   Vec       Xa = th->Xa, Va = th->Va;
74b07a2398SLisandro Dalcin   Vec       X0 = th->X0, V0 = th->V0;
75b07a2398SLisandro Dalcin   PetscReal dt      = ts->time_step;
76b07a2398SLisandro Dalcin   PetscReal Alpha_m = th->Alpha_m;
77b07a2398SLisandro Dalcin   PetscReal Alpha_f = th->Alpha_f;
78b07a2398SLisandro Dalcin   PetscReal Gamma   = th->Gamma;
79b07a2398SLisandro Dalcin 
80b07a2398SLisandro Dalcin   PetscFunctionBegin;
81b07a2398SLisandro Dalcin   /* V1 = 1/(Gamma*dT)*(X1-X0) + (1-1/Gamma)*V0 */
829566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(V1, -1.0, X0, X1));
839566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(V1, 1 - 1 / Gamma, 1 / (Gamma * dt), V0));
84b07a2398SLisandro Dalcin   /* Xa = X0 + Alpha_f*(X1-X0) */
859566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Xa, -1.0, X0, X1));
869566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Xa, Alpha_f, X0));
87b07a2398SLisandro Dalcin   /* Va = V0 + Alpha_m*(V1-V0) */
889566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Va, -1.0, V0, V1));
899566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Va, Alpha_m, V0));
903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
91b07a2398SLisandro Dalcin }
92b07a2398SLisandro Dalcin 
93d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_SNESSolve(TS ts, Vec b, Vec x)
94d71ae5a4SJacob Faibussowitsch {
95b07a2398SLisandro Dalcin   PetscInt nits, lits;
96b07a2398SLisandro Dalcin 
97b07a2398SLisandro Dalcin   PetscFunctionBegin;
989566063dSJacob Faibussowitsch   PetscCall(SNESSolve(ts->snes, b, x));
999566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(ts->snes, &nits));
1009566063dSJacob Faibussowitsch   PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits));
1019371c9d4SSatish Balay   ts->snes_its += nits;
1029371c9d4SSatish Balay   ts->ksp_its += lits;
1033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
104b07a2398SLisandro Dalcin }
105b07a2398SLisandro Dalcin 
106b07a2398SLisandro Dalcin /*
107b07a2398SLisandro Dalcin   Compute a consistent initial state for the generalized-alpha method.
108b07a2398SLisandro Dalcin   - Solve two successive backward Euler steps with halved time step.
109b07a2398SLisandro Dalcin   - Compute the initial time derivative using backward differences.
110b07a2398SLisandro Dalcin   - If using adaptivity, estimate the LTE of the initial step.
111b07a2398SLisandro Dalcin */
112d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_Restart(TS ts, PetscBool *initok)
113d71ae5a4SJacob Faibussowitsch {
114b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
115b07a2398SLisandro Dalcin   PetscReal time_step;
116b07a2398SLisandro Dalcin   PetscReal alpha_m, alpha_f, gamma;
117b07a2398SLisandro Dalcin   Vec       X0 = ts->vec_sol, X1, X2 = th->X1;
118b07a2398SLisandro Dalcin   PetscBool stageok;
119b07a2398SLisandro Dalcin 
120b07a2398SLisandro Dalcin   PetscFunctionBegin;
1219566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X0, &X1));
122b07a2398SLisandro Dalcin 
123b07a2398SLisandro Dalcin   /* Setup backward Euler with halved time step */
1249566063dSJacob Faibussowitsch   PetscCall(TSAlphaGetParams(ts, &alpha_m, &alpha_f, &gamma));
1259566063dSJacob Faibussowitsch   PetscCall(TSAlphaSetParams(ts, 1, 1, 1));
1269566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &time_step));
127b07a2398SLisandro Dalcin   ts->time_step = time_step / 2;
1289566063dSJacob Faibussowitsch   PetscCall(TSAlpha_StageTime(ts));
129b07a2398SLisandro Dalcin   th->stage_time = ts->ptime;
1309566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(th->V0));
131b07a2398SLisandro Dalcin 
132b07a2398SLisandro Dalcin   /* First BE step, (t0,X0) -> (t1,X1) */
133b07a2398SLisandro Dalcin   th->stage_time += ts->time_step;
1349566063dSJacob Faibussowitsch   PetscCall(VecCopy(X0, th->X0));
1359566063dSJacob Faibussowitsch   PetscCall(TSPreStage(ts, th->stage_time));
1369566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->X0, X1));
1379566063dSJacob Faibussowitsch   PetscCall(TSAlpha_SNESSolve(ts, NULL, X1));
1389566063dSJacob Faibussowitsch   PetscCall(TSPostStage(ts, th->stage_time, 0, &X1));
1399566063dSJacob Faibussowitsch   PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, X1, &stageok));
140b07a2398SLisandro Dalcin   if (!stageok) goto finally;
141b07a2398SLisandro Dalcin 
142b07a2398SLisandro Dalcin   /* Second BE step, (t1,X1) -> (t2,X2) */
143b07a2398SLisandro Dalcin   th->stage_time += ts->time_step;
1449566063dSJacob Faibussowitsch   PetscCall(VecCopy(X1, th->X0));
1459566063dSJacob Faibussowitsch   PetscCall(TSPreStage(ts, th->stage_time));
1469566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->X0, X2));
1479566063dSJacob Faibussowitsch   PetscCall(TSAlpha_SNESSolve(ts, NULL, X2));
1489566063dSJacob Faibussowitsch   PetscCall(TSPostStage(ts, th->stage_time, 0, &X2));
1499566063dSJacob Faibussowitsch   PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, X2, &stageok));
150b07a2398SLisandro Dalcin   if (!stageok) goto finally;
151b07a2398SLisandro Dalcin 
152b07a2398SLisandro Dalcin   /* Compute V0 ~ dX/dt at t0 with backward differences */
1539566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(th->V0));
1549566063dSJacob Faibussowitsch   PetscCall(VecAXPY(th->V0, -3 / ts->time_step, X0));
1559566063dSJacob Faibussowitsch   PetscCall(VecAXPY(th->V0, +4 / ts->time_step, X1));
1569566063dSJacob Faibussowitsch   PetscCall(VecAXPY(th->V0, -1 / ts->time_step, X2));
157b07a2398SLisandro Dalcin 
158b07a2398SLisandro Dalcin   /* Rough, lower-order estimate LTE of the initial step */
1592ffb9264SLisandro Dalcin   if (th->vec_lte_work) {
1609566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(th->vec_lte_work));
1619566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work, +2, X2));
1629566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work, -4, X1));
1639566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work, +2, X0));
164b07a2398SLisandro Dalcin   }
165b07a2398SLisandro Dalcin 
166b07a2398SLisandro Dalcin finally:
167b07a2398SLisandro Dalcin   /* Revert TSAlpha to the initial state (t0,X0) */
168b07a2398SLisandro Dalcin   if (initok) *initok = stageok;
1699566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, time_step));
1709566063dSJacob Faibussowitsch   PetscCall(TSAlphaSetParams(ts, alpha_m, alpha_f, gamma));
1719566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol, th->X0));
172b07a2398SLisandro Dalcin 
1739566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X1));
1743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
175b07a2398SLisandro Dalcin }
176b07a2398SLisandro Dalcin 
177d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Alpha(TS ts)
178d71ae5a4SJacob Faibussowitsch {
179b07a2398SLisandro Dalcin   TS_Alpha *th         = (TS_Alpha *)ts->data;
180b07a2398SLisandro Dalcin   PetscInt  rejections = 0;
181b07a2398SLisandro Dalcin   PetscBool stageok, accept = PETSC_TRUE;
182b07a2398SLisandro Dalcin   PetscReal next_time_step = ts->time_step;
183b07a2398SLisandro Dalcin 
184b07a2398SLisandro Dalcin   PetscFunctionBegin;
1859566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(citation, &cited));
186b07a2398SLisandro Dalcin 
187b07a2398SLisandro Dalcin   if (!ts->steprollback) {
1889566063dSJacob Faibussowitsch     if (th->vec_sol_prev) PetscCall(VecCopy(th->X0, th->vec_sol_prev));
1899566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_sol, th->X0));
1909566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->V1, th->V0));
191b07a2398SLisandro Dalcin   }
192b07a2398SLisandro Dalcin 
1931566a47fSLisandro Dalcin   th->status = TS_STEP_INCOMPLETE;
194b07a2398SLisandro Dalcin   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
195fecfb714SLisandro Dalcin     if (ts->steprestart) {
1969566063dSJacob Faibussowitsch       PetscCall(TSAlpha_Restart(ts, &stageok));
197fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
198b07a2398SLisandro Dalcin     }
199b07a2398SLisandro Dalcin 
2009566063dSJacob Faibussowitsch     PetscCall(TSAlpha_StageTime(ts));
2019566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->X0, th->X1));
2029566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, th->stage_time));
2039566063dSJacob Faibussowitsch     PetscCall(TSAlpha_SNESSolve(ts, NULL, th->X1));
2049566063dSJacob Faibussowitsch     PetscCall(TSPostStage(ts, th->stage_time, 0, &th->Xa));
2059566063dSJacob Faibussowitsch     PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->Xa, &stageok));
206fecfb714SLisandro Dalcin     if (!stageok) goto reject_step;
207b07a2398SLisandro Dalcin 
2081566a47fSLisandro Dalcin     th->status = TS_STEP_PENDING;
2099566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->X1, ts->vec_sol));
2109566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
2111566a47fSLisandro Dalcin     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
212b07a2398SLisandro Dalcin     if (!accept) {
2139566063dSJacob Faibussowitsch       PetscCall(VecCopy(th->X0, ts->vec_sol));
214be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
215be5899b3SLisandro Dalcin       goto reject_step;
216b07a2398SLisandro Dalcin     }
217b07a2398SLisandro Dalcin 
218b07a2398SLisandro Dalcin     ts->ptime += ts->time_step;
219b07a2398SLisandro Dalcin     ts->time_step = next_time_step;
220b07a2398SLisandro Dalcin     break;
221b07a2398SLisandro Dalcin 
222b07a2398SLisandro Dalcin   reject_step:
2239371c9d4SSatish Balay     ts->reject++;
2249371c9d4SSatish Balay     accept = PETSC_FALSE;
225b07a2398SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
226b07a2398SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
22763a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
228b07a2398SLisandro Dalcin     }
229b07a2398SLisandro Dalcin   }
2303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
231b07a2398SLisandro Dalcin }
232b07a2398SLisandro Dalcin 
233d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateWLTE_Alpha(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte)
234d71ae5a4SJacob Faibussowitsch {
235b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
2369808bdc1SLisandro Dalcin   Vec       X  = th->X1;           /* X = solution */
2371566a47fSLisandro Dalcin   Vec       Y  = th->vec_lte_work; /* Y = X + LTE  */
2387453f775SEmil Constantinescu   PetscReal wltea, wlter;
239b07a2398SLisandro Dalcin 
240b07a2398SLisandro Dalcin   PetscFunctionBegin;
2419371c9d4SSatish Balay   if (!th->vec_sol_prev) {
2429371c9d4SSatish Balay     *wlte = -1;
2433ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2449371c9d4SSatish Balay   }
2459371c9d4SSatish Balay   if (!th->vec_lte_work) {
2469371c9d4SSatish Balay     *wlte = -1;
2473ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2489371c9d4SSatish Balay   }
249fecfb714SLisandro Dalcin   if (ts->steprestart) {
250fecfb714SLisandro Dalcin     /* th->vec_lte_work is set to the LTE in TSAlpha_Restart() */
2519566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Y, 1, X));
252b07a2398SLisandro Dalcin   } else {
253b07a2398SLisandro Dalcin     /* Compute LTE using backward differences with non-constant time step */
254be5899b3SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
255be5899b3SLisandro Dalcin     PetscReal   a = 1 + h_prev / h;
2569371c9d4SSatish Balay     PetscScalar scal[3];
2579371c9d4SSatish Balay     Vec         vecs[3];
2589371c9d4SSatish Balay     scal[0] = +1 / a;
2599371c9d4SSatish Balay     scal[1] = -1 / (a - 1);
2609371c9d4SSatish Balay     scal[2] = +1 / (a * (a - 1));
2619371c9d4SSatish Balay     vecs[0] = th->X1;
2629371c9d4SSatish Balay     vecs[1] = th->X0;
2639371c9d4SSatish Balay     vecs[2] = th->vec_sol_prev;
2649566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Y));
2659566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(Y, 3, scal, vecs));
266b07a2398SLisandro Dalcin   }
2679566063dSJacob Faibussowitsch   PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, wlte, &wltea, &wlter));
2689808bdc1SLisandro Dalcin   if (order) *order = 2;
2693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
270b07a2398SLisandro Dalcin }
271b07a2398SLisandro Dalcin 
272d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_Alpha(TS ts)
273d71ae5a4SJacob Faibussowitsch {
274b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
275b07a2398SLisandro Dalcin 
276b07a2398SLisandro Dalcin   PetscFunctionBegin;
2779566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->X0, ts->vec_sol));
2783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
279b07a2398SLisandro Dalcin }
280b07a2398SLisandro Dalcin 
281d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_Alpha(TS ts, PetscReal t, Vec X)
282d71ae5a4SJacob Faibussowitsch {
283b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
284b07a2398SLisandro Dalcin   PetscReal dt = t - ts->ptime;
285b07a2398SLisandro Dalcin 
286b07a2398SLisandro Dalcin   PetscFunctionBegin;
2879566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol, X));
2889566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, th->Gamma * dt, th->V1));
2899566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, (1 - th->Gamma) * dt, th->V0));
2903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
291b07a2398SLisandro Dalcin }
292b07a2398SLisandro Dalcin 
293d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes, Vec X, Vec F, TS ts)
294d71ae5a4SJacob Faibussowitsch {
295b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
296b07a2398SLisandro Dalcin   PetscReal ta = th->stage_time;
297b07a2398SLisandro Dalcin   Vec       Xa = th->Xa, Va = th->Va;
298b07a2398SLisandro Dalcin 
299b07a2398SLisandro Dalcin   PetscFunctionBegin;
3009566063dSJacob Faibussowitsch   PetscCall(TSAlpha_StageVecs(ts, X));
301b07a2398SLisandro Dalcin   /* F = Function(ta,Xa,Va) */
3029566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, ta, Xa, Va, F, PETSC_FALSE));
3039566063dSJacob Faibussowitsch   PetscCall(VecScale(F, th->scale_F));
3043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
305b07a2398SLisandro Dalcin }
306b07a2398SLisandro Dalcin 
307d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes, PETSC_UNUSED Vec X, Mat J, Mat P, TS ts)
308d71ae5a4SJacob Faibussowitsch {
309b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
310b07a2398SLisandro Dalcin   PetscReal ta = th->stage_time;
311b07a2398SLisandro Dalcin   Vec       Xa = th->Xa, Va = th->Va;
312b07a2398SLisandro Dalcin   PetscReal dVdX = th->shift_V;
313b07a2398SLisandro Dalcin 
314b07a2398SLisandro Dalcin   PetscFunctionBegin;
315b07a2398SLisandro Dalcin   /* J,P = Jacobian(ta,Xa,Va) */
3169566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, ta, Xa, Va, dVdX, J, P, PETSC_FALSE));
3173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
318b07a2398SLisandro Dalcin }
319b07a2398SLisandro Dalcin 
320d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Alpha(TS ts)
321d71ae5a4SJacob Faibussowitsch {
322b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
323b07a2398SLisandro Dalcin 
324b07a2398SLisandro Dalcin   PetscFunctionBegin;
3259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X0));
3269566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Xa));
3279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X1));
3289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->V0));
3299566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Va));
3309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->V1));
3319566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_sol_prev));
3329566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_lte_work));
3333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
334b07a2398SLisandro Dalcin }
335b07a2398SLisandro Dalcin 
336d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Alpha(TS ts)
337d71ae5a4SJacob Faibussowitsch {
338b07a2398SLisandro Dalcin   PetscFunctionBegin;
3399566063dSJacob Faibussowitsch   PetscCall(TSReset_Alpha(ts));
3409566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
341b07a2398SLisandro Dalcin 
3429566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetRadius_C", NULL));
3439566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetParams_C", NULL));
3449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaGetParams_C", NULL));
3453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
346b07a2398SLisandro Dalcin }
347b07a2398SLisandro Dalcin 
348d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Alpha(TS ts)
349d71ae5a4SJacob Faibussowitsch {
350b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
3512ffb9264SLisandro Dalcin   PetscBool match;
352b07a2398SLisandro Dalcin 
353b07a2398SLisandro Dalcin   PetscFunctionBegin;
3548ec9177eSStefano Zampini   if (!th->X0) PetscCall(VecDuplicate(ts->vec_sol, &th->X0));
3559566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->Xa));
3569566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->X1));
3579566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->V0));
3589566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->Va));
3599566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->V1));
3601566a47fSLisandro Dalcin 
3619566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &ts->adapt));
3629566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt));
3639566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match));
3642ffb9264SLisandro Dalcin   if (!match) {
3659566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev));
3669566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work));
367b07a2398SLisandro Dalcin   }
3681566a47fSLisandro Dalcin 
3699566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &ts->snes));
3703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
371b07a2398SLisandro Dalcin }
372b07a2398SLisandro Dalcin 
373d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_Alpha(TS ts, PetscOptionItems *PetscOptionsObject)
374d71ae5a4SJacob Faibussowitsch {
375b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
376b07a2398SLisandro Dalcin 
377b07a2398SLisandro Dalcin   PetscFunctionBegin;
378d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Generalized-Alpha ODE solver options");
379b07a2398SLisandro Dalcin   {
380b07a2398SLisandro Dalcin     PetscBool flg;
381b07a2398SLisandro Dalcin     PetscReal radius = 1;
3829566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_radius", "Spectral radius (high-frequency dissipation)", "TSAlphaSetRadius", radius, &radius, &flg));
3839566063dSJacob Faibussowitsch     if (flg) PetscCall(TSAlphaSetRadius(ts, radius));
3849566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_alpha_m", "Algorithmic parameter alpha_m", "TSAlphaSetParams", th->Alpha_m, &th->Alpha_m, NULL));
3859566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_alpha_f", "Algorithmic parameter alpha_f", "TSAlphaSetParams", th->Alpha_f, &th->Alpha_f, NULL));
3869566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_gamma", "Algorithmic parameter gamma", "TSAlphaSetParams", th->Gamma, &th->Gamma, NULL));
3879566063dSJacob Faibussowitsch     PetscCall(TSAlphaSetParams(ts, th->Alpha_m, th->Alpha_f, th->Gamma));
388b07a2398SLisandro Dalcin   }
389d0609cedSBarry Smith   PetscOptionsHeadEnd();
3903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
391b07a2398SLisandro Dalcin }
392b07a2398SLisandro Dalcin 
393d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Alpha(TS ts, PetscViewer viewer)
394d71ae5a4SJacob Faibussowitsch {
395b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
3969c334d8fSLisandro Dalcin   PetscBool iascii;
397b07a2398SLisandro Dalcin 
398b07a2398SLisandro Dalcin   PetscFunctionBegin;
3999566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
40048a46eb9SPierre 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));
4013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
402b07a2398SLisandro Dalcin }
403b07a2398SLisandro Dalcin 
404d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlphaSetRadius_Alpha(TS ts, PetscReal radius)
405d71ae5a4SJacob Faibussowitsch {
406b07a2398SLisandro Dalcin   PetscReal alpha_m, alpha_f, gamma;
407b07a2398SLisandro Dalcin 
408b07a2398SLisandro Dalcin   PetscFunctionBegin;
409cad9d221SBarry Smith   PetscCheck(radius >= 0 && radius <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius);
410b07a2398SLisandro Dalcin   alpha_m = (PetscReal)0.5 * (3 - radius) / (1 + radius);
411b07a2398SLisandro Dalcin   alpha_f = 1 / (1 + radius);
412b07a2398SLisandro Dalcin   gamma   = (PetscReal)0.5 + alpha_m - alpha_f;
4139566063dSJacob Faibussowitsch   PetscCall(TSAlphaSetParams(ts, alpha_m, alpha_f, gamma));
4143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
415b07a2398SLisandro Dalcin }
416b07a2398SLisandro Dalcin 
417d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlphaSetParams_Alpha(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma)
418d71ae5a4SJacob Faibussowitsch {
419b07a2398SLisandro Dalcin   TS_Alpha *th  = (TS_Alpha *)ts->data;
420b07a2398SLisandro Dalcin   PetscReal tol = 100 * PETSC_MACHINE_EPSILON;
421b07a2398SLisandro Dalcin   PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma;
422b07a2398SLisandro Dalcin 
423b07a2398SLisandro Dalcin   PetscFunctionBegin;
424b07a2398SLisandro Dalcin   th->Alpha_m = alpha_m;
425b07a2398SLisandro Dalcin   th->Alpha_f = alpha_f;
426b07a2398SLisandro Dalcin   th->Gamma   = gamma;
427b07a2398SLisandro Dalcin   th->order   = (PetscAbsReal(res) < tol) ? 2 : 1;
4283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
429b07a2398SLisandro Dalcin }
430b07a2398SLisandro Dalcin 
431d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlphaGetParams_Alpha(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma)
432d71ae5a4SJacob Faibussowitsch {
433b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
434b07a2398SLisandro Dalcin 
435b07a2398SLisandro Dalcin   PetscFunctionBegin;
436b07a2398SLisandro Dalcin   if (alpha_m) *alpha_m = th->Alpha_m;
437b07a2398SLisandro Dalcin   if (alpha_f) *alpha_f = th->Alpha_f;
438b07a2398SLisandro Dalcin   if (gamma) *gamma = th->Gamma;
4393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
440b07a2398SLisandro Dalcin }
441b07a2398SLisandro Dalcin 
442b07a2398SLisandro Dalcin /*MC
44314d0ab18SJacob Faibussowitsch   TSALPHA - ODE/DAE solver using the implicit Generalized-Alpha method for first-order systems
444b07a2398SLisandro Dalcin 
445b07a2398SLisandro Dalcin   Level: beginner
446b07a2398SLisandro Dalcin 
447b07a2398SLisandro Dalcin   References:
448606c0280SSatish Balay + * - K.E. Jansen, C.H. Whiting, G.M. Hulber, "A generalized-alpha
449b07a2398SLisandro Dalcin   method for integrating the filtered Navier-Stokes equations with a
450b07a2398SLisandro Dalcin   stabilized finite element method", Computer Methods in Applied
451b07a2398SLisandro Dalcin   Mechanics and Engineering, 190, 305-319, 2000.
452b07a2398SLisandro Dalcin   DOI: 10.1016/S0045-7825(00)00203-6.
453606c0280SSatish Balay - * -  J. Chung, G.M.Hubert. "A Time Integration Algorithm for Structural
454b07a2398SLisandro Dalcin   Dynamics with Improved Numerical Dissipation: The Generalized-alpha
455b07a2398SLisandro Dalcin   Method" ASME Journal of Applied Mechanics, 60, 371:375, 1993.
456b07a2398SLisandro Dalcin 
4571cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSCreate()`, `TSSetType()`, `TSAlphaSetRadius()`, `TSAlphaSetParams()`
458b07a2398SLisandro Dalcin M*/
459d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Alpha(TS ts)
460d71ae5a4SJacob Faibussowitsch {
461b07a2398SLisandro Dalcin   TS_Alpha *th;
462b07a2398SLisandro Dalcin 
463b07a2398SLisandro Dalcin   PetscFunctionBegin;
464b07a2398SLisandro Dalcin   ts->ops->reset          = TSReset_Alpha;
465b07a2398SLisandro Dalcin   ts->ops->destroy        = TSDestroy_Alpha;
466b07a2398SLisandro Dalcin   ts->ops->view           = TSView_Alpha;
467b07a2398SLisandro Dalcin   ts->ops->setup          = TSSetUp_Alpha;
468b07a2398SLisandro Dalcin   ts->ops->setfromoptions = TSSetFromOptions_Alpha;
469b07a2398SLisandro Dalcin   ts->ops->step           = TSStep_Alpha;
4709808bdc1SLisandro Dalcin   ts->ops->evaluatewlte   = TSEvaluateWLTE_Alpha;
471b07a2398SLisandro Dalcin   ts->ops->rollback       = TSRollBack_Alpha;
472b07a2398SLisandro Dalcin   ts->ops->interpolate    = TSInterpolate_Alpha;
4738ec9177eSStefano Zampini   ts->ops->resizeregister = TSResizeRegister_Alpha;
474b07a2398SLisandro Dalcin   ts->ops->snesfunction   = SNESTSFormFunction_Alpha;
475b07a2398SLisandro Dalcin   ts->ops->snesjacobian   = SNESTSFormJacobian_Alpha;
4762ffb9264SLisandro Dalcin   ts->default_adapt_type  = TSADAPTNONE;
477b07a2398SLisandro Dalcin 
478efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
479efd4aadfSBarry Smith 
4804dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&th));
481b07a2398SLisandro Dalcin   ts->data = (void *)th;
482b07a2398SLisandro Dalcin 
483b07a2398SLisandro Dalcin   th->Alpha_m = 0.5;
484b07a2398SLisandro Dalcin   th->Alpha_f = 0.5;
485b07a2398SLisandro Dalcin   th->Gamma   = 0.5;
486b07a2398SLisandro Dalcin   th->order   = 2;
487b07a2398SLisandro Dalcin 
4889566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetRadius_C", TSAlphaSetRadius_Alpha));
4899566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetParams_C", TSAlphaSetParams_Alpha));
4909566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaGetParams_C", TSAlphaGetParams_Alpha));
4913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
492b07a2398SLisandro Dalcin }
493b07a2398SLisandro Dalcin 
494b07a2398SLisandro Dalcin /*@
495bcf0153eSBarry Smith   TSAlphaSetRadius - sets the desired spectral radius of the method for `TSALPHA`
496b07a2398SLisandro Dalcin   (i.e. high-frequency numerical damping)
497b07a2398SLisandro Dalcin 
498c3339decSBarry Smith   Logically Collective
499b07a2398SLisandro Dalcin 
500d8d19677SJose E. Roman   Input Parameters:
501b07a2398SLisandro Dalcin + ts     - timestepping context
502b07a2398SLisandro Dalcin - radius - the desired spectral radius
503b07a2398SLisandro Dalcin 
504bcf0153eSBarry Smith   Options Database Key:
50567b8a455SSatish Balay . -ts_alpha_radius <radius> - set alpha radius
506b07a2398SLisandro Dalcin 
507b07a2398SLisandro Dalcin   Level: intermediate
508b07a2398SLisandro Dalcin 
50914d0ab18SJacob Faibussowitsch   Notes:
51014d0ab18SJacob Faibussowitsch   The algorithmic parameters $\alpha_m$ and $\alpha_f$ of the generalized-$\alpha$ method can
51114d0ab18SJacob Faibussowitsch   be computed in terms of a specified spectral radius $\rho$ in [0, 1] for infinite time step
51214d0ab18SJacob Faibussowitsch   in order to control high-frequency numerical damping\:
513*562efe2eSBarry Smith 
51414d0ab18SJacob Faibussowitsch   $$
515*562efe2eSBarry Smith   \begin{align*}
516*562efe2eSBarry Smith   \alpha_m = 0.5*(3-\rho)/(1+\rho) \\
51714d0ab18SJacob Faibussowitsch   \alpha_f = 1/(1+\rho)
518*562efe2eSBarry Smith   \end{align*}
51914d0ab18SJacob Faibussowitsch   $$
52014d0ab18SJacob Faibussowitsch 
5211cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA`, `TSAlphaSetParams()`, `TSAlphaGetParams()`
522b07a2398SLisandro Dalcin @*/
523d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlphaSetRadius(TS ts, PetscReal radius)
524d71ae5a4SJacob Faibussowitsch {
525b07a2398SLisandro Dalcin   PetscFunctionBegin;
526b07a2398SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
527b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, radius, 2);
528cad9d221SBarry Smith   PetscCheck(radius >= 0 && radius <= 1, ((PetscObject)ts)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius);
529cac4c232SBarry Smith   PetscTryMethod(ts, "TSAlphaSetRadius_C", (TS, PetscReal), (ts, radius));
5303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
531b07a2398SLisandro Dalcin }
532b07a2398SLisandro Dalcin 
533b07a2398SLisandro Dalcin /*@
534bcf0153eSBarry Smith   TSAlphaSetParams - sets the algorithmic parameters for `TSALPHA`
535b07a2398SLisandro Dalcin 
536c3339decSBarry Smith   Logically Collective
537b07a2398SLisandro Dalcin 
538d8d19677SJose E. Roman   Input Parameters:
539b07a2398SLisandro Dalcin + ts      - timestepping context
5402fe279fdSBarry Smith . alpha_m - algorithmic parameter
5412fe279fdSBarry Smith . alpha_f - algorithmic parameter
5422fe279fdSBarry Smith - gamma   - algorithmic parameter
543b07a2398SLisandro Dalcin 
544bcf0153eSBarry Smith   Options Database Keys:
54567b8a455SSatish Balay + -ts_alpha_alpha_m <alpha_m> - set alpha_m
54667b8a455SSatish Balay . -ts_alpha_alpha_f <alpha_f> - set alpha_f
54767b8a455SSatish Balay - -ts_alpha_gamma   <gamma>   - set gamma
548b07a2398SLisandro Dalcin 
549bcf0153eSBarry Smith   Level: advanced
550bcf0153eSBarry Smith 
551b07a2398SLisandro Dalcin   Note:
552*562efe2eSBarry Smith   Second-order accuracy can be obtained so long as\:  $\gamma = 0.5 + \alpha_m - \alpha_f$
55314d0ab18SJacob Faibussowitsch 
554*562efe2eSBarry Smith   Unconditional stability requires\: $\alpha_m >= \alpha_f >= 0.5$
55514d0ab18SJacob Faibussowitsch 
556*562efe2eSBarry Smith   Backward Euler method is recovered with\: $\alpha_m = \alpha_f = \gamma = 1$
55714d0ab18SJacob Faibussowitsch 
55814d0ab18SJacob Faibussowitsch   Use of this function is normally only required to hack `TSALPHA` to use a modified
55914d0ab18SJacob Faibussowitsch   integration scheme. Users should call `TSAlphaSetRadius()` to set the desired spectral radius
56014d0ab18SJacob Faibussowitsch   of the methods (i.e. high-frequency damping) in order so select optimal values for these
56114d0ab18SJacob Faibussowitsch   parameters.
562b07a2398SLisandro Dalcin 
5631cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA`, `TSAlphaSetRadius()`, `TSAlphaGetParams()`
564b07a2398SLisandro Dalcin @*/
565d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlphaSetParams(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma)
566d71ae5a4SJacob Faibussowitsch {
567b07a2398SLisandro Dalcin   PetscFunctionBegin;
568b07a2398SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
569b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, alpha_m, 2);
570b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, alpha_f, 3);
571b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, gamma, 4);
572cac4c232SBarry Smith   PetscTryMethod(ts, "TSAlphaSetParams_C", (TS, PetscReal, PetscReal, PetscReal), (ts, alpha_m, alpha_f, gamma));
5733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
574b07a2398SLisandro Dalcin }
575b07a2398SLisandro Dalcin 
576b07a2398SLisandro Dalcin /*@
577bcf0153eSBarry Smith   TSAlphaGetParams - gets the algorithmic parameters for `TSALPHA`
578b07a2398SLisandro Dalcin 
579b07a2398SLisandro Dalcin   Not Collective
580b07a2398SLisandro Dalcin 
581b07a2398SLisandro Dalcin   Input Parameter:
582b07a2398SLisandro Dalcin . ts - timestepping context
583b07a2398SLisandro Dalcin 
584b07a2398SLisandro Dalcin   Output Parameters:
5852fe279fdSBarry Smith + alpha_m - algorithmic parameter
5862fe279fdSBarry Smith . alpha_f - algorithmic parameter
5872fe279fdSBarry Smith - gamma   - algorithmic parameter
588b07a2398SLisandro Dalcin 
589bcf0153eSBarry Smith   Level: advanced
590bcf0153eSBarry Smith 
591b07a2398SLisandro Dalcin   Note:
59214d0ab18SJacob Faibussowitsch   Use of this function is normally only required to hack `TSALPHA` to use a modified
59314d0ab18SJacob Faibussowitsch   integration scheme. Users should call `TSAlphaSetRadius()` to set the high-frequency damping
59414d0ab18SJacob Faibussowitsch   (i.e. spectral radius of the method) in order so select optimal values for these parameters.
595b07a2398SLisandro Dalcin 
5961cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA`, `TSAlphaSetRadius()`, `TSAlphaSetParams()`
597b07a2398SLisandro Dalcin @*/
598d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlphaGetParams(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma)
599d71ae5a4SJacob Faibussowitsch {
600b07a2398SLisandro Dalcin   PetscFunctionBegin;
601b07a2398SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
6024f572ea9SToby Isaac   if (alpha_m) PetscAssertPointer(alpha_m, 2);
6034f572ea9SToby Isaac   if (alpha_f) PetscAssertPointer(alpha_f, 3);
6044f572ea9SToby Isaac   if (gamma) PetscAssertPointer(gamma, 4);
605cac4c232SBarry Smith   PetscUseMethod(ts, "TSAlphaGetParams_C", (TS, PetscReal *, PetscReal *, PetscReal *), (ts, alpha_m, alpha_f, gamma));
6063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
607b07a2398SLisandro Dalcin }
608