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