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