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