xref: /petsc/src/ts/impls/implicit/alpha/alpha2.c (revision 562efe2ef922487c6beae96ba39e19afd4eefbe6)
1818efac9SLisandro Dalcin /*
2818efac9SLisandro Dalcin   Code for timestepping with implicit generalized-\alpha method
3818efac9SLisandro Dalcin   for second order systems.
4818efac9SLisandro Dalcin */
5818efac9SLisandro Dalcin #include <petsc/private/tsimpl.h> /*I   "petscts.h"   I*/
6818efac9SLisandro Dalcin 
7818efac9SLisandro Dalcin static PetscBool  cited      = PETSC_FALSE;
89371c9d4SSatish Balay static const char citation[] = "@article{Chung1993,\n"
9818efac9SLisandro Dalcin                                "  title   = {A Time Integration Algorithm for Structural Dynamics with Improved Numerical Dissipation: The Generalized-$\\alpha$ Method},\n"
10818efac9SLisandro Dalcin                                "  author  = {J. Chung, G. M. Hubert},\n"
11818efac9SLisandro Dalcin                                "  journal = {ASME Journal of Applied Mechanics},\n"
12818efac9SLisandro Dalcin                                "  volume  = {60},\n"
13818efac9SLisandro Dalcin                                "  number  = {2},\n"
14818efac9SLisandro Dalcin                                "  pages   = {371--375},\n"
15818efac9SLisandro Dalcin                                "  year    = {1993},\n"
16818efac9SLisandro Dalcin                                "  issn    = {0021-8936},\n"
17818efac9SLisandro Dalcin                                "  doi     = {http://dx.doi.org/10.1115/1.2900803}\n}\n";
18818efac9SLisandro Dalcin 
19818efac9SLisandro Dalcin typedef struct {
20818efac9SLisandro Dalcin   PetscReal stage_time;
21818efac9SLisandro Dalcin   PetscReal shift_V;
22818efac9SLisandro Dalcin   PetscReal shift_A;
23818efac9SLisandro Dalcin   PetscReal scale_F;
24818efac9SLisandro Dalcin   Vec       X0, Xa, X1;
25818efac9SLisandro Dalcin   Vec       V0, Va, V1;
26818efac9SLisandro Dalcin   Vec       A0, Aa, A1;
27818efac9SLisandro Dalcin 
28818efac9SLisandro Dalcin   Vec vec_dot;
29818efac9SLisandro Dalcin 
30818efac9SLisandro Dalcin   PetscReal Alpha_m;
31818efac9SLisandro Dalcin   PetscReal Alpha_f;
32818efac9SLisandro Dalcin   PetscReal Gamma;
33818efac9SLisandro Dalcin   PetscReal Beta;
34818efac9SLisandro Dalcin   PetscInt  order;
35818efac9SLisandro Dalcin 
36818efac9SLisandro Dalcin   Vec vec_sol_prev;
37818efac9SLisandro Dalcin   Vec vec_dot_prev;
38818efac9SLisandro Dalcin   Vec vec_lte_work[2];
39818efac9SLisandro Dalcin 
40818efac9SLisandro Dalcin   TSStepStatus status;
41818efac9SLisandro Dalcin } TS_Alpha;
42818efac9SLisandro Dalcin 
43d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_StageTime(TS ts)
44d71ae5a4SJacob Faibussowitsch {
45818efac9SLisandro Dalcin   TS_Alpha *th      = (TS_Alpha *)ts->data;
46818efac9SLisandro Dalcin   PetscReal t       = ts->ptime;
47818efac9SLisandro Dalcin   PetscReal dt      = ts->time_step;
48818efac9SLisandro Dalcin   PetscReal Alpha_m = th->Alpha_m;
49818efac9SLisandro Dalcin   PetscReal Alpha_f = th->Alpha_f;
50818efac9SLisandro Dalcin   PetscReal Gamma   = th->Gamma;
51818efac9SLisandro Dalcin   PetscReal Beta    = th->Beta;
52818efac9SLisandro Dalcin 
53818efac9SLisandro Dalcin   PetscFunctionBegin;
54818efac9SLisandro Dalcin   th->stage_time = t + Alpha_f * dt;
55818efac9SLisandro Dalcin   th->shift_V    = Gamma / (dt * Beta);
56818efac9SLisandro Dalcin   th->shift_A    = Alpha_m / (Alpha_f * dt * dt * Beta);
57818efac9SLisandro Dalcin   th->scale_F    = 1 / Alpha_f;
583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
59818efac9SLisandro Dalcin }
60818efac9SLisandro Dalcin 
61d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_StageVecs(TS ts, Vec X)
62d71ae5a4SJacob Faibussowitsch {
63818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
64818efac9SLisandro Dalcin   Vec       X1 = X, V1 = th->V1, A1 = th->A1;
65818efac9SLisandro Dalcin   Vec       Xa = th->Xa, Va = th->Va, Aa = th->Aa;
66818efac9SLisandro Dalcin   Vec       X0 = th->X0, V0 = th->V0, A0 = th->A0;
67818efac9SLisandro Dalcin   PetscReal dt      = ts->time_step;
68818efac9SLisandro Dalcin   PetscReal Alpha_m = th->Alpha_m;
69818efac9SLisandro Dalcin   PetscReal Alpha_f = th->Alpha_f;
70818efac9SLisandro Dalcin   PetscReal Gamma   = th->Gamma;
71818efac9SLisandro Dalcin   PetscReal Beta    = th->Beta;
72818efac9SLisandro Dalcin 
73818efac9SLisandro Dalcin   PetscFunctionBegin;
74818efac9SLisandro Dalcin   /* A1 = ... */
759566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(A1, -1.0, X0, X1));
769566063dSJacob Faibussowitsch   PetscCall(VecAXPY(A1, -dt, V0));
779566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(A1, -(1 - 2 * Beta) / (2 * Beta), 1 / (dt * dt * Beta), A0));
78818efac9SLisandro Dalcin   /* V1 = ... */
799566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(V1, (1.0 - Gamma) / Gamma, A0, A1));
809566063dSJacob Faibussowitsch   PetscCall(VecAYPX(V1, dt * Gamma, V0));
81818efac9SLisandro Dalcin   /* Xa = X0 + Alpha_f*(X1-X0) */
829566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Xa, -1.0, X0, X1));
839566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Xa, Alpha_f, X0));
84818efac9SLisandro Dalcin   /* Va = V0 + Alpha_f*(V1-V0) */
859566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Va, -1.0, V0, V1));
869566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Va, Alpha_f, V0));
87818efac9SLisandro Dalcin   /* Aa = A0 + Alpha_m*(A1-A0) */
889566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Aa, -1.0, A0, A1));
899566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Aa, Alpha_m, A0));
903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
91818efac9SLisandro Dalcin }
92818efac9SLisandro Dalcin 
93d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_SNESSolve(TS ts, Vec b, Vec x)
94d71ae5a4SJacob Faibussowitsch {
95818efac9SLisandro Dalcin   PetscInt nits, lits;
96818efac9SLisandro Dalcin 
97818efac9SLisandro Dalcin   PetscFunctionBegin;
989566063dSJacob Faibussowitsch   PetscCall(SNESSolve(ts->snes, b, x));
999566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(ts->snes, &nits));
1009566063dSJacob Faibussowitsch   PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits));
1019371c9d4SSatish Balay   ts->snes_its += nits;
1029371c9d4SSatish Balay   ts->ksp_its += lits;
1033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
104818efac9SLisandro Dalcin }
105818efac9SLisandro Dalcin 
106818efac9SLisandro Dalcin /*
107818efac9SLisandro Dalcin   Compute a consistent initial state for the generalized-alpha method.
108818efac9SLisandro Dalcin   - Solve two successive backward Euler steps with halved time step.
109818efac9SLisandro Dalcin   - Compute the initial second time derivative using backward differences.
110818efac9SLisandro Dalcin   - If using adaptivity, estimate the LTE of the initial step.
111818efac9SLisandro Dalcin */
112d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha_Restart(TS ts, PetscBool *initok)
113d71ae5a4SJacob Faibussowitsch {
114818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
115818efac9SLisandro Dalcin   PetscReal time_step;
116818efac9SLisandro Dalcin   PetscReal alpha_m, alpha_f, gamma, beta;
117818efac9SLisandro Dalcin   Vec       X0 = ts->vec_sol, X1, X2 = th->X1;
118818efac9SLisandro Dalcin   Vec       V0 = ts->vec_dot, V1, V2 = th->V1;
119818efac9SLisandro Dalcin   PetscBool stageok;
120818efac9SLisandro Dalcin 
121818efac9SLisandro Dalcin   PetscFunctionBegin;
1229566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X0, &X1));
1239566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(V0, &V1));
124818efac9SLisandro Dalcin 
125818efac9SLisandro Dalcin   /* Setup backward Euler with halved time step */
1269566063dSJacob Faibussowitsch   PetscCall(TSAlpha2GetParams(ts, &alpha_m, &alpha_f, &gamma, &beta));
1279566063dSJacob Faibussowitsch   PetscCall(TSAlpha2SetParams(ts, 1, 1, 1, 0.5));
1289566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &time_step));
129818efac9SLisandro Dalcin   ts->time_step = time_step / 2;
1309566063dSJacob Faibussowitsch   PetscCall(TSAlpha_StageTime(ts));
131818efac9SLisandro Dalcin   th->stage_time = ts->ptime;
1329566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(th->A0));
133818efac9SLisandro Dalcin 
134818efac9SLisandro Dalcin   /* First BE step, (t0,X0,V0) -> (t1,X1,V1) */
135818efac9SLisandro Dalcin   th->stage_time += ts->time_step;
1369566063dSJacob Faibussowitsch   PetscCall(VecCopy(X0, th->X0));
1379566063dSJacob Faibussowitsch   PetscCall(VecCopy(V0, th->V0));
1389566063dSJacob Faibussowitsch   PetscCall(TSPreStage(ts, th->stage_time));
1399566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->X0, X1));
1409566063dSJacob Faibussowitsch   PetscCall(TSAlpha_SNESSolve(ts, NULL, X1));
1419566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->V1, V1));
1429566063dSJacob Faibussowitsch   PetscCall(TSPostStage(ts, th->stage_time, 0, &X1));
1439566063dSJacob Faibussowitsch   PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, X1, &stageok));
144818efac9SLisandro Dalcin   if (!stageok) goto finally;
145818efac9SLisandro Dalcin 
146818efac9SLisandro Dalcin   /* Second BE step, (t1,X1,V1) -> (t2,X2,V2) */
147818efac9SLisandro Dalcin   th->stage_time += ts->time_step;
1489566063dSJacob Faibussowitsch   PetscCall(VecCopy(X1, th->X0));
1499566063dSJacob Faibussowitsch   PetscCall(VecCopy(V1, th->V0));
1509566063dSJacob Faibussowitsch   PetscCall(TSPreStage(ts, th->stage_time));
1519566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->X0, X2));
1529566063dSJacob Faibussowitsch   PetscCall(TSAlpha_SNESSolve(ts, NULL, X2));
1539566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->V1, V2));
1549566063dSJacob Faibussowitsch   PetscCall(TSPostStage(ts, th->stage_time, 0, &X2));
1559566063dSJacob Faibussowitsch   PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, X1, &stageok));
156818efac9SLisandro Dalcin   if (!stageok) goto finally;
157818efac9SLisandro Dalcin 
158818efac9SLisandro Dalcin   /* Compute A0 ~ dV/dt at t0 with backward differences */
1599566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(th->A0));
1609566063dSJacob Faibussowitsch   PetscCall(VecAXPY(th->A0, -3 / ts->time_step, V0));
1619566063dSJacob Faibussowitsch   PetscCall(VecAXPY(th->A0, +4 / ts->time_step, V1));
1629566063dSJacob Faibussowitsch   PetscCall(VecAXPY(th->A0, -1 / ts->time_step, V2));
163818efac9SLisandro Dalcin 
164818efac9SLisandro Dalcin   /* Rough, lower-order estimate LTE of the initial step */
1652ffb9264SLisandro Dalcin   if (th->vec_lte_work[0]) {
1669566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(th->vec_lte_work[0]));
1679566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work[0], +2, X2));
1689566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work[0], -4, X1));
1699566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work[0], +2, X0));
170818efac9SLisandro Dalcin   }
1712ffb9264SLisandro Dalcin   if (th->vec_lte_work[1]) {
1729566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(th->vec_lte_work[1]));
1739566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work[1], +2, V2));
1749566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work[1], -4, V1));
1759566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work[1], +2, V0));
176818efac9SLisandro Dalcin   }
177818efac9SLisandro Dalcin 
178818efac9SLisandro Dalcin finally:
179818efac9SLisandro Dalcin   /* Revert TSAlpha to the initial state (t0,X0,V0) */
180818efac9SLisandro Dalcin   if (initok) *initok = stageok;
1819566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, time_step));
1829566063dSJacob Faibussowitsch   PetscCall(TSAlpha2SetParams(ts, alpha_m, alpha_f, gamma, beta));
1839566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol, th->X0));
1849566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_dot, th->V0));
185818efac9SLisandro Dalcin 
1869566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X1));
1879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&V1));
1883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
189818efac9SLisandro Dalcin }
190818efac9SLisandro Dalcin 
191d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_Alpha(TS ts)
192d71ae5a4SJacob Faibussowitsch {
193818efac9SLisandro Dalcin   TS_Alpha *th         = (TS_Alpha *)ts->data;
194818efac9SLisandro Dalcin   PetscInt  rejections = 0;
195818efac9SLisandro Dalcin   PetscBool stageok, accept = PETSC_TRUE;
196818efac9SLisandro Dalcin   PetscReal next_time_step = ts->time_step;
197818efac9SLisandro Dalcin 
198818efac9SLisandro Dalcin   PetscFunctionBegin;
1999566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(citation, &cited));
200818efac9SLisandro Dalcin 
201818efac9SLisandro Dalcin   if (!ts->steprollback) {
2029566063dSJacob Faibussowitsch     if (th->vec_sol_prev) PetscCall(VecCopy(th->X0, th->vec_sol_prev));
2039566063dSJacob Faibussowitsch     if (th->vec_dot_prev) PetscCall(VecCopy(th->V0, th->vec_dot_prev));
2049566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_sol, th->X0));
2059566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_dot, th->V0));
2069566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->A1, th->A0));
207818efac9SLisandro Dalcin   }
208818efac9SLisandro Dalcin 
209818efac9SLisandro Dalcin   th->status = TS_STEP_INCOMPLETE;
210818efac9SLisandro Dalcin   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
211818efac9SLisandro Dalcin     if (ts->steprestart) {
2129566063dSJacob Faibussowitsch       PetscCall(TSAlpha_Restart(ts, &stageok));
213818efac9SLisandro Dalcin       if (!stageok) goto reject_step;
214818efac9SLisandro Dalcin     }
215818efac9SLisandro Dalcin 
2169566063dSJacob Faibussowitsch     PetscCall(TSAlpha_StageTime(ts));
2179566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->X0, th->X1));
2189566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, th->stage_time));
2199566063dSJacob Faibussowitsch     PetscCall(TSAlpha_SNESSolve(ts, NULL, th->X1));
2209566063dSJacob Faibussowitsch     PetscCall(TSPostStage(ts, th->stage_time, 0, &th->Xa));
2219566063dSJacob Faibussowitsch     PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->Xa, &stageok));
222818efac9SLisandro Dalcin     if (!stageok) goto reject_step;
223818efac9SLisandro Dalcin 
224818efac9SLisandro Dalcin     th->status = TS_STEP_PENDING;
2259566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->X1, ts->vec_sol));
2269566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->V1, ts->vec_dot));
2279566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
228818efac9SLisandro Dalcin     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
229818efac9SLisandro Dalcin     if (!accept) {
2309566063dSJacob Faibussowitsch       PetscCall(VecCopy(th->X0, ts->vec_sol));
2319566063dSJacob Faibussowitsch       PetscCall(VecCopy(th->V0, ts->vec_dot));
232818efac9SLisandro Dalcin       ts->time_step = next_time_step;
233818efac9SLisandro Dalcin       goto reject_step;
234818efac9SLisandro Dalcin     }
235818efac9SLisandro Dalcin 
236818efac9SLisandro Dalcin     ts->ptime += ts->time_step;
237818efac9SLisandro Dalcin     ts->time_step = next_time_step;
238818efac9SLisandro Dalcin     break;
239818efac9SLisandro Dalcin 
240818efac9SLisandro Dalcin   reject_step:
2419371c9d4SSatish Balay     ts->reject++;
2429371c9d4SSatish Balay     accept = PETSC_FALSE;
243818efac9SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
244818efac9SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
24563a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
246818efac9SLisandro Dalcin     }
247818efac9SLisandro Dalcin   }
2483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
249818efac9SLisandro Dalcin }
250818efac9SLisandro Dalcin 
251d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateWLTE_Alpha(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte)
252d71ae5a4SJacob Faibussowitsch {
253818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
254818efac9SLisandro Dalcin   Vec       X  = th->X1;              /* X = solution */
255818efac9SLisandro Dalcin   Vec       V  = th->V1;              /* V = solution */
256818efac9SLisandro Dalcin   Vec       Y  = th->vec_lte_work[0]; /* Y = X + LTE  */
257818efac9SLisandro Dalcin   Vec       Z  = th->vec_lte_work[1]; /* Z = V + LTE  */
2587453f775SEmil Constantinescu   PetscReal enormX, enormV, enormXa, enormVa, enormXr, enormVr;
259818efac9SLisandro Dalcin 
260818efac9SLisandro Dalcin   PetscFunctionBegin;
2619371c9d4SSatish Balay   if (!th->vec_sol_prev) {
2629371c9d4SSatish Balay     *wlte = -1;
2633ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2649371c9d4SSatish Balay   }
2659371c9d4SSatish Balay   if (!th->vec_dot_prev) {
2669371c9d4SSatish Balay     *wlte = -1;
2673ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2689371c9d4SSatish Balay   }
2699371c9d4SSatish Balay   if (!th->vec_lte_work[0]) {
2709371c9d4SSatish Balay     *wlte = -1;
2713ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2729371c9d4SSatish Balay   }
2739371c9d4SSatish Balay   if (!th->vec_lte_work[1]) {
2749371c9d4SSatish Balay     *wlte = -1;
2753ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2769371c9d4SSatish Balay   }
277818efac9SLisandro Dalcin   if (ts->steprestart) {
2782ffb9264SLisandro Dalcin     /* th->vec_lte_prev is set to the LTE in TSAlpha_Restart() */
2799566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Y, 1, X));
2809566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Z, 1, V));
281818efac9SLisandro Dalcin   } else {
282818efac9SLisandro Dalcin     /* Compute LTE using backward differences with non-constant time step */
283818efac9SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
284818efac9SLisandro Dalcin     PetscReal   a = 1 + h_prev / h;
2859371c9d4SSatish Balay     PetscScalar scal[3];
2869371c9d4SSatish Balay     Vec         vecX[3], vecV[3];
2879371c9d4SSatish Balay     scal[0] = +1 / a;
2889371c9d4SSatish Balay     scal[1] = -1 / (a - 1);
2899371c9d4SSatish Balay     scal[2] = +1 / (a * (a - 1));
2909371c9d4SSatish Balay     vecX[0] = th->X1;
2919371c9d4SSatish Balay     vecX[1] = th->X0;
2929371c9d4SSatish Balay     vecX[2] = th->vec_sol_prev;
2939371c9d4SSatish Balay     vecV[0] = th->V1;
2949371c9d4SSatish Balay     vecV[1] = th->V0;
2959371c9d4SSatish Balay     vecV[2] = th->vec_dot_prev;
2969566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Y));
2979566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(Y, 3, scal, vecX));
2989566063dSJacob Faibussowitsch     PetscCall(VecCopy(V, Z));
2999566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(Z, 3, scal, vecV));
300818efac9SLisandro Dalcin   }
301818efac9SLisandro Dalcin   /* XXX ts->atol and ts->vatol are not appropriate for computing enormV */
3029566063dSJacob Faibussowitsch   PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, &enormX, &enormXa, &enormXr));
3039566063dSJacob Faibussowitsch   PetscCall(TSErrorWeightedNorm(ts, V, Z, wnormtype, &enormV, &enormVa, &enormVr));
3049371c9d4SSatish Balay   if (wnormtype == NORM_2) *wlte = PetscSqrtReal(PetscSqr(enormX) / 2 + PetscSqr(enormV) / 2);
3059371c9d4SSatish Balay   else *wlte = PetscMax(enormX, enormV);
306818efac9SLisandro Dalcin   if (order) *order = 2;
3073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
308818efac9SLisandro Dalcin }
309818efac9SLisandro Dalcin 
310d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_Alpha(TS ts)
311d71ae5a4SJacob Faibussowitsch {
312818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
313818efac9SLisandro Dalcin 
314818efac9SLisandro Dalcin   PetscFunctionBegin;
3159566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->X0, ts->vec_sol));
3169566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->V0, ts->vec_dot));
3173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
318818efac9SLisandro Dalcin }
319818efac9SLisandro Dalcin 
320818efac9SLisandro Dalcin /*
321818efac9SLisandro Dalcin static PetscErrorCode TSInterpolate_Alpha(TS ts,PetscReal t,Vec X,Vec V)
322818efac9SLisandro Dalcin {
323818efac9SLisandro Dalcin   TS_Alpha       *th = (TS_Alpha*)ts->data;
324818efac9SLisandro Dalcin   PetscReal      dt  = t - ts->ptime;
325818efac9SLisandro Dalcin 
326818efac9SLisandro Dalcin   PetscFunctionBegin;
3279566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_dot,V));
3289566063dSJacob Faibussowitsch   PetscCall(VecAXPY(V,dt*(1-th->Gamma),th->A0));
3299566063dSJacob Faibussowitsch   PetscCall(VecAXPY(V,dt*th->Gamma,th->A1));
3309566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol,X));
3319566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X,dt,V));
3329566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X,dt*dt*((PetscReal)0.5-th->Beta),th->A0));
3339566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X,dt*dt*th->Beta,th->A1));
3343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
335818efac9SLisandro Dalcin }
336818efac9SLisandro Dalcin */
337818efac9SLisandro Dalcin 
338d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes, Vec X, Vec F, TS ts)
339d71ae5a4SJacob Faibussowitsch {
340818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
341818efac9SLisandro Dalcin   PetscReal ta = th->stage_time;
342818efac9SLisandro Dalcin   Vec       Xa = th->Xa, Va = th->Va, Aa = th->Aa;
343818efac9SLisandro Dalcin 
344818efac9SLisandro Dalcin   PetscFunctionBegin;
3459566063dSJacob Faibussowitsch   PetscCall(TSAlpha_StageVecs(ts, X));
346818efac9SLisandro Dalcin   /* F = Function(ta,Xa,Va,Aa) */
3479566063dSJacob Faibussowitsch   PetscCall(TSComputeI2Function(ts, ta, Xa, Va, Aa, F));
3489566063dSJacob Faibussowitsch   PetscCall(VecScale(F, th->scale_F));
3493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
350818efac9SLisandro Dalcin }
351818efac9SLisandro Dalcin 
352d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes, PETSC_UNUSED Vec X, Mat J, Mat P, TS ts)
353d71ae5a4SJacob Faibussowitsch {
354818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
355818efac9SLisandro Dalcin   PetscReal ta = th->stage_time;
356818efac9SLisandro Dalcin   Vec       Xa = th->Xa, Va = th->Va, Aa = th->Aa;
357818efac9SLisandro Dalcin   PetscReal dVdX = th->shift_V, dAdX = th->shift_A;
358818efac9SLisandro Dalcin 
359818efac9SLisandro Dalcin   PetscFunctionBegin;
360818efac9SLisandro Dalcin   /* J,P = Jacobian(ta,Xa,Va,Aa) */
3619566063dSJacob Faibussowitsch   PetscCall(TSComputeI2Jacobian(ts, ta, Xa, Va, Aa, dVdX, dAdX, J, P));
3623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
363818efac9SLisandro Dalcin }
364818efac9SLisandro Dalcin 
365d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_Alpha(TS ts)
366d71ae5a4SJacob Faibussowitsch {
367818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
368818efac9SLisandro Dalcin 
369818efac9SLisandro Dalcin   PetscFunctionBegin;
3709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X0));
3719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Xa));
3729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X1));
3739566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->V0));
3749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Va));
3759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->V1));
3769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->A0));
3779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Aa));
3789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->A1));
3799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_sol_prev));
3809566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_dot_prev));
3819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_lte_work[0]));
3829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_lte_work[1]));
3833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
384818efac9SLisandro Dalcin }
385818efac9SLisandro Dalcin 
386d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_Alpha(TS ts)
387d71ae5a4SJacob Faibussowitsch {
388818efac9SLisandro Dalcin   PetscFunctionBegin;
3899566063dSJacob Faibussowitsch   PetscCall(TSReset_Alpha(ts));
3909566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
391818efac9SLisandro Dalcin 
3929566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetRadius_C", NULL));
3939566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetParams_C", NULL));
3949566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2GetParams_C", NULL));
3953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
396818efac9SLisandro Dalcin }
397818efac9SLisandro Dalcin 
398d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_Alpha(TS ts)
399d71ae5a4SJacob Faibussowitsch {
400818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
4012ffb9264SLisandro Dalcin   PetscBool match;
402818efac9SLisandro Dalcin 
403818efac9SLisandro Dalcin   PetscFunctionBegin;
4049566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->X0));
4059566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->Xa));
4069566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->X1));
4079566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->V0));
4089566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->Va));
4099566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->V1));
4109566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->A0));
4119566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->Aa));
4129566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->A1));
413818efac9SLisandro Dalcin 
4149566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &ts->adapt));
4159566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt));
4169566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match));
4172ffb9264SLisandro Dalcin   if (!match) {
4189566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev));
4199566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_dot_prev));
4209566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work[0]));
4219566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work[1]));
422818efac9SLisandro Dalcin   }
423818efac9SLisandro Dalcin 
4249566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &ts->snes));
4253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
426818efac9SLisandro Dalcin }
427818efac9SLisandro Dalcin 
428d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_Alpha(TS ts, PetscOptionItems *PetscOptionsObject)
429d71ae5a4SJacob Faibussowitsch {
430818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
431818efac9SLisandro Dalcin 
432818efac9SLisandro Dalcin   PetscFunctionBegin;
433d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Generalized-Alpha ODE solver options");
434818efac9SLisandro Dalcin   {
435818efac9SLisandro Dalcin     PetscBool flg;
436818efac9SLisandro Dalcin     PetscReal radius = 1;
4379566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_radius", "Spectral radius (high-frequency dissipation)", "TSAlpha2SetRadius", radius, &radius, &flg));
4389566063dSJacob Faibussowitsch     if (flg) PetscCall(TSAlpha2SetRadius(ts, radius));
4399566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_alpha_m", "Algorithmic parameter alpha_m", "TSAlpha2SetParams", th->Alpha_m, &th->Alpha_m, NULL));
4409566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_alpha_f", "Algorithmic parameter alpha_f", "TSAlpha2SetParams", th->Alpha_f, &th->Alpha_f, NULL));
4419566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_gamma", "Algorithmic parameter gamma", "TSAlpha2SetParams", th->Gamma, &th->Gamma, NULL));
4429566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_beta", "Algorithmic parameter beta", "TSAlpha2SetParams", th->Beta, &th->Beta, NULL));
4439566063dSJacob Faibussowitsch     PetscCall(TSAlpha2SetParams(ts, th->Alpha_m, th->Alpha_f, th->Gamma, th->Beta));
444818efac9SLisandro Dalcin   }
445d0609cedSBarry Smith   PetscOptionsHeadEnd();
4463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
447818efac9SLisandro Dalcin }
448818efac9SLisandro Dalcin 
449d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_Alpha(TS ts, PetscViewer viewer)
450d71ae5a4SJacob Faibussowitsch {
451818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
452818efac9SLisandro Dalcin   PetscBool iascii;
453818efac9SLisandro Dalcin 
454818efac9SLisandro Dalcin   PetscFunctionBegin;
4559566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
45648a46eb9SPierre Jolivet   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));
4573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
458818efac9SLisandro Dalcin }
459818efac9SLisandro Dalcin 
460d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha2SetRadius_Alpha(TS ts, PetscReal radius)
461d71ae5a4SJacob Faibussowitsch {
462818efac9SLisandro Dalcin   PetscReal alpha_m, alpha_f, gamma, beta;
463818efac9SLisandro Dalcin 
464818efac9SLisandro Dalcin   PetscFunctionBegin;
465cad9d221SBarry Smith   PetscCheck(radius >= 0 && radius <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius);
466818efac9SLisandro Dalcin   alpha_m = (2 - radius) / (1 + radius);
467818efac9SLisandro Dalcin   alpha_f = 1 / (1 + radius);
468818efac9SLisandro Dalcin   gamma   = (PetscReal)0.5 + alpha_m - alpha_f;
4699371c9d4SSatish Balay   beta    = (PetscReal)0.5 * (1 + alpha_m - alpha_f);
4709371c9d4SSatish Balay   beta *= beta;
4719566063dSJacob Faibussowitsch   PetscCall(TSAlpha2SetParams(ts, alpha_m, alpha_f, gamma, beta));
4723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
473818efac9SLisandro Dalcin }
474818efac9SLisandro Dalcin 
475d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha2SetParams_Alpha(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma, PetscReal beta)
476d71ae5a4SJacob Faibussowitsch {
477818efac9SLisandro Dalcin   TS_Alpha *th  = (TS_Alpha *)ts->data;
478818efac9SLisandro Dalcin   PetscReal tol = 100 * PETSC_MACHINE_EPSILON;
479818efac9SLisandro Dalcin   PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma;
480818efac9SLisandro Dalcin 
481818efac9SLisandro Dalcin   PetscFunctionBegin;
482818efac9SLisandro Dalcin   th->Alpha_m = alpha_m;
483818efac9SLisandro Dalcin   th->Alpha_f = alpha_f;
484818efac9SLisandro Dalcin   th->Gamma   = gamma;
485818efac9SLisandro Dalcin   th->Beta    = beta;
486818efac9SLisandro Dalcin   th->order   = (PetscAbsReal(res) < tol) ? 2 : 1;
4873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
488818efac9SLisandro Dalcin }
489818efac9SLisandro Dalcin 
490d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAlpha2GetParams_Alpha(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma, PetscReal *beta)
491d71ae5a4SJacob Faibussowitsch {
492818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
493818efac9SLisandro Dalcin 
494818efac9SLisandro Dalcin   PetscFunctionBegin;
495818efac9SLisandro Dalcin   if (alpha_m) *alpha_m = th->Alpha_m;
496818efac9SLisandro Dalcin   if (alpha_f) *alpha_f = th->Alpha_f;
497818efac9SLisandro Dalcin   if (gamma) *gamma = th->Gamma;
498818efac9SLisandro Dalcin   if (beta) *beta = th->Beta;
4993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
500818efac9SLisandro Dalcin }
501818efac9SLisandro Dalcin 
502818efac9SLisandro Dalcin /*MC
50314d0ab18SJacob Faibussowitsch   TSALPHA2 - ODE/DAE solver using the implicit Generalized-Alpha method for second-order systems
504818efac9SLisandro Dalcin 
505818efac9SLisandro Dalcin   Level: beginner
506818efac9SLisandro Dalcin 
507818efac9SLisandro Dalcin   References:
508606c0280SSatish Balay . * - J. Chung, G.M.Hubert. "A Time Integration Algorithm for Structural
509818efac9SLisandro Dalcin   Dynamics with Improved Numerical Dissipation: The Generalized-alpha
510818efac9SLisandro Dalcin   Method" ASME Journal of Applied Mechanics, 60, 371:375, 1993.
511818efac9SLisandro Dalcin 
5121cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSCreate()`, `TSSetType()`, `TSAlpha2SetRadius()`, `TSAlpha2SetParams()`
513818efac9SLisandro Dalcin M*/
514d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_Alpha2(TS ts)
515d71ae5a4SJacob Faibussowitsch {
516818efac9SLisandro Dalcin   TS_Alpha *th;
517818efac9SLisandro Dalcin 
518818efac9SLisandro Dalcin   PetscFunctionBegin;
519818efac9SLisandro Dalcin   ts->ops->reset          = TSReset_Alpha;
520818efac9SLisandro Dalcin   ts->ops->destroy        = TSDestroy_Alpha;
521818efac9SLisandro Dalcin   ts->ops->view           = TSView_Alpha;
522818efac9SLisandro Dalcin   ts->ops->setup          = TSSetUp_Alpha;
523818efac9SLisandro Dalcin   ts->ops->setfromoptions = TSSetFromOptions_Alpha;
524818efac9SLisandro Dalcin   ts->ops->step           = TSStep_Alpha;
525818efac9SLisandro Dalcin   ts->ops->evaluatewlte   = TSEvaluateWLTE_Alpha;
526818efac9SLisandro Dalcin   ts->ops->rollback       = TSRollBack_Alpha;
527818efac9SLisandro Dalcin   /*ts->ops->interpolate  = TSInterpolate_Alpha;*/
528818efac9SLisandro Dalcin   ts->ops->snesfunction  = SNESTSFormFunction_Alpha;
529818efac9SLisandro Dalcin   ts->ops->snesjacobian  = SNESTSFormJacobian_Alpha;
5302ffb9264SLisandro Dalcin   ts->default_adapt_type = TSADAPTNONE;
531818efac9SLisandro Dalcin 
532efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
533efd4aadfSBarry Smith 
5344dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&th));
535818efac9SLisandro Dalcin   ts->data = (void *)th;
536818efac9SLisandro Dalcin 
537818efac9SLisandro Dalcin   th->Alpha_m = 0.5;
538818efac9SLisandro Dalcin   th->Alpha_f = 0.5;
539818efac9SLisandro Dalcin   th->Gamma   = 0.5;
540818efac9SLisandro Dalcin   th->Beta    = 0.25;
541818efac9SLisandro Dalcin   th->order   = 2;
542818efac9SLisandro Dalcin 
5439566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetRadius_C", TSAlpha2SetRadius_Alpha));
5449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetParams_C", TSAlpha2SetParams_Alpha));
5459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2GetParams_C", TSAlpha2GetParams_Alpha));
5463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
547818efac9SLisandro Dalcin }
548818efac9SLisandro Dalcin 
549818efac9SLisandro Dalcin /*@
550bcf0153eSBarry Smith   TSAlpha2SetRadius - sets the desired spectral radius of the method for `TSALPHA2`
551818efac9SLisandro Dalcin   (i.e. high-frequency numerical damping)
552818efac9SLisandro Dalcin 
553c3339decSBarry Smith   Logically Collective
554818efac9SLisandro Dalcin 
555d8d19677SJose E. Roman   Input Parameters:
556818efac9SLisandro Dalcin + ts     - timestepping context
557818efac9SLisandro Dalcin - radius - the desired spectral radius
558818efac9SLisandro Dalcin 
559bcf0153eSBarry Smith   Options Database Key:
56067b8a455SSatish Balay . -ts_alpha_radius <radius> - set the desired spectral radius
561818efac9SLisandro Dalcin 
562818efac9SLisandro Dalcin   Level: intermediate
563818efac9SLisandro Dalcin 
56414d0ab18SJacob Faibussowitsch   Notes:
56514d0ab18SJacob Faibussowitsch 
56614d0ab18SJacob Faibussowitsch   The algorithmic parameters $\alpha_m$ and $\alpha_f$ of the generalized-$\alpha$ method can
56714d0ab18SJacob Faibussowitsch   be computed in terms of a specified spectral radius $\rho$ in `[0, 1]` for infinite time step
56814d0ab18SJacob Faibussowitsch   in order to control high-frequency numerical damping\:
569*562efe2eSBarry Smith 
57014d0ab18SJacob Faibussowitsch   $$
571*562efe2eSBarry Smith   \begin{align*}
572*562efe2eSBarry Smith   \alpha_m = (2-\rho)/(1+\rho) \\
57314d0ab18SJacob Faibussowitsch   \alpha_f = 1/(1+\rho)
574*562efe2eSBarry Smith   \end{align*}
57514d0ab18SJacob Faibussowitsch   $$
57614d0ab18SJacob Faibussowitsch 
5771cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA2`, `TSAlpha2SetParams()`, `TSAlpha2GetParams()`
578818efac9SLisandro Dalcin @*/
579d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlpha2SetRadius(TS ts, PetscReal radius)
580d71ae5a4SJacob Faibussowitsch {
581818efac9SLisandro Dalcin   PetscFunctionBegin;
582818efac9SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
583818efac9SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, radius, 2);
584cad9d221SBarry Smith   PetscCheck(radius >= 0 && radius <= 1, ((PetscObject)ts)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius);
585cac4c232SBarry Smith   PetscTryMethod(ts, "TSAlpha2SetRadius_C", (TS, PetscReal), (ts, radius));
5863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
587818efac9SLisandro Dalcin }
588818efac9SLisandro Dalcin 
589818efac9SLisandro Dalcin /*@
590bcf0153eSBarry Smith   TSAlpha2SetParams - sets the algorithmic parameters for `TSALPHA2`
591818efac9SLisandro Dalcin 
592c3339decSBarry Smith   Logically Collective
593818efac9SLisandro Dalcin 
594d8d19677SJose E. Roman   Input Parameters:
595818efac9SLisandro Dalcin + ts      - timestepping context
5962fe279fdSBarry Smith . alpha_m - algorithmic parameter
5972fe279fdSBarry Smith . alpha_f - algorithmic parameter
5982fe279fdSBarry Smith . gamma   - algorithmic parameter
5992fe279fdSBarry Smith - beta    - algorithmic parameter
600818efac9SLisandro Dalcin 
601bcf0153eSBarry Smith   Options Database Keys:
60267b8a455SSatish Balay + -ts_alpha_alpha_m <alpha_m> - set alpha_m
60367b8a455SSatish Balay . -ts_alpha_alpha_f <alpha_f> - set alpha_f
60467b8a455SSatish Balay . -ts_alpha_gamma   <gamma>   - set gamma
60567b8a455SSatish Balay - -ts_alpha_beta    <beta>    - set beta
606818efac9SLisandro Dalcin 
607bcf0153eSBarry Smith   Level: advanced
608bcf0153eSBarry Smith 
60914d0ab18SJacob Faibussowitsch   Notes:
61014d0ab18SJacob Faibussowitsch   Second-order accuracy can be obtained so long as\:
611*562efe2eSBarry Smith 
61214d0ab18SJacob Faibussowitsch   $$
613*562efe2eSBarry Smith   \begin{align*}
614*562efe2eSBarry Smith   \gamma = 1/2 + \alpha_m - \alpha_f \\
615*562efe2eSBarry Smith   \beta  = 1/4 (1 + \alpha_m - \alpha_f)^2.
616*562efe2eSBarry Smith   \end{align*}
61714d0ab18SJacob Faibussowitsch   $$
61814d0ab18SJacob Faibussowitsch 
61914d0ab18SJacob Faibussowitsch   Unconditional stability requires\:
62014d0ab18SJacob Faibussowitsch   $$
621*562efe2eSBarry Smith   \alpha_m >= \alpha_f >= 1/2.
62214d0ab18SJacob Faibussowitsch   $$
62314d0ab18SJacob Faibussowitsch 
62414d0ab18SJacob Faibussowitsch   Use of this function is normally only required to hack `TSALPHA2` to use a modified
62514d0ab18SJacob Faibussowitsch   integration scheme. Users should call `TSAlpha2SetRadius()` to set the desired spectral
62614d0ab18SJacob Faibussowitsch   radius of the methods (i.e. high-frequency damping) in order so select optimal values for
627818efac9SLisandro Dalcin   these parameters.
628818efac9SLisandro Dalcin 
6291cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA2`, `TSAlpha2SetRadius()`, `TSAlpha2GetParams()`
630818efac9SLisandro Dalcin @*/
631d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlpha2SetParams(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma, PetscReal beta)
632d71ae5a4SJacob Faibussowitsch {
633818efac9SLisandro Dalcin   PetscFunctionBegin;
634818efac9SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
635818efac9SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, alpha_m, 2);
636818efac9SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, alpha_f, 3);
637818efac9SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, gamma, 4);
638818efac9SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, beta, 5);
639cac4c232SBarry Smith   PetscTryMethod(ts, "TSAlpha2SetParams_C", (TS, PetscReal, PetscReal, PetscReal, PetscReal), (ts, alpha_m, alpha_f, gamma, beta));
6403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
641818efac9SLisandro Dalcin }
642818efac9SLisandro Dalcin 
643818efac9SLisandro Dalcin /*@
644bcf0153eSBarry Smith   TSAlpha2GetParams - gets the algorithmic parameters for `TSALPHA2`
645818efac9SLisandro Dalcin 
646818efac9SLisandro Dalcin   Not Collective
647818efac9SLisandro Dalcin 
648818efac9SLisandro Dalcin   Input Parameter:
649818efac9SLisandro Dalcin . ts - timestepping context
650818efac9SLisandro Dalcin 
651818efac9SLisandro Dalcin   Output Parameters:
6522fe279fdSBarry Smith + alpha_m - algorithmic parameter
6532fe279fdSBarry Smith . alpha_f - algorithmic parameter
6542fe279fdSBarry Smith . gamma   - algorithmic parameter
6552fe279fdSBarry Smith - beta    - algorithmic parameter
656818efac9SLisandro Dalcin 
657bcf0153eSBarry Smith   Level: advanced
658bcf0153eSBarry Smith 
659818efac9SLisandro Dalcin   Note:
66014d0ab18SJacob Faibussowitsch   Use of this function is normally only required to hack `TSALPHA2` to use a modified
66114d0ab18SJacob Faibussowitsch   integration scheme. Users should call `TSAlpha2SetRadius()` to set the high-frequency damping
66214d0ab18SJacob Faibussowitsch   (i.e. spectral radius of the method) in order so select optimal values for these parameters.
663818efac9SLisandro Dalcin 
6641cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSALPHA2`, `TSAlpha2SetRadius()`, `TSAlpha2SetParams()`
665818efac9SLisandro Dalcin @*/
666d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAlpha2GetParams(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma, PetscReal *beta)
667d71ae5a4SJacob Faibussowitsch {
668818efac9SLisandro Dalcin   PetscFunctionBegin;
669818efac9SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
6704f572ea9SToby Isaac   if (alpha_m) PetscAssertPointer(alpha_m, 2);
6714f572ea9SToby Isaac   if (alpha_f) PetscAssertPointer(alpha_f, 3);
6724f572ea9SToby Isaac   if (gamma) PetscAssertPointer(gamma, 4);
6734f572ea9SToby Isaac   if (beta) PetscAssertPointer(beta, 5);
674cac4c232SBarry Smith   PetscUseMethod(ts, "TSAlpha2GetParams_C", (TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *), (ts, alpha_m, alpha_f, gamma, beta));
6753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
676818efac9SLisandro Dalcin }
677