Actual source code: sundials.c

petsc-3.3-p2 2012-07-13
  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 sundials.h

 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,
 17:                     booleantype jok,booleantype *jcurPtr,
 18:                     realtype _gamma,void *P_data,
 19:                     N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3)
 20: {
 21:   TS             ts = (TS) P_data;
 22:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
 23:   PC             pc;
 25:   Mat            J,P;
 26:   Vec            yy = cvode->w1,yydot = cvode->ydot;
 27:   PetscReal      gm = (PetscReal)_gamma;
 28:   MatStructure   str = DIFFERENT_NONZERO_PATTERN;
 29:   PetscScalar    *y_data;

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

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

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

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

 76: /*
 77:         TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
 78: */
 81: int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx)
 82: {
 83:   TS              ts = (TS) ctx;
 84:   MPI_Comm        comm = ((PetscObject)ts)->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;
 88:   PetscErrorCode  ierr;

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

 97:   /* now compute the right hand side function */
 98:   if (!ts->userops->ifunction) {
 99:     TSComputeRHSFunction(ts,t,yy,yyd);
100:   } else {                      /* If rhsfunction is also set, this computes both parts and shifts them to the right */
101:     VecZeroEntries(yydot);
102:     TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE); CHKERRABORT(comm,ierr);
103:     VecScale(yyd,-1.);
104:   }
105:   VecResetArray(yy); CHKERRABORT(comm,ierr);
106:   VecResetArray(yyd); CHKERRABORT(comm,ierr);
107:   return(0);
108: }

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

126:   mem  = cvode->mem;
127:   tout = ts->max_time;
128:   VecGetArray(ts->vec_sol,&y_data);
129:   N_VSetArrayPointer((realtype *)y_data,cvode->y);
130:   VecRestoreArray(ts->vec_sol,PETSC_NULL);

132:   TSPreStep(ts);

134:   if (cvode->monitorstep) {
135:     flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
136:   } else {
137:     flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL);
138:   }

140:   if (flag){ /* display error message */
141:     switch (flag){
142:       case CV_ILL_INPUT:
143:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT");
144:         break;
145:       case CV_TOO_CLOSE:
146:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE");
147:         break;
148:       case CV_TOO_MUCH_WORK: {
149:         PetscReal      tcur;
150:         CVodeGetNumSteps(mem,&nsteps);
151:         CVodeGetCurrentTime(mem,&tcur);
152:         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()",tcur,nsteps,ts->max_steps);
153:       } break;
154:       case CV_TOO_MUCH_ACC:
155:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC");
156:         break;
157:       case CV_ERR_FAILURE:
158:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE");
159:         break;
160:       case CV_CONV_FAILURE:
161:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE");
162:         break;
163:       case CV_LINIT_FAIL:
164:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL");
165:         break;
166:       case CV_LSETUP_FAIL:
167:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL");
168:         break;
169:       case CV_LSOLVE_FAIL:
170:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL");
171:         break;
172:       case CV_RHSFUNC_FAIL:
173:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL");
174:         break;
175:       case CV_FIRST_RHSFUNC_ERR:
176:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR");
177:         break;
178:       case CV_REPTD_RHSFUNC_ERR:
179:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR");
180:         break;
181:       case CV_UNREC_RHSFUNC_ERR:
182:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR");
183:         break;
184:       case CV_RTFUNC_FAIL:
185:         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL");
186:         break;
187:       default:
188:         SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag);
189:     }
190:   }

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);
197:   CVodeGetNumNonlinSolvIters(mem,&its);
198:   CVSpilsGetNumLinIters(mem, &its);
199:   ts->snes_its = its; ts->ksp_its = its;

201:   ts->time_step = t - ts->ptime;
202:   ts->ptime     = t;
203:   ts->steps++;

205:   CVodeGetNumSteps(mem,&nsteps);
206:   if (!cvode->monitorstep) ts->steps = nsteps;
207:   return(0);
208: }

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


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

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

230:   VecGetArray(X,&x_data);
231:   N_VSetArrayPointer((realtype *)x_data,y);
232:   CVodeGetDky(cvode->mem,t,0,y);
233:   VecRestoreArray(X,&x_data);

235:   return(0);
236: }

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

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

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

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

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

292:   /* get the vector size */
293:   VecGetSize(ts->vec_sol,&glosize);
294:   VecGetLocalSize(ts->vec_sol,&locsize);

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

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

