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