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