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