Actual source code: sundials.c

  1: #define PETSCTS_DLL

  3: /*
  4:     Provides a PETSc interface to SUNDIALS/CVODE solver.
  5:     The interface to PVODE (old version of CVODE) was originally contributed 
  6:     by Liyang Xu. It has been redone by Hong Zhang and Dinesh Kaushik.

  8:     Reference: sundials-2.4.0/examples/cvode/parallel/cvDiurnal_kry_p.c
  9: */
 10:  #include sundials.h

 12: /*
 13:       TSPrecond_Sundials - function that we provide to SUNDIALS to
 14:                         evaluate the preconditioner.
 15: */
 18: PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy,
 19:                     booleantype jok,booleantype *jcurPtr,
 20:                     realtype _gamma,void *P_data,
 21:                     N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3)
 22: {
 23:   TS             ts = (TS) P_data;
 24:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
 25:   PC             pc = cvode->pc;
 27:   Mat            Jac = ts->B;
 28:   Vec            yy = cvode->w1;
 29:   PetscScalar    one = 1.0,gm;
 30:   MatStructure   str = DIFFERENT_NONZERO_PATTERN;
 31:   PetscScalar    *y_data;
 32: 
 34:   /* This allows us to construct preconditioners in-place if we like */
 35:   MatSetUnfactored(Jac);
 36: 
 37:   /* jok - TRUE means reuse current Jacobian else recompute Jacobian */
 38:   if (jok) {
 39:     MatCopy(cvode->pmat,Jac,str);
 40:     *jcurPtr = FALSE;
 41:   } else {
 42:     /* make PETSc vector yy point to SUNDIALS vector y */
 43:     y_data = (PetscScalar *) N_VGetArrayPointer(y);
 44:     VecPlaceArray(yy,y_data);

 46:     /* compute the Jacobian */
 47:     TSComputeRHSJacobian(ts,ts->ptime,yy,&Jac,&Jac,&str);
 48:     VecResetArray(yy);

 50:     /* copy the Jacobian matrix */
 51:     if (!cvode->pmat) {
 52:       MatDuplicate(Jac,MAT_COPY_VALUES,&cvode->pmat);
 53:       PetscLogObjectParent(ts,cvode->pmat);
 54:     } else {
 55:       MatCopy(Jac,cvode->pmat,str);
 56:     }
 57:     *jcurPtr = TRUE;
 58:   }
 59: 
 60:   /* construct I-gamma*Jac  */
 61:   gm   = -_gamma;
 62:   MatScale(Jac,gm);
 63:   MatShift(Jac,one);
 64: 
 65:   PCSetOperators(pc,Jac,Jac,str);
 66:   return(0);
 67: }

 69: /*
 70:      TSPSolve_Sundials -  routine that we provide to Sundials that applies the preconditioner.
 71: */
 74: PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z,
 75:                                  realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp)
 76: {
 77:   TS              ts = (TS) P_data;
 78:   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
 79:   PC              pc = cvode->pc;
 80:   Vec             rr = cvode->w1,zz = cvode->w2;
 81:   PetscErrorCode  ierr;
 82:   PetscScalar     *r_data,*z_data;

 85:   /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
 86:   r_data  = (PetscScalar *) N_VGetArrayPointer(r);
 87:   z_data  = (PetscScalar *) N_VGetArrayPointer(z);
 88:   VecPlaceArray(rr,r_data);
 89:   VecPlaceArray(zz,z_data);

 91:   /* Solve the Px=r and put the result in zz */
 92:   PCApply(pc,rr,zz);
 93:   VecResetArray(rr);
 94:   VecResetArray(zz);
 95:   return(0);
 96: }

 98: /*
 99:         TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
100: */
103: int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx)
104: {
105:   TS              ts = (TS) ctx;
106:   MPI_Comm        comm = ((PetscObject)ts)->comm;
107:   TS_Sundials     *cvode = (TS_Sundials*)ts->data;
108:   Vec             yy = cvode->w1,yyd = cvode->w2;
109:   PetscScalar     *y_data,*ydot_data;
110:   PetscErrorCode  ierr;

113:   /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
114:   y_data     = (PetscScalar *) N_VGetArrayPointer(y);
115:   ydot_data  = (PetscScalar *) N_VGetArrayPointer(ydot);
116:   VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr);
117:   VecPlaceArray(yyd,ydot_data); CHKERRABORT(comm,ierr);

