xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision f9426fe092dba0ba2fdf65dfec8d938c4b10a31c)
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,booleantype jok,booleantype *jcurPtr,
17                                   realtype _gamma,void *P_data,N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3)
18 {
19   TS             ts     = (TS) P_data;
20   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
21   PC             pc;
22   PetscErrorCode ierr;
23   Mat            J,P;
24   Vec            yy  = cvode->w1,yydot = cvode->ydot;
25   PetscReal      gm  = (PetscReal)_gamma;
26   MatStructure   str = DIFFERENT_NONZERO_PATTERN;
27   PetscScalar    *y_data;
28 
29   PetscFunctionBegin;
30   ierr   = TSGetIJacobian(ts,&J,&P,NULL,NULL);CHKERRQ(ierr);
31   y_data = (PetscScalar*) N_VGetArrayPointer(y);
32   ierr   = VecPlaceArray(yy,y_data);CHKERRQ(ierr);
33   ierr   = VecZeroEntries(yydot);CHKERRQ(ierr); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */
34   /* compute the shifted Jacobian   (1/gm)*I + Jrest */
35   ierr     = TSComputeIJacobian(ts,ts->ptime,yy,yydot,1/gm,&J,&P,&str,PETSC_FALSE);CHKERRQ(ierr);
36   ierr     = VecResetArray(yy);CHKERRQ(ierr);
37   ierr     = MatScale(P,gm);CHKERRQ(ierr); /* turn into I-gm*Jrest, J is not used by Sundials  */
38   *jcurPtr = TRUE;
39   ierr     = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr);
40   ierr     = PCSetOperators(pc,J,P,str);CHKERRQ(ierr);
41   PetscFunctionReturn(0);
42 }
43 
44 /*
45      TSPSolve_Sundials -  routine that we provide to Sundials that applies the preconditioner.
46 */
47 #undef __FUNCT__
48 #define __FUNCT__ "TSPSolve_Sundials"
49 PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z,
50                                  realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp)
51 {
52   TS             ts     = (TS) P_data;
53   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
54   PC             pc;
55   Vec            rr = cvode->w1,zz = cvode->w2;
56   PetscErrorCode ierr;
57   PetscScalar    *r_data,*z_data;
58 
59   PetscFunctionBegin;
60   /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
61   r_data = (PetscScalar*) N_VGetArrayPointer(r);
62   z_data = (PetscScalar*) N_VGetArrayPointer(z);
63   ierr   = VecPlaceArray(rr,r_data);CHKERRQ(ierr);
64   ierr   = VecPlaceArray(zz,z_data);CHKERRQ(ierr);
65 
66   /* Solve the Px=r and put the result in zz */
67   ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr);
68   ierr = PCApply(pc,rr,zz);CHKERRQ(ierr);
69   ierr = VecResetArray(rr);CHKERRQ(ierr);
70   ierr = VecResetArray(zz);CHKERRQ(ierr);
71   PetscFunctionReturn(0);
72 }
73 
74 /*
75         TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
76 */
77 #undef __FUNCT__
78 #define __FUNCT__ "TSFunction_Sundials"
79 int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx)
80 {
81   TS             ts = (TS) ctx;
82   DM             dm;
83   DMTS           tsdm;
84   TSIFunction    ifunction;
85   MPI_Comm       comm;
86   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
87   Vec            yy     = cvode->w1,yyd = cvode->w2,yydot = cvode->ydot;
88   PetscScalar    *y_data,*ydot_data;
89   PetscErrorCode ierr;
90 
91   PetscFunctionBegin;
92   ierr = PetscObjectGetComm((PetscObject)ts,&comm);CHKERRQ(ierr);
93   /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
94   y_data    = (PetscScalar*) N_VGetArrayPointer(y);
95   ydot_data = (PetscScalar*) N_VGetArrayPointer(ydot);
96   ierr      = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr);
97   ierr      = VecPlaceArray(yyd,ydot_data);CHKERRABORT(comm,ierr);
98 
99   /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */
100   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
101   ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr);
102   ierr = DMTSGetIFunction(dm,&ifunction,NULL);CHKERRQ(ierr);
103   if (!ifunction) {
104     ierr = TSComputeRHSFunction(ts,t,yy,yyd);CHKERRQ(ierr);
105   } else {                      /* If rhsfunction is also set, this computes both parts and shifts them to the right */
106     ierr = VecZeroEntries(yydot);CHKERRQ(ierr);
107     ierr = TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE);CHKERRABORT(comm,ierr);
108     ierr = VecScale(yyd,-1.);CHKERRQ(ierr);
109   }
110   ierr = VecResetArray(yy);CHKERRABORT(comm,ierr);
111   ierr = VecResetArray(yyd);CHKERRABORT(comm,ierr);
112   PetscFunctionReturn(0);
113 }
114 
115 /*
116        TSStep_Sundials - Calls Sundials to integrate the ODE.
117 */
118 #undef __FUNCT__
119 #define __FUNCT__ "TSStep_Sundials"
120 PetscErrorCode TSStep_Sundials(TS ts)
121 {
122   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
123   PetscErrorCode ierr;
124   PetscInt       flag;
125   long int       its,nsteps;
126   realtype       t,tout;
127   PetscScalar    *y_data;
128   void           *mem;
129 
130   PetscFunctionBegin;
131   mem  = cvode->mem;
132   tout = ts->max_time;
133   ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr);
134   N_VSetArrayPointer((realtype*)y_data,cvode->y);
135   ierr = VecRestoreArray(ts->vec_sol,NULL);CHKERRQ(ierr);
136 
137   ierr = TSPreStep(ts);CHKERRQ(ierr);
138 
139   /* We would like to call TSPreStep() when starting each step (including rejections) and TSPreStage() before each
140    * stage solve, but CVode does not appear to support this. */
141   if (cvode->monitorstep) flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
142   else flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL);
143 
144   if (flag) { /* display error message */
145     switch (flag) {
146       case CV_ILL_INPUT:
147         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT");
148         break;
149       case CV_TOO_CLOSE:
150         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE");
151         break;
152       case CV_TOO_MUCH_WORK: {
153         PetscReal      tcur;
154         ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
155         ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr);
156         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);
157       } break;
158       case CV_TOO_MUCH_ACC:
159         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC");
160         break;
161       case CV_ERR_FAILURE:
162         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE");
163         break;
164       case CV_CONV_FAILURE:
165         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE");
166         break;
167       case CV_LINIT_FAIL:
168         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL");
169         break;
170       case CV_LSETUP_FAIL:
171         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL");
172         break;
173       case CV_LSOLVE_FAIL:
174         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL");
175         break;
176       case CV_RHSFUNC_FAIL:
177         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL");
178         break;
179       case CV_FIRST_RHSFUNC_ERR:
180         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR");
181         break;
182       case CV_REPTD_RHSFUNC_ERR:
183         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR");
184         break;
185       case CV_UNREC_RHSFUNC_ERR:
186         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR");
187         break;
188       case CV_RTFUNC_FAIL:
189         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL");
190         break;
191       default:
192         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag);
193     }
194   }
195 
196   /* copy the solution from cvode->y to cvode->update and sol */
197   ierr = VecPlaceArray(cvode->w1,y_data);CHKERRQ(ierr);
198   ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr);
199   ierr = VecResetArray(cvode->w1);CHKERRQ(ierr);
200   ierr = VecCopy(cvode->update,ts->vec_sol);CHKERRQ(ierr);
201   ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr);
202   ierr = CVSpilsGetNumLinIters(mem, &its);
203   ts->snes_its = its; ts->ksp_its = its;
204 
205   ts->time_step = t - ts->ptime;
206   ts->ptime     = t;
207   ts->steps++;
208 
209   ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
210   if (!cvode->monitorstep) ts->steps = nsteps;
211   PetscFunctionReturn(0);
212 }
213 
214 #undef __FUNCT__
215 #define __FUNCT__ "TSInterpolate_Sundials"
216 static PetscErrorCode TSInterpolate_Sundials(TS ts,PetscReal t,Vec X)
217 {
218   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
219   N_Vector       y;
220   PetscErrorCode ierr;
221   PetscScalar    *x_data;
222   PetscInt       glosize,locsize;
223 
224   PetscFunctionBegin;
225   /* get the vector size */
226   ierr = VecGetSize(X,&glosize);CHKERRQ(ierr);
227   ierr = VecGetLocalSize(X,&locsize);CHKERRQ(ierr);
228 
229   /* allocate the memory for N_Vec y */
230   y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
231   if (!y) SETERRQ(PETSC_COMM_SELF,1,"Interpolated y is not allocated");
232 
233   ierr = VecGetArray(X,&x_data);CHKERRQ(ierr);
234   N_VSetArrayPointer((realtype*)x_data,y);
235   ierr = CVodeGetDky(cvode->mem,t,0,y);CHKERRQ(ierr);
236   ierr = VecRestoreArray(X,&x_data);CHKERRQ(ierr);
237   PetscFunctionReturn(0);
238 }
239 
240 #undef __FUNCT__
241 #define __FUNCT__ "TSReset_Sundials"
242 PetscErrorCode TSReset_Sundials(TS ts)
243 {
244   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
245   PetscErrorCode ierr;
246 
247   PetscFunctionBegin;
248   ierr = VecDestroy(&cvode->update);CHKERRQ(ierr);
249   ierr = VecDestroy(&cvode->ydot);CHKERRQ(ierr);
250   ierr = VecDestroy(&cvode->w1);CHKERRQ(ierr);
251   ierr = VecDestroy(&cvode->w2);CHKERRQ(ierr);
252   if (cvode->mem) CVodeFree(&cvode->mem);
253   PetscFunctionReturn(0);
254 }
255 
256 #undef __FUNCT__
257 #define __FUNCT__ "TSDestroy_Sundials"
258 PetscErrorCode TSDestroy_Sundials(TS ts)
259 {
260   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
261   PetscErrorCode ierr;
262 
263   PetscFunctionBegin;
264   ierr = TSReset_Sundials(ts);CHKERRQ(ierr);
265   ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr);
266   ierr = PetscFree(ts->data);CHKERRQ(ierr);
267   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",NULL);CHKERRQ(ierr);
268   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",NULL);CHKERRQ(ierr);
269   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",NULL);CHKERRQ(ierr);
270   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",NULL);CHKERRQ(ierr);
271   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",NULL);CHKERRQ(ierr);
272   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",NULL);CHKERRQ(ierr);
273   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",NULL);CHKERRQ(ierr);
274   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",NULL);CHKERRQ(ierr);
275   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",NULL);CHKERRQ(ierr);
276   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",NULL);CHKERRQ(ierr);
277   PetscFunctionReturn(0);
278 }
279 
280 #undef __FUNCT__
281 #define __FUNCT__ "TSSetUp_Sundials"
282 PetscErrorCode TSSetUp_Sundials(TS ts)
283 {
284   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
285   PetscErrorCode ierr;
286   PetscInt       glosize,locsize,i,flag;
287   PetscScalar    *y_data,*parray;
288   void           *mem;
289   PC             pc;
290   PCType         pctype;
291   PetscBool      pcnone;
292 
293   PetscFunctionBegin;
294   /* get the vector size */
295   ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr);
296   ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr);
297 
298   /* allocate the memory for N_Vec y */
299   cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
300   if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated");
301 
302   /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
303   ierr   = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr);
304   y_data = (PetscScalar*) N_VGetArrayPointer(cvode->y);
305   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
306   ierr = VecRestoreArray(ts->vec_sol,NULL);CHKERRQ(ierr);
307 
308   ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr);
309   ierr = VecDuplicate(ts->vec_sol,&cvode->ydot);CHKERRQ(ierr);
310   ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->update);CHKERRQ(ierr);
311   ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->ydot);CHKERRQ(ierr);
312 
313   /*
314     Create work vectors for the TSPSolve_Sundials() routine. Note these are
315     allocated with zero space arrays because the actual array space is provided
316     by Sundials and set using VecPlaceArray().
317   */
318   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr);
319   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr);
320   ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w1);CHKERRQ(ierr);
321   ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w2);CHKERRQ(ierr);
322 
323   /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
324   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
325   if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails");
326   cvode->mem = mem;
327 
328   /* Set the pointer to user-defined data */
329   flag = CVodeSetUserData(mem, ts);
330   if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails");
331 
332   /* Sundials may choose to use a smaller initial step, but will never use a larger step. */
333   flag = CVodeSetInitStep(mem,(realtype)ts->time_step);
334   if (flag) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetInitStep() failed");
335   if (cvode->mindt > 0) {
336     flag = CVodeSetMinStep(mem,(realtype)cvode->mindt);
337     if (flag) {
338       if (flag == CV_MEM_NULL) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL");
339       else if (flag == CV_ILL_INPUT) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size");
340       else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed");
341     }
342   }
343   if (cvode->maxdt > 0) {
344     flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt);
345     if (flag) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMaxStep() failed");
346   }
347 
348   /* Call CVodeInit to initialize the integrator memory and specify the
349    * user's right hand side function in u'=f(t,u), the inital time T0, and
350    * the initial dependent variable vector cvode->y */
351   flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y);
352   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag);
353 
354   /* specifies scalar relative and absolute tolerances */
355   flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol);
356   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag);
357 
358   /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */
359   flag = CVodeSetMaxNumSteps(mem,ts->max_steps);
360 
361   /* call CVSpgmr to use GMRES as the linear solver.        */
362   /* setup the ode integrator with the given preconditioner */
363   ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr);
364   ierr = PCGetType(pc,&pctype);CHKERRQ(ierr);
365   ierr = PetscObjectTypeCompare((PetscObject)pc,PCNONE,&pcnone);CHKERRQ(ierr);
366   if (pcnone) {
367     flag = CVSpgmr(mem,PREC_NONE,0);
368     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
369   } else {
370     flag = CVSpgmr(mem,PREC_LEFT,cvode->maxl);
371     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
372 
373     /* Set preconditioner and solve routines Precond and PSolve,
374      and the pointer to the user-defined block data */
375     flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials);
376     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag);
377   }
378 
379   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
380   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag);
381   PetscFunctionReturn(0);
382 }
383 
384 /* type of CVODE linear multistep method */
385 const char *const TSSundialsLmmTypes[] = {"","ADAMS","BDF","TSSundialsLmmType","SUNDIALS_",0};
386 /* type of G-S orthogonalization used by CVODE linear solver */
387 const char *const TSSundialsGramSchmidtTypes[] = {"","MODIFIED","CLASSICAL","TSSundialsGramSchmidtType","SUNDIALS_",0};
388 
389 #undef __FUNCT__
390 #define __FUNCT__ "TSSetFromOptions_Sundials"
391 PetscErrorCode TSSetFromOptions_Sundials(TS ts)
392 {
393   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
394   PetscErrorCode ierr;
395   int            indx;
396   PetscBool      flag;
397   PC             pc;
398 
399   PetscFunctionBegin;
400   ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr);
401   ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr);
402   if (flag) {
403     ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr);
404   }
405   ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr);
406   if (flag) {
407     ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr);
408   }
409   ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,NULL);CHKERRQ(ierr);
410   ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,NULL);CHKERRQ(ierr);
411   ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,NULL);CHKERRQ(ierr);
412   ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,NULL);CHKERRQ(ierr);
413   ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr);
414   ierr = PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,&flag);CHKERRQ(ierr);
415   ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,NULL);CHKERRQ(ierr);
416   ierr = PetscOptionsTail();CHKERRQ(ierr);
417   ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr);
418   ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
419   PetscFunctionReturn(0);
420 }
421 
422 #undef __FUNCT__
423 #define __FUNCT__ "TSView_Sundials"
424 PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
425 {
426   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
427   PetscErrorCode ierr;
428   char           *type;
429   char           atype[] = "Adams";
430   char           btype[] = "BDF: backward differentiation formula";
431   PetscBool      iascii,isstring;
432   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
433   PetscInt       qlast,qcur;
434   PetscReal      hinused,hlast,hcur,tcur,tolsfac;
435   PC             pc;
436 
437   PetscFunctionBegin;
438   if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype;
439   else                                     type = btype;
440 
441   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
442   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
443   if (iascii) {
444     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr);
445     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr);
446     ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr);
447     ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr);
448     ierr = PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %D\n",cvode->maxl);CHKERRQ(ierr);
449     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
450       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
451     } else {
452       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
453     }
454     if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);}
455     if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);}
456 
457     /* Outputs from CVODE, CVSPILS */
458     ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr);
459     ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr);
460     ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
461                                    &nlinsetups,&nfails,&qlast,&qcur,
462                                    &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr);
463     ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr);
464     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr);
465     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr);
466     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr);
467 
468     ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr);
469     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr);
470     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr);
471 
472     ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */
473     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr);
474     ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr);
475     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr);
476 
477     ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr);
478     ierr = PCView(pc,viewer);CHKERRQ(ierr);
479     ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr);
480     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr);
481     ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr);
482     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr);
483 
484     ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr);
485     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr);
486     ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr);
487     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr);
488   } else if (isstring) {
489     ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr);
490   }
491   PetscFunctionReturn(0);
492 }
493 
494 
495 /* --------------------------------------------------------------------------*/
496 #undef __FUNCT__
497 #define __FUNCT__ "TSSundialsSetType_Sundials"
498 PetscErrorCode  TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
499 {
500   TS_Sundials *cvode = (TS_Sundials*)ts->data;
501 
502   PetscFunctionBegin;
503   cvode->cvode_type = type;
504   PetscFunctionReturn(0);
505 }
506 
507 #undef __FUNCT__
508 #define __FUNCT__ "TSSundialsSetMaxl_Sundials"
509 PetscErrorCode  TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl)
510 {
511   TS_Sundials *cvode = (TS_Sundials*)ts->data;
512 
513   PetscFunctionBegin;
514   cvode->maxl = maxl;
515   PetscFunctionReturn(0);
516 }
517 
518 #undef __FUNCT__
519 #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials"
520 PetscErrorCode  TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
521 {
522   TS_Sundials *cvode = (TS_Sundials*)ts->data;
523 
524   PetscFunctionBegin;
525   cvode->linear_tol = tol;
526   PetscFunctionReturn(0);
527 }
528 
529 #undef __FUNCT__
530 #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials"
531 PetscErrorCode  TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
532 {
533   TS_Sundials *cvode = (TS_Sundials*)ts->data;
534 
535   PetscFunctionBegin;
536   cvode->gtype = type;
537   PetscFunctionReturn(0);
538 }
539 
540 #undef __FUNCT__
541 #define __FUNCT__ "TSSundialsSetTolerance_Sundials"
542 PetscErrorCode  TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
543 {
544   TS_Sundials *cvode = (TS_Sundials*)ts->data;
545 
546   PetscFunctionBegin;
547   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
548   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
549   PetscFunctionReturn(0);
550 }
551 
552 #undef __FUNCT__
553 #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials"
554 PetscErrorCode  TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)
555 {
556   TS_Sundials *cvode = (TS_Sundials*)ts->data;
557 
558   PetscFunctionBegin;
559   cvode->mindt = mindt;
560   PetscFunctionReturn(0);
561 }
562 
563 #undef __FUNCT__
564 #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials"
565 PetscErrorCode  TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)
566 {
567   TS_Sundials *cvode = (TS_Sundials*)ts->data;
568 
569   PetscFunctionBegin;
570   cvode->maxdt = maxdt;
571   PetscFunctionReturn(0);
572 }
573 #undef __FUNCT__
574 #define __FUNCT__ "TSSundialsGetPC_Sundials"
575 PetscErrorCode  TSSundialsGetPC_Sundials(TS ts,PC *pc)
576 {
577   SNES           snes;
578   KSP            ksp;
579   PetscErrorCode ierr;
580 
581   PetscFunctionBegin;
582   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
583   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
584   ierr = KSPGetPC(ksp,pc);CHKERRQ(ierr);
585   PetscFunctionReturn(0);
586 }
587 
588 #undef __FUNCT__
589 #define __FUNCT__ "TSSundialsGetIterations_Sundials"
590 PetscErrorCode  TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
591 {
592   PetscFunctionBegin;
593   if (nonlin) *nonlin = ts->snes_its;
594   if (lin)    *lin    = ts->ksp_its;
595   PetscFunctionReturn(0);
596 }
597 
598 #undef __FUNCT__
599 #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials"
600 PetscErrorCode  TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s)
601 {
602   TS_Sundials *cvode = (TS_Sundials*)ts->data;
603 
604   PetscFunctionBegin;
605   cvode->monitorstep = s;
606   PetscFunctionReturn(0);
607 }
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(), TSSundialsSetMaxl(),
632           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
633           TSSundialsGetIterations(), TSSundialsSetType(),
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(),  TSSundialsSetMaxl(),
662           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
663           TSSundialsGetIterations(), TSSundialsSetType(),
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__ "TSSundialsSetMaxl"
678 /*@
679    TSSundialsSetMaxl - 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 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 -    maxl - number of direction vectors (the dimension of Krylov subspace).
688 
689    Level: advanced
690 
691 .keywords: GMRES
692 
693 .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
694           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
695           TSSundialsGetIterations(), TSSundialsSetType(),
696           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
697           TSSetExactFinalTime()
698 
699 @*/
700 PetscErrorCode  TSSundialsSetMaxl(TS ts,PetscInt maxl)
701 {
702   PetscErrorCode ierr;
703 
704   PetscFunctionBegin;
705   PetscValidLogicalCollectiveInt(ts,maxl,2);
706   ierr = PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));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(), TSSundialsSetMaxl(),
728           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
729           TSSundialsGetIterations(), TSSundialsSetType(),
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(), TSSundialsSetMaxl(),
761           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
762           TSSundialsGetIterations(), TSSundialsSetType(),
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(), TSSundialsSetGMRESMaxl(),
797           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
798           TSSundialsGetIterations(), TSSundialsSetType(),
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(), TSSundialsSetMaxl(),
826           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
827           TSSundialsGetIterations(), TSSundialsSetType(),
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(), TSSundialsSetMaxl(),
899           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
900           TSSundialsGetIterations(), TSSundialsSetType(),
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_maxl <maxl> - Max dimension of the Krylov subspace
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(), TSSundialsSetMaxl(), TSSundialsSetLinearTolerance(),
931            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime()
932 
933 M*/
934 #undef __FUNCT__
935 #define __FUNCT__ "TSCreate_Sundials"
936 PETSC_EXTERN PetscErrorCode TSCreate_Sundials(TS ts)
937 {
938   TS_Sundials    *cvode;
939   PetscErrorCode ierr;
940   PC             pc;
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->step           = TSStep_Sundials;
948   ts->ops->interpolate    = TSInterpolate_Sundials;
949   ts->ops->setfromoptions = TSSetFromOptions_Sundials;
950 
951   ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr);
952 
953   ts->data           = (void*)cvode;
954   cvode->cvode_type  = SUNDIALS_BDF;
955   cvode->gtype       = SUNDIALS_CLASSICAL_GS;
956   cvode->maxl        = 5;
957   cvode->linear_tol  = .05;
958   cvode->monitorstep = PETSC_TRUE;
959 
960   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)ts),&(cvode->comm_sundials));CHKERRQ(ierr);
961 
962   cvode->mindt = -1.;
963   cvode->maxdt = -1.;
964 
965   /* set tolerance for Sundials */
966   cvode->reltol = 1e-6;
967   cvode->abstol = 1e-6;
968 
969   /* set PCNONE as default pctype */
970   ierr = TSSundialsGetPC_Sundials(ts,&pc);CHKERRQ(ierr);
971   ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
972 
973   if (ts->exact_final_time != TS_EXACTFINALTIME_STEPOVER) ts->exact_final_time = TS_EXACTFINALTIME_STEPOVER;
974 
975   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",TSSundialsSetType_Sundials);CHKERRQ(ierr);
976   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",TSSundialsSetMaxl_Sundials);CHKERRQ(ierr);
977   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
978   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
979   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
980   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr);
981   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr);
982   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",TSSundialsGetPC_Sundials);CHKERRQ(ierr);
983   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
984   ierr = PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr);
985   PetscFunctionReturn(0);
986 }
987