Actual source code: sundials.c

petsc-3.12.5 2020-03-29
Report Typos and Errors
  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.

  6:     Reference: sundials-2.4.0/examples/cvode/parallel/cvDiurnal_kry_p.c
  7: */
  8:  #include <../src/ts/impls/implicit/sundials/sundials.h>

 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;
 21:   Mat            J,P;
 22:   Vec            yy  = cvode->w1,yydot = cvode->ydot;
 23:   PetscReal      gm  = (PetscReal)_gamma;
 24:   PetscScalar    *y_data;

 27:   TSGetIJacobian(ts,&J,&P,NULL,NULL);
 28:   y_data = (PetscScalar*) N_VGetArrayPointer(y);
 29:   VecPlaceArray(yy,y_data);
 30:   VecZeroEntries(yydot); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */
 31:   /* compute the shifted Jacobian   (1/gm)*I + Jrest */
 32:   TSComputeIJacobian(ts,ts->ptime,yy,yydot,1/gm,J,P,PETSC_FALSE);
 33:   VecResetArray(yy);
 34:   MatScale(P,gm); /* turn into I-gm*Jrest, J is not used by Sundials  */
 35:   *jcurPtr = TRUE;
 36:   TSSundialsGetPC(ts,&pc);
 37:   PCSetOperators(pc,J,P);
 38:   return(0);
 39: }

 41: /*
 42:      TSPSolve_Sundials -  routine that we provide to Sundials that applies the preconditioner.
 43: */
 44: PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z,
 45:                                  realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp)
 46: {
 47:   TS             ts     = (TS) P_data;
 48:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
 49:   PC             pc;
 50:   Vec            rr = cvode->w1,zz = cvode->w2;
 52:   PetscScalar    *r_data,*z_data;

 55:   /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
 56:   r_data = (PetscScalar*) N_VGetArrayPointer(r);
 57:   z_data = (PetscScalar*) N_VGetArrayPointer(z);
 58:   VecPlaceArray(rr,r_data);
 59:   VecPlaceArray(zz,z_data);

 61:   /* Solve the Px=r and put the result in zz */
 62:   TSSundialsGetPC(ts,&pc);
 63:   PCApply(pc,rr,zz);
 64:   VecResetArray(rr);
 65:   VecResetArray(zz);
 66:   return(0);
 67: }

 69: /*
 70:         TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
 71: */
 72: int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx)
 73: {
 74:   TS             ts = (TS) ctx;
 75:   DM             dm;
 76:   DMTS           tsdm;
 77:   TSIFunction    ifunction;
 78:   MPI_Comm       comm;
 79:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
 80:   Vec            yy     = cvode->w1,yyd = cvode->w2,yydot = cvode->ydot;
 81:   PetscScalar    *y_data,*ydot_data;

 85:   PetscObjectGetComm((PetscObject)ts,&comm);
 86:   /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
 87:   y_data    = (PetscScalar*) N_VGetArrayPointer(y);
 88:   ydot_data = (PetscScalar*) N_VGetArrayPointer(ydot);
 89:   VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr);
 90:   VecPlaceArray(yyd,ydot_data);CHKERRABORT(comm,ierr);

 92:   /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */
 93:   TSGetDM(ts,&dm);
 94:   DMGetDMTS(dm,&tsdm);
 95:   DMTSGetIFunction(dm,&ifunction,NULL);
 96:   if (!ifunction) {
 97:     TSComputeRHSFunction(ts,t,yy,yyd);
 98:   } else {                      /* If rhsfunction is also set, this computes both parts and shifts them to the right */
 99:     VecZeroEntries(yydot);
100:     TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE);CHKERRABORT(comm,ierr);
101:     VecScale(yyd,-1.);
102:   }
103:   VecResetArray(yy);CHKERRABORT(comm,ierr);
104:   VecResetArray(yyd);CHKERRABORT(comm,ierr);
105:   return(0);
106: }

108: /*
109:        TSStep_Sundials - Calls Sundials to integrate the ODE.
110: */
111: PetscErrorCode TSStep_Sundials(TS ts)
112: {
113:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
115:   PetscInt       flag;
116:   long int       nits,lits,nsteps;
117:   realtype       t,tout;
118:   PetscScalar    *y_data;
119:   void           *mem;

122:   mem  = cvode->mem;
123:   tout = ts->max_time;
124:   VecGetArray(ts->vec_sol,&y_data);
125:   N_VSetArrayPointer((realtype*)y_data,cvode->y);
126:   VecRestoreArray(ts->vec_sol,NULL);

128:   /* We would like to TSPreStage() and TSPostStage()
129:    * before each stage solve but CVode does not appear to support this. */
130:   if (cvode->monitorstep)
131:     flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
132:   else
133:     flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL);

