xref: /petsc/src/ts/impls/implicit/alpha/alpha1.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1b07a2398SLisandro Dalcin /*
2b07a2398SLisandro Dalcin   Code for timestepping with implicit generalized-\alpha method
3b07a2398SLisandro Dalcin   for first order systems.
4b07a2398SLisandro Dalcin */
5b07a2398SLisandro Dalcin #include <petsc/private/tsimpl.h> /*I   "petscts.h"   I*/
6b07a2398SLisandro Dalcin 
7b07a2398SLisandro Dalcin static PetscBool  cited      = PETSC_FALSE;
89371c9d4SSatish Balay static const char citation[] = "@article{Jansen2000,\n"
9b07a2398SLisandro Dalcin                                "  title   = {A generalized-$\\alpha$ method for integrating the filtered {N}avier--{S}tokes equations with a stabilized finite element method},\n"
10b07a2398SLisandro Dalcin                                "  author  = {Kenneth E. Jansen and Christian H. Whiting and Gregory M. Hulbert},\n"
11b07a2398SLisandro Dalcin                                "  journal = {Computer Methods in Applied Mechanics and Engineering},\n"
12b07a2398SLisandro Dalcin                                "  volume  = {190},\n"
13b07a2398SLisandro Dalcin                                "  number  = {3--4},\n"
14b07a2398SLisandro Dalcin                                "  pages   = {305--319},\n"
15b07a2398SLisandro Dalcin                                "  year    = {2000},\n"
16b07a2398SLisandro Dalcin                                "  issn    = {0045-7825},\n"
17b07a2398SLisandro Dalcin                                "  doi     = {http://dx.doi.org/10.1016/S0045-7825(00)00203-6}\n}\n";
18b07a2398SLisandro Dalcin 
19b07a2398SLisandro Dalcin typedef struct {
20b07a2398SLisandro Dalcin   PetscReal stage_time;
21b07a2398SLisandro Dalcin   PetscReal shift_V;
22b07a2398SLisandro Dalcin   PetscReal scale_F;
23b07a2398SLisandro Dalcin   Vec       X0, Xa, X1;
24b07a2398SLisandro Dalcin   Vec       V0, Va, V1;
25b07a2398SLisandro Dalcin 
26b07a2398SLisandro Dalcin   PetscReal Alpha_m;
27b07a2398SLisandro Dalcin   PetscReal Alpha_f;
28b07a2398SLisandro Dalcin   PetscReal Gamma;
29b07a2398SLisandro Dalcin   PetscInt  order;
30b07a2398SLisandro Dalcin 
31b07a2398SLisandro Dalcin   Vec vec_sol_prev;
321566a47fSLisandro Dalcin   Vec vec_lte_work;
33b07a2398SLisandro Dalcin 
34b07a2398SLisandro Dalcin   TSStepStatus status;
35b07a2398SLisandro Dalcin } TS_Alpha;
36b07a2398SLisandro Dalcin 
379371c9d4SSatish Balay static PetscErrorCode TSAlpha_StageTime(TS ts) {
38b07a2398SLisandro Dalcin   TS_Alpha *th      = (TS_Alpha *)ts->data;
39b07a2398SLisandro Dalcin   PetscReal t       = ts->ptime;
40b07a2398SLisandro Dalcin   PetscReal dt      = ts->time_step;
41b07a2398SLisandro Dalcin   PetscReal Alpha_m = th->Alpha_m;
42b07a2398SLisandro Dalcin   PetscReal Alpha_f = th->Alpha_f;
43b07a2398SLisandro Dalcin   PetscReal Gamma   = th->Gamma;
44b07a2398SLisandro Dalcin 
45b07a2398SLisandro Dalcin   PetscFunctionBegin;
46b07a2398SLisandro Dalcin   th->stage_time = t + Alpha_f * dt;
47b07a2398SLisandro Dalcin   th->shift_V    = Alpha_m / (Alpha_f * Gamma * dt);
48b07a2398SLisandro Dalcin   th->scale_F    = 1 / Alpha_f;
49b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
50b07a2398SLisandro Dalcin }
51b07a2398SLisandro Dalcin 
529371c9d4SSatish Balay static PetscErrorCode TSAlpha_StageVecs(TS ts, Vec X) {
53b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
54b07a2398SLisandro Dalcin   Vec       X1 = X, V1 = th->V1;
55b07a2398SLisandro Dalcin   Vec       Xa = th->Xa, Va = th->Va;
56b07a2398SLisandro Dalcin   Vec       X0 = th->X0, V0 = th->V0;
57b07a2398SLisandro Dalcin   PetscReal dt      = ts->time_step;
58b07a2398SLisandro Dalcin   PetscReal Alpha_m = th->Alpha_m;
59b07a2398SLisandro Dalcin   PetscReal Alpha_f = th->Alpha_f;
60b07a2398SLisandro Dalcin   PetscReal Gamma   = th->Gamma;
61b07a2398SLisandro Dalcin 
62b07a2398SLisandro Dalcin   PetscFunctionBegin;
63b07a2398SLisandro Dalcin   /* V1 = 1/(Gamma*dT)*(X1-X0) + (1-1/Gamma)*V0 */
649566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(V1, -1.0, X0, X1));
659566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(V1, 1 - 1 / Gamma, 1 / (Gamma * dt), V0));
66b07a2398SLisandro Dalcin   /* Xa = X0 + Alpha_f*(X1-X0) */
679566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Xa, -1.0, X0, X1));
689566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Xa, Alpha_f, X0));
69b07a2398SLisandro Dalcin   /* Va = V0 + Alpha_m*(V1-V0) */
709566063dSJacob Faibussowitsch   PetscCall(VecWAXPY(Va, -1.0, V0, V1));
719566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Va, Alpha_m, V0));
72b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
73b07a2398SLisandro Dalcin }
74b07a2398SLisandro Dalcin 
759371c9d4SSatish Balay static PetscErrorCode TSAlpha_SNESSolve(TS ts, Vec b, Vec x) {
76b07a2398SLisandro Dalcin   PetscInt nits, lits;
77b07a2398SLisandro Dalcin 
78b07a2398SLisandro Dalcin   PetscFunctionBegin;
799566063dSJacob Faibussowitsch   PetscCall(SNESSolve(ts->snes, b, x));
809566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(ts->snes, &nits));
819566063dSJacob Faibussowitsch   PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits));
829371c9d4SSatish Balay   ts->snes_its += nits;
839371c9d4SSatish Balay   ts->ksp_its += lits;
84b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
85b07a2398SLisandro Dalcin }
86b07a2398SLisandro Dalcin 
87b07a2398SLisandro Dalcin /*
88b07a2398SLisandro Dalcin   Compute a consistent initial state for the generalized-alpha method.
89b07a2398SLisandro Dalcin   - Solve two successive backward Euler steps with halved time step.
90b07a2398SLisandro Dalcin   - Compute the initial time derivative using backward differences.
91b07a2398SLisandro Dalcin   - If using adaptivity, estimate the LTE of the initial step.
92b07a2398SLisandro Dalcin */
939371c9d4SSatish Balay static PetscErrorCode TSAlpha_Restart(TS ts, PetscBool *initok) {
94b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
95b07a2398SLisandro Dalcin   PetscReal time_step;
96b07a2398SLisandro Dalcin   PetscReal alpha_m, alpha_f, gamma;
97b07a2398SLisandro Dalcin   Vec       X0 = ts->vec_sol, X1, X2 = th->X1;
98b07a2398SLisandro Dalcin   PetscBool stageok;
99b07a2398SLisandro Dalcin 
100b07a2398SLisandro Dalcin   PetscFunctionBegin;
1019566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X0, &X1));
102b07a2398SLisandro Dalcin 
103b07a2398SLisandro Dalcin   /* Setup backward Euler with halved time step */
1049566063dSJacob Faibussowitsch   PetscCall(TSAlphaGetParams(ts, &alpha_m, &alpha_f, &gamma));
1059566063dSJacob Faibussowitsch   PetscCall(TSAlphaSetParams(ts, 1, 1, 1));
1069566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &time_step));
107b07a2398SLisandro Dalcin   ts->time_step = time_step / 2;
1089566063dSJacob Faibussowitsch   PetscCall(TSAlpha_StageTime(ts));
109b07a2398SLisandro Dalcin   th->stage_time = ts->ptime;
1109566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(th->V0));
111b07a2398SLisandro Dalcin 
112b07a2398SLisandro Dalcin   /* First BE step, (t0,X0) -> (t1,X1) */
113b07a2398SLisandro Dalcin   th->stage_time += ts->time_step;
1149566063dSJacob Faibussowitsch   PetscCall(VecCopy(X0, th->X0));
1159566063dSJacob Faibussowitsch   PetscCall(TSPreStage(ts, th->stage_time));
1169566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->X0, X1));
1179566063dSJacob Faibussowitsch   PetscCall(TSAlpha_SNESSolve(ts, NULL, X1));
1189566063dSJacob Faibussowitsch   PetscCall(TSPostStage(ts, th->stage_time, 0, &X1));
1199566063dSJacob Faibussowitsch   PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, X1, &stageok));
120b07a2398SLisandro Dalcin   if (!stageok) goto finally;
121b07a2398SLisandro Dalcin 
122b07a2398SLisandro Dalcin   /* Second BE step, (t1,X1) -> (t2,X2) */
123b07a2398SLisandro Dalcin   th->stage_time += ts->time_step;
1249566063dSJacob Faibussowitsch   PetscCall(VecCopy(X1, th->X0));
1259566063dSJacob Faibussowitsch   PetscCall(TSPreStage(ts, th->stage_time));
1269566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->X0, X2));
1279566063dSJacob Faibussowitsch   PetscCall(TSAlpha_SNESSolve(ts, NULL, X2));
1289566063dSJacob Faibussowitsch   PetscCall(TSPostStage(ts, th->stage_time, 0, &X2));
1299566063dSJacob Faibussowitsch   PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, X2, &stageok));
130b07a2398SLisandro Dalcin   if (!stageok) goto finally;
131b07a2398SLisandro Dalcin 
132b07a2398SLisandro Dalcin   /* Compute V0 ~ dX/dt at t0 with backward differences */
1339566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(th->V0));
1349566063dSJacob Faibussowitsch   PetscCall(VecAXPY(th->V0, -3 / ts->time_step, X0));
1359566063dSJacob Faibussowitsch   PetscCall(VecAXPY(th->V0, +4 / ts->time_step, X1));
1369566063dSJacob Faibussowitsch   PetscCall(VecAXPY(th->V0, -1 / ts->time_step, X2));
137b07a2398SLisandro Dalcin 
138b07a2398SLisandro Dalcin   /* Rough, lower-order estimate LTE of the initial step */
1392ffb9264SLisandro Dalcin   if (th->vec_lte_work) {
1409566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(th->vec_lte_work));
1419566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work, +2, X2));
1429566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work, -4, X1));
1439566063dSJacob Faibussowitsch     PetscCall(VecAXPY(th->vec_lte_work, +2, X0));
144b07a2398SLisandro Dalcin   }
145b07a2398SLisandro Dalcin 
146b07a2398SLisandro Dalcin finally:
147b07a2398SLisandro Dalcin   /* Revert TSAlpha to the initial state (t0,X0) */
148b07a2398SLisandro Dalcin   if (initok) *initok = stageok;
1499566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, time_step));
1509566063dSJacob Faibussowitsch   PetscCall(TSAlphaSetParams(ts, alpha_m, alpha_f, gamma));
1519566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol, th->X0));
152b07a2398SLisandro Dalcin 
1539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X1));
154b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
155b07a2398SLisandro Dalcin }
156b07a2398SLisandro Dalcin 
1579371c9d4SSatish Balay static PetscErrorCode TSStep_Alpha(TS ts) {
158b07a2398SLisandro Dalcin   TS_Alpha *th         = (TS_Alpha *)ts->data;
159b07a2398SLisandro Dalcin   PetscInt  rejections = 0;
160b07a2398SLisandro Dalcin   PetscBool stageok, accept = PETSC_TRUE;
161b07a2398SLisandro Dalcin   PetscReal next_time_step = ts->time_step;
162b07a2398SLisandro Dalcin 
163b07a2398SLisandro Dalcin   PetscFunctionBegin;
1649566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(citation, &cited));
165b07a2398SLisandro Dalcin 
166b07a2398SLisandro Dalcin   if (!ts->steprollback) {
1679566063dSJacob Faibussowitsch     if (th->vec_sol_prev) PetscCall(VecCopy(th->X0, th->vec_sol_prev));
1689566063dSJacob Faibussowitsch     PetscCall(VecCopy(ts->vec_sol, th->X0));
1699566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->V1, th->V0));
170b07a2398SLisandro Dalcin   }
171b07a2398SLisandro Dalcin 
1721566a47fSLisandro Dalcin   th->status = TS_STEP_INCOMPLETE;
173b07a2398SLisandro Dalcin   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
174fecfb714SLisandro Dalcin     if (ts->steprestart) {
1759566063dSJacob Faibussowitsch       PetscCall(TSAlpha_Restart(ts, &stageok));
176fecfb714SLisandro Dalcin       if (!stageok) goto reject_step;
177b07a2398SLisandro Dalcin     }
178b07a2398SLisandro Dalcin 
1799566063dSJacob Faibussowitsch     PetscCall(TSAlpha_StageTime(ts));
1809566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->X0, th->X1));
1819566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, th->stage_time));
1829566063dSJacob Faibussowitsch     PetscCall(TSAlpha_SNESSolve(ts, NULL, th->X1));
1839566063dSJacob Faibussowitsch     PetscCall(TSPostStage(ts, th->stage_time, 0, &th->Xa));
1849566063dSJacob Faibussowitsch     PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->Xa, &stageok));
185fecfb714SLisandro Dalcin     if (!stageok) goto reject_step;
186b07a2398SLisandro Dalcin 
1871566a47fSLisandro Dalcin     th->status = TS_STEP_PENDING;
1889566063dSJacob Faibussowitsch     PetscCall(VecCopy(th->X1, ts->vec_sol));
1899566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
1901566a47fSLisandro Dalcin     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
191b07a2398SLisandro Dalcin     if (!accept) {
1929566063dSJacob Faibussowitsch       PetscCall(VecCopy(th->X0, ts->vec_sol));
193be5899b3SLisandro Dalcin       ts->time_step = next_time_step;
194be5899b3SLisandro Dalcin       goto reject_step;
195b07a2398SLisandro Dalcin     }
196b07a2398SLisandro Dalcin 
197b07a2398SLisandro Dalcin     ts->ptime += ts->time_step;
198b07a2398SLisandro Dalcin     ts->time_step = next_time_step;
199b07a2398SLisandro Dalcin     break;
200b07a2398SLisandro Dalcin 
201b07a2398SLisandro Dalcin   reject_step:
2029371c9d4SSatish Balay     ts->reject++;
2039371c9d4SSatish Balay     accept = PETSC_FALSE;
204b07a2398SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
205b07a2398SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
20663a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
207b07a2398SLisandro Dalcin     }
208b07a2398SLisandro Dalcin   }
209b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
210b07a2398SLisandro Dalcin }
211b07a2398SLisandro Dalcin 
2129371c9d4SSatish Balay static PetscErrorCode TSEvaluateWLTE_Alpha(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte) {
213b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
2149808bdc1SLisandro Dalcin   Vec       X  = th->X1;           /* X = solution */
2151566a47fSLisandro Dalcin   Vec       Y  = th->vec_lte_work; /* Y = X + LTE  */
2167453f775SEmil Constantinescu   PetscReal wltea, wlter;
217b07a2398SLisandro Dalcin 
218b07a2398SLisandro Dalcin   PetscFunctionBegin;
2199371c9d4SSatish Balay   if (!th->vec_sol_prev) {
2209371c9d4SSatish Balay     *wlte = -1;
2219371c9d4SSatish Balay     PetscFunctionReturn(0);
2229371c9d4SSatish Balay   }
2239371c9d4SSatish Balay   if (!th->vec_lte_work) {
2249371c9d4SSatish Balay     *wlte = -1;
2259371c9d4SSatish Balay     PetscFunctionReturn(0);
2269371c9d4SSatish Balay   }
227fecfb714SLisandro Dalcin   if (ts->steprestart) {
228fecfb714SLisandro Dalcin     /* th->vec_lte_work is set to the LTE in TSAlpha_Restart() */
2299566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Y, 1, X));
230b07a2398SLisandro Dalcin   } else {
231b07a2398SLisandro Dalcin     /* Compute LTE using backward differences with non-constant time step */
232be5899b3SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
233be5899b3SLisandro Dalcin     PetscReal   a = 1 + h_prev / h;
2349371c9d4SSatish Balay     PetscScalar scal[3];
2359371c9d4SSatish Balay     Vec         vecs[3];
2369371c9d4SSatish Balay     scal[0] = +1 / a;
2379371c9d4SSatish Balay     scal[1] = -1 / (a - 1);
2389371c9d4SSatish Balay     scal[2] = +1 / (a * (a - 1));
2399371c9d4SSatish Balay     vecs[0] = th->X1;
2409371c9d4SSatish Balay     vecs[1] = th->X0;
2419371c9d4SSatish Balay     vecs[2] = th->vec_sol_prev;
2429566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Y));
2439566063dSJacob Faibussowitsch     PetscCall(VecMAXPY(Y, 3, scal, vecs));
244b07a2398SLisandro Dalcin   }
2459566063dSJacob Faibussowitsch   PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, wlte, &wltea, &wlter));
2469808bdc1SLisandro Dalcin   if (order) *order = 2;
247b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
248b07a2398SLisandro Dalcin }
249b07a2398SLisandro Dalcin 
2509371c9d4SSatish Balay static PetscErrorCode TSRollBack_Alpha(TS ts) {
251b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
252b07a2398SLisandro Dalcin 
253b07a2398SLisandro Dalcin   PetscFunctionBegin;
2549566063dSJacob Faibussowitsch   PetscCall(VecCopy(th->X0, ts->vec_sol));
255b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
256b07a2398SLisandro Dalcin }
257b07a2398SLisandro Dalcin 
2589371c9d4SSatish Balay static PetscErrorCode TSInterpolate_Alpha(TS ts, PetscReal t, Vec X) {
259b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
260b07a2398SLisandro Dalcin   PetscReal dt = t - ts->ptime;
261b07a2398SLisandro Dalcin 
262b07a2398SLisandro Dalcin   PetscFunctionBegin;
2639566063dSJacob Faibussowitsch   PetscCall(VecCopy(ts->vec_sol, X));
2649566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, th->Gamma * dt, th->V1));
2659566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, (1 - th->Gamma) * dt, th->V0));
266b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
267b07a2398SLisandro Dalcin }
268b07a2398SLisandro Dalcin 
2699371c9d4SSatish Balay static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes, Vec X, Vec F, TS ts) {
270b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
271b07a2398SLisandro Dalcin   PetscReal ta = th->stage_time;
272b07a2398SLisandro Dalcin   Vec       Xa = th->Xa, Va = th->Va;
273b07a2398SLisandro Dalcin 
274b07a2398SLisandro Dalcin   PetscFunctionBegin;
2759566063dSJacob Faibussowitsch   PetscCall(TSAlpha_StageVecs(ts, X));
276b07a2398SLisandro Dalcin   /* F = Function(ta,Xa,Va) */
2779566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, ta, Xa, Va, F, PETSC_FALSE));
2789566063dSJacob Faibussowitsch   PetscCall(VecScale(F, th->scale_F));
279b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
280b07a2398SLisandro Dalcin }
281b07a2398SLisandro Dalcin 
2829371c9d4SSatish Balay static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes, PETSC_UNUSED Vec X, Mat J, Mat P, TS ts) {
283b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
284b07a2398SLisandro Dalcin   PetscReal ta = th->stage_time;
285b07a2398SLisandro Dalcin   Vec       Xa = th->Xa, Va = th->Va;
286b07a2398SLisandro Dalcin   PetscReal dVdX = th->shift_V;
287b07a2398SLisandro Dalcin 
288b07a2398SLisandro Dalcin   PetscFunctionBegin;
289b07a2398SLisandro Dalcin   /* J,P = Jacobian(ta,Xa,Va) */
2909566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, ta, Xa, Va, dVdX, J, P, PETSC_FALSE));
291b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
292b07a2398SLisandro Dalcin }
293b07a2398SLisandro Dalcin 
2949371c9d4SSatish Balay static PetscErrorCode TSReset_Alpha(TS ts) {
295b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
296b07a2398SLisandro Dalcin 
297b07a2398SLisandro Dalcin   PetscFunctionBegin;
2989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X0));
2999566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Xa));
3009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->X1));
3019566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->V0));
3029566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->Va));
3039566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->V1));
3049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_sol_prev));
3059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&th->vec_lte_work));
306b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
307b07a2398SLisandro Dalcin }
308b07a2398SLisandro Dalcin 
3099371c9d4SSatish Balay static PetscErrorCode TSDestroy_Alpha(TS ts) {
310b07a2398SLisandro Dalcin   PetscFunctionBegin;
3119566063dSJacob Faibussowitsch   PetscCall(TSReset_Alpha(ts));
3129566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
313b07a2398SLisandro Dalcin 
3149566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetRadius_C", NULL));
3159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetParams_C", NULL));
3169566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaGetParams_C", NULL));
317b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
318b07a2398SLisandro Dalcin }
319b07a2398SLisandro Dalcin 
3209371c9d4SSatish Balay static PetscErrorCode TSSetUp_Alpha(TS ts) {
321b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
3222ffb9264SLisandro Dalcin   PetscBool match;
323b07a2398SLisandro Dalcin 
324b07a2398SLisandro Dalcin   PetscFunctionBegin;
3259566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->X0));
3269566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->Xa));
3279566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->X1));
3289566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->V0));
3299566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->Va));
3309566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol, &th->V1));
3311566a47fSLisandro Dalcin 
3329566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts, &ts->adapt));
3339566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt));
3349566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match));
3352ffb9264SLisandro Dalcin   if (!match) {
3369566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev));
3379566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work));
338b07a2398SLisandro Dalcin   }
3391566a47fSLisandro Dalcin 
3409566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &ts->snes));
341b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
342b07a2398SLisandro Dalcin }
343b07a2398SLisandro Dalcin 
3449371c9d4SSatish Balay static PetscErrorCode TSSetFromOptions_Alpha(TS ts, PetscOptionItems *PetscOptionsObject) {
345b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
346b07a2398SLisandro Dalcin 
347b07a2398SLisandro Dalcin   PetscFunctionBegin;
348d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Generalized-Alpha ODE solver options");
349b07a2398SLisandro Dalcin   {
350b07a2398SLisandro Dalcin     PetscBool flg;
351b07a2398SLisandro Dalcin     PetscReal radius = 1;
3529566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_radius", "Spectral radius (high-frequency dissipation)", "TSAlphaSetRadius", radius, &radius, &flg));
3539566063dSJacob Faibussowitsch     if (flg) PetscCall(TSAlphaSetRadius(ts, radius));
3549566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_alpha_m", "Algorithmic parameter alpha_m", "TSAlphaSetParams", th->Alpha_m, &th->Alpha_m, NULL));
3559566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_alpha_f", "Algorithmic parameter alpha_f", "TSAlphaSetParams", th->Alpha_f, &th->Alpha_f, NULL));
3569566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-ts_alpha_gamma", "Algorithmic parameter gamma", "TSAlphaSetParams", th->Gamma, &th->Gamma, NULL));
3579566063dSJacob Faibussowitsch     PetscCall(TSAlphaSetParams(ts, th->Alpha_m, th->Alpha_f, th->Gamma));
358b07a2398SLisandro Dalcin   }
359d0609cedSBarry Smith   PetscOptionsHeadEnd();
360b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
361b07a2398SLisandro Dalcin }
362b07a2398SLisandro Dalcin 
3639371c9d4SSatish Balay static PetscErrorCode TSView_Alpha(TS ts, PetscViewer viewer) {
364b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
3659c334d8fSLisandro Dalcin   PetscBool iascii;
366b07a2398SLisandro Dalcin 
367b07a2398SLisandro Dalcin   PetscFunctionBegin;
3689566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
369*48a46eb9SPierre Jolivet   if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, "  Alpha_m=%g, Alpha_f=%g, Gamma=%g\n", (double)th->Alpha_m, (double)th->Alpha_f, (double)th->Gamma));
370b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
371b07a2398SLisandro Dalcin }
372b07a2398SLisandro Dalcin 
3739371c9d4SSatish Balay static PetscErrorCode TSAlphaSetRadius_Alpha(TS ts, PetscReal radius) {
374b07a2398SLisandro Dalcin   PetscReal alpha_m, alpha_f, gamma;
375b07a2398SLisandro Dalcin 
376b07a2398SLisandro Dalcin   PetscFunctionBegin;
377cad9d221SBarry Smith   PetscCheck(radius >= 0 && radius <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius);
378b07a2398SLisandro Dalcin   alpha_m = (PetscReal)0.5 * (3 - radius) / (1 + radius);
379b07a2398SLisandro Dalcin   alpha_f = 1 / (1 + radius);
380b07a2398SLisandro Dalcin   gamma   = (PetscReal)0.5 + alpha_m - alpha_f;
3819566063dSJacob Faibussowitsch   PetscCall(TSAlphaSetParams(ts, alpha_m, alpha_f, gamma));
382b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
383b07a2398SLisandro Dalcin }
384b07a2398SLisandro Dalcin 
3859371c9d4SSatish Balay static PetscErrorCode TSAlphaSetParams_Alpha(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma) {
386b07a2398SLisandro Dalcin   TS_Alpha *th  = (TS_Alpha *)ts->data;
387b07a2398SLisandro Dalcin   PetscReal tol = 100 * PETSC_MACHINE_EPSILON;
388b07a2398SLisandro Dalcin   PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma;
389b07a2398SLisandro Dalcin 
390b07a2398SLisandro Dalcin   PetscFunctionBegin;
391b07a2398SLisandro Dalcin   th->Alpha_m = alpha_m;
392b07a2398SLisandro Dalcin   th->Alpha_f = alpha_f;
393b07a2398SLisandro Dalcin   th->Gamma   = gamma;
394b07a2398SLisandro Dalcin   th->order   = (PetscAbsReal(res) < tol) ? 2 : 1;
395b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
396b07a2398SLisandro Dalcin }
397b07a2398SLisandro Dalcin 
3989371c9d4SSatish Balay static PetscErrorCode TSAlphaGetParams_Alpha(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma) {
399b07a2398SLisandro Dalcin   TS_Alpha *th = (TS_Alpha *)ts->data;
400b07a2398SLisandro Dalcin 
401b07a2398SLisandro Dalcin   PetscFunctionBegin;
402b07a2398SLisandro Dalcin   if (alpha_m) *alpha_m = th->Alpha_m;
403b07a2398SLisandro Dalcin   if (alpha_f) *alpha_f = th->Alpha_f;
404b07a2398SLisandro Dalcin   if (gamma) *gamma = th->Gamma;
405b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
406b07a2398SLisandro Dalcin }
407b07a2398SLisandro Dalcin 
408b07a2398SLisandro Dalcin /*MC
409b07a2398SLisandro Dalcin       TSALPHA - ODE/DAE solver using the implicit Generalized-Alpha method
410b07a2398SLisandro Dalcin                 for first-order systems
411b07a2398SLisandro Dalcin 
412b07a2398SLisandro Dalcin   Level: beginner
413b07a2398SLisandro Dalcin 
414b07a2398SLisandro Dalcin   References:
415606c0280SSatish Balay + * - K.E. Jansen, C.H. Whiting, G.M. Hulber, "A generalized-alpha
416b07a2398SLisandro Dalcin   method for integrating the filtered Navier-Stokes equations with a
417b07a2398SLisandro Dalcin   stabilized finite element method", Computer Methods in Applied
418b07a2398SLisandro Dalcin   Mechanics and Engineering, 190, 305-319, 2000.
419b07a2398SLisandro Dalcin   DOI: 10.1016/S0045-7825(00)00203-6.
420606c0280SSatish Balay - * -  J. Chung, G.M.Hubert. "A Time Integration Algorithm for Structural
421b07a2398SLisandro Dalcin   Dynamics with Improved Numerical Dissipation: The Generalized-alpha
422b07a2398SLisandro Dalcin   Method" ASME Journal of Applied Mechanics, 60, 371:375, 1993.
423b07a2398SLisandro Dalcin 
424db781477SPatrick Sanan .seealso: `TS`, `TSCreate()`, `TSSetType()`, `TSAlphaSetRadius()`, `TSAlphaSetParams()`
425b07a2398SLisandro Dalcin M*/
4269371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSCreate_Alpha(TS ts) {
427b07a2398SLisandro Dalcin   TS_Alpha *th;
428b07a2398SLisandro Dalcin 
429b07a2398SLisandro Dalcin   PetscFunctionBegin;
430b07a2398SLisandro Dalcin   ts->ops->reset          = TSReset_Alpha;
431b07a2398SLisandro Dalcin   ts->ops->destroy        = TSDestroy_Alpha;
432b07a2398SLisandro Dalcin   ts->ops->view           = TSView_Alpha;
433b07a2398SLisandro Dalcin   ts->ops->setup          = TSSetUp_Alpha;
434b07a2398SLisandro Dalcin   ts->ops->setfromoptions = TSSetFromOptions_Alpha;
435b07a2398SLisandro Dalcin   ts->ops->step           = TSStep_Alpha;
4369808bdc1SLisandro Dalcin   ts->ops->evaluatewlte   = TSEvaluateWLTE_Alpha;
437b07a2398SLisandro Dalcin   ts->ops->rollback       = TSRollBack_Alpha;
438b07a2398SLisandro Dalcin   ts->ops->interpolate    = TSInterpolate_Alpha;
439b07a2398SLisandro Dalcin   ts->ops->snesfunction   = SNESTSFormFunction_Alpha;
440b07a2398SLisandro Dalcin   ts->ops->snesjacobian   = SNESTSFormJacobian_Alpha;
4412ffb9264SLisandro Dalcin   ts->default_adapt_type  = TSADAPTNONE;
442b07a2398SLisandro Dalcin 
443efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
444efd4aadfSBarry Smith 
4459566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ts, &th));
446b07a2398SLisandro Dalcin   ts->data = (void *)th;
447b07a2398SLisandro Dalcin 
448b07a2398SLisandro Dalcin   th->Alpha_m = 0.5;
449b07a2398SLisandro Dalcin   th->Alpha_f = 0.5;
450b07a2398SLisandro Dalcin   th->Gamma   = 0.5;
451b07a2398SLisandro Dalcin   th->order   = 2;
452b07a2398SLisandro Dalcin 
4539566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetRadius_C", TSAlphaSetRadius_Alpha));
4549566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaSetParams_C", TSAlphaSetParams_Alpha));
4559566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSAlphaGetParams_C", TSAlphaGetParams_Alpha));
456b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
457b07a2398SLisandro Dalcin }
458b07a2398SLisandro Dalcin 
459b07a2398SLisandro Dalcin /*@
460b07a2398SLisandro Dalcin   TSAlphaSetRadius - sets the desired spectral radius of the method
461b07a2398SLisandro Dalcin                      (i.e. high-frequency numerical damping)
462b07a2398SLisandro Dalcin 
463b07a2398SLisandro Dalcin   Logically Collective on TS
464b07a2398SLisandro Dalcin 
465b07a2398SLisandro Dalcin   The algorithmic parameters \alpha_m and \alpha_f of the
466b07a2398SLisandro Dalcin   generalized-\alpha method can be computed in terms of a specified
467b07a2398SLisandro Dalcin   spectral radius \rho in [0,1] for infinite time step in order to
468b07a2398SLisandro Dalcin   control high-frequency numerical damping:
469b07a2398SLisandro Dalcin     \alpha_m = 0.5*(3-\rho)/(1+\rho)
470b07a2398SLisandro Dalcin     \alpha_f = 1/(1+\rho)
471b07a2398SLisandro Dalcin 
472d8d19677SJose E. Roman   Input Parameters:
473b07a2398SLisandro Dalcin +  ts - timestepping context
474b07a2398SLisandro Dalcin -  radius - the desired spectral radius
475b07a2398SLisandro Dalcin 
476b07a2398SLisandro Dalcin   Options Database:
47767b8a455SSatish Balay .  -ts_alpha_radius <radius> - set alpha radius
478b07a2398SLisandro Dalcin 
479b07a2398SLisandro Dalcin   Level: intermediate
480b07a2398SLisandro Dalcin 
481db781477SPatrick Sanan .seealso: `TSAlphaSetParams()`, `TSAlphaGetParams()`
482b07a2398SLisandro Dalcin @*/
4839371c9d4SSatish Balay PetscErrorCode TSAlphaSetRadius(TS ts, PetscReal radius) {
484b07a2398SLisandro Dalcin   PetscFunctionBegin;
485b07a2398SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
486b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, radius, 2);
487cad9d221SBarry Smith   PetscCheck(radius >= 0 && radius <= 1, ((PetscObject)ts)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Radius %g not in range [0,1]", (double)radius);
488cac4c232SBarry Smith   PetscTryMethod(ts, "TSAlphaSetRadius_C", (TS, PetscReal), (ts, radius));
489b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
490b07a2398SLisandro Dalcin }
491b07a2398SLisandro Dalcin 
492b07a2398SLisandro Dalcin /*@
493b07a2398SLisandro Dalcin   TSAlphaSetParams - sets the algorithmic parameters for TSALPHA
494b07a2398SLisandro Dalcin 
495b07a2398SLisandro Dalcin   Logically Collective on TS
496b07a2398SLisandro Dalcin 
497b07a2398SLisandro Dalcin   Second-order accuracy can be obtained so long as:
498b07a2398SLisandro Dalcin     \gamma = 0.5 + alpha_m - alpha_f
499b07a2398SLisandro Dalcin 
500b07a2398SLisandro Dalcin   Unconditional stability requires:
501b07a2398SLisandro Dalcin     \alpha_m >= \alpha_f >= 0.5
502b07a2398SLisandro Dalcin 
503b07a2398SLisandro Dalcin   Backward Euler method is recovered with:
504b07a2398SLisandro Dalcin     \alpha_m = \alpha_f = gamma = 1
505b07a2398SLisandro Dalcin 
506d8d19677SJose E. Roman   Input Parameters:
507b07a2398SLisandro Dalcin +  ts - timestepping context
508a5b23f4aSJose E. Roman .  \alpha_m - algorithmic parameter
509a5b23f4aSJose E. Roman .  \alpha_f - algorithmic parameter
510a5b23f4aSJose E. Roman -  \gamma   - algorithmic parameter
511b07a2398SLisandro Dalcin 
512b07a2398SLisandro Dalcin    Options Database:
51367b8a455SSatish Balay +  -ts_alpha_alpha_m <alpha_m> - set alpha_m
51467b8a455SSatish Balay .  -ts_alpha_alpha_f <alpha_f> - set alpha_f
51567b8a455SSatish Balay -  -ts_alpha_gamma   <gamma> - set gamma
516b07a2398SLisandro Dalcin 
517b07a2398SLisandro Dalcin   Note:
518b07a2398SLisandro Dalcin   Use of this function is normally only required to hack TSALPHA to
519b07a2398SLisandro Dalcin   use a modified integration scheme. Users should call
520b07a2398SLisandro Dalcin   TSAlphaSetRadius() to set the desired spectral radius of the methods
521b07a2398SLisandro Dalcin   (i.e. high-frequency damping) in order so select optimal values for
522b07a2398SLisandro Dalcin   these parameters.
523b07a2398SLisandro Dalcin 
524b07a2398SLisandro Dalcin   Level: advanced
525b07a2398SLisandro Dalcin 
526db781477SPatrick Sanan .seealso: `TSAlphaSetRadius()`, `TSAlphaGetParams()`
527b07a2398SLisandro Dalcin @*/
5289371c9d4SSatish Balay PetscErrorCode TSAlphaSetParams(TS ts, PetscReal alpha_m, PetscReal alpha_f, PetscReal gamma) {
529b07a2398SLisandro Dalcin   PetscFunctionBegin;
530b07a2398SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
531b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, alpha_m, 2);
532b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, alpha_f, 3);
533b07a2398SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts, gamma, 4);
534cac4c232SBarry Smith   PetscTryMethod(ts, "TSAlphaSetParams_C", (TS, PetscReal, PetscReal, PetscReal), (ts, alpha_m, alpha_f, gamma));
535b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
536b07a2398SLisandro Dalcin }
537b07a2398SLisandro Dalcin 
538b07a2398SLisandro Dalcin /*@
539b07a2398SLisandro Dalcin   TSAlphaGetParams - gets the algorithmic parameters for TSALPHA
540b07a2398SLisandro Dalcin 
541b07a2398SLisandro Dalcin   Not Collective
542b07a2398SLisandro Dalcin 
543b07a2398SLisandro Dalcin   Input Parameter:
544b07a2398SLisandro Dalcin .  ts - timestepping context
545b07a2398SLisandro Dalcin 
546b07a2398SLisandro Dalcin   Output Parameters:
547b07a2398SLisandro Dalcin +  \alpha_m - algorithmic parameter
548b07a2398SLisandro Dalcin .  \alpha_f - algorithmic parameter
549b07a2398SLisandro Dalcin -  \gamma   - algorithmic parameter
550b07a2398SLisandro Dalcin 
551b07a2398SLisandro Dalcin   Note:
552b07a2398SLisandro Dalcin   Use of this function is normally only required to hack TSALPHA to
553b07a2398SLisandro Dalcin   use a modified integration scheme. Users should call
554b07a2398SLisandro Dalcin   TSAlphaSetRadius() to set the high-frequency damping (i.e. spectral
555b07a2398SLisandro Dalcin   radius of the method) in order so select optimal values for these
556b07a2398SLisandro Dalcin   parameters.
557b07a2398SLisandro Dalcin 
558b07a2398SLisandro Dalcin   Level: advanced
559b07a2398SLisandro Dalcin 
560db781477SPatrick Sanan .seealso: `TSAlphaSetRadius()`, `TSAlphaSetParams()`
561b07a2398SLisandro Dalcin @*/
5629371c9d4SSatish Balay PetscErrorCode TSAlphaGetParams(TS ts, PetscReal *alpha_m, PetscReal *alpha_f, PetscReal *gamma) {
563b07a2398SLisandro Dalcin   PetscFunctionBegin;
564b07a2398SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
565b07a2398SLisandro Dalcin   if (alpha_m) PetscValidRealPointer(alpha_m, 2);
566b07a2398SLisandro Dalcin   if (alpha_f) PetscValidRealPointer(alpha_f, 3);
567b07a2398SLisandro Dalcin   if (gamma) PetscValidRealPointer(gamma, 4);
568cac4c232SBarry Smith   PetscUseMethod(ts, "TSAlphaGetParams_C", (TS, PetscReal *, PetscReal *, PetscReal *), (ts, alpha_m, alpha_f, gamma));
569b07a2398SLisandro Dalcin   PetscFunctionReturn(0);
570b07a2398SLisandro Dalcin }
571