Actual source code: petscpvode.c

  2: /*
  3:     Provides a PETSc interface to PVODE. Alan Hindmarsh's parallel ODE
  4:    solver.
  5: */

  7: #include "src/ts/impls/implicit/pvode/petscpvode.h"  /*I "petscts.h" I*/    

 11: /*@C
 12:    TSPVodeGetParameters - Extract "iopt" and "ropt" PVODE parameters

 14:    Input Parameter:
 15: .    ts - the time-step context

 17:    Output Parameters:
 18: +   opt_size - size of the parameter arrays (will be set to PVODE value OPT_SIZE)
 19: .   iopt - the integer parameters
 20: -   ropt - the double paramters


 23:    Level: advanced
 24: g
 25:    Notes: You may pass PETSC_NULL for any value you do not desire

 27:           PETSc initializes these array with the default PVODE values of 0, you may 
 28:           change them before calling the TS solver.

 30:           See the PVODE include file cvode.h for the definitions of the fields

 32:     Suggested by: Timothy William Chevalier

 34: .seealso: TSPVodeGetIterations(), TSPVodeSetType(), TSPVodeSetGMRESRestart(),
 35:           TSPVodeSetLinearTolerance(), TSPVodeSetGramSchmidtType(), TSPVodeSetTolerance(),
 36:           TSPVodeGetIterations(), TSPVodeSetType(), TSPVodeSetGMRESRestart(),
 37:           TSPVodeSetLinearTolerance(), TSPVodeSetTolerance()
 38: @*/
 39: PetscErrorCode TSPVodeGetParameters(TS ts,int *opt_size,long int *iopt[],double *ropt[])
 40: {
 41:   TS_PVode     *cvode = (TS_PVode*)ts->data;

 44:   if (opt_size) *opt_size = OPT_SIZE;
 45:   if (iopt)     *iopt     = cvode->iopt;
 46:   if (ropt)     *ropt     = cvode->ropt;
 47:   return(0);
 48: }

 50: /*
 51:       TSPrecond_PVode - function that we provide to PVODE to
 52:                         evaluate the preconditioner.

 54:     Contributed by: Liyang Xu

 56: */
 59: PetscErrorCode TSPrecond_PVode(integertype N,realtype tn,N_Vector y,N_Vector fy,booleantype jok,
 60:                     booleantype *jcurPtr,realtype _gamma,N_Vector ewt,realtype h,
 61:                     realtype uround,long int *nfePtr,void *P_data,
 62:                     N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3)
 63: {
 64:   TS           ts = (TS) P_data;
 65:   TS_PVode     *cvode = (TS_PVode*)ts->data;
 66:   PC           pc = cvode->pc;
 68:   Mat          Jac = ts->B;
 69:   Vec          tmpy = cvode->w1;
 70:   PetscScalar  one = 1.0,gm;
 71:   MatStructure str = DIFFERENT_NONZERO_PATTERN;
 72: 
 74:   /* This allows us to construct preconditioners in-place if we like */
 75:   MatSetUnfactored(Jac);

 77:   /*
 78:        jok - TRUE means reuse current Jacobian else recompute Jacobian
 79:   */
 80:   if (jok) {
 81:     MatCopy(cvode->pmat,Jac,SAME_NONZERO_PATTERN);
 82:     str      = SAME_NONZERO_PATTERN;
 83:     *jcurPtr = FALSE;
 84:   } else {
 85:     /* make PETSc vector tmpy point to PVODE vector y */
 86:     VecPlaceArray(tmpy,N_VGetData(y));

 88:     /* compute the Jacobian */
 89:     TSComputeRHSJacobian(ts,ts->ptime,tmpy,&Jac,&Jac,&str);

 91:     /* copy the Jacobian matrix */
 92:     if (!cvode->pmat) {
 93:       MatDuplicate(Jac,MAT_COPY_VALUES,&cvode->pmat);
 94:       PetscLogObjectParent(ts,cvode->pmat);
 95:     }
 96:     MatCopy(Jac,cvode->pmat,SAME_NONZERO_PATTERN);

 98:     *jcurPtr = TRUE;
 99:   }

101:   /* construct I-gamma*Jac  */
102:   gm   = -_gamma;
103:   MatScale(&gm,Jac);
104:   MatShift(&one,Jac);
105: 
106:   PCSetOperators(pc,Jac,Jac,str);
107:   return(0);
108: }

