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