119:   /* now compute the right hand side function */
120:   TSComputeRHSFunction(ts,t,yy,yyd); CHKERRABORT(comm,ierr);
121:   VecResetArray(yy); CHKERRABORT(comm,ierr);
122:   VecResetArray(yyd); CHKERRABORT(comm,ierr);
123:   return(0);
124: }

126: /*
127:        TSStep_Sundials_Nonlinear - Calls Sundials to integrate the ODE.
128: */
131: PetscErrorCode TSStep_Sundials_Nonlinear(TS ts,int *steps,double *time)
132: {
133:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
134:   Vec            sol = ts->vec_sol;
136:   PetscInt       i,max_steps = ts->max_steps,flag;
137:   long int       its;
138:   realtype       t,tout;
139:   PetscScalar    *y_data;
140:   void           *mem;
141: 
143:   mem  = cvode->mem;
144:   tout = ts->max_time;
145:   VecGetArray(ts->vec_sol,&y_data);
146:   N_VSetArrayPointer((realtype *)y_data,cvode->y);
147:   VecRestoreArray(ts->vec_sol,PETSC_NULL);
148:   for (i = 0; i < max_steps; i++) {
149:     if (ts->ptime >= ts->max_time) break;
150:     TSPreStep(ts);
151:     if (cvode->monitorstep){
152:       flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
153:     } else {
154:       flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL);
155:     }
156:     if (flag)SETERRQ1(PETSC_ERR_LIB,"CVode() fails, flag %d",flag);
157:     if (t > ts->max_time && cvode->exact_final_time) {
158:       /* interpolate to final requested time */
159:       CVodeGetDky(mem,tout,0,cvode->y);
160:       t = tout;
161:     }
162:     ts->time_step = t - ts->ptime;
163:     ts->ptime     = t;

165:     /* copy the solution from cvode->y to cvode->update and sol */
166:     VecPlaceArray(cvode->w1,y_data);
167:     VecCopy(cvode->w1,cvode->update);
168:     VecResetArray(cvode->w1);
169:     VecCopy(cvode->update,sol);
170:     CVodeGetNumNonlinSolvIters(mem,&its);
171:     ts->nonlinear_its = its;
172:     CVSpilsGetNumLinIters(mem, &its);
173:     ts->linear_its = its;
174:     ts->steps++;
175:     TSPostStep(ts);
176:     TSMonitor(ts,ts->steps,t,sol);
177:   }
178:   *steps += ts->steps;
179:   *time   = t;
180:   return(0);
181: }

185: PetscErrorCode TSDestroy_Sundials(TS ts)
186: {
187:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;

191:   if (cvode->pmat)   {MatDestroy(cvode->pmat);}
192:   if (cvode->pc)     {PCDestroy(cvode->pc);}
193:   if (cvode->update) {VecDestroy(cvode->update);}
194:   if (cvode->func)   {VecDestroy(cvode->func);}
195:   if (cvode->rhs)    {VecDestroy(cvode->rhs);}
196:   if (cvode->w1)     {VecDestroy(cvode->w1);}
197:   if (cvode->w2)     {VecDestroy(cvode->w2);}
198:   MPI_Comm_free(&(cvode->comm_sundials));
199:   if (cvode->mem) {CVodeFree(&cvode->mem);}
200:   PetscFree(cvode);
201:   return(0);
202: }

206: PetscErrorCode TSSetUp_Sundials_Nonlinear(TS ts)
207: {
208:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
210:   PetscInt       glosize,locsize,i,flag;
211:   PetscScalar    *y_data,*parray;
212:   void           *mem;
213:   const PCType   pctype;
214:   PetscTruth     pcnone;
215:   Vec            sol = ts->vec_sol;

218:   PCSetFromOptions(cvode->pc);
219:   /* get the vector size */
220:   VecGetSize(ts->vec_sol,&glosize);
221:   VecGetLocalSize(ts->vec_sol,&locsize);

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

227:   /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
228:   VecGetArray(ts->vec_sol,&parray);
229:   y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y);
230:   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
231:   /*PetscMemcpy(y_data,parray,locsize*sizeof(PETSC_SCALAR)); */
232:   VecRestoreArray(ts->vec_sol,PETSC_NULL);
233:   VecDuplicate(ts->vec_sol,&cvode->update);
234:   VecDuplicate(ts->vec_sol,&cvode->func);
235:   PetscLogObjectParent(ts,cvode->update);
236:   PetscLogObjectParent(ts,cvode->func);

