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