xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision 829b6ff0376e8d5b8b7c6c54aef4e99b13aa00d1)
1 #define PETSCTS_DLL
2 
3 /*
4     Provides a PETSc interface to SUNDIALS/CVODE solver.
5     The interface to PVODE (old version of CVODE) was originally contributed
6     by Liyang Xu. It has been redone by Hong Zhang and Dinesh Kaushik.
7 
8     Reference: sundials-2.4.0/examples/cvode/parallel/cvDiurnal_kry_p.c
9 */
10 #include "sundials.h"  /*I "petscts.h" I*/
11 
12 /*
13       TSPrecond_Sundials - function that we provide to SUNDIALS to
14                         evaluate the preconditioner.
15 */
16 #undef __FUNCT__
17 #define __FUNCT__ "TSPrecond_Sundials"
18 PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy,
19                     booleantype jok,booleantype *jcurPtr,
20                     realtype _gamma,void *P_data,
21                     N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3)
22 {
23   TS             ts = (TS) P_data;
24   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
25   PC             pc = cvode->pc;
26   PetscErrorCode ierr;
27   Mat            Jac = ts->B;
28   Vec            yy = cvode->w1;
29   PetscScalar    one = 1.0,gm;
30   MatStructure   str = DIFFERENT_NONZERO_PATTERN;
31   PetscScalar    *y_data;
32 
33   PetscFunctionBegin;
34   /* This allows us to construct preconditioners in-place if we like */
35   ierr = MatSetUnfactored(Jac);CHKERRQ(ierr);
36 
37   /* jok - TRUE means reuse current Jacobian else recompute Jacobian */
38   if (jok) {
39     ierr     = MatCopy(cvode->pmat,Jac,str);CHKERRQ(ierr);
40     *jcurPtr = FALSE;
41   } else {
42     /* make PETSc vector yy point to SUNDIALS vector y */
43     y_data = (PetscScalar *) N_VGetArrayPointer(y);
44     ierr   = VecPlaceArray(yy,y_data); CHKERRQ(ierr);
45 
46     /* compute the Jacobian */
47     ierr = TSComputeRHSJacobian(ts,ts->ptime,yy,&Jac,&Jac,&str);CHKERRQ(ierr);
48     ierr = VecResetArray(yy); CHKERRQ(ierr);
49 
50     /* copy the Jacobian matrix */
51     if (!cvode->pmat) {
52       ierr = MatDuplicate(Jac,MAT_COPY_VALUES,&cvode->pmat);CHKERRQ(ierr);
53       ierr = PetscLogObjectParent(ts,cvode->pmat);CHKERRQ(ierr);
54     } else {
55       ierr = MatCopy(Jac,cvode->pmat,str);CHKERRQ(ierr);
56     }
57     *jcurPtr = TRUE;
58   }
59 
60   /* construct I-gamma*Jac  */
61   gm   = -_gamma;
62   ierr = MatScale(Jac,gm);CHKERRQ(ierr);
63   ierr = MatShift(Jac,one);CHKERRQ(ierr);
64 
65   ierr = PCSetOperators(pc,Jac,Jac,str);CHKERRQ(ierr);
66   PetscFunctionReturn(0);
67 }
68 
69 /*
70      TSPSolve_Sundials -  routine that we provide to Sundials that applies the preconditioner.
71 */
72 #undef __FUNCT__
73 #define __FUNCT__ "TSPSolve_Sundials"
74 PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z,
75                                  realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp)
76 {
77   TS              ts = (TS) P_data;
78   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
79   PC              pc = cvode->pc;
80   Vec             rr = cvode->w1,zz = cvode->w2;
81   PetscErrorCode  ierr;
82   PetscScalar     *r_data,*z_data;
83 
84   PetscFunctionBegin;
85   /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
86   r_data  = (PetscScalar *) N_VGetArrayPointer(r);
87   z_data  = (PetscScalar *) N_VGetArrayPointer(z);
88   ierr = VecPlaceArray(rr,r_data); CHKERRQ(ierr);
89   ierr = VecPlaceArray(zz,z_data); CHKERRQ(ierr);
90 
91   /* Solve the Px=r and put the result in zz */
92   ierr = PCApply(pc,rr,zz); CHKERRQ(ierr);
93   ierr = VecResetArray(rr); CHKERRQ(ierr);
94   ierr = VecResetArray(zz); CHKERRQ(ierr);
95   PetscFunctionReturn(0);
96 }
97 
98 /*
99         TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
100 */
101 #undef __FUNCT__
102 #define __FUNCT__ "TSFunction_Sundials"
103 int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx)
104 {
105   TS              ts = (TS) ctx;
106   MPI_Comm        comm = ((PetscObject)ts)->comm;
107   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
108   Vec             yy = cvode->w1,yyd = cvode->w2;
109   PetscScalar     *y_data,*ydot_data;
110   PetscErrorCode  ierr;
111 
112   PetscFunctionBegin;
113   /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
114   y_data     = (PetscScalar *) N_VGetArrayPointer(y);
115   ydot_data  = (PetscScalar *) N_VGetArrayPointer(ydot);
116   ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr);
117   ierr = VecPlaceArray(yyd,ydot_data); CHKERRABORT(comm,ierr);
118 
119   /* now compute the right hand side function */
120   ierr = TSComputeRHSFunction(ts,t,yy,yyd); CHKERRABORT(comm,ierr);
121   ierr = VecResetArray(yy); CHKERRABORT(comm,ierr);
122   ierr = VecResetArray(yyd); CHKERRABORT(comm,ierr);
123   PetscFunctionReturn(0);
124 }
125 
126 /*
127        TSStep_Sundials_Nonlinear - Calls Sundials to integrate the ODE.
128 */
129 #undef __FUNCT__
130 #define __FUNCT__ "TSStep_Sundials_Nonlinear"
131 PetscErrorCode TSStep_Sundials_Nonlinear(TS ts,int *steps,double *time)
132 {
133   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
134   Vec            sol = ts->vec_sol;
135   PetscErrorCode ierr;
136   PetscInt       i,max_steps = ts->max_steps,flag;
137   long int       its;
138   realtype       t,tout;
139   PetscScalar    *y_data;
140   void           *mem;
141 
142   PetscFunctionBegin;
143   mem  = cvode->mem;
144   tout = ts->max_time;
145   ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr);
146   N_VSetArrayPointer((realtype *)y_data,cvode->y);
147   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
148   for (i = 0; i < max_steps; i++) {
149     if (ts->ptime >= ts->max_time) break;
150     if (cvode->monitorstep){
151       flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
152     } else {
153       flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL);
154     }
155     if (flag)SETERRQ1(PETSC_ERR_LIB,"CVode() fails, flag %d",flag);
156     if (t > ts->max_time && cvode->exact_final_time) {
157       /* interpolate to final requested time */
158       ierr = CVodeGetDky(mem,tout,0,cvode->y);CHKERRQ(ierr);
159       t = tout;
160     }
161     ts->time_step = t - ts->ptime;
162     ts->ptime     = t;
163 
164     /* copy the solution from cvode->y to cvode->update and sol */
165     ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr);
166     ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr);
167     ierr = VecResetArray(cvode->w1); CHKERRQ(ierr);
168     ierr = VecCopy(cvode->update,sol);CHKERRQ(ierr);
169     ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr);
170     ts->nonlinear_its = its;
171     ierr = CVSpilsGetNumLinIters(mem, &its);
172     ts->linear_its = its;
173     ts->steps++;
174     ierr = TSMonitor(ts,ts->steps,t,sol);CHKERRQ(ierr);
175   }
176   *steps += ts->steps;
177   *time   = t;
178   PetscFunctionReturn(0);
179 }
180 
181 #undef __FUNCT__
182 #define __FUNCT__ "TSDestroy_Sundials"
183 PetscErrorCode TSDestroy_Sundials(TS ts)
184 {
185   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
186   PetscErrorCode ierr;
187 
188   PetscFunctionBegin;
189   if (cvode->pmat)   {ierr = MatDestroy(cvode->pmat);CHKERRQ(ierr);}
190   if (cvode->pc)     {ierr = PCDestroy(cvode->pc);CHKERRQ(ierr);}
191   if (cvode->update) {ierr = VecDestroy(cvode->update);CHKERRQ(ierr);}
192   if (cvode->func)   {ierr = VecDestroy(cvode->func);CHKERRQ(ierr);}
193   if (cvode->rhs)    {ierr = VecDestroy(cvode->rhs);CHKERRQ(ierr);}
194   if (cvode->w1)     {ierr = VecDestroy(cvode->w1);CHKERRQ(ierr);}
195   if (cvode->w2)     {ierr = VecDestroy(cvode->w2);CHKERRQ(ierr);}
196   ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr);
197   if (cvode->mem) {CVodeFree(&cvode->mem);}
198   ierr = PetscFree(cvode);CHKERRQ(ierr);
199   PetscFunctionReturn(0);
200 }
201 
202 #undef __FUNCT__
203 #define __FUNCT__ "TSSetUp_Sundials_Nonlinear"
204 PetscErrorCode TSSetUp_Sundials_Nonlinear(TS ts)
205 {
206   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
207   PetscErrorCode ierr;
208   PetscInt       glosize,locsize,i,flag;
209   PetscScalar    *y_data,*parray;
210   void           *mem;
211   const PCType   pctype;
212   PetscTruth     pcnone;
213   Vec            sol = ts->vec_sol;
214 
215   PetscFunctionBegin;
216   ierr = PCSetFromOptions(cvode->pc);CHKERRQ(ierr);
217   /* get the vector size */
218   ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr);
219   ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr);
220 
221   /* allocate the memory for N_Vec y */
222   cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
223   if (!cvode->y) SETERRQ(1,"cvode->y is not allocated");
224 
225   /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
226   ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr);
227   y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y);
228   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
229   /*ierr = PetscMemcpy(y_data,parray,locsize*sizeof(PETSC_SCALAR)); CHKERRQ(ierr);*/
230   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
231   ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr);
232   ierr = VecDuplicate(ts->vec_sol,&cvode->func);CHKERRQ(ierr);
233   ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr);
234   ierr = PetscLogObjectParent(ts,cvode->func);CHKERRQ(ierr);
235 
236   /*
237     Create work vectors for the TSPSolve_Sundials() routine. Note these are
238     allocated with zero space arrays because the actual array space is provided
239     by Sundials and set using VecPlaceArray().
240   */
241   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr);
242   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr);
243   ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr);
244   ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr);
245 
246   /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
247   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
248   if (!mem) SETERRQ(PETSC_ERR_MEM,"CVodeCreate() fails");
249   cvode->mem = mem;
250 
251   /* Set the pointer to user-defined data */
252   flag = CVodeSetUserData(mem, ts);
253   if (flag) SETERRQ(PETSC_ERR_LIB,"CVodeSetUserData() fails");
254 
255   /* Call CVodeInit to initialize the integrator memory and specify the
256    * user's right hand side function in u'=f(t,u), the inital time T0, and
257    * the initial dependent variable vector cvode->y */
258   flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y);
259   if (flag){
260     SETERRQ1(PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag);
261   }
262 
263   flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol);
264   if (flag){
265     SETERRQ1(PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag);
266   }
267 
268   /* initialize the number of steps */
269   ierr   = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
270 
271   /* call CVSpgmr to use GMRES as the linear solver.        */
272   /* setup the ode integrator with the given preconditioner */
273   ierr = PCGetType(cvode->pc,&pctype);CHKERRQ(ierr);
274   ierr = PetscTypeCompare((PetscObject)cvode->pc,PCNONE,&pcnone);CHKERRQ(ierr);
275   if (pcnone){
276     flag  = CVSpgmr(mem,PREC_NONE,0);
277     if (flag) SETERRQ1(PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
278   } else {
279     flag  = CVSpgmr(mem,PREC_LEFT,0);
280     if (flag) SETERRQ1(PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
281 
282     /* Set preconditioner and solve routines Precond and PSolve,
283      and the pointer to the user-defined block data */
284     flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials);
285     if (flag) SETERRQ1(PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag);
286   }
287 
288   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
289   if (flag) SETERRQ1(PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag);
290   PetscFunctionReturn(0);
291 }
292 
293 /* type of CVODE linear multistep method */
294 const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0};
295 /* type of G-S orthogonalization used by CVODE linear solver */
296 const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0};
297 
298 #undef __FUNCT__
299 #define __FUNCT__ "TSSetFromOptions_Sundials_Nonlinear"
300 PetscErrorCode TSSetFromOptions_Sundials_Nonlinear(TS ts)
301 {
302   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
303   PetscErrorCode ierr;
304   int            indx;
305   PetscTruth     flag;
306 
307   PetscFunctionBegin;
308   ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr);
309     ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr);
310     if (flag) {
311       ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr);
312     }
313     ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr);
314     if (flag) {
315       ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr);
316     }
317     ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr);
318     ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr);
319     ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr);
320     ierr = PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);CHKERRQ(ierr);
321     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);
322     ierr = PetscOptionsTruth("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);CHKERRQ(ierr);
323   ierr = PetscOptionsTail();CHKERRQ(ierr);
324   PetscFunctionReturn(0);
325 }
326 
327 #undef __FUNCT__
328 #define __FUNCT__ "TSView_Sundials"
329 PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
330 {
331   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
332   PetscErrorCode ierr;
333   char           *type;
334   char           atype[] = "Adams";
335   char           btype[] = "BDF: backward differentiation formula";
336   PetscTruth     iascii,isstring;
337   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
338   PetscInt       qlast,qcur;
339   PetscReal      hinused,hlast,hcur,tcur,tolsfac;
340 
341   PetscFunctionBegin;
342   if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;}
343   else                                     {type = btype;}
344 
345   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
346   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
347   if (iascii) {
348     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr);
349     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr);
350     ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr);
351     ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr);
352     ierr = PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);CHKERRQ(ierr);
353     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
354       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
355     } else {
356       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
357     }
358 
359     /* Outputs from CVODE, CVSPILS */
360     ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr);
361     ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr);
362     ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
363                                    &nlinsetups,&nfails,&qlast,&qcur,
364                                    &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr);
365     ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr);
366     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr);
367     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr);
368     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr);
369 
370     ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr);
371     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr);
372     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr);
373 
374     ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */
375     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr);
376     ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr);
377     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr);
378     ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr);
379     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr);
380     ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr);
381     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr);
382     ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr);
383     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr);
384     ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr);
385     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr);
386   } else if (isstring) {
387     ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr);
388   } else {
389     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by TS Sundials",((PetscObject)viewer)->type_name);
390   }
391   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
392   ierr = PCView(cvode->pc,viewer);CHKERRQ(ierr);
393   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
394   PetscFunctionReturn(0);
395 }
396 
397 
398 /* --------------------------------------------------------------------------*/
399 EXTERN_C_BEGIN
400 #undef __FUNCT__
401 #define __FUNCT__ "TSSundialsSetType_Sundials"
402 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
403 {
404   TS_Sundials *cvode = (TS_Sundials*)ts->data;
405 
406   PetscFunctionBegin;
407   cvode->cvode_type = type;
408   PetscFunctionReturn(0);
409 }
410 EXTERN_C_END
411 
412 EXTERN_C_BEGIN
413 #undef __FUNCT__
414 #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials"
415 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart_Sundials(TS ts,int restart)
416 {
417   TS_Sundials *cvode = (TS_Sundials*)ts->data;
418 
419   PetscFunctionBegin;
420   cvode->restart = restart;
421   PetscFunctionReturn(0);
422 }
423 EXTERN_C_END
424 
425 EXTERN_C_BEGIN
426 #undef __FUNCT__
427 #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials"
428 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
429 {
430   TS_Sundials *cvode = (TS_Sundials*)ts->data;
431 
432   PetscFunctionBegin;
433   cvode->linear_tol = tol;
434   PetscFunctionReturn(0);
435 }
436 EXTERN_C_END
437 
438 EXTERN_C_BEGIN
439 #undef __FUNCT__
440 #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials"
441 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
442 {
443   TS_Sundials *cvode = (TS_Sundials*)ts->data;
444 
445   PetscFunctionBegin;
446   cvode->gtype = type;
447   PetscFunctionReturn(0);
448 }
449 EXTERN_C_END
450 
451 EXTERN_C_BEGIN
452 #undef __FUNCT__
453 #define __FUNCT__ "TSSundialsSetTolerance_Sundials"
454 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
455 {
456   TS_Sundials *cvode = (TS_Sundials*)ts->data;
457 
458   PetscFunctionBegin;
459   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
460   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
461   PetscFunctionReturn(0);
462 }
463 EXTERN_C_END
464 
465 EXTERN_C_BEGIN
466 #undef __FUNCT__
467 #define __FUNCT__ "TSSundialsGetPC_Sundials"
468 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC_Sundials(TS ts,PC *pc)
469 {
470   TS_Sundials *cvode = (TS_Sundials*)ts->data;
471 
472   PetscFunctionBegin;
473   *pc = cvode->pc;
474   PetscFunctionReturn(0);
475 }
476 EXTERN_C_END
477 
478 EXTERN_C_BEGIN
479 #undef __FUNCT__
480 #define __FUNCT__ "TSSundialsGetIterations_Sundials"
481 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
482 {
483   PetscFunctionBegin;
484   if (nonlin) *nonlin = ts->nonlinear_its;
485   if (lin)    *lin    = ts->linear_its;
486   PetscFunctionReturn(0);
487 }
488 EXTERN_C_END
489 
490 EXTERN_C_BEGIN
491 #undef __FUNCT__
492 #define __FUNCT__ "TSSundialsSetExactFinalTime_Sundials"
493 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime_Sundials(TS ts,PetscTruth s)
494 {
495   TS_Sundials *cvode = (TS_Sundials*)ts->data;
496 
497   PetscFunctionBegin;
498   cvode->exact_final_time = s;
499   PetscFunctionReturn(0);
500 }
501 EXTERN_C_END
502 
503 EXTERN_C_BEGIN
504 #undef __FUNCT__
505 #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials"
506 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscTruth s)
507 {
508   TS_Sundials *cvode = (TS_Sundials*)ts->data;
509 
510   PetscFunctionBegin;
511   cvode->monitorstep = s;
512   PetscFunctionReturn(0);
513 }
514 EXTERN_C_END
515 /* -------------------------------------------------------------------------------------------*/
516 
517 #undef __FUNCT__
518 #define __FUNCT__ "TSSundialsGetIterations"
519 /*@C
520    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
521 
522    Not Collective
523 
524    Input parameters:
525 .    ts     - the time-step context
526 
527    Output Parameters:
528 +   nonlin - number of nonlinear iterations
529 -   lin    - number of linear iterations
530 
531    Level: advanced
532 
533    Notes:
534     These return the number since the creation of the TS object
535 
536 .keywords: non-linear iterations, linear iterations
537 
538 .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(),
539           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
540           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
541           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
542           TSSundialsSetExactFinalTime()
543 
544 @*/
545 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
546 {
547   PetscErrorCode ierr,(*f)(TS,int*,int*);
548 
549   PetscFunctionBegin;
550   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetIterations_C",(void (**)(void))&f);CHKERRQ(ierr);
551   if (f) {
552     ierr = (*f)(ts,nonlin,lin);CHKERRQ(ierr);
553   }
554   PetscFunctionReturn(0);
555 }
556 
557 #undef __FUNCT__
558 #define __FUNCT__ "TSSundialsSetType"
559 /*@
560    TSSundialsSetType - Sets the method that Sundials will use for integration.
561 
562    Collective on TS
563 
564    Input parameters:
565 +    ts     - the time-step context
566 -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
567 
568    Level: intermediate
569 
570 .keywords: Adams, backward differentiation formula
571 
572 .seealso: TSSundialsGetIterations(),  TSSundialsSetGMRESRestart(),
573           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
574           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
575           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
576           TSSundialsSetExactFinalTime()
577 @*/
578 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType(TS ts,TSSundialsLmmType type)
579 {
580   PetscErrorCode ierr,(*f)(TS,TSSundialsLmmType);
581 
582   PetscFunctionBegin;
583   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
584   if (f) {
585     ierr = (*f)(ts,type);CHKERRQ(ierr);
586   }
587   PetscFunctionReturn(0);
588 }
589 
590 #undef __FUNCT__
591 #define __FUNCT__ "TSSundialsSetGMRESRestart"
592 /*@
593    TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by
594        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
595        this is ALSO the maximum number of GMRES steps that will be used.
596 
597    Collective on TS
598 
599    Input parameters:
600 +    ts      - the time-step context
601 -    restart - number of direction vectors (the restart size).
602 
603    Level: advanced
604 
605 .keywords: GMRES, restart
606 
607 .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
608           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
609           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
610           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
611           TSSundialsSetExactFinalTime()
612 
613 @*/
614 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart(TS ts,int restart)
615 {
616   PetscErrorCode ierr,(*f)(TS,int);
617 
618   PetscFunctionBegin;
619   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGMRESRestart_C",(void (**)(void))&f);CHKERRQ(ierr);
620   if (f) {
621     ierr = (*f)(ts,restart);CHKERRQ(ierr);
622   }
623   PetscFunctionReturn(0);
624 }
625 
626 #undef __FUNCT__
627 #define __FUNCT__ "TSSundialsSetLinearTolerance"
628 /*@
629    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
630        system by SUNDIALS.
631 
632    Collective on TS
633 
634    Input parameters:
635 +    ts     - the time-step context
636 -    tol    - the factor by which the tolerance on the nonlinear solver is
637              multiplied to get the tolerance on the linear solver, .05 by default.
638 
639    Level: advanced
640 
641 .keywords: GMRES, linear convergence tolerance, SUNDIALS
642 
643 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
644           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
645           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
646           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
647           TSSundialsSetExactFinalTime()
648 
649 @*/
650 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance(TS ts,double tol)
651 {
652   PetscErrorCode ierr,(*f)(TS,double);
653 
654   PetscFunctionBegin;
655   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
656   if (f) {
657     ierr = (*f)(ts,tol);CHKERRQ(ierr);
658   }
659   PetscFunctionReturn(0);
660 }
661 
662 #undef __FUNCT__
663 #define __FUNCT__ "TSSundialsSetGramSchmidtType"
664 /*@
665    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
666         in GMRES method by SUNDIALS linear solver.
667 
668    Collective on TS
669 
670    Input parameters:
671 +    ts  - the time-step context
672 -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
673 
674    Level: advanced
675 
676 .keywords: Sundials, orthogonalization
677 
678 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
679           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
680           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
681           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
682           TSSundialsSetExactFinalTime()
683 
684 @*/
685 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
686 {
687   PetscErrorCode ierr,(*f)(TS,TSSundialsGramSchmidtType);
688 
689   PetscFunctionBegin;
690   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",(void (**)(void))&f);CHKERRQ(ierr);
691   if (f) {
692     ierr = (*f)(ts,type);CHKERRQ(ierr);
693   }
694   PetscFunctionReturn(0);
695 }
696 
697 #undef __FUNCT__
698 #define __FUNCT__ "TSSundialsSetTolerance"
699 /*@
700    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
701                          Sundials for error control.
702 
703    Collective on TS
704 
705    Input parameters:
706 +    ts  - the time-step context
707 .    aabs - the absolute tolerance
708 -    rel - the relative tolerance
709 
710      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
711     these regulate the size of the error for a SINGLE timestep.
712 
713    Level: intermediate
714 
715 .keywords: Sundials, tolerance
716 
717 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
718           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
719           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
720           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
721           TSSundialsSetExactFinalTime()
722 
723 @*/
724 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance(TS ts,double aabs,double rel)
725 {
726   PetscErrorCode ierr,(*f)(TS,double,double);
727 
728   PetscFunctionBegin;
729   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
730   if (f) {
731     ierr = (*f)(ts,aabs,rel);CHKERRQ(ierr);
732   }
733   PetscFunctionReturn(0);
734 }
735 
736 #undef __FUNCT__
737 #define __FUNCT__ "TSSundialsGetPC"
738 /*@
739    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
740 
741    Input Parameter:
742 .    ts - the time-step context
743 
744    Output Parameter:
745 .    pc - the preconditioner context
746 
747    Level: advanced
748 
749 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
750           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
751           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
752           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
753 @*/
754 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC(TS ts,PC *pc)
755 {
756   PetscErrorCode ierr,(*f)(TS,PC *);
757 
758   PetscFunctionBegin;
759   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetPC_C",(void (**)(void))&f);CHKERRQ(ierr);
760   if (f) {
761     ierr = (*f)(ts,pc);CHKERRQ(ierr);
762   } else {
763     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"TS must be of Sundials type to extract the PC");
764   }
765 
766   PetscFunctionReturn(0);
767 }
768 
769 #undef __FUNCT__
770 #define __FUNCT__ "TSSundialsSetExactFinalTime"
771 /*@
772    TSSundialsSetExactFinalTime - Determines if Sundials interpolates solution to the
773       exact final time requested by the user or just returns it at the final time
774       it computed. (Defaults to true).
775 
776    Input Parameter:
777 +   ts - the time-step context
778 -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE
779 
780    Level: beginner
781 
782 .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
783           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
784           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
785           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
786 @*/
787 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime(TS ts,PetscTruth ft)
788 {
789   PetscErrorCode ierr,(*f)(TS,PetscTruth);
790 
791   PetscFunctionBegin;
792   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetExactFinalTime_C",(void (**)(void))&f);CHKERRQ(ierr);
793   if (f) {
794     ierr = (*f)(ts,ft);CHKERRQ(ierr);
795   }
796   PetscFunctionReturn(0);
797 }
798 
799 #undef __FUNCT__
800 #define __FUNCT__ "TSSundialsMonitorInternalSteps"
801 /*@
802    TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false).
803 
804    Input Parameter:
805 +   ts - the time-step context
806 -   ft - PETSC_TRUE if monitor, else PETSC_FALSE
807 
808    Level: beginner
809 
810 .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
811           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
812           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
813           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
814 @*/
815 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsMonitorInternalSteps(TS ts,PetscTruth ft)
816 {
817   PetscErrorCode ierr,(*f)(TS,PetscTruth);
818 
819   PetscFunctionBegin;
820   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",(void (**)(void))&f);CHKERRQ(ierr);
821   if (f) {
822     ierr = (*f)(ts,ft);CHKERRQ(ierr);
823   }
824   PetscFunctionReturn(0);
825 }
826 /* -------------------------------------------------------------------------------------------*/
827 /*MC
828       TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
829 
830    Options Database:
831 +    -ts_sundials_type <bdf,adams>
832 .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
833 .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
834 .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
835 .    -ts_sundials_linear_tolerance <tol>
836 .    -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions
837 .    -ts_sundials_exact_final_time - Interpolate output to stop exactly at the final time
838 -    -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps
839 
840 
841     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
842            only PETSc PC options
843 
844     Level: beginner
845 
846 .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(),
847            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSundialsSetExactFinalTime()
848 
849 M*/
850 EXTERN_C_BEGIN
851 #undef __FUNCT__
852 #define __FUNCT__ "TSCreate_Sundials"
853 PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Sundials(TS ts)
854 {
855   TS_Sundials *cvode;
856   PetscErrorCode ierr;
857 
858   PetscFunctionBegin;
859   if (ts->problem_type != TS_NONLINEAR) {
860     SETERRQ(PETSC_ERR_SUP,"Only support for nonlinear problems");
861   }
862   ts->ops->destroy        = TSDestroy_Sundials;
863   ts->ops->view           = TSView_Sundials;
864   ts->ops->setup          = TSSetUp_Sundials_Nonlinear;
865   ts->ops->step           = TSStep_Sundials_Nonlinear;
866   ts->ops->setfromoptions = TSSetFromOptions_Sundials_Nonlinear;
867 
868   ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr);
869   ierr = PCCreate(((PetscObject)ts)->comm,&cvode->pc);CHKERRQ(ierr);
870   ierr = PetscLogObjectParent(ts,cvode->pc);CHKERRQ(ierr);
871   ts->data          = (void*)cvode;
872   cvode->cvode_type = SUNDIALS_BDF;
873   cvode->gtype      = SUNDIALS_CLASSICAL_GS;
874   cvode->restart    = 5;
875   cvode->linear_tol = .05;
876 
877   cvode->exact_final_time = PETSC_TRUE;
878   cvode->monitorstep      = PETSC_FALSE;
879 
880   ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr);
881   /* set tolerance for Sundials */
882   cvode->reltol = 1e-6;
883   cvode->abstol = 1e-6;
884 
885   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
886                     TSSundialsSetType_Sundials);CHKERRQ(ierr);
887   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C",
888                     "TSSundialsSetGMRESRestart_Sundials",
889                     TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr);
890   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
891                     "TSSundialsSetLinearTolerance_Sundials",
892                      TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
893   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
894                     "TSSundialsSetGramSchmidtType_Sundials",
895                      TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
896   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
897                     "TSSundialsSetTolerance_Sundials",
898                      TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
899   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
900                     "TSSundialsGetPC_Sundials",
901                      TSSundialsGetPC_Sundials);CHKERRQ(ierr);
902   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
903                     "TSSundialsGetIterations_Sundials",
904                      TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
905   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C",
906                     "TSSundialsSetExactFinalTime_Sundials",
907                      TSSundialsSetExactFinalTime_Sundials);CHKERRQ(ierr);
908   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",
909                     "TSSundialsMonitorInternalSteps_Sundials",
910                      TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr);
911 
912   PetscFunctionReturn(0);
913 }
914 EXTERN_C_END
915 
916 
917 
918 
919 
920 
921 
922 
923 
924 
925