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