xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision 117ce3f2ee664074fd8b517fcd34f2f53338bcd0)
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   PetscFunctionReturn(0);
215 }
216 
217 #undef __FUNCT__
218 #define __FUNCT__ "TSReset_Sundials"
219 PetscErrorCode TSReset_Sundials(TS ts)
220 {
221   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
222   PetscErrorCode ierr;
223 
224   PetscFunctionBegin;
225   ierr = VecDestroy(&cvode->update);CHKERRQ(ierr);
226   ierr = VecDestroy(&cvode->ydot);CHKERRQ(ierr);
227   ierr = VecDestroy(&cvode->w1);CHKERRQ(ierr);
228   ierr = VecDestroy(&cvode->w2);CHKERRQ(ierr);
229   if (cvode->mem)    {CVodeFree(&cvode->mem);}
230   PetscFunctionReturn(0);
231 }
232 
233 #undef __FUNCT__
234 #define __FUNCT__ "TSDestroy_Sundials"
235 PetscErrorCode TSDestroy_Sundials(TS ts)
236 {
237   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
238   PetscErrorCode ierr;
239 
240   PetscFunctionBegin;
241   ierr = TSReset_Sundials(ts);CHKERRQ(ierr);
242   ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr);
243   ierr = PetscFree(ts->data);CHKERRQ(ierr);
244   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","",PETSC_NULL);CHKERRQ(ierr);
245   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C","",PETSC_NULL);CHKERRQ(ierr);
246   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C","",PETSC_NULL);CHKERRQ(ierr);
247   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C","",PETSC_NULL);CHKERRQ(ierr);
248   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C","",PETSC_NULL);CHKERRQ(ierr);
249   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
250   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
251   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C","",PETSC_NULL);CHKERRQ(ierr);
252   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C","",PETSC_NULL);CHKERRQ(ierr);
253   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C","",PETSC_NULL);CHKERRQ(ierr);
254   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C","",PETSC_NULL);CHKERRQ(ierr);
255   PetscFunctionReturn(0);
256 }
257 
258 #undef __FUNCT__
259 #define __FUNCT__ "TSSetUp_Sundials"
260 PetscErrorCode TSSetUp_Sundials(TS ts)
261 {
262   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
263   PetscErrorCode ierr;
264   PetscInt       glosize,locsize,i,flag;
265   PetscScalar    *y_data,*parray;
266   void           *mem;
267   PC             pc;
268   const PCType   pctype;
269   PetscBool      pcnone;
270 
271   PetscFunctionBegin;
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->ydot);CHKERRQ(ierr);
289   ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr);
290   ierr = PetscLogObjectParent(ts,cvode->ydot);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 
348   /* call CVSpgmr to use GMRES as the linear solver.        */
349   /* setup the ode integrator with the given preconditioner */
350   ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
351   ierr = PCGetType(pc,&pctype);CHKERRQ(ierr);
352   ierr = PetscTypeCompare((PetscObject)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"
378 PetscErrorCode TSSetFromOptions_Sundials(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 = 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);
402     ierr = PetscOptionsBool("-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 = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
475   PetscFunctionReturn(0);
476 }
477 
478 
479 /* --------------------------------------------------------------------------*/
480 EXTERN_C_BEGIN
481 #undef __FUNCT__
482 #define __FUNCT__ "TSSundialsSetType_Sundials"
483 PetscErrorCode  TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
484 {
485   TS_Sundials *cvode = (TS_Sundials*)ts->data;
486 
487   PetscFunctionBegin;
488   cvode->cvode_type = type;
489   PetscFunctionReturn(0);
490 }
491 EXTERN_C_END
492 
493 EXTERN_C_BEGIN
494 #undef __FUNCT__
495 #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials"
496 PetscErrorCode  TSSundialsSetGMRESRestart_Sundials(TS ts,int restart)
497 {
498   TS_Sundials *cvode = (TS_Sundials*)ts->data;
499 
500   PetscFunctionBegin;
501   cvode->restart = restart;
502   PetscFunctionReturn(0);
503 }
504 EXTERN_C_END
505 
506 EXTERN_C_BEGIN
507 #undef __FUNCT__
508 #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials"
509 PetscErrorCode  TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
510 {
511   TS_Sundials *cvode = (TS_Sundials*)ts->data;
512 
513   PetscFunctionBegin;
514   cvode->linear_tol = tol;
515   PetscFunctionReturn(0);
516 }
517 EXTERN_C_END
518 
519 EXTERN_C_BEGIN
520 #undef __FUNCT__
521 #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials"
522 PetscErrorCode  TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
523 {
524   TS_Sundials *cvode = (TS_Sundials*)ts->data;
525 
526   PetscFunctionBegin;
527   cvode->gtype = type;
528   PetscFunctionReturn(0);
529 }
530 EXTERN_C_END
531 
532 EXTERN_C_BEGIN
533 #undef __FUNCT__
534 #define __FUNCT__ "TSSundialsSetTolerance_Sundials"
535 PetscErrorCode  TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
536 {
537   TS_Sundials *cvode = (TS_Sundials*)ts->data;
538 
539   PetscFunctionBegin;
540   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
541   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
542   PetscFunctionReturn(0);
543 }
544 EXTERN_C_END
545 
546 EXTERN_C_BEGIN
547 #undef __FUNCT__
548 #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials"
549 PetscErrorCode  TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)
550 {
551   TS_Sundials *cvode = (TS_Sundials*)ts->data;
552 
553   PetscFunctionBegin;
554   cvode->mindt = mindt;
555   PetscFunctionReturn(0);
556 }
557 EXTERN_C_END
558 
559 EXTERN_C_BEGIN
560 #undef __FUNCT__
561 #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials"
562 PetscErrorCode  TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)
563 {
564   TS_Sundials *cvode = (TS_Sundials*)ts->data;
565 
566   PetscFunctionBegin;
567   cvode->maxdt = maxdt;
568   PetscFunctionReturn(0);
569 }
570 EXTERN_C_END
571 
572 EXTERN_C_BEGIN
573 #undef __FUNCT__
574 #define __FUNCT__ "TSSundialsGetPC_Sundials"
575 PetscErrorCode  TSSundialsGetPC_Sundials(TS ts,PC *pc)
576 {
577   SNES            snes;
578   KSP             ksp;
579   PetscErrorCode  ierr;
580 
581   PetscFunctionBegin;
582   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
583   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
584   ierr = KSPGetPC(ksp,pc); CHKERRQ(ierr);
585   PetscFunctionReturn(0);
586 }
587 EXTERN_C_END
588 
589 EXTERN_C_BEGIN
590 #undef __FUNCT__
591 #define __FUNCT__ "TSSundialsGetIterations_Sundials"
592 PetscErrorCode  TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
593 {
594   PetscFunctionBegin;
595   if (nonlin) *nonlin = ts->nonlinear_its;
596   if (lin)    *lin    = ts->linear_its;
597   PetscFunctionReturn(0);
598 }
599 EXTERN_C_END
600 
601 EXTERN_C_BEGIN
602 #undef __FUNCT__
603 #define __FUNCT__ "TSSundialsSetExactFinalTime_Sundials"
604 PetscErrorCode  TSSundialsSetExactFinalTime_Sundials(TS ts,PetscBool  s)
605 {
606   TS_Sundials *cvode = (TS_Sundials*)ts->data;
607 
608   PetscFunctionBegin;
609   cvode->exact_final_time = s;
610   PetscFunctionReturn(0);
611 }
612 EXTERN_C_END
613 
614 EXTERN_C_BEGIN
615 #undef __FUNCT__
616 #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials"
617 PetscErrorCode  TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool  s)
618 {
619   TS_Sundials *cvode = (TS_Sundials*)ts->data;
620 
621   PetscFunctionBegin;
622   cvode->monitorstep = s;
623   PetscFunctionReturn(0);
624 }
625 EXTERN_C_END
626 /* -------------------------------------------------------------------------------------------*/
627 
628 #undef __FUNCT__
629 #define __FUNCT__ "TSSundialsGetIterations"
630 /*@C
631    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
632 
633    Not Collective
634 
635    Input parameters:
636 .    ts     - the time-step context
637 
638    Output Parameters:
639 +   nonlin - number of nonlinear iterations
640 -   lin    - number of linear iterations
641 
642    Level: advanced
643 
644    Notes:
645     These return the number since the creation of the TS object
646 
647 .keywords: non-linear iterations, linear iterations
648 
649 .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(),
650           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
651           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
652           TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSundialsSetExactFinalTime()
653 
654 @*/
655 PetscErrorCode  TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
656 {
657   PetscErrorCode ierr;
658 
659   PetscFunctionBegin;
660   ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr);
661   PetscFunctionReturn(0);
662 }
663 
664 #undef __FUNCT__
665 #define __FUNCT__ "TSSundialsSetType"
666 /*@
667    TSSundialsSetType - Sets the method that Sundials will use for integration.
668 
669    Logically Collective on TS
670 
671    Input parameters:
672 +    ts     - the time-step context
673 -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
674 
675    Level: intermediate
676 
677 .keywords: Adams, backward differentiation formula
678 
679 .seealso: TSSundialsGetIterations(),  TSSundialsSetGMRESRestart(),
680           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
681           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
682           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
683           TSSundialsSetExactFinalTime()
684 @*/
685 PetscErrorCode  TSSundialsSetType(TS ts,TSSundialsLmmType type)
686 {
687   PetscErrorCode ierr;
688 
689   PetscFunctionBegin;
690   ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr);
691   PetscFunctionReturn(0);
692 }
693 
694 #undef __FUNCT__
695 #define __FUNCT__ "TSSundialsSetGMRESRestart"
696 /*@
697    TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by
698        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
699        this is ALSO the maximum number of GMRES steps that will be used.
700 
701    Logically Collective on TS
702 
703    Input parameters:
704 +    ts      - the time-step context
705 -    restart - number of direction vectors (the restart size).
706 
707    Level: advanced
708 
709 .keywords: GMRES, restart
710 
711 .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
712           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
713           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
714           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
715           TSSundialsSetExactFinalTime()
716 
717 @*/
718 PetscErrorCode  TSSundialsSetGMRESRestart(TS ts,int restart)
719 {
720   PetscErrorCode ierr;
721 
722   PetscFunctionBegin;
723   PetscValidLogicalCollectiveInt(ts,restart,2);
724   ierr = PetscTryMethod(ts,"TSSundialsSetGMRESRestart_C",(TS,int),(ts,restart));CHKERRQ(ierr);
725   PetscFunctionReturn(0);
726 }
727 
728 #undef __FUNCT__
729 #define __FUNCT__ "TSSundialsSetLinearTolerance"
730 /*@
731    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
732        system by SUNDIALS.
733 
734    Logically Collective on TS
735 
736    Input parameters:
737 +    ts     - the time-step context
738 -    tol    - the factor by which the tolerance on the nonlinear solver is
739              multiplied to get the tolerance on the linear solver, .05 by default.
740 
741    Level: advanced
742 
743 .keywords: GMRES, linear convergence tolerance, SUNDIALS
744 
745 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
746           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
747           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
748           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
749           TSSundialsSetExactFinalTime()
750 
751 @*/
752 PetscErrorCode  TSSundialsSetLinearTolerance(TS ts,double tol)
753 {
754   PetscErrorCode ierr;
755 
756   PetscFunctionBegin;
757   PetscValidLogicalCollectiveReal(ts,tol,2);
758   ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr);
759   PetscFunctionReturn(0);
760 }
761 
762 #undef __FUNCT__
763 #define __FUNCT__ "TSSundialsSetGramSchmidtType"
764 /*@
765    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
766         in GMRES method by SUNDIALS linear solver.
767 
768    Logically Collective on TS
769 
770    Input parameters:
771 +    ts  - the time-step context
772 -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
773 
774    Level: advanced
775 
776 .keywords: Sundials, orthogonalization
777 
778 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
779           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
780           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
781           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
782           TSSundialsSetExactFinalTime()
783 
784 @*/
785 PetscErrorCode  TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
786 {
787   PetscErrorCode ierr;
788 
789   PetscFunctionBegin;
790   ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr);
791   PetscFunctionReturn(0);
792 }
793 
794 #undef __FUNCT__
795 #define __FUNCT__ "TSSundialsSetTolerance"
796 /*@
797    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
798                          Sundials for error control.
799 
800    Logically Collective on TS
801 
802    Input parameters:
803 +    ts  - the time-step context
804 .    aabs - the absolute tolerance
805 -    rel - the relative tolerance
806 
807      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
808     these regulate the size of the error for a SINGLE timestep.
809 
810    Level: intermediate
811 
812 .keywords: Sundials, tolerance
813 
814 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
815           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
816           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
817           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
818           TSSundialsSetExactFinalTime()
819 
820 @*/
821 PetscErrorCode  TSSundialsSetTolerance(TS ts,double aabs,double rel)
822 {
823   PetscErrorCode ierr;
824 
825   PetscFunctionBegin;
826   ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr);
827   PetscFunctionReturn(0);
828 }
829 
830 #undef __FUNCT__
831 #define __FUNCT__ "TSSundialsGetPC"
832 /*@
833    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
834 
835    Input Parameter:
836 .    ts - the time-step context
837 
838    Output Parameter:
839 .    pc - the preconditioner context
840 
841    Level: advanced
842 
843 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
844           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
845           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
846           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
847 @*/
848 PetscErrorCode  TSSundialsGetPC(TS ts,PC *pc)
849 {
850   PetscErrorCode ierr;
851 
852   PetscFunctionBegin;
853   ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC *),(ts,pc));CHKERRQ(ierr);
854   PetscFunctionReturn(0);
855 }
856 
857 #undef __FUNCT__
858 #define __FUNCT__ "TSSundialsSetExactFinalTime"
859 /*@
860    TSSundialsSetExactFinalTime - Determines if Sundials interpolates solution to the
861       exact final time requested by the user or just returns it at the final time
862       it computed. (Defaults to true).
863 
864    Input Parameter:
865 +   ts - the time-step context
866 -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE
867 
868    Level: beginner
869 
870 .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
871           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
872           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
873           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
874 @*/
875 PetscErrorCode  TSSundialsSetExactFinalTime(TS ts,PetscBool  ft)
876 {
877   PetscErrorCode ierr;
878 
879   PetscFunctionBegin;
880   ierr = PetscTryMethod(ts,"TSSundialsSetExactFinalTime_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr);
881   PetscFunctionReturn(0);
882 }
883 
884 #undef __FUNCT__
885 #define __FUNCT__ "TSSundialsSetMinTimeStep"
886 /*@
887    TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller.
888 
889    Input Parameter:
890 +   ts - the time-step context
891 -   mindt - lowest time step if positive, negative to deactivate
892 
893    Note:
894    Sundials will error if it is not possible to keep the estimated truncation error below
895    the tolerance set with TSSundialsSetTolerance() without going below this step size.
896 
897    Level: beginner
898 
899 .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
900 @*/
901 PetscErrorCode  TSSundialsSetMinTimeStep(TS ts,PetscReal mindt)
902 {
903   PetscErrorCode ierr;
904 
905   PetscFunctionBegin;
906   ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr);
907   PetscFunctionReturn(0);
908 }
909 
910 #undef __FUNCT__
911 #define __FUNCT__ "TSSundialsSetMaxTimeStep"
912 /*@
913    TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller.
914 
915    Input Parameter:
916 +   ts - the time-step context
917 -   maxdt - lowest time step if positive, negative to deactivate
918 
919    Level: beginner
920 
921 .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
922 @*/
923 PetscErrorCode  TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt)
924 {
925   PetscErrorCode ierr;
926 
927   PetscFunctionBegin;
928   ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr);
929   PetscFunctionReturn(0);
930 }
931 
932 #undef __FUNCT__
933 #define __FUNCT__ "TSSundialsMonitorInternalSteps"
934 /*@
935    TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false).
936 
937    Input Parameter:
938 +   ts - the time-step context
939 -   ft - PETSC_TRUE if monitor, else PETSC_FALSE
940 
941    Level: beginner
942 
943 .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
944           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
945           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
946           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
947 @*/
948 PetscErrorCode  TSSundialsMonitorInternalSteps(TS ts,PetscBool  ft)
949 {
950   PetscErrorCode ierr;
951 
952   PetscFunctionBegin;
953   ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr);
954   PetscFunctionReturn(0);
955 }
956 /* -------------------------------------------------------------------------------------------*/
957 /*MC
958       TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
959 
960    Options Database:
961 +    -ts_sundials_type <bdf,adams>
962 .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
963 .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
964 .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
965 .    -ts_sundials_linear_tolerance <tol>
966 .    -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions
967 .    -ts_sundials_exact_final_time - Interpolate output to stop exactly at the final time
968 -    -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps
969 
970 
971     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
972            only PETSc PC options
973 
974     Level: beginner
975 
976 .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(),
977            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSundialsSetExactFinalTime()
978 
979 M*/
980 EXTERN_C_BEGIN
981 #undef __FUNCT__
982 #define __FUNCT__ "TSCreate_Sundials"
983 PetscErrorCode  TSCreate_Sundials(TS ts)
984 {
985   TS_Sundials *cvode;
986   PetscErrorCode ierr;
987 
988   PetscFunctionBegin;
989   ts->ops->reset          = TSReset_Sundials;
990   ts->ops->destroy        = TSDestroy_Sundials;
991   ts->ops->view           = TSView_Sundials;
992   ts->ops->setup          = TSSetUp_Sundials;
993   ts->ops->solve          = TSSolve_Sundials;
994   ts->ops->setfromoptions = TSSetFromOptions_Sundials;
995 
996   ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr);
997   ts->data          = (void*)cvode;
998   cvode->cvode_type = SUNDIALS_BDF;
999   cvode->gtype      = SUNDIALS_CLASSICAL_GS;
1000   cvode->restart    = 5;
1001   cvode->linear_tol = .05;
1002 
1003   cvode->exact_final_time = PETSC_TRUE;
1004   cvode->monitorstep      = PETSC_FALSE;
1005 
1006   ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr);
1007 
1008   cvode->mindt = -1.;
1009   cvode->maxdt = -1.;
1010 
1011   /* set tolerance for Sundials */
1012   cvode->reltol = 1e-6;
1013   cvode->abstol = 1e-6;
1014 
1015   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
1016                     TSSundialsSetType_Sundials);CHKERRQ(ierr);
1017   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C",
1018                     "TSSundialsSetGMRESRestart_Sundials",
1019                     TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr);
1020   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
1021                     "TSSundialsSetLinearTolerance_Sundials",
1022                      TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
1023   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
1024                     "TSSundialsSetGramSchmidtType_Sundials",
1025                      TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
1026   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
1027                     "TSSundialsSetTolerance_Sundials",
1028                      TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
1029   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C",
1030                     "TSSundialsSetMinTimeStep_Sundials",
1031                      TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr);
1032   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",
1033                     "TSSundialsSetMaxTimeStep_Sundials",
1034                      TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr);
1035   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
1036                     "TSSundialsGetPC_Sundials",
1037                      TSSundialsGetPC_Sundials);CHKERRQ(ierr);
1038   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
1039                     "TSSundialsGetIterations_Sundials",
1040                      TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
1041   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C",
1042                     "TSSundialsSetExactFinalTime_Sundials",
1043                      TSSundialsSetExactFinalTime_Sundials);CHKERRQ(ierr);
1044   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",
1045                     "TSSundialsMonitorInternalSteps_Sundials",
1046                      TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr);
1047 
1048   PetscFunctionReturn(0);
1049 }
1050 EXTERN_C_END
1051