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.
  7:     Reference: sundials-2.3.0/examples/cvode/parallel/cvkryx_p.c
  8: */
 9:  #include sundials.h

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

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

 49:     /* copy the Jacobian matrix */
 50:     if (!cvode->pmat) {
 51:       MatDuplicate(Jac,MAT_COPY_VALUES,&cvode->pmat);
 52:       PetscLogObjectParent(ts,cvode->pmat);
 53:     }
 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:   cvode->linear_solves++;
 96:   return(0);
 97: }

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

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

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

127: /*
128:        TSStep_Sundials_Nonlinear - Calls Sundials to integrate the ODE.
129: */
132: /* 
133:     TSStep_Sundials_Nonlinear - 
134:   
135:    steps - number of time steps
136:    time - time that integrater is  terminated. 
137: */
138: PetscErrorCode TSStep_Sundials_Nonlinear(TS ts,int *steps,double *time)
139: {
140:   TS_Sundials  *cvode = (TS_Sundials*)ts->data;
141:   Vec          sol = ts->vec_sol;
143:   int          i,max_steps = ts->max_steps,flag;
144:   long int     its;
145:   realtype     t,tout;
146:   PetscScalar  *y_data;
147:   void         *mem;

150:   /* Call CVodeCreate to create the solver memory */
151:   mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
152:   if (!mem) SETERRQ(1,"CVodeCreate() fails");
153:   flag = CVodeSetFdata(mem,ts);
154:   if (flag) SETERRQ(1,"CVodeSetFdata() fails");

156:   /* 
157:      Call CVodeMalloc to initialize the integrator memory: 
158:      mem is the pointer to the integrator memory returned by CVodeCreate
159:      f       is the user's right hand side function in y'=f(t,y)
160:      T0      is the initial time
161:      u       is the initial dependent variable vector
162:      CV_SS   specifies scalar relative and absolute tolerances
163:      reltol  is the relative tolerance
164:      &abstol is a pointer to the scalar absolute tolerance
165:   */
166:   flag = CVodeMalloc(mem,TSFunction_Sundials,ts->ptime,cvode->y,CV_SS,cvode->reltol,&cvode->abstol);
167:   if (flag) SETERRQ(1,"CVodeMalloc() fails");

169:   /* initialize the number of steps */
170:   *steps = -ts->steps;
171:   TSMonitor(ts,ts->steps,ts->ptime,sol);

173:   /* call CVSpgmr to use GMRES as the linear solver.        */
174:   /* setup the ode integrator with the given preconditioner */
175:   flag  = CVSpgmr(mem,PREC_LEFT,0);
176:   if (flag) SETERRQ(1,"CVSpgmr() fails");
177: 
178:   flag = CVSpilsSetGSType(mem, MODIFIED_GS);
179:   if (flag) SETERRQ(1,"CVSpgmrSetGSType() fails");

181:   /* Set preconditioner setup and solve routines Precond and PSolve, 
182:      and the pointer to the user-defined block data */
183:   flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials,ts);
184:   if (flag) SETERRQ(1,"CVSpgmrSetPreconditioner() fails");

186:   tout = ts->max_time;
187:   VecGetArray(ts->vec_sol,&y_data);
188:   N_VSetArrayPointer((realtype *)y_data,cvode->y);
189:   VecRestoreArray(ts->vec_sol,PETSC_NULL);
190:   for (i = 0; i < max_steps; i++) {
191:     if (ts->ptime >= ts->max_time) break;
192:     CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
193:     CVodeGetNumNonlinSolvIters(mem,&its);
194:     cvode->nonlinear_solves += its;

196:     if (t > ts->max_time && cvode->exact_final_time) {
197:       /* interpolate to final requested time */
198:       CVodeGetDky(mem,tout,0,cvode->y);
199:       t = tout;
200:     }
201:     ts->time_step = t - ts->ptime;
202:     ts->ptime     = t;

204:     /* copy the solution from cvode->y to cvode->update and sol */
205:     VecPlaceArray(cvode->w1,y_data);
206:     VecCopy(cvode->w1,cvode->update);
207:     VecResetArray(cvode->w1);
208:     VecCopy(cvode->update,sol);
209:     CVodeGetNumNonlinSolvIters(mem,&its);
210:     ts->nonlinear_its = its;
211:     CVSpilsGetNumLinIters(mem, &its);
212:     ts->linear_its = its;
213:     ts->steps++;
214:     TSMonitor(ts,ts->steps,t,sol);
215:   }
216:   CVodeFree(&mem);
217:   *steps += ts->steps;
218:   *time   = t;
219:   return(0);
220: }

