xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision a201590f66d78a0ee8e1bf4e6f9fcce39f32a809)
1 /*
2     Provides a PETSc interface to SUNDIALS/CVODE solver.
3     The interface to PVODE (old version of CVODE) was originally contributed
4     by Liyang Xu. It has been redone by Hong Zhang and Dinesh Kaushik.
5 
6     Reference: sundials-2.4.0/examples/cvode/parallel/cvDiurnal_kry_p.c
7 */
8 #include "sundials.h"  /*I "petscts.h" I*/
9 
10 /*
11       TSPrecond_Sundials - function that we provide to SUNDIALS to
12                         evaluate the preconditioner.
13 */
14 #undef __FUNCT__
15 #define __FUNCT__ "TSPrecond_Sundials"
16 PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy,
17                     booleantype jok,booleantype *jcurPtr,
18                     realtype _gamma,void *P_data,
19                     N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3)
20 {
21   TS             ts = (TS) P_data;
22   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
23   PC             pc;
24   PetscErrorCode ierr;
25   Mat            J,P;
26   Vec            yy = cvode->w1,yydot = cvode->ydot;
27   PetscReal      gm = (PetscReal)_gamma;
28   MatStructure   str = DIFFERENT_NONZERO_PATTERN;
29   PetscScalar    *y_data;
30 
31   PetscFunctionBegin;
32   ierr = TSGetIJacobian(ts,&J,&P,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
33   y_data = (PetscScalar *) N_VGetArrayPointer(y);
34   ierr = VecPlaceArray(yy,y_data); CHKERRQ(ierr);
35   ierr = VecZeroEntries(yydot);CHKERRQ(ierr); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */
36   /* compute the shifted Jacobian   (1/gm)*I + Jrest */
37   ierr = TSComputeIJacobian(ts,ts->ptime,yy,yydot,1/gm,&J,&P,&str,PETSC_FALSE);CHKERRQ(ierr);
38   ierr = VecResetArray(yy); CHKERRQ(ierr);
39   ierr = MatScale(P,gm);CHKERRQ(ierr);  /* turn into I-gm*Jrest, J is not used by Sundials  */
40   *jcurPtr = TRUE;
41   ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
42   ierr = PCSetOperators(pc,J,P,str);CHKERRQ(ierr);
43   PetscFunctionReturn(0);
44 }
45 
46 /*
47      TSPSolve_Sundials -  routine that we provide to Sundials that applies the preconditioner.
48 */
49 #undef __FUNCT__
50 #define __FUNCT__ "TSPSolve_Sundials"
51 PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z,
52                                  realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp)
53 {
54   TS              ts = (TS) P_data;
55   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
56   PC              pc;
57   Vec             rr = cvode->w1,zz = cvode->w2;
58   PetscErrorCode  ierr;
59   PetscScalar     *r_data,*z_data;
60 
61   PetscFunctionBegin;
62   /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
63   r_data  = (PetscScalar *) N_VGetArrayPointer(r);
64   z_data  = (PetscScalar *) N_VGetArrayPointer(z);
65   ierr = VecPlaceArray(rr,r_data); CHKERRQ(ierr);
66   ierr = VecPlaceArray(zz,z_data); CHKERRQ(ierr);
67 
68   /* Solve the Px=r and put the result in zz */
69   ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
70   ierr = PCApply(pc,rr,zz); CHKERRQ(ierr);
71   ierr = VecResetArray(rr); CHKERRQ(ierr);
72   ierr = VecResetArray(zz); CHKERRQ(ierr);
73   PetscFunctionReturn(0);
74 }
75 
76 /*
77         TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
78 */
79 #undef __FUNCT__
80 #define __FUNCT__ "TSFunction_Sundials"
81 int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx)
82 {
83   TS              ts = (TS) ctx;
84   MPI_Comm        comm = ((PetscObject)ts)->comm;
85   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
86   Vec             yy = cvode->w1,yyd = cvode->w2,yydot = cvode->ydot;
87   PetscScalar     *y_data,*ydot_data;
88   PetscErrorCode  ierr;
89 
90   PetscFunctionBegin;
91   /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
92   y_data     = (PetscScalar *) N_VGetArrayPointer(y);
93   ydot_data  = (PetscScalar *) N_VGetArrayPointer(ydot);
94   ierr = VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr);
95   ierr = VecPlaceArray(yyd,ydot_data); CHKERRABORT(comm,ierr);
96 
97   /* now compute the right hand side function */
98   if (!ts->userops->ifunction) {
99     ierr = TSComputeRHSFunction(ts,t,yy,yyd);CHKERRQ(ierr);
100   } else {                      /* If rhsfunction is also set, this computes both parts and shifts them to the right */
101     ierr = VecZeroEntries(yydot);CHKERRQ(ierr);
102     ierr = TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE); CHKERRABORT(comm,ierr);
103     ierr = VecScale(yyd,-1.);CHKERRQ(ierr);
104   }
105   ierr = VecResetArray(yy); CHKERRABORT(comm,ierr);
106   ierr = VecResetArray(yyd); CHKERRABORT(comm,ierr);
107   PetscFunctionReturn(0);
108 }
109 
110 /*
111        TSStep_Sundials - Calls Sundials to integrate the ODE.
112 */
113 #undef __FUNCT__
114 #define __FUNCT__ "TSStep_Sundials"
115 PetscErrorCode TSStep_Sundials(TS ts)
116 {
117   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
118   PetscErrorCode ierr;
119   PetscInt       flag;
120   long int       its,nsteps;
121   realtype       t,tout;
122   PetscScalar    *y_data;
123   void           *mem;
124 
125   PetscFunctionBegin;
126   mem  = cvode->mem;
127   tout = ts->max_time;
128   ierr = VecGetArray(ts->vec_sol,&y_data);CHKERRQ(ierr);
129   N_VSetArrayPointer((realtype *)y_data,cvode->y);
130   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
131 
132   ierr = TSPreStep(ts);CHKERRQ(ierr);
133 
134   if (cvode->monitorstep) {
135     flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
136   } else {
137     flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL);
138   }
139 
140   if (flag){ /* display error message */
141     switch (flag){
142       case CV_ILL_INPUT:
143         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT");
144         break;
145       case CV_TOO_CLOSE:
146         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE");
147         break;
148       case CV_TOO_MUCH_WORK: {
149         PetscReal      tcur;
150         ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
151         ierr = CVodeGetCurrentTime(mem,&tcur);CHKERRQ(ierr);
152         SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_WORK. At t=%G, nsteps %D exceeds mxstep %D. Increase '-ts_max_steps <>' or modify TSSetDuration()",tcur,nsteps,ts->max_steps);
153       } break;
154       case CV_TOO_MUCH_ACC:
155         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC");
156         break;
157       case CV_ERR_FAILURE:
158         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE");
159         break;
160       case CV_CONV_FAILURE:
161         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE");
162         break;
163       case CV_LINIT_FAIL:
164         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL");
165         break;
166       case CV_LSETUP_FAIL:
167         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL");
168         break;
169       case CV_LSOLVE_FAIL:
170         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL");
171         break;
172       case CV_RHSFUNC_FAIL:
173         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL");
174         break;
175       case CV_FIRST_RHSFUNC_ERR:
176         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR");
177         break;
178       case CV_REPTD_RHSFUNC_ERR:
179         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR");
180         break;
181       case CV_UNREC_RHSFUNC_ERR:
182         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR");
183         break;
184       case CV_RTFUNC_FAIL:
185         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL");
186         break;
187       default:
188         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag);
189     }
190   }
191 
192   /* copy the solution from cvode->y to cvode->update and sol */
193   ierr = VecPlaceArray(cvode->w1,y_data); CHKERRQ(ierr);
194   ierr = VecCopy(cvode->w1,cvode->update);CHKERRQ(ierr);
195   ierr = VecResetArray(cvode->w1); CHKERRQ(ierr);
196   ierr = VecCopy(cvode->update,ts->vec_sol);CHKERRQ(ierr);
197   ierr = CVodeGetNumNonlinSolvIters(mem,&its);CHKERRQ(ierr);
198   ierr = CVSpilsGetNumLinIters(mem, &its);
199   ts->snes_its = its; ts->ksp_its = its;
200 
201   ts->time_step = t - ts->ptime;
202   ts->ptime     = t;
203   ts->steps++;
204 
205   ierr = CVodeGetNumSteps(mem,&nsteps);CHKERRQ(ierr);
206   if (!cvode->monitorstep) ts->steps = nsteps;
207   PetscFunctionReturn(0);
208 }
209 
210 #undef __FUNCT__
211 #define __FUNCT__ "TSInterpolate_Sundials"
212 static PetscErrorCode TSInterpolate_Sundials(TS ts,PetscReal t,Vec X)
213 {
214   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
215   N_Vector        y;
216   PetscErrorCode  ierr;
217   PetscScalar     *x_data;
218   PetscInt        glosize,locsize;
219 
220   PetscFunctionBegin;
221 
222   /* get the vector size */
223   ierr = VecGetSize(X,&glosize);CHKERRQ(ierr);
224   ierr = VecGetLocalSize(X,&locsize);CHKERRQ(ierr);
225 
226   /* allocate the memory for N_Vec y */
227   y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
228   if (!y) SETERRQ(PETSC_COMM_SELF,1,"Interpolated y is not allocated");
229 
230   ierr = VecGetArray(X,&x_data);CHKERRQ(ierr);
231   N_VSetArrayPointer((realtype *)x_data,y);
232   ierr = CVodeGetDky(cvode->mem,t,0,y);CHKERRQ(ierr);
233   ierr = VecRestoreArray(X,&x_data);CHKERRQ(ierr);
234 
235   PetscFunctionReturn(0);
236 }
237 
238 #undef __FUNCT__
239 #define __FUNCT__ "TSReset_Sundials"
240 PetscErrorCode TSReset_Sundials(TS ts)
241 {
242   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
243   PetscErrorCode ierr;
244 
245   PetscFunctionBegin;
246   ierr = VecDestroy(&cvode->update);CHKERRQ(ierr);
247   ierr = VecDestroy(&cvode->ydot);CHKERRQ(ierr);
248   ierr = VecDestroy(&cvode->w1);CHKERRQ(ierr);
249   ierr = VecDestroy(&cvode->w2);CHKERRQ(ierr);
250   if (cvode->mem)    {CVodeFree(&cvode->mem);}
251   PetscFunctionReturn(0);
252 }
253 
254 #undef __FUNCT__
255 #define __FUNCT__ "TSDestroy_Sundials"
256 PetscErrorCode TSDestroy_Sundials(TS ts)
257 {
258   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
259   PetscErrorCode ierr;
260 
261   PetscFunctionBegin;
262   ierr = TSReset_Sundials(ts);CHKERRQ(ierr);
263   ierr = MPI_Comm_free(&(cvode->comm_sundials));CHKERRQ(ierr);
264   ierr = PetscFree(ts->data);CHKERRQ(ierr);
265   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","",PETSC_NULL);CHKERRQ(ierr);
266   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxl_C","",PETSC_NULL);CHKERRQ(ierr);
267   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C","",PETSC_NULL);CHKERRQ(ierr);
268   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C","",PETSC_NULL);CHKERRQ(ierr);
269   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C","",PETSC_NULL);CHKERRQ(ierr);
270   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
271   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C","",PETSC_NULL);CHKERRQ(ierr);
272   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C","",PETSC_NULL);CHKERRQ(ierr);
273   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C","",PETSC_NULL);CHKERRQ(ierr);
274   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C","",PETSC_NULL);CHKERRQ(ierr);
275   PetscFunctionReturn(0);
276 }
277 
278 #undef __FUNCT__
279 #define __FUNCT__ "TSSetUp_Sundials"
280 PetscErrorCode TSSetUp_Sundials(TS ts)
281 {
282   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
283   PetscErrorCode ierr;
284   PetscInt       glosize,locsize,i,flag;
285   PetscScalar    *y_data,*parray;
286   void           *mem;
287   PC             pc;
288   const PCType   pctype;
289   PetscBool      pcnone;
290 
291   PetscFunctionBegin;
292   /* get the vector size */
293   ierr = VecGetSize(ts->vec_sol,&glosize);CHKERRQ(ierr);
294   ierr = VecGetLocalSize(ts->vec_sol,&locsize);CHKERRQ(ierr);
295 
296   /* allocate the memory for N_Vec y */
297   cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
298   if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated");
299 
300   /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
301   ierr = VecGetArray(ts->vec_sol,&parray);CHKERRQ(ierr);
302   y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y);
303   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
304   ierr = VecRestoreArray(ts->vec_sol,PETSC_NULL);CHKERRQ(ierr);
305 
306   ierr = VecDuplicate(ts->vec_sol,&cvode->update);CHKERRQ(ierr);
307   ierr = VecDuplicate(ts->vec_sol,&cvode->ydot);CHKERRQ(ierr);
308   ierr = PetscLogObjectParent(ts,cvode->update);CHKERRQ(ierr);
309   ierr = PetscLogObjectParent(ts,cvode->ydot);CHKERRQ(ierr);
310 
311   /*
312     Create work vectors for the TSPSolve_Sundials() routine. Note these are
313     allocated with zero space arrays because the actual array space is provided
314     by Sundials and set using VecPlaceArray().
315   */
316   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,1,locsize,PETSC_DECIDE,0,&cvode->w1);CHKERRQ(ierr);
317   ierr = VecCreateMPIWithArray(((PetscObject)ts)->comm,1,locsize,PETSC_DECIDE,0,&cvode->w2);CHKERRQ(ierr);
318   ierr = PetscLogObjectParent(ts,cvode->w1);CHKERRQ(ierr);
319   ierr = PetscLogObjectParent(ts,cvode->w2);CHKERRQ(ierr);
320 
321   /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
322   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
323   if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails");
324   cvode->mem = mem;
325 
326   /* Set the pointer to user-defined data */
327   flag = CVodeSetUserData(mem, ts);
328   if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails");
329 
330   /* Sundials may choose to use a smaller initial step, but will never use a larger step. */
331   flag = CVodeSetInitStep(mem,(realtype)ts->time_step);
332   if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetInitStep() failed");
333   if (cvode->mindt > 0) {
334     flag = CVodeSetMinStep(mem,(realtype)cvode->mindt);
335     if (flag){
336       if (flag == CV_MEM_NULL){
337         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL");
338       } else if (flag == CV_ILL_INPUT){
339         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size");
340       } else {
341         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed");
342       }
343     }
344   }
345   if (cvode->maxdt > 0) {
346     flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt);
347     if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMaxStep() failed");
348   }
349 
350   /* Call CVodeInit to initialize the integrator memory and specify the
351    * user's right hand side function in u'=f(t,u), the inital time T0, and
352    * the initial dependent variable vector cvode->y */
353   flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y);
354   if (flag){
355     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag);
356   }
357 
358   /* specifies scalar relative and absolute tolerances */
359   flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol);
360   if (flag){
361     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag);
362   }
363 
364   /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */
365   flag = CVodeSetMaxNumSteps(mem,ts->max_steps);
366 
367   /* call CVSpgmr to use GMRES as the linear solver.        */
368   /* setup the ode integrator with the given preconditioner */
369   ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
370   ierr = PCGetType(pc,&pctype);CHKERRQ(ierr);
371   ierr = PetscObjectTypeCompare((PetscObject)pc,PCNONE,&pcnone);CHKERRQ(ierr);
372   if (pcnone){
373     flag  = CVSpgmr(mem,PREC_NONE,0);
374     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
375   } else {
376     flag  = CVSpgmr(mem,PREC_LEFT,cvode->maxl);
377     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
378 
379     /* Set preconditioner and solve routines Precond and PSolve,
380      and the pointer to the user-defined block data */
381     flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials);
382     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag);
383   }
384 
385   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
386   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag);
387   PetscFunctionReturn(0);
388 }
389 
390 /* type of CVODE linear multistep method */
391 const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0};
392 /* type of G-S orthogonalization used by CVODE linear solver */
393 const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0};
394 
395 #undef __FUNCT__
396 #define __FUNCT__ "TSSetFromOptions_Sundials"
397 PetscErrorCode TSSetFromOptions_Sundials(TS ts)
398 {
399   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
400   PetscErrorCode ierr;
401   int            indx;
402   PetscBool      flag;
403   PC             pc;
404 
405   PetscFunctionBegin;
406   ierr = PetscOptionsHead("SUNDIALS ODE solver options");CHKERRQ(ierr);
407     ierr = PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);CHKERRQ(ierr);
408     if (flag) {
409       ierr = TSSundialsSetType(ts,(TSSundialsLmmType)indx);CHKERRQ(ierr);
410     }
411     ierr = PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);CHKERRQ(ierr);
412     if (flag) {
413       ierr = TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);CHKERRQ(ierr);
414     }
415     ierr = PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);CHKERRQ(ierr);
416     ierr = PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);CHKERRQ(ierr);
417     ierr = PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,PETSC_NULL);CHKERRQ(ierr);
418     ierr = PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,PETSC_NULL);CHKERRQ(ierr);
419     ierr = PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);CHKERRQ(ierr);
420     ierr = PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,&flag);CHKERRQ(ierr);
421     ierr = PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);CHKERRQ(ierr);
422   ierr = PetscOptionsTail();CHKERRQ(ierr);
423   ierr = TSSundialsGetPC(ts,&pc);CHKERRQ(ierr);
424   ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
425   PetscFunctionReturn(0);
426 }
427 
428 #undef __FUNCT__
429 #define __FUNCT__ "TSView_Sundials"
430 PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
431 {
432   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
433   PetscErrorCode ierr;
434   char           *type;
435   char           atype[] = "Adams";
436   char           btype[] = "BDF: backward differentiation formula";
437   PetscBool      iascii,isstring;
438   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
439   PetscInt       qlast,qcur;
440   PetscReal      hinused,hlast,hcur,tcur,tolsfac;
441   PC             pc;
442 
443   PetscFunctionBegin;
444   if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;}
445   else                                     {type = btype;}
446 
447   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
448   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
449   if (iascii) {
450     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");CHKERRQ(ierr);
451     ierr = PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);CHKERRQ(ierr);
452     ierr = PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);CHKERRQ(ierr);
453     ierr = PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);CHKERRQ(ierr);
454     ierr = PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %D\n",cvode->maxl);CHKERRQ(ierr);
455     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
456       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
457     } else {
458       ierr = PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");CHKERRQ(ierr);
459     }
460     if (cvode->mindt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);CHKERRQ(ierr);}
461     if (cvode->maxdt > 0) {ierr = PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);CHKERRQ(ierr);}
462 
463     /* Outputs from CVODE, CVSPILS */
464     ierr = CVodeGetTolScaleFactor(cvode->mem,&tolsfac);CHKERRQ(ierr);
465     ierr = PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);CHKERRQ(ierr);
466     ierr = CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
467                                    &nlinsetups,&nfails,&qlast,&qcur,
468                                    &hinused,&hlast,&hcur,&tcur);CHKERRQ(ierr);
469     ierr = PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);CHKERRQ(ierr);
470     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);CHKERRQ(ierr);
471     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);CHKERRQ(ierr);
472     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);CHKERRQ(ierr);
473 
474     ierr = CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);CHKERRQ(ierr);
475     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);CHKERRQ(ierr);
476     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);CHKERRQ(ierr);
477 
478     ierr = CVSpilsGetNumLinIters(cvode->mem, &its);CHKERRQ(ierr); /* its = no. of calls to TSPrecond_Sundials() */
479     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);CHKERRQ(ierr);
480     ierr = CVSpilsGetNumConvFails(cvode->mem,&itmp);CHKERRQ(ierr);
481     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);CHKERRQ(ierr);
482 
483     ierr = TSSundialsGetPC(ts,&pc); CHKERRQ(ierr);
484     ierr = PCView(pc,viewer);CHKERRQ(ierr);
485     ierr = CVSpilsGetNumPrecEvals(cvode->mem,&itmp);CHKERRQ(ierr);
486     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);CHKERRQ(ierr);
487     ierr = CVSpilsGetNumPrecSolves(cvode->mem,&itmp);CHKERRQ(ierr);
488     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);CHKERRQ(ierr);
489 
490     ierr = CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);CHKERRQ(ierr);
491     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);CHKERRQ(ierr);
492     ierr = CVSpilsGetNumRhsEvals(cvode->mem,&itmp);CHKERRQ(ierr);
493     ierr = PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);CHKERRQ(ierr);
494   } else if (isstring) {
495     ierr = PetscViewerStringSPrintf(viewer,"Sundials type %s",type);CHKERRQ(ierr);
496   }
497   PetscFunctionReturn(0);
498 }
499 
500 
501 /* --------------------------------------------------------------------------*/
502 EXTERN_C_BEGIN
503 #undef __FUNCT__
504 #define __FUNCT__ "TSSundialsSetType_Sundials"
505 PetscErrorCode  TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
506 {
507   TS_Sundials *cvode = (TS_Sundials*)ts->data;
508 
509   PetscFunctionBegin;
510   cvode->cvode_type = type;
511   PetscFunctionReturn(0);
512 }
513 EXTERN_C_END
514 
515 EXTERN_C_BEGIN
516 #undef __FUNCT__
517 #define __FUNCT__ "TSSundialsSetMaxl_Sundials"
518 PetscErrorCode  TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl)
519 {
520   TS_Sundials *cvode = (TS_Sundials*)ts->data;
521 
522   PetscFunctionBegin;
523   cvode->maxl = maxl;
524   PetscFunctionReturn(0);
525 }
526 EXTERN_C_END
527 
528 EXTERN_C_BEGIN
529 #undef __FUNCT__
530 #define __FUNCT__ "TSSundialsSetLinearTolerance_Sundials"
531 PetscErrorCode  TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
532 {
533   TS_Sundials *cvode = (TS_Sundials*)ts->data;
534 
535   PetscFunctionBegin;
536   cvode->linear_tol = tol;
537   PetscFunctionReturn(0);
538 }
539 EXTERN_C_END
540 
541 EXTERN_C_BEGIN
542 #undef __FUNCT__
543 #define __FUNCT__ "TSSundialsSetGramSchmidtType_Sundials"
544 PetscErrorCode  TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
545 {
546   TS_Sundials *cvode = (TS_Sundials*)ts->data;
547 
548   PetscFunctionBegin;
549   cvode->gtype = type;
550   PetscFunctionReturn(0);
551 }
552 EXTERN_C_END
553 
554 EXTERN_C_BEGIN
555 #undef __FUNCT__
556 #define __FUNCT__ "TSSundialsSetTolerance_Sundials"
557 PetscErrorCode  TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
558 {
559   TS_Sundials *cvode = (TS_Sundials*)ts->data;
560 
561   PetscFunctionBegin;
562   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
563   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
564   PetscFunctionReturn(0);
565 }
566 EXTERN_C_END
567 
568 EXTERN_C_BEGIN
569 #undef __FUNCT__
570 #define __FUNCT__ "TSSundialsSetMinTimeStep_Sundials"
571 PetscErrorCode  TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)
572 {
573   TS_Sundials *cvode = (TS_Sundials*)ts->data;
574 
575   PetscFunctionBegin;
576   cvode->mindt = mindt;
577   PetscFunctionReturn(0);
578 }
579 EXTERN_C_END
580 
581 EXTERN_C_BEGIN
582 #undef __FUNCT__
583 #define __FUNCT__ "TSSundialsSetMaxTimeStep_Sundials"
584 PetscErrorCode  TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)
585 {
586   TS_Sundials *cvode = (TS_Sundials*)ts->data;
587 
588   PetscFunctionBegin;
589   cvode->maxdt = maxdt;
590   PetscFunctionReturn(0);
591 }
592 EXTERN_C_END
593 
594 EXTERN_C_BEGIN
595 #undef __FUNCT__
596 #define __FUNCT__ "TSSundialsGetPC_Sundials"
597 PetscErrorCode  TSSundialsGetPC_Sundials(TS ts,PC *pc)
598 {
599   SNES            snes;
600   KSP             ksp;
601   PetscErrorCode  ierr;
602 
603   PetscFunctionBegin;
604   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
605   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
606   ierr = KSPGetPC(ksp,pc); CHKERRQ(ierr);
607   PetscFunctionReturn(0);
608 }
609 EXTERN_C_END
610 
611 EXTERN_C_BEGIN
612 #undef __FUNCT__
613 #define __FUNCT__ "TSSundialsGetIterations_Sundials"
614 PetscErrorCode  TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
615 {
616   PetscFunctionBegin;
617   if (nonlin) *nonlin = ts->snes_its;
618   if (lin)    *lin    = ts->ksp_its;
619   PetscFunctionReturn(0);
620 }
621 EXTERN_C_END
622 
623 EXTERN_C_BEGIN
624 #undef __FUNCT__
625 #define __FUNCT__ "TSSundialsMonitorInternalSteps_Sundials"
626 PetscErrorCode  TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool  s)
627 {
628   TS_Sundials *cvode = (TS_Sundials*)ts->data;
629 
630   PetscFunctionBegin;
631   cvode->monitorstep = s;
632   PetscFunctionReturn(0);
633 }
634 EXTERN_C_END
635 /* -------------------------------------------------------------------------------------------*/
636 
637 #undef __FUNCT__
638 #define __FUNCT__ "TSSundialsGetIterations"
639 /*@C
640    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
641 
642    Not Collective
643 
644    Input parameters:
645 .    ts     - the time-step context
646 
647    Output Parameters:
648 +   nonlin - number of nonlinear iterations
649 -   lin    - number of linear iterations
650 
651    Level: advanced
652 
653    Notes:
654     These return the number since the creation of the TS object
655 
656 .keywords: non-linear iterations, linear iterations
657 
658 .seealso: TSSundialsSetType(), TSSundialsSetMaxl(),
659           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
660           TSSundialsGetIterations(), TSSundialsSetType(),
661           TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime()
662 
663 @*/
664 PetscErrorCode  TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
665 {
666   PetscErrorCode ierr;
667 
668   PetscFunctionBegin;
669   ierr = PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));CHKERRQ(ierr);
670   PetscFunctionReturn(0);
671 }
672 
673 #undef __FUNCT__
674 #define __FUNCT__ "TSSundialsSetType"
675 /*@
676    TSSundialsSetType - Sets the method that Sundials will use for integration.
677 
678    Logically Collective on TS
679 
680    Input parameters:
681 +    ts     - the time-step context
682 -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
683 
684    Level: intermediate
685 
686 .keywords: Adams, backward differentiation formula
687 
688 .seealso: TSSundialsGetIterations(),  TSSundialsSetMaxl(),
689           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
690           TSSundialsGetIterations(), TSSundialsSetType(),
691           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
692           TSSetExactFinalTime()
693 @*/
694 PetscErrorCode  TSSundialsSetType(TS ts,TSSundialsLmmType type)
695 {
696   PetscErrorCode ierr;
697 
698   PetscFunctionBegin;
699   ierr = PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));CHKERRQ(ierr);
700   PetscFunctionReturn(0);
701 }
702 
703 #undef __FUNCT__
704 #define __FUNCT__ "TSSundialsSetMaxl"
705 /*@
706    TSSundialsSetMaxl - Sets the dimension of the Krylov space used by
707        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
708        this is the maximum number of GMRES steps that will be used.
709 
710    Logically Collective on TS
711 
712    Input parameters:
713 +    ts      - the time-step context
714 -    maxl - number of direction vectors (the dimension of Krylov subspace).
715 
716    Level: advanced
717 
718 .keywords: GMRES
719 
720 .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
721           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
722           TSSundialsGetIterations(), TSSundialsSetType(),
723           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
724           TSSetExactFinalTime()
725 
726 @*/
727 PetscErrorCode  TSSundialsSetMaxl(TS ts,PetscInt maxl)
728 {
729   PetscErrorCode ierr;
730 
731   PetscFunctionBegin;
732   PetscValidLogicalCollectiveInt(ts,maxl,2);
733   ierr = PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));CHKERRQ(ierr);
734   PetscFunctionReturn(0);
735 }
736 
737 #undef __FUNCT__
738 #define __FUNCT__ "TSSundialsSetLinearTolerance"
739 /*@
740    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
741        system by SUNDIALS.
742 
743    Logically Collective on TS
744 
745    Input parameters:
746 +    ts     - the time-step context
747 -    tol    - the factor by which the tolerance on the nonlinear solver is
748              multiplied to get the tolerance on the linear solver, .05 by default.
749 
750    Level: advanced
751 
752 .keywords: GMRES, linear convergence tolerance, SUNDIALS
753 
754 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
755           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
756           TSSundialsGetIterations(), TSSundialsSetType(),
757           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
758           TSSetExactFinalTime()
759 
760 @*/
761 PetscErrorCode  TSSundialsSetLinearTolerance(TS ts,double tol)
762 {
763   PetscErrorCode ierr;
764 
765   PetscFunctionBegin;
766   PetscValidLogicalCollectiveReal(ts,tol,2);
767   ierr = PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));CHKERRQ(ierr);
768   PetscFunctionReturn(0);
769 }
770 
771 #undef __FUNCT__
772 #define __FUNCT__ "TSSundialsSetGramSchmidtType"
773 /*@
774    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
775         in GMRES method by SUNDIALS linear solver.
776 
777    Logically Collective on TS
778 
779    Input parameters:
780 +    ts  - the time-step context
781 -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
782 
783    Level: advanced
784 
785 .keywords: Sundials, orthogonalization
786 
787 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
788           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
789           TSSundialsGetIterations(), TSSundialsSetType(),
790           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
791           TSSetExactFinalTime()
792 
793 @*/
794 PetscErrorCode  TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
795 {
796   PetscErrorCode ierr;
797 
798   PetscFunctionBegin;
799   ierr = PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));CHKERRQ(ierr);
800   PetscFunctionReturn(0);
801 }
802 
803 #undef __FUNCT__
804 #define __FUNCT__ "TSSundialsSetTolerance"
805 /*@
806    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
807                          Sundials for error control.
808 
809    Logically Collective on TS
810 
811    Input parameters:
812 +    ts  - the time-step context
813 .    aabs - the absolute tolerance
814 -    rel - the relative tolerance
815 
816      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
817     these regulate the size of the error for a SINGLE timestep.
818 
819    Level: intermediate
820 
821 .keywords: Sundials, tolerance
822 
823 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESMaxl(),
824           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
825           TSSundialsGetIterations(), TSSundialsSetType(),
826           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
827           TSSetExactFinalTime()
828 
829 @*/
830 PetscErrorCode  TSSundialsSetTolerance(TS ts,double aabs,double rel)
831 {
832   PetscErrorCode ierr;
833 
834   PetscFunctionBegin;
835   ierr = PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));CHKERRQ(ierr);
836   PetscFunctionReturn(0);
837 }
838 
839 #undef __FUNCT__
840 #define __FUNCT__ "TSSundialsGetPC"
841 /*@
842    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
843 
844    Input Parameter:
845 .    ts - the time-step context
846 
847    Output Parameter:
848 .    pc - the preconditioner context
849 
850    Level: advanced
851 
852 .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
853           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
854           TSSundialsGetIterations(), TSSundialsSetType(),
855           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
856 @*/
857 PetscErrorCode  TSSundialsGetPC(TS ts,PC *pc)
858 {
859   PetscErrorCode ierr;
860 
861   PetscFunctionBegin;
862   ierr = PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC *),(ts,pc));CHKERRQ(ierr);
863   PetscFunctionReturn(0);
864 }
865 
866 #undef __FUNCT__
867 #define __FUNCT__ "TSSundialsSetMinTimeStep"
868 /*@
869    TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller.
870 
871    Input Parameter:
872 +   ts - the time-step context
873 -   mindt - lowest time step if positive, negative to deactivate
874 
875    Note:
876    Sundials will error if it is not possible to keep the estimated truncation error below
877    the tolerance set with TSSundialsSetTolerance() without going below this step size.
878 
879    Level: beginner
880 
881 .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
882 @*/
883 PetscErrorCode  TSSundialsSetMinTimeStep(TS ts,PetscReal mindt)
884 {
885   PetscErrorCode ierr;
886 
887   PetscFunctionBegin;
888   ierr = PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));CHKERRQ(ierr);
889   PetscFunctionReturn(0);
890 }
891 
892 #undef __FUNCT__
893 #define __FUNCT__ "TSSundialsSetMaxTimeStep"
894 /*@
895    TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller.
896 
897    Input Parameter:
898 +   ts - the time-step context
899 -   maxdt - lowest time step if positive, negative to deactivate
900 
901    Level: beginner
902 
903 .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
904 @*/
905 PetscErrorCode  TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt)
906 {
907   PetscErrorCode ierr;
908 
909   PetscFunctionBegin;
910   ierr = PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));CHKERRQ(ierr);
911   PetscFunctionReturn(0);
912 }
913 
914 #undef __FUNCT__
915 #define __FUNCT__ "TSSundialsMonitorInternalSteps"
916 /*@
917    TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false).
918 
919    Input Parameter:
920 +   ts - the time-step context
921 -   ft - PETSC_TRUE if monitor, else PETSC_FALSE
922 
923    Level: beginner
924 
925 .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
926           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
927           TSSundialsGetIterations(), TSSundialsSetType(),
928           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
929 @*/
930 PetscErrorCode  TSSundialsMonitorInternalSteps(TS ts,PetscBool  ft)
931 {
932   PetscErrorCode ierr;
933 
934   PetscFunctionBegin;
935   ierr = PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));CHKERRQ(ierr);
936   PetscFunctionReturn(0);
937 }
938 /* -------------------------------------------------------------------------------------------*/
939 /*MC
940       TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
941 
942    Options Database:
943 +    -ts_sundials_type <bdf,adams>
944 .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
945 .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
946 .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
947 .    -ts_sundials_linear_tolerance <tol>
948 .    -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace
949 -    -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps
950 
951 
952     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
953            only PETSc PC options
954 
955     Level: beginner
956 
957 .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetMaxl(), TSSundialsSetLinearTolerance(),
958            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime()
959 
960 M*/
961 EXTERN_C_BEGIN
962 #undef __FUNCT__
963 #define __FUNCT__ "TSCreate_Sundials"
964 PetscErrorCode  TSCreate_Sundials(TS ts)
965 {
966   TS_Sundials    *cvode;
967   PetscErrorCode ierr;
968   PC             pc;
969 
970   PetscFunctionBegin;
971   ts->ops->reset          = TSReset_Sundials;
972   ts->ops->destroy        = TSDestroy_Sundials;
973   ts->ops->view           = TSView_Sundials;
974   ts->ops->setup          = TSSetUp_Sundials;
975   ts->ops->step           = TSStep_Sundials;
976   ts->ops->interpolate    = TSInterpolate_Sundials;
977   ts->ops->setfromoptions = TSSetFromOptions_Sundials;
978 
979   ierr = PetscNewLog(ts,TS_Sundials,&cvode);CHKERRQ(ierr);
980   ts->data                = (void*)cvode;
981   cvode->cvode_type       = SUNDIALS_BDF;
982   cvode->gtype            = SUNDIALS_CLASSICAL_GS;
983   cvode->maxl             = 5;
984   cvode->linear_tol       = .05;
985 
986   cvode->monitorstep      = PETSC_TRUE;
987 
988   ierr = MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));CHKERRQ(ierr);
989 
990   cvode->mindt = -1.;
991   cvode->maxdt = -1.;
992 
993   /* set tolerance for Sundials */
994   cvode->reltol = 1e-6;
995   cvode->abstol = 1e-6;
996 
997   /* set PCNONE as default pctype */
998   ierr = TSSundialsGetPC_Sundials(ts,&pc);CHKERRQ(ierr);
999   ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1000 
1001   if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_FALSE;
1002 
1003   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
1004                     TSSundialsSetType_Sundials);CHKERRQ(ierr);
1005   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxl_C",
1006                     "TSSundialsSetMaxl_Sundials",
1007                     TSSundialsSetMaxl_Sundials);CHKERRQ(ierr);
1008   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
1009                     "TSSundialsSetLinearTolerance_Sundials",
1010                      TSSundialsSetLinearTolerance_Sundials);CHKERRQ(ierr);
1011   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
1012                     "TSSundialsSetGramSchmidtType_Sundials",
1013                      TSSundialsSetGramSchmidtType_Sundials);CHKERRQ(ierr);
1014   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
1015                     "TSSundialsSetTolerance_Sundials",
1016                      TSSundialsSetTolerance_Sundials);CHKERRQ(ierr);
1017   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C",
1018                     "TSSundialsSetMinTimeStep_Sundials",
1019                      TSSundialsSetMinTimeStep_Sundials);CHKERRQ(ierr);
1020   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",
1021                     "TSSundialsSetMaxTimeStep_Sundials",
1022                      TSSundialsSetMaxTimeStep_Sundials);CHKERRQ(ierr);
1023   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
1024                     "TSSundialsGetPC_Sundials",
1025                      TSSundialsGetPC_Sundials);CHKERRQ(ierr);
1026   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
1027                     "TSSundialsGetIterations_Sundials",
1028                      TSSundialsGetIterations_Sundials);CHKERRQ(ierr);
1029   ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",
1030                     "TSSundialsMonitorInternalSteps_Sundials",
1031                      TSSundialsMonitorInternalSteps_Sundials);CHKERRQ(ierr);
1032 
1033   PetscFunctionReturn(0);
1034 }
1035 EXTERN_C_END
1036