xref: /petsc/src/ts/impls/implicit/alpha/alpha1.c (revision bcda9346efad4e5ba2d553af84eb238771ba1e25)
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 
TSResizeRegister_Alpha(TS ts,PetscBool reg)378ec9177eSStefano Zampini static PetscErrorCode TSResizeRegister_Alpha(TS ts, PetscBool reg)
388ec9177eSStefano Zampini {
398ec9177eSStefano Zampini   TS_Alpha *th = (TS_Alpha *)ts->data;
408ec9177eSStefano Zampini 
418ec9177eSStefano Zampini   PetscFunctionBegin;
42ecc87898SStefano Zampini   if (reg) {
43ecc87898SStefano Zampini     PetscCall(TSResizeRegisterVec(ts, "ts:theta:sol_prev", th->vec_sol_prev));
44ecc87898SStefano Zampini     PetscCall(TSResizeRegisterVec(ts, "ts:theta:X0", th->X0));
45ecc87898SStefano Zampini   } else {
46ecc87898SStefano Zampini     PetscCall(TSResizeRetrieveVec(ts, "ts:theta:sol_prev", &th->vec_sol_prev));
47ecc87898SStefano Zampini     PetscCall(PetscObjectReference((PetscObject)th->vec_sol_prev));
48ecc87898SStefano Zampini     PetscCall(TSResizeRetrieveVec(ts, "ts:theta:X0", &th->X0));
498ec9177eSStefano Zampini     PetscCall(PetscObjectReference((PetscObject)th->X0));
508ec9177eSStefano Zampini   }
518ec9177eSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
528ec9177eSStefano Zampini }
538ec9177eSStefano Zampini 
TSAlpha_StageTime(TS ts)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 
TSAlpha_StageVecs(TS ts,Vec X)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 
TSAlpha_SNESSolve(TS ts,Vec b,Vec x)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 */
TSAlpha_Restart(TS ts,PetscBool * initok)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));
155070f04f8SDavid Kamensky   PetscCall(VecAXPY(th->V0, -3 / time_step, X0));
156070f04f8SDavid Kamensky   PetscCall(VecAXPY(th->V0, +4 / time_step, X1));
157070f04f8SDavid Kamensky   PetscCall(VecAXPY(th->V0, -1 / 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:
1684e797053SDavid Kamensky   /* Revert TSAlpha to the initial state (t0,X0), but retain
1694e797053SDavid Kamensky      potential time step reduction by factor after failed solve. */
170b07a2398SLisandro Dalcin   if (initok) *initok = stageok;
1714e797053SDavid Kamensky   PetscCall(TSSetTimeStep(ts, 2 * ts->time_step));
1729566063dSJacob Faibussowitsch   PetscCall(TSAlphaSetParams(ts, alpha_m, alpha_f, gamma));
1739566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol, th->X0));
174b07a2398SLisandro Dalcin 
1759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X1));
1763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
177b07a2398SLisandro Dalcin }
178b07a2398SLisandro Dalcin 
TSStep_Alpha(TS ts)179d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Alpha(TS ts)
180d71ae5a4SJacob Faibussowitsch {
181b07a2398SLisandro Dalcin   TS_Alpha *th         = (TS_Alpha *)ts->data;
182b07a2398SLisandro Dalcin   PetscInt  rejections = 0;
183b07a2398SLisandro Dalcin   PetscBool stageok, accept = PETSC_TRUE;
184b07a2398SLisandro Dalcin   PetscReal next_time_step = ts->time_step;
185b07a2398SLisandro Dalcin 
186b07a2398SLisandro Dalcin   PetscFunctionBegin;
1879566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(citation, &cited));
188b07a2398SLisandro Dalcin 
189b07a2398SLisandro Dalcin   if (!ts->steprollback) {
1909566063dSJacob Faibussowitsch     if (th->vec_sol_prev) PetscCall(VecCopy(th->X0, th->vec_sol_prev));
1919566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_sol, th->X0));
1929566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->V1, th->V0));
193b07a2398SLisandro Dalcin   }
194b07a2398SLisandro Dalcin 
1951566a47fSLisandro Dalcin   th->status = TS_STEP_INCOMPLETE;
196b07a2398SLisandro Dalcin   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
197fecfb714SLisandro Dalcin     if (ts->steprestart) {
1989566063dSJacob Faibussowitsch       PetscCall(TSAlpha_Restart(ts, &stageok));
199fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
200b07a2398SLisandro Dalcin     }
201b07a2398SLisandro Dalcin 
2029566063dSJacob Faibussowitsch     PetscCall(TSAlpha_StageTime(ts));
2039566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->X0, th->X1));
2049566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, th->stage_time));
2059566063dSJacob Faibussowitsch     PetscCall(TSAlpha_SNESSolve(ts, NULL, th->X1));
2069566063dSJacob Faibussowitsch     PetscCall(TSPostStage(ts, th->stage_time, 0, &th->Xa));
2079566063dSJacob Faibussowitsch     PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->Xa, &stageok));
208fecfb714SLisandro Dalcin     if (!stageok) goto reject_step;
209b07a2398SLisandro Dalcin 
2101566a47fSLisandro Dalcin     th->status = TS_STEP_PENDING;
2119566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->X1, ts->vec_sol));
2129566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
2131566a47fSLisandro Dalcin     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
214b07a2398SLisandro Dalcin     if (!accept) {
2159566063dSJacob Faibussowitsch       PetscCall(VecCopy(th->X0, ts->vec_sol));
216be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
217be5899b3SLisandro Dalcin       goto reject_step;
218b07a2398SLisandro Dalcin     }
219b07a2398SLisandro Dalcin 
220b07a2398SLisandro Dalcin     ts->ptime += ts->time_step;
221b07a2398SLisandro Dalcin     ts->time_step = next_time_step;
222b07a2398SLisandro Dalcin     break;
223b07a2398SLisandro Dalcin 
224b07a2398SLisandro Dalcin   reject_step:
2259371c9d4SSatish Balay     ts->reject++;
2269371c9d4SSatish Balay     accept = PETSC_FALSE;
227b07a2398SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
228b07a2398SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
22963a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
230b07a2398SLisandro Dalcin     }
231b07a2398SLisandro Dalcin   }
2323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
233b07a2398SLisandro Dalcin }
234b07a2398SLisandro Dalcin 
TSEvaluateWLTE_Alpha(TS ts,NormType wnormtype,PetscInt * order,PetscReal * wlte)235d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateWLTE_Alpha(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte)
236d71ae5a4SJacob Faibussowitsch {
237b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
2389808bdc1SLisandro Dalcin   Vec       X  = th->X1;           /* X = solution */
2391566a47fSLisandro Dalcin   Vec       Y  = th->vec_lte_work; /* Y = X + LTE  */
2407453f775SEmil Constantinescu   PetscReal wltea, wlter;
241b07a2398SLisandro Dalcin 
242b07a2398SLisandro Dalcin   PetscFunctionBegin;
2439371c9d4SSatish Balay   if (!th->vec_sol_prev) {
2449371c9d4SSatish Balay     *wlte = -1;
2453ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2469371c9d4SSatish Balay   }
2479371c9d4SSatish Balay   if (!th->vec_lte_work) {
2489371c9d4SSatish Balay     *wlte = -1;
2493ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2509371c9d4SSatish Balay   }
251fecfb714SLisandro Dalcin   if (ts->steprestart) {
252fecfb714SLisandro Dalcin     /* th->vec_lte_work is set to the LTE in TSAlpha_Restart() */
2539566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Y, 1, X));
254b07a2398SLisandro Dalcin   } else {
255b07a2398SLisandro Dalcin     /* Compute LTE using backward differences with non-constant time step */
256be5899b3SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
257be5899b3SLisandro Dalcin     PetscReal   a = 1 + h_prev / h;
2589371c9d4SSatish Balay     PetscScalar scal[3];
2599371c9d4SSatish Balay     Vec         vecs[3];
2609371c9d4SSatish Balay     scal[0] = +1 / a;
2619371c9d4SSatish Balay     scal[1] = -1 / (a - 1);
2629371c9d4SSatish Balay     scal[2] = +1 / (a * (a - 1));
2639371c9d4SSatish Balay     vecs[0] = th->X1;
2649371c9d4SSatish Balay     vecs[1] = th->X0;
2659371c9d4SSatish Balay     vecs[2] = th->vec_sol_prev;
2669566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Y));
2679566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(Y, 3, scal, vecs));
268b07a2398SLisandro Dalcin   }
2699566063dSJacob Faibussowitsch   PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, wlte, &wltea, &wlter));
2709808bdc1SLisandro Dalcin   if (order) *order = 2;
2713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
272b07a2398SLisandro Dalcin }
273b07a2398SLisandro Dalcin 
TSInterpolate_Alpha(TS ts,PetscReal t,Vec X)274d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_Alpha(TS ts, PetscReal t, Vec X)
275d71ae5a4SJacob Faibussowitsch {
276b07a2398SLisandro Dalcin   TS_Alpha *th    = (TS_Alpha *)ts->data;
277b07a2398SLisandro Dalcin   PetscReal dt    = t - ts->ptime;
27851adc54cSHong Zhang   PetscReal Gamma = th->Gamma;
279b07a2398SLisandro Dalcin 
280b07a2398SLisandro Dalcin   PetscFunctionBegin;
28151adc54cSHong Zhang   PetscCall(VecWAXPY(th->V1, -1.0, th->X0, ts->vec_sol));
28251adc54cSHong Zhang   PetscCall(VecAXPBY(th->V1, 1 - 1 / Gamma, 1 / (Gamma * ts->time_step), th->V0));
2839566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol, X));
28451adc54cSHong Zhang   /* X = X + Gamma*dT*V1 */
2859566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, th->Gamma * dt, th->V1));
28651adc54cSHong Zhang   /* X = X + (1-Gamma)*dT*V0 */
2879566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, (1 - th->Gamma) * dt, th->V0));
2883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
289b07a2398SLisandro Dalcin }
290b07a2398SLisandro Dalcin 
SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes,Vec X,Vec F,TS ts)291d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes, Vec X, Vec F, 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 
297b07a2398SLisandro Dalcin   PetscFunctionBegin;
2989566063dSJacob Faibussowitsch   PetscCall(TSAlpha_StageVecs(ts, X));
299b07a2398SLisandro Dalcin   /* F = Function(ta,Xa,Va) */
3009566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, ta, Xa, Va, F, PETSC_FALSE));
3019566063dSJacob Faibussowitsch   PetscCall(VecScale(F, th->scale_F));
3023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
303b07a2398SLisandro Dalcin }
304b07a2398SLisandro Dalcin 
SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes,PETSC_UNUSED Vec X,Mat J,Mat P,TS ts)305d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes, PETSC_UNUSED Vec X, Mat J, Mat P, TS ts)
306d71ae5a4SJacob Faibussowitsch {
307b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
308b07a2398SLisandro Dalcin   PetscReal ta = th->stage_time;
309b07a2398SLisandro Dalcin   Vec       Xa = th->Xa, Va = th->Va;
310b07a2398SLisandro Dalcin   PetscReal dVdX = th->shift_V;
311b07a2398SLisandro Dalcin 
312b07a2398SLisandro Dalcin   PetscFunctionBegin;
313b07a2398SLisandro Dalcin   /* J,P = Jacobian(ta,Xa,Va) */
3149566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, ta, Xa, Va, dVdX, J, P, PETSC_FALSE));
3153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
316b07a2398SLisandro Dalcin }
317b07a2398SLisandro Dalcin 
TSReset_Alpha(TS ts)318d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Alpha(TS ts)
319d71ae5a4SJacob Faibussowitsch {
320b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
321b07a2398SLisandro Dalcin 
322b07a2398SLisandro Dalcin   PetscFunctionBegin;
3239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X0));
3249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Xa));
3259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X1));
3269566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->V0));
3279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Va));
3289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->V1));
3299566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_sol_prev));
3309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_lte_work));
3313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
332b07a2398SLisandro Dalcin }
333b07a2398SLisandro Dalcin 
TSDestroy_Alpha(TS ts)334d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Alpha(TS ts)
335d71ae5a4SJacob Faibussowitsch {
336b07a2398SLisandro Dalcin   PetscFunctionBegin;
3379566063dSJacob Faibussowitsch   PetscCall(TSReset_Alpha(ts));
3389566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
339b07a2398SLisandro Dalcin 
3409566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetRadius_C", NULL));
3419566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetParams_C", NULL));
3429566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaGetParams_C", NULL));
3433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
344b07a2398SLisandro Dalcin }
345b07a2398SLisandro Dalcin 
TSSetUp_Alpha(TS ts)346d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Alpha(TS ts)
347d71ae5a4SJacob Faibussowitsch {
348b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
3492ffb9264SLisandro Dalcin   PetscBool match;
350b07a2398SLisandro Dalcin 
351b07a2398SLisandro Dalcin   PetscFunctionBegin;
3528ec9177eSStefano Zampini   if (!th->X0) PetscCall(VecDuplicate(ts->vec_sol, &th->X0));
3539566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->Xa));
3549566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->X1));
3559566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->V0));
3569566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->Va));
3579566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->V1));
3581566a47fSLisandro Dalcin 
3599566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &ts->adapt));
3609566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt));
3619566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match));
3622ffb9264SLisandro Dalcin   if (!match) {
363ecc87898SStefano Zampini     if (!th->vec_sol_prev) PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev));
364ecc87898SStefano Zampini     if (!th->vec_lte_work) PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work));
365b07a2398SLisandro Dalcin   }
3661566a47fSLisandro Dalcin 
3679566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &ts->snes));
3683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
369b07a2398SLisandro Dalcin }
370b07a2398SLisandro Dalcin 
TSSetFromOptions_Alpha(TS ts,PetscOptionItems PetscOptionsObject)371ce78bad3SBarry Smith static PetscErrorCode TSSetFromOptions_Alpha(TS ts, PetscOptionItems PetscOptionsObject)
372d71ae5a4SJacob Faibussowitsch {
373b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
374b07a2398SLisandro Dalcin 
375b07a2398SLisandro Dalcin   PetscFunctionBegin;
376d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Generalized-Alpha ODE solver options");
377b07a2398SLisandro Dalcin   {
378b07a2398SLisandro Dalcin     PetscBool flg;
379b07a2398SLisandro Dalcin     PetscReal radius = 1;
3809566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_radius", "Spectral radius (high-frequency dissipation)", "TSAlphaSetRadius", radius, &radius, &flg));
3819566063dSJacob Faibussowitsch     if (flg) PetscCall(TSAlphaSetRadius(ts, radius));
3829566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_alpha_m", "Algorithmic parameter alpha_m", "TSAlphaSetParams", th->Alpha_m, &th->Alpha_m, NULL));
3839566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_alpha_f", "Algorithmic parameter alpha_f", "TSAlphaSetParams", th->Alpha_f, &th->Alpha_f, NULL));
3849566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_gamma", "Algorithmic parameter gamma", "TSAlphaSetParams", th->Gamma, &th->Gamma, NULL));
3859566063dSJacob Faibussowitsch     PetscCall(TSAlphaSetParams(ts, th->Alpha_m, th->Alpha_f, th->Gamma));
386b07a2398SLisandro Dalcin   }
387d0609cedSBarry Smith   PetscOptionsHeadEnd();
3883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
389b07a2398SLisandro Dalcin }
390b07a2398SLisandro Dalcin 
TSView_Alpha(TS ts,PetscViewer viewer)391d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Alpha(TS ts, PetscViewer viewer)
392d71ae5a4SJacob Faibussowitsch {
393b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
394*9f196a02SMartin Diehl   PetscBool isascii;
395b07a2398SLisandro Dalcin 
396b07a2398SLisandro Dalcin   PetscFunctionBegin;
397*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
398*9f196a02SMartin Diehl   if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "  Alpha_m=%g, Alpha_f=%g, Gamma=%g\n", (double)th->Alpha_m, (double)th->Alpha_f, (double)th->Gamma));
3993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
400b07a2398SLisandro Dalcin }
401b07a2398SLisandro Dalcin 
TSAlphaSetRadius_Alpha(TS ts,PetscReal radius)402d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlphaSetRadius_Alpha(TS ts, PetscReal radius)
403d71ae5a4SJacob Faibussowitsch {
404b07a2398SLisandro Dalcin   PetscReal alpha_m, alpha_f, gamma;
405b07a2398SLisandro Dalcin 
406b07a2398SLisandro Dalcin   PetscFunctionBegin;
407cad9d221SBarry Smith   PetscCheck(radius >= 0 && radius <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius);
408b07a2398SLisandro Dalcin   alpha_m = (PetscReal)0.5 * (3 - radius) / (1 + radius);
409b07a2398SLisandro Dalcin   alpha_f = 1 / (1 + radius);
410b07a2398SLisandro Dalcin   gamma   = (PetscReal)0.5 + alpha_m - alpha_f;
4119566063dSJacob Faibussowitsch   PetscCall(TSAlphaSetParams(ts, alpha_m, alpha_f, gamma));
4123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
413b07a2398SLisandro Dalcin }
414b07a2398SLisandro Dalcin 
TSAlphaSetParams_Alpha(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma)415d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlphaSetParams_Alpha(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma)
416d71ae5a4SJacob Faibussowitsch {
417b07a2398SLisandro Dalcin   TS_Alpha *th  = (TS_Alpha *)ts->data;
418b07a2398SLisandro Dalcin   PetscReal tol = 100 * PETSC_MACHINE_EPSILON;
419b07a2398SLisandro Dalcin   PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma;
420b07a2398SLisandro Dalcin 
421b07a2398SLisandro Dalcin   PetscFunctionBegin;
422b07a2398SLisandro Dalcin   th->Alpha_m = alpha_m;
423b07a2398SLisandro Dalcin   th->Alpha_f = alpha_f;
424b07a2398SLisandro Dalcin   th->Gamma   = gamma;
425b07a2398SLisandro Dalcin   th->order   = (PetscAbsReal(res) < tol) ? 2 : 1;
4263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
427b07a2398SLisandro Dalcin }
428b07a2398SLisandro Dalcin 
TSAlphaGetParams_Alpha(TS ts,PetscReal * alpha_m,PetscReal * alpha_f,PetscReal * gamma)429d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlphaGetParams_Alpha(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma)
430d71ae5a4SJacob Faibussowitsch {
431b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
432b07a2398SLisandro Dalcin 
433b07a2398SLisandro Dalcin   PetscFunctionBegin;
434b07a2398SLisandro Dalcin   if (alpha_m) *alpha_m = th->Alpha_m;
435b07a2398SLisandro Dalcin   if (alpha_f) *alpha_f = th->Alpha_f;
436b07a2398SLisandro Dalcin   if (gamma) *gamma = th->Gamma;
4373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
438b07a2398SLisandro Dalcin }
439b07a2398SLisandro Dalcin 
440b07a2398SLisandro Dalcin /*MC
4411d27aa22SBarry Smith   TSALPHA - ODE/DAE solver using the implicit Generalized-Alpha method {cite}`jansen_2000` {cite}`chung1993` for first-order systems
442b07a2398SLisandro Dalcin 
443b07a2398SLisandro Dalcin   Level: beginner
444b07a2398SLisandro Dalcin 
4451cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSCreate()`, `TSSetType()`, `TSAlphaSetRadius()`, `TSAlphaSetParams()`
446b07a2398SLisandro Dalcin M*/
TSCreate_Alpha(TS ts)447d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Alpha(TS ts)
448d71ae5a4SJacob Faibussowitsch {
449b07a2398SLisandro Dalcin   TS_Alpha *th;
450b07a2398SLisandro Dalcin 
451b07a2398SLisandro Dalcin   PetscFunctionBegin;
452b07a2398SLisandro Dalcin   ts->ops->reset          = TSReset_Alpha;
453b07a2398SLisandro Dalcin   ts->ops->destroy        = TSDestroy_Alpha;
454b07a2398SLisandro Dalcin   ts->ops->view           = TSView_Alpha;
455b07a2398SLisandro Dalcin   ts->ops->setup          = TSSetUp_Alpha;
456b07a2398SLisandro Dalcin   ts->ops->setfromoptions = TSSetFromOptions_Alpha;
457b07a2398SLisandro Dalcin   ts->ops->step           = TSStep_Alpha;
4589808bdc1SLisandro Dalcin   ts->ops->evaluatewlte   = TSEvaluateWLTE_Alpha;
459b07a2398SLisandro Dalcin   ts->ops->interpolate    = TSInterpolate_Alpha;
4608ec9177eSStefano Zampini   ts->ops->resizeregister = TSResizeRegister_Alpha;
461b07a2398SLisandro Dalcin   ts->ops->snesfunction   = SNESTSFormFunction_Alpha;
462b07a2398SLisandro Dalcin   ts->ops->snesjacobian   = SNESTSFormJacobian_Alpha;
4632ffb9264SLisandro Dalcin   ts->default_adapt_type  = TSADAPTNONE;
464b07a2398SLisandro Dalcin 
465efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
466efd4aadfSBarry Smith 
4674dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&th));
468b07a2398SLisandro Dalcin   ts->data = (void *)th;
469b07a2398SLisandro Dalcin 
470b07a2398SLisandro Dalcin   th->Alpha_m = 0.5;
471b07a2398SLisandro Dalcin   th->Alpha_f = 0.5;
472b07a2398SLisandro Dalcin   th->Gamma   = 0.5;
473b07a2398SLisandro Dalcin   th->order   = 2;
474b07a2398SLisandro Dalcin 
4759566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetRadius_C", TSAlphaSetRadius_Alpha));
4769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetParams_C", TSAlphaSetParams_Alpha));
4779566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaGetParams_C", TSAlphaGetParams_Alpha));
4783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
479b07a2398SLisandro Dalcin }
480b07a2398SLisandro Dalcin 
481b07a2398SLisandro Dalcin /*@
482bcf0153eSBarry Smith   TSAlphaSetRadius - sets the desired spectral radius of the method for `TSALPHA`
483b07a2398SLisandro Dalcin   (i.e. high-frequency numerical damping)
484b07a2398SLisandro Dalcin 
485c3339decSBarry Smith   Logically Collective
486b07a2398SLisandro Dalcin 
487d8d19677SJose E. Roman   Input Parameters:
488b07a2398SLisandro Dalcin + ts     - timestepping context
489b07a2398SLisandro Dalcin - radius - the desired spectral radius
490b07a2398SLisandro Dalcin 
491bcf0153eSBarry Smith   Options Database Key:
49267b8a455SSatish Balay . -ts_alpha_radius <radius> - set alpha radius
493b07a2398SLisandro Dalcin 
494b07a2398SLisandro Dalcin   Level: intermediate
495b07a2398SLisandro Dalcin 
49614d0ab18SJacob Faibussowitsch   Notes:
49714d0ab18SJacob Faibussowitsch   The algorithmic parameters $\alpha_m$ and $\alpha_f$ of the generalized-$\alpha$ method can
49814d0ab18SJacob Faibussowitsch   be computed in terms of a specified spectral radius $\rho$ in [0, 1] for infinite time step
49914d0ab18SJacob Faibussowitsch   in order to control high-frequency numerical damping\:
500562efe2eSBarry Smith 
50114d0ab18SJacob Faibussowitsch   $$
502562efe2eSBarry Smith   \begin{align*}
503562efe2eSBarry Smith   \alpha_m = 0.5*(3-\rho)/(1+\rho) \\
50414d0ab18SJacob Faibussowitsch   \alpha_f = 1/(1+\rho)
505562efe2eSBarry Smith   \end{align*}
50614d0ab18SJacob Faibussowitsch   $$
50714d0ab18SJacob Faibussowitsch 
5081cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA`, `TSAlphaSetParams()`, `TSAlphaGetParams()`
509b07a2398SLisandro Dalcin @*/
TSAlphaSetRadius(TS ts,PetscReal radius)510d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlphaSetRadius(TS ts, PetscReal radius)
511d71ae5a4SJacob Faibussowitsch {
512b07a2398SLisandro Dalcin   PetscFunctionBegin;
513b07a2398SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
514b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, radius, 2);
515cad9d221SBarry Smith   PetscCheck(radius >= 0 && radius <= 1, ((PetscObject)ts)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius);
516cac4c232SBarry Smith   PetscTryMethod(ts, "TSAlphaSetRadius_C", (TS, PetscReal), (ts, radius));
5173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
518b07a2398SLisandro Dalcin }
519b07a2398SLisandro Dalcin 
520b07a2398SLisandro Dalcin /*@
521bcf0153eSBarry Smith   TSAlphaSetParams - sets the algorithmic parameters for `TSALPHA`
522b07a2398SLisandro Dalcin 
523c3339decSBarry Smith   Logically Collective
524b07a2398SLisandro Dalcin 
525d8d19677SJose E. Roman   Input Parameters:
526b07a2398SLisandro Dalcin + ts      - timestepping context
5272fe279fdSBarry Smith . alpha_m - algorithmic parameter
5282fe279fdSBarry Smith . alpha_f - algorithmic parameter
5292fe279fdSBarry Smith - gamma   - algorithmic parameter
530b07a2398SLisandro Dalcin 
531bcf0153eSBarry Smith   Options Database Keys:
53267b8a455SSatish Balay + -ts_alpha_alpha_m <alpha_m> - set alpha_m
53367b8a455SSatish Balay . -ts_alpha_alpha_f <alpha_f> - set alpha_f
53467b8a455SSatish Balay - -ts_alpha_gamma   <gamma>   - set gamma
535b07a2398SLisandro Dalcin 
536bcf0153eSBarry Smith   Level: advanced
537bcf0153eSBarry Smith 
538b07a2398SLisandro Dalcin   Note:
539562efe2eSBarry Smith   Second-order accuracy can be obtained so long as\:  $\gamma = 0.5 + \alpha_m - \alpha_f$
54014d0ab18SJacob Faibussowitsch 
541562efe2eSBarry Smith   Unconditional stability requires\: $\alpha_m >= \alpha_f >= 0.5$
54214d0ab18SJacob Faibussowitsch 
543562efe2eSBarry Smith   Backward Euler method is recovered with\: $\alpha_m = \alpha_f = \gamma = 1$
54414d0ab18SJacob Faibussowitsch 
54514d0ab18SJacob Faibussowitsch   Use of this function is normally only required to hack `TSALPHA` to use a modified
54614d0ab18SJacob Faibussowitsch   integration scheme. Users should call `TSAlphaSetRadius()` to set the desired spectral radius
54714d0ab18SJacob Faibussowitsch   of the methods (i.e. high-frequency damping) in order so select optimal values for these
54814d0ab18SJacob Faibussowitsch   parameters.
549b07a2398SLisandro Dalcin 
5501cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA`, `TSAlphaSetRadius()`, `TSAlphaGetParams()`
551b07a2398SLisandro Dalcin @*/
TSAlphaSetParams(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma)552d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlphaSetParams(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma)
553d71ae5a4SJacob Faibussowitsch {
554b07a2398SLisandro Dalcin   PetscFunctionBegin;
555b07a2398SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
556b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, alpha_m, 2);
557b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, alpha_f, 3);
558b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, gamma, 4);
559cac4c232SBarry Smith   PetscTryMethod(ts, "TSAlphaSetParams_C", (TS, PetscReal, PetscReal, PetscReal), (ts, alpha_m, alpha_f, gamma));
5603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
561b07a2398SLisandro Dalcin }
562b07a2398SLisandro Dalcin 
563b07a2398SLisandro Dalcin /*@
564bcf0153eSBarry Smith   TSAlphaGetParams - gets the algorithmic parameters for `TSALPHA`
565b07a2398SLisandro Dalcin 
566b07a2398SLisandro Dalcin   Not Collective
567b07a2398SLisandro Dalcin 
568b07a2398SLisandro Dalcin   Input Parameter:
569b07a2398SLisandro Dalcin . ts - timestepping context
570b07a2398SLisandro Dalcin 
571b07a2398SLisandro Dalcin   Output Parameters:
5722fe279fdSBarry Smith + alpha_m - algorithmic parameter
5732fe279fdSBarry Smith . alpha_f - algorithmic parameter
5742fe279fdSBarry Smith - gamma   - algorithmic parameter
575b07a2398SLisandro Dalcin 
576bcf0153eSBarry Smith   Level: advanced
577bcf0153eSBarry Smith 
578b07a2398SLisandro Dalcin   Note:
57914d0ab18SJacob Faibussowitsch   Use of this function is normally only required to hack `TSALPHA` to use a modified
58014d0ab18SJacob Faibussowitsch   integration scheme. Users should call `TSAlphaSetRadius()` to set the high-frequency damping
58114d0ab18SJacob Faibussowitsch   (i.e. spectral radius of the method) in order so select optimal values for these parameters.
582b07a2398SLisandro Dalcin 
5831cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA`, `TSAlphaSetRadius()`, `TSAlphaSetParams()`
584b07a2398SLisandro Dalcin @*/
TSAlphaGetParams(TS ts,PetscReal * alpha_m,PetscReal * alpha_f,PetscReal * gamma)585d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlphaGetParams(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma)
586d71ae5a4SJacob Faibussowitsch {
587b07a2398SLisandro Dalcin   PetscFunctionBegin;
588b07a2398SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
5894f572ea9SToby Isaac   if (alpha_m) PetscAssertPointer(alpha_m, 2);
5904f572ea9SToby Isaac   if (alpha_f) PetscAssertPointer(alpha_f, 3);
5914f572ea9SToby Isaac   if (gamma) PetscAssertPointer(gamma, 4);
592cac4c232SBarry Smith   PetscUseMethod(ts, "TSAlphaGetParams_C", (TS, PetscReal *, PetscReal *, PetscReal *), (ts, alpha_m, alpha_f, gamma));
5933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
594b07a2398SLisandro Dalcin }
595