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