110: /*
111:      TSPSolve_PVode -  routine that we provide to PVode that applies the preconditioner.
112:       
113:     Contributed by: Liyang Xu

115: */
118: PetscErrorCode TSPSolve_PVode(integertype N,realtype tn,N_Vector y,N_Vector fy,N_Vector vtemp,
119:                    realtype _gamma,N_Vector ewt,realtype delta,long int *nfePtr,
120:                    N_Vector r,int lr,void *P_data,N_Vector z)
121: {
122:   TS       ts = (TS) P_data;
123:   TS_PVode *cvode = (TS_PVode*)ts->data;
124:   PC       pc = cvode->pc;
125:   Vec      rr = cvode->w1,xx = cvode->w2;

129:   /*
130:       Make the PETSc work vectors rr and xx point to the arrays in the PVODE vectors 
131:   */
132:   VecPlaceArray(rr,N_VGetData(r));
133:   VecPlaceArray(xx,N_VGetData(z));

135:   /* 
136:       Solve the Px=r and put the result in xx 
137:   */
138:   PCApply(pc,rr,xx);
139:   cvode->linear_solves++;


142:   return(0);
143: }

145: /*
146:         TSFunction_PVode - routine that we provide to PVode that applies the right hand side.
147:       
148:     Contributed by: Liyang Xu
149: */
152: void TSFunction_PVode(int N,double t,N_Vector y,N_Vector ydot,void *ctx)
153: {
154:   TS        ts = (TS) ctx;
155:   TS_PVode *cvode = (TS_PVode*)ts->data;
156:   Vec       tmpx = cvode->w1,tmpy = cvode->w2;

160:   /*
161:       Make the PETSc work vectors tmpx and tmpy point to the arrays in the PVODE vectors 
162:   */
163:   VecPlaceArray(tmpx,N_VGetData(y));
164:   if (ierr) {
165:     (*PetscErrorPrintf)("TSFunction_PVode:Could not place array. Error code %d",(int)ierr);
166:   }
167:   VecPlaceArray(tmpy,N_VGetData(ydot));
168:   if (ierr) {
169:     (*PetscErrorPrintf)("TSFunction_PVode:Could not place array. Error code %d",(int)ierr);
170:   }

172:   /* now compute the right hand side function */
173:   TSComputeRHSFunction(ts,t,tmpx,tmpy);
174:   if (ierr) {
175:     (*PetscErrorPrintf)("TSFunction_PVode:Could not compute RHS function. Error code %d",(int)ierr);
176:   }
177: }

179: /*
180:        TSStep_PVode_Nonlinear - Calls PVode to integrate the ODE.

182:     Contributed by: Liyang Xu
183: */
186: /* 
187:     TSStep_PVode_Nonlinear - 
188:   
189:    steps - number of time steps
190:    time - time that integrater is  terminated. 

192: */
193: PetscErrorCode TSStep_PVode_Nonlinear(TS ts,int *steps,double *time)
194: {
195:   TS_PVode  *cvode = (TS_PVode*)ts->data;
196:   Vec       sol = ts->vec_sol;
197:   int       ierr;
198:   int i,max_steps = ts->max_steps,flag;
199:   double    t,tout;
200:   realtype  *tmp;

203:   /* initialize the number of steps */
204:   *steps = -ts->steps;
205:   TSMonitor(ts,ts->steps,ts->ptime,sol);

207:   /* call CVSpgmr to use GMRES as the linear solver. */
208:   /* setup the ode integrator with the given preconditioner */
209:   CVSpgmr(cvode->mem,LEFT,cvode->gtype,cvode->restart,cvode->linear_tol,TSPrecond_PVode,TSPSolve_PVode,ts,0,0);

211:   tout = ts->max_time;
212:   for (i=0; i<max_steps; i++) {
213:     if (ts->ptime >= tout) break;
214:     VecGetArray(ts->vec_sol,&tmp);
215:     N_VSetData(tmp,cvode->y);
216:     flag = CVode(cvode->mem,tout,cvode->y,&t,ONE_STEP);
217:     cvode->nonlinear_solves += cvode->iopt[NNI];
218:     VecRestoreArray(ts->vec_sol,PETSC_NULL);
219:     if (flag != SUCCESS) SETERRQ(PETSC_ERR_LIB,"PVODE failed");

221:     if (t > tout && cvode->exact_final_time) {
222:       /* interpolate to final requested time */
223:       VecGetArray(ts->vec_sol,&tmp);
224:       N_VSetData(tmp,cvode->y);
225:       flag = CVodeDky(cvode->mem,tout,0,cvode->y);
226:       VecRestoreArray(ts->vec_sol,PETSC_NULL);
227:       if (flag != SUCCESS) SETERRQ(PETSC_ERR_LIB,"PVODE interpolation to final time failed");
228:       t = tout;
229:     }

231:     ts->time_step = t - ts->ptime;
232:     ts->ptime     = t;

234:     /*
235:        copy the solution from cvode->y to cvode->update and sol 
236:     */
237:     VecPlaceArray(cvode->w1,N_VGetData(cvode->y));
238:     VecCopy(cvode->w1,cvode->update);
239:     VecCopy(cvode->update,sol);
240: 
241:     ts->steps++;
242:     TSMonitor(ts,ts->steps,t,sol);
243:     ts->nonlinear_its = cvode->iopt[NNI];
244:     ts->linear_its    = cvode->iopt[SPGMR_NLI];
245:   }

247:   *steps           += ts->steps;
248:   *time             = t;

250:   return(0);
251: }

