xref: /petsc/src/ts/impls/implicit/sundials/sundials.c (revision 4ffacfe27a72f4cdf51b68a3bbb6aed96040fb2f)
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   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetUseDense_C",NULL));
265   PetscFunctionReturn(0);
266 }
267 
268 PetscErrorCode TSSetUp_Sundials(TS ts)
269 {
270   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
271   PetscInt       glosize,locsize,i,flag;
272   PetscScalar    *y_data,*parray;
273   void           *mem;
274   PC             pc;
275   PCType         pctype;
276   PetscBool      pcnone;
277 
278   PetscFunctionBegin;
279   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");
280 
281   /* get the vector size */
282   PetscCall(VecGetSize(ts->vec_sol,&glosize));
283   PetscCall(VecGetLocalSize(ts->vec_sol,&locsize));
284 
285   /* allocate the memory for N_Vec y */
286   if (cvode->use_dense) {
287     PetscMPIInt size;
288 
289     PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
290     PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"TSSUNDIALS only supports a dense solve in the serial case");
291     cvode->y = N_VNew_Serial(locsize);
292   } else {
293     cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
294   }
295   PetscCheck(cvode->y,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cvode->y is not allocated");
296 
297   /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
298   PetscCall(VecGetArray(ts->vec_sol,&parray));
299   y_data = (PetscScalar*) N_VGetArrayPointer(cvode->y);
300   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
301   PetscCall(VecRestoreArray(ts->vec_sol,NULL));
302 
303   PetscCall(VecDuplicate(ts->vec_sol,&cvode->update));
304   PetscCall(VecDuplicate(ts->vec_sol,&cvode->ydot));
305   PetscCall(PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->update));
306   PetscCall(PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->ydot));
307 
308   /*
309     Create work vectors for the TSPSolve_Sundials() routine. Note these are
310     allocated with zero space arrays because the actual array space is provided
311     by Sundials and set using VecPlaceArray().
312   */
313   PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,NULL,&cvode->w1));
314   PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,NULL,&cvode->w2));
315   PetscCall(PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w1));
316   PetscCall(PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w2));
317 
318   /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
319   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
320   PetscCheck(mem,PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails");
321   cvode->mem = mem;
322 
323   /* Set the pointer to user-defined data */
324   flag = CVodeSetUserData(mem, ts);
325   PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails");
326 
327   /* Sundials may choose to use a smaller initial step, but will never use a larger step. */
328   flag = CVodeSetInitStep(mem,(realtype)ts->time_step);
329   PetscCheck(!flag,PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetInitStep() failed");
330   if (cvode->mindt > 0) {
331     flag = CVodeSetMinStep(mem,(realtype)cvode->mindt);
332     if (flag) {
333       PetscCheck(flag != CV_MEM_NULL,PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL");
334       PetscCheck(flag != CV_ILL_INPUT,PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size");
335       SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed");
336     }
337   }
338   if (cvode->maxdt > 0) {
339     flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt);
340     PetscCheck(!flag,PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMaxStep() failed");
341   }
342 
343   /* Call CVodeInit to initialize the integrator memory and specify the
344    * user's right hand side function in u'=f(t,u), the initial time T0, and
345    * the initial dependent variable vector cvode->y */
346   flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y);
347   PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag);
348 
349   /* specifies scalar relative and absolute tolerances */
350   flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol);
351   PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag);
352 
353   /* Specify max order of BDF / ADAMS method */
354   if (cvode->maxord != PETSC_DEFAULT) {
355     flag = CVodeSetMaxOrd(mem,cvode->maxord);
356     PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetMaxOrd() fails, flag %d",flag);
357   }
358 
359   /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */
360   flag = CVodeSetMaxNumSteps(mem,ts->max_steps);
361   PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetMaxNumSteps() fails, flag %d",flag);
362 
363   if (cvode->use_dense) {
364     /* call CVDense to use a dense linear solver. */
365     PetscMPIInt size;
366 
367     PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
368     PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"TSSUNDIALS only supports a dense solve in the serial case");
369     flag = CVDense(mem,locsize);
370     PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVDense() fails, flag %d",flag);
371   } else {
372     /* call CVSpgmr to use GMRES as the linear solver.        */
373     /* setup the ode integrator with the given preconditioner */
374     PetscCall(TSSundialsGetPC(ts,&pc));
375     PetscCall(PCGetType(pc,&pctype));
376     PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCNONE,&pcnone));
377     if (pcnone) {
378       flag = CVSpgmr(mem,PREC_NONE,0);
379       PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
380     } else {
381       flag = CVSpgmr(mem,PREC_LEFT,cvode->maxl);
382       PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
383 
384       /* Set preconditioner and solve routines Precond and PSolve,
385          and the pointer to the user-defined block data */
386       flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials);
387       PetscCheck(!flag,PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag);
388     }
389   }
390   PetscFunctionReturn(0);
391 }
392 
393 /* type of CVODE linear multistep method */
394 const char *const TSSundialsLmmTypes[] = {"","ADAMS","BDF","TSSundialsLmmType","SUNDIALS_",NULL};
395 /* type of G-S orthogonalization used by CVODE linear solver */
396 const char *const TSSundialsGramSchmidtTypes[] = {"","MODIFIED","CLASSICAL","TSSundialsGramSchmidtType","SUNDIALS_",NULL};
397 
398 PetscErrorCode TSSetFromOptions_Sundials(PetscOptionItems *PetscOptionsObject,TS ts)
399 {
400   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
401   int            indx;
402   PetscBool      flag;
403   PC             pc;
404 
405   PetscFunctionBegin;
406   PetscOptionsHeadBegin(PetscOptionsObject,"SUNDIALS ODE solver options");
407   PetscCall(PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag));
408   if (flag) PetscCall(TSSundialsSetType(ts,(TSSundialsLmmType)indx));
409   PetscCall(PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag));
410   if (flag) PetscCall(TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx));
411   PetscCall(PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,NULL));
412   PetscCall(PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,NULL));
413   PetscCall(PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,NULL));
414   PetscCall(PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,NULL));
415   PetscCall(PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,NULL));
416   PetscCall(PetscOptionsInt("-ts_sundials_maxord","Max Order for BDF/Adams method","TSSundialsSetMaxOrd",cvode->maxord,&cvode->maxord,NULL));
417   PetscCall(PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,NULL));
418   PetscCall(PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internal steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,NULL));
419   PetscCall(PetscOptionsBool("-ts_sundials_use_dense","Use dense internal solver in SUNDIALS (serial only)","TSSundialsSetUseDense",cvode->use_dense,&cvode->use_dense,NULL));
420   PetscOptionsHeadEnd();
421   PetscCall(TSSundialsGetPC(ts,&pc));
422   PetscCall(PCSetFromOptions(pc));
423   PetscFunctionReturn(0);
424 }
425 
426 PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
427 {
428   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
429   char           *type;
430   char           atype[] = "Adams";
431   char           btype[] = "BDF: backward differentiation formula";
432   PetscBool      iascii,isstring;
433   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
434   PetscInt       qlast,qcur;
435   PetscReal      hinused,hlast,hcur,tcur,tolsfac;
436   PC             pc;
437 
438   PetscFunctionBegin;
439   if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype;
440   else                                     type = btype;
441 
442   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
443   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring));
444   if (iascii) {
445     PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials integrator does not use SNES!\n"));
446     PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials integrator type %s\n",type));
447     PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials maxord %" PetscInt_FMT "\n",cvode->maxord));
448     PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",(double)cvode->abstol,(double)cvode->reltol));
449     if (cvode->use_dense) {
450       PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials integrator using a dense linear solve\n"));
451     } else {
452       PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",(double)cvode->linear_tol));
453       PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %" PetscInt_FMT "\n",cvode->maxl));
454       if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
455         PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n"));
456       } else {
457         PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n"));
458       }
459     }
460     if (cvode->mindt > 0) PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",(double)cvode->mindt));
461     if (cvode->maxdt > 0) PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",(double)cvode->maxdt));
462 
463     /* Outputs from CVODE, CVSPILS */
464     PetscCall(CVodeGetTolScaleFactor(cvode->mem,&tolsfac));
465     PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac));
466     PetscCall(CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,&nlinsetups,&nfails,&qlast,&qcur,&hinused,&hlast,&hcur,&tcur));
467     PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %ld\n",nsteps));
468     PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %ld\n",nfevals));
469     PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %ld\n",nlinsetups));
470     PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %ld\n",nfails));
471 
472     PetscCall(CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails));
473     PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %ld\n",its));
474     PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %ld\n",nfails));
475     if (!cvode->use_dense) {
476       PetscCall(CVSpilsGetNumLinIters(cvode->mem, &its)); /* its = no. of calls to TSPrecond_Sundials() */
477       PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %ld\n",its));
478       PetscCall(CVSpilsGetNumConvFails(cvode->mem,&itmp));
479       PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %ld\n",itmp));
480 
481       PetscCall(TSSundialsGetPC(ts,&pc));
482       PetscCall(PCView(pc,viewer));
483       PetscCall(CVSpilsGetNumPrecEvals(cvode->mem,&itmp));
484       PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %ld\n",itmp));
485       PetscCall(CVSpilsGetNumPrecSolves(cvode->mem,&itmp));
486       PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %ld\n",itmp));
487     }
488     PetscCall(CVSpilsGetNumJtimesEvals(cvode->mem,&itmp));
489     PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %ld\n",itmp));
490     PetscCall(CVSpilsGetNumRhsEvals(cvode->mem,&itmp));
491     PetscCall(PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %ld\n",itmp));
492   } else if (isstring) {
493     PetscCall(PetscViewerStringSPrintf(viewer,"Sundials type %s",type));
494   }
495   PetscFunctionReturn(0);
496 }
497 
498 /* --------------------------------------------------------------------------*/
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 
508 PetscErrorCode  TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl)
509 {
510   TS_Sundials *cvode = (TS_Sundials*)ts->data;
511 
512   PetscFunctionBegin;
513   cvode->maxl = maxl;
514   PetscFunctionReturn(0);
515 }
516 
517 PetscErrorCode  TSSundialsSetLinearTolerance_Sundials(TS ts,PetscReal tol)
518 {
519   TS_Sundials *cvode = (TS_Sundials*)ts->data;
520 
521   PetscFunctionBegin;
522   cvode->linear_tol = tol;
523   PetscFunctionReturn(0);
524 }
525 
526 PetscErrorCode  TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
527 {
528   TS_Sundials *cvode = (TS_Sundials*)ts->data;
529 
530   PetscFunctionBegin;
531   cvode->gtype = type;
532   PetscFunctionReturn(0);
533 }
534 
535 PetscErrorCode  TSSundialsSetTolerance_Sundials(TS ts,PetscReal aabs,PetscReal rel)
536 {
537   TS_Sundials *cvode = (TS_Sundials*)ts->data;
538 
539   PetscFunctionBegin;
540   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
541   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
542   PetscFunctionReturn(0);
543 }
544 
545 PetscErrorCode  TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)
546 {
547   TS_Sundials *cvode = (TS_Sundials*)ts->data;
548 
549   PetscFunctionBegin;
550   cvode->mindt = mindt;
551   PetscFunctionReturn(0);
552 }
553 
554 PetscErrorCode  TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)
555 {
556   TS_Sundials *cvode = (TS_Sundials*)ts->data;
557 
558   PetscFunctionBegin;
559   cvode->maxdt = maxdt;
560   PetscFunctionReturn(0);
561 }
562 
563 PetscErrorCode  TSSundialsSetUseDense_Sundials(TS ts,PetscBool use_dense)
564 {
565   TS_Sundials *cvode = (TS_Sundials*)ts->data;
566 
567   PetscFunctionBegin;
568   cvode->use_dense = use_dense;
569   PetscFunctionReturn(0);
570 }
571 
572 PetscErrorCode  TSSundialsGetPC_Sundials(TS ts,PC *pc)
573 {
574   SNES           snes;
575   KSP            ksp;
576 
577   PetscFunctionBegin;
578   PetscCall(TSGetSNES(ts,&snes));
579   PetscCall(SNESGetKSP(snes,&ksp));
580   PetscCall(KSPGetPC(ksp,pc));
581   PetscFunctionReturn(0);
582 }
583 
584 PetscErrorCode  TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
585 {
586   PetscFunctionBegin;
587   if (nonlin) *nonlin = ts->snes_its;
588   if (lin)    *lin    = ts->ksp_its;
589   PetscFunctionReturn(0);
590 }
591 
592 PetscErrorCode  TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s)
593 {
594   TS_Sundials *cvode = (TS_Sundials*)ts->data;
595 
596   PetscFunctionBegin;
597   cvode->monitorstep = s;
598   PetscFunctionReturn(0);
599 }
600 /* -------------------------------------------------------------------------------------------*/
601 
602 /*@C
603    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
604 
605    Not Collective
606 
607    Input Parameter:
608 .    ts     - the time-step context
609 
610    Output Parameters:
611 +   nonlin - number of nonlinear iterations
612 -   lin    - number of linear iterations
613 
614    Level: advanced
615 
616    Notes:
617     These return the number since the creation of the TS object
618 
619 .seealso: `TSSundialsSetType()`, `TSSundialsSetMaxl()`,
620           `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`,
621           `TSSundialsGetIterations()`, `TSSundialsSetType()`,
622           `TSSundialsSetLinearTolerance()`, `TSSundialsGetPC()`, `TSSetExactFinalTime()`
623 
624 @*/
625 PetscErrorCode  TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
626 {
627   PetscFunctionBegin;
628   PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));
629   PetscFunctionReturn(0);
630 }
631 
632 /*@
633    TSSundialsSetType - Sets the method that Sundials will use for integration.
634 
635    Logically Collective on TS
636 
637    Input Parameters:
638 +    ts     - the time-step context
639 -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF
640 
641    Level: intermediate
642 
643 .seealso: `TSSundialsGetIterations()`, `TSSundialsSetMaxl()`,
644           `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`,
645           `TSSundialsGetIterations()`, `TSSundialsSetType()`,
646           `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`,
647           `TSSetExactFinalTime()`
648 @*/
649 PetscErrorCode  TSSundialsSetType(TS ts,TSSundialsLmmType type)
650 {
651   PetscFunctionBegin;
652   PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));
653   PetscFunctionReturn(0);
654 }
655 
656 /*@
657    TSSundialsSetMaxord - Sets the maximum order for BDF/Adams method used by SUNDIALS.
658 
659    Logically Collective on TS
660 
661    Input Parameters:
662 +    ts      - the time-step context
663 -    maxord  - maximum order of BDF / Adams method
664 
665    Level: advanced
666 
667 .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`,
668           `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`,
669           `TSSundialsGetIterations()`, `TSSundialsSetType()`,
670           `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`,
671           `TSSetExactFinalTime()`
672 
673 @*/
674 PetscErrorCode  TSSundialsSetMaxord(TS ts,PetscInt maxord)
675 {
676   PetscFunctionBegin;
677   PetscValidLogicalCollectiveInt(ts,maxord,2);
678   PetscTryMethod(ts,"TSSundialsSetMaxOrd_C",(TS,PetscInt),(ts,maxord));
679   PetscFunctionReturn(0);
680 }
681 
682 /*@
683    TSSundialsSetMaxl - Sets the dimension of the Krylov space used by
684        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
685        this is the maximum number of GMRES steps that will be used.
686 
687    Logically Collective on TS
688 
689    Input Parameters:
690 +    ts      - the time-step context
691 -    maxl - number of direction vectors (the dimension of Krylov subspace).
692 
693    Level: advanced
694 
695 .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`,
696           `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`,
697           `TSSundialsGetIterations()`, `TSSundialsSetType()`,
698           `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`,
699           `TSSetExactFinalTime()`
700 
701 @*/
702 PetscErrorCode  TSSundialsSetMaxl(TS ts,PetscInt maxl)
703 {
704   PetscFunctionBegin;
705   PetscValidLogicalCollectiveInt(ts,maxl,2);
706   PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));
707   PetscFunctionReturn(0);
708 }
709 
710 /*@
711    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
712        system by SUNDIALS.
713 
714    Logically Collective on TS
715 
716    Input Parameters:
717 +    ts     - the time-step context
718 -    tol    - the factor by which the tolerance on the nonlinear solver is
719              multiplied to get the tolerance on the linear solver, .05 by default.
720 
721    Level: advanced
722 
723 .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`,
724           `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`,
725           `TSSundialsGetIterations()`, `TSSundialsSetType()`,
726           `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`,
727           `TSSetExactFinalTime()`
728 
729 @*/
730 PetscErrorCode  TSSundialsSetLinearTolerance(TS ts,PetscReal tol)
731 {
732   PetscFunctionBegin;
733   PetscValidLogicalCollectiveReal(ts,tol,2);
734   PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,PetscReal),(ts,tol));
735   PetscFunctionReturn(0);
736 }
737 
738 /*@
739    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
740         in GMRES method by SUNDIALS linear solver.
741 
742    Logically Collective on TS
743 
744    Input Parameters:
745 +    ts  - the time-step context
746 -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
747 
748    Level: advanced
749 
750 .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`,
751           `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`,
752           `TSSundialsGetIterations()`, `TSSundialsSetType()`,
753           `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`,
754           `TSSetExactFinalTime()`
755 
756 @*/
757 PetscErrorCode  TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
758 {
759   PetscFunctionBegin;
760   PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));
761   PetscFunctionReturn(0);
762 }
763 
764 /*@
765    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
766                          Sundials for error control.
767 
768    Logically Collective on TS
769 
770    Input Parameters:
771 +    ts  - the time-step context
772 .    aabs - the absolute tolerance
773 -    rel - the relative tolerance
774 
775      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
776     these regulate the size of the error for a SINGLE timestep.
777 
778    Level: intermediate
779 
780 .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetGMRESMaxl()`,
781           `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`,
782           `TSSundialsGetIterations()`, `TSSundialsSetType()`,
783           `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`,
784           `TSSetExactFinalTime()`
785 
786 @*/
787 PetscErrorCode  TSSundialsSetTolerance(TS ts,PetscReal aabs,PetscReal rel)
788 {
789   PetscFunctionBegin;
790   PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,PetscReal,PetscReal),(ts,aabs,rel));
791   PetscFunctionReturn(0);
792 }
793 
794 /*@
795    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
796 
797    Input Parameter:
798 .    ts - the time-step context
799 
800    Output Parameter:
801 .    pc - the preconditioner context
802 
803    Level: advanced
804 
805 .seealso: `TSSundialsGetIterations()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`,
806           `TSSundialsSetLinearTolerance()`, `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`,
807           `TSSundialsGetIterations()`, `TSSundialsSetType()`,
808           `TSSundialsSetLinearTolerance()`, `TSSundialsSetTolerance()`
809 @*/
810 PetscErrorCode  TSSundialsGetPC(TS ts,PC *pc)
811 {
812   PetscFunctionBegin;
813   PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC*),(ts,pc));
814   PetscFunctionReturn(0);
815 }
816 
817 /*@
818    TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller.
819 
820    Input Parameters:
821 +   ts - the time-step context
822 -   mindt - lowest time step if positive, negative to deactivate
823 
824    Note:
825    Sundials will error if it is not possible to keep the estimated truncation error below
826    the tolerance set with TSSundialsSetTolerance() without going below this step size.
827 
828    Level: beginner
829 
830 .seealso: `TSSundialsSetType()`, `TSSundialsSetTolerance()`,
831 @*/
832 PetscErrorCode  TSSundialsSetMinTimeStep(TS ts,PetscReal mindt)
833 {
834   PetscFunctionBegin;
835   PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));
836   PetscFunctionReturn(0);
837 }
838 
839 /*@
840    TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller.
841 
842    Input Parameters:
843 +   ts - the time-step context
844 -   maxdt - lowest time step if positive, negative to deactivate
845 
846    Level: beginner
847 
848 .seealso: `TSSundialsSetType()`, `TSSundialsSetTolerance()`,
849 @*/
850 PetscErrorCode  TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt)
851 {
852   PetscFunctionBegin;
853   PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));
854   PetscFunctionReturn(0);
855 }
856 
857 /*@
858    TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false).
859 
860    Input Parameters:
861 +   ts - the time-step context
862 -   ft - PETSC_TRUE if monitor, else PETSC_FALSE
863 
864    Level: beginner
865 
866 .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
867           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
868           TSSundialsGetIterations(), TSSundialsSetType(),
869           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
870 @*/
871 PetscErrorCode  TSSundialsMonitorInternalSteps(TS ts,PetscBool ft)
872 {
873   PetscFunctionBegin;
874   PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));
875   PetscFunctionReturn(0);
876 }
877 
878 /*@
879    TSSundialsSetUseDense - Set a flag to use a dense linear solver in SUNDIALS (serial only)
880 
881    Logically Collective
882 
883    Input Parameters:
884 +    ts         - the time-step context
885 -    use_dense  - PETSC_TRUE to use the dense solver
886 
887    Level: advanced
888 
889 .seealso: `TSSUNDIALS`
890 
891 @*/
892 PetscErrorCode  TSSundialsSetUseDense(TS ts,PetscBool use_dense)
893 {
894   PetscFunctionBegin;
895   PetscValidLogicalCollectiveInt(ts,use_dense,2);
896   PetscTryMethod(ts,"TSSundialsSetUseDense_C",(TS,PetscBool),(ts,use_dense));
897   PetscFunctionReturn(0);
898 }
899 
900 /* -------------------------------------------------------------------------------------------*/
901 /*MC
902       TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
903 
904    Options Database:
905 +    -ts_sundials_type <bdf,adams> -
906 .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
907 .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
908 .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
909 .    -ts_sundials_linear_tolerance <tol> -
910 .    -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace
911 .    -ts_sundials_monitor_steps - Monitor SUNDIALS internal steps
912 -    -ts_sundials_use_dense - Use a dense linear solver within CVODE (serial only)
913 
914     Notes:
915     This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply,
916            only PETSc PC options.
917 
918     Level: beginner
919 
920 .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSSundialsSetType()`, `TSSundialsSetMaxl()`, `TSSundialsSetLinearTolerance()`,
921           `TSSundialsSetGramSchmidtType()`, `TSSundialsSetTolerance()`, `TSSundialsGetPC()`, `TSSundialsGetIterations()`, `TSSetExactFinalTime()`
922 
923 M*/
924 PETSC_EXTERN PetscErrorCode TSCreate_Sundials(TS ts)
925 {
926   TS_Sundials    *cvode;
927   PC             pc;
928 
929   PetscFunctionBegin;
930   ts->ops->reset          = TSReset_Sundials;
931   ts->ops->destroy        = TSDestroy_Sundials;
932   ts->ops->view           = TSView_Sundials;
933   ts->ops->setup          = TSSetUp_Sundials;
934   ts->ops->step           = TSStep_Sundials;
935   ts->ops->interpolate    = TSInterpolate_Sundials;
936   ts->ops->setfromoptions = TSSetFromOptions_Sundials;
937   ts->default_adapt_type  = TSADAPTNONE;
938 
939   PetscCall(PetscNewLog(ts,&cvode));
940 
941   ts->usessnes = PETSC_TRUE;
942 
943   ts->data           = (void*)cvode;
944   cvode->cvode_type  = SUNDIALS_BDF;
945   cvode->gtype       = SUNDIALS_CLASSICAL_GS;
946   cvode->maxl        = 5;
947   cvode->maxord      = PETSC_DEFAULT;
948   cvode->linear_tol  = .05;
949   cvode->monitorstep = PETSC_TRUE;
950   cvode->use_dense   = PETSC_FALSE;
951 
952   PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)ts),&(cvode->comm_sundials)));
953 
954   cvode->mindt = -1.;
955   cvode->maxdt = -1.;
956 
957   /* set tolerance for Sundials */
958   cvode->reltol = 1e-6;
959   cvode->abstol = 1e-6;
960 
961   /* set PCNONE as default pctype */
962   PetscCall(TSSundialsGetPC_Sundials(ts,&pc));
963   PetscCall(PCSetType(pc,PCNONE));
964 
965   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",TSSundialsSetType_Sundials));
966   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",TSSundialsSetMaxl_Sundials));
967   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",TSSundialsSetLinearTolerance_Sundials));
968   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",TSSundialsSetGramSchmidtType_Sundials));
969   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",TSSundialsSetTolerance_Sundials));
970   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",TSSundialsSetMinTimeStep_Sundials));
971   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",TSSundialsSetMaxTimeStep_Sundials));
972   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",TSSundialsGetPC_Sundials));
973   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",TSSundialsGetIterations_Sundials));
974   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",TSSundialsMonitorInternalSteps_Sundials));
975   PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetUseDense_C",TSSundialsSetUseDense_Sundials));
976   PetscFunctionReturn(0);
977 }
978