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