253: /*

255:     Contributed by: Liyang Xu
256: */
259: PetscErrorCode TSDestroy_PVode(TS ts)
260: {
261:   TS_PVode *cvode = (TS_PVode*)ts->data;

265:   if (cvode->pmat)   {MatDestroy(cvode->pmat);}
266:   if (cvode->pc)     {PCDestroy(cvode->pc);}
267:   if (cvode->update) {VecDestroy(cvode->update);}
268:   if (cvode->func)   {VecDestroy(cvode->func);}
269:   if (cvode->rhs)    {VecDestroy(cvode->rhs);}
270:   if (cvode->w1)     {VecDestroy(cvode->w1);}
271:   if (cvode->w2)     {VecDestroy(cvode->w2);}
272:   MPI_Comm_free(&(cvode->comm_pvode));
273:   PetscFree(cvode);
274:   return(0);
275: }

277: /*

279:     Contributed by: Liyang Xu
280: */
283: PetscErrorCode TSSetUp_PVode_Nonlinear(TS ts)
284: {
285:   TS_PVode    *cvode = (TS_PVode*)ts->data;
287:   int M,locsize;
288:   M_Env       machEnv;
289:   realtype    *tmp;

292:   PCSetFromOptions(cvode->pc);
293:   /* get the vector size */
294:   VecGetSize(ts->vec_sol,&M);
295:   VecGetLocalSize(ts->vec_sol,&locsize);

297:   /* allocate the memory for machEnv */
298:   /* machEnv = PVInitMPI(cvode->>comm_pvode,locsize,M);  */
299:   machEnv = M_EnvInit_Parallel(cvode->comm_pvode, locsize, M, 0, 0);


302:   /* allocate the memory for N_Vec y */
303:   cvode->y         = N_VNew(M,machEnv);
304:   VecGetArray(ts->vec_sol,&tmp);
305:   N_VSetData(tmp,cvode->y);
306:   VecRestoreArray(ts->vec_sol,PETSC_NULL);

308:   /* initializing vector update and func */
309:   VecDuplicate(ts->vec_sol,&cvode->update);
310:   VecDuplicate(ts->vec_sol,&cvode->func);
311:   PetscLogObjectParent(ts,cvode->update);
312:   PetscLogObjectParent(ts,cvode->func);

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

324:   /* allocate memory for PVode */
325:   VecGetArray(ts->vec_sol,&tmp);
326:   N_VSetData(tmp,cvode->y);
327:   cvode->mem = CVodeMalloc(M,TSFunction_PVode,ts->ptime,cvode->y,cvode->cvode_type,
328:                            NEWTON,SS,&cvode->reltol,&cvode->abstol,ts,NULL,TRUE,cvode->iopt,
329:                            cvode->ropt,machEnv);
330:   VecRestoreArray(ts->vec_sol,PETSC_NULL);
331:   return(0);
332: }

