xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision a04cf4d80e220366494dd5711ec4e8408d1a483f)
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     }
54     else {
55       ierr = MatCopy(Jac,cvode->pmat,str);CHKERRQ(ierr);
56     }
57     *jcurPtr = TRUE;
58   }
59 
60   /* construct I-gamma*Jac  */
61   gm   = -_gamma;
62   ierr = MatScale(Jac,gm);CHKERRQ(ierr);
63   ierr = MatShift(Jac,one);CHKERRQ(ierr);
64 
65   ierr = PCSetOperators(pc,Jac,Jac,str);CHKERRQ(ierr);
66   PetscFunctionReturn(0);
67 }
68 
69 /*
70      TSPSolve_Sundials -  routine that we provide to Sundials that applies the preconditioner.
71 */
72 #undef __FUNCT__
73 #define __FUNCT__ "TSPSolve_Sundials"
74 PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z,
75                                  realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp)
76 {
77   TS              ts = (TS) P_data;
78   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
79   PC              pc = cvode->pc;
80   Vec             rr = cvode->w1,zz = cvode->w2;
81   PetscErrorCode  ierr;
82   PetscScalar     *r_data,*z_data;
83 
84   PetscFunctionBegin;
85   /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
86   r_data  = (PetscScalar *) N_VGetArrayPointer(r);
87   z_data  = (PetscScalar *) N_VGetArrayPointer(z);
88   ierr = VecPlaceArray(rr,r_data); CHKERRQ(ierr);
89   ierr = VecPlaceArray(zz,z_data); CHKERRQ(ierr);
90 
91   /* Solve the Px=r and put the result in zz */
92   ierr = PCApply(pc,rr,zz); CHKERRQ(ierr);
93   ierr = VecResetArray(rr); CHKERRQ(ierr);
94   ierr = VecResetArray(zz); CHKERRQ(ierr);
95   cvode->linear_solves++;
96   PetscFunctionReturn(0);
97 }
98 
99 /*
100         TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
101 */
102 #undef __FUNCT__
103 #define __FUNCT__ "TSFunction_Sundials"
104 int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx)
105 {
106   TS              ts = (TS) ctx;
107   MPI_Comm        comm = ts->comm;
108   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
109   Vec             yy = cvode->w1,yyd = cvode->w2;
110   PetscScalar     *y_data,*ydot_data;
111   PetscErrorCode  ierr;
112 
113   PetscFunctionBegin;
114   /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
115   y_data     = (PetscScalar *) N_VGetArrayPointer(y);
116   ydot_data  = (PetscScalar *) N_VGetArrayPointer(ydot);
117   ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr);
118   ierr = VecPlaceArray(yyd,ydot_data); CHKERRABORT(comm,ierr);
119 
120   /* now compute the right hand side function */
121   ierr = TSComputeRHSFunction(ts,t,yy,yyd); CHKERRABORT(comm,ierr);
122   ierr = VecResetArray(yy); CHKERRABORT(comm,ierr);
123   ierr = VecResetArray(yyd); CHKERRABORT(comm,ierr);
124   PetscFunctionReturn(0);
125 }
126 
127 /*
128        TSStep_Sundials_Nonlinear - Calls Sundials to integrate the ODE.
129 */
130 #undef __FUNCT__
131 #define __FUNCT__ "TSStep_Sundials_Nonlinear"
132 /*
133     TSStep_Sundials_Nonlinear -
134 
135    steps - number of time steps
136    time - time that integrater is  terminated.
137 */
138 PetscErrorCode TSStep_Sundials_Nonlinear(TS ts,int *steps,double *time)
139 {
140   TS_Sundials  *cvode = (TS_Sundials*)ts->data;
141   Vec          sol = ts->vec_sol;
142   PetscErrorCode ierr;
143   int          i,max_steps = ts->max_steps,flag;
144   long int     its;
145   realtype     t,tout;
146   PetscScalar  *y_data;
147   void         *mem;
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   flag  = CVSpgmr(mem,PREC_LEFT,0);
176   if (flag) SETERRQ(1,"CVSpgmr() fails");
177 
178   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
179   if (flag) SETERRQ(1,"CVSpgmrSetGSType() fails");
180 
181   /* Set preconditioner setup and solve routines Precond and PSolve,
182      and the pointer to the user-defined block data */
183   flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials,ts);
184   if (flag) SETERRQ(1,"CVSpgmrSetPreconditioner() fails");
185 
186   tout = ts->max_time;
187   ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr);
188   N_VSetArrayPointer((realtype *)y_data,cvode->y);
189   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
190   for (i = 0; i < max_steps; i++) {
191     if (ts->ptime >= ts->max_time) break;
192     ierr = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);CHKERRQ(ierr);
193     ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr);
194     cvode->nonlinear_solves += its;
195 
196     if (t > ts->max_time && cvode->exact_final_time) {
197       /* interpolate to final requested time */
198       ierr = CVodeGetDky(mem,tout,0,cvode->y);CHKERRQ(ierr);
199       t = tout;
200     }
201     ts->time_step = t - ts->ptime;
202     ts->ptime     = t;
203 
204     /* copy the solution from cvode->y to cvode->update and sol */
205     ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr);
206     ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr);
207     ierr = VecResetArray(cvode->w1); CHKERRQ(ierr);
208     ierr = VecCopy(cvode->update,sol);CHKERRQ(ierr);
209     ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr);
210     ts->nonlinear_its = its;
211     ierr = CVSpilsGetNumLinIters(mem, &its);
212     ts->linear_its = its;
213     ts->steps++;
214     ierr = TSMonitor(ts,ts->steps,t,sol);CHKERRQ(ierr);
215   }
216   CVodeFree(&mem);
217   *steps += ts->steps;
218   *time   = t;
219   PetscFunctionReturn(0);
220 }
221 
222 #undef __FUNCT__
223 #define __FUNCT__ "TSDestroy_Sundials"
224 PetscErrorCode TSDestroy_Sundials(TS ts)
225 {
226   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
227   PetscErrorCode ierr;
228 
229   PetscFunctionBegin;
230   if (cvode->pmat)   {ierr = MatDestroy(cvode->pmat);CHKERRQ(ierr);}
231   if (cvode->pc)     {ierr = PCDestroy(cvode->pc);CHKERRQ(ierr);}
232   if (cvode->update) {ierr = VecDestroy(cvode->update);CHKERRQ(ierr);}
233   if (cvode->func)   {ierr = VecDestroy(cvode->func);CHKERRQ(ierr);}
234   if (cvode->rhs)    {ierr = VecDestroy(cvode->rhs);CHKERRQ(ierr);}
235   if (cvode->w1)     {ierr = VecDestroy(cvode->w1);CHKERRQ(ierr);}
236   if (cvode->w2)     {ierr = VecDestroy(cvode->w2);CHKERRQ(ierr);}
237   ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr);
238   ierr = PetscFree(cvode);CHKERRQ(ierr);
239   PetscFunctionReturn(0);
240 }
241 
242 #undef __FUNCT__
243 #define __FUNCT__ "TSSetUp_Sundials_Nonlinear"
244 PetscErrorCode TSSetUp_Sundials_Nonlinear(TS ts)
245 {
246   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
247   PetscErrorCode ierr;
248   int            glosize,locsize,i;
249   PetscScalar    *y_data,*parray;
250 
251   PetscFunctionBegin;
252   ierr = PCSetFromOptions(cvode->pc);CHKERRQ(ierr);
253   /* get the vector size */
254   ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr);
255   ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr);
256 
257   /* allocate the memory for N_Vec y */
258   cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
259   if (!cvode->y) SETERRQ(1,"cvode->y is not allocated");
260 
261   /* initialize N_Vec y */
262   ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr);
263   y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y);
264   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
265   /*ierr = PetscMemcpy(y_data,parray,locsize*sizeof(PETSC_SCALAR)); CHKERRQ(ierr);*/
266   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
267   ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr);
268   ierr = VecDuplicate(ts->vec_sol,&cvode->func);CHKERRQ(ierr);
269   ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr);
270   ierr = PetscLogObjectParent(ts,cvode->func);CHKERRQ(ierr);
271 
272   /*
273       Create work vectors for the TSPSolve_Sundials() routine. Note these are
274     allocated with zero space arrays because the actual array space is provided
275     by Sundials and set using VecPlaceArray().
276   */
277   ierr = VecCreateMPIWithArray(ts->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr);
278   ierr = VecCreateMPIWithArray(ts->comm,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr);
279   ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr);
280   ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr);
281   PetscFunctionReturn(0);
282 }
283 
284 const char *TSSundialsTypes[] = {"Undefined","adams","bdf","TSSundialsType","SUNDIALS_",0};
285 const char *TSSundialsTypes[] = {"modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0};
286 
287 #undef __FUNCT__
288 #define __FUNCT__ "TSSetFromOptions_Sundials_Nonlinear"
289 PetscErrorCode TSSetFromOptions_Sundials_Nonlinear(TS ts)
290 {
291   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
292   PetscErrorCode ierr;
293   int            indx;
294   PetscTruth     flag;
295 
296   PetscFunctionBegin;
297   ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr);
298     ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsTypes,3,TSSundialsTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr);
299     if (flag) {
300       ierr = TSSundialsSetType(ts,(TSSundialsType)indx);CHKERRQ(ierr);
301     }
302     ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,2,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr);
303     if (flag) {
304       ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr);
305     }
306     ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr);
307     ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr);
308     ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr);
309     ierr = PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);CHKERRQ(ierr);
310     ierr = PetscOptionsName("-ts_sundials_exact_final_time","Allow SUNDIALS to stop near the final time, not exactly on it","TSSundialsSetExactFinalTime",&flag);CHKERRQ(ierr);
311     if (flag) cvode->exact_final_time = PETSC_TRUE;
312   ierr = PetscOptionsTail();CHKERRQ(ierr);
313   PetscFunctionReturn(0);
314 }
315 
316 #undef __FUNCT__
317 #define __FUNCT__ "TSView_Sundials"
318 PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
319 {
320   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
321   PetscErrorCode ierr;
322   char           *type;
323   char            atype[] = "Adams";
324   char            btype[] = "BDF: backward differentiation formula";
325   PetscTruth     iascii,isstring;
326 
327   PetscFunctionBegin;
328   if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;}
329   else                                     {type = btype;}
330 
331   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
332   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
333   if (iascii) {
334     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr);
335     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr);
336     ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr);
337     ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr);
338     ierr = PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);CHKERRQ(ierr);
339     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
340       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
341     } else {
342       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
343     }
344   } else if (isstring) {
345     ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr);
346   } else {
347     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by TS Sundials",((PetscObject)viewer)->type_name);
348   }
349   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
350   ierr = PCView(cvode->pc,viewer);CHKERRQ(ierr);
351   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
352   PetscFunctionReturn(0);
353 }
354 
355 
356 /* --------------------------------------------------------------------------*/
357 EXTERN_C_BEGIN
358 #undef __FUNCT__
359 #define __FUNCT__ "TSSundialsSetType_Sundials"
360 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType_Sundials(TS ts,TSSundialsType type)
361 {
362   TS_Sundials *cvode = (TS_Sundials*)ts->data;
363 
364   PetscFunctionBegin;
365   cvode->cvode_type = type;
366   PetscFunctionReturn(0);
367 }
368 EXTERN_C_END
369 
370 EXTERN_C_BEGIN
371 #undef __FUNCT__
372 #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials"
373 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart_Sundials(TS ts,int restart)
374 {
375   TS_Sundials *cvode = (TS_Sundials*)ts->data;
376 
377   PetscFunctionBegin;
378   cvode->restart = restart;
379   PetscFunctionReturn(0);
380 }
381 EXTERN_C_END
382 
383 EXTERN_C_BEGIN
384 #undef __FUNCT__
385 #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials"
386 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
387 {
388   TS_Sundials *cvode = (TS_Sundials*)ts->data;
389 
390   PetscFunctionBegin;
391   cvode->linear_tol = tol;
392   PetscFunctionReturn(0);
393 }
394 EXTERN_C_END
395 
396 EXTERN_C_BEGIN
397 #undef __FUNCT__
398 #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials"
399 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
400 {
401   TS_Sundials *cvode = (TS_Sundials*)ts->data;
402 
403   PetscFunctionBegin;
404   cvode->gtype = type;
405   PetscFunctionReturn(0);
406 }
407 EXTERN_C_END
408 
409 EXTERN_C_BEGIN
410 #undef __FUNCT__
411 #define __FUNCT__ "TSSundialsSetTolerance_Sundials"
412 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
413 {
414   TS_Sundials *cvode = (TS_Sundials*)ts->data;
415 
416   PetscFunctionBegin;
417   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
418   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
419   PetscFunctionReturn(0);
420 }
421 EXTERN_C_END
422 
423 EXTERN_C_BEGIN
424 #undef __FUNCT__
425 #define __FUNCT__ "TSSundialsGetPC_Sundials"
426 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC_Sundials(TS ts,PC *pc)
427 {
428   TS_Sundials *cvode = (TS_Sundials*)ts->data;
429 
430   PetscFunctionBegin;
431   *pc = cvode->pc;
432   PetscFunctionReturn(0);
433 }
434 EXTERN_C_END
435 
436 EXTERN_C_BEGIN
437 #undef __FUNCT__
438 #define __FUNCT__ "TSSundialsGetIterations_Sundials"
439 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
440 {
441   TS_Sundials *cvode = (TS_Sundials*)ts->data;
442 
443   PetscFunctionBegin;
444   if (nonlin) *nonlin = cvode->nonlinear_solves;
445   if (lin)    *lin    = cvode->linear_solves;
446   PetscFunctionReturn(0);
447 }
448 EXTERN_C_END
449 
450 EXTERN_C_BEGIN
451 #undef __FUNCT__
452 #define __FUNCT__ "TSSundialsSetExactFinalTime_Sundials"
453 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime_Sundials(TS ts,PetscTruth s)
454 {
455   TS_Sundials *cvode = (TS_Sundials*)ts->data;
456 
457   PetscFunctionBegin;
458   cvode->exact_final_time = s;
459   PetscFunctionReturn(0);
460 }
461 EXTERN_C_END
462 /* -------------------------------------------------------------------------------------------*/
463 
464 #undef __FUNCT__
465 #define __FUNCT__ "TSSundialsGetIterations"
466 /*@C
467    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
468 
469    Not Collective
470 
471    Input parameters:
472 .    ts     - the time-step context
473 
474    Output Parameters:
475 +   nonlin - number of nonlinear iterations
476 -   lin    - number of linear iterations
477 
478    Level: advanced
479 
480    Notes:
481     These return the number since the creation of the TS object
482 
483 .keywords: non-linear iterations, linear iterations
484 
485 .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(),
486           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
487           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
488           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
489           TSSundialsSetExactFinalTime()
490 
491 @*/
492 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
493 {
494   PetscErrorCode ierr,(*f)(TS,int*,int*);
495 
496   PetscFunctionBegin;
497   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetIterations_C",(void (**)(void))&f);CHKERRQ(ierr);
498   if (f) {
499     ierr = (*f)(ts,nonlin,lin);CHKERRQ(ierr);
500   }
501   PetscFunctionReturn(0);
502 }
503 
504 #undef __FUNCT__
505 #define __FUNCT__ "TSSundialsSetType"
506 /*@
507    TSSundialsSetType - Sets the method that Sundials will use for integration.
508 
509    Collective on TS
510 
511    Input parameters:
512 +    ts     - the time-step context
513 -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
514 
515    Level: intermediate
516 
517 .keywords: Adams, backward differentiation formula
518 
519 .seealso: TSSundialsGetIterations(),  TSSundialsSetGMRESRestart(),
520           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
521           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
522           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
523           TSSundialsSetExactFinalTime()
524 @*/
525 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType(TS ts,TSSundialsType type)
526 {
527   PetscErrorCode ierr,(*f)(TS,TSSundialsType);
528 
529   PetscFunctionBegin;
530   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
531   if (f) {
532     ierr = (*f)(ts,type);CHKERRQ(ierr);
533   }
534   PetscFunctionReturn(0);
535 }
536 
537 #undef __FUNCT__
538 #define __FUNCT__ "TSSundialsSetGMRESRestart"
539 /*@
540    TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by
541        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
542        this is ALSO the maximum number of GMRES steps that will be used.
543 
544    Collective on TS
545 
546    Input parameters:
547 +    ts      - the time-step context
548 -    restart - number of direction vectors (the restart size).
549 
550    Level: advanced
551 
552 .keywords: GMRES, restart
553 
554 .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
555           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
556           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
557           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
558           TSSundialsSetExactFinalTime()
559 
560 @*/
561 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart(TS ts,int restart)
562 {
563   PetscErrorCode ierr,(*f)(TS,int);
564 
565   PetscFunctionBegin;
566   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGMRESRestart_C",(void (**)(void))&f);CHKERRQ(ierr);
567   if (f) {
568     ierr = (*f)(ts,restart);CHKERRQ(ierr);
569   }
570   PetscFunctionReturn(0);
571 }
572 
573 #undef __FUNCT__
574 #define __FUNCT__ "TSSundialsSetLinearTolerance"
575 /*@
576    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
577        system by SUNDIALS.
578 
579    Collective on TS
580 
581    Input parameters:
582 +    ts     - the time-step context
583 -    tol    - the factor by which the tolerance on the nonlinear solver is
584              multiplied to get the tolerance on the linear solver, .05 by default.
585 
586    Level: advanced
587 
588 .keywords: GMRES, linear convergence tolerance, SUNDIALS
589 
590 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
591           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
592           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
593           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
594           TSSundialsSetExactFinalTime()
595 
596 @*/
597 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance(TS ts,double tol)
598 {
599   PetscErrorCode ierr,(*f)(TS,double);
600 
601   PetscFunctionBegin;
602   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
603   if (f) {
604     ierr = (*f)(ts,tol);CHKERRQ(ierr);
605   }
606   PetscFunctionReturn(0);
607 }
608 
609 #undef __FUNCT__
610 #define __FUNCT__ "TSSundialsSetGramSchmidtType"
611 /*@
612    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
613         in GMRES method by SUNDIALS linear solver.
614 
615    Collective on TS
616 
617    Input parameters:
618 +    ts  - the time-step context
619 -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
620 
621    Level: advanced
622 
623 .keywords: Sundials, orthogonalization
624 
625 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
626           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
627           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
628           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
629           TSSundialsSetExactFinalTime()
630 
631 @*/
632 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
633 {
634   PetscErrorCode ierr,(*f)(TS,TSSundialsGramSchmidtType);
635 
636   PetscFunctionBegin;
637   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",(void (**)(void))&f);CHKERRQ(ierr);
638   if (f) {
639     ierr = (*f)(ts,type);CHKERRQ(ierr);
640   }
641   PetscFunctionReturn(0);
642 }
643 
644 #undef __FUNCT__
645 #define __FUNCT__ "TSSundialsSetTolerance"
646 /*@
647    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
648                          Sundials for error control.
649 
650    Collective on TS
651 
652    Input parameters:
653 +    ts  - the time-step context
654 .    aabs - the absolute tolerance
655 -    rel - the relative tolerance
656 
657      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
658     these regulate the size of the error for a SINGLE timestep.
659 
660    Level: intermediate
661 
662 .keywords: Sundials, tolerance
663 
664 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
665           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
666           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
667           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
668           TSSundialsSetExactFinalTime()
669 
670 @*/
671 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance(TS ts,double aabs,double rel)
672 {
673   PetscErrorCode ierr,(*f)(TS,double,double);
674 
675   PetscFunctionBegin;
676   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
677   if (f) {
678     ierr = (*f)(ts,aabs,rel);CHKERRQ(ierr);
679   }
680   PetscFunctionReturn(0);
681 }
682 
683 #undef __FUNCT__
684 #define __FUNCT__ "TSSundialsGetPC"
685 /*@
686    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
687 
688    Input Parameter:
689 .    ts - the time-step context
690 
691    Output Parameter:
692 .    pc - the preconditioner context
693 
694    Level: advanced
695 
696 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
697           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
698           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
699           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
700 @*/
701 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC(TS ts,PC *pc)
702 {
703   PetscErrorCode ierr,(*f)(TS,PC *);
704 
705   PetscFunctionBegin;
706   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetPC_C",(void (**)(void))&f);CHKERRQ(ierr);
707   if (f) {
708     ierr = (*f)(ts,pc);CHKERRQ(ierr);
709   } else {
710     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"TS must be of Sundials type to extract the PC");
711   }
712 
713   PetscFunctionReturn(0);
714 }
715 
716 #undef __FUNCT__
717 #define __FUNCT__ "TSSundialsSetExactFinalTime"
718 /*@
719    TSSundialsSetExactFinalTime - Determines if Sundials interpolates solution to the
720       exact final time requested by the user or just returns it at the final time
721       it computed. (Defaults to true).
722 
723    Input Parameter:
724 +   ts - the time-step context
725 -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE
726 
727    Level: beginner
728 
729 .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
730           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
731           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
732           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
733 @*/
734 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime(TS ts,PetscTruth ft)
735 {
736   PetscErrorCode ierr,(*f)(TS,PetscTruth);
737 
738   PetscFunctionBegin;
739   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetExactFinalTime_C",(void (**)(void))&f);CHKERRQ(ierr);
740   if (f) {
741     ierr = (*f)(ts,ft);CHKERRQ(ierr);
742   }
743   PetscFunctionReturn(0);
744 }
745 
746 /* -------------------------------------------------------------------------------------------*/
747 /*MC
748       TS_Sundials - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
749 
750    Options Database:
751 +    -ts_sundials_type <bdf,adams>
752 .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
753 .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
754 .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
755 .    -ts_sundials_linear_tolerance <tol>
756 .    -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions
757 -    -ts_sundials_not_exact_final_time -Allow SUNDIALS to stop near the final time, not exactly on it
758 
759     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
760            only PETSc PC options
761 
762     Level: beginner
763 
764 .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(),
765            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSundialsSetExactFinalTime()
766 
767 M*/
768 EXTERN_C_BEGIN
769 #undef __FUNCT__
770 #define __FUNCT__ "TSCreate_Sundials"
771 PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Sundials(TS ts)
772 {
773   TS_Sundials *cvode;
774   PetscErrorCode ierr;
775 
776   PetscFunctionBegin;
777   ts->ops->destroy         = TSDestroy_Sundials;
778   ts->ops->view            = TSView_Sundials;
779 
780   if (ts->problem_type != TS_NONLINEAR) {
781     SETERRQ(PETSC_ERR_SUP,"Only support for nonlinear problems");
782   }
783   ts->ops->setup           = TSSetUp_Sundials_Nonlinear;
784   ts->ops->step            = TSStep_Sundials_Nonlinear;
785   ts->ops->setfromoptions  = TSSetFromOptions_Sundials_Nonlinear;
786 
787   ierr = PetscNew(TS_Sundials,&cvode);CHKERRQ(ierr);
788   ierr = PCCreate(ts->comm,&cvode->pc);CHKERRQ(ierr);
789   ierr = PetscLogObjectParent(ts,cvode->pc);CHKERRQ(ierr);
790   ts->data          = (void*)cvode;
791   cvode->cvode_type = CV_BDF;
792   cvode->gtype      = SUNDIALS_UNMODIFIED_GS;
793   cvode->restart    = 5;
794   cvode->linear_tol = .05;
795 
796   cvode->exact_final_time = PETSC_FALSE;
797 
798   ierr = MPI_Comm_dup(ts->comm,&(cvode->comm_sundials));CHKERRQ(ierr);
799   /* set tolerance for Sundials */
800   cvode->abstol = 1e-6;
801   cvode->reltol = 1e-6;
802 
803   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
804                     TSSundialsSetType_Sundials);CHKERRQ(ierr);
805   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C",
806                     "TSSundialsSetGMRESRestart_Sundials",
807                     TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr);
808   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
809                     "TSSundialsSetLinearTolerance_Sundials",
810                      TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
811   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
812                     "TSSundialsSetGramSchmidtType_Sundials",
813                      TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
814   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
815                     "TSSundialsSetTolerance_Sundials",
816                      TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
817   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
818                     "TSSundialsGetPC_Sundials",
819                      TSSundialsGetPC_Sundials);CHKERRQ(ierr);
820   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
821                     "TSSundialsGetIterations_Sundials",
822                      TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
823   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C",
824                     "TSSundialsSetExactFinalTime_Sundials",
825                      TSSundialsSetExactFinalTime_Sundials);CHKERRQ(ierr);
826   PetscFunctionReturn(0);
827 }
828 EXTERN_C_END
829 
830 
831 
832 
833 
834 
835 
836 
837 
838 
839