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