xref: /petsc/src/ts/impls/implicit/glle/glle.c (revision 4ad8454beace47809662cdae21ee081016eaa39a)
1 #include <../src/ts/impls/implicit/glle/glle.h> /*I   "petscts.h"   I*/
2 #include <petscdm.h>
3 #include <petscblaslapack.h>
4 
5 static const char       *TSGLLEErrorDirections[] = {"FORWARD", "BACKWARD", "TSGLLEErrorDirection", "TSGLLEERROR_", NULL};
6 static PetscFunctionList TSGLLEList;
7 static PetscFunctionList TSGLLEAcceptList;
8 static PetscBool         TSGLLEPackageInitialized;
9 static PetscBool         TSGLLERegisterAllCalled;
10 
11 /* This function is pure */
12 static PetscScalar Factorial(PetscInt n)
13 {
14   PetscInt i;
15   if (n < 12) { /* Can compute with 32-bit integers */
16     PetscInt f = 1;
17     for (i = 2; i <= n; i++) f *= i;
18     return (PetscScalar)f;
19   } else {
20     PetscScalar f = 1.;
21     for (i = 2; i <= n; i++) f *= (PetscScalar)i;
22     return f;
23   }
24 }
25 
26 /* This function is pure */
27 static PetscScalar CPowF(PetscScalar c, PetscInt p)
28 {
29   return PetscPowRealInt(PetscRealPart(c), p) / Factorial(p);
30 }
31 
32 static PetscErrorCode TSGLLEGetVecs(TS ts, DM dm, Vec *Z, Vec *Ydotstage)
33 {
34   TS_GLLE *gl = (TS_GLLE *)ts->data;
35 
36   PetscFunctionBegin;
37   if (Z) {
38     if (dm && dm != ts->dm) {
39       PetscCall(DMGetNamedGlobalVector(dm, "TSGLLE_Z", Z));
40     } else *Z = gl->Z;
41   }
42   if (Ydotstage) {
43     if (dm && dm != ts->dm) {
44       PetscCall(DMGetNamedGlobalVector(dm, "TSGLLE_Ydot", Ydotstage));
45     } else *Ydotstage = gl->Ydot[gl->stage];
46   }
47   PetscFunctionReturn(PETSC_SUCCESS);
48 }
49 
50 static PetscErrorCode TSGLLERestoreVecs(TS ts, DM dm, Vec *Z, Vec *Ydotstage)
51 {
52   PetscFunctionBegin;
53   if (Z) {
54     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSGLLE_Z", Z));
55   }
56   if (Ydotstage) {
57     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSGLLE_Ydot", Ydotstage));
58   }
59   PetscFunctionReturn(PETSC_SUCCESS);
60 }
61 
62 static PetscErrorCode DMCoarsenHook_TSGLLE(DM fine, DM coarse, void *ctx)
63 {
64   PetscFunctionBegin;
65   PetscFunctionReturn(PETSC_SUCCESS);
66 }
67 
68 static PetscErrorCode DMRestrictHook_TSGLLE(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx)
69 {
70   TS  ts = (TS)ctx;
71   Vec Ydot, Ydot_c;
72 
73   PetscFunctionBegin;
74   PetscCall(TSGLLEGetVecs(ts, fine, NULL, &Ydot));
75   PetscCall(TSGLLEGetVecs(ts, coarse, NULL, &Ydot_c));
76   PetscCall(MatRestrict(restrct, Ydot, Ydot_c));
77   PetscCall(VecPointwiseMult(Ydot_c, rscale, Ydot_c));
78   PetscCall(TSGLLERestoreVecs(ts, fine, NULL, &Ydot));
79   PetscCall(TSGLLERestoreVecs(ts, coarse, NULL, &Ydot_c));
80   PetscFunctionReturn(PETSC_SUCCESS);
81 }
82 
83 static PetscErrorCode DMSubDomainHook_TSGLLE(DM dm, DM subdm, void *ctx)
84 {
85   PetscFunctionBegin;
86   PetscFunctionReturn(PETSC_SUCCESS);
87 }
88 
89 static PetscErrorCode DMSubDomainRestrictHook_TSGLLE(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
90 {
91   TS  ts = (TS)ctx;
92   Vec Ydot, Ydot_s;
93 
94   PetscFunctionBegin;
95   PetscCall(TSGLLEGetVecs(ts, dm, NULL, &Ydot));
96   PetscCall(TSGLLEGetVecs(ts, subdm, NULL, &Ydot_s));
97 
98   PetscCall(VecScatterBegin(gscat, Ydot, Ydot_s, INSERT_VALUES, SCATTER_FORWARD));
99   PetscCall(VecScatterEnd(gscat, Ydot, Ydot_s, INSERT_VALUES, SCATTER_FORWARD));
100 
101   PetscCall(TSGLLERestoreVecs(ts, dm, NULL, &Ydot));
102   PetscCall(TSGLLERestoreVecs(ts, subdm, NULL, &Ydot_s));
103   PetscFunctionReturn(PETSC_SUCCESS);
104 }
105 
106 static PetscErrorCode TSGLLESchemeCreate(PetscInt p, PetscInt q, PetscInt r, PetscInt s, const PetscScalar *c, const PetscScalar *a, const PetscScalar *b, const PetscScalar *u, const PetscScalar *v, TSGLLEScheme *inscheme)
107 {
108   TSGLLEScheme scheme;
109   PetscInt     j;
110 
111   PetscFunctionBegin;
112   PetscCheck(p >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Scheme order must be positive");
113   PetscCheck(r >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "At least one item must be carried between steps");
114   PetscCheck(s >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "At least one stage is required");
115   PetscAssertPointer(inscheme, 10);
116   *inscheme = NULL;
117   PetscCall(PetscNew(&scheme));
118   scheme->p = p;
119   scheme->q = q;
120   scheme->r = r;
121   scheme->s = s;
122 
123   PetscCall(PetscMalloc5(s, &scheme->c, s * s, &scheme->a, r * s, &scheme->b, r * s, &scheme->u, r * r, &scheme->v));
124   PetscCall(PetscArraycpy(scheme->c, c, s));
125   for (j = 0; j < s * s; j++) scheme->a[j] = (PetscAbsScalar(a[j]) < 1e-12) ? 0 : a[j];
126   for (j = 0; j < r * s; j++) scheme->b[j] = (PetscAbsScalar(b[j]) < 1e-12) ? 0 : b[j];
127   for (j = 0; j < s * r; j++) scheme->u[j] = (PetscAbsScalar(u[j]) < 1e-12) ? 0 : u[j];
128   for (j = 0; j < r * r; j++) scheme->v[j] = (PetscAbsScalar(v[j]) < 1e-12) ? 0 : v[j];
129 
130   PetscCall(PetscMalloc6(r, &scheme->alpha, r, &scheme->beta, r, &scheme->gamma, 3 * s, &scheme->phi, 3 * r, &scheme->psi, r, &scheme->stage_error));
131   {
132     PetscInt     i, j, k, ss = s + 2;
133     PetscBLASInt m, n, one = 1, *ipiv, lwork = 4 * ((s + 3) * 3 + 3), info, ldb;
134     PetscReal    rcond, *sing, *workreal;
135     PetscScalar *ImV, *H, *bmat, *workscalar, *c = scheme->c, *a = scheme->a, *b = scheme->b, *u = scheme->u, *v = scheme->v;
136     PetscBLASInt rank;
137     PetscCall(PetscMalloc7(PetscSqr(r), &ImV, 3 * s, &H, 3 * ss, &bmat, lwork, &workscalar, 5 * (3 + r), &workreal, r + s, &sing, r + s, &ipiv));
138 
139     /* column-major input */
140     for (i = 0; i < r - 1; i++) {
141       for (j = 0; j < r - 1; j++) ImV[i + j * r] = 1.0 * (i == j) - v[(i + 1) * r + j + 1];
142     }
143     /* Build right hand side for alpha (tp - glm.B(2:end,:)*(glm.c.^(p)./factorial(p))) */
144     for (i = 1; i < r; i++) {
145       scheme->alpha[i] = 1. / Factorial(p + 1 - i);
146       for (j = 0; j < s; j++) scheme->alpha[i] -= b[i * s + j] * CPowF(c[j], p);
147     }
148     PetscCall(PetscBLASIntCast(r - 1, &m));
149     PetscCall(PetscBLASIntCast(r, &n));
150     PetscCallBLAS("LAPACKgesv", LAPACKgesv_(&m, &one, ImV, &n, ipiv, scheme->alpha + 1, &n, &info));
151     PetscCheck(info >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Bad argument to GESV");
152     PetscCheck(info <= 0, PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Bad LU factorization");
153 
154     /* Build right hand side for beta (tp1 - glm.B(2:end,:)*(glm.c.^(p+1)./factorial(p+1)) - e.alpha) */
155     for (i = 1; i < r; i++) {
156       scheme->beta[i] = 1. / Factorial(p + 2 - i) - scheme->alpha[i];
157       for (j = 0; j < s; j++) scheme->beta[i] -= b[i * s + j] * CPowF(c[j], p + 1);
158     }
159     PetscCallBLAS("LAPACKgetrs", LAPACKgetrs_("No transpose", &m, &one, ImV, &n, ipiv, scheme->beta + 1, &n, &info));
160     PetscCheck(info >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Bad argument to GETRS");
161     PetscCheck(info <= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Should not happen");
162 
163     /* Build stage_error vector
164            xi = glm.c.^(p+1)/factorial(p+1) - glm.A*glm.c.^p/factorial(p) + glm.U(:,2:end)*e.alpha;
165     */
166     for (i = 0; i < s; i++) {
167       scheme->stage_error[i] = CPowF(c[i], p + 1);
168       for (j = 0; j < s; j++) scheme->stage_error[i] -= a[i * s + j] * CPowF(c[j], p);
169       for (j = 1; j < r; j++) scheme->stage_error[i] += u[i * r + j] * scheme->alpha[j];
170     }
171 
172     /* alpha[0] (epsilon in B,J,W 2007)
173            epsilon = 1/factorial(p+1) - B(1,:)*c.^p/factorial(p) + V(1,2:end)*e.alpha;
174     */
175     scheme->alpha[0] = 1. / Factorial(p + 1);
176     for (j = 0; j < s; j++) scheme->alpha[0] -= b[0 * s + j] * CPowF(c[j], p);
177     for (j = 1; j < r; j++) scheme->alpha[0] += v[0 * r + j] * scheme->alpha[j];
178 
179     /* right hand side for gamma (glm.B(2:end,:)*e.xi - e.epsilon*eye(s-1,1)) */
180     for (i = 1; i < r; i++) {
181       scheme->gamma[i] = (i == 1 ? -1. : 0) * scheme->alpha[0];
182       for (j = 0; j < s; j++) scheme->gamma[i] += b[i * s + j] * scheme->stage_error[j];
183     }
184     PetscCallBLAS("LAPACKgetrs", LAPACKgetrs_("No transpose", &m, &one, ImV, &n, ipiv, scheme->gamma + 1, &n, &info));
185     PetscCheck(info >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Bad argument to GETRS");
186     PetscCheck(info <= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Should not happen");
187 
188     /* beta[0] (rho in B,J,W 2007)
189         e.rho = 1/factorial(p+2) - glm.B(1,:)*glm.c.^(p+1)/factorial(p+1) ...
190             + glm.V(1,2:end)*e.beta;% - e.epsilon;
191     % Note: The paper (B,J,W 2007) includes the last term in their definition
192     * */
193     scheme->beta[0] = 1. / Factorial(p + 2);
194     for (j = 0; j < s; j++) scheme->beta[0] -= b[0 * s + j] * CPowF(c[j], p + 1);
195     for (j = 1; j < r; j++) scheme->beta[0] += v[0 * r + j] * scheme->beta[j];
196 
197     /* gamma[0] (sigma in B,J,W 2007)
198     *   e.sigma = glm.B(1,:)*e.xi + glm.V(1,2:end)*e.gamma;
199     * */
200     scheme->gamma[0] = 0.0;
201     for (j = 0; j < s; j++) scheme->gamma[0] += b[0 * s + j] * scheme->stage_error[j];
202     for (j = 1; j < r; j++) scheme->gamma[0] += v[0 * s + j] * scheme->gamma[j];
203 
204     /* Assemble H
205     *    % " PetscInt_FMT "etermine the error estimators phi
206        H = [[cpow(glm.c,p) + C*e.alpha] [cpow(glm.c,p+1) + C*e.beta] ...
207                [e.xi - C*(e.gamma + 0*e.epsilon*eye(s-1,1))]]';
208     % Paper has formula above without the 0, but that term must be left
209     % out to satisfy the conditions they propose and to make the
210     % example schemes work
211     e.H = H;
212     e.phi = (H \ [1 0 0;1 1 0;0 0 -1])';
213     e.psi = -e.phi*C;
214     * */
215     for (j = 0; j < s; j++) {
216       H[0 + j * 3] = CPowF(c[j], p);
217       H[1 + j * 3] = CPowF(c[j], p + 1);
218       H[2 + j * 3] = scheme->stage_error[j];
219       for (k = 1; k < r; k++) {
220         H[0 + j * 3] += CPowF(c[j], k - 1) * scheme->alpha[k];
221         H[1 + j * 3] += CPowF(c[j], k - 1) * scheme->beta[k];
222         H[2 + j * 3] -= CPowF(c[j], k - 1) * scheme->gamma[k];
223       }
224     }
225     bmat[0 + 0 * ss] = 1.;
226     bmat[0 + 1 * ss] = 0.;
227     bmat[0 + 2 * ss] = 0.;
228     bmat[1 + 0 * ss] = 1.;
229     bmat[1 + 1 * ss] = 1.;
230     bmat[1 + 2 * ss] = 0.;
231     bmat[2 + 0 * ss] = 0.;
232     bmat[2 + 1 * ss] = 0.;
233     bmat[2 + 2 * ss] = -1.;
234     m                = 3;
235     PetscCall(PetscBLASIntCast(s, &n));
236     PetscCall(PetscBLASIntCast(ss, &ldb));
237     rcond = 1e-12;
238 #if defined(PETSC_USE_COMPLEX)
239     /* ZGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, RWORK, INFO) */
240     PetscCallBLAS("LAPACKgelss", LAPACKgelss_(&m, &n, &m, H, &m, bmat, &ldb, sing, &rcond, &rank, workscalar, &lwork, workreal, &info));
241 #else
242     /* DGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, INFO) */
243     PetscCallBLAS("LAPACKgelss", LAPACKgelss_(&m, &n, &m, H, &m, bmat, &ldb, sing, &rcond, &rank, workscalar, &lwork, &info));
244 #endif
245     PetscCheck(info >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Bad argument to GELSS");
246     PetscCheck(info <= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "SVD failed to converge");
247 
248     for (j = 0; j < 3; j++) {
249       for (k = 0; k < s; k++) scheme->phi[k + j * s] = bmat[k + j * ss];
250     }
251 
252     /* the other part of the error estimator, psi in B,J,W 2007 */
253     scheme->psi[0 * r + 0] = 0.;
254     scheme->psi[1 * r + 0] = 0.;
255     scheme->psi[2 * r + 0] = 0.;
256     for (j = 1; j < r; j++) {
257       scheme->psi[0 * r + j] = 0.;
258       scheme->psi[1 * r + j] = 0.;
259       scheme->psi[2 * r + j] = 0.;
260       for (k = 0; k < s; k++) {
261         scheme->psi[0 * r + j] -= CPowF(c[k], j - 1) * scheme->phi[0 * s + k];
262         scheme->psi[1 * r + j] -= CPowF(c[k], j - 1) * scheme->phi[1 * s + k];
263         scheme->psi[2 * r + j] -= CPowF(c[k], j - 1) * scheme->phi[2 * s + k];
264       }
265     }
266     PetscCall(PetscFree7(ImV, H, bmat, workscalar, workreal, sing, ipiv));
267   }
268   /* Check which properties are satisfied */
269   scheme->stiffly_accurate = PETSC_TRUE;
270   if (scheme->c[s - 1] != 1.) scheme->stiffly_accurate = PETSC_FALSE;
271   for (j = 0; j < s; j++)
272     if (a[(s - 1) * s + j] != b[j]) scheme->stiffly_accurate = PETSC_FALSE;
273   for (j = 0; j < r; j++)
274     if (u[(s - 1) * r + j] != v[j]) scheme->stiffly_accurate = PETSC_FALSE;
275   scheme->fsal = scheme->stiffly_accurate; /* FSAL is stronger */
276   for (j = 0; j < s - 1; j++)
277     if (r > 1 && b[1 * s + j] != 0.) scheme->fsal = PETSC_FALSE;
278   if (b[1 * s + r - 1] != 1.) scheme->fsal = PETSC_FALSE;
279   for (j = 0; j < r; j++)
280     if (r > 1 && v[1 * r + j] != 0.) scheme->fsal = PETSC_FALSE;
281 
282   *inscheme = scheme;
283   PetscFunctionReturn(PETSC_SUCCESS);
284 }
285 
286 static PetscErrorCode TSGLLESchemeDestroy(TSGLLEScheme sc)
287 {
288   PetscFunctionBegin;
289   PetscCall(PetscFree5(sc->c, sc->a, sc->b, sc->u, sc->v));
290   PetscCall(PetscFree6(sc->alpha, sc->beta, sc->gamma, sc->phi, sc->psi, sc->stage_error));
291   PetscCall(PetscFree(sc));
292   PetscFunctionReturn(PETSC_SUCCESS);
293 }
294 
295 static PetscErrorCode TSGLLEDestroy_Default(TS_GLLE *gl)
296 {
297   PetscInt i;
298 
299   PetscFunctionBegin;
300   for (i = 0; i < gl->nschemes; i++) {
301     if (gl->schemes[i]) PetscCall(TSGLLESchemeDestroy(gl->schemes[i]));
302   }
303   PetscCall(PetscFree(gl->schemes));
304   gl->nschemes = 0;
305   PetscCall(PetscMemzero(gl->type_name, sizeof(gl->type_name)));
306   PetscFunctionReturn(PETSC_SUCCESS);
307 }
308 
309 static PetscErrorCode TSGLLEViewTable_Private(PetscViewer viewer, PetscInt m, PetscInt n, const PetscScalar a[], const char name[])
310 {
311   PetscBool iascii;
312   PetscInt  i, j;
313 
314   PetscFunctionBegin;
315   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
316   if (iascii) {
317     PetscCall(PetscViewerASCIIPrintf(viewer, "%30s = [", name));
318     for (i = 0; i < m; i++) {
319       if (i) PetscCall(PetscViewerASCIIPrintf(viewer, "%30s   [", ""));
320       PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
321       for (j = 0; j < n; j++) PetscCall(PetscViewerASCIIPrintf(viewer, " %12.8g", (double)PetscRealPart(a[i * n + j])));
322       PetscCall(PetscViewerASCIIPrintf(viewer, "]\n"));
323       PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
324     }
325   }
326   PetscFunctionReturn(PETSC_SUCCESS);
327 }
328 
329 static PetscErrorCode TSGLLESchemeView(TSGLLEScheme sc, PetscBool view_details, PetscViewer viewer)
330 {
331   PetscBool iascii;
332 
333   PetscFunctionBegin;
334   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
335   if (iascii) {
336     PetscCall(PetscViewerASCIIPrintf(viewer, "GL scheme p,q,r,s = %" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "\n", sc->p, sc->q, sc->r, sc->s));
337     PetscCall(PetscViewerASCIIPushTab(viewer));
338     PetscCall(PetscViewerASCIIPrintf(viewer, "Stiffly accurate: %s,  FSAL: %s\n", sc->stiffly_accurate ? "yes" : "no", sc->fsal ? "yes" : "no"));
339     PetscCall(PetscViewerASCIIPrintf(viewer, "Leading error constants: %10.3e  %10.3e  %10.3e\n", (double)PetscRealPart(sc->alpha[0]), (double)PetscRealPart(sc->beta[0]), (double)PetscRealPart(sc->gamma[0])));
340     PetscCall(TSGLLEViewTable_Private(viewer, 1, sc->s, sc->c, "Abscissas c"));
341     if (view_details) {
342       PetscCall(TSGLLEViewTable_Private(viewer, sc->s, sc->s, sc->a, "A"));
343       PetscCall(TSGLLEViewTable_Private(viewer, sc->r, sc->s, sc->b, "B"));
344       PetscCall(TSGLLEViewTable_Private(viewer, sc->s, sc->r, sc->u, "U"));
345       PetscCall(TSGLLEViewTable_Private(viewer, sc->r, sc->r, sc->v, "V"));
346 
347       PetscCall(TSGLLEViewTable_Private(viewer, 3, sc->s, sc->phi, "Error estimate phi"));
348       PetscCall(TSGLLEViewTable_Private(viewer, 3, sc->r, sc->psi, "Error estimate psi"));
349       PetscCall(TSGLLEViewTable_Private(viewer, 1, sc->r, sc->alpha, "Modify alpha"));
350       PetscCall(TSGLLEViewTable_Private(viewer, 1, sc->r, sc->beta, "Modify beta"));
351       PetscCall(TSGLLEViewTable_Private(viewer, 1, sc->r, sc->gamma, "Modify gamma"));
352       PetscCall(TSGLLEViewTable_Private(viewer, 1, sc->s, sc->stage_error, "Stage error xi"));
353     }
354     PetscCall(PetscViewerASCIIPopTab(viewer));
355   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Viewer type %s not supported", ((PetscObject)viewer)->type_name);
356   PetscFunctionReturn(PETSC_SUCCESS);
357 }
358 
359 static PetscErrorCode TSGLLEEstimateHigherMoments_Default(TSGLLEScheme sc, PetscReal h, Vec Ydot[], Vec Xold[], Vec hm[])
360 {
361   PetscInt i;
362 
363   PetscFunctionBegin;
364   PetscCheck(sc->r <= 64 && sc->s <= 64, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Ridiculous number of stages or items passed between stages");
365   /* build error vectors*/
366   for (i = 0; i < 3; i++) {
367     PetscScalar phih[64];
368     PetscInt    j;
369     for (j = 0; j < sc->s; j++) phih[j] = sc->phi[i * sc->s + j] * h;
370     PetscCall(VecZeroEntries(hm[i]));
371     PetscCall(VecMAXPY(hm[i], sc->s, phih, Ydot));
372     PetscCall(VecMAXPY(hm[i], sc->r, &sc->psi[i * sc->r], Xold));
373   }
374   PetscFunctionReturn(PETSC_SUCCESS);
375 }
376 
377 static PetscErrorCode TSGLLECompleteStep_Rescale(TSGLLEScheme sc, PetscReal h, TSGLLEScheme next_sc, PetscReal next_h, Vec Ydot[], Vec Xold[], Vec X[])
378 {
379   PetscScalar brow[32], vrow[32];
380   PetscInt    i, j, r, s;
381 
382   PetscFunctionBegin;
383   /* Build the new solution from (X,Ydot) */
384   r = sc->r;
385   s = sc->s;
386   for (i = 0; i < r; i++) {
387     PetscCall(VecZeroEntries(X[i]));
388     for (j = 0; j < s; j++) brow[j] = h * sc->b[i * s + j];
389     PetscCall(VecMAXPY(X[i], s, brow, Ydot));
390     for (j = 0; j < r; j++) vrow[j] = sc->v[i * r + j];
391     PetscCall(VecMAXPY(X[i], r, vrow, Xold));
392   }
393   PetscFunctionReturn(PETSC_SUCCESS);
394 }
395 
396 static PetscErrorCode TSGLLECompleteStep_RescaleAndModify(TSGLLEScheme sc, PetscReal h, TSGLLEScheme next_sc, PetscReal next_h, Vec Ydot[], Vec Xold[], Vec X[])
397 {
398   PetscScalar brow[32], vrow[32];
399   PetscReal   ratio;
400   PetscInt    i, j, p, r, s;
401 
402   PetscFunctionBegin;
403   /* Build the new solution from (X,Ydot) */
404   p     = sc->p;
405   r     = sc->r;
406   s     = sc->s;
407   ratio = next_h / h;
408   for (i = 0; i < r; i++) {
409     PetscCall(VecZeroEntries(X[i]));
410     for (j = 0; j < s; j++) {
411       brow[j] = h * (PetscPowRealInt(ratio, i) * sc->b[i * s + j] + (PetscPowRealInt(ratio, i) - PetscPowRealInt(ratio, p + 1)) * (+sc->alpha[i] * sc->phi[0 * s + j]) + (PetscPowRealInt(ratio, i) - PetscPowRealInt(ratio, p + 2)) * (+sc->beta[i] * sc->phi[1 * s + j] + sc->gamma[i] * sc->phi[2 * s + j]));
412     }
413     PetscCall(VecMAXPY(X[i], s, brow, Ydot));
414     for (j = 0; j < r; j++) {
415       vrow[j] = (PetscPowRealInt(ratio, i) * sc->v[i * r + j] + (PetscPowRealInt(ratio, i) - PetscPowRealInt(ratio, p + 1)) * (+sc->alpha[i] * sc->psi[0 * r + j]) + (PetscPowRealInt(ratio, i) - PetscPowRealInt(ratio, p + 2)) * (+sc->beta[i] * sc->psi[1 * r + j] + sc->gamma[i] * sc->psi[2 * r + j]));
416     }
417     PetscCall(VecMAXPY(X[i], r, vrow, Xold));
418   }
419   if (r < next_sc->r) {
420     PetscCheck(r + 1 == next_sc->r, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot accommodate jump in r greater than 1");
421     PetscCall(VecZeroEntries(X[r]));
422     for (j = 0; j < s; j++) brow[j] = h * PetscPowRealInt(ratio, p + 1) * sc->phi[0 * s + j];
423     PetscCall(VecMAXPY(X[r], s, brow, Ydot));
424     for (j = 0; j < r; j++) vrow[j] = PetscPowRealInt(ratio, p + 1) * sc->psi[0 * r + j];
425     PetscCall(VecMAXPY(X[r], r, vrow, Xold));
426   }
427   PetscFunctionReturn(PETSC_SUCCESS);
428 }
429 
430 static PetscErrorCode TSGLLECreate_IRKS(TS ts)
431 {
432   TS_GLLE *gl = (TS_GLLE *)ts->data;
433 
434   PetscFunctionBegin;
435   gl->Destroy               = TSGLLEDestroy_Default;
436   gl->EstimateHigherMoments = TSGLLEEstimateHigherMoments_Default;
437   gl->CompleteStep          = TSGLLECompleteStep_RescaleAndModify;
438   PetscCall(PetscMalloc1(10, &gl->schemes));
439   gl->nschemes = 0;
440 
441   {
442     /* p=1,q=1, r=s=2, A- and L-stable with error estimates of order 2 and 3
443     * Listed in Butcher & Podhaisky 2006. On error estimation in general linear methods for stiff ODE.
444     * irks(0.3,0,[.3,1],[1],1)
445     * Note: can be made to have classical order (not stage order) 2 by replacing 0.3 with 1-sqrt(1/2)
446     * but doing so would sacrifice the error estimator.
447     */
448     const PetscScalar c[2]    = {3. / 10., 1.};
449     const PetscScalar a[2][2] = {
450       {3. / 10., 0       },
451       {7. / 10., 3. / 10.}
452     };
453     const PetscScalar b[2][2] = {
454       {7. / 10., 3. / 10.},
455       {0,        1       }
456     };
457     const PetscScalar u[2][2] = {
458       {1, 0},
459       {1, 0}
460     };
461     const PetscScalar v[2][2] = {
462       {1, 0},
463       {0, 0}
464     };
465     PetscCall(TSGLLESchemeCreate(1, 1, 2, 2, c, *a, *b, *u, *v, &gl->schemes[gl->nschemes++]));
466   }
467 
468   {
469     /* p=q=2, r=s=3: irks(4/9,0,[1:3]/3,[0.33852],1) */
470     /* http://www.math.auckland.ac.nz/~hpod/atlas/i2a.html */
471     const PetscScalar c[3]    = {1. / 3., 2. / 3., 1};
472     const PetscScalar a[3][3] = {
473       {4. / 9.,              0,                     0      },
474       {1.03750643704090e+00, 4. / 9.,               0      },
475       {7.67024779410304e-01, -3.81140216918943e-01, 4. / 9.}
476     };
477     const PetscScalar b[3][3] = {
478       {0.767024779410304,  -0.381140216918943, 4. / 9.          },
479       {0.000000000000000,  0.000000000000000,  1.000000000000000},
480       {-2.075048385225385, 0.621728385225383,  1.277197204924873}
481     };
482     const PetscScalar u[3][3] = {
483       {1.0000000000000000, -0.1111111111111109, -0.0925925925925922},
484       {1.0000000000000000, -0.8152842148186744, -0.4199095530877056},
485       {1.0000000000000000, 0.1696709930641948,  0.0539741070314165 }
486     };
487     const PetscScalar v[3][3] = {
488       {1.0000000000000000, 0.1696709930641948, 0.0539741070314165},
489       {0.000000000000000,  0.000000000000000,  0.000000000000000 },
490       {0.000000000000000,  0.176122795075129,  0.000000000000000 }
491     };
492     PetscCall(TSGLLESchemeCreate(2, 2, 3, 3, c, *a, *b, *u, *v, &gl->schemes[gl->nschemes++]));
493   }
494   {
495     /* p=q=3, r=s=4: irks(9/40,0,[1:4]/4,[0.3312 1.0050],[0.49541 1;1 0]) */
496     const PetscScalar c[4]    = {0.25, 0.5, 0.75, 1.0};
497     const PetscScalar a[4][4] = {
498       {9. / 40.,             0,                     0,                 0       },
499       {2.11286958887701e-01, 9. / 40.,              0,                 0       },
500       {9.46338294287584e-01, -3.42942861246094e-01, 9. / 40.,          0       },
501       {0.521490453970721,    -0.662474225622980,    0.490476425459734, 9. / 40.}
502     };
503     const PetscScalar b[4][4] = {
504       {0.521490453970721,  -0.662474225622980, 0.490476425459734,  9. / 40.         },
505       {0.000000000000000,  0.000000000000000,  0.000000000000000,  1.000000000000000},
506       {-0.084677029310348, 1.390757514776085,  -1.568157386206001, 2.023192696767826},
507       {0.465383797936408,  1.478273530625148,  -1.930836081010182, 1.644872111193354}
508     };
509     const PetscScalar u[4][4] = {
510       {1.00000000000000000, 0.02500000000001035,  -0.02499999999999053, -0.00442708333332865},
511       {1.00000000000000000, 0.06371304111232945,  -0.04032173972189845, -0.01389438413189452},
512       {1.00000000000000000, -0.07839543304147778, 0.04738685705116663,  0.02032603595928376 },
513       {1.00000000000000000, 0.42550734619251651,  0.10800718022400080,  -0.01726712647760034}
514     };
515     const PetscScalar v[4][4] = {
516       {1.00000000000000000, 0.42550734619251651, 0.10800718022400080, -0.01726712647760034},
517       {0.000000000000000,   0.000000000000000,   0.000000000000000,   0.000000000000000   },
518       {0.000000000000000,   -1.761115796027561,  -0.521284157173780,  0.258249384305463   },
519       {0.000000000000000,   -1.657693358744728,  -1.052227765232394,  0.521284157173780   }
520     };
521     PetscCall(TSGLLESchemeCreate(3, 3, 4, 4, c, *a, *b, *u, *v, &gl->schemes[gl->nschemes++]));
522   }
523   {
524     /* p=q=4, r=s=5:
525           irks(3/11,0,[1:5]/5, [0.1715   -0.1238    0.6617],...
526           [ -0.0812    0.4079    1.0000
527              1.0000         0         0
528              0.8270    1.0000         0])
529     */
530     const PetscScalar c[5]    = {0.2, 0.4, 0.6, 0.8, 1.0};
531     const PetscScalar a[5][5] = {
532       {2.72727272727352e-01,  0.00000000000000e+00,  0.00000000000000e+00,  0.00000000000000e+00, 0.00000000000000e+00},
533       {-1.03980153733431e-01, 2.72727272727405e-01,  0.00000000000000e+00,  0.00000000000000e+00, 0.00000000000000e+00},
534       {-1.58615400341492e+00, 7.44168951881122e-01,  2.72727272727309e-01,  0.00000000000000e+00, 0.00000000000000e+00},
535       {-8.73658042865628e-01, 5.37884671894595e-01,  -1.63298538799523e-01, 2.72727272726996e-01, 0.00000000000000e+00},
536       {2.95489397443992e-01,  -1.18481693910097e+00, -6.68029812659953e-01, 1.00716687860943e+00, 2.72727272727288e-01}
537     };
538     const PetscScalar b[5][5] = {
539       {2.95489397443992e-01,  -1.18481693910097e+00, -6.68029812659953e-01, 1.00716687860943e+00,  2.72727272727288e-01},
540       {0.00000000000000e+00,  1.11022302462516e-16,  -2.22044604925031e-16, 0.00000000000000e+00,  1.00000000000000e+00},
541       {-4.05882503986005e+00, -4.00924006567769e+00, -1.38930610972481e+00, 4.45223930308488e+00,  6.32331093108427e-01},
542       {8.35690179937017e+00,  -2.26640927349732e+00, 6.86647884973826e+00,  -5.22595158025740e+00, 4.50893068837431e+00},
543       {1.27656267027479e+01,  2.80882153840821e+00,  8.91173096522890e+00,  -1.07936444078906e+01, 4.82534148988854e+00}
544     };
545     const PetscScalar u[5][5] = {
546       {1.00000000000000e+00, -7.27272727273551e-02, -3.45454545454419e-02, -4.12121212119565e-03, -2.96969696964014e-04},
547       {1.00000000000000e+00, 2.31252881006154e-01,  -8.29487834416481e-03, -9.07191207681020e-03, -1.70378403743473e-03},
548       {1.00000000000000e+00, 1.16925777880663e+00,  3.59268562942635e-02,  -4.09013451730615e-02, -1.02411119670164e-02},
549       {1.00000000000000e+00, 1.02634463704356e+00,  1.59375044913405e-01,  1.89673015035370e-03,  -4.89987231897569e-03},
550       {1.00000000000000e+00, 1.27746320298021e+00,  2.37186008132728e-01,  -8.28694373940065e-02, -5.34396510196430e-02}
551     };
552     const PetscScalar v[5][5] = {
553       {1.00000000000000e+00, 1.27746320298021e+00,  2.37186008132728e-01,  -8.28694373940065e-02, -5.34396510196430e-02},
554       {0.00000000000000e+00, -1.77635683940025e-15, -1.99840144432528e-15, -9.99200722162641e-16, -3.33066907387547e-16},
555       {0.00000000000000e+00, 4.37280081906924e+00,  5.49221645016377e-02,  -8.88913177394943e-02, 1.12879077989154e-01 },
556       {0.00000000000000e+00, -1.22399504837280e+01, -5.21287338448645e+00, -8.03952325565291e-01, 4.60298678047147e-01 },
557       {0.00000000000000e+00, -1.85178762883829e+01, -5.21411849862624e+00, -1.04283436528809e+00, 7.49030161063651e-01 }
558     };
559     PetscCall(TSGLLESchemeCreate(4, 4, 5, 5, c, *a, *b, *u, *v, &gl->schemes[gl->nschemes++]));
560   }
561   {
562     /* p=q=5, r=s=6;
563        irks(1/3,0,[1:6]/6,...
564           [-0.0489    0.4228   -0.8814    0.9021],...
565           [-0.3474   -0.6617    0.6294    0.2129
566             0.0044   -0.4256   -0.1427   -0.8936
567            -0.8267    0.4821    0.1371   -0.2557
568            -0.4426   -0.3855   -0.7514    0.3014])
569     */
570     const PetscScalar c[6]    = {1. / 6, 2. / 6, 3. / 6, 4. / 6, 5. / 6, 1.};
571     const PetscScalar a[6][6] = {
572       {3.33333333333940e-01,  0,                     0,                     0,                     0,                    0                   },
573       {-8.64423857333350e-02, 3.33333333332888e-01,  0,                     0,                     0,                    0                   },
574       {-2.16850174258252e+00, -2.23619072028839e+00, 3.33333333335204e-01,  0,                     0,                    0                   },
575       {-4.73160970138997e+00, -3.89265344629268e+00, -2.76318716520933e-01, 3.33333333335759e-01,  0,                    0                   },
576       {-6.75187540297338e+00, -7.90756533769377e+00, 7.90245051802259e-01,  -4.48352364517632e-01, 3.33333333328483e-01, 0                   },
577       {-4.26488287921548e+00, -1.19320395589302e+01, 3.38924509887755e+00,  -2.23969848002481e+00, 6.62807710124007e-01, 3.33333333335440e-01}
578     };
579     const PetscScalar b[6][6] = {
580       {-4.26488287921548e+00, -1.19320395589302e+01, 3.38924509887755e+00,  -2.23969848002481e+00, 6.62807710124007e-01,  3.33333333335440e-01 },
581       {-8.88178419700125e-16, 4.44089209850063e-16,  -1.54737334057131e-15, -8.88178419700125e-16, 0.00000000000000e+00,  1.00000000000001e+00 },
582       {-2.87780425770651e+01, -1.13520448264971e+01, 2.62002318943161e+01,  2.56943874812797e+01,  -3.06702268304488e+01, 6.68067773510103e+00 },
583       {5.47971245256474e+01,  6.80366875868284e+01,  -6.50952588861999e+01, -8.28643975339097e+01, 8.17416943896414e+01,  -1.17819043489036e+01},
584       {-2.33332114788869e+02, 6.12942539462634e+01,  -4.91850135865944e+01, 1.82716844135480e+02,  -1.29788173979395e+02, 3.09968095651099e+01 },
585       {-1.72049132343751e+02, 8.60194713593999e+00,  7.98154219170200e-01,  1.50371386053218e+02,  -1.18515423962066e+02, 2.50898277784663e+01 }
586     };
587     const PetscScalar u[6][6] = {
588       {1.00000000000000e+00, -1.66666666666870e-01, -4.16666666664335e-02, -3.85802469124815e-03, -2.25051440302250e-04, -9.64506172339142e-06},
589       {1.00000000000000e+00, 8.64423857327162e-02,  -4.11484912671353e-02, -1.11450903217645e-02, -1.47651050487126e-03, -1.34395070766826e-04},
590       {1.00000000000000e+00, 4.57135912953434e+00,  1.06514719719137e+00,  1.33517564218007e-01,  1.11365952968659e-02,  6.12382756769504e-04 },
591       {1.00000000000000e+00, 9.23391519753404e+00,  2.22431212392095e+00,  2.91823807741891e-01,  2.52058456411084e-02,  1.22800542949647e-03 },
592       {1.00000000000000e+00, 1.48175480533865e+01,  3.73439117461835e+00,  5.14648336541804e-01,  4.76430038853402e-02,  2.56798515502156e-03 },
593       {1.00000000000000e+00, 1.50512347758335e+01,  4.10099701165164e+00,  5.66039141003603e-01,  3.91213893800891e-02,  -2.99136269067853e-03}
594     };
595     const PetscScalar v[6][6] = {
596       {1.00000000000000e+00, 1.50512347758335e+01,  4.10099701165164e+00,  5.66039141003603e-01,  3.91213893800891e-02,  -2.99136269067853e-03},
597       {0.00000000000000e+00, -4.88498130835069e-15, -6.43929354282591e-15, -3.55271367880050e-15, -1.22124532708767e-15, -3.12250225675825e-16},
598       {0.00000000000000e+00, 1.22250171233141e+01,  -1.77150760606169e+00, 3.54516769879390e-01,  6.22298845883398e-01,  2.31647447450276e-01 },
599       {0.00000000000000e+00, -4.48339457331040e+01, -3.57363126641880e-01, 5.18750173123425e-01,  6.55727990241799e-02,  1.63175368287079e-01 },
600       {0.00000000000000e+00, 1.37297394708005e+02,  -1.60145272991317e+00, -5.05319555199441e+00, 1.55328940390990e-01,  9.16629423682464e-01 },
601       {0.00000000000000e+00, 1.05703241119022e+02,  -1.16610260983038e+00, -2.99767252773859e+00, -1.13472315553890e-01, 1.09742849254729e+00 }
602     };
603     PetscCall(TSGLLESchemeCreate(5, 5, 6, 6, c, *a, *b, *u, *v, &gl->schemes[gl->nschemes++]));
604   }
605   PetscFunctionReturn(PETSC_SUCCESS);
606 }
607 
608 /*@C
609   TSGLLESetType - sets the class of general linear method, `TSGLLE` to use for time-stepping
610 
611   Collective
612 
613   Input Parameters:
614 + ts   - the `TS` context
615 - type - a method
616 
617   Options Database Key:
618 . -ts_gl_type <type> - sets the method, use -help for a list of available method (e.g. irks)
619 
620   Level: intermediate
621 
622   Notes:
623   See "petsc/include/petscts.h" for available methods (for instance)
624 .    TSGLLE_IRKS - Diagonally implicit methods with inherent Runge-Kutta stability (for stiff problems)
625 
626   Normally, it is best to use the `TSSetFromOptions()` command and then set the `TSGLLE` type
627   from the options database rather than by using this routine.  Using the options database
628   provides the user with maximum flexibility in evaluating the many different solvers.  The
629   `TSGLLESetType()` routine is provided for those situations where it is necessary to set the
630   timestepping solver independently of the command line or options database.  This might be the
631   case, for example, when the choice of solver changes during the execution of the program, and
632   the user's application is taking responsibility for choosing the appropriate method.
633 
634 .seealso: [](ch_ts), `TS`, `TSGLLEType`, `TSGLLE`
635 @*/
636 PetscErrorCode TSGLLESetType(TS ts, TSGLLEType type)
637 {
638   PetscFunctionBegin;
639   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
640   PetscAssertPointer(type, 2);
641   PetscTryMethod(ts, "TSGLLESetType_C", (TS, TSGLLEType), (ts, type));
642   PetscFunctionReturn(PETSC_SUCCESS);
643 }
644 
645 /*@C
646   TSGLLESetAcceptType - sets the acceptance test for `TSGLLE`
647 
648   Logically Collective
649 
650   Input Parameters:
651 + ts   - the `TS` context
652 - type - the type
653 
654   Options Database Key:
655 . -ts_gl_accept_type <type> - sets the method used to determine whether to accept or reject a step
656 
657   Level: intermediate
658 
659   Notes:
660   Time integrators that need to control error must have the option to reject a time step based
661   on local error estimates. This function allows different schemes to be set.
662 
663 .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLEAcceptRegister()`, `TSGLLEAdapt`
664 @*/
665 PetscErrorCode TSGLLESetAcceptType(TS ts, TSGLLEAcceptType type)
666 {
667   PetscFunctionBegin;
668   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
669   PetscAssertPointer(type, 2);
670   PetscTryMethod(ts, "TSGLLESetAcceptType_C", (TS, TSGLLEAcceptType), (ts, type));
671   PetscFunctionReturn(PETSC_SUCCESS);
672 }
673 
674 /*@C
675   TSGLLEGetAdapt - gets the `TSGLLEAdapt` object from the `TS`
676 
677   Not Collective
678 
679   Input Parameter:
680 . ts - the `TS` context
681 
682   Output Parameter:
683 . adapt - the `TSGLLEAdapt` context
684 
685   Level: advanced
686 
687   Note:
688   This allows the user set options on the `TSGLLEAdapt` object. Usually it is better to do this
689   using the options database, so this function is rarely needed.
690 
691 .seealso: [](ch_ts), `TS`, `TSGLLE`, `TSGLLEAdapt`, `TSGLLEAdaptRegister()`
692 @*/
693 PetscErrorCode TSGLLEGetAdapt(TS ts, TSGLLEAdapt *adapt)
694 {
695   PetscFunctionBegin;
696   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
697   PetscAssertPointer(adapt, 2);
698   PetscUseMethod(ts, "TSGLLEGetAdapt_C", (TS, TSGLLEAdapt *), (ts, adapt));
699   PetscFunctionReturn(PETSC_SUCCESS);
700 }
701 
702 static PetscErrorCode TSGLLEAccept_Always(TS ts, PetscReal tleft, PetscReal h, const PetscReal enorms[], PetscBool *accept)
703 {
704   PetscFunctionBegin;
705   *accept = PETSC_TRUE;
706   PetscFunctionReturn(PETSC_SUCCESS);
707 }
708 
709 static PetscErrorCode TSGLLEUpdateWRMS(TS ts)
710 {
711   TS_GLLE     *gl = (TS_GLLE *)ts->data;
712   PetscScalar *x, *w;
713   PetscInt     n, i;
714 
715   PetscFunctionBegin;
716   PetscCall(VecGetArray(gl->X[0], &x));
717   PetscCall(VecGetArray(gl->W, &w));
718   PetscCall(VecGetLocalSize(gl->W, &n));
719   for (i = 0; i < n; i++) w[i] = 1. / (gl->wrms_atol + gl->wrms_rtol * PetscAbsScalar(x[i]));
720   PetscCall(VecRestoreArray(gl->X[0], &x));
721   PetscCall(VecRestoreArray(gl->W, &w));
722   PetscFunctionReturn(PETSC_SUCCESS);
723 }
724 
725 static PetscErrorCode TSGLLEVecNormWRMS(TS ts, Vec X, PetscReal *nrm)
726 {
727   TS_GLLE     *gl = (TS_GLLE *)ts->data;
728   PetscScalar *x, *w;
729   PetscReal    sum = 0.0, gsum;
730   PetscInt     n, N, i;
731 
732   PetscFunctionBegin;
733   PetscCall(VecGetArray(X, &x));
734   PetscCall(VecGetArray(gl->W, &w));
735   PetscCall(VecGetLocalSize(gl->W, &n));
736   for (i = 0; i < n; i++) sum += PetscAbsScalar(PetscSqr(x[i] * w[i]));
737   PetscCall(VecRestoreArray(X, &x));
738   PetscCall(VecRestoreArray(gl->W, &w));
739   PetscCall(MPIU_Allreduce(&sum, &gsum, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)ts)));
740   PetscCall(VecGetSize(gl->W, &N));
741   *nrm = PetscSqrtReal(gsum / (1. * N));
742   PetscFunctionReturn(PETSC_SUCCESS);
743 }
744 
745 static PetscErrorCode TSGLLESetType_GLLE(TS ts, TSGLLEType type)
746 {
747   PetscBool same;
748   TS_GLLE  *gl = (TS_GLLE *)ts->data;
749   PetscErrorCode (*r)(TS);
750 
751   PetscFunctionBegin;
752   if (gl->type_name[0]) {
753     PetscCall(PetscStrcmp(gl->type_name, type, &same));
754     if (same) PetscFunctionReturn(PETSC_SUCCESS);
755     PetscCall((*gl->Destroy)(gl));
756   }
757 
758   PetscCall(PetscFunctionListFind(TSGLLEList, type, &r));
759   PetscCheck(r, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSGLLE type \"%s\" given", type);
760   PetscCall((*r)(ts));
761   PetscCall(PetscStrncpy(gl->type_name, type, sizeof(gl->type_name)));
762   PetscFunctionReturn(PETSC_SUCCESS);
763 }
764 
765 static PetscErrorCode TSGLLESetAcceptType_GLLE(TS ts, TSGLLEAcceptType type)
766 {
767   TSGLLEAcceptFn *r;
768   TS_GLLE        *gl = (TS_GLLE *)ts->data;
769 
770   PetscFunctionBegin;
771   PetscCall(PetscFunctionListFind(TSGLLEAcceptList, type, &r));
772   PetscCheck(r, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown TSGLLEAccept type \"%s\" given", type);
773   gl->Accept = r;
774   PetscCall(PetscStrncpy(gl->accept_name, type, sizeof(gl->accept_name)));
775   PetscFunctionReturn(PETSC_SUCCESS);
776 }
777 
778 static PetscErrorCode TSGLLEGetAdapt_GLLE(TS ts, TSGLLEAdapt *adapt)
779 {
780   TS_GLLE *gl = (TS_GLLE *)ts->data;
781 
782   PetscFunctionBegin;
783   if (!gl->adapt) {
784     PetscCall(TSGLLEAdaptCreate(PetscObjectComm((PetscObject)ts), &gl->adapt));
785     PetscCall(PetscObjectIncrementTabLevel((PetscObject)gl->adapt, (PetscObject)ts, 1));
786   }
787   *adapt = gl->adapt;
788   PetscFunctionReturn(PETSC_SUCCESS);
789 }
790 
791 static PetscErrorCode TSGLLEChooseNextScheme(TS ts, PetscReal h, const PetscReal hmnorm[], PetscInt *next_scheme, PetscReal *next_h, PetscBool *finish)
792 {
793   TS_GLLE  *gl = (TS_GLLE *)ts->data;
794   PetscInt  i, n, cur_p, cur, next_sc, candidates[64], orders[64];
795   PetscReal errors[64], costs[64], tleft;
796 
797   PetscFunctionBegin;
798   cur   = -1;
799   cur_p = gl->schemes[gl->current_scheme]->p;
800   tleft = ts->max_time - (ts->ptime + ts->time_step);
801   for (i = 0, n = 0; i < gl->nschemes; i++) {
802     TSGLLEScheme sc = gl->schemes[i];
803     if (sc->p < gl->min_order || gl->max_order < sc->p) continue;
804     if (sc->p == cur_p - 1) errors[n] = PetscAbsScalar(sc->alpha[0]) * hmnorm[0];
805     else if (sc->p == cur_p) errors[n] = PetscAbsScalar(sc->alpha[0]) * hmnorm[1];
806     else if (sc->p == cur_p + 1) errors[n] = PetscAbsScalar(sc->alpha[0]) * (hmnorm[2] + hmnorm[3]);
807     else continue;
808     candidates[n] = i;
809     orders[n]     = PetscMin(sc->p, sc->q); /* order of global truncation error */
810     costs[n]      = sc->s;                  /* estimate the cost as the number of stages */
811     if (i == gl->current_scheme) cur = n;
812     n++;
813   }
814   PetscCheck(cur >= 0 && gl->nschemes > cur, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Current scheme not found in scheme list");
815   PetscCall(TSGLLEAdaptChoose(gl->adapt, n, orders, errors, costs, cur, h, tleft, &next_sc, next_h, finish));
816   *next_scheme = candidates[next_sc];
817   PetscCall(PetscInfo(ts, "Adapt chose scheme %" PetscInt_FMT " (%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ") with step size %6.2e, finish=%s\n", *next_scheme, gl->schemes[*next_scheme]->p, gl->schemes[*next_scheme]->q,
818                       gl->schemes[*next_scheme]->r, gl->schemes[*next_scheme]->s, (double)*next_h, PetscBools[*finish]));
819   PetscFunctionReturn(PETSC_SUCCESS);
820 }
821 
822 static PetscErrorCode TSGLLEGetMaxSizes(TS ts, PetscInt *max_r, PetscInt *max_s)
823 {
824   TS_GLLE *gl = (TS_GLLE *)ts->data;
825 
826   PetscFunctionBegin;
827   *max_r = gl->schemes[gl->nschemes - 1]->r;
828   *max_s = gl->schemes[gl->nschemes - 1]->s;
829   PetscFunctionReturn(PETSC_SUCCESS);
830 }
831 
832 static PetscErrorCode TSSolve_GLLE(TS ts)
833 {
834   TS_GLLE            *gl = (TS_GLLE *)ts->data;
835   PetscInt            i, k, its, lits, max_r, max_s;
836   PetscBool           final_step, finish;
837   SNESConvergedReason snesreason;
838 
839   PetscFunctionBegin;
840   PetscCall(TSMonitor(ts, ts->steps, ts->ptime, ts->vec_sol));
841 
842   PetscCall(TSGLLEGetMaxSizes(ts, &max_r, &max_s));
843   PetscCall(VecCopy(ts->vec_sol, gl->X[0]));
844   for (i = 1; i < max_r; i++) PetscCall(VecZeroEntries(gl->X[i]));
845   PetscCall(TSGLLEUpdateWRMS(ts));
846 
847   if (0) {
848     /* Find consistent initial data for DAE */
849     gl->stage_time = ts->ptime + ts->time_step;
850     gl->scoeff     = 1.;
851     gl->stage      = 0;
852 
853     PetscCall(VecCopy(ts->vec_sol, gl->Z));
854     PetscCall(VecCopy(ts->vec_sol, gl->Y));
855     PetscCall(SNESSolve(ts->snes, NULL, gl->Y));
856     PetscCall(SNESGetIterationNumber(ts->snes, &its));
857     PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits));
858     PetscCall(SNESGetConvergedReason(ts->snes, &snesreason));
859 
860     ts->snes_its += its;
861     ts->ksp_its += lits;
862     if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) {
863       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
864       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", nonlinear solve failures %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, ts->num_snes_failures));
865       PetscFunctionReturn(PETSC_SUCCESS);
866     }
867   }
868 
869   PetscCheck(gl->current_scheme >= 0, PETSC_COMM_SELF, PETSC_ERR_ORDER, "A starting scheme has not been provided");
870 
871   for (k = 0, final_step = PETSC_FALSE, finish = PETSC_FALSE; k < ts->max_steps && !finish; k++) {
872     PetscInt           j, r, s, next_scheme = 0, rejections;
873     PetscReal          h, hmnorm[4], enorm[3], next_h;
874     PetscBool          accept;
875     const PetscScalar *c, *a, *u;
876     Vec               *X, *Ydot, Y;
877     TSGLLEScheme       scheme = gl->schemes[gl->current_scheme];
878 
879     r    = scheme->r;
880     s    = scheme->s;
881     c    = scheme->c;
882     a    = scheme->a;
883     u    = scheme->u;
884     h    = ts->time_step;
885     X    = gl->X;
886     Ydot = gl->Ydot;
887     Y    = gl->Y;
888 
889     if (ts->ptime > ts->max_time) break;
890 
891     /*
892       We only call PreStep at the start of each STEP, not each STAGE.  This is because it is
893       possible to fail (have to restart a step) after multiple stages.
894     */
895     PetscCall(TSPreStep(ts));
896 
897     rejections = 0;
898     while (1) {
899       for (i = 0; i < s; i++) {
900         PetscScalar shift;
901         gl->scoeff     = 1. / PetscRealPart(a[i * s + i]);
902         shift          = gl->scoeff / ts->time_step;
903         gl->stage      = i;
904         gl->stage_time = ts->ptime + PetscRealPart(c[i]) * h;
905 
906         /*
907         * Stage equation: Y = h A Y' + U X
908         * We assume that A is lower-triangular so that we can solve the stages (Y,Y') sequentially
909         * Build the affine vector z_i = -[1/(h a_ii)](h sum_j a_ij y'_j + sum_j u_ij x_j)
910         * Then y'_i = z + 1/(h a_ii) y_i
911         */
912         PetscCall(VecZeroEntries(gl->Z));
913         for (j = 0; j < r; j++) PetscCall(VecAXPY(gl->Z, -shift * u[i * r + j], X[j]));
914         for (j = 0; j < i; j++) PetscCall(VecAXPY(gl->Z, -shift * h * a[i * s + j], Ydot[j]));
915         /* Note: Z is used within function evaluation, Ydot = Z + shift*Y */
916 
917         /* Compute an estimate of Y to start Newton iteration */
918         if (gl->extrapolate) {
919           if (i == 0) {
920             /* Linear extrapolation on the first stage */
921             PetscCall(VecWAXPY(Y, c[i] * h, X[1], X[0]));
922           } else {
923             /* Linear extrapolation from the last stage */
924             PetscCall(VecAXPY(Y, (c[i] - c[i - 1]) * h, Ydot[i - 1]));
925           }
926         } else if (i == 0) { /* Directly use solution from the last step, otherwise reuse the last stage (do nothing) */
927           PetscCall(VecCopy(X[0], Y));
928         }
929 
930         /* Solve this stage (Ydot[i] is computed during function evaluation) */
931         PetscCall(SNESSolve(ts->snes, NULL, Y));
932         PetscCall(SNESGetIterationNumber(ts->snes, &its));
933         PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits));
934         PetscCall(SNESGetConvergedReason(ts->snes, &snesreason));
935         ts->snes_its += its;
936         ts->ksp_its += lits;
937         if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) {
938           ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
939           PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", nonlinear solve failures %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, ts->num_snes_failures));
940           PetscFunctionReturn(PETSC_SUCCESS);
941         }
942       }
943 
944       gl->stage_time = ts->ptime + ts->time_step;
945 
946       PetscCall((*gl->EstimateHigherMoments)(scheme, h, Ydot, gl->X, gl->himom));
947       /* hmnorm[i] = h^{p+i}x^{(p+i)} with i=0,1,2; hmnorm[3] = h^{p+2}(dx'/dx) x^{(p+1)} */
948       for (i = 0; i < 3; i++) PetscCall(TSGLLEVecNormWRMS(ts, gl->himom[i], &hmnorm[i + 1]));
949       enorm[0] = PetscRealPart(scheme->alpha[0]) * hmnorm[1];
950       enorm[1] = PetscRealPart(scheme->beta[0]) * hmnorm[2];
951       enorm[2] = PetscRealPart(scheme->gamma[0]) * hmnorm[3];
952       PetscCall((*gl->Accept)(ts, ts->max_time - gl->stage_time, h, enorm, &accept));
953       if (accept) goto accepted;
954       rejections++;
955       PetscCall(PetscInfo(ts, "Step %" PetscInt_FMT " (t=%g) not accepted, rejections=%" PetscInt_FMT "\n", k, (double)gl->stage_time, rejections));
956       if (rejections > gl->max_step_rejections) break;
957       /*
958         There are lots of reasons why a step might be rejected, including solvers not converging and other factors that
959         TSGLLEChooseNextScheme does not support.  Additionally, the error estimates may be very screwed up, so I'm not
960         convinced that it's safe to just compute a new error estimate using the same interface as the current adaptor
961         (the adaptor interface probably has to change).  Here we make an arbitrary and naive choice.  This assumes that
962         steps were written in Nordsieck form.  The "correct" method would be to re-complete the previous time step with
963         the correct "next" step size.  It is unclear to me whether the present ad-hoc method of rescaling X is stable.
964       */
965       h *= 0.5;
966       for (i = 1; i < scheme->r; i++) PetscCall(VecScale(X[i], PetscPowRealInt(0.5, i)));
967     }
968     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_CONV_FAILED, "Time step %" PetscInt_FMT " (t=%g) not accepted after %" PetscInt_FMT " failures", k, (double)gl->stage_time, rejections);
969 
970   accepted:
971     /* This term is not error, but it *would* be the leading term for a lower order method */
972     PetscCall(TSGLLEVecNormWRMS(ts, gl->X[scheme->r - 1], &hmnorm[0]));
973     /* Correct scaling so that these are equivalent to norms of the Nordsieck vectors */
974 
975     PetscCall(PetscInfo(ts, "Last moment norm %10.2e, estimated error norms %10.2e %10.2e %10.2e\n", (double)hmnorm[0], (double)enorm[0], (double)enorm[1], (double)enorm[2]));
976     if (!final_step) {
977       PetscCall(TSGLLEChooseNextScheme(ts, h, hmnorm, &next_scheme, &next_h, &final_step));
978     } else {
979       /* Dummy values to complete the current step in a consistent manner */
980       next_scheme = gl->current_scheme;
981       next_h      = h;
982       finish      = PETSC_TRUE;
983     }
984 
985     X        = gl->Xold;
986     gl->Xold = gl->X;
987     gl->X    = X;
988     PetscCall((*gl->CompleteStep)(scheme, h, gl->schemes[next_scheme], next_h, Ydot, gl->Xold, gl->X));
989 
990     PetscCall(TSGLLEUpdateWRMS(ts));
991 
992     /* Post the solution for the user, we could avoid this copy with a small bit of cleverness */
993     PetscCall(VecCopy(gl->X[0], ts->vec_sol));
994     ts->ptime += h;
995     ts->steps++;
996 
997     PetscCall(TSPostEvaluate(ts));
998     PetscCall(TSPostStep(ts));
999     PetscCall(TSMonitor(ts, ts->steps, ts->ptime, ts->vec_sol));
1000 
1001     gl->current_scheme = next_scheme;
1002     ts->time_step      = next_h;
1003   }
1004   PetscFunctionReturn(PETSC_SUCCESS);
1005 }
1006 
1007 /*------------------------------------------------------------*/
1008 
1009 static PetscErrorCode TSReset_GLLE(TS ts)
1010 {
1011   TS_GLLE *gl = (TS_GLLE *)ts->data;
1012   PetscInt max_r, max_s;
1013 
1014   PetscFunctionBegin;
1015   if (gl->setupcalled) {
1016     PetscCall(TSGLLEGetMaxSizes(ts, &max_r, &max_s));
1017     PetscCall(VecDestroyVecs(max_r, &gl->Xold));
1018     PetscCall(VecDestroyVecs(max_r, &gl->X));
1019     PetscCall(VecDestroyVecs(max_s, &gl->Ydot));
1020     PetscCall(VecDestroyVecs(3, &gl->himom));
1021     PetscCall(VecDestroy(&gl->W));
1022     PetscCall(VecDestroy(&gl->Y));
1023     PetscCall(VecDestroy(&gl->Z));
1024   }
1025   gl->setupcalled = PETSC_FALSE;
1026   PetscFunctionReturn(PETSC_SUCCESS);
1027 }
1028 
1029 static PetscErrorCode TSDestroy_GLLE(TS ts)
1030 {
1031   TS_GLLE *gl = (TS_GLLE *)ts->data;
1032 
1033   PetscFunctionBegin;
1034   PetscCall(TSReset_GLLE(ts));
1035   if (ts->dm) {
1036     PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSGLLE, DMRestrictHook_TSGLLE, ts));
1037     PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSGLLE, DMSubDomainRestrictHook_TSGLLE, ts));
1038   }
1039   if (gl->adapt) PetscCall(TSGLLEAdaptDestroy(&gl->adapt));
1040   if (gl->Destroy) PetscCall((*gl->Destroy)(gl));
1041   PetscCall(PetscFree(ts->data));
1042   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLESetType_C", NULL));
1043   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLESetAcceptType_C", NULL));
1044   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLEGetAdapt_C", NULL));
1045   PetscFunctionReturn(PETSC_SUCCESS);
1046 }
1047 
1048 /*
1049     This defines the nonlinear equation that is to be solved with SNES
1050     g(x) = f(t,x,z+shift*x) = 0
1051 */
1052 static PetscErrorCode SNESTSFormFunction_GLLE(SNES snes, Vec x, Vec f, TS ts)
1053 {
1054   TS_GLLE *gl = (TS_GLLE *)ts->data;
1055   Vec      Z, Ydot;
1056   DM       dm, dmsave;
1057 
1058   PetscFunctionBegin;
1059   PetscCall(SNESGetDM(snes, &dm));
1060   PetscCall(TSGLLEGetVecs(ts, dm, &Z, &Ydot));
1061   PetscCall(VecWAXPY(Ydot, gl->scoeff / ts->time_step, x, Z));
1062   dmsave = ts->dm;
1063   ts->dm = dm;
1064   PetscCall(TSComputeIFunction(ts, gl->stage_time, x, Ydot, f, PETSC_FALSE));
1065   ts->dm = dmsave;
1066   PetscCall(TSGLLERestoreVecs(ts, dm, &Z, &Ydot));
1067   PetscFunctionReturn(PETSC_SUCCESS);
1068 }
1069 
1070 static PetscErrorCode SNESTSFormJacobian_GLLE(SNES snes, Vec x, Mat A, Mat B, TS ts)
1071 {
1072   TS_GLLE *gl = (TS_GLLE *)ts->data;
1073   Vec      Z, Ydot;
1074   DM       dm, dmsave;
1075 
1076   PetscFunctionBegin;
1077   PetscCall(SNESGetDM(snes, &dm));
1078   PetscCall(TSGLLEGetVecs(ts, dm, &Z, &Ydot));
1079   dmsave = ts->dm;
1080   ts->dm = dm;
1081   /* gl->Xdot will have already been computed in SNESTSFormFunction_GLLE */
1082   PetscCall(TSComputeIJacobian(ts, gl->stage_time, x, gl->Ydot[gl->stage], gl->scoeff / ts->time_step, A, B, PETSC_FALSE));
1083   ts->dm = dmsave;
1084   PetscCall(TSGLLERestoreVecs(ts, dm, &Z, &Ydot));
1085   PetscFunctionReturn(PETSC_SUCCESS);
1086 }
1087 
1088 static PetscErrorCode TSSetUp_GLLE(TS ts)
1089 {
1090   TS_GLLE *gl = (TS_GLLE *)ts->data;
1091   PetscInt max_r, max_s;
1092   DM       dm;
1093 
1094   PetscFunctionBegin;
1095   gl->setupcalled = PETSC_TRUE;
1096   PetscCall(TSGLLEGetMaxSizes(ts, &max_r, &max_s));
1097   PetscCall(VecDuplicateVecs(ts->vec_sol, max_r, &gl->X));
1098   PetscCall(VecDuplicateVecs(ts->vec_sol, max_r, &gl->Xold));
1099   PetscCall(VecDuplicateVecs(ts->vec_sol, max_s, &gl->Ydot));
1100   PetscCall(VecDuplicateVecs(ts->vec_sol, 3, &gl->himom));
1101   PetscCall(VecDuplicate(ts->vec_sol, &gl->W));
1102   PetscCall(VecDuplicate(ts->vec_sol, &gl->Y));
1103   PetscCall(VecDuplicate(ts->vec_sol, &gl->Z));
1104 
1105   /* Default acceptance tests and adaptivity */
1106   if (!gl->Accept) PetscCall(TSGLLESetAcceptType(ts, TSGLLEACCEPT_ALWAYS));
1107   if (!gl->adapt) PetscCall(TSGLLEGetAdapt(ts, &gl->adapt));
1108 
1109   if (gl->current_scheme < 0) {
1110     PetscInt i;
1111     for (i = 0;; i++) {
1112       if (gl->schemes[i]->p == gl->start_order) break;
1113       PetscCheck(i + 1 != gl->nschemes, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No schemes available with requested start order %" PetscInt_FMT, i);
1114     }
1115     gl->current_scheme = i;
1116   }
1117   PetscCall(TSGetDM(ts, &dm));
1118   PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSGLLE, DMRestrictHook_TSGLLE, ts));
1119   PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSGLLE, DMSubDomainRestrictHook_TSGLLE, ts));
1120   PetscFunctionReturn(PETSC_SUCCESS);
1121 }
1122 /*------------------------------------------------------------*/
1123 
1124 static PetscErrorCode TSSetFromOptions_GLLE(TS ts, PetscOptionItems *PetscOptionsObject)
1125 {
1126   TS_GLLE *gl         = (TS_GLLE *)ts->data;
1127   char     tname[256] = TSGLLE_IRKS, completef[256] = "rescale-and-modify";
1128 
1129   PetscFunctionBegin;
1130   PetscOptionsHeadBegin(PetscOptionsObject, "General Linear ODE solver options");
1131   {
1132     PetscBool flg;
1133     PetscCall(PetscOptionsFList("-ts_gl_type", "Type of GL method", "TSGLLESetType", TSGLLEList, gl->type_name[0] ? gl->type_name : tname, tname, sizeof(tname), &flg));
1134     if (flg || !gl->type_name[0]) PetscCall(TSGLLESetType(ts, tname));
1135     PetscCall(PetscOptionsInt("-ts_gl_max_step_rejections", "Maximum number of times to attempt a step", "None", gl->max_step_rejections, &gl->max_step_rejections, NULL));
1136     PetscCall(PetscOptionsInt("-ts_gl_max_order", "Maximum order to try", "TSGLLESetMaxOrder", gl->max_order, &gl->max_order, NULL));
1137     PetscCall(PetscOptionsInt("-ts_gl_min_order", "Minimum order to try", "TSGLLESetMinOrder", gl->min_order, &gl->min_order, NULL));
1138     PetscCall(PetscOptionsInt("-ts_gl_start_order", "Initial order to try", "TSGLLESetMinOrder", gl->start_order, &gl->start_order, NULL));
1139     PetscCall(PetscOptionsEnum("-ts_gl_error_direction", "Which direction to look when estimating error", "TSGLLESetErrorDirection", TSGLLEErrorDirections, (PetscEnum)gl->error_direction, (PetscEnum *)&gl->error_direction, NULL));
1140     PetscCall(PetscOptionsBool("-ts_gl_extrapolate", "Extrapolate stage solution from previous solution (sometimes unstable)", "TSGLLESetExtrapolate", gl->extrapolate, &gl->extrapolate, NULL));
1141     PetscCall(PetscOptionsReal("-ts_gl_atol", "Absolute tolerance", "TSGLLESetTolerances", gl->wrms_atol, &gl->wrms_atol, NULL));
1142     PetscCall(PetscOptionsReal("-ts_gl_rtol", "Relative tolerance", "TSGLLESetTolerances", gl->wrms_rtol, &gl->wrms_rtol, NULL));
1143     PetscCall(PetscOptionsString("-ts_gl_complete", "Method to use for completing the step", "none", completef, completef, sizeof(completef), &flg));
1144     if (flg) {
1145       PetscBool match1, match2;
1146       PetscCall(PetscStrcmp(completef, "rescale", &match1));
1147       PetscCall(PetscStrcmp(completef, "rescale-and-modify", &match2));
1148       if (match1) gl->CompleteStep = TSGLLECompleteStep_Rescale;
1149       else if (match2) gl->CompleteStep = TSGLLECompleteStep_RescaleAndModify;
1150       else SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "%s", completef);
1151     }
1152     {
1153       char type[256] = TSGLLEACCEPT_ALWAYS;
1154       PetscCall(PetscOptionsFList("-ts_gl_accept_type", "Method to use for determining whether to accept a step", "TSGLLESetAcceptType", TSGLLEAcceptList, gl->accept_name[0] ? gl->accept_name : type, type, sizeof(type), &flg));
1155       if (flg || !gl->accept_name[0]) PetscCall(TSGLLESetAcceptType(ts, type));
1156     }
1157     {
1158       TSGLLEAdapt adapt;
1159       PetscCall(TSGLLEGetAdapt(ts, &adapt));
1160       PetscCall(TSGLLEAdaptSetFromOptions(adapt, PetscOptionsObject));
1161     }
1162   }
1163   PetscOptionsHeadEnd();
1164   PetscFunctionReturn(PETSC_SUCCESS);
1165 }
1166 
1167 static PetscErrorCode TSView_GLLE(TS ts, PetscViewer viewer)
1168 {
1169   TS_GLLE  *gl = (TS_GLLE *)ts->data;
1170   PetscInt  i;
1171   PetscBool iascii, details;
1172 
1173   PetscFunctionBegin;
1174   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1175   if (iascii) {
1176     PetscCall(PetscViewerASCIIPrintf(viewer, "  min order %" PetscInt_FMT ", max order %" PetscInt_FMT ", current order %" PetscInt_FMT "\n", gl->min_order, gl->max_order, gl->schemes[gl->current_scheme]->p));
1177     PetscCall(PetscViewerASCIIPrintf(viewer, "  Error estimation: %s\n", TSGLLEErrorDirections[gl->error_direction]));
1178     PetscCall(PetscViewerASCIIPrintf(viewer, "  Extrapolation: %s\n", gl->extrapolate ? "yes" : "no"));
1179     PetscCall(PetscViewerASCIIPrintf(viewer, "  Acceptance test: %s\n", gl->accept_name[0] ? gl->accept_name : "(not yet set)"));
1180     PetscCall(PetscViewerASCIIPushTab(viewer));
1181     PetscCall(TSGLLEAdaptView(gl->adapt, viewer));
1182     PetscCall(PetscViewerASCIIPopTab(viewer));
1183     PetscCall(PetscViewerASCIIPrintf(viewer, "  type: %s\n", gl->type_name[0] ? gl->type_name : "(not yet set)"));
1184     PetscCall(PetscViewerASCIIPrintf(viewer, "Schemes within family (%" PetscInt_FMT "):\n", gl->nschemes));
1185     details = PETSC_FALSE;
1186     PetscCall(PetscOptionsGetBool(((PetscObject)ts)->options, ((PetscObject)ts)->prefix, "-ts_gl_view_detailed", &details, NULL));
1187     PetscCall(PetscViewerASCIIPushTab(viewer));
1188     for (i = 0; i < gl->nschemes; i++) PetscCall(TSGLLESchemeView(gl->schemes[i], details, viewer));
1189     if (gl->View) PetscCall((*gl->View)(gl, viewer));
1190     PetscCall(PetscViewerASCIIPopTab(viewer));
1191   }
1192   PetscFunctionReturn(PETSC_SUCCESS);
1193 }
1194 
1195 /*@C
1196   TSGLLERegister -  adds a `TSGLLE` implementation
1197 
1198   Not Collective
1199 
1200   Input Parameters:
1201 + sname    - name of user-defined general linear scheme
1202 - function - routine to create method context
1203 
1204   Level: advanced
1205 
1206   Note:
1207   `TSGLLERegister()` may be called multiple times to add several user-defined families.
1208 
1209   Example Usage:
1210 .vb
1211   TSGLLERegister("my_scheme", MySchemeCreate);
1212 .ve
1213 
1214   Then, your scheme can be chosen with the procedural interface via
1215 $ TSGLLESetType(ts, "my_scheme")
1216   or at runtime via the option
1217 $ -ts_gl_type my_scheme
1218 
1219 .seealso: [](ch_ts), `TSGLLE`, `TSGLLEType`, `TSGLLERegisterAll()`
1220 @*/
1221 PetscErrorCode TSGLLERegister(const char sname[], PetscErrorCode (*function)(TS))
1222 {
1223   PetscFunctionBegin;
1224   PetscCall(TSGLLEInitializePackage());
1225   PetscCall(PetscFunctionListAdd(&TSGLLEList, sname, function));
1226   PetscFunctionReturn(PETSC_SUCCESS);
1227 }
1228 
1229 /*@C
1230   TSGLLEAcceptRegister -  adds a `TSGLLE` acceptance scheme
1231 
1232   Not Collective
1233 
1234   Input Parameters:
1235 + sname    - name of user-defined acceptance scheme
1236 - function - routine to create method context, see `TSGLLEAcceptFn` for the calling sequence
1237 
1238   Level: advanced
1239 
1240   Note:
1241   `TSGLLEAcceptRegister()` may be called multiple times to add several user-defined families.
1242 
1243   Example Usage:
1244 .vb
1245   TSGLLEAcceptRegister("my_scheme", MySchemeCreate);
1246 .ve
1247 
1248   Then, your scheme can be chosen with the procedural interface via
1249 .vb
1250   TSGLLESetAcceptType(ts, "my_scheme")
1251 .ve
1252   or at runtime via the option `-ts_gl_accept_type my_scheme`
1253 
1254 .seealso: [](ch_ts), `TSGLLE`, `TSGLLEType`, `TSGLLERegisterAll()`, `TSGLLEAcceptFn`
1255 @*/
1256 PetscErrorCode TSGLLEAcceptRegister(const char sname[], TSGLLEAcceptFn *function)
1257 {
1258   PetscFunctionBegin;
1259   PetscCall(PetscFunctionListAdd(&TSGLLEAcceptList, sname, function));
1260   PetscFunctionReturn(PETSC_SUCCESS);
1261 }
1262 
1263 /*@C
1264   TSGLLERegisterAll - Registers all of the general linear methods in `TSGLLE`
1265 
1266   Not Collective
1267 
1268   Level: advanced
1269 
1270 .seealso: [](ch_ts), `TSGLLE`, `TSGLLERegisterDestroy()`
1271 @*/
1272 PetscErrorCode TSGLLERegisterAll(void)
1273 {
1274   PetscFunctionBegin;
1275   if (TSGLLERegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
1276   TSGLLERegisterAllCalled = PETSC_TRUE;
1277 
1278   PetscCall(TSGLLERegister(TSGLLE_IRKS, TSGLLECreate_IRKS));
1279   PetscCall(TSGLLEAcceptRegister(TSGLLEACCEPT_ALWAYS, TSGLLEAccept_Always));
1280   PetscFunctionReturn(PETSC_SUCCESS);
1281 }
1282 
1283 /*@C
1284   TSGLLEInitializePackage - This function initializes everything in the `TSGLLE` package. It is called
1285   from `TSInitializePackage()`.
1286 
1287   Level: developer
1288 
1289 .seealso: [](ch_ts), `PetscInitialize()`, `TSInitializePackage()`, `TSGLLEFinalizePackage()`
1290 @*/
1291 PetscErrorCode TSGLLEInitializePackage(void)
1292 {
1293   PetscFunctionBegin;
1294   if (TSGLLEPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
1295   TSGLLEPackageInitialized = PETSC_TRUE;
1296   PetscCall(TSGLLERegisterAll());
1297   PetscCall(PetscRegisterFinalize(TSGLLEFinalizePackage));
1298   PetscFunctionReturn(PETSC_SUCCESS);
1299 }
1300 
1301 /*@C
1302   TSGLLEFinalizePackage - This function destroys everything in the `TSGLLE` package. It is
1303   called from `PetscFinalize()`.
1304 
1305   Level: developer
1306 
1307 .seealso: [](ch_ts), `PetscFinalize()`, `TSGLLEInitializePackage()`, `TSInitializePackage()`
1308 @*/
1309 PetscErrorCode TSGLLEFinalizePackage(void)
1310 {
1311   PetscFunctionBegin;
1312   PetscCall(PetscFunctionListDestroy(&TSGLLEList));
1313   PetscCall(PetscFunctionListDestroy(&TSGLLEAcceptList));
1314   TSGLLEPackageInitialized = PETSC_FALSE;
1315   TSGLLERegisterAllCalled  = PETSC_FALSE;
1316   PetscFunctionReturn(PETSC_SUCCESS);
1317 }
1318 
1319 /* ------------------------------------------------------------ */
1320 /*MC
1321   TSGLLE - DAE solver using implicit General Linear methods {cite}`butcher_2007` {cite}`butcher2016numerical`
1322 
1323   Options Database Keys:
1324 +  -ts_gl_type <type> - the class of general linear method (irks)
1325 .  -ts_gl_rtol <tol>  - relative error
1326 .  -ts_gl_atol <tol>  - absolute error
1327 .  -ts_gl_min_order <p> - minimum order method to consider (default=1)
1328 .  -ts_gl_max_order <p> - maximum order method to consider (default=3)
1329 .  -ts_gl_start_order <p> - order of starting method (default=1)
1330 .  -ts_gl_complete <method> - method to use for completing the step (rescale-and-modify or rescale)
1331 -  -ts_adapt_type <method> - adaptive controller to use (none step both)
1332 
1333   Level: beginner
1334 
1335   Notes:
1336   These methods contain Runge-Kutta and multistep schemes as special cases. These special cases
1337   have some fundamental limitations. For example, diagonally implicit Runge-Kutta cannot have
1338   stage order greater than 1 which limits their applicability to very stiff systems.
1339   Meanwhile, multistep methods cannot be A-stable for order greater than 2 and BDF are not
1340   0-stable for order greater than 6. GL methods can be A- and L-stable with arbitrarily high
1341   stage order and reliable error estimates for both 1 and 2 orders higher to facilitate
1342   adaptive step sizes and adaptive order schemes. All this is possible while preserving a
1343   singly diagonally implicit structure.
1344 
1345   This integrator can be applied to DAE.
1346 
1347   Diagonally implicit general linear (DIGL) methods are a generalization of diagonally implicit
1348   Runge-Kutta (DIRK). They are represented by the tableau
1349 
1350 .vb
1351   A  |  U
1352   -------
1353   B  |  V
1354 .ve
1355 
1356   combined with a vector c of abscissa. "Diagonally implicit" means that $A$ is lower
1357   triangular. A step of the general method reads
1358 
1359   $$
1360   \begin{align*}
1361   [ Y ] = [A  U] [  Y'   ] \\
1362   [X^k] = [B  V] [X^{k-1}]
1363   \end{align*}
1364   $$
1365 
1366   where Y is the multivector of stage values, $Y'$ is the multivector of stage derivatives, $X^k$
1367   is the Nordsieck vector of the solution at step $k$. The Nordsieck vector consists of the first
1368   $r$ moments of the solution, given by
1369 
1370   $$
1371   X = [x_0,x_1,...,x_{r-1}] = [x, h x', h^2 x'', ..., h^{r-1} x^{(r-1)} ]
1372   $$
1373 
1374   If $A$ is lower triangular, we can solve the stages $(Y, Y')$ sequentially
1375 
1376   $$
1377   y_i = h \sum_{j=0}^{s-1} (a_{ij} y'_j) + \sum_{j=0}^{r-1} u_{ij} x_j, \, \,   i=0,...,{s-1}
1378   $$
1379 
1380   and then construct the pieces to carry to the next step
1381 
1382   $$
1383   xx_i = h \sum_{j=0}^{s-1} b_{ij} y'_j  + \sum_{j=0}^{r-1} v_{ij} x_j,  \, \,  i=0,...,{r-1}
1384   $$
1385 
1386   Note that when the equations are cast in implicit form, we are using the stage equation to
1387   define $y'_i$ in terms of $y_i$ and known stuff ($y_j$ for $j<i$ and $x_j$ for all $j$).
1388 
1389   Error estimation
1390 
1391   At present, the most attractive GL methods for stiff problems are singly diagonally implicit
1392   schemes which posses Inherent Runge-Kutta Stability (`TSIRKS`).  These methods have $r=s$, the
1393   number of items passed between steps is equal to the number of stages.  The order and
1394   stage-order are one less than the number of stages.  We use the error estimates in the 2007
1395   paper which provide the following estimates
1396 
1397   $$
1398   \begin{align*}
1399   h^{p+1} X^{(p+1)}          = \phi_0^T Y' + [0 \psi_0^T] Xold \\
1400   h^{p+2} X^{(p+2)}          = \phi_1^T Y' + [0 \psi_1^T] Xold \\
1401   h^{p+2} (dx'/dx) X^{(p+1)} = \phi_2^T Y' + [0 \psi_2^T] Xold
1402   \end{align*}
1403   $$
1404 
1405   These estimates are accurate to $ O(h^{p+3})$.
1406 
1407   Changing the step size
1408 
1409   Uses the generalized "rescale and modify" scheme, see equation (4.5) of {cite}`butcher_2007`.
1410 
1411 .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSType`
1412 M*/
1413 PETSC_EXTERN PetscErrorCode TSCreate_GLLE(TS ts)
1414 {
1415   TS_GLLE *gl;
1416 
1417   PetscFunctionBegin;
1418   PetscCall(TSGLLEInitializePackage());
1419 
1420   PetscCall(PetscNew(&gl));
1421   ts->data = (void *)gl;
1422 
1423   ts->ops->reset          = TSReset_GLLE;
1424   ts->ops->destroy        = TSDestroy_GLLE;
1425   ts->ops->view           = TSView_GLLE;
1426   ts->ops->setup          = TSSetUp_GLLE;
1427   ts->ops->solve          = TSSolve_GLLE;
1428   ts->ops->setfromoptions = TSSetFromOptions_GLLE;
1429   ts->ops->snesfunction   = SNESTSFormFunction_GLLE;
1430   ts->ops->snesjacobian   = SNESTSFormJacobian_GLLE;
1431 
1432   ts->usessnes = PETSC_TRUE;
1433 
1434   gl->max_step_rejections = 1;
1435   gl->min_order           = 1;
1436   gl->max_order           = 3;
1437   gl->start_order         = 1;
1438   gl->current_scheme      = -1;
1439   gl->extrapolate         = PETSC_FALSE;
1440 
1441   gl->wrms_atol = 1e-8;
1442   gl->wrms_rtol = 1e-5;
1443 
1444   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLESetType_C", &TSGLLESetType_GLLE));
1445   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLESetAcceptType_C", &TSGLLESetAcceptType_GLLE));
1446   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSGLLEGetAdapt_C", &TSGLLEGetAdapt_GLLE));
1447   PetscFunctionReturn(PETSC_SUCCESS);
1448 }
1449