xref: /petsc/src/ts/impls/bdf/bdf.c (revision cac4c232dc4f93991e342196e27ef7b0655dac7b)
1211a84d6SLisandro Dalcin /*
2211a84d6SLisandro Dalcin   Code for timestepping with BDF methods
3211a84d6SLisandro Dalcin */
4211a84d6SLisandro Dalcin #include <petsc/private/tsimpl.h>  /*I "petscts.h" I*/
51117961dSLisandro Dalcin #include <petscdm.h>
6211a84d6SLisandro Dalcin 
7211a84d6SLisandro Dalcin static PetscBool  cited = PETSC_FALSE;
8211a84d6SLisandro Dalcin static const char citation[] =
9211a84d6SLisandro Dalcin   "@book{Brenan1995,\n"
10211a84d6SLisandro Dalcin   "  title     = {Numerical Solution of Initial-Value Problems in Differential-Algebraic Equations},\n"
11211a84d6SLisandro Dalcin   "  author    = {Brenan, K. and Campbell, S. and Petzold, L.},\n"
12211a84d6SLisandro Dalcin   "  publisher = {Society for Industrial and Applied Mathematics},\n"
13211a84d6SLisandro Dalcin   "  year      = {1995},\n"
14211a84d6SLisandro Dalcin   "  doi       = {10.1137/1.9781611971224},\n}\n";
15211a84d6SLisandro Dalcin 
16211a84d6SLisandro Dalcin typedef struct {
17211a84d6SLisandro Dalcin   PetscInt  k,n;
18211a84d6SLisandro Dalcin   PetscReal time[6+2];
19211a84d6SLisandro Dalcin   Vec       work[6+2];
20e3c11fc1SJed Brown   Vec       tvwork[6+2];
21211a84d6SLisandro Dalcin   PetscReal shift;
22e3c11fc1SJed Brown   Vec       vec_dot;            /* Xdot when !transientvar, else Cdot where C(X) is the transient variable. */
231117961dSLisandro Dalcin   Vec       vec_wrk;
24211a84d6SLisandro Dalcin   Vec       vec_lte;
25211a84d6SLisandro Dalcin 
26e3c11fc1SJed Brown   PetscBool    transientvar;
27211a84d6SLisandro Dalcin   PetscInt     order;
28211a84d6SLisandro Dalcin   TSStepStatus status;
29211a84d6SLisandro Dalcin } TS_BDF;
30211a84d6SLisandro Dalcin 
31e3c11fc1SJed Brown /* Compute Lagrange polynomials on T[:n] evaluated at t.
32e3c11fc1SJed Brown  * If one has data (T[i], Y[i]), then the interpolation/extrapolation is f(t) = \sum_i L[i]*Y[i].
33e3c11fc1SJed Brown  */
349fbee547SJacob Faibussowitsch static inline void LagrangeBasisVals(PetscInt n,PetscReal t,const PetscReal T[],PetscScalar L[])
35211a84d6SLisandro Dalcin {
36211a84d6SLisandro Dalcin   PetscInt k,j;
37211a84d6SLisandro Dalcin   for (k=0; k<n; k++)
38211a84d6SLisandro Dalcin     for (L[k]=1, j=0; j<n; j++)
39211a84d6SLisandro Dalcin       if (j != k)
40211a84d6SLisandro Dalcin         L[k] *= (t - T[j])/(T[k] - T[j]);
41211a84d6SLisandro Dalcin }
42211a84d6SLisandro Dalcin 
439fbee547SJacob Faibussowitsch static inline void LagrangeBasisDers(PetscInt n,PetscReal t,const PetscReal T[],PetscScalar dL[])
44211a84d6SLisandro Dalcin {
45211a84d6SLisandro Dalcin   PetscInt  k,j,i;
46211a84d6SLisandro Dalcin   for (k=0; k<n; k++)
47211a84d6SLisandro Dalcin     for (dL[k]=0, j=0; j<n; j++)
48211a84d6SLisandro Dalcin       if (j != k) {
49211a84d6SLisandro Dalcin         PetscReal L = 1/(T[k] - T[j]);
50211a84d6SLisandro Dalcin         for (i=0; i<n; i++)
51211a84d6SLisandro Dalcin           if (i != j && i != k)
52211a84d6SLisandro Dalcin             L *= (t - T[i])/(T[k] - T[i]);
53211a84d6SLisandro Dalcin         dL[k] += L;
54211a84d6SLisandro Dalcin       }
55211a84d6SLisandro Dalcin }
56211a84d6SLisandro Dalcin 
571117961dSLisandro Dalcin static PetscErrorCode TSBDF_GetVecs(TS ts,DM dm,Vec *Xdot,Vec *Ydot)
581117961dSLisandro Dalcin {
591117961dSLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
601117961dSLisandro Dalcin 
611117961dSLisandro Dalcin   PetscFunctionBegin;
621117961dSLisandro Dalcin   if (dm && dm != ts->dm) {
639566063dSJacob Faibussowitsch     PetscCall(DMGetNamedGlobalVector(dm,"TSBDF_Vec_Xdot",Xdot));
649566063dSJacob Faibussowitsch     PetscCall(DMGetNamedGlobalVector(dm,"TSBDF_Vec_Ydot",Ydot));
651117961dSLisandro Dalcin   } else {
661117961dSLisandro Dalcin     *Xdot = bdf->vec_dot;
671117961dSLisandro Dalcin     *Ydot = bdf->vec_wrk;
681117961dSLisandro Dalcin   }
691117961dSLisandro Dalcin   PetscFunctionReturn(0);
701117961dSLisandro Dalcin }
711117961dSLisandro Dalcin 
721117961dSLisandro Dalcin static PetscErrorCode TSBDF_RestoreVecs(TS ts,DM dm,Vec *Xdot,Vec *Ydot)
731117961dSLisandro Dalcin {
741117961dSLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
751117961dSLisandro Dalcin 
761117961dSLisandro Dalcin   PetscFunctionBegin;
771117961dSLisandro Dalcin   if (dm && dm != ts->dm) {
789566063dSJacob Faibussowitsch     PetscCall(DMRestoreNamedGlobalVector(dm,"TSBDF_Vec_Xdot",Xdot));
799566063dSJacob Faibussowitsch     PetscCall(DMRestoreNamedGlobalVector(dm,"TSBDF_Vec_Ydot",Ydot));
801117961dSLisandro Dalcin   } else {
813c633725SBarry Smith     PetscCheck(*Xdot == bdf->vec_dot,PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"Vec does not match the cache");
823c633725SBarry Smith     PetscCheck(*Ydot == bdf->vec_wrk,PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"Vec does not match the cache");
831117961dSLisandro Dalcin     *Xdot = NULL;
841117961dSLisandro Dalcin     *Ydot = NULL;
851117961dSLisandro Dalcin   }
861117961dSLisandro Dalcin   PetscFunctionReturn(0);
871117961dSLisandro Dalcin }
881117961dSLisandro Dalcin 
891117961dSLisandro Dalcin static PetscErrorCode DMCoarsenHook_TSBDF(DM fine,DM coarse,void *ctx)
901117961dSLisandro Dalcin {
911117961dSLisandro Dalcin   PetscFunctionBegin;
921117961dSLisandro Dalcin   PetscFunctionReturn(0);
931117961dSLisandro Dalcin }
941117961dSLisandro Dalcin 
951117961dSLisandro Dalcin static PetscErrorCode DMRestrictHook_TSBDF(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
961117961dSLisandro Dalcin {
971117961dSLisandro Dalcin   TS             ts = (TS)ctx;
981117961dSLisandro Dalcin   Vec            Ydot,Ydot_c;
991117961dSLisandro Dalcin   Vec            Xdot,Xdot_c;
1001117961dSLisandro Dalcin 
1011117961dSLisandro Dalcin   PetscFunctionBegin;
1029566063dSJacob Faibussowitsch   PetscCall(TSBDF_GetVecs(ts,fine,&Xdot,&Ydot));
1039566063dSJacob Faibussowitsch   PetscCall(TSBDF_GetVecs(ts,coarse,&Xdot_c,&Ydot_c));
1041117961dSLisandro Dalcin 
1059566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct,Ydot,Ydot_c));
1069566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(Ydot_c,rscale,Ydot_c));
1071117961dSLisandro Dalcin 
1089566063dSJacob Faibussowitsch   PetscCall(TSBDF_RestoreVecs(ts,fine,&Xdot,&Ydot));
1099566063dSJacob Faibussowitsch   PetscCall(TSBDF_RestoreVecs(ts,coarse,&Xdot_c,&Ydot_c));
1101117961dSLisandro Dalcin   PetscFunctionReturn(0);
1111117961dSLisandro Dalcin }
1121117961dSLisandro Dalcin 
113211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_Advance(TS ts,PetscReal t,Vec X)
114211a84d6SLisandro Dalcin {
115211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
116211a84d6SLisandro Dalcin   PetscInt       i,n = (PetscInt)(sizeof(bdf->work)/sizeof(Vec));
117e3c11fc1SJed Brown   Vec            tail = bdf->work[n-1],tvtail = bdf->tvwork[n-1];
118211a84d6SLisandro Dalcin 
119211a84d6SLisandro Dalcin   PetscFunctionBegin;
120211a84d6SLisandro Dalcin   for (i=n-1; i>=2; i--) {
121211a84d6SLisandro Dalcin     bdf->time[i] = bdf->time[i-1];
122211a84d6SLisandro Dalcin     bdf->work[i] = bdf->work[i-1];
123e3c11fc1SJed Brown     bdf->tvwork[i] = bdf->tvwork[i-1];
124211a84d6SLisandro Dalcin   }
125211a84d6SLisandro Dalcin   bdf->n       = PetscMin(bdf->n+1,n-1);
126211a84d6SLisandro Dalcin   bdf->time[1] = t;
127211a84d6SLisandro Dalcin   bdf->work[1] = tail;
128e3c11fc1SJed Brown   bdf->tvwork[1] = tvtail;
1299566063dSJacob Faibussowitsch   PetscCall(VecCopy(X,tail));
1309566063dSJacob Faibussowitsch   PetscCall(TSComputeTransientVariable(ts,tail,tvtail));
131211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
132211a84d6SLisandro Dalcin }
133211a84d6SLisandro Dalcin 
134211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_VecLTE(TS ts,PetscInt order,Vec lte)
135211a84d6SLisandro Dalcin {
136211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
137211a84d6SLisandro Dalcin   PetscInt       i,n = order+1;
138211a84d6SLisandro Dalcin   PetscReal      *time = bdf->time;
139211a84d6SLisandro Dalcin   Vec            *vecs = bdf->work;
140211a84d6SLisandro Dalcin   PetscScalar    a[8],b[8],alpha[8];
141211a84d6SLisandro Dalcin 
142211a84d6SLisandro Dalcin   PetscFunctionBegin;
143211a84d6SLisandro Dalcin   LagrangeBasisDers(n+0,time[0],time,a); a[n] =0;
144211a84d6SLisandro Dalcin   LagrangeBasisDers(n+1,time[0],time,b);
145211a84d6SLisandro Dalcin   for (i=0; i<n+1; i++) alpha[i] = (a[i]-b[i])/a[0];
1469566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(lte));
1479566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(lte,n+1,alpha,vecs));
148211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
149211a84d6SLisandro Dalcin }
150211a84d6SLisandro Dalcin 
151211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_Extrapolate(TS ts,PetscInt order,PetscReal t,Vec X)
152211a84d6SLisandro Dalcin {
153211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
154211a84d6SLisandro Dalcin   PetscInt       n = order+1;
155211a84d6SLisandro Dalcin   PetscReal      *time = bdf->time+1;
156211a84d6SLisandro Dalcin   Vec            *vecs = bdf->work+1;
157211a84d6SLisandro Dalcin   PetscScalar    alpha[7];
158211a84d6SLisandro Dalcin 
159211a84d6SLisandro Dalcin   PetscFunctionBegin;
160211a84d6SLisandro Dalcin   n = PetscMin(n,bdf->n);
161211a84d6SLisandro Dalcin   LagrangeBasisVals(n,t,time,alpha);
1629566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(X));
1639566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(X,n,alpha,vecs));
164211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
165211a84d6SLisandro Dalcin }
166211a84d6SLisandro Dalcin 
167211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_Interpolate(TS ts,PetscInt order,PetscReal t,Vec X)
168211a84d6SLisandro Dalcin {
169211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
170211a84d6SLisandro Dalcin   PetscInt       n = order+1;
171211a84d6SLisandro Dalcin   PetscReal      *time = bdf->time;
172211a84d6SLisandro Dalcin   Vec            *vecs = bdf->work;
173211a84d6SLisandro Dalcin   PetscScalar    alpha[7];
174211a84d6SLisandro Dalcin 
175211a84d6SLisandro Dalcin   PetscFunctionBegin;
176211a84d6SLisandro Dalcin   LagrangeBasisVals(n,t,time,alpha);
1779566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(X));
1789566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(X,n,alpha,vecs));
179211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
180211a84d6SLisandro Dalcin }
181211a84d6SLisandro Dalcin 
182e3c11fc1SJed Brown /* Compute the affine term V0 such that Xdot = shift*X + V0.
183e3c11fc1SJed Brown  *
184e3c11fc1SJed Brown  * When using transient variables, we're computing Cdot = shift*C(X) + V0, and thus choose a linear combination of tvwork.
185e3c11fc1SJed Brown  */
1861117961dSLisandro Dalcin static PetscErrorCode TSBDF_PreSolve(TS ts)
1871117961dSLisandro Dalcin {
1881117961dSLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
1891117961dSLisandro Dalcin   PetscInt       i,n = PetscMax(bdf->k,1) + 1;
1901117961dSLisandro Dalcin   Vec            V,V0;
1911117961dSLisandro Dalcin   Vec            vecs[7];
1921117961dSLisandro Dalcin   PetscScalar    alpha[7];
1931117961dSLisandro Dalcin 
1941117961dSLisandro Dalcin   PetscFunctionBegin;
1959566063dSJacob Faibussowitsch   PetscCall(TSBDF_GetVecs(ts,NULL,&V,&V0));
1961117961dSLisandro Dalcin   LagrangeBasisDers(n,bdf->time[0],bdf->time,alpha);
197e3c11fc1SJed Brown   for (i=1; i<n; i++) {
198e3c11fc1SJed Brown     vecs[i] = bdf->transientvar ? bdf->tvwork[i] : bdf->work[i];
199e3c11fc1SJed Brown   }
2009566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(V0));
2019566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(V0,n-1,alpha+1,vecs+1));
2021117961dSLisandro Dalcin   bdf->shift = PetscRealPart(alpha[0]);
2039566063dSJacob Faibussowitsch   PetscCall(TSBDF_RestoreVecs(ts,NULL,&V,&V0));
2041117961dSLisandro Dalcin   PetscFunctionReturn(0);
2051117961dSLisandro Dalcin }
2061117961dSLisandro Dalcin 
207470880abSPatrick Sanan static PetscErrorCode TSBDF_SNESSolve(TS ts,Vec b,Vec x)
208211a84d6SLisandro Dalcin {
209211a84d6SLisandro Dalcin   PetscInt       nits,lits;
210211a84d6SLisandro Dalcin 
211211a84d6SLisandro Dalcin   PetscFunctionBegin;
2129566063dSJacob Faibussowitsch   PetscCall(TSBDF_PreSolve(ts));
2139566063dSJacob Faibussowitsch   PetscCall(SNESSolve(ts->snes,b,x));
2149566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(ts->snes,&nits));
2159566063dSJacob Faibussowitsch   PetscCall(SNESGetLinearSolveIterations(ts->snes,&lits));
216211a84d6SLisandro Dalcin   ts->snes_its += nits; ts->ksp_its += lits;
217211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
218211a84d6SLisandro Dalcin }
219211a84d6SLisandro Dalcin 
220211a84d6SLisandro Dalcin static PetscErrorCode TSBDF_Restart(TS ts,PetscBool *accept)
221211a84d6SLisandro Dalcin {
222211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
223211a84d6SLisandro Dalcin 
224211a84d6SLisandro Dalcin   PetscFunctionBegin;
225211a84d6SLisandro Dalcin   bdf->k = 1; bdf->n = 0;
2269566063dSJacob Faibussowitsch   PetscCall(TSBDF_Advance(ts,ts->ptime,ts->vec_sol));
227211a84d6SLisandro Dalcin 
228211a84d6SLisandro Dalcin   bdf->time[0] = ts->ptime + ts->time_step/2;
2299566063dSJacob Faibussowitsch   PetscCall(VecCopy(bdf->work[1],bdf->work[0]));
2309566063dSJacob Faibussowitsch   PetscCall(TSPreStage(ts,bdf->time[0]));
2319566063dSJacob Faibussowitsch   PetscCall(TSBDF_SNESSolve(ts,NULL,bdf->work[0]));
2329566063dSJacob Faibussowitsch   PetscCall(TSPostStage(ts,bdf->time[0],0,&bdf->work[0]));
2339566063dSJacob Faibussowitsch   PetscCall(TSAdaptCheckStage(ts->adapt,ts,bdf->time[0],bdf->work[0],accept));
234211a84d6SLisandro Dalcin   if (!*accept) PetscFunctionReturn(0);
235211a84d6SLisandro Dalcin 
236211a84d6SLisandro Dalcin   bdf->k = PetscMin(2,bdf->order); bdf->n++;
2379566063dSJacob Faibussowitsch   PetscCall(VecCopy(bdf->work[0],bdf->work[2]));
238211a84d6SLisandro Dalcin   bdf->time[2] = bdf->time[0];
2399566063dSJacob Faibussowitsch   PetscCall(TSComputeTransientVariable(ts,bdf->work[2],bdf->tvwork[2]));
240211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
241211a84d6SLisandro Dalcin }
242211a84d6SLisandro Dalcin 
243211a84d6SLisandro Dalcin static const char *const BDF_SchemeName[] = {"", "1", "2", "3", "4", "5", "6"};
244211a84d6SLisandro Dalcin 
245211a84d6SLisandro Dalcin static PetscErrorCode TSStep_BDF(TS ts)
246211a84d6SLisandro Dalcin {
247211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
248211a84d6SLisandro Dalcin   PetscInt       rejections = 0;
249211a84d6SLisandro Dalcin   PetscBool      stageok,accept = PETSC_TRUE;
250211a84d6SLisandro Dalcin   PetscReal      next_time_step = ts->time_step;
251211a84d6SLisandro Dalcin 
252211a84d6SLisandro Dalcin   PetscFunctionBegin;
2539566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(citation,&cited));
254211a84d6SLisandro Dalcin 
255211a84d6SLisandro Dalcin   if (!ts->steprollback && !ts->steprestart) {
256211a84d6SLisandro Dalcin     bdf->k = PetscMin(bdf->k+1,bdf->order);
2579566063dSJacob Faibussowitsch     PetscCall(TSBDF_Advance(ts,ts->ptime,ts->vec_sol));
258211a84d6SLisandro Dalcin   }
259211a84d6SLisandro Dalcin 
260211a84d6SLisandro Dalcin   bdf->status = TS_STEP_INCOMPLETE;
261211a84d6SLisandro Dalcin   while (!ts->reason && bdf->status != TS_STEP_COMPLETE) {
262211a84d6SLisandro Dalcin 
263211a84d6SLisandro Dalcin     if (ts->steprestart) {
2649566063dSJacob Faibussowitsch       PetscCall(TSBDF_Restart(ts,&stageok));
265211a84d6SLisandro Dalcin       if (!stageok) goto reject_step;
266211a84d6SLisandro Dalcin     }
267211a84d6SLisandro Dalcin 
268211a84d6SLisandro Dalcin     bdf->time[0] = ts->ptime + ts->time_step;
2699566063dSJacob Faibussowitsch     PetscCall(TSBDF_Extrapolate(ts,bdf->k-(accept?0:1),bdf->time[0],bdf->work[0]));
2709566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts,bdf->time[0]));
2719566063dSJacob Faibussowitsch     PetscCall(TSBDF_SNESSolve(ts,NULL,bdf->work[0]));
2729566063dSJacob Faibussowitsch     PetscCall(TSPostStage(ts,bdf->time[0],0,&bdf->work[0]));
2739566063dSJacob Faibussowitsch     PetscCall(TSAdaptCheckStage(ts->adapt,ts,bdf->time[0],bdf->work[0],&stageok));
274211a84d6SLisandro Dalcin     if (!stageok) goto reject_step;
275211a84d6SLisandro Dalcin 
276211a84d6SLisandro Dalcin     bdf->status = TS_STEP_PENDING;
2779566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidatesClear(ts->adapt));
2789566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidateAdd(ts->adapt,BDF_SchemeName[bdf->k],bdf->k,1,1.0,1.0,PETSC_TRUE));
2799566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept));
280211a84d6SLisandro Dalcin     bdf->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
281211a84d6SLisandro Dalcin     if (!accept) { ts->time_step = next_time_step; goto reject_step; }
282211a84d6SLisandro Dalcin 
2839566063dSJacob Faibussowitsch     PetscCall(VecCopy(bdf->work[0],ts->vec_sol));
284211a84d6SLisandro Dalcin     ts->ptime += ts->time_step;
285211a84d6SLisandro Dalcin     ts->time_step = next_time_step;
286211a84d6SLisandro Dalcin     break;
287211a84d6SLisandro Dalcin 
288211a84d6SLisandro Dalcin   reject_step:
289211a84d6SLisandro Dalcin     ts->reject++; accept = PETSC_FALSE;
290211a84d6SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
2919566063dSJacob Faibussowitsch       PetscCall(PetscInfo(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections));
292211a84d6SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
293211a84d6SLisandro Dalcin     }
294211a84d6SLisandro Dalcin   }
295211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
296211a84d6SLisandro Dalcin }
297211a84d6SLisandro Dalcin 
298211a84d6SLisandro Dalcin static PetscErrorCode TSInterpolate_BDF(TS ts,PetscReal t,Vec X)
299211a84d6SLisandro Dalcin {
300211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
301211a84d6SLisandro Dalcin 
302211a84d6SLisandro Dalcin   PetscFunctionBegin;
3039566063dSJacob Faibussowitsch   PetscCall(TSBDF_Interpolate(ts,bdf->k,t,X));
304211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
305211a84d6SLisandro Dalcin }
306211a84d6SLisandro Dalcin 
307211a84d6SLisandro Dalcin static PetscErrorCode TSEvaluateWLTE_BDF(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte)
308211a84d6SLisandro Dalcin {
309211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
310211a84d6SLisandro Dalcin   PetscInt       k = bdf->k;
3117453f775SEmil Constantinescu   PetscReal      wltea,wlter;
312211a84d6SLisandro Dalcin   Vec            X = bdf->work[0], Y = bdf->vec_lte;
313211a84d6SLisandro Dalcin 
314211a84d6SLisandro Dalcin   PetscFunctionBegin;
315211a84d6SLisandro Dalcin   k = PetscMin(k,bdf->n-1);
3169566063dSJacob Faibussowitsch   PetscCall(TSBDF_VecLTE(ts,k,Y));
3179566063dSJacob Faibussowitsch   PetscCall(VecAXPY(Y,1,X));
3189566063dSJacob Faibussowitsch   PetscCall(TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter));
319211a84d6SLisandro Dalcin   if (order) *order = k + 1;
320211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
321211a84d6SLisandro Dalcin }
322211a84d6SLisandro Dalcin 
323211a84d6SLisandro Dalcin static PetscErrorCode TSRollBack_BDF(TS ts)
324211a84d6SLisandro Dalcin {
325211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
326211a84d6SLisandro Dalcin 
327211a84d6SLisandro Dalcin   PetscFunctionBegin;
3289566063dSJacob Faibussowitsch   PetscCall(VecCopy(bdf->work[1],ts->vec_sol));
329211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
330211a84d6SLisandro Dalcin }
331211a84d6SLisandro Dalcin 
3321117961dSLisandro Dalcin static PetscErrorCode SNESTSFormFunction_BDF(SNES snes,Vec X,Vec F,TS ts)
333211a84d6SLisandro Dalcin {
334211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
3351117961dSLisandro Dalcin   DM             dm, dmsave = ts->dm;
336211a84d6SLisandro Dalcin   PetscReal      t = bdf->time[0];
3371117961dSLisandro Dalcin   PetscReal      shift = bdf->shift;
3381117961dSLisandro Dalcin   Vec            V,V0;
339211a84d6SLisandro Dalcin 
340211a84d6SLisandro Dalcin   PetscFunctionBegin;
3419566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
3429566063dSJacob Faibussowitsch   PetscCall(TSBDF_GetVecs(ts,dm,&V,&V0));
343e3c11fc1SJed Brown   if (bdf->transientvar) {      /* shift*C(X) + V0 */
3449566063dSJacob Faibussowitsch     PetscCall(TSComputeTransientVariable(ts,X,V));
3459566063dSJacob Faibussowitsch     PetscCall(VecAYPX(V,shift,V0));
346e3c11fc1SJed Brown   } else {                      /* shift*X + V0 */
3479566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(V,shift,X,V0));
348e3c11fc1SJed Brown   }
3491117961dSLisandro Dalcin 
350211a84d6SLisandro Dalcin   /* F = Function(t,X,V) */
3511117961dSLisandro Dalcin   ts->dm = dm;
3529566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts,t,X,V,F,PETSC_FALSE));
3531117961dSLisandro Dalcin   ts->dm = dmsave;
3541117961dSLisandro Dalcin 
3559566063dSJacob Faibussowitsch   PetscCall(TSBDF_RestoreVecs(ts,dm,&V,&V0));
356211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
357211a84d6SLisandro Dalcin }
358211a84d6SLisandro Dalcin 
3591117961dSLisandro Dalcin static PetscErrorCode SNESTSFormJacobian_BDF(SNES snes,Vec X,Mat J,Mat P,TS ts)
360211a84d6SLisandro Dalcin {
361211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
3621117961dSLisandro Dalcin   DM             dm, dmsave = ts->dm;
363211a84d6SLisandro Dalcin   PetscReal      t = bdf->time[0];
3641117961dSLisandro Dalcin   PetscReal      shift = bdf->shift;
3651117961dSLisandro Dalcin   Vec            V,V0;
366211a84d6SLisandro Dalcin 
367211a84d6SLisandro Dalcin   PetscFunctionBegin;
3689566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
3699566063dSJacob Faibussowitsch   PetscCall(TSBDF_GetVecs(ts,dm,&V,&V0));
3701117961dSLisandro Dalcin 
371211a84d6SLisandro Dalcin   /* J,P = Jacobian(t,X,V) */
3721117961dSLisandro Dalcin   ts->dm = dm;
3739566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts,t,X,V,shift,J,P,PETSC_FALSE));
3741117961dSLisandro Dalcin   ts->dm = dmsave;
3751117961dSLisandro Dalcin 
3769566063dSJacob Faibussowitsch   PetscCall(TSBDF_RestoreVecs(ts,dm,&V,&V0));
377211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
378211a84d6SLisandro Dalcin }
379211a84d6SLisandro Dalcin 
380211a84d6SLisandro Dalcin static PetscErrorCode TSReset_BDF(TS ts)
381211a84d6SLisandro Dalcin {
382211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
383211a84d6SLisandro Dalcin   size_t         i,n = sizeof(bdf->work)/sizeof(Vec);
384211a84d6SLisandro Dalcin 
385211a84d6SLisandro Dalcin   PetscFunctionBegin;
3861117961dSLisandro Dalcin   bdf->k = bdf->n = 0;
387e3c11fc1SJed Brown   for (i=0; i<n; i++) {
3889566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bdf->work[i]));
3899566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bdf->tvwork[i]));
390e3c11fc1SJed Brown   }
3919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bdf->vec_dot));
3929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bdf->vec_wrk));
3939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bdf->vec_lte));
3949566063dSJacob Faibussowitsch   if (ts->dm) PetscCall(DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSBDF,DMRestrictHook_TSBDF,ts));
395211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
396211a84d6SLisandro Dalcin }
397211a84d6SLisandro Dalcin 
398211a84d6SLisandro Dalcin static PetscErrorCode TSDestroy_BDF(TS ts)
399211a84d6SLisandro Dalcin {
400211a84d6SLisandro Dalcin   PetscFunctionBegin;
4019566063dSJacob Faibussowitsch   PetscCall(TSReset_BDF(ts));
4029566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
4039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSBDFSetOrder_C",NULL));
4049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSBDFGetOrder_C",NULL));
405211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
406211a84d6SLisandro Dalcin }
407211a84d6SLisandro Dalcin 
408211a84d6SLisandro Dalcin static PetscErrorCode TSSetUp_BDF(TS ts)
409211a84d6SLisandro Dalcin {
410211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
411211a84d6SLisandro Dalcin   size_t         i,n = sizeof(bdf->work)/sizeof(Vec);
4122ffb9264SLisandro Dalcin   PetscReal      low,high,two = 2;
413211a84d6SLisandro Dalcin 
414211a84d6SLisandro Dalcin   PetscFunctionBegin;
4159566063dSJacob Faibussowitsch   PetscCall(TSHasTransientVariable(ts,&bdf->transientvar));
416211a84d6SLisandro Dalcin   bdf->k = bdf->n = 0;
417e3c11fc1SJed Brown   for (i=0; i<n; i++) {
4189566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol,&bdf->work[i]));
419e3c11fc1SJed Brown     if (i && bdf->transientvar) {
4209566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(ts->vec_sol,&bdf->tvwork[i]));
421e3c11fc1SJed Brown     }
422e3c11fc1SJed Brown   }
4239566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol,&bdf->vec_dot));
4249566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol,&bdf->vec_wrk));
4259566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(ts->vec_sol,&bdf->vec_lte));
4269566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&ts->dm));
4279566063dSJacob Faibussowitsch   PetscCall(DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSBDF,DMRestrictHook_TSBDF,ts));
428211a84d6SLisandro Dalcin 
4299566063dSJacob Faibussowitsch   PetscCall(TSGetAdapt(ts,&ts->adapt));
4309566063dSJacob Faibussowitsch   PetscCall(TSAdaptCandidatesClear(ts->adapt));
4319566063dSJacob Faibussowitsch   PetscCall(TSAdaptGetClip(ts->adapt,&low,&high));
4329566063dSJacob Faibussowitsch   PetscCall(TSAdaptSetClip(ts->adapt,low,PetscMin(high,two)));
433211a84d6SLisandro Dalcin 
4349566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts,&ts->snes));
435211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
436211a84d6SLisandro Dalcin }
437211a84d6SLisandro Dalcin 
438211a84d6SLisandro Dalcin static PetscErrorCode TSSetFromOptions_BDF(PetscOptionItems *PetscOptionsObject,TS ts)
439211a84d6SLisandro Dalcin {
440211a84d6SLisandro Dalcin   PetscFunctionBegin;
4419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHead(PetscOptionsObject,"BDF ODE solver options"));
442211a84d6SLisandro Dalcin   {
443211a84d6SLisandro Dalcin     PetscBool flg;
444e5b8ffdfSLisandro Dalcin     PetscInt  order;
4459566063dSJacob Faibussowitsch     PetscCall(TSBDFGetOrder(ts,&order));
4469566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-ts_bdf_order","Order of the BDF method","TSBDFSetOrder",order,&order,&flg));
4479566063dSJacob Faibussowitsch     if (flg) PetscCall(TSBDFSetOrder(ts,order));
448211a84d6SLisandro Dalcin   }
4499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsTail());
450211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
451211a84d6SLisandro Dalcin }
452211a84d6SLisandro Dalcin 
453211a84d6SLisandro Dalcin static PetscErrorCode TSView_BDF(TS ts,PetscViewer viewer)
454211a84d6SLisandro Dalcin {
455211a84d6SLisandro Dalcin   TS_BDF         *bdf = (TS_BDF*)ts->data;
456211a84d6SLisandro Dalcin   PetscBool      iascii;
457211a84d6SLisandro Dalcin 
458211a84d6SLisandro Dalcin   PetscFunctionBegin;
4599566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
460211a84d6SLisandro Dalcin   if (iascii) {
4619566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  Order=%D\n",bdf->order));
462211a84d6SLisandro Dalcin   }
463211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
464211a84d6SLisandro Dalcin }
465211a84d6SLisandro Dalcin 
466211a84d6SLisandro Dalcin /* ------------------------------------------------------------ */
467211a84d6SLisandro Dalcin 
468211a84d6SLisandro Dalcin static PetscErrorCode TSBDFSetOrder_BDF(TS ts,PetscInt order)
469211a84d6SLisandro Dalcin {
470211a84d6SLisandro Dalcin   TS_BDF *bdf = (TS_BDF*)ts->data;
471211a84d6SLisandro Dalcin 
472211a84d6SLisandro Dalcin   PetscFunctionBegin;
473211a84d6SLisandro Dalcin   if (order == bdf->order) PetscFunctionReturn(0);
4743c633725SBarry Smith   PetscCheck(order >= 1 && order <= 6,PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"BDF Order %D not implemented",order);
475211a84d6SLisandro Dalcin   bdf->order = order;
476211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
477211a84d6SLisandro Dalcin }
478211a84d6SLisandro Dalcin 
479211a84d6SLisandro Dalcin static PetscErrorCode TSBDFGetOrder_BDF(TS ts,PetscInt *order)
480211a84d6SLisandro Dalcin {
481211a84d6SLisandro Dalcin   TS_BDF *bdf = (TS_BDF*)ts->data;
482211a84d6SLisandro Dalcin 
483211a84d6SLisandro Dalcin   PetscFunctionBegin;
484211a84d6SLisandro Dalcin   *order = bdf->order;
485211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
486211a84d6SLisandro Dalcin }
487211a84d6SLisandro Dalcin 
488211a84d6SLisandro Dalcin /* ------------------------------------------------------------ */
489211a84d6SLisandro Dalcin 
490211a84d6SLisandro Dalcin /*MC
491211a84d6SLisandro Dalcin       TSBDF - DAE solver using BDF methods
492211a84d6SLisandro Dalcin 
493211a84d6SLisandro Dalcin   Level: beginner
494211a84d6SLisandro Dalcin 
495211a84d6SLisandro Dalcin .seealso:  TS, TSCreate(), TSSetType()
496211a84d6SLisandro Dalcin M*/
497211a84d6SLisandro Dalcin PETSC_EXTERN PetscErrorCode TSCreate_BDF(TS ts)
498211a84d6SLisandro Dalcin {
499211a84d6SLisandro Dalcin   TS_BDF         *bdf;
500211a84d6SLisandro Dalcin 
501211a84d6SLisandro Dalcin   PetscFunctionBegin;
502211a84d6SLisandro Dalcin   ts->ops->reset          = TSReset_BDF;
503211a84d6SLisandro Dalcin   ts->ops->destroy        = TSDestroy_BDF;
504211a84d6SLisandro Dalcin   ts->ops->view           = TSView_BDF;
505211a84d6SLisandro Dalcin   ts->ops->setup          = TSSetUp_BDF;
506211a84d6SLisandro Dalcin   ts->ops->setfromoptions = TSSetFromOptions_BDF;
507211a84d6SLisandro Dalcin   ts->ops->step           = TSStep_BDF;
508211a84d6SLisandro Dalcin   ts->ops->evaluatewlte   = TSEvaluateWLTE_BDF;
509211a84d6SLisandro Dalcin   ts->ops->rollback       = TSRollBack_BDF;
510211a84d6SLisandro Dalcin   ts->ops->interpolate    = TSInterpolate_BDF;
511211a84d6SLisandro Dalcin   ts->ops->snesfunction   = SNESTSFormFunction_BDF;
512211a84d6SLisandro Dalcin   ts->ops->snesjacobian   = SNESTSFormJacobian_BDF;
5132ffb9264SLisandro Dalcin   ts->default_adapt_type  = TSADAPTBASIC;
514211a84d6SLisandro Dalcin 
515efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
516efd4aadfSBarry Smith 
5179566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(ts,&bdf));
518211a84d6SLisandro Dalcin   ts->data = (void*)bdf;
519211a84d6SLisandro Dalcin 
520211a84d6SLisandro Dalcin   bdf->status = TS_STEP_COMPLETE;
521211a84d6SLisandro Dalcin 
5229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSBDFSetOrder_C",TSBDFSetOrder_BDF));
5239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSBDFGetOrder_C",TSBDFGetOrder_BDF));
5249566063dSJacob Faibussowitsch   PetscCall(TSBDFSetOrder(ts,2));
525211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
526211a84d6SLisandro Dalcin }
527211a84d6SLisandro Dalcin 
528211a84d6SLisandro Dalcin /* ------------------------------------------------------------ */
529211a84d6SLisandro Dalcin 
530211a84d6SLisandro Dalcin /*@
531211a84d6SLisandro Dalcin   TSBDFSetOrder - Set the order of the BDF method
532211a84d6SLisandro Dalcin 
533211a84d6SLisandro Dalcin   Logically Collective on TS
534211a84d6SLisandro Dalcin 
535d8d19677SJose E. Roman   Input Parameters:
536211a84d6SLisandro Dalcin +  ts - timestepping context
537211a84d6SLisandro Dalcin -  order - order of the method
538211a84d6SLisandro Dalcin 
539211a84d6SLisandro Dalcin   Options Database:
540147403d9SBarry Smith .  -ts_bdf_order <order> - select the order
541211a84d6SLisandro Dalcin 
542211a84d6SLisandro Dalcin   Level: intermediate
543211a84d6SLisandro Dalcin 
544211a84d6SLisandro Dalcin @*/
545211a84d6SLisandro Dalcin PetscErrorCode TSBDFSetOrder(TS ts,PetscInt order)
546211a84d6SLisandro Dalcin {
547211a84d6SLisandro Dalcin   PetscFunctionBegin;
548211a84d6SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
549211a84d6SLisandro Dalcin   PetscValidLogicalCollectiveInt(ts,order,2);
550*cac4c232SBarry Smith   PetscTryMethod(ts,"TSBDFSetOrder_C",(TS,PetscInt),(ts,order));
551211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
552211a84d6SLisandro Dalcin }
553211a84d6SLisandro Dalcin 
554211a84d6SLisandro Dalcin /*@
555211a84d6SLisandro Dalcin   TSBDFGetOrder - Get the order of the BDF method
556211a84d6SLisandro Dalcin 
557211a84d6SLisandro Dalcin   Not Collective
558211a84d6SLisandro Dalcin 
559211a84d6SLisandro Dalcin   Input Parameter:
560211a84d6SLisandro Dalcin .  ts - timestepping context
561211a84d6SLisandro Dalcin 
562211a84d6SLisandro Dalcin   Output Parameter:
563211a84d6SLisandro Dalcin .  order - order of the method
564211a84d6SLisandro Dalcin 
565211a84d6SLisandro Dalcin   Level: intermediate
566211a84d6SLisandro Dalcin 
567211a84d6SLisandro Dalcin @*/
568211a84d6SLisandro Dalcin PetscErrorCode TSBDFGetOrder(TS ts,PetscInt *order)
569211a84d6SLisandro Dalcin {
570211a84d6SLisandro Dalcin   PetscFunctionBegin;
571211a84d6SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
572211a84d6SLisandro Dalcin   PetscValidIntPointer(order,2);
573*cac4c232SBarry Smith   PetscUseMethod(ts,"TSBDFGetOrder_C",(TS,PetscInt*),(ts,order));
574211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
575211a84d6SLisandro Dalcin }
576