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