135:   if (flag) { /* display error message */
136:     switch (flag) {
137:       case CV_ILL_INPUT:
138:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT");
139:         break;
140:       case CV_TOO_CLOSE:
141:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE");
142:         break;
143:       case CV_TOO_MUCH_WORK: {
144:         PetscReal tcur;
145:         CVodeGetNumSteps(mem,&nsteps);
146:         CVodeGetCurrentTime(mem,&tcur);
147:         SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_WORK. At t=%g, nsteps %D exceeds maxstep %D. Increase '-ts_max_steps <>' or modify TSSetMaxSteps()",(double)tcur,nsteps,ts->max_steps);
148:       } break;
149:       case CV_TOO_MUCH_ACC:
150:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC");
151:         break;
152:       case CV_ERR_FAILURE:
153:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE");
154:         break;
155:       case CV_CONV_FAILURE:
156:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE");
157:         break;
158:       case CV_LINIT_FAIL:
159:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL");
160:         break;
161:       case CV_LSETUP_FAIL:
162:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL");
163:         break;
164:       case CV_LSOLVE_FAIL:
165:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL");
166:         break;
167:       case CV_RHSFUNC_FAIL:
168:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL");
169:         break;
170:       case CV_FIRST_RHSFUNC_ERR:
171:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR");
172:         break;
173:       case CV_REPTD_RHSFUNC_ERR:
174:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR");
175:         break;
176:       case CV_UNREC_RHSFUNC_ERR:
177:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR");
178:         break;
179:       case CV_RTFUNC_FAIL:
180:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL");
181:         break;
182:       default:
183:         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag);
184:     }
185:   }

187:   /* log inner nonlinear and linear iterations */
188:   CVodeGetNumNonlinSolvIters(mem,&nits);
189:   CVSpilsGetNumLinIters(mem,&lits);
190:   ts->snes_its += nits; ts->ksp_its = lits;

192:   /* copy the solution from cvode->y to cvode->update and sol */
193:   VecPlaceArray(cvode->w1,y_data);
194:   VecCopy(cvode->w1,cvode->update);
195:   VecResetArray(cvode->w1);
196:   VecCopy(cvode->update,ts->vec_sol);

198:   ts->time_step = t - ts->ptime;
199:   ts->ptime = t;

201:   CVodeGetNumSteps(mem,&nsteps);
202:   if (!cvode->monitorstep) ts->steps += nsteps - 1; /* TSStep() increments the step counter by one */
203:   return(0);
204: }

206: static PetscErrorCode TSInterpolate_Sundials(TS ts,PetscReal t,Vec X)
207: {
208:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
209:   N_Vector       y;
211:   PetscScalar    *x_data;
212:   PetscInt       glosize,locsize;

215:   /* get the vector size */
216:   VecGetSize(X,&glosize);
217:   VecGetLocalSize(X,&locsize);
218:   VecGetArray(X,&x_data);

220:   /* Initialize N_Vec y with x_data */
221:   y = N_VMake_Parallel(cvode->comm_sundials,locsize,glosize,(realtype*)x_data);
222:   if (!y) SETERRQ(PETSC_COMM_SELF,1,"Interpolated y is not allocated");

224:   CVodeGetDky(cvode->mem,t,0,y);
225:   VecRestoreArray(X,&x_data);
226:   return(0);
227: }

229: PetscErrorCode TSReset_Sundials(TS ts)
230: {
231:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;

235:   VecDestroy(&cvode->update);
236:   VecDestroy(&cvode->ydot);
237:   VecDestroy(&cvode->w1);
238:   VecDestroy(&cvode->w2);
239:   if (cvode->mem) CVodeFree(&cvode->mem);
240:   return(0);
241: }

