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