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