243: PetscErrorCode TSDestroy_Sundials(TS ts)
244: {
245:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;

249:   TSReset_Sundials(ts);
250:   MPI_Comm_free(&(cvode->comm_sundials));
251:   PetscFree(ts->data);
252:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",NULL);
253:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",NULL);
254:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",NULL);
255:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",NULL);
256:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",NULL);
257:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",NULL);
258:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",NULL);
259:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",NULL);
260:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",NULL);
261:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",NULL);
262:   return(0);
263: }

265: PetscErrorCode TSSetUp_Sundials(TS ts)
266: {
267:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
269:   PetscInt       glosize,locsize,i,flag;
270:   PetscScalar    *y_data,*parray;
271:   void           *mem;
272:   PC             pc;
273:   PCType         pctype;
274:   PetscBool      pcnone;

277:   if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for exact final time option 'MATCHSTEP' when using Sundials");

279:   /* get the vector size */
280:   VecGetSize(ts->vec_sol,&glosize);
281:   VecGetLocalSize(ts->vec_sol,&locsize);

283:   /* allocate the memory for N_Vec y */
284:   cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
285:   if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated");

287:   /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
288:   VecGetArray(ts->vec_sol,&parray);
289:   y_data = (PetscScalar*) N_VGetArrayPointer(cvode->y);
290:   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
291:   VecRestoreArray(ts->vec_sol,NULL);

293:   VecDuplicate(ts->vec_sol,&cvode->update);
294:   VecDuplicate(ts->vec_sol,&cvode->ydot);
295:   PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->update);
296:   PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->ydot);

298:   /*
299:     Create work vectors for the TSPSolve_Sundials() routine. Note these are
300:     allocated with zero space arrays because the actual array space is provided
301:     by Sundials and set using VecPlaceArray().
302:   */
303:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,0,&cvode->w1);
304:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,0,&cvode->w2);
305:   PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w1);
306:   PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w2);

308:   /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
309:   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
310:   if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails");
311:   cvode->mem = mem;

313:   /* Set the pointer to user-defined data */
314:   flag = CVodeSetUserData(mem, ts);
315:   if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails");

317:   /* Sundials may choose to use a smaller initial step, but will never use a larger step. */
318:   flag = CVodeSetInitStep(mem,(realtype)ts->time_step);
319:   if (flag) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetInitStep() failed");
320:   if (cvode->mindt > 0) {
321:     flag = CVodeSetMinStep(mem,(realtype)cvode->mindt);
322:     if (flag) {
323:       if (flag == CV_MEM_NULL) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL");
324:       else if (flag == CV_ILL_INPUT) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size");
325:       else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed");
326:     }
327:   }
328:   if (cvode->maxdt > 0) {
329:     flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt);
330:     if (flag) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMaxStep() failed");
331:   }

333:   /* Call CVodeInit to initialize the integrator memory and specify the
334:    * user's right hand side function in u'=f(t,u), the inital time T0, and
335:    * the initial dependent variable vector cvode->y */
336:   flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y);
337:   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag);

339:   /* specifies scalar relative and absolute tolerances */
340:   flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol);
341:   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag);

343:   /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */
344:   flag = CVodeSetMaxNumSteps(mem,ts->max_steps);
345:   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetMaxNumSteps() fails, flag %d",flag);

347:   /* call CVSpgmr to use GMRES as the linear solver.        */
348:   /* setup the ode integrator with the given preconditioner */
349:   TSSundialsGetPC(ts,&pc);
350:   PCGetType(pc,&pctype);
351:   PetscObjectTypeCompare((PetscObject)pc,PCNONE,&pcnone);
352:   if (pcnone) {
353:     flag = CVSpgmr(mem,PREC_NONE,0);
354:     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
355:   } else {
356:     flag = CVSpgmr(mem,PREC_LEFT,cvode->maxl);
357:     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);

359:     /* Set preconditioner and solve routines Precond and PSolve,
360:      and the pointer to the user-defined block data */
361:     flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials);
362:     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag);
363:   }

365:   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
366:   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag);
367:   return(0);
368: }

370: /* type of CVODE linear multistep method */
371: const char *const TSSundialsLmmTypes[] = {"","ADAMS","BDF","TSSundialsLmmType","SUNDIALS_",0};
372: /* type of G-S orthogonalization used by CVODE linear solver */
373: const char *const TSSundialsGramSchmidtTypes[] = {"","MODIFIED","CLASSICAL","TSSundialsGramSchmidtType","SUNDIALS_",0};

