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