xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision ec7bbb8b40366f23f2cdc05f6cb757f9ad47bb42)
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;
24   PetscErrorCode ierr;
25   Mat            J,P;
26   Vec            yy = cvode->w1,yydot = cvode->ydot;
27   PetscReal      gm = (PetscReal)_gamma;
28   MatStructure   str = DIFFERENT_NONZERO_PATTERN;
29   PetscScalar    *y_data;
30 
31   PetscFunctionBegin;
32   ierr = TSGetIJacobian(ts,&J,&P,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
33   y_data = (PetscScalar *) N_VGetArrayPointer(y);
34   ierr = VecPlaceArray(yy,y_data); CHKERRQ(ierr);
35   ierr = VecZeroEntries(yydot);CHKERRQ(ierr); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */
36   /* compute the shifted Jacobian   (1/gm)*I + Jrest */
37   ierr = TSComputeIJacobian(ts,ts->ptime,yy,yydot,1/gm,&J,&P,&str,PETSC_FALSE);CHKERRQ(ierr);
38   ierr = VecResetArray(yy); CHKERRQ(ierr);
39   ierr = MatScale(P,gm);CHKERRQ(ierr);  /* turn into I-gm*Jrest, J is not used by Sundials  */
40   *jcurPtr = TRUE;
41   ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
42   ierr = PCSetOperators(pc,J,P,str);CHKERRQ(ierr);
43   PetscFunctionReturn(0);
44 }
45 
46 /*
47      TSPSolve_Sundials -  routine that we provide to Sundials that applies the preconditioner.
48 */
49 #undef __FUNCT__
50 #define __FUNCT__ "TSPSolve_Sundials"
51 PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z,
52                                  realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp)
53 {
54   TS              ts = (TS) P_data;
55   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
56   PC              pc;
57   Vec             rr = cvode->w1,zz = cvode->w2;
58   PetscErrorCode  ierr;
59   PetscScalar     *r_data,*z_data;
60 
61   PetscFunctionBegin;
62   /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
63   r_data  = (PetscScalar *) N_VGetArrayPointer(r);
64   z_data  = (PetscScalar *) N_VGetArrayPointer(z);
65   ierr = VecPlaceArray(rr,r_data); CHKERRQ(ierr);
66   ierr = VecPlaceArray(zz,z_data); CHKERRQ(ierr);
67 
68   /* Solve the Px=r and put the result in zz */
69   ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
70   ierr = PCApply(pc,rr,zz); CHKERRQ(ierr);
71   ierr = VecResetArray(rr); CHKERRQ(ierr);
72   ierr = VecResetArray(zz); CHKERRQ(ierr);
73   PetscFunctionReturn(0);
74 }
75 
76 /*
77         TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
78 */
79 #undef __FUNCT__
80 #define __FUNCT__ "TSFunction_Sundials"
81 int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx)
82 {
83   TS              ts = (TS) ctx;
84   MPI_Comm        comm = ((PetscObject)ts)->comm;
85   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
86   Vec             yy = cvode->w1,yyd = cvode->w2,yydot = cvode->ydot;
87   PetscScalar     *y_data,*ydot_data;
88   PetscErrorCode  ierr;
89 
90   PetscFunctionBegin;
91   /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
92   y_data     = (PetscScalar *) N_VGetArrayPointer(y);
93   ydot_data  = (PetscScalar *) N_VGetArrayPointer(ydot);
94   ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr);
95   ierr = VecPlaceArray(yyd,ydot_data); CHKERRABORT(comm,ierr);
96   ierr = VecZeroEntries(yydot);CHKERRQ(ierr);
97 
98   /* now compute the right hand side function */
99   ierr = TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE); CHKERRABORT(comm,ierr);
100   ierr = VecScale(yyd,-1.);CHKERRQ(ierr);
101   ierr = VecResetArray(yy); CHKERRABORT(comm,ierr);
102   ierr = VecResetArray(yyd); CHKERRABORT(comm,ierr);
103   PetscFunctionReturn(0);
104 }
105 
106 /*
107        TSSolve_Sundials - Calls Sundials to integrate the ODE.
108 */
109 #undef __FUNCT__
110 #define __FUNCT__ "TSSolve_Sundials"
111 PetscErrorCode TSSolve_Sundials(TS ts)
112 {
113   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
114   Vec            sol = ts->vec_sol;
115   PetscErrorCode ierr;
116   PetscInt       i,flag;
117   long int       its,nsteps;
118   realtype       t,tout;
119   PetscScalar    *y_data;
120   void           *mem;
121 
122   PetscFunctionBegin;
123   mem  = cvode->mem;
124   tout = ts->max_time;
125   ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr);
126   N_VSetArrayPointer((realtype *)y_data,cvode->y);
127   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
128 
129   for (i = 0; i < ts->max_steps; i++) {
130     if (ts->ptime >= ts->max_time) break;
131     ierr = TSPreStep(ts);CHKERRQ(ierr);
132 
133     if (cvode->monitorstep) {
134       flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
135     } else {
136       flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL);
137     }
138 
139     if (flag){ /* display error message */
140       switch (flag){
141       case CV_ILL_INPUT:
142         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT");
143         break;
144       case CV_TOO_CLOSE:
145         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE");
146         break;
147       case CV_TOO_MUCH_WORK: {
148         PetscReal      tcur;
149         ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
150         ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr);
151         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);
152       } break;
153       case CV_TOO_MUCH_ACC:
154         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC");
155         break;
156       case CV_ERR_FAILURE:
157         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE");
158         break;
159       case CV_CONV_FAILURE:
160         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE");
161         break;
162       case CV_LINIT_FAIL:
163         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL");
164         break;
165       case CV_LSETUP_FAIL:
166         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL");
167         break;
168       case CV_LSOLVE_FAIL:
169         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL");
170         break;
171       case CV_RHSFUNC_FAIL:
172         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL");
173         break;
174       case CV_FIRST_RHSFUNC_ERR:
175         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR");
176         break;
177       case CV_REPTD_RHSFUNC_ERR:
178         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR");
179         break;
180       case CV_UNREC_RHSFUNC_ERR:
181         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR");
182         break;
183       case CV_RTFUNC_FAIL:
184         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL");
185         break;
186       default:
187         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag);
188       }
189     }
190 
191     if (t > ts->max_time && cvode->exact_final_time) {
192       /* interpolate to final requested time */
193       ierr = CVodeGetDky(mem,tout,0,cvode->y);CHKERRQ(ierr);
194       t = tout;
195     }
196 
197     /* copy the solution from cvode->y to cvode->update and sol */
198     ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr);
199     ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr);
200     ierr = VecResetArray(cvode->w1); CHKERRQ(ierr);
201     ierr = VecCopy(cvode->update,sol);CHKERRQ(ierr);
202     ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr);
203     ierr = CVSpilsGetNumLinIters(mem, &its);
204     ts->nonlinear_its = its; ts->linear_its = its;
205 
206     ts->time_step = t - ts->ptime;
207     ts->ptime     = t;
208     ts->steps++;
209 
210     ierr = TSPostStep(ts);CHKERRQ(ierr);
211     ierr = TSMonitor(ts,ts->steps,t,sol);CHKERRQ(ierr);
212   }
213   ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
214   ts->steps = nsteps;
215   PetscFunctionReturn(0);
216 }
217 
218 #undef __FUNCT__
219 #define __FUNCT__ "TSReset_Sundials"
220 PetscErrorCode TSReset_Sundials(TS ts)
221 {
222   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
223   PetscErrorCode ierr;
224 
225   PetscFunctionBegin;
226   ierr = VecDestroy(&cvode->update);CHKERRQ(ierr);
227   ierr = VecDestroy(&cvode->ydot);CHKERRQ(ierr);
228   ierr = VecDestroy(&cvode->w1);CHKERRQ(ierr);
229   ierr = VecDestroy(&cvode->w2);CHKERRQ(ierr);
230   if (cvode->mem)    {CVodeFree(&cvode->mem);}
231   PetscFunctionReturn(0);
232 }
233 
234 #undef __FUNCT__
235 #define __FUNCT__ "TSDestroy_Sundials"
236 PetscErrorCode TSDestroy_Sundials(TS ts)
237 {
238   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
239   PetscErrorCode ierr;
240 
241   PetscFunctionBegin;
242   ierr = TSReset_Sundials(ts);CHKERRQ(ierr);
243   ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr);
244   ierr = PetscFree(ts->data);CHKERRQ(ierr);
245   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","",PETSC_NULL);CHKERRQ(ierr);
246   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C","",PETSC_NULL);CHKERRQ(ierr);
247   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C","",PETSC_NULL);CHKERRQ(ierr);
248   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C","",PETSC_NULL);CHKERRQ(ierr);
249   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C","",PETSC_NULL);CHKERRQ(ierr);
250   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
251   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
252   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C","",PETSC_NULL);CHKERRQ(ierr);
253   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C","",PETSC_NULL);CHKERRQ(ierr);
254   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C","",PETSC_NULL);CHKERRQ(ierr);
255   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C","",PETSC_NULL);CHKERRQ(ierr);
256   PetscFunctionReturn(0);
257 }
258 
259 #undef __FUNCT__
260 #define __FUNCT__ "TSSetUp_Sundials"
261 PetscErrorCode TSSetUp_Sundials(TS ts)
262 {
263   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
264   PetscErrorCode ierr;
265   PetscInt       glosize,locsize,i,flag;
266   PetscScalar    *y_data,*parray;
267   void           *mem;
268   PC             pc;
269   const PCType   pctype;
270   PetscBool      pcnone;
271 
272   PetscFunctionBegin;
273 
274   /* get the vector size */
275   ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr);
276   ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr);
277 
278   /* allocate the memory for N_Vec y */
279   cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
280   if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated");
281 
282   /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
283   ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr);
284   y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y);
285   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
286   /*ierr = PetscMemcpy(y_data,parray,locsize*sizeof(PETSC_SCALAR)); CHKERRQ(ierr);*/
287   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
288   ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr);
289   ierr = VecDuplicate(ts->vec_sol,&cvode->ydot);CHKERRQ(ierr);
290   ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr);
291   ierr = PetscLogObjectParent(ts,cvode->ydot);CHKERRQ(ierr);
292 
293   /*
294     Create work vectors for the TSPSolve_Sundials() routine. Note these are
295     allocated with zero space arrays because the actual array space is provided
296     by Sundials and set using VecPlaceArray().
297   */
298   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr);
299   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr);
300   ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr);
301   ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr);
302 
303   /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
304   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
305   if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails");
306   cvode->mem = mem;
307 
308   /* Set the pointer to user-defined data */
309   flag = CVodeSetUserData(mem, ts);
310   if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails");
311 
312   /* Sundials may choose to use a smaller initial step, but will never use a larger step. */
313   flag = CVodeSetInitStep(mem,(realtype)ts->initial_time_step);
314   if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetInitStep() failed");
315   if (cvode->mindt > 0) {
316     flag = CVodeSetMinStep(mem,(realtype)cvode->mindt);
317     if (flag){
318       if (flag == CV_MEM_NULL){
319         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL");
320       } else if (flag == CV_ILL_INPUT){
321         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size");
322       } else {
323         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed");
324       }
325     }
326   }
327   if (cvode->maxdt > 0) {
328     flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt);
329     if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMaxStep() failed");
330   }
331 
332   /* Call CVodeInit to initialize the integrator memory and specify the
333    * user's right hand side function in u'=f(t,u), the inital time T0, and
334    * the initial dependent variable vector cvode->y */
335   flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y);
336   if (flag){
337     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag);
338   }
339 
340   /* specifies scalar relative and absolute tolerances */
341   flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol);
342   if (flag){
343     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag);
344   }
345 
346   /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */
347   flag = CVodeSetMaxNumSteps(mem,ts->max_steps);
348 
349   /* call CVSpgmr to use GMRES as the linear solver.        */
350   /* setup the ode integrator with the given preconditioner */
351   ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
352   ierr = PCGetType(pc,&pctype);CHKERRQ(ierr);
353   ierr = PetscTypeCompare((PetscObject)pc,PCNONE,&pcnone);CHKERRQ(ierr);
354   if (pcnone){
355     flag  = CVSpgmr(mem,PREC_NONE,0);
356     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
357   } else {
358     flag  = CVSpgmr(mem,PREC_LEFT,0);
359     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
360 
361     /* Set preconditioner and solve routines Precond and PSolve,
362      and the pointer to the user-defined block data */
363     flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials);
364     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag);
365   }
366 
367   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
368   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag);
369   PetscFunctionReturn(0);
370 }
371 
372 /* type of CVODE linear multistep method */
373 const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0};
374 /* type of G-S orthogonalization used by CVODE linear solver */
375 const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0};
376 
377 #undef __FUNCT__
378 #define __FUNCT__ "TSSetFromOptions_Sundials"
379 PetscErrorCode TSSetFromOptions_Sundials(TS ts)
380 {
381   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
382   PetscErrorCode ierr;
383   int            indx;
384   PetscBool      flag;
385 
386   PetscFunctionBegin;
387   ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr);
388     ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr);
389     if (flag) {
390       ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr);
391     }
392     ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr);
393     if (flag) {
394       ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr);
395     }
396     ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr);
397     ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr);
398     ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,PETSC_NULL);CHKERRQ(ierr);
399     ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,PETSC_NULL);CHKERRQ(ierr);
400     ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr);
401     ierr = PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);CHKERRQ(ierr);
402     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);
403     ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);CHKERRQ(ierr);
404   ierr = PetscOptionsTail();CHKERRQ(ierr);
405   PetscFunctionReturn(0);
406 }
407 
408 #undef __FUNCT__
409 #define __FUNCT__ "TSView_Sundials"
410 PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
411 {
412   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
413   PetscErrorCode ierr;
414   char           *type;
415   char           atype[] = "Adams";
416   char           btype[] = "BDF: backward differentiation formula";
417   PetscBool      iascii,isstring;
418   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
419   PetscInt       qlast,qcur;
420   PetscReal      hinused,hlast,hcur,tcur,tolsfac;
421 
422   PetscFunctionBegin;
423   if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;}
424   else                                     {type = btype;}
425 
426   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
427   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
428   if (iascii) {
429     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr);
430     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr);
431     ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr);
432     ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr);
433     ierr = PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);CHKERRQ(ierr);
434     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
435       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
436     } else {
437       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
438     }
439     if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);}
440     if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);}
441 
442     /* Outputs from CVODE, CVSPILS */
443     ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr);
444     ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr);
445     ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
446                                    &nlinsetups,&nfails,&qlast,&qcur,
447                                    &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr);
448     ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr);
449     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr);
450     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr);
451     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr);
452 
453     ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr);
454     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr);
455     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr);
456 
457     ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */
458     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr);
459     ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr);
460     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr);
461     ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr);
462     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr);
463     ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr);
464     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr);
465     ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr);
466     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr);
467     ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr);
468     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr);
469   } else if (isstring) {
470     ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr);
471   }
472   PetscFunctionReturn(0);
473 }
474 
475 
476 /* --------------------------------------------------------------------------*/
477 EXTERN_C_BEGIN
478 #undef __FUNCT__
479 #define __FUNCT__ "TSSundialsSetType_Sundials"
480 PetscErrorCode  TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
481 {
482   TS_Sundials *cvode = (TS_Sundials*)ts->data;
483 
484   PetscFunctionBegin;
485   cvode->cvode_type = type;
486   PetscFunctionReturn(0);
487 }
488 EXTERN_C_END
489 
490 EXTERN_C_BEGIN
491 #undef __FUNCT__
492 #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials"
493 PetscErrorCode  TSSundialsSetGMRESRestart_Sundials(TS ts,int restart)
494 {
495   TS_Sundials *cvode = (TS_Sundials*)ts->data;
496 
497   PetscFunctionBegin;
498   cvode->restart = restart;
499   PetscFunctionReturn(0);
500 }
501 EXTERN_C_END
502 
503 EXTERN_C_BEGIN
504 #undef __FUNCT__
505 #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials"
506 PetscErrorCode  TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
507 {
508   TS_Sundials *cvode = (TS_Sundials*)ts->data;
509 
510   PetscFunctionBegin;
511   cvode->linear_tol = tol;
512   PetscFunctionReturn(0);
513 }
514 EXTERN_C_END
515 
516 EXTERN_C_BEGIN
517 #undef __FUNCT__
518 #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials"
519 PetscErrorCode  TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
520 {
521   TS_Sundials *cvode = (TS_Sundials*)ts->data;
522 
523   PetscFunctionBegin;
524   cvode->gtype = type;
525   PetscFunctionReturn(0);
526 }
527 EXTERN_C_END
528 
529 EXTERN_C_BEGIN
530 #undef __FUNCT__
531 #define __FUNCT__ "TSSundialsSetTolerance_Sundials"
532 PetscErrorCode  TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
533 {
534   TS_Sundials *cvode = (TS_Sundials*)ts->data;
535 
536   PetscFunctionBegin;
537   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
538   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
539   PetscFunctionReturn(0);
540 }
541 EXTERN_C_END
542 
543 EXTERN_C_BEGIN
544 #undef __FUNCT__
545 #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials"
546 PetscErrorCode  TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)
547 {
548   TS_Sundials *cvode = (TS_Sundials*)ts->data;
549 
550   PetscFunctionBegin;
551   cvode->mindt = mindt;
552   PetscFunctionReturn(0);
553 }
554 EXTERN_C_END
555 
556 EXTERN_C_BEGIN
557 #undef __FUNCT__
558 #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials"
559 PetscErrorCode  TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)
560 {
561   TS_Sundials *cvode = (TS_Sundials*)ts->data;
562 
563   PetscFunctionBegin;
564   cvode->maxdt = maxdt;
565   PetscFunctionReturn(0);
566 }
567 EXTERN_C_END
568 
569 EXTERN_C_BEGIN
570 #undef __FUNCT__
571 #define __FUNCT__ "TSSundialsGetPC_Sundials"
572 PetscErrorCode  TSSundialsGetPC_Sundials(TS ts,PC *pc)
573 {
574   SNES            snes;
575   KSP             ksp;
576   PetscErrorCode  ierr;
577 
578   PetscFunctionBegin;
579   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
580   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
581   ierr = KSPGetPC(ksp,pc); CHKERRQ(ierr);
582   PetscFunctionReturn(0);
583 }
584 EXTERN_C_END
585 
586 EXTERN_C_BEGIN
587 #undef __FUNCT__
588 #define __FUNCT__ "TSSundialsGetIterations_Sundials"
589 PetscErrorCode  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  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  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  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  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  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  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  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  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  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  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  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  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  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  TSCreate_Sundials(TS ts)
981 {
982   TS_Sundials *cvode;
983   PetscErrorCode ierr;
984 
985   PetscFunctionBegin;
986   ts->ops->reset          = TSReset_Sundials;
987   ts->ops->destroy        = TSDestroy_Sundials;
988   ts->ops->view           = TSView_Sundials;
989   ts->ops->setup          = TSSetUp_Sundials;
990   ts->ops->solve          = TSSolve_Sundials;
991   ts->ops->setfromoptions = TSSetFromOptions_Sundials;
992 
993   ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr);
994   ts->data          = (void*)cvode;
995   cvode->cvode_type = SUNDIALS_BDF;
996   cvode->gtype      = SUNDIALS_CLASSICAL_GS;
997   cvode->restart    = 5;
998   cvode->linear_tol = .05;
999 
1000   cvode->exact_final_time = PETSC_TRUE;
1001   cvode->monitorstep      = PETSC_FALSE;
1002 
1003   ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr);
1004 
1005   cvode->mindt = -1.;
1006   cvode->maxdt = -1.;
1007 
1008   /* set tolerance for Sundials */
1009   cvode->reltol = 1e-6;
1010   cvode->abstol = 1e-6;
1011 
1012   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
1013                     TSSundialsSetType_Sundials);CHKERRQ(ierr);
1014   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C",
1015                     "TSSundialsSetGMRESRestart_Sundials",
1016                     TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr);
1017   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
1018                     "TSSundialsSetLinearTolerance_Sundials",
1019                      TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
1020   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
1021                     "TSSundialsSetGramSchmidtType_Sundials",
1022                      TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
1023   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
1024                     "TSSundialsSetTolerance_Sundials",
1025                      TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
1026   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C",
1027                     "TSSundialsSetMinTimeStep_Sundials",
1028                      TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr);
1029   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",
1030                     "TSSundialsSetMaxTimeStep_Sundials",
1031                      TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr);
1032   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
1033                     "TSSundialsGetPC_Sundials",
1034                      TSSundialsGetPC_Sundials);CHKERRQ(ierr);
1035   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
1036                     "TSSundialsGetIterations_Sundials",
1037                      TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
1038   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C",
1039                     "TSSundialsSetExactFinalTime_Sundials",
1040                      TSSundialsSetExactFinalTime_Sundials);CHKERRQ(ierr);
1041   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",
1042                     "TSSundialsMonitorInternalSteps_Sundials",
1043                      TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr);
1044 
1045   PetscFunctionReturn(0);
1046 }
1047 EXTERN_C_END
1048