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