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