Actual source code: sundials.c

petsc-3.6.1 2015-08-06
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>  /*I "petscts.h" I*/

 10: /*
 11:       TSPrecond_Sundials - function that we provide to SUNDIALS to
 12:                         evaluate the preconditioner.
 13: */
 16: PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy,booleantype jok,booleantype *jcurPtr,
 17:                                   realtype _gamma,void *P_data,N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3)
 18: {
 19:   TS             ts     = (TS) P_data;
 20:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
 21:   PC             pc;
 23:   Mat            J,P;
 24:   Vec            yy  = cvode->w1,yydot = cvode->ydot;
 25:   PetscReal      gm  = (PetscReal)_gamma;
 26:   PetscScalar    *y_data;

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

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

 59:   /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
 60:   r_data = (PetscScalar*) N_VGetArrayPointer(r);
 61:   z_data = (PetscScalar*) N_VGetArrayPointer(z);
 62:   VecPlaceArray(rr,r_data);
 63:   VecPlaceArray(zz,z_data);

 65:   /* Solve the Px=r and put the result in zz */
 66:   TSSundialsGetPC(ts,&pc);
 67:   PCApply(pc,rr,zz);
 68:   VecResetArray(rr);
 69:   VecResetArray(zz);
 70:   return(0);
 71: }

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

 91:   PetscObjectGetComm((PetscObject)ts,&comm);
 92:   /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
 93:   y_data    = (PetscScalar*) N_VGetArrayPointer(y);
 94:   ydot_data = (PetscScalar*) N_VGetArrayPointer(ydot);
 95:   VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr);
 96:   VecPlaceArray(yyd,ydot_data);CHKERRABORT(comm,ierr);

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

114: /*
115:        TSStep_Sundials - Calls Sundials to integrate the ODE.
116: */
119: PetscErrorCode TSStep_Sundials(TS ts)
120: {
121:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
123:   PetscInt       flag;
124:   long int       its,nsteps;
125:   realtype       t,tout;
126:   PetscScalar    *y_data;
127:   void           *mem;

130:   mem  = cvode->mem;
131:   tout = ts->max_time;
132:   VecGetArray(ts->vec_sol,&y_data);
133:   N_VSetArrayPointer((realtype*)y_data,cvode->y);
134:   VecRestoreArray(ts->vec_sol,NULL);

136:   TSPreStep(ts);

138:   /* We would like to call TSPreStep() when starting each step (including rejections), TSPreStage(),
139:    * and TSPostStage() before each stage solve, but CVode does not appear to support this. */
140:   if (cvode->monitorstep) flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
141:   else flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL);

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

195:   /* copy the solution from cvode->y to cvode->update and sol */
196:   VecPlaceArray(cvode->w1,y_data);
197:   VecCopy(cvode->w1,cvode->update);
198:   VecResetArray(cvode->w1);
199:   VecCopy(cvode->update,ts->vec_sol);
200:   CVodeGetNumNonlinSolvIters(mem,&its);
201:   CVSpilsGetNumLinIters(mem, &its);
202:   ts->snes_its = its; ts->ksp_its = its;

204:   ts->time_step = t - ts->ptime;
205:   ts->ptime     = t;
206:   ts->steps++;

208:   CVodeGetNumSteps(mem,&nsteps);
209:   if (!cvode->monitorstep) ts->steps = nsteps;
210:   return(0);
211: }

215: static PetscErrorCode TSInterpolate_Sundials(TS ts,PetscReal t,Vec X)
216: {
217:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
218:   N_Vector       y;
220:   PetscScalar    *x_data;
221:   PetscInt       glosize,locsize;

224:   /* get the vector size */
225:   VecGetSize(X,&glosize);
226:   VecGetLocalSize(X,&locsize);

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

232:   VecGetArray(X,&x_data);
233:   N_VSetArrayPointer((realtype*)x_data,y);
234:   CVodeGetDky(cvode->mem,t,0,y);
235:   VecRestoreArray(X,&x_data);
236:   return(0);
237: }

241: PetscErrorCode TSReset_Sundials(TS ts)
242: {
243:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;

247:   VecDestroy(&cvode->update);
248:   VecDestroy(&cvode->ydot);
249:   VecDestroy(&cvode->w1);
250:   VecDestroy(&cvode->w2);
251:   if (cvode->mem) CVodeFree(&cvode->mem);
252:   return(0);
253: }

257: PetscErrorCode TSDestroy_Sundials(TS ts)
258: {
259:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;

263:   TSReset_Sundials(ts);
264:   MPI_Comm_free(&(cvode->comm_sundials));
265:   PetscFree(ts->data);
266:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",NULL);
267:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",NULL);
268:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",NULL);
269:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",NULL);
270:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",NULL);
271:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",NULL);
272:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",NULL);
273:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",NULL);
274:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",NULL);
275:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",NULL);
276:   return(0);
277: }

