xref: /petsc/src/ts/impls/bdf/bdf.c (revision efd4aadf157bf1ba2d80c2be092fcf4247860003)
1211a84d6SLisandro Dalcin /*
2211a84d6SLisandro Dalcin   Code for timestepping with BDF methods
3211a84d6SLisandro Dalcin */
4211a84d6SLisandro Dalcin #include <petsc/private/tsimpl.h>  /*I "petscts.h" I*/
5211a84d6SLisandro Dalcin 
6211a84d6SLisandro Dalcin static PetscBool  cited = PETSC_FALSE;
7211a84d6SLisandro Dalcin static const char citation[] =
8211a84d6SLisandro Dalcin   "@book{Brenan1995,\n"
9211a84d6SLisandro Dalcin   "  title     = {Numerical Solution of Initial-Value Problems in Differential-Algebraic Equations},\n"
10211a84d6SLisandro Dalcin   "  author    = {Brenan, K. and Campbell, S. and Petzold, L.},\n"
11211a84d6SLisandro Dalcin   "  publisher = {Society for Industrial and Applied Mathematics},\n"
12211a84d6SLisandro Dalcin   "  year      = {1995},\n"
13211a84d6SLisandro Dalcin   "  doi       = {10.1137/1.9781611971224},\n}\n";
14211a84d6SLisandro Dalcin 
15211a84d6SLisandro Dalcin typedef struct {
16211a84d6SLisandro Dalcin   PetscInt  k,n;
17211a84d6SLisandro Dalcin   PetscReal time[6+2];
18211a84d6SLisandro Dalcin   Vec       work[6+2];
19211a84d6SLisandro Dalcin   PetscReal shift;
20211a84d6SLisandro Dalcin   Vec       vec_dot;
21211a84d6SLisandro Dalcin   Vec       vec_lte;
22211a84d6SLisandro Dalcin 
23211a84d6SLisandro Dalcin   PetscInt     order;
24211a84d6SLisandro Dalcin   TSStepStatus status;
25211a84d6SLisandro Dalcin } TS_BDF;
26211a84d6SLisandro Dalcin 
27211a84d6SLisandro Dalcin 
28211a84d6SLisandro Dalcin PETSC_STATIC_INLINE void LagrangeBasisVals(PetscInt n,PetscReal t,const PetscReal T[],PetscScalar L[])
29211a84d6SLisandro Dalcin {
30211a84d6SLisandro Dalcin   PetscInt k,j;
31211a84d6SLisandro Dalcin   for (k=0; k<n; k++)
32211a84d6SLisandro Dalcin     for (L[k]=1, j=0; j<n; j++)
33211a84d6SLisandro Dalcin       if (j != k)
34211a84d6SLisandro Dalcin         L[k] *= (t - T[j])/(T[k] - T[j]);
35211a84d6SLisandro Dalcin }
36211a84d6SLisandro Dalcin 
37211a84d6SLisandro Dalcin PETSC_STATIC_INLINE void LagrangeBasisDers(PetscInt n,PetscReal t,const PetscReal T[],PetscScalar dL[])
38211a84d6SLisandro Dalcin {
39211a84d6SLisandro Dalcin   PetscInt  k,j,i;
40211a84d6SLisandro Dalcin   for (k=0; k<n; k++)
41211a84d6SLisandro Dalcin     for (dL[k]=0, j=0; j<n; j++)
42211a84d6SLisandro Dalcin       if (j != k) {
43211a84d6SLisandro Dalcin         PetscReal L = 1/(T[k] - T[j]);
44211a84d6SLisandro Dalcin         for (i=0; i<n; i++)
45211a84d6SLisandro Dalcin           if (i != j && i != k)
46211a84d6SLisandro Dalcin             L *= (t - T[i])/(T[k] - T[i]);
47211a84d6SLisandro Dalcin         dL[k] += L;
48211a84d6SLisandro Dalcin       }
49211a84d6SLisandro Dalcin }
50211a84d6SLisandro Dalcin 
51211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_Advance(TS ts,PetscReal t,Vec X)
52211a84d6SLisandro Dalcin {
53211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
54211a84d6SLisandro Dalcin   PetscInt       i,n = (PetscInt)(sizeof(bdf->work)/sizeof(Vec));
55211a84d6SLisandro Dalcin   Vec            tail = bdf->work[n-1];
56211a84d6SLisandro Dalcin   PetscErrorCode ierr;
57211a84d6SLisandro Dalcin 
58211a84d6SLisandro Dalcin   PetscFunctionBegin;
59211a84d6SLisandro Dalcin   for (i=n-1; i>=2; i--) {
60211a84d6SLisandro Dalcin     bdf->time[i] = bdf->time[i-1];
61211a84d6SLisandro Dalcin     bdf->work[i] = bdf->work[i-1];
62211a84d6SLisandro Dalcin   }
63211a84d6SLisandro Dalcin   bdf->n       = PetscMin(bdf->n+1,n-1);
64211a84d6SLisandro Dalcin   bdf->time[1] = t;
65211a84d6SLisandro Dalcin   bdf->work[1] = tail;
66211a84d6SLisandro Dalcin   ierr = VecCopy(X,tail);CHKERRQ(ierr);
67211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
68211a84d6SLisandro Dalcin }
69211a84d6SLisandro Dalcin 
70211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_VecDot(TS ts,PetscInt order,PetscReal t,Vec X,Vec Xdot,PetscReal *shift)
71211a84d6SLisandro Dalcin {
72211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
73211a84d6SLisandro Dalcin   PetscInt       i,n = order+1;
74211a84d6SLisandro Dalcin   PetscReal      time[7];
75211a84d6SLisandro Dalcin   Vec            vecs[7];
76211a84d6SLisandro Dalcin   PetscScalar    alpha[7];
77211a84d6SLisandro Dalcin   PetscErrorCode ierr;
78211a84d6SLisandro Dalcin 
79211a84d6SLisandro Dalcin   PetscFunctionBegin;
80211a84d6SLisandro Dalcin   n = PetscMax(n,2);
81211a84d6SLisandro Dalcin   time[0] = t; for (i=1; i<n; i++) time[i] = bdf->time[i];
82211a84d6SLisandro Dalcin   vecs[0] = X; for (i=1; i<n; i++) vecs[i] = bdf->work[i];
83211a84d6SLisandro Dalcin   LagrangeBasisDers(n,t,time,alpha);
84211a84d6SLisandro Dalcin   ierr = VecZeroEntries(Xdot);CHKERRQ(ierr);
85211a84d6SLisandro Dalcin   ierr = VecMAXPY(Xdot,n,alpha,vecs);CHKERRQ(ierr);
86211a84d6SLisandro Dalcin   if (shift) *shift = PetscRealPart(alpha[0]);
87211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
88211a84d6SLisandro Dalcin }
89211a84d6SLisandro Dalcin 
90211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_VecLTE(TS ts,PetscInt order,Vec lte)
91211a84d6SLisandro Dalcin {
92211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
93211a84d6SLisandro Dalcin   PetscInt       i,n = order+1;
94211a84d6SLisandro Dalcin   PetscReal      *time = bdf->time;
95211a84d6SLisandro Dalcin   Vec            *vecs = bdf->work;
96211a84d6SLisandro Dalcin   PetscScalar    a[8],b[8],alpha[8];
97211a84d6SLisandro Dalcin   PetscErrorCode ierr;
98211a84d6SLisandro Dalcin 
99211a84d6SLisandro Dalcin   PetscFunctionBegin;
100211a84d6SLisandro Dalcin   LagrangeBasisDers(n+0,time[0],time,a); a[n] =0;
101211a84d6SLisandro Dalcin   LagrangeBasisDers(n+1,time[0],time,b);
102211a84d6SLisandro Dalcin   for (i=0; i<n+1; i++) alpha[i] = (a[i]-b[i])/a[0];
103211a84d6SLisandro Dalcin   ierr = VecZeroEntries(lte);CHKERRQ(ierr);
104211a84d6SLisandro Dalcin   ierr = VecMAXPY(lte,n+1,alpha,vecs);CHKERRQ(ierr);
105211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
106211a84d6SLisandro Dalcin }
107211a84d6SLisandro Dalcin 
108211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_Extrapolate(TS ts,PetscInt order,PetscReal t,Vec X)
109211a84d6SLisandro Dalcin {
110211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
111211a84d6SLisandro Dalcin   PetscInt       n = order+1;
112211a84d6SLisandro Dalcin   PetscReal      *time = bdf->time+1;
113211a84d6SLisandro Dalcin   Vec            *vecs = bdf->work+1;
114211a84d6SLisandro Dalcin   PetscScalar    alpha[7];
115211a84d6SLisandro Dalcin   PetscErrorCode ierr;
116211a84d6SLisandro Dalcin 
117211a84d6SLisandro Dalcin   PetscFunctionBegin;
118211a84d6SLisandro Dalcin   n = PetscMin(n,bdf->n);
119211a84d6SLisandro Dalcin   LagrangeBasisVals(n,t,time,alpha);
120211a84d6SLisandro Dalcin   ierr = VecZeroEntries(X);CHKERRQ(ierr);
121211a84d6SLisandro Dalcin   ierr = VecMAXPY(X,n,alpha,vecs);CHKERRQ(ierr);
122211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
123211a84d6SLisandro Dalcin }
124211a84d6SLisandro Dalcin 
125211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_Interpolate(TS ts,PetscInt order,PetscReal t,Vec X)
126211a84d6SLisandro Dalcin {
127211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
128211a84d6SLisandro Dalcin   PetscInt       n = order+1;
129211a84d6SLisandro Dalcin   PetscReal      *time = bdf->time;
130211a84d6SLisandro Dalcin   Vec            *vecs = bdf->work;
131211a84d6SLisandro Dalcin   PetscScalar    alpha[7];
132211a84d6SLisandro Dalcin   PetscErrorCode ierr;
133211a84d6SLisandro Dalcin 
134211a84d6SLisandro Dalcin   PetscFunctionBegin;
135211a84d6SLisandro Dalcin   LagrangeBasisVals(n,t,time,alpha);
136211a84d6SLisandro Dalcin   ierr = VecZeroEntries(X);CHKERRQ(ierr);
137211a84d6SLisandro Dalcin   ierr = VecMAXPY(X,n,alpha,vecs);CHKERRQ(ierr);
138211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
139211a84d6SLisandro Dalcin }
140211a84d6SLisandro Dalcin 
141211a84d6SLisandro Dalcin static PetscErrorCode TS_SNESSolve(TS ts,Vec b,Vec x)
142211a84d6SLisandro Dalcin {
143211a84d6SLisandro Dalcin   PetscInt       nits,lits;
144211a84d6SLisandro Dalcin   PetscErrorCode ierr;
145211a84d6SLisandro Dalcin 
146211a84d6SLisandro Dalcin   PetscFunctionBegin;
147211a84d6SLisandro Dalcin   ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr);
148211a84d6SLisandro Dalcin   ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr);
149211a84d6SLisandro Dalcin   ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr);
150211a84d6SLisandro Dalcin   ts->snes_its += nits; ts->ksp_its += lits;
151211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
152211a84d6SLisandro Dalcin }
153211a84d6SLisandro Dalcin 
154211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_Restart(TS ts,PetscBool *accept)
155211a84d6SLisandro Dalcin {
156211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
157211a84d6SLisandro Dalcin   PetscErrorCode ierr;
158211a84d6SLisandro Dalcin 
159211a84d6SLisandro Dalcin   PetscFunctionBegin;
160211a84d6SLisandro Dalcin   bdf->k = 1; bdf->n = 0;
161211a84d6SLisandro Dalcin   ierr = TSBDF_Advance(ts,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
162211a84d6SLisandro Dalcin 
163211a84d6SLisandro Dalcin   bdf->time[0] = ts->ptime + ts->time_step/2;
164211a84d6SLisandro Dalcin   ierr = VecCopy(bdf->work[1],bdf->work[0]);CHKERRQ(ierr);
165211a84d6SLisandro Dalcin   ierr = TSPreStage(ts,bdf->time[0]);CHKERRQ(ierr);
166211a84d6SLisandro Dalcin   ierr = TS_SNESSolve(ts,NULL,bdf->work[0]);CHKERRQ(ierr);
167211a84d6SLisandro Dalcin   ierr = TSPostStage(ts,bdf->time[0],0,&bdf->work[0]);CHKERRQ(ierr);
168211a84d6SLisandro Dalcin   ierr = TSAdaptCheckStage(ts->adapt,ts,bdf->time[0],bdf->work[0],accept);CHKERRQ(ierr);
169211a84d6SLisandro Dalcin   if (!*accept) PetscFunctionReturn(0);
170211a84d6SLisandro Dalcin 
171211a84d6SLisandro Dalcin   bdf->k = PetscMin(2,bdf->order); bdf->n++;
172211a84d6SLisandro Dalcin   ierr = VecCopy(bdf->work[0],bdf->work[2]);CHKERRQ(ierr);
173211a84d6SLisandro Dalcin   bdf->time[2] = bdf->time[0];
174211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
175211a84d6SLisandro Dalcin }
176211a84d6SLisandro Dalcin 
177211a84d6SLisandro Dalcin static const char *const BDF_SchemeName[] = {"", "1", "2", "3", "4", "5", "6"};
178211a84d6SLisandro Dalcin 
179211a84d6SLisandro Dalcin static PetscErrorCode TSStep_BDF(TS ts)
180211a84d6SLisandro Dalcin {
181211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
182211a84d6SLisandro Dalcin   PetscInt       rejections = 0;
183211a84d6SLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
184211a84d6SLisandro Dalcin   PetscReal      next_time_step = ts->time_step;
185211a84d6SLisandro Dalcin   PetscErrorCode ierr;
186211a84d6SLisandro Dalcin 
187211a84d6SLisandro Dalcin   PetscFunctionBegin;
188211a84d6SLisandro Dalcin   ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr);
189211a84d6SLisandro Dalcin 
190211a84d6SLisandro Dalcin   if (!ts->steprollback && !ts->steprestart) {
191211a84d6SLisandro Dalcin     bdf->k = PetscMin(bdf->k+1,bdf->order);
192211a84d6SLisandro Dalcin     ierr = TSBDF_Advance(ts,ts->ptime,ts->vec_sol);CHKERRQ(ierr);
193211a84d6SLisandro Dalcin   }
194211a84d6SLisandro Dalcin 
195211a84d6SLisandro Dalcin   bdf->status = TS_STEP_INCOMPLETE;
196211a84d6SLisandro Dalcin   while (!ts->reason && bdf->status != TS_STEP_COMPLETE) {
197211a84d6SLisandro Dalcin 
198211a84d6SLisandro Dalcin     if (ts->steprestart) {
199211a84d6SLisandro Dalcin       ierr = TSBDF_Restart(ts,&stageok);CHKERRQ(ierr);
200211a84d6SLisandro Dalcin       if (!stageok) goto reject_step;
201211a84d6SLisandro Dalcin     }
202211a84d6SLisandro Dalcin 
203211a84d6SLisandro Dalcin     bdf->time[0] = ts->ptime + ts->time_step;
204211a84d6SLisandro Dalcin     ierr = TSBDF_Extrapolate(ts,bdf->k-(accept?0:1),bdf->time[0],bdf->work[0]);CHKERRQ(ierr);
205211a84d6SLisandro Dalcin     ierr = TSPreStage(ts,bdf->time[0]);CHKERRQ(ierr);
206211a84d6SLisandro Dalcin     ierr = TS_SNESSolve(ts,NULL,bdf->work[0]);CHKERRQ(ierr);
207211a84d6SLisandro Dalcin     ierr = TSPostStage(ts,bdf->time[0],0,&bdf->work[0]);CHKERRQ(ierr);
208211a84d6SLisandro Dalcin     ierr = TSAdaptCheckStage(ts->adapt,ts,bdf->time[0],bdf->work[0],&stageok);CHKERRQ(ierr);
209211a84d6SLisandro Dalcin     if (!stageok) goto reject_step;
210211a84d6SLisandro Dalcin 
211211a84d6SLisandro Dalcin     bdf->status = TS_STEP_PENDING;
212211a84d6SLisandro Dalcin     ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
213211a84d6SLisandro Dalcin     ierr = TSAdaptCandidateAdd(ts->adapt,BDF_SchemeName[bdf->k],bdf->k,1,1.0,1.0,PETSC_TRUE);CHKERRQ(ierr);
214211a84d6SLisandro Dalcin     ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr);
215211a84d6SLisandro Dalcin     bdf->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
216211a84d6SLisandro Dalcin     if (!accept) { ts->time_step = next_time_step; goto reject_step; }
217211a84d6SLisandro Dalcin 
218211a84d6SLisandro Dalcin     ierr = VecCopy(bdf->work[0],ts->vec_sol);CHKERRQ(ierr);
219211a84d6SLisandro Dalcin     ts->ptime += ts->time_step;
220211a84d6SLisandro Dalcin     ts->time_step = next_time_step;
221211a84d6SLisandro Dalcin     break;
222211a84d6SLisandro Dalcin 
223211a84d6SLisandro Dalcin   reject_step:
224211a84d6SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
225211a84d6SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
226211a84d6SLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr);
227211a84d6SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
228211a84d6SLisandro Dalcin     }
229211a84d6SLisandro Dalcin   }
230211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
231211a84d6SLisandro Dalcin }
232211a84d6SLisandro Dalcin 
233211a84d6SLisandro Dalcin static PetscErrorCode TSInterpolate_BDF(TS ts,PetscReal t,Vec X)
234211a84d6SLisandro Dalcin {
235211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
236211a84d6SLisandro Dalcin   PetscErrorCode ierr;
237211a84d6SLisandro Dalcin 
238211a84d6SLisandro Dalcin   PetscFunctionBegin;
239211a84d6SLisandro Dalcin   ierr = TSBDF_Interpolate(ts,bdf->k,t,X);CHKERRQ(ierr);
240211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
241211a84d6SLisandro Dalcin }
242211a84d6SLisandro Dalcin 
243211a84d6SLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_BDF(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
244211a84d6SLisandro Dalcin {
245211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
246211a84d6SLisandro Dalcin   PetscInt       k = bdf->k;
2477453f775SEmil Constantinescu   PetscReal      wltea,wlter;
248211a84d6SLisandro Dalcin   Vec            X = bdf->work[0], Y = bdf->vec_lte;
249211a84d6SLisandro Dalcin   PetscErrorCode ierr;
250211a84d6SLisandro Dalcin 
251211a84d6SLisandro Dalcin   PetscFunctionBegin;
252211a84d6SLisandro Dalcin   k = PetscMin(k,bdf->n-1);
253211a84d6SLisandro Dalcin   ierr = TSBDF_VecLTE(ts,k,Y);CHKERRQ(ierr);
254211a84d6SLisandro Dalcin   ierr = VecAXPY(Y,1,X);CHKERRQ(ierr);
2557453f775SEmil Constantinescu   ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr);
256211a84d6SLisandro Dalcin   if (order) *order = k + 1;
257211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
258211a84d6SLisandro Dalcin }
259211a84d6SLisandro Dalcin 
260211a84d6SLisandro Dalcin static PetscErrorCode TSRollBack_BDF(TS ts)
261211a84d6SLisandro Dalcin {
262211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
263211a84d6SLisandro Dalcin   PetscErrorCode ierr;
264211a84d6SLisandro Dalcin 
265211a84d6SLisandro Dalcin   PetscFunctionBegin;
266211a84d6SLisandro Dalcin   ierr = VecCopy(bdf->work[1],ts->vec_sol);CHKERRQ(ierr);
267211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
268211a84d6SLisandro Dalcin }
269211a84d6SLisandro Dalcin 
270211a84d6SLisandro Dalcin static PetscErrorCode SNESTSFormFunction_BDF(PETSC_UNUSED SNES snes,Vec X,Vec F,TS ts)
271211a84d6SLisandro Dalcin {
272211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
273211a84d6SLisandro Dalcin   PetscInt       k = bdf->k;
274211a84d6SLisandro Dalcin   PetscReal      t = bdf->time[0];
275211a84d6SLisandro Dalcin   Vec            V = bdf->vec_dot;
276211a84d6SLisandro Dalcin   PetscErrorCode ierr;
277211a84d6SLisandro Dalcin 
278211a84d6SLisandro Dalcin   PetscFunctionBegin;
279211a84d6SLisandro Dalcin   ierr = TSBDF_VecDot(ts,k,t,X,V,&bdf->shift);CHKERRQ(ierr);
280211a84d6SLisandro Dalcin   /* F = Function(t,X,V) */
281211a84d6SLisandro Dalcin   ierr = TSComputeIFunction(ts,t,X,V,F,PETSC_FALSE);CHKERRQ(ierr);
282211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
283211a84d6SLisandro Dalcin }
284211a84d6SLisandro Dalcin 
285211a84d6SLisandro Dalcin static PetscErrorCode SNESTSFormJacobian_BDF(PETSC_UNUSED SNES snes,
286211a84d6SLisandro Dalcin                                              PETSC_UNUSED Vec X,
287211a84d6SLisandro Dalcin                                              Mat J,Mat P,
288211a84d6SLisandro Dalcin                                              TS ts)
289211a84d6SLisandro Dalcin {
290211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
291211a84d6SLisandro Dalcin   PetscReal      t = bdf->time[0];
292211a84d6SLisandro Dalcin   Vec            V = bdf->vec_dot;
293211a84d6SLisandro Dalcin   PetscReal      dVdX = bdf->shift;
294211a84d6SLisandro Dalcin   PetscErrorCode ierr;
295211a84d6SLisandro Dalcin 
296211a84d6SLisandro Dalcin   PetscFunctionBegin;
297211a84d6SLisandro Dalcin   /* J,P = Jacobian(t,X,V) */
298211a84d6SLisandro Dalcin   ierr = TSComputeIJacobian(ts,t,X,V,dVdX,J,P,PETSC_FALSE);CHKERRQ(ierr);
299211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
300211a84d6SLisandro Dalcin }
301211a84d6SLisandro Dalcin 
302211a84d6SLisandro Dalcin static PetscErrorCode TSReset_BDF(TS ts)
303211a84d6SLisandro Dalcin {
304211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
305211a84d6SLisandro Dalcin   size_t         i,n = sizeof(bdf->work)/sizeof(Vec);
306211a84d6SLisandro Dalcin   PetscErrorCode ierr;
307211a84d6SLisandro Dalcin 
308211a84d6SLisandro Dalcin   PetscFunctionBegin;
309211a84d6SLisandro Dalcin   for (i=0; i<n; i++) {ierr = VecDestroy(&bdf->work[i]);CHKERRQ(ierr);}
310211a84d6SLisandro Dalcin   ierr = VecDestroy(&bdf->vec_dot);CHKERRQ(ierr);
311211a84d6SLisandro Dalcin   ierr = VecDestroy(&bdf->vec_lte);CHKERRQ(ierr);
312211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
313211a84d6SLisandro Dalcin }
314211a84d6SLisandro Dalcin 
315211a84d6SLisandro Dalcin static PetscErrorCode TSDestroy_BDF(TS ts)
316211a84d6SLisandro Dalcin {
317211a84d6SLisandro Dalcin   PetscErrorCode ierr;
318211a84d6SLisandro Dalcin 
319211a84d6SLisandro Dalcin   PetscFunctionBegin;
320211a84d6SLisandro Dalcin   ierr = TSReset_BDF(ts);CHKERRQ(ierr);
321211a84d6SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
322211a84d6SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBDFSetOrder_C",NULL);CHKERRQ(ierr);
323211a84d6SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBDFGetOrder_C",NULL);CHKERRQ(ierr);
324211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
325211a84d6SLisandro Dalcin }
326211a84d6SLisandro Dalcin 
327211a84d6SLisandro Dalcin static PetscErrorCode TSSetUp_BDF(TS ts)
328211a84d6SLisandro Dalcin {
329211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
330211a84d6SLisandro Dalcin   size_t         i,n = sizeof(bdf->work)/sizeof(Vec);
3312ffb9264SLisandro Dalcin   PetscReal      low,high,two = 2;
332211a84d6SLisandro Dalcin   PetscErrorCode ierr;
333211a84d6SLisandro Dalcin 
334211a84d6SLisandro Dalcin   PetscFunctionBegin;
335211a84d6SLisandro Dalcin   bdf->k = bdf->n = 0;
336211a84d6SLisandro Dalcin   for (i=0; i<n; i++) {ierr = VecDuplicate(ts->vec_sol,&bdf->work[i]);CHKERRQ(ierr);}
337211a84d6SLisandro Dalcin   ierr = VecDuplicate(ts->vec_sol,&bdf->vec_dot);CHKERRQ(ierr);
338211a84d6SLisandro Dalcin   ierr = VecDuplicate(ts->vec_sol,&bdf->vec_lte);CHKERRQ(ierr);
339211a84d6SLisandro Dalcin 
340211a84d6SLisandro Dalcin   ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr);
341211a84d6SLisandro Dalcin   ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr);
3421917a363SLisandro Dalcin   ierr = TSAdaptGetClip(ts->adapt,&low,&high);CHKERRQ(ierr);
3432ffb9264SLisandro Dalcin   ierr = TSAdaptSetClip(ts->adapt,low,PetscMin(high,two));CHKERRQ(ierr);
344211a84d6SLisandro Dalcin 
345211a84d6SLisandro Dalcin   ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr);
346211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
347211a84d6SLisandro Dalcin }
348211a84d6SLisandro Dalcin 
349211a84d6SLisandro Dalcin static PetscErrorCode TSSetFromOptions_BDF(PetscOptionItems *PetscOptionsObject,TS ts)
350211a84d6SLisandro Dalcin {
351211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
352211a84d6SLisandro Dalcin   PetscErrorCode ierr;
353211a84d6SLisandro Dalcin 
354211a84d6SLisandro Dalcin   PetscFunctionBegin;
355211a84d6SLisandro Dalcin   ierr = PetscOptionsHead(PetscOptionsObject,"BDF ODE solver options");CHKERRQ(ierr);
356211a84d6SLisandro Dalcin   {
357211a84d6SLisandro Dalcin     PetscBool flg;
358211a84d6SLisandro Dalcin     PetscInt  order = bdf->order;
359211a84d6SLisandro Dalcin     ierr = PetscOptionsInt("-ts_bdf_order","Order of the BDF method","TSBDFSetOrder",order,&order,&flg);CHKERRQ(ierr);
360211a84d6SLisandro Dalcin     if (flg) {ierr = TSBDFSetOrder(ts,order);CHKERRQ(ierr);}
361211a84d6SLisandro Dalcin   }
362211a84d6SLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
363211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
364211a84d6SLisandro Dalcin }
365211a84d6SLisandro Dalcin 
366211a84d6SLisandro Dalcin static PetscErrorCode TSView_BDF(TS ts,PetscViewer viewer)
367211a84d6SLisandro Dalcin {
368211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
369211a84d6SLisandro Dalcin   PetscBool      iascii;
370211a84d6SLisandro Dalcin   PetscErrorCode ierr;
371211a84d6SLisandro Dalcin 
372211a84d6SLisandro Dalcin   PetscFunctionBegin;
373211a84d6SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
374211a84d6SLisandro Dalcin   if (iascii) {
375211a84d6SLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer,"  Order=%D\n",bdf->order);CHKERRQ(ierr);
376211a84d6SLisandro Dalcin   }
377211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
378211a84d6SLisandro Dalcin }
379211a84d6SLisandro Dalcin 
380211a84d6SLisandro Dalcin /* ------------------------------------------------------------ */
381211a84d6SLisandro Dalcin 
382211a84d6SLisandro Dalcin static PetscErrorCode TSBDFSetOrder_BDF(TS ts,PetscInt order)
383211a84d6SLisandro Dalcin {
384211a84d6SLisandro Dalcin   TS_BDF *bdf = (TS_BDF*)ts->data;
385211a84d6SLisandro Dalcin 
386211a84d6SLisandro Dalcin   PetscFunctionBegin;
387211a84d6SLisandro Dalcin   if (order == bdf->order) PetscFunctionReturn(0);
388211a84d6SLisandro Dalcin   if (order < 1 || order > 6) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"BDF Order %D not implemented",order);
389211a84d6SLisandro Dalcin   bdf->order = order;
390211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
391211a84d6SLisandro Dalcin }
392211a84d6SLisandro Dalcin 
393211a84d6SLisandro Dalcin static PetscErrorCode TSBDFGetOrder_BDF(TS ts,PetscInt *order)
394211a84d6SLisandro Dalcin {
395211a84d6SLisandro Dalcin   TS_BDF *bdf = (TS_BDF*)ts->data;
396211a84d6SLisandro Dalcin 
397211a84d6SLisandro Dalcin   PetscFunctionBegin;
398211a84d6SLisandro Dalcin   *order = bdf->order;
399211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
400211a84d6SLisandro Dalcin }
401211a84d6SLisandro Dalcin 
402211a84d6SLisandro Dalcin /* ------------------------------------------------------------ */
403211a84d6SLisandro Dalcin 
404211a84d6SLisandro Dalcin /*MC
405211a84d6SLisandro Dalcin       TSBDF - DAE solver using BDF methods
406211a84d6SLisandro Dalcin 
407211a84d6SLisandro Dalcin   Level: beginner
408211a84d6SLisandro Dalcin 
409211a84d6SLisandro Dalcin .seealso:  TS, TSCreate(), TSSetType()
410211a84d6SLisandro Dalcin M*/
411211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSCreate_BDF(TS ts)
412211a84d6SLisandro Dalcin {
413211a84d6SLisandro Dalcin   TS_BDF         *bdf;
414211a84d6SLisandro Dalcin   PetscErrorCode ierr;
415211a84d6SLisandro Dalcin 
416211a84d6SLisandro Dalcin   PetscFunctionBegin;
417211a84d6SLisandro Dalcin   ts->ops->reset          = TSReset_BDF;
418211a84d6SLisandro Dalcin   ts->ops->destroy        = TSDestroy_BDF;
419211a84d6SLisandro Dalcin   ts->ops->view           = TSView_BDF;
420211a84d6SLisandro Dalcin   ts->ops->setup          = TSSetUp_BDF;
421211a84d6SLisandro Dalcin   ts->ops->setfromoptions = TSSetFromOptions_BDF;
422211a84d6SLisandro Dalcin   ts->ops->step           = TSStep_BDF;
423211a84d6SLisandro Dalcin   ts->ops->evaluatewlte   = TSEvaluateWLTE_BDF;
424211a84d6SLisandro Dalcin   ts->ops->rollback       = TSRollBack_BDF;
425211a84d6SLisandro Dalcin   ts->ops->interpolate    = TSInterpolate_BDF;
426211a84d6SLisandro Dalcin   ts->ops->snesfunction   = SNESTSFormFunction_BDF;
427211a84d6SLisandro Dalcin   ts->ops->snesjacobian   = SNESTSFormJacobian_BDF;
4282ffb9264SLisandro Dalcin   ts->default_adapt_type  = TSADAPTBASIC;
429211a84d6SLisandro Dalcin 
430*efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
431*efd4aadfSBarry Smith 
432211a84d6SLisandro Dalcin   ierr = PetscNewLog(ts,&bdf);CHKERRQ(ierr);
433211a84d6SLisandro Dalcin   ts->data = (void*)bdf;
434211a84d6SLisandro Dalcin 
435211a84d6SLisandro Dalcin   bdf->order  = 2;
436211a84d6SLisandro Dalcin   bdf->status = TS_STEP_COMPLETE;
437211a84d6SLisandro Dalcin 
438211a84d6SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBDFSetOrder_C",TSBDFSetOrder_BDF);CHKERRQ(ierr);
439211a84d6SLisandro Dalcin   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSBDFGetOrder_C",TSBDFGetOrder_BDF);CHKERRQ(ierr);
440211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
441211a84d6SLisandro Dalcin }
442211a84d6SLisandro Dalcin 
443211a84d6SLisandro Dalcin /* ------------------------------------------------------------ */
444211a84d6SLisandro Dalcin 
445211a84d6SLisandro Dalcin /*@
446211a84d6SLisandro Dalcin   TSBDFSetOrder - Set the order of the BDF method
447211a84d6SLisandro Dalcin 
448211a84d6SLisandro Dalcin   Logically Collective on TS
449211a84d6SLisandro Dalcin 
450211a84d6SLisandro Dalcin   Input Parameter:
451211a84d6SLisandro Dalcin +  ts - timestepping context
452211a84d6SLisandro Dalcin -  order - order of the method
453211a84d6SLisandro Dalcin 
454211a84d6SLisandro Dalcin   Options Database:
455211a84d6SLisandro Dalcin .  -ts_bdf_order <order>
456211a84d6SLisandro Dalcin 
457211a84d6SLisandro Dalcin   Level: intermediate
458211a84d6SLisandro Dalcin 
459211a84d6SLisandro Dalcin @*/
460211a84d6SLisandro Dalcin PetscErrorCode TSBDFSetOrder(TS ts,PetscInt order)
461211a84d6SLisandro Dalcin {
462211a84d6SLisandro Dalcin   PetscErrorCode ierr;
463211a84d6SLisandro Dalcin 
464211a84d6SLisandro Dalcin   PetscFunctionBegin;
465211a84d6SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
466211a84d6SLisandro Dalcin   PetscValidLogicalCollectiveInt(ts,order,2);
467211a84d6SLisandro Dalcin   ierr = PetscTryMethod(ts,"TSBDFSetOrder_C",(TS,PetscInt),(ts,order));CHKERRQ(ierr);
468211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
469211a84d6SLisandro Dalcin }
470211a84d6SLisandro Dalcin 
471211a84d6SLisandro Dalcin /*@
472211a84d6SLisandro Dalcin   TSBDFGetOrder - Get the order of the BDF method
473211a84d6SLisandro Dalcin 
474211a84d6SLisandro Dalcin   Not Collective
475211a84d6SLisandro Dalcin 
476211a84d6SLisandro Dalcin   Input Parameter:
477211a84d6SLisandro Dalcin .  ts - timestepping context
478211a84d6SLisandro Dalcin 
479211a84d6SLisandro Dalcin   Output Parameter:
480211a84d6SLisandro Dalcin .  order - order of the method
481211a84d6SLisandro Dalcin 
482211a84d6SLisandro Dalcin   Level: intermediate
483211a84d6SLisandro Dalcin 
484211a84d6SLisandro Dalcin @*/
485211a84d6SLisandro Dalcin PetscErrorCode TSBDFGetOrder(TS ts,PetscInt *order)
486211a84d6SLisandro Dalcin {
487211a84d6SLisandro Dalcin   PetscErrorCode ierr;
488211a84d6SLisandro Dalcin 
489211a84d6SLisandro Dalcin   PetscFunctionBegin;
490211a84d6SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
491211a84d6SLisandro Dalcin   PetscValidIntPointer(order,2);
492211a84d6SLisandro Dalcin   ierr = PetscUseMethod(ts,"TSBDFGetOrder_C",(TS,PetscInt*),(ts,order));CHKERRQ(ierr);
493211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
494211a84d6SLisandro Dalcin }
495