xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision b4eba00b70283cae48e61570c0129be6b6b148d1)
17cb94ee6SHong Zhang /*
2db268bfaSHong Zhang     Provides a PETSc interface to SUNDIALS/CVODE solver.
37cb94ee6SHong Zhang     The interface to PVODE (old version of CVODE) was originally contributed
47cb94ee6SHong Zhang     by Liyang Xu. It has been redone by Hong Zhang and Dinesh Kaushik.
528adb3f7SHong Zhang 
628adb3f7SHong Zhang     Reference: sundials-2.4.0/examples/cvode/parallel/cvDiurnal_kry_p.c
77cb94ee6SHong Zhang */
856a740aaSHong Zhang #include "sundials.h"  /*I "petscts.h" I*/
97cb94ee6SHong Zhang 
107cb94ee6SHong Zhang /*
117cb94ee6SHong Zhang       TSPrecond_Sundials - function that we provide to SUNDIALS to
127cb94ee6SHong Zhang                         evaluate the preconditioner.
137cb94ee6SHong Zhang */
147cb94ee6SHong Zhang #undef __FUNCT__
157cb94ee6SHong Zhang #define __FUNCT__ "TSPrecond_Sundials"
167cb94ee6SHong Zhang PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy,
177cb94ee6SHong Zhang                     booleantype jok,booleantype *jcurPtr,
187cb94ee6SHong Zhang                     realtype _gamma,void *P_data,
197cb94ee6SHong Zhang                     N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3)
207cb94ee6SHong Zhang {
217cb94ee6SHong Zhang   TS             ts = (TS) P_data;
227cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
23dcbc6d53SSean Farley   PC             pc;
247cb94ee6SHong Zhang   PetscErrorCode ierr;
250679a0aeSJed Brown   Mat            J,P;
260679a0aeSJed Brown   Vec            yy = cvode->w1,yydot = cvode->ydot;
270679a0aeSJed Brown   PetscReal      gm = (PetscReal)_gamma;
287cb94ee6SHong Zhang   MatStructure   str = DIFFERENT_NONZERO_PATTERN;
297cb94ee6SHong Zhang   PetscScalar    *y_data;
307cb94ee6SHong Zhang 
317cb94ee6SHong Zhang   PetscFunctionBegin;
320679a0aeSJed Brown   ierr = TSGetIJacobian(ts,&J,&P,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
337cb94ee6SHong Zhang   y_data = (PetscScalar *) N_VGetArrayPointer(y);
347cb94ee6SHong Zhang   ierr = VecPlaceArray(yy,y_data); CHKERRQ(ierr);
350679a0aeSJed Brown   ierr = VecZeroEntries(yydot);CHKERRQ(ierr); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */
360679a0aeSJed Brown   /* compute the shifted Jacobian   (1/gm)*I + Jrest */
370679a0aeSJed Brown   ierr = TSComputeIJacobian(ts,ts->ptime,yy,yydot,1/gm,&J,&P,&str,PETSC_FALSE);CHKERRQ(ierr);
387cb94ee6SHong Zhang   ierr = VecResetArray(yy); CHKERRQ(ierr);
390679a0aeSJed Brown   ierr = MatScale(P,gm);CHKERRQ(ierr);  /* turn into I-gm*Jrest, J is not used by Sundials  */
407cb94ee6SHong Zhang   *jcurPtr = TRUE;
41dcbc6d53SSean Farley   ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
420679a0aeSJed Brown   ierr = PCSetOperators(pc,J,P,str);CHKERRQ(ierr);
437cb94ee6SHong Zhang   PetscFunctionReturn(0);
447cb94ee6SHong Zhang }
457cb94ee6SHong Zhang 
467cb94ee6SHong Zhang /*
477cb94ee6SHong Zhang      TSPSolve_Sundials -  routine that we provide to Sundials that applies the preconditioner.
487cb94ee6SHong Zhang */
497cb94ee6SHong Zhang #undef __FUNCT__
507cb94ee6SHong Zhang #define __FUNCT__ "TSPSolve_Sundials"
514ac7836bSHong Zhang PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z,
524ac7836bSHong Zhang                                  realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp)
537cb94ee6SHong Zhang {
547cb94ee6SHong Zhang   TS              ts = (TS) P_data;
557cb94ee6SHong Zhang   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
56dcbc6d53SSean Farley   PC              pc;
577cb94ee6SHong Zhang   Vec             rr = cvode->w1,zz = cvode->w2;
587cb94ee6SHong Zhang   PetscErrorCode  ierr;
597cb94ee6SHong Zhang   PetscScalar     *r_data,*z_data;
607cb94ee6SHong Zhang 
617cb94ee6SHong Zhang   PetscFunctionBegin;
627cb94ee6SHong Zhang   /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
637cb94ee6SHong Zhang   r_data  = (PetscScalar *) N_VGetArrayPointer(r);
647cb94ee6SHong Zhang   z_data  = (PetscScalar *) N_VGetArrayPointer(z);
657cb94ee6SHong Zhang   ierr = VecPlaceArray(rr,r_data); CHKERRQ(ierr);
667cb94ee6SHong Zhang   ierr = VecPlaceArray(zz,z_data); CHKERRQ(ierr);
674ac7836bSHong Zhang 
687cb94ee6SHong Zhang   /* Solve the Px=r and put the result in zz */
69dcbc6d53SSean Farley   ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
707cb94ee6SHong Zhang   ierr = PCApply(pc,rr,zz); CHKERRQ(ierr);
717cb94ee6SHong Zhang   ierr = VecResetArray(rr); CHKERRQ(ierr);
727cb94ee6SHong Zhang   ierr = VecResetArray(zz); CHKERRQ(ierr);
737cb94ee6SHong Zhang   PetscFunctionReturn(0);
747cb94ee6SHong Zhang }
757cb94ee6SHong Zhang 
767cb94ee6SHong Zhang /*
777cb94ee6SHong Zhang         TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
787cb94ee6SHong Zhang */
797cb94ee6SHong Zhang #undef __FUNCT__
807cb94ee6SHong Zhang #define __FUNCT__ "TSFunction_Sundials"
814ac7836bSHong Zhang int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx)
827cb94ee6SHong Zhang {
837cb94ee6SHong Zhang   TS              ts = (TS) ctx;
847adad957SLisandro Dalcin   MPI_Comm        comm = ((PetscObject)ts)->comm;
857cb94ee6SHong Zhang   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
860679a0aeSJed Brown   Vec             yy = cvode->w1,yyd = cvode->w2,yydot = cvode->ydot;
877cb94ee6SHong Zhang   PetscScalar     *y_data,*ydot_data;
887cb94ee6SHong Zhang   PetscErrorCode  ierr;
897cb94ee6SHong Zhang 
907cb94ee6SHong Zhang   PetscFunctionBegin;
915ff03cf7SDinesh Kaushik   /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
927cb94ee6SHong Zhang   y_data     = (PetscScalar *) N_VGetArrayPointer(y);
937cb94ee6SHong Zhang   ydot_data  = (PetscScalar *) N_VGetArrayPointer(ydot);
944ac7836bSHong Zhang   ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr);
954ac7836bSHong Zhang   ierr = VecPlaceArray(yyd,ydot_data); CHKERRABORT(comm,ierr);
960679a0aeSJed Brown   ierr = VecZeroEntries(yydot);CHKERRQ(ierr);
974ac7836bSHong Zhang 
987cb94ee6SHong Zhang   /* now compute the right hand side function */
990679a0aeSJed Brown   ierr = TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE); CHKERRABORT(comm,ierr);
1000679a0aeSJed Brown   ierr = VecScale(yyd,-1.);CHKERRQ(ierr);
1017cb94ee6SHong Zhang   ierr = VecResetArray(yy); CHKERRABORT(comm,ierr);
1027cb94ee6SHong Zhang   ierr = VecResetArray(yyd); CHKERRABORT(comm,ierr);
1034ac7836bSHong Zhang   PetscFunctionReturn(0);
1047cb94ee6SHong Zhang }
1057cb94ee6SHong Zhang 
1067cb94ee6SHong Zhang /*
107*b4eba00bSSean Farley        TSStep_Sundials - Calls Sundials to integrate the ODE.
1087cb94ee6SHong Zhang */
1097cb94ee6SHong Zhang #undef __FUNCT__
110*b4eba00bSSean Farley #define __FUNCT__ "TSStep_Sundials"
111*b4eba00bSSean Farley PetscErrorCode TSStep_Sundials(TS ts)
1127cb94ee6SHong Zhang {
1137cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
1147cb94ee6SHong Zhang   Vec            sol = ts->vec_sol;
1157cb94ee6SHong Zhang   PetscErrorCode ierr;
116*b4eba00bSSean Farley   PetscInt       flag;
1179f94935aSHong Zhang   long int       its,nsteps;
1187cb94ee6SHong Zhang   realtype       t,tout;
1197cb94ee6SHong Zhang   PetscScalar    *y_data;
1207cb94ee6SHong Zhang   void           *mem;
121abe29043SSean Farley   SNES           snes;
122abe29043SSean Farley   Vec            res; /* This, together with snes, will check if the SNES vec_func has been set */
1237cb94ee6SHong Zhang 
1247cb94ee6SHong Zhang   PetscFunctionBegin;
12516016685SHong Zhang   mem  = cvode->mem;
1267cb94ee6SHong Zhang   tout = ts->max_time;
1277cb94ee6SHong Zhang   ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr);
1287cb94ee6SHong Zhang   N_VSetArrayPointer((realtype *)y_data,cvode->y);
1297cb94ee6SHong Zhang   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
130186e87acSLisandro Dalcin 
131*b4eba00bSSean Farley   /* Should think about moving this outside the loop */
132abe29043SSean Farley   ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr);
133abe29043SSean Farley   ierr = SNESGetFunction(snes, &res, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
134abe29043SSean Farley   if (!res) {
135abe29043SSean Farley     ierr = TSSetIFunction(ts, sol, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
136abe29043SSean Farley   }
137abe29043SSean Farley 
1383f2090d5SJed Brown   ierr = TSPreStep(ts);CHKERRQ(ierr);
139186e87acSLisandro Dalcin 
1402bfc04deSHong Zhang   if (cvode->monitorstep) {
14128adb3f7SHong Zhang     flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
1422bfc04deSHong Zhang   } else {
1432bfc04deSHong Zhang     flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL);
1442bfc04deSHong Zhang   }
1459f94935aSHong Zhang 
1469f94935aSHong Zhang   if (flag){ /* display error message */
1479f94935aSHong Zhang     switch (flag){
1489f94935aSHong Zhang       case CV_ILL_INPUT:
1499f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT");
1509f94935aSHong Zhang         break;
1519f94935aSHong Zhang       case CV_TOO_CLOSE:
1529f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE");
1539f94935aSHong Zhang         break;
1543c7fefeeSJed Brown       case CV_TOO_MUCH_WORK: {
1559f94935aSHong Zhang         PetscReal      tcur;
1569f94935aSHong Zhang         ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
1579f94935aSHong Zhang         ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr);
1589f94935aSHong Zhang         SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_WORK. At t=%G, nsteps %D exceeds mxstep %D. Increase '-ts_max_steps <>' or modify TSSetDuration()",tcur,nsteps,ts->max_steps);
1593c7fefeeSJed Brown       } break;
1609f94935aSHong Zhang       case CV_TOO_MUCH_ACC:
1619f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC");
1629f94935aSHong Zhang         break;
1639f94935aSHong Zhang       case CV_ERR_FAILURE:
1649f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE");
1659f94935aSHong Zhang         break;
1669f94935aSHong Zhang       case CV_CONV_FAILURE:
1679f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE");
1689f94935aSHong Zhang         break;
1699f94935aSHong Zhang       case CV_LINIT_FAIL:
1709f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL");
1719f94935aSHong Zhang         break;
1729f94935aSHong Zhang       case CV_LSETUP_FAIL:
1739f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL");
1749f94935aSHong Zhang         break;
1759f94935aSHong Zhang       case CV_LSOLVE_FAIL:
1769f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL");
1779f94935aSHong Zhang         break;
1789f94935aSHong Zhang       case CV_RHSFUNC_FAIL:
1799f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL");
1809f94935aSHong Zhang         break;
1819f94935aSHong Zhang       case CV_FIRST_RHSFUNC_ERR:
1829f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR");
1839f94935aSHong Zhang         break;
1849f94935aSHong Zhang       case CV_REPTD_RHSFUNC_ERR:
1859f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR");
1869f94935aSHong Zhang         break;
1879f94935aSHong Zhang       case CV_UNREC_RHSFUNC_ERR:
1889f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR");
1899f94935aSHong Zhang         break;
1909f94935aSHong Zhang       case CV_RTFUNC_FAIL:
1919f94935aSHong Zhang         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL");
1929f94935aSHong Zhang         break;
1939f94935aSHong Zhang       default:
1949f94935aSHong Zhang         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag);
1959f94935aSHong Zhang     }
1969f94935aSHong Zhang   }
1979f94935aSHong Zhang 
1987cb94ee6SHong Zhang   /* copy the solution from cvode->y to cvode->update and sol */
1997cb94ee6SHong Zhang   ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr);
2007cb94ee6SHong Zhang   ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr);
2017cb94ee6SHong Zhang   ierr = VecResetArray(cvode->w1); CHKERRQ(ierr);
2027cb94ee6SHong Zhang   ierr = VecCopy(cvode->update,sol);CHKERRQ(ierr);
2037cb94ee6SHong Zhang   ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr);
2044ac7836bSHong Zhang   ierr = CVSpilsGetNumLinIters(mem, &its);
205186e87acSLisandro Dalcin   ts->nonlinear_its = its; ts->linear_its = its;
206186e87acSLisandro Dalcin 
207186e87acSLisandro Dalcin   ts->time_step = t - ts->ptime;
208186e87acSLisandro Dalcin   ts->ptime     = t;
2097cb94ee6SHong Zhang   ts->steps++;
210186e87acSLisandro Dalcin 
2119f94935aSHong Zhang   ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
212*b4eba00bSSean Farley   if (!cvode->monitorstep) ts->steps = nsteps;
213*b4eba00bSSean Farley   PetscFunctionReturn(0);
214*b4eba00bSSean Farley }
215*b4eba00bSSean Farley 
216*b4eba00bSSean Farley #undef __FUNCT__
217*b4eba00bSSean Farley #define __FUNCT__ "TSInterpolate_Sundials"
218*b4eba00bSSean Farley static PetscErrorCode TSInterpolate_Sundials(TS ts,PetscReal t,Vec X)
219*b4eba00bSSean Farley {
220*b4eba00bSSean Farley   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
221*b4eba00bSSean Farley   N_Vector        y;
222*b4eba00bSSean Farley   PetscErrorCode  ierr;
223*b4eba00bSSean Farley   PetscScalar     *x_data;
224*b4eba00bSSean Farley   PetscInt        glosize,locsize;
225*b4eba00bSSean Farley 
226*b4eba00bSSean Farley   PetscFunctionBegin;
227*b4eba00bSSean Farley 
228*b4eba00bSSean Farley   /* get the vector size */
229*b4eba00bSSean Farley   ierr = VecGetSize(X,&glosize);CHKERRQ(ierr);
230*b4eba00bSSean Farley   ierr = VecGetLocalSize(X,&locsize);CHKERRQ(ierr);
231*b4eba00bSSean Farley 
232*b4eba00bSSean Farley   /* allocate the memory for N_Vec y */
233*b4eba00bSSean Farley   y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
234*b4eba00bSSean Farley   if (!y) SETERRQ(PETSC_COMM_SELF,1,"Interpolated y is not allocated");
235*b4eba00bSSean Farley 
236*b4eba00bSSean Farley   ierr = VecGetArray(X,&x_data);CHKERRQ(ierr);
237*b4eba00bSSean Farley   N_VSetArrayPointer((realtype *)x_data,y);
238*b4eba00bSSean Farley   ierr = CVodeGetDky(cvode->mem,t,0,y);CHKERRQ(ierr);
239*b4eba00bSSean Farley   ierr = VecRestoreArray(X,&x_data);CHKERRQ(ierr);
240*b4eba00bSSean Farley 
2417cb94ee6SHong Zhang   PetscFunctionReturn(0);
2427cb94ee6SHong Zhang }
2437cb94ee6SHong Zhang 
2447cb94ee6SHong Zhang #undef __FUNCT__
245277b19d0SLisandro Dalcin #define __FUNCT__ "TSReset_Sundials"
246277b19d0SLisandro Dalcin PetscErrorCode TSReset_Sundials(TS ts)
247277b19d0SLisandro Dalcin {
248277b19d0SLisandro Dalcin   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
249277b19d0SLisandro Dalcin   PetscErrorCode ierr;
250277b19d0SLisandro Dalcin 
251277b19d0SLisandro Dalcin   PetscFunctionBegin;
2525c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->update);CHKERRQ(ierr);
2530679a0aeSJed Brown   ierr = VecDestroy(&cvode->ydot);CHKERRQ(ierr);
2545c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->w1);CHKERRQ(ierr);
2555c6b4a3dSJed Brown   ierr = VecDestroy(&cvode->w2);CHKERRQ(ierr);
256277b19d0SLisandro Dalcin   if (cvode->mem)    {CVodeFree(&cvode->mem);}
257277b19d0SLisandro Dalcin   PetscFunctionReturn(0);
258277b19d0SLisandro Dalcin }
259277b19d0SLisandro Dalcin 
260277b19d0SLisandro Dalcin #undef __FUNCT__
2617cb94ee6SHong Zhang #define __FUNCT__ "TSDestroy_Sundials"
2627cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts)
2637cb94ee6SHong Zhang {
2647cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2657cb94ee6SHong Zhang   PetscErrorCode ierr;
2667cb94ee6SHong Zhang 
2677cb94ee6SHong Zhang   PetscFunctionBegin;
268277b19d0SLisandro Dalcin   ierr = TSReset_Sundials(ts);CHKERRQ(ierr);
269a07356b8SSatish Balay   ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr);
270277b19d0SLisandro Dalcin   ierr = PetscFree(ts->data);CHKERRQ(ierr);
271335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","",PETSC_NULL);CHKERRQ(ierr);
272335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C","",PETSC_NULL);CHKERRQ(ierr);
273335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C","",PETSC_NULL);CHKERRQ(ierr);
274335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C","",PETSC_NULL);CHKERRQ(ierr);
275335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C","",PETSC_NULL);CHKERRQ(ierr);
276335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
277335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
278335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C","",PETSC_NULL);CHKERRQ(ierr);
279335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C","",PETSC_NULL);CHKERRQ(ierr);
280335f802eSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C","",PETSC_NULL);CHKERRQ(ierr);
2817cb94ee6SHong Zhang   PetscFunctionReturn(0);
2827cb94ee6SHong Zhang }
2837cb94ee6SHong Zhang 
2847cb94ee6SHong Zhang #undef __FUNCT__
285214bc6a2SJed Brown #define __FUNCT__ "TSSetUp_Sundials"
286214bc6a2SJed Brown PetscErrorCode TSSetUp_Sundials(TS ts)
2877cb94ee6SHong Zhang {
2887cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2897cb94ee6SHong Zhang   PetscErrorCode ierr;
29016016685SHong Zhang   PetscInt       glosize,locsize,i,flag;
2917cb94ee6SHong Zhang   PetscScalar    *y_data,*parray;
29216016685SHong Zhang   void           *mem;
293dcbc6d53SSean Farley   PC             pc;
29416016685SHong Zhang   const PCType   pctype;
295ace3abfcSBarry Smith   PetscBool      pcnone;
2967cb94ee6SHong Zhang 
2977cb94ee6SHong Zhang   PetscFunctionBegin;
2989f94935aSHong Zhang 
2997cb94ee6SHong Zhang   /* get the vector size */
3007cb94ee6SHong Zhang   ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr);
3017cb94ee6SHong Zhang   ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr);
3027cb94ee6SHong Zhang 
3037cb94ee6SHong Zhang   /* allocate the memory for N_Vec y */
304a07356b8SSatish Balay   cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
305e32f2f54SBarry Smith   if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated");
3067cb94ee6SHong Zhang 
30728adb3f7SHong Zhang   /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
3087cb94ee6SHong Zhang   ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr);
3097cb94ee6SHong Zhang   y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y);
3107cb94ee6SHong Zhang   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
3117cb94ee6SHong Zhang   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
31242528757SHong Zhang 
3137cb94ee6SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr);
3140679a0aeSJed Brown   ierr = VecDuplicate(ts->vec_sol,&cvode->ydot);CHKERRQ(ierr);
3157cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr);
3160679a0aeSJed Brown   ierr = PetscLogObjectParent(ts,cvode->ydot);CHKERRQ(ierr);
3177cb94ee6SHong Zhang 
3187cb94ee6SHong Zhang   /*
3197cb94ee6SHong Zhang     Create work vectors for the TSPSolve_Sundials() routine. Note these are
3207cb94ee6SHong Zhang     allocated with zero space arrays because the actual array space is provided
3217cb94ee6SHong Zhang     by Sundials and set using VecPlaceArray().
3227cb94ee6SHong Zhang   */
3237adad957SLisandro Dalcin   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr);
3247adad957SLisandro Dalcin   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr);
3257cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr);
3267cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr);
32716016685SHong Zhang 
32816016685SHong Zhang   /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
32916016685SHong Zhang   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
330e32f2f54SBarry Smith   if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails");
33116016685SHong Zhang   cvode->mem = mem;
33216016685SHong Zhang 
33316016685SHong Zhang   /* Set the pointer to user-defined data */
33416016685SHong Zhang   flag = CVodeSetUserData(mem, ts);
335e32f2f54SBarry Smith   if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails");
33616016685SHong Zhang 
337fc6b9e64SJed Brown   /* Sundials may choose to use a smaller initial step, but will never use a larger step. */
338fc6b9e64SJed Brown   flag = CVodeSetInitStep(mem,(realtype)ts->initial_time_step);
339fc6b9e64SJed Brown   if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetInitStep() failed");
340f1cd61daSJed Brown   if (cvode->mindt > 0) {
341f1cd61daSJed Brown     flag = CVodeSetMinStep(mem,(realtype)cvode->mindt);
3429f94935aSHong Zhang     if (flag){
3439f94935aSHong Zhang       if (flag == CV_MEM_NULL){
3449f94935aSHong Zhang         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL");
3459f94935aSHong Zhang       } else if (flag == CV_ILL_INPUT){
3469f94935aSHong Zhang         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size");
3479f94935aSHong Zhang       } else {
3489f94935aSHong Zhang         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed");
3499f94935aSHong Zhang       }
3509f94935aSHong Zhang     }
351f1cd61daSJed Brown   }
352f1cd61daSJed Brown   if (cvode->maxdt > 0) {
353f1cd61daSJed Brown     flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt);
354f1cd61daSJed Brown     if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMaxStep() failed");
355f1cd61daSJed Brown   }
356f1cd61daSJed Brown 
35716016685SHong Zhang   /* Call CVodeInit to initialize the integrator memory and specify the
35816016685SHong Zhang    * user's right hand side function in u'=f(t,u), the inital time T0, and
35916016685SHong Zhang    * the initial dependent variable vector cvode->y */
36016016685SHong Zhang   flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y);
36116016685SHong Zhang   if (flag){
362e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag);
36316016685SHong Zhang   }
36416016685SHong Zhang 
3659f94935aSHong Zhang   /* specifies scalar relative and absolute tolerances */
36616016685SHong Zhang   flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol);
36716016685SHong Zhang   if (flag){
368e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag);
36916016685SHong Zhang   }
37016016685SHong Zhang 
3719f94935aSHong Zhang   /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */
3729f94935aSHong Zhang   flag = CVodeSetMaxNumSteps(mem,ts->max_steps);
37316016685SHong Zhang 
37416016685SHong Zhang   /* call CVSpgmr to use GMRES as the linear solver.        */
37516016685SHong Zhang   /* setup the ode integrator with the given preconditioner */
376dcbc6d53SSean Farley   ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
377dcbc6d53SSean Farley   ierr = PCGetType(pc,&pctype);CHKERRQ(ierr);
378dcbc6d53SSean Farley   ierr = PetscTypeCompare((PetscObject)pc,PCNONE,&pcnone);CHKERRQ(ierr);
37916016685SHong Zhang   if (pcnone){
38016016685SHong Zhang     flag  = CVSpgmr(mem,PREC_NONE,0);
381e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
38216016685SHong Zhang   } else {
38316016685SHong Zhang     flag  = CVSpgmr(mem,PREC_LEFT,0);
384e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
38516016685SHong Zhang 
38616016685SHong Zhang     /* Set preconditioner and solve routines Precond and PSolve,
38716016685SHong Zhang      and the pointer to the user-defined block data */
38816016685SHong Zhang     flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials);
389e32f2f54SBarry Smith     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag);
39016016685SHong Zhang   }
39116016685SHong Zhang 
39216016685SHong Zhang   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
393e32f2f54SBarry Smith   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag);
3947cb94ee6SHong Zhang   PetscFunctionReturn(0);
3957cb94ee6SHong Zhang }
3967cb94ee6SHong Zhang 
3976fadb2cdSHong Zhang /* type of CVODE linear multistep method */
398dfcd32c8SHong Zhang const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0};
3996fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */
400dfcd32c8SHong Zhang const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0};
401a04cf4d8SBarry Smith 
4027cb94ee6SHong Zhang #undef __FUNCT__
403214bc6a2SJed Brown #define __FUNCT__ "TSSetFromOptions_Sundials"
404214bc6a2SJed Brown PetscErrorCode TSSetFromOptions_Sundials(TS ts)
4057cb94ee6SHong Zhang {
4067cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
4077cb94ee6SHong Zhang   PetscErrorCode ierr;
4087cb94ee6SHong Zhang   int            indx;
409ace3abfcSBarry Smith   PetscBool      flag;
4107cb94ee6SHong Zhang 
4117cb94ee6SHong Zhang   PetscFunctionBegin;
4127cb94ee6SHong Zhang   ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr);
4136fadb2cdSHong Zhang     ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr);
4147cb94ee6SHong Zhang     if (flag) {
4156fadb2cdSHong Zhang       ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr);
4167cb94ee6SHong Zhang     }
4176fadb2cdSHong Zhang     ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr);
4187cb94ee6SHong Zhang     if (flag) {
4197cb94ee6SHong Zhang       ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr);
4207cb94ee6SHong Zhang     }
4217cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr);
4227cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr);
423f1cd61daSJed Brown     ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,PETSC_NULL);CHKERRQ(ierr);
424f1cd61daSJed Brown     ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,PETSC_NULL);CHKERRQ(ierr);
4257cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr);
4267cb94ee6SHong Zhang     ierr = PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);CHKERRQ(ierr);
427acfcf0e5SJed Brown     ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);CHKERRQ(ierr);
4287cb94ee6SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
4297cb94ee6SHong Zhang   PetscFunctionReturn(0);
4307cb94ee6SHong Zhang }
4317cb94ee6SHong Zhang 
4327cb94ee6SHong Zhang #undef __FUNCT__
4337cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials"
4347cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
4357cb94ee6SHong Zhang {
4367cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
4377cb94ee6SHong Zhang   PetscErrorCode ierr;
4387cb94ee6SHong Zhang   char           *type;
4397cb94ee6SHong Zhang   char           atype[] = "Adams";
4407cb94ee6SHong Zhang   char           btype[] = "BDF: backward differentiation formula";
441ace3abfcSBarry Smith   PetscBool      iascii,isstring;
4422c823083SHong Zhang   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
4432c823083SHong Zhang   PetscInt       qlast,qcur;
4445d47aee6SHong Zhang   PetscReal      hinused,hlast,hcur,tcur,tolsfac;
44542528757SHong Zhang   PC             pc;
4467cb94ee6SHong Zhang 
4477cb94ee6SHong Zhang   PetscFunctionBegin;
4487cb94ee6SHong Zhang   if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;}
4497cb94ee6SHong Zhang   else                                     {type = btype;}
4507cb94ee6SHong Zhang 
4512692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4522692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
4537cb94ee6SHong Zhang   if (iascii) {
4547cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr);
4557cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr);
4567cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr);
4577cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr);
4587cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);CHKERRQ(ierr);
4597cb94ee6SHong Zhang     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
4607cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
4617cb94ee6SHong Zhang     } else {
4627cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
4637cb94ee6SHong Zhang     }
464f1cd61daSJed Brown     if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);}
465f1cd61daSJed Brown     if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);}
4662c823083SHong Zhang 
4675d47aee6SHong Zhang     /* Outputs from CVODE, CVSPILS */
4685d47aee6SHong Zhang     ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr);
4695d47aee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr);
4702c823083SHong Zhang     ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
4712c823083SHong Zhang                                    &nlinsetups,&nfails,&qlast,&qcur,
4722c823083SHong Zhang                                    &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr);
4732c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr);
4742c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr);
4752c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr);
4762c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr);
4772c823083SHong Zhang 
4782c823083SHong Zhang     ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr);
4792c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr);
4802c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr);
4812c823083SHong Zhang 
4822c823083SHong Zhang     ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */
4832c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr);
4842c823083SHong Zhang     ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr);
4852c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr);
48642528757SHong Zhang 
48742528757SHong Zhang     ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
48842528757SHong Zhang     ierr = PCView(pc,viewer);CHKERRQ(ierr);
4892c823083SHong Zhang     ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4902c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr);
4912c823083SHong Zhang     ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr);
4922c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr);
49342528757SHong Zhang 
4942c823083SHong Zhang     ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4952c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr);
4962c823083SHong Zhang     ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr);
4972c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr);
4987cb94ee6SHong Zhang   } else if (isstring) {
4997cb94ee6SHong Zhang     ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr);
5007cb94ee6SHong Zhang   }
5017cb94ee6SHong Zhang   PetscFunctionReturn(0);
5027cb94ee6SHong Zhang }
5037cb94ee6SHong Zhang 
5047cb94ee6SHong Zhang 
5057cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/
5067cb94ee6SHong Zhang EXTERN_C_BEGIN
5077cb94ee6SHong Zhang #undef __FUNCT__
5087cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials"
5097087cfbeSBarry Smith PetscErrorCode  TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
5107cb94ee6SHong Zhang {
5117cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5127cb94ee6SHong Zhang 
5137cb94ee6SHong Zhang   PetscFunctionBegin;
5147cb94ee6SHong Zhang   cvode->cvode_type = type;
5157cb94ee6SHong Zhang   PetscFunctionReturn(0);
5167cb94ee6SHong Zhang }
5177cb94ee6SHong Zhang EXTERN_C_END
5187cb94ee6SHong Zhang 
5197cb94ee6SHong Zhang EXTERN_C_BEGIN
5207cb94ee6SHong Zhang #undef __FUNCT__
5217cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials"
5227087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGMRESRestart_Sundials(TS ts,int restart)
5237cb94ee6SHong Zhang {
5247cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5257cb94ee6SHong Zhang 
5267cb94ee6SHong Zhang   PetscFunctionBegin;
5277cb94ee6SHong Zhang   cvode->restart = restart;
5287cb94ee6SHong Zhang   PetscFunctionReturn(0);
5297cb94ee6SHong Zhang }
5307cb94ee6SHong Zhang EXTERN_C_END
5317cb94ee6SHong Zhang 
5327cb94ee6SHong Zhang EXTERN_C_BEGIN
5337cb94ee6SHong Zhang #undef __FUNCT__
5347cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials"
5357087cfbeSBarry Smith PetscErrorCode  TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
5367cb94ee6SHong Zhang {
5377cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5387cb94ee6SHong Zhang 
5397cb94ee6SHong Zhang   PetscFunctionBegin;
5407cb94ee6SHong Zhang   cvode->linear_tol = tol;
5417cb94ee6SHong Zhang   PetscFunctionReturn(0);
5427cb94ee6SHong Zhang }
5437cb94ee6SHong Zhang EXTERN_C_END
5447cb94ee6SHong Zhang 
5457cb94ee6SHong Zhang EXTERN_C_BEGIN
5467cb94ee6SHong Zhang #undef __FUNCT__
5477cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials"
5487087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
5497cb94ee6SHong Zhang {
5507cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5517cb94ee6SHong Zhang 
5527cb94ee6SHong Zhang   PetscFunctionBegin;
5537cb94ee6SHong Zhang   cvode->gtype = type;
5547cb94ee6SHong Zhang   PetscFunctionReturn(0);
5557cb94ee6SHong Zhang }
5567cb94ee6SHong Zhang EXTERN_C_END
5577cb94ee6SHong Zhang 
5587cb94ee6SHong Zhang EXTERN_C_BEGIN
5597cb94ee6SHong Zhang #undef __FUNCT__
5607cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials"
5617087cfbeSBarry Smith PetscErrorCode  TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
5627cb94ee6SHong Zhang {
5637cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5647cb94ee6SHong Zhang 
5657cb94ee6SHong Zhang   PetscFunctionBegin;
5667cb94ee6SHong Zhang   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
5677cb94ee6SHong Zhang   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
5687cb94ee6SHong Zhang   PetscFunctionReturn(0);
5697cb94ee6SHong Zhang }
5707cb94ee6SHong Zhang EXTERN_C_END
5717cb94ee6SHong Zhang 
5727cb94ee6SHong Zhang EXTERN_C_BEGIN
5737cb94ee6SHong Zhang #undef __FUNCT__
574f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials"
5757087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)
576f1cd61daSJed Brown {
577f1cd61daSJed Brown   TS_Sundials *cvode = (TS_Sundials*)ts->data;
578f1cd61daSJed Brown 
579f1cd61daSJed Brown   PetscFunctionBegin;
580f1cd61daSJed Brown   cvode->mindt = mindt;
581f1cd61daSJed Brown   PetscFunctionReturn(0);
582f1cd61daSJed Brown }
583f1cd61daSJed Brown EXTERN_C_END
584f1cd61daSJed Brown 
585f1cd61daSJed Brown EXTERN_C_BEGIN
586f1cd61daSJed Brown #undef __FUNCT__
587f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials"
5887087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)
589f1cd61daSJed Brown {
590f1cd61daSJed Brown   TS_Sundials *cvode = (TS_Sundials*)ts->data;
591f1cd61daSJed Brown 
592f1cd61daSJed Brown   PetscFunctionBegin;
593f1cd61daSJed Brown   cvode->maxdt = maxdt;
594f1cd61daSJed Brown   PetscFunctionReturn(0);
595f1cd61daSJed Brown }
596f1cd61daSJed Brown EXTERN_C_END
597f1cd61daSJed Brown 
598f1cd61daSJed Brown EXTERN_C_BEGIN
599f1cd61daSJed Brown #undef __FUNCT__
6007cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials"
6017087cfbeSBarry Smith PetscErrorCode  TSSundialsGetPC_Sundials(TS ts,PC *pc)
6027cb94ee6SHong Zhang {
603dcbc6d53SSean Farley   SNES            snes;
604dcbc6d53SSean Farley   KSP             ksp;
605dcbc6d53SSean Farley   PetscErrorCode  ierr;
6067cb94ee6SHong Zhang 
6077cb94ee6SHong Zhang   PetscFunctionBegin;
608dcbc6d53SSean Farley   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
609dcbc6d53SSean Farley   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
610dcbc6d53SSean Farley   ierr = KSPGetPC(ksp,pc); CHKERRQ(ierr);
6117cb94ee6SHong Zhang   PetscFunctionReturn(0);
6127cb94ee6SHong Zhang }
6137cb94ee6SHong Zhang EXTERN_C_END
6147cb94ee6SHong Zhang 
6157cb94ee6SHong Zhang EXTERN_C_BEGIN
6167cb94ee6SHong Zhang #undef __FUNCT__
6177cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials"
6187087cfbeSBarry Smith PetscErrorCode  TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
6197cb94ee6SHong Zhang {
6207cb94ee6SHong Zhang   PetscFunctionBegin;
6212c823083SHong Zhang   if (nonlin) *nonlin = ts->nonlinear_its;
6222c823083SHong Zhang   if (lin)    *lin    = ts->linear_its;
6237cb94ee6SHong Zhang   PetscFunctionReturn(0);
6247cb94ee6SHong Zhang }
6257cb94ee6SHong Zhang EXTERN_C_END
6267cb94ee6SHong Zhang 
6277cb94ee6SHong Zhang EXTERN_C_BEGIN
6287cb94ee6SHong Zhang #undef __FUNCT__
6292bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials"
6307087cfbeSBarry Smith PetscErrorCode  TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool  s)
6312bfc04deSHong Zhang {
6322bfc04deSHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
6332bfc04deSHong Zhang 
6342bfc04deSHong Zhang   PetscFunctionBegin;
6352bfc04deSHong Zhang   cvode->monitorstep = s;
6362bfc04deSHong Zhang   PetscFunctionReturn(0);
6372bfc04deSHong Zhang }
6382bfc04deSHong Zhang EXTERN_C_END
6397cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
6407cb94ee6SHong Zhang 
6417cb94ee6SHong Zhang #undef __FUNCT__
6427cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations"
6437cb94ee6SHong Zhang /*@C
6447cb94ee6SHong Zhang    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
6457cb94ee6SHong Zhang 
6467cb94ee6SHong Zhang    Not Collective
6477cb94ee6SHong Zhang 
6487cb94ee6SHong Zhang    Input parameters:
6497cb94ee6SHong Zhang .    ts     - the time-step context
6507cb94ee6SHong Zhang 
6517cb94ee6SHong Zhang    Output Parameters:
6527cb94ee6SHong Zhang +   nonlin - number of nonlinear iterations
6537cb94ee6SHong Zhang -   lin    - number of linear iterations
6547cb94ee6SHong Zhang 
6557cb94ee6SHong Zhang    Level: advanced
6567cb94ee6SHong Zhang 
6577cb94ee6SHong Zhang    Notes:
6587cb94ee6SHong Zhang     These return the number since the creation of the TS object
6597cb94ee6SHong Zhang 
6607cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations
6617cb94ee6SHong Zhang 
6627cb94ee6SHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6637cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
6647cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
665a43b19c4SJed Brown           TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime()
6667cb94ee6SHong Zhang 
6677cb94ee6SHong Zhang @*/
6687087cfbeSBarry Smith PetscErrorCode  TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
6697cb94ee6SHong Zhang {
6704ac538c5SBarry Smith   PetscErrorCode ierr;
6717cb94ee6SHong Zhang 
6727cb94ee6SHong Zhang   PetscFunctionBegin;
6734ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr);
6747cb94ee6SHong Zhang   PetscFunctionReturn(0);
6757cb94ee6SHong Zhang }
6767cb94ee6SHong Zhang 
6777cb94ee6SHong Zhang #undef __FUNCT__
6787cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType"
6797cb94ee6SHong Zhang /*@
6807cb94ee6SHong Zhang    TSSundialsSetType - Sets the method that Sundials will use for integration.
6817cb94ee6SHong Zhang 
6823f9fe445SBarry Smith    Logically Collective on TS
6837cb94ee6SHong Zhang 
6847cb94ee6SHong Zhang    Input parameters:
6857cb94ee6SHong Zhang +    ts     - the time-step context
6867cb94ee6SHong Zhang -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
6877cb94ee6SHong Zhang 
6887cb94ee6SHong Zhang    Level: intermediate
6897cb94ee6SHong Zhang 
6907cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula
6917cb94ee6SHong Zhang 
6927cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(),  TSSundialsSetGMRESRestart(),
6937cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
6947cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6957cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
696a43b19c4SJed Brown           TSSetExactFinalTime()
6977cb94ee6SHong Zhang @*/
6987087cfbeSBarry Smith PetscErrorCode  TSSundialsSetType(TS ts,TSSundialsLmmType type)
6997cb94ee6SHong Zhang {
7004ac538c5SBarry Smith   PetscErrorCode ierr;
7017cb94ee6SHong Zhang 
7027cb94ee6SHong Zhang   PetscFunctionBegin;
7034ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr);
7047cb94ee6SHong Zhang   PetscFunctionReturn(0);
7057cb94ee6SHong Zhang }
7067cb94ee6SHong Zhang 
7077cb94ee6SHong Zhang #undef __FUNCT__
7087cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart"
7097cb94ee6SHong Zhang /*@
7107cb94ee6SHong Zhang    TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by
7117cb94ee6SHong Zhang        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
7127cb94ee6SHong Zhang        this is ALSO the maximum number of GMRES steps that will be used.
7137cb94ee6SHong Zhang 
7143f9fe445SBarry Smith    Logically Collective on TS
7157cb94ee6SHong Zhang 
7167cb94ee6SHong Zhang    Input parameters:
7177cb94ee6SHong Zhang +    ts      - the time-step context
7187cb94ee6SHong Zhang -    restart - number of direction vectors (the restart size).
7197cb94ee6SHong Zhang 
7207cb94ee6SHong Zhang    Level: advanced
7217cb94ee6SHong Zhang 
7227cb94ee6SHong Zhang .keywords: GMRES, restart
7237cb94ee6SHong Zhang 
7247cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
7257cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
7267cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7277cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
728a43b19c4SJed Brown           TSSetExactFinalTime()
7297cb94ee6SHong Zhang 
7307cb94ee6SHong Zhang @*/
7317087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGMRESRestart(TS ts,int restart)
7327cb94ee6SHong Zhang {
7334ac538c5SBarry Smith   PetscErrorCode ierr;
7347cb94ee6SHong Zhang 
7357cb94ee6SHong Zhang   PetscFunctionBegin;
736c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(ts,restart,2);
7374ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetGMRESRestart_C",(TS,int),(ts,restart));CHKERRQ(ierr);
7387cb94ee6SHong Zhang   PetscFunctionReturn(0);
7397cb94ee6SHong Zhang }
7407cb94ee6SHong Zhang 
7417cb94ee6SHong Zhang #undef __FUNCT__
7427cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance"
7437cb94ee6SHong Zhang /*@
7447cb94ee6SHong Zhang    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
7457cb94ee6SHong Zhang        system by SUNDIALS.
7467cb94ee6SHong Zhang 
7473f9fe445SBarry Smith    Logically Collective on TS
7487cb94ee6SHong Zhang 
7497cb94ee6SHong Zhang    Input parameters:
7507cb94ee6SHong Zhang +    ts     - the time-step context
7517cb94ee6SHong Zhang -    tol    - the factor by which the tolerance on the nonlinear solver is
7527cb94ee6SHong Zhang              multiplied to get the tolerance on the linear solver, .05 by default.
7537cb94ee6SHong Zhang 
7547cb94ee6SHong Zhang    Level: advanced
7557cb94ee6SHong Zhang 
7567cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS
7577cb94ee6SHong Zhang 
7587cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7597cb94ee6SHong Zhang           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
7607cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7617cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
762a43b19c4SJed Brown           TSSetExactFinalTime()
7637cb94ee6SHong Zhang 
7647cb94ee6SHong Zhang @*/
7657087cfbeSBarry Smith PetscErrorCode  TSSundialsSetLinearTolerance(TS ts,double tol)
7667cb94ee6SHong Zhang {
7674ac538c5SBarry Smith   PetscErrorCode ierr;
7687cb94ee6SHong Zhang 
7697cb94ee6SHong Zhang   PetscFunctionBegin;
770c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(ts,tol,2);
7714ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr);
7727cb94ee6SHong Zhang   PetscFunctionReturn(0);
7737cb94ee6SHong Zhang }
7747cb94ee6SHong Zhang 
7757cb94ee6SHong Zhang #undef __FUNCT__
7767cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType"
7777cb94ee6SHong Zhang /*@
7787cb94ee6SHong Zhang    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
7797cb94ee6SHong Zhang         in GMRES method by SUNDIALS linear solver.
7807cb94ee6SHong Zhang 
7813f9fe445SBarry Smith    Logically Collective on TS
7827cb94ee6SHong Zhang 
7837cb94ee6SHong Zhang    Input parameters:
7847cb94ee6SHong Zhang +    ts  - the time-step context
7857cb94ee6SHong Zhang -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
7867cb94ee6SHong Zhang 
7877cb94ee6SHong Zhang    Level: advanced
7887cb94ee6SHong Zhang 
7897cb94ee6SHong Zhang .keywords: Sundials, orthogonalization
7907cb94ee6SHong Zhang 
7917cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7927cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
7937cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7947cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
795a43b19c4SJed Brown           TSSetExactFinalTime()
7967cb94ee6SHong Zhang 
7977cb94ee6SHong Zhang @*/
7987087cfbeSBarry Smith PetscErrorCode  TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
7997cb94ee6SHong Zhang {
8004ac538c5SBarry Smith   PetscErrorCode ierr;
8017cb94ee6SHong Zhang 
8027cb94ee6SHong Zhang   PetscFunctionBegin;
8034ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr);
8047cb94ee6SHong Zhang   PetscFunctionReturn(0);
8057cb94ee6SHong Zhang }
8067cb94ee6SHong Zhang 
8077cb94ee6SHong Zhang #undef __FUNCT__
8087cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance"
8097cb94ee6SHong Zhang /*@
8107cb94ee6SHong Zhang    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
8117cb94ee6SHong Zhang                          Sundials for error control.
8127cb94ee6SHong Zhang 
8133f9fe445SBarry Smith    Logically Collective on TS
8147cb94ee6SHong Zhang 
8157cb94ee6SHong Zhang    Input parameters:
8167cb94ee6SHong Zhang +    ts  - the time-step context
8177cb94ee6SHong Zhang .    aabs - the absolute tolerance
8187cb94ee6SHong Zhang -    rel - the relative tolerance
8197cb94ee6SHong Zhang 
8207cb94ee6SHong Zhang      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
8217cb94ee6SHong Zhang     these regulate the size of the error for a SINGLE timestep.
8227cb94ee6SHong Zhang 
8237cb94ee6SHong Zhang    Level: intermediate
8247cb94ee6SHong Zhang 
8257cb94ee6SHong Zhang .keywords: Sundials, tolerance
8267cb94ee6SHong Zhang 
8277cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8287cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
8297cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8307cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
831a43b19c4SJed Brown           TSSetExactFinalTime()
8327cb94ee6SHong Zhang 
8337cb94ee6SHong Zhang @*/
8347087cfbeSBarry Smith PetscErrorCode  TSSundialsSetTolerance(TS ts,double aabs,double rel)
8357cb94ee6SHong Zhang {
8364ac538c5SBarry Smith   PetscErrorCode ierr;
8377cb94ee6SHong Zhang 
8387cb94ee6SHong Zhang   PetscFunctionBegin;
8394ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr);
8407cb94ee6SHong Zhang   PetscFunctionReturn(0);
8417cb94ee6SHong Zhang }
8427cb94ee6SHong Zhang 
8437cb94ee6SHong Zhang #undef __FUNCT__
8447cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC"
8457cb94ee6SHong Zhang /*@
8467cb94ee6SHong Zhang    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
8477cb94ee6SHong Zhang 
8487cb94ee6SHong Zhang    Input Parameter:
8497cb94ee6SHong Zhang .    ts - the time-step context
8507cb94ee6SHong Zhang 
8517cb94ee6SHong Zhang    Output Parameter:
8527cb94ee6SHong Zhang .    pc - the preconditioner context
8537cb94ee6SHong Zhang 
8547cb94ee6SHong Zhang    Level: advanced
8557cb94ee6SHong Zhang 
8567cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8577cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
8587cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8597cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
8607cb94ee6SHong Zhang @*/
8617087cfbeSBarry Smith PetscErrorCode  TSSundialsGetPC(TS ts,PC *pc)
8627cb94ee6SHong Zhang {
8634ac538c5SBarry Smith   PetscErrorCode ierr;
8647cb94ee6SHong Zhang 
8657cb94ee6SHong Zhang   PetscFunctionBegin;
8664ac538c5SBarry Smith   ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC *),(ts,pc));CHKERRQ(ierr);
8677cb94ee6SHong Zhang   PetscFunctionReturn(0);
8687cb94ee6SHong Zhang }
8697cb94ee6SHong Zhang 
8707cb94ee6SHong Zhang #undef __FUNCT__
871f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMinTimeStep"
872f1cd61daSJed Brown /*@
873f1cd61daSJed Brown    TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller.
874f1cd61daSJed Brown 
875f1cd61daSJed Brown    Input Parameter:
876f1cd61daSJed Brown +   ts - the time-step context
877f1cd61daSJed Brown -   mindt - lowest time step if positive, negative to deactivate
878f1cd61daSJed Brown 
879fc6b9e64SJed Brown    Note:
880fc6b9e64SJed Brown    Sundials will error if it is not possible to keep the estimated truncation error below
881fc6b9e64SJed Brown    the tolerance set with TSSundialsSetTolerance() without going below this step size.
882fc6b9e64SJed Brown 
883f1cd61daSJed Brown    Level: beginner
884f1cd61daSJed Brown 
885f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
886f1cd61daSJed Brown @*/
8877087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMinTimeStep(TS ts,PetscReal mindt)
888f1cd61daSJed Brown {
8894ac538c5SBarry Smith   PetscErrorCode ierr;
890f1cd61daSJed Brown 
891f1cd61daSJed Brown   PetscFunctionBegin;
8924ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr);
893f1cd61daSJed Brown   PetscFunctionReturn(0);
894f1cd61daSJed Brown }
895f1cd61daSJed Brown 
896f1cd61daSJed Brown #undef __FUNCT__
897f1cd61daSJed Brown #define __FUNCT__ "TSSundialsSetMaxTimeStep"
898f1cd61daSJed Brown /*@
899f1cd61daSJed Brown    TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller.
900f1cd61daSJed Brown 
901f1cd61daSJed Brown    Input Parameter:
902f1cd61daSJed Brown +   ts - the time-step context
903f1cd61daSJed Brown -   maxdt - lowest time step if positive, negative to deactivate
904f1cd61daSJed Brown 
905f1cd61daSJed Brown    Level: beginner
906f1cd61daSJed Brown 
907f1cd61daSJed Brown .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
908f1cd61daSJed Brown @*/
9097087cfbeSBarry Smith PetscErrorCode  TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt)
910f1cd61daSJed Brown {
9114ac538c5SBarry Smith   PetscErrorCode ierr;
912f1cd61daSJed Brown 
913f1cd61daSJed Brown   PetscFunctionBegin;
9144ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr);
915f1cd61daSJed Brown   PetscFunctionReturn(0);
916f1cd61daSJed Brown }
917f1cd61daSJed Brown 
918f1cd61daSJed Brown #undef __FUNCT__
9192bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps"
9202bfc04deSHong Zhang /*@
9212bfc04deSHong Zhang    TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false).
9222bfc04deSHong Zhang 
9232bfc04deSHong Zhang    Input Parameter:
9242bfc04deSHong Zhang +   ts - the time-step context
9252bfc04deSHong Zhang -   ft - PETSC_TRUE if monitor, else PETSC_FALSE
9262bfc04deSHong Zhang 
9272bfc04deSHong Zhang    Level: beginner
9282bfc04deSHong Zhang 
9292bfc04deSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
9302bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
9312bfc04deSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
9322bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
9332bfc04deSHong Zhang @*/
9347087cfbeSBarry Smith PetscErrorCode  TSSundialsMonitorInternalSteps(TS ts,PetscBool  ft)
9352bfc04deSHong Zhang {
9364ac538c5SBarry Smith   PetscErrorCode ierr;
9372bfc04deSHong Zhang 
9382bfc04deSHong Zhang   PetscFunctionBegin;
9394ac538c5SBarry Smith   ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr);
9402bfc04deSHong Zhang   PetscFunctionReturn(0);
9412bfc04deSHong Zhang }
9427cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
9437cb94ee6SHong Zhang /*MC
94496f5712cSJed Brown       TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
9457cb94ee6SHong Zhang 
9467cb94ee6SHong Zhang    Options Database:
9477cb94ee6SHong Zhang +    -ts_sundials_type <bdf,adams>
9487cb94ee6SHong Zhang .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
9497cb94ee6SHong Zhang .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
9507cb94ee6SHong Zhang .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
9517cb94ee6SHong Zhang .    -ts_sundials_linear_tolerance <tol>
9527cb94ee6SHong Zhang .    -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions
95316016685SHong Zhang -    -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps
95416016685SHong Zhang 
9557cb94ee6SHong Zhang 
9567cb94ee6SHong Zhang     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
9577cb94ee6SHong Zhang            only PETSc PC options
9587cb94ee6SHong Zhang 
9597cb94ee6SHong Zhang     Level: beginner
9607cb94ee6SHong Zhang 
9617cb94ee6SHong Zhang .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(),
962a43b19c4SJed Brown            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime()
9637cb94ee6SHong Zhang 
9647cb94ee6SHong Zhang M*/
9657cb94ee6SHong Zhang EXTERN_C_BEGIN
9667cb94ee6SHong Zhang #undef __FUNCT__
9677cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials"
9687087cfbeSBarry Smith PetscErrorCode  TSCreate_Sundials(TS ts)
9697cb94ee6SHong Zhang {
9707cb94ee6SHong Zhang   TS_Sundials    *cvode;
9717cb94ee6SHong Zhang   PetscErrorCode ierr;
97242528757SHong Zhang   PC             pc;
9737cb94ee6SHong Zhang 
9747cb94ee6SHong Zhang   PetscFunctionBegin;
975277b19d0SLisandro Dalcin   ts->ops->reset          = TSReset_Sundials;
97628adb3f7SHong Zhang   ts->ops->destroy        = TSDestroy_Sundials;
97728adb3f7SHong Zhang   ts->ops->view           = TSView_Sundials;
978214bc6a2SJed Brown   ts->ops->setup          = TSSetUp_Sundials;
979*b4eba00bSSean Farley   ts->ops->step           = TSStep_Sundials;
980*b4eba00bSSean Farley   ts->ops->interpolate    = TSInterpolate_Sundials;
981214bc6a2SJed Brown   ts->ops->setfromoptions = TSSetFromOptions_Sundials;
9827cb94ee6SHong Zhang 
98338f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr);
9847cb94ee6SHong Zhang   ts->data                = (void*)cvode;
9856fadb2cdSHong Zhang   cvode->cvode_type       = SUNDIALS_BDF;
9866fadb2cdSHong Zhang   cvode->gtype            = SUNDIALS_CLASSICAL_GS;
9877cb94ee6SHong Zhang   cvode->restart          = 5;
9887cb94ee6SHong Zhang   cvode->linear_tol       = .05;
9897cb94ee6SHong Zhang 
990*b4eba00bSSean Farley   cvode->monitorstep      = PETSC_TRUE;
9917cb94ee6SHong Zhang 
992a07356b8SSatish Balay   ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr);
993f1cd61daSJed Brown 
994f1cd61daSJed Brown   cvode->mindt = -1.;
995f1cd61daSJed Brown   cvode->maxdt = -1.;
996f1cd61daSJed Brown 
9977cb94ee6SHong Zhang   /* set tolerance for Sundials */
9987cb94ee6SHong Zhang   cvode->reltol = 1e-6;
9992c823083SHong Zhang   cvode->abstol = 1e-6;
10007cb94ee6SHong Zhang 
100142528757SHong Zhang   /* set PCNONE as default pctype */
100242528757SHong Zhang   ierr = TSSundialsGetPC_Sundials(ts,&pc);CHKERRQ(ierr);
100342528757SHong Zhang   ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
100442528757SHong Zhang 
1005*b4eba00bSSean Farley   if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_FALSE;
1006a43b19c4SJed Brown 
10077cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
10087cb94ee6SHong Zhang                     TSSundialsSetType_Sundials);CHKERRQ(ierr);
10097cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C",
10107cb94ee6SHong Zhang                     "TSSundialsSetGMRESRestart_Sundials",
10117cb94ee6SHong Zhang                     TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr);
10127cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
10137cb94ee6SHong Zhang                     "TSSundialsSetLinearTolerance_Sundials",
10147cb94ee6SHong Zhang                      TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
10157cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
10167cb94ee6SHong Zhang                     "TSSundialsSetGramSchmidtType_Sundials",
10177cb94ee6SHong Zhang                      TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
10187cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
10197cb94ee6SHong Zhang                     "TSSundialsSetTolerance_Sundials",
10207cb94ee6SHong Zhang                      TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
1021f1cd61daSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C",
1022f1cd61daSJed Brown                     "TSSundialsSetMinTimeStep_Sundials",
1023f1cd61daSJed Brown                      TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr);
1024f1cd61daSJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",
1025f1cd61daSJed Brown                     "TSSundialsSetMaxTimeStep_Sundials",
1026f1cd61daSJed Brown                      TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr);
10277cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
10287cb94ee6SHong Zhang                     "TSSundialsGetPC_Sundials",
10297cb94ee6SHong Zhang                      TSSundialsGetPC_Sundials);CHKERRQ(ierr);
10307cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
10317cb94ee6SHong Zhang                     "TSSundialsGetIterations_Sundials",
10327cb94ee6SHong Zhang                      TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
10332bfc04deSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",
10342bfc04deSHong Zhang                     "TSSundialsMonitorInternalSteps_Sundials",
10352bfc04deSHong Zhang                      TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr);
10362bfc04deSHong Zhang 
10377cb94ee6SHong Zhang   PetscFunctionReturn(0);
10387cb94ee6SHong Zhang }
10397cb94ee6SHong Zhang EXTERN_C_END
1040