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