375: PetscErrorCode TSSetFromOptions_Sundials(PetscOptionItems *PetscOptionsObject,TS ts)
376: {
377:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
379:   int            indx;
380:   PetscBool      flag;
381:   PC             pc;

384:   PetscOptionsHead(PetscOptionsObject,"SUNDIALS ODE solver options");
385:   PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);
386:   if (flag) {
387:     TSSundialsSetType(ts,(TSSundialsLmmType)indx);
388:   }
389:   PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);
390:   if (flag) {
391:     TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);
392:   }
393:   PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,NULL);
394:   PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,NULL);
395:   PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,NULL);
396:   PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,NULL);
397:   PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,NULL);
398:   PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,NULL);
399:   PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internal steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,NULL);
400:   PetscOptionsTail();
401:   TSSundialsGetPC(ts,&pc);
402:   PCSetFromOptions(pc);
403:   return(0);
404: }

406: PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
407: {
408:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
410:   char           *type;
411:   char           atype[] = "Adams";
412:   char           btype[] = "BDF: backward differentiation formula";
413:   PetscBool      iascii,isstring;
414:   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
415:   PetscInt       qlast,qcur;
416:   PetscReal      hinused,hlast,hcur,tcur,tolsfac;
417:   PC             pc;

420:   if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype;
421:   else                                     type = btype;

423:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
424:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
425:   if (iascii) {
426:     PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");
427:     PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);
428:     PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);
429:     PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);
430:     PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %D\n",cvode->maxl);
431:     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
432:       PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");
433:     } else {
434:       PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");
435:     }
436:     if (cvode->mindt > 0) {PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);}
437:     if (cvode->maxdt > 0) {PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);}

439:     /* Outputs from CVODE, CVSPILS */
440:     CVodeGetTolScaleFactor(cvode->mem,&tolsfac);
441:     PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);
442:     CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
443:                                    &nlinsetups,&nfails,&qlast,&qcur,
444:                                    &hinused,&hlast,&hcur,&tcur);
445:     PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);
446:     PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);
447:     PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);
448:     PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);

450:     CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);
451:     PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);
452:     PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);

454:     CVSpilsGetNumLinIters(cvode->mem, &its); /* its = no. of calls to TSPrecond_Sundials() */
455:     PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);
456:     CVSpilsGetNumConvFails(cvode->mem,&itmp);
457:     PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);

459:     TSSundialsGetPC(ts,&pc);
460:     PCView(pc,viewer);
461:     CVSpilsGetNumPrecEvals(cvode->mem,&itmp);
462:     PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);
463:     CVSpilsGetNumPrecSolves(cvode->mem,&itmp);
464:     PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);

466:     CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);
467:     PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);
468:     CVSpilsGetNumRhsEvals(cvode->mem,&itmp);
469:     PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);
470:   } else if (isstring) {
471:     PetscViewerStringSPrintf(viewer,"Sundials type %s",type);
472:   }
473:   return(0);
474: }


477: /* --------------------------------------------------------------------------*/
478: PetscErrorCode  TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
479: {
480:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

483:   cvode->cvode_type = type;
484:   return(0);
485: }

487: PetscErrorCode  TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl)
488: {
489:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

492:   cvode->maxl = maxl;
493:   return(0);
494: }

496: PetscErrorCode  TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
497: {
498:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

501:   cvode->linear_tol = tol;
502:   return(0);
503: }

505: PetscErrorCode  TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
506: {
507:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

510:   cvode->gtype = type;
511:   return(0);
512: }

514: PetscErrorCode  TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
515: {
516:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

519:   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
520:   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
521:   return(0);
522: }

524: PetscErrorCode  TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)
525: {
526:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

529:   cvode->mindt = mindt;
530:   return(0);
531: }

533: PetscErrorCode  TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)
534: {
535:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

538:   cvode->maxdt = maxdt;
539:   return(0);
540: }
541: PetscErrorCode  TSSundialsGetPC_Sundials(TS ts,PC *pc)
542: {
543:   SNES           snes;
544:   KSP            ksp;

548:   TSGetSNES(ts,&snes);
549:   SNESGetKSP(snes,&ksp);
550:   KSPGetPC(ksp,pc);
551:   return(0);
552: }

554: PetscErrorCode  TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
555: {
557:   if (nonlin) *nonlin = ts->snes_its;
558:   if (lin)    *lin    = ts->ksp_its;
559:   return(0);
560: }