224: PetscErrorCode TSDestroy_Sundials(TS ts)
225: {
226:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;

230:   if (cvode->pmat)   {MatDestroy(cvode->pmat);}
231:   if (cvode->pc)     {PCDestroy(cvode->pc);}
232:   if (cvode->update) {VecDestroy(cvode->update);}
233:   if (cvode->func)   {VecDestroy(cvode->func);}
234:   if (cvode->rhs)    {VecDestroy(cvode->rhs);}
235:   if (cvode->w1)     {VecDestroy(cvode->w1);}
236:   if (cvode->w2)     {VecDestroy(cvode->w2);}
237:   MPI_Comm_free(&(cvode->comm_sundials));
238:   PetscFree(cvode);
239:   return(0);
240: }

244: PetscErrorCode TSSetUp_Sundials_Nonlinear(TS ts)
245: {
246:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
248:   int            glosize,locsize,i;
249:   PetscScalar    *y_data,*parray;

252:   PCSetFromOptions(cvode->pc);
253:   /* get the vector size */
254:   VecGetSize(ts->vec_sol,&glosize);
255:   VecGetLocalSize(ts->vec_sol,&locsize);

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

261:   /* initialize N_Vec y */
262:   VecGetArray(ts->vec_sol,&parray);
263:   y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y);
264:   for (i = 0; i < locsize; i++) y_data[i] = parray[i];
265:   /*PetscMemcpy(y_data,parray,locsize*sizeof(PETSC_SCALAR)); */
266:   VecRestoreArray(ts->vec_sol,PETSC_NULL);
267:   VecDuplicate(ts->vec_sol,&cvode->update);
268:   VecDuplicate(ts->vec_sol,&cvode->func);
269:   PetscLogObjectParent(ts,cvode->update);
270:   PetscLogObjectParent(ts,cvode->func);

272:   /* 
273:       Create work vectors for the TSPSolve_Sundials() routine. Note these are
274:     allocated with zero space arrays because the actual array space is provided 
275:     by Sundials and set using VecPlaceArray().
276:   */
277:   VecCreateMPIWithArray(ts->comm,locsize,PETSC_DECIDE,0,&cvode->w1);
278:   VecCreateMPIWithArray(ts->comm,locsize,PETSC_DECIDE,0,&cvode->w2);
279:   PetscLogObjectParent(ts,cvode->w1);
280:   PetscLogObjectParent(ts,cvode->w2);
281:   return(0);
282: }

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

291: PetscErrorCode TSSetFromOptions_Sundials_Nonlinear(TS ts)
292: {
293:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
295:   int            indx;
296:   PetscTruth     flag;

299:   PetscOptionsHead("SUNDIALS ODE solver options");
300:     PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);
301:     if (flag) {
302:       TSSundialsSetType(ts,(TSSundialsLmmType)indx);
303:     }
304:     PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);
305:     if (flag) {
306:       TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);
307:     }
308:     PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);
309:     PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);
310:     PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);
311:     PetscOptionsInt("-ts_sundials_gmres_restart","Number of GMRES orthogonalization directions","TSSundialsSetGMRESRestart",cvode->restart,&cvode->restart,&flag);
312:     PetscOptionsName("-ts_sundials_exact_final_time","Allow SUNDIALS to stop near the final time, not exactly on it","TSSundialsSetExactFinalTime",&flag);
313:     if (flag) cvode->exact_final_time = PETSC_TRUE;
314:   PetscOptionsTail();
315:   return(0);
316: }

