xref: /petsc/src/ts/impls/eimex/eimex.c (revision a336c15037c72f93cd561f5a5e11e93175f2efd9)
1 #include <petsc/private/tsimpl.h> /*I   "petscts.h"   I*/
2 #include <petscdm.h>
3 
4 static const PetscInt TSEIMEXDefault = 3;
5 
6 typedef struct {
7   PetscInt     row_ind;    /* Return the term T[row_ind][col_ind] */
8   PetscInt     col_ind;    /* Return the term T[row_ind][col_ind] */
9   PetscInt     nstages;    /* Numbers of stages in current scheme */
10   PetscInt     max_rows;   /* Maximum number of rows */
11   PetscInt    *N;          /* Harmonic sequence N[max_rows] */
12   Vec          Y;          /* States computed during the step, used to complete the step */
13   Vec          Z;          /* For shift*(Y-Z) */
14   Vec         *T;          /* Working table, size determined by nstages */
15   Vec          YdotRHS;    /* g(x) Work vector holding YdotRHS during residual evaluation */
16   Vec          YdotI;      /* xdot-f(x) Work vector holding YdotI = F(t,x,xdot) when xdot =0 */
17   Vec          Ydot;       /* f(x)+g(x) Work vector */
18   Vec          VecSolPrev; /* Work vector holding the solution from the previous step (used for interpolation) */
19   PetscReal    shift;
20   PetscReal    ctime;
21   PetscBool    recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */
22   PetscBool    ord_adapt;          /* order adapativity */
23   TSStepStatus status;
24 } TS_EIMEX;
25 
26 /* This function is pure */
27 static PetscInt Map(PetscInt i, PetscInt j, PetscInt s)
28 {
29   return (2 * s - j + 1) * j / 2 + i - j;
30 }
31 
32 static PetscErrorCode TSEvaluateStep_EIMEX(TS ts, PetscInt order, Vec X, PetscBool *done)
33 {
34   TS_EIMEX      *ext = (TS_EIMEX *)ts->data;
35   const PetscInt ns  = ext->nstages;
36 
37   PetscFunctionBegin;
38   PetscCall(VecCopy(ext->T[Map(ext->row_ind, ext->col_ind, ns)], X));
39   PetscFunctionReturn(PETSC_SUCCESS);
40 }
41 
42 static PetscErrorCode TSStage_EIMEX(TS ts, PetscInt istage)
43 {
44   TS_EIMEX *ext = (TS_EIMEX *)ts->data;
45   PetscReal h;
46   Vec       Y = ext->Y, Z = ext->Z;
47   SNES      snes;
48   TSAdapt   adapt;
49   PetscInt  i, its, lits;
50   PetscBool accept;
51 
52   PetscFunctionBegin;
53   PetscCall(TSGetSNES(ts, &snes));
54   h          = ts->time_step / ext->N[istage]; /* step size for the istage-th stage */
55   ext->shift = 1. / h;
56   PetscCall(SNESSetLagJacobian(snes, -2)); /* Recompute the Jacobian on this solve, but not again */
57   PetscCall(VecCopy(ext->VecSolPrev, Y));  /* Take the previous solution as initial step */
58 
59   for (i = 0; i < ext->N[istage]; i++) {
60     ext->ctime = ts->ptime + h * i;
61     PetscCall(VecCopy(Y, Z)); /* Save the solution of the previous substep */
62     PetscCall(SNESSolve(snes, NULL, Y));
63     PetscCall(SNESGetIterationNumber(snes, &its));
64     PetscCall(SNESGetLinearSolveIterations(snes, &lits));
65     ts->snes_its += its;
66     ts->ksp_its += lits;
67     PetscCall(TSGetAdapt(ts, &adapt));
68     PetscCall(TSAdaptCheckStage(adapt, ts, ext->ctime, Y, &accept));
69   }
70   PetscFunctionReturn(PETSC_SUCCESS);
71 }
72 
73 static PetscErrorCode TSStep_EIMEX(TS ts)
74 {
75   TS_EIMEX      *ext = (TS_EIMEX *)ts->data;
76   const PetscInt ns  = ext->nstages;
77   Vec           *T = ext->T, Y = ext->Y;
78   SNES           snes;
79   PetscInt       i, j;
80   PetscBool      accept = PETSC_FALSE;
81   PetscReal      alpha, local_error, local_error_a, local_error_r;
82 
83   PetscFunctionBegin;
84   PetscCall(TSGetSNES(ts, &snes));
85   PetscCall(SNESSetType(snes, "ksponly"));
86   ext->status = TS_STEP_INCOMPLETE;
87 
88   PetscCall(VecCopy(ts->vec_sol, ext->VecSolPrev));
89 
90   /* Apply n_j steps of the base method to obtain solutions of T(j,1),1<=j<=s */
91   for (j = 0; j < ns; j++) {
92     PetscCall(TSStage_EIMEX(ts, j));
93     PetscCall(VecCopy(Y, T[j]));
94   }
95 
96   for (i = 1; i < ns; i++) {
97     for (j = i; j < ns; j++) {
98       alpha = -(PetscReal)ext->N[j] / ext->N[j - i];
99       PetscCall(VecAXPBYPCZ(T[Map(j, i, ns)], alpha, 1.0, 0, T[Map(j, i - 1, ns)], T[Map(j - 1, i - 1, ns)])); /* T[j][i]=alpha*T[j][i-1]+T[j-1][i-1] */
100       alpha = 1.0 / (1.0 + alpha);
101       PetscCall(VecScale(T[Map(j, i, ns)], alpha));
102     }
103   }
104 
105   PetscCall(TSEvaluateStep(ts, ns, ts->vec_sol, NULL)); /*update ts solution */
106 
107   if (ext->ord_adapt && ext->nstages < ext->max_rows) {
108     accept = PETSC_FALSE;
109     while (!accept && ext->nstages < ext->max_rows) {
110       PetscCall(TSErrorWeightedNorm(ts, ts->vec_sol, T[Map(ext->nstages - 1, ext->nstages - 2, ext->nstages)], ts->adapt->wnormtype, &local_error, &local_error_a, &local_error_r));
111       accept = (local_error < 1.0) ? PETSC_TRUE : PETSC_FALSE;
112 
113       if (!accept) { /* add one more stage*/
114         PetscCall(TSStage_EIMEX(ts, ext->nstages));
115         ext->nstages++;
116         ext->row_ind++;
117         ext->col_ind++;
118         /*T table need to be recycled*/
119         PetscCall(VecDuplicateVecs(ts->vec_sol, (1 + ext->nstages) * ext->nstages / 2, &ext->T));
120         for (i = 0; i < ext->nstages - 1; i++) {
121           for (j = 0; j <= i; j++) PetscCall(VecCopy(T[Map(i, j, ext->nstages - 1)], ext->T[Map(i, j, ext->nstages)]));
122         }
123         PetscCall(VecDestroyVecs(ext->nstages * (ext->nstages - 1) / 2, &T));
124         T = ext->T; /*reset the pointer*/
125         /*recycling finished, store the new solution*/
126         PetscCall(VecCopy(Y, T[ext->nstages - 1]));
127         /*extrapolation for the newly added stage*/
128         for (i = 1; i < ext->nstages; i++) {
129           alpha = -(PetscReal)ext->N[ext->nstages - 1] / ext->N[ext->nstages - 1 - i];
130           PetscCall(VecAXPBYPCZ(T[Map(ext->nstages - 1, i, ext->nstages)], alpha, 1.0, 0, T[Map(ext->nstages - 1, i - 1, ext->nstages)], T[Map(ext->nstages - 1 - 1, i - 1, ext->nstages)])); /*T[ext->nstages-1][i]=alpha*T[ext->nstages-1][i-1]+T[ext->nstages-1-1][i-1]*/
131           alpha = 1.0 / (1.0 + alpha);
132           PetscCall(VecScale(T[Map(ext->nstages - 1, i, ext->nstages)], alpha));
133         }
134         /*update ts solution */
135         PetscCall(TSEvaluateStep(ts, ext->nstages, ts->vec_sol, NULL));
136       } /*end if !accept*/
137     } /*end while*/
138 
139     if (ext->nstages == ext->max_rows) PetscCall(PetscInfo(ts, "Max number of rows has been used\n"));
140   } /*end if ext->ord_adapt*/
141   ts->ptime += ts->time_step;
142   ext->status = TS_STEP_COMPLETE;
143 
144   if (ext->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
145   PetscFunctionReturn(PETSC_SUCCESS);
146 }
147 
148 /* cubic Hermit spline */
149 static PetscErrorCode TSInterpolate_EIMEX(TS ts, PetscReal itime, Vec X)
150 {
151   TS_EIMEX       *ext = (TS_EIMEX *)ts->data;
152   PetscReal       t, a, b;
153   Vec             Y0 = ext->VecSolPrev, Y1 = ext->Y, Ydot = ext->Ydot, YdotI = ext->YdotI;
154   const PetscReal h = ts->ptime - ts->ptime_prev;
155 
156   PetscFunctionBegin;
157   t = (itime - ts->ptime + h) / h;
158   /* YdotI = -f(x)-g(x) */
159 
160   PetscCall(VecZeroEntries(Ydot));
161   PetscCall(TSComputeIFunction(ts, ts->ptime - h, Y0, Ydot, YdotI, PETSC_FALSE));
162 
163   a = 2.0 * t * t * t - 3.0 * t * t + 1.0;
164   b = -(t * t * t - 2.0 * t * t + t) * h;
165   PetscCall(VecAXPBYPCZ(X, a, b, 0.0, Y0, YdotI));
166 
167   PetscCall(TSComputeIFunction(ts, ts->ptime, Y1, Ydot, YdotI, PETSC_FALSE));
168   a = -2.0 * t * t * t + 3.0 * t * t;
169   b = -(t * t * t - t * t) * h;
170   PetscCall(VecAXPBYPCZ(X, a, b, 1.0, Y1, YdotI));
171   PetscFunctionReturn(PETSC_SUCCESS);
172 }
173 
174 static PetscErrorCode TSReset_EIMEX(TS ts)
175 {
176   TS_EIMEX *ext = (TS_EIMEX *)ts->data;
177   PetscInt  ns;
178 
179   PetscFunctionBegin;
180   ns = ext->nstages;
181   PetscCall(VecDestroyVecs((1 + ns) * ns / 2, &ext->T));
182   PetscCall(VecDestroy(&ext->Y));
183   PetscCall(VecDestroy(&ext->Z));
184   PetscCall(VecDestroy(&ext->YdotRHS));
185   PetscCall(VecDestroy(&ext->YdotI));
186   PetscCall(VecDestroy(&ext->Ydot));
187   PetscCall(VecDestroy(&ext->VecSolPrev));
188   PetscCall(PetscFree(ext->N));
189   PetscFunctionReturn(PETSC_SUCCESS);
190 }
191 
192 static PetscErrorCode TSDestroy_EIMEX(TS ts)
193 {
194   PetscFunctionBegin;
195   PetscCall(TSReset_EIMEX(ts));
196   PetscCall(PetscFree(ts->data));
197   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSEIMEXSetMaxRows_C", NULL));
198   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSEIMEXSetRowCol_C", NULL));
199   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSEIMEXSetOrdAdapt_C", NULL));
200   PetscFunctionReturn(PETSC_SUCCESS);
201 }
202 
203 static PetscErrorCode TSEIMEXGetVecs(TS ts, DM dm, Vec *Z, Vec *Ydot, Vec *YdotI, Vec *YdotRHS)
204 {
205   TS_EIMEX *ext = (TS_EIMEX *)ts->data;
206 
207   PetscFunctionBegin;
208   if (Z) {
209     if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSEIMEX_Z", Z));
210     else *Z = ext->Z;
211   }
212   if (Ydot) {
213     if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSEIMEX_Ydot", Ydot));
214     else *Ydot = ext->Ydot;
215   }
216   if (YdotI) {
217     if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSEIMEX_YdotI", YdotI));
218     else *YdotI = ext->YdotI;
219   }
220   if (YdotRHS) {
221     if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSEIMEX_YdotRHS", YdotRHS));
222     else *YdotRHS = ext->YdotRHS;
223   }
224   PetscFunctionReturn(PETSC_SUCCESS);
225 }
226 
227 static PetscErrorCode TSEIMEXRestoreVecs(TS ts, DM dm, Vec *Z, Vec *Ydot, Vec *YdotI, Vec *YdotRHS)
228 {
229   PetscFunctionBegin;
230   if (Z) {
231     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSEIMEX_Z", Z));
232   }
233   if (Ydot) {
234     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSEIMEX_Ydot", Ydot));
235   }
236   if (YdotI) {
237     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSEIMEX_YdotI", YdotI));
238   }
239   if (YdotRHS) {
240     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSEIMEX_YdotRHS", YdotRHS));
241   }
242   PetscFunctionReturn(PETSC_SUCCESS);
243 }
244 
245 /*
246   This defines the nonlinear equation that is to be solved with SNES
247   Fn[t0+Theta*dt, U, (U-U0)*shift] = 0
248   In the case of Backward Euler, Fn = (U-U0)/h-g(t1,U))
249   Since FormIFunction calculates G = ydot - g(t,y), ydot will be set to (U-U0)/h
250 */
251 static PetscErrorCode SNESTSFormFunction_EIMEX(SNES snes, Vec X, Vec G, TS ts)
252 {
253   TS_EIMEX *ext = (TS_EIMEX *)ts->data;
254   Vec       Ydot, Z;
255   DM        dm, dmsave;
256 
257   PetscFunctionBegin;
258   PetscCall(VecZeroEntries(G));
259 
260   PetscCall(SNESGetDM(snes, &dm));
261   PetscCall(TSEIMEXGetVecs(ts, dm, &Z, &Ydot, NULL, NULL));
262   PetscCall(VecZeroEntries(Ydot));
263   dmsave = ts->dm;
264   ts->dm = dm;
265   PetscCall(TSComputeIFunction(ts, ext->ctime, X, Ydot, G, PETSC_FALSE));
266   /* PETSC_FALSE indicates non-imex, adding explicit RHS to the implicit I function.  */
267   PetscCall(VecCopy(G, Ydot));
268   ts->dm = dmsave;
269   PetscCall(TSEIMEXRestoreVecs(ts, dm, &Z, &Ydot, NULL, NULL));
270   PetscFunctionReturn(PETSC_SUCCESS);
271 }
272 
273 /*
274  This defined the Jacobian matrix for SNES. Jn = (I/h-g'(t,y))
275  */
276 static PetscErrorCode SNESTSFormJacobian_EIMEX(SNES snes, Vec X, Mat A, Mat B, TS ts)
277 {
278   TS_EIMEX *ext = (TS_EIMEX *)ts->data;
279   Vec       Ydot;
280   DM        dm, dmsave;
281 
282   PetscFunctionBegin;
283   PetscCall(SNESGetDM(snes, &dm));
284   PetscCall(TSEIMEXGetVecs(ts, dm, NULL, &Ydot, NULL, NULL));
285   /*  PetscCall(VecZeroEntries(Ydot)); */
286   /* ext->Ydot have already been computed in SNESTSFormFunction_EIMEX (SNES guarantees this) */
287   dmsave = ts->dm;
288   ts->dm = dm;
289   PetscCall(TSComputeIJacobian(ts, ts->ptime, X, Ydot, ext->shift, A, B, PETSC_TRUE));
290   ts->dm = dmsave;
291   PetscCall(TSEIMEXRestoreVecs(ts, dm, NULL, &Ydot, NULL, NULL));
292   PetscFunctionReturn(PETSC_SUCCESS);
293 }
294 
295 static PetscErrorCode DMCoarsenHook_TSEIMEX(DM fine, DM coarse, PetscCtx ctx)
296 {
297   PetscFunctionBegin;
298   PetscFunctionReturn(PETSC_SUCCESS);
299 }
300 
301 static PetscErrorCode DMRestrictHook_TSEIMEX(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, PetscCtx ctx)
302 {
303   TS  ts = (TS)ctx;
304   Vec Z, Z_c;
305 
306   PetscFunctionBegin;
307   PetscCall(TSEIMEXGetVecs(ts, fine, &Z, NULL, NULL, NULL));
308   PetscCall(TSEIMEXGetVecs(ts, coarse, &Z_c, NULL, NULL, NULL));
309   PetscCall(MatRestrict(restrct, Z, Z_c));
310   PetscCall(VecPointwiseMult(Z_c, rscale, Z_c));
311   PetscCall(TSEIMEXRestoreVecs(ts, fine, &Z, NULL, NULL, NULL));
312   PetscCall(TSEIMEXRestoreVecs(ts, coarse, &Z_c, NULL, NULL, NULL));
313   PetscFunctionReturn(PETSC_SUCCESS);
314 }
315 
316 static PetscErrorCode TSSetUp_EIMEX(TS ts)
317 {
318   TS_EIMEX *ext = (TS_EIMEX *)ts->data;
319   DM        dm;
320 
321   PetscFunctionBegin;
322   if (!ext->N) { /* ext->max_rows not set */
323     PetscCall(TSEIMEXSetMaxRows(ts, TSEIMEXDefault));
324   }
325   if (-1 == ext->row_ind && -1 == ext->col_ind) {
326     PetscCall(TSEIMEXSetRowCol(ts, ext->max_rows, ext->max_rows));
327   } else { /* ext->row_ind and col_ind already set */
328     if (ext->ord_adapt) PetscCall(PetscInfo(ts, "Order adaptivity is enabled and TSEIMEXSetRowCol or -ts_eimex_row_col option will take no effect\n"));
329   }
330 
331   if (ext->ord_adapt) {
332     ext->nstages = 2; /* Start with the 2-stage scheme */
333     PetscCall(TSEIMEXSetRowCol(ts, ext->nstages, ext->nstages));
334   } else {
335     ext->nstages = ext->max_rows; /* by default nstages is the same as max_rows, this can be changed by setting order adaptivity */
336   }
337 
338   PetscCall(TSGetAdapt(ts, &ts->adapt));
339 
340   PetscCall(VecDuplicateVecs(ts->vec_sol, (1 + ext->nstages) * ext->nstages / 2, &ext->T)); /* full T table */
341   PetscCall(VecDuplicate(ts->vec_sol, &ext->YdotI));
342   PetscCall(VecDuplicate(ts->vec_sol, &ext->YdotRHS));
343   PetscCall(VecDuplicate(ts->vec_sol, &ext->Ydot));
344   PetscCall(VecDuplicate(ts->vec_sol, &ext->VecSolPrev));
345   PetscCall(VecDuplicate(ts->vec_sol, &ext->Y));
346   PetscCall(VecDuplicate(ts->vec_sol, &ext->Z));
347   PetscCall(TSGetDM(ts, &dm));
348   if (dm) PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSEIMEX, DMRestrictHook_TSEIMEX, ts));
349   PetscFunctionReturn(PETSC_SUCCESS);
350 }
351 
352 static PetscErrorCode TSSetFromOptions_EIMEX(TS ts, PetscOptionItems PetscOptionsObject)
353 {
354   TS_EIMEX *ext = (TS_EIMEX *)ts->data;
355   PetscInt  tindex[2];
356   PetscInt  np = 2, nrows = TSEIMEXDefault;
357 
358   PetscFunctionBegin;
359   tindex[0] = TSEIMEXDefault;
360   tindex[1] = TSEIMEXDefault;
361   PetscOptionsHeadBegin(PetscOptionsObject, "EIMEX ODE solver options");
362   {
363     PetscBool flg;
364     PetscCall(PetscOptionsInt("-ts_eimex_max_rows", "Define the maximum number of rows used", "TSEIMEXSetMaxRows", nrows, &nrows, &flg)); /* default value 3 */
365     if (flg) PetscCall(TSEIMEXSetMaxRows(ts, nrows));
366     PetscCall(PetscOptionsIntArray("-ts_eimex_row_col", "Return the specific term in the T table", "TSEIMEXSetRowCol", tindex, &np, &flg));
367     if (flg) PetscCall(TSEIMEXSetRowCol(ts, tindex[0], tindex[1]));
368     PetscCall(PetscOptionsBool("-ts_eimex_order_adapt", "Solve the problem with adaptive order", "TSEIMEXSetOrdAdapt", ext->ord_adapt, &ext->ord_adapt, NULL));
369   }
370   PetscOptionsHeadEnd();
371   PetscFunctionReturn(PETSC_SUCCESS);
372 }
373 
374 static PetscErrorCode TSView_EIMEX(TS ts, PetscViewer viewer)
375 {
376   PetscFunctionBegin;
377   PetscFunctionReturn(PETSC_SUCCESS);
378 }
379 
380 /*@
381   TSEIMEXSetMaxRows - Set the maximum number of rows for `TSEIMEX` schemes
382 
383   Logically Collective
384 
385   Input Parameters:
386 + ts    - timestepping context
387 - nrows - maximum number of rows
388 
389   Level: intermediate
390 
391 .seealso: [](ch_ts), `TSEIMEXSetRowCol()`, `TSEIMEXSetOrdAdapt()`, `TSEIMEX`
392 @*/
393 PetscErrorCode TSEIMEXSetMaxRows(TS ts, PetscInt nrows)
394 {
395   PetscFunctionBegin;
396   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
397   PetscTryMethod(ts, "TSEIMEXSetMaxRows_C", (TS, PetscInt), (ts, nrows));
398   PetscFunctionReturn(PETSC_SUCCESS);
399 }
400 
401 /*@
402   TSEIMEXSetRowCol - Set the number of rows and the number of columns for the tableau that represents the T solution in the `TSEIMEX` scheme
403 
404   Logically Collective
405 
406   Input Parameters:
407 + ts  - timestepping context
408 . row - the row
409 - col - the column
410 
411   Level: intermediate
412 
413 .seealso: [](ch_ts), `TSEIMEXSetMaxRows()`, `TSEIMEXSetOrdAdapt()`, `TSEIMEX`
414 @*/
415 PetscErrorCode TSEIMEXSetRowCol(TS ts, PetscInt row, PetscInt col)
416 {
417   PetscFunctionBegin;
418   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
419   PetscTryMethod(ts, "TSEIMEXSetRowCol_C", (TS, PetscInt, PetscInt), (ts, row, col));
420   PetscFunctionReturn(PETSC_SUCCESS);
421 }
422 
423 /*@
424   TSEIMEXSetOrdAdapt - Set the order adaptativity for the `TSEIMEX` schemes
425 
426   Logically Collective
427 
428   Input Parameters:
429 + ts  - timestepping context
430 - flg - index in the T table
431 
432   Level: intermediate
433 
434 .seealso: [](ch_ts), `TSEIMEXSetRowCol()`, `TSEIMEX`
435 @*/
436 PetscErrorCode TSEIMEXSetOrdAdapt(TS ts, PetscBool flg)
437 {
438   PetscFunctionBegin;
439   PetscValidHeaderSpecific(ts, TS_CLASSID, 1);
440   PetscTryMethod(ts, "TSEIMEXSetOrdAdapt_C", (TS, PetscBool), (ts, flg));
441   PetscFunctionReturn(PETSC_SUCCESS);
442 }
443 
444 static PetscErrorCode TSEIMEXSetMaxRows_EIMEX(TS ts, PetscInt nrows)
445 {
446   TS_EIMEX *ext = (TS_EIMEX *)ts->data;
447   PetscInt  i;
448 
449   PetscFunctionBegin;
450   PetscCheck(nrows >= 0 && nrows <= 100, ((PetscObject)ts)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Max number of rows (current value %" PetscInt_FMT ") should be an integer number between 1 and 100", nrows);
451   PetscCall(PetscFree(ext->N));
452   ext->max_rows = nrows;
453   PetscCall(PetscMalloc1(nrows, &ext->N));
454   for (i = 0; i < nrows; i++) ext->N[i] = i + 1;
455   PetscFunctionReturn(PETSC_SUCCESS);
456 }
457 
458 static PetscErrorCode TSEIMEXSetRowCol_EIMEX(TS ts, PetscInt row, PetscInt col)
459 {
460   TS_EIMEX *ext = (TS_EIMEX *)ts->data;
461 
462   PetscFunctionBegin;
463   PetscCheck(row >= 1 && col >= 1, ((PetscObject)ts)->comm, PETSC_ERR_ARG_OUTOFRANGE, "The row or column index (current value %" PetscInt_FMT ",%" PetscInt_FMT ") should not be less than 1 ", row, col);
464   PetscCheck(row <= ext->max_rows && col <= ext->max_rows, ((PetscObject)ts)->comm, PETSC_ERR_ARG_OUTOFRANGE, "The row or column index (current value %" PetscInt_FMT ",%" PetscInt_FMT ") exceeds the maximum number of rows %" PetscInt_FMT, row, col,
465              ext->max_rows);
466   PetscCheck(col <= row, ((PetscObject)ts)->comm, PETSC_ERR_ARG_OUTOFRANGE, "The column index (%" PetscInt_FMT ") exceeds the row index (%" PetscInt_FMT ")", col, row);
467 
468   ext->row_ind = row - 1;
469   ext->col_ind = col - 1; /* Array index in C starts from 0 */
470   PetscFunctionReturn(PETSC_SUCCESS);
471 }
472 
473 static PetscErrorCode TSEIMEXSetOrdAdapt_EIMEX(TS ts, PetscBool flg)
474 {
475   TS_EIMEX *ext = (TS_EIMEX *)ts->data;
476 
477   PetscFunctionBegin;
478   ext->ord_adapt = flg;
479   PetscFunctionReturn(PETSC_SUCCESS);
480 }
481 
482 /*MC
483    TSEIMEX - Time stepping with Extrapolated W-IMEX methods {cite}`constantinescu_a2010a`.
484 
485    These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly nonlinear such that it
486    is expensive to solve with a fully implicit method. The user should provide the stiff part of the equation using `TSSetIFunction()` and the
487    non-stiff part with `TSSetRHSFunction()`.
488 
489       Level: beginner
490 
491   Notes:
492   The default is a 3-stage scheme, it can be changed with `TSEIMEXSetMaxRows()` or -ts_eimex_max_rows
493 
494   This method currently only works with ODEs, for which the stiff part $ F(t,X,Xdot) $  has the form $ Xdot + Fhat(t,X)$.
495 
496   The general system is written as
497 
498   $$
499   F(t,X,Xdot) = G(t,X)
500   $$
501 
502   where F represents the stiff part and G represents the non-stiff part. The user should provide the stiff part
503   of the equation using TSSetIFunction() and the non-stiff part with `TSSetRHSFunction()`.
504   This method is designed to be linearly implicit on G and can use an approximate and lagged Jacobian.
505 
506   Another common form for the system is
507 
508   $$
509   y'=f(x)+g(x)
510   $$
511 
512   The relationship between F,G and f,g is
513 
514   $$
515   F = y'-f(x), G = g(x)
516   $$
517 
518 .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSEIMEXSetMaxRows()`, `TSEIMEXSetRowCol()`, `TSEIMEXSetOrdAdapt()`, `TSType`
519  M*/
520 PETSC_EXTERN PetscErrorCode TSCreate_EIMEX(TS ts)
521 {
522   TS_EIMEX *ext;
523 
524   PetscFunctionBegin;
525   ts->ops->reset          = TSReset_EIMEX;
526   ts->ops->destroy        = TSDestroy_EIMEX;
527   ts->ops->view           = TSView_EIMEX;
528   ts->ops->setup          = TSSetUp_EIMEX;
529   ts->ops->step           = TSStep_EIMEX;
530   ts->ops->interpolate    = TSInterpolate_EIMEX;
531   ts->ops->evaluatestep   = TSEvaluateStep_EIMEX;
532   ts->ops->setfromoptions = TSSetFromOptions_EIMEX;
533   ts->ops->snesfunction   = SNESTSFormFunction_EIMEX;
534   ts->ops->snesjacobian   = SNESTSFormJacobian_EIMEX;
535   ts->default_adapt_type  = TSADAPTNONE;
536 
537   ts->usessnes = PETSC_TRUE;
538 
539   PetscCall(PetscNew(&ext));
540   ts->data = (void *)ext;
541 
542   ext->ord_adapt = PETSC_FALSE; /* By default, no order adapativity */
543   ext->row_ind   = -1;
544   ext->col_ind   = -1;
545   ext->max_rows  = TSEIMEXDefault;
546   ext->nstages   = TSEIMEXDefault;
547 
548   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSEIMEXSetMaxRows_C", TSEIMEXSetMaxRows_EIMEX));
549   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSEIMEXSetRowCol_C", TSEIMEXSetRowCol_EIMEX));
550   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSEIMEXSetOrdAdapt_C", TSEIMEXSetOrdAdapt_EIMEX));
551   PetscFunctionReturn(PETSC_SUCCESS);
552 }
553