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