320: PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
321: {
322:   TS_Sundials    *cvode = (TS_Sundials*)ts->data;
324:   char           *type;
325:   char            atype[] = "Adams";
326:   char            btype[] = "BDF: backward differentiation formula";
327:   PetscTruth     iascii,isstring;

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

333:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
334:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
335:   if (iascii) {
336:     PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");
337:     PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);
338:     PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);
339:     PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);
340:     PetscViewerASCIIPrintf(viewer,"Sundials GMRES max iterations (same as restart in SUNDIALS) %D\n",cvode->restart);
341:     if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
342:       PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");
343:     } else {
344:       PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");
345:     }
346:   } else if (isstring) {
347:     PetscViewerStringSPrintf(viewer,"Sundials type %s",type);
348:   } else {
349:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by TS Sundials",((PetscObject)viewer)->type_name);
350:   }
351:   PetscViewerASCIIPushTab(viewer);
352:   PCView(cvode->pc,viewer);
353:   PetscViewerASCIIPopTab(viewer);
354:   return(0);
355: }


358: /* --------------------------------------------------------------------------*/
362: PetscErrorCode  TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
363: {
364:   TS_Sundials *cvode = (TS_Sundials*)ts->data;
365: 
367:   cvode->cvode_type = type;
368:   return(0);
369: }

375: PetscErrorCode  TSSundialsSetGMRESRestart_Sundials(TS ts,int restart)
376: {
377:   TS_Sundials *cvode = (TS_Sundials*)ts->data;
378: 
380:   cvode->restart = restart;
381:   return(0);
382: }

388: PetscErrorCode  TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
389: {
390:   TS_Sundials *cvode = (TS_Sundials*)ts->data;
391: 
393:   cvode->linear_tol = tol;
394:   return(0);
395: }

401: PetscErrorCode  TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
402: {
403:   TS_Sundials *cvode = (TS_Sundials*)ts->data;
404: 
406:   cvode->gtype = type;
407:   return(0);
408: }

414: PetscErrorCode  TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
415: {
416:   TS_Sundials *cvode = (TS_Sundials*)ts->data;
417: 
419:   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
420:   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
421:   return(0);
422: }

428: PetscErrorCode  TSSundialsGetPC_Sundials(TS ts,PC *pc)
429: {
430:   TS_Sundials *cvode = (TS_Sundials*)ts->data;

433:   *pc = cvode->pc;
434:   return(0);
435: }

441: PetscErrorCode  TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
442: {
443:   TS_Sundials *cvode = (TS_Sundials*)ts->data;
444: 
446:   if (nonlin) *nonlin = cvode->nonlinear_solves;
447:   if (lin)    *lin    = cvode->linear_solves;
448:   return(0);
449: }
451: 
455: PetscErrorCode  TSSundialsSetExactFinalTime_Sundials(TS ts,PetscTruth s)
456: {
457:   TS_Sundials *cvode = (TS_Sundials*)ts->data;
458: 
460:   cvode->exact_final_time = s;
461:   return(0);
462: }
464: /* -------------------------------------------------------------------------------------------*/

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

471:    Not Collective

473:    Input parameters:
474: .    ts     - the time-step context

476:    Output Parameters:
477: +   nonlin - number of nonlinear iterations
478: -   lin    - number of linear iterations

480:    Level: advanced

482:    Notes:
483:     These return the number since the creation of the TS object

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

487: .seealso: TSSundialsSetType(), TSSundialsSetGMRESRestart(),
488:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
489:           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
490:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
491:           TSSundialsSetExactFinalTime()

493: @*/
494: PetscErrorCode  TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
495: {
496:   PetscErrorCode ierr,(*f)(TS,int*,int*);
497: 
499:   PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetIterations_C",(void (**)(void))&f);
500:   if (f) {
501:     (*f)(ts,nonlin,lin);
502:   }
503:   return(0);
504: }

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

511:    Collective on TS

513:    Input parameters:
514: +    ts     - the time-step context
515: -    type   - one of  SUNDIALS_ADAMS or SUNDIALS_BDF

517:    Level: intermediate

519: .keywords: Adams, backward differentiation formula