238:   /* 
239:     Create work vectors for the TSPSolve_Sundials() routine. Note these are
240:     allocated with zero space arrays because the actual array space is provided 
241:     by Sundials and set using VecPlaceArray().
242:   */
243:   VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);
244:   VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w2);
245:   PetscLogObjectParent(ts,cvode->w1);
246:   PetscLogObjectParent(ts,cvode->w2);

248:   /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
249:   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
250:   if (!mem) SETERRQ(PETSC_ERR_MEM,"CVodeCreate() fails");
251:   cvode->mem = mem;

253:   /* Set the pointer to user-defined data */
254:   flag = CVodeSetUserData(mem, ts);
255:   if (flag) SETERRQ(PETSC_ERR_LIB,"CVodeSetUserData() fails");

257:   /* Call CVodeInit to initialize the integrator memory and specify the
258:    * user's right hand side function in u'=f(t,u), the inital time T0, and
259:    * the initial dependent variable vector cvode->y */
260:   flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y);
261:   if (flag){
262:     SETERRQ1(PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag);
263:   }

265:   flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol);
266:   if (flag){
267:     SETERRQ1(PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag);
268:   }

270:   /* initialize the number of steps */
271:   TSMonitor(ts,ts->steps,ts->ptime,sol);

273:   /* call CVSpgmr to use GMRES as the linear solver.        */
274:   /* setup the ode integrator with the given preconditioner */
275:   PCGetType(cvode->pc,&pctype);
276:   PetscTypeCompare((PetscObject)cvode->pc,PCNONE,&pcnone);
277:   if (pcnone){
278:     flag  = CVSpgmr(mem,PREC_NONE,0);
279:     if (flag) SETERRQ1(PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
280:   } else {
281:     flag  = CVSpgmr(mem,PREC_LEFT,0);
282:     if (flag) SETERRQ1(PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);

284:     /* Set preconditioner and solve routines Precond and PSolve, 
285:      and the pointer to the user-defined block data */
286:     flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials);
287:     if (flag) SETERRQ1(PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag);
288:   }

290:   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
291:   if (flag) SETERRQ1(PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag);
292:   return(0);
293: }

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

302: PetscErrorCode TSSetFromOptions_Sundials_Nonlinear(TS ts)
303: {
304:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
306:   int            indx;
307:   PetscTruth     flag;

310:   PetscOptionsHead("SUNDIALS ODE solver options");
311:     PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);
312:     if (flag) {
313:       TSSundialsSetType(ts,(TSSundialsLmmType)indx);
314:     }
315:     PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);
316:     if (flag) {
317:       TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);
318:     }
319:     PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);
320:     PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);
321:     PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);
322:     PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);
323:     PetscOptionsTruth("-ts_sundials_exact_final_time","Interpolate output to stop exactly at the final time","TSSundialsSetExactFinalTime",cvode->exact_final_time,&cvode->exact_final_time,PETSC_NULL);
324:     PetscOptionsTruth("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);
325:   PetscOptionsTail();
326:   return(0);
327: }

331: PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
332: {
333:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
335:   char           *type;
336:   char           atype[] = "Adams";
337:   char           btype[] = "BDF: backward differentiation formula";
338:   PetscTruth     iascii,isstring;
339:   long int       nsteps,its,nfevals,nlinsetups,nfails,itmp;
340:   PetscInt       qlast,qcur;
341:   PetscReal      hinused,hlast,hcur,tcur,tolsfac;

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

347:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
348:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
349:   if (iascii) {
350:     PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");
351:     PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);
352:     PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);
353:     PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);
354:     PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);
355:     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
356:       PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");
357:     } else {
358:       PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");
359:     }
360: 
361:     /* Outputs from CVODE, CVSPILS */
362:     CVodeGetTolScaleFactor(cvode->mem,&tolsfac);
363:     PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);
364:     CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
365:                                    &nlinsetups,&nfails,&qlast,&qcur,
366:                                    &hinused,&hlast,&hcur,&tcur);
367:     PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);
368:     PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);
369:     PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);
370:     PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);

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

