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