334: /*

336:     Contributed by: Liyang Xu
337: */
340: PetscErrorCode TSSetFromOptions_PVode_Nonlinear(TS ts)
341: {
342:   TS_PVode   *cvode = (TS_PVode*)ts->data;
344:   int indx;
345:   const char *btype[] = {"bdf","adams"},*otype[] = {"modified","unmodified"};
346:   PetscTruth flag;

349:   PetscOptionsHead("PVODE ODE solver options");
350:     PetscOptionsEList("-ts_pvode_type","Scheme","TSPVodeSetType",btype,2,"bdf",&indx,&flag);
351:     if (flag) {
352:       TSPVodeSetType(ts,(TSPVodeType)indx);
353:     }
354:     PetscOptionsEList("-ts_pvode_gramschmidt_type","Type of orthogonalization","TSPVodeSetGramSchmidtType",otype,2,"unmodified",&indx,&flag);
355:     if (flag) {
356:       TSPVodeSetGramSchmidtType(ts,(TSPVodeGramSchmidtType)indx);
357:     }
358:     PetscOptionsReal("-ts_pvode_atol","Absolute tolerance for convergence","TSPVodeSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);
359:     PetscOptionsReal("-ts_pvode_rtol","Relative tolerance for convergence","TSPVodeSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);
360:     PetscOptionsReal("-ts_pvode_linear_tolerance","Convergence tolerance for linear solve","TSPVodeSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);
361:     PetscOptionsInt("-ts_pvode_gmres_restart","Number of GMRES orthogonalization directions","TSPVodeSetGMRESRestart",cvode->restart,&cvode->restart,&flag);
362:     PetscOptionsName("-ts_pvode_not_exact_final_time","Allow PVODE to stop near the final time, not exactly on it","TSPVodeSetExactFinalTime",&cvode->exact_final_time);
363:   PetscOptionsTail();

365:   return(0);
366: }

368: /*

370:     Contributed by: Liyang Xu
371: */
374: PetscErrorCode TSPrintHelp_PVode(TS ts,char *p)
375: {

379:   (*PetscHelpPrintf)(ts->comm," Options for TSPVODE integrater:\n");
380:   (*PetscHelpPrintf)(ts->comm," -ts_pvode_type <bdf,adams>: integration approach\n",p);
381:   (*PetscHelpPrintf)(ts->comm," -ts_pvode_atol aabs: absolute tolerance of ODE solution\n",p);
382:   (*PetscHelpPrintf)(ts->comm," -ts_pvode_rtol rel: relative tolerance of ODE solution\n",p);
383:   (*PetscHelpPrintf)(ts->comm," -ts_pvode_gramschmidt_type <unmodified,modified>\n");
384:   (*PetscHelpPrintf)(ts->comm," -ts_pvode_gmres_restart <restart_size> (also max. GMRES its)\n");
385:   (*PetscHelpPrintf)(ts->comm," -ts_pvode_linear_tolerance <tol>\n");
386:   (*PetscHelpPrintf)(ts->comm," -ts_pvode_not_exact_final_time\n");

388:   return(0);
389: }

391: /*

393:     Contributed by: Liyang Xu
394: */
397: PetscErrorCode TSView_PVode(TS ts,PetscViewer viewer)
398: {
399:   TS_PVode   *cvode = (TS_PVode*)ts->data;
401:   char       *type;
402:   PetscTruth iascii,isstring;

405:   if (cvode->cvode_type == PVODE_ADAMS) {type = "Adams";}
406:   else                                  {type = "BDF: backward differentiation formula";}

408:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
409:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);
410:   if (iascii) {
411:     PetscViewerASCIIPrintf(viewer,"PVode integrater does not use SNES!\n");
412:     PetscViewerASCIIPrintf(viewer,"PVode integrater type %s\n",type);
413:     PetscViewerASCIIPrintf(viewer,"PVode abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);
414:     PetscViewerASCIIPrintf(viewer,"PVode linear solver tolerance factor %g\n",cvode->linear_tol);
415:     PetscViewerASCIIPrintf(viewer,"PVode GMRES max iterations (same as restart in PVODE) %D\n",cvode->restart);
416:     if (cvode->gtype == PVODE_MODIFIED_GS) {
417:       PetscViewerASCIIPrintf(viewer,"PVode using modified Gram-Schmidt for orthogonalization in GMRES\n");
418:     } else {
419:       PetscViewerASCIIPrintf(viewer,"PVode using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");
420:     }
421:   } else if (isstring) {
422:     PetscViewerStringSPrintf(viewer,"Pvode type %s",type);
423:   } else {
424:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by TS PVode",((PetscObject)viewer)->type_name);
425:   }
426:   PetscViewerASCIIPushTab(viewer);
427:   PCView(cvode->pc,viewer);
428:   PetscViewerASCIIPopTab(viewer);

430:   return(0);
431: }


434: /* --------------------------------------------------------------------------*/
438: PetscErrorCode TSPVodeSetType_PVode(TS ts,TSPVodeType type)
439: {
440:   TS_PVode *cvode = (TS_PVode*)ts->data;
441: 
443:   cvode->cvode_type = type;
444:   return(0);
445: }

451: PetscErrorCode TSPVodeSetGMRESRestart_PVode(TS ts,int restart)
452: {
453:   TS_PVode *cvode = (TS_PVode*)ts->data;
454: 
456:   cvode->restart = restart;
457:   return(0);
458: }

464: PetscErrorCode TSPVodeSetLinearTolerance_PVode(TS ts,double tol)
465: {
466:   TS_PVode *cvode = (TS_PVode*)ts->data;
467: 
469:   cvode->linear_tol = tol;
470:   return(0);
471: }

477: PetscErrorCode TSPVodeSetGramSchmidtType_PVode(TS ts,TSPVodeGramSchmidtType type)
478: {
479:   TS_PVode *cvode = (TS_PVode*)ts->data;
480: 
482:   cvode->gtype = type;

484:   return(0);
485: }

491: PetscErrorCode TSPVodeSetTolerance_PVode(TS ts,double aabs,double rel)
492: {
493:   TS_PVode *cvode = (TS_PVode*)ts->data;
494: 
496:   if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
497:   if (rel != PETSC_DECIDE)  cvode->reltol = rel;
498:   return(0);
499: }

505: PetscErrorCode TSPVodeGetPC_PVode(TS ts,PC *pc)
506: {
507:   TS_PVode *cvode = (TS_PVode*)ts->data;

510:   *pc = cvode->pc;

512:   return(0);
513: }

519: PetscErrorCode TSPVodeGetIterations_PVode(TS ts,int *nonlin,int *lin)
520: {
521:   TS_PVode *cvode = (TS_PVode*)ts->data;
522: 
524:   if (nonlin) *nonlin = cvode->nonlinear_solves;
525:   if (lin)    *lin    = cvode->linear_solves;
526:   return(0);
527: }
529: 
533: PetscErrorCode TSPVodeSetExactFinalTime_PVode(TS ts,PetscTruth s)
534: {
535:   TS_PVode *cvode = (TS_PVode*)ts->data;
536: 
538:   cvode->exact_final_time = s;
539:   return(0);
540: }
542: /* -------------------------------------------------------------------------------------------*/

546: /*@C
547:    TSPVodeGetIterations - Gets the number of nonlinear and linear iterations used so far by PVode.

549:    Not Collective

551:    Input parameters:
552: .    ts     - the time-step context

554:    Output Parameters:
555: +   nonlin - number of nonlinear iterations
556: -   lin    - number of linear iterations

558:    Level: advanced

560:    Notes:
561:     These return the number since the creation of the TS object

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

565: .seealso: TSPVodeSetType(), TSPVodeSetGMRESRestart(),
566:           TSPVodeSetLinearTolerance(), TSPVodeSetGramSchmidtType(), TSPVodeSetTolerance(),
567:           TSPVodeGetIterations(), TSPVodeSetType(), TSPVodeSetGMRESRestart(),
568:           TSPVodeSetLinearTolerance(), TSPVodeSetTolerance(), TSPVodeGetPC(),
569:           TSPVodeSetExactFinalTime()

571: @*/
572: PetscErrorCode TSPVodeGetIterations(TS ts,int *nonlin,int *lin)
573: {
574:   PetscErrorCode ierr,(*f)(TS,int*,int*);
575: 
577:   PetscObjectQueryFunction((PetscObject)ts,"TSPVodeGetIterations_C",(void (**)(void))&f);
578:   if (f) {
579:     (*f)(ts,nonlin,lin);
580:   }
581:   return(0);
582: }

586: /*@
587:    TSPVodeSetType - Sets the method that PVode will use for integration.

589:    Collective on TS

591:    Input parameters:
592: +    ts     - the time-step context
593: -    type   - one of  PVODE_ADAMS or PVODE_BDF

595:     Contributed by: Liyang Xu

597:    Level: intermediate

599: .keywords: Adams, backward differentiation formula

601: .seealso: TSPVodeGetIterations(),  TSPVodeSetGMRESRestart(),
602:           TSPVodeSetLinearTolerance(), TSPVodeSetGramSchmidtType(), TSPVodeSetTolerance(),
603:           TSPVodeGetIterations(), TSPVodeSetType(), TSPVodeSetGMRESRestart(),
604:           TSPVodeSetLinearTolerance(), TSPVodeSetTolerance(), TSPVodeGetPC(),
605:           TSPVodeSetExactFinalTime()
606: @*/
607: PetscErrorCode TSPVodeSetType(TS ts,TSPVodeType type)
608: {
609:   PetscErrorCode ierr,(*f)(TS,TSPVodeType);
610: 
612:   PetscObjectQueryFunction((PetscObject)ts,"TSPVodeSetType_C",(void (**)(void))&f);
613:   if (f) {
614:     (*f)(ts,type);
615:   }
616:   return(0);
617: }

621: /*@
622:    TSPVodeSetGMRESRestart - Sets the dimension of the Krylov space used by 
623:        GMRES in the linear solver in PVODE. PVODE DOES NOT use restarted GMRES so
624:        this is ALSO the maximum number of GMRES steps that will be used.

626:    Collective on TS

628:    Input parameters:
629: +    ts      - the time-step context
630: -    restart - number of direction vectors (the restart size).

632:    Level: advanced

634: .keywords: GMRES, restart

636: .seealso: TSPVodeGetIterations(), TSPVodeSetType(), 
637:           TSPVodeSetLinearTolerance(), TSPVodeSetGramSchmidtType(), TSPVodeSetTolerance(),
638:           TSPVodeGetIterations(), TSPVodeSetType(), TSPVodeSetGMRESRestart(),
639:           TSPVodeSetLinearTolerance(), TSPVodeSetTolerance(), TSPVodeGetPC(),
640:           TSPVodeSetExactFinalTime()

642: @*/
643: PetscErrorCode TSPVodeSetGMRESRestart(TS ts,int restart)
644: {
645:   PetscErrorCode ierr,(*f)(TS,int);

648:   PetscObjectQueryFunction((PetscObject)ts,"TSPVodeSetGMRESRestart_C",(void (**)(void))&f);
649:   if (f) {
650:     (*f)(ts,restart);
651:   }

653:   return(0);
654: }

658: /*@
659:    TSPVodeSetLinearTolerance - Sets the tolerance used to solve the linear
660:        system by PVODE.

662:    Collective on TS

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

669:    Level: advanced

671: .keywords: GMRES, linear convergence tolerance, PVODE

673: .seealso: TSPVodeGetIterations(), TSPVodeSetType(), TSPVodeSetGMRESRestart(),
674:           TSPVodeSetGramSchmidtType(), TSPVodeSetTolerance(),
675:           TSPVodeGetIterations(), TSPVodeSetType(), TSPVodeSetGMRESRestart(),
676:           TSPVodeSetLinearTolerance(), TSPVodeSetTolerance(), TSPVodeGetPC(),
677:           TSPVodeSetExactFinalTime()

679: @*/
680: PetscErrorCode TSPVodeSetLinearTolerance(TS ts,double tol)
681: {
682:   PetscErrorCode ierr,(*f)(TS,double);
683: 
685:   PetscObjectQueryFunction((PetscObject)ts,"TSPVodeSetLinearTolerance_C",(void (**)(void))&f);
686:   if (f) {
687:     (*f)(ts,tol);
688:   }
689:   return(0);
690: }

694: /*@
695:    TSPVodeSetGramSchmidtType - Sets type of orthogonalization used
696:         in GMRES method by PVODE linear solver.

698:    Collective on TS

700:    Input parameters:
701: +    ts  - the time-step context
702: -    type - either PVODE_MODIFIED_GS or PVODE_CLASSICAL_GS

704:    Level: advanced

706: .keywords: PVode, orthogonalization

708: .seealso: TSPVodeGetIterations(), TSPVodeSetType(), TSPVodeSetGMRESRestart(),
709:           TSPVodeSetLinearTolerance(),  TSPVodeSetTolerance(),
710:           TSPVodeGetIterations(), TSPVodeSetType(), TSPVodeSetGMRESRestart(),
711:           TSPVodeSetLinearTolerance(), TSPVodeSetTolerance(), TSPVodeGetPC(),
712:           TSPVodeSetExactFinalTime()

714: @*/
715: PetscErrorCode TSPVodeSetGramSchmidtType(TS ts,TSPVodeGramSchmidtType type)
716: {
717:   PetscErrorCode ierr,(*f)(TS,TSPVodeGramSchmidtType);
718: 
720:   PetscObjectQueryFunction((PetscObject)ts,"TSPVodeSetGramSchmidtType_C",(void (**)(void))&f);
721:   if (f) {
722:     (*f)(ts,type);
723:   }
724:   return(0);
725: }

729: /*@
730:    TSPVodeSetTolerance - Sets the absolute and relative tolerance used by 
731:                          PVode for error control.

733:    Collective on TS

735:    Input parameters:
736: +    ts  - the time-step context
737: .    aabs - the absolute tolerance  
738: -    rel - the relative tolerance

740:     Contributed by: Liyang Xu

742:      See the Cvode/Pvode users manual for exact details on these parameters. Essentially
743:     these regulate the size of the error for a SINGLE timestep.

745:    Level: intermediate

747: .keywords: PVode, tolerance

749: .seealso: TSPVodeGetIterations(), TSPVodeSetType(), TSPVodeSetGMRESRestart(),
750:           TSPVodeSetLinearTolerance(), TSPVodeSetGramSchmidtType(), 
751:           TSPVodeGetIterations(), TSPVodeSetType(), TSPVodeSetGMRESRestart(),
752:           TSPVodeSetLinearTolerance(), TSPVodeSetTolerance(), TSPVodeGetPC(),
753:           TSPVodeSetExactFinalTime()

755: @*/
756: PetscErrorCode TSPVodeSetTolerance(TS ts,double aabs,double rel)
757: {
758:   PetscErrorCode ierr,(*f)(TS,double,double);
759: 
761:   PetscObjectQueryFunction((PetscObject)ts,"TSPVodeSetTolerance_C",(void (**)(void))&f);
762:   if (f) {
763:     (*f)(ts,aabs,rel);
764:   }
765:   return(0);
766: }

770: /*@
771:    TSPVodeGetPC - Extract the PC context from a time-step context for PVode.

773:    Input Parameter:
774: .    ts - the time-step context

776:    Output Parameter:
777: .    pc - the preconditioner context

779:    Level: advanced

781:     Contributed by: Liyang Xu

783: .seealso: TSPVodeGetIterations(), TSPVodeSetType(), TSPVodeSetGMRESRestart(),
784:           TSPVodeSetLinearTolerance(), TSPVodeSetGramSchmidtType(), TSPVodeSetTolerance(),
785:           TSPVodeGetIterations(), TSPVodeSetType(), TSPVodeSetGMRESRestart(),
786:           TSPVodeSetLinearTolerance(), TSPVodeSetTolerance()
787: @*/
788: PetscErrorCode TSPVodeGetPC(TS ts,PC *pc)
789: {
790:   PetscErrorCode ierr,(*f)(TS,PC *);

793:   PetscObjectQueryFunction((PetscObject)ts,"TSPVodeGetPC_C",(void (**)(void))&f);
794:   if (f) {
795:     (*f)(ts,pc);
796:   } else {
797:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"TS must be of PVode type to extract the PC");
798:   }

800:   return(0);
801: }

805: /*@
806:    TSPVodeSetExactFinalTime - Determines if PVode interpolates solution to the 
807:       exact final time requested by the user or just returns it at the final time
808:       it computed. (Defaults to true).

810:    Input Parameter:
811: +   ts - the time-step context
812: -   ft - PETSC_TRUE if interpolates, else PETSC_FALSE

814:    Level: beginner

816: .seealso:TSPVodeGetIterations(), TSPVodeSetType(), TSPVodeSetGMRESRestart(),
817:           TSPVodeSetLinearTolerance(), TSPVodeSetGramSchmidtType(), TSPVodeSetTolerance(),
818:           TSPVodeGetIterations(), TSPVodeSetType(), TSPVodeSetGMRESRestart(),
819:           TSPVodeSetLinearTolerance(), TSPVodeSetTolerance(), TSPVodeGetPC() 
820: @*/
821: PetscErrorCode TSPVodeSetExactFinalTime(TS ts,PetscTruth ft)
822: {
823:   PetscErrorCode ierr,(*f)(TS,PetscTruth);

826:   PetscObjectQueryFunction((PetscObject)ts,"TSPVodeSetExactFinalTime_C",(void (**)(void))&f);
827:   if (f) {
828:     (*f)(ts,ft);
829:   }

831:   return(0);
832: }

834: /* -------------------------------------------------------------------------------------------*/
835: /*MC
836:       TS_PVode - ODE solver using the LLNL CVODE/PVODE package (now called SUNDIALS)

838:    Options Database:
839: +    -ts_pvode_type <bdf,adams>
840: .    -ts_pvode_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
841: .    -ts_pvode_atol <tol> - Absolute tolerance for convergence
842: .    -ts_pvode_rtol <tol> - Relative tolerance for convergence
843: .    -ts_pvode_linear_tolerance <tol> 
844: .    -ts_pvode_gmres_restart <restart> - Number of GMRES orthogonalization directions
845: -    -ts_pvode_not_exact_final_time -Allow PVODE to stop near the final time, not exactly on it

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

850:     Contributed by: Liyang Xu

852:     Level: beginner

854: .seealso:  TSCreate(), TS, TSSetType(), TSPVodeSetType(), TSPVodeSetGMRESRestart(), TSPVodeSetLinearTolerance(),
855:            TSPVodeSetGramSchmidtType(), TSPVodeSetTolerance(), TSPVodeGetPC(), TSPVodeGetIterations(), TSPVodeSetExactFinalTime()

857: M*/
861: PetscErrorCode TSCreate_PVode(TS ts)
862: {
863:   TS_PVode *cvode;

867:   ts->ops->destroy         = TSDestroy_PVode;
868:   ts->ops->view            = TSView_PVode;

870:   if (ts->problem_type != TS_NONLINEAR) {
871:     SETERRQ(PETSC_ERR_SUP,"Only support for nonlinear problems");
872:   }
873:   ts->ops->setup           = TSSetUp_PVode_Nonlinear;
874:   ts->ops->step            = TSStep_PVode_Nonlinear;
875:   ts->ops->setfromoptions  = TSSetFromOptions_PVode_Nonlinear;

877:   PetscNew(TS_PVode,&cvode);
878:   PetscMemzero(cvode,sizeof(TS_PVode));
879:   PCCreate(ts->comm,&cvode->pc);
880:   PetscLogObjectParent(ts,cvode->pc);
881:   ts->data          = (void*)cvode;
882:   cvode->cvode_type = BDF;
883:   cvode->gtype      = PVODE_UNMODIFIED_GS;
884:   cvode->restart    = 5;
885:   cvode->linear_tol = .05;

887:   cvode->exact_final_time = PETSC_TRUE;

889:   MPI_Comm_dup(ts->comm,&(cvode->comm_pvode));
890:   /* set tolerance for PVode */
891:   cvode->abstol = 1e-6;
892:   cvode->reltol = 1e-6;

894:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPVodeSetType_C","TSPVodeSetType_PVode",
895:                     TSPVodeSetType_PVode);
896:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPVodeSetGMRESRestart_C",
897:                     "TSPVodeSetGMRESRestart_PVode",
898:                     TSPVodeSetGMRESRestart_PVode);
899:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPVodeSetLinearTolerance_C",
900:                     "TSPVodeSetLinearTolerance_PVode",
901:                      TSPVodeSetLinearTolerance_PVode);
902:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPVodeSetGramSchmidtType_C",
903:                     "TSPVodeSetGramSchmidtType_PVode",
904:                      TSPVodeSetGramSchmidtType_PVode);
905:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPVodeSetTolerance_C",
906:                     "TSPVodeSetTolerance_PVode",
907:                      TSPVodeSetTolerance_PVode);
908:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPVodeGetPC_C",
909:                     "TSPVodeGetPC_PVode",
910:                      TSPVodeGetPC_PVode);
911:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPVodeGetIterations_C",
912:                     "TSPVodeGetIterations_PVode",
913:                      TSPVodeGetIterations_PVode);
914:   PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSPVodeSetExactFinalTime_C",
915:                     "TSPVodeSetExactFinalTime_PVode",
916:                      TSPVodeSetExactFinalTime_PVode);
917:   return(0);
918: }