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