xref: /petsc/src/ts/impls/implicit/alpha/alpha2.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
439371c9d4SSatish Balay static PetscErrorCode TSAlpha_StageTime(TS ts) {
44818efac9SLisandro Dalcin   TS_Alpha *th      = (TS_Alpha *)ts->data;
45818efac9SLisandro Dalcin   PetscReal t       = ts->ptime;
46818efac9SLisandro Dalcin   PetscReal dt      = ts->time_step;
47818efac9SLisandro Dalcin   PetscReal Alpha_m = th->Alpha_m;
48818efac9SLisandro Dalcin   PetscReal Alpha_f = th->Alpha_f;
49818efac9SLisandro Dalcin   PetscReal Gamma   = th->Gamma;
50818efac9SLisandro Dalcin   PetscReal Beta    = th->Beta;
51818efac9SLisandro Dalcin 
52818efac9SLisandro Dalcin   PetscFunctionBegin;
53818efac9SLisandro Dalcin   th->stage_time = t + Alpha_f * dt;
54818efac9SLisandro Dalcin   th->shift_V    = Gamma / (dt * Beta);
55818efac9SLisandro Dalcin   th->shift_A    = Alpha_m / (Alpha_f * dt * dt * Beta);
56818efac9SLisandro Dalcin   th->scale_F    = 1 / Alpha_f;
57818efac9SLisandro Dalcin   PetscFunctionReturn(0);
58818efac9SLisandro Dalcin }
59818efac9SLisandro Dalcin 
609371c9d4SSatish Balay static PetscErrorCode TSAlpha_StageVecs(TS ts, Vec X) {
61818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
62818efac9SLisandro Dalcin   Vec       X1 = X, V1 = th->V1, A1 = th->A1;
63818efac9SLisandro Dalcin   Vec       Xa = th->Xa, Va = th->Va, Aa = th->Aa;
64818efac9SLisandro Dalcin   Vec       X0 = th->X0, V0 = th->V0, A0 = th->A0;
65818efac9SLisandro Dalcin   PetscReal dt      = ts->time_step;
66818efac9SLisandro Dalcin   PetscReal Alpha_m = th->Alpha_m;
67818efac9SLisandro Dalcin   PetscReal Alpha_f = th->Alpha_f;
68818efac9SLisandro Dalcin   PetscReal Gamma   = th->Gamma;
69818efac9SLisandro Dalcin   PetscReal Beta    = th->Beta;
70818efac9SLisandro Dalcin 
71818efac9SLisandro Dalcin   PetscFunctionBegin;
72818efac9SLisandro Dalcin   /* A1 = ... */
739566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(A1, -1.0, X0, X1));
749566063dSJacob Faibussowitsch   PetscCall(VecAXPY(A1, -dt, V0));
759566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(A1, -(1 - 2 * Beta) / (2 * Beta), 1 / (dt * dt * Beta), A0));
76818efac9SLisandro Dalcin   /* V1 = ... */
779566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(V1, (1.0 - Gamma) / Gamma, A0, A1));
789566063dSJacob Faibussowitsch   PetscCall(VecAYPX(V1, dt * Gamma, V0));
79818efac9SLisandro Dalcin   /* Xa = X0 + Alpha_f*(X1-X0) */
809566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Xa, -1.0, X0, X1));
819566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Xa, Alpha_f, X0));
82818efac9SLisandro Dalcin   /* Va = V0 + Alpha_f*(V1-V0) */
839566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Va, -1.0, V0, V1));
849566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Va, Alpha_f, V0));
85818efac9SLisandro Dalcin   /* Aa = A0 + Alpha_m*(A1-A0) */
869566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Aa, -1.0, A0, A1));
879566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Aa, Alpha_m, A0));
88818efac9SLisandro Dalcin   PetscFunctionReturn(0);
89818efac9SLisandro Dalcin }
90818efac9SLisandro Dalcin 
919371c9d4SSatish Balay static PetscErrorCode TSAlpha_SNESSolve(TS ts, Vec b, Vec x) {
92818efac9SLisandro Dalcin   PetscInt nits, lits;
93818efac9SLisandro Dalcin 
94818efac9SLisandro Dalcin   PetscFunctionBegin;
959566063dSJacob Faibussowitsch   PetscCall(SNESSolve(ts->snes, b, x));
969566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(ts->snes, &nits));
979566063dSJacob Faibussowitsch   PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits));
989371c9d4SSatish Balay   ts->snes_its += nits;
999371c9d4SSatish Balay   ts->ksp_its += lits;
100818efac9SLisandro Dalcin   PetscFunctionReturn(0);
101818efac9SLisandro Dalcin }
102818efac9SLisandro Dalcin 
103818efac9SLisandro Dalcin /*
104818efac9SLisandro Dalcin   Compute a consistent initial state for the generalized-alpha method.
105818efac9SLisandro Dalcin   - Solve two successive backward Euler steps with halved time step.
106818efac9SLisandro Dalcin   - Compute the initial second time derivative using backward differences.
107818efac9SLisandro Dalcin   - If using adaptivity, estimate the LTE of the initial step.
108818efac9SLisandro Dalcin */
1099371c9d4SSatish Balay static PetscErrorCode TSAlpha_Restart(TS ts, PetscBool *initok) {
110818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
111818efac9SLisandro Dalcin   PetscReal time_step;
112818efac9SLisandro Dalcin   PetscReal alpha_m, alpha_f, gamma, beta;
113818efac9SLisandro Dalcin   Vec       X0 = ts->vec_sol, X1, X2 = th->X1;
114818efac9SLisandro Dalcin   Vec       V0 = ts->vec_dot, V1, V2 = th->V1;
115818efac9SLisandro Dalcin   PetscBool stageok;
116818efac9SLisandro Dalcin 
117818efac9SLisandro Dalcin   PetscFunctionBegin;
1189566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X0, &X1));
1199566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(V0, &V1));
120818efac9SLisandro Dalcin 
121818efac9SLisandro Dalcin   /* Setup backward Euler with halved time step */
1229566063dSJacob Faibussowitsch   PetscCall(TSAlpha2GetParams(ts, &alpha_m, &alpha_f, &gamma, &beta));
1239566063dSJacob Faibussowitsch   PetscCall(TSAlpha2SetParams(ts, 1, 1, 1, 0.5));
1249566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &time_step));
125818efac9SLisandro Dalcin   ts->time_step = time_step / 2;
1269566063dSJacob Faibussowitsch   PetscCall(TSAlpha_StageTime(ts));
127818efac9SLisandro Dalcin   th->stage_time = ts->ptime;
1289566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(th->A0));
129818efac9SLisandro Dalcin 
130818efac9SLisandro Dalcin   /* First BE step, (t0,X0,V0) -> (t1,X1,V1) */
131818efac9SLisandro Dalcin   th->stage_time += ts->time_step;
1329566063dSJacob Faibussowitsch   PetscCall(VecCopy(X0, th->X0));
1339566063dSJacob Faibussowitsch   PetscCall(VecCopy(V0, th->V0));
1349566063dSJacob Faibussowitsch   PetscCall(TSPreStage(ts, th->stage_time));
1359566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->X0, X1));
1369566063dSJacob Faibussowitsch   PetscCall(TSAlpha_SNESSolve(ts, NULL, X1));
1379566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->V1, V1));
1389566063dSJacob Faibussowitsch   PetscCall(TSPostStage(ts, th->stage_time, 0, &X1));
1399566063dSJacob Faibussowitsch   PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, X1, &stageok));
140818efac9SLisandro Dalcin   if (!stageok) goto finally;
141818efac9SLisandro Dalcin 
142818efac9SLisandro Dalcin   /* Second BE step, (t1,X1,V1) -> (t2,X2,V2) */
143818efac9SLisandro Dalcin   th->stage_time += ts->time_step;
1449566063dSJacob Faibussowitsch   PetscCall(VecCopy(X1, th->X0));
1459566063dSJacob Faibussowitsch   PetscCall(VecCopy(V1, th->V0));
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(VecCopy(th->V1, V2));
1509566063dSJacob Faibussowitsch   PetscCall(TSPostStage(ts, th->stage_time, 0, &X2));
1519566063dSJacob Faibussowitsch   PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, X1, &stageok));
152818efac9SLisandro Dalcin   if (!stageok) goto finally;
153818efac9SLisandro Dalcin 
154818efac9SLisandro Dalcin   /* Compute A0 ~ dV/dt at t0 with backward differences */
1559566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(th->A0));
1569566063dSJacob Faibussowitsch   PetscCall(VecAXPY(th->A0, -3 / ts->time_step, V0));
1579566063dSJacob Faibussowitsch   PetscCall(VecAXPY(th->A0, +4 / ts->time_step, V1));
1589566063dSJacob Faibussowitsch   PetscCall(VecAXPY(th->A0, -1 / ts->time_step, V2));
159818efac9SLisandro Dalcin 
160818efac9SLisandro Dalcin   /* Rough, lower-order estimate LTE of the initial step */
1612ffb9264SLisandro Dalcin   if (th->vec_lte_work[0]) {
1629566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(th->vec_lte_work[0]));
1639566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work[0], +2, X2));
1649566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work[0], -4, X1));
1659566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work[0], +2, X0));
166818efac9SLisandro Dalcin   }
1672ffb9264SLisandro Dalcin   if (th->vec_lte_work[1]) {
1689566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(th->vec_lte_work[1]));
1699566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work[1], +2, V2));
1709566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work[1], -4, V1));
1719566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work[1], +2, V0));
172818efac9SLisandro Dalcin   }
173818efac9SLisandro Dalcin 
174818efac9SLisandro Dalcin finally:
175818efac9SLisandro Dalcin   /* Revert TSAlpha to the initial state (t0,X0,V0) */
176818efac9SLisandro Dalcin   if (initok) *initok = stageok;
1779566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, time_step));
1789566063dSJacob Faibussowitsch   PetscCall(TSAlpha2SetParams(ts, alpha_m, alpha_f, gamma, beta));
1799566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol, th->X0));
1809566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_dot, th->V0));
181818efac9SLisandro Dalcin 
1829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X1));
1839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&V1));
184818efac9SLisandro Dalcin   PetscFunctionReturn(0);
185818efac9SLisandro Dalcin }
186818efac9SLisandro Dalcin 
1879371c9d4SSatish Balay static PetscErrorCode TSStep_Alpha(TS ts) {
188818efac9SLisandro Dalcin   TS_Alpha *th         = (TS_Alpha *)ts->data;
189818efac9SLisandro Dalcin   PetscInt  rejections = 0;
190818efac9SLisandro Dalcin   PetscBool stageok, accept = PETSC_TRUE;
191818efac9SLisandro Dalcin   PetscReal next_time_step = ts->time_step;
192818efac9SLisandro Dalcin 
193818efac9SLisandro Dalcin   PetscFunctionBegin;
1949566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(citation, &cited));
195818efac9SLisandro Dalcin 
196818efac9SLisandro Dalcin   if (!ts->steprollback) {
1979566063dSJacob Faibussowitsch     if (th->vec_sol_prev) PetscCall(VecCopy(th->X0, th->vec_sol_prev));
1989566063dSJacob Faibussowitsch     if (th->vec_dot_prev) PetscCall(VecCopy(th->V0, th->vec_dot_prev));
1999566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_sol, th->X0));
2009566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_dot, th->V0));
2019566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->A1, th->A0));
202818efac9SLisandro Dalcin   }
203818efac9SLisandro Dalcin 
204818efac9SLisandro Dalcin   th->status = TS_STEP_INCOMPLETE;
205818efac9SLisandro Dalcin   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
206818efac9SLisandro Dalcin     if (ts->steprestart) {
2079566063dSJacob Faibussowitsch       PetscCall(TSAlpha_Restart(ts, &stageok));
208818efac9SLisandro Dalcin       if (!stageok) goto reject_step;
209818efac9SLisandro Dalcin     }
210818efac9SLisandro Dalcin 
2119566063dSJacob Faibussowitsch     PetscCall(TSAlpha_StageTime(ts));
2129566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->X0, th->X1));
2139566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, th->stage_time));
2149566063dSJacob Faibussowitsch     PetscCall(TSAlpha_SNESSolve(ts, NULL, th->X1));
2159566063dSJacob Faibussowitsch     PetscCall(TSPostStage(ts, th->stage_time, 0, &th->Xa));
2169566063dSJacob Faibussowitsch     PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->Xa, &stageok));
217818efac9SLisandro Dalcin     if (!stageok) goto reject_step;
218818efac9SLisandro Dalcin 
219818efac9SLisandro Dalcin     th->status = TS_STEP_PENDING;
2209566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->X1, ts->vec_sol));
2219566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->V1, ts->vec_dot));
2229566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
223818efac9SLisandro Dalcin     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
224818efac9SLisandro Dalcin     if (!accept) {
2259566063dSJacob Faibussowitsch       PetscCall(VecCopy(th->X0, ts->vec_sol));
2269566063dSJacob Faibussowitsch       PetscCall(VecCopy(th->V0, ts->vec_dot));
227818efac9SLisandro Dalcin       ts->time_step = next_time_step;
228818efac9SLisandro Dalcin       goto reject_step;
229818efac9SLisandro Dalcin     }
230818efac9SLisandro Dalcin 
231818efac9SLisandro Dalcin     ts->ptime += ts->time_step;
232818efac9SLisandro Dalcin     ts->time_step = next_time_step;
233818efac9SLisandro Dalcin     break;
234818efac9SLisandro Dalcin 
235818efac9SLisandro Dalcin   reject_step:
2369371c9d4SSatish Balay     ts->reject++;
2379371c9d4SSatish Balay     accept = PETSC_FALSE;
238818efac9SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
239818efac9SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
24063a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
241818efac9SLisandro Dalcin     }
242818efac9SLisandro Dalcin   }
243818efac9SLisandro Dalcin   PetscFunctionReturn(0);
244818efac9SLisandro Dalcin }
245818efac9SLisandro Dalcin 
2469371c9d4SSatish Balay static PetscErrorCode TSEvaluateWLTE_Alpha(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte) {
247818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
248818efac9SLisandro Dalcin   Vec       X  = th->X1;              /* X = solution */
249818efac9SLisandro Dalcin   Vec       V  = th->V1;              /* V = solution */
250818efac9SLisandro Dalcin   Vec       Y  = th->vec_lte_work[0]; /* Y = X + LTE  */
251818efac9SLisandro Dalcin   Vec       Z  = th->vec_lte_work[1]; /* Z = V + LTE  */
2527453f775SEmil Constantinescu   PetscReal enormX, enormV, enormXa, enormVa, enormXr, enormVr;
253818efac9SLisandro Dalcin 
254818efac9SLisandro Dalcin   PetscFunctionBegin;
2559371c9d4SSatish Balay   if (!th->vec_sol_prev) {
2569371c9d4SSatish Balay     *wlte = -1;
2579371c9d4SSatish Balay     PetscFunctionReturn(0);
2589371c9d4SSatish Balay   }
2599371c9d4SSatish Balay   if (!th->vec_dot_prev) {
2609371c9d4SSatish Balay     *wlte = -1;
2619371c9d4SSatish Balay     PetscFunctionReturn(0);
2629371c9d4SSatish Balay   }
2639371c9d4SSatish Balay   if (!th->vec_lte_work[0]) {
2649371c9d4SSatish Balay     *wlte = -1;
2659371c9d4SSatish Balay     PetscFunctionReturn(0);
2669371c9d4SSatish Balay   }
2679371c9d4SSatish Balay   if (!th->vec_lte_work[1]) {
2689371c9d4SSatish Balay     *wlte = -1;
2699371c9d4SSatish Balay     PetscFunctionReturn(0);
2709371c9d4SSatish Balay   }
271818efac9SLisandro Dalcin   if (ts->steprestart) {
2722ffb9264SLisandro Dalcin     /* th->vec_lte_prev is set to the LTE in TSAlpha_Restart() */
2739566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Y, 1, X));
2749566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Z, 1, V));
275818efac9SLisandro Dalcin   } else {
276818efac9SLisandro Dalcin     /* Compute LTE using backward differences with non-constant time step */
277818efac9SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
278818efac9SLisandro Dalcin     PetscReal   a = 1 + h_prev / h;
2799371c9d4SSatish Balay     PetscScalar scal[3];
2809371c9d4SSatish Balay     Vec         vecX[3], vecV[3];
2819371c9d4SSatish Balay     scal[0] = +1 / a;
2829371c9d4SSatish Balay     scal[1] = -1 / (a - 1);
2839371c9d4SSatish Balay     scal[2] = +1 / (a * (a - 1));
2849371c9d4SSatish Balay     vecX[0] = th->X1;
2859371c9d4SSatish Balay     vecX[1] = th->X0;
2869371c9d4SSatish Balay     vecX[2] = th->vec_sol_prev;
2879371c9d4SSatish Balay     vecV[0] = th->V1;
2889371c9d4SSatish Balay     vecV[1] = th->V0;
2899371c9d4SSatish Balay     vecV[2] = th->vec_dot_prev;
2909566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Y));
2919566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(Y, 3, scal, vecX));
2929566063dSJacob Faibussowitsch     PetscCall(VecCopy(V, Z));
2939566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(Z, 3, scal, vecV));
294818efac9SLisandro Dalcin   }
295818efac9SLisandro Dalcin   /* XXX ts->atol and ts->vatol are not appropriate for computing enormV */
2969566063dSJacob Faibussowitsch   PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, &enormX, &enormXa, &enormXr));
2979566063dSJacob Faibussowitsch   PetscCall(TSErrorWeightedNorm(ts, V, Z, wnormtype, &enormV, &enormVa, &enormVr));
2989371c9d4SSatish Balay   if (wnormtype == NORM_2) *wlte = PetscSqrtReal(PetscSqr(enormX) / 2 + PetscSqr(enormV) / 2);
2999371c9d4SSatish Balay   else *wlte = PetscMax(enormX, enormV);
300818efac9SLisandro Dalcin   if (order) *order = 2;
301818efac9SLisandro Dalcin   PetscFunctionReturn(0);
302818efac9SLisandro Dalcin }
303818efac9SLisandro Dalcin 
3049371c9d4SSatish Balay static PetscErrorCode TSRollBack_Alpha(TS ts) {
305818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
306818efac9SLisandro Dalcin 
307818efac9SLisandro Dalcin   PetscFunctionBegin;
3089566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->X0, ts->vec_sol));
3099566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->V0, ts->vec_dot));
310818efac9SLisandro Dalcin   PetscFunctionReturn(0);
311818efac9SLisandro Dalcin }
312818efac9SLisandro Dalcin 
313818efac9SLisandro Dalcin /*
314818efac9SLisandro Dalcin static PetscErrorCode TSInterpolate_Alpha(TS ts,PetscReal t,Vec X,Vec V)
315818efac9SLisandro Dalcin {
316818efac9SLisandro Dalcin   TS_Alpha       *th = (TS_Alpha*)ts->data;
317818efac9SLisandro Dalcin   PetscReal      dt  = t - ts->ptime;
318818efac9SLisandro Dalcin 
319818efac9SLisandro Dalcin   PetscFunctionBegin;
3209566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_dot,V));
3219566063dSJacob Faibussowitsch   PetscCall(VecAXPY(V,dt*(1-th->Gamma),th->A0));
3229566063dSJacob Faibussowitsch   PetscCall(VecAXPY(V,dt*th->Gamma,th->A1));
3239566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol,X));
3249566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X,dt,V));
3259566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X,dt*dt*((PetscReal)0.5-th->Beta),th->A0));
3269566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X,dt*dt*th->Beta,th->A1));
327818efac9SLisandro Dalcin   PetscFunctionReturn(0);
328818efac9SLisandro Dalcin }
329818efac9SLisandro Dalcin */
330818efac9SLisandro Dalcin 
3319371c9d4SSatish Balay static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes, Vec X, Vec F, TS ts) {
332818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
333818efac9SLisandro Dalcin   PetscReal ta = th->stage_time;
334818efac9SLisandro Dalcin   Vec       Xa = th->Xa, Va = th->Va, Aa = th->Aa;
335818efac9SLisandro Dalcin 
336818efac9SLisandro Dalcin   PetscFunctionBegin;
3379566063dSJacob Faibussowitsch   PetscCall(TSAlpha_StageVecs(ts, X));
338818efac9SLisandro Dalcin   /* F = Function(ta,Xa,Va,Aa) */
3399566063dSJacob Faibussowitsch   PetscCall(TSComputeI2Function(ts, ta, Xa, Va, Aa, F));
3409566063dSJacob Faibussowitsch   PetscCall(VecScale(F, th->scale_F));
341818efac9SLisandro Dalcin   PetscFunctionReturn(0);
342818efac9SLisandro Dalcin }
343818efac9SLisandro Dalcin 
3449371c9d4SSatish Balay static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes, PETSC_UNUSED Vec X, Mat J, Mat P, TS ts) {
345818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
346818efac9SLisandro Dalcin   PetscReal ta = th->stage_time;
347818efac9SLisandro Dalcin   Vec       Xa = th->Xa, Va = th->Va, Aa = th->Aa;
348818efac9SLisandro Dalcin   PetscReal dVdX = th->shift_V, dAdX = th->shift_A;
349818efac9SLisandro Dalcin 
350818efac9SLisandro Dalcin   PetscFunctionBegin;
351818efac9SLisandro Dalcin   /* J,P = Jacobian(ta,Xa,Va,Aa) */
3529566063dSJacob Faibussowitsch   PetscCall(TSComputeI2Jacobian(ts, ta, Xa, Va, Aa, dVdX, dAdX, J, P));
353818efac9SLisandro Dalcin   PetscFunctionReturn(0);
354818efac9SLisandro Dalcin }
355818efac9SLisandro Dalcin 
3569371c9d4SSatish Balay static PetscErrorCode TSReset_Alpha(TS ts) {
357818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
358818efac9SLisandro Dalcin 
359818efac9SLisandro Dalcin   PetscFunctionBegin;
3609566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X0));
3619566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Xa));
3629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X1));
3639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->V0));
3649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Va));
3659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->V1));
3669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->A0));
3679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Aa));
3689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->A1));
3699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_sol_prev));
3709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_dot_prev));
3719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_lte_work[0]));
3729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_lte_work[1]));
373818efac9SLisandro Dalcin   PetscFunctionReturn(0);
374818efac9SLisandro Dalcin }
375818efac9SLisandro Dalcin 
3769371c9d4SSatish Balay static PetscErrorCode TSDestroy_Alpha(TS ts) {
377818efac9SLisandro Dalcin   PetscFunctionBegin;
3789566063dSJacob Faibussowitsch   PetscCall(TSReset_Alpha(ts));
3799566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
380818efac9SLisandro Dalcin 
3819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetRadius_C", NULL));
3829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetParams_C", NULL));
3839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2GetParams_C", NULL));
384818efac9SLisandro Dalcin   PetscFunctionReturn(0);
385818efac9SLisandro Dalcin }
386818efac9SLisandro Dalcin 
3879371c9d4SSatish Balay static PetscErrorCode TSSetUp_Alpha(TS ts) {
388818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
3892ffb9264SLisandro Dalcin   PetscBool match;
390818efac9SLisandro Dalcin 
391818efac9SLisandro Dalcin   PetscFunctionBegin;
3929566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->X0));
3939566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->Xa));
3949566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->X1));
3959566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->V0));
3969566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->Va));
3979566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->V1));
3989566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->A0));
3999566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->Aa));
4009566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->A1));
401818efac9SLisandro Dalcin 
4029566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &ts->adapt));
4039566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt));
4049566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match));
4052ffb9264SLisandro Dalcin   if (!match) {
4069566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev));
4079566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_dot_prev));
4089566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work[0]));
4099566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work[1]));
410818efac9SLisandro Dalcin   }
411818efac9SLisandro Dalcin 
4129566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &ts->snes));
413818efac9SLisandro Dalcin   PetscFunctionReturn(0);
414818efac9SLisandro Dalcin }
415818efac9SLisandro Dalcin 
4169371c9d4SSatish Balay static PetscErrorCode TSSetFromOptions_Alpha(TS ts, PetscOptionItems *PetscOptionsObject) {
417818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
418818efac9SLisandro Dalcin 
419818efac9SLisandro Dalcin   PetscFunctionBegin;
420d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Generalized-Alpha ODE solver options");
421818efac9SLisandro Dalcin   {
422818efac9SLisandro Dalcin     PetscBool flg;
423818efac9SLisandro Dalcin     PetscReal radius = 1;
4249566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_radius", "Spectral radius (high-frequency dissipation)", "TSAlpha2SetRadius", radius, &radius, &flg));
4259566063dSJacob Faibussowitsch     if (flg) PetscCall(TSAlpha2SetRadius(ts, radius));
4269566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_alpha_m", "Algorithmic parameter alpha_m", "TSAlpha2SetParams", th->Alpha_m, &th->Alpha_m, NULL));
4279566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_alpha_f", "Algorithmic parameter alpha_f", "TSAlpha2SetParams", th->Alpha_f, &th->Alpha_f, NULL));
4289566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_gamma", "Algorithmic parameter gamma", "TSAlpha2SetParams", th->Gamma, &th->Gamma, NULL));
4299566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_beta", "Algorithmic parameter beta", "TSAlpha2SetParams", th->Beta, &th->Beta, NULL));
4309566063dSJacob Faibussowitsch     PetscCall(TSAlpha2SetParams(ts, th->Alpha_m, th->Alpha_f, th->Gamma, th->Beta));
431818efac9SLisandro Dalcin   }
432d0609cedSBarry Smith   PetscOptionsHeadEnd();
433818efac9SLisandro Dalcin   PetscFunctionReturn(0);
434818efac9SLisandro Dalcin }
435818efac9SLisandro Dalcin 
4369371c9d4SSatish Balay static PetscErrorCode TSView_Alpha(TS ts, PetscViewer viewer) {
437818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
438818efac9SLisandro Dalcin   PetscBool iascii;
439818efac9SLisandro Dalcin 
440818efac9SLisandro Dalcin   PetscFunctionBegin;
4419566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
442*48a46eb9SPierre 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));
443818efac9SLisandro Dalcin   PetscFunctionReturn(0);
444818efac9SLisandro Dalcin }
445818efac9SLisandro Dalcin 
4469371c9d4SSatish Balay static PetscErrorCode TSAlpha2SetRadius_Alpha(TS ts, PetscReal radius) {
447818efac9SLisandro Dalcin   PetscReal alpha_m, alpha_f, gamma, beta;
448818efac9SLisandro Dalcin 
449818efac9SLisandro Dalcin   PetscFunctionBegin;
450cad9d221SBarry Smith   PetscCheck(radius >= 0 && radius <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius);
451818efac9SLisandro Dalcin   alpha_m = (2 - radius) / (1 + radius);
452818efac9SLisandro Dalcin   alpha_f = 1 / (1 + radius);
453818efac9SLisandro Dalcin   gamma   = (PetscReal)0.5 + alpha_m - alpha_f;
4549371c9d4SSatish Balay   beta    = (PetscReal)0.5 * (1 + alpha_m - alpha_f);
4559371c9d4SSatish Balay   beta *= beta;
4569566063dSJacob Faibussowitsch   PetscCall(TSAlpha2SetParams(ts, alpha_m, alpha_f, gamma, beta));
457818efac9SLisandro Dalcin   PetscFunctionReturn(0);
458818efac9SLisandro Dalcin }
459818efac9SLisandro Dalcin 
4609371c9d4SSatish Balay static PetscErrorCode TSAlpha2SetParams_Alpha(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma, PetscReal beta) {
461818efac9SLisandro Dalcin   TS_Alpha *th  = (TS_Alpha *)ts->data;
462818efac9SLisandro Dalcin   PetscReal tol = 100 * PETSC_MACHINE_EPSILON;
463818efac9SLisandro Dalcin   PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma;
464818efac9SLisandro Dalcin 
465818efac9SLisandro Dalcin   PetscFunctionBegin;
466818efac9SLisandro Dalcin   th->Alpha_m = alpha_m;
467818efac9SLisandro Dalcin   th->Alpha_f = alpha_f;
468818efac9SLisandro Dalcin   th->Gamma   = gamma;
469818efac9SLisandro Dalcin   th->Beta    = beta;
470818efac9SLisandro Dalcin   th->order   = (PetscAbsReal(res) < tol) ? 2 : 1;
471818efac9SLisandro Dalcin   PetscFunctionReturn(0);
472818efac9SLisandro Dalcin }
473818efac9SLisandro Dalcin 
4749371c9d4SSatish Balay static PetscErrorCode TSAlpha2GetParams_Alpha(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma, PetscReal *beta) {
475818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
476818efac9SLisandro Dalcin 
477818efac9SLisandro Dalcin   PetscFunctionBegin;
478818efac9SLisandro Dalcin   if (alpha_m) *alpha_m = th->Alpha_m;
479818efac9SLisandro Dalcin   if (alpha_f) *alpha_f = th->Alpha_f;
480818efac9SLisandro Dalcin   if (gamma) *gamma = th->Gamma;
481818efac9SLisandro Dalcin   if (beta) *beta = th->Beta;
482818efac9SLisandro Dalcin   PetscFunctionReturn(0);
483818efac9SLisandro Dalcin }
484818efac9SLisandro Dalcin 
485818efac9SLisandro Dalcin /*MC
486818efac9SLisandro Dalcin       TSALPHA2 - ODE/DAE solver using the implicit Generalized-Alpha method
487818efac9SLisandro Dalcin                  for second-order systems
488818efac9SLisandro Dalcin 
489818efac9SLisandro Dalcin   Level: beginner
490818efac9SLisandro Dalcin 
491818efac9SLisandro Dalcin   References:
492606c0280SSatish Balay . * - J. Chung, G.M.Hubert. "A Time Integration Algorithm for Structural
493818efac9SLisandro Dalcin   Dynamics with Improved Numerical Dissipation: The Generalized-alpha
494818efac9SLisandro Dalcin   Method" ASME Journal of Applied Mechanics, 60, 371:375, 1993.
495818efac9SLisandro Dalcin 
496db781477SPatrick Sanan .seealso: `TS`, `TSCreate()`, `TSSetType()`, `TSAlpha2SetRadius()`, `TSAlpha2SetParams()`
497818efac9SLisandro Dalcin M*/
4989371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSCreate_Alpha2(TS ts) {
499818efac9SLisandro Dalcin   TS_Alpha *th;
500818efac9SLisandro Dalcin 
501818efac9SLisandro Dalcin   PetscFunctionBegin;
502818efac9SLisandro Dalcin   ts->ops->reset          = TSReset_Alpha;
503818efac9SLisandro Dalcin   ts->ops->destroy        = TSDestroy_Alpha;
504818efac9SLisandro Dalcin   ts->ops->view           = TSView_Alpha;
505818efac9SLisandro Dalcin   ts->ops->setup          = TSSetUp_Alpha;
506818efac9SLisandro Dalcin   ts->ops->setfromoptions = TSSetFromOptions_Alpha;
507818efac9SLisandro Dalcin   ts->ops->step           = TSStep_Alpha;
508818efac9SLisandro Dalcin   ts->ops->evaluatewlte   = TSEvaluateWLTE_Alpha;
509818efac9SLisandro Dalcin   ts->ops->rollback       = TSRollBack_Alpha;
510818efac9SLisandro Dalcin   /*ts->ops->interpolate  = TSInterpolate_Alpha;*/
511818efac9SLisandro Dalcin   ts->ops->snesfunction   = SNESTSFormFunction_Alpha;
512818efac9SLisandro Dalcin   ts->ops->snesjacobian   = SNESTSFormJacobian_Alpha;
5132ffb9264SLisandro Dalcin   ts->default_adapt_type  = TSADAPTNONE;
514818efac9SLisandro Dalcin 
515efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
516efd4aadfSBarry Smith 
5179566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ts, &th));
518818efac9SLisandro Dalcin   ts->data = (void *)th;
519818efac9SLisandro Dalcin 
520818efac9SLisandro Dalcin   th->Alpha_m = 0.5;
521818efac9SLisandro Dalcin   th->Alpha_f = 0.5;
522818efac9SLisandro Dalcin   th->Gamma   = 0.5;
523818efac9SLisandro Dalcin   th->Beta    = 0.25;
524818efac9SLisandro Dalcin   th->order   = 2;
525818efac9SLisandro Dalcin 
5269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetRadius_C", TSAlpha2SetRadius_Alpha));
5279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2SetParams_C", TSAlpha2SetParams_Alpha));
5289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlpha2GetParams_C", TSAlpha2GetParams_Alpha));
529818efac9SLisandro Dalcin   PetscFunctionReturn(0);
530818efac9SLisandro Dalcin }
531818efac9SLisandro Dalcin 
532818efac9SLisandro Dalcin /*@
533818efac9SLisandro Dalcin   TSAlpha2SetRadius - sets the desired spectral radius of the method
534818efac9SLisandro Dalcin                       (i.e. high-frequency numerical damping)
535818efac9SLisandro Dalcin 
536818efac9SLisandro Dalcin   Logically Collective on TS
537818efac9SLisandro Dalcin 
538818efac9SLisandro Dalcin   The algorithmic parameters \alpha_m and \alpha_f of the
539818efac9SLisandro Dalcin   generalized-\alpha method can be computed in terms of a specified
540818efac9SLisandro Dalcin   spectral radius \rho in [0,1] for infinite time step in order to
541818efac9SLisandro Dalcin   control high-frequency numerical damping:
542818efac9SLisandro Dalcin     \alpha_m = (2-\rho)/(1+\rho)
543818efac9SLisandro Dalcin     \alpha_f = 1/(1+\rho)
544818efac9SLisandro Dalcin 
545d8d19677SJose E. Roman   Input Parameters:
546818efac9SLisandro Dalcin +  ts - timestepping context
547818efac9SLisandro Dalcin -  radius - the desired spectral radius
548818efac9SLisandro Dalcin 
549818efac9SLisandro Dalcin   Options Database:
55067b8a455SSatish Balay .  -ts_alpha_radius <radius> - set the desired spectral radius
551818efac9SLisandro Dalcin 
552818efac9SLisandro Dalcin   Level: intermediate
553818efac9SLisandro Dalcin 
554db781477SPatrick Sanan .seealso: `TSAlpha2SetParams()`, `TSAlpha2GetParams()`
555818efac9SLisandro Dalcin @*/
5569371c9d4SSatish Balay PetscErrorCode TSAlpha2SetRadius(TS ts, PetscReal radius) {
557818efac9SLisandro Dalcin   PetscFunctionBegin;
558818efac9SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
559818efac9SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, radius, 2);
560cad9d221SBarry Smith   PetscCheck(radius >= 0 && radius <= 1, ((PetscObject)ts)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius);
561cac4c232SBarry Smith   PetscTryMethod(ts, "TSAlpha2SetRadius_C", (TS, PetscReal), (ts, radius));
562818efac9SLisandro Dalcin   PetscFunctionReturn(0);
563818efac9SLisandro Dalcin }
564818efac9SLisandro Dalcin 
565818efac9SLisandro Dalcin /*@
566818efac9SLisandro Dalcin   TSAlpha2SetParams - sets the algorithmic parameters for TSALPHA2
567818efac9SLisandro Dalcin 
568818efac9SLisandro Dalcin   Logically Collective on TS
569818efac9SLisandro Dalcin 
570818efac9SLisandro Dalcin   Second-order accuracy can be obtained so long as:
571818efac9SLisandro Dalcin     \gamma = 1/2 + alpha_m - alpha_f
572818efac9SLisandro Dalcin     \beta  = 1/4 (1 + alpha_m - alpha_f)^2
573818efac9SLisandro Dalcin 
574818efac9SLisandro Dalcin   Unconditional stability requires:
575818efac9SLisandro Dalcin     \alpha_m >= \alpha_f >= 1/2
576818efac9SLisandro Dalcin 
577d8d19677SJose E. Roman   Input Parameters:
578818efac9SLisandro Dalcin + ts - timestepping context
579a5b23f4aSJose E. Roman . \alpha_m - algorithmic parameter
580a5b23f4aSJose E. Roman . \alpha_f - algorithmic parameter
581a5b23f4aSJose E. Roman . \gamma   - algorithmic parameter
582a5b23f4aSJose E. Roman - \beta    - algorithmic parameter
583818efac9SLisandro Dalcin 
584818efac9SLisandro Dalcin   Options Database:
58567b8a455SSatish Balay + -ts_alpha_alpha_m <alpha_m> - set alpha_m
58667b8a455SSatish Balay . -ts_alpha_alpha_f <alpha_f> - set alpha_f
58767b8a455SSatish Balay . -ts_alpha_gamma   <gamma> - set gamma
58867b8a455SSatish Balay - -ts_alpha_beta    <beta> - set beta
589818efac9SLisandro Dalcin 
590818efac9SLisandro Dalcin   Note:
591818efac9SLisandro Dalcin   Use of this function is normally only required to hack TSALPHA2 to
592818efac9SLisandro Dalcin   use a modified integration scheme. Users should call
593818efac9SLisandro Dalcin   TSAlpha2SetRadius() to set the desired spectral radius of the methods
594818efac9SLisandro Dalcin   (i.e. high-frequency damping) in order so select optimal values for
595818efac9SLisandro Dalcin   these parameters.
596818efac9SLisandro Dalcin 
597818efac9SLisandro Dalcin   Level: advanced
598818efac9SLisandro Dalcin 
599db781477SPatrick Sanan .seealso: `TSAlpha2SetRadius()`, `TSAlpha2GetParams()`
600818efac9SLisandro Dalcin @*/
6019371c9d4SSatish Balay PetscErrorCode TSAlpha2SetParams(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma, PetscReal beta) {
602818efac9SLisandro Dalcin   PetscFunctionBegin;
603818efac9SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
604818efac9SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, alpha_m, 2);
605818efac9SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, alpha_f, 3);
606818efac9SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, gamma, 4);
607818efac9SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, beta, 5);
608cac4c232SBarry Smith   PetscTryMethod(ts, "TSAlpha2SetParams_C", (TS, PetscReal, PetscReal, PetscReal, PetscReal), (ts, alpha_m, alpha_f, gamma, beta));
609818efac9SLisandro Dalcin   PetscFunctionReturn(0);
610818efac9SLisandro Dalcin }
611818efac9SLisandro Dalcin 
612818efac9SLisandro Dalcin /*@
613818efac9SLisandro Dalcin   TSAlpha2GetParams - gets the algorithmic parameters for TSALPHA2
614818efac9SLisandro Dalcin 
615818efac9SLisandro Dalcin   Not Collective
616818efac9SLisandro Dalcin 
617818efac9SLisandro Dalcin   Input Parameter:
618818efac9SLisandro Dalcin . ts - timestepping context
619818efac9SLisandro Dalcin 
620818efac9SLisandro Dalcin   Output Parameters:
621818efac9SLisandro Dalcin + \alpha_m - algorithmic parameter
622818efac9SLisandro Dalcin . \alpha_f - algorithmic parameter
623818efac9SLisandro Dalcin . \gamma   - algorithmic parameter
624818efac9SLisandro Dalcin - \beta    - algorithmic parameter
625818efac9SLisandro Dalcin 
626818efac9SLisandro Dalcin   Note:
627818efac9SLisandro Dalcin   Use of this function is normally only required to hack TSALPHA2 to
628818efac9SLisandro Dalcin   use a modified integration scheme. Users should call
629818efac9SLisandro Dalcin   TSAlpha2SetRadius() to set the high-frequency damping (i.e. spectral
630818efac9SLisandro Dalcin   radius of the method) in order so select optimal values for these
631818efac9SLisandro Dalcin   parameters.
632818efac9SLisandro Dalcin 
633818efac9SLisandro Dalcin   Level: advanced
634818efac9SLisandro Dalcin 
635db781477SPatrick Sanan .seealso: `TSAlpha2SetRadius()`, `TSAlpha2SetParams()`
636818efac9SLisandro Dalcin @*/
6379371c9d4SSatish Balay PetscErrorCode TSAlpha2GetParams(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma, PetscReal *beta) {
638818efac9SLisandro Dalcin   PetscFunctionBegin;
639818efac9SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
640818efac9SLisandro Dalcin   if (alpha_m) PetscValidRealPointer(alpha_m, 2);
641818efac9SLisandro Dalcin   if (alpha_f) PetscValidRealPointer(alpha_f, 3);
642818efac9SLisandro Dalcin   if (gamma) PetscValidRealPointer(gamma, 4);
643818efac9SLisandro Dalcin   if (beta) PetscValidRealPointer(beta, 5);
644cac4c232SBarry Smith   PetscUseMethod(ts, "TSAlpha2GetParams_C", (TS, PetscReal *, PetscReal *, PetscReal *, PetscReal *), (ts, alpha_m, alpha_f, gamma, beta));
645818efac9SLisandro Dalcin   PetscFunctionReturn(0);
646818efac9SLisandro Dalcin }
647