xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision 90d69ab76acc68be38016f1844fbdbca1424683c)
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.
7db268bfaSHong Zhang     Reference: sundials-2.3.0/examples/cvode/parallel/cvkryx_p.c
87cb94ee6SHong Zhang */
956a740aaSHong Zhang #include "sundials.h"  /*I "petscts.h" I*/
107cb94ee6SHong Zhang 
117cb94ee6SHong Zhang /*
127cb94ee6SHong Zhang       TSPrecond_Sundials - function that we provide to SUNDIALS to
137cb94ee6SHong Zhang                         evaluate the preconditioner.
147cb94ee6SHong Zhang */
157cb94ee6SHong Zhang #undef __FUNCT__
167cb94ee6SHong Zhang #define __FUNCT__ "TSPrecond_Sundials"
177cb94ee6SHong Zhang PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy,
187cb94ee6SHong Zhang                     booleantype jok,booleantype *jcurPtr,
197cb94ee6SHong Zhang                     realtype _gamma,void *P_data,
207cb94ee6SHong Zhang                     N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3)
217cb94ee6SHong Zhang {
227cb94ee6SHong Zhang   TS             ts = (TS) P_data;
237cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
247cb94ee6SHong Zhang   PC             pc = cvode->pc;
257cb94ee6SHong Zhang   PetscErrorCode ierr;
267cb94ee6SHong Zhang   Mat            Jac = ts->B;
277cb94ee6SHong Zhang   Vec            yy = cvode->w1;
287cb94ee6SHong Zhang   PetscScalar    one = 1.0,gm;
297cb94ee6SHong Zhang   MatStructure   str = DIFFERENT_NONZERO_PATTERN;
307cb94ee6SHong Zhang   PetscScalar    *y_data;
317cb94ee6SHong Zhang 
327cb94ee6SHong Zhang   PetscFunctionBegin;
337cb94ee6SHong Zhang   /* This allows us to construct preconditioners in-place if we like */
347cb94ee6SHong Zhang   ierr = MatSetUnfactored(Jac);CHKERRQ(ierr);
357cb94ee6SHong Zhang 
367cb94ee6SHong Zhang   /* jok - TRUE means reuse current Jacobian else recompute Jacobian */
377cb94ee6SHong Zhang   if (jok) {
387cb94ee6SHong Zhang     ierr     = MatCopy(cvode->pmat,Jac,str);CHKERRQ(ierr);
397cb94ee6SHong Zhang     *jcurPtr = FALSE;
407cb94ee6SHong Zhang   } else {
417cb94ee6SHong Zhang     /* make PETSc vector yy point to SUNDIALS vector y */
427cb94ee6SHong Zhang     y_data = (PetscScalar *) N_VGetArrayPointer(y);
437cb94ee6SHong Zhang     ierr   = VecPlaceArray(yy,y_data); CHKERRQ(ierr);
444ac7836bSHong Zhang 
457cb94ee6SHong Zhang     /* compute the Jacobian */
467cb94ee6SHong Zhang     ierr = TSComputeRHSJacobian(ts,ts->ptime,yy,&Jac,&Jac,&str);CHKERRQ(ierr);
477cb94ee6SHong Zhang     ierr = VecResetArray(yy); CHKERRQ(ierr);
484ac7836bSHong Zhang 
497cb94ee6SHong Zhang     /* copy the Jacobian matrix */
507cb94ee6SHong Zhang     if (!cvode->pmat) {
517cb94ee6SHong Zhang       ierr = MatDuplicate(Jac,MAT_COPY_VALUES,&cvode->pmat);CHKERRQ(ierr);
527cb94ee6SHong Zhang       ierr = PetscLogObjectParent(ts,cvode->pmat);CHKERRQ(ierr);
532c823083SHong Zhang     } else {
547cb94ee6SHong Zhang       ierr = MatCopy(Jac,cvode->pmat,str);CHKERRQ(ierr);
557cb94ee6SHong Zhang     }
567cb94ee6SHong Zhang     *jcurPtr = TRUE;
577cb94ee6SHong Zhang   }
587cb94ee6SHong Zhang 
597cb94ee6SHong Zhang   /* construct I-gamma*Jac  */
607cb94ee6SHong Zhang   gm   = -_gamma;
617cb94ee6SHong Zhang   ierr = MatScale(Jac,gm);CHKERRQ(ierr);
627cb94ee6SHong Zhang   ierr = MatShift(Jac,one);CHKERRQ(ierr);
637cb94ee6SHong Zhang 
647cb94ee6SHong Zhang   ierr = PCSetOperators(pc,Jac,Jac,str);CHKERRQ(ierr);
657cb94ee6SHong Zhang   PetscFunctionReturn(0);
667cb94ee6SHong Zhang }
677cb94ee6SHong Zhang 
687cb94ee6SHong Zhang /*
697cb94ee6SHong Zhang      TSPSolve_Sundials -  routine that we provide to Sundials that applies the preconditioner.
707cb94ee6SHong Zhang */
717cb94ee6SHong Zhang #undef __FUNCT__
727cb94ee6SHong Zhang #define __FUNCT__ "TSPSolve_Sundials"
734ac7836bSHong Zhang PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z,
744ac7836bSHong Zhang                                  realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp)
757cb94ee6SHong Zhang {
767cb94ee6SHong Zhang   TS              ts = (TS) P_data;
777cb94ee6SHong Zhang   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
787cb94ee6SHong Zhang   PC              pc = cvode->pc;
797cb94ee6SHong Zhang   Vec             rr = cvode->w1,zz = cvode->w2;
807cb94ee6SHong Zhang   PetscErrorCode  ierr;
817cb94ee6SHong Zhang   PetscScalar     *r_data,*z_data;
827cb94ee6SHong Zhang 
837cb94ee6SHong Zhang   PetscFunctionBegin;
847cb94ee6SHong Zhang   /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
857cb94ee6SHong Zhang   r_data  = (PetscScalar *) N_VGetArrayPointer(r);
867cb94ee6SHong Zhang   z_data  = (PetscScalar *) N_VGetArrayPointer(z);
877cb94ee6SHong Zhang   ierr = VecPlaceArray(rr,r_data); CHKERRQ(ierr);
887cb94ee6SHong Zhang   ierr = VecPlaceArray(zz,z_data); CHKERRQ(ierr);
894ac7836bSHong Zhang 
907cb94ee6SHong Zhang   /* Solve the Px=r and put the result in zz */
917cb94ee6SHong Zhang   ierr = PCApply(pc,rr,zz); CHKERRQ(ierr);
927cb94ee6SHong Zhang   ierr = VecResetArray(rr); CHKERRQ(ierr);
937cb94ee6SHong Zhang   ierr = VecResetArray(zz); CHKERRQ(ierr);
947cb94ee6SHong Zhang   PetscFunctionReturn(0);
957cb94ee6SHong Zhang }
967cb94ee6SHong Zhang 
977cb94ee6SHong Zhang /*
987cb94ee6SHong Zhang         TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
997cb94ee6SHong Zhang */
1007cb94ee6SHong Zhang #undef __FUNCT__
1017cb94ee6SHong Zhang #define __FUNCT__ "TSFunction_Sundials"
1024ac7836bSHong Zhang int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx)
1037cb94ee6SHong Zhang {
1047cb94ee6SHong Zhang   TS              ts = (TS) ctx;
1057adad957SLisandro Dalcin   MPI_Comm        comm = ((PetscObject)ts)->comm;
1067cb94ee6SHong Zhang   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
1077cb94ee6SHong Zhang   Vec             yy = cvode->w1,yyd = cvode->w2;
1087cb94ee6SHong Zhang   PetscScalar     *y_data,*ydot_data;
1097cb94ee6SHong Zhang   PetscErrorCode  ierr;
1107cb94ee6SHong Zhang 
1117cb94ee6SHong Zhang   PetscFunctionBegin;
1125ff03cf7SDinesh Kaushik   /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
1137cb94ee6SHong Zhang   y_data     = (PetscScalar *) N_VGetArrayPointer(y);
1147cb94ee6SHong Zhang   ydot_data  = (PetscScalar *) N_VGetArrayPointer(ydot);
1154ac7836bSHong Zhang   ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr);
1164ac7836bSHong Zhang   ierr = VecPlaceArray(yyd,ydot_data); CHKERRABORT(comm,ierr);
1174ac7836bSHong Zhang 
1187cb94ee6SHong Zhang   /* now compute the right hand side function */
1197cb94ee6SHong Zhang   ierr = TSComputeRHSFunction(ts,t,yy,yyd); CHKERRABORT(comm,ierr);
1207cb94ee6SHong Zhang   ierr = VecResetArray(yy); CHKERRABORT(comm,ierr);
1217cb94ee6SHong Zhang   ierr = VecResetArray(yyd); CHKERRABORT(comm,ierr);
1224ac7836bSHong Zhang   PetscFunctionReturn(0);
1237cb94ee6SHong Zhang }
1247cb94ee6SHong Zhang 
1257cb94ee6SHong Zhang /*
1267cb94ee6SHong Zhang        TSStep_Sundials_Nonlinear - Calls Sundials to integrate the ODE.
1277cb94ee6SHong Zhang */
1287cb94ee6SHong Zhang #undef __FUNCT__
1297cb94ee6SHong Zhang #define __FUNCT__ "TSStep_Sundials_Nonlinear"
1307cb94ee6SHong Zhang /*
1317cb94ee6SHong Zhang     TSStep_Sundials_Nonlinear -
1327cb94ee6SHong Zhang 
1337cb94ee6SHong Zhang    steps - number of time steps
1347cb94ee6SHong Zhang    time - time that integrater is  terminated.
1357cb94ee6SHong Zhang */
1367cb94ee6SHong Zhang PetscErrorCode TSStep_Sundials_Nonlinear(TS ts,int *steps,double *time)
1377cb94ee6SHong Zhang {
1387cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
1397cb94ee6SHong Zhang   Vec            sol = ts->vec_sol;
1407cb94ee6SHong Zhang   PetscErrorCode ierr;
1415d47aee6SHong Zhang   PetscInt       i,max_steps = ts->max_steps,flag;
1427cb94ee6SHong Zhang   long int       its;
1437cb94ee6SHong Zhang   realtype       t,tout;
1447cb94ee6SHong Zhang   PetscScalar    *y_data;
1457cb94ee6SHong Zhang   void           *mem;
1465d47aee6SHong Zhang   const PCType   pctype;
1475d47aee6SHong Zhang   PetscTruth     pcnone;
1487cb94ee6SHong Zhang 
1497cb94ee6SHong Zhang   PetscFunctionBegin;
150db268bfaSHong Zhang   /* Call CVodeCreate to create the solver memory */
1517cb94ee6SHong Zhang   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
1527cb94ee6SHong Zhang   if (!mem) SETERRQ(1,"CVodeCreate() fails");
1537cb94ee6SHong Zhang   flag = CVodeSetFdata(mem,ts);
1547cb94ee6SHong Zhang   if (flag) SETERRQ(1,"CVodeSetFdata() fails");
1557cb94ee6SHong Zhang 
1567cb94ee6SHong Zhang   /*
1577cb94ee6SHong Zhang      Call CVodeMalloc to initialize the integrator memory:
1587cb94ee6SHong Zhang      mem is the pointer to the integrator memory returned by CVodeCreate
1597cb94ee6SHong Zhang      f       is the user's right hand side function in y'=f(t,y)
1607cb94ee6SHong Zhang      T0      is the initial time
1617cb94ee6SHong Zhang      u       is the initial dependent variable vector
1627cb94ee6SHong Zhang      CV_SS   specifies scalar relative and absolute tolerances
1637cb94ee6SHong Zhang      reltol  is the relative tolerance
1647cb94ee6SHong Zhang      &abstol is a pointer to the scalar absolute tolerance
1657cb94ee6SHong Zhang   */
1667cb94ee6SHong Zhang   flag = CVodeMalloc(mem,TSFunction_Sundials,ts->ptime,cvode->y,CV_SS,cvode->reltol,&cvode->abstol);
1677cb94ee6SHong Zhang   if (flag) SETERRQ(1,"CVodeMalloc() fails");
1687cb94ee6SHong Zhang 
1697cb94ee6SHong Zhang   /* initialize the number of steps */
1707cb94ee6SHong Zhang   *steps = -ts->steps;
1717cb94ee6SHong Zhang   ierr   = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
1727cb94ee6SHong Zhang 
1737cb94ee6SHong Zhang   /* call CVSpgmr to use GMRES as the linear solver.        */
1747cb94ee6SHong Zhang   /* setup the ode integrator with the given preconditioner */
1755d47aee6SHong Zhang   ierr = PCGetType(cvode->pc,&pctype);CHKERRQ(ierr);
1765d47aee6SHong Zhang   ierr = PetscTypeCompare((PetscObject)cvode->pc,PCNONE,&pcnone);CHKERRQ(ierr);
1775d47aee6SHong Zhang   if (pcnone){
1785d47aee6SHong Zhang     flag  = CVSpgmr(mem,PREC_NONE,0);
1795d47aee6SHong Zhang   } else {
1807cb94ee6SHong Zhang     flag  = CVSpgmr(mem,PREC_LEFT,0);
1815d47aee6SHong Zhang   }
1827cb94ee6SHong Zhang   if (flag) SETERRQ(1,"CVSpgmr() fails");
1834ac7836bSHong Zhang 
1847cb94ee6SHong Zhang   /* Set preconditioner setup and solve routines Precond and PSolve,
1857cb94ee6SHong Zhang      and the pointer to the user-defined block data */
1864ac7836bSHong Zhang   flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials,ts);
1877cb94ee6SHong Zhang   if (flag) SETERRQ(1,"CVSpgmrSetPreconditioner() fails");
1887cb94ee6SHong Zhang 
1895d47aee6SHong Zhang   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
1905d47aee6SHong Zhang   if (flag) SETERRQ(1,"CVSpgmrSetGSType() fails");
1915d47aee6SHong Zhang 
1927cb94ee6SHong Zhang   tout = ts->max_time;
1937cb94ee6SHong Zhang   ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr);
1947cb94ee6SHong Zhang   N_VSetArrayPointer((realtype *)y_data,cvode->y);
1957cb94ee6SHong Zhang   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
1967cb94ee6SHong Zhang   for (i = 0; i < max_steps; i++) {
1977cb94ee6SHong Zhang     if (ts->ptime >= ts->max_time) break;
1987cb94ee6SHong Zhang     ierr = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);CHKERRQ(ierr);
199e44a3fc5SHong Zhang     if (ts->ops->postupdate){
200e44a3fc5SHong Zhang       ierr = (*ts->ops->postupdate)(ts,ts->ptime,PETSC_NULL);CHKERRQ(ierr);
201e44a3fc5SHong Zhang     }
2027cb94ee6SHong Zhang     if (t > ts->max_time && cvode->exact_final_time) {
2037cb94ee6SHong Zhang       /* interpolate to final requested time */
2047cb94ee6SHong Zhang       ierr = CVodeGetDky(mem,tout,0,cvode->y);CHKERRQ(ierr);
2057cb94ee6SHong Zhang       t = tout;
2067cb94ee6SHong Zhang     }
2077cb94ee6SHong Zhang     ts->time_step = t - ts->ptime;
2087cb94ee6SHong Zhang     ts->ptime     = t;
2097cb94ee6SHong Zhang 
2107cb94ee6SHong Zhang     /* copy the solution from cvode->y to cvode->update and sol */
2117cb94ee6SHong Zhang     ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr);
2127cb94ee6SHong Zhang     ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr);
2137cb94ee6SHong Zhang     ierr = VecResetArray(cvode->w1); CHKERRQ(ierr);
2147cb94ee6SHong Zhang     ierr = VecCopy(cvode->update,sol);CHKERRQ(ierr);
2157cb94ee6SHong Zhang     ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr);
2167cb94ee6SHong Zhang     ts->nonlinear_its = its;
2174ac7836bSHong Zhang     ierr = CVSpilsGetNumLinIters(mem, &its);
2187cb94ee6SHong Zhang     ts->linear_its = its;
2197cb94ee6SHong Zhang     ts->steps++;
2207cb94ee6SHong Zhang     ierr = TSMonitor(ts,ts->steps,t,sol);CHKERRQ(ierr);
2217cb94ee6SHong Zhang   }
2222c823083SHong Zhang   cvode->mem = mem;
2237cb94ee6SHong Zhang   *steps    += ts->steps;
2247cb94ee6SHong Zhang   *time      = t;
2257cb94ee6SHong Zhang   PetscFunctionReturn(0);
2267cb94ee6SHong Zhang }
2277cb94ee6SHong Zhang 
2287cb94ee6SHong Zhang #undef __FUNCT__
2297cb94ee6SHong Zhang #define __FUNCT__ "TSDestroy_Sundials"
2307cb94ee6SHong Zhang PetscErrorCode TSDestroy_Sundials(TS ts)
2317cb94ee6SHong Zhang {
2327cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2337cb94ee6SHong Zhang   PetscErrorCode ierr;
2347cb94ee6SHong Zhang 
2357cb94ee6SHong Zhang   PetscFunctionBegin;
2367cb94ee6SHong Zhang   if (cvode->pmat)   {ierr = MatDestroy(cvode->pmat);CHKERRQ(ierr);}
2377cb94ee6SHong Zhang   if (cvode->pc)     {ierr = PCDestroy(cvode->pc);CHKERRQ(ierr);}
2387cb94ee6SHong Zhang   if (cvode->update) {ierr = VecDestroy(cvode->update);CHKERRQ(ierr);}
2397cb94ee6SHong Zhang   if (cvode->func)   {ierr = VecDestroy(cvode->func);CHKERRQ(ierr);}
2407cb94ee6SHong Zhang   if (cvode->rhs)    {ierr = VecDestroy(cvode->rhs);CHKERRQ(ierr);}
2417cb94ee6SHong Zhang   if (cvode->w1)     {ierr = VecDestroy(cvode->w1);CHKERRQ(ierr);}
2427cb94ee6SHong Zhang   if (cvode->w2)     {ierr = VecDestroy(cvode->w2);CHKERRQ(ierr);}
243a07356b8SSatish Balay   ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr);
2442c823083SHong Zhang   if (cvode->mem) {CVodeFree(&cvode->mem);}
2457cb94ee6SHong Zhang   ierr = PetscFree(cvode);CHKERRQ(ierr);
2467cb94ee6SHong Zhang   PetscFunctionReturn(0);
2477cb94ee6SHong Zhang }
2487cb94ee6SHong Zhang 
2497cb94ee6SHong Zhang #undef __FUNCT__
2507cb94ee6SHong Zhang #define __FUNCT__ "TSSetUp_Sundials_Nonlinear"
2517cb94ee6SHong Zhang PetscErrorCode TSSetUp_Sundials_Nonlinear(TS ts)
2527cb94ee6SHong Zhang {
2537cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
2547cb94ee6SHong Zhang   PetscErrorCode ierr;
2557cb94ee6SHong Zhang   int            glosize,locsize,i;
2567cb94ee6SHong Zhang   PetscScalar    *y_data,*parray;
2577cb94ee6SHong Zhang 
2587cb94ee6SHong Zhang   PetscFunctionBegin;
2597cb94ee6SHong Zhang   ierr = PCSetFromOptions(cvode->pc);CHKERRQ(ierr);
2607cb94ee6SHong Zhang   /* get the vector size */
2617cb94ee6SHong Zhang   ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr);
2627cb94ee6SHong Zhang   ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr);
2637cb94ee6SHong Zhang 
2647cb94ee6SHong Zhang   /* allocate the memory for N_Vec y */
265a07356b8SSatish Balay   cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
2667cb94ee6SHong Zhang   if (!cvode->y) SETERRQ(1,"cvode->y is not allocated");
2677cb94ee6SHong Zhang 
2687cb94ee6SHong Zhang   /* initialize N_Vec y */
2697cb94ee6SHong Zhang   ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr);
2707cb94ee6SHong Zhang   y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y);
2717cb94ee6SHong Zhang   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
2727cb94ee6SHong Zhang   /*ierr = PetscMemcpy(y_data,parray,locsize*sizeof(PETSC_SCALAR)); CHKERRQ(ierr);*/
2737cb94ee6SHong Zhang   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
2747cb94ee6SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr);
2757cb94ee6SHong Zhang   ierr = VecDuplicate(ts->vec_sol,&cvode->func);CHKERRQ(ierr);
2767cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr);
2777cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->func);CHKERRQ(ierr);
2787cb94ee6SHong Zhang 
2797cb94ee6SHong Zhang   /*
2807cb94ee6SHong Zhang       Create work vectors for the TSPSolve_Sundials() routine. Note these are
2817cb94ee6SHong Zhang     allocated with zero space arrays because the actual array space is provided
2827cb94ee6SHong Zhang     by Sundials and set using VecPlaceArray().
2837cb94ee6SHong Zhang   */
2847adad957SLisandro Dalcin   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr);
2857adad957SLisandro Dalcin   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr);
2867cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr);
2877cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr);
2887cb94ee6SHong Zhang   PetscFunctionReturn(0);
2897cb94ee6SHong Zhang }
2907cb94ee6SHong Zhang 
2916fadb2cdSHong Zhang /* type of CVODE linear multistep method */
292dfcd32c8SHong Zhang const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0};
2936fadb2cdSHong Zhang /* type of G-S orthogonalization used by CVODE linear solver */
294dfcd32c8SHong Zhang const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0};
295a04cf4d8SBarry Smith 
2967cb94ee6SHong Zhang #undef __FUNCT__
2977cb94ee6SHong Zhang #define __FUNCT__ "TSSetFromOptions_Sundials_Nonlinear"
2987cb94ee6SHong Zhang PetscErrorCode TSSetFromOptions_Sundials_Nonlinear(TS ts)
2997cb94ee6SHong Zhang {
3007cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
3017cb94ee6SHong Zhang   PetscErrorCode ierr;
3027cb94ee6SHong Zhang   int            indx;
3037cb94ee6SHong Zhang   PetscTruth     flag;
3047cb94ee6SHong Zhang 
3057cb94ee6SHong Zhang   PetscFunctionBegin;
3067cb94ee6SHong Zhang   ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr);
3076fadb2cdSHong Zhang     ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr);
3087cb94ee6SHong Zhang     if (flag) {
3096fadb2cdSHong Zhang       ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr);
3107cb94ee6SHong Zhang     }
3116fadb2cdSHong Zhang     ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr);
3127cb94ee6SHong Zhang     if (flag) {
3137cb94ee6SHong Zhang       ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr);
3147cb94ee6SHong Zhang     }
3157cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr);
3167cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr);
3177cb94ee6SHong Zhang     ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr);
3187cb94ee6SHong Zhang     ierr = PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);CHKERRQ(ierr);
319*90d69ab7SBarry Smith     ierr = PetscOptionsTruth("-ts_sundials_exact_final_time","Allow SUNDIALS to stop near the final time, not exactly on it","TSSundialsSetExactFinalTime",cvode->exact_final_time,&cvode->exact_final_time,PETSC_NULL);CHKERRQ(ierr);
3207cb94ee6SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
3217cb94ee6SHong Zhang   PetscFunctionReturn(0);
3227cb94ee6SHong Zhang }
3237cb94ee6SHong Zhang 
3247cb94ee6SHong Zhang #undef __FUNCT__
3257cb94ee6SHong Zhang #define __FUNCT__ "TSView_Sundials"
3267cb94ee6SHong Zhang PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
3277cb94ee6SHong Zhang {
3287cb94ee6SHong Zhang   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
3297cb94ee6SHong Zhang   PetscErrorCode ierr;
3307cb94ee6SHong Zhang   char           *type;
3317cb94ee6SHong Zhang   char           atype[] = "Adams";
3327cb94ee6SHong Zhang   char           btype[] = "BDF: backward differentiation formula";
3337cb94ee6SHong Zhang   PetscTruth     iascii,isstring;
3342c823083SHong Zhang   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
3352c823083SHong Zhang   PetscInt       qlast,qcur;
3365d47aee6SHong Zhang   PetscReal      hinused,hlast,hcur,tcur,tolsfac;
3377cb94ee6SHong Zhang 
3387cb94ee6SHong Zhang   PetscFunctionBegin;
3397cb94ee6SHong Zhang   if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;}
3407cb94ee6SHong Zhang   else                                     {type = btype;}
3417cb94ee6SHong Zhang 
3427cb94ee6SHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
3437cb94ee6SHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
3447cb94ee6SHong Zhang   if (iascii) {
3457cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr);
3467cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr);
3477cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr);
3487cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr);
3497cb94ee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);CHKERRQ(ierr);
3507cb94ee6SHong Zhang     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
3517cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
3527cb94ee6SHong Zhang     } else {
3537cb94ee6SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
3547cb94ee6SHong Zhang     }
3552c823083SHong Zhang 
3565d47aee6SHong Zhang     /* Outputs from CVODE, CVSPILS */
3575d47aee6SHong Zhang     ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr);
3585d47aee6SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr);
3592c823083SHong Zhang     ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
3602c823083SHong Zhang                                    &nlinsetups,&nfails,&qlast,&qcur,
3612c823083SHong Zhang                                    &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr);
3622c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr);
3632c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr);
3642c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr);
3652c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr);
3662c823083SHong Zhang 
3672c823083SHong Zhang     ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr);
3682c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr);
3692c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr);
3702c823083SHong Zhang 
3712c823083SHong Zhang     ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */
3722c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr);
3732c823083SHong Zhang     ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr);
3742c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr);
3752c823083SHong Zhang     ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr);
3762c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr);
3772c823083SHong Zhang     ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr);
3782c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr);
3792c823083SHong Zhang     ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr);
3802c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr);
3812c823083SHong Zhang     ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr);
3822c823083SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr);
3837cb94ee6SHong Zhang   } else if (isstring) {
3847cb94ee6SHong Zhang     ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr);
3857cb94ee6SHong Zhang   } else {
3867cb94ee6SHong Zhang     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by TS Sundials",((PetscObject)viewer)->type_name);
3877cb94ee6SHong Zhang   }
3887cb94ee6SHong Zhang   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
3897cb94ee6SHong Zhang   ierr = PCView(cvode->pc,viewer);CHKERRQ(ierr);
3907cb94ee6SHong Zhang   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
3917cb94ee6SHong Zhang   PetscFunctionReturn(0);
3927cb94ee6SHong Zhang }
3937cb94ee6SHong Zhang 
3947cb94ee6SHong Zhang 
3957cb94ee6SHong Zhang /* --------------------------------------------------------------------------*/
3967cb94ee6SHong Zhang EXTERN_C_BEGIN
3977cb94ee6SHong Zhang #undef __FUNCT__
3987cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType_Sundials"
3996fadb2cdSHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
4007cb94ee6SHong Zhang {
4017cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4027cb94ee6SHong Zhang 
4037cb94ee6SHong Zhang   PetscFunctionBegin;
4047cb94ee6SHong Zhang   cvode->cvode_type = type;
4057cb94ee6SHong Zhang   PetscFunctionReturn(0);
4067cb94ee6SHong Zhang }
4077cb94ee6SHong Zhang EXTERN_C_END
4087cb94ee6SHong Zhang 
4097cb94ee6SHong Zhang EXTERN_C_BEGIN
4107cb94ee6SHong Zhang #undef __FUNCT__
4117cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials"
4127cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart_Sundials(TS ts,int restart)
4137cb94ee6SHong Zhang {
4147cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4157cb94ee6SHong Zhang 
4167cb94ee6SHong Zhang   PetscFunctionBegin;
4177cb94ee6SHong Zhang   cvode->restart = restart;
4187cb94ee6SHong Zhang   PetscFunctionReturn(0);
4197cb94ee6SHong Zhang }
4207cb94ee6SHong Zhang EXTERN_C_END
4217cb94ee6SHong Zhang 
4227cb94ee6SHong Zhang EXTERN_C_BEGIN
4237cb94ee6SHong Zhang #undef __FUNCT__
4247cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials"
4257cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
4267cb94ee6SHong Zhang {
4277cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4287cb94ee6SHong Zhang 
4297cb94ee6SHong Zhang   PetscFunctionBegin;
4307cb94ee6SHong Zhang   cvode->linear_tol = tol;
4317cb94ee6SHong Zhang   PetscFunctionReturn(0);
4327cb94ee6SHong Zhang }
4337cb94ee6SHong Zhang EXTERN_C_END
4347cb94ee6SHong Zhang 
4357cb94ee6SHong Zhang EXTERN_C_BEGIN
4367cb94ee6SHong Zhang #undef __FUNCT__
4377cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials"
4387cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
4397cb94ee6SHong Zhang {
4407cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4417cb94ee6SHong Zhang 
4427cb94ee6SHong Zhang   PetscFunctionBegin;
4437cb94ee6SHong Zhang   cvode->gtype = type;
4447cb94ee6SHong Zhang   PetscFunctionReturn(0);
4457cb94ee6SHong Zhang }
4467cb94ee6SHong Zhang EXTERN_C_END
4477cb94ee6SHong Zhang 
4487cb94ee6SHong Zhang EXTERN_C_BEGIN
4497cb94ee6SHong Zhang #undef __FUNCT__
4507cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance_Sundials"
4517cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
4527cb94ee6SHong Zhang {
4537cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4547cb94ee6SHong Zhang 
4557cb94ee6SHong Zhang   PetscFunctionBegin;
4567cb94ee6SHong Zhang   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
4577cb94ee6SHong Zhang   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
4587cb94ee6SHong Zhang   PetscFunctionReturn(0);
4597cb94ee6SHong Zhang }
4607cb94ee6SHong Zhang EXTERN_C_END
4617cb94ee6SHong Zhang 
4627cb94ee6SHong Zhang EXTERN_C_BEGIN
4637cb94ee6SHong Zhang #undef __FUNCT__
4647cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC_Sundials"
4657cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC_Sundials(TS ts,PC *pc)
4667cb94ee6SHong Zhang {
4677cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4687cb94ee6SHong Zhang 
4697cb94ee6SHong Zhang   PetscFunctionBegin;
4707cb94ee6SHong Zhang   *pc = cvode->pc;
4717cb94ee6SHong Zhang   PetscFunctionReturn(0);
4727cb94ee6SHong Zhang }
4737cb94ee6SHong Zhang EXTERN_C_END
4747cb94ee6SHong Zhang 
4757cb94ee6SHong Zhang EXTERN_C_BEGIN
4767cb94ee6SHong Zhang #undef __FUNCT__
4777cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations_Sundials"
4787cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
4797cb94ee6SHong Zhang {
4807cb94ee6SHong Zhang   PetscFunctionBegin;
4812c823083SHong Zhang   if (nonlin) *nonlin = ts->nonlinear_its;
4822c823083SHong Zhang   if (lin)    *lin    = ts->linear_its;
4837cb94ee6SHong Zhang   PetscFunctionReturn(0);
4847cb94ee6SHong Zhang }
4857cb94ee6SHong Zhang EXTERN_C_END
4867cb94ee6SHong Zhang 
4877cb94ee6SHong Zhang EXTERN_C_BEGIN
4887cb94ee6SHong Zhang #undef __FUNCT__
4897cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime_Sundials"
4907cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime_Sundials(TS ts,PetscTruth s)
4917cb94ee6SHong Zhang {
4927cb94ee6SHong Zhang   TS_Sundials *cvode = (TS_Sundials*)ts->data;
4937cb94ee6SHong Zhang 
4947cb94ee6SHong Zhang   PetscFunctionBegin;
4957cb94ee6SHong Zhang   cvode->exact_final_time = s;
4967cb94ee6SHong Zhang   PetscFunctionReturn(0);
4977cb94ee6SHong Zhang }
4987cb94ee6SHong Zhang EXTERN_C_END
4997cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
5007cb94ee6SHong Zhang 
5017cb94ee6SHong Zhang #undef __FUNCT__
5027cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetIterations"
5037cb94ee6SHong Zhang /*@C
5047cb94ee6SHong Zhang    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
5057cb94ee6SHong Zhang 
5067cb94ee6SHong Zhang    Not Collective
5077cb94ee6SHong Zhang 
5087cb94ee6SHong Zhang    Input parameters:
5097cb94ee6SHong Zhang .    ts     - the time-step context
5107cb94ee6SHong Zhang 
5117cb94ee6SHong Zhang    Output Parameters:
5127cb94ee6SHong Zhang +   nonlin - number of nonlinear iterations
5137cb94ee6SHong Zhang -   lin    - number of linear iterations
5147cb94ee6SHong Zhang 
5157cb94ee6SHong Zhang    Level: advanced
5167cb94ee6SHong Zhang 
5177cb94ee6SHong Zhang    Notes:
5187cb94ee6SHong Zhang     These return the number since the creation of the TS object
5197cb94ee6SHong Zhang 
5207cb94ee6SHong Zhang .keywords: non-linear iterations, linear iterations
5217cb94ee6SHong Zhang 
5227cb94ee6SHong Zhang .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(),
5237cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
5247cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
5257cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
5267cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
5277cb94ee6SHong Zhang 
5287cb94ee6SHong Zhang @*/
5297cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
5307cb94ee6SHong Zhang {
5317cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,int*,int*);
5327cb94ee6SHong Zhang 
5337cb94ee6SHong Zhang   PetscFunctionBegin;
5347cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetIterations_C",(void (**)(void))&f);CHKERRQ(ierr);
5357cb94ee6SHong Zhang   if (f) {
5367cb94ee6SHong Zhang     ierr = (*f)(ts,nonlin,lin);CHKERRQ(ierr);
5377cb94ee6SHong Zhang   }
5387cb94ee6SHong Zhang   PetscFunctionReturn(0);
5397cb94ee6SHong Zhang }
5407cb94ee6SHong Zhang 
5417cb94ee6SHong Zhang #undef __FUNCT__
5427cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetType"
5437cb94ee6SHong Zhang /*@
5447cb94ee6SHong Zhang    TSSundialsSetType - Sets the method that Sundials will use for integration.
5457cb94ee6SHong Zhang 
5467cb94ee6SHong Zhang    Collective on TS
5477cb94ee6SHong Zhang 
5487cb94ee6SHong Zhang    Input parameters:
5497cb94ee6SHong Zhang +    ts     - the time-step context
5507cb94ee6SHong Zhang -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
5517cb94ee6SHong Zhang 
5527cb94ee6SHong Zhang    Level: intermediate
5537cb94ee6SHong Zhang 
5547cb94ee6SHong Zhang .keywords: Adams, backward differentiation formula
5557cb94ee6SHong Zhang 
5567cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(),  TSSundialsSetGMRESRestart(),
5577cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
5587cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
5597cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
5607cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
5617cb94ee6SHong Zhang @*/
5626fadb2cdSHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType(TS ts,TSSundialsLmmType type)
5637cb94ee6SHong Zhang {
5646fadb2cdSHong Zhang   PetscErrorCode ierr,(*f)(TS,TSSundialsLmmType);
5657cb94ee6SHong Zhang 
5667cb94ee6SHong Zhang   PetscFunctionBegin;
5677cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
5687cb94ee6SHong Zhang   if (f) {
5697cb94ee6SHong Zhang     ierr = (*f)(ts,type);CHKERRQ(ierr);
5707cb94ee6SHong Zhang   }
5717cb94ee6SHong Zhang   PetscFunctionReturn(0);
5727cb94ee6SHong Zhang }
5737cb94ee6SHong Zhang 
5747cb94ee6SHong Zhang #undef __FUNCT__
5757cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGMRESRestart"
5767cb94ee6SHong Zhang /*@
5777cb94ee6SHong Zhang    TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by
5787cb94ee6SHong Zhang        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
5797cb94ee6SHong Zhang        this is ALSO the maximum number of GMRES steps that will be used.
5807cb94ee6SHong Zhang 
5817cb94ee6SHong Zhang    Collective on TS
5827cb94ee6SHong Zhang 
5837cb94ee6SHong Zhang    Input parameters:
5847cb94ee6SHong Zhang +    ts      - the time-step context
5857cb94ee6SHong Zhang -    restart - number of direction vectors (the restart size).
5867cb94ee6SHong Zhang 
5877cb94ee6SHong Zhang    Level: advanced
5887cb94ee6SHong Zhang 
5897cb94ee6SHong Zhang .keywords: GMRES, restart
5907cb94ee6SHong Zhang 
5917cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
5927cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
5937cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
5947cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
5957cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
5967cb94ee6SHong Zhang 
5977cb94ee6SHong Zhang @*/
5987cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart(TS ts,int restart)
5997cb94ee6SHong Zhang {
6007cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,int);
6017cb94ee6SHong Zhang 
6027cb94ee6SHong Zhang   PetscFunctionBegin;
6037cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGMRESRestart_C",(void (**)(void))&f);CHKERRQ(ierr);
6047cb94ee6SHong Zhang   if (f) {
6057cb94ee6SHong Zhang     ierr = (*f)(ts,restart);CHKERRQ(ierr);
6067cb94ee6SHong Zhang   }
6077cb94ee6SHong Zhang   PetscFunctionReturn(0);
6087cb94ee6SHong Zhang }
6097cb94ee6SHong Zhang 
6107cb94ee6SHong Zhang #undef __FUNCT__
6117cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetLinearTolerance"
6127cb94ee6SHong Zhang /*@
6137cb94ee6SHong Zhang    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
6147cb94ee6SHong Zhang        system by SUNDIALS.
6157cb94ee6SHong Zhang 
6167cb94ee6SHong Zhang    Collective on TS
6177cb94ee6SHong Zhang 
6187cb94ee6SHong Zhang    Input parameters:
6197cb94ee6SHong Zhang +    ts     - the time-step context
6207cb94ee6SHong Zhang -    tol    - the factor by which the tolerance on the nonlinear solver is
6217cb94ee6SHong Zhang              multiplied to get the tolerance on the linear solver, .05 by default.
6227cb94ee6SHong Zhang 
6237cb94ee6SHong Zhang    Level: advanced
6247cb94ee6SHong Zhang 
6257cb94ee6SHong Zhang .keywords: GMRES, linear convergence tolerance, SUNDIALS
6267cb94ee6SHong Zhang 
6277cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6287cb94ee6SHong Zhang           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
6297cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6307cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
6317cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
6327cb94ee6SHong Zhang 
6337cb94ee6SHong Zhang @*/
6347cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance(TS ts,double tol)
6357cb94ee6SHong Zhang {
6367cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,double);
6377cb94ee6SHong Zhang 
6387cb94ee6SHong Zhang   PetscFunctionBegin;
6397cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
6407cb94ee6SHong Zhang   if (f) {
6417cb94ee6SHong Zhang     ierr = (*f)(ts,tol);CHKERRQ(ierr);
6427cb94ee6SHong Zhang   }
6437cb94ee6SHong Zhang   PetscFunctionReturn(0);
6447cb94ee6SHong Zhang }
6457cb94ee6SHong Zhang 
6467cb94ee6SHong Zhang #undef __FUNCT__
6477cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetGramSchmidtType"
6487cb94ee6SHong Zhang /*@
6497cb94ee6SHong Zhang    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
6507cb94ee6SHong Zhang         in GMRES method by SUNDIALS linear solver.
6517cb94ee6SHong Zhang 
6527cb94ee6SHong Zhang    Collective on TS
6537cb94ee6SHong Zhang 
6547cb94ee6SHong Zhang    Input parameters:
6557cb94ee6SHong Zhang +    ts  - the time-step context
6567cb94ee6SHong Zhang -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
6577cb94ee6SHong Zhang 
6587cb94ee6SHong Zhang    Level: advanced
6597cb94ee6SHong Zhang 
6607cb94ee6SHong Zhang .keywords: Sundials, orthogonalization
6617cb94ee6SHong Zhang 
6627cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6637cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
6647cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
6657cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
6667cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
6677cb94ee6SHong Zhang 
6687cb94ee6SHong Zhang @*/
6697cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
6707cb94ee6SHong Zhang {
6717cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,TSSundialsGramSchmidtType);
6727cb94ee6SHong Zhang 
6737cb94ee6SHong Zhang   PetscFunctionBegin;
6747cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",(void (**)(void))&f);CHKERRQ(ierr);
6757cb94ee6SHong Zhang   if (f) {
6767cb94ee6SHong Zhang     ierr = (*f)(ts,type);CHKERRQ(ierr);
6777cb94ee6SHong Zhang   }
6787cb94ee6SHong Zhang   PetscFunctionReturn(0);
6797cb94ee6SHong Zhang }
6807cb94ee6SHong Zhang 
6817cb94ee6SHong Zhang #undef __FUNCT__
6827cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetTolerance"
6837cb94ee6SHong Zhang /*@
6847cb94ee6SHong Zhang    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
6857cb94ee6SHong Zhang                          Sundials for error control.
6867cb94ee6SHong Zhang 
6877cb94ee6SHong Zhang    Collective on TS
6887cb94ee6SHong Zhang 
6897cb94ee6SHong Zhang    Input parameters:
6907cb94ee6SHong Zhang +    ts  - the time-step context
6917cb94ee6SHong Zhang .    aabs - the absolute tolerance
6927cb94ee6SHong Zhang -    rel - the relative tolerance
6937cb94ee6SHong Zhang 
6947cb94ee6SHong Zhang      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
6957cb94ee6SHong Zhang     these regulate the size of the error for a SINGLE timestep.
6967cb94ee6SHong Zhang 
6977cb94ee6SHong Zhang    Level: intermediate
6987cb94ee6SHong Zhang 
6997cb94ee6SHong Zhang .keywords: Sundials, tolerance
7007cb94ee6SHong Zhang 
7017cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7027cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
7037cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7047cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
7057cb94ee6SHong Zhang           TSSundialsSetExactFinalTime()
7067cb94ee6SHong Zhang 
7077cb94ee6SHong Zhang @*/
7087cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance(TS ts,double aabs,double rel)
7097cb94ee6SHong Zhang {
7107cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,double,double);
7117cb94ee6SHong Zhang 
7127cb94ee6SHong Zhang   PetscFunctionBegin;
7137cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
7147cb94ee6SHong Zhang   if (f) {
7157cb94ee6SHong Zhang     ierr = (*f)(ts,aabs,rel);CHKERRQ(ierr);
7167cb94ee6SHong Zhang   }
7177cb94ee6SHong Zhang   PetscFunctionReturn(0);
7187cb94ee6SHong Zhang }
7197cb94ee6SHong Zhang 
7207cb94ee6SHong Zhang #undef __FUNCT__
7217cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsGetPC"
7227cb94ee6SHong Zhang /*@
7237cb94ee6SHong Zhang    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
7247cb94ee6SHong Zhang 
7257cb94ee6SHong Zhang    Input Parameter:
7267cb94ee6SHong Zhang .    ts - the time-step context
7277cb94ee6SHong Zhang 
7287cb94ee6SHong Zhang    Output Parameter:
7297cb94ee6SHong Zhang .    pc - the preconditioner context
7307cb94ee6SHong Zhang 
7317cb94ee6SHong Zhang    Level: advanced
7327cb94ee6SHong Zhang 
7337cb94ee6SHong Zhang .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7347cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
7357cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7367cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
7377cb94ee6SHong Zhang @*/
7387cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC(TS ts,PC *pc)
7397cb94ee6SHong Zhang {
7407cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,PC *);
7417cb94ee6SHong Zhang 
7427cb94ee6SHong Zhang   PetscFunctionBegin;
7437cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetPC_C",(void (**)(void))&f);CHKERRQ(ierr);
7447cb94ee6SHong Zhang   if (f) {
7457cb94ee6SHong Zhang     ierr = (*f)(ts,pc);CHKERRQ(ierr);
7467cb94ee6SHong Zhang   } else {
7477cb94ee6SHong Zhang     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"TS must be of Sundials type to extract the PC");
7487cb94ee6SHong Zhang   }
7497cb94ee6SHong Zhang 
7507cb94ee6SHong Zhang   PetscFunctionReturn(0);
7517cb94ee6SHong Zhang }
7527cb94ee6SHong Zhang 
7537cb94ee6SHong Zhang #undef __FUNCT__
7547cb94ee6SHong Zhang #define __FUNCT__ "TSSundialsSetExactFinalTime"
7557cb94ee6SHong Zhang /*@
7567cb94ee6SHong Zhang    TSSundialsSetExactFinalTime - Determines if Sundials interpolates solution to the
7577cb94ee6SHong Zhang       exact final time requested by the user or just returns it at the final time
7587cb94ee6SHong Zhang       it computed. (Defaults to true).
7597cb94ee6SHong Zhang 
7607cb94ee6SHong Zhang    Input Parameter:
7617cb94ee6SHong Zhang +   ts - the time-step context
7627cb94ee6SHong Zhang -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE
7637cb94ee6SHong Zhang 
7647cb94ee6SHong Zhang    Level: beginner
7657cb94ee6SHong Zhang 
7667cb94ee6SHong Zhang .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7677cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
7687cb94ee6SHong Zhang           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
7697cb94ee6SHong Zhang           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
7707cb94ee6SHong Zhang @*/
7717cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime(TS ts,PetscTruth ft)
7727cb94ee6SHong Zhang {
7737cb94ee6SHong Zhang   PetscErrorCode ierr,(*f)(TS,PetscTruth);
7747cb94ee6SHong Zhang 
7757cb94ee6SHong Zhang   PetscFunctionBegin;
7767cb94ee6SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetExactFinalTime_C",(void (**)(void))&f);CHKERRQ(ierr);
7777cb94ee6SHong Zhang   if (f) {
7787cb94ee6SHong Zhang     ierr = (*f)(ts,ft);CHKERRQ(ierr);
7797cb94ee6SHong Zhang   }
7807cb94ee6SHong Zhang   PetscFunctionReturn(0);
7817cb94ee6SHong Zhang }
7827cb94ee6SHong Zhang 
7837cb94ee6SHong Zhang /* -------------------------------------------------------------------------------------------*/
7847cb94ee6SHong Zhang /*MC
7857cb94ee6SHong Zhang       TS_Sundials - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
7867cb94ee6SHong Zhang 
7877cb94ee6SHong Zhang    Options Database:
7887cb94ee6SHong Zhang +    -ts_sundials_type <bdf,adams>
7897cb94ee6SHong Zhang .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
7907cb94ee6SHong Zhang .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
7917cb94ee6SHong Zhang .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
7927cb94ee6SHong Zhang .    -ts_sundials_linear_tolerance <tol>
7937cb94ee6SHong Zhang .    -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions
7947cb94ee6SHong Zhang -    -ts_sundials_not_exact_final_time -Allow SUNDIALS to stop near the final time, not exactly on it
7957cb94ee6SHong Zhang 
7967cb94ee6SHong Zhang     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
7977cb94ee6SHong Zhang            only PETSc PC options
7987cb94ee6SHong Zhang 
7997cb94ee6SHong Zhang     Level: beginner
8007cb94ee6SHong Zhang 
8017cb94ee6SHong Zhang .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(),
8027cb94ee6SHong Zhang            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSundialsSetExactFinalTime()
8037cb94ee6SHong Zhang 
8047cb94ee6SHong Zhang M*/
8057cb94ee6SHong Zhang EXTERN_C_BEGIN
8067cb94ee6SHong Zhang #undef __FUNCT__
8077cb94ee6SHong Zhang #define __FUNCT__ "TSCreate_Sundials"
8087cb94ee6SHong Zhang PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Sundials(TS ts)
8097cb94ee6SHong Zhang {
8107cb94ee6SHong Zhang   TS_Sundials *cvode;
8117cb94ee6SHong Zhang   PetscErrorCode ierr;
8127cb94ee6SHong Zhang 
8137cb94ee6SHong Zhang   PetscFunctionBegin;
8147cb94ee6SHong Zhang   ts->ops->destroy         = TSDestroy_Sundials;
8157cb94ee6SHong Zhang   ts->ops->view            = TSView_Sundials;
8167cb94ee6SHong Zhang 
8177cb94ee6SHong Zhang   if (ts->problem_type != TS_NONLINEAR) {
8187cb94ee6SHong Zhang     SETERRQ(PETSC_ERR_SUP,"Only support for nonlinear problems");
8197cb94ee6SHong Zhang   }
8207cb94ee6SHong Zhang   ts->ops->setup           = TSSetUp_Sundials_Nonlinear;
8217cb94ee6SHong Zhang   ts->ops->step            = TSStep_Sundials_Nonlinear;
8227cb94ee6SHong Zhang   ts->ops->setfromoptions  = TSSetFromOptions_Sundials_Nonlinear;
8237cb94ee6SHong Zhang 
82438f2d2fdSLisandro Dalcin   ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr);
8257adad957SLisandro Dalcin   ierr = PCCreate(((PetscObject)ts)->comm,&cvode->pc);CHKERRQ(ierr);
8267cb94ee6SHong Zhang   ierr = PetscLogObjectParent(ts,cvode->pc);CHKERRQ(ierr);
8277cb94ee6SHong Zhang   ts->data          = (void*)cvode;
8286fadb2cdSHong Zhang   cvode->cvode_type = SUNDIALS_BDF;
8296fadb2cdSHong Zhang   cvode->gtype      = SUNDIALS_CLASSICAL_GS;
8307cb94ee6SHong Zhang   cvode->restart    = 5;
8317cb94ee6SHong Zhang   cvode->linear_tol = .05;
8327cb94ee6SHong Zhang 
8337cb94ee6SHong Zhang   cvode->exact_final_time = PETSC_FALSE;
8347cb94ee6SHong Zhang 
835a07356b8SSatish Balay   ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr);
8367cb94ee6SHong Zhang   /* set tolerance for Sundials */
8377cb94ee6SHong Zhang   cvode->reltol = 1e-6;
8382c823083SHong Zhang   cvode->abstol = 1e-6;
8397cb94ee6SHong Zhang 
8407cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
8417cb94ee6SHong Zhang                     TSSundialsSetType_Sundials);CHKERRQ(ierr);
8427cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C",
8437cb94ee6SHong Zhang                     "TSSundialsSetGMRESRestart_Sundials",
8447cb94ee6SHong Zhang                     TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr);
8457cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
8467cb94ee6SHong Zhang                     "TSSundialsSetLinearTolerance_Sundials",
8477cb94ee6SHong Zhang                      TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
8487cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
8497cb94ee6SHong Zhang                     "TSSundialsSetGramSchmidtType_Sundials",
8507cb94ee6SHong Zhang                      TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
8517cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
8527cb94ee6SHong Zhang                     "TSSundialsSetTolerance_Sundials",
8537cb94ee6SHong Zhang                      TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
8547cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
8557cb94ee6SHong Zhang                     "TSSundialsGetPC_Sundials",
8567cb94ee6SHong Zhang                      TSSundialsGetPC_Sundials);CHKERRQ(ierr);
8577cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
8587cb94ee6SHong Zhang                     "TSSundialsGetIterations_Sundials",
8597cb94ee6SHong Zhang                      TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
8607cb94ee6SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C",
8617cb94ee6SHong Zhang                     "TSSundialsSetExactFinalTime_Sundials",
8627cb94ee6SHong Zhang                      TSSundialsSetExactFinalTime_Sundials);CHKERRQ(ierr);
8637cb94ee6SHong Zhang   PetscFunctionReturn(0);
8647cb94ee6SHong Zhang }
8657cb94ee6SHong Zhang EXTERN_C_END
8667cb94ee6SHong Zhang 
8677cb94ee6SHong Zhang 
8687cb94ee6SHong Zhang 
8697cb94ee6SHong Zhang 
8707cb94ee6SHong Zhang 
8717cb94ee6SHong Zhang 
8727cb94ee6SHong Zhang 
8737cb94ee6SHong Zhang 
8747cb94ee6SHong Zhang 
8757cb94ee6SHong Zhang 
876