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