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