xref: /petsc/src/ts/impls/implicit/alpha/alpha2.c (revision 220f924abfb7c127e3e73e9ee375e6c993b7366b)
1 /*
2   Code for timestepping with implicit generalized-\alpha method
3   for second order systems.
4 */
5 #include <petsc/private/tsimpl.h> /*I   "petscts.h"   I*/
6 
7 static PetscBool  cited      = PETSC_FALSE;
8 static const char citation[] = "@article{Chung1993,\n"
9                                "  title   = {A Time Integration Algorithm for Structural Dynamics with Improved Numerical Dissipation: The Generalized-$\\alpha$ Method},\n"
10                                "  author  = {J. Chung, G. M. Hubert},\n"
11                                "  journal = {ASME Journal of Applied Mechanics},\n"
12                                "  volume  = {60},\n"
13                                "  number  = {2},\n"
14                                "  pages   = {371--375},\n"
15                                "  year    = {1993},\n"
16                                "  issn    = {0021-8936},\n"
17                                "  doi     = {http://dx.doi.org/10.1115/1.2900803}\n}\n";
18 
19 typedef struct {
20   PetscReal stage_time;
21   PetscReal shift_V;
22   PetscReal shift_A;
23   PetscReal scale_F;
24   Vec       X0, Xa, X1;
25   Vec       V0, Va, V1;
26   Vec       A0, Aa, A1;
27 
28   Vec vec_dot;
29 
30   PetscReal Alpha_m;
31   PetscReal Alpha_f;
32   PetscReal Gamma;
33   PetscReal Beta;
34   PetscInt  order;
35 
36   Vec vec_sol_prev;
37   Vec vec_dot_prev;
38   Vec vec_lte_work[2];
39 
40   TSStepStatus status;
41 
42   TSAlpha2Predictor predictor;
43   void             *predictor_ctx;
44 } TS_Alpha;
45 
46 /*@C
47   TSAlpha2SetPredictor - sets the callback for computing a predictor (i.e., initial guess
48   for the nonlinear solver).
49 
50   Input Parameters:
51 + ts        - timestepping context
52 . predictor - callback to set the predictor in each step
53 - ctx       - the application context, which may be set to `NULL` if not used
54 
55   Level: intermediate
56 
57   Notes:
58 
59   If this function is never called, a same-state-vector predictor will be used, i.e.,
60   the initial guess will be the converged solution from the previous time step, without regard
61   for the previous velocity or acceleration.
62 
63 .seealso: [](ch_ts), `TS`, `TSALPHA2`, `TSAlpha2Predictor`
64 @*/
65 PetscErrorCode TSAlpha2SetPredictor(TS ts, TSAlpha2Predictor predictor, void *ctx)
66 {
67   TS_Alpha *th = (TS_Alpha *)(ts->data);
68 
69   PetscFunctionBegin;
70   th->predictor     = predictor;
71   th->predictor_ctx = ctx;
72   PetscFunctionReturn(PETSC_SUCCESS);
73 }
74 
75 static PetscErrorCode TSAlpha_ApplyPredictor(TS ts, Vec X1)
76 {
77   /* Apply a custom predictor if set, or default to same-displacement. */
78   TS_Alpha *th = (TS_Alpha *)(ts->data);
79 
80   PetscFunctionBegin;
81   if (th->predictor) PetscCall(th->predictor(ts, th->X0, th->V0, th->A0, X1, th->predictor_ctx));
82   else PetscCall(VecCopy(th->X0, X1));
83   PetscFunctionReturn(PETSC_SUCCESS);
84 }
85 
86 static PetscErrorCode TSAlpha_StageTime(TS ts)
87 {
88   TS_Alpha *th      = (TS_Alpha *)ts->data;
89   PetscReal t       = ts->ptime;
90   PetscReal dt      = ts->time_step;
91   PetscReal Alpha_m = th->Alpha_m;
92   PetscReal Alpha_f = th->Alpha_f;
93   PetscReal Gamma   = th->Gamma;
94   PetscReal Beta    = th->Beta;
95 
96   PetscFunctionBegin;
97   th->stage_time = t + Alpha_f * dt;
98   th->shift_V    = Gamma / (dt * Beta);
99   th->shift_A    = Alpha_m / (Alpha_f * dt * dt * Beta);
100   th->scale_F    = 1 / Alpha_f;
101   PetscFunctionReturn(PETSC_SUCCESS);
102 }
103 
104 static PetscErrorCode TSAlpha_StageVecs(TS ts, Vec X)
105 {
106   TS_Alpha *th = (TS_Alpha *)ts->data;
107   Vec       X1 = X, V1 = th->V1, A1 = th->A1;
108   Vec       Xa = th->Xa, Va = th->Va, Aa = th->Aa;
109   Vec       X0 = th->X0, V0 = th->V0, A0 = th->A0;
110   PetscReal dt      = ts->time_step;
111   PetscReal Alpha_m = th->Alpha_m;
112   PetscReal Alpha_f = th->Alpha_f;
113   PetscReal Gamma   = th->Gamma;
114   PetscReal Beta    = th->Beta;
115 
116   PetscFunctionBegin;
117   /* A1 = ... */
118   PetscCall(VecWAXPY(A1, -1.0, X0, X1));
119   PetscCall(VecAXPY(A1, -dt, V0));
120   PetscCall(VecAXPBY(A1, -(1 - 2 * Beta) / (2 * Beta), 1 / (dt * dt * Beta), A0));
121   /* V1 = ... */
122   PetscCall(VecWAXPY(V1, (1.0 - Gamma) / Gamma, A0, A1));
123   PetscCall(VecAYPX(V1, dt * Gamma, V0));
124   /* Xa = X0 + Alpha_f*(X1-X0) */
125   PetscCall(VecWAXPY(Xa, -1.0, X0, X1));
126   PetscCall(VecAYPX(Xa, Alpha_f, X0));
127   /* Va = V0 + Alpha_f*(V1-V0) */
128   PetscCall(VecWAXPY(Va, -1.0, V0, V1));
129   PetscCall(VecAYPX(Va, Alpha_f, V0));
130   /* Aa = A0 + Alpha_m*(A1-A0) */
131   PetscCall(VecWAXPY(Aa, -1.0, A0, A1));
132   PetscCall(VecAYPX(Aa, Alpha_m, A0));
133   PetscFunctionReturn(PETSC_SUCCESS);
134 }
135 
136 static PetscErrorCode TSAlpha_SNESSolve(TS ts, Vec b, Vec x)
137 {
138   PetscInt nits, lits;
139 
140   PetscFunctionBegin;
141   PetscCall(SNESSolve(ts->snes, b, x));
142   PetscCall(SNESGetIterationNumber(ts->snes, &nits));
143   PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits));
144   ts->snes_its += nits;
145   ts->ksp_its += lits;
146   PetscFunctionReturn(PETSC_SUCCESS);
147 }
148 
149 /*
150   Compute a consistent initial state for the generalized-alpha method.
151   - Solve two successive backward Euler steps with halved time step.
152   - Compute the initial second time derivative using backward differences.
153   - If using adaptivity, estimate the LTE of the initial step.
154 */
155 static PetscErrorCode TSAlpha_Restart(TS ts, PetscBool *initok)
156 {
157   TS_Alpha *th = (TS_Alpha *)ts->data;
158   PetscReal time_step;
159   PetscReal alpha_m, alpha_f, gamma, beta;
160   Vec       X0 = ts->vec_sol, X1, X2 = th->X1;
161   Vec       V0 = ts->vec_dot, V1, V2 = th->V1;
162   PetscBool stageok;
163 
164   PetscFunctionBegin;
165   PetscCall(VecDuplicate(X0, &X1));
166   PetscCall(VecDuplicate(V0, &V1));
167 
168   /* Setup backward Euler with halved time step */
169   PetscCall(TSAlpha2GetParams(ts, &alpha_m, &alpha_f, &gamma, &beta));
170   PetscCall(TSAlpha2SetParams(ts, 1, 1, 1, 0.5));
171   PetscCall(TSGetTimeStep(ts, &time_step));
172   ts->time_step = time_step / 2;
173   PetscCall(TSAlpha_StageTime(ts));
174   th->stage_time = ts->ptime;
175   PetscCall(VecZeroEntries(th->A0));
176 
177   /* First BE step, (t0,X0,V0) -> (t1,X1,V1) */
178   th->stage_time += ts->time_step;
179   PetscCall(VecCopy(X0, th->X0));
180   PetscCall(VecCopy(V0, th->V0));
181   PetscCall(TSPreStage(ts, th->stage_time));
182   PetscCall(TSAlpha_ApplyPredictor(ts, X1));
183   PetscCall(TSAlpha_SNESSolve(ts, NULL, X1));
184   PetscCall(VecCopy(th->V1, V1));
185   PetscCall(TSPostStage(ts, th->stage_time, 0, &X1));
186   PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, X1, &stageok));
187   if (!stageok) goto finally;
188 
189   /* Second BE step, (t1,X1,V1) -> (t2,X2,V2) */
190   th->stage_time += ts->time_step;
191   PetscCall(VecCopy(X1, th->X0));
192   PetscCall(VecCopy(V1, th->V0));
193   PetscCall(TSPreStage(ts, th->stage_time));
194   PetscCall(TSAlpha_ApplyPredictor(ts, X2));
195   PetscCall(TSAlpha_SNESSolve(ts, NULL, X2));
196   PetscCall(VecCopy(th->V1, V2));
197   PetscCall(TSPostStage(ts, th->stage_time, 0, &X2));
198   PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, X1, &stageok));
199   if (!stageok) goto finally;
200 
201   /* Compute A0 ~ dV/dt at t0 with backward differences */
202   PetscCall(VecZeroEntries(th->A0));
203   PetscCall(VecAXPY(th->A0, -3 / ts->time_step, V0));
204   PetscCall(VecAXPY(th->A0, +4 / ts->time_step, V1));
205   PetscCall(VecAXPY(th->A0, -1 / ts->time_step, V2));
206 
207   /* Rough, lower-order estimate LTE of the initial step */
208   if (th->vec_lte_work[0]) {
209     PetscCall(VecZeroEntries(th->vec_lte_work[0]));
210     PetscCall(VecAXPY(th->vec_lte_work[0], +2, X2));
211     PetscCall(VecAXPY(th->vec_lte_work[0], -4, X1));
212     PetscCall(VecAXPY(th->vec_lte_work[0], +2, X0));
213   }
214   if (th->vec_lte_work[1]) {
215     PetscCall(VecZeroEntries(th->vec_lte_work[1]));
216     PetscCall(VecAXPY(th->vec_lte_work[1], +2, V2));
217     PetscCall(VecAXPY(th->vec_lte_work[1], -4, V1));
218     PetscCall(VecAXPY(th->vec_lte_work[1], +2, V0));
219   }
220 
221 finally:
222   /* Revert TSAlpha to the initial state (t0,X0,V0) */
223   if (initok) *initok = stageok;
224   PetscCall(TSSetTimeStep(ts, time_step));
225   PetscCall(TSAlpha2SetParams(ts, alpha_m, alpha_f, gamma, beta));
226   PetscCall(VecCopy(ts->vec_sol, th->X0));
227   PetscCall(VecCopy(ts->vec_dot, th->V0));
228 
229   PetscCall(VecDestroy(&X1));
230   PetscCall(VecDestroy(&V1));
231   PetscFunctionReturn(PETSC_SUCCESS);
232 }
233 
234 static PetscErrorCode TSStep_Alpha(TS ts)
235 {
236   TS_Alpha *th         = (TS_Alpha *)ts->data;
237   PetscInt  rejections = 0;
238   PetscBool stageok, accept = PETSC_TRUE;
239   PetscReal next_time_step = ts->time_step;
240 
241   PetscFunctionBegin;
242   PetscCall(PetscCitationsRegister(citation, &cited));
243 
244   if (!ts->steprollback) {
245     if (th->vec_sol_prev) PetscCall(VecCopy(th->X0, th->vec_sol_prev));
246     if (th->vec_dot_prev) PetscCall(VecCopy(th->V0, th->vec_dot_prev));
247     PetscCall(VecCopy(ts->vec_sol, th->X0));
248     PetscCall(VecCopy(ts->vec_dot, th->V0));
249     PetscCall(VecCopy(th->A1, th->A0));
250   }
251 
252   th->status = TS_STEP_INCOMPLETE;
253   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
254     if (ts->steprestart) {
255       PetscCall(TSAlpha_Restart(ts, &stageok));
256       if (!stageok) goto reject_step;
257     }
258 
259     PetscCall(TSAlpha_StageTime(ts));
260     PetscCall(TSAlpha_ApplyPredictor(ts, th->X1));
261     PetscCall(TSPreStage(ts, th->stage_time));
262     PetscCall(TSAlpha_SNESSolve(ts, NULL, th->X1));
263     PetscCall(TSPostStage(ts, th->stage_time, 0, &th->Xa));
264     PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->Xa, &stageok));
265     if (!stageok) goto reject_step;
266 
267     th->status = TS_STEP_PENDING;
268     PetscCall(VecCopy(th->X1, ts->vec_sol));
269     PetscCall(VecCopy(th->V1, ts->vec_dot));
270     PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
271     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
272     if (!accept) {
273       PetscCall(VecCopy(th->X0, ts->vec_sol));
274       PetscCall(VecCopy(th->V0, ts->vec_dot));
275       ts->time_step = next_time_step;
276       goto reject_step;
277     }
278 
279     ts->ptime += ts->time_step;
280     ts->time_step = next_time_step;
281     break;
282 
283   reject_step:
284     ts->reject++;
285     accept = PETSC_FALSE;
286     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
287       ts->reason = TS_DIVERGED_STEP_REJECTED;
288       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
289     }
290   }
291   PetscFunctionReturn(PETSC_SUCCESS);
292 }
293 
294 static PetscErrorCode TSEvaluateWLTE_Alpha(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte)
295 {
296   TS_Alpha *th = (TS_Alpha *)ts->data;
297   Vec       X  = th->X1;              /* X = solution */
298   Vec       V  = th->V1;              /* V = solution */
299   Vec       Y  = th->vec_lte_work[0]; /* Y = X + LTE  */
300   Vec       Z  = th->vec_lte_work[1]; /* Z = V + LTE  */
301   PetscReal enormX, enormV, enormXa, enormVa, enormXr, enormVr;
302 
303   PetscFunctionBegin;
304   if (!th->vec_sol_prev) {
305     *wlte = -1;
306     PetscFunctionReturn(PETSC_SUCCESS);
307   }
308   if (!th->vec_dot_prev) {
309     *wlte = -1;
310     PetscFunctionReturn(PETSC_SUCCESS);
311   }
312   if (!th->vec_lte_work[0]) {
313     *wlte = -1;
314     PetscFunctionReturn(PETSC_SUCCESS);
315   }
316   if (!th->vec_lte_work[1]) {
317     *wlte = -1;
318     PetscFunctionReturn(PETSC_SUCCESS);
319   }
320   if (ts->steprestart) {
321     /* th->vec_lte_prev is set to the LTE in TSAlpha_Restart() */
322     PetscCall(VecAXPY(Y, 1, X));
323     PetscCall(VecAXPY(Z, 1, V));
324   } else {
325     /* Compute LTE using backward differences with non-constant time step */
326     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
327     PetscReal   a = 1 + h_prev / h;
328     PetscScalar scal[3];
329     Vec         vecX[3], vecV[3];
330     scal[0] = +1 / a;
331     scal[1] = -1 / (a - 1);
332     scal[2] = +1 / (a * (a - 1));
333     vecX[0] = th->X1;
334     vecX[1] = th->X0;
335     vecX[2] = th->vec_sol_prev;
336     vecV[0] = th->V1;
337     vecV[1] = th->V0;
338     vecV[2] = th->vec_dot_prev;
339     PetscCall(VecCopy(X, Y));
340     PetscCall(VecMAXPY(Y, 3, scal, vecX));
341     PetscCall(VecCopy(V, Z));
342     PetscCall(VecMAXPY(Z, 3, scal, vecV));
343   }
344   /* XXX ts->atol and ts->vatol are not appropriate for computing enormV */
345   PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, &enormX, &enormXa, &enormXr));
346   PetscCall(TSErrorWeightedNorm(ts, V, Z, wnormtype, &enormV, &enormVa, &enormVr));
347   if (wnormtype == NORM_2) *wlte = PetscSqrtReal(PetscSqr(enormX) / 2 + PetscSqr(enormV) / 2);
348   else *wlte = PetscMax(enormX, enormV);
349   if (order) *order = 2;
350   PetscFunctionReturn(PETSC_SUCCESS);
351 }
352 
353 static PetscErrorCode TSRollBack_Alpha(TS ts)
354 {
355   TS_Alpha *th = (TS_Alpha *)ts->data;
356 
357   PetscFunctionBegin;
358   PetscCall(VecCopy(th->X0, ts->vec_sol));
359   PetscCall(VecCopy(th->V0, ts->vec_dot));
360   PetscFunctionReturn(PETSC_SUCCESS);
361 }
362 
363 /*
364 static PetscErrorCode TSInterpolate_Alpha(TS ts,PetscReal t,Vec X,Vec V)
365 {
366   TS_Alpha       *th = (TS_Alpha*)ts->data;
367   PetscReal      dt  = t - ts->ptime;
368 
369   PetscFunctionBegin;
370   PetscCall(VecCopy(ts->vec_dot,V));
371   PetscCall(VecAXPY(V,dt*(1-th->Gamma),th->A0));
372   PetscCall(VecAXPY(V,dt*th->Gamma,th->A1));
373   PetscCall(VecCopy(ts->vec_sol,X));
374   PetscCall(VecAXPY(X,dt,V));
375   PetscCall(VecAXPY(X,dt*dt*((PetscReal)0.5-th->Beta),th->A0));
376   PetscCall(VecAXPY(X,dt*dt*th->Beta,th->A1));
377   PetscFunctionReturn(PETSC_SUCCESS);
378 }
379 */
380 
381 static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes, Vec X, Vec F, TS ts)
382 {
383   TS_Alpha *th = (TS_Alpha *)ts->data;
384   PetscReal ta = th->stage_time;
385   Vec       Xa = th->Xa, Va = th->Va, Aa = th->Aa;
386 
387   PetscFunctionBegin;
388   PetscCall(TSAlpha_StageVecs(ts, X));
389   /* F = Function(ta,Xa,Va,Aa) */
390   PetscCall(TSComputeI2Function(ts, ta, Xa, Va, Aa, F));
391   PetscCall(VecScale(F, th->scale_F));
392   PetscFunctionReturn(PETSC_SUCCESS);
393 }
394 
395 static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes, PETSC_UNUSED Vec X, Mat J, Mat P, TS ts)
396 {
397   TS_Alpha *th = (TS_Alpha *)ts->data;
398   PetscReal ta = th->stage_time;
399   Vec       Xa = th->Xa, Va = th->Va, Aa = th->Aa;
400   PetscReal dVdX = th->shift_V, dAdX = th->shift_A;
401 
402   PetscFunctionBegin;
403   /* J,P = Jacobian(ta,Xa,Va,Aa) */
404   PetscCall(TSComputeI2Jacobian(ts, ta, Xa, Va, Aa, dVdX, dAdX, J, P));
405   PetscFunctionReturn(PETSC_SUCCESS);
406 }
407 
408 static PetscErrorCode TSReset_Alpha(TS ts)
409 {
410   TS_Alpha *th = (TS_Alpha *)ts->data;
411 
412   PetscFunctionBegin;
413   PetscCall(VecDestroy(&th->X0));
414   PetscCall(VecDestroy(&th->Xa));
415   PetscCall(VecDestroy(&th->X1));
416   PetscCall(VecDestroy(&th->V0));
417   PetscCall(VecDestroy(&th->Va));
418   PetscCall(VecDestroy(&th->V1));
419   PetscCall(VecDestroy(&th->A0));
420   PetscCall(VecDestroy(&th->Aa));
421   PetscCall(VecDestroy(&th->A1));
422   PetscCall(VecDestroy(&th->vec_sol_prev));
423   PetscCall(VecDestroy(&th->vec_dot_prev));
424   PetscCall(VecDestroy(&th->vec_lte_work[0]));
425   PetscCall(VecDestroy(&th->vec_lte_work[1]));
426   PetscFunctionReturn(PETSC_SUCCESS);
427 }
428 
429 static PetscErrorCode TSDestroy_Alpha(TS ts)
430 {
431   PetscFunctionBegin;
432   PetscCall(TSReset_Alpha(ts));
433   PetscCall(PetscFree(ts->data));
434 
435   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetRadius_C", NULL));
436   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetParams_C", NULL));
437   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2GetParams_C", NULL));
438   PetscFunctionReturn(PETSC_SUCCESS);
439 }
440 
441 static PetscErrorCode TSSetUp_Alpha(TS ts)
442 {
443   TS_Alpha *th = (TS_Alpha *)ts->data;
444   PetscBool match;
445 
446   PetscFunctionBegin;
447   PetscCall(VecDuplicate(ts->vec_sol, &th->X0));
448   PetscCall(VecDuplicate(ts->vec_sol, &th->Xa));
449   PetscCall(VecDuplicate(ts->vec_sol, &th->X1));
450   PetscCall(VecDuplicate(ts->vec_sol, &th->V0));
451   PetscCall(VecDuplicate(ts->vec_sol, &th->Va));
452   PetscCall(VecDuplicate(ts->vec_sol, &th->V1));
453   PetscCall(VecDuplicate(ts->vec_sol, &th->A0));
454   PetscCall(VecDuplicate(ts->vec_sol, &th->Aa));
455   PetscCall(VecDuplicate(ts->vec_sol, &th->A1));
456 
457   PetscCall(TSGetAdapt(ts, &ts->adapt));
458   PetscCall(TSAdaptCandidatesClear(ts->adapt));
459   PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match));
460   if (!match) {
461     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev));
462     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_dot_prev));
463     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work[0]));
464     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work[1]));
465   }
466 
467   PetscCall(TSGetSNES(ts, &ts->snes));
468   PetscFunctionReturn(PETSC_SUCCESS);
469 }
470 
471 static PetscErrorCode TSSetFromOptions_Alpha(TS ts, PetscOptionItems *PetscOptionsObject)
472 {
473   TS_Alpha *th = (TS_Alpha *)ts->data;
474 
475   PetscFunctionBegin;
476   PetscOptionsHeadBegin(PetscOptionsObject, "Generalized-Alpha ODE solver options");
477   {
478     PetscBool flg;
479     PetscReal radius = 1;
480     PetscCall(PetscOptionsReal("-ts_alpha_radius", "Spectral radius (high-frequency dissipation)", "TSAlpha2SetRadius", radius, &radius, &flg));
481     if (flg) PetscCall(TSAlpha2SetRadius(ts, radius));
482     PetscCall(PetscOptionsReal("-ts_alpha_alpha_m", "Algorithmic parameter alpha_m", "TSAlpha2SetParams", th->Alpha_m, &th->Alpha_m, NULL));
483     PetscCall(PetscOptionsReal("-ts_alpha_alpha_f", "Algorithmic parameter alpha_f", "TSAlpha2SetParams", th->Alpha_f, &th->Alpha_f, NULL));
484     PetscCall(PetscOptionsReal("-ts_alpha_gamma", "Algorithmic parameter gamma", "TSAlpha2SetParams", th->Gamma, &th->Gamma, NULL));
485     PetscCall(PetscOptionsReal("-ts_alpha_beta", "Algorithmic parameter beta", "TSAlpha2SetParams", th->Beta, &th->Beta, NULL));
486     PetscCall(TSAlpha2SetParams(ts, th->Alpha_m, th->Alpha_f, th->Gamma, th->Beta));
487   }
488   PetscOptionsHeadEnd();
489   PetscFunctionReturn(PETSC_SUCCESS);
490 }
491 
492 static PetscErrorCode TSView_Alpha(TS ts, PetscViewer viewer)
493 {
494   TS_Alpha *th = (TS_Alpha *)ts->data;
495   PetscBool iascii;
496 
497   PetscFunctionBegin;
498   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
499   if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, "  Alpha_m=%g, Alpha_f=%g, Gamma=%g, Beta=%g\n", (double)th->Alpha_m, (double)th->Alpha_f, (double)th->Gamma, (double)th->Beta));
500   PetscFunctionReturn(PETSC_SUCCESS);
501 }
502 
503 static PetscErrorCode TSAlpha2SetRadius_Alpha(TS ts, PetscReal radius)
504 {
505   PetscReal alpha_m, alpha_f, gamma, beta;
506 
507   PetscFunctionBegin;
508   PetscCheck(radius >= 0 && radius <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius);
509   alpha_m = (2 - radius) / (1 + radius);
510   alpha_f = 1 / (1 + radius);
511   gamma   = (PetscReal)0.5 + alpha_m - alpha_f;
512   beta    = (PetscReal)0.5 * (1 + alpha_m - alpha_f);
513   beta *= beta;
514   PetscCall(TSAlpha2SetParams(ts, alpha_m, alpha_f, gamma, beta));
515   PetscFunctionReturn(PETSC_SUCCESS);
516 }
517 
518 static PetscErrorCode TSAlpha2SetParams_Alpha(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma, PetscReal beta)
519 {
520   TS_Alpha *th  = (TS_Alpha *)ts->data;
521   PetscReal tol = 100 * PETSC_MACHINE_EPSILON;
522   PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma;
523 
524   PetscFunctionBegin;
525   th->Alpha_m = alpha_m;
526   th->Alpha_f = alpha_f;
527   th->Gamma   = gamma;
528   th->Beta    = beta;
529   th->order   = (PetscAbsReal(res) < tol) ? 2 : 1;
530   PetscFunctionReturn(PETSC_SUCCESS);
531 }
532 
533 static PetscErrorCode TSAlpha2GetParams_Alpha(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma, PetscReal *beta)
534 {
535   TS_Alpha *th = (TS_Alpha *)ts->data;
536 
537   PetscFunctionBegin;
538   if (alpha_m) *alpha_m = th->Alpha_m;
539   if (alpha_f) *alpha_f = th->Alpha_f;
540   if (gamma) *gamma = th->Gamma;
541   if (beta) *beta = th->Beta;
542   PetscFunctionReturn(PETSC_SUCCESS);
543 }
544 
545 /*MC
546   TSALPHA2 - ODE/DAE solver using the implicit Generalized-Alpha method for second-order systems
547 
548   Level: beginner
549 
550   References:
551 . * - J. Chung, G.M.Hubert. "A Time Integration Algorithm for Structural
552   Dynamics with Improved Numerical Dissipation: The Generalized-alpha
553   Method" ASME Journal of Applied Mechanics, 60, 371:375, 1993.
554 
555 .seealso: [](ch_ts), `TS`, `TSCreate()`, `TSSetType()`, `TSAlpha2SetRadius()`, `TSAlpha2SetParams()`
556 M*/
557 PETSC_EXTERN PetscErrorCode TSCreate_Alpha2(TS ts)
558 {
559   TS_Alpha *th;
560 
561   PetscFunctionBegin;
562   ts->ops->reset          = TSReset_Alpha;
563   ts->ops->destroy        = TSDestroy_Alpha;
564   ts->ops->view           = TSView_Alpha;
565   ts->ops->setup          = TSSetUp_Alpha;
566   ts->ops->setfromoptions = TSSetFromOptions_Alpha;
567   ts->ops->step           = TSStep_Alpha;
568   ts->ops->evaluatewlte   = TSEvaluateWLTE_Alpha;
569   ts->ops->rollback       = TSRollBack_Alpha;
570   /*ts->ops->interpolate  = TSInterpolate_Alpha;*/
571   ts->ops->snesfunction  = SNESTSFormFunction_Alpha;
572   ts->ops->snesjacobian  = SNESTSFormJacobian_Alpha;
573   ts->default_adapt_type = TSADAPTNONE;
574 
575   ts->usessnes = PETSC_TRUE;
576 
577   PetscCall(PetscNew(&th));
578   ts->data = (void *)th;
579 
580   th->Alpha_m = 0.5;
581   th->Alpha_f = 0.5;
582   th->Gamma   = 0.5;
583   th->Beta    = 0.25;
584   th->order   = 2;
585 
586   th->predictor     = NULL;
587   th->predictor_ctx = NULL;
588 
589   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetRadius_C", TSAlpha2SetRadius_Alpha));
590   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetParams_C", TSAlpha2SetParams_Alpha));
591   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2GetParams_C", TSAlpha2GetParams_Alpha));
592   PetscFunctionReturn(PETSC_SUCCESS);
593 }
594 
595 /*@
596   TSAlpha2SetRadius - sets the desired spectral radius of the method for `TSALPHA2`
597   (i.e. high-frequency numerical damping)
598 
599   Logically Collective
600 
601   Input Parameters:
602 + ts     - timestepping context
603 - radius - the desired spectral radius
604 
605   Options Database Key:
606 . -ts_alpha_radius <radius> - set the desired spectral radius
607 
608   Level: intermediate
609 
610   Notes:
611 
612   The algorithmic parameters $\alpha_m$ and $\alpha_f$ of the generalized-$\alpha$ method can
613   be computed in terms of a specified spectral radius $\rho$ in `[0, 1]` for infinite time step
614   in order to control high-frequency numerical damping\:
615   $$
616   \alpha_m = (2-\rho)/(1+\rho)
617   \alpha_f = 1/(1+\rho)
618   $$
619 
620 .seealso: [](ch_ts), `TS`, `TSALPHA2`, `TSAlpha2SetParams()`, `TSAlpha2GetParams()`
621 @*/
622 PetscErrorCode TSAlpha2SetRadius(TS ts, PetscReal radius)
623 {
624   PetscFunctionBegin;
625   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
626   PetscValidLogicalCollectiveReal(ts, radius, 2);
627   PetscCheck(radius >= 0 && radius <= 1, ((PetscObject)ts)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius);
628   PetscTryMethod(ts, "TSAlpha2SetRadius_C", (TS, PetscReal), (ts, radius));
629   PetscFunctionReturn(PETSC_SUCCESS);
630 }
631 
632 /*@
633   TSAlpha2SetParams - sets the algorithmic parameters for `TSALPHA2`
634 
635   Logically Collective
636 
637   Input Parameters:
638 + ts      - timestepping context
639 . alpha_m - algorithmic parameter
640 . alpha_f - algorithmic parameter
641 . gamma   - algorithmic parameter
642 - beta    - algorithmic parameter
643 
644   Options Database Keys:
645 + -ts_alpha_alpha_m <alpha_m> - set alpha_m
646 . -ts_alpha_alpha_f <alpha_f> - set alpha_f
647 . -ts_alpha_gamma   <gamma>   - set gamma
648 - -ts_alpha_beta    <beta>    - set beta
649 
650   Level: advanced
651 
652   Notes:
653   Second-order accuracy can be obtained so long as\:
654   $$
655   \gamma = 1/2 + alpha_m - alpha_f
656   \beta  = 1/4 (1 + alpha_m - alpha_f)^2
657   $$
658 
659   Unconditional stability requires\:
660   $$
661   \alpha_m >= \alpha_f >= 1/2
662   $$
663 
664   Use of this function is normally only required to hack `TSALPHA2` to use a modified
665   integration scheme. Users should call `TSAlpha2SetRadius()` to set the desired spectral
666   radius of the methods (i.e. high-frequency damping) in order so select optimal values for
667   these parameters.
668 
669 .seealso: [](ch_ts), `TS`, `TSALPHA2`, `TSAlpha2SetRadius()`, `TSAlpha2GetParams()`
670 @*/
671 PetscErrorCode TSAlpha2SetParams(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma, PetscReal beta)
672 {
673   PetscFunctionBegin;
674   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
675   PetscValidLogicalCollectiveReal(ts, alpha_m, 2);
676   PetscValidLogicalCollectiveReal(ts, alpha_f, 3);
677   PetscValidLogicalCollectiveReal(ts, gamma, 4);
678   PetscValidLogicalCollectiveReal(ts, beta, 5);
679   PetscTryMethod(ts, "TSAlpha2SetParams_C", (TS, PetscReal, PetscReal, PetscReal, PetscReal), (ts, alpha_m, alpha_f, gamma, beta));
680   PetscFunctionReturn(PETSC_SUCCESS);
681 }
682 
683 /*@
684   TSAlpha2GetParams - gets the algorithmic parameters for `TSALPHA2`
685 
686   Not Collective
687 
688   Input Parameter:
689 . ts - timestepping context
690 
691   Output Parameters:
692 + alpha_m - algorithmic parameter
693 . alpha_f - algorithmic parameter
694 . gamma   - algorithmic parameter
695 - beta    - algorithmic parameter
696 
697   Level: advanced
698 
699   Note:
700   Use of this function is normally only required to hack `TSALPHA2` to use a modified
701   integration scheme. Users should call `TSAlpha2SetRadius()` to set the high-frequency damping
702   (i.e. spectral radius of the method) in order so select optimal values for these parameters.
703 
704 .seealso: [](ch_ts), `TS`, `TSALPHA2`, `TSAlpha2SetRadius()`, `TSAlpha2SetParams()`
705 @*/
706 PetscErrorCode TSAlpha2GetParams(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma, PetscReal *beta)
707 {
708   PetscFunctionBegin;
709   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
710   if (alpha_m) PetscAssertPointer(alpha_m, 2);
711   if (alpha_f) PetscAssertPointer(alpha_f, 3);
712   if (gamma) PetscAssertPointer(gamma, 4);
713   if (beta) PetscAssertPointer(beta, 5);
714   PetscUseMethod(ts, "TSAlpha2GetParams_C", (TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *), (ts, alpha_m, alpha_f, gamma, beta));
715   PetscFunctionReturn(PETSC_SUCCESS);
716 }
717