306:   VecDuplicate(ts->vec_sol,&cvode->update);
307:   VecDuplicate(ts->vec_sol,&cvode->ydot);
308:   PetscLogObjectParent(ts,cvode->update);
309:   PetscLogObjectParent(ts,cvode->ydot);

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

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

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

330:   /* Sundials may choose to use a smaller initial step, but will never use a larger step. */
331:   flag = CVodeSetInitStep(mem,(realtype)ts->time_step);
332:   if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetInitStep() failed");
333:   if (cvode->mindt > 0) {
334:     flag = CVodeSetMinStep(mem,(realtype)cvode->mindt);
335:     if (flag){
336:       if (flag == CV_MEM_NULL){
337:         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL");
338:       } else if (flag == CV_ILL_INPUT){
339:         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size");
340:       } else {
341:         SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed");
342:       }
343:     }
344:   }
345:   if (cvode->maxdt > 0) {
346:     flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt);
347:     if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMaxStep() failed");
348:   }

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

358:   /* specifies scalar relative and absolute tolerances */
359:   flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol);
360:   if (flag){
361:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag);
362:   }

364:   /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */
365:   flag = CVodeSetMaxNumSteps(mem,ts->max_steps);

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

379:     /* Set preconditioner and solve routines Precond and PSolve,
380:      and the pointer to the user-defined block data */
381:     flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials);
382:     if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag);
383:   }

385:   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
386:   if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag);
387:   return(0);
388: }

390: /* type of CVODE linear multistep method */
391: const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0};
392: /* type of G-S orthogonalization used by CVODE linear solver */
393: const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0};

397: PetscErrorCode TSSetFromOptions_Sundials(TS ts)
398: {
399:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
401:   int            indx;
402:   PetscBool      flag;
403:   PC             pc;

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

430: PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
431: {
432:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
434:   char           *type;
435:   char           atype[] = "Adams";
436:   char           btype[] = "BDF: backward differentiation formula";
437:   PetscBool      iascii,isstring;
438:   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
439:   PetscInt       qlast,qcur;
440:   PetscReal      hinused,hlast,hcur,tcur,tolsfac;
441:   PC             pc;

444:   if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;}
445:   else                                     {type = btype;}

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

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

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

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

483:     TSSundialsGetPC(ts,&pc);
484:     PCView(pc,viewer);
485:     CVSpilsGetNumPrecEvals(cvode->mem,&itmp);
486:     PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);
487:     CVSpilsGetNumPrecSolves(cvode->mem,&itmp);
488:     PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);

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


501: /* --------------------------------------------------------------------------*/
502: EXTERN_C_BEGIN
505: PetscErrorCode  TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
506: {
507:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

510:   cvode->cvode_type = type;
511:   return(0);
512: }
513: EXTERN_C_END

515: EXTERN_C_BEGIN
518: PetscErrorCode  TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl)
519: {
520:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

523:   cvode->maxl = maxl;
524:   return(0);
525: }
526: EXTERN_C_END

528: EXTERN_C_BEGIN
531: PetscErrorCode  TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
532: {
533:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

536:   cvode->linear_tol = tol;
537:   return(0);
538: }
539: EXTERN_C_END

541: EXTERN_C_BEGIN
544: PetscErrorCode  TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
545: {
546:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

549:   cvode->gtype = type;
550:   return(0);
551: }
552: EXTERN_C_END

554: EXTERN_C_BEGIN
557: PetscErrorCode  TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
558: {
559:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

562:   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
563:   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
564:   return(0);
565: }
566: EXTERN_C_END

568: EXTERN_C_BEGIN
571: PetscErrorCode  TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)
572: {
573:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

576:   cvode->mindt = mindt;
577:   return(0);
578: }
579: EXTERN_C_END

581: EXTERN_C_BEGIN
584: PetscErrorCode  TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)
585: {
586:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

589:   cvode->maxdt = maxdt;
590:   return(0);
591: }
592: EXTERN_C_END

594: EXTERN_C_BEGIN
597: PetscErrorCode  TSSundialsGetPC_Sundials(TS ts,PC *pc)
598: {
599:   SNES            snes;
600:   KSP             ksp;
601:   PetscErrorCode  ierr;

604:   TSGetSNES(ts,&snes);
605:   SNESGetKSP(snes,&ksp);
606:   KSPGetPC(ksp,pc);
607:   return(0);
608: }
609: EXTERN_C_END

