xref: /petsc/src/ts/impls/implicit/alpha/alpha1.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
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 
37d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_StageTime(TS ts)
38d71ae5a4SJacob Faibussowitsch {
39b07a2398SLisandro Dalcin   TS_Alpha *th      = (TS_Alpha *)ts->data;
40b07a2398SLisandro Dalcin   PetscReal t       = ts->ptime;
41b07a2398SLisandro Dalcin   PetscReal dt      = ts->time_step;
42b07a2398SLisandro Dalcin   PetscReal Alpha_m = th->Alpha_m;
43b07a2398SLisandro Dalcin   PetscReal Alpha_f = th->Alpha_f;
44b07a2398SLisandro Dalcin   PetscReal Gamma   = th->Gamma;
45b07a2398SLisandro Dalcin 
46b07a2398SLisandro Dalcin   PetscFunctionBegin;
47b07a2398SLisandro Dalcin   th->stage_time = t + Alpha_f * dt;
48b07a2398SLisandro Dalcin   th->shift_V    = Alpha_m / (Alpha_f * Gamma * dt);
49b07a2398SLisandro Dalcin   th->scale_F    = 1 / Alpha_f;
503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
51b07a2398SLisandro Dalcin }
52b07a2398SLisandro Dalcin 
53d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_StageVecs(TS ts, Vec X)
54d71ae5a4SJacob Faibussowitsch {
55b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
56b07a2398SLisandro Dalcin   Vec       X1 = X, V1 = th->V1;
57b07a2398SLisandro Dalcin   Vec       Xa = th->Xa, Va = th->Va;
58b07a2398SLisandro Dalcin   Vec       X0 = th->X0, V0 = th->V0;
59b07a2398SLisandro Dalcin   PetscReal dt      = ts->time_step;
60b07a2398SLisandro Dalcin   PetscReal Alpha_m = th->Alpha_m;
61b07a2398SLisandro Dalcin   PetscReal Alpha_f = th->Alpha_f;
62b07a2398SLisandro Dalcin   PetscReal Gamma   = th->Gamma;
63b07a2398SLisandro Dalcin 
64b07a2398SLisandro Dalcin   PetscFunctionBegin;
65b07a2398SLisandro Dalcin   /* V1 = 1/(Gamma*dT)*(X1-X0) + (1-1/Gamma)*V0 */
669566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(V1, -1.0, X0, X1));
679566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(V1, 1 - 1 / Gamma, 1 / (Gamma * dt), V0));
68b07a2398SLisandro Dalcin   /* Xa = X0 + Alpha_f*(X1-X0) */
699566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Xa, -1.0, X0, X1));
709566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Xa, Alpha_f, X0));
71b07a2398SLisandro Dalcin   /* Va = V0 + Alpha_m*(V1-V0) */
729566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Va, -1.0, V0, V1));
739566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Va, Alpha_m, V0));
743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
75b07a2398SLisandro Dalcin }
76b07a2398SLisandro Dalcin 
77d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_SNESSolve(TS ts, Vec b, Vec x)
78d71ae5a4SJacob Faibussowitsch {
79b07a2398SLisandro Dalcin   PetscInt nits, lits;
80b07a2398SLisandro Dalcin 
81b07a2398SLisandro Dalcin   PetscFunctionBegin;
829566063dSJacob Faibussowitsch   PetscCall(SNESSolve(ts->snes, b, x));
839566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(ts->snes, &nits));
849566063dSJacob Faibussowitsch   PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits));
859371c9d4SSatish Balay   ts->snes_its += nits;
869371c9d4SSatish Balay   ts->ksp_its += lits;
873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
88b07a2398SLisandro Dalcin }
89b07a2398SLisandro Dalcin 
90b07a2398SLisandro Dalcin /*
91b07a2398SLisandro Dalcin   Compute a consistent initial state for the generalized-alpha method.
92b07a2398SLisandro Dalcin   - Solve two successive backward Euler steps with halved time step.
93b07a2398SLisandro Dalcin   - Compute the initial time derivative using backward differences.
94b07a2398SLisandro Dalcin   - If using adaptivity, estimate the LTE of the initial step.
95b07a2398SLisandro Dalcin */
96d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_Restart(TS ts, PetscBool *initok)
97d71ae5a4SJacob Faibussowitsch {
98b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
99b07a2398SLisandro Dalcin   PetscReal time_step;
100b07a2398SLisandro Dalcin   PetscReal alpha_m, alpha_f, gamma;
101b07a2398SLisandro Dalcin   Vec       X0 = ts->vec_sol, X1, X2 = th->X1;
102b07a2398SLisandro Dalcin   PetscBool stageok;
103b07a2398SLisandro Dalcin 
104b07a2398SLisandro Dalcin   PetscFunctionBegin;
1059566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X0, &X1));
106b07a2398SLisandro Dalcin 
107b07a2398SLisandro Dalcin   /* Setup backward Euler with halved time step */
1089566063dSJacob Faibussowitsch   PetscCall(TSAlphaGetParams(ts, &alpha_m, &alpha_f, &gamma));
1099566063dSJacob Faibussowitsch   PetscCall(TSAlphaSetParams(ts, 1, 1, 1));
1109566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &time_step));
111b07a2398SLisandro Dalcin   ts->time_step = time_step / 2;
1129566063dSJacob Faibussowitsch   PetscCall(TSAlpha_StageTime(ts));
113b07a2398SLisandro Dalcin   th->stage_time = ts->ptime;
1149566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(th->V0));
115b07a2398SLisandro Dalcin 
116b07a2398SLisandro Dalcin   /* First BE step, (t0,X0) -> (t1,X1) */
117b07a2398SLisandro Dalcin   th->stage_time += ts->time_step;
1189566063dSJacob Faibussowitsch   PetscCall(VecCopy(X0, th->X0));
1199566063dSJacob Faibussowitsch   PetscCall(TSPreStage(ts, th->stage_time));
1209566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->X0, X1));
1219566063dSJacob Faibussowitsch   PetscCall(TSAlpha_SNESSolve(ts, NULL, X1));
1229566063dSJacob Faibussowitsch   PetscCall(TSPostStage(ts, th->stage_time, 0, &X1));
1239566063dSJacob Faibussowitsch   PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, X1, &stageok));
124b07a2398SLisandro Dalcin   if (!stageok) goto finally;
125b07a2398SLisandro Dalcin 
126b07a2398SLisandro Dalcin   /* Second BE step, (t1,X1) -> (t2,X2) */
127b07a2398SLisandro Dalcin   th->stage_time += ts->time_step;
1289566063dSJacob Faibussowitsch   PetscCall(VecCopy(X1, th->X0));
1299566063dSJacob Faibussowitsch   PetscCall(TSPreStage(ts, th->stage_time));
1309566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->X0, X2));
1319566063dSJacob Faibussowitsch   PetscCall(TSAlpha_SNESSolve(ts, NULL, X2));
1329566063dSJacob Faibussowitsch   PetscCall(TSPostStage(ts, th->stage_time, 0, &X2));
1339566063dSJacob Faibussowitsch   PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, X2, &stageok));
134b07a2398SLisandro Dalcin   if (!stageok) goto finally;
135b07a2398SLisandro Dalcin 
136b07a2398SLisandro Dalcin   /* Compute V0 ~ dX/dt at t0 with backward differences */
1379566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(th->V0));
1389566063dSJacob Faibussowitsch   PetscCall(VecAXPY(th->V0, -3 / ts->time_step, X0));
1399566063dSJacob Faibussowitsch   PetscCall(VecAXPY(th->V0, +4 / ts->time_step, X1));
1409566063dSJacob Faibussowitsch   PetscCall(VecAXPY(th->V0, -1 / ts->time_step, X2));
141b07a2398SLisandro Dalcin 
142b07a2398SLisandro Dalcin   /* Rough, lower-order estimate LTE of the initial step */
1432ffb9264SLisandro Dalcin   if (th->vec_lte_work) {
1449566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(th->vec_lte_work));
1459566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work, +2, X2));
1469566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work, -4, X1));
1479566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work, +2, X0));
148b07a2398SLisandro Dalcin   }
149b07a2398SLisandro Dalcin 
150b07a2398SLisandro Dalcin finally:
151b07a2398SLisandro Dalcin   /* Revert TSAlpha to the initial state (t0,X0) */
152b07a2398SLisandro Dalcin   if (initok) *initok = stageok;
1539566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, time_step));
1549566063dSJacob Faibussowitsch   PetscCall(TSAlphaSetParams(ts, alpha_m, alpha_f, gamma));
1559566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol, th->X0));
156b07a2398SLisandro Dalcin 
1579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X1));
1583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
159b07a2398SLisandro Dalcin }
160b07a2398SLisandro Dalcin 
161d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Alpha(TS ts)
162d71ae5a4SJacob Faibussowitsch {
163b07a2398SLisandro Dalcin   TS_Alpha *th         = (TS_Alpha *)ts->data;
164b07a2398SLisandro Dalcin   PetscInt  rejections = 0;
165b07a2398SLisandro Dalcin   PetscBool stageok, accept = PETSC_TRUE;
166b07a2398SLisandro Dalcin   PetscReal next_time_step = ts->time_step;
167b07a2398SLisandro Dalcin 
168b07a2398SLisandro Dalcin   PetscFunctionBegin;
1699566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(citation, &cited));
170b07a2398SLisandro Dalcin 
171b07a2398SLisandro Dalcin   if (!ts->steprollback) {
1729566063dSJacob Faibussowitsch     if (th->vec_sol_prev) PetscCall(VecCopy(th->X0, th->vec_sol_prev));
1739566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_sol, th->X0));
1749566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->V1, th->V0));
175b07a2398SLisandro Dalcin   }
176b07a2398SLisandro Dalcin 
1771566a47fSLisandro Dalcin   th->status = TS_STEP_INCOMPLETE;
178b07a2398SLisandro Dalcin   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
179fecfb714SLisandro Dalcin     if (ts->steprestart) {
1809566063dSJacob Faibussowitsch       PetscCall(TSAlpha_Restart(ts, &stageok));
181fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
182b07a2398SLisandro Dalcin     }
183b07a2398SLisandro Dalcin 
1849566063dSJacob Faibussowitsch     PetscCall(TSAlpha_StageTime(ts));
1859566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->X0, th->X1));
1869566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, th->stage_time));
1879566063dSJacob Faibussowitsch     PetscCall(TSAlpha_SNESSolve(ts, NULL, th->X1));
1889566063dSJacob Faibussowitsch     PetscCall(TSPostStage(ts, th->stage_time, 0, &th->Xa));
1899566063dSJacob Faibussowitsch     PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->Xa, &stageok));
190fecfb714SLisandro Dalcin     if (!stageok) goto reject_step;
191b07a2398SLisandro Dalcin 
1921566a47fSLisandro Dalcin     th->status = TS_STEP_PENDING;
1939566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->X1, ts->vec_sol));
1949566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
1951566a47fSLisandro Dalcin     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
196b07a2398SLisandro Dalcin     if (!accept) {
1979566063dSJacob Faibussowitsch       PetscCall(VecCopy(th->X0, ts->vec_sol));
198be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
199be5899b3SLisandro Dalcin       goto reject_step;
200b07a2398SLisandro Dalcin     }
201b07a2398SLisandro Dalcin 
202b07a2398SLisandro Dalcin     ts->ptime += ts->time_step;
203b07a2398SLisandro Dalcin     ts->time_step = next_time_step;
204b07a2398SLisandro Dalcin     break;
205b07a2398SLisandro Dalcin 
206b07a2398SLisandro Dalcin   reject_step:
2079371c9d4SSatish Balay     ts->reject++;
2089371c9d4SSatish Balay     accept = PETSC_FALSE;
209b07a2398SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
210b07a2398SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
21163a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
212b07a2398SLisandro Dalcin     }
213b07a2398SLisandro Dalcin   }
2143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
215b07a2398SLisandro Dalcin }
216b07a2398SLisandro Dalcin 
217d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateWLTE_Alpha(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte)
218d71ae5a4SJacob Faibussowitsch {
219b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
2209808bdc1SLisandro Dalcin   Vec       X  = th->X1;           /* X = solution */
2211566a47fSLisandro Dalcin   Vec       Y  = th->vec_lte_work; /* Y = X + LTE  */
2227453f775SEmil Constantinescu   PetscReal wltea, wlter;
223b07a2398SLisandro Dalcin 
224b07a2398SLisandro Dalcin   PetscFunctionBegin;
2259371c9d4SSatish Balay   if (!th->vec_sol_prev) {
2269371c9d4SSatish Balay     *wlte = -1;
2273ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2289371c9d4SSatish Balay   }
2299371c9d4SSatish Balay   if (!th->vec_lte_work) {
2309371c9d4SSatish Balay     *wlte = -1;
2313ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2329371c9d4SSatish Balay   }
233fecfb714SLisandro Dalcin   if (ts->steprestart) {
234fecfb714SLisandro Dalcin     /* th->vec_lte_work is set to the LTE in TSAlpha_Restart() */
2359566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Y, 1, X));
236b07a2398SLisandro Dalcin   } else {
237b07a2398SLisandro Dalcin     /* Compute LTE using backward differences with non-constant time step */
238be5899b3SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
239be5899b3SLisandro Dalcin     PetscReal   a = 1 + h_prev / h;
2409371c9d4SSatish Balay     PetscScalar scal[3];
2419371c9d4SSatish Balay     Vec         vecs[3];
2429371c9d4SSatish Balay     scal[0] = +1 / a;
2439371c9d4SSatish Balay     scal[1] = -1 / (a - 1);
2449371c9d4SSatish Balay     scal[2] = +1 / (a * (a - 1));
2459371c9d4SSatish Balay     vecs[0] = th->X1;
2469371c9d4SSatish Balay     vecs[1] = th->X0;
2479371c9d4SSatish Balay     vecs[2] = th->vec_sol_prev;
2489566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Y));
2499566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(Y, 3, scal, vecs));
250b07a2398SLisandro Dalcin   }
2519566063dSJacob Faibussowitsch   PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, wlte, &wltea, &wlter));
2529808bdc1SLisandro Dalcin   if (order) *order = 2;
2533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
254b07a2398SLisandro Dalcin }
255b07a2398SLisandro Dalcin 
256d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_Alpha(TS ts)
257d71ae5a4SJacob Faibussowitsch {
258b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
259b07a2398SLisandro Dalcin 
260b07a2398SLisandro Dalcin   PetscFunctionBegin;
2619566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->X0, ts->vec_sol));
2623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
263b07a2398SLisandro Dalcin }
264b07a2398SLisandro Dalcin 
265d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_Alpha(TS ts, PetscReal t, Vec X)
266d71ae5a4SJacob Faibussowitsch {
267b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
268b07a2398SLisandro Dalcin   PetscReal dt = t - ts->ptime;
269b07a2398SLisandro Dalcin 
270b07a2398SLisandro Dalcin   PetscFunctionBegin;
2719566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol, X));
2729566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, th->Gamma * dt, th->V1));
2739566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, (1 - th->Gamma) * dt, th->V0));
2743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
275b07a2398SLisandro Dalcin }
276b07a2398SLisandro Dalcin 
277d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes, Vec X, Vec F, TS ts)
278d71ae5a4SJacob Faibussowitsch {
279b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
280b07a2398SLisandro Dalcin   PetscReal ta = th->stage_time;
281b07a2398SLisandro Dalcin   Vec       Xa = th->Xa, Va = th->Va;
282b07a2398SLisandro Dalcin 
283b07a2398SLisandro Dalcin   PetscFunctionBegin;
2849566063dSJacob Faibussowitsch   PetscCall(TSAlpha_StageVecs(ts, X));
285b07a2398SLisandro Dalcin   /* F = Function(ta,Xa,Va) */
2869566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, ta, Xa, Va, F, PETSC_FALSE));
2879566063dSJacob Faibussowitsch   PetscCall(VecScale(F, th->scale_F));
2883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
289b07a2398SLisandro Dalcin }
290b07a2398SLisandro Dalcin 
291d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes, PETSC_UNUSED Vec X, Mat J, Mat P, TS ts)
292d71ae5a4SJacob Faibussowitsch {
293b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
294b07a2398SLisandro Dalcin   PetscReal ta = th->stage_time;
295b07a2398SLisandro Dalcin   Vec       Xa = th->Xa, Va = th->Va;
296b07a2398SLisandro Dalcin   PetscReal dVdX = th->shift_V;
297b07a2398SLisandro Dalcin 
298b07a2398SLisandro Dalcin   PetscFunctionBegin;
299b07a2398SLisandro Dalcin   /* J,P = Jacobian(ta,Xa,Va) */
3009566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, ta, Xa, Va, dVdX, J, P, PETSC_FALSE));
3013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
302b07a2398SLisandro Dalcin }
303b07a2398SLisandro Dalcin 
304d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Alpha(TS ts)
305d71ae5a4SJacob Faibussowitsch {
306b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
307b07a2398SLisandro Dalcin 
308b07a2398SLisandro Dalcin   PetscFunctionBegin;
3099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X0));
3109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Xa));
3119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X1));
3129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->V0));
3139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Va));
3149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->V1));
3159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_sol_prev));
3169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_lte_work));
3173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
318b07a2398SLisandro Dalcin }
319b07a2398SLisandro Dalcin 
320d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Alpha(TS ts)
321d71ae5a4SJacob Faibussowitsch {
322b07a2398SLisandro Dalcin   PetscFunctionBegin;
3239566063dSJacob Faibussowitsch   PetscCall(TSReset_Alpha(ts));
3249566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
325b07a2398SLisandro Dalcin 
3269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetRadius_C", NULL));
3279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetParams_C", NULL));
3289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaGetParams_C", NULL));
3293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
330b07a2398SLisandro Dalcin }
331b07a2398SLisandro Dalcin 
332d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Alpha(TS ts)
333d71ae5a4SJacob Faibussowitsch {
334b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
3352ffb9264SLisandro Dalcin   PetscBool match;
336b07a2398SLisandro Dalcin 
337b07a2398SLisandro Dalcin   PetscFunctionBegin;
3389566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->X0));
3399566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->Xa));
3409566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->X1));
3419566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->V0));
3429566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->Va));
3439566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->V1));
3441566a47fSLisandro Dalcin 
3459566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &ts->adapt));
3469566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt));
3479566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match));
3482ffb9264SLisandro Dalcin   if (!match) {
3499566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev));
3509566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work));
351b07a2398SLisandro Dalcin   }
3521566a47fSLisandro Dalcin 
3539566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &ts->snes));
3543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
355b07a2398SLisandro Dalcin }
356b07a2398SLisandro Dalcin 
357d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_Alpha(TS ts, PetscOptionItems *PetscOptionsObject)
358d71ae5a4SJacob Faibussowitsch {
359b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
360b07a2398SLisandro Dalcin 
361b07a2398SLisandro Dalcin   PetscFunctionBegin;
362d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Generalized-Alpha ODE solver options");
363b07a2398SLisandro Dalcin   {
364b07a2398SLisandro Dalcin     PetscBool flg;
365b07a2398SLisandro Dalcin     PetscReal radius = 1;
3669566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_radius", "Spectral radius (high-frequency dissipation)", "TSAlphaSetRadius", radius, &radius, &flg));
3679566063dSJacob Faibussowitsch     if (flg) PetscCall(TSAlphaSetRadius(ts, radius));
3689566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_alpha_m", "Algorithmic parameter alpha_m", "TSAlphaSetParams", th->Alpha_m, &th->Alpha_m, NULL));
3699566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_alpha_f", "Algorithmic parameter alpha_f", "TSAlphaSetParams", th->Alpha_f, &th->Alpha_f, NULL));
3709566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_gamma", "Algorithmic parameter gamma", "TSAlphaSetParams", th->Gamma, &th->Gamma, NULL));
3719566063dSJacob Faibussowitsch     PetscCall(TSAlphaSetParams(ts, th->Alpha_m, th->Alpha_f, th->Gamma));
372b07a2398SLisandro Dalcin   }
373d0609cedSBarry Smith   PetscOptionsHeadEnd();
3743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
375b07a2398SLisandro Dalcin }
376b07a2398SLisandro Dalcin 
377d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Alpha(TS ts, PetscViewer viewer)
378d71ae5a4SJacob Faibussowitsch {
379b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
3809c334d8fSLisandro Dalcin   PetscBool iascii;
381b07a2398SLisandro Dalcin 
382b07a2398SLisandro Dalcin   PetscFunctionBegin;
3839566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
38448a46eb9SPierre 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));
3853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
386b07a2398SLisandro Dalcin }
387b07a2398SLisandro Dalcin 
388d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlphaSetRadius_Alpha(TS ts, PetscReal radius)
389d71ae5a4SJacob Faibussowitsch {
390b07a2398SLisandro Dalcin   PetscReal alpha_m, alpha_f, gamma;
391b07a2398SLisandro Dalcin 
392b07a2398SLisandro Dalcin   PetscFunctionBegin;
393cad9d221SBarry Smith   PetscCheck(radius >= 0 && radius <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius);
394b07a2398SLisandro Dalcin   alpha_m = (PetscReal)0.5 * (3 - radius) / (1 + radius);
395b07a2398SLisandro Dalcin   alpha_f = 1 / (1 + radius);
396b07a2398SLisandro Dalcin   gamma   = (PetscReal)0.5 + alpha_m - alpha_f;
3979566063dSJacob Faibussowitsch   PetscCall(TSAlphaSetParams(ts, alpha_m, alpha_f, gamma));
3983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
399b07a2398SLisandro Dalcin }
400b07a2398SLisandro Dalcin 
401d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlphaSetParams_Alpha(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma)
402d71ae5a4SJacob Faibussowitsch {
403b07a2398SLisandro Dalcin   TS_Alpha *th  = (TS_Alpha *)ts->data;
404b07a2398SLisandro Dalcin   PetscReal tol = 100 * PETSC_MACHINE_EPSILON;
405b07a2398SLisandro Dalcin   PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma;
406b07a2398SLisandro Dalcin 
407b07a2398SLisandro Dalcin   PetscFunctionBegin;
408b07a2398SLisandro Dalcin   th->Alpha_m = alpha_m;
409b07a2398SLisandro Dalcin   th->Alpha_f = alpha_f;
410b07a2398SLisandro Dalcin   th->Gamma   = gamma;
411b07a2398SLisandro Dalcin   th->order   = (PetscAbsReal(res) < tol) ? 2 : 1;
4123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
413b07a2398SLisandro Dalcin }
414b07a2398SLisandro Dalcin 
415d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlphaGetParams_Alpha(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma)
416d71ae5a4SJacob Faibussowitsch {
417b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
418b07a2398SLisandro Dalcin 
419b07a2398SLisandro Dalcin   PetscFunctionBegin;
420b07a2398SLisandro Dalcin   if (alpha_m) *alpha_m = th->Alpha_m;
421b07a2398SLisandro Dalcin   if (alpha_f) *alpha_f = th->Alpha_f;
422b07a2398SLisandro Dalcin   if (gamma) *gamma = th->Gamma;
4233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
424b07a2398SLisandro Dalcin }
425b07a2398SLisandro Dalcin 
426b07a2398SLisandro Dalcin /*MC
427b07a2398SLisandro Dalcin       TSALPHA - ODE/DAE solver using the implicit Generalized-Alpha method
428b07a2398SLisandro Dalcin                 for first-order systems
429b07a2398SLisandro Dalcin 
430b07a2398SLisandro Dalcin   Level: beginner
431b07a2398SLisandro Dalcin 
432b07a2398SLisandro Dalcin   References:
433606c0280SSatish Balay + * - K.E. Jansen, C.H. Whiting, G.M. Hulber, "A generalized-alpha
434b07a2398SLisandro Dalcin   method for integrating the filtered Navier-Stokes equations with a
435b07a2398SLisandro Dalcin   stabilized finite element method", Computer Methods in Applied
436b07a2398SLisandro Dalcin   Mechanics and Engineering, 190, 305-319, 2000.
437b07a2398SLisandro Dalcin   DOI: 10.1016/S0045-7825(00)00203-6.
438606c0280SSatish Balay - * -  J. Chung, G.M.Hubert. "A Time Integration Algorithm for Structural
439b07a2398SLisandro Dalcin   Dynamics with Improved Numerical Dissipation: The Generalized-alpha
440b07a2398SLisandro Dalcin   Method" ASME Journal of Applied Mechanics, 60, 371:375, 1993.
441b07a2398SLisandro Dalcin 
442*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSCreate()`, `TSSetType()`, `TSAlphaSetRadius()`, `TSAlphaSetParams()`
443b07a2398SLisandro Dalcin M*/
444d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Alpha(TS ts)
445d71ae5a4SJacob Faibussowitsch {
446b07a2398SLisandro Dalcin   TS_Alpha *th;
447b07a2398SLisandro Dalcin 
448b07a2398SLisandro Dalcin   PetscFunctionBegin;
449b07a2398SLisandro Dalcin   ts->ops->reset          = TSReset_Alpha;
450b07a2398SLisandro Dalcin   ts->ops->destroy        = TSDestroy_Alpha;
451b07a2398SLisandro Dalcin   ts->ops->view           = TSView_Alpha;
452b07a2398SLisandro Dalcin   ts->ops->setup          = TSSetUp_Alpha;
453b07a2398SLisandro Dalcin   ts->ops->setfromoptions = TSSetFromOptions_Alpha;
454b07a2398SLisandro Dalcin   ts->ops->step           = TSStep_Alpha;
4559808bdc1SLisandro Dalcin   ts->ops->evaluatewlte   = TSEvaluateWLTE_Alpha;
456b07a2398SLisandro Dalcin   ts->ops->rollback       = TSRollBack_Alpha;
457b07a2398SLisandro Dalcin   ts->ops->interpolate    = TSInterpolate_Alpha;
458b07a2398SLisandro Dalcin   ts->ops->snesfunction   = SNESTSFormFunction_Alpha;
459b07a2398SLisandro Dalcin   ts->ops->snesjacobian   = SNESTSFormJacobian_Alpha;
4602ffb9264SLisandro Dalcin   ts->default_adapt_type  = TSADAPTNONE;
461b07a2398SLisandro Dalcin 
462efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
463efd4aadfSBarry Smith 
4644dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&th));
465b07a2398SLisandro Dalcin   ts->data = (void *)th;
466b07a2398SLisandro Dalcin 
467b07a2398SLisandro Dalcin   th->Alpha_m = 0.5;
468b07a2398SLisandro Dalcin   th->Alpha_f = 0.5;
469b07a2398SLisandro Dalcin   th->Gamma   = 0.5;
470b07a2398SLisandro Dalcin   th->order   = 2;
471b07a2398SLisandro Dalcin 
4729566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetRadius_C", TSAlphaSetRadius_Alpha));
4739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetParams_C", TSAlphaSetParams_Alpha));
4749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaGetParams_C", TSAlphaGetParams_Alpha));
4753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
476b07a2398SLisandro Dalcin }
477b07a2398SLisandro Dalcin 
478b07a2398SLisandro Dalcin /*@
479bcf0153eSBarry Smith   TSAlphaSetRadius - sets the desired spectral radius of the method for `TSALPHA`
480b07a2398SLisandro Dalcin                      (i.e. high-frequency numerical damping)
481b07a2398SLisandro Dalcin 
482c3339decSBarry Smith   Logically Collective
483b07a2398SLisandro Dalcin 
484b07a2398SLisandro Dalcin   The algorithmic parameters \alpha_m and \alpha_f of the
485b07a2398SLisandro Dalcin   generalized-\alpha method can be computed in terms of a specified
486b07a2398SLisandro Dalcin   spectral radius \rho in [0,1] for infinite time step in order to
487b07a2398SLisandro Dalcin   control high-frequency numerical damping:
488b07a2398SLisandro Dalcin     \alpha_m = 0.5*(3-\rho)/(1+\rho)
489b07a2398SLisandro Dalcin     \alpha_f = 1/(1+\rho)
490b07a2398SLisandro Dalcin 
491d8d19677SJose E. Roman   Input Parameters:
492b07a2398SLisandro Dalcin +  ts - timestepping context
493b07a2398SLisandro Dalcin -  radius - the desired spectral radius
494b07a2398SLisandro Dalcin 
495bcf0153eSBarry Smith   Options Database Key:
49667b8a455SSatish Balay .  -ts_alpha_radius <radius> - set alpha radius
497b07a2398SLisandro Dalcin 
498b07a2398SLisandro Dalcin   Level: intermediate
499b07a2398SLisandro Dalcin 
500*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA`, `TSAlphaSetParams()`, `TSAlphaGetParams()`
501b07a2398SLisandro Dalcin @*/
502d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlphaSetRadius(TS ts, PetscReal radius)
503d71ae5a4SJacob Faibussowitsch {
504b07a2398SLisandro Dalcin   PetscFunctionBegin;
505b07a2398SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
506b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, radius, 2);
507cad9d221SBarry Smith   PetscCheck(radius >= 0 && radius <= 1, ((PetscObject)ts)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius);
508cac4c232SBarry Smith   PetscTryMethod(ts, "TSAlphaSetRadius_C", (TS, PetscReal), (ts, radius));
5093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
510b07a2398SLisandro Dalcin }
511b07a2398SLisandro Dalcin 
512b07a2398SLisandro Dalcin /*@
513bcf0153eSBarry Smith   TSAlphaSetParams - sets the algorithmic parameters for `TSALPHA`
514b07a2398SLisandro Dalcin 
515c3339decSBarry Smith   Logically Collective
516b07a2398SLisandro Dalcin 
517b07a2398SLisandro Dalcin   Second-order accuracy can be obtained so long as:
518b07a2398SLisandro Dalcin     \gamma = 0.5 + alpha_m - alpha_f
519b07a2398SLisandro Dalcin 
520b07a2398SLisandro Dalcin   Unconditional stability requires:
521b07a2398SLisandro Dalcin     \alpha_m >= \alpha_f >= 0.5
522b07a2398SLisandro Dalcin 
523b07a2398SLisandro Dalcin   Backward Euler method is recovered with:
524b07a2398SLisandro Dalcin     \alpha_m = \alpha_f = gamma = 1
525b07a2398SLisandro Dalcin 
526d8d19677SJose E. Roman   Input Parameters:
527b07a2398SLisandro Dalcin +  ts - timestepping context
5282fe279fdSBarry Smith .  alpha_m - algorithmic parameter
5292fe279fdSBarry Smith .  alpha_f - algorithmic parameter
5302fe279fdSBarry Smith -  gamma   - algorithmic parameter
531b07a2398SLisandro Dalcin 
532bcf0153eSBarry Smith    Options Database Keys:
53367b8a455SSatish Balay +  -ts_alpha_alpha_m <alpha_m> - set alpha_m
53467b8a455SSatish Balay .  -ts_alpha_alpha_f <alpha_f> - set alpha_f
53567b8a455SSatish Balay -  -ts_alpha_gamma   <gamma> - set gamma
536b07a2398SLisandro Dalcin 
537bcf0153eSBarry Smith   Level: advanced
538bcf0153eSBarry Smith 
539b07a2398SLisandro Dalcin   Note:
5402fe279fdSBarry Smith   Use of this function is normally only required to hack `TSALPHA` to
541b07a2398SLisandro Dalcin   use a modified integration scheme. Users should call
542bcf0153eSBarry Smith   `TSAlphaSetRadius()` to set the desired spectral radius of the methods
543b07a2398SLisandro Dalcin   (i.e. high-frequency damping) in order so select optimal values for
544b07a2398SLisandro Dalcin   these parameters.
545b07a2398SLisandro Dalcin 
546*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA`, `TSAlphaSetRadius()`, `TSAlphaGetParams()`
547b07a2398SLisandro Dalcin @*/
548d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlphaSetParams(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma)
549d71ae5a4SJacob Faibussowitsch {
550b07a2398SLisandro Dalcin   PetscFunctionBegin;
551b07a2398SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
552b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, alpha_m, 2);
553b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, alpha_f, 3);
554b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, gamma, 4);
555cac4c232SBarry Smith   PetscTryMethod(ts, "TSAlphaSetParams_C", (TS, PetscReal, PetscReal, PetscReal), (ts, alpha_m, alpha_f, gamma));
5563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
557b07a2398SLisandro Dalcin }
558b07a2398SLisandro Dalcin 
559b07a2398SLisandro Dalcin /*@
560bcf0153eSBarry Smith   TSAlphaGetParams - gets the algorithmic parameters for `TSALPHA`
561b07a2398SLisandro Dalcin 
562b07a2398SLisandro Dalcin   Not Collective
563b07a2398SLisandro Dalcin 
564b07a2398SLisandro Dalcin   Input Parameter:
565b07a2398SLisandro Dalcin .  ts - timestepping context
566b07a2398SLisandro Dalcin 
567b07a2398SLisandro Dalcin   Output Parameters:
5682fe279fdSBarry Smith +  alpha_m - algorithmic parameter
5692fe279fdSBarry Smith .  alpha_f - algorithmic parameter
5702fe279fdSBarry Smith -  gamma   - algorithmic parameter
571b07a2398SLisandro Dalcin 
572bcf0153eSBarry Smith   Level: advanced
573bcf0153eSBarry Smith 
574b07a2398SLisandro Dalcin   Note:
575bcf0153eSBarry Smith   Use of this function is normally only required to hack `TSALPHA` to
576b07a2398SLisandro Dalcin   use a modified integration scheme. Users should call
577bcf0153eSBarry Smith   `TSAlphaSetRadius()` to set the high-frequency damping (i.e. spectral
578b07a2398SLisandro Dalcin   radius of the method) in order so select optimal values for these
579b07a2398SLisandro Dalcin   parameters.
580b07a2398SLisandro Dalcin 
581*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA`, `TSAlphaSetRadius()`, `TSAlphaSetParams()`
582b07a2398SLisandro Dalcin @*/
583d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlphaGetParams(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma)
584d71ae5a4SJacob Faibussowitsch {
585b07a2398SLisandro Dalcin   PetscFunctionBegin;
586b07a2398SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
587b07a2398SLisandro Dalcin   if (alpha_m) PetscValidRealPointer(alpha_m, 2);
588b07a2398SLisandro Dalcin   if (alpha_f) PetscValidRealPointer(alpha_f, 3);
589b07a2398SLisandro Dalcin   if (gamma) PetscValidRealPointer(gamma, 4);
590cac4c232SBarry Smith   PetscUseMethod(ts, "TSAlphaGetParams_C", (TS, PetscReal *, PetscReal *, PetscReal *), (ts, alpha_m, alpha_f, gamma));
5913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
592b07a2398SLisandro Dalcin }
593