xref: /petsc/src/ts/impls/implicit/alpha/alpha1.c (revision ecc878985a1a123bb75746e4d26307f3bec57277)
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;
42*ecc87898SStefano Zampini   if (reg) {
43*ecc87898SStefano Zampini     PetscCall(TSResizeRegisterVec(ts, "ts:theta:sol_prev", th->vec_sol_prev));
44*ecc87898SStefano Zampini     PetscCall(TSResizeRegisterVec(ts, "ts:theta:X0", th->X0));
45*ecc87898SStefano Zampini   } else {
46*ecc87898SStefano Zampini     PetscCall(TSResizeRetrieveVec(ts, "ts:theta:sol_prev", &th->vec_sol_prev));
47*ecc87898SStefano Zampini     PetscCall(PetscObjectReference((PetscObject)th->vec_sol_prev));
48*ecc87898SStefano 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;
277b07a2398SLisandro Dalcin 
278b07a2398SLisandro Dalcin   PetscFunctionBegin;
2799566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol, X));
2809566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, th->Gamma * dt, th->V1));
2819566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, (1 - th->Gamma) * dt, th->V0));
2823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
283b07a2398SLisandro Dalcin }
284b07a2398SLisandro Dalcin 
285d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes, Vec X, Vec F, TS ts)
286d71ae5a4SJacob Faibussowitsch {
287b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
288b07a2398SLisandro Dalcin   PetscReal ta = th->stage_time;
289b07a2398SLisandro Dalcin   Vec       Xa = th->Xa, Va = th->Va;
290b07a2398SLisandro Dalcin 
291b07a2398SLisandro Dalcin   PetscFunctionBegin;
2929566063dSJacob Faibussowitsch   PetscCall(TSAlpha_StageVecs(ts, X));
293b07a2398SLisandro Dalcin   /* F = Function(ta,Xa,Va) */
2949566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, ta, Xa, Va, F, PETSC_FALSE));
2959566063dSJacob Faibussowitsch   PetscCall(VecScale(F, th->scale_F));
2963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
297b07a2398SLisandro Dalcin }
298b07a2398SLisandro Dalcin 
299d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes, PETSC_UNUSED Vec X, Mat J, Mat P, TS ts)
300d71ae5a4SJacob Faibussowitsch {
301b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
302b07a2398SLisandro Dalcin   PetscReal ta = th->stage_time;
303b07a2398SLisandro Dalcin   Vec       Xa = th->Xa, Va = th->Va;
304b07a2398SLisandro Dalcin   PetscReal dVdX = th->shift_V;
305b07a2398SLisandro Dalcin 
306b07a2398SLisandro Dalcin   PetscFunctionBegin;
307b07a2398SLisandro Dalcin   /* J,P = Jacobian(ta,Xa,Va) */
3089566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, ta, Xa, Va, dVdX, J, P, PETSC_FALSE));
3093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
310b07a2398SLisandro Dalcin }
311b07a2398SLisandro Dalcin 
312d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Alpha(TS ts)
313d71ae5a4SJacob Faibussowitsch {
314b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
315b07a2398SLisandro Dalcin 
316b07a2398SLisandro Dalcin   PetscFunctionBegin;
3179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X0));
3189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Xa));
3199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X1));
3209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->V0));
3219566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Va));
3229566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->V1));
3239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_sol_prev));
3249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_lte_work));
3253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
326b07a2398SLisandro Dalcin }
327b07a2398SLisandro Dalcin 
328d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Alpha(TS ts)
329d71ae5a4SJacob Faibussowitsch {
330b07a2398SLisandro Dalcin   PetscFunctionBegin;
3319566063dSJacob Faibussowitsch   PetscCall(TSReset_Alpha(ts));
3329566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
333b07a2398SLisandro Dalcin 
3349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetRadius_C", NULL));
3359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetParams_C", NULL));
3369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaGetParams_C", NULL));
3373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
338b07a2398SLisandro Dalcin }
339b07a2398SLisandro Dalcin 
340d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Alpha(TS ts)
341d71ae5a4SJacob Faibussowitsch {
342b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
3432ffb9264SLisandro Dalcin   PetscBool match;
344b07a2398SLisandro Dalcin 
345b07a2398SLisandro Dalcin   PetscFunctionBegin;
3468ec9177eSStefano Zampini   if (!th->X0) PetscCall(VecDuplicate(ts->vec_sol, &th->X0));
3479566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->Xa));
3489566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->X1));
3499566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->V0));
3509566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->Va));
3519566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->V1));
3521566a47fSLisandro Dalcin 
3539566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &ts->adapt));
3549566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt));
3559566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match));
3562ffb9264SLisandro Dalcin   if (!match) {
357*ecc87898SStefano Zampini     if (!th->vec_sol_prev) PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev));
358*ecc87898SStefano Zampini     if (!th->vec_lte_work) PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work));
359b07a2398SLisandro Dalcin   }
3601566a47fSLisandro Dalcin 
3619566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &ts->snes));
3623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
363b07a2398SLisandro Dalcin }
364b07a2398SLisandro Dalcin 
365d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_Alpha(TS ts, PetscOptionItems *PetscOptionsObject)
366d71ae5a4SJacob Faibussowitsch {
367b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
368b07a2398SLisandro Dalcin 
369b07a2398SLisandro Dalcin   PetscFunctionBegin;
370d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Generalized-Alpha ODE solver options");
371b07a2398SLisandro Dalcin   {
372b07a2398SLisandro Dalcin     PetscBool flg;
373b07a2398SLisandro Dalcin     PetscReal radius = 1;
3749566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_radius", "Spectral radius (high-frequency dissipation)", "TSAlphaSetRadius", radius, &radius, &flg));
3759566063dSJacob Faibussowitsch     if (flg) PetscCall(TSAlphaSetRadius(ts, radius));
3769566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_alpha_m", "Algorithmic parameter alpha_m", "TSAlphaSetParams", th->Alpha_m, &th->Alpha_m, NULL));
3779566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_alpha_f", "Algorithmic parameter alpha_f", "TSAlphaSetParams", th->Alpha_f, &th->Alpha_f, NULL));
3789566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_gamma", "Algorithmic parameter gamma", "TSAlphaSetParams", th->Gamma, &th->Gamma, NULL));
3799566063dSJacob Faibussowitsch     PetscCall(TSAlphaSetParams(ts, th->Alpha_m, th->Alpha_f, th->Gamma));
380b07a2398SLisandro Dalcin   }
381d0609cedSBarry Smith   PetscOptionsHeadEnd();
3823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
383b07a2398SLisandro Dalcin }
384b07a2398SLisandro Dalcin 
385d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Alpha(TS ts, PetscViewer viewer)
386d71ae5a4SJacob Faibussowitsch {
387b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
3889c334d8fSLisandro Dalcin   PetscBool iascii;
389b07a2398SLisandro Dalcin 
390b07a2398SLisandro Dalcin   PetscFunctionBegin;
3919566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
39248a46eb9SPierre 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));
3933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
394b07a2398SLisandro Dalcin }
395b07a2398SLisandro Dalcin 
396d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlphaSetRadius_Alpha(TS ts, PetscReal radius)
397d71ae5a4SJacob Faibussowitsch {
398b07a2398SLisandro Dalcin   PetscReal alpha_m, alpha_f, gamma;
399b07a2398SLisandro Dalcin 
400b07a2398SLisandro Dalcin   PetscFunctionBegin;
401cad9d221SBarry Smith   PetscCheck(radius >= 0 && radius <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius);
402b07a2398SLisandro Dalcin   alpha_m = (PetscReal)0.5 * (3 - radius) / (1 + radius);
403b07a2398SLisandro Dalcin   alpha_f = 1 / (1 + radius);
404b07a2398SLisandro Dalcin   gamma   = (PetscReal)0.5 + alpha_m - alpha_f;
4059566063dSJacob Faibussowitsch   PetscCall(TSAlphaSetParams(ts, alpha_m, alpha_f, gamma));
4063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
407b07a2398SLisandro Dalcin }
408b07a2398SLisandro Dalcin 
409d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlphaSetParams_Alpha(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma)
410d71ae5a4SJacob Faibussowitsch {
411b07a2398SLisandro Dalcin   TS_Alpha *th  = (TS_Alpha *)ts->data;
412b07a2398SLisandro Dalcin   PetscReal tol = 100 * PETSC_MACHINE_EPSILON;
413b07a2398SLisandro Dalcin   PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma;
414b07a2398SLisandro Dalcin 
415b07a2398SLisandro Dalcin   PetscFunctionBegin;
416b07a2398SLisandro Dalcin   th->Alpha_m = alpha_m;
417b07a2398SLisandro Dalcin   th->Alpha_f = alpha_f;
418b07a2398SLisandro Dalcin   th->Gamma   = gamma;
419b07a2398SLisandro Dalcin   th->order   = (PetscAbsReal(res) < tol) ? 2 : 1;
4203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
421b07a2398SLisandro Dalcin }
422b07a2398SLisandro Dalcin 
423d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlphaGetParams_Alpha(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma)
424d71ae5a4SJacob Faibussowitsch {
425b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
426b07a2398SLisandro Dalcin 
427b07a2398SLisandro Dalcin   PetscFunctionBegin;
428b07a2398SLisandro Dalcin   if (alpha_m) *alpha_m = th->Alpha_m;
429b07a2398SLisandro Dalcin   if (alpha_f) *alpha_f = th->Alpha_f;
430b07a2398SLisandro Dalcin   if (gamma) *gamma = th->Gamma;
4313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
432b07a2398SLisandro Dalcin }
433b07a2398SLisandro Dalcin 
434b07a2398SLisandro Dalcin /*MC
4351d27aa22SBarry Smith   TSALPHA - ODE/DAE solver using the implicit Generalized-Alpha method {cite}`jansen_2000` {cite}`chung1993` for first-order systems
436b07a2398SLisandro Dalcin 
437b07a2398SLisandro Dalcin   Level: beginner
438b07a2398SLisandro Dalcin 
4391cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSCreate()`, `TSSetType()`, `TSAlphaSetRadius()`, `TSAlphaSetParams()`
440b07a2398SLisandro Dalcin M*/
441d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Alpha(TS ts)
442d71ae5a4SJacob Faibussowitsch {
443b07a2398SLisandro Dalcin   TS_Alpha *th;
444b07a2398SLisandro Dalcin 
445b07a2398SLisandro Dalcin   PetscFunctionBegin;
446b07a2398SLisandro Dalcin   ts->ops->reset          = TSReset_Alpha;
447b07a2398SLisandro Dalcin   ts->ops->destroy        = TSDestroy_Alpha;
448b07a2398SLisandro Dalcin   ts->ops->view           = TSView_Alpha;
449b07a2398SLisandro Dalcin   ts->ops->setup          = TSSetUp_Alpha;
450b07a2398SLisandro Dalcin   ts->ops->setfromoptions = TSSetFromOptions_Alpha;
451b07a2398SLisandro Dalcin   ts->ops->step           = TSStep_Alpha;
4529808bdc1SLisandro Dalcin   ts->ops->evaluatewlte   = TSEvaluateWLTE_Alpha;
453b07a2398SLisandro Dalcin   ts->ops->interpolate    = TSInterpolate_Alpha;
4548ec9177eSStefano Zampini   ts->ops->resizeregister = TSResizeRegister_Alpha;
455b07a2398SLisandro Dalcin   ts->ops->snesfunction   = SNESTSFormFunction_Alpha;
456b07a2398SLisandro Dalcin   ts->ops->snesjacobian   = SNESTSFormJacobian_Alpha;
4572ffb9264SLisandro Dalcin   ts->default_adapt_type  = TSADAPTNONE;
458b07a2398SLisandro Dalcin 
459efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
460efd4aadfSBarry Smith 
4614dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&th));
462b07a2398SLisandro Dalcin   ts->data = (void *)th;
463b07a2398SLisandro Dalcin 
464b07a2398SLisandro Dalcin   th->Alpha_m = 0.5;
465b07a2398SLisandro Dalcin   th->Alpha_f = 0.5;
466b07a2398SLisandro Dalcin   th->Gamma   = 0.5;
467b07a2398SLisandro Dalcin   th->order   = 2;
468b07a2398SLisandro Dalcin 
4699566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetRadius_C", TSAlphaSetRadius_Alpha));
4709566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetParams_C", TSAlphaSetParams_Alpha));
4719566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaGetParams_C", TSAlphaGetParams_Alpha));
4723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
473b07a2398SLisandro Dalcin }
474b07a2398SLisandro Dalcin 
475b07a2398SLisandro Dalcin /*@
476bcf0153eSBarry Smith   TSAlphaSetRadius - sets the desired spectral radius of the method for `TSALPHA`
477b07a2398SLisandro Dalcin   (i.e. high-frequency numerical damping)
478b07a2398SLisandro Dalcin 
479c3339decSBarry Smith   Logically Collective
480b07a2398SLisandro Dalcin 
481d8d19677SJose E. Roman   Input Parameters:
482b07a2398SLisandro Dalcin + ts     - timestepping context
483b07a2398SLisandro Dalcin - radius - the desired spectral radius
484b07a2398SLisandro Dalcin 
485bcf0153eSBarry Smith   Options Database Key:
48667b8a455SSatish Balay . -ts_alpha_radius <radius> - set alpha radius
487b07a2398SLisandro Dalcin 
488b07a2398SLisandro Dalcin   Level: intermediate
489b07a2398SLisandro Dalcin 
49014d0ab18SJacob Faibussowitsch   Notes:
49114d0ab18SJacob Faibussowitsch   The algorithmic parameters $\alpha_m$ and $\alpha_f$ of the generalized-$\alpha$ method can
49214d0ab18SJacob Faibussowitsch   be computed in terms of a specified spectral radius $\rho$ in [0, 1] for infinite time step
49314d0ab18SJacob Faibussowitsch   in order to control high-frequency numerical damping\:
494562efe2eSBarry Smith 
49514d0ab18SJacob Faibussowitsch   $$
496562efe2eSBarry Smith   \begin{align*}
497562efe2eSBarry Smith   \alpha_m = 0.5*(3-\rho)/(1+\rho) \\
49814d0ab18SJacob Faibussowitsch   \alpha_f = 1/(1+\rho)
499562efe2eSBarry Smith   \end{align*}
50014d0ab18SJacob Faibussowitsch   $$
50114d0ab18SJacob Faibussowitsch 
5021cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA`, `TSAlphaSetParams()`, `TSAlphaGetParams()`
503b07a2398SLisandro Dalcin @*/
504d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlphaSetRadius(TS ts, PetscReal radius)
505d71ae5a4SJacob Faibussowitsch {
506b07a2398SLisandro Dalcin   PetscFunctionBegin;
507b07a2398SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
508b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, radius, 2);
509cad9d221SBarry Smith   PetscCheck(radius >= 0 && radius <= 1, ((PetscObject)ts)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius);
510cac4c232SBarry Smith   PetscTryMethod(ts, "TSAlphaSetRadius_C", (TS, PetscReal), (ts, radius));
5113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
512b07a2398SLisandro Dalcin }
513b07a2398SLisandro Dalcin 
514b07a2398SLisandro Dalcin /*@
515bcf0153eSBarry Smith   TSAlphaSetParams - sets the algorithmic parameters for `TSALPHA`
516b07a2398SLisandro Dalcin 
517c3339decSBarry Smith   Logically Collective
518b07a2398SLisandro Dalcin 
519d8d19677SJose E. Roman   Input Parameters:
520b07a2398SLisandro Dalcin + ts      - timestepping context
5212fe279fdSBarry Smith . alpha_m - algorithmic parameter
5222fe279fdSBarry Smith . alpha_f - algorithmic parameter
5232fe279fdSBarry Smith - gamma   - algorithmic parameter
524b07a2398SLisandro Dalcin 
525bcf0153eSBarry Smith   Options Database Keys:
52667b8a455SSatish Balay + -ts_alpha_alpha_m <alpha_m> - set alpha_m
52767b8a455SSatish Balay . -ts_alpha_alpha_f <alpha_f> - set alpha_f
52867b8a455SSatish Balay - -ts_alpha_gamma   <gamma>   - set gamma
529b07a2398SLisandro Dalcin 
530bcf0153eSBarry Smith   Level: advanced
531bcf0153eSBarry Smith 
532b07a2398SLisandro Dalcin   Note:
533562efe2eSBarry Smith   Second-order accuracy can be obtained so long as\:  $\gamma = 0.5 + \alpha_m - \alpha_f$
53414d0ab18SJacob Faibussowitsch 
535562efe2eSBarry Smith   Unconditional stability requires\: $\alpha_m >= \alpha_f >= 0.5$
53614d0ab18SJacob Faibussowitsch 
537562efe2eSBarry Smith   Backward Euler method is recovered with\: $\alpha_m = \alpha_f = \gamma = 1$
53814d0ab18SJacob Faibussowitsch 
53914d0ab18SJacob Faibussowitsch   Use of this function is normally only required to hack `TSALPHA` to use a modified
54014d0ab18SJacob Faibussowitsch   integration scheme. Users should call `TSAlphaSetRadius()` to set the desired spectral radius
54114d0ab18SJacob Faibussowitsch   of the methods (i.e. high-frequency damping) in order so select optimal values for these
54214d0ab18SJacob Faibussowitsch   parameters.
543b07a2398SLisandro Dalcin 
5441cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA`, `TSAlphaSetRadius()`, `TSAlphaGetParams()`
545b07a2398SLisandro Dalcin @*/
546d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlphaSetParams(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma)
547d71ae5a4SJacob Faibussowitsch {
548b07a2398SLisandro Dalcin   PetscFunctionBegin;
549b07a2398SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
550b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, alpha_m, 2);
551b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, alpha_f, 3);
552b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, gamma, 4);
553cac4c232SBarry Smith   PetscTryMethod(ts, "TSAlphaSetParams_C", (TS, PetscReal, PetscReal, PetscReal), (ts, alpha_m, alpha_f, gamma));
5543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
555b07a2398SLisandro Dalcin }
556b07a2398SLisandro Dalcin 
557b07a2398SLisandro Dalcin /*@
558bcf0153eSBarry Smith   TSAlphaGetParams - gets the algorithmic parameters for `TSALPHA`
559b07a2398SLisandro Dalcin 
560b07a2398SLisandro Dalcin   Not Collective
561b07a2398SLisandro Dalcin 
562b07a2398SLisandro Dalcin   Input Parameter:
563b07a2398SLisandro Dalcin . ts - timestepping context
564b07a2398SLisandro Dalcin 
565b07a2398SLisandro Dalcin   Output Parameters:
5662fe279fdSBarry Smith + alpha_m - algorithmic parameter
5672fe279fdSBarry Smith . alpha_f - algorithmic parameter
5682fe279fdSBarry Smith - gamma   - algorithmic parameter
569b07a2398SLisandro Dalcin 
570bcf0153eSBarry Smith   Level: advanced
571bcf0153eSBarry Smith 
572b07a2398SLisandro Dalcin   Note:
57314d0ab18SJacob Faibussowitsch   Use of this function is normally only required to hack `TSALPHA` to use a modified
57414d0ab18SJacob Faibussowitsch   integration scheme. Users should call `TSAlphaSetRadius()` to set the high-frequency damping
57514d0ab18SJacob Faibussowitsch   (i.e. spectral radius of the method) in order so select optimal values for these parameters.
576b07a2398SLisandro Dalcin 
5771cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA`, `TSAlphaSetRadius()`, `TSAlphaSetParams()`
578b07a2398SLisandro Dalcin @*/
579d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlphaGetParams(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma)
580d71ae5a4SJacob Faibussowitsch {
581b07a2398SLisandro Dalcin   PetscFunctionBegin;
582b07a2398SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5834f572ea9SToby Isaac   if (alpha_m) PetscAssertPointer(alpha_m, 2);
5844f572ea9SToby Isaac   if (alpha_f) PetscAssertPointer(alpha_f, 3);
5854f572ea9SToby Isaac   if (gamma) PetscAssertPointer(gamma, 4);
586cac4c232SBarry Smith   PetscUseMethod(ts, "TSAlphaGetParams_C", (TS, PetscReal *, PetscReal *, PetscReal *), (ts, alpha_m, alpha_f, gamma));
5873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
588b07a2398SLisandro Dalcin }
589