Actual source code: eimex.c

  1: #include <petsc/private/tsimpl.h>
  2: #include <petscdm.h>

  4: static const PetscInt TSEIMEXDefault = 3;

  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;    /* f(x) Work vector holding YdotRHS during residual evaluation */
 16:   Vec          YdotI;      /* xdot-g(x) Work vector holding YdotI = G(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;

 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: }

 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;

 37:   PetscFunctionBegin;
 38:   PetscCall(VecCopy(ext->T[Map(ext->row_ind, ext->col_ind, ns)], X));
 39:   PetscFunctionReturn(PETSC_SUCCESS);
 40: }

 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;

 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 */

 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: }

 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;

 83:   PetscFunctionBegin;
 84:   PetscCall(TSGetSNES(ts, &snes));
 85:   PetscCall(SNESSetType(snes, "ksponly"));
 86:   ext->status = TS_STEP_INCOMPLETE;

 88:   PetscCall(VecCopy(ts->vec_sol, ext->VecSolPrev));

 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:   }

 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:   }

105:   PetscCall(TSEvaluateStep(ts, ns, ts->vec_sol, NULL)); /*update ts solution */

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;

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*/

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;

144:   if (ext->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
145:   PetscFunctionReturn(PETSC_SUCCESS);
146: }

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;

156:   PetscFunctionBegin;
157:   t = (itime - ts->ptime + h) / h;
158:   /* YdotI = -f(x)-g(x) */

160:   PetscCall(VecZeroEntries(Ydot));
161:   PetscCall(TSComputeIFunction(ts, ts->ptime - h, Y0, Ydot, YdotI, PETSC_FALSE));

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));

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: }

174: static PetscErrorCode TSReset_EIMEX(TS ts)
175: {
176:   TS_EIMEX *ext = (TS_EIMEX *)ts->data;
177:   PetscInt  ns;

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: }

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: }

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;

207:   PetscFunctionBegin;
208:   if (Z) {
209:     if (dm && dm != ts->dm) {
210:       PetscCall(DMGetNamedGlobalVector(dm, "TSEIMEX_Z", Z));
211:     } else *Z = ext->Z;
212:   }
213:   if (Ydot) {
214:     if (dm && dm != ts->dm) {
215:       PetscCall(DMGetNamedGlobalVector(dm, "TSEIMEX_Ydot", Ydot));
216:     } else *Ydot = ext->Ydot;
217:   }
218:   if (YdotI) {
219:     if (dm && dm != ts->dm) {
220:       PetscCall(DMGetNamedGlobalVector(dm, "TSEIMEX_YdotI", YdotI));
221:     } else *YdotI = ext->YdotI;
222:   }
223:   if (YdotRHS) {
224:     if (dm && dm != ts->dm) {
225:       PetscCall(DMGetNamedGlobalVector(dm, "TSEIMEX_YdotRHS", YdotRHS));
226:     } else *YdotRHS = ext->YdotRHS;
227:   }
228:   PetscFunctionReturn(PETSC_SUCCESS);
229: }

231: static PetscErrorCode TSEIMEXRestoreVecs(TS ts, DM dm, Vec *Z, Vec *Ydot, Vec *YdotI, Vec *YdotRHS)
232: {
233:   PetscFunctionBegin;
234:   if (Z) {
235:     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSEIMEX_Z", Z));
236:   }
237:   if (Ydot) {
238:     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSEIMEX_Ydot", Ydot));
239:   }
240:   if (YdotI) {
241:     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSEIMEX_YdotI", YdotI));
242:   }
243:   if (YdotRHS) {
244:     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSEIMEX_YdotRHS", YdotRHS));
245:   }
246:   PetscFunctionReturn(PETSC_SUCCESS);
247: }

249: /*
250:   This defines the nonlinear equation that is to be solved with SNES
251:   Fn[t0+Theta*dt, U, (U-U0)*shift] = 0
252:   In the case of Backward Euler, Fn = (U-U0)/h-g(t1,U))
253:   Since FormIFunction calculates G = ydot - g(t,y), ydot will be set to (U-U0)/h
254: */
255: static PetscErrorCode SNESTSFormFunction_EIMEX(SNES snes, Vec X, Vec G, TS ts)
256: {
257:   TS_EIMEX *ext = (TS_EIMEX *)ts->data;
258:   Vec       Ydot, Z;
259:   DM        dm, dmsave;

261:   PetscFunctionBegin;
262:   PetscCall(VecZeroEntries(G));

264:   PetscCall(SNESGetDM(snes, &dm));
265:   PetscCall(TSEIMEXGetVecs(ts, dm, &Z, &Ydot, NULL, NULL));
266:   PetscCall(VecZeroEntries(Ydot));
267:   dmsave = ts->dm;
268:   ts->dm = dm;
269:   PetscCall(TSComputeIFunction(ts, ext->ctime, X, Ydot, G, PETSC_FALSE));
270:   /* PETSC_FALSE indicates non-imex, adding explicit RHS to the implicit I function.  */
271:   PetscCall(VecCopy(G, Ydot));
272:   ts->dm = dmsave;
273:   PetscCall(TSEIMEXRestoreVecs(ts, dm, &Z, &Ydot, NULL, NULL));
274:   PetscFunctionReturn(PETSC_SUCCESS);
275: }

