xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision 780a600e80509eba75ebd5b520fae1ebfeb16014)
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 = ((PetscObject)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(&(((PetscObject)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(((PetscObject)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(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr);
278   ierr = VecCreateMPIWithArray(((PetscObject)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 /* type of CVODE linear multistep method */
285 const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0};
286 /* type of G-S orthogonalization used by CVODE linear solver */
287 const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0};
288 
289 #undef __FUNCT__
290 #define __FUNCT__ "TSSetFromOptions_Sundials_Nonlinear"
291 PetscErrorCode TSSetFromOptions_Sundials_Nonlinear(TS ts)
292 {
293   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
294   PetscErrorCode ierr;
295   int            indx;
296   PetscTruth     flag;
297 
298   PetscFunctionBegin;
299   ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr);
300     ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr);
301     if (flag) {
302       ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr);
303     }
304     ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr);
305     if (flag) {
306       ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr);
307     }
308     ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr);
309     ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr);
310     ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr);
311     ierr = PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);CHKERRQ(ierr);
312     ierr = PetscOptionsName("-ts_sundials_exact_final_time","Allow SUNDIALS to stop near the final time, not exactly on it","TSSundialsSetExactFinalTime",&flag);CHKERRQ(ierr);
313     if (flag) cvode->exact_final_time = PETSC_TRUE;
314   ierr = PetscOptionsTail();CHKERRQ(ierr);
315   PetscFunctionReturn(0);
316 }
317 
318 #undef __FUNCT__
319 #define __FUNCT__ "TSView_Sundials"
320 PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
321 {
322   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
323   PetscErrorCode ierr;
324   char           *type;
325   char            atype[] = "Adams";
326   char            btype[] = "BDF: backward differentiation formula";
327   PetscTruth     iascii,isstring;
328 
329   PetscFunctionBegin;
330   if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;}
331   else                                     {type = btype;}
332 
333   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
334   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
335   if (iascii) {
336     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr);
337     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr);
338     ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr);
339     ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr);
340     ierr = PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);CHKERRQ(ierr);
341     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
342       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
343     } else {
344       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
345     }
346   } else if (isstring) {
347     ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr);
348   } else {
349     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by TS Sundials",((PetscObject)viewer)->type_name);
350   }
351   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
352   ierr = PCView(cvode->pc,viewer);CHKERRQ(ierr);
353   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
354   PetscFunctionReturn(0);
355 }
356 
357 
358 /* --------------------------------------------------------------------------*/
359 EXTERN_C_BEGIN
360 #undef __FUNCT__
361 #define __FUNCT__ "TSSundialsSetType_Sundials"
362 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
363 {
364   TS_Sundials *cvode = (TS_Sundials*)ts->data;
365 
366   PetscFunctionBegin;
367   cvode->cvode_type = type;
368   PetscFunctionReturn(0);
369 }
370 EXTERN_C_END
371 
372 EXTERN_C_BEGIN
373 #undef __FUNCT__
374 #define __FUNCT__ "TSSundialsSetGMRESRestart_Sundials"
375 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart_Sundials(TS ts,int restart)
376 {
377   TS_Sundials *cvode = (TS_Sundials*)ts->data;
378 
379   PetscFunctionBegin;
380   cvode->restart = restart;
381   PetscFunctionReturn(0);
382 }
383 EXTERN_C_END
384 
385 EXTERN_C_BEGIN
386 #undef __FUNCT__
387 #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials"
388 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
389 {
390   TS_Sundials *cvode = (TS_Sundials*)ts->data;
391 
392   PetscFunctionBegin;
393   cvode->linear_tol = tol;
394   PetscFunctionReturn(0);
395 }
396 EXTERN_C_END
397 
398 EXTERN_C_BEGIN
399 #undef __FUNCT__
400 #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials"
401 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
402 {
403   TS_Sundials *cvode = (TS_Sundials*)ts->data;
404 
405   PetscFunctionBegin;
406   cvode->gtype = type;
407   PetscFunctionReturn(0);
408 }
409 EXTERN_C_END
410 
411 EXTERN_C_BEGIN
412 #undef __FUNCT__
413 #define __FUNCT__ "TSSundialsSetTolerance_Sundials"
414 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
415 {
416   TS_Sundials *cvode = (TS_Sundials*)ts->data;
417 
418   PetscFunctionBegin;
419   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
420   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
421   PetscFunctionReturn(0);
422 }
423 EXTERN_C_END
424 
425 EXTERN_C_BEGIN
426 #undef __FUNCT__
427 #define __FUNCT__ "TSSundialsGetPC_Sundials"
428 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC_Sundials(TS ts,PC *pc)
429 {
430   TS_Sundials *cvode = (TS_Sundials*)ts->data;
431 
432   PetscFunctionBegin;
433   *pc = cvode->pc;
434   PetscFunctionReturn(0);
435 }
436 EXTERN_C_END
437 
438 EXTERN_C_BEGIN
439 #undef __FUNCT__
440 #define __FUNCT__ "TSSundialsGetIterations_Sundials"
441 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
442 {
443   TS_Sundials *cvode = (TS_Sundials*)ts->data;
444 
445   PetscFunctionBegin;
446   if (nonlin) *nonlin = cvode->nonlinear_solves;
447   if (lin)    *lin    = cvode->linear_solves;
448   PetscFunctionReturn(0);
449 }
450 EXTERN_C_END
451 
452 EXTERN_C_BEGIN
453 #undef __FUNCT__
454 #define __FUNCT__ "TSSundialsSetExactFinalTime_Sundials"
455 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime_Sundials(TS ts,PetscTruth s)
456 {
457   TS_Sundials *cvode = (TS_Sundials*)ts->data;
458 
459   PetscFunctionBegin;
460   cvode->exact_final_time = s;
461   PetscFunctionReturn(0);
462 }
463 EXTERN_C_END
464 /* -------------------------------------------------------------------------------------------*/
465 
466 #undef __FUNCT__
467 #define __FUNCT__ "TSSundialsGetIterations"
468 /*@C
469    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
470 
471    Not Collective
472 
473    Input parameters:
474 .    ts     - the time-step context
475 
476    Output Parameters:
477 +   nonlin - number of nonlinear iterations
478 -   lin    - number of linear iterations
479 
480    Level: advanced
481 
482    Notes:
483     These return the number since the creation of the TS object
484 
485 .keywords: non-linear iterations, linear iterations
486 
487 .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(),
488           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
489           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
490           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
491           TSSundialsSetExactFinalTime()
492 
493 @*/
494 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
495 {
496   PetscErrorCode ierr,(*f)(TS,int*,int*);
497 
498   PetscFunctionBegin;
499   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetIterations_C",(void (**)(void))&f);CHKERRQ(ierr);
500   if (f) {
501     ierr = (*f)(ts,nonlin,lin);CHKERRQ(ierr);
502   }
503   PetscFunctionReturn(0);
504 }
505 
506 #undef __FUNCT__
507 #define __FUNCT__ "TSSundialsSetType"
508 /*@
509    TSSundialsSetType - Sets the method that Sundials will use for integration.
510 
511    Collective on TS
512 
513    Input parameters:
514 +    ts     - the time-step context
515 -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
516 
517    Level: intermediate
518 
519 .keywords: Adams, backward differentiation formula
520 
521 .seealso: TSSundialsGetIterations(),  TSSundialsSetGMRESRestart(),
522           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
523           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
524           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
525           TSSundialsSetExactFinalTime()
526 @*/
527 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetType(TS ts,TSSundialsLmmType type)
528 {
529   PetscErrorCode ierr,(*f)(TS,TSSundialsLmmType);
530 
531   PetscFunctionBegin;
532   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
533   if (f) {
534     ierr = (*f)(ts,type);CHKERRQ(ierr);
535   }
536   PetscFunctionReturn(0);
537 }
538 
539 #undef __FUNCT__
540 #define __FUNCT__ "TSSundialsSetGMRESRestart"
541 /*@
542    TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by
543        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
544        this is ALSO the maximum number of GMRES steps that will be used.
545 
546    Collective on TS
547 
548    Input parameters:
549 +    ts      - the time-step context
550 -    restart - number of direction vectors (the restart size).
551 
552    Level: advanced
553 
554 .keywords: GMRES, restart
555 
556 .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
557           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
558           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
559           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
560           TSSundialsSetExactFinalTime()
561 
562 @*/
563 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGMRESRestart(TS ts,int restart)
564 {
565   PetscErrorCode ierr,(*f)(TS,int);
566 
567   PetscFunctionBegin;
568   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGMRESRestart_C",(void (**)(void))&f);CHKERRQ(ierr);
569   if (f) {
570     ierr = (*f)(ts,restart);CHKERRQ(ierr);
571   }
572   PetscFunctionReturn(0);
573 }
574 
575 #undef __FUNCT__
576 #define __FUNCT__ "TSSundialsSetLinearTolerance"
577 /*@
578    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
579        system by SUNDIALS.
580 
581    Collective on TS
582 
583    Input parameters:
584 +    ts     - the time-step context
585 -    tol    - the factor by which the tolerance on the nonlinear solver is
586              multiplied to get the tolerance on the linear solver, .05 by default.
587 
588    Level: advanced
589 
590 .keywords: GMRES, linear convergence tolerance, SUNDIALS
591 
592 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
593           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
594           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
595           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
596           TSSundialsSetExactFinalTime()
597 
598 @*/
599 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetLinearTolerance(TS ts,double tol)
600 {
601   PetscErrorCode ierr,(*f)(TS,double);
602 
603   PetscFunctionBegin;
604   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
605   if (f) {
606     ierr = (*f)(ts,tol);CHKERRQ(ierr);
607   }
608   PetscFunctionReturn(0);
609 }
610 
611 #undef __FUNCT__
612 #define __FUNCT__ "TSSundialsSetGramSchmidtType"
613 /*@
614    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
615         in GMRES method by SUNDIALS linear solver.
616 
617    Collective on TS
618 
619    Input parameters:
620 +    ts  - the time-step context
621 -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
622 
623    Level: advanced
624 
625 .keywords: Sundials, orthogonalization
626 
627 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
628           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
629           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
630           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
631           TSSundialsSetExactFinalTime()
632 
633 @*/
634 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
635 {
636   PetscErrorCode ierr,(*f)(TS,TSSundialsGramSchmidtType);
637 
638   PetscFunctionBegin;
639   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",(void (**)(void))&f);CHKERRQ(ierr);
640   if (f) {
641     ierr = (*f)(ts,type);CHKERRQ(ierr);
642   }
643   PetscFunctionReturn(0);
644 }
645 
646 #undef __FUNCT__
647 #define __FUNCT__ "TSSundialsSetTolerance"
648 /*@
649    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
650                          Sundials for error control.
651 
652    Collective on TS
653 
654    Input parameters:
655 +    ts  - the time-step context
656 .    aabs - the absolute tolerance
657 -    rel - the relative tolerance
658 
659      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
660     these regulate the size of the error for a SINGLE timestep.
661 
662    Level: intermediate
663 
664 .keywords: Sundials, tolerance
665 
666 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
667           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
668           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
669           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
670           TSSundialsSetExactFinalTime()
671 
672 @*/
673 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetTolerance(TS ts,double aabs,double rel)
674 {
675   PetscErrorCode ierr,(*f)(TS,double,double);
676 
677   PetscFunctionBegin;
678   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetTolerance_C",(void (**)(void))&f);CHKERRQ(ierr);
679   if (f) {
680     ierr = (*f)(ts,aabs,rel);CHKERRQ(ierr);
681   }
682   PetscFunctionReturn(0);
683 }
684 
685 #undef __FUNCT__
686 #define __FUNCT__ "TSSundialsGetPC"
687 /*@
688    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
689 
690    Input Parameter:
691 .    ts - the time-step context
692 
693    Output Parameter:
694 .    pc - the preconditioner context
695 
696    Level: advanced
697 
698 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
699           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
700           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
701           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
702 @*/
703 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsGetPC(TS ts,PC *pc)
704 {
705   PetscErrorCode ierr,(*f)(TS,PC *);
706 
707   PetscFunctionBegin;
708   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetPC_C",(void (**)(void))&f);CHKERRQ(ierr);
709   if (f) {
710     ierr = (*f)(ts,pc);CHKERRQ(ierr);
711   } else {
712     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"TS must be of Sundials type to extract the PC");
713   }
714 
715   PetscFunctionReturn(0);
716 }
717 
718 #undef __FUNCT__
719 #define __FUNCT__ "TSSundialsSetExactFinalTime"
720 /*@
721    TSSundialsSetExactFinalTime - Determines if Sundials interpolates solution to the
722       exact final time requested by the user or just returns it at the final time
723       it computed. (Defaults to true).
724 
725    Input Parameter:
726 +   ts - the time-step context
727 -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE
728 
729    Level: beginner
730 
731 .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
732           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
733           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
734           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
735 @*/
736 PetscErrorCode PETSCTS_DLLEXPORT TSSundialsSetExactFinalTime(TS ts,PetscTruth ft)
737 {
738   PetscErrorCode ierr,(*f)(TS,PetscTruth);
739 
740   PetscFunctionBegin;
741   ierr = PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetExactFinalTime_C",(void (**)(void))&f);CHKERRQ(ierr);
742   if (f) {
743     ierr = (*f)(ts,ft);CHKERRQ(ierr);
744   }
745   PetscFunctionReturn(0);
746 }
747 
748 /* -------------------------------------------------------------------------------------------*/
749 /*MC
750       TS_Sundials - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
751 
752    Options Database:
753 +    -ts_sundials_type <bdf,adams>
754 .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
755 .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
756 .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
757 .    -ts_sundials_linear_tolerance <tol>
758 .    -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions
759 -    -ts_sundials_not_exact_final_time -Allow SUNDIALS to stop near the final time, not exactly on it
760 
761     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
762            only PETSc PC options
763 
764     Level: beginner
765 
766 .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(),
767            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSundialsSetExactFinalTime()
768 
769 M*/
770 EXTERN_C_BEGIN
771 #undef __FUNCT__
772 #define __FUNCT__ "TSCreate_Sundials"
773 PetscErrorCode PETSCTS_DLLEXPORT TSCreate_Sundials(TS ts)
774 {
775   TS_Sundials *cvode;
776   PetscErrorCode ierr;
777 
778   PetscFunctionBegin;
779   ts->ops->destroy         = TSDestroy_Sundials;
780   ts->ops->view            = TSView_Sundials;
781 
782   if (ts->problem_type != TS_NONLINEAR) {
783     SETERRQ(PETSC_ERR_SUP,"Only support for nonlinear problems");
784   }
785   ts->ops->setup           = TSSetUp_Sundials_Nonlinear;
786   ts->ops->step            = TSStep_Sundials_Nonlinear;
787   ts->ops->setfromoptions  = TSSetFromOptions_Sundials_Nonlinear;
788 
789   ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr);
790   ierr = PCCreate(((PetscObject)ts)->comm,&cvode->pc);CHKERRQ(ierr);
791   ierr = PetscLogObjectParent(ts,cvode->pc);CHKERRQ(ierr);
792   ts->data          = (void*)cvode;
793   cvode->cvode_type = SUNDIALS_BDF;
794   cvode->gtype      = SUNDIALS_CLASSICAL_GS;
795   cvode->restart    = 5;
796   cvode->linear_tol = .05;
797 
798   cvode->exact_final_time = PETSC_FALSE;
799 
800   ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(((PetscObject)cvode)->comm_sundials));CHKERRQ(ierr);
801   /* set tolerance for Sundials */
802   cvode->abstol = 1e-6;
803   cvode->reltol = 1e-6;
804 
805   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
806                     TSSundialsSetType_Sundials);CHKERRQ(ierr);
807   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C",
808                     "TSSundialsSetGMRESRestart_Sundials",
809                     TSSundialsSetGMRESRestart_Sundials);CHKERRQ(ierr);
810   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
811                     "TSSundialsSetLinearTolerance_Sundials",
812                      TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
813   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
814                     "TSSundialsSetGramSchmidtType_Sundials",
815                      TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
816   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
817                     "TSSundialsSetTolerance_Sundials",
818                      TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
819   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
820                     "TSSundialsGetPC_Sundials",
821                      TSSundialsGetPC_Sundials);CHKERRQ(ierr);
822   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
823                     "TSSundialsGetIterations_Sundials",
824                      TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
825   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C",
826                     "TSSundialsSetExactFinalTime_Sundials",
827                      TSSundialsSetExactFinalTime_Sundials);CHKERRQ(ierr);
828   PetscFunctionReturn(0);
829 }
830 EXTERN_C_END
831 
832 
833 
834 
835 
836 
837 
838 
839 
840 
841