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