xref: /petsc/src/ts/impls/eimex/eimex.c (revision 679f3ef1d3d88f2e9fcf07e91e863c07d2617150)
17d5697caSHong Zhang 
2af0996ceSBarry Smith #include <petsc/private/tsimpl.h>                /*I   "petscts.h"   I*/
34042b796SJed Brown #include <petscdm.h>
47d5697caSHong Zhang 
57d5697caSHong Zhang static const PetscInt TSEIMEXDefault = 3;
67d5697caSHong Zhang 
77d5697caSHong Zhang typedef struct {
87d5697caSHong Zhang   PetscInt     row_ind;         /* Return the term T[row_ind][col_ind] */
97d5697caSHong Zhang   PetscInt     col_ind;         /* Return the term T[row_ind][col_ind] */
107d5697caSHong Zhang   PetscInt     nstages;         /* Numbers of stages in current scheme */
117d5697caSHong Zhang   PetscInt     max_rows;        /* Maximum number of rows */
127d5697caSHong Zhang   PetscInt     *N;              /* Harmonic sequence N[max_rows] */
137d5697caSHong Zhang   Vec          Y;               /* States computed during the step, used to complete the step */
147d5697caSHong Zhang   Vec          Z;               /* For shift*(Y-Z) */
157d5697caSHong Zhang   Vec          *T;              /* Working table, size determined by nstages */
167d5697caSHong Zhang   Vec          YdotRHS;         /* f(x) Work vector holding YdotRHS during residual evaluation */
177d5697caSHong Zhang   Vec          YdotI;           /* xdot-g(x) Work vector holding YdotI = G(t,x,xdot) when xdot =0 */
187d5697caSHong Zhang   Vec          Ydot;            /* f(x)+g(x) Work vector */
197d5697caSHong Zhang   Vec          VecSolPrev;      /* Work vector holding the solution from the previous step (used for interpolation) */
207d5697caSHong Zhang   PetscReal    shift;
217d5697caSHong Zhang   PetscReal    ctime;
227d5697caSHong Zhang   PetscBool    recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */
237d5697caSHong Zhang   PetscBool    ord_adapt;       /* order adapativity */
247d5697caSHong Zhang   TSStepStatus status;
257d5697caSHong Zhang } TS_EIMEX;
267d5697caSHong Zhang 
277d5697caSHong Zhang /* This function is pure */
287d5697caSHong Zhang static PetscInt Map(PetscInt i, PetscInt j, PetscInt s)
297d5697caSHong Zhang {
307d5697caSHong Zhang   return ((2*s-j+1)*j/2+i-j);
317d5697caSHong Zhang }
327d5697caSHong Zhang 
337d5697caSHong Zhang 
347d5697caSHong Zhang #undef __FUNCT__
357d5697caSHong Zhang #define __FUNCT__ "TSEvaluateStep_EIMEX"
367d5697caSHong Zhang static PetscErrorCode TSEvaluateStep_EIMEX(TS ts,PetscInt order,Vec X,PetscBool *done)
377d5697caSHong Zhang {
387d5697caSHong Zhang   TS_EIMEX        *ext = (TS_EIMEX*)ts->data;
397d5697caSHong Zhang   const PetscInt  ns = ext->nstages;
407d5697caSHong Zhang   PetscErrorCode ierr;
417d5697caSHong Zhang   PetscFunctionBegin;
427d5697caSHong Zhang   ierr = VecCopy(ext->T[Map(ext->row_ind,ext->col_ind,ns)],X);CHKERRQ(ierr);
437d5697caSHong Zhang   PetscFunctionReturn(0);
447d5697caSHong Zhang }
457d5697caSHong Zhang 
467d5697caSHong Zhang 
477d5697caSHong Zhang #undef __FUNCT__
487d5697caSHong Zhang #define __FUNCT__ "TSStage_EIMEX"
497d5697caSHong Zhang static PetscErrorCode TSStage_EIMEX(TS ts,PetscInt istage)
507d5697caSHong Zhang {
517d5697caSHong Zhang   TS_EIMEX        *ext = (TS_EIMEX*)ts->data;
527d5697caSHong Zhang   PetscReal       h;
537d5697caSHong Zhang   Vec             Y=ext->Y, Z=ext->Z;
547d5697caSHong Zhang   SNES            snes;
557d5697caSHong Zhang   TSAdapt         adapt;
567d5697caSHong Zhang   PetscInt        i,its,lits;
577d5697caSHong Zhang   PetscBool       accept;
587d5697caSHong Zhang   PetscErrorCode  ierr;
597d5697caSHong Zhang 
607d5697caSHong Zhang   PetscFunctionBegin;
617d5697caSHong Zhang   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
627d5697caSHong Zhang   h = ts->time_step/ext->N[istage];/* step size for the istage-th stage */
637d5697caSHong Zhang   ext->shift = 1./h;
647d5697caSHong Zhang   ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */
657d5697caSHong Zhang   ierr = VecCopy(ext->VecSolPrev,Y);CHKERRQ(ierr); /* Take the previous solution as intial step */
667d5697caSHong Zhang 
677d5697caSHong Zhang   for(i=0; i<ext->N[istage]; i++){
687d5697caSHong Zhang     ext->ctime = ts->ptime + h*i;
697d5697caSHong Zhang     ierr = VecCopy(Y,Z);CHKERRQ(ierr);/* Save the solution of the previous substep */
7016dc5bacSEmil Constantinescu     ierr = SNESSolve(snes,NULL,Y);CHKERRQ(ierr);
717d5697caSHong Zhang     ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
727d5697caSHong Zhang     ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr);
737d5697caSHong Zhang     ts->snes_its += its; ts->ksp_its += lits;
74552698daSJed Brown     ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr);
75b295832fSPierre Barbier de Reuille     ierr = TSAdaptCheckStage(adapt,ts,ext->ctime,Y,&accept);CHKERRQ(ierr);
767d5697caSHong Zhang   }
777d5697caSHong Zhang 
787d5697caSHong Zhang   PetscFunctionReturn(0);
797d5697caSHong Zhang }
807d5697caSHong Zhang 
817d5697caSHong Zhang 
827d5697caSHong Zhang #undef __FUNCT__
837d5697caSHong Zhang #define __FUNCT__ "TSStep_EIMEX"
847d5697caSHong Zhang static PetscErrorCode TSStep_EIMEX(TS ts)
857d5697caSHong Zhang {
867d5697caSHong Zhang   TS_EIMEX        *ext = (TS_EIMEX*)ts->data;
877d5697caSHong Zhang   const PetscInt  ns = ext->nstages;
887d5697caSHong Zhang   Vec             *T=ext->T, Y=ext->Y;
897d5697caSHong Zhang 
907d5697caSHong Zhang   SNES            snes;
917d5697caSHong Zhang   PetscInt        i,j;
927d5697caSHong Zhang   PetscBool       accept = PETSC_FALSE;
937d5697caSHong Zhang   PetscErrorCode  ierr;
947d5697caSHong Zhang   PetscReal       alpha,local_error;
957d5697caSHong Zhang   PetscFunctionBegin;
967d5697caSHong Zhang 
977d5697caSHong Zhang   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
987d5697caSHong Zhang   ierr = SNESSetType(snes,"ksponly"); CHKERRQ(ierr);
997d5697caSHong Zhang   ext->status = TS_STEP_INCOMPLETE;
1007d5697caSHong Zhang 
1017d5697caSHong Zhang   ierr = VecCopy(ts->vec_sol,ext->VecSolPrev);CHKERRQ(ierr);
1027d5697caSHong Zhang 
1037d5697caSHong Zhang   /* Apply n_j steps of the base method to obtain solutions of T(j,1),1<=j<=s */
1047d5697caSHong Zhang   for(j=0; j<ns; j++){
1057d5697caSHong Zhang         ierr = TSStage_EIMEX(ts,j);CHKERRQ(ierr);
1067d5697caSHong Zhang         ierr = VecCopy(Y,T[j]); CHKERRQ(ierr);
1077d5697caSHong Zhang   }
1087d5697caSHong Zhang 
1097d5697caSHong Zhang   for(i=1;i<ns;i++){
1107d5697caSHong Zhang     for(j=i;j<ns;j++){
1117d5697caSHong Zhang       alpha = -(PetscReal)ext->N[j]/ext->N[j-i];
112302440fdSBarry Smith       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);
1137d5697caSHong Zhang       alpha = 1.0/(1.0 + alpha);
1147d5697caSHong Zhang       ierr  = VecScale(T[Map(j,i,ns)],alpha);CHKERRQ(ierr);
1157d5697caSHong Zhang     }
1167d5697caSHong Zhang   }
1177d5697caSHong Zhang 
11816dc5bacSEmil Constantinescu   ierr = TSEvaluateStep(ts,ns,ts->vec_sol,NULL);CHKERRQ(ierr);/*update ts solution */
1197d5697caSHong Zhang 
1207d5697caSHong Zhang   if(ext->ord_adapt && ext->nstages < ext->max_rows){
1217d5697caSHong Zhang 	accept = PETSC_FALSE;
1227d5697caSHong Zhang 	while(!accept && ext->nstages < ext->max_rows){
123a4868fbcSLisandro Dalcin 	  ierr = TSErrorWeightedNorm(ts,ts->vec_sol,T[Map(ext->nstages-1,ext->nstages-2,ext->nstages)],ts->adapt->wnormtype,&local_error);CHKERRQ(ierr);
1247d5697caSHong Zhang 	  accept = (local_error < 1.0)? PETSC_TRUE : PETSC_FALSE;
1257d5697caSHong Zhang 
1267d5697caSHong Zhang 	  if(!accept){/* add one more stage*/
1277d5697caSHong Zhang 	    ierr = TSStage_EIMEX(ts,ext->nstages);CHKERRQ(ierr);
1287d5697caSHong Zhang 	    ext->nstages++; ext->row_ind++; ext->col_ind++;
1297d5697caSHong Zhang 	    /*T table need to be recycled*/
1307d5697caSHong Zhang 	    ierr = VecDuplicateVecs(ts->vec_sol,(1+ext->nstages)*ext->nstages/2,&ext->T);CHKERRQ(ierr);
1317d5697caSHong Zhang 	    for(i=0; i<ext->nstages-1; i++){
1327d5697caSHong Zhang 	      for(j=0; j<=i; j++){
1337d5697caSHong Zhang 	        ierr = VecCopy(T[Map(i,j,ext->nstages-1)],ext->T[Map(i,j,ext->nstages)]);CHKERRQ(ierr);
1347d5697caSHong Zhang 	      }
1357d5697caSHong Zhang 	    }
1367d5697caSHong Zhang 	    ierr = VecDestroyVecs(ext->nstages*(ext->nstages-1)/2,&T);CHKERRQ(ierr);
1377d5697caSHong Zhang 	    T = ext->T; /*reset the pointer*/
1387d5697caSHong Zhang 	    /*recycling finished, store the new solution*/
1397d5697caSHong Zhang 	    ierr = VecCopy(Y,T[ext->nstages-1]); CHKERRQ(ierr);
1407d5697caSHong Zhang 	    /*extrapolation for the newly added stage*/
1417d5697caSHong Zhang 	    for(i=1;i<ext->nstages;i++){
1427d5697caSHong Zhang 	      alpha = -(PetscReal)ext->N[ext->nstages-1]/ext->N[ext->nstages-1-i];
143302440fdSBarry Smith 	      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);
1447d5697caSHong Zhang 	      alpha = 1.0/(1.0 + alpha);
1457d5697caSHong Zhang 	      ierr  = VecScale(T[Map(ext->nstages-1,i,ext->nstages)],alpha);CHKERRQ(ierr);
1467d5697caSHong Zhang 	    }
1477d5697caSHong Zhang 	    /*update ts solution */
14816dc5bacSEmil Constantinescu 	    ierr = TSEvaluateStep(ts,ext->nstages,ts->vec_sol,NULL);CHKERRQ(ierr);
1497d5697caSHong Zhang 	  }/*end if !accept*/
1507d5697caSHong Zhang 	}/*end while*/
1517d5697caSHong Zhang 
1527d5697caSHong Zhang 	if(ext->nstages == ext->max_rows){
1537d5697caSHong Zhang 		ierr = PetscInfo(ts,"Max number of rows has been used\n");CHKERRQ(ierr);
1547d5697caSHong Zhang 	}
1557d5697caSHong Zhang   }/*end if ext->ord_adapt*/
1567d5697caSHong Zhang   ts->ptime += ts->time_step;
1577d5697caSHong Zhang   ext->status = TS_STEP_COMPLETE;
1587d5697caSHong Zhang 
1597d5697caSHong Zhang   if (ext->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED;
1607d5697caSHong Zhang   PetscFunctionReturn(0);
1617d5697caSHong Zhang }
1627d5697caSHong Zhang 
1637d5697caSHong Zhang /* cubic Hermit spline */
1647d5697caSHong Zhang #undef __FUNCT__
1657d5697caSHong Zhang #define __FUNCT__ "TSInterpolate_EIMEX"
1667d5697caSHong Zhang static PetscErrorCode TSInterpolate_EIMEX(TS ts,PetscReal itime,Vec X)
1677d5697caSHong Zhang {
1687d5697caSHong Zhang   TS_EIMEX       *ext = (TS_EIMEX*)ts->data;
1697d5697caSHong Zhang   PetscReal      t,a,b;
170be5899b3SLisandro Dalcin   Vec            Y0=ext->VecSolPrev,Y1=ext->Y,Ydot=ext->Ydot,YdotI=ext->YdotI;
171be5899b3SLisandro Dalcin   const PetscReal h = ts->ptime - ts->ptime_prev;
1727d5697caSHong Zhang   PetscErrorCode ierr;
1737d5697caSHong Zhang   PetscFunctionBegin;
1747d5697caSHong Zhang   t = (itime -ts->ptime + h)/h;
1757d5697caSHong Zhang   /* YdotI = -f(x)-g(x) */
1767d5697caSHong Zhang 
1777d5697caSHong Zhang   ierr = VecZeroEntries(Ydot);CHKERRQ(ierr);
1787d5697caSHong Zhang   ierr = TSComputeIFunction(ts,ts->ptime-h,Y0,Ydot,YdotI,PETSC_FALSE);CHKERRQ(ierr);
1797d5697caSHong Zhang 
18016dc5bacSEmil Constantinescu   a    = 2.0*t*t*t - 3.0*t*t + 1.0;
18116dc5bacSEmil Constantinescu   b    = -(t*t*t - 2.0*t*t + t)*h;
1827d5697caSHong Zhang   ierr = VecAXPBYPCZ(X,a,b,0.0,Y0,YdotI);CHKERRQ(ierr);
1837d5697caSHong Zhang 
1847d5697caSHong Zhang   ierr = TSComputeIFunction(ts,ts->ptime,Y1,Ydot,YdotI,PETSC_FALSE);CHKERRQ(ierr);
18516dc5bacSEmil Constantinescu   a    = -2.0*t*t*t+3.0*t*t;
18616dc5bacSEmil Constantinescu   b    = -(t*t*t - t*t)*h;
1877d5697caSHong Zhang   ierr = VecAXPBYPCZ(X,a,b,1.0,Y1,YdotI);CHKERRQ(ierr);
1887d5697caSHong Zhang 
1897d5697caSHong Zhang   PetscFunctionReturn(0);
1907d5697caSHong Zhang }
1917d5697caSHong Zhang 
1927d5697caSHong Zhang 
1937d5697caSHong Zhang #undef __FUNCT__
1947d5697caSHong Zhang #define __FUNCT__ "TSReset_EIMEX"
1957d5697caSHong Zhang static PetscErrorCode TSReset_EIMEX(TS ts)
1967d5697caSHong Zhang {
1977d5697caSHong Zhang   TS_EIMEX        *ext = (TS_EIMEX*)ts->data;
1987d5697caSHong Zhang   PetscInt        ns;
1997d5697caSHong Zhang   PetscErrorCode  ierr;
2007d5697caSHong Zhang 
2017d5697caSHong Zhang   PetscFunctionBegin;
2027d5697caSHong Zhang   ns = ext->nstages;
2037d5697caSHong Zhang   ierr = VecDestroyVecs((1+ns)*ns/2,&ext->T);CHKERRQ(ierr);
2047d5697caSHong Zhang   ierr = VecDestroy(&ext->Y);CHKERRQ(ierr);
2057d5697caSHong Zhang   ierr = VecDestroy(&ext->Z);CHKERRQ(ierr);
2067d5697caSHong Zhang   ierr = VecDestroy(&ext->YdotRHS);CHKERRQ(ierr);
2077d5697caSHong Zhang   ierr = VecDestroy(&ext->YdotI);CHKERRQ(ierr);
2087d5697caSHong Zhang   ierr = VecDestroy(&ext->Ydot);CHKERRQ(ierr);
2097d5697caSHong Zhang   ierr = VecDestroy(&ext->VecSolPrev);CHKERRQ(ierr);
2107d5697caSHong Zhang   ierr = PetscFree(ext->N);CHKERRQ(ierr);
2117d5697caSHong Zhang   PetscFunctionReturn(0);
2127d5697caSHong Zhang }
2137d5697caSHong Zhang 
2147d5697caSHong Zhang #undef __FUNCT__
2157d5697caSHong Zhang #define __FUNCT__ "TSDestroy_EIMEX"
2167d5697caSHong Zhang static PetscErrorCode TSDestroy_EIMEX(TS ts)
2177d5697caSHong Zhang {
2187d5697caSHong Zhang   PetscErrorCode  ierr;
2197d5697caSHong Zhang 
2207d5697caSHong Zhang   PetscFunctionBegin;
2217d5697caSHong Zhang   ierr = TSReset_EIMEX(ts);CHKERRQ(ierr);
2227d5697caSHong Zhang   ierr = PetscFree(ts->data);CHKERRQ(ierr);
223e1d27e54SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetMaxRows_C",NULL);CHKERRQ(ierr);
224e1d27e54SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetRowCol_C",NULL);CHKERRQ(ierr);
225e1d27e54SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetOrdAdapt_C",NULL);CHKERRQ(ierr);
2267d5697caSHong Zhang 
2277d5697caSHong Zhang   PetscFunctionReturn(0);
2287d5697caSHong Zhang }
2297d5697caSHong Zhang 
2307d5697caSHong Zhang 
2317d5697caSHong Zhang #undef __FUNCT__
2327d5697caSHong Zhang #define __FUNCT__ "TSEIMEXGetVecs"
2337d5697caSHong Zhang static PetscErrorCode TSEIMEXGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydot,Vec *YdotI, Vec *YdotRHS)
2347d5697caSHong Zhang {
2357d5697caSHong Zhang   TS_EIMEX       *ext = (TS_EIMEX*)ts->data;
2367d5697caSHong Zhang   PetscErrorCode ierr;
2377d5697caSHong Zhang 
2387d5697caSHong Zhang   PetscFunctionBegin;
2397d5697caSHong Zhang   if (Z) {
2407d5697caSHong Zhang     if (dm && dm != ts->dm) {
2417d5697caSHong Zhang       ierr = DMGetNamedGlobalVector(dm,"TSEIMEX_Z",Z);CHKERRQ(ierr);
2427d5697caSHong Zhang     } else *Z = ext->Z;
2437d5697caSHong Zhang   }
2447d5697caSHong Zhang   if (Ydot) {
2457d5697caSHong Zhang     if (dm && dm != ts->dm) {
2467d5697caSHong Zhang       ierr = DMGetNamedGlobalVector(dm,"TSEIMEX_Ydot",Ydot);CHKERRQ(ierr);
2477d5697caSHong Zhang     } else *Ydot = ext->Ydot;
2487d5697caSHong Zhang   }
2497d5697caSHong Zhang   if (YdotI) {
2507d5697caSHong Zhang     if (dm && dm != ts->dm) {
2517d5697caSHong Zhang       ierr = DMGetNamedGlobalVector(dm,"TSEIMEX_YdotI",YdotI);CHKERRQ(ierr);
2527d5697caSHong Zhang     } else *YdotI = ext->YdotI;
2537d5697caSHong Zhang   }
2547d5697caSHong Zhang   if (YdotRHS) {
2557d5697caSHong Zhang     if (dm && dm != ts->dm) {
2567d5697caSHong Zhang       ierr = DMGetNamedGlobalVector(dm,"TSEIMEX_YdotRHS",YdotRHS);CHKERRQ(ierr);
2577d5697caSHong Zhang     } else *YdotRHS = ext->YdotRHS;
2587d5697caSHong Zhang   }
2597d5697caSHong Zhang   PetscFunctionReturn(0);
2607d5697caSHong Zhang }
2617d5697caSHong Zhang 
2627d5697caSHong Zhang 
2637d5697caSHong Zhang #undef __FUNCT__
2647d5697caSHong Zhang #define __FUNCT__ "TSEIMEXRestoreVecs"
2657d5697caSHong Zhang static PetscErrorCode TSEIMEXRestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydot,Vec *YdotI,Vec *YdotRHS)
2667d5697caSHong Zhang {
2677d5697caSHong Zhang   PetscErrorCode ierr;
2687d5697caSHong Zhang 
2697d5697caSHong Zhang   PetscFunctionBegin;
2707d5697caSHong Zhang   if (Z) {
2717d5697caSHong Zhang     if (dm && dm != ts->dm) {
2727d5697caSHong Zhang       ierr = DMRestoreNamedGlobalVector(dm,"TSEIMEX_Z",Z);CHKERRQ(ierr);
2737d5697caSHong Zhang     }
2747d5697caSHong Zhang   }
2757d5697caSHong Zhang   if (Ydot) {
2767d5697caSHong Zhang     if (dm && dm != ts->dm) {
2777d5697caSHong Zhang       ierr = DMRestoreNamedGlobalVector(dm,"TSEIMEX_Ydot",Ydot);CHKERRQ(ierr);
2787d5697caSHong Zhang     }
2797d5697caSHong Zhang   }
2807d5697caSHong Zhang   if (YdotI) {
2817d5697caSHong Zhang     if (dm && dm != ts->dm) {
2827d5697caSHong Zhang       ierr = DMRestoreNamedGlobalVector(dm,"TSEIMEX_YdotI",YdotI);CHKERRQ(ierr);
2837d5697caSHong Zhang     }
2847d5697caSHong Zhang   }
2857d5697caSHong Zhang   if (YdotRHS) {
2867d5697caSHong Zhang     if (dm && dm != ts->dm) {
2877d5697caSHong Zhang       ierr = DMRestoreNamedGlobalVector(dm,"TSEIMEX_YdotRHS",YdotRHS);CHKERRQ(ierr);
2887d5697caSHong Zhang     }
2897d5697caSHong Zhang   }
2907d5697caSHong Zhang   PetscFunctionReturn(0);
2917d5697caSHong Zhang }
2927d5697caSHong Zhang 
2937d5697caSHong Zhang 
2947d5697caSHong Zhang /*
2957d5697caSHong Zhang   This defines the nonlinear equation that is to be solved with SNES
2967d5697caSHong Zhang   Fn[t0+Theta*dt, U, (U-U0)*shift] = 0
2977d5697caSHong Zhang   In the case of Backward Euler, Fn = (U-U0)/h-g(t1,U))
2987d5697caSHong Zhang   Since FormIFunction calculates G = ydot - g(t,y), ydot will be set to (U-U0)/h
2997d5697caSHong Zhang */
3007d5697caSHong Zhang #undef __FUNCT__
3017d5697caSHong Zhang #define __FUNCT__ "SNESTSFormFunction_EIMEX"
3027d5697caSHong Zhang static PetscErrorCode SNESTSFormFunction_EIMEX(SNES snes,Vec X,Vec G,TS ts)
3037d5697caSHong Zhang {
3047d5697caSHong Zhang   TS_EIMEX        *ext = (TS_EIMEX*)ts->data;
3057d5697caSHong Zhang   PetscErrorCode  ierr;
3067d5697caSHong Zhang   Vec             Ydot,Z;
3077d5697caSHong Zhang   DM              dm,dmsave;
3087d5697caSHong Zhang 
3097d5697caSHong Zhang   PetscFunctionBegin;
3107d5697caSHong Zhang   ierr = VecZeroEntries(G);CHKERRQ(ierr);
3117d5697caSHong Zhang 
3127d5697caSHong Zhang   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
31316dc5bacSEmil Constantinescu   ierr = TSEIMEXGetVecs(ts,dm,&Z,&Ydot,NULL,NULL);CHKERRQ(ierr);
3147d5697caSHong Zhang   ierr = VecZeroEntries(Ydot);CHKERRQ(ierr);
3157d5697caSHong Zhang   dmsave = ts->dm;
3167d5697caSHong Zhang   ts->dm = dm;
3177d5697caSHong Zhang   ierr = TSComputeIFunction(ts,ext->ctime,X,Ydot,G,PETSC_FALSE);CHKERRQ(ierr);
3187d5697caSHong Zhang   /* PETSC_FALSE indicates non-imex, adding explicit RHS to the implicit I function.  */
3197d5697caSHong Zhang   ierr = VecCopy(G,Ydot);CHKERRQ(ierr);
3207d5697caSHong Zhang   ts->dm = dmsave;
32116dc5bacSEmil Constantinescu   ierr = TSEIMEXRestoreVecs(ts,dm,&Z,&Ydot,NULL,NULL);CHKERRQ(ierr);
3227d5697caSHong Zhang 
3237d5697caSHong Zhang   PetscFunctionReturn(0);
3247d5697caSHong Zhang }
3257d5697caSHong Zhang 
3267d5697caSHong Zhang /*
3277d5697caSHong Zhang  This defined the Jacobian matrix for SNES. Jn = (I/h-g'(t,y))
3287d5697caSHong Zhang  */
3297d5697caSHong Zhang #undef __FUNCT__
3307d5697caSHong Zhang #define __FUNCT__ "SNESTSFormJacobian_EIMEX"
331d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_EIMEX(SNES snes,Vec X,Mat A,Mat B,TS ts)
3327d5697caSHong Zhang {
3337d5697caSHong Zhang   TS_EIMEX        *ext = (TS_EIMEX*)ts->data;
3347d5697caSHong Zhang   Vec             Ydot;
3357d5697caSHong Zhang   PetscErrorCode  ierr;
3367d5697caSHong Zhang   DM              dm,dmsave;
3377d5697caSHong Zhang   PetscFunctionBegin;
3387d5697caSHong Zhang   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
33916dc5bacSEmil Constantinescu   ierr = TSEIMEXGetVecs(ts,dm,NULL,&Ydot,NULL,NULL);CHKERRQ(ierr);
3407d5697caSHong Zhang   /*  ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); */
3417d5697caSHong Zhang   /* ext->Ydot have already been computed in SNESTSFormFunction_EIMEX (SNES guarantees this) */
3427d5697caSHong Zhang   dmsave = ts->dm;
3437d5697caSHong Zhang   ts->dm = dm;
344d1e9a80fSBarry Smith   ierr = TSComputeIJacobian(ts,ts->ptime,X,Ydot,ext->shift,A,B,PETSC_TRUE);CHKERRQ(ierr);
3457d5697caSHong Zhang   ts->dm = dmsave;
34616dc5bacSEmil Constantinescu   ierr = TSEIMEXRestoreVecs(ts,dm,NULL,&Ydot,NULL,NULL);CHKERRQ(ierr);
3477d5697caSHong Zhang   PetscFunctionReturn(0);
3487d5697caSHong Zhang }
3497d5697caSHong Zhang 
3507d5697caSHong Zhang #undef __FUNCT__
3517d5697caSHong Zhang #define __FUNCT__ "DMCoarsenHook_TSEIMEX"
3527d5697caSHong Zhang static PetscErrorCode DMCoarsenHook_TSEIMEX(DM fine,DM coarse,void *ctx)
3537d5697caSHong Zhang {
3547d5697caSHong Zhang 
3557d5697caSHong Zhang   PetscFunctionBegin;
3567d5697caSHong Zhang   PetscFunctionReturn(0);
3577d5697caSHong Zhang }
3587d5697caSHong Zhang 
3597d5697caSHong Zhang #undef __FUNCT__
3607d5697caSHong Zhang #define __FUNCT__ "DMRestrictHook_TSEIMEX"
3617d5697caSHong Zhang static PetscErrorCode DMRestrictHook_TSEIMEX(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx)
3627d5697caSHong Zhang {
3637d5697caSHong Zhang   TS ts = (TS)ctx;
3647d5697caSHong Zhang   PetscErrorCode ierr;
3657d5697caSHong Zhang   Vec Z,Z_c;
3667d5697caSHong Zhang 
3677d5697caSHong Zhang   PetscFunctionBegin;
36816dc5bacSEmil Constantinescu   ierr = TSEIMEXGetVecs(ts,fine,&Z,NULL,NULL,NULL);CHKERRQ(ierr);
36916dc5bacSEmil Constantinescu   ierr = TSEIMEXGetVecs(ts,coarse,&Z_c,NULL,NULL,NULL);CHKERRQ(ierr);
3707d5697caSHong Zhang   ierr = MatRestrict(restrct,Z,Z_c);CHKERRQ(ierr);
3717d5697caSHong Zhang   ierr = VecPointwiseMult(Z_c,rscale,Z_c);CHKERRQ(ierr);
37216dc5bacSEmil Constantinescu   ierr = TSEIMEXRestoreVecs(ts,fine,&Z,NULL,NULL,NULL);CHKERRQ(ierr);
37316dc5bacSEmil Constantinescu   ierr = TSEIMEXRestoreVecs(ts,coarse,&Z_c,NULL,NULL,NULL);CHKERRQ(ierr);
3747d5697caSHong Zhang   PetscFunctionReturn(0);
3757d5697caSHong Zhang }
3767d5697caSHong Zhang 
3777d5697caSHong Zhang 
3787d5697caSHong Zhang #undef __FUNCT__
3797d5697caSHong Zhang #define __FUNCT__ "TSSetUp_EIMEX"
3807d5697caSHong Zhang static PetscErrorCode TSSetUp_EIMEX(TS ts)
3817d5697caSHong Zhang {
3827d5697caSHong Zhang   TS_EIMEX       *ext = (TS_EIMEX*)ts->data;
3837d5697caSHong Zhang   PetscErrorCode ierr;
3847d5697caSHong Zhang   DM             dm;
3857d5697caSHong Zhang 
386955c1f14SBarry Smith   PetscFunctionBegin;
3877d5697caSHong Zhang   if (!ext->N){ /* ext->max_rows not set */
3887d5697caSHong Zhang     ierr = TSEIMEXSetMaxRows(ts,TSEIMEXDefault);CHKERRQ(ierr);
3897d5697caSHong Zhang   }
3907d5697caSHong Zhang   if(-1 == ext->row_ind && -1 == ext->col_ind){
3917d5697caSHong Zhang         ierr = TSEIMEXSetRowCol(ts,ext->max_rows,ext->max_rows);CHKERRQ(ierr);
3927d5697caSHong Zhang   } else{/* ext->row_ind and col_ind already set */
3937d5697caSHong Zhang     if (ext->ord_adapt){
3947d5697caSHong Zhang       ierr = PetscInfo(ts,"Order adaptivity is enabled and TSEIMEXSetRowCol or -ts_eimex_row_col option will take no effect\n");CHKERRQ(ierr);
3957d5697caSHong Zhang     }
3967d5697caSHong Zhang   }
3977d5697caSHong Zhang 
3987d5697caSHong Zhang   if(ext->ord_adapt){
3997d5697caSHong Zhang     ext->nstages = 2; /* Start with the 2-stage scheme */
4007d5697caSHong Zhang     ierr = TSEIMEXSetRowCol(ts,ext->nstages,ext->nstages);CHKERRQ(ierr);
4017d5697caSHong Zhang   } else{
4027d5697caSHong Zhang     ext->nstages = ext->max_rows; /* by default nstages is the same as max_rows, this can be changed by setting order adaptivity */
4037d5697caSHong Zhang   }
4047d5697caSHong Zhang 
4057d5697caSHong Zhang   ierr = VecDuplicateVecs(ts->vec_sol,(1+ext->nstages)*ext->nstages/2,&ext->T);CHKERRQ(ierr);/* full T table */
4067d5697caSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&ext->YdotI);CHKERRQ(ierr);
4077d5697caSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&ext->YdotRHS);CHKERRQ(ierr);
4087d5697caSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&ext->Ydot);CHKERRQ(ierr);
4097d5697caSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&ext->VecSolPrev);CHKERRQ(ierr);
4107d5697caSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&ext->Y);CHKERRQ(ierr);
4117d5697caSHong Zhang   ierr = VecDuplicate(ts->vec_sol,&ext->Z);CHKERRQ(ierr);
4127d5697caSHong Zhang   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
4137d5697caSHong Zhang   if (dm) {
4147d5697caSHong Zhang     ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSEIMEX,DMRestrictHook_TSEIMEX,ts);CHKERRQ(ierr);
4157d5697caSHong Zhang   }
4167d5697caSHong Zhang   PetscFunctionReturn(0);
4177d5697caSHong Zhang }
4187d5697caSHong Zhang 
4197d5697caSHong Zhang #undef __FUNCT__
4207d5697caSHong Zhang #define __FUNCT__ "TSSetFromOptions_EIMEX"
4214416b707SBarry Smith static PetscErrorCode TSSetFromOptions_EIMEX(PetscOptionItems *PetscOptionsObject,TS ts)
4227d5697caSHong Zhang {
4237d5697caSHong Zhang   TS_EIMEX       *ext = (TS_EIMEX*)ts->data;
4247d5697caSHong Zhang   PetscErrorCode ierr;
4259200755eSBarry Smith   PetscInt       tindex[2];
4267d5697caSHong Zhang   PetscInt       np = 2, nrows=TSEIMEXDefault;
4279200755eSBarry Smith 
4287d5697caSHong Zhang   PetscFunctionBegin;
4299200755eSBarry Smith   tindex[0] = TSEIMEXDefault;
4309200755eSBarry Smith   tindex[1] = TSEIMEXDefault;
431e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"EIMEX ODE solver options");CHKERRQ(ierr);
4327d5697caSHong Zhang   {
4337d5697caSHong Zhang     PetscBool flg;
4347d5697caSHong Zhang     ierr = PetscOptionsInt("-ts_eimex_max_rows","Define the maximum number of rows used","TSEIMEXSetMaxRows",nrows,&nrows,&flg);CHKERRQ(ierr); /* default value 3 */
4357d5697caSHong Zhang     if(flg){
4367d5697caSHong Zhang       ierr = TSEIMEXSetMaxRows(ts,nrows);CHKERRQ(ierr);
4377d5697caSHong Zhang     }
4387d5697caSHong Zhang     ierr = PetscOptionsIntArray("-ts_eimex_row_col","Return the specific term in the T table","TSEIMEXSetRowCol",tindex,&np,&flg);CHKERRQ(ierr);
4397d5697caSHong Zhang     if(flg){
4407d5697caSHong Zhang       ierr = TSEIMEXSetRowCol(ts,tindex[0],tindex[1]);CHKERRQ(ierr);
4417d5697caSHong Zhang     }
44216dc5bacSEmil Constantinescu     ierr = PetscOptionsBool("-ts_eimex_order_adapt","Solve the problem with adaptive order","TSEIMEXSetOrdAdapt",ext->ord_adapt,&ext->ord_adapt,NULL);CHKERRQ(ierr);
4437d5697caSHong Zhang   }
4447d5697caSHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
4457d5697caSHong Zhang   PetscFunctionReturn(0);
4467d5697caSHong Zhang }
4477d5697caSHong Zhang 
4487d5697caSHong Zhang #undef __FUNCT__
4497d5697caSHong Zhang #define __FUNCT__ "TSView_EIMEX"
4507d5697caSHong Zhang static PetscErrorCode TSView_EIMEX(TS ts,PetscViewer viewer)
4517d5697caSHong Zhang {
4527d5697caSHong Zhang   /*  TS_EIMEX         *ext = (TS_EIMEX*)ts->data; */
4537d5697caSHong Zhang   PetscBool        iascii;
4547d5697caSHong Zhang   PetscErrorCode   ierr;
4557d5697caSHong Zhang 
4567d5697caSHong Zhang   PetscFunctionBegin;
4577d5697caSHong Zhang   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4587d5697caSHong Zhang   if (iascii) {
4597d5697caSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  EIMEX\n");CHKERRQ(ierr);
4607d5697caSHong Zhang   }
4617d5697caSHong Zhang   ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr);
4627d5697caSHong Zhang   PetscFunctionReturn(0);
4637d5697caSHong Zhang }
4647d5697caSHong Zhang 
4657d5697caSHong Zhang 
4667d5697caSHong Zhang #undef __FUNCT__
4677d5697caSHong Zhang #define __FUNCT__ "TSEIMEXSetMaxRows"
4687d5697caSHong Zhang /*@C
4697d5697caSHong Zhang   TSEIMEXSetMaxRows - Set the maximum number of rows for EIMEX schemes
4707d5697caSHong Zhang 
4717d5697caSHong Zhang   Logically collective
4727d5697caSHong Zhang 
4737d5697caSHong Zhang   Input Parameter:
4747d5697caSHong Zhang +  ts - timestepping context
4757d5697caSHong Zhang -  nrows - maximum number of rows
4767d5697caSHong Zhang 
4777d5697caSHong Zhang   Level: intermediate
4787d5697caSHong Zhang 
4797d5697caSHong Zhang .seealso: TSEIMEXSetRowCol(), TSEIMEXSetOrdAdapt(), TSEIMEX
4807d5697caSHong Zhang @*/
4817d5697caSHong Zhang PetscErrorCode TSEIMEXSetMaxRows(TS ts, PetscInt nrows)
4827d5697caSHong Zhang {
4837d5697caSHong Zhang   PetscErrorCode ierr;
4847d5697caSHong Zhang   PetscFunctionBegin;
4857d5697caSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
4867d5697caSHong Zhang   ierr = PetscTryMethod(ts,"TSEIMEXSetMaxRows_C",(TS,PetscInt),(ts,nrows));CHKERRQ(ierr);
4877d5697caSHong Zhang   PetscFunctionReturn(0);
4887d5697caSHong Zhang }
4897d5697caSHong Zhang 
4907d5697caSHong Zhang 
4917d5697caSHong Zhang #undef __FUNCT__
4927d5697caSHong Zhang #define __FUNCT__ "TSEIMEXSetRowCol"
4937d5697caSHong Zhang /*@C
4947d5697caSHong Zhang   TSEIMEXSetRowCol - Set the type index in the T table for the return value
4957d5697caSHong Zhang 
4967d5697caSHong Zhang   Logically collective
4977d5697caSHong Zhang 
4987d5697caSHong Zhang   Input Parameter:
4997d5697caSHong Zhang +  ts - timestepping context
5007d5697caSHong Zhang -  tindex - index in the T table
5017d5697caSHong Zhang 
5027d5697caSHong Zhang   Level: intermediate
5037d5697caSHong Zhang 
5047d5697caSHong Zhang .seealso: TSEIMEXSetMaxRows(), TSEIMEXSetOrdAdapt(), TSEIMEX
5057d5697caSHong Zhang @*/
5067d5697caSHong Zhang PetscErrorCode TSEIMEXSetRowCol(TS ts, PetscInt row, PetscInt col)
5077d5697caSHong Zhang {
5087d5697caSHong Zhang   PetscErrorCode ierr;
5097d5697caSHong Zhang   PetscFunctionBegin;
5107d5697caSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5117d5697caSHong Zhang   ierr = PetscTryMethod(ts,"TSEIMEXSetRowCol_C",(TS,PetscInt, PetscInt),(ts,row,col));CHKERRQ(ierr);
5127d5697caSHong Zhang   PetscFunctionReturn(0);
5137d5697caSHong Zhang }
5147d5697caSHong Zhang 
5157d5697caSHong Zhang 
5167d5697caSHong Zhang #undef __FUNCT__
5177d5697caSHong Zhang #define __FUNCT__ "TSEIMEXSetOrdAdapt"
5187d5697caSHong Zhang /*@C
5197d5697caSHong Zhang   TSEIMEXSetOrdAdapt - Set the order adaptativity
5207d5697caSHong Zhang 
5217d5697caSHong Zhang   Logically collective
5227d5697caSHong Zhang 
5237d5697caSHong Zhang   Input Parameter:
5247d5697caSHong Zhang +  ts - timestepping context
5257d5697caSHong Zhang -  tindex - index in the T table
5267d5697caSHong Zhang 
5277d5697caSHong Zhang   Level: intermediate
5287d5697caSHong Zhang 
5297d5697caSHong Zhang .seealso: TSEIMEXSetRowCol(), TSEIMEXSetOrdAdapt(), TSEIMEX
5307d5697caSHong Zhang @*/
5317d5697caSHong Zhang PetscErrorCode TSEIMEXSetOrdAdapt(TS ts, PetscBool flg)
5327d5697caSHong Zhang {
5337d5697caSHong Zhang   PetscErrorCode ierr;
5347d5697caSHong Zhang   PetscFunctionBegin;
5357d5697caSHong Zhang   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
5367d5697caSHong Zhang   ierr = PetscTryMethod(ts,"TSEIMEXSetOrdAdapt_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr);
5377d5697caSHong Zhang   PetscFunctionReturn(0);
5387d5697caSHong Zhang }
5397d5697caSHong Zhang 
5407d5697caSHong Zhang 
5417d5697caSHong Zhang #undef __FUNCT__
5427d5697caSHong Zhang #define __FUNCT__ "TSEIMEXSetMaxRows_EIMEX"
5434042b796SJed Brown static PetscErrorCode TSEIMEXSetMaxRows_EIMEX(TS ts,PetscInt nrows)
5447d5697caSHong Zhang {
5457d5697caSHong Zhang   TS_EIMEX *ext = (TS_EIMEX*)ts->data;
5467d5697caSHong Zhang   PetscErrorCode ierr;
5477d5697caSHong Zhang   PetscInt       i;
5487d5697caSHong Zhang 
5497d5697caSHong Zhang   PetscFunctionBegin;
5507d5697caSHong 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);
5517d5697caSHong Zhang   ierr = PetscFree(ext->N);CHKERRQ(ierr);
5527d5697caSHong Zhang   ext->max_rows = nrows;
553785e854fSJed Brown   ierr = PetscMalloc1(nrows,&ext->N);CHKERRQ(ierr);
5547d5697caSHong Zhang   for(i=0;i<nrows;i++) ext->N[i]=i+1;
5557d5697caSHong Zhang   PetscFunctionReturn(0);
5567d5697caSHong Zhang }
5577d5697caSHong Zhang 
5587d5697caSHong Zhang #undef __FUNCT__
5597d5697caSHong Zhang #define __FUNCT__ "TSEIMEXSetRowCol_EIMEX"
5604042b796SJed Brown static PetscErrorCode TSEIMEXSetRowCol_EIMEX(TS ts,PetscInt row,PetscInt col)
5617d5697caSHong Zhang {
5627d5697caSHong Zhang   TS_EIMEX *ext = (TS_EIMEX*)ts->data;
5637d5697caSHong Zhang 
5647d5697caSHong Zhang   PetscFunctionBegin;
5657d5697caSHong 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);
5667d5697caSHong 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);
5677d5697caSHong Zhang   if (col > row) SETERRQ2(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"The column index (%d) exceeds the row index (%d)\n",col,row);
5687d5697caSHong Zhang 
5697d5697caSHong Zhang   ext->row_ind = row - 1;
5707d5697caSHong Zhang   ext->col_ind = col - 1; /* Array index in C starts from 0 */
5717d5697caSHong Zhang   PetscFunctionReturn(0);
5727d5697caSHong Zhang }
5737d5697caSHong Zhang 
5747d5697caSHong Zhang #undef __FUNCT__
5757d5697caSHong Zhang #define __FUNCT__ "TSEIMEXSetOrdAdapt_EIMEX"
5764042b796SJed Brown static PetscErrorCode TSEIMEXSetOrdAdapt_EIMEX(TS ts,PetscBool flg)
5777d5697caSHong Zhang {
5787d5697caSHong Zhang   TS_EIMEX *ext = (TS_EIMEX*)ts->data;
5797d5697caSHong Zhang   PetscFunctionBegin;
5807d5697caSHong Zhang   ext->ord_adapt = flg;
5817d5697caSHong Zhang   PetscFunctionReturn(0);
5827d5697caSHong Zhang }
5837d5697caSHong Zhang 
5847d5697caSHong Zhang /*MC
585*679f3ef1SBarry Smith       TSEIMEX - Time stepping with Extrapolated IMEX methods.
586*679f3ef1SBarry Smith 
587*679f3ef1SBarry Smith    These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly nonlinear such that it
588*679f3ef1SBarry Smith    is expensive to solve with a fully implicit method. The user should provide the stiff part of the equation using TSSetIFunction() and the
589*679f3ef1SBarry Smith    non-stiff part with TSSetRHSFunction().
5907d5697caSHong Zhang 
5917d5697caSHong Zhang    Notes:
5927d5697caSHong Zhang   The default is a 3-stage scheme, it can be changed with TSEIMEXSetMaxRows() or -ts_eimex_max_rows
5937d5697caSHong Zhang 
5947d5697caSHong Zhang   This method currently only works with ODE, for which the stiff part G(t,X,Xdot) has the form Xdot + Ghat(t,X).
5957d5697caSHong Zhang 
596*679f3ef1SBarry Smith   The general system is written as
597*679f3ef1SBarry Smith 
598*679f3ef1SBarry Smith   G(t,X,Xdot) = F(t,X)
599*679f3ef1SBarry Smith 
600*679f3ef1SBarry Smith   where G represents the stiff part and F represents the non-stiff part. The user should provide the stiff part
601*679f3ef1SBarry Smith   of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction().
602*679f3ef1SBarry Smith   This method is designed to be linearly implicit on G and can use an approximate and lagged Jacobian.
603*679f3ef1SBarry Smith 
604*679f3ef1SBarry Smith   Another common form for the system is
605*679f3ef1SBarry Smith 
606*679f3ef1SBarry Smith   y'=f(x)+g(x)
607*679f3ef1SBarry Smith 
608*679f3ef1SBarry Smith   The relationship between F,G and f,g is
609*679f3ef1SBarry Smith 
610*679f3ef1SBarry Smith   G = y'-g(x), F = f(x)
611*679f3ef1SBarry Smith 
612*679f3ef1SBarry Smith  References
613*679f3ef1SBarry Smith   E. Constantinescu and A. Sandu, Extrapolated implicit-explicit time stepping, SIAM Journal on Scientific
614*679f3ef1SBarry Smith Computing, 31 (2010), pp. 4452-4477.
615*679f3ef1SBarry Smith 
6167d5697caSHong Zhang       Level: beginner
6177d5697caSHong Zhang 
618*679f3ef1SBarry Smith .seealso:  TSCreate(), TS, TSSetType(), TSEIMEXSetMaxRows(), TSEIMEXSetRowCol(), TSEIMEXSetOrdAdapt()
619*679f3ef1SBarry Smith 
6207d5697caSHong Zhang  M*/
6217d5697caSHong Zhang #undef __FUNCT__
6227d5697caSHong Zhang #define __FUNCT__ "TSCreate_EIMEX"
6234042b796SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_EIMEX(TS ts)
6247d5697caSHong Zhang {
6257d5697caSHong Zhang   TS_EIMEX       *ext;
6267d5697caSHong Zhang   PetscErrorCode ierr;
6277d5697caSHong Zhang 
6287d5697caSHong Zhang   PetscFunctionBegin;
6297d5697caSHong Zhang 
6307d5697caSHong Zhang   ts->ops->reset          = TSReset_EIMEX;
6317d5697caSHong Zhang   ts->ops->destroy        = TSDestroy_EIMEX;
6327d5697caSHong Zhang   ts->ops->view           = TSView_EIMEX;
6337d5697caSHong Zhang   ts->ops->setup          = TSSetUp_EIMEX;
6347d5697caSHong Zhang   ts->ops->step           = TSStep_EIMEX;
6357d5697caSHong Zhang   ts->ops->interpolate    = TSInterpolate_EIMEX;
6367d5697caSHong Zhang   ts->ops->evaluatestep   = TSEvaluateStep_EIMEX;
6377d5697caSHong Zhang   ts->ops->setfromoptions = TSSetFromOptions_EIMEX;
6387d5697caSHong Zhang   ts->ops->snesfunction   = SNESTSFormFunction_EIMEX;
6397d5697caSHong Zhang   ts->ops->snesjacobian   = SNESTSFormJacobian_EIMEX;
6407d5697caSHong Zhang 
641b00a9115SJed Brown   ierr = PetscNewLog(ts,&ext);CHKERRQ(ierr);
6427d5697caSHong Zhang   ts->data = (void*)ext;
6437d5697caSHong Zhang 
6447d5697caSHong Zhang   ext->ord_adapt = PETSC_FALSE; /* By default, no order adapativity */
6457d5697caSHong Zhang   ext->row_ind   = -1;
6467d5697caSHong Zhang   ext->col_ind   = -1;
6477d5697caSHong Zhang   ext->max_rows  = TSEIMEXDefault;
6487d5697caSHong Zhang   ext->nstages   = TSEIMEXDefault;
6497d5697caSHong Zhang 
650e1d27e54SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetMaxRows_C", TSEIMEXSetMaxRows_EIMEX);CHKERRQ(ierr);
651e1d27e54SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetRowCol_C",  TSEIMEXSetRowCol_EIMEX);CHKERRQ(ierr);
652e1d27e54SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetOrdAdapt_C",TSEIMEXSetOrdAdapt_EIMEX);CHKERRQ(ierr);
6537d5697caSHong Zhang   PetscFunctionReturn(0);
6547d5697caSHong Zhang }
655