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