611: EXTERN_C_BEGIN
614: PetscErrorCode  TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
615: {
617:   if (nonlin) *nonlin = ts->snes_its;
618:   if (lin)    *lin    = ts->ksp_its;
619:   return(0);
620: }
621: EXTERN_C_END

623: EXTERN_C_BEGIN
626: PetscErrorCode  TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool  s)
627: {
628:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

631:   cvode->monitorstep = s;
632:   return(0);
633: }
634: EXTERN_C_END
635: /* -------------------------------------------------------------------------------------------*/

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

642:    Not Collective

644:    Input parameters:
645: .    ts     - the time-step context

647:    Output Parameters:
648: +   nonlin - number of nonlinear iterations
649: -   lin    - number of linear iterations

651:    Level: advanced

653:    Notes:
654:     These return the number since the creation of the TS object

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

658: .seealso: TSSundialsSetType(), TSSundialsSetMaxl(),
659:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
660:           TSSundialsGetIterations(), TSSundialsSetType(),
661:           TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime()

663: @*/
664: PetscErrorCode  TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
665: {

669:   PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));
670:   return(0);
671: }

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

678:    Logically Collective on TS

680:    Input parameters:
681: +    ts     - the time-step context
682: -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF

684:    Level: intermediate

686: .keywords: Adams, backward differentiation formula

688: .seealso: TSSundialsGetIterations(),  TSSundialsSetMaxl(),
689:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
690:           TSSundialsGetIterations(), TSSundialsSetType(), 
691:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
692:           TSSetExactFinalTime()
693: @*/
694: PetscErrorCode  TSSundialsSetType(TS ts,TSSundialsLmmType type)
695: {

699:   PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));
700:   return(0);
701: }

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

710:    Logically Collective on TS

712:    Input parameters:
713: +    ts      - the time-step context
714: -    maxl - number of direction vectors (the dimension of Krylov subspace).

716:    Level: advanced

718: .keywords: GMRES

720: .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
721:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
722:           TSSundialsGetIterations(), TSSundialsSetType(),
723:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
724:           TSSetExactFinalTime()

726: @*/
727: PetscErrorCode  TSSundialsSetMaxl(TS ts,PetscInt maxl)
728: {

733:   PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));
734:   return(0);
735: }

739: /*@
740:    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
741:        system by SUNDIALS.

743:    Logically Collective on TS

745:    Input parameters:
746: +    ts     - the time-step context
747: -    tol    - the factor by which the tolerance on the nonlinear solver is
748:              multiplied to get the tolerance on the linear solver, .05 by default.

750:    Level: advanced

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

754: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
755:           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
756:           TSSundialsGetIterations(), TSSundialsSetType(),
757:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
758:           TSSetExactFinalTime()

760: @*/
761: PetscErrorCode  TSSundialsSetLinearTolerance(TS ts,double tol)
762: {

767:   PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));
768:   return(0);
769: }

773: /*@
774:    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
775:         in GMRES method by SUNDIALS linear solver.

777:    Logically Collective on TS

779:    Input parameters:
780: +    ts  - the time-step context
781: -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS

783:    Level: advanced

785: .keywords: Sundials, orthogonalization

787: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
788:           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
789:           TSSundialsGetIterations(), TSSundialsSetType(), 
790:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
791:           TSSetExactFinalTime()

793: @*/
794: PetscErrorCode  TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
795: {

799:   PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));
800:   return(0);
801: }

805: /*@
806:    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
807:                          Sundials for error control.

809:    Logically Collective on TS

811:    Input parameters:
812: +    ts  - the time-step context
813: .    aabs - the absolute tolerance
814: -    rel - the relative tolerance

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

819:    Level: intermediate

821: .keywords: Sundials, tolerance

823: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESMaxl(),
824:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
825:           TSSundialsGetIterations(), TSSundialsSetType(),
826:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
827:           TSSetExactFinalTime()

829: @*/
830: PetscErrorCode  TSSundialsSetTolerance(TS ts,double aabs,double rel)
831: {

835:   PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));
836:   return(0);
837: }

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

844:    Input Parameter:
845: .    ts - the time-step context

847:    Output Parameter:
848: .    pc - the preconditioner context

850:    Level: advanced

852: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
853:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
854:           TSSundialsGetIterations(), TSSundialsSetType(), 
855:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
856: @*/
857: PetscErrorCode  TSSundialsGetPC(TS ts,PC *pc)
858: {

862:   PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC *),(ts,pc));
863:   return(0);
864: }

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

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

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

