xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision ee1cef2c4f92b4b45b8aa3ceac59cb944062fab9)
1 #define PETSCTS_DLL
2 
3 /*
4     Provides a PETSc interface to SUNDIALS/CVODE solver.
5     The interface to PVODE (old version of CVODE) was originally contributed
6     by Liyang Xu. It has been redone by Hong Zhang and Dinesh Kaushik.
7     Reference: sundials-2.3.0/examples/cvode/parallel/cvkryx_p.c
8 */
9 #include "sundials.h"  /*I "petscts.h" I*/
10 
11 /*
12       TSPrecond_Sundials - function that we provide to SUNDIALS to
13                         evaluate the preconditioner.
14 */
15 #undef __FUNCT__
16 #define __FUNCT__ "TSPrecond_Sundials"
17 PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy,
18                     booleantype jok,booleantype *jcurPtr,
19                     realtype _gamma,void *P_data,
20                     N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3)
21 {
22   TS             ts = (TS) P_data;
23   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
24   PC             pc = cvode->pc;
25   PetscErrorCode ierr;
26   Mat            Jac = ts->B;
27   Vec            yy = cvode->w1;
28   PetscScalar    one = 1.0,gm;
29   MatStructure   str = DIFFERENT_NONZERO_PATTERN;
30   PetscScalar    *y_data;
31 
32   PetscFunctionBegin;
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 /*
131     TSStep_Sundials_Nonlinear -
132 
133    steps - number of time steps
134    time - time that integrater is  terminated.
135 */
136 PetscErrorCode TSStep_Sundials_Nonlinear(TS ts,int *steps,double *time)
137 {
138   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
139   Vec            sol = ts->vec_sol;
140   PetscErrorCode ierr;
141   PetscInt       i,max_steps = ts->max_steps,flag;
142   long int       its;
143   realtype       t,tout;
144   PetscScalar    *y_data;
145   void           *mem;
146   const PCType   pctype;
147   PetscTruth     pcnone;
148 
149   PetscFunctionBegin;
150   /* Call CVodeCreate to create the solver memory */
151   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
152   if (!mem) SETERRQ(1,"CVodeCreate() fails");
153   flag = CVodeSetFdata(mem,ts);
154   if (flag) SETERRQ(1,"CVodeSetFdata() fails");
155 
156   /*
157      Call CVodeMalloc to initialize the integrator memory:
158      mem is the pointer to the integrator memory returned by CVodeCreate
159      f       is the user's right hand side function in y'=f(t,y)
160      T0      is the initial time
161      u       is the initial dependent variable vector
162      CV_SS   specifies scalar relative and absolute tolerances
163      reltol  is the relative tolerance
164      &abstol is a pointer to the scalar absolute tolerance
165   */
166   flag = CVodeMalloc(mem,TSFunction_Sundials,ts->ptime,cvode->y,CV_SS,cvode->reltol,&cvode->abstol);
167   if (flag) SETERRQ(1,"CVodeMalloc() fails");
168 
169   /* initialize the number of steps */
170   *steps = -ts->steps;
171   ierr   = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr);
172 
173   /* call CVSpgmr to use GMRES as the linear solver.        */
174   /* setup the ode integrator with the given preconditioner */
175   ierr = PCGetType(cvode->pc,&pctype);CHKERRQ(ierr);
176   ierr = PetscTypeCompare((PetscObject)cvode->pc,PCNONE,&pcnone);CHKERRQ(ierr);
177   if (pcnone){
178     flag  = CVSpgmr(mem,PREC_NONE,0);
179   } else {
180     flag  = CVSpgmr(mem,PREC_LEFT,0);
181   }
182   if (flag) SETERRQ(1,"CVSpgmr() fails");
183 
184   /* Set preconditioner setup and solve routines Precond and PSolve,
185      and the pointer to the user-defined block data */
186   flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials,ts);
187   if (flag) SETERRQ(1,"CVSpgmrSetPreconditioner() fails");
188 
189   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
190   if (flag) SETERRQ(1,"CVSpgmrSetGSType() fails");
191 
192   tout = ts->max_time;
193   ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr);
194   N_VSetArrayPointer((realtype *)y_data,cvode->y);
195   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
196   for (i = 0; i < max_steps; i++) {
197     if (ts->ptime >= ts->max_time) break;
198     ierr = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);CHKERRQ(ierr);
199     if (t > ts->max_time && cvode->exact_final_time) {
200       /* interpolate to final requested time */
201       ierr = CVodeGetDky(mem,tout,0,cvode->y);CHKERRQ(ierr);
202       t = tout;
203     }
204     ts->time_step = t - ts->ptime;
205     ts->ptime     = t;
206 
207     /* copy the solution from cvode->y to cvode->update and sol */
208     ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr);
209     ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr);
210     ierr = VecResetArray(cvode->w1); CHKERRQ(ierr);
211     ierr = VecCopy(cvode->update,sol);CHKERRQ(ierr);
212     ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr);
213     ts->nonlinear_its = its;
214     ierr = CVSpilsGetNumLinIters(mem, &its);
215     ts->linear_its = its;
216     ts->steps++;
217     ierr = TSMonitor(ts,ts->steps,t,sol);CHKERRQ(ierr);
218   }
219   cvode->mem = mem;
220   *steps    += ts->steps;
221   *time      = t;
222   PetscFunctionReturn(0);
223 }
224 
225 #undef __FUNCT__
226 #define __FUNCT__ "TSDestroy_Sundials"
227 PetscErrorCode TSDestroy_Sundials(TS ts)
228 {
229   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
230   PetscErrorCode ierr;
231 
232   PetscFunctionBegin;
233   if (cvode->pmat)   {ierr = MatDestroy(cvode->pmat);CHKERRQ(ierr);}
234   if (cvode->pc)     {ierr = PCDestroy(cvode->pc);CHKERRQ(ierr);}
235   if (cvode->update) {ierr = VecDestroy(cvode->update);CHKERRQ(ierr);}
236   if (cvode->func)   {ierr = VecDestroy(cvode->func);CHKERRQ(ierr);}
237   if (cvode->rhs)    {ierr = VecDestroy(cvode->rhs);CHKERRQ(ierr);}
238   if (cvode->w1)     {ierr = VecDestroy(cvode->w1);CHKERRQ(ierr);}
239   if (cvode->w2)     {ierr = VecDestroy(cvode->w2);CHKERRQ(ierr);}
240   ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr);
241   if (cvode->mem) {CVodeFree(&cvode->mem);}
242   ierr = PetscFree(cvode);CHKERRQ(ierr);
243   PetscFunctionReturn(0);
244 }
245 
246 #undef __FUNCT__
247 #define __FUNCT__ "TSSetUp_Sundials_Nonlinear"
248 PetscErrorCode TSSetUp_Sundials_Nonlinear(TS ts)
249 {
250   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
251   PetscErrorCode ierr;
252   int            glosize,locsize,i;
253   PetscScalar    *y_data,*parray;
254 
255   PetscFunctionBegin;
256   ierr = PCSetFromOptions(cvode->pc);CHKERRQ(ierr);
257   /* get the vector size */
258   ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr);
259   ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr);
260 
261   /* allocate the memory for N_Vec y */
262   cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
263   if (!cvode->y) SETERRQ(1,"cvode->y is not allocated");
264 
265   /* initialize N_Vec y */
266   ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr);
267   y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y);
268   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
269   /*ierr = PetscMemcpy(y_data,parray,locsize*sizeof(PETSC_SCALAR)); CHKERRQ(ierr);*/
270   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
271   ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr);
272   ierr = VecDuplicate(ts->vec_sol,&cvode->func);CHKERRQ(ierr);
273   ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr);
274   ierr = PetscLogObjectParent(ts,cvode->func);CHKERRQ(ierr);
275 
276   /*
277       Create work vectors for the TSPSolve_Sundials() routine. Note these are
278     allocated with zero space arrays because the actual array space is provided
279     by Sundials and set using VecPlaceArray().
280   */
281   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr);
282   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr);
283   ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr);
284   ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr);
285   PetscFunctionReturn(0);
286 }
287 
288 /* type of CVODE linear multistep method */
289 const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0};
290 /* type of G-S orthogonalization used by CVODE linear solver */
291 const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0};
292 
293 #undef __FUNCT__
294 #define __FUNCT__ "TSSetFromOptions_Sundials_Nonlinear"
295 PetscErrorCode TSSetFromOptions_Sundials_Nonlinear(TS ts)
296 {
297   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
298   PetscErrorCode ierr;
299   int            indx;
300   PetscTruth     flag;
301 
302   PetscFunctionBegin;
303   ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr);
304     ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr);
305     if (flag) {
306       ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr);
307     }
308     ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr);
309     if (flag) {
310       ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr);
311     }
312     ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr);
313     ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr);
314     ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr);
315     ierr = PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);CHKERRQ(ierr);
316     ierr = PetscOptionsTruth("-ts_sundials_exact_final_time","Allow SUNDIALS to stop near the final time, not exactly on it","TSSundialsSetExactFinalTime",cvode->exact_final_time,&cvode->exact_final_time,PETSC_NULL);CHKERRQ(ierr);
317   ierr = PetscOptionsTail();CHKERRQ(ierr);
318   PetscFunctionReturn(0);
319 }
320 
321 #undef __FUNCT__
322 #define __FUNCT__ "TSView_Sundials"
323 PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
324 {
325   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
326   PetscErrorCode ierr;
327   char           *type;
328   char           atype[] = "Adams";
329   char           btype[] = "BDF: backward differentiation formula";
330   PetscTruth     iascii,isstring;
331   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
332   PetscInt       qlast,qcur;
333   PetscReal      hinused,hlast,hcur,tcur,tolsfac;
334 
335   PetscFunctionBegin;
336   if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;}
337   else                                     {type = btype;}
338 
339   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
340   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
341   if (iascii) {
342     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr);
343     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr);
344     ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr);
345     ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr);
346     ierr = PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);CHKERRQ(ierr);
347     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
348       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
349     } else {
350       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
351     }
352 
353     /* Outputs from CVODE, CVSPILS */
354     ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr);
355     ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr);
356     ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
357                                    &nlinsetups,&nfails,&qlast,&qcur,
358                                    &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr);
359     ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr);
360     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr);
361     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr);
362     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr);
363 
364     ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr);
365     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr);
366     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr);
367 
368     ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */
369     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr);
370     ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr);
371     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr);
372     ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr);
373     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr);
374     ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr);
375     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr);
376     ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr);
377     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr);
378     ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr);
379     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr);
380   } else if (isstring) {
381     ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr);
382   } else {
383     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by TS Sundials",((PetscObject)viewer)->type_name);
384   }
385   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
386   ierr = PCView(cvode->pc,viewer);CHKERRQ(ierr);
387   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
388   PetscFunctionReturn(0);
389 }
390 
391 
392 /* --------------------------------------------------------------------------*/
393 EXTERN_C_BEGIN
394 #undef __FUNCT__
395 #define __FUNCT__ "TSSundialsSetType_Sundials"
396 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
397 {
398   TS_Sundials *cvode = (TS_Sundials*)ts->data;
399 
400   PetscFunctionBegin;
401   cvode->cvode_type = type;
402   PetscFunctionReturn(0);
403 }
404 EXTERN_C_END
405 
406 EXTERN_C_BEGIN
407 #undef __FUNCT__
408 #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials"
409 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart_Sundials(TS ts,int restart)
410 {
411   TS_Sundials *cvode = (TS_Sundials*)ts->data;
412 
413   PetscFunctionBegin;
414   cvode->restart = restart;
415   PetscFunctionReturn(0);
416 }
417 EXTERN_C_END
418 
419 EXTERN_C_BEGIN
420 #undef __FUNCT__
421 #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials"
422 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
423 {
424   TS_Sundials *cvode = (TS_Sundials*)ts->data;
425 
426   PetscFunctionBegin;
427   cvode->linear_tol = tol;
428   PetscFunctionReturn(0);
429 }
430 EXTERN_C_END
431 
432 EXTERN_C_BEGIN
433 #undef __FUNCT__
434 #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials"
435 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
436 {
437   TS_Sundials *cvode = (TS_Sundials*)ts->data;
438 
439   PetscFunctionBegin;
440   cvode->gtype = type;
441   PetscFunctionReturn(0);
442 }
443 EXTERN_C_END
444 
445 EXTERN_C_BEGIN
446 #undef __FUNCT__
447 #define __FUNCT__ "TSSundialsSetTolerance_Sundials"
448 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
449 {
450   TS_Sundials *cvode = (TS_Sundials*)ts->data;
451 
452   PetscFunctionBegin;
453   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
454   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
455   PetscFunctionReturn(0);
456 }
457 EXTERN_C_END
458 
459 EXTERN_C_BEGIN
460 #undef __FUNCT__
461 #define __FUNCT__ "TSSundialsGetPC_Sundials"
462 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC_Sundials(TS ts,PC *pc)
463 {
464   TS_Sundials *cvode = (TS_Sundials*)ts->data;
465 
466   PetscFunctionBegin;
467   *pc = cvode->pc;
468   PetscFunctionReturn(0);
469 }
470 EXTERN_C_END
471 
472 EXTERN_C_BEGIN
473 #undef __FUNCT__
474 #define __FUNCT__ "TSSundialsGetIterations_Sundials"
475 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
476 {
477   PetscFunctionBegin;
478   if (nonlin) *nonlin = ts->nonlinear_its;
479   if (lin)    *lin    = ts->linear_its;
480   PetscFunctionReturn(0);
481 }
482 EXTERN_C_END
483 
484 EXTERN_C_BEGIN
485 #undef __FUNCT__
486 #define __FUNCT__ "TSSundialsSetExactFinalTime_Sundials"
487 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime_Sundials(TS ts,PetscTruth s)
488 {
489   TS_Sundials *cvode = (TS_Sundials*)ts->data;
490 
491   PetscFunctionBegin;
492   cvode->exact_final_time = s;
493   PetscFunctionReturn(0);
494 }
495 EXTERN_C_END
496 /* -------------------------------------------------------------------------------------------*/
497 
498 #undef __FUNCT__
499 #define __FUNCT__ "TSSundialsGetIterations"
500 /*@C
501    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
502 
503    Not Collective
504 
505    Input parameters:
506 .    ts     - the time-step context
507 
508    Output Parameters:
509 +   nonlin - number of nonlinear iterations
510 -   lin    - number of linear iterations
511 
512    Level: advanced
513 
514    Notes:
515     These return the number since the creation of the TS object
516 
517 .keywords: non-linear iterations, linear iterations
518 
519 .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(),
520           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
521           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
522           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
523           TSSundialsSetExactFinalTime()
524 
525 @*/
526 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
527 {
528   PetscErrorCode ierr,(*f)(TS,int*,int*);
529 
530   PetscFunctionBegin;
531   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetIterations_C",(void (**)(void))&f);CHKERRQ(ierr);
532   if (f) {
533     ierr = (*f)(ts,nonlin,lin);CHKERRQ(ierr);
534   }
535   PetscFunctionReturn(0);
536 }
537 
538 #undef __FUNCT__
539 #define __FUNCT__ "TSSundialsSetType"
540 /*@
541    TSSundialsSetType - Sets the method that Sundials will use for integration.
542 
543    Collective on TS
544 
545    Input parameters:
546 +    ts     - the time-step context
547 -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
548 
549    Level: intermediate
550 
551 .keywords: Adams, backward differentiation formula
552 
553 .seealso: TSSundialsGetIterations(),  TSSundialsSetGMRESRestart(),
554           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
555           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
556           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
557           TSSundialsSetExactFinalTime()
558 @*/
559 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType(TS ts,TSSundialsLmmType type)
560 {
561   PetscErrorCode ierr,(*f)(TS,TSSundialsLmmType);
562 
563   PetscFunctionBegin;
564   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
565   if (f) {
566     ierr = (*f)(ts,type);CHKERRQ(ierr);
567   }
568   PetscFunctionReturn(0);
569 }
570 
571 #undef __FUNCT__
572 #define __FUNCT__ "TSSundialsSetGMRESRestart"
573 /*@
574    TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by
575        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
576        this is ALSO the maximum number of GMRES steps that will be used.
577 
578    Collective on TS
579 
580    Input parameters:
581 +    ts      - the time-step context
582 -    restart - number of direction vectors (the restart size).
583 
584    Level: advanced
585 
586 .keywords: GMRES, restart
587 
588 .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
589           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
590           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
591           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
592           TSSundialsSetExactFinalTime()
593 
594 @*/
595 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart(TS ts,int restart)
596 {
597   PetscErrorCode ierr,(*f)(TS,int);
598 
599   PetscFunctionBegin;
600   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGMRESRestart_C",(void (**)(void))&f);CHKERRQ(ierr);
601   if (f) {
602     ierr = (*f)(ts,restart);CHKERRQ(ierr);
603   }
604   PetscFunctionReturn(0);
605 }
606 
607 #undef __FUNCT__
608 #define __FUNCT__ "TSSundialsSetLinearTolerance"
609 /*@
610    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
611        system by SUNDIALS.
612 
613    Collective on TS
614 
615    Input parameters:
616 +    ts     - the time-step context
617 -    tol    - the factor by which the tolerance on the nonlinear solver is
618              multiplied to get the tolerance on the linear solver, .05 by default.
619 
620    Level: advanced
621 
622 .keywords: GMRES, linear convergence tolerance, SUNDIALS
623 
624 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
625           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
626           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
627           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
628           TSSundialsSetExactFinalTime()
629 
630 @*/
631 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance(TS ts,double tol)
632 {
633   PetscErrorCode ierr,(*f)(TS,double);
634 
635   PetscFunctionBegin;
636   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
637   if (f) {
638     ierr = (*f)(ts,tol);CHKERRQ(ierr);
639   }
640   PetscFunctionReturn(0);
641 }
642 
643 #undef __FUNCT__
644 #define __FUNCT__ "TSSundialsSetGramSchmidtType"
645 /*@
646    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
647         in GMRES method by SUNDIALS linear solver.
648 
649    Collective on TS
650 
651    Input parameters:
652 +    ts  - the time-step context
653 -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
654 
655    Level: advanced
656 
657 .keywords: Sundials, orthogonalization
658 
659 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
660           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
661           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
662           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
663           TSSundialsSetExactFinalTime()
664 
665 @*/
666 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
667 {
668   PetscErrorCode ierr,(*f)(TS,TSSundialsGramSchmidtType);
669 
670   PetscFunctionBegin;
671   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",(void (**)(void))&f);CHKERRQ(ierr);
672   if (f) {
673     ierr = (*f)(ts,type);CHKERRQ(ierr);
674   }
675   PetscFunctionReturn(0);
676 }
677 
678 #undef __FUNCT__
679 #define __FUNCT__ "TSSundialsSetTolerance"
680 /*@
681    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
682                          Sundials for error control.
683 
684    Collective on TS
685 
686    Input parameters:
687 +    ts  - the time-step context
688 .    aabs - the absolute tolerance
689 -    rel - the relative tolerance
690 
691      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
692     these regulate the size of the error for a SINGLE timestep.
693 
694    Level: intermediate
695 
696 .keywords: Sundials, tolerance
697 
698 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
699           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
700           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
701           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
702           TSSundialsSetExactFinalTime()
703 
704 @*/
705 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance(TS ts,double aabs,double rel)
706 {
707   PetscErrorCode ierr,(*f)(TS,double,double);
708 
709   PetscFunctionBegin;
710   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
711   if (f) {
712     ierr = (*f)(ts,aabs,rel);CHKERRQ(ierr);
713   }
714   PetscFunctionReturn(0);
715 }
716 
717 #undef __FUNCT__
718 #define __FUNCT__ "TSSundialsGetPC"
719 /*@
720    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
721 
722    Input Parameter:
723 .    ts - the time-step context
724 
725    Output Parameter:
726 .    pc - the preconditioner context
727 
728    Level: advanced
729 
730 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
731           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
732           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
733           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
734 @*/
735 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC(TS ts,PC *pc)
736 {
737   PetscErrorCode ierr,(*f)(TS,PC *);
738 
739   PetscFunctionBegin;
740   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetPC_C",(void (**)(void))&f);CHKERRQ(ierr);
741   if (f) {
742     ierr = (*f)(ts,pc);CHKERRQ(ierr);
743   } else {
744     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"TS must be of Sundials type to extract the PC");
745   }
746 
747   PetscFunctionReturn(0);
748 }
749 
750 #undef __FUNCT__
751 #define __FUNCT__ "TSSundialsSetExactFinalTime"
752 /*@
753    TSSundialsSetExactFinalTime - Determines if Sundials interpolates solution to the
754       exact final time requested by the user or just returns it at the final time
755       it computed. (Defaults to true).
756 
757    Input Parameter:
758 +   ts - the time-step context
759 -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE
760 
761    Level: beginner
762 
763 .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
764           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
765           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
766           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
767 @*/
768 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime(TS ts,PetscTruth ft)
769 {
770   PetscErrorCode ierr,(*f)(TS,PetscTruth);
771 
772   PetscFunctionBegin;
773   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetExactFinalTime_C",(void (**)(void))&f);CHKERRQ(ierr);
774   if (f) {
775     ierr = (*f)(ts,ft);CHKERRQ(ierr);
776   }
777   PetscFunctionReturn(0);
778 }
779 
780 /* -------------------------------------------------------------------------------------------*/
781 /*MC
782       TS_Sundials - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
783 
784    Options Database:
785 +    -ts_sundials_type <bdf,adams>
786 .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
787 .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
788 .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
789 .    -ts_sundials_linear_tolerance <tol>
790 .    -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions
791 -    -ts_sundials_not_exact_final_time -Allow SUNDIALS to stop near the final time, not exactly on it
792 
793     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
794            only PETSc PC options
795 
796     Level: beginner
797 
798 .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(),
799            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSundialsSetExactFinalTime()
800 
801 M*/
802 EXTERN_C_BEGIN
803 #undef __FUNCT__
804 #define __FUNCT__ "TSCreate_Sundials"
805 PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Sundials(TS ts)
806 {
807   TS_Sundials *cvode;
808   PetscErrorCode ierr;
809 
810   PetscFunctionBegin;
811   ts->ops->destroy         = TSDestroy_Sundials;
812   ts->ops->view            = TSView_Sundials;
813 
814   if (ts->problem_type != TS_NONLINEAR) {
815     SETERRQ(PETSC_ERR_SUP,"Only support for nonlinear problems");
816   }
817   ts->ops->setup           = TSSetUp_Sundials_Nonlinear;
818   ts->ops->step            = TSStep_Sundials_Nonlinear;
819   ts->ops->setfromoptions  = TSSetFromOptions_Sundials_Nonlinear;
820 
821   ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr);
822   ierr = PCCreate(((PetscObject)ts)->comm,&cvode->pc);CHKERRQ(ierr);
823   ierr = PetscLogObjectParent(ts,cvode->pc);CHKERRQ(ierr);
824   ts->data          = (void*)cvode;
825   cvode->cvode_type = SUNDIALS_BDF;
826   cvode->gtype      = SUNDIALS_CLASSICAL_GS;
827   cvode->restart    = 5;
828   cvode->linear_tol = .05;
829 
830   cvode->exact_final_time = PETSC_FALSE;
831 
832   ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr);
833   /* set tolerance for Sundials */
834   cvode->reltol = 1e-6;
835   cvode->abstol = 1e-6;
836 
837   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
838                     TSSundialsSetType_Sundials);CHKERRQ(ierr);
839   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C",
840                     "TSSundialsSetGMRESRestart_Sundials",
841                     TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr);
842   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
843                     "TSSundialsSetLinearTolerance_Sundials",
844                      TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
845   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
846                     "TSSundialsSetGramSchmidtType_Sundials",
847                      TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
848   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
849                     "TSSundialsSetTolerance_Sundials",
850                      TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
851   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
852                     "TSSundialsGetPC_Sundials",
853                      TSSundialsGetPC_Sundials);CHKERRQ(ierr);
854   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
855                     "TSSundialsGetIterations_Sundials",
856                      TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
857   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C",
858                     "TSSundialsSetExactFinalTime_Sundials",
859                      TSSundialsSetExactFinalTime_Sundials);CHKERRQ(ierr);
860   PetscFunctionReturn(0);
861 }
862 EXTERN_C_END
863 
864 
865 
866 
867 
868 
869 
870 
871 
872 
873