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