Actual source code: sundials.c
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);
96: VecZeroEntries(yydot);
98: /* now compute the right hand side function */
99: TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE); CHKERRABORT(comm,ierr);
100: VecScale(yyd,-1.);
101: VecResetArray(yy); CHKERRABORT(comm,ierr);
102: VecResetArray(yyd); CHKERRABORT(comm,ierr);
103: return(0);
104: }
106: /*
107: TSStep_Sundials - Calls Sundials to integrate the ODE.
108: */
111: PetscErrorCode TSStep_Sundials(TS ts)
112: {
113: TS_Sundials *cvode = (TS_Sundials*)ts->data;
114: Vec sol = ts->vec_sol;
116: PetscInt flag;
117: long int its,nsteps;
118: realtype t,tout;
119: PetscScalar *y_data;
120: void *mem;
121: SNES snes;
122: Vec res; /* This, together with snes, will check if the SNES vec_func has been set */
125: mem = cvode->mem;
126: tout = ts->max_time;
127: VecGetArray(ts->vec_sol,&y_data);
128: N_VSetArrayPointer((realtype *)y_data,cvode->y);
129: VecRestoreArray(ts->vec_sol,PETSC_NULL);
131: /* Should think about moving this outside the loop */
132: TSGetSNES(ts, &snes);
133: SNESGetFunction(snes, &res, PETSC_NULL, PETSC_NULL);
134: if (!res) {
135: TSSetIFunction(ts, sol, PETSC_NULL, PETSC_NULL);
136: }
138: TSPreStep(ts);
140: if (cvode->monitorstep) {
141: flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
142: } else {
143: flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL);
144: }
146: if (flag){ /* display error message */
147: switch (flag){
148: case CV_ILL_INPUT:
149: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT");
150: break;
151: case CV_TOO_CLOSE:
152: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE");
153: break;
154: case CV_TOO_MUCH_WORK: {
155: PetscReal tcur;
156: CVodeGetNumSteps(mem,&nsteps);
157: CVodeGetCurrentTime(mem,&tcur);
158: 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);
159: } break;
160: case CV_TOO_MUCH_ACC:
161: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC");
162: break;
163: case CV_ERR_FAILURE:
164: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE");
165: break;
166: case CV_CONV_FAILURE:
167: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE");
168: break;
169: case CV_LINIT_FAIL:
170: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL");
171: break;
172: case CV_LSETUP_FAIL:
173: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL");
174: break;
175: case CV_LSOLVE_FAIL:
176: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL");
177: break;
178: case CV_RHSFUNC_FAIL:
179: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL");
180: break;
181: case CV_FIRST_RHSFUNC_ERR:
182: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR");
183: break;
184: case CV_REPTD_RHSFUNC_ERR:
185: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR");
186: break;
187: case CV_UNREC_RHSFUNC_ERR:
188: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR");
189: break;
190: case CV_RTFUNC_FAIL:
191: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL");
192: break;
193: default:
194: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag);
195: }
196: }
198: /* copy the solution from cvode->y to cvode->update and sol */
199: VecPlaceArray(cvode->w1,y_data);
200: VecCopy(cvode->w1,cvode->update);
201: VecResetArray(cvode->w1);
202: VecCopy(cvode->update,sol);
203: CVodeGetNumNonlinSolvIters(mem,&its);
204: CVSpilsGetNumLinIters(mem, &its);
205: ts->nonlinear_its = its; ts->linear_its = its;
207: ts->time_step = t - ts->ptime;
208: ts->ptime = t;
209: ts->steps++;
211: CVodeGetNumSteps(mem,&nsteps);
212: if (!cvode->monitorstep) ts->steps = nsteps;
213: return(0);
214: }
218: static PetscErrorCode TSInterpolate_Sundials(TS ts,PetscReal t,Vec X)
219: {
220: TS_Sundials *cvode = (TS_Sundials*)ts->data;
221: N_Vector y;
222: PetscErrorCode ierr;
223: PetscScalar *x_data;
224: PetscInt glosize,locsize;
228: /* get the vector size */
229: VecGetSize(X,&glosize);
230: VecGetLocalSize(X,&locsize);
232: /* allocate the memory for N_Vec y */
233: y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
234: if (!y) SETERRQ(PETSC_COMM_SELF,1,"Interpolated y is not allocated");
236: VecGetArray(X,&x_data);
237: N_VSetArrayPointer((realtype *)x_data,y);
238: CVodeGetDky(cvode->mem,t,0,y);
239: VecRestoreArray(X,&x_data);
241: return(0);
242: }
246: PetscErrorCode TSReset_Sundials(TS ts)
247: {
248: TS_Sundials *cvode = (TS_Sundials*)ts->data;
252: VecDestroy(&cvode->update);
253: VecDestroy(&cvode->ydot);
254: VecDestroy(&cvode->w1);
255: VecDestroy(&cvode->w2);
256: if (cvode->mem) {CVodeFree(&cvode->mem);}
257: return(0);
258: }
262: PetscErrorCode TSDestroy_Sundials(TS ts)
263: {
264: TS_Sundials *cvode = (TS_Sundials*)ts->data;
268: TSReset_Sundials(ts);
269: MPI_Comm_free(&(cvode->comm_sundials));
270: PetscFree(ts->data);
271: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","",PETSC_NULL);
272: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxl_C","",PETSC_NULL);
273: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C","",PETSC_NULL);
274: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C","",PETSC_NULL);
275: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C","",PETSC_NULL);
276: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C","",PETSC_NULL);
277: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C","",PETSC_NULL);
278: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C","",PETSC_NULL);
279: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C","",PETSC_NULL);
280: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C","",PETSC_NULL);
281: return(0);
282: }
286: PetscErrorCode TSSetUp_Sundials(TS ts)
287: {
288: TS_Sundials *cvode = (TS_Sundials*)ts->data;
290: PetscInt glosize,locsize,i,flag;
291: PetscScalar *y_data,*parray;
292: void *mem;
293: PC pc;
294: const PCType pctype;
295: PetscBool pcnone;
298: /* get the vector size */
299: VecGetSize(ts->vec_sol,&glosize);
300: VecGetLocalSize(ts->vec_sol,&locsize);
302: /* allocate the memory for N_Vec y */
303: cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
304: if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated");
306: /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
307: VecGetArray(ts->vec_sol,&parray);
308: y_data = (PetscScalar *) N_VGetArrayPointer(cvode->y);
309: for (i = 0; i < locsize; i++) y_data[i] = parray[i];
310: VecRestoreArray(ts->vec_sol,PETSC_NULL);
312: VecDuplicate(ts->vec_sol,&cvode->update);
313: VecDuplicate(ts->vec_sol,&cvode->ydot);
314: PetscLogObjectParent(ts,cvode->update);
315: PetscLogObjectParent(ts,cvode->ydot);
317: /*
318: Create work vectors for the TSPSolve_Sundials() routine. Note these are
319: allocated with zero space arrays because the actual array space is provided
320: by Sundials and set using VecPlaceArray().
321: */
322: VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w1);
323: VecCreateMPIWithArray(((PetscObject)ts)->comm,locsize,PETSC_DECIDE,0,&cvode->w2);
324: PetscLogObjectParent(ts,cvode->w1);
325: PetscLogObjectParent(ts,cvode->w2);
327: /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
328: mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
329: if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails");
330: cvode->mem = mem;
332: /* Set the pointer to user-defined data */
333: flag = CVodeSetUserData(mem, ts);
334: if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails");
336: /* Sundials may choose to use a smaller initial step, but will never use a larger step. */
337: flag = CVodeSetInitStep(mem,(realtype)ts->time_step);
338: if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetInitStep() failed");
339: if (cvode->mindt > 0) {
340: flag = CVodeSetMinStep(mem,(realtype)cvode->mindt);
341: if (flag){
342: if (flag == CV_MEM_NULL){
343: SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL");
344: } else if (flag == CV_ILL_INPUT){
345: SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size");
346: } else {
347: SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMinStep() failed");
348: }
349: }
350: }
351: if (cvode->maxdt > 0) {
352: flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt);
353: if (flag) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_LIB,"CVodeSetMaxStep() failed");
354: }
356: /* Call CVodeInit to initialize the integrator memory and specify the
357: * user's right hand side function in u'=f(t,u), the inital time T0, and
358: * the initial dependent variable vector cvode->y */
359: flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y);
360: if (flag){
361: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag);
362: }
364: /* specifies scalar relative and absolute tolerances */
365: flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol);
366: if (flag){
367: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag);
368: }
370: /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */
371: flag = CVodeSetMaxNumSteps(mem,ts->max_steps);
373: /* call CVSpgmr to use GMRES as the linear solver. */
374: /* setup the ode integrator with the given preconditioner */
375: TSSundialsGetPC(ts,&pc);
376: PCGetType(pc,&pctype);
377: PetscTypeCompare((PetscObject)pc,PCNONE,&pcnone);
378: if (pcnone){
379: flag = CVSpgmr(mem,PREC_NONE,0);
380: if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
381: } else {
382: flag = CVSpgmr(mem,PREC_LEFT,cvode->maxl);
383: if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
385: /* Set preconditioner and solve routines Precond and PSolve,
386: and the pointer to the user-defined block data */
387: flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials);
388: if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag);
389: }
391: flag = CVSpilsSetGSType(mem, MODIFIED_GS);
392: if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag);
393: return(0);
394: }
396: /* type of CVODE linear multistep method */
397: const char *TSSundialsLmmTypes[] = {"","adams","bdf","TSSundialsLmmType","SUNDIALS_",0};
398: /* type of G-S orthogonalization used by CVODE linear solver */
399: const char *TSSundialsGramSchmidtTypes[] = {"","modified","classical","TSSundialsGramSchmidtType","SUNDIALS_",0};
403: PetscErrorCode TSSetFromOptions_Sundials(TS ts)
404: {
405: TS_Sundials *cvode = (TS_Sundials*)ts->data;
407: int indx;
408: PetscBool flag;
409: PC pc;
412: PetscOptionsHead("SUNDIALS ODE solver options");
413: PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);
414: if (flag) {
415: TSSundialsSetType(ts,(TSSundialsLmmType)indx);
416: }
417: PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);
418: if (flag) {
419: TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);
420: }
421: PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,PETSC_NULL);
422: PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,PETSC_NULL);
423: PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,PETSC_NULL);
424: PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,PETSC_NULL);
425: PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,&flag);
426: PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,&flag);
427: PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internel steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,PETSC_NULL);
428: PetscOptionsTail();
429: TSSundialsGetPC(ts,&pc);
430: PCSetFromOptions(pc);
431: return(0);
432: }
436: PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
437: {
438: TS_Sundials *cvode = (TS_Sundials*)ts->data;
440: char *type;
441: char atype[] = "Adams";
442: char btype[] = "BDF: backward differentiation formula";
443: PetscBool iascii,isstring;
444: long int nsteps,its,nfevals,nlinsetups,nfails,itmp;
445: PetscInt qlast,qcur;
446: PetscReal hinused,hlast,hcur,tcur,tolsfac;
447: PC pc;
450: if (cvode->cvode_type == SUNDIALS_ADAMS) {type = atype;}
451: else {type = btype;}
453: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
454: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
455: if (iascii) {
456: PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");
457: PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);
458: PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);
459: PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);
460: PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %D\n",cvode->maxl);
461: if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
462: PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");
463: } else {
464: PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");
465: }
466: if (cvode->mindt > 0) {PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);}
467: if (cvode->maxdt > 0) {PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);}
469: /* Outputs from CVODE, CVSPILS */
470: CVodeGetTolScaleFactor(cvode->mem,&tolsfac);
471: PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);
472: CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
473: &nlinsetups,&nfails,&qlast,&qcur,
474: &hinused,&hlast,&hcur,&tcur);
475: PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);
476: PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);
477: PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);
478: PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);
480: CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);
481: PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);
482: PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);
484: CVSpilsGetNumLinIters(cvode->mem, &its); /* its = no. of calls to TSPrecond_Sundials() */
485: PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);
486: CVSpilsGetNumConvFails(cvode->mem,&itmp);
487: PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);
489: TSSundialsGetPC(ts,&pc);
490: PCView(pc,viewer);
491: CVSpilsGetNumPrecEvals(cvode->mem,&itmp);
492: PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);
493: CVSpilsGetNumPrecSolves(cvode->mem,&itmp);
494: PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);
496: CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);
497: PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);
498: CVSpilsGetNumRhsEvals(cvode->mem,&itmp);
499: PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);
500: } else if (isstring) {
501: PetscViewerStringSPrintf(viewer,"Sundials type %s",type);
502: }
503: return(0);
504: }
507: /* --------------------------------------------------------------------------*/
511: PetscErrorCode TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
512: {
513: TS_Sundials *cvode = (TS_Sundials*)ts->data;
516: cvode->cvode_type = type;
517: return(0);
518: }
524: PetscErrorCode TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl)
525: {
526: TS_Sundials *cvode = (TS_Sundials*)ts->data;
529: cvode->maxl = maxl;
530: return(0);
531: }
537: PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
538: {
539: TS_Sundials *cvode = (TS_Sundials*)ts->data;
542: cvode->linear_tol = tol;
543: return(0);
544: }
550: PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
551: {
552: TS_Sundials *cvode = (TS_Sundials*)ts->data;
555: cvode->gtype = type;
556: return(0);
557: }
563: PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
564: {
565: TS_Sundials *cvode = (TS_Sundials*)ts->data;
568: if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
569: if (rel != PETSC_DECIDE) cvode->reltol = rel;
570: return(0);
571: }
577: PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)
578: {
579: TS_Sundials *cvode = (TS_Sundials*)ts->data;
582: cvode->mindt = mindt;
583: return(0);
584: }
590: PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)
591: {
592: TS_Sundials *cvode = (TS_Sundials*)ts->data;
595: cvode->maxdt = maxdt;
596: return(0);
597: }
603: PetscErrorCode TSSundialsGetPC_Sundials(TS ts,PC *pc)
604: {
605: SNES snes;
606: KSP ksp;
607: PetscErrorCode ierr;
610: TSGetSNES(ts,&snes);
611: SNESGetKSP(snes,&ksp);
612: KSPGetPC(ksp,pc);
613: return(0);
614: }
620: PetscErrorCode TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
621: {
623: if (nonlin) *nonlin = ts->nonlinear_its;
624: if (lin) *lin = ts->linear_its;
625: return(0);
626: }
632: PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s)
633: {
634: TS_Sundials *cvode = (TS_Sundials*)ts->data;
637: cvode->monitorstep = s;
638: return(0);
639: }
641: /* -------------------------------------------------------------------------------------------*/
645: /*@C
646: TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
648: Not Collective
650: Input parameters:
651: . ts - the time-step context
653: Output Parameters:
654: + nonlin - number of nonlinear iterations
655: - lin - number of linear iterations
657: Level: advanced
659: Notes:
660: These return the number since the creation of the TS object
662: .keywords: non-linear iterations, linear iterations
664: .seealso: TSSundialsSetType(), TSSundialsSetMaxl(),
665: TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
666: TSSundialsGetIterations(), TSSundialsSetType(),
667: TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime()
669: @*/
670: PetscErrorCode TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
671: {
675: PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));
676: return(0);
677: }
681: /*@
682: TSSundialsSetType - Sets the method that Sundials will use for integration.
684: Logically Collective on TS
686: Input parameters:
687: + ts - the time-step context
688: - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF
690: Level: intermediate
692: .keywords: Adams, backward differentiation formula
694: .seealso: TSSundialsGetIterations(), TSSundialsSetMaxl(),
695: TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
696: TSSundialsGetIterations(), TSSundialsSetType(),
697: TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
698: TSSetExactFinalTime()
699: @*/
700: PetscErrorCode TSSundialsSetType(TS ts,TSSundialsLmmType type)
701: {
705: PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));
706: return(0);
707: }
711: /*@
712: TSSundialsSetMaxl - Sets the dimension of the Krylov space used by
713: GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
714: this is the maximum number of GMRES steps that will be used.
716: Logically Collective on TS
718: Input parameters:
719: + ts - the time-step context
720: - maxl - number of direction vectors (the dimension of Krylov subspace).
722: Level: advanced
724: .keywords: GMRES
726: .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
727: TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
728: TSSundialsGetIterations(), TSSundialsSetType(),
729: TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
730: TSSetExactFinalTime()
732: @*/
733: PetscErrorCode TSSundialsSetMaxl(TS ts,PetscInt maxl)
734: {
739: PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));
740: return(0);
741: }
745: /*@
746: TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
747: system by SUNDIALS.
749: Logically Collective on TS
751: Input parameters:
752: + ts - the time-step context
753: - tol - the factor by which the tolerance on the nonlinear solver is
754: multiplied to get the tolerance on the linear solver, .05 by default.
756: Level: advanced
758: .keywords: GMRES, linear convergence tolerance, SUNDIALS
760: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
761: TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
762: TSSundialsGetIterations(), TSSundialsSetType(),
763: TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
764: TSSetExactFinalTime()
766: @*/
767: PetscErrorCode TSSundialsSetLinearTolerance(TS ts,double tol)
768: {
773: PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));
774: return(0);
775: }
779: /*@
780: TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
781: in GMRES method by SUNDIALS linear solver.
783: Logically Collective on TS
785: Input parameters:
786: + ts - the time-step context
787: - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
789: Level: advanced
791: .keywords: Sundials, orthogonalization
793: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
794: TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(),
795: TSSundialsGetIterations(), TSSundialsSetType(),
796: TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
797: TSSetExactFinalTime()
799: @*/
800: PetscErrorCode TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
801: {
805: PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));
806: return(0);
807: }
811: /*@
812: TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
813: Sundials for error control.
815: Logically Collective on TS
817: Input parameters:
818: + ts - the time-step context
819: . aabs - the absolute tolerance
820: - rel - the relative tolerance
822: See the Cvode/Sundials users manual for exact details on these parameters. Essentially
823: these regulate the size of the error for a SINGLE timestep.
825: Level: intermediate
827: .keywords: Sundials, tolerance
829: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESMaxl(),
830: TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
831: TSSundialsGetIterations(), TSSundialsSetType(),
832: TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
833: TSSetExactFinalTime()
835: @*/
836: PetscErrorCode TSSundialsSetTolerance(TS ts,double aabs,double rel)
837: {
841: PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));
842: return(0);
843: }
847: /*@
848: TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
850: Input Parameter:
851: . ts - the time-step context
853: Output Parameter:
854: . pc - the preconditioner context
856: Level: advanced
858: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
859: TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
860: TSSundialsGetIterations(), TSSundialsSetType(),
861: TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
862: @*/
863: PetscErrorCode TSSundialsGetPC(TS ts,PC *pc)
864: {
868: PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC *),(ts,pc));
869: return(0);
870: }
874: /*@
875: TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller.
877: Input Parameter:
878: + ts - the time-step context
879: - mindt - lowest time step if positive, negative to deactivate
881: Note:
882: Sundials will error if it is not possible to keep the estimated truncation error below
883: the tolerance set with TSSundialsSetTolerance() without going below this step size.
885: Level: beginner
887: .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
888: @*/
889: PetscErrorCode TSSundialsSetMinTimeStep(TS ts,PetscReal mindt)
890: {
894: PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));
895: return(0);
896: }
900: /*@
901: TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller.
903: Input Parameter:
904: + ts - the time-step context
905: - maxdt - lowest time step if positive, negative to deactivate
907: Level: beginner
909: .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
910: @*/
911: PetscErrorCode TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt)
912: {
916: PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));
917: return(0);
918: }
922: /*@
923: TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false).
925: Input Parameter:
926: + ts - the time-step context
927: - ft - PETSC_TRUE if monitor, else PETSC_FALSE
929: Level: beginner
931: .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
932: TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
933: TSSundialsGetIterations(), TSSundialsSetType(),
934: TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
935: @*/
936: PetscErrorCode TSSundialsMonitorInternalSteps(TS ts,PetscBool ft)
937: {
941: PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));
942: return(0);
943: }
944: /* -------------------------------------------------------------------------------------------*/
945: /*MC
946: TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
948: Options Database:
949: + -ts_sundials_type <bdf,adams>
950: . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
951: . -ts_sundials_atol <tol> - Absolute tolerance for convergence
952: . -ts_sundials_rtol <tol> - Relative tolerance for convergence
953: . -ts_sundials_linear_tolerance <tol>
954: . -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace
955: - -ts_sundials_monitor_steps - Monitor SUNDIALS internel steps
958: Notes: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply
959: only PETSc PC options
961: Level: beginner
963: .seealso: TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetMaxl(), TSSundialsSetLinearTolerance(),
964: TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime()
966: M*/
970: PetscErrorCode TSCreate_Sundials(TS ts)
971: {
972: TS_Sundials *cvode;
974: PC pc;
977: ts->ops->reset = TSReset_Sundials;
978: ts->ops->destroy = TSDestroy_Sundials;
979: ts->ops->view = TSView_Sundials;
980: ts->ops->setup = TSSetUp_Sundials;
981: ts->ops->step = TSStep_Sundials;
982: ts->ops->interpolate = TSInterpolate_Sundials;
983: ts->ops->setfromoptions = TSSetFromOptions_Sundials;
985: PetscNewLog(ts,TS_Sundials,&cvode);
986: ts->data = (void*)cvode;
987: cvode->cvode_type = SUNDIALS_BDF;
988: cvode->gtype = SUNDIALS_CLASSICAL_GS;
989: cvode->maxl = 5;
990: cvode->linear_tol = .05;
992: cvode->monitorstep = PETSC_TRUE;
994: MPI_Comm_dup(((PetscObject)ts)->comm,&(cvode->comm_sundials));
996: cvode->mindt = -1.;
997: cvode->maxdt = -1.;
999: /* set tolerance for Sundials */
1000: cvode->reltol = 1e-6;
1001: cvode->abstol = 1e-6;
1003: /* set PCNONE as default pctype */
1004: TSSundialsGetPC_Sundials(ts,&pc);
1005: PCSetType(pc,PCNONE);
1007: if (ts->exact_final_time == PETSC_DECIDE) ts->exact_final_time = PETSC_FALSE;
1009: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetType_C","TSSundialsSetType_Sundials",
1010: TSSundialsSetType_Sundials);
1011: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxl_C",
1012: "TSSundialsSetMaxl_Sundials",
1013: TSSundialsSetMaxl_Sundials);
1014: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetLinearTolerance_C",
1015: "TSSundialsSetLinearTolerance_Sundials",
1016: TSSundialsSetLinearTolerance_Sundials);
1017: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",
1018: "TSSundialsSetGramSchmidtType_Sundials",
1019: TSSundialsSetGramSchmidtType_Sundials);
1020: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetTolerance_C",
1021: "TSSundialsSetTolerance_Sundials",
1022: TSSundialsSetTolerance_Sundials);
1023: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMinTimeStep_C",
1024: "TSSundialsSetMinTimeStep_Sundials",
1025: TSSundialsSetMinTimeStep_Sundials);
1026: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",
1027: "TSSundialsSetMaxTimeStep_Sundials",
1028: TSSundialsSetMaxTimeStep_Sundials);
1029: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetPC_C",
1030: "TSSundialsGetPC_Sundials",
1031: TSSundialsGetPC_Sundials);
1032: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsGetIterations_C",
1033: "TSSundialsGetIterations_Sundials",
1034: TSSundialsGetIterations_Sundials);
1035: PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",
1036: "TSSundialsMonitorInternalSteps_Sundials",
1037: TSSundialsMonitorInternalSteps_Sundials);
1039: return(0);
1040: }