376:     CVSpilsGetNumLinIters(cvode->mem, &its); /* its = no. of calls to TSPrecond_Sundials() */
377:     PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);
378:     CVSpilsGetNumConvFails(cvode->mem,&itmp);
379:     PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);
380:     CVSpilsGetNumPrecEvals(cvode->mem,&itmp);
381:     PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);
382:     CVSpilsGetNumPrecSolves(cvode->mem,&itmp);
383:     PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);
384:     CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);
385:     PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);
386:     CVSpilsGetNumRhsEvals(cvode->mem,&itmp);
387:     PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);
388:   } else if (isstring) {
389:     PetscViewerStringSPrintf(viewer,"Sundials type %s",type);
390:   } else {
391:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by TS Sundials",((PetscObject)viewer)->type_name);
392:   }
393:   PetscViewerASCIIPushTab(viewer);
394:   PCView(cvode->pc,viewer);
395:   PetscViewerASCIIPopTab(viewer);
396:   return(0);
397: }


400: /* --------------------------------------------------------------------------*/
404: PetscErrorCode  TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
405: {
406:   TS_Sundials *cvode = (TS_Sundials*)ts->data;
407: 
409:   cvode->cvode_type = type;
410:   return(0);
411: }

417: PetscErrorCode  TSSundialsSetGMRESRestart_Sundials(TS ts,int restart)
418: {
419:   TS_Sundials *cvode = (TS_Sundials*)ts->data;
420: 
422:   cvode->restart = restart;
423:   return(0);
424: }

430: PetscErrorCode  TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
431: {
432:   TS_Sundials *cvode = (TS_Sundials*)ts->data;
433: 
435:   cvode->linear_tol = tol;
436:   return(0);
437: }

443: PetscErrorCode  TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
444: {
445:   TS_Sundials *cvode = (TS_Sundials*)ts->data;
446: 
448:   cvode->gtype = type;
449:   return(0);
450: }

456: PetscErrorCode  TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
457: {
458:   TS_Sundials *cvode = (TS_Sundials*)ts->data;
459: 
461:   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
462:   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
463:   return(0);
464: }

470: PetscErrorCode  TSSundialsGetPC_Sundials(TS ts,PC *pc)
471: {
472:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

475:   *pc = cvode->pc;
476:   return(0);
477: }

483: PetscErrorCode  TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
484: {
486:   if (nonlin) *nonlin = ts->nonlinear_its;
487:   if (lin)    *lin    = ts->linear_its;
488:   return(0);
489: }
491: 
495: PetscErrorCode  TSSundialsSetExactFinalTime_Sundials(TS ts,PetscTruth s)
496: {
497:   TS_Sundials *cvode = (TS_Sundials*)ts->data;
498: 
500:   cvode->exact_final_time = s;
501:   return(0);
502: }

508: PetscErrorCode  TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscTruth s)
509: {
510:   TS_Sundials *cvode = (TS_Sundials*)ts->data;
511: 
513:   cvode->monitorstep = s;
514:   return(0);
515: }
517: /* -------------------------------------------------------------------------------------------*/

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

524:    Not Collective

526:    Input parameters:
527: .    ts     - the time-step context

529:    Output Parameters:
530: +   nonlin - number of nonlinear iterations
531: -   lin    - number of linear iterations

533:    Level: advanced

535:    Notes:
536:     These return the number since the creation of the TS object

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

540: .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(),
541:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
542:           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
543:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
544:           TSSundialsSetExactFinalTime()

546: @*/
547: PetscErrorCode  TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
548: {
549:   PetscErrorCode ierr,(*f)(TS,int*,int*);
550: 
552:   PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetIterations_C",(void (**)(void))&f);
553:   if (f) {
554:     (*f)(ts,nonlin,lin);
555:   }
556:   return(0);
557: }

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

564:    Collective on TS

566:    Input parameters:
567: +    ts     - the time-step context
568: -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF

570:    Level: intermediate

572: .keywords: Adams, backward differentiation formula