281: PetscErrorCode TSSetUp_Sundials(TS ts)
282: {
283:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
285:   PetscInt       glosize,locsize,i,flag;
286:   PetscScalar    *y_data,*parray;
287:   void           *mem;
288:   PC             pc;
289:   PCType         pctype;
290:   PetscBool      pcnone;

293:   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");

295:   /* get the vector size */
296:   VecGetSize(ts->vec_sol,&glosize);
297:   VecGetLocalSize(ts->vec_sol,&locsize);

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

303:   /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
304:   VecGetArray(ts->vec_sol,&parray);
305:   y_data = (PetscScalar*) N_VGetArrayPointer(cvode->y);
306:   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
307:   VecRestoreArray(ts->vec_sol,NULL);

309:   VecDuplicate(ts->vec_sol,&cvode->update);
310:   VecDuplicate(ts->vec_sol,&cvode->ydot);
311:   PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->update);
312:   PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->ydot);

314:   /*
315:     Create work vectors for the TSPSolve_Sundials() routine. Note these are
316:     allocated with zero space arrays because the actual array space is provided
317:     by Sundials and set using VecPlaceArray().
318:   */
319:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,0,&cvode->w1);
320:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,0,&cvode->w2);
321:   PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w1);
322:   PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w2);

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

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

333:   /* Sundials may choose to use a smaller initial step, but will never use a larger step. */
334:   flag = CVodeSetInitStep(mem,(realtype)ts->time_step);
335:   if (flag) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetInitStep() failed");
336:   if (cvode->mindt > 0) {
337:     flag = CVodeSetMinStep(mem,(realtype)cvode->mindt);
338:     if (flag) {
339:       if (flag == CV_MEM_NULL) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL");
340:       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");
341:       else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed");
342:     }
343:   }
344:   if (cvode->maxdt > 0) {
345:     flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt);
346:     if (flag) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMaxStep() failed");
347:   }

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

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

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);

362:   /* call CVSpgmr to use GMRES as the linear solver.        */
363:   /* setup the ode integrator with the given preconditioner */
364:   TSSundialsGetPC(ts,&pc);
365:   PCGetType(pc,&pctype);
366:   PetscObjectTypeCompare((PetscObject)pc,PCNONE,&pcnone);
367:   if (pcnone) {
368:     flag = CVSpgmr(mem,PREC_NONE,0);
369:     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
370:   } else {
371:     flag = CVSpgmr(mem,PREC_LEFT,cvode->maxl);
372:     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);

374:     /* Set preconditioner and solve routines Precond and PSolve,
375:      and the pointer to the user-defined block data */
376:     flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials);
377:     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag);
378:   }

380:   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
381:   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag);
382:   return(0);
383: }

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

392: PetscErrorCode TSSetFromOptions_Sundials(PetscOptions *PetscOptionsObject,TS ts)
393: {
394:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
396:   int            indx;
397:   PetscBool      flag;
398:   PC             pc;

401:   PetscOptionsHead(PetscOptionsObject,"SUNDIALS ODE solver options");
402:   PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);
403:   if (flag) {
404:     TSSundialsSetType(ts,(TSSundialsLmmType)indx);
405:   }
406:   PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);
407:   if (flag) {
408:     TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);
409:   }
410:   PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,NULL);
411:   PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,NULL);
412:   PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,NULL);
413:   PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,NULL);
414:   PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,NULL);
415:   PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,NULL);
416:   PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,NULL);
417:   PetscOptionsTail();
418:   TSSundialsGetPC(ts,&pc);
419:   PCSetFromOptions(pc);
420:   return(0);
421: }

425: PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
426: {
427:   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;

439:   if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype;
440:   else                                     type = btype;

442:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
443:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
444:   if (iascii) {
445:     PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");
446:     PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);
447:     PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);
448:     PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);
449:     PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %D\n",cvode->maxl);
450:     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
451:       PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");
452:     } else {
453:       PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");
454:     }
455:     if (cvode->mindt > 0) {PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);}
456:     if (cvode->maxdt > 0) {PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);}

458:     /* Outputs from CVODE, CVSPILS */
459:     CVodeGetTolScaleFactor(cvode->mem,&tolsfac);
460:     PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);
461:     CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
462:                                    &nlinsetups,&nfails,&qlast,&qcur,
463:                                    &hinused,&hlast,&hcur,&tcur);
464:     PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);
465:     PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);
466:     PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);
467:     PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);

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

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

478:     TSSundialsGetPC(ts,&pc);
479:     PCView(pc,viewer);
480:     CVSpilsGetNumPrecEvals(cvode->mem,&itmp);
481:     PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);
482:     CVSpilsGetNumPrecSolves(cvode->mem,&itmp);
483:     PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);

