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