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