485:     CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);
486:     PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);
487:     CVSpilsGetNumRhsEvals(cvode->mem,&itmp);
488:     PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);
489:   } else if (isstring) {
490:     PetscViewerStringSPrintf(viewer,"Sundials type %s",type);
491:   }
492:   return(0);
493: }


496: /* --------------------------------------------------------------------------*/
499: PetscErrorCode  TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
500: {
501:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

504:   cvode->cvode_type = type;
505:   return(0);
506: }

510: PetscErrorCode  TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl)
511: {
512:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

515:   cvode->maxl = maxl;
516:   return(0);
517: }

521: PetscErrorCode  TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
522: {
523:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

526:   cvode->linear_tol = tol;
527:   return(0);
528: }

532: PetscErrorCode  TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
533: {
534:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

537:   cvode->gtype = type;
538:   return(0);
539: }

543: PetscErrorCode  TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
544: {
545:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

548:   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
549:   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
550:   return(0);
551: }

555: PetscErrorCode  TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)
556: {
557:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

560:   cvode->mindt = mindt;
561:   return(0);
562: }

566: PetscErrorCode  TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)
567: {
568:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

571:   cvode->maxdt = maxdt;
572:   return(0);
573: }
576: PetscErrorCode  TSSundialsGetPC_Sundials(TS ts,PC *pc)
577: {
578:   SNES           snes;
579:   KSP            ksp;

583:   TSGetSNES(ts,&snes);
584:   SNESGetKSP(snes,&ksp);
585:   KSPGetPC(ksp,pc);
586:   return(0);
587: }

591: PetscErrorCode  TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
592: {
594:   if (nonlin) *nonlin = ts->snes_its;
595:   if (lin)    *lin    = ts->ksp_its;
596:   return(0);
597: }

601: PetscErrorCode  TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s)
602: {
603:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

606:   cvode->monitorstep = s;
607:   return(0);
608: }
609: /* -------------------------------------------------------------------------------------------*/

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

616:    Not Collective

618:    Input parameters:
619: .    ts     - the time-step context

621:    Output Parameters:
622: +   nonlin - number of nonlinear iterations
623: -   lin    - number of linear iterations

625:    Level: advanced

627:    Notes:
628:     These return the number since the creation of the TS object

630: .keywords: non-linear iterations, linear iterations

632: .seealso: TSSundialsSetType(), TSSundialsSetMaxl(),
633:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
634:           TSSundialsGetIterations(), TSSundialsSetType(),
635:           TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime()

637: @*/
638: PetscErrorCode  TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
639: {

643:   PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));
644:   return(0);
645: }

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

652:    Logically Collective on TS

654:    Input parameters:
655: +    ts     - the time-step context
656: -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF

658:    Level: intermediate

660: .keywords: Adams, backward differentiation formula

662: .seealso: TSSundialsGetIterations(),  TSSundialsSetMaxl(),
663:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
664:           TSSundialsGetIterations(), TSSundialsSetType(),
665:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
666:           TSSetExactFinalTime()
667: @*/
668: PetscErrorCode  TSSundialsSetType(TS ts,TSSundialsLmmType type)
669: {

673:   PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));
674:   return(0);
675: }

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

684:    Logically Collective on TS

686:    Input parameters:
687: +    ts      - the time-step context
688: -    maxl - number of direction vectors (the dimension of Krylov subspace).

690:    Level: advanced

692: .keywords: GMRES

694: .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
695:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
696:           TSSundialsGetIterations(), TSSundialsSetType(),
697:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
698:           TSSetExactFinalTime()

700: @*/
701: PetscErrorCode  TSSundialsSetMaxl(TS ts,PetscInt maxl)
702: {

707:   PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));
708:   return(0);
709: }

713: /*@
714:    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
715:        system by SUNDIALS.

717:    Logically Collective on TS

719:    Input parameters:
720: +    ts     - the time-step context
721: -    tol    - the factor by which the tolerance on the nonlinear solver is
722:              multiplied to get the tolerance on the linear solver, .05 by default.

724:    Level: advanced

726: .keywords: GMRES, linear convergence tolerance, SUNDIALS

728: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
729:           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
730:           TSSundialsGetIterations(), TSSundialsSetType(),
731:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
732:           TSSetExactFinalTime()

734: @*/
735: PetscErrorCode  TSSundialsSetLinearTolerance(TS ts,double tol)
736: {

741:   PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));
742:   return(0);
743: }

747: /*@
748:    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
749:         in GMRES method by SUNDIALS linear solver.

751:    Logically Collective on TS

753:    Input parameters:
754: +    ts  - the time-step context
755: -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS

757:    Level: advanced

759: .keywords: Sundials, orthogonalization

761: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
762:           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
763:           TSSundialsGetIterations(), TSSundialsSetType(),
764:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
765:           TSSetExactFinalTime()

767: @*/
768: PetscErrorCode  TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
769: {

773:   PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));
774:   return(0);
775: }

