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