xref: /petsc/src/ts/impls/implicit/alpha/alpha2.c (revision 7453f77509fc006cbcc7de2dd772f60dc49feac5)
1818efac9SLisandro Dalcin /*
2818efac9SLisandro Dalcin   Code for timestepping with implicit generalized-\alpha method
3818efac9SLisandro Dalcin   for second order systems.
4818efac9SLisandro Dalcin */
5818efac9SLisandro Dalcin #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
6818efac9SLisandro Dalcin 
7818efac9SLisandro Dalcin static PetscBool  cited = PETSC_FALSE;
8818efac9SLisandro Dalcin static const char citation[] =
9818efac9SLisandro Dalcin   "@article{Chung1993,\n"
10818efac9SLisandro Dalcin   "  title   = {A Time Integration Algorithm for Structural Dynamics with Improved Numerical Dissipation: The Generalized-$\\alpha$ Method},\n"
11818efac9SLisandro Dalcin   "  author  = {J. Chung, G. M. Hubert},\n"
12818efac9SLisandro Dalcin   "  journal = {ASME Journal of Applied Mechanics},\n"
13818efac9SLisandro Dalcin   "  volume  = {60},\n"
14818efac9SLisandro Dalcin   "  number  = {2},\n"
15818efac9SLisandro Dalcin   "  pages   = {371--375},\n"
16818efac9SLisandro Dalcin   "  year    = {1993},\n"
17818efac9SLisandro Dalcin   "  issn    = {0021-8936},\n"
18818efac9SLisandro Dalcin   "  doi     = {http://dx.doi.org/10.1115/1.2900803}\n}\n";
19818efac9SLisandro Dalcin 
20818efac9SLisandro Dalcin typedef struct {
21818efac9SLisandro Dalcin   PetscReal stage_time;
22818efac9SLisandro Dalcin   PetscReal shift_V;
23818efac9SLisandro Dalcin   PetscReal shift_A;
24818efac9SLisandro Dalcin   PetscReal scale_F;
25818efac9SLisandro Dalcin   Vec       X0,Xa,X1;
26818efac9SLisandro Dalcin   Vec       V0,Va,V1;
27818efac9SLisandro Dalcin   Vec       A0,Aa,A1;
28818efac9SLisandro Dalcin 
29818efac9SLisandro Dalcin   Vec       vec_dot;
30818efac9SLisandro Dalcin 
31818efac9SLisandro Dalcin   PetscReal Alpha_m;
32818efac9SLisandro Dalcin   PetscReal Alpha_f;
33818efac9SLisandro Dalcin   PetscReal Gamma;
34818efac9SLisandro Dalcin   PetscReal Beta;
35818efac9SLisandro Dalcin   PetscInt  order;
36818efac9SLisandro Dalcin 
37818efac9SLisandro Dalcin   PetscBool adapt;
38818efac9SLisandro Dalcin   Vec       vec_sol_prev;
39818efac9SLisandro Dalcin   Vec       vec_dot_prev;
40818efac9SLisandro Dalcin   Vec       vec_lte_work[2];
41818efac9SLisandro Dalcin 
42818efac9SLisandro Dalcin   TSStepStatus status;
43818efac9SLisandro Dalcin } TS_Alpha;
44818efac9SLisandro Dalcin 
45818efac9SLisandro Dalcin #undef __FUNCT__
46818efac9SLisandro Dalcin #define __FUNCT__ "TSAlpha_StageTime"
47818efac9SLisandro Dalcin static PetscErrorCode TSAlpha_StageTime(TS ts)
48818efac9SLisandro Dalcin {
49818efac9SLisandro Dalcin   TS_Alpha  *th = (TS_Alpha*)ts->data;
50818efac9SLisandro Dalcin   PetscReal t  = ts->ptime;
51818efac9SLisandro Dalcin   PetscReal dt = ts->time_step;
52818efac9SLisandro Dalcin   PetscReal Alpha_m = th->Alpha_m;
53818efac9SLisandro Dalcin   PetscReal Alpha_f = th->Alpha_f;
54818efac9SLisandro Dalcin   PetscReal Gamma   = th->Gamma;
55818efac9SLisandro Dalcin   PetscReal Beta    = th->Beta;
56818efac9SLisandro Dalcin 
57818efac9SLisandro Dalcin   PetscFunctionBegin;
58818efac9SLisandro Dalcin   th->stage_time = t + Alpha_f*dt;
59818efac9SLisandro Dalcin   th->shift_V = Gamma/(dt*Beta);
60818efac9SLisandro Dalcin   th->shift_A = Alpha_m/(Alpha_f*dt*dt*Beta);
61818efac9SLisandro Dalcin   th->scale_F = 1/Alpha_f;
62818efac9SLisandro Dalcin   PetscFunctionReturn(0);
63818efac9SLisandro Dalcin }
64818efac9SLisandro Dalcin 
65818efac9SLisandro Dalcin #undef __FUNCT__
66818efac9SLisandro Dalcin #define __FUNCT__ "TSAlpha_StageVecs"
67818efac9SLisandro Dalcin static PetscErrorCode TSAlpha_StageVecs(TS ts,Vec X)
68818efac9SLisandro Dalcin {
69818efac9SLisandro Dalcin   TS_Alpha       *th = (TS_Alpha*)ts->data;
70818efac9SLisandro Dalcin   Vec            X1 = X,      V1 = th->V1, A1 = th->A1;
71818efac9SLisandro Dalcin   Vec            Xa = th->Xa, Va = th->Va, Aa = th->Aa;
72818efac9SLisandro Dalcin   Vec            X0 = th->X0, V0 = th->V0, A0 = th->A0;
73818efac9SLisandro Dalcin   PetscReal      dt = ts->time_step;
74818efac9SLisandro Dalcin   PetscReal      Alpha_m = th->Alpha_m;
75818efac9SLisandro Dalcin   PetscReal      Alpha_f = th->Alpha_f;
76818efac9SLisandro Dalcin   PetscReal      Gamma   = th->Gamma;
77818efac9SLisandro Dalcin   PetscReal      Beta    = th->Beta;
78818efac9SLisandro Dalcin   PetscErrorCode ierr;
79818efac9SLisandro Dalcin 
80818efac9SLisandro Dalcin   PetscFunctionBegin;
81818efac9SLisandro Dalcin   /* A1 = ... */
82818efac9SLisandro Dalcin   ierr = VecWAXPY(A1,-1.0,X0,X1);CHKERRQ(ierr);
83818efac9SLisandro Dalcin   ierr = VecAXPY (A1,-dt,V0);CHKERRQ(ierr);
84818efac9SLisandro Dalcin   ierr = VecAXPBY(A1,-(1-2*Beta)/(2*Beta),1/(dt*dt*Beta),A0);CHKERRQ(ierr);
85818efac9SLisandro Dalcin   /* V1 = ... */
86818efac9SLisandro Dalcin   ierr = VecWAXPY(V1,(1.0-Gamma)/Gamma,A0,A1);CHKERRQ(ierr);
87818efac9SLisandro Dalcin   ierr = VecAYPX (V1,dt*Gamma,V0);CHKERRQ(ierr);
88818efac9SLisandro Dalcin   /* Xa = X0 + Alpha_f*(X1-X0) */
89818efac9SLisandro Dalcin   ierr = VecWAXPY(Xa,-1.0,X0,X1);CHKERRQ(ierr);
90818efac9SLisandro Dalcin   ierr = VecAYPX (Xa,Alpha_f,X0);CHKERRQ(ierr);
91818efac9SLisandro Dalcin   /* Va = V0 + Alpha_f*(V1-V0) */
92818efac9SLisandro Dalcin   ierr = VecWAXPY(Va,-1.0,V0,V1);CHKERRQ(ierr);
93818efac9SLisandro Dalcin   ierr = VecAYPX (Va,Alpha_f,V0);CHKERRQ(ierr);
94818efac9SLisandro Dalcin   /* Aa = A0 + Alpha_m*(A1-A0) */
95818efac9SLisandro Dalcin   ierr = VecWAXPY(Aa,-1.0,A0,A1);CHKERRQ(ierr);
96818efac9SLisandro Dalcin   ierr = VecAYPX (Aa,Alpha_m,A0);CHKERRQ(ierr);
97818efac9SLisandro Dalcin   PetscFunctionReturn(0);
98818efac9SLisandro Dalcin }
99818efac9SLisandro Dalcin 
100818efac9SLisandro Dalcin #undef __FUNCT__
101818efac9SLisandro Dalcin #define __FUNCT__ "TS_SNESSolve"
102818efac9SLisandro Dalcin static PetscErrorCode TS_SNESSolve(TS ts,Vec b,Vec x)
103818efac9SLisandro Dalcin {
104818efac9SLisandro Dalcin   PetscInt       nits,lits;
105818efac9SLisandro Dalcin   PetscErrorCode ierr;
106818efac9SLisandro Dalcin 
107818efac9SLisandro Dalcin   PetscFunctionBegin;
108818efac9SLisandro Dalcin   ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr);
109818efac9SLisandro Dalcin   ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr);
110818efac9SLisandro Dalcin   ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
111818efac9SLisandro Dalcin   ts->snes_its += nits; ts->ksp_its += lits;
112818efac9SLisandro Dalcin   PetscFunctionReturn(0);
113818efac9SLisandro Dalcin }
114818efac9SLisandro Dalcin 
115818efac9SLisandro Dalcin /*
116818efac9SLisandro Dalcin   Compute a consistent initial state for the generalized-alpha method.
117818efac9SLisandro Dalcin   - Solve two successive backward Euler steps with halved time step.
118818efac9SLisandro Dalcin   - Compute the initial second time derivative using backward differences.
119818efac9SLisandro Dalcin   - If using adaptivity, estimate the LTE of the initial step.
120818efac9SLisandro Dalcin */
121818efac9SLisandro Dalcin #undef __FUNCT__
122818efac9SLisandro Dalcin #define __FUNCT__ "TSAlpha_Restart"
123818efac9SLisandro Dalcin static PetscErrorCode TSAlpha_Restart(TS ts,PetscBool *initok)
124818efac9SLisandro Dalcin {
125818efac9SLisandro Dalcin   TS_Alpha       *th = (TS_Alpha*)ts->data;
126818efac9SLisandro Dalcin   PetscReal      time_step;
127818efac9SLisandro Dalcin   PetscReal      alpha_m,alpha_f,gamma,beta;
128818efac9SLisandro Dalcin   Vec            X0 = ts->vec_sol, X1, X2 = th->X1;
129818efac9SLisandro Dalcin   Vec            V0 = ts->vec_dot, V1, V2 = th->V1;
130818efac9SLisandro Dalcin   PetscBool      stageok;
131818efac9SLisandro Dalcin   PetscErrorCode ierr;
132818efac9SLisandro Dalcin 
133818efac9SLisandro Dalcin   PetscFunctionBegin;
134818efac9SLisandro Dalcin   ierr = VecDuplicate(X0,&X1);CHKERRQ(ierr);
135818efac9SLisandro Dalcin   ierr = VecDuplicate(V0,&V1);CHKERRQ(ierr);
136818efac9SLisandro Dalcin 
137818efac9SLisandro Dalcin   /* Setup backward Euler with halved time step */
138818efac9SLisandro Dalcin   ierr = TSAlpha2GetParams(ts,&alpha_m,&alpha_f,&gamma,&beta);CHKERRQ(ierr);
139818efac9SLisandro Dalcin   ierr = TSAlpha2SetParams(ts,1,1,1,0.5);CHKERRQ(ierr);
140818efac9SLisandro Dalcin   ierr = TSGetTimeStep(ts,&time_step);CHKERRQ(ierr);
141818efac9SLisandro Dalcin   ts->time_step = time_step/2;
142818efac9SLisandro Dalcin   ierr = TSAlpha_StageTime(ts);CHKERRQ(ierr);
143818efac9SLisandro Dalcin   th->stage_time = ts->ptime;
144818efac9SLisandro Dalcin   ierr = VecZeroEntries(th->A0);CHKERRQ(ierr);
145818efac9SLisandro Dalcin 
146818efac9SLisandro Dalcin   /* First BE step, (t0,X0,V0) -> (t1,X1,V1) */
147818efac9SLisandro Dalcin   th->stage_time += ts->time_step;
148818efac9SLisandro Dalcin   ierr = VecCopy(X0,th->X0);CHKERRQ(ierr);
149818efac9SLisandro Dalcin   ierr = VecCopy(V0,th->V0);CHKERRQ(ierr);
150818efac9SLisandro Dalcin   ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
151818efac9SLisandro Dalcin   ierr = VecCopy(th->X0,X1);CHKERRQ(ierr);
152818efac9SLisandro Dalcin   ierr = TS_SNESSolve(ts,NULL,X1);CHKERRQ(ierr);
153818efac9SLisandro Dalcin   ierr = VecCopy(th->V1,V1);CHKERRQ(ierr);
154818efac9SLisandro Dalcin   ierr = TSPostStage(ts,th->stage_time,0,&X1);CHKERRQ(ierr);
155818efac9SLisandro Dalcin   ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,X1,&stageok);CHKERRQ(ierr);
156818efac9SLisandro Dalcin   if (!stageok) goto finally;
157818efac9SLisandro Dalcin 
158818efac9SLisandro Dalcin   /* Second BE step, (t1,X1,V1) -> (t2,X2,V2) */
159818efac9SLisandro Dalcin   th->stage_time += ts->time_step;
160818efac9SLisandro Dalcin   ierr = VecCopy(X1,th->X0);CHKERRQ(ierr);
161818efac9SLisandro Dalcin   ierr = VecCopy(V1,th->V0);CHKERRQ(ierr);
162818efac9SLisandro Dalcin   ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
163818efac9SLisandro Dalcin   ierr = VecCopy(th->X0,X2);CHKERRQ(ierr);
164818efac9SLisandro Dalcin   ierr = TS_SNESSolve(ts,NULL,X2);CHKERRQ(ierr);
165818efac9SLisandro Dalcin   ierr = VecCopy(th->V1,V2);CHKERRQ(ierr);
166818efac9SLisandro Dalcin   ierr = TSPostStage(ts,th->stage_time,0,&X2);CHKERRQ(ierr);
167818efac9SLisandro Dalcin   ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,X1,&stageok);CHKERRQ(ierr);
168818efac9SLisandro Dalcin   if (!stageok) goto finally;
169818efac9SLisandro Dalcin 
170818efac9SLisandro Dalcin   /* Compute A0 ~ dV/dt at t0 with backward differences */
171818efac9SLisandro Dalcin   ierr = VecZeroEntries(th->A0);CHKERRQ(ierr);
172818efac9SLisandro Dalcin   ierr = VecAXPY(th->A0,-3/ts->time_step,V0);CHKERRQ(ierr);
173818efac9SLisandro Dalcin   ierr = VecAXPY(th->A0,+4/ts->time_step,V1);CHKERRQ(ierr);
174818efac9SLisandro Dalcin   ierr = VecAXPY(th->A0,-1/ts->time_step,V2);CHKERRQ(ierr);
175818efac9SLisandro Dalcin 
176818efac9SLisandro Dalcin   /* Rough, lower-order estimate LTE of the initial step */
177818efac9SLisandro Dalcin   if (th->adapt) {
178818efac9SLisandro Dalcin     ierr = VecZeroEntries(th->vec_lte_work[0]);CHKERRQ(ierr);
179818efac9SLisandro Dalcin     ierr = VecAXPY(th->vec_lte_work[0],+2,X2);CHKERRQ(ierr);
180818efac9SLisandro Dalcin     ierr = VecAXPY(th->vec_lte_work[0],-4,X1);CHKERRQ(ierr);
181818efac9SLisandro Dalcin     ierr = VecAXPY(th->vec_lte_work[0],+2,X0);CHKERRQ(ierr);
182818efac9SLisandro Dalcin   }
183818efac9SLisandro Dalcin   if (th->adapt) {
184818efac9SLisandro Dalcin     ierr = VecZeroEntries(th->vec_lte_work[1]);CHKERRQ(ierr);
185818efac9SLisandro Dalcin     ierr = VecAXPY(th->vec_lte_work[1],+2,V2);CHKERRQ(ierr);
186818efac9SLisandro Dalcin     ierr = VecAXPY(th->vec_lte_work[1],-4,V1);CHKERRQ(ierr);
187818efac9SLisandro Dalcin     ierr = VecAXPY(th->vec_lte_work[1],+2,V0);CHKERRQ(ierr);
188818efac9SLisandro Dalcin   }
189818efac9SLisandro Dalcin 
190818efac9SLisandro Dalcin  finally:
191818efac9SLisandro Dalcin   /* Revert TSAlpha to the initial state (t0,X0,V0) */
192818efac9SLisandro Dalcin   if (initok) *initok = stageok;
193818efac9SLisandro Dalcin   ierr = TSSetTimeStep(ts,time_step);CHKERRQ(ierr);
194818efac9SLisandro Dalcin   ierr = TSAlpha2SetParams(ts,alpha_m,alpha_f,gamma,beta);CHKERRQ(ierr);
195818efac9SLisandro Dalcin   ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
196818efac9SLisandro Dalcin   ierr = VecCopy(ts->vec_dot,th->V0);CHKERRQ(ierr);
197818efac9SLisandro Dalcin 
198818efac9SLisandro Dalcin   ierr = VecDestroy(&X1);CHKERRQ(ierr);
199818efac9SLisandro Dalcin   ierr = VecDestroy(&V1);CHKERRQ(ierr);
200818efac9SLisandro Dalcin   PetscFunctionReturn(0);
201818efac9SLisandro Dalcin }
202818efac9SLisandro Dalcin 
203818efac9SLisandro Dalcin #undef __FUNCT__
204818efac9SLisandro Dalcin #define __FUNCT__ "TSStep_Alpha"
205818efac9SLisandro Dalcin static PetscErrorCode TSStep_Alpha(TS ts)
206818efac9SLisandro Dalcin {
207818efac9SLisandro Dalcin   TS_Alpha       *th = (TS_Alpha*)ts->data;
208818efac9SLisandro Dalcin   PetscInt       rejections = 0;
209818efac9SLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
210818efac9SLisandro Dalcin   PetscReal      next_time_step = ts->time_step;
211818efac9SLisandro Dalcin   PetscErrorCode ierr;
212818efac9SLisandro Dalcin 
213818efac9SLisandro Dalcin   PetscFunctionBegin;
214818efac9SLisandro Dalcin   ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr);
215818efac9SLisandro Dalcin 
216818efac9SLisandro Dalcin   if (!ts->steprollback) {
217818efac9SLisandro Dalcin     if (th->adapt) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); }
218818efac9SLisandro Dalcin     if (th->adapt) { ierr = VecCopy(th->V0,th->vec_dot_prev);CHKERRQ(ierr); }
219818efac9SLisandro Dalcin     ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr);
220818efac9SLisandro Dalcin     ierr = VecCopy(ts->vec_dot,th->V0);CHKERRQ(ierr);
221818efac9SLisandro Dalcin     ierr = VecCopy(th->A1,th->A0);CHKERRQ(ierr);
222818efac9SLisandro Dalcin   }
223818efac9SLisandro Dalcin 
224818efac9SLisandro Dalcin   th->status = TS_STEP_INCOMPLETE;
225818efac9SLisandro Dalcin   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
226818efac9SLisandro Dalcin 
227818efac9SLisandro Dalcin     if (ts->steprestart) {
228818efac9SLisandro Dalcin       ierr = TSAlpha_Restart(ts,&stageok);CHKERRQ(ierr);
229818efac9SLisandro Dalcin       if (!stageok) goto reject_step;
230818efac9SLisandro Dalcin     }
231818efac9SLisandro Dalcin 
232818efac9SLisandro Dalcin     ierr = TSAlpha_StageTime(ts);CHKERRQ(ierr);
233818efac9SLisandro Dalcin     ierr = VecCopy(th->X0,th->X1);CHKERRQ(ierr);
234818efac9SLisandro Dalcin     ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr);
235818efac9SLisandro Dalcin     ierr = TS_SNESSolve(ts,NULL,th->X1);CHKERRQ(ierr);
236818efac9SLisandro Dalcin     ierr = TSPostStage(ts,th->stage_time,0,&th->Xa);CHKERRQ(ierr);
237818efac9SLisandro Dalcin     ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->Xa,&stageok);CHKERRQ(ierr);
238818efac9SLisandro Dalcin     if (!stageok) goto reject_step;
239818efac9SLisandro Dalcin 
240818efac9SLisandro Dalcin     th->status = TS_STEP_PENDING;
241818efac9SLisandro Dalcin     ierr = VecCopy(th->X1,ts->vec_sol);CHKERRQ(ierr);
242818efac9SLisandro Dalcin     ierr = VecCopy(th->V1,ts->vec_dot);CHKERRQ(ierr);
243818efac9SLisandro Dalcin     ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
244818efac9SLisandro Dalcin     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
245818efac9SLisandro Dalcin     if (!accept) {
246818efac9SLisandro Dalcin       ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
247818efac9SLisandro Dalcin       ierr = VecCopy(th->V0,ts->vec_dot);CHKERRQ(ierr);
248818efac9SLisandro Dalcin       ts->time_step = next_time_step;
249818efac9SLisandro Dalcin       goto reject_step;
250818efac9SLisandro Dalcin     }
251818efac9SLisandro Dalcin 
252818efac9SLisandro Dalcin     ts->ptime += ts->time_step;
253818efac9SLisandro Dalcin     ts->time_step = next_time_step;
254818efac9SLisandro Dalcin     break;
255818efac9SLisandro Dalcin 
256818efac9SLisandro Dalcin   reject_step:
257818efac9SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
258818efac9SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
259818efac9SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
260818efac9SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
261818efac9SLisandro Dalcin     }
262818efac9SLisandro Dalcin 
263818efac9SLisandro Dalcin   }
264818efac9SLisandro Dalcin   PetscFunctionReturn(0);
265818efac9SLisandro Dalcin }
266818efac9SLisandro Dalcin 
267818efac9SLisandro Dalcin #undef __FUNCT__
268818efac9SLisandro Dalcin #define __FUNCT__ "TSEvaluateWLTE_Alpha"
269818efac9SLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_Alpha(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
270818efac9SLisandro Dalcin {
271818efac9SLisandro Dalcin   TS_Alpha       *th = (TS_Alpha*)ts->data;
272818efac9SLisandro Dalcin   Vec            X = th->X1;              /* X = solution */
273818efac9SLisandro Dalcin   Vec            V = th->V1;              /* V = solution */
274818efac9SLisandro Dalcin   Vec            Y = th->vec_lte_work[0]; /* Y = X + LTE  */
275818efac9SLisandro Dalcin   Vec            Z = th->vec_lte_work[1]; /* Z = V + LTE  */
276*7453f775SEmil Constantinescu   PetscReal      enormX,enormV,enormXa,enormVa,enormXr,enormVr;
277818efac9SLisandro Dalcin   PetscErrorCode ierr;
278818efac9SLisandro Dalcin 
279818efac9SLisandro Dalcin   PetscFunctionBegin;
280818efac9SLisandro Dalcin   if (ts->steprestart) {
281818efac9SLisandro Dalcin     /* th->vec_{sol|dot}_prev is set to the LTE in TSAlpha_Restart() */
282818efac9SLisandro Dalcin     ierr = VecAXPY(Y,1,X);CHKERRQ(ierr);
283818efac9SLisandro Dalcin     ierr = VecAXPY(Z,1,V);CHKERRQ(ierr);
284818efac9SLisandro Dalcin   } else {
285818efac9SLisandro Dalcin     /* Compute LTE using backward differences with non-constant time step */
286818efac9SLisandro Dalcin     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
287818efac9SLisandro Dalcin     PetscReal   a = 1 + h_prev/h;
288818efac9SLisandro Dalcin     PetscScalar scal[3]; Vec vecX[3],vecV[3];
289818efac9SLisandro Dalcin     scal[0] = +1/a;   scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1));
290818efac9SLisandro Dalcin     vecX[0] = th->X1; vecX[1] = th->X0;   vecX[2] = th->vec_sol_prev;
291818efac9SLisandro Dalcin     vecV[0] = th->V1; vecV[1] = th->V0;   vecV[2] = th->vec_dot_prev;
292818efac9SLisandro Dalcin     ierr = VecCopy(X,Y);CHKERRQ(ierr);
293818efac9SLisandro Dalcin     ierr = VecMAXPY(Y,3,scal,vecX);CHKERRQ(ierr);
294818efac9SLisandro Dalcin     ierr = VecCopy(V,Z);CHKERRQ(ierr);
295818efac9SLisandro Dalcin     ierr = VecMAXPY(Z,3,scal,vecV);CHKERRQ(ierr);
296818efac9SLisandro Dalcin   }
297818efac9SLisandro Dalcin   /* XXX ts->atol and ts->vatol are not appropriate for computing enormV */
298*7453f775SEmil Constantinescu   ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,&enormX,&enormXa,&enormXr);CHKERRQ(ierr);
299*7453f775SEmil Constantinescu   ierr = TSErrorWeightedNorm(ts,V,Z,wnormtype,&enormV,&enormVa,&enormVr);CHKERRQ(ierr);
300818efac9SLisandro Dalcin   if (wnormtype == NORM_2)
301818efac9SLisandro Dalcin     *wlte = PetscSqrtReal(PetscSqr(enormX)/2 + PetscSqr(enormV)/2);
302818efac9SLisandro Dalcin   else
303818efac9SLisandro Dalcin     *wlte = PetscMax(enormX,enormV);
304818efac9SLisandro Dalcin   if (order) *order = 2;
305818efac9SLisandro Dalcin   PetscFunctionReturn(0);
306818efac9SLisandro Dalcin }
307818efac9SLisandro Dalcin 
308818efac9SLisandro Dalcin #undef __FUNCT__
309818efac9SLisandro Dalcin #define __FUNCT__ "TSRollBack_Alpha"
310818efac9SLisandro Dalcin static PetscErrorCode TSRollBack_Alpha(TS ts)
311818efac9SLisandro Dalcin {
312818efac9SLisandro Dalcin   TS_Alpha       *th = (TS_Alpha*)ts->data;
313818efac9SLisandro Dalcin   PetscErrorCode ierr;
314818efac9SLisandro Dalcin 
315818efac9SLisandro Dalcin   PetscFunctionBegin;
316818efac9SLisandro Dalcin   ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr);
317818efac9SLisandro Dalcin   ierr = VecCopy(th->V0,ts->vec_dot);CHKERRQ(ierr);
318818efac9SLisandro Dalcin   PetscFunctionReturn(0);
319818efac9SLisandro Dalcin }
320818efac9SLisandro Dalcin 
321818efac9SLisandro Dalcin /*
322818efac9SLisandro Dalcin #undef __FUNCT__
323818efac9SLisandro Dalcin #define __FUNCT__ "TSInterpolate_Alpha"
324818efac9SLisandro Dalcin static PetscErrorCode TSInterpolate_Alpha(TS ts,PetscReal t,Vec X,Vec V)
325818efac9SLisandro Dalcin {
326818efac9SLisandro Dalcin   TS_Alpha       *th = (TS_Alpha*)ts->data;
327818efac9SLisandro Dalcin   PetscReal      dt  = t - ts->ptime;
328818efac9SLisandro Dalcin   PetscErrorCode ierr;
329818efac9SLisandro Dalcin 
330818efac9SLisandro Dalcin   PetscFunctionBegin;
331818efac9SLisandro Dalcin   ierr = VecCopy(ts->vec_dot,V);CHKERRQ(ierr);
332818efac9SLisandro Dalcin   ierr = VecAXPY(V,dt*(1-th->Gamma),th->A0);CHKERRQ(ierr);
333818efac9SLisandro Dalcin   ierr = VecAXPY(V,dt*th->Gamma,th->A1);CHKERRQ(ierr);
334818efac9SLisandro Dalcin   ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr);
335818efac9SLisandro Dalcin   ierr = VecAXPY(X,dt,V);CHKERRQ(ierr);
336818efac9SLisandro Dalcin   ierr = VecAXPY(X,dt*dt*((PetscReal)0.5-th->Beta),th->A0);CHKERRQ(ierr);
337818efac9SLisandro Dalcin   ierr = VecAXPY(X,dt*dt*th->Beta,th->A1);CHKERRQ(ierr);
338818efac9SLisandro Dalcin   PetscFunctionReturn(0);
339818efac9SLisandro Dalcin }
340818efac9SLisandro Dalcin */
341818efac9SLisandro Dalcin 
342818efac9SLisandro Dalcin #undef __FUNCT__
343818efac9SLisandro Dalcin #define __FUNCT__ "SNESTSFormFunction_Alpha"
344818efac9SLisandro Dalcin static PetscErrorCode SNESTSFormFunction_Alpha(PETSC_UNUSED SNES snes,Vec X,Vec F,TS ts)
345818efac9SLisandro Dalcin {
346818efac9SLisandro Dalcin   TS_Alpha       *th = (TS_Alpha*)ts->data;
347818efac9SLisandro Dalcin   PetscReal      ta = th->stage_time;
348818efac9SLisandro Dalcin   Vec            Xa = th->Xa, Va = th->Va, Aa = th->Aa;
349818efac9SLisandro Dalcin   PetscErrorCode ierr;
350818efac9SLisandro Dalcin 
351818efac9SLisandro Dalcin   PetscFunctionBegin;
352818efac9SLisandro Dalcin   ierr = TSAlpha_StageVecs(ts,X);CHKERRQ(ierr);
353818efac9SLisandro Dalcin   /* F = Function(ta,Xa,Va,Aa) */
354818efac9SLisandro Dalcin   ierr = TSComputeI2Function(ts,ta,Xa,Va,Aa,F);CHKERRQ(ierr);
355818efac9SLisandro Dalcin   ierr = VecScale(F,th->scale_F);CHKERRQ(ierr);
356818efac9SLisandro Dalcin   PetscFunctionReturn(0);
357818efac9SLisandro Dalcin }
358818efac9SLisandro Dalcin 
359818efac9SLisandro Dalcin #undef __FUNCT__
360818efac9SLisandro Dalcin #define __FUNCT__ "SNESTSFormJacobian_Alpha"
361818efac9SLisandro Dalcin static PetscErrorCode SNESTSFormJacobian_Alpha(PETSC_UNUSED SNES snes,PETSC_UNUSED Vec X,Mat J,Mat P,TS ts)
362818efac9SLisandro Dalcin {
363818efac9SLisandro Dalcin   TS_Alpha       *th = (TS_Alpha*)ts->data;
364818efac9SLisandro Dalcin   PetscReal      ta = th->stage_time;
365818efac9SLisandro Dalcin   Vec            Xa = th->Xa, Va = th->Va, Aa = th->Aa;
366818efac9SLisandro Dalcin   PetscReal      dVdX = th->shift_V, dAdX = th->shift_A;
367818efac9SLisandro Dalcin   PetscErrorCode ierr;
368818efac9SLisandro Dalcin 
369818efac9SLisandro Dalcin   PetscFunctionBegin;
370818efac9SLisandro Dalcin   /* J,P = Jacobian(ta,Xa,Va,Aa) */
371818efac9SLisandro Dalcin   ierr = TSComputeI2Jacobian(ts,ta,Xa,Va,Aa,dVdX,dAdX,J,P);CHKERRQ(ierr);
372818efac9SLisandro Dalcin   PetscFunctionReturn(0);
373818efac9SLisandro Dalcin }
374818efac9SLisandro Dalcin 
375818efac9SLisandro Dalcin #undef __FUNCT__
376818efac9SLisandro Dalcin #define __FUNCT__ "TSReset_Alpha"
377818efac9SLisandro Dalcin static PetscErrorCode TSReset_Alpha(TS ts)
378818efac9SLisandro Dalcin {
379818efac9SLisandro Dalcin   TS_Alpha       *th = (TS_Alpha*)ts->data;
380818efac9SLisandro Dalcin   PetscErrorCode ierr;
381818efac9SLisandro Dalcin 
382818efac9SLisandro Dalcin   PetscFunctionBegin;
383818efac9SLisandro Dalcin   ierr = VecDestroy(&th->X0);CHKERRQ(ierr);
384818efac9SLisandro Dalcin   ierr = VecDestroy(&th->Xa);CHKERRQ(ierr);
385818efac9SLisandro Dalcin   ierr = VecDestroy(&th->X1);CHKERRQ(ierr);
386818efac9SLisandro Dalcin   ierr = VecDestroy(&th->V0);CHKERRQ(ierr);
387818efac9SLisandro Dalcin   ierr = VecDestroy(&th->Va);CHKERRQ(ierr);
388818efac9SLisandro Dalcin   ierr = VecDestroy(&th->V1);CHKERRQ(ierr);
389818efac9SLisandro Dalcin   ierr = VecDestroy(&th->A0);CHKERRQ(ierr);
390818efac9SLisandro Dalcin   ierr = VecDestroy(&th->Aa);CHKERRQ(ierr);
391818efac9SLisandro Dalcin   ierr = VecDestroy(&th->A1);CHKERRQ(ierr);
392818efac9SLisandro Dalcin   ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr);
393818efac9SLisandro Dalcin   ierr = VecDestroy(&th->vec_dot_prev);CHKERRQ(ierr);
394818efac9SLisandro Dalcin   ierr = VecDestroy(&th->vec_lte_work[0]);CHKERRQ(ierr);
395818efac9SLisandro Dalcin   ierr = VecDestroy(&th->vec_lte_work[1]);CHKERRQ(ierr);
396818efac9SLisandro Dalcin   PetscFunctionReturn(0);
397818efac9SLisandro Dalcin }
398818efac9SLisandro Dalcin 
399818efac9SLisandro Dalcin #undef __FUNCT__
400818efac9SLisandro Dalcin #define __FUNCT__ "TSDestroy_Alpha"
401818efac9SLisandro Dalcin static PetscErrorCode TSDestroy_Alpha(TS ts)
402818efac9SLisandro Dalcin {
403818efac9SLisandro Dalcin   PetscErrorCode ierr;
404818efac9SLisandro Dalcin 
405818efac9SLisandro Dalcin   PetscFunctionBegin;
406818efac9SLisandro Dalcin   ierr = TSReset_Alpha(ts);CHKERRQ(ierr);
407818efac9SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
408818efac9SLisandro Dalcin 
409818efac9SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2UseAdapt_C",NULL);CHKERRQ(ierr);
410818efac9SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2SetRadius_C",NULL);CHKERRQ(ierr);
411818efac9SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2SetParams_C",NULL);CHKERRQ(ierr);
412818efac9SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2GetParams_C",NULL);CHKERRQ(ierr);
413818efac9SLisandro Dalcin   PetscFunctionReturn(0);
414818efac9SLisandro Dalcin }
415818efac9SLisandro Dalcin 
416818efac9SLisandro Dalcin #undef __FUNCT__
417818efac9SLisandro Dalcin #define __FUNCT__ "TSSetUp_Alpha"
418818efac9SLisandro Dalcin static PetscErrorCode TSSetUp_Alpha(TS ts)
419818efac9SLisandro Dalcin {
420818efac9SLisandro Dalcin   TS_Alpha       *th = (TS_Alpha*)ts->data;
421818efac9SLisandro Dalcin   PetscErrorCode ierr;
422818efac9SLisandro Dalcin 
423818efac9SLisandro Dalcin   PetscFunctionBegin;
424818efac9SLisandro Dalcin   ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr);
425818efac9SLisandro Dalcin   ierr = VecDuplicate(ts->vec_sol,&th->Xa);CHKERRQ(ierr);
426818efac9SLisandro Dalcin   ierr = VecDuplicate(ts->vec_sol,&th->X1);CHKERRQ(ierr);
427818efac9SLisandro Dalcin   ierr = VecDuplicate(ts->vec_sol,&th->V0);CHKERRQ(ierr);
428818efac9SLisandro Dalcin   ierr = VecDuplicate(ts->vec_sol,&th->Va);CHKERRQ(ierr);
429818efac9SLisandro Dalcin   ierr = VecDuplicate(ts->vec_sol,&th->V1);CHKERRQ(ierr);
430818efac9SLisandro Dalcin   ierr = VecDuplicate(ts->vec_sol,&th->A0);CHKERRQ(ierr);
431818efac9SLisandro Dalcin   ierr = VecDuplicate(ts->vec_sol,&th->Aa);CHKERRQ(ierr);
432818efac9SLisandro Dalcin   ierr = VecDuplicate(ts->vec_sol,&th->A1);CHKERRQ(ierr);
433818efac9SLisandro Dalcin 
434818efac9SLisandro Dalcin   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
435818efac9SLisandro Dalcin   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
436818efac9SLisandro Dalcin   if (!th->adapt) {
437818efac9SLisandro Dalcin     ierr = TSAdaptSetType(ts->adapt,TSADAPTNONE);CHKERRQ(ierr);
438818efac9SLisandro Dalcin   } else {
439818efac9SLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr);
440818efac9SLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_dot_prev);CHKERRQ(ierr);
441818efac9SLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work[0]);CHKERRQ(ierr);
442818efac9SLisandro Dalcin     ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work[1]);CHKERRQ(ierr);
443818efac9SLisandro Dalcin     if (ts->exact_final_time == TS_EXACTFINALTIME_UNSPECIFIED)
444818efac9SLisandro Dalcin       ts->exact_final_time = TS_EXACTFINALTIME_MATCHSTEP;
445818efac9SLisandro Dalcin   }
446818efac9SLisandro Dalcin 
447818efac9SLisandro Dalcin   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
448818efac9SLisandro Dalcin   PetscFunctionReturn(0);
449818efac9SLisandro Dalcin }
450818efac9SLisandro Dalcin 
451818efac9SLisandro Dalcin #undef __FUNCT__
452818efac9SLisandro Dalcin #define __FUNCT__ "TSSetFromOptions_Alpha"
453818efac9SLisandro Dalcin static PetscErrorCode TSSetFromOptions_Alpha(PetscOptionItems *PetscOptionsObject,TS ts)
454818efac9SLisandro Dalcin {
455818efac9SLisandro Dalcin   TS_Alpha       *th = (TS_Alpha*)ts->data;
456818efac9SLisandro Dalcin   PetscErrorCode ierr;
457818efac9SLisandro Dalcin 
458818efac9SLisandro Dalcin   PetscFunctionBegin;
459818efac9SLisandro Dalcin   ierr = PetscOptionsHead(PetscOptionsObject,"Generalized-Alpha ODE solver options");CHKERRQ(ierr);
460818efac9SLisandro Dalcin   {
461818efac9SLisandro Dalcin     PetscBool flg;
462818efac9SLisandro Dalcin     PetscReal radius = 1;
463818efac9SLisandro Dalcin     PetscBool adapt  = th->adapt;
464818efac9SLisandro Dalcin     ierr = PetscOptionsReal("-ts_alpha_radius","Spectral radius (high-frequency dissipation)","TSAlpha2SetRadius",radius,&radius,&flg);CHKERRQ(ierr);
465818efac9SLisandro Dalcin     if (flg) {ierr = TSAlpha2SetRadius(ts,radius);CHKERRQ(ierr);}
466818efac9SLisandro Dalcin     ierr = PetscOptionsReal("-ts_alpha_alpha_m","Algoritmic parameter alpha_m","TSAlpha2SetParams",th->Alpha_m,&th->Alpha_m,NULL);CHKERRQ(ierr);
467818efac9SLisandro Dalcin     ierr = PetscOptionsReal("-ts_alpha_alpha_f","Algoritmic parameter alpha_f","TSAlpha2SetParams",th->Alpha_f,&th->Alpha_f,NULL);CHKERRQ(ierr);
468818efac9SLisandro Dalcin     ierr = PetscOptionsReal("-ts_alpha_gamma","Algoritmic parameter gamma","TSAlpha2SetParams",th->Gamma,&th->Gamma,NULL);CHKERRQ(ierr);
469818efac9SLisandro Dalcin     ierr = PetscOptionsReal("-ts_alpha_beta","Algoritmic parameter beta","TSAlpha2SetParams",th->Beta,&th->Beta,NULL);CHKERRQ(ierr);
470818efac9SLisandro Dalcin     ierr = TSAlpha2SetParams(ts,th->Alpha_m,th->Alpha_f,th->Gamma,th->Beta);CHKERRQ(ierr);
471818efac9SLisandro Dalcin     ierr = PetscOptionsBool("-ts_alpha_adapt","Use time-step adaptivity with the Alpha method","TSAlpha2UseAdapt",adapt,&adapt,&flg);CHKERRQ(ierr);
472818efac9SLisandro Dalcin     if (flg) {ierr = TSAlpha2UseAdapt(ts,adapt);CHKERRQ(ierr);}
473818efac9SLisandro Dalcin   }
474818efac9SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
475818efac9SLisandro Dalcin   PetscFunctionReturn(0);
476818efac9SLisandro Dalcin }
477818efac9SLisandro Dalcin 
478818efac9SLisandro Dalcin #undef __FUNCT__
479818efac9SLisandro Dalcin #define __FUNCT__ "TSView_Alpha"
480818efac9SLisandro Dalcin static PetscErrorCode TSView_Alpha(TS ts,PetscViewer viewer)
481818efac9SLisandro Dalcin {
482818efac9SLisandro Dalcin   TS_Alpha       *th = (TS_Alpha*)ts->data;
483818efac9SLisandro Dalcin   PetscBool      iascii;
484818efac9SLisandro Dalcin   PetscErrorCode ierr;
485818efac9SLisandro Dalcin 
486818efac9SLisandro Dalcin   PetscFunctionBegin;
487818efac9SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
488818efac9SLisandro Dalcin   if (iascii) {
489818efac9SLisandro 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);
490818efac9SLisandro Dalcin   }
491818efac9SLisandro Dalcin   if (ts->adapt) {ierr = TSAdaptView(ts->adapt,viewer);CHKERRQ(ierr);}
492818efac9SLisandro Dalcin   if (ts->snes)  {ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);}
493818efac9SLisandro Dalcin   PetscFunctionReturn(0);
494818efac9SLisandro Dalcin }
495818efac9SLisandro Dalcin 
496818efac9SLisandro Dalcin #undef __FUNCT__
497818efac9SLisandro Dalcin #define __FUNCT__ "TSAlpha2UseAdapt_Alpha"
498818efac9SLisandro Dalcin static PetscErrorCode TSAlpha2UseAdapt_Alpha(TS ts,PetscBool use)
499818efac9SLisandro Dalcin {
500818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha*)ts->data;
501818efac9SLisandro Dalcin 
502818efac9SLisandro Dalcin   PetscFunctionBegin;
503818efac9SLisandro Dalcin   if (use == th->adapt) PetscFunctionReturn(0);
504818efac9SLisandro Dalcin   if (ts->setupcalled) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ORDER,"Cannot change adaptivity after TSSetUp()");
505818efac9SLisandro Dalcin   th->adapt = use;
506818efac9SLisandro Dalcin   PetscFunctionReturn(0);
507818efac9SLisandro Dalcin }
508818efac9SLisandro Dalcin 
509818efac9SLisandro Dalcin #undef __FUNCT__
510818efac9SLisandro Dalcin #define __FUNCT__ "TSAlpha2SetRadius_Alpha"
511818efac9SLisandro Dalcin static PetscErrorCode TSAlpha2SetRadius_Alpha(TS ts,PetscReal radius)
512818efac9SLisandro Dalcin {
513818efac9SLisandro Dalcin   PetscReal      alpha_m,alpha_f,gamma,beta;
514818efac9SLisandro Dalcin   PetscErrorCode ierr;
515818efac9SLisandro Dalcin 
516818efac9SLisandro Dalcin   PetscFunctionBegin;
517818efac9SLisandro Dalcin   if (radius < 0 || radius > 1) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Radius %g not in range [0,1]",(double)radius);
518818efac9SLisandro Dalcin   alpha_m = (2-radius)/(1+radius);
519818efac9SLisandro Dalcin   alpha_f = 1/(1+radius);
520818efac9SLisandro Dalcin   gamma   = (PetscReal)0.5 + alpha_m - alpha_f;
521818efac9SLisandro Dalcin   beta    = (PetscReal)0.5 * (1 + alpha_m - alpha_f); beta *= beta;
522818efac9SLisandro Dalcin   ierr = TSAlpha2SetParams(ts,alpha_m,alpha_f,gamma,beta);CHKERRQ(ierr);
523818efac9SLisandro Dalcin   PetscFunctionReturn(0);
524818efac9SLisandro Dalcin }
525818efac9SLisandro Dalcin 
526818efac9SLisandro Dalcin #undef __FUNCT__
527818efac9SLisandro Dalcin #define __FUNCT__ "TSAlpha2SetParams_Alpha"
528818efac9SLisandro Dalcin static PetscErrorCode TSAlpha2SetParams_Alpha(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma,PetscReal beta)
529818efac9SLisandro Dalcin {
530818efac9SLisandro Dalcin   TS_Alpha  *th = (TS_Alpha*)ts->data;
531818efac9SLisandro Dalcin   PetscReal tol = 100*PETSC_MACHINE_EPSILON;
532818efac9SLisandro Dalcin   PetscReal res = ((PetscReal)0.5 + alpha_m - alpha_f) - gamma;
533818efac9SLisandro Dalcin 
534818efac9SLisandro Dalcin   PetscFunctionBegin;
535818efac9SLisandro Dalcin   th->Alpha_m = alpha_m;
536818efac9SLisandro Dalcin   th->Alpha_f = alpha_f;
537818efac9SLisandro Dalcin   th->Gamma   = gamma;
538818efac9SLisandro Dalcin   th->Beta    = beta;
539818efac9SLisandro Dalcin   th->order   = (PetscAbsReal(res) < tol) ? 2 : 1;
540818efac9SLisandro Dalcin   PetscFunctionReturn(0);
541818efac9SLisandro Dalcin }
542818efac9SLisandro Dalcin 
543818efac9SLisandro Dalcin #undef __FUNCT__
544818efac9SLisandro Dalcin #define __FUNCT__ "TSAlpha2GetParams_Alpha"
545818efac9SLisandro Dalcin static PetscErrorCode TSAlpha2GetParams_Alpha(TS ts,PetscReal *alpha_m,PetscReal *alpha_f,PetscReal *gamma,PetscReal *beta)
546818efac9SLisandro Dalcin {
547818efac9SLisandro Dalcin   TS_Alpha *th = (TS_Alpha*)ts->data;
548818efac9SLisandro Dalcin 
549818efac9SLisandro Dalcin   PetscFunctionBegin;
550818efac9SLisandro Dalcin   if (alpha_m) *alpha_m = th->Alpha_m;
551818efac9SLisandro Dalcin   if (alpha_f) *alpha_f = th->Alpha_f;
552818efac9SLisandro Dalcin   if (gamma)   *gamma   = th->Gamma;
553818efac9SLisandro Dalcin   if (beta)    *beta    = th->Beta;
554818efac9SLisandro Dalcin   PetscFunctionReturn(0);
555818efac9SLisandro Dalcin }
556818efac9SLisandro Dalcin 
557818efac9SLisandro Dalcin /*MC
558818efac9SLisandro Dalcin       TSALPHA2 - ODE/DAE solver using the implicit Generalized-Alpha method
559818efac9SLisandro Dalcin                  for second-order systems
560818efac9SLisandro Dalcin 
561818efac9SLisandro Dalcin   Level: beginner
562818efac9SLisandro Dalcin 
563818efac9SLisandro Dalcin   References:
564818efac9SLisandro Dalcin   J. Chung, G.M.Hubert. "A Time Integration Algorithm for Structural
565818efac9SLisandro Dalcin   Dynamics with Improved Numerical Dissipation: The Generalized-alpha
566818efac9SLisandro Dalcin   Method" ASME Journal of Applied Mechanics, 60, 371:375, 1993.
567818efac9SLisandro Dalcin 
568818efac9SLisandro Dalcin .seealso:  TS, TSCreate(), TSSetType(), TSAlpha2SetRadius(), TSAlpha2SetParams()
569818efac9SLisandro Dalcin M*/
570818efac9SLisandro Dalcin #undef __FUNCT__
571818efac9SLisandro Dalcin #define __FUNCT__ "TSCreate_Alpha2"
572818efac9SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSCreate_Alpha2(TS ts)
573818efac9SLisandro Dalcin {
574818efac9SLisandro Dalcin   TS_Alpha       *th;
575818efac9SLisandro Dalcin   PetscErrorCode ierr;
576818efac9SLisandro Dalcin 
577818efac9SLisandro Dalcin   PetscFunctionBegin;
578818efac9SLisandro Dalcin   ts->ops->reset          = TSReset_Alpha;
579818efac9SLisandro Dalcin   ts->ops->destroy        = TSDestroy_Alpha;
580818efac9SLisandro Dalcin   ts->ops->view           = TSView_Alpha;
581818efac9SLisandro Dalcin   ts->ops->setup          = TSSetUp_Alpha;
582818efac9SLisandro Dalcin   ts->ops->setfromoptions = TSSetFromOptions_Alpha;
583818efac9SLisandro Dalcin   ts->ops->step           = TSStep_Alpha;
584818efac9SLisandro Dalcin   ts->ops->evaluatewlte   = TSEvaluateWLTE_Alpha;
585818efac9SLisandro Dalcin   ts->ops->rollback       = TSRollBack_Alpha;
586818efac9SLisandro Dalcin   /*ts->ops->interpolate  = TSInterpolate_Alpha;*/
587818efac9SLisandro Dalcin   ts->ops->snesfunction   = SNESTSFormFunction_Alpha;
588818efac9SLisandro Dalcin   ts->ops->snesjacobian   = SNESTSFormJacobian_Alpha;
589818efac9SLisandro Dalcin 
590818efac9SLisandro Dalcin   ierr = PetscNewLog(ts,&th);CHKERRQ(ierr);
591818efac9SLisandro Dalcin   ts->data = (void*)th;
592818efac9SLisandro Dalcin 
593818efac9SLisandro Dalcin   th->Alpha_m = 0.5;
594818efac9SLisandro Dalcin   th->Alpha_f = 0.5;
595818efac9SLisandro Dalcin   th->Gamma   = 0.5;
596818efac9SLisandro Dalcin   th->Beta    = 0.25;
597818efac9SLisandro Dalcin   th->order   = 2;
598818efac9SLisandro Dalcin 
599818efac9SLisandro Dalcin   th->adapt = PETSC_FALSE;
600818efac9SLisandro Dalcin 
601818efac9SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2UseAdapt_C",TSAlpha2UseAdapt_Alpha);CHKERRQ(ierr);
602818efac9SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2SetRadius_C",TSAlpha2SetRadius_Alpha);CHKERRQ(ierr);
603818efac9SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2SetParams_C",TSAlpha2SetParams_Alpha);CHKERRQ(ierr);
604818efac9SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSAlpha2GetParams_C",TSAlpha2GetParams_Alpha);CHKERRQ(ierr);
605818efac9SLisandro Dalcin   PetscFunctionReturn(0);
606818efac9SLisandro Dalcin }
607818efac9SLisandro Dalcin 
608818efac9SLisandro Dalcin #undef __FUNCT__
609818efac9SLisandro Dalcin #define __FUNCT__ "TSAlpha2UseAdapt"
610818efac9SLisandro Dalcin /*@
611818efac9SLisandro Dalcin   TSAlpha2UseAdapt - Use time-step adaptivity with the Alpha method
612818efac9SLisandro Dalcin 
613818efac9SLisandro Dalcin   Logically Collective on TS
614818efac9SLisandro Dalcin 
615818efac9SLisandro Dalcin   Input Parameter:
616818efac9SLisandro Dalcin +  ts - timestepping context
617818efac9SLisandro Dalcin -  use - flag to use adaptivity
618818efac9SLisandro Dalcin 
619818efac9SLisandro Dalcin   Options Database:
620818efac9SLisandro Dalcin .  -ts_alpha_adapt
621818efac9SLisandro Dalcin 
622818efac9SLisandro Dalcin   Level: intermediate
623818efac9SLisandro Dalcin 
624818efac9SLisandro Dalcin .seealso: TSAdapt, TSADAPTBASIC
625818efac9SLisandro Dalcin @*/
626818efac9SLisandro Dalcin PetscErrorCode TSAlpha2UseAdapt(TS ts,PetscBool use)
627818efac9SLisandro Dalcin {
628818efac9SLisandro Dalcin   PetscErrorCode ierr;
629818efac9SLisandro Dalcin 
630818efac9SLisandro Dalcin   PetscFunctionBegin;
631818efac9SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
632818efac9SLisandro Dalcin   PetscValidLogicalCollectiveBool(ts,use,2);
633818efac9SLisandro Dalcin   ierr = PetscTryMethod(ts,"TSAlpha2UseAdapt_C",(TS,PetscBool),(ts,use));CHKERRQ(ierr);
634818efac9SLisandro Dalcin   PetscFunctionReturn(0);
635818efac9SLisandro Dalcin }
636818efac9SLisandro Dalcin 
637818efac9SLisandro Dalcin #undef __FUNCT__
638818efac9SLisandro Dalcin #define __FUNCT__ "TSAlpha2SetRadius"
639818efac9SLisandro Dalcin /*@
640818efac9SLisandro Dalcin   TSAlpha2SetRadius - sets the desired spectral radius of the method
641818efac9SLisandro Dalcin                       (i.e. high-frequency numerical damping)
642818efac9SLisandro Dalcin 
643818efac9SLisandro Dalcin   Logically Collective on TS
644818efac9SLisandro Dalcin 
645818efac9SLisandro Dalcin   The algorithmic parameters \alpha_m and \alpha_f of the
646818efac9SLisandro Dalcin   generalized-\alpha method can be computed in terms of a specified
647818efac9SLisandro Dalcin   spectral radius \rho in [0,1] for infinite time step in order to
648818efac9SLisandro Dalcin   control high-frequency numerical damping:
649818efac9SLisandro Dalcin     \alpha_m = (2-\rho)/(1+\rho)
650818efac9SLisandro Dalcin     \alpha_f = 1/(1+\rho)
651818efac9SLisandro Dalcin 
652818efac9SLisandro Dalcin   Input Parameter:
653818efac9SLisandro Dalcin +  ts - timestepping context
654818efac9SLisandro Dalcin -  radius - the desired spectral radius
655818efac9SLisandro Dalcin 
656818efac9SLisandro Dalcin   Options Database:
657818efac9SLisandro Dalcin .  -ts_alpha_radius <radius>
658818efac9SLisandro Dalcin 
659818efac9SLisandro Dalcin   Level: intermediate
660818efac9SLisandro Dalcin 
661818efac9SLisandro Dalcin .seealso: TSAlpha2SetParams(), TSAlpha2GetParams()
662818efac9SLisandro Dalcin @*/
663818efac9SLisandro Dalcin PetscErrorCode TSAlpha2SetRadius(TS ts,PetscReal radius)
664818efac9SLisandro Dalcin {
665818efac9SLisandro Dalcin   PetscErrorCode ierr;
666818efac9SLisandro Dalcin 
667818efac9SLisandro Dalcin   PetscFunctionBegin;
668818efac9SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
669818efac9SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts,radius,2);
670818efac9SLisandro Dalcin   if (radius < 0 || radius > 1) SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Radius %g not in range [0,1]",(double)radius);
671818efac9SLisandro Dalcin   ierr = PetscTryMethod(ts,"TSAlpha2SetRadius_C",(TS,PetscReal),(ts,radius));CHKERRQ(ierr);
672818efac9SLisandro Dalcin   PetscFunctionReturn(0);
673818efac9SLisandro Dalcin }
674818efac9SLisandro Dalcin 
675818efac9SLisandro Dalcin #undef __FUNCT__
676818efac9SLisandro Dalcin #define __FUNCT__ "TSAlpha2SetParams"
677818efac9SLisandro Dalcin /*@
678818efac9SLisandro Dalcin   TSAlpha2SetParams - sets the algorithmic parameters for TSALPHA2
679818efac9SLisandro Dalcin 
680818efac9SLisandro Dalcin   Logically Collective on TS
681818efac9SLisandro Dalcin 
682818efac9SLisandro Dalcin   Second-order accuracy can be obtained so long as:
683818efac9SLisandro Dalcin     \gamma = 1/2 + alpha_m - alpha_f
684818efac9SLisandro Dalcin     \beta  = 1/4 (1 + alpha_m - alpha_f)^2
685818efac9SLisandro Dalcin 
686818efac9SLisandro Dalcin   Unconditional stability requires:
687818efac9SLisandro Dalcin     \alpha_m >= \alpha_f >= 1/2
688818efac9SLisandro Dalcin 
689818efac9SLisandro Dalcin 
690818efac9SLisandro Dalcin   Input Parameter:
691818efac9SLisandro Dalcin + ts - timestepping context
692818efac9SLisandro Dalcin . \alpha_m - algorithmic paramenter
693818efac9SLisandro Dalcin . \alpha_f - algorithmic paramenter
694818efac9SLisandro Dalcin . \gamma   - algorithmic paramenter
695818efac9SLisandro Dalcin - \beta    - algorithmic paramenter
696818efac9SLisandro Dalcin 
697818efac9SLisandro Dalcin   Options Database:
698818efac9SLisandro Dalcin + -ts_alpha_alpha_m <alpha_m>
699818efac9SLisandro Dalcin . -ts_alpha_alpha_f <alpha_f>
700818efac9SLisandro Dalcin . -ts_alpha_gamma   <gamma>
701818efac9SLisandro Dalcin - -ts_alpha_beta    <beta>
702818efac9SLisandro Dalcin 
703818efac9SLisandro Dalcin   Note:
704818efac9SLisandro Dalcin   Use of this function is normally only required to hack TSALPHA2 to
705818efac9SLisandro Dalcin   use a modified integration scheme. Users should call
706818efac9SLisandro Dalcin   TSAlpha2SetRadius() to set the desired spectral radius of the methods
707818efac9SLisandro Dalcin   (i.e. high-frequency damping) in order so select optimal values for
708818efac9SLisandro Dalcin   these parameters.
709818efac9SLisandro Dalcin 
710818efac9SLisandro Dalcin   Level: advanced
711818efac9SLisandro Dalcin 
712818efac9SLisandro Dalcin .seealso: TSAlpha2SetRadius(), TSAlpha2GetParams()
713818efac9SLisandro Dalcin @*/
714818efac9SLisandro Dalcin PetscErrorCode TSAlpha2SetParams(TS ts,PetscReal alpha_m,PetscReal alpha_f,PetscReal gamma,PetscReal beta)
715818efac9SLisandro Dalcin {
716818efac9SLisandro Dalcin   PetscErrorCode ierr;
717818efac9SLisandro Dalcin 
718818efac9SLisandro Dalcin   PetscFunctionBegin;
719818efac9SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
720818efac9SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts,alpha_m,2);
721818efac9SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts,alpha_f,3);
722818efac9SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts,gamma,4);
723818efac9SLisandro Dalcin   PetscValidLogicalCollectiveReal(ts,beta,5);
724818efac9SLisandro Dalcin   ierr = PetscTryMethod(ts,"TSAlpha2SetParams_C",(TS,PetscReal,PetscReal,PetscReal,PetscReal),(ts,alpha_m,alpha_f,gamma,beta));CHKERRQ(ierr);
725818efac9SLisandro Dalcin   PetscFunctionReturn(0);
726818efac9SLisandro Dalcin }
727818efac9SLisandro Dalcin 
728818efac9SLisandro Dalcin #undef __FUNCT__
729818efac9SLisandro Dalcin #define __FUNCT__ "TSAlpha2GetParams"
730818efac9SLisandro Dalcin /*@
731818efac9SLisandro Dalcin   TSAlpha2GetParams - gets the algorithmic parameters for TSALPHA2
732818efac9SLisandro Dalcin 
733818efac9SLisandro Dalcin   Not Collective
734818efac9SLisandro Dalcin 
735818efac9SLisandro Dalcin   Input Parameter:
736818efac9SLisandro Dalcin . ts - timestepping context
737818efac9SLisandro Dalcin 
738818efac9SLisandro Dalcin   Output Parameters:
739818efac9SLisandro Dalcin + \alpha_m - algorithmic parameter
740818efac9SLisandro Dalcin . \alpha_f - algorithmic parameter
741818efac9SLisandro Dalcin . \gamma   - algorithmic parameter
742818efac9SLisandro Dalcin - \beta    - algorithmic parameter
743818efac9SLisandro Dalcin 
744818efac9SLisandro Dalcin   Note:
745818efac9SLisandro Dalcin   Use of this function is normally only required to hack TSALPHA2 to
746818efac9SLisandro Dalcin   use a modified integration scheme. Users should call
747818efac9SLisandro Dalcin   TSAlpha2SetRadius() to set the high-frequency damping (i.e. spectral
748818efac9SLisandro Dalcin   radius of the method) in order so select optimal values for these
749818efac9SLisandro Dalcin   parameters.
750818efac9SLisandro Dalcin 
751818efac9SLisandro Dalcin   Level: advanced
752818efac9SLisandro Dalcin 
753818efac9SLisandro Dalcin .seealso: TSAlpha2SetRadius(), TSAlpha2SetParams()
754818efac9SLisandro Dalcin @*/
755818efac9SLisandro Dalcin PetscErrorCode TSAlpha2GetParams(TS ts,PetscReal *alpha_m,PetscReal *alpha_f,PetscReal *gamma,PetscReal *beta)
756818efac9SLisandro Dalcin {
757818efac9SLisandro Dalcin   PetscErrorCode ierr;
758818efac9SLisandro Dalcin 
759818efac9SLisandro Dalcin   PetscFunctionBegin;
760818efac9SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
761818efac9SLisandro Dalcin   if (alpha_m) PetscValidRealPointer(alpha_m,2);
762818efac9SLisandro Dalcin   if (alpha_f) PetscValidRealPointer(alpha_f,3);
763818efac9SLisandro Dalcin   if (gamma)   PetscValidRealPointer(gamma,4);
764818efac9SLisandro Dalcin   if (beta)    PetscValidRealPointer(beta,5);
765818efac9SLisandro Dalcin   ierr = PetscUseMethod(ts,"TSAlpha2GetParams_C",(TS,PetscReal*,PetscReal*,PetscReal*,PetscReal*),(ts,alpha_m,alpha_f,gamma,beta));CHKERRQ(ierr);
766818efac9SLisandro Dalcin   PetscFunctionReturn(0);
767818efac9SLisandro Dalcin }
768