xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision 156232b76c11958d2c0daf032f3a45a36a39c32d)
1 /*
2     Provides a PETSc interface to SUNDIALS/CVODE solver.
3     The interface to PVODE (old version of CVODE) was originally contributed
4     by Liyang Xu. It has been redone by Hong Zhang and Dinesh Kaushik.
5 
6     Reference: sundials-2.4.0/examples/cvode/parallel/cvDiurnal_kry_p.c
7 */
8 #include "sundials.h"  /*I "petscts.h" I*/
9 
10 /*
11       TSPrecond_Sundials - function that we provide to SUNDIALS to
12                         evaluate the preconditioner.
13 */
14 #undef __FUNCT__
15 #define __FUNCT__ "TSPrecond_Sundials"
16 PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy,
17                     booleantype jok,booleantype *jcurPtr,
18                     realtype _gamma,void *P_data,
19                     N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3)
20 {
21   TS             ts = (TS) P_data;
22   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
23   PC             pc;
24   PetscErrorCode ierr;
25   Mat            J,P;
26   Vec            yy = cvode->w1,yydot = cvode->ydot;
27   PetscReal      gm = (PetscReal)_gamma;
28   MatStructure   str = DIFFERENT_NONZERO_PATTERN;
29   PetscScalar    *y_data;
30 
31   PetscFunctionBegin;
32   ierr = TSGetIJacobian(ts,&J,&P,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
33   y_data = (PetscScalar *) N_VGetArrayPointer(y);
34   ierr = VecPlaceArray(yy,y_data); CHKERRQ(ierr);
35   ierr = VecZeroEntries(yydot);CHKERRQ(ierr); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */
36   /* compute the shifted Jacobian   (1/gm)*I + Jrest */
37   ierr = TSComputeIJacobian(ts,ts->ptime,yy,yydot,1/gm,&J,&P,&str,PETSC_FALSE);CHKERRQ(ierr);
38   ierr = VecResetArray(yy); CHKERRQ(ierr);
39   ierr = MatScale(P,gm);CHKERRQ(ierr);  /* turn into I-gm*Jrest, J is not used by Sundials  */
40   *jcurPtr = TRUE;
41   ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
42   ierr = PCSetOperators(pc,J,P,str);CHKERRQ(ierr);
43   PetscFunctionReturn(0);
44 }
45 
46 /*
47      TSPSolve_Sundials -  routine that we provide to Sundials that applies the preconditioner.
48 */
49 #undef __FUNCT__
50 #define __FUNCT__ "TSPSolve_Sundials"
51 PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z,
52                                  realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp)
53 {
54   TS              ts = (TS) P_data;
55   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
56   PC              pc;
57   Vec             rr = cvode->w1,zz = cvode->w2;
58   PetscErrorCode  ierr;
59   PetscScalar     *r_data,*z_data;
60 
61   PetscFunctionBegin;
62   /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
63   r_data  = (PetscScalar *) N_VGetArrayPointer(r);
64   z_data  = (PetscScalar *) N_VGetArrayPointer(z);
65   ierr = VecPlaceArray(rr,r_data); CHKERRQ(ierr);
66   ierr = VecPlaceArray(zz,z_data); CHKERRQ(ierr);
67 
68   /* Solve the Px=r and put the result in zz */
69   ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
70   ierr = PCApply(pc,rr,zz); CHKERRQ(ierr);
71   ierr = VecResetArray(rr); CHKERRQ(ierr);
72   ierr = VecResetArray(zz); CHKERRQ(ierr);
73   PetscFunctionReturn(0);
74 }
75 
76 /*
77         TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
78 */
79 #undef __FUNCT__
80 #define __FUNCT__ "TSFunction_Sundials"
81 int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx)
82 {
83   TS              ts = (TS) ctx;
84   MPI_Comm        comm = ((PetscObject)ts)->comm;
85   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
86   Vec             yy = cvode->w1,yyd = cvode->w2,yydot = cvode->ydot;
87   PetscScalar     *y_data,*ydot_data;
88   PetscErrorCode  ierr;
89 
90   PetscFunctionBegin;
91   /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
92   y_data     = (PetscScalar *) N_VGetArrayPointer(y);
93   ydot_data  = (PetscScalar *) N_VGetArrayPointer(ydot);
94   ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr);
95   ierr = VecPlaceArray(yyd,ydot_data); CHKERRABORT(comm,ierr);
96   ierr = VecZeroEntries(yydot);CHKERRQ(ierr);
97 
98   /* now compute the right hand side function */
99   ierr = TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE); CHKERRABORT(comm,ierr);
100   ierr = VecScale(yyd,-1.);CHKERRQ(ierr);
101   ierr = VecResetArray(yy); CHKERRABORT(comm,ierr);
102   ierr = VecResetArray(yyd); CHKERRABORT(comm,ierr);
103   PetscFunctionReturn(0);
104 }
105 
106 /*
107        TSStep_Sundials - Calls Sundials to integrate the ODE.
108 */
109 #undef __FUNCT__
110 #define __FUNCT__ "TSStep_Sundials"
111 PetscErrorCode TSStep_Sundials(TS ts)
112 {
113   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
114   Vec            sol = ts->vec_sol;
115   PetscErrorCode ierr;
116   PetscInt       flag;
117   long int       its,nsteps;
118   realtype       t,tout;
119   PetscScalar    *y_data;
120   void           *mem;
121   SNES           snes;
122   Vec            res; /* This, together with snes, will check if the SNES vec_func has been set */
123 
124   PetscFunctionBegin;
125   mem  = cvode->mem;
126   tout = ts->max_time;
127   ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr);
128   N_VSetArrayPointer((realtype *)y_data,cvode->y);
129   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
130 
131   /* Should think about moving this outside the loop */
132   ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr);
133   ierr = SNESGetFunction(snes, &res, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
134   if (!res) {
135     ierr = TSSetIFunction(ts, sol, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
136   }
137 
138   ierr = TSPreStep(ts);CHKERRQ(ierr);
139 
140   if (cvode->monitorstep) {
141     flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
142   } else {
143     flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL);
144   }
145 
146   if (flag){ /* display error message */
147     switch (flag){
148       case CV_ILL_INPUT:
149         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT");
150         break;
151       case CV_TOO_CLOSE:
152         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE");
153         break;
154       case CV_TOO_MUCH_WORK: {
155         PetscReal      tcur;
156         ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
157         ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr);
158         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);
159       } break;
160       case CV_TOO_MUCH_ACC:
161         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC");
162         break;
163       case CV_ERR_FAILURE:
164         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE");
165         break;
166       case CV_CONV_FAILURE:
167         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE");
168         break;
169       case CV_LINIT_FAIL:
170         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL");
171         break;
172       case CV_LSETUP_FAIL:
173         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL");
174         break;
175       case CV_LSOLVE_FAIL:
176         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL");
177         break;
178       case CV_RHSFUNC_FAIL:
179         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL");
180         break;
181       case CV_FIRST_RHSFUNC_ERR:
182         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR");
183         break;
184       case CV_REPTD_RHSFUNC_ERR:
185         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR");
186         break;
187       case CV_UNREC_RHSFUNC_ERR:
188         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR");
189         break;
190       case CV_RTFUNC_FAIL:
191         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL");
192         break;
193       default:
194         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag);
195     }
196   }
197 
198   /* copy the solution from cvode->y to cvode->update and sol */
199   ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr);
200   ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr);
201   ierr = VecResetArray(cvode->w1); CHKERRQ(ierr);
202   ierr = VecCopy(cvode->update,sol);CHKERRQ(ierr);
203   ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr);
204   ierr = CVSpilsGetNumLinIters(mem, &its);
205   ts->nonlinear_its = its; ts->linear_its = its;
206 
207   ts->time_step = t - ts->ptime;
208   ts->ptime     = t;
209   ts->steps++;
210 
211   ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
212   if (!cvode->monitorstep) ts->steps = nsteps;
213   PetscFunctionReturn(0);
214 }
215 
216 #undef __FUNCT__
217 #define __FUNCT__ "TSInterpolate_Sundials"
218 static PetscErrorCode TSInterpolate_Sundials(TS ts,PetscReal t,Vec X)
219 {
220   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
221   N_Vector        y;
222   PetscErrorCode  ierr;
223   PetscScalar     *x_data;
224   PetscInt        glosize,locsize;
225 
226   PetscFunctionBegin;
227 
228   /* get the vector size */
229   ierr = VecGetSize(X,&glosize);CHKERRQ(ierr);
230   ierr = VecGetLocalSize(X,&locsize);CHKERRQ(ierr);
231 
232   /* allocate the memory for N_Vec y */
233   y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
234   if (!y) SETERRQ(PETSC_COMM_SELF,1,"Interpolated y is not allocated");
235 
236   ierr = VecGetArray(X,&x_data);CHKERRQ(ierr);
237   N_VSetArrayPointer((realtype *)x_data,y);
238   ierr = CVodeGetDky(cvode->mem,t,0,y);CHKERRQ(ierr);
239   ierr = VecRestoreArray(X,&x_data);CHKERRQ(ierr);
240 
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   const 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,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr);
323   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,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->initial_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 = PetscTypeCompare((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 *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0};
398 /* type of G-S orthogonalization used by CVODE linear solver */
399 const char *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 
410   PetscFunctionBegin;
411   ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr);
412     ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr);
413     if (flag) {
414       ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr);
415     }
416     ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr);
417     if (flag) {
418       ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr);
419     }
420     ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr);
421     ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr);
422     ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,PETSC_NULL);CHKERRQ(ierr);
423     ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,PETSC_NULL);CHKERRQ(ierr);
424     ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr);
425     ierr = PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,&flag);CHKERRQ(ierr);
426     ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);CHKERRQ(ierr);
427   ierr = PetscOptionsTail();CHKERRQ(ierr);
428   PetscFunctionReturn(0);
429 }
430 
431 #undef __FUNCT__
432 #define __FUNCT__ "TSView_Sundials"
433 PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
434 {
435   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
436   PetscErrorCode ierr;
437   char           *type;
438   char           atype[] = "Adams";
439   char           btype[] = "BDF: backward differentiation formula";
440   PetscBool      iascii,isstring;
441   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
442   PetscInt       qlast,qcur;
443   PetscReal      hinused,hlast,hcur,tcur,tolsfac;
444   PC             pc;
445 
446   PetscFunctionBegin;
447   if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;}
448   else                                     {type = btype;}
449 
450   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
451   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
452   if (iascii) {
453     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr);
454     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr);
455     ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr);
456     ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr);
457     ierr = PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %D\n",cvode->maxl);CHKERRQ(ierr);
458     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
459       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
460     } else {
461       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
462     }
463     if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);}
464     if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);}
465 
466     /* Outputs from CVODE, CVSPILS */
467     ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr);
468     ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr);
469     ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
470                                    &nlinsetups,&nfails,&qlast,&qcur,
471                                    &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr);
472     ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr);
473     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr);
474     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr);
475     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr);
476 
477     ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr);
478     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr);
479     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr);
480 
481     ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */
482     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr);
483     ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr);
484     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr);
485 
486     ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
487     ierr = PCView(pc,viewer);CHKERRQ(ierr);
488     ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr);
489     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr);
490     ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr);
491     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr);
492 
493     ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr);
494     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr);
495     ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr);
496     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr);
497   } else if (isstring) {
498     ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr);
499   }
500   PetscFunctionReturn(0);
501 }
502 
503 
504 /* --------------------------------------------------------------------------*/
505 EXTERN_C_BEGIN
506 #undef __FUNCT__
507 #define __FUNCT__ "TSSundialsSetType_Sundials"
508 PetscErrorCode  TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
509 {
510   TS_Sundials *cvode = (TS_Sundials*)ts->data;
511 
512   PetscFunctionBegin;
513   cvode->cvode_type = type;
514   PetscFunctionReturn(0);
515 }
516 EXTERN_C_END
517 
518 EXTERN_C_BEGIN
519 #undef __FUNCT__
520 #define __FUNCT__ "TSSundialsSetMaxl_Sundials"
521 PetscErrorCode  TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl)
522 {
523   TS_Sundials *cvode = (TS_Sundials*)ts->data;
524 
525   PetscFunctionBegin;
526   cvode->maxl = maxl;
527   PetscFunctionReturn(0);
528 }
529 EXTERN_C_END
530 
531 EXTERN_C_BEGIN
532 #undef __FUNCT__
533 #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials"
534 PetscErrorCode  TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
535 {
536   TS_Sundials *cvode = (TS_Sundials*)ts->data;
537 
538   PetscFunctionBegin;
539   cvode->linear_tol = tol;
540   PetscFunctionReturn(0);
541 }
542 EXTERN_C_END
543 
544 EXTERN_C_BEGIN
545 #undef __FUNCT__
546 #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials"
547 PetscErrorCode  TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
548 {
549   TS_Sundials *cvode = (TS_Sundials*)ts->data;
550 
551   PetscFunctionBegin;
552   cvode->gtype = type;
553   PetscFunctionReturn(0);
554 }
555 EXTERN_C_END
556 
557 EXTERN_C_BEGIN
558 #undef __FUNCT__
559 #define __FUNCT__ "TSSundialsSetTolerance_Sundials"
560 PetscErrorCode  TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
561 {
562   TS_Sundials *cvode = (TS_Sundials*)ts->data;
563 
564   PetscFunctionBegin;
565   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
566   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
567   PetscFunctionReturn(0);
568 }
569 EXTERN_C_END
570 
571 EXTERN_C_BEGIN
572 #undef __FUNCT__
573 #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials"
574 PetscErrorCode  TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)
575 {
576   TS_Sundials *cvode = (TS_Sundials*)ts->data;
577 
578   PetscFunctionBegin;
579   cvode->mindt = mindt;
580   PetscFunctionReturn(0);
581 }
582 EXTERN_C_END
583 
584 EXTERN_C_BEGIN
585 #undef __FUNCT__
586 #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials"
587 PetscErrorCode  TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)
588 {
589   TS_Sundials *cvode = (TS_Sundials*)ts->data;
590 
591   PetscFunctionBegin;
592   cvode->maxdt = maxdt;
593   PetscFunctionReturn(0);
594 }
595 EXTERN_C_END
596 
597 EXTERN_C_BEGIN
598 #undef __FUNCT__
599 #define __FUNCT__ "TSSundialsGetPC_Sundials"
600 PetscErrorCode  TSSundialsGetPC_Sundials(TS ts,PC *pc)
601 {
602   SNES            snes;
603   KSP             ksp;
604   PetscErrorCode  ierr;
605 
606   PetscFunctionBegin;
607   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
608   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
609   ierr = KSPGetPC(ksp,pc); CHKERRQ(ierr);
610   PetscFunctionReturn(0);
611 }
612 EXTERN_C_END
613 
614 EXTERN_C_BEGIN
615 #undef __FUNCT__
616 #define __FUNCT__ "TSSundialsGetIterations_Sundials"
617 PetscErrorCode  TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
618 {
619   PetscFunctionBegin;
620   if (nonlin) *nonlin = ts->nonlinear_its;
621   if (lin)    *lin    = ts->linear_its;
622   PetscFunctionReturn(0);
623 }
624 EXTERN_C_END
625 
626 EXTERN_C_BEGIN
627 #undef __FUNCT__
628 #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials"
629 PetscErrorCode  TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool  s)
630 {
631   TS_Sundials *cvode = (TS_Sundials*)ts->data;
632 
633   PetscFunctionBegin;
634   cvode->monitorstep = s;
635   PetscFunctionReturn(0);
636 }
637 EXTERN_C_END
638 /* -------------------------------------------------------------------------------------------*/
639 
640 #undef __FUNCT__
641 #define __FUNCT__ "TSSundialsGetIterations"
642 /*@C
643    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
644 
645    Not Collective
646 
647    Input parameters:
648 .    ts     - the time-step context
649 
650    Output Parameters:
651 +   nonlin - number of nonlinear iterations
652 -   lin    - number of linear iterations
653 
654    Level: advanced
655 
656    Notes:
657     These return the number since the creation of the TS object
658 
659 .keywords: non-linear iterations, linear iterations
660 
661 .seealso: TSSundialsSetType(), TSSundialsSetMaxl(),
662           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
663           TSSundialsGetIterations(), TSSundialsSetType(),
664           TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime()
665 
666 @*/
667 PetscErrorCode  TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
668 {
669   PetscErrorCode ierr;
670 
671   PetscFunctionBegin;
672   ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr);
673   PetscFunctionReturn(0);
674 }
675 
676 #undef __FUNCT__
677 #define __FUNCT__ "TSSundialsSetType"
678 /*@
679    TSSundialsSetType - Sets the method that Sundials will use for integration.
680 
681    Logically Collective on TS
682 
683    Input parameters:
684 +    ts     - the time-step context
685 -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
686 
687    Level: intermediate
688 
689 .keywords: Adams, backward differentiation formula
690 
691 .seealso: TSSundialsGetIterations(),  TSSundialsSetMaxl(),
692           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
693           TSSundialsGetIterations(), TSSundialsSetType(),
694           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
695           TSSetExactFinalTime()
696 @*/
697 PetscErrorCode  TSSundialsSetType(TS ts,TSSundialsLmmType type)
698 {
699   PetscErrorCode ierr;
700 
701   PetscFunctionBegin;
702   ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr);
703   PetscFunctionReturn(0);
704 }
705 
706 #undef __FUNCT__
707 #define __FUNCT__ "TSSundialsSetMaxl"
708 /*@
709    TSSundialsSetMaxl - Sets the dimension of the Krylov space used by
710        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
711        this is the maximum number of GMRES steps that will be used.
712 
713    Logically Collective on TS
714 
715    Input parameters:
716 +    ts      - the time-step context
717 -    maxl - number of direction vectors (the dimension of Krylov subspace).
718 
719    Level: advanced
720 
721 .keywords: GMRES
722 
723 .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
724           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
725           TSSundialsGetIterations(), TSSundialsSetType(),
726           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
727           TSSetExactFinalTime()
728 
729 @*/
730 PetscErrorCode  TSSundialsSetMaxl(TS ts,PetscInt maxl)
731 {
732   PetscErrorCode ierr;
733 
734   PetscFunctionBegin;
735   PetscValidLogicalCollectiveInt(ts,maxl,2);
736   ierr = PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));CHKERRQ(ierr);
737   PetscFunctionReturn(0);
738 }
739 
740 #undef __FUNCT__
741 #define __FUNCT__ "TSSundialsSetLinearTolerance"
742 /*@
743    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
744        system by SUNDIALS.
745 
746    Logically Collective on TS
747 
748    Input parameters:
749 +    ts     - the time-step context
750 -    tol    - the factor by which the tolerance on the nonlinear solver is
751              multiplied to get the tolerance on the linear solver, .05 by default.
752 
753    Level: advanced
754 
755 .keywords: GMRES, linear convergence tolerance, SUNDIALS
756 
757 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
758           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
759           TSSundialsGetIterations(), TSSundialsSetType(),
760           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
761           TSSetExactFinalTime()
762 
763 @*/
764 PetscErrorCode  TSSundialsSetLinearTolerance(TS ts,double tol)
765 {
766   PetscErrorCode ierr;
767 
768   PetscFunctionBegin;
769   PetscValidLogicalCollectiveReal(ts,tol,2);
770   ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr);
771   PetscFunctionReturn(0);
772 }
773 
774 #undef __FUNCT__
775 #define __FUNCT__ "TSSundialsSetGramSchmidtType"
776 /*@
777    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
778         in GMRES method by SUNDIALS linear solver.
779 
780    Logically Collective on TS
781 
782    Input parameters:
783 +    ts  - the time-step context
784 -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
785 
786    Level: advanced
787 
788 .keywords: Sundials, orthogonalization
789 
790 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
791           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
792           TSSundialsGetIterations(), TSSundialsSetType(),
793           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
794           TSSetExactFinalTime()
795 
796 @*/
797 PetscErrorCode  TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
798 {
799   PetscErrorCode ierr;
800 
801   PetscFunctionBegin;
802   ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr);
803   PetscFunctionReturn(0);
804 }
805 
806 #undef __FUNCT__
807 #define __FUNCT__ "TSSundialsSetTolerance"
808 /*@
809    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
810                          Sundials for error control.
811 
812    Logically Collective on TS
813 
814    Input parameters:
815 +    ts  - the time-step context
816 .    aabs - the absolute tolerance
817 -    rel - the relative tolerance
818 
819      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
820     these regulate the size of the error for a SINGLE timestep.
821 
822    Level: intermediate
823 
824 .keywords: Sundials, tolerance
825 
826 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESMaxl(),
827           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
828           TSSundialsGetIterations(), TSSundialsSetType(),
829           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
830           TSSetExactFinalTime()
831 
832 @*/
833 PetscErrorCode  TSSundialsSetTolerance(TS ts,double aabs,double rel)
834 {
835   PetscErrorCode ierr;
836 
837   PetscFunctionBegin;
838   ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr);
839   PetscFunctionReturn(0);
840 }
841 
842 #undef __FUNCT__
843 #define __FUNCT__ "TSSundialsGetPC"
844 /*@
845    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
846 
847    Input Parameter:
848 .    ts - the time-step context
849 
850    Output Parameter:
851 .    pc - the preconditioner context
852 
853    Level: advanced
854 
855 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
856           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
857           TSSundialsGetIterations(), TSSundialsSetType(),
858           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
859 @*/
860 PetscErrorCode  TSSundialsGetPC(TS ts,PC *pc)
861 {
862   PetscErrorCode ierr;
863 
864   PetscFunctionBegin;
865   ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC *),(ts,pc));CHKERRQ(ierr);
866   PetscFunctionReturn(0);
867 }
868 
869 #undef __FUNCT__
870 #define __FUNCT__ "TSSundialsSetMinTimeStep"
871 /*@
872    TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller.
873 
874    Input Parameter:
875 +   ts - the time-step context
876 -   mindt - lowest time step if positive, negative to deactivate
877 
878    Note:
879    Sundials will error if it is not possible to keep the estimated truncation error below
880    the tolerance set with TSSundialsSetTolerance() without going below this step size.
881 
882    Level: beginner
883 
884 .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
885 @*/
886 PetscErrorCode  TSSundialsSetMinTimeStep(TS ts,PetscReal mindt)
887 {
888   PetscErrorCode ierr;
889 
890   PetscFunctionBegin;
891   ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr);
892   PetscFunctionReturn(0);
893 }
894 
895 #undef __FUNCT__
896 #define __FUNCT__ "TSSundialsSetMaxTimeStep"
897 /*@
898    TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller.
899 
900    Input Parameter:
901 +   ts - the time-step context
902 -   maxdt - lowest time step if positive, negative to deactivate
903 
904    Level: beginner
905 
906 .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
907 @*/
908 PetscErrorCode  TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt)
909 {
910   PetscErrorCode ierr;
911 
912   PetscFunctionBegin;
913   ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr);
914   PetscFunctionReturn(0);
915 }
916 
917 #undef __FUNCT__
918 #define __FUNCT__ "TSSundialsMonitorInternalSteps"
919 /*@
920    TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false).
921 
922    Input Parameter:
923 +   ts - the time-step context
924 -   ft - PETSC_TRUE if monitor, else PETSC_FALSE
925 
926    Level: beginner
927 
928 .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
929           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
930           TSSundialsGetIterations(), TSSundialsSetType(),
931           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
932 @*/
933 PetscErrorCode  TSSundialsMonitorInternalSteps(TS ts,PetscBool  ft)
934 {
935   PetscErrorCode ierr;
936 
937   PetscFunctionBegin;
938   ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr);
939   PetscFunctionReturn(0);
940 }
941 /* -------------------------------------------------------------------------------------------*/
942 /*MC
943       TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
944 
945    Options Database:
946 +    -ts_sundials_type <bdf,adams>
947 .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
948 .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
949 .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
950 .    -ts_sundials_linear_tolerance <tol>
951 .    -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace
952 -    -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps
953 
954 
955     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
956            only PETSc PC options
957 
958     Level: beginner
959 
960 .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetMaxl(), TSSundialsSetLinearTolerance(),
961            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime()
962 
963 M*/
964 EXTERN_C_BEGIN
965 #undef __FUNCT__
966 #define __FUNCT__ "TSCreate_Sundials"
967 PetscErrorCode  TSCreate_Sundials(TS ts)
968 {
969   TS_Sundials    *cvode;
970   PetscErrorCode ierr;
971   PC             pc;
972 
973   PetscFunctionBegin;
974   ts->ops->reset          = TSReset_Sundials;
975   ts->ops->destroy        = TSDestroy_Sundials;
976   ts->ops->view           = TSView_Sundials;
977   ts->ops->setup          = TSSetUp_Sundials;
978   ts->ops->step           = TSStep_Sundials;
979   ts->ops->interpolate    = TSInterpolate_Sundials;
980   ts->ops->setfromoptions = TSSetFromOptions_Sundials;
981 
982   ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr);
983   ts->data                = (void*)cvode;
984   cvode->cvode_type       = SUNDIALS_BDF;
985   cvode->gtype            = SUNDIALS_CLASSICAL_GS;
986   cvode->maxl             = 5;
987   cvode->linear_tol       = .05;
988 
989   cvode->monitorstep      = PETSC_TRUE;
990 
991   ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr);
992 
993   cvode->mindt = -1.;
994   cvode->maxdt = -1.;
995 
996   /* set tolerance for Sundials */
997   cvode->reltol = 1e-6;
998   cvode->abstol = 1e-6;
999 
1000   /* set PCNONE as default pctype */
1001   ierr = TSSundialsGetPC_Sundials(ts,&pc);CHKERRQ(ierr);
1002   ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1003 
1004   if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_FALSE;
1005 
1006   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
1007                     TSSundialsSetType_Sundials);CHKERRQ(ierr);
1008   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxl_C",
1009                     "TSSundialsSetMaxl_Sundials",
1010                     TSSundialsSetMaxl_Sundials);CHKERRQ(ierr);
1011   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
1012                     "TSSundialsSetLinearTolerance_Sundials",
1013                      TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
1014   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
1015                     "TSSundialsSetGramSchmidtType_Sundials",
1016                      TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
1017   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
1018                     "TSSundialsSetTolerance_Sundials",
1019                      TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
1020   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C",
1021                     "TSSundialsSetMinTimeStep_Sundials",
1022                      TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr);
1023   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",
1024                     "TSSundialsSetMaxTimeStep_Sundials",
1025                      TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr);
1026   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
1027                     "TSSundialsGetPC_Sundials",
1028                      TSSundialsGetPC_Sundials);CHKERRQ(ierr);
1029   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
1030                     "TSSundialsGetIterations_Sundials",
1031                      TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
1032   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",
1033                     "TSSundialsMonitorInternalSteps_Sundials",
1034                      TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr);
1035 
1036   PetscFunctionReturn(0);
1037 }
1038 EXTERN_C_END
1039