562: PetscErrorCode  TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s)
563: {
564:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

567:   cvode->monitorstep = s;
568:   return(0);
569: }
570: /* -------------------------------------------------------------------------------------------*/

572: /*@C
573:    TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.

575:    Not Collective

577:    Input parameters:
578: .    ts     - the time-step context

580:    Output Parameters:
581: +   nonlin - number of nonlinear iterations
582: -   lin    - number of linear iterations

584:    Level: advanced

586:    Notes:
587:     These return the number since the creation of the TS object

589: .seealso: TSSundialsSetType(), TSSundialsSetMaxl(),
590:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
591:           TSSundialsGetIterations(), TSSundialsSetType(),
592:           TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime()

594: @*/
595: PetscErrorCode  TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
596: {

600:   PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));
601:   return(0);
602: }

604: /*@
605:    TSSundialsSetType - Sets the method that Sundials will use for integration.

607:    Logically Collective on TS

609:    Input parameters:
610: +    ts     - the time-step context
611: -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF

613:    Level: intermediate

615: .seealso: TSSundialsGetIterations(),  TSSundialsSetMaxl(),
616:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
617:           TSSundialsGetIterations(), TSSundialsSetType(),
618:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
619:           TSSetExactFinalTime()
620: @*/
621: PetscErrorCode  TSSundialsSetType(TS ts,TSSundialsLmmType type)
622: {

626:   PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));
627:   return(0);
628: }

630: /*@
631:    TSSundialsSetMaxl - Sets the dimension of the Krylov space used by
632:        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
633:        this is the maximum number of GMRES steps that will be used.

635:    Logically Collective on TS

637:    Input parameters:
638: +    ts      - the time-step context
639: -    maxl - number of direction vectors (the dimension of Krylov subspace).

641:    Level: advanced

643: .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
644:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
645:           TSSundialsGetIterations(), TSSundialsSetType(),
646:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
647:           TSSetExactFinalTime()

649: @*/
650: PetscErrorCode  TSSundialsSetMaxl(TS ts,PetscInt maxl)
651: {

656:   PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));
657:   return(0);
658: }

660: /*@
661:    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
662:        system by SUNDIALS.

664:    Logically Collective on TS

666:    Input parameters:
667: +    ts     - the time-step context
668: -    tol    - the factor by which the tolerance on the nonlinear solver is
669:              multiplied to get the tolerance on the linear solver, .05 by default.

671:    Level: advanced

673: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
674:           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
675:           TSSundialsGetIterations(), TSSundialsSetType(),
676:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
677:           TSSetExactFinalTime()

679: @*/
680: PetscErrorCode  TSSundialsSetLinearTolerance(TS ts,double tol)
681: {

686:   PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));
687:   return(0);
688: }

690: /*@
691:    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
692:         in GMRES method by SUNDIALS linear solver.

694:    Logically Collective on TS

696:    Input parameters:
697: +    ts  - the time-step context
698: -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS

700:    Level: advanced

702: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
703:           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
704:           TSSundialsGetIterations(), TSSundialsSetType(),
705:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
706:           TSSetExactFinalTime()

708: @*/
709: PetscErrorCode  TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
710: {

714:   PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));
715:   return(0);
716: }

718: /*@
719:    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
720:                          Sundials for error control.

722:    Logically Collective on TS

724:    Input parameters:
725: +    ts  - the time-step context
726: .    aabs - the absolute tolerance
727: -    rel - the relative tolerance

729:      See the Cvode/Sundials users manual for exact details on these parameters. Essentially
730:     these regulate the size of the error for a SINGLE timestep.

732:    Level: intermediate

734: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESMaxl(),
735:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
736:           TSSundialsGetIterations(), TSSundialsSetType(),
737:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
738:           TSSetExactFinalTime()

740: @*/
741: PetscErrorCode  TSSundialsSetTolerance(TS ts,double aabs,double rel)
742: {

746:   PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));
747:   return(0);
748: }

750: /*@
751:    TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.

753:    Input Parameter:
754: .    ts - the time-step context

756:    Output Parameter:
757: .    pc - the preconditioner context

759:    Level: advanced

761: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
762:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
763:           TSSundialsGetIterations(), TSSundialsSetType(),
764:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
765: @*/
766: PetscErrorCode  TSSundialsGetPC(TS ts,PC *pc)
767: {

771:   PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC*),(ts,pc));
772:   return(0);
773: }