574: .seealso: TSSundialsGetIterations(),  TSSundialsSetGMRESRestart(),
575:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
576:           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
577:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
578:           TSSundialsSetExactFinalTime()
579: @*/
580: PetscErrorCode  TSSundialsSetType(TS ts,TSSundialsLmmType type)
581: {
582:   PetscErrorCode ierr,(*f)(TS,TSSundialsLmmType);
583: 
585:   PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetType_C",(void (**)(void))&f);
586:   if (f) {
587:     (*f)(ts,type);
588:   }
589:   return(0);
590: }

594: /*@
595:    TSSundialsSetGMRESRestart - Sets the dimension of the Krylov space used by 
596:        GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
597:        this is ALSO the maximum number of GMRES steps that will be used.

599:    Collective on TS

601:    Input parameters:
602: +    ts      - the time-step context
603: -    restart - number of direction vectors (the restart size).

605:    Level: advanced

607: .keywords: GMRES, restart

609: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 
610:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
611:           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
612:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
613:           TSSundialsSetExactFinalTime()

615: @*/
616: PetscErrorCode  TSSundialsSetGMRESRestart(TS ts,int restart)
617: {
618:   PetscErrorCode ierr,(*f)(TS,int);

621:   PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGMRESRestart_C",(void (**)(void))&f);
622:   if (f) {
623:     (*f)(ts,restart);
624:   }
625:   return(0);
626: }

630: /*@
631:    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
632:        system by SUNDIALS.

634:    Collective on TS

636:    Input parameters:
637: +    ts     - the time-step context
638: -    tol    - the factor by which the tolerance on the nonlinear solver is
639:              multiplied to get the tolerance on the linear solver, .05 by default.

641:    Level: advanced

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

645: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
646:           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
647:           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
648:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
649:           TSSundialsSetExactFinalTime()

651: @*/
652: PetscErrorCode  TSSundialsSetLinearTolerance(TS ts,double tol)
653: {
654:   PetscErrorCode ierr,(*f)(TS,double);
655: 
657:   PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",(void (**)(void))&f);
658:   if (f) {
659:     (*f)(ts,tol);
660:   }
661:   return(0);
662: }

666: /*@
667:    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
668:         in GMRES method by SUNDIALS linear solver.

670:    Collective on TS

672:    Input parameters:
673: +    ts  - the time-step context
674: -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS

676:    Level: advanced

678: .keywords: Sundials, orthogonalization

680: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
681:           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
682:           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
683:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
684:           TSSundialsSetExactFinalTime()

686: @*/
687: PetscErrorCode  TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
688: {
689:   PetscErrorCode ierr,(*f)(TS,TSSundialsGramSchmidtType);
690: 
692:   PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",(void (**)(void))&f);
693:   if (f) {
694:     (*f)(ts,type);
695:   }
696:   return(0);
697: }

701: /*@
702:    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 
703:                          Sundials for error control.

705:    Collective on TS

707:    Input parameters:
708: +    ts  - the time-step context
709: .    aabs - the absolute tolerance  
710: -    rel - the relative tolerance

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

715:    Level: intermediate

717: .keywords: Sundials, tolerance

719: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
720:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), 
721:           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
722:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
723:           TSSundialsSetExactFinalTime()

725: @*/
726: PetscErrorCode  TSSundialsSetTolerance(TS ts,double aabs,double rel)
727: {
728:   PetscErrorCode ierr,(*f)(TS,double,double);
729: 
731:   PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetTolerance_C",(void (**)(void))&f);
732:   if (f) {
733:     (*f)(ts,aabs,rel);
734:   }
735:   return(0);
736: }

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

743:    Input Parameter:
744: .    ts - the time-step context

746:    Output Parameter:
747: .    pc - the preconditioner context

749:    Level: advanced

751: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
752:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
753:           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
754:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
755: @*/
756: PetscErrorCode  TSSundialsGetPC(TS ts,PC *pc)
757: {
758:   PetscErrorCode ierr,(*f)(TS,PC *);

761:   PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetPC_C",(void (**)(void))&f);
762:   if (f) {
763:     (*f)(ts,pc);
764:   } else {
765:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"TS must be of Sundials type to extract the PC");
766:   }

768:   return(0);
769: }

