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