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