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