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 static PetscErrorCode TSEvaluateStep_EIMEX(TS ts,PetscInt order,Vec X,PetscBool *done) 357d5697caSHong Zhang { 367d5697caSHong Zhang TS_EIMEX *ext = (TS_EIMEX*)ts->data; 377d5697caSHong Zhang const PetscInt ns = ext->nstages; 387d5697caSHong Zhang PetscErrorCode ierr; 397d5697caSHong Zhang PetscFunctionBegin; 407d5697caSHong Zhang ierr = VecCopy(ext->T[Map(ext->row_ind,ext->col_ind,ns)],X);CHKERRQ(ierr); 417d5697caSHong Zhang PetscFunctionReturn(0); 427d5697caSHong Zhang } 437d5697caSHong Zhang 447d5697caSHong Zhang 457d5697caSHong Zhang static PetscErrorCode TSStage_EIMEX(TS ts,PetscInt istage) 467d5697caSHong Zhang { 477d5697caSHong Zhang TS_EIMEX *ext = (TS_EIMEX*)ts->data; 487d5697caSHong Zhang PetscReal h; 497d5697caSHong Zhang Vec Y=ext->Y, Z=ext->Z; 507d5697caSHong Zhang SNES snes; 517d5697caSHong Zhang TSAdapt adapt; 527d5697caSHong Zhang PetscInt i,its,lits; 537d5697caSHong Zhang PetscBool accept; 547d5697caSHong Zhang PetscErrorCode ierr; 557d5697caSHong Zhang 567d5697caSHong Zhang PetscFunctionBegin; 577d5697caSHong Zhang ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 587d5697caSHong Zhang h = ts->time_step/ext->N[istage];/* step size for the istage-th stage */ 597d5697caSHong Zhang ext->shift = 1./h; 607d5697caSHong Zhang ierr = SNESSetLagJacobian(snes,-2);CHKERRQ(ierr); /* Recompute the Jacobian on this solve, but not again */ 617d5697caSHong Zhang ierr = VecCopy(ext->VecSolPrev,Y);CHKERRQ(ierr); /* Take the previous solution as intial step */ 627d5697caSHong Zhang 637d5697caSHong Zhang for(i=0; i<ext->N[istage]; i++){ 647d5697caSHong Zhang ext->ctime = ts->ptime + h*i; 657d5697caSHong Zhang ierr = VecCopy(Y,Z);CHKERRQ(ierr);/* Save the solution of the previous substep */ 6616dc5bacSEmil Constantinescu ierr = SNESSolve(snes,NULL,Y);CHKERRQ(ierr); 677d5697caSHong Zhang ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 687d5697caSHong Zhang ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 697d5697caSHong Zhang ts->snes_its += its; ts->ksp_its += lits; 70552698daSJed Brown ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 71b295832fSPierre Barbier de Reuille ierr = TSAdaptCheckStage(adapt,ts,ext->ctime,Y,&accept);CHKERRQ(ierr); 727d5697caSHong Zhang } 737d5697caSHong Zhang 747d5697caSHong Zhang PetscFunctionReturn(0); 757d5697caSHong Zhang } 767d5697caSHong Zhang 777d5697caSHong Zhang 787d5697caSHong Zhang static PetscErrorCode TSStep_EIMEX(TS ts) 797d5697caSHong Zhang { 807d5697caSHong Zhang TS_EIMEX *ext = (TS_EIMEX*)ts->data; 817d5697caSHong Zhang const PetscInt ns = ext->nstages; 827d5697caSHong Zhang Vec *T=ext->T, Y=ext->Y; 837d5697caSHong Zhang 847d5697caSHong Zhang SNES snes; 857d5697caSHong Zhang PetscInt i,j; 867d5697caSHong Zhang PetscBool accept = PETSC_FALSE; 877d5697caSHong Zhang PetscErrorCode ierr; 887453f775SEmil Constantinescu PetscReal alpha,local_error,local_error_a,local_error_r; 897d5697caSHong Zhang PetscFunctionBegin; 907d5697caSHong Zhang 917d5697caSHong Zhang ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 927d5697caSHong Zhang ierr = SNESSetType(snes,"ksponly"); CHKERRQ(ierr); 937d5697caSHong Zhang ext->status = TS_STEP_INCOMPLETE; 947d5697caSHong Zhang 957d5697caSHong Zhang ierr = VecCopy(ts->vec_sol,ext->VecSolPrev);CHKERRQ(ierr); 967d5697caSHong Zhang 977d5697caSHong Zhang /* Apply n_j steps of the base method to obtain solutions of T(j,1),1<=j<=s */ 987d5697caSHong Zhang for(j=0; j<ns; j++){ 997d5697caSHong Zhang ierr = TSStage_EIMEX(ts,j);CHKERRQ(ierr); 1007d5697caSHong Zhang ierr = VecCopy(Y,T[j]); CHKERRQ(ierr); 1017d5697caSHong Zhang } 1027d5697caSHong Zhang 1037d5697caSHong Zhang for(i=1;i<ns;i++){ 1047d5697caSHong Zhang for(j=i;j<ns;j++){ 1057d5697caSHong Zhang alpha = -(PetscReal)ext->N[j]/ext->N[j-i]; 106302440fdSBarry 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); 1077d5697caSHong Zhang alpha = 1.0/(1.0 + alpha); 1087d5697caSHong Zhang ierr = VecScale(T[Map(j,i,ns)],alpha);CHKERRQ(ierr); 1097d5697caSHong Zhang } 1107d5697caSHong Zhang } 1117d5697caSHong Zhang 11216dc5bacSEmil Constantinescu ierr = TSEvaluateStep(ts,ns,ts->vec_sol,NULL);CHKERRQ(ierr);/*update ts solution */ 1137d5697caSHong Zhang 1147d5697caSHong Zhang if(ext->ord_adapt && ext->nstages < ext->max_rows){ 1157d5697caSHong Zhang accept = PETSC_FALSE; 1167d5697caSHong Zhang while(!accept && ext->nstages < ext->max_rows){ 1177453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,ts->vec_sol,T[Map(ext->nstages-1,ext->nstages-2,ext->nstages)],ts->adapt->wnormtype,&local_error,&local_error_a,&local_error_r);CHKERRQ(ierr); 1187d5697caSHong Zhang accept = (local_error < 1.0)? PETSC_TRUE : PETSC_FALSE; 1197d5697caSHong Zhang 1207d5697caSHong Zhang if(!accept){/* add one more stage*/ 1217d5697caSHong Zhang ierr = TSStage_EIMEX(ts,ext->nstages);CHKERRQ(ierr); 1227d5697caSHong Zhang ext->nstages++; ext->row_ind++; ext->col_ind++; 1237d5697caSHong Zhang /*T table need to be recycled*/ 1247d5697caSHong Zhang ierr = VecDuplicateVecs(ts->vec_sol,(1+ext->nstages)*ext->nstages/2,&ext->T);CHKERRQ(ierr); 1257d5697caSHong Zhang for(i=0; i<ext->nstages-1; i++){ 1267d5697caSHong Zhang for(j=0; j<=i; j++){ 1277d5697caSHong Zhang ierr = VecCopy(T[Map(i,j,ext->nstages-1)],ext->T[Map(i,j,ext->nstages)]);CHKERRQ(ierr); 1287d5697caSHong Zhang } 1297d5697caSHong Zhang } 1307d5697caSHong Zhang ierr = VecDestroyVecs(ext->nstages*(ext->nstages-1)/2,&T);CHKERRQ(ierr); 1317d5697caSHong Zhang T = ext->T; /*reset the pointer*/ 1327d5697caSHong Zhang /*recycling finished, store the new solution*/ 1337d5697caSHong Zhang ierr = VecCopy(Y,T[ext->nstages-1]); CHKERRQ(ierr); 1347d5697caSHong Zhang /*extrapolation for the newly added stage*/ 1357d5697caSHong Zhang for(i=1;i<ext->nstages;i++){ 1367d5697caSHong Zhang alpha = -(PetscReal)ext->N[ext->nstages-1]/ext->N[ext->nstages-1-i]; 137302440fdSBarry 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); 1387d5697caSHong Zhang alpha = 1.0/(1.0 + alpha); 1397d5697caSHong Zhang ierr = VecScale(T[Map(ext->nstages-1,i,ext->nstages)],alpha);CHKERRQ(ierr); 1407d5697caSHong Zhang } 1417d5697caSHong Zhang /*update ts solution */ 14216dc5bacSEmil Constantinescu ierr = TSEvaluateStep(ts,ext->nstages,ts->vec_sol,NULL);CHKERRQ(ierr); 1437d5697caSHong Zhang }/*end if !accept*/ 1447d5697caSHong Zhang }/*end while*/ 1457d5697caSHong Zhang 1467d5697caSHong Zhang if(ext->nstages == ext->max_rows){ 1477d5697caSHong Zhang ierr = PetscInfo(ts,"Max number of rows has been used\n");CHKERRQ(ierr); 1487d5697caSHong Zhang } 1497d5697caSHong Zhang }/*end if ext->ord_adapt*/ 1507d5697caSHong Zhang ts->ptime += ts->time_step; 1517d5697caSHong Zhang ext->status = TS_STEP_COMPLETE; 1527d5697caSHong Zhang 1537d5697caSHong Zhang if (ext->status != TS_STEP_COMPLETE && !ts->reason) ts->reason = TS_DIVERGED_STEP_REJECTED; 1547d5697caSHong Zhang PetscFunctionReturn(0); 1557d5697caSHong Zhang } 1567d5697caSHong Zhang 1577d5697caSHong Zhang /* cubic Hermit spline */ 1587d5697caSHong Zhang static PetscErrorCode TSInterpolate_EIMEX(TS ts,PetscReal itime,Vec X) 1597d5697caSHong Zhang { 1607d5697caSHong Zhang TS_EIMEX *ext = (TS_EIMEX*)ts->data; 1617d5697caSHong Zhang PetscReal t,a,b; 162be5899b3SLisandro Dalcin Vec Y0=ext->VecSolPrev,Y1=ext->Y,Ydot=ext->Ydot,YdotI=ext->YdotI; 163be5899b3SLisandro Dalcin const PetscReal h = ts->ptime - ts->ptime_prev; 1647d5697caSHong Zhang PetscErrorCode ierr; 1657d5697caSHong Zhang PetscFunctionBegin; 1667d5697caSHong Zhang t = (itime -ts->ptime + h)/h; 1677d5697caSHong Zhang /* YdotI = -f(x)-g(x) */ 1687d5697caSHong Zhang 1697d5697caSHong Zhang ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); 1707d5697caSHong Zhang ierr = TSComputeIFunction(ts,ts->ptime-h,Y0,Ydot,YdotI,PETSC_FALSE);CHKERRQ(ierr); 1717d5697caSHong Zhang 17216dc5bacSEmil Constantinescu a = 2.0*t*t*t - 3.0*t*t + 1.0; 17316dc5bacSEmil Constantinescu b = -(t*t*t - 2.0*t*t + t)*h; 1747d5697caSHong Zhang ierr = VecAXPBYPCZ(X,a,b,0.0,Y0,YdotI);CHKERRQ(ierr); 1757d5697caSHong Zhang 1767d5697caSHong Zhang ierr = TSComputeIFunction(ts,ts->ptime,Y1,Ydot,YdotI,PETSC_FALSE);CHKERRQ(ierr); 17716dc5bacSEmil Constantinescu a = -2.0*t*t*t+3.0*t*t; 17816dc5bacSEmil Constantinescu b = -(t*t*t - t*t)*h; 1797d5697caSHong Zhang ierr = VecAXPBYPCZ(X,a,b,1.0,Y1,YdotI);CHKERRQ(ierr); 1807d5697caSHong Zhang 1817d5697caSHong Zhang PetscFunctionReturn(0); 1827d5697caSHong Zhang } 1837d5697caSHong Zhang 1847d5697caSHong Zhang 1857d5697caSHong Zhang static PetscErrorCode TSReset_EIMEX(TS ts) 1867d5697caSHong Zhang { 1877d5697caSHong Zhang TS_EIMEX *ext = (TS_EIMEX*)ts->data; 1887d5697caSHong Zhang PetscInt ns; 1897d5697caSHong Zhang PetscErrorCode ierr; 1907d5697caSHong Zhang 1917d5697caSHong Zhang PetscFunctionBegin; 1927d5697caSHong Zhang ns = ext->nstages; 1937d5697caSHong Zhang ierr = VecDestroyVecs((1+ns)*ns/2,&ext->T);CHKERRQ(ierr); 1947d5697caSHong Zhang ierr = VecDestroy(&ext->Y);CHKERRQ(ierr); 1957d5697caSHong Zhang ierr = VecDestroy(&ext->Z);CHKERRQ(ierr); 1967d5697caSHong Zhang ierr = VecDestroy(&ext->YdotRHS);CHKERRQ(ierr); 1977d5697caSHong Zhang ierr = VecDestroy(&ext->YdotI);CHKERRQ(ierr); 1987d5697caSHong Zhang ierr = VecDestroy(&ext->Ydot);CHKERRQ(ierr); 1997d5697caSHong Zhang ierr = VecDestroy(&ext->VecSolPrev);CHKERRQ(ierr); 2007d5697caSHong Zhang ierr = PetscFree(ext->N);CHKERRQ(ierr); 2017d5697caSHong Zhang PetscFunctionReturn(0); 2027d5697caSHong Zhang } 2037d5697caSHong Zhang 2047d5697caSHong Zhang static PetscErrorCode TSDestroy_EIMEX(TS ts) 2057d5697caSHong Zhang { 2067d5697caSHong Zhang PetscErrorCode ierr; 2077d5697caSHong Zhang 2087d5697caSHong Zhang PetscFunctionBegin; 2097d5697caSHong Zhang ierr = TSReset_EIMEX(ts);CHKERRQ(ierr); 2107d5697caSHong Zhang ierr = PetscFree(ts->data);CHKERRQ(ierr); 211e1d27e54SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetMaxRows_C",NULL);CHKERRQ(ierr); 212e1d27e54SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetRowCol_C",NULL);CHKERRQ(ierr); 213e1d27e54SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetOrdAdapt_C",NULL);CHKERRQ(ierr); 2147d5697caSHong Zhang 2157d5697caSHong Zhang PetscFunctionReturn(0); 2167d5697caSHong Zhang } 2177d5697caSHong Zhang 2187d5697caSHong Zhang 2197d5697caSHong Zhang static PetscErrorCode TSEIMEXGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydot,Vec *YdotI, Vec *YdotRHS) 2207d5697caSHong Zhang { 2217d5697caSHong Zhang TS_EIMEX *ext = (TS_EIMEX*)ts->data; 2227d5697caSHong Zhang PetscErrorCode ierr; 2237d5697caSHong Zhang 2247d5697caSHong Zhang PetscFunctionBegin; 2257d5697caSHong Zhang if (Z) { 2267d5697caSHong Zhang if (dm && dm != ts->dm) { 2277d5697caSHong Zhang ierr = DMGetNamedGlobalVector(dm,"TSEIMEX_Z",Z);CHKERRQ(ierr); 2287d5697caSHong Zhang } else *Z = ext->Z; 2297d5697caSHong Zhang } 2307d5697caSHong Zhang if (Ydot) { 2317d5697caSHong Zhang if (dm && dm != ts->dm) { 2327d5697caSHong Zhang ierr = DMGetNamedGlobalVector(dm,"TSEIMEX_Ydot",Ydot);CHKERRQ(ierr); 2337d5697caSHong Zhang } else *Ydot = ext->Ydot; 2347d5697caSHong Zhang } 2357d5697caSHong Zhang if (YdotI) { 2367d5697caSHong Zhang if (dm && dm != ts->dm) { 2377d5697caSHong Zhang ierr = DMGetNamedGlobalVector(dm,"TSEIMEX_YdotI",YdotI);CHKERRQ(ierr); 2387d5697caSHong Zhang } else *YdotI = ext->YdotI; 2397d5697caSHong Zhang } 2407d5697caSHong Zhang if (YdotRHS) { 2417d5697caSHong Zhang if (dm && dm != ts->dm) { 2427d5697caSHong Zhang ierr = DMGetNamedGlobalVector(dm,"TSEIMEX_YdotRHS",YdotRHS);CHKERRQ(ierr); 2437d5697caSHong Zhang } else *YdotRHS = ext->YdotRHS; 2447d5697caSHong Zhang } 2457d5697caSHong Zhang PetscFunctionReturn(0); 2467d5697caSHong Zhang } 2477d5697caSHong Zhang 2487d5697caSHong Zhang 2497d5697caSHong Zhang static PetscErrorCode TSEIMEXRestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydot,Vec *YdotI,Vec *YdotRHS) 2507d5697caSHong Zhang { 2517d5697caSHong Zhang PetscErrorCode ierr; 2527d5697caSHong Zhang 2537d5697caSHong Zhang PetscFunctionBegin; 2547d5697caSHong Zhang if (Z) { 2557d5697caSHong Zhang if (dm && dm != ts->dm) { 2567d5697caSHong Zhang ierr = DMRestoreNamedGlobalVector(dm,"TSEIMEX_Z",Z);CHKERRQ(ierr); 2577d5697caSHong Zhang } 2587d5697caSHong Zhang } 2597d5697caSHong Zhang if (Ydot) { 2607d5697caSHong Zhang if (dm && dm != ts->dm) { 2617d5697caSHong Zhang ierr = DMRestoreNamedGlobalVector(dm,"TSEIMEX_Ydot",Ydot);CHKERRQ(ierr); 2627d5697caSHong Zhang } 2637d5697caSHong Zhang } 2647d5697caSHong Zhang if (YdotI) { 2657d5697caSHong Zhang if (dm && dm != ts->dm) { 2667d5697caSHong Zhang ierr = DMRestoreNamedGlobalVector(dm,"TSEIMEX_YdotI",YdotI);CHKERRQ(ierr); 2677d5697caSHong Zhang } 2687d5697caSHong Zhang } 2697d5697caSHong Zhang if (YdotRHS) { 2707d5697caSHong Zhang if (dm && dm != ts->dm) { 2717d5697caSHong Zhang ierr = DMRestoreNamedGlobalVector(dm,"TSEIMEX_YdotRHS",YdotRHS);CHKERRQ(ierr); 2727d5697caSHong Zhang } 2737d5697caSHong Zhang } 2747d5697caSHong Zhang PetscFunctionReturn(0); 2757d5697caSHong Zhang } 2767d5697caSHong Zhang 2777d5697caSHong Zhang 2787d5697caSHong Zhang /* 2797d5697caSHong Zhang This defines the nonlinear equation that is to be solved with SNES 2807d5697caSHong Zhang Fn[t0+Theta*dt, U, (U-U0)*shift] = 0 2817d5697caSHong Zhang In the case of Backward Euler, Fn = (U-U0)/h-g(t1,U)) 2827d5697caSHong Zhang Since FormIFunction calculates G = ydot - g(t,y), ydot will be set to (U-U0)/h 2837d5697caSHong Zhang */ 2847d5697caSHong Zhang static PetscErrorCode SNESTSFormFunction_EIMEX(SNES snes,Vec X,Vec G,TS ts) 2857d5697caSHong Zhang { 2867d5697caSHong Zhang TS_EIMEX *ext = (TS_EIMEX*)ts->data; 2877d5697caSHong Zhang PetscErrorCode ierr; 2887d5697caSHong Zhang Vec Ydot,Z; 2897d5697caSHong Zhang DM dm,dmsave; 2907d5697caSHong Zhang 2917d5697caSHong Zhang PetscFunctionBegin; 2927d5697caSHong Zhang ierr = VecZeroEntries(G);CHKERRQ(ierr); 2937d5697caSHong Zhang 2947d5697caSHong Zhang ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 29516dc5bacSEmil Constantinescu ierr = TSEIMEXGetVecs(ts,dm,&Z,&Ydot,NULL,NULL);CHKERRQ(ierr); 2967d5697caSHong Zhang ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); 2977d5697caSHong Zhang dmsave = ts->dm; 2987d5697caSHong Zhang ts->dm = dm; 2997d5697caSHong Zhang ierr = TSComputeIFunction(ts,ext->ctime,X,Ydot,G,PETSC_FALSE);CHKERRQ(ierr); 3007d5697caSHong Zhang /* PETSC_FALSE indicates non-imex, adding explicit RHS to the implicit I function. */ 3017d5697caSHong Zhang ierr = VecCopy(G,Ydot);CHKERRQ(ierr); 3027d5697caSHong Zhang ts->dm = dmsave; 30316dc5bacSEmil Constantinescu ierr = TSEIMEXRestoreVecs(ts,dm,&Z,&Ydot,NULL,NULL);CHKERRQ(ierr); 3047d5697caSHong Zhang 3057d5697caSHong Zhang PetscFunctionReturn(0); 3067d5697caSHong Zhang } 3077d5697caSHong Zhang 3087d5697caSHong Zhang /* 3097d5697caSHong Zhang This defined the Jacobian matrix for SNES. Jn = (I/h-g'(t,y)) 3107d5697caSHong Zhang */ 311d1e9a80fSBarry Smith static PetscErrorCode SNESTSFormJacobian_EIMEX(SNES snes,Vec X,Mat A,Mat B,TS ts) 3127d5697caSHong Zhang { 3137d5697caSHong Zhang TS_EIMEX *ext = (TS_EIMEX*)ts->data; 3147d5697caSHong Zhang Vec Ydot; 3157d5697caSHong Zhang PetscErrorCode ierr; 3167d5697caSHong Zhang DM dm,dmsave; 3177d5697caSHong Zhang PetscFunctionBegin; 3187d5697caSHong Zhang ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 31916dc5bacSEmil Constantinescu ierr = TSEIMEXGetVecs(ts,dm,NULL,&Ydot,NULL,NULL);CHKERRQ(ierr); 3207d5697caSHong Zhang /* ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); */ 3217d5697caSHong Zhang /* ext->Ydot have already been computed in SNESTSFormFunction_EIMEX (SNES guarantees this) */ 3227d5697caSHong Zhang dmsave = ts->dm; 3237d5697caSHong Zhang ts->dm = dm; 324d1e9a80fSBarry Smith ierr = TSComputeIJacobian(ts,ts->ptime,X,Ydot,ext->shift,A,B,PETSC_TRUE);CHKERRQ(ierr); 3257d5697caSHong Zhang ts->dm = dmsave; 32616dc5bacSEmil Constantinescu ierr = TSEIMEXRestoreVecs(ts,dm,NULL,&Ydot,NULL,NULL);CHKERRQ(ierr); 3277d5697caSHong Zhang PetscFunctionReturn(0); 3287d5697caSHong Zhang } 3297d5697caSHong Zhang 3307d5697caSHong Zhang static PetscErrorCode DMCoarsenHook_TSEIMEX(DM fine,DM coarse,void *ctx) 3317d5697caSHong Zhang { 3327d5697caSHong Zhang 3337d5697caSHong Zhang PetscFunctionBegin; 3347d5697caSHong Zhang PetscFunctionReturn(0); 3357d5697caSHong Zhang } 3367d5697caSHong Zhang 3377d5697caSHong Zhang static PetscErrorCode DMRestrictHook_TSEIMEX(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 3387d5697caSHong Zhang { 3397d5697caSHong Zhang TS ts = (TS)ctx; 3407d5697caSHong Zhang PetscErrorCode ierr; 3417d5697caSHong Zhang Vec Z,Z_c; 3427d5697caSHong Zhang 3437d5697caSHong Zhang PetscFunctionBegin; 34416dc5bacSEmil Constantinescu ierr = TSEIMEXGetVecs(ts,fine,&Z,NULL,NULL,NULL);CHKERRQ(ierr); 34516dc5bacSEmil Constantinescu ierr = TSEIMEXGetVecs(ts,coarse,&Z_c,NULL,NULL,NULL);CHKERRQ(ierr); 3467d5697caSHong Zhang ierr = MatRestrict(restrct,Z,Z_c);CHKERRQ(ierr); 3477d5697caSHong Zhang ierr = VecPointwiseMult(Z_c,rscale,Z_c);CHKERRQ(ierr); 34816dc5bacSEmil Constantinescu ierr = TSEIMEXRestoreVecs(ts,fine,&Z,NULL,NULL,NULL);CHKERRQ(ierr); 34916dc5bacSEmil Constantinescu ierr = TSEIMEXRestoreVecs(ts,coarse,&Z_c,NULL,NULL,NULL);CHKERRQ(ierr); 3507d5697caSHong Zhang PetscFunctionReturn(0); 3517d5697caSHong Zhang } 3527d5697caSHong Zhang 3537d5697caSHong Zhang 3547d5697caSHong Zhang static PetscErrorCode TSSetUp_EIMEX(TS ts) 3557d5697caSHong Zhang { 3567d5697caSHong Zhang TS_EIMEX *ext = (TS_EIMEX*)ts->data; 3577d5697caSHong Zhang PetscErrorCode ierr; 3587d5697caSHong Zhang DM dm; 3597d5697caSHong Zhang 360955c1f14SBarry Smith PetscFunctionBegin; 3617d5697caSHong Zhang if (!ext->N){ /* ext->max_rows not set */ 3627d5697caSHong Zhang ierr = TSEIMEXSetMaxRows(ts,TSEIMEXDefault);CHKERRQ(ierr); 3637d5697caSHong Zhang } 3647d5697caSHong Zhang if(-1 == ext->row_ind && -1 == ext->col_ind){ 3657d5697caSHong Zhang ierr = TSEIMEXSetRowCol(ts,ext->max_rows,ext->max_rows);CHKERRQ(ierr); 3667d5697caSHong Zhang } else{/* ext->row_ind and col_ind already set */ 3677d5697caSHong Zhang if (ext->ord_adapt){ 3687d5697caSHong Zhang ierr = PetscInfo(ts,"Order adaptivity is enabled and TSEIMEXSetRowCol or -ts_eimex_row_col option will take no effect\n");CHKERRQ(ierr); 3697d5697caSHong Zhang } 3707d5697caSHong Zhang } 3717d5697caSHong Zhang 3727d5697caSHong Zhang if(ext->ord_adapt){ 3737d5697caSHong Zhang ext->nstages = 2; /* Start with the 2-stage scheme */ 3747d5697caSHong Zhang ierr = TSEIMEXSetRowCol(ts,ext->nstages,ext->nstages);CHKERRQ(ierr); 3757d5697caSHong Zhang } else{ 3767d5697caSHong Zhang ext->nstages = ext->max_rows; /* by default nstages is the same as max_rows, this can be changed by setting order adaptivity */ 3777d5697caSHong Zhang } 3787d5697caSHong Zhang 3792ffb9264SLisandro Dalcin ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 3802ffb9264SLisandro Dalcin 3817d5697caSHong Zhang ierr = VecDuplicateVecs(ts->vec_sol,(1+ext->nstages)*ext->nstages/2,&ext->T);CHKERRQ(ierr);/* full T table */ 3827d5697caSHong Zhang ierr = VecDuplicate(ts->vec_sol,&ext->YdotI);CHKERRQ(ierr); 3837d5697caSHong Zhang ierr = VecDuplicate(ts->vec_sol,&ext->YdotRHS);CHKERRQ(ierr); 3847d5697caSHong Zhang ierr = VecDuplicate(ts->vec_sol,&ext->Ydot);CHKERRQ(ierr); 3857d5697caSHong Zhang ierr = VecDuplicate(ts->vec_sol,&ext->VecSolPrev);CHKERRQ(ierr); 3867d5697caSHong Zhang ierr = VecDuplicate(ts->vec_sol,&ext->Y);CHKERRQ(ierr); 3877d5697caSHong Zhang ierr = VecDuplicate(ts->vec_sol,&ext->Z);CHKERRQ(ierr); 3887d5697caSHong Zhang ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 3897d5697caSHong Zhang if (dm) { 3907d5697caSHong Zhang ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSEIMEX,DMRestrictHook_TSEIMEX,ts);CHKERRQ(ierr); 3917d5697caSHong Zhang } 3927d5697caSHong Zhang PetscFunctionReturn(0); 3937d5697caSHong Zhang } 3947d5697caSHong Zhang 3954416b707SBarry Smith static PetscErrorCode TSSetFromOptions_EIMEX(PetscOptionItems *PetscOptionsObject,TS ts) 3967d5697caSHong Zhang { 3977d5697caSHong Zhang TS_EIMEX *ext = (TS_EIMEX*)ts->data; 3987d5697caSHong Zhang PetscErrorCode ierr; 3999200755eSBarry Smith PetscInt tindex[2]; 4007d5697caSHong Zhang PetscInt np = 2, nrows=TSEIMEXDefault; 4019200755eSBarry Smith 4027d5697caSHong Zhang PetscFunctionBegin; 4039200755eSBarry Smith tindex[0] = TSEIMEXDefault; 4049200755eSBarry Smith tindex[1] = TSEIMEXDefault; 405e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"EIMEX ODE solver options");CHKERRQ(ierr); 4067d5697caSHong Zhang { 4077d5697caSHong Zhang PetscBool flg; 4087d5697caSHong Zhang ierr = PetscOptionsInt("-ts_eimex_max_rows","Define the maximum number of rows used","TSEIMEXSetMaxRows",nrows,&nrows,&flg);CHKERRQ(ierr); /* default value 3 */ 4097d5697caSHong Zhang if(flg){ 4107d5697caSHong Zhang ierr = TSEIMEXSetMaxRows(ts,nrows);CHKERRQ(ierr); 4117d5697caSHong Zhang } 4127d5697caSHong Zhang ierr = PetscOptionsIntArray("-ts_eimex_row_col","Return the specific term in the T table","TSEIMEXSetRowCol",tindex,&np,&flg);CHKERRQ(ierr); 4137d5697caSHong Zhang if(flg){ 4147d5697caSHong Zhang ierr = TSEIMEXSetRowCol(ts,tindex[0],tindex[1]);CHKERRQ(ierr); 4157d5697caSHong Zhang } 41616dc5bacSEmil Constantinescu ierr = PetscOptionsBool("-ts_eimex_order_adapt","Solve the problem with adaptive order","TSEIMEXSetOrdAdapt",ext->ord_adapt,&ext->ord_adapt,NULL);CHKERRQ(ierr); 4177d5697caSHong Zhang } 4187d5697caSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4197d5697caSHong Zhang PetscFunctionReturn(0); 4207d5697caSHong Zhang } 4217d5697caSHong Zhang 4227d5697caSHong Zhang static PetscErrorCode TSView_EIMEX(TS ts,PetscViewer viewer) 4237d5697caSHong Zhang { 4247d5697caSHong Zhang PetscFunctionBegin; 4257d5697caSHong Zhang PetscFunctionReturn(0); 4267d5697caSHong Zhang } 4277d5697caSHong Zhang 4287d5697caSHong Zhang 4297d5697caSHong Zhang /*@C 4307d5697caSHong Zhang TSEIMEXSetMaxRows - Set the maximum number of rows for EIMEX schemes 4317d5697caSHong Zhang 4327d5697caSHong Zhang Logically collective 4337d5697caSHong Zhang 4347d5697caSHong Zhang Input Parameter: 4357d5697caSHong Zhang + ts - timestepping context 4367d5697caSHong Zhang - nrows - maximum number of rows 4377d5697caSHong Zhang 4387d5697caSHong Zhang Level: intermediate 4397d5697caSHong Zhang 4407d5697caSHong Zhang .seealso: TSEIMEXSetRowCol(), TSEIMEXSetOrdAdapt(), TSEIMEX 4417d5697caSHong Zhang @*/ 4427d5697caSHong Zhang PetscErrorCode TSEIMEXSetMaxRows(TS ts, PetscInt nrows) 4437d5697caSHong Zhang { 4447d5697caSHong Zhang PetscErrorCode ierr; 4457d5697caSHong Zhang PetscFunctionBegin; 4467d5697caSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4477d5697caSHong Zhang ierr = PetscTryMethod(ts,"TSEIMEXSetMaxRows_C",(TS,PetscInt),(ts,nrows));CHKERRQ(ierr); 4487d5697caSHong Zhang PetscFunctionReturn(0); 4497d5697caSHong Zhang } 4507d5697caSHong Zhang 4517d5697caSHong Zhang 4527d5697caSHong Zhang /*@C 4537d5697caSHong Zhang TSEIMEXSetRowCol - Set the type index in the T table for the return value 4547d5697caSHong Zhang 4557d5697caSHong Zhang Logically collective 4567d5697caSHong Zhang 4577d5697caSHong Zhang Input Parameter: 4587d5697caSHong Zhang + ts - timestepping context 4597d5697caSHong Zhang - tindex - index in the T table 4607d5697caSHong Zhang 4617d5697caSHong Zhang Level: intermediate 4627d5697caSHong Zhang 4637d5697caSHong Zhang .seealso: TSEIMEXSetMaxRows(), TSEIMEXSetOrdAdapt(), TSEIMEX 4647d5697caSHong Zhang @*/ 4657d5697caSHong Zhang PetscErrorCode TSEIMEXSetRowCol(TS ts, PetscInt row, PetscInt col) 4667d5697caSHong Zhang { 4677d5697caSHong Zhang PetscErrorCode ierr; 4687d5697caSHong Zhang PetscFunctionBegin; 4697d5697caSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4707d5697caSHong Zhang ierr = PetscTryMethod(ts,"TSEIMEXSetRowCol_C",(TS,PetscInt, PetscInt),(ts,row,col));CHKERRQ(ierr); 4717d5697caSHong Zhang PetscFunctionReturn(0); 4727d5697caSHong Zhang } 4737d5697caSHong Zhang 4747d5697caSHong Zhang 4757d5697caSHong Zhang /*@C 4767d5697caSHong Zhang TSEIMEXSetOrdAdapt - Set the order adaptativity 4777d5697caSHong Zhang 4787d5697caSHong Zhang Logically collective 4797d5697caSHong Zhang 4807d5697caSHong Zhang Input Parameter: 4817d5697caSHong Zhang + ts - timestepping context 4827d5697caSHong Zhang - tindex - index in the T table 4837d5697caSHong Zhang 4847d5697caSHong Zhang Level: intermediate 4857d5697caSHong Zhang 4867d5697caSHong Zhang .seealso: TSEIMEXSetRowCol(), TSEIMEXSetOrdAdapt(), TSEIMEX 4877d5697caSHong Zhang @*/ 4887d5697caSHong Zhang PetscErrorCode TSEIMEXSetOrdAdapt(TS ts, PetscBool flg) 4897d5697caSHong Zhang { 4907d5697caSHong Zhang PetscErrorCode ierr; 4917d5697caSHong Zhang PetscFunctionBegin; 4927d5697caSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 4937d5697caSHong Zhang ierr = PetscTryMethod(ts,"TSEIMEXSetOrdAdapt_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 4947d5697caSHong Zhang PetscFunctionReturn(0); 4957d5697caSHong Zhang } 4967d5697caSHong Zhang 4977d5697caSHong Zhang 4984042b796SJed Brown static PetscErrorCode TSEIMEXSetMaxRows_EIMEX(TS ts,PetscInt nrows) 4997d5697caSHong Zhang { 5007d5697caSHong Zhang TS_EIMEX *ext = (TS_EIMEX*)ts->data; 5017d5697caSHong Zhang PetscErrorCode ierr; 5027d5697caSHong Zhang PetscInt i; 5037d5697caSHong Zhang 5047d5697caSHong Zhang PetscFunctionBegin; 5057d5697caSHong 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); 5067d5697caSHong Zhang ierr = PetscFree(ext->N);CHKERRQ(ierr); 5077d5697caSHong Zhang ext->max_rows = nrows; 508785e854fSJed Brown ierr = PetscMalloc1(nrows,&ext->N);CHKERRQ(ierr); 5097d5697caSHong Zhang for(i=0;i<nrows;i++) ext->N[i]=i+1; 5107d5697caSHong Zhang PetscFunctionReturn(0); 5117d5697caSHong Zhang } 5127d5697caSHong Zhang 5134042b796SJed Brown static PetscErrorCode TSEIMEXSetRowCol_EIMEX(TS ts,PetscInt row,PetscInt col) 5147d5697caSHong Zhang { 5157d5697caSHong Zhang TS_EIMEX *ext = (TS_EIMEX*)ts->data; 5167d5697caSHong Zhang 5177d5697caSHong Zhang PetscFunctionBegin; 5187d5697caSHong 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); 5197d5697caSHong 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); 5207d5697caSHong Zhang if (col > row) SETERRQ2(((PetscObject)ts)->comm,PETSC_ERR_ARG_OUTOFRANGE,"The column index (%d) exceeds the row index (%d)\n",col,row); 5217d5697caSHong Zhang 5227d5697caSHong Zhang ext->row_ind = row - 1; 5237d5697caSHong Zhang ext->col_ind = col - 1; /* Array index in C starts from 0 */ 5247d5697caSHong Zhang PetscFunctionReturn(0); 5257d5697caSHong Zhang } 5267d5697caSHong Zhang 5274042b796SJed Brown static PetscErrorCode TSEIMEXSetOrdAdapt_EIMEX(TS ts,PetscBool flg) 5287d5697caSHong Zhang { 5297d5697caSHong Zhang TS_EIMEX *ext = (TS_EIMEX*)ts->data; 5307d5697caSHong Zhang PetscFunctionBegin; 5317d5697caSHong Zhang ext->ord_adapt = flg; 5327d5697caSHong Zhang PetscFunctionReturn(0); 5337d5697caSHong Zhang } 5347d5697caSHong Zhang 5357d5697caSHong Zhang /*MC 536679f3ef1SBarry Smith TSEIMEX - Time stepping with Extrapolated IMEX methods. 537679f3ef1SBarry Smith 538679f3ef1SBarry Smith These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly nonlinear such that it 539679f3ef1SBarry Smith is expensive to solve with a fully implicit method. The user should provide the stiff part of the equation using TSSetIFunction() and the 540679f3ef1SBarry Smith non-stiff part with TSSetRHSFunction(). 5417d5697caSHong Zhang 5427d5697caSHong Zhang Notes: 5437d5697caSHong Zhang The default is a 3-stage scheme, it can be changed with TSEIMEXSetMaxRows() or -ts_eimex_max_rows 5447d5697caSHong Zhang 5457d5697caSHong Zhang This method currently only works with ODE, for which the stiff part G(t,X,Xdot) has the form Xdot + Ghat(t,X). 5467d5697caSHong Zhang 547679f3ef1SBarry Smith The general system is written as 548679f3ef1SBarry Smith 549679f3ef1SBarry Smith G(t,X,Xdot) = F(t,X) 550679f3ef1SBarry Smith 551679f3ef1SBarry Smith where G represents the stiff part and F represents the non-stiff part. The user should provide the stiff part 552679f3ef1SBarry Smith of the equation using TSSetIFunction() and the non-stiff part with TSSetRHSFunction(). 553679f3ef1SBarry Smith This method is designed to be linearly implicit on G and can use an approximate and lagged Jacobian. 554679f3ef1SBarry Smith 555679f3ef1SBarry Smith Another common form for the system is 556679f3ef1SBarry Smith 557679f3ef1SBarry Smith y'=f(x)+g(x) 558679f3ef1SBarry Smith 559679f3ef1SBarry Smith The relationship between F,G and f,g is 560679f3ef1SBarry Smith 561679f3ef1SBarry Smith G = y'-g(x), F = f(x) 562679f3ef1SBarry Smith 563679f3ef1SBarry Smith References 564679f3ef1SBarry Smith E. Constantinescu and A. Sandu, Extrapolated implicit-explicit time stepping, SIAM Journal on Scientific 565679f3ef1SBarry Smith Computing, 31 (2010), pp. 4452-4477. 566679f3ef1SBarry Smith 5677d5697caSHong Zhang Level: beginner 5687d5697caSHong Zhang 569679f3ef1SBarry Smith .seealso: TSCreate(), TS, TSSetType(), TSEIMEXSetMaxRows(), TSEIMEXSetRowCol(), TSEIMEXSetOrdAdapt() 570679f3ef1SBarry Smith 5717d5697caSHong Zhang M*/ 5724042b796SJed Brown PETSC_EXTERN PetscErrorCode TSCreate_EIMEX(TS ts) 5737d5697caSHong Zhang { 5747d5697caSHong Zhang TS_EIMEX *ext; 5757d5697caSHong Zhang PetscErrorCode ierr; 5767d5697caSHong Zhang 5777d5697caSHong Zhang PetscFunctionBegin; 5787d5697caSHong Zhang 5797d5697caSHong Zhang ts->ops->reset = TSReset_EIMEX; 5807d5697caSHong Zhang ts->ops->destroy = TSDestroy_EIMEX; 5817d5697caSHong Zhang ts->ops->view = TSView_EIMEX; 5827d5697caSHong Zhang ts->ops->setup = TSSetUp_EIMEX; 5837d5697caSHong Zhang ts->ops->step = TSStep_EIMEX; 5847d5697caSHong Zhang ts->ops->interpolate = TSInterpolate_EIMEX; 5857d5697caSHong Zhang ts->ops->evaluatestep = TSEvaluateStep_EIMEX; 5867d5697caSHong Zhang ts->ops->setfromoptions = TSSetFromOptions_EIMEX; 5877d5697caSHong Zhang ts->ops->snesfunction = SNESTSFormFunction_EIMEX; 5887d5697caSHong Zhang ts->ops->snesjacobian = SNESTSFormJacobian_EIMEX; 5892ffb9264SLisandro Dalcin ts->default_adapt_type = TSADAPTNONE; 5907d5697caSHong Zhang 591*efd4aadfSBarry Smith ts->usessnes = PETSC_TRUE; 592*efd4aadfSBarry Smith 593b00a9115SJed Brown ierr = PetscNewLog(ts,&ext);CHKERRQ(ierr); 5947d5697caSHong Zhang ts->data = (void*)ext; 5957d5697caSHong Zhang 5967d5697caSHong Zhang ext->ord_adapt = PETSC_FALSE; /* By default, no order adapativity */ 5977d5697caSHong Zhang ext->row_ind = -1; 5987d5697caSHong Zhang ext->col_ind = -1; 5997d5697caSHong Zhang ext->max_rows = TSEIMEXDefault; 6007d5697caSHong Zhang ext->nstages = TSEIMEXDefault; 6017d5697caSHong Zhang 602e1d27e54SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetMaxRows_C", TSEIMEXSetMaxRows_EIMEX);CHKERRQ(ierr); 603e1d27e54SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetRowCol_C", TSEIMEXSetRowCol_EIMEX);CHKERRQ(ierr); 604e1d27e54SJed Brown ierr = PetscObjectComposeFunction((PetscObject)ts,"TSEIMEXSetOrdAdapt_C",TSEIMEXSetOrdAdapt_EIMEX);CHKERRQ(ierr); 6057d5697caSHong Zhang PetscFunctionReturn(0); 6067d5697caSHong Zhang } 607