277: /*
278:  This defined the Jacobian matrix for SNES. Jn = (I/h-g'(t,y))
279:  */
280: static PetscErrorCode SNESTSFormJacobian_EIMEX(SNES snes, Vec X, Mat A, Mat B, TS ts)
281: {
282:   TS_EIMEX *ext = (TS_EIMEX *)ts->data;
283:   Vec       Ydot;
284:   DM        dm, dmsave;

286:   PetscFunctionBegin;
287:   PetscCall(SNESGetDM(snes, &dm));
288:   PetscCall(TSEIMEXGetVecs(ts, dm, NULL, &Ydot, NULL, NULL));
289:   /*  PetscCall(VecZeroEntries(Ydot)); */
290:   /* ext->Ydot have already been computed in SNESTSFormFunction_EIMEX (SNES guarantees this) */
291:   dmsave = ts->dm;
292:   ts->dm = dm;
293:   PetscCall(TSComputeIJacobian(ts, ts->ptime, X, Ydot, ext->shift, A, B, PETSC_TRUE));
294:   ts->dm = dmsave;
295:   PetscCall(TSEIMEXRestoreVecs(ts, dm, NULL, &Ydot, NULL, NULL));
296:   PetscFunctionReturn(PETSC_SUCCESS);
297: }

299: static PetscErrorCode DMCoarsenHook_TSEIMEX(DM fine, DM coarse, void *ctx)
300: {
301:   PetscFunctionBegin;
302:   PetscFunctionReturn(PETSC_SUCCESS);
303: }

305: static PetscErrorCode DMRestrictHook_TSEIMEX(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx)
306: {
307:   TS  ts = (TS)ctx;
308:   Vec Z, Z_c;

310:   PetscFunctionBegin;
311:   PetscCall(TSEIMEXGetVecs(ts, fine, &Z, NULL, NULL, NULL));
312:   PetscCall(TSEIMEXGetVecs(ts, coarse, &Z_c, NULL, NULL, NULL));
313:   PetscCall(MatRestrict(restrct, Z, Z_c));
314:   PetscCall(VecPointwiseMult(Z_c, rscale, Z_c));
315:   PetscCall(TSEIMEXRestoreVecs(ts, fine, &Z, NULL, NULL, NULL));
316:   PetscCall(TSEIMEXRestoreVecs(ts, coarse, &Z_c, NULL, NULL, NULL));
317:   PetscFunctionReturn(PETSC_SUCCESS);
318: }

320: static PetscErrorCode TSSetUp_EIMEX(TS ts)
321: {
322:   TS_EIMEX *ext = (TS_EIMEX *)ts->data;
323:   DM        dm;

325:   PetscFunctionBegin;
326:   if (!ext->N) { /* ext->max_rows not set */
327:     PetscCall(TSEIMEXSetMaxRows(ts, TSEIMEXDefault));
328:   }
329:   if (-1 == ext->row_ind && -1 == ext->col_ind) {
330:     PetscCall(TSEIMEXSetRowCol(ts, ext->max_rows, ext->max_rows));
331:   } else { /* ext->row_ind and col_ind already set */
332:     if (ext->ord_adapt) PetscCall(PetscInfo(ts, "Order adaptivity is enabled and TSEIMEXSetRowCol or -ts_eimex_row_col option will take no effect\n"));
333:   }

335:   if (ext->ord_adapt) {
336:     ext->nstages = 2; /* Start with the 2-stage scheme */
337:     PetscCall(TSEIMEXSetRowCol(ts, ext->nstages, ext->nstages));
338:   } else {
339:     ext->nstages = ext->max_rows; /* by default nstages is the same as max_rows, this can be changed by setting order adaptivity */
340:   }

342:   PetscCall(TSGetAdapt(ts, &ts->adapt));

344:   PetscCall(VecDuplicateVecs(ts->vec_sol, (1 + ext->nstages) * ext->nstages / 2, &ext->T)); /* full T table */
345:   PetscCall(VecDuplicate(ts->vec_sol, &ext->YdotI));
346:   PetscCall(VecDuplicate(ts->vec_sol, &ext->YdotRHS));
347:   PetscCall(VecDuplicate(ts->vec_sol, &ext->Ydot));
348:   PetscCall(VecDuplicate(ts->vec_sol, &ext->VecSolPrev));
349:   PetscCall(VecDuplicate(ts->vec_sol, &ext->Y));
350:   PetscCall(VecDuplicate(ts->vec_sol, &ext->Z));
351:   PetscCall(TSGetDM(ts, &dm));
352:   if (dm) PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSEIMEX, DMRestrictHook_TSEIMEX, ts));
353:   PetscFunctionReturn(PETSC_SUCCESS);
354: }

