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