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