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