779: /*@
780:    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
781:                          Sundials for error control.

783:    Logically Collective on TS

785:    Input parameters:
786: +    ts  - the time-step context
787: .    aabs - the absolute tolerance
788: -    rel - the relative tolerance

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

793:    Level: intermediate

795: .keywords: Sundials, tolerance

797: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESMaxl(),
798:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
799:           TSSundialsGetIterations(), TSSundialsSetType(),
800:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
801:           TSSetExactFinalTime()

803: @*/
804: PetscErrorCode  TSSundialsSetTolerance(TS ts,double aabs,double rel)
805: {

809:   PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));
810:   return(0);
811: }

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

818:    Input Parameter:
819: .    ts - the time-step context

821:    Output Parameter:
822: .    pc - the preconditioner context

824:    Level: advanced

826: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
827:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
828:           TSSundialsGetIterations(), TSSundialsSetType(),
829:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
830: @*/
831: PetscErrorCode  TSSundialsGetPC(TS ts,PC *pc)
832: {

836:   PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC*),(ts,pc));
837:   return(0);
838: }

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

845:    Input Parameter:
846: +   ts - the time-step context
847: -   mindt - lowest time step if positive, negative to deactivate

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

853:    Level: beginner

855: .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
856: @*/
857: PetscErrorCode  TSSundialsSetMinTimeStep(TS ts,PetscReal mindt)
858: {

862:   PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));
863:   return(0);
864: }

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

871:    Input Parameter:
872: +   ts - the time-step context
873: -   maxdt - lowest time step if positive, negative to deactivate

875:    Level: beginner

877: .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
878: @*/
879: PetscErrorCode  TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt)
880: {

884:   PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));
885:   return(0);
886: }

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

893:    Input Parameter:
894: +   ts - the time-step context
895: -   ft - PETSC_TRUE if monitor, else PETSC_FALSE

897:    Level: beginner

899: .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
900:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
901:           TSSundialsGetIterations(), TSSundialsSetType(),
902:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
903: @*/
904: PetscErrorCode  TSSundialsMonitorInternalSteps(TS ts,PetscBool ft)
905: {

909:   PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));
910:   return(0);
911: }
912: /* -------------------------------------------------------------------------------------------*/
913: /*MC
914:       TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)

916:    Options Database:
917: +    -ts_sundials_type <bdf,adams>
918: .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
919: .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
920: .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
921: .    -ts_sundials_linear_tolerance <tol>
922: .    -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace
923: -    -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps


926:     Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
927:            only PETSc PC options

929:     Level: beginner

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

934: M*/
937: PETSC_EXTERN PetscErrorCode TSCreate_Sundials(TS ts)
938: {
939:   TS_Sundials    *cvode;
941:   PC             pc;

944:   ts->ops->reset          = TSReset_Sundials;
945:   ts->ops->destroy        = TSDestroy_Sundials;
946:   ts->ops->view           = TSView_Sundials;
947:   ts->ops->setup          = TSSetUp_Sundials;
948:   ts->ops->step           = TSStep_Sundials;
949:   ts->ops->interpolate    = TSInterpolate_Sundials;
950:   ts->ops->setfromoptions = TSSetFromOptions_Sundials;

952:   PetscNewLog(ts,&cvode);

954:   ts->data           = (void*)cvode;
955:   cvode->cvode_type  = SUNDIALS_BDF;
956:   cvode->gtype       = SUNDIALS_CLASSICAL_GS;
957:   cvode->maxl        = 5;
958:   cvode->linear_tol  = .05;
959:   cvode->monitorstep = PETSC_TRUE;

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

963:   cvode->mindt = -1.;
964:   cvode->maxdt = -1.;

966:   /* set tolerance for Sundials */
967:   cvode->reltol = 1e-6;
968:   cvode->abstol = 1e-6;

970:   /* set PCNONE as default pctype */
971:   TSSundialsGetPC_Sundials(ts,&pc);
972:   PCSetType(pc,PCNONE);

974:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",TSSundialsSetType_Sundials);
975:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",TSSundialsSetMaxl_Sundials);
976:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",TSSundialsSetLinearTolerance_Sundials);
977:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",TSSundialsSetGramSchmidtType_Sundials);
978:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",TSSundialsSetTolerance_Sundials);
979:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",TSSundialsSetMinTimeStep_Sundials);
980:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",TSSundialsSetMaxTimeStep_Sundials);
981:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",TSSundialsGetPC_Sundials);
982:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",TSSundialsGetIterations_Sundials);
983:   PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",TSSundialsMonitorInternalSteps_Sundials);
984:   return(0);
985: }