xref: /petsc/src/ts/impls/bdf/bdf.c (revision bcf0153e883cfed9568ef4557dcc209048fb58f7)
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;
89371c9d4SSatish Balay static const char citation[] = "@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];
19e3c11fc1SJed Brown   Vec       tvwork[6 + 2];
20211a84d6SLisandro Dalcin   PetscReal shift;
21e3c11fc1SJed Brown   Vec       vec_dot; /* Xdot when !transientvar, else Cdot where C(X) is the transient variable. */
221117961dSLisandro Dalcin   Vec       vec_wrk;
23211a84d6SLisandro Dalcin   Vec       vec_lte;
24211a84d6SLisandro Dalcin 
25e3c11fc1SJed Brown   PetscBool    transientvar;
26211a84d6SLisandro Dalcin   PetscInt     order;
27211a84d6SLisandro Dalcin   TSStepStatus status;
28211a84d6SLisandro Dalcin } TS_BDF;
29211a84d6SLisandro Dalcin 
30e3c11fc1SJed Brown /* Compute Lagrange polynomials on T[:n] evaluated at t.
31e3c11fc1SJed Brown  * If one has data (T[i], Y[i]), then the interpolation/extrapolation is f(t) = \sum_i L[i]*Y[i].
32e3c11fc1SJed Brown  */
33d71ae5a4SJacob Faibussowitsch static inline void LagrangeBasisVals(PetscInt n, PetscReal t, const PetscReal T[], PetscScalar L[])
34d71ae5a4SJacob Faibussowitsch {
35211a84d6SLisandro Dalcin   PetscInt k, j;
36211a84d6SLisandro Dalcin   for (k = 0; k < n; k++)
37211a84d6SLisandro Dalcin     for (L[k] = 1, j = 0; j < n; j++)
389371c9d4SSatish Balay       if (j != k) L[k] *= (t - T[j]) / (T[k] - T[j]);
39211a84d6SLisandro Dalcin }
40211a84d6SLisandro Dalcin 
41d71ae5a4SJacob Faibussowitsch static inline void LagrangeBasisDers(PetscInt n, PetscReal t, const PetscReal T[], PetscScalar dL[])
42d71ae5a4SJacob Faibussowitsch {
43211a84d6SLisandro Dalcin   PetscInt k, j, i;
44211a84d6SLisandro Dalcin   for (k = 0; k < n; k++)
45211a84d6SLisandro Dalcin     for (dL[k] = 0, j = 0; j < n; j++)
46211a84d6SLisandro Dalcin       if (j != k) {
47211a84d6SLisandro Dalcin         PetscReal L = 1 / (T[k] - T[j]);
48211a84d6SLisandro Dalcin         for (i = 0; i < n; i++)
499371c9d4SSatish Balay           if (i != j && i != k) L *= (t - T[i]) / (T[k] - T[i]);
50211a84d6SLisandro Dalcin         dL[k] += L;
51211a84d6SLisandro Dalcin       }
52211a84d6SLisandro Dalcin }
53211a84d6SLisandro Dalcin 
54d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_GetVecs(TS ts, DM dm, Vec *Xdot, Vec *Ydot)
55d71ae5a4SJacob Faibussowitsch {
561117961dSLisandro Dalcin   TS_BDF *bdf = (TS_BDF *)ts->data;
571117961dSLisandro Dalcin 
581117961dSLisandro Dalcin   PetscFunctionBegin;
591117961dSLisandro Dalcin   if (dm && dm != ts->dm) {
609566063dSJacob Faibussowitsch     PetscCall(DMGetNamedGlobalVector(dm, "TSBDF_Vec_Xdot", Xdot));
619566063dSJacob Faibussowitsch     PetscCall(DMGetNamedGlobalVector(dm, "TSBDF_Vec_Ydot", Ydot));
621117961dSLisandro Dalcin   } else {
631117961dSLisandro Dalcin     *Xdot = bdf->vec_dot;
641117961dSLisandro Dalcin     *Ydot = bdf->vec_wrk;
651117961dSLisandro Dalcin   }
661117961dSLisandro Dalcin   PetscFunctionReturn(0);
671117961dSLisandro Dalcin }
681117961dSLisandro Dalcin 
69d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_RestoreVecs(TS ts, DM dm, Vec *Xdot, Vec *Ydot)
70d71ae5a4SJacob Faibussowitsch {
711117961dSLisandro Dalcin   TS_BDF *bdf = (TS_BDF *)ts->data;
721117961dSLisandro Dalcin 
731117961dSLisandro Dalcin   PetscFunctionBegin;
741117961dSLisandro Dalcin   if (dm && dm != ts->dm) {
759566063dSJacob Faibussowitsch     PetscCall(DMRestoreNamedGlobalVector(dm, "TSBDF_Vec_Xdot", Xdot));
769566063dSJacob Faibussowitsch     PetscCall(DMRestoreNamedGlobalVector(dm, "TSBDF_Vec_Ydot", Ydot));
771117961dSLisandro Dalcin   } else {
783c633725SBarry Smith     PetscCheck(*Xdot == bdf->vec_dot, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_INCOMP, "Vec does not match the cache");
793c633725SBarry Smith     PetscCheck(*Ydot == bdf->vec_wrk, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_INCOMP, "Vec does not match the cache");
801117961dSLisandro Dalcin     *Xdot = NULL;
811117961dSLisandro Dalcin     *Ydot = NULL;
821117961dSLisandro Dalcin   }
831117961dSLisandro Dalcin   PetscFunctionReturn(0);
841117961dSLisandro Dalcin }
851117961dSLisandro Dalcin 
86d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCoarsenHook_TSBDF(DM fine, DM coarse, void *ctx)
87d71ae5a4SJacob Faibussowitsch {
881117961dSLisandro Dalcin   PetscFunctionBegin;
891117961dSLisandro Dalcin   PetscFunctionReturn(0);
901117961dSLisandro Dalcin }
911117961dSLisandro Dalcin 
92d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMRestrictHook_TSBDF(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx)
93d71ae5a4SJacob Faibussowitsch {
941117961dSLisandro Dalcin   TS  ts = (TS)ctx;
951117961dSLisandro Dalcin   Vec Ydot, Ydot_c;
961117961dSLisandro Dalcin   Vec Xdot, Xdot_c;
971117961dSLisandro Dalcin 
981117961dSLisandro Dalcin   PetscFunctionBegin;
999566063dSJacob Faibussowitsch   PetscCall(TSBDF_GetVecs(ts, fine, &Xdot, &Ydot));
1009566063dSJacob Faibussowitsch   PetscCall(TSBDF_GetVecs(ts, coarse, &Xdot_c, &Ydot_c));
1011117961dSLisandro Dalcin 
1029566063dSJacob Faibussowitsch   PetscCall(MatRestrict(restrct, Ydot, Ydot_c));
1039566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMult(Ydot_c, rscale, Ydot_c));
1041117961dSLisandro Dalcin 
1059566063dSJacob Faibussowitsch   PetscCall(TSBDF_RestoreVecs(ts, fine, &Xdot, &Ydot));
1069566063dSJacob Faibussowitsch   PetscCall(TSBDF_RestoreVecs(ts, coarse, &Xdot_c, &Ydot_c));
1071117961dSLisandro Dalcin   PetscFunctionReturn(0);
1081117961dSLisandro Dalcin }
1091117961dSLisandro Dalcin 
110d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_Advance(TS ts, PetscReal t, Vec X)
111d71ae5a4SJacob Faibussowitsch {
112211a84d6SLisandro Dalcin   TS_BDF  *bdf = (TS_BDF *)ts->data;
113211a84d6SLisandro Dalcin   PetscInt i, n = (PetscInt)(sizeof(bdf->work) / sizeof(Vec));
114e3c11fc1SJed Brown   Vec      tail = bdf->work[n - 1], tvtail = bdf->tvwork[n - 1];
115211a84d6SLisandro Dalcin 
116211a84d6SLisandro Dalcin   PetscFunctionBegin;
117211a84d6SLisandro Dalcin   for (i = n - 1; i >= 2; i--) {
118211a84d6SLisandro Dalcin     bdf->time[i]   = bdf->time[i - 1];
119211a84d6SLisandro Dalcin     bdf->work[i]   = bdf->work[i - 1];
120e3c11fc1SJed Brown     bdf->tvwork[i] = bdf->tvwork[i - 1];
121211a84d6SLisandro Dalcin   }
122211a84d6SLisandro Dalcin   bdf->n         = PetscMin(bdf->n + 1, n - 1);
123211a84d6SLisandro Dalcin   bdf->time[1]   = t;
124211a84d6SLisandro Dalcin   bdf->work[1]   = tail;
125e3c11fc1SJed Brown   bdf->tvwork[1] = tvtail;
1269566063dSJacob Faibussowitsch   PetscCall(VecCopy(X, tail));
1279566063dSJacob Faibussowitsch   PetscCall(TSComputeTransientVariable(ts, tail, tvtail));
128211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
129211a84d6SLisandro Dalcin }
130211a84d6SLisandro Dalcin 
131d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_VecLTE(TS ts, PetscInt order, Vec lte)
132d71ae5a4SJacob Faibussowitsch {
133211a84d6SLisandro Dalcin   TS_BDF     *bdf = (TS_BDF *)ts->data;
134211a84d6SLisandro Dalcin   PetscInt    i, n = order + 1;
135211a84d6SLisandro Dalcin   PetscReal  *time = bdf->time;
136211a84d6SLisandro Dalcin   Vec        *vecs = bdf->work;
137211a84d6SLisandro Dalcin   PetscScalar a[8], b[8], alpha[8];
138211a84d6SLisandro Dalcin 
139211a84d6SLisandro Dalcin   PetscFunctionBegin;
1409371c9d4SSatish Balay   LagrangeBasisDers(n + 0, time[0], time, a);
1419371c9d4SSatish Balay   a[n] = 0;
142211a84d6SLisandro Dalcin   LagrangeBasisDers(n + 1, time[0], time, b);
143211a84d6SLisandro Dalcin   for (i = 0; i < n + 1; i++) alpha[i] = (a[i] - b[i]) / a[0];
1449566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(lte));
1459566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(lte, n + 1, alpha, vecs));
146211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
147211a84d6SLisandro Dalcin }
148211a84d6SLisandro Dalcin 
149d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_Extrapolate(TS ts, PetscInt order, PetscReal t, Vec X)
150d71ae5a4SJacob Faibussowitsch {
151211a84d6SLisandro Dalcin   TS_BDF     *bdf  = (TS_BDF *)ts->data;
152211a84d6SLisandro Dalcin   PetscInt    n    = order + 1;
153211a84d6SLisandro Dalcin   PetscReal  *time = bdf->time + 1;
154211a84d6SLisandro Dalcin   Vec        *vecs = bdf->work + 1;
155211a84d6SLisandro Dalcin   PetscScalar alpha[7];
156211a84d6SLisandro Dalcin 
157211a84d6SLisandro Dalcin   PetscFunctionBegin;
158211a84d6SLisandro Dalcin   n = PetscMin(n, bdf->n);
159211a84d6SLisandro Dalcin   LagrangeBasisVals(n, t, time, alpha);
1609566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(X));
1619566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(X, n, alpha, vecs));
162211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
163211a84d6SLisandro Dalcin }
164211a84d6SLisandro Dalcin 
165d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_Interpolate(TS ts, PetscInt order, PetscReal t, Vec X)
166d71ae5a4SJacob Faibussowitsch {
167211a84d6SLisandro Dalcin   TS_BDF     *bdf  = (TS_BDF *)ts->data;
168211a84d6SLisandro Dalcin   PetscInt    n    = order + 1;
169211a84d6SLisandro Dalcin   PetscReal  *time = bdf->time;
170211a84d6SLisandro Dalcin   Vec        *vecs = bdf->work;
171211a84d6SLisandro Dalcin   PetscScalar alpha[7];
172211a84d6SLisandro Dalcin 
173211a84d6SLisandro Dalcin   PetscFunctionBegin;
174211a84d6SLisandro Dalcin   LagrangeBasisVals(n, t, time, alpha);
1759566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(X));
1769566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(X, n, alpha, vecs));
177211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
178211a84d6SLisandro Dalcin }
179211a84d6SLisandro Dalcin 
180e3c11fc1SJed Brown /* Compute the affine term V0 such that Xdot = shift*X + V0.
181e3c11fc1SJed Brown  *
182e3c11fc1SJed Brown  * When using transient variables, we're computing Cdot = shift*C(X) + V0, and thus choose a linear combination of tvwork.
183e3c11fc1SJed Brown  */
184d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_PreSolve(TS ts)
185d71ae5a4SJacob Faibussowitsch {
1861117961dSLisandro Dalcin   TS_BDF     *bdf = (TS_BDF *)ts->data;
1871117961dSLisandro Dalcin   PetscInt    i, n = PetscMax(bdf->k, 1) + 1;
1881117961dSLisandro Dalcin   Vec         V, V0;
1891117961dSLisandro Dalcin   Vec         vecs[7];
1901117961dSLisandro Dalcin   PetscScalar alpha[7];
1911117961dSLisandro Dalcin 
1921117961dSLisandro Dalcin   PetscFunctionBegin;
1939566063dSJacob Faibussowitsch   PetscCall(TSBDF_GetVecs(ts, NULL, &V, &V0));
1941117961dSLisandro Dalcin   LagrangeBasisDers(n, bdf->time[0], bdf->time, alpha);
195ad540459SPierre Jolivet   for (i = 1; i < n; i++) vecs[i] = bdf->transientvar ? bdf->tvwork[i] : bdf->work[i];
1969566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(V0));
1979566063dSJacob Faibussowitsch   PetscCall(VecMAXPY(V0, n - 1, alpha + 1, vecs + 1));
1981117961dSLisandro Dalcin   bdf->shift = PetscRealPart(alpha[0]);
1999566063dSJacob Faibussowitsch   PetscCall(TSBDF_RestoreVecs(ts, NULL, &V, &V0));
2001117961dSLisandro Dalcin   PetscFunctionReturn(0);
2011117961dSLisandro Dalcin }
2021117961dSLisandro Dalcin 
203d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_SNESSolve(TS ts, Vec b, Vec x)
204d71ae5a4SJacob Faibussowitsch {
205211a84d6SLisandro Dalcin   PetscInt nits, lits;
206211a84d6SLisandro Dalcin 
207211a84d6SLisandro Dalcin   PetscFunctionBegin;
2089566063dSJacob Faibussowitsch   PetscCall(TSBDF_PreSolve(ts));
2099566063dSJacob Faibussowitsch   PetscCall(SNESSolve(ts->snes, b, x));
2109566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(ts->snes, &nits));
2119566063dSJacob Faibussowitsch   PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits));
2129371c9d4SSatish Balay   ts->snes_its += nits;
2139371c9d4SSatish Balay   ts->ksp_its += lits;
214211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
215211a84d6SLisandro Dalcin }
216211a84d6SLisandro Dalcin 
217d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDF_Restart(TS ts, PetscBool *accept)
218d71ae5a4SJacob Faibussowitsch {
219211a84d6SLisandro Dalcin   TS_BDF *bdf = (TS_BDF *)ts->data;
220211a84d6SLisandro Dalcin 
221211a84d6SLisandro Dalcin   PetscFunctionBegin;
2229371c9d4SSatish Balay   bdf->k = 1;
2239371c9d4SSatish Balay   bdf->n = 0;
2249566063dSJacob Faibussowitsch   PetscCall(TSBDF_Advance(ts, ts->ptime, ts->vec_sol));
225211a84d6SLisandro Dalcin 
226211a84d6SLisandro Dalcin   bdf->time[0] = ts->ptime + ts->time_step / 2;
2279566063dSJacob Faibussowitsch   PetscCall(VecCopy(bdf->work[1], bdf->work[0]));
2289566063dSJacob Faibussowitsch   PetscCall(TSPreStage(ts, bdf->time[0]));
2299566063dSJacob Faibussowitsch   PetscCall(TSBDF_SNESSolve(ts, NULL, bdf->work[0]));
2309566063dSJacob Faibussowitsch   PetscCall(TSPostStage(ts, bdf->time[0], 0, &bdf->work[0]));
2319566063dSJacob Faibussowitsch   PetscCall(TSAdaptCheckStage(ts->adapt, ts, bdf->time[0], bdf->work[0], accept));
232211a84d6SLisandro Dalcin   if (!*accept) PetscFunctionReturn(0);
233211a84d6SLisandro Dalcin 
2349371c9d4SSatish Balay   bdf->k = PetscMin(2, bdf->order);
2359371c9d4SSatish Balay   bdf->n++;
2369566063dSJacob Faibussowitsch   PetscCall(VecCopy(bdf->work[0], bdf->work[2]));
237211a84d6SLisandro Dalcin   bdf->time[2] = bdf->time[0];
2389566063dSJacob Faibussowitsch   PetscCall(TSComputeTransientVariable(ts, bdf->work[2], bdf->tvwork[2]));
239211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
240211a84d6SLisandro Dalcin }
241211a84d6SLisandro Dalcin 
242211a84d6SLisandro Dalcin static const char *const BDF_SchemeName[] = {"", "1", "2", "3", "4", "5", "6"};
243211a84d6SLisandro Dalcin 
244d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSStep_BDF(TS ts)
245d71ae5a4SJacob Faibussowitsch {
246211a84d6SLisandro Dalcin   TS_BDF   *bdf        = (TS_BDF *)ts->data;
247211a84d6SLisandro Dalcin   PetscInt  rejections = 0;
248211a84d6SLisandro Dalcin   PetscBool stageok, accept = PETSC_TRUE;
249211a84d6SLisandro Dalcin   PetscReal next_time_step = ts->time_step;
250211a84d6SLisandro Dalcin 
251211a84d6SLisandro Dalcin   PetscFunctionBegin;
2529566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(citation, &cited));
253211a84d6SLisandro Dalcin 
254211a84d6SLisandro Dalcin   if (!ts->steprollback && !ts->steprestart) {
255211a84d6SLisandro Dalcin     bdf->k = PetscMin(bdf->k + 1, bdf->order);
2569566063dSJacob Faibussowitsch     PetscCall(TSBDF_Advance(ts, ts->ptime, ts->vec_sol));
257211a84d6SLisandro Dalcin   }
258211a84d6SLisandro Dalcin 
259211a84d6SLisandro Dalcin   bdf->status = TS_STEP_INCOMPLETE;
260211a84d6SLisandro Dalcin   while (!ts->reason && bdf->status != TS_STEP_COMPLETE) {
261211a84d6SLisandro Dalcin     if (ts->steprestart) {
2629566063dSJacob Faibussowitsch       PetscCall(TSBDF_Restart(ts, &stageok));
263211a84d6SLisandro Dalcin       if (!stageok) goto reject_step;
264211a84d6SLisandro Dalcin     }
265211a84d6SLisandro Dalcin 
266211a84d6SLisandro Dalcin     bdf->time[0] = ts->ptime + ts->time_step;
2679566063dSJacob Faibussowitsch     PetscCall(TSBDF_Extrapolate(ts, bdf->k - (accept ? 0 : 1), bdf->time[0], bdf->work[0]));
2689566063dSJacob Faibussowitsch     PetscCall(TSPreStage(ts, bdf->time[0]));
2699566063dSJacob Faibussowitsch     PetscCall(TSBDF_SNESSolve(ts, NULL, bdf->work[0]));
2709566063dSJacob Faibussowitsch     PetscCall(TSPostStage(ts, bdf->time[0], 0, &bdf->work[0]));
2719566063dSJacob Faibussowitsch     PetscCall(TSAdaptCheckStage(ts->adapt, ts, bdf->time[0], bdf->work[0], &stageok));
272211a84d6SLisandro Dalcin     if (!stageok) goto reject_step;
273211a84d6SLisandro Dalcin 
274211a84d6SLisandro Dalcin     bdf->status = TS_STEP_PENDING;
2759566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidatesClear(ts->adapt));
2769566063dSJacob Faibussowitsch     PetscCall(TSAdaptCandidateAdd(ts->adapt, BDF_SchemeName[bdf->k], bdf->k, 1, 1.0, 1.0, PETSC_TRUE));
2779566063dSJacob Faibussowitsch     PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
278211a84d6SLisandro Dalcin     bdf->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
2799371c9d4SSatish Balay     if (!accept) {
2809371c9d4SSatish Balay       ts->time_step = next_time_step;
2819371c9d4SSatish Balay       goto reject_step;
2829371c9d4SSatish Balay     }
283211a84d6SLisandro Dalcin 
2849566063dSJacob Faibussowitsch     PetscCall(VecCopy(bdf->work[0], ts->vec_sol));
285211a84d6SLisandro Dalcin     ts->ptime += ts->time_step;
286211a84d6SLisandro Dalcin     ts->time_step = next_time_step;
287211a84d6SLisandro Dalcin     break;
288211a84d6SLisandro Dalcin 
289211a84d6SLisandro Dalcin   reject_step:
2909371c9d4SSatish Balay     ts->reject++;
2919371c9d4SSatish Balay     accept = PETSC_FALSE;
292211a84d6SLisandro Dalcin     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
29363a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
294211a84d6SLisandro Dalcin       ts->reason = TS_DIVERGED_STEP_REJECTED;
295211a84d6SLisandro Dalcin     }
296211a84d6SLisandro Dalcin   }
297211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
298211a84d6SLisandro Dalcin }
299211a84d6SLisandro Dalcin 
300d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSInterpolate_BDF(TS ts, PetscReal t, Vec X)
301d71ae5a4SJacob Faibussowitsch {
302211a84d6SLisandro Dalcin   TS_BDF *bdf = (TS_BDF *)ts->data;
303211a84d6SLisandro Dalcin 
304211a84d6SLisandro Dalcin   PetscFunctionBegin;
3059566063dSJacob Faibussowitsch   PetscCall(TSBDF_Interpolate(ts, bdf->k, t, X));
306211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
307211a84d6SLisandro Dalcin }
308211a84d6SLisandro Dalcin 
309d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSEvaluateWLTE_BDF(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte)
310d71ae5a4SJacob Faibussowitsch {
311211a84d6SLisandro Dalcin   TS_BDF   *bdf = (TS_BDF *)ts->data;
312211a84d6SLisandro Dalcin   PetscInt  k   = bdf->k;
3137453f775SEmil Constantinescu   PetscReal wltea, wlter;
314211a84d6SLisandro Dalcin   Vec       X = bdf->work[0], Y = bdf->vec_lte;
315211a84d6SLisandro Dalcin 
316211a84d6SLisandro Dalcin   PetscFunctionBegin;
317211a84d6SLisandro Dalcin   k = PetscMin(k, bdf->n - 1);
3189566063dSJacob Faibussowitsch   PetscCall(TSBDF_VecLTE(ts, k, Y));
3199566063dSJacob Faibussowitsch   PetscCall(VecAXPY(Y, 1, X));
3209566063dSJacob Faibussowitsch   PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, wlte, &wltea, &wlter));
321211a84d6SLisandro Dalcin   if (order) *order = k + 1;
322211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
323211a84d6SLisandro Dalcin }
324211a84d6SLisandro Dalcin 
325d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSRollBack_BDF(TS ts)
326d71ae5a4SJacob Faibussowitsch {
327211a84d6SLisandro Dalcin   TS_BDF *bdf = (TS_BDF *)ts->data;
328211a84d6SLisandro Dalcin 
329211a84d6SLisandro Dalcin   PetscFunctionBegin;
3309566063dSJacob Faibussowitsch   PetscCall(VecCopy(bdf->work[1], ts->vec_sol));
331211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
332211a84d6SLisandro Dalcin }
333211a84d6SLisandro Dalcin 
334d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormFunction_BDF(SNES snes, Vec X, Vec F, TS ts)
335d71ae5a4SJacob Faibussowitsch {
336211a84d6SLisandro Dalcin   TS_BDF   *bdf = (TS_BDF *)ts->data;
3371117961dSLisandro Dalcin   DM        dm, dmsave = ts->dm;
338211a84d6SLisandro Dalcin   PetscReal t     = bdf->time[0];
3391117961dSLisandro Dalcin   PetscReal shift = bdf->shift;
3401117961dSLisandro Dalcin   Vec       V, V0;
341211a84d6SLisandro Dalcin 
342211a84d6SLisandro Dalcin   PetscFunctionBegin;
3439566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
3449566063dSJacob Faibussowitsch   PetscCall(TSBDF_GetVecs(ts, dm, &V, &V0));
345e3c11fc1SJed Brown   if (bdf->transientvar) { /* shift*C(X) + V0 */
3469566063dSJacob Faibussowitsch     PetscCall(TSComputeTransientVariable(ts, X, V));
3479566063dSJacob Faibussowitsch     PetscCall(VecAYPX(V, shift, V0));
348e3c11fc1SJed Brown   } else { /* shift*X + V0 */
3499566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(V, shift, X, V0));
350e3c11fc1SJed Brown   }
3511117961dSLisandro Dalcin 
352211a84d6SLisandro Dalcin   /* F = Function(t,X,V) */
3531117961dSLisandro Dalcin   ts->dm = dm;
3549566063dSJacob Faibussowitsch   PetscCall(TSComputeIFunction(ts, t, X, V, F, PETSC_FALSE));
3551117961dSLisandro Dalcin   ts->dm = dmsave;
3561117961dSLisandro Dalcin 
3579566063dSJacob Faibussowitsch   PetscCall(TSBDF_RestoreVecs(ts, dm, &V, &V0));
358211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
359211a84d6SLisandro Dalcin }
360211a84d6SLisandro Dalcin 
361d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTSFormJacobian_BDF(SNES snes, Vec X, Mat J, Mat P, TS ts)
362d71ae5a4SJacob Faibussowitsch {
363211a84d6SLisandro Dalcin   TS_BDF   *bdf = (TS_BDF *)ts->data;
3641117961dSLisandro Dalcin   DM        dm, dmsave = ts->dm;
365211a84d6SLisandro Dalcin   PetscReal t     = bdf->time[0];
3661117961dSLisandro Dalcin   PetscReal shift = bdf->shift;
3671117961dSLisandro Dalcin   Vec       V, V0;
368211a84d6SLisandro Dalcin 
369211a84d6SLisandro Dalcin   PetscFunctionBegin;
3709566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
3719566063dSJacob Faibussowitsch   PetscCall(TSBDF_GetVecs(ts, dm, &V, &V0));
3721117961dSLisandro Dalcin 
373211a84d6SLisandro Dalcin   /* J,P = Jacobian(t,X,V) */
3741117961dSLisandro Dalcin   ts->dm = dm;
3759566063dSJacob Faibussowitsch   PetscCall(TSComputeIJacobian(ts, t, X, V, shift, J, P, PETSC_FALSE));
3761117961dSLisandro Dalcin   ts->dm = dmsave;
3771117961dSLisandro Dalcin 
3789566063dSJacob Faibussowitsch   PetscCall(TSBDF_RestoreVecs(ts, dm, &V, &V0));
379211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
380211a84d6SLisandro Dalcin }
381211a84d6SLisandro Dalcin 
382d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSReset_BDF(TS ts)
383d71ae5a4SJacob Faibussowitsch {
384211a84d6SLisandro Dalcin   TS_BDF *bdf = (TS_BDF *)ts->data;
385211a84d6SLisandro Dalcin   size_t  i, n = sizeof(bdf->work) / sizeof(Vec);
386211a84d6SLisandro Dalcin 
387211a84d6SLisandro Dalcin   PetscFunctionBegin;
3881117961dSLisandro Dalcin   bdf->k = bdf->n = 0;
389e3c11fc1SJed Brown   for (i = 0; i < n; i++) {
3909566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bdf->work[i]));
3919566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bdf->tvwork[i]));
392e3c11fc1SJed Brown   }
3939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bdf->vec_dot));
3949566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bdf->vec_wrk));
3959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&bdf->vec_lte));
3969566063dSJacob Faibussowitsch   if (ts->dm) PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSBDF, DMRestrictHook_TSBDF, ts));
397211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
398211a84d6SLisandro Dalcin }
399211a84d6SLisandro Dalcin 
400d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSDestroy_BDF(TS ts)
401d71ae5a4SJacob Faibussowitsch {
402211a84d6SLisandro Dalcin   PetscFunctionBegin;
4039566063dSJacob Faibussowitsch   PetscCall(TSReset_BDF(ts));
4049566063dSJacob Faibussowitsch   PetscCall(PetscFree(ts->data));
4059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBDFSetOrder_C", NULL));
4069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBDFGetOrder_C", NULL));
407211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
408211a84d6SLisandro Dalcin }
409211a84d6SLisandro Dalcin 
410d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetUp_BDF(TS ts)
411d71ae5a4SJacob Faibussowitsch {
412211a84d6SLisandro Dalcin   TS_BDF   *bdf = (TS_BDF *)ts->data;
413211a84d6SLisandro Dalcin   size_t    i, n = sizeof(bdf->work) / sizeof(Vec);
4142ffb9264SLisandro Dalcin   PetscReal low, high, two = 2;
415211a84d6SLisandro Dalcin 
416211a84d6SLisandro Dalcin   PetscFunctionBegin;
4179566063dSJacob Faibussowitsch   PetscCall(TSHasTransientVariable(ts, &bdf->transientvar));
418211a84d6SLisandro Dalcin   bdf->k = bdf->n = 0;
419e3c11fc1SJed Brown   for (i = 0; i < n; i++) {
4209566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(ts->vec_sol, &bdf->work[i]));
42148a46eb9SPierre Jolivet     if (i && bdf->transientvar) PetscCall(VecDuplicate(ts->vec_sol, &bdf->tvwork[i]));
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 
438d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSSetFromOptions_BDF(TS ts, PetscOptionItems *PetscOptionsObject)
439d71ae5a4SJacob Faibussowitsch {
440211a84d6SLisandro Dalcin   PetscFunctionBegin;
441d0609cedSBarry Smith   PetscOptionsHeadBegin(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   }
449d0609cedSBarry Smith   PetscOptionsHeadEnd();
450211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
451211a84d6SLisandro Dalcin }
452211a84d6SLisandro Dalcin 
453d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSView_BDF(TS ts, PetscViewer viewer)
454d71ae5a4SJacob Faibussowitsch {
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));
46048a46eb9SPierre Jolivet   if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, "  Order=%" PetscInt_FMT "\n", bdf->order));
461211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
462211a84d6SLisandro Dalcin }
463211a84d6SLisandro Dalcin 
464211a84d6SLisandro Dalcin /* ------------------------------------------------------------ */
465211a84d6SLisandro Dalcin 
466d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDFSetOrder_BDF(TS ts, PetscInt order)
467d71ae5a4SJacob Faibussowitsch {
468211a84d6SLisandro Dalcin   TS_BDF *bdf = (TS_BDF *)ts->data;
469211a84d6SLisandro Dalcin 
470211a84d6SLisandro Dalcin   PetscFunctionBegin;
471211a84d6SLisandro Dalcin   if (order == bdf->order) PetscFunctionReturn(0);
47263a3b9bcSJacob Faibussowitsch   PetscCheck(order >= 1 && order <= 6, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "BDF Order %" PetscInt_FMT " not implemented", order);
473211a84d6SLisandro Dalcin   bdf->order = order;
474211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
475211a84d6SLisandro Dalcin }
476211a84d6SLisandro Dalcin 
477d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSBDFGetOrder_BDF(TS ts, PetscInt *order)
478d71ae5a4SJacob Faibussowitsch {
479211a84d6SLisandro Dalcin   TS_BDF *bdf = (TS_BDF *)ts->data;
480211a84d6SLisandro Dalcin 
481211a84d6SLisandro Dalcin   PetscFunctionBegin;
482211a84d6SLisandro Dalcin   *order = bdf->order;
483211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
484211a84d6SLisandro Dalcin }
485211a84d6SLisandro Dalcin 
486211a84d6SLisandro Dalcin /* ------------------------------------------------------------ */
487211a84d6SLisandro Dalcin 
488211a84d6SLisandro Dalcin /*MC
489211a84d6SLisandro Dalcin       TSBDF - DAE solver using BDF methods
490211a84d6SLisandro Dalcin 
491211a84d6SLisandro Dalcin   Level: beginner
492211a84d6SLisandro Dalcin 
493*bcf0153eSBarry Smith .seealso: [](chapter_ts), `TS`, `TSCreate()`, `TSSetType()`, `TSType`
494211a84d6SLisandro Dalcin M*/
495d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSCreate_BDF(TS ts)
496d71ae5a4SJacob Faibussowitsch {
497211a84d6SLisandro Dalcin   TS_BDF *bdf;
498211a84d6SLisandro Dalcin 
499211a84d6SLisandro Dalcin   PetscFunctionBegin;
500211a84d6SLisandro Dalcin   ts->ops->reset          = TSReset_BDF;
501211a84d6SLisandro Dalcin   ts->ops->destroy        = TSDestroy_BDF;
502211a84d6SLisandro Dalcin   ts->ops->view           = TSView_BDF;
503211a84d6SLisandro Dalcin   ts->ops->setup          = TSSetUp_BDF;
504211a84d6SLisandro Dalcin   ts->ops->setfromoptions = TSSetFromOptions_BDF;
505211a84d6SLisandro Dalcin   ts->ops->step           = TSStep_BDF;
506211a84d6SLisandro Dalcin   ts->ops->evaluatewlte   = TSEvaluateWLTE_BDF;
507211a84d6SLisandro Dalcin   ts->ops->rollback       = TSRollBack_BDF;
508211a84d6SLisandro Dalcin   ts->ops->interpolate    = TSInterpolate_BDF;
509211a84d6SLisandro Dalcin   ts->ops->snesfunction   = SNESTSFormFunction_BDF;
510211a84d6SLisandro Dalcin   ts->ops->snesjacobian   = SNESTSFormJacobian_BDF;
5112ffb9264SLisandro Dalcin   ts->default_adapt_type  = TSADAPTBASIC;
512211a84d6SLisandro Dalcin 
513efd4aadfSBarry Smith   ts->usessnes = PETSC_TRUE;
514efd4aadfSBarry Smith 
5154dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&bdf));
516211a84d6SLisandro Dalcin   ts->data = (void *)bdf;
517211a84d6SLisandro Dalcin 
518211a84d6SLisandro Dalcin   bdf->status = TS_STEP_COMPLETE;
519211a84d6SLisandro Dalcin 
5209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBDFSetOrder_C", TSBDFSetOrder_BDF));
5219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSBDFGetOrder_C", TSBDFGetOrder_BDF));
5229566063dSJacob Faibussowitsch   PetscCall(TSBDFSetOrder(ts, 2));
523211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
524211a84d6SLisandro Dalcin }
525211a84d6SLisandro Dalcin 
526211a84d6SLisandro Dalcin /* ------------------------------------------------------------ */
527211a84d6SLisandro Dalcin 
528211a84d6SLisandro Dalcin /*@
529*bcf0153eSBarry Smith   TSBDFSetOrder - Set the order of the `TSBDF` method
530211a84d6SLisandro Dalcin 
531*bcf0153eSBarry Smith   Logically Collective on ts
532211a84d6SLisandro Dalcin 
533d8d19677SJose E. Roman   Input Parameters:
534211a84d6SLisandro Dalcin +  ts - timestepping context
535211a84d6SLisandro Dalcin -  order - order of the method
536211a84d6SLisandro Dalcin 
537*bcf0153eSBarry Smith   Options Database Key:
538147403d9SBarry Smith .  -ts_bdf_order <order> - select the order
539211a84d6SLisandro Dalcin 
540211a84d6SLisandro Dalcin   Level: intermediate
541211a84d6SLisandro Dalcin 
542*bcf0153eSBarry Smith .seealso: `TSBDFGetOrder()`, `TS`, `TSBDF`
543211a84d6SLisandro Dalcin @*/
544d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBDFSetOrder(TS ts, PetscInt order)
545d71ae5a4SJacob Faibussowitsch {
546211a84d6SLisandro Dalcin   PetscFunctionBegin;
547211a84d6SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
548211a84d6SLisandro Dalcin   PetscValidLogicalCollectiveInt(ts, order, 2);
549cac4c232SBarry Smith   PetscTryMethod(ts, "TSBDFSetOrder_C", (TS, PetscInt), (ts, order));
550211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
551211a84d6SLisandro Dalcin }
552211a84d6SLisandro Dalcin 
553211a84d6SLisandro Dalcin /*@
554*bcf0153eSBarry Smith   TSBDFGetOrder - Get the order of the `TSBDF` method
555211a84d6SLisandro Dalcin 
556211a84d6SLisandro Dalcin   Not Collective
557211a84d6SLisandro Dalcin 
558211a84d6SLisandro Dalcin   Input Parameter:
559211a84d6SLisandro Dalcin .  ts - timestepping context
560211a84d6SLisandro Dalcin 
561211a84d6SLisandro Dalcin   Output Parameter:
562211a84d6SLisandro Dalcin .  order - order of the method
563211a84d6SLisandro Dalcin 
564211a84d6SLisandro Dalcin   Level: intermediate
565211a84d6SLisandro Dalcin 
566*bcf0153eSBarry Smith .seealso: `TSBDFSetOrder()`, `TS`, `TSBDF`
567211a84d6SLisandro Dalcin @*/
568d71ae5a4SJacob Faibussowitsch PetscErrorCode TSBDFGetOrder(TS ts, PetscInt *order)
569d71ae5a4SJacob Faibussowitsch {
570211a84d6SLisandro Dalcin   PetscFunctionBegin;
571211a84d6SLisandro Dalcin   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
572211a84d6SLisandro Dalcin   PetscValidIntPointer(order, 2);
573cac4c232SBarry Smith   PetscUseMethod(ts, "TSBDFGetOrder_C", (TS, PetscInt *), (ts, order));
574211a84d6SLisandro Dalcin   PetscFunctionReturn(0);
575211a84d6SLisandro Dalcin }
576