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