521: .seealso: TSSundialsGetIterations(),  TSSundialsSetGMRESRestart(),
522:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
523:           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
524:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
525:           TSSundialsSetExactFinalTime()
526: @*/
527: PetscErrorCode  TSSundialsSetType(TS ts,TSSundialsLmmType type)
528: {
529:   PetscErrorCode ierr,(*f)(TS,TSSundialsLmmType);
530: 
532:   PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetType_C",(void (**)(void))&f);
533:   if (f) {
534:     (*f)(ts,type);
535:   }
536:   return(0);
537: }

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

546:    Collective on TS

548:    Input parameters:
549: +    ts      - the time-step context
550: -    restart - number of direction vectors (the restart size).

552:    Level: advanced

554: .keywords: GMRES, restart

556: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), 
557:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
558:           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
559:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
560:           TSSundialsSetExactFinalTime()

562: @*/
563: PetscErrorCode  TSSundialsSetGMRESRestart(TS ts,int restart)
564: {
565:   PetscErrorCode ierr,(*f)(TS,int);

568:   PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGMRESRestart_C",(void (**)(void))&f);
569:   if (f) {
570:     (*f)(ts,restart);
571:   }
572:   return(0);
573: }

577: /*@
578:    TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
579:        system by SUNDIALS.

581:    Collective on TS

583:    Input parameters:
584: +    ts     - the time-step context
585: -    tol    - the factor by which the tolerance on the nonlinear solver is
586:              multiplied to get the tolerance on the linear solver, .05 by default.

588:    Level: advanced

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

592: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
593:           TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
594:           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
595:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
596:           TSSundialsSetExactFinalTime()

598: @*/
599: PetscErrorCode  TSSundialsSetLinearTolerance(TS ts,double tol)
600: {
601:   PetscErrorCode ierr,(*f)(TS,double);
602: 
604:   PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",(void (**)(void))&f);
605:   if (f) {
606:     (*f)(ts,tol);
607:   }
608:   return(0);
609: }

613: /*@
614:    TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
615:         in GMRES method by SUNDIALS linear solver.

617:    Collective on TS

619:    Input parameters:
620: +    ts  - the time-step context
621: -    type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS

623:    Level: advanced

625: .keywords: Sundials, orthogonalization

627: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
628:           TSSundialsSetLinearTolerance(),  TSSundialsSetTolerance(),
629:           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
630:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
631:           TSSundialsSetExactFinalTime()

633: @*/
634: PetscErrorCode  TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
635: {
636:   PetscErrorCode ierr,(*f)(TS,TSSundialsGramSchmidtType);
637: 
639:   PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",(void (**)(void))&f);
640:   if (f) {
641:     (*f)(ts,type);
642:   }
643:   return(0);
644: }

648: /*@
649:    TSSundialsSetTolerance - Sets the absolute and relative tolerance used by 
650:                          Sundials for error control.

652:    Collective on TS

654:    Input parameters:
655: +    ts  - the time-step context
656: .    aabs - the absolute tolerance  
657: -    rel - the relative tolerance

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

662:    Level: intermediate

664: .keywords: Sundials, tolerance

666: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
667:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), 
668:           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
669:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
670:           TSSundialsSetExactFinalTime()

672: @*/
673: PetscErrorCode  TSSundialsSetTolerance(TS ts,double aabs,double rel)
674: {
675:   PetscErrorCode ierr,(*f)(TS,double,double);
676: 
678:   PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetTolerance_C",(void (**)(void))&f);
679:   if (f) {
680:     (*f)(ts,aabs,rel);
681:   }
682:   return(0);
683: }

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

690:    Input Parameter:
691: .    ts - the time-step context

693:    Output Parameter:
694: .    pc - the preconditioner context

696:    Level: advanced

698: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
699:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
700:           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
701:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
702: @*/
703: PetscErrorCode  TSSundialsGetPC(TS ts,PC *pc)
704: {
705:   PetscErrorCode ierr,(*f)(TS,PC *);

708:   PetscObjectQueryFunction((PetscObject)ts,"TSSundialsGetPC_C",(void (**)(void))&f);
709:   if (f) {
710:     (*f)(ts,pc);
711:   } else {
712:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"TS must be of Sundials type to extract the PC");
713:   }

