xref: /petsc/src/ts/impls/implicit/alpha/alpha1.c (revision 2f6eced2a19e978d64f27de66fbc6c26cc5c7934)
1 /*
2   Code for timestepping with implicit generalized-\alpha method
3   for first order systems.
4 */
5 #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
6 
7 static PetscBool  cited = PETSC_FALSE;
8 static const char citation[] =
9   "@article{Jansen2000,\n"
10   "  title   = {A generalized-$\\alpha$ method for integrating the filtered {N}avier--{S}tokes equations with a stabilized finite element method},\n"
11   "  author  = {Kenneth E. Jansen and Christian H. Whiting and Gregory M. Hulbert},\n"
12   "  journal = {Computer Methods in Applied Mechanics and Engineering},\n"
13   "  volume  = {190},\n"
14   "  number  = {3--4},\n"
15   "  pages   = {305--319},\n"
16   "  year    = {2000},\n"
17   "  issn    = {0045-7825},\n"
18   "  doi     = {http://dx.doi.org/10.1016/S0045-7825(00)00203-6}\n}\n";
19 
20 typedef struct {
21   PetscReal stage_time;
22   PetscReal shift_V;
23   PetscReal scale_F;
24   Vec       X0,Xa,X1;
25   Vec       V0,Va,V1;
26 
27   PetscReal Alpha_m;
28   PetscReal Alpha_f;
29   PetscReal Gamma;
30   PetscInt  order;
31 
32   PetscBool adapt;
33   Vec       vec_sol_prev;
34   Vec       vec_lte_work;
35 
36   TSStepStatus status;
37 } TS_Alpha;
38 
39 static PetscErrorCode TSAlpha_StageTime(TS ts)
40 {
41   TS_Alpha  *th = (TS_Alpha*)ts->data;
42   PetscReal t  = ts->ptime;
43   PetscReal dt = ts->time_step;
44   PetscReal Alpha_m = th->Alpha_m;
45   PetscReal Alpha_f = th->Alpha_f;
46   PetscReal Gamma   = th->Gamma;
47 
48   PetscFunctionBegin;
49   th->stage_time = t + Alpha_f*dt;
50   th->shift_V = Alpha_m/(Alpha_f*Gamma*dt);
51   th->scale_F = 1/Alpha_f;
52   PetscFunctionReturn(0);
53 }
54 
55 static PetscErrorCode TSAlpha_StageVecs(TS ts,Vec X)
56 {
57   TS_Alpha       *th = (TS_Alpha*)ts->data;
58   Vec            X1 = X,      V1 = th->V1;
59   Vec            Xa = th->Xa, Va = th->Va;
60   Vec            X0 = th->X0, V0 = th->V0;
61   PetscReal      dt = ts->time_step;
62   PetscReal      Alpha_m = th->Alpha_m;
63   PetscReal      Alpha_f = th->Alpha_f;
64   PetscReal      Gamma   = th->Gamma;
65   PetscErrorCode ierr;
66 
67   PetscFunctionBegin;
68   /* V1 = 1/(Gamma*dT)*(X1-X0) + (1-1/Gamma)*V0 */
69   ierr = VecWAXPY(V1,-1.0,X0,X1);CHKERRQ(ierr);
70   ierr = VecAXPBY(V1,1-1/Gamma,1/(Gamma*dt),V0);CHKERRQ(ierr);
71   /* Xa = X0 + Alpha_f*(X1-X0) */
72   ierr = VecWAXPY(Xa,-1.0,X0,X1);CHKERRQ(ierr);
73   ierr = VecAYPX(Xa,Alpha_f,X0);CHKERRQ(ierr);
74   /* Va = V0 + Alpha_m*(V1-V0) */
75   ierr = VecWAXPY(Va,-1.0,V0,V1);CHKERRQ(ierr);
76   ierr = VecAYPX(Va,Alpha_m,V0);CHKERRQ(ierr);
77   PetscFunctionReturn(0);
78 }
79 
80 static PetscErrorCode TS_SNESSolve(TS ts,Vec b,Vec x)
81 {
82   PetscInt       nits,lits;
83   PetscErrorCode ierr;
84 
85   PetscFunctionBegin;
86   ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr);
87   ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr);
88   ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
89   ts->snes_its += nits; ts->ksp_its += lits;
90   PetscFunctionReturn(0);
91 }
92 
93 /*
94   Compute a consistent initial state for the generalized-alpha method.
95   - Solve two successive backward Euler steps with halved time step.
96   - Compute the initial time derivative using backward differences.
97   - If using adaptivity, estimate the LTE of the initial step.
98 */
99 static PetscErrorCode TSAlpha_Restart(TS ts,PetscBool *initok)
100 {
101   TS_Alpha       *th = (TS_Alpha*)ts->data;
102   PetscReal      time_step;
103   PetscReal      alpha_m,alpha_f,gamma;
104   Vec            X0 = ts->vec_sol, X1, X2 = th->X1;
105   PetscBool      stageok;
106   PetscErrorCode ierr;
107 
108   PetscFunctionBegin;
109   ierr = VecDuplicate(X0,&X1);CHKERRQ(ierr);
110 
111   /* Setup backward Euler with halved time step */
112   ierr = TSAlphaGetParams(ts,&alpha_m,&alpha_f,&gamma);CHKERRQ(ierr);
113   ierr = TSAlphaSetParams(ts,1,1,1);CHKERRQ(ierr);
114   ierr = TSGetTimeStep(ts,&time_step);CHKERRQ(ierr);
115   ts->time_step = time_step/2;
116   ierr = TSAlpha_StageTime(ts);CHKERRQ(ierr);
117   th->stage_time = ts->ptime;
118   ierr = VecZeroEntries(th->V0);CHKERRQ(ierr);
119 
120   /* First BE step, (t0,X0) -> (t1,X1) */
121   th->stage_time += ts->time_step;
122   ierr = VecCopy(X0,th->X0);CHKERRQ(ierr);
123   ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
124   ierr = VecCopy(th->X0,X1);CHKERRQ(ierr);
125   ierr = TS_SNESSolve(ts,NULL,X1);CHKERRQ(ierr);
126   ierr = TSPostStage(ts,th->stage_time,0,&X1);CHKERRQ(ierr);
127   ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,X1,&stageok);CHKERRQ(ierr);
128   if (!stageok) goto finally;
129 
130   /* Second BE step, (t1,X1) -> (t2,X2) */
131   th->stage_time += ts->time_step;
132   ierr = VecCopy(X1,th->X0);CHKERRQ(ierr);
133   ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
134   ierr = VecCopy(th->X0,X2);CHKERRQ(ierr);
135   ierr = TS_SNESSolve(ts,NULL,X2);CHKERRQ(ierr);
136   ierr = TSPostStage(ts,th->stage_time,0,&X2);CHKERRQ(ierr);
137   ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,X2,&stageok);CHKERRQ(ierr);
138   if (!stageok) goto finally;
139 
140   /* Compute V0 ~ dX/dt at t0 with backward differences */
141   ierr = VecZeroEntries(th->V0);CHKERRQ(ierr);
142   ierr = VecAXPY(th->V0,-3/ts->time_step,X0);CHKERRQ(ierr);
143   ierr = VecAXPY(th->V0,+4/ts->time_step,X1);CHKERRQ(ierr);
144   ierr = VecAXPY(th->V0,-1/ts->time_step,X2);CHKERRQ(ierr);
145 
146   /* Rough, lower-order estimate LTE of the initial step */
147   if (th->adapt) {
148     ierr = VecZeroEntries(th->vec_lte_work);CHKERRQ(ierr);
149     ierr = VecAXPY(th->vec_lte_work,+2,X2);CHKERRQ(ierr);
150     ierr = VecAXPY(th->vec_lte_work,-4,X1);CHKERRQ(ierr);
151     ierr = VecAXPY(th->vec_lte_work,+2,X0);CHKERRQ(ierr);
152   }
153 
154  finally:
155   /* Revert TSAlpha to the initial state (t0,X0) */
156   if (initok) *initok = stageok;
157   ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr);
158   ierr = TSAlphaSetParams(ts,alpha_m,alpha_f,gamma);CHKERRQ(ierr);
159   ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
160 
161   ierr = VecDestroy(&X1);CHKERRQ(ierr);
162   PetscFunctionReturn(0);
163 }
164 
165 static PetscErrorCode TSStep_Alpha(TS ts)
166 {
167   TS_Alpha       *th = (TS_Alpha*)ts->data;
168   PetscInt       rejections = 0;
169   PetscBool      stageok,accept = PETSC_TRUE;
170   PetscReal      next_time_step = ts->time_step;
171   PetscErrorCode ierr;
172 
173   PetscFunctionBegin;
174   ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr);
175 
176   if (!ts->steprollback) {
177     if (th->adapt) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); }
178     ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
179     ierr = VecCopy(th->V1,th->V0);CHKERRQ(ierr);
180   }
181 
182   th->status = TS_STEP_INCOMPLETE;
183   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
184 
185     if (ts->steprestart) {
186       ierr = TSAlpha_Restart(ts,&stageok);CHKERRQ(ierr);
187       if (!stageok) goto reject_step;
188     }
189 
190     ierr = TSAlpha_StageTime(ts);CHKERRQ(ierr);
191     ierr = VecCopy(th->X0,th->X1);CHKERRQ(ierr);
192     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
193     ierr = TS_SNESSolve(ts,NULL,th->X1);CHKERRQ(ierr);
194     ierr = TSPostStage(ts,th->stage_time,0,&th->Xa);CHKERRQ(ierr);
195     ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->Xa,&stageok);CHKERRQ(ierr);
196     if (!stageok) goto reject_step;
197 
198     th->status = TS_STEP_PENDING;
199     ierr = VecCopy(th->X1,ts->vec_sol);CHKERRQ(ierr);
200     ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
201     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
202     if (!accept) {
203       ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
204       ts->time_step = next_time_step;
205       goto reject_step;
206     }
207 
208     ts->ptime += ts->time_step;
209     ts->time_step = next_time_step;
210     break;
211 
212   reject_step:
213     ts->reject++; accept = PETSC_FALSE;
214     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
215       ts->reason = TS_DIVERGED_STEP_REJECTED;
216       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
217     }
218 
219   }
220   PetscFunctionReturn(0);
221 }
222 
223 static PetscErrorCode TSEvaluateWLTE_Alpha(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
224 {
225   TS_Alpha       *th = (TS_Alpha*)ts->data;
226   Vec            X = th->X1;           /* X = solution */
227   Vec            Y = th->vec_lte_work; /* Y = X + LTE  */
228   PetscReal      wltea,wlter;
229   PetscErrorCode ierr;
230 
231   PetscFunctionBegin;
232   if (ts->steprestart) {
233     /* th->vec_lte_work is set to the LTE in TSAlpha_Restart() */
234     ierr = VecAXPY(Y,1,X);CHKERRQ(ierr);
235   } else {
236     /* Compute LTE using backward differences with non-constant time step */
237     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
238     PetscReal   a = 1 + h_prev/h;
239     PetscScalar scal[3]; Vec vecs[3];
240     scal[0] = +1/a;   scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1));
241     vecs[0] = th->X1; vecs[1] = th->X0;   vecs[2] = th->vec_sol_prev;
242     ierr = VecCopy(X,Y);CHKERRQ(ierr);
243     ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr);
244   }
245   ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr);
246   if (order) *order = 2;
247   PetscFunctionReturn(0);
248 }
249 
250 static PetscErrorCode TSRollBack_Alpha(TS ts)
251 {
252   TS_Alpha       *th = (TS_Alpha*)ts->data;
253   PetscErrorCode ierr;
254 
255   PetscFunctionBegin;
256   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
257   PetscFunctionReturn(0);
258 }
259 
260 static PetscErrorCode TSInterpolate_Alpha(TS ts,PetscReal t,Vec X)
261 {
262   TS_Alpha       *th = (TS_Alpha*)ts->data;
263   PetscReal      dt  = t - ts->ptime;
264   PetscErrorCode ierr;
265 
266   PetscFunctionBegin;
267   ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
268   ierr = VecAXPY(X,th->Gamma*dt,th->V1);CHKERRQ(ierr);
269   ierr = VecAXPY(X,(1-th->Gamma)*dt,th->V0);CHKERRQ(ierr);
270   PetscFunctionReturn(0);
271 }
272 
273 static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes,Vec X,Vec F,TS ts)
274 {
275   TS_Alpha       *th = (TS_Alpha*)ts->data;
276   PetscReal      ta = th->stage_time;
277   Vec            Xa = th->Xa, Va = th->Va;
278   PetscErrorCode ierr;
279 
280   PetscFunctionBegin;
281   ierr = TSAlpha_StageVecs(ts,X);CHKERRQ(ierr);
282   /* F = Function(ta,Xa,Va) */
283   ierr = TSComputeIFunction(ts,ta,Xa,Va,F,PETSC_FALSE);CHKERRQ(ierr);
284   ierr = VecScale(F,th->scale_F);CHKERRQ(ierr);
285   PetscFunctionReturn(0);
286 }
287 
288 static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes,PETSC_UNUSED Vec X,Mat J,Mat P,TS ts)
289 {
290   TS_Alpha       *th = (TS_Alpha*)ts->data;
291   PetscReal      ta = th->stage_time;
292   Vec            Xa = th->Xa, Va = th->Va;
293   PetscReal      dVdX = th->shift_V;
294   PetscErrorCode ierr;
295 
296   PetscFunctionBegin;
297   /* J,P = Jacobian(ta,Xa,Va) */
298   ierr = TSComputeIJacobian(ts,ta,Xa,Va,dVdX,J,P,PETSC_FALSE);CHKERRQ(ierr);
299   PetscFunctionReturn(0);
300 }
301 
302 static PetscErrorCode TSReset_Alpha(TS ts)
303 {
304   TS_Alpha       *th = (TS_Alpha*)ts->data;
305   PetscErrorCode ierr;
306 
307   PetscFunctionBegin;
308   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
309   ierr = VecDestroy(&th->Xa);CHKERRQ(ierr);
310   ierr = VecDestroy(&th->X1);CHKERRQ(ierr);
311   ierr = VecDestroy(&th->V0);CHKERRQ(ierr);
312   ierr = VecDestroy(&th->Va);CHKERRQ(ierr);
313   ierr = VecDestroy(&th->V1);CHKERRQ(ierr);
314   ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr);
315   ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr);
316   PetscFunctionReturn(0);
317 }
318 
319 static PetscErrorCode TSDestroy_Alpha(TS ts)
320 {
321   PetscErrorCode ierr;
322 
323   PetscFunctionBegin;
324   ierr = TSReset_Alpha(ts);CHKERRQ(ierr);
325   ierr = PetscFree(ts->data);CHKERRQ(ierr);
326 
327   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaUseAdapt_C",NULL);CHKERRQ(ierr);
328   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetRadius_C",NULL);CHKERRQ(ierr);
329   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetParams_C",NULL);CHKERRQ(ierr);
330   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaGetParams_C",NULL);CHKERRQ(ierr);
331   PetscFunctionReturn(0);
332 }
333 
334 static PetscErrorCode TSSetUp_Alpha(TS ts)
335 {
336   TS_Alpha       *th = (TS_Alpha*)ts->data;
337   PetscErrorCode ierr;
338 
339   PetscFunctionBegin;
340   ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
341   ierr = VecDuplicate(ts->vec_sol,&th->Xa);CHKERRQ(ierr);
342   ierr = VecDuplicate(ts->vec_sol,&th->X1);CHKERRQ(ierr);
343   ierr = VecDuplicate(ts->vec_sol,&th->V0);CHKERRQ(ierr);
344   ierr = VecDuplicate(ts->vec_sol,&th->Va);CHKERRQ(ierr);
345   ierr = VecDuplicate(ts->vec_sol,&th->V1);CHKERRQ(ierr);
346 
347   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
348   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
349   if (!th->adapt) {
350     ierr = TSAdaptSetType(ts->adapt,TSADAPTNONE);CHKERRQ(ierr);
351   } else {
352     ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr);
353     ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr);
354     if (ts->exact_final_time == TS_EXACTFINALTIME_UNSPECIFIED)
355       ts->exact_final_time = TS_EXACTFINALTIME_MATCHSTEP;
356   }
357 
358   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
359   PetscFunctionReturn(0);
360 }
361 
362 static PetscErrorCode TSSetFromOptions_Alpha(PetscOptionItems *PetscOptionsObject,TS ts)
363 {
364   TS_Alpha       *th = (TS_Alpha*)ts->data;
365   PetscErrorCode ierr;
366 
367   PetscFunctionBegin;
368   ierr = PetscOptionsHead(PetscOptionsObject,"Generalized-Alpha ODE solver options");CHKERRQ(ierr);
369   {
370     PetscBool flg;
371     PetscReal radius = 1;
372     PetscBool adapt  = th->adapt;
373     ierr = PetscOptionsReal("-ts_alpha_radius","Spectral radius (high-frequency dissipation)","TSAlphaSetRadius",radius,&radius,&flg);CHKERRQ(ierr);
374     if (flg) {ierr = TSAlphaSetRadius(ts,radius);CHKERRQ(ierr);}
375     ierr = PetscOptionsReal("-ts_alpha_alpha_m","Algoritmic parameter alpha_m","TSAlphaSetParams",th->Alpha_m,&th->Alpha_m,NULL);CHKERRQ(ierr);
376     ierr = PetscOptionsReal("-ts_alpha_alpha_f","Algoritmic parameter alpha_f","TSAlphaSetParams",th->Alpha_f,&th->Alpha_f,NULL);CHKERRQ(ierr);
377     ierr = PetscOptionsReal("-ts_alpha_gamma","Algoritmic parameter gamma","TSAlphaSetParams",th->Gamma,&th->Gamma,NULL);CHKERRQ(ierr);
378     ierr = TSAlphaSetParams(ts,th->Alpha_m,th->Alpha_f,th->Gamma);CHKERRQ(ierr);
379     ierr = PetscOptionsBool("-ts_alpha_adapt","Use time-step adaptivity with the Alpha method","TSAlpha2UseAdapt",adapt,&adapt,&flg);CHKERRQ(ierr);
380     if (flg) {ierr = TSAlphaUseAdapt(ts,adapt);CHKERRQ(ierr);}
381   }
382   ierr = PetscOptionsTail();CHKERRQ(ierr);
383   PetscFunctionReturn(0);
384 }
385 
386 static PetscErrorCode TSView_Alpha(TS ts,PetscViewer viewer)
387 {
388   TS_Alpha       *th = (TS_Alpha*)ts->data;
389   PetscBool      iascii;
390   PetscErrorCode ierr;
391 
392   PetscFunctionBegin;
393   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
394   if (iascii) {
395     ierr = PetscViewerASCIIPrintf(viewer,"  Alpha_m=%g, Alpha_f=%g, Gamma=%g\n",(double)th->Alpha_m,(double)th->Alpha_f,(double)th->Gamma);CHKERRQ(ierr);
396   }
397   if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);}
398   if (ts->snes)  {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
399   PetscFunctionReturn(0);
400 }
401 
402 static PetscErrorCode TSAlphaUseAdapt_Alpha(TS ts,PetscBool use)
403 {
404   TS_Alpha *th = (TS_Alpha*)ts->data;
405 
406   PetscFunctionBegin;
407   if (use == th->adapt) PetscFunctionReturn(0);
408   if (ts->setupcalled) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ORDER,"Cannot change adaptivity after TSSetUp()");
409   th->adapt = use;
410   PetscFunctionReturn(0);
411 }
412 
413 static PetscErrorCode TSAlphaSetRadius_Alpha(TS ts,PetscReal radius)
414 {
415   PetscReal      alpha_m,alpha_f,gamma;
416   PetscErrorCode ierr;
417 
418   PetscFunctionBegin;
419   if (radius < 0 || radius > 1) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Radius %g not in range [0,1]",(double)radius);
420   alpha_m = (PetscReal)0.5*(3-radius)/(1+radius);
421   alpha_f = 1/(1+radius);
422   gamma   = (PetscReal)0.5 + alpha_m - alpha_f;
423   ierr = TSAlphaSetParams(ts,alpha_m,alpha_f,gamma);CHKERRQ(ierr);
424   PetscFunctionReturn(0);
425 }
426 
427 static PetscErrorCode TSAlphaSetParams_Alpha(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma)
428 {
429   TS_Alpha  *th = (TS_Alpha*)ts->data;
430   PetscReal tol = 100*PETSC_MACHINE_EPSILON;
431   PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma;
432 
433   PetscFunctionBegin;
434   th->Alpha_m = alpha_m;
435   th->Alpha_f = alpha_f;
436   th->Gamma   = gamma;
437   th->order   = (PetscAbsReal(res) < tol) ? 2 : 1;
438   PetscFunctionReturn(0);
439 }
440 
441 static PetscErrorCode TSAlphaGetParams_Alpha(TS ts,PetscReal *alpha_m,PetscReal *alpha_f,PetscReal *gamma)
442 {
443   TS_Alpha *th = (TS_Alpha*)ts->data;
444 
445   PetscFunctionBegin;
446   if (alpha_m) *alpha_m = th->Alpha_m;
447   if (alpha_f) *alpha_f = th->Alpha_f;
448   if (gamma)   *gamma   = th->Gamma;
449   PetscFunctionReturn(0);
450 }
451 
452 /*MC
453       TSALPHA - ODE/DAE solver using the implicit Generalized-Alpha method
454                 for first-order systems
455 
456   Level: beginner
457 
458   References:
459   K.E. Jansen, C.H. Whiting, G.M. Hulber, "A generalized-alpha
460   method for integrating the filtered Navier-Stokes equations with a
461   stabilized finite element method", Computer Methods in Applied
462   Mechanics and Engineering, 190, 305-319, 2000.
463   DOI: 10.1016/S0045-7825(00)00203-6.
464 
465   J. Chung, G.M.Hubert. "A Time Integration Algorithm for Structural
466   Dynamics with Improved Numerical Dissipation: The Generalized-alpha
467   Method" ASME Journal of Applied Mechanics, 60, 371:375, 1993.
468 
469 .seealso:  TS, TSCreate(), TSSetType(), TSAlphaSetRadius(), TSAlphaSetParams()
470 M*/
471 PETSC_EXTERN PetscErrorCode TSCreate_Alpha(TS ts)
472 {
473   TS_Alpha       *th;
474   PetscErrorCode ierr;
475 
476   PetscFunctionBegin;
477   ts->ops->reset          = TSReset_Alpha;
478   ts->ops->destroy        = TSDestroy_Alpha;
479   ts->ops->view           = TSView_Alpha;
480   ts->ops->setup          = TSSetUp_Alpha;
481   ts->ops->setfromoptions = TSSetFromOptions_Alpha;
482   ts->ops->step           = TSStep_Alpha;
483   ts->ops->evaluatewlte   = TSEvaluateWLTE_Alpha;
484   ts->ops->rollback       = TSRollBack_Alpha;
485   ts->ops->interpolate    = TSInterpolate_Alpha;
486   ts->ops->snesfunction   = SNESTSFormFunction_Alpha;
487   ts->ops->snesjacobian   = SNESTSFormJacobian_Alpha;
488 
489   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
490   ts->data = (void*)th;
491 
492   th->Alpha_m = 0.5;
493   th->Alpha_f = 0.5;
494   th->Gamma   = 0.5;
495   th->order   = 2;
496 
497   th->adapt = PETSC_FALSE;
498 
499   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaUseAdapt_C",TSAlphaUseAdapt_Alpha);CHKERRQ(ierr);
500   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetRadius_C",TSAlphaSetRadius_Alpha);CHKERRQ(ierr);
501   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaSetParams_C",TSAlphaSetParams_Alpha);CHKERRQ(ierr);
502   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlphaGetParams_C",TSAlphaGetParams_Alpha);CHKERRQ(ierr);
503   PetscFunctionReturn(0);
504 }
505 
506 /*@
507   TSAlphaUseAdapt - Use time-step adaptivity with the Alpha method
508 
509   Logically Collective on TS
510 
511   Input Parameter:
512 +  ts - timestepping context
513 -  use - flag to use adaptivity
514 
515   Options Database:
516 .  -ts_alpha_adapt
517 
518   Level: intermediate
519 
520 .seealso: TSAdapt, TSADAPTBASIC
521 @*/
522 PetscErrorCode TSAlphaUseAdapt(TS ts,PetscBool use)
523 {
524   PetscErrorCode ierr;
525 
526   PetscFunctionBegin;
527   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
528   PetscValidLogicalCollectiveBool(ts,use,2);
529   ierr = PetscTryMethod(ts,"TSAlphaUseAdapt_C",(TS,PetscBool),(ts,use));CHKERRQ(ierr);
530   PetscFunctionReturn(0);
531 }
532 
533 /*@
534   TSAlphaSetRadius - sets the desired spectral radius of the method
535                      (i.e. high-frequency numerical damping)
536 
537   Logically Collective on TS
538 
539   The algorithmic parameters \alpha_m and \alpha_f of the
540   generalized-\alpha method can be computed in terms of a specified
541   spectral radius \rho in [0,1] for infinite time step in order to
542   control high-frequency numerical damping:
543     \alpha_m = 0.5*(3-\rho)/(1+\rho)
544     \alpha_f = 1/(1+\rho)
545 
546   Input Parameter:
547 +  ts - timestepping context
548 -  radius - the desired spectral radius
549 
550   Options Database:
551 .  -ts_alpha_radius <radius>
552 
553   Level: intermediate
554 
555 .seealso: TSAlphaSetParams(), TSAlphaGetParams()
556 @*/
557 PetscErrorCode TSAlphaSetRadius(TS ts,PetscReal radius)
558 {
559   PetscErrorCode ierr;
560 
561   PetscFunctionBegin;
562   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
563   PetscValidLogicalCollectiveReal(ts,radius,2);
564   if (radius < 0 || radius > 1) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Radius %g not in range [0,1]",(double)radius);
565   ierr = PetscTryMethod(ts,"TSAlphaSetRadius_C",(TS,PetscReal),(ts,radius));CHKERRQ(ierr);
566   PetscFunctionReturn(0);
567 }
568 
569 /*@
570   TSAlphaSetParams - sets the algorithmic parameters for TSALPHA
571 
572   Logically Collective on TS
573 
574   Second-order accuracy can be obtained so long as:
575     \gamma = 0.5 + alpha_m - alpha_f
576 
577   Unconditional stability requires:
578     \alpha_m >= \alpha_f >= 0.5
579 
580   Backward Euler method is recovered with:
581     \alpha_m = \alpha_f = gamma = 1
582 
583   Input Parameter:
584 +  ts - timestepping context
585 .  \alpha_m - algorithmic paramenter
586 .  \alpha_f - algorithmic paramenter
587 -  \gamma   - algorithmic paramenter
588 
589    Options Database:
590 +  -ts_alpha_alpha_m <alpha_m>
591 .  -ts_alpha_alpha_f <alpha_f>
592 -  -ts_alpha_gamma   <gamma>
593 
594   Note:
595   Use of this function is normally only required to hack TSALPHA to
596   use a modified integration scheme. Users should call
597   TSAlphaSetRadius() to set the desired spectral radius of the methods
598   (i.e. high-frequency damping) in order so select optimal values for
599   these parameters.
600 
601   Level: advanced
602 
603 .seealso: TSAlphaSetRadius(), TSAlphaGetParams()
604 @*/
605 PetscErrorCode TSAlphaSetParams(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma)
606 {
607   PetscErrorCode ierr;
608 
609   PetscFunctionBegin;
610   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
611   PetscValidLogicalCollectiveReal(ts,alpha_m,2);
612   PetscValidLogicalCollectiveReal(ts,alpha_f,3);
613   PetscValidLogicalCollectiveReal(ts,gamma,4);
614   ierr = PetscTryMethod(ts,"TSAlphaSetParams_C",(TS,PetscReal,PetscReal,PetscReal),(ts,alpha_m,alpha_f,gamma));CHKERRQ(ierr);
615   PetscFunctionReturn(0);
616 }
617 
618 /*@
619   TSAlphaGetParams - gets the algorithmic parameters for TSALPHA
620 
621   Not Collective
622 
623   Input Parameter:
624 .  ts - timestepping context
625 
626   Output Parameters:
627 +  \alpha_m - algorithmic parameter
628 .  \alpha_f - algorithmic parameter
629 -  \gamma   - algorithmic parameter
630 
631   Note:
632   Use of this function is normally only required to hack TSALPHA to
633   use a modified integration scheme. Users should call
634   TSAlphaSetRadius() to set the high-frequency damping (i.e. spectral
635   radius of the method) in order so select optimal values for these
636   parameters.
637 
638   Level: advanced
639 
640 .seealso: TSAlphaSetRadius(), TSAlphaSetParams()
641 @*/
642 PetscErrorCode TSAlphaGetParams(TS ts,PetscReal *alpha_m,PetscReal *alpha_f,PetscReal *gamma)
643 {
644   PetscErrorCode ierr;
645 
646   PetscFunctionBegin;
647   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
648   if (alpha_m) PetscValidRealPointer(alpha_m,2);
649   if (alpha_f) PetscValidRealPointer(alpha_f,3);
650   if (gamma)   PetscValidRealPointer(gamma,4);
651   ierr = PetscUseMethod(ts,"TSAlphaGetParams_C",(TS,PetscReal*,PetscReal*,PetscReal*),(ts,alpha_m,alpha_f,gamma));CHKERRQ(ierr);
652   PetscFunctionReturn(0);
653 }
654