879:    Level: beginner

881: .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
882: @*/
883: PetscErrorCode  TSSundialsSetMinTimeStep(TS ts,PetscReal mindt)
884: {

888:   PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));
889:   return(0);
890: }

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

897:    Input Parameter:
898: +   ts - the time-step context
899: -   maxdt - lowest time step if positive, negative to deactivate

901:    Level: beginner

903: .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
904: @*/
905: PetscErrorCode  TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt)
906: {

910:   PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));
911:   return(0);
912: }

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

919:    Input Parameter:
920: +   ts - the time-step context
921: -   ft - PETSC_TRUE if monitor, else PETSC_FALSE

923:    Level: beginner

925: .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
926:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
927:           TSSundialsGetIterations(), TSSundialsSetType(), 
928:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
929: @*/
930: PetscErrorCode  TSSundialsMonitorInternalSteps(TS ts,PetscBool  ft)
931: {

935:   PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));
936:   return(0);
937: }
938: /* -------------------------------------------------------------------------------------------*/
939: /*MC
940:       TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)

942:    Options Database:
943: +    -ts_sundials_type <bdf,adams>
944: .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
945: .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
946: .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
947: .    -ts_sundials_linear_tolerance <tol>
948: .    -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace
949: -    -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps


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

955:     Level: beginner

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

960: M*/
961: EXTERN_C_BEGIN
964: PetscErrorCode  TSCreate_Sundials(TS ts)
965: {
966:   TS_Sundials    *cvode;
968:   PC             pc;

971:   ts->ops->reset          = TSReset_Sundials;
972:   ts->ops->destroy        = TSDestroy_Sundials;
973:   ts->ops->view           = TSView_Sundials;
974:   ts->ops->setup          = TSSetUp_Sundials;
975:   ts->ops->step           = TSStep_Sundials;
976:   ts->ops->interpolate    = TSInterpolate_Sundials;
977:   ts->ops->setfromoptions = TSSetFromOptions_Sundials;

979:   PetscNewLog(ts,TS_Sundials,&cvode);
980:   ts->data                = (void*)cvode;
981:   cvode->cvode_type       = SUNDIALS_BDF;
982:   cvode->gtype            = SUNDIALS_CLASSICAL_GS;
983:   cvode->maxl             = 5;
984:   cvode->linear_tol       = .05;

986:   cvode->monitorstep      = PETSC_TRUE;

988:   MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));

990:   cvode->mindt = -1.;
991:   cvode->maxdt = -1.;

993:   /* set tolerance for Sundials */
994:   cvode->reltol = 1e-6;
995:   cvode->abstol = 1e-6;

997:   /* set PCNONE as default pctype */
998:   TSSundialsGetPC_Sundials(ts,&pc);
999:   PCSetType(pc,PCNONE);

1001:   if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_FALSE;

1003:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
1004:                     TSSundialsSetType_Sundials);
1005:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxl_C",
1006:                     "TSSundialsSetMaxl_Sundials",
1007:                     TSSundialsSetMaxl_Sundials);
1008:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
1009:                     "TSSundialsSetLinearTolerance_Sundials",
1010:                      TSSundialsSetLinearTolerance_Sundials);
1011:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
1012:                     "TSSundialsSetGramSchmidtType_Sundials",
1013:                      TSSundialsSetGramSchmidtType_Sundials);
1014:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
1015:                     "TSSundialsSetTolerance_Sundials",
1016:                      TSSundialsSetTolerance_Sundials);
1017:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C",
1018:                     "TSSundialsSetMinTimeStep_Sundials",
1019:                      TSSundialsSetMinTimeStep_Sundials);
1020:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",
1021:                     "TSSundialsSetMaxTimeStep_Sundials",
1022:                      TSSundialsSetMaxTimeStep_Sundials);
1023:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
1024:                     "TSSundialsGetPC_Sundials",
1025:                      TSSundialsGetPC_Sundials);
1026:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
1027:                     "TSSundialsGetIterations_Sundials",
1028:                      TSSundialsGetIterations_Sundials);
1029:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",
1030:                     "TSSundialsMonitorInternalSteps_Sundials",
1031:                      TSSundialsMonitorInternalSteps_Sundials);

1033:   return(0);
1034: }
1035: EXTERN_C_END