715:   return(0);
716: }

720: /*@
721:    TSSundialsSetExactFinalTime - Determines if Sundials interpolates solution to the 
722:       exact final time requested by the user or just returns it at the final time
723:       it computed. (Defaults to true).

725:    Input Parameter:
726: +   ts - the time-step context
727: -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE

729:    Level: beginner

731: .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
732:           TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
733:           TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESRestart(),
734:           TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC() 
735: @*/
736: PetscErrorCode  TSSundialsSetExactFinalTime(TS ts,PetscTruth ft)
737: {
738:   PetscErrorCode ierr,(*f)(TS,PetscTruth);

741:   PetscObjectQueryFunction((PetscObject)ts,"TSSundialsSetExactFinalTime_C",(void (**)(void))&f);
742:   if (f) {
743:     (*f)(ts,ft);
744:   }
745:   return(0);
746: }

748: /* -------------------------------------------------------------------------------------------*/
749: /*MC
750:       TS_Sundials - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)

752:    Options Database:
753: +    -ts_sundials_type <bdf,adams>
754: .    -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
755: .    -ts_sundials_atol <tol> - Absolute tolerance for convergence
756: .    -ts_sundials_rtol <tol> - Relative tolerance for convergence
757: .    -ts_sundials_linear_tolerance <tol> 
758: .    -ts_sundials_gmres_restart <restart> - Number of GMRES orthogonalization directions
759: -    -ts_sundials_not_exact_final_time -Allow SUNDIALS to stop near the final time, not exactly on it

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

764:     Level: beginner

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

769: M*/
773: PetscErrorCode  TSCreate_Sundials(TS ts)
774: {
775:   TS_Sundials *cvode;

779:   ts->ops->destroy         = TSDestroy_Sundials;
780:   ts->ops->view            = TSView_Sundials;

782:   if (ts->problem_type != TS_NONLINEAR) {
783:     SETERRQ(PETSC_ERR_SUP,"Only support for nonlinear problems");
784:   }
785:   ts->ops->setup           = TSSetUp_Sundials_Nonlinear;
786:   ts->ops->step            = TSStep_Sundials_Nonlinear;
787:   ts->ops->setfromoptions  = TSSetFromOptions_Sundials_Nonlinear;

789:   PetscNew(TS_Sundials,&cvode);
790:   PCCreate(ts->comm,&cvode->pc);
791:   PetscLogObjectParent(ts,cvode->pc);
792:   ts->data          = (void*)cvode;
793:   cvode->cvode_type = SUNDIALS_BDF;
794:   cvode->gtype      = SUNDIALS_CLASSICAL_GS;
795:   cvode->restart    = 5;
796:   cvode->linear_tol = .05;

798:   cvode->exact_final_time = PETSC_FALSE;

800:   MPI_Comm_dup(ts->comm,&(cvode->comm_sundials));
801:   /* set tolerance for Sundials */
802:   cvode->abstol = 1e-6;
803:   cvode->reltol = 1e-6;

805:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
806:                     TSSundialsSetType_Sundials);
807:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGMRESRestart_C",
808:                     "TSSundialsSetGMRESRestart_Sundials",
809:                     TSSundialsSetGMRESRestart_Sundials);
810:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
811:                     "TSSundialsSetLinearTolerance_Sundials",
812:                      TSSundialsSetLinearTolerance_Sundials);
813:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
814:                     "TSSundialsSetGramSchmidtType_Sundials",
815:                      TSSundialsSetGramSchmidtType_Sundials);
816:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
817:                     "TSSundialsSetTolerance_Sundials",
818:                      TSSundialsSetTolerance_Sundials);
819:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
820:                     "TSSundialsGetPC_Sundials",
821:                      TSSundialsGetPC_Sundials);
822:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
823:                     "TSSundialsGetIterations_Sundials",
824:                      TSSundialsGetIterations_Sundials);
825:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetExactFinalTime_C",
826:                     "TSSundialsSetExactFinalTime_Sundials",
827:                      TSSundialsSetExactFinalTime_Sundials);
828:   return(0);
829: }