775: /*@
776:    TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller.

778:    Input Parameter:
779: +   ts - the time-step context
780: -   mindt - lowest time step if positive, negative to deactivate

782:    Note:
783:    Sundials will error if it is not possible to keep the estimated truncation error below
784:    the tolerance set with TSSundialsSetTolerance() without going below this step size.

786:    Level: beginner

788: .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
789: @*/
790: PetscErrorCode  TSSundialsSetMinTimeStep(TS ts,PetscReal mindt)
791: {

795:   PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));
796:   return(0);
797: }

799: /*@
800:    TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller.

802:    Input Parameter:
803: +   ts - the time-step context
804: -   maxdt - lowest time step if positive, negative to deactivate

806:    Level: beginner

808: .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
809: @*/
810: PetscErrorCode  TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt)
811: {

815:   PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));
816:   return(0);
817: }

819: /*@
820:    TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false).

822:    Input Parameter:
823: +   ts - the time-step context
824: -   ft - PETSC_TRUE if monitor, else PETSC_FALSE

826:    Level: beginner

828: .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
829:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
830:           TSSundialsGetIterations(), TSSundialsSetType(),
831:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
832: @*/
833: PetscErrorCode  TSSundialsMonitorInternalSteps(TS ts,PetscBool ft)
834: {

838:   PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));
839:   return(0);
840: }
841: /* -------------------------------------------------------------------------------------------*/
842: /*MC
843:       TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)

845:    Options Database:
846: +    -ts_sundials_type <bdf,adams> -
847: .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
848: .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
849: .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
850: .    -ts_sundials_linear_tolerance <tol> - 
851: .    -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace
852: -    -ts_sundials_monitor_steps - Monitor SUNDIALS internal steps


855:     Notes:
856:     This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply,
857:            only PETSc PC options.

859:     Level: beginner

861: .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetMaxl(), TSSundialsSetLinearTolerance(),
862:            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime()

864: M*/
865: PETSC_EXTERN PetscErrorCode TSCreate_Sundials(TS ts)
866: {
867:   TS_Sundials    *cvode;
869:   PC             pc;

872:   ts->ops->reset          = TSReset_Sundials;
873:   ts->ops->destroy        = TSDestroy_Sundials;
874:   ts->ops->view           = TSView_Sundials;
875:   ts->ops->setup          = TSSetUp_Sundials;
876:   ts->ops->step           = TSStep_Sundials;
877:   ts->ops->interpolate    = TSInterpolate_Sundials;
878:   ts->ops->setfromoptions = TSSetFromOptions_Sundials;
879:   ts->default_adapt_type  = TSADAPTNONE;

881:   PetscNewLog(ts,&cvode);

883:   ts->usessnes = PETSC_TRUE;

885:   ts->data           = (void*)cvode;
886:   cvode->cvode_type  = SUNDIALS_BDF;
887:   cvode->gtype       = SUNDIALS_CLASSICAL_GS;
888:   cvode->maxl        = 5;
889:   cvode->linear_tol  = .05;
890:   cvode->monitorstep = PETSC_TRUE;

892:   MPI_Comm_dup(PetscObjectComm((PetscObject)ts),&(cvode->comm_sundials));

894:   cvode->mindt = -1.;
895:   cvode->maxdt = -1.;

897:   /* set tolerance for Sundials */
898:   cvode->reltol = 1e-6;
899:   cvode->abstol = 1e-6;

901:   /* set PCNONE as default pctype */
902:   TSSundialsGetPC_Sundials(ts,&pc);
903:   PCSetType(pc,PCNONE);

905:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",TSSundialsSetType_Sundials);
906:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",TSSundialsSetMaxl_Sundials);
907:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",TSSundialsSetLinearTolerance_Sundials);
908:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",TSSundialsSetGramSchmidtType_Sundials);
909:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",TSSundialsSetTolerance_Sundials);
910:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",TSSundialsSetMinTimeStep_Sundials);
911:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",TSSundialsSetMaxTimeStep_Sundials);
912:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",TSSundialsGetPC_Sundials);
913:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",TSSundialsGetIterations_Sundials);
914:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",TSSundialsMonitorInternalSteps_Sundials);
915:   return(0);
916: }