xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision 96f5712c7cba09c0c28604a0f5a443f4ccbe148b)
17cb94ee6SHong Zhang #define PETSCTS_DLL
27cb94ee6SHong Zhang 
37cb94ee6SHong Zhang /*
4db268bfaSHong Zhang     Provides a PETSc interface to SUNDIALS/CVODE solver.
57cb94ee6SHong Zhang     The interface to PVODE (old version of CVODE) was originally contributed
67cb94ee6SHong Zhang     by Liyang Xu. It has been redone by Hong Zhang and Dinesh Kaushik.
728adb3f7SHong Zhang 
828adb3f7SHong Zhang     Reference: sundials-2.4.0/examples/cvode/parallel/cvDiurnal_kry_p.c
97cb94ee6SHong Zhang */
1056a740aaSHong Zhang #include "sundials.h"  /*I "petscts.h" I*/
117cb94ee6SHong Zhang 
127cb94ee6SHong Zhang /*
137cb94ee6SHong Zhang       TSPrecond_Sundials - function that we provide to SUNDIALS to
147cb94ee6SHong Zhang                         evaluate the preconditioner.
157cb94ee6SHong Zhang */
167cb94ee6SHong Zhang #undef __FUNCT__
177cb94ee6SHong Zhang #define __FUNCT__ "TSPrecond_Sundials"
187cb94ee6SHong Zhang PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy,
197cb94ee6SHong Zhang                     booleantype jok,booleantype *jcurPtr,
207cb94ee6SHong Zhang                     realtype _gamma,void *P_data,
217cb94ee6SHong Zhang                     N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3)
227cb94ee6SHong Zhang {
237cb94ee6SHong Zhang   TS             ts = (TS) P_data;
247cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
257cb94ee6SHong Zhang   PC             pc = cvode->pc;
267cb94ee6SHong Zhang   PetscErrorCode ierr;
277cb94ee6SHong Zhang   Mat            Jac = ts->B;
287cb94ee6SHong Zhang   Vec            yy = cvode->w1;
297cb94ee6SHong Zhang   PetscScalar    one = 1.0,gm;
307cb94ee6SHong Zhang   MatStructure   str = DIFFERENT_NONZERO_PATTERN;
317cb94ee6SHong Zhang   PetscScalar    *y_data;
327cb94ee6SHong Zhang 
337cb94ee6SHong Zhang   PetscFunctionBegin;
347cb94ee6SHong Zhang   /* This allows us to construct preconditioners in-place if we like */
357cb94ee6SHong Zhang   ierr = MatSetUnfactored(Jac);CHKERRQ(ierr);
367cb94ee6SHong Zhang 
377cb94ee6SHong Zhang   /* jok - TRUE means reuse current Jacobian else recompute Jacobian */
387cb94ee6SHong Zhang   if (jok) {
397cb94ee6SHong Zhang     ierr     = MatCopy(cvode->pmat,Jac,str);CHKERRQ(ierr);
407cb94ee6SHong Zhang     *jcurPtr = FALSE;
417cb94ee6SHong Zhang   } else {
427cb94ee6SHong Zhang     /* make PETSc vector yy point to SUNDIALS vector y */
437cb94ee6SHong Zhang     y_data = (PetscScalar *) N_VGetArrayPointer(y);
447cb94ee6SHong Zhang     ierr   = VecPlaceArray(yy,y_data); CHKERRQ(ierr);
454ac7836bSHong Zhang 
467cb94ee6SHong Zhang     /* compute the Jacobian */
477cb94ee6SHong Zhang     ierr = TSComputeRHSJacobian(ts,ts->ptime,yy,&Jac,&Jac,&str);CHKERRQ(ierr);
487cb94ee6SHong Zhang     ierr = VecResetArray(yy); CHKERRQ(ierr);
494ac7836bSHong Zhang 
507cb94ee6SHong Zhang     /* copy the Jacobian matrix */
517cb94ee6SHong Zhang     if (!cvode->pmat) {
527cb94ee6SHong Zhang       ierr = MatDuplicate(Jac,MAT_COPY_VALUES,&cvode->pmat);CHKERRQ(ierr);
537cb94ee6SHong Zhang       ierr = PetscLogObjectParent(ts,cvode->pmat);CHKERRQ(ierr);
542c823083SHong Zhang     } else {
557cb94ee6SHong Zhang       ierr = MatCopy(Jac,cvode->pmat,str);CHKERRQ(ierr);
567cb94ee6SHong Zhang     }
577cb94ee6SHong Zhang     *jcurPtr = TRUE;
587cb94ee6SHong Zhang   }
597cb94ee6SHong Zhang 
607cb94ee6SHong Zhang   /* construct I-gamma*Jac  */
617cb94ee6SHong Zhang   gm   = -_gamma;
627cb94ee6SHong Zhang   ierr = MatScale(Jac,gm);CHKERRQ(ierr);
637cb94ee6SHong Zhang   ierr = MatShift(Jac,one);CHKERRQ(ierr);
647cb94ee6SHong Zhang 
657cb94ee6SHong Zhang   ierr = PCSetOperators(pc,Jac,Jac,str);CHKERRQ(ierr);
667cb94ee6SHong Zhang   PetscFunctionReturn(0);
677cb94ee6SHong Zhang }
687cb94ee6SHong Zhang 
697cb94ee6SHong Zhang /*
707cb94ee6SHong Zhang      TSPSolve_Sundials -  routine that we provide to Sundials that applies the preconditioner.
717cb94ee6SHong Zhang */
727cb94ee6SHong Zhang #undef __FUNCT__
737cb94ee6SHong Zhang #define __FUNCT__ "TSPSolve_Sundials"
744ac7836bSHong Zhang PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z,
754ac7836bSHong Zhang                                  realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp)
767cb94ee6SHong Zhang {
777cb94ee6SHong Zhang   TS              ts = (TS) P_data;
787cb94ee6SHong Zhang   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
797cb94ee6SHong Zhang   PC              pc = cvode->pc;
807cb94ee6SHong Zhang   Vec             rr = cvode->w1,zz = cvode->w2;
817cb94ee6SHong Zhang   PetscErrorCode  ierr;
827cb94ee6SHong Zhang   PetscScalar     *r_data,*z_data;
837cb94ee6SHong Zhang 
847cb94ee6SHong Zhang   PetscFunctionBegin;
857cb94ee6SHong Zhang   /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
867cb94ee6SHong Zhang   r_data  = (PetscScalar *) N_VGetArrayPointer(r);
877cb94ee6SHong Zhang   z_data  = (PetscScalar *) N_VGetArrayPointer(z);
887cb94ee6SHong Zhang   ierr = VecPlaceArray(rr,r_data); CHKERRQ(ierr);
897cb94ee6SHong Zhang   ierr = VecPlaceArray(zz,z_data); CHKERRQ(ierr);
904ac7836bSHong Zhang 
917cb94ee6SHong Zhang   /* Solve the Px=r and put the result in zz */
927cb94ee6SHong Zhang   ierr = PCApply(pc,rr,zz); CHKERRQ(ierr);
937cb94ee6SHong Zhang   ierr = VecResetArray(rr); CHKERRQ(ierr);
947cb94ee6SHong Zhang   ierr = VecResetArray(zz); CHKERRQ(ierr);
957cb94ee6SHong Zhang   PetscFunctionReturn(0);
967cb94ee6SHong Zhang }
977cb94ee6SHong Zhang 
987cb94ee6SHong Zhang /*
997cb94ee6SHong Zhang         TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
1007cb94ee6SHong Zhang */
1017cb94ee6SHong Zhang #undef __FUNCT__
1027cb94ee6SHong Zhang #define __FUNCT__ "TSFunction_Sundials"
1034ac7836bSHong Zhang int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx)
1047cb94ee6SHong Zhang {
1057cb94ee6SHong Zhang   TS              ts = (TS) ctx;
1067adad957SLisandro Dalcin   MPI_Comm        comm = ((PetscObject)ts)->comm;
1077cb94ee6SHong Zhang   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
1087cb94ee6SHong Zhang   Vec             yy = cvode->w1,yyd = cvode->w2;
1097cb94ee6SHong Zhang   PetscScalar     *y_data,*ydot_data;
1107cb94ee6SHong Zhang   PetscErrorCode  ierr;
1117cb94ee6SHong Zhang 
1127cb94ee6SHong Zhang   PetscFunctionBegin;
1135ff03cf7SDinesh Kaushik   /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
1147cb94ee6SHong Zhang   y_data     = (PetscScalar *) N_VGetArrayPointer(y);
1157cb94ee6SHong Zhang   ydot_data  = (PetscScalar *) N_VGetArrayPointer(ydot);
1164ac7836bSHong Zhang   ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr);
1174ac7836bSHong Zhang   ierr = VecPlaceArray(yyd,ydot_data); CHKERRABORT(comm,ierr);
1184ac7836bSHong Zhang 
1197cb94ee6SHong Zhang   /* now compute the right hand side function */
1207cb94ee6SHong Zhang   ierr = TSComputeRHSFunction(ts,t,yy,yyd); CHKERRABORT(comm,ierr);
1217cb94ee6SHong Zhang   ierr = VecResetArray(yy); CHKERRABORT(comm,ierr);
1227cb94ee6SHong Zhang   ierr = VecResetArray(yyd); CHKERRABORT(comm,ierr);
1234ac7836bSHong Zhang   PetscFunctionReturn(0);
1247cb94ee6SHong Zhang }
1257cb94ee6SHong Zhang 
1267cb94ee6SHong Zhang /*
1277cb94ee6SHong Zhang        TSStep_Sundials_Nonlinear - Calls Sundials to integrate the ODE.
1287cb94ee6SHong Zhang */
1297cb94ee6SHong Zhang #undef __FUNCT__
1307cb94ee6SHong Zhang #define __FUNCT__ "TSStep_Sundials_Nonlinear"
1317cb94ee6SHong Zhang PetscErrorCode TSStep_Sundials_Nonlinear(TS ts,int *steps,double *time)
1327cb94ee6SHong Zhang {
1337cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
1347cb94ee6SHong Zhang   Vec            sol = ts->vec_sol;
1357cb94ee6SHong Zhang   PetscErrorCode ierr;
1365d47aee6SHong Zhang   PetscInt       i,max_steps = ts->max_steps,flag;
1377cb94ee6SHong Zhang   long int       its;
1387cb94ee6SHong Zhang   realtype       t,tout;
1397cb94ee6SHong Zhang   PetscScalar    *y_data;
1407cb94ee6SHong Zhang   void           *mem;
1417cb94ee6SHong Zhang 
1427cb94ee6SHong Zhang   PetscFunctionBegin;
14316016685SHong Zhang   mem  = cvode->mem;
1447cb94ee6SHong Zhang   tout = ts->max_time;
1457cb94ee6SHong Zhang   ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr);
1467cb94ee6SHong Zhang   N_VSetArrayPointer((realtype *)y_data,cvode->y);
1477cb94ee6SHong Zhang   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
1487cb94ee6SHong Zhang   for (i = 0; i < max_steps; i++) {
1497cb94ee6SHong Zhang     if (ts->ptime >= ts->max_time) break;
1502bfc04deSHong Zhang     if (cvode->monitorstep){
15128adb3f7SHong Zhang       flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
1522bfc04deSHong Zhang     } else {
1532bfc04deSHong Zhang       flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL);
1542bfc04deSHong Zhang     }
155ce8aa302SHong Zhang     if (flag)SETERRQ1(PETSC_ERR_LIB,"CVode() fails, flag %d",flag);
1567cb94ee6SHong Zhang     if (t > ts->max_time && cvode->exact_final_time) {
1577cb94ee6SHong Zhang       /* interpolate to final requested time */
1587cb94ee6SHong Zhang       ierr = CVodeGetDky(mem,tout,0,cvode->y);CHKERRQ(ierr);
1597cb94ee6SHong Zhang       t = tout;
1607cb94ee6SHong Zhang     }
1617cb94ee6SHong Zhang     ts->time_step = t - ts->ptime;
1627cb94ee6SHong Zhang     ts->ptime     = t;
1637cb94ee6SHong Zhang 
1647cb94ee6SHong Zhang     /* copy the solution from cvode->y to cvode->update and sol */
1657cb94ee6SHong Zhang     ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr);
1667cb94ee6SHong Zhang     ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr);
1677cb94ee6SHong Zhang     ierr = VecResetArray(cvode->w1); CHKERRQ(ierr);
1687cb94ee6SHong Zhang     ierr = VecCopy(cvode->update,sol);CHKERRQ(ierr);
1697cb94ee6SHong Zhang     ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr);
1707cb94ee6SHong Zhang     ts->nonlinear_its = its;
1714ac7836bSHong Zhang     ierr = CVSpilsGetNumLinIters(mem, &its);
1727cb94ee6SHong Zhang     ts->linear_its = its;
1737cb94ee6SHong Zhang     ts->steps++;
1747cb94ee6SHong Zhang     ierr = TSMonitor(ts,ts->steps,t,sol);CHKERRQ(ierr);
1757cb94ee6SHong Zhang   }
1767cb94ee6SHong Zhang   *steps += ts->steps;
1777cb94ee6SHong Zhang   *time   = t;
1787cb94ee6SHong Zhang   PetscFunctionReturn(0);
1797cb94ee6SHong Zhang }
1807cb94ee6SHong Zhang 
1817cb94ee6SHong Zhang #undef __FUNCT__
1827cb94ee6SHong Zhang #define __FUNCT__ "TSDestroy_Sundials"
1837cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts)
1847cb94ee6SHong Zhang {
1857cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
1867cb94ee6SHong Zhang   PetscErrorCode ierr;
1877cb94ee6SHong Zhang 
1887cb94ee6SHong Zhang   PetscFunctionBegin;
1897cb94ee6SHong Zhang   if (cvode->pmat)   {ierr = MatDestroy(cvode->pmat);CHKERRQ(ierr);}
1907cb94ee6SHong Zhang   if (cvode->pc)     {ierr = PCDestroy(cvode->pc);CHKERRQ(ierr);}
1917cb94ee6SHong Zhang   if (cvode->update) {ierr = VecDestroy(cvode->update);CHKERRQ(ierr);}
1927cb94ee6SHong Zhang   if (cvode->func)   {ierr = VecDestroy(cvode->func);CHKERRQ(ierr);}
1937cb94ee6SHong Zhang   if (cvode->rhs)    {ierr = VecDestroy(cvode->rhs);CHKERRQ(ierr);}
1947cb94ee6SHong Zhang   if (cvode->w1)     {ierr = VecDestroy(cvode->w1);CHKERRQ(ierr);}
1957cb94ee6SHong Zhang   if (cvode->w2)     {ierr = VecDestroy(cvode->w2);CHKERRQ(ierr);}
196a07356b8SSatish Balay   ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr);
1972c823083SHong Zhang   if (cvode->mem) {CVodeFree(&cvode->mem);}
1987cb94ee6SHong Zhang   ierr = PetscFree(cvode);CHKERRQ(ierr);
1997cb94ee6SHong Zhang   PetscFunctionReturn(0);
2007cb94ee6SHong Zhang }
2017cb94ee6SHong Zhang 
2027cb94ee6SHong Zhang #undef __FUNCT__
2037cb94ee6SHong Zhang #define __FUNCT__ "TSSetUp_Sundials_Nonlinear"
2047cb94ee6SHong Zhang PetscErrorCode TSSetUp_Sundials_Nonlinear(TS ts)
2057cb94ee6SHong Zhang {
2067cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2077cb94ee6SHong Zhang   PetscErrorCode ierr;
20816016685SHong Zhang   PetscInt       glosize,locsize,i,flag;
2097cb94ee6SHong Zhang   PetscScalar    *y_data,*parray;
21016016685SHong Zhang   void           *mem;
21116016685SHong Zhang   const PCType   pctype;
21216016685SHong Zhang   PetscTruth     pcnone;
21316016685SHong Zhang   Vec            sol = ts->vec_sol;
2147cb94ee6SHong Zhang 
2157cb94ee6SHong Zhang   PetscFunctionBegin;
2167cb94ee6SHong Zhang   ierr = PCSetFromOptions(cvode->pc);CHKERRQ(ierr);
2177cb94ee6SHong Zhang   /* get the vector size */
2187cb94ee6SHong Zhang   ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr);
2197cb94ee6SHong Zhang   ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr);
2207cb94ee6SHong Zhang 
2217cb94ee6SHong Zhang   /* allocate the memory for N_Vec y */
222a07356b8SSatish Balay   cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
2237cb94ee6SHong Zhang   if (!cvode->y) SETERRQ(1,"cvode->y is not allocated");
2247cb94ee6SHong Zhang 
22528adb3f7SHong Zhang   /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
2267cb94ee6SHong Zhang   ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr);
2277cb94ee6SHong Zhang   y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y);
2287cb94ee6SHong Zhang   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
2297cb94ee6SHong Zhang   /*ierr = PetscMemcpy(y_data,parray,locsize*sizeof(PETSC_SCALAR)); CHKERRQ(ierr);*/
2307cb94ee6SHong Zhang   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
2317cb94ee6SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr);
2327cb94ee6SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&cvode->func);CHKERRQ(ierr);
2337cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr);
2347cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->func);CHKERRQ(ierr);
2357cb94ee6SHong Zhang 
2367cb94ee6SHong Zhang   /*
2377cb94ee6SHong Zhang     Create work vectors for the TSPSolve_Sundials() routine. Note these are
2387cb94ee6SHong Zhang     allocated with zero space arrays because the actual array space is provided
2397cb94ee6SHong Zhang     by Sundials and set using VecPlaceArray().
2407cb94ee6SHong Zhang   */
2417adad957SLisandro Dalcin   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr);
2427adad957SLisandro Dalcin   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr);
2437cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr);
2447cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr);
24516016685SHong Zhang 
24616016685SHong Zhang   /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
24716016685SHong Zhang   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
24816016685SHong Zhang   if (!mem) SETERRQ(PETSC_ERR_MEM,"CVodeCreate() fails");
24916016685SHong Zhang   cvode->mem = mem;
25016016685SHong Zhang 
25116016685SHong Zhang   /* Set the pointer to user-defined data */
25216016685SHong Zhang   flag = CVodeSetUserData(mem, ts);
253ce8aa302SHong Zhang   if (flag) SETERRQ(PETSC_ERR_LIB,"CVodeSetUserData() fails");
25416016685SHong Zhang 
25516016685SHong Zhang   /* Call CVodeInit to initialize the integrator memory and specify the
25616016685SHong Zhang    * user's right hand side function in u'=f(t,u), the inital time T0, and
25716016685SHong Zhang    * the initial dependent variable vector cvode->y */
25816016685SHong Zhang   flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y);
25916016685SHong Zhang   if (flag){
260ce8aa302SHong Zhang     SETERRQ1(PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag);
26116016685SHong Zhang   }
26216016685SHong Zhang 
26316016685SHong Zhang   flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol);
26416016685SHong Zhang   if (flag){
265ce8aa302SHong Zhang     SETERRQ1(PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag);
26616016685SHong Zhang   }
26716016685SHong Zhang 
26816016685SHong Zhang   /* initialize the number of steps */
26916016685SHong Zhang   ierr   = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
27016016685SHong Zhang 
27116016685SHong Zhang   /* call CVSpgmr to use GMRES as the linear solver.        */
27216016685SHong Zhang   /* setup the ode integrator with the given preconditioner */
27316016685SHong Zhang   ierr = PCGetType(cvode->pc,&pctype);CHKERRQ(ierr);
27416016685SHong Zhang   ierr = PetscTypeCompare((PetscObject)cvode->pc,PCNONE,&pcnone);CHKERRQ(ierr);
27516016685SHong Zhang   if (pcnone){
27616016685SHong Zhang     flag  = CVSpgmr(mem,PREC_NONE,0);
277ce8aa302SHong Zhang     if (flag) SETERRQ1(PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
27816016685SHong Zhang   } else {
27916016685SHong Zhang     flag  = CVSpgmr(mem,PREC_LEFT,0);
280ce8aa302SHong Zhang     if (flag) SETERRQ1(PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
28116016685SHong Zhang 
28216016685SHong Zhang     /* Set preconditioner and solve routines Precond and PSolve,
28316016685SHong Zhang      and the pointer to the user-defined block data */
28416016685SHong Zhang     flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials);
285ce8aa302SHong Zhang     if (flag) SETERRQ1(PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag);
28616016685SHong Zhang   }
28716016685SHong Zhang 
28816016685SHong Zhang   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
289ce8aa302SHong Zhang   if (flag) SETERRQ1(PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag);
2907cb94ee6SHong Zhang   PetscFunctionReturn(0);
2917cb94ee6SHong Zhang }
2927cb94ee6SHong Zhang 
2936fadb2cdSHong Zhang /* type of CVODE linear multistep method */
294dfcd32c8SHong Zhang const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0};
2956fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */
296dfcd32c8SHong Zhang const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0};
297a04cf4d8SBarry Smith 
2987cb94ee6SHong Zhang #undef __FUNCT__
2997cb94ee6SHong Zhang #define __FUNCT__ "TSSetFromOptions_Sundials_Nonlinear"
3007cb94ee6SHong Zhang PetscErrorCode TSSetFromOptions_Sundials_Nonlinear(TS ts)
3017cb94ee6SHong Zhang {
3027cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
3037cb94ee6SHong Zhang   PetscErrorCode ierr;
3047cb94ee6SHong Zhang   int            indx;
3057cb94ee6SHong Zhang   PetscTruth     flag;
3067cb94ee6SHong Zhang 
3077cb94ee6SHong Zhang   PetscFunctionBegin;
3087cb94ee6SHong Zhang   ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr);
3096fadb2cdSHong Zhang     ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr);
3107cb94ee6SHong Zhang     if (flag) {
3116fadb2cdSHong Zhang       ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr);
3127cb94ee6SHong Zhang     }
3136fadb2cdSHong Zhang     ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr);
3147cb94ee6SHong Zhang     if (flag) {
3157cb94ee6SHong Zhang       ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr);
3167cb94ee6SHong Zhang     }
3177cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr);
3187cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr);
3197cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr);
3207cb94ee6SHong Zhang     ierr = PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);CHKERRQ(ierr);
32116016685SHong Zhang     ierr = PetscOptionsTruth("-ts_sundials_exact_final_time","Interpolate output to stop exactly at the final time","TSSundialsSetExactFinalTime",cvode->exact_final_time,&cvode->exact_final_time,PETSC_NULL);CHKERRQ(ierr);
3222bfc04deSHong Zhang     ierr = PetscOptionsTruth("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);CHKERRQ(ierr);
3237cb94ee6SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
3247cb94ee6SHong Zhang   PetscFunctionReturn(0);
3257cb94ee6SHong Zhang }
3267cb94ee6SHong Zhang 
3277cb94ee6SHong Zhang #undef __FUNCT__
3287cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials"
3297cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
3307cb94ee6SHong Zhang {
3317cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
3327cb94ee6SHong Zhang   PetscErrorCode ierr;
3337cb94ee6SHong Zhang   char           *type;
3347cb94ee6SHong Zhang   char           atype[] = "Adams";
3357cb94ee6SHong Zhang   char           btype[] = "BDF: backward differentiation formula";
3367cb94ee6SHong Zhang   PetscTruth     iascii,isstring;
3372c823083SHong Zhang   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
3382c823083SHong Zhang   PetscInt       qlast,qcur;
3395d47aee6SHong Zhang   PetscReal      hinused,hlast,hcur,tcur,tolsfac;
3407cb94ee6SHong Zhang 
3417cb94ee6SHong Zhang   PetscFunctionBegin;
3427cb94ee6SHong Zhang   if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;}
3437cb94ee6SHong Zhang   else                                     {type = btype;}
3447cb94ee6SHong Zhang 
3457cb94ee6SHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
3467cb94ee6SHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
3477cb94ee6SHong Zhang   if (iascii) {
3487cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr);
3497cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr);
3507cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr);
3517cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr);
3527cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);CHKERRQ(ierr);
3537cb94ee6SHong Zhang     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
3547cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
3557cb94ee6SHong Zhang     } else {
3567cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
3577cb94ee6SHong Zhang     }
3582c823083SHong Zhang 
3595d47aee6SHong Zhang     /* Outputs from CVODE, CVSPILS */
3605d47aee6SHong Zhang     ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr);
3615d47aee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr);
3622c823083SHong Zhang     ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
3632c823083SHong Zhang                                    &nlinsetups,&nfails,&qlast,&qcur,
3642c823083SHong Zhang                                    &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr);
3652c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr);
3662c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr);
3672c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr);
3682c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr);
3692c823083SHong Zhang 
3702c823083SHong Zhang     ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr);
3712c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr);
3722c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr);
3732c823083SHong Zhang 
3742c823083SHong Zhang     ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */
3752c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr);
3762c823083SHong Zhang     ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr);
3772c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr);
3782c823083SHong Zhang     ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr);
3792c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr);
3802c823083SHong Zhang     ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr);
3812c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr);
3822c823083SHong Zhang     ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr);
3832c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr);
3842c823083SHong Zhang     ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr);
3852c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr);
3867cb94ee6SHong Zhang   } else if (isstring) {
3877cb94ee6SHong Zhang     ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr);
3887cb94ee6SHong Zhang   } else {
3897cb94ee6SHong Zhang     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by TS Sundials",((PetscObject)viewer)->type_name);
3907cb94ee6SHong Zhang   }
3917cb94ee6SHong Zhang   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
3927cb94ee6SHong Zhang   ierr = PCView(cvode->pc,viewer);CHKERRQ(ierr);
3937cb94ee6SHong Zhang   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
3947cb94ee6SHong Zhang   PetscFunctionReturn(0);
3957cb94ee6SHong Zhang }
3967cb94ee6SHong Zhang 
3977cb94ee6SHong Zhang 
3987cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/
3997cb94ee6SHong Zhang EXTERN_C_BEGIN
4007cb94ee6SHong Zhang #undef __FUNCT__
4017cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials"
4026fadb2cdSHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
4037cb94ee6SHong Zhang {
4047cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4057cb94ee6SHong Zhang 
4067cb94ee6SHong Zhang   PetscFunctionBegin;
4077cb94ee6SHong Zhang   cvode->cvode_type = type;
4087cb94ee6SHong Zhang   PetscFunctionReturn(0);
4097cb94ee6SHong Zhang }
4107cb94ee6SHong Zhang EXTERN_C_END
4117cb94ee6SHong Zhang 
4127cb94ee6SHong Zhang EXTERN_C_BEGIN
4137cb94ee6SHong Zhang #undef __FUNCT__
4147cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials"
4157cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart_Sundials(TS ts,int restart)
4167cb94ee6SHong Zhang {
4177cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4187cb94ee6SHong Zhang 
4197cb94ee6SHong Zhang   PetscFunctionBegin;
4207cb94ee6SHong Zhang   cvode->restart = restart;
4217cb94ee6SHong Zhang   PetscFunctionReturn(0);
4227cb94ee6SHong Zhang }
4237cb94ee6SHong Zhang EXTERN_C_END
4247cb94ee6SHong Zhang 
4257cb94ee6SHong Zhang EXTERN_C_BEGIN
4267cb94ee6SHong Zhang #undef __FUNCT__
4277cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials"
4287cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
4297cb94ee6SHong Zhang {
4307cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4317cb94ee6SHong Zhang 
4327cb94ee6SHong Zhang   PetscFunctionBegin;
4337cb94ee6SHong Zhang   cvode->linear_tol = tol;
4347cb94ee6SHong Zhang   PetscFunctionReturn(0);
4357cb94ee6SHong Zhang }
4367cb94ee6SHong Zhang EXTERN_C_END
4377cb94ee6SHong Zhang 
4387cb94ee6SHong Zhang EXTERN_C_BEGIN
4397cb94ee6SHong Zhang #undef __FUNCT__
4407cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials"
4417cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
4427cb94ee6SHong Zhang {
4437cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4447cb94ee6SHong Zhang 
4457cb94ee6SHong Zhang   PetscFunctionBegin;
4467cb94ee6SHong Zhang   cvode->gtype = type;
4477cb94ee6SHong Zhang   PetscFunctionReturn(0);
4487cb94ee6SHong Zhang }
4497cb94ee6SHong Zhang EXTERN_C_END
4507cb94ee6SHong Zhang 
4517cb94ee6SHong Zhang EXTERN_C_BEGIN
4527cb94ee6SHong Zhang #undef __FUNCT__
4537cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials"
4547cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
4557cb94ee6SHong Zhang {
4567cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4577cb94ee6SHong Zhang 
4587cb94ee6SHong Zhang   PetscFunctionBegin;
4597cb94ee6SHong Zhang   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
4607cb94ee6SHong Zhang   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
4617cb94ee6SHong Zhang   PetscFunctionReturn(0);
4627cb94ee6SHong Zhang }
4637cb94ee6SHong Zhang EXTERN_C_END
4647cb94ee6SHong Zhang 
4657cb94ee6SHong Zhang EXTERN_C_BEGIN
4667cb94ee6SHong Zhang #undef __FUNCT__
4677cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials"
4687cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC_Sundials(TS ts,PC *pc)
4697cb94ee6SHong Zhang {
4707cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4717cb94ee6SHong Zhang 
4727cb94ee6SHong Zhang   PetscFunctionBegin;
4737cb94ee6SHong Zhang   *pc = cvode->pc;
4747cb94ee6SHong Zhang   PetscFunctionReturn(0);
4757cb94ee6SHong Zhang }
4767cb94ee6SHong Zhang EXTERN_C_END
4777cb94ee6SHong Zhang 
4787cb94ee6SHong Zhang EXTERN_C_BEGIN
4797cb94ee6SHong Zhang #undef __FUNCT__
4807cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials"
4817cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
4827cb94ee6SHong Zhang {
4837cb94ee6SHong Zhang   PetscFunctionBegin;
4842c823083SHong Zhang   if (nonlin) *nonlin = ts->nonlinear_its;
4852c823083SHong Zhang   if (lin)    *lin    = ts->linear_its;
4867cb94ee6SHong Zhang   PetscFunctionReturn(0);
4877cb94ee6SHong Zhang }
4887cb94ee6SHong Zhang EXTERN_C_END
4897cb94ee6SHong Zhang 
4907cb94ee6SHong Zhang EXTERN_C_BEGIN
4917cb94ee6SHong Zhang #undef __FUNCT__
4927cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime_Sundials"
4937cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime_Sundials(TS ts,PetscTruth s)
4947cb94ee6SHong Zhang {
4957cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4967cb94ee6SHong Zhang 
4977cb94ee6SHong Zhang   PetscFunctionBegin;
4987cb94ee6SHong Zhang   cvode->exact_final_time = s;
4997cb94ee6SHong Zhang   PetscFunctionReturn(0);
5007cb94ee6SHong Zhang }
5017cb94ee6SHong Zhang EXTERN_C_END
5022bfc04deSHong Zhang 
5032bfc04deSHong Zhang EXTERN_C_BEGIN
5042bfc04deSHong Zhang #undef __FUNCT__
5052bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials"
5062bfc04deSHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscTruth s)
5072bfc04deSHong Zhang {
5082bfc04deSHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
5092bfc04deSHong Zhang 
5102bfc04deSHong Zhang   PetscFunctionBegin;
5112bfc04deSHong Zhang   cvode->monitorstep = s;
5122bfc04deSHong Zhang   PetscFunctionReturn(0);
5132bfc04deSHong Zhang }
5142bfc04deSHong Zhang EXTERN_C_END
5157cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
5167cb94ee6SHong Zhang 
5177cb94ee6SHong Zhang #undef __FUNCT__
5187cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations"
5197cb94ee6SHong Zhang /*@C
5207cb94ee6SHong Zhang    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
5217cb94ee6SHong Zhang 
5227cb94ee6SHong Zhang    Not Collective
5237cb94ee6SHong Zhang 
5247cb94ee6SHong Zhang    Input parameters:
5257cb94ee6SHong Zhang .    ts     - the time-step context
5267cb94ee6SHong Zhang 
5277cb94ee6SHong Zhang    Output Parameters:
5287cb94ee6SHong Zhang +   nonlin - number of nonlinear iterations
5297cb94ee6SHong Zhang -   lin    - number of linear iterations
5307cb94ee6SHong Zhang 
5317cb94ee6SHong Zhang    Level: advanced
5327cb94ee6SHong Zhang 
5337cb94ee6SHong Zhang    Notes:
5347cb94ee6SHong Zhang     These return the number since the creation of the TS object
5357cb94ee6SHong Zhang 
5367cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations
5377cb94ee6SHong Zhang 
5387cb94ee6SHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(),
5397cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
5407cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
5417cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
5427cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
5437cb94ee6SHong Zhang 
5447cb94ee6SHong Zhang @*/
5457cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
5467cb94ee6SHong Zhang {
5477cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,int*,int*);
5487cb94ee6SHong Zhang 
5497cb94ee6SHong Zhang   PetscFunctionBegin;
5507cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetIterations_C",(void (**)(void))&f);CHKERRQ(ierr);
5517cb94ee6SHong Zhang   if (f) {
5527cb94ee6SHong Zhang     ierr = (*f)(ts,nonlin,lin);CHKERRQ(ierr);
5537cb94ee6SHong Zhang   }
5547cb94ee6SHong Zhang   PetscFunctionReturn(0);
5557cb94ee6SHong Zhang }
5567cb94ee6SHong Zhang 
5577cb94ee6SHong Zhang #undef __FUNCT__
5587cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType"
5597cb94ee6SHong Zhang /*@
5607cb94ee6SHong Zhang    TSSundialsSetType - Sets the method that Sundials will use for integration.
5617cb94ee6SHong Zhang 
5627cb94ee6SHong Zhang    Collective on TS
5637cb94ee6SHong Zhang 
5647cb94ee6SHong Zhang    Input parameters:
5657cb94ee6SHong Zhang +    ts     - the time-step context
5667cb94ee6SHong Zhang -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
5677cb94ee6SHong Zhang 
5687cb94ee6SHong Zhang    Level: intermediate
5697cb94ee6SHong Zhang 
5707cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula
5717cb94ee6SHong Zhang 
5727cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(),  TSSundialsSetGMRESRestart(),
5737cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
5747cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
5757cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
5767cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
5777cb94ee6SHong Zhang @*/
5786fadb2cdSHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType(TS ts,TSSundialsLmmType type)
5797cb94ee6SHong Zhang {
5806fadb2cdSHong Zhang   PetscErrorCode ierr,(*f)(TS,TSSundialsLmmType);
5817cb94ee6SHong Zhang 
5827cb94ee6SHong Zhang   PetscFunctionBegin;
5837cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
5847cb94ee6SHong Zhang   if (f) {
5857cb94ee6SHong Zhang     ierr = (*f)(ts,type);CHKERRQ(ierr);
5867cb94ee6SHong Zhang   }
5877cb94ee6SHong Zhang   PetscFunctionReturn(0);
5887cb94ee6SHong Zhang }
5897cb94ee6SHong Zhang 
5907cb94ee6SHong Zhang #undef __FUNCT__
5917cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart"
5927cb94ee6SHong Zhang /*@
5937cb94ee6SHong Zhang    TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by
5947cb94ee6SHong Zhang        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
5957cb94ee6SHong Zhang        this is ALSO the maximum number of GMRES steps that will be used.
5967cb94ee6SHong Zhang 
5977cb94ee6SHong Zhang    Collective on TS
5987cb94ee6SHong Zhang 
5997cb94ee6SHong Zhang    Input parameters:
6007cb94ee6SHong Zhang +    ts      - the time-step context
6017cb94ee6SHong Zhang -    restart - number of direction vectors (the restart size).
6027cb94ee6SHong Zhang 
6037cb94ee6SHong Zhang    Level: advanced
6047cb94ee6SHong Zhang 
6057cb94ee6SHong Zhang .keywords: GMRES, restart
6067cb94ee6SHong Zhang 
6077cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
6087cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
6097cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6107cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
6117cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
6127cb94ee6SHong Zhang 
6137cb94ee6SHong Zhang @*/
6147cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart(TS ts,int restart)
6157cb94ee6SHong Zhang {
6167cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,int);
6177cb94ee6SHong Zhang 
6187cb94ee6SHong Zhang   PetscFunctionBegin;
6197cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGMRESRestart_C",(void (**)(void))&f);CHKERRQ(ierr);
6207cb94ee6SHong Zhang   if (f) {
6217cb94ee6SHong Zhang     ierr = (*f)(ts,restart);CHKERRQ(ierr);
6227cb94ee6SHong Zhang   }
6237cb94ee6SHong Zhang   PetscFunctionReturn(0);
6247cb94ee6SHong Zhang }
6257cb94ee6SHong Zhang 
6267cb94ee6SHong Zhang #undef __FUNCT__
6277cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance"
6287cb94ee6SHong Zhang /*@
6297cb94ee6SHong Zhang    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
6307cb94ee6SHong Zhang        system by SUNDIALS.
6317cb94ee6SHong Zhang 
6327cb94ee6SHong Zhang    Collective on TS
6337cb94ee6SHong Zhang 
6347cb94ee6SHong Zhang    Input parameters:
6357cb94ee6SHong Zhang +    ts     - the time-step context
6367cb94ee6SHong Zhang -    tol    - the factor by which the tolerance on the nonlinear solver is
6377cb94ee6SHong Zhang              multiplied to get the tolerance on the linear solver, .05 by default.
6387cb94ee6SHong Zhang 
6397cb94ee6SHong Zhang    Level: advanced
6407cb94ee6SHong Zhang 
6417cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS
6427cb94ee6SHong Zhang 
6437cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6447cb94ee6SHong Zhang           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
6457cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6467cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
6477cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
6487cb94ee6SHong Zhang 
6497cb94ee6SHong Zhang @*/
6507cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance(TS ts,double tol)
6517cb94ee6SHong Zhang {
6527cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,double);
6537cb94ee6SHong Zhang 
6547cb94ee6SHong Zhang   PetscFunctionBegin;
6557cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
6567cb94ee6SHong Zhang   if (f) {
6577cb94ee6SHong Zhang     ierr = (*f)(ts,tol);CHKERRQ(ierr);
6587cb94ee6SHong Zhang   }
6597cb94ee6SHong Zhang   PetscFunctionReturn(0);
6607cb94ee6SHong Zhang }
6617cb94ee6SHong Zhang 
6627cb94ee6SHong Zhang #undef __FUNCT__
6637cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType"
6647cb94ee6SHong Zhang /*@
6657cb94ee6SHong Zhang    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
6667cb94ee6SHong Zhang         in GMRES method by SUNDIALS linear solver.
6677cb94ee6SHong Zhang 
6687cb94ee6SHong Zhang    Collective on TS
6697cb94ee6SHong Zhang 
6707cb94ee6SHong Zhang    Input parameters:
6717cb94ee6SHong Zhang +    ts  - the time-step context
6727cb94ee6SHong Zhang -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
6737cb94ee6SHong Zhang 
6747cb94ee6SHong Zhang    Level: advanced
6757cb94ee6SHong Zhang 
6767cb94ee6SHong Zhang .keywords: Sundials, orthogonalization
6777cb94ee6SHong Zhang 
6787cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6797cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
6807cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6817cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
6827cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
6837cb94ee6SHong Zhang 
6847cb94ee6SHong Zhang @*/
6857cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
6867cb94ee6SHong Zhang {
6877cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,TSSundialsGramSchmidtType);
6887cb94ee6SHong Zhang 
6897cb94ee6SHong Zhang   PetscFunctionBegin;
6907cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",(void (**)(void))&f);CHKERRQ(ierr);
6917cb94ee6SHong Zhang   if (f) {
6927cb94ee6SHong Zhang     ierr = (*f)(ts,type);CHKERRQ(ierr);
6937cb94ee6SHong Zhang   }
6947cb94ee6SHong Zhang   PetscFunctionReturn(0);
6957cb94ee6SHong Zhang }
6967cb94ee6SHong Zhang 
6977cb94ee6SHong Zhang #undef __FUNCT__
6987cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance"
6997cb94ee6SHong Zhang /*@
7007cb94ee6SHong Zhang    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
7017cb94ee6SHong Zhang                          Sundials for error control.
7027cb94ee6SHong Zhang 
7037cb94ee6SHong Zhang    Collective on TS
7047cb94ee6SHong Zhang 
7057cb94ee6SHong Zhang    Input parameters:
7067cb94ee6SHong Zhang +    ts  - the time-step context
7077cb94ee6SHong Zhang .    aabs - the absolute tolerance
7087cb94ee6SHong Zhang -    rel - the relative tolerance
7097cb94ee6SHong Zhang 
7107cb94ee6SHong Zhang      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
7117cb94ee6SHong Zhang     these regulate the size of the error for a SINGLE timestep.
7127cb94ee6SHong Zhang 
7137cb94ee6SHong Zhang    Level: intermediate
7147cb94ee6SHong Zhang 
7157cb94ee6SHong Zhang .keywords: Sundials, tolerance
7167cb94ee6SHong Zhang 
7177cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7187cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
7197cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7207cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
7217cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
7227cb94ee6SHong Zhang 
7237cb94ee6SHong Zhang @*/
7247cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance(TS ts,double aabs,double rel)
7257cb94ee6SHong Zhang {
7267cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,double,double);
7277cb94ee6SHong Zhang 
7287cb94ee6SHong Zhang   PetscFunctionBegin;
7297cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
7307cb94ee6SHong Zhang   if (f) {
7317cb94ee6SHong Zhang     ierr = (*f)(ts,aabs,rel);CHKERRQ(ierr);
7327cb94ee6SHong Zhang   }
7337cb94ee6SHong Zhang   PetscFunctionReturn(0);
7347cb94ee6SHong Zhang }
7357cb94ee6SHong Zhang 
7367cb94ee6SHong Zhang #undef __FUNCT__
7377cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC"
7387cb94ee6SHong Zhang /*@
7397cb94ee6SHong Zhang    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
7407cb94ee6SHong Zhang 
7417cb94ee6SHong Zhang    Input Parameter:
7427cb94ee6SHong Zhang .    ts - the time-step context
7437cb94ee6SHong Zhang 
7447cb94ee6SHong Zhang    Output Parameter:
7457cb94ee6SHong Zhang .    pc - the preconditioner context
7467cb94ee6SHong Zhang 
7477cb94ee6SHong Zhang    Level: advanced
7487cb94ee6SHong Zhang 
7497cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7507cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
7517cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7527cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
7537cb94ee6SHong Zhang @*/
7547cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC(TS ts,PC *pc)
7557cb94ee6SHong Zhang {
7567cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,PC *);
7577cb94ee6SHong Zhang 
7587cb94ee6SHong Zhang   PetscFunctionBegin;
7597cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetPC_C",(void (**)(void))&f);CHKERRQ(ierr);
7607cb94ee6SHong Zhang   if (f) {
7617cb94ee6SHong Zhang     ierr = (*f)(ts,pc);CHKERRQ(ierr);
7627cb94ee6SHong Zhang   } else {
7637cb94ee6SHong Zhang     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"TS must be of Sundials type to extract the PC");
7647cb94ee6SHong Zhang   }
7657cb94ee6SHong Zhang 
7667cb94ee6SHong Zhang   PetscFunctionReturn(0);
7677cb94ee6SHong Zhang }
7687cb94ee6SHong Zhang 
7697cb94ee6SHong Zhang #undef __FUNCT__
7707cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime"
7717cb94ee6SHong Zhang /*@
7727cb94ee6SHong Zhang    TSSundialsSetExactFinalTime - Determines if Sundials interpolates solution to the
7737cb94ee6SHong Zhang       exact final time requested by the user or just returns it at the final time
7747cb94ee6SHong Zhang       it computed. (Defaults to true).
7757cb94ee6SHong Zhang 
7767cb94ee6SHong Zhang    Input Parameter:
7777cb94ee6SHong Zhang +   ts - the time-step context
7787cb94ee6SHong Zhang -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE
7797cb94ee6SHong Zhang 
7807cb94ee6SHong Zhang    Level: beginner
7817cb94ee6SHong Zhang 
7827cb94ee6SHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7837cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
7847cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7857cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
7867cb94ee6SHong Zhang @*/
7877cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime(TS ts,PetscTruth ft)
7887cb94ee6SHong Zhang {
7897cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,PetscTruth);
7907cb94ee6SHong Zhang 
7917cb94ee6SHong Zhang   PetscFunctionBegin;
7927cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetExactFinalTime_C",(void (**)(void))&f);CHKERRQ(ierr);
7937cb94ee6SHong Zhang   if (f) {
7947cb94ee6SHong Zhang     ierr = (*f)(ts,ft);CHKERRQ(ierr);
7957cb94ee6SHong Zhang   }
7967cb94ee6SHong Zhang   PetscFunctionReturn(0);
7977cb94ee6SHong Zhang }
7987cb94ee6SHong Zhang 
7992bfc04deSHong Zhang #undef __FUNCT__
8002bfc04deSHong Zhang #define __FUNCT__ "TSSundialsMonitorInternalSteps"
8012bfc04deSHong Zhang /*@
8022bfc04deSHong Zhang    TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false).
8032bfc04deSHong Zhang 
8042bfc04deSHong Zhang    Input Parameter:
8052bfc04deSHong Zhang +   ts - the time-step context
8062bfc04deSHong Zhang -   ft - PETSC_TRUE if monitor, else PETSC_FALSE
8072bfc04deSHong Zhang 
8082bfc04deSHong Zhang    Level: beginner
8092bfc04deSHong Zhang 
8102bfc04deSHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8112bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
8122bfc04deSHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
8132bfc04deSHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
8142bfc04deSHong Zhang @*/
8152bfc04deSHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsMonitorInternalSteps(TS ts,PetscTruth ft)
8162bfc04deSHong Zhang {
8172bfc04deSHong Zhang   PetscErrorCode ierr,(*f)(TS,PetscTruth);
8182bfc04deSHong Zhang 
8192bfc04deSHong Zhang   PetscFunctionBegin;
8202bfc04deSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",(void (**)(void))&f);CHKERRQ(ierr);
8212bfc04deSHong Zhang   if (f) {
8222bfc04deSHong Zhang     ierr = (*f)(ts,ft);CHKERRQ(ierr);
8232bfc04deSHong Zhang   }
8242bfc04deSHong Zhang   PetscFunctionReturn(0);
8252bfc04deSHong Zhang }
8267cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
8277cb94ee6SHong Zhang /*MC
828*96f5712cSJed Brown       TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
8297cb94ee6SHong Zhang 
8307cb94ee6SHong Zhang    Options Database:
8317cb94ee6SHong Zhang +    -ts_sundials_type <bdf,adams>
8327cb94ee6SHong Zhang .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
8337cb94ee6SHong Zhang .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
8347cb94ee6SHong Zhang .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
8357cb94ee6SHong Zhang .    -ts_sundials_linear_tolerance <tol>
8367cb94ee6SHong Zhang .    -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions
83716016685SHong Zhang .    -ts_sundials_exact_final_time - Interpolate output to stop exactly at the final time
83816016685SHong Zhang -    -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps
83916016685SHong Zhang 
8407cb94ee6SHong Zhang 
8417cb94ee6SHong Zhang     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
8427cb94ee6SHong Zhang            only PETSc PC options
8437cb94ee6SHong Zhang 
8447cb94ee6SHong Zhang     Level: beginner
8457cb94ee6SHong Zhang 
8467cb94ee6SHong Zhang .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(),
8477cb94ee6SHong Zhang            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSundialsSetExactFinalTime()
8487cb94ee6SHong Zhang 
8497cb94ee6SHong Zhang M*/
8507cb94ee6SHong Zhang EXTERN_C_BEGIN
8517cb94ee6SHong Zhang #undef __FUNCT__
8527cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials"
8537cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Sundials(TS ts)
8547cb94ee6SHong Zhang {
8557cb94ee6SHong Zhang   TS_Sundials *cvode;
8567cb94ee6SHong Zhang   PetscErrorCode ierr;
8577cb94ee6SHong Zhang 
8587cb94ee6SHong Zhang   PetscFunctionBegin;
8597cb94ee6SHong Zhang   if (ts->problem_type != TS_NONLINEAR) {
8607cb94ee6SHong Zhang     SETERRQ(PETSC_ERR_SUP,"Only support for nonlinear problems");
8617cb94ee6SHong Zhang   }
86228adb3f7SHong Zhang   ts->ops->destroy        = TSDestroy_Sundials;
86328adb3f7SHong Zhang   ts->ops->view           = TSView_Sundials;
8647cb94ee6SHong Zhang   ts->ops->setup          = TSSetUp_Sundials_Nonlinear;
8657cb94ee6SHong Zhang   ts->ops->step           = TSStep_Sundials_Nonlinear;
8667cb94ee6SHong Zhang   ts->ops->setfromoptions = TSSetFromOptions_Sundials_Nonlinear;
8677cb94ee6SHong Zhang 
86838f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr);
8697adad957SLisandro Dalcin   ierr = PCCreate(((PetscObject)ts)->comm,&cvode->pc);CHKERRQ(ierr);
8707cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->pc);CHKERRQ(ierr);
8717cb94ee6SHong Zhang   ts->data          = (void*)cvode;
8726fadb2cdSHong Zhang   cvode->cvode_type = SUNDIALS_BDF;
8736fadb2cdSHong Zhang   cvode->gtype      = SUNDIALS_CLASSICAL_GS;
8747cb94ee6SHong Zhang   cvode->restart    = 5;
8757cb94ee6SHong Zhang   cvode->linear_tol = .05;
8767cb94ee6SHong Zhang 
8772bfc04deSHong Zhang   cvode->exact_final_time = PETSC_TRUE;
8782bfc04deSHong Zhang   cvode->monitorstep      = PETSC_FALSE;
8797cb94ee6SHong Zhang 
880a07356b8SSatish Balay   ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr);
8817cb94ee6SHong Zhang   /* set tolerance for Sundials */
8827cb94ee6SHong Zhang   cvode->reltol = 1e-6;
8832c823083SHong Zhang   cvode->abstol = 1e-6;
8847cb94ee6SHong Zhang 
8857cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
8867cb94ee6SHong Zhang                     TSSundialsSetType_Sundials);CHKERRQ(ierr);
8877cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C",
8887cb94ee6SHong Zhang                     "TSSundialsSetGMRESRestart_Sundials",
8897cb94ee6SHong Zhang                     TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr);
8907cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
8917cb94ee6SHong Zhang                     "TSSundialsSetLinearTolerance_Sundials",
8927cb94ee6SHong Zhang                      TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
8937cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
8947cb94ee6SHong Zhang                     "TSSundialsSetGramSchmidtType_Sundials",
8957cb94ee6SHong Zhang                      TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
8967cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
8977cb94ee6SHong Zhang                     "TSSundialsSetTolerance_Sundials",
8987cb94ee6SHong Zhang                      TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
8997cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
9007cb94ee6SHong Zhang                     "TSSundialsGetPC_Sundials",
9017cb94ee6SHong Zhang                      TSSundialsGetPC_Sundials);CHKERRQ(ierr);
9027cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
9037cb94ee6SHong Zhang                     "TSSundialsGetIterations_Sundials",
9047cb94ee6SHong Zhang                      TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
9057cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C",
9067cb94ee6SHong Zhang                     "TSSundialsSetExactFinalTime_Sundials",
9077cb94ee6SHong Zhang                      TSSundialsSetExactFinalTime_Sundials);CHKERRQ(ierr);
9082bfc04deSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",
9092bfc04deSHong Zhang                     "TSSundialsMonitorInternalSteps_Sundials",
9102bfc04deSHong Zhang                      TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr);
9112bfc04deSHong Zhang 
9127cb94ee6SHong Zhang   PetscFunctionReturn(0);
9137cb94ee6SHong Zhang }
9147cb94ee6SHong Zhang EXTERN_C_END
9157cb94ee6SHong Zhang 
9167cb94ee6SHong Zhang 
9177cb94ee6SHong Zhang 
9187cb94ee6SHong Zhang 
9197cb94ee6SHong Zhang 
9207cb94ee6SHong Zhang 
9217cb94ee6SHong Zhang 
9227cb94ee6SHong Zhang 
9237cb94ee6SHong Zhang 
9247cb94ee6SHong Zhang 
925