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