773: /*@
774:    TSSundialsSetExactFinalTime - Determines if Sundials interpolates solution to the 
775:       exact final time requested by the user or just returns it at the final time
776:       it computed. (Defaults to true).

778:    Input Parameter:
779: +   ts - the time-step context
780: -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE

782:    Level: beginner

784: .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
785:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
786:           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
787:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 
788: @*/
789: PetscErrorCode  TSSundialsSetExactFinalTime(TS ts,PetscTruth ft)
790: {
791:   PetscErrorCode ierr,(*f)(TS,PetscTruth);

794:   PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetExactFinalTime_C",(void (**)(void))&f);
795:   if (f) {
796:     (*f)(ts,ft);
797:   }
798:   return(0);
799: }

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

806:    Input Parameter:
807: +   ts - the time-step context
808: -   ft - PETSC_TRUE if monitor, else PETSC_FALSE

810:    Level: beginner

812: .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
813:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
814:           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
815:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 
816: @*/
817: PetscErrorCode  TSSundialsMonitorInternalSteps(TS ts,PetscTruth ft)
818: {
819:   PetscErrorCode ierr,(*f)(TS,PetscTruth);

822:   PetscObjectQueryFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",(void (**)(void))&f);
823:   if (f) {
824:     (*f)(ts,ft);
825:   }
826:   return(0);
827: }
828: /* -------------------------------------------------------------------------------------------*/
829: /*MC
830:       TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)

832:    Options Database:
833: +    -ts_sundials_type <bdf,adams>
834: .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
835: .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
836: .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
837: .    -ts_sundials_linear_tolerance <tol> 
838: .    -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions
839: .    -ts_sundials_exact_final_time - Interpolate output to stop exactly at the final time
840: -    -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps


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

846:     Level: beginner

848: .seealso:  TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetGMRESRestart(), TSSundialsSetLinearTolerance(),
849:            TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSundialsSetExactFinalTime()

851: M*/
855: PetscErrorCode  TSCreate_Sundials(TS ts)
856: {
857:   TS_Sundials *cvode;

861:   if (ts->problem_type != TS_NONLINEAR) {
862:     SETERRQ(PETSC_ERR_SUP,"Only support for nonlinear problems");
863:   }
864:   ts->ops->destroy        = TSDestroy_Sundials;
865:   ts->ops->view           = TSView_Sundials;
866:   ts->ops->setup          = TSSetUp_Sundials_Nonlinear;
867:   ts->ops->step           = TSStep_Sundials_Nonlinear;
868:   ts->ops->setfromoptions = TSSetFromOptions_Sundials_Nonlinear;

870:   PetscNewLog(ts,TS_Sundials,&cvode);
871:   PCCreate(((PetscObject)ts)->comm,&cvode->pc);
872:   PetscLogObjectParent(ts,cvode->pc);
873:   ts->data          = (void*)cvode;
874:   cvode->cvode_type = SUNDIALS_BDF;
875:   cvode->gtype      = SUNDIALS_CLASSICAL_GS;
876:   cvode->restart    = 5;
877:   cvode->linear_tol = .05;

879:   cvode->exact_final_time = PETSC_TRUE;
880:   cvode->monitorstep      = PETSC_FALSE;

882:   MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));
883:   /* set tolerance for Sundials */
884:   cvode->reltol = 1e-6;
885:   cvode->abstol = 1e-6;

887:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
888:                     TSSundialsSetType_Sundials);
889:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C",
890:                     "TSSundialsSetGMRESRestart_Sundials",
891:                     TSSundialsSetGMRESRestart_Sundials);
892:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
893:                     "TSSundialsSetLinearTolerance_Sundials",
894:                      TSSundialsSetLinearTolerance_Sundials);
895:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
896:                     "TSSundialsSetGramSchmidtType_Sundials",
897:                      TSSundialsSetGramSchmidtType_Sundials);
898:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
899:                     "TSSundialsSetTolerance_Sundials",
900:                      TSSundialsSetTolerance_Sundials);
901:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
902:                     "TSSundialsGetPC_Sundials",
903:                      TSSundialsGetPC_Sundials);
904:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
905:                     "TSSundialsGetIterations_Sundials",
906:                      TSSundialsGetIterations_Sundials);
907:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C",
908:                     "TSSundialsSetExactFinalTime_Sundials",
909:                      TSSundialsSetExactFinalTime_Sundials);
910:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",
911:                     "TSSundialsMonitorInternalSteps_Sundials",
912:                      TSSundialsMonitorInternalSteps_Sundials);

914:   return(0);
915: }