356: static PetscErrorCode TSSetFromOptions_EIMEX(TS ts, PetscOptionItems *PetscOptionsObject)
357: {
358:   TS_EIMEX *ext = (TS_EIMEX *)ts->data;
359:   PetscInt  tindex[2];
360:   PetscInt  np = 2, nrows = TSEIMEXDefault;

362:   PetscFunctionBegin;
363:   tindex[0] = TSEIMEXDefault;
364:   tindex[1] = TSEIMEXDefault;
365:   PetscOptionsHeadBegin(PetscOptionsObject, "EIMEX ODE solver options");
366:   {
367:     PetscBool flg;
368:     PetscCall(PetscOptionsInt("-ts_eimex_max_rows", "Define the maximum number of rows used", "TSEIMEXSetMaxRows", nrows, &nrows, &flg)); /* default value 3 */
369:     if (flg) PetscCall(TSEIMEXSetMaxRows(ts, nrows));
370:     PetscCall(PetscOptionsIntArray("-ts_eimex_row_col", "Return the specific term in the T table", "TSEIMEXSetRowCol", tindex, &np, &flg));
371:     if (flg) PetscCall(TSEIMEXSetRowCol(ts, tindex[0], tindex[1]));
372:     PetscCall(PetscOptionsBool("-ts_eimex_order_adapt", "Solve the problem with adaptive order", "TSEIMEXSetOrdAdapt", ext->ord_adapt, &ext->ord_adapt, NULL));
373:   }
374:   PetscOptionsHeadEnd();
375:   PetscFunctionReturn(PETSC_SUCCESS);
376: }

378: static PetscErrorCode TSView_EIMEX(TS ts, PetscViewer viewer)
379: {
380:   PetscFunctionBegin;
381:   PetscFunctionReturn(PETSC_SUCCESS);
382: }

384: /*@C
385:   TSEIMEXSetMaxRows - Set the maximum number of rows for `TSEIMEX` schemes

387:   Logically Collective

389:   Input Parameters:
390: + ts    - timestepping context
391: - nrows - maximum number of rows

393:   Level: intermediate

395: .seealso: [](ch_ts), `TSEIMEXSetRowCol()`, `TSEIMEXSetOrdAdapt()`, `TSEIMEX`
396: @*/
397: PetscErrorCode TSEIMEXSetMaxRows(TS ts, PetscInt nrows)
398: {
399:   PetscFunctionBegin;
401:   PetscTryMethod(ts, "TSEIMEXSetMaxRows_C", (TS, PetscInt), (ts, nrows));
402:   PetscFunctionReturn(PETSC_SUCCESS);
403: }

405: /*@C
406:   TSEIMEXSetRowCol - Set the number of rows and the number of columns for the tableau that represents the T solution in the `TSEIMEX` scheme

408:   Logically Collective

410:   Input Parameters:
411: + ts  - timestepping context
412: . row - the row
413: - col - the column

415:   Level: intermediate

417: .seealso: [](ch_ts), `TSEIMEXSetMaxRows()`, `TSEIMEXSetOrdAdapt()`, `TSEIMEX`
418: @*/
419: PetscErrorCode TSEIMEXSetRowCol(TS ts, PetscInt row, PetscInt col)
420: {
421:   PetscFunctionBegin;
423:   PetscTryMethod(ts, "TSEIMEXSetRowCol_C", (TS, PetscInt, PetscInt), (ts, row, col));
424:   PetscFunctionReturn(PETSC_SUCCESS);
425: }

427: /*@C
428:   TSEIMEXSetOrdAdapt - Set the order adaptativity for the `TSEIMEX` schemes

430:   Logically Collective

432:   Input Parameters:
433: + ts  - timestepping context
434: - flg - index in the T table

436:   Level: intermediate

438: .seealso: [](ch_ts), `TSEIMEXSetRowCol()`, `TSEIMEX`
439: @*/
440: PetscErrorCode TSEIMEXSetOrdAdapt(TS ts, PetscBool flg)
441: {
442:   PetscFunctionBegin;
444:   PetscTryMethod(ts, "TSEIMEXSetOrdAdapt_C", (TS, PetscBool), (ts, flg));
445:   PetscFunctionReturn(PETSC_SUCCESS);
446: }

448: static PetscErrorCode TSEIMEXSetMaxRows_EIMEX(TS ts, PetscInt nrows)
449: {
450:   TS_EIMEX *ext = (TS_EIMEX *)ts->data;
451:   PetscInt  i;

453:   PetscFunctionBegin;
454:   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);
455:   PetscCall(PetscFree(ext->N));
456:   ext->max_rows = nrows;
457:   PetscCall(PetscMalloc1(nrows, &ext->N));
458:   for (i = 0; i < nrows; i++) ext->N[i] = i + 1;
459:   PetscFunctionReturn(PETSC_SUCCESS);
460: }

462: static PetscErrorCode TSEIMEXSetRowCol_EIMEX(TS ts, PetscInt row, PetscInt col)
463: {
464:   TS_EIMEX *ext = (TS_EIMEX *)ts->data;

466:   PetscFunctionBegin;
467:   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);
468:   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,
469:              ext->max_rows);
470:   PetscCheck(col <= row, ((PetscObject)ts)->comm, PETSC_ERR_ARG_OUTOFRANGE, "The column index (%" PetscInt_FMT ") exceeds the row index (%" PetscInt_FMT ")", col, row);

472:   ext->row_ind = row - 1;
473:   ext->col_ind = col - 1; /* Array index in C starts from 0 */
474:   PetscFunctionReturn(PETSC_SUCCESS);
475: }

477: static PetscErrorCode TSEIMEXSetOrdAdapt_EIMEX(TS ts, PetscBool flg)
478: {
479:   TS_EIMEX *ext = (TS_EIMEX *)ts->data;

481:   PetscFunctionBegin;
482:   ext->ord_adapt = flg;
483:   PetscFunctionReturn(PETSC_SUCCESS);
484: }

486: /*MC
487:    TSEIMEX - Time stepping with Extrapolated IMEX methods {cite}`constantinescu_a2010a`.

489:    These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly nonlinear such that it
490:    is expensive to solve with a fully implicit method. The user should provide the stiff part of the equation using `TSSetIFunction()` and the
491:    non-stiff part with `TSSetRHSFunction()`.

493:       Level: beginner

495:   Notes:
496:   The default is a 3-stage scheme, it can be changed with `TSEIMEXSetMaxRows()` or -ts_eimex_max_rows

498:   This method currently only works with ODE, for which the stiff part $ G(t,X,Xdot) $  has the form $ Xdot + Ghat(t,X)$.

500:   The general system is written as

502:   $$
503:   G(t,X,Xdot) = F(t,X)
504:   $$

506:   where G represents the stiff part and F represents the non-stiff part. The user should provide the stiff part
507:   of the equation using TSSetIFunction() and the non-stiff part with `TSSetRHSFunction()`.
508:   This method is designed to be linearly implicit on G and can use an approximate and lagged Jacobian.

510:   Another common form for the system is

512:   $$
513:   y'=f(x)+g(x)
514:   $$

516:   The relationship between F,G and f,g is

518:   $$
519:   G = y'-g(x), F = f(x)
520:   $$

522: .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSEIMEXSetMaxRows()`, `TSEIMEXSetRowCol()`, `TSEIMEXSetOrdAdapt()`, `TSType`
523:  M*/
524: PETSC_EXTERN PetscErrorCode TSCreate_EIMEX(TS ts)
525: {
526:   TS_EIMEX *ext;

528:   PetscFunctionBegin;
529:   ts->ops->reset          = TSReset_EIMEX;
530:   ts->ops->destroy        = TSDestroy_EIMEX;
531:   ts->ops->view           = TSView_EIMEX;
532:   ts->ops->setup          = TSSetUp_EIMEX;
533:   ts->ops->step           = TSStep_EIMEX;
534:   ts->ops->interpolate    = TSInterpolate_EIMEX;
535:   ts->ops->evaluatestep   = TSEvaluateStep_EIMEX;
536:   ts->ops->setfromoptions = TSSetFromOptions_EIMEX;
537:   ts->ops->snesfunction   = SNESTSFormFunction_EIMEX;
538:   ts->ops->snesjacobian   = SNESTSFormJacobian_EIMEX;
539:   ts->default_adapt_type  = TSADAPTNONE;

541:   ts->usessnes = PETSC_TRUE;

543:   PetscCall(PetscNew(&ext));
544:   ts->data = (void *)ext;

546:   ext->ord_adapt = PETSC_FALSE; /* By default, no order adapativity */
547:   ext->row_ind   = -1;
548:   ext->col_ind   = -1;
549:   ext->max_rows  = TSEIMEXDefault;
550:   ext->nstages   = TSEIMEXDefault;

552:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSEIMEXSetMaxRows_C", TSEIMEXSetMaxRows_EIMEX));
553:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSEIMEXSetRowCol_C", TSEIMEXSetRowCol_EIMEX));
554:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSEIMEXSetOrdAdapt_C", TSEIMEXSetOrdAdapt_EIMEX));
555:   PetscFunctionReturn(PETSC_SUCCESS);
556: }