Actual source code: sundials.c
petsc-3.11.4 2019-09-28
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 <../src/ts/impls/implicit/sundials/sundials.h>
10: /*
11: TSPrecond_Sundials - function that we provide to SUNDIALS to
12: evaluate the preconditioner.
13: */
14: PetscErrorCode TSPrecond_Sundials(realtype tn,N_Vector y,N_Vector fy,booleantype jok,booleantype *jcurPtr,
15: realtype _gamma,void *P_data,N_Vector vtemp1,N_Vector vtemp2,N_Vector vtemp3)
16: {
17: TS ts = (TS) P_data;
18: TS_Sundials *cvode = (TS_Sundials*)ts->data;
19: PC pc;
21: Mat J,P;
22: Vec yy = cvode->w1,yydot = cvode->ydot;
23: PetscReal gm = (PetscReal)_gamma;
24: PetscScalar *y_data;
27: TSGetIJacobian(ts,&J,&P,NULL,NULL);
28: y_data = (PetscScalar*) N_VGetArrayPointer(y);
29: VecPlaceArray(yy,y_data);
30: VecZeroEntries(yydot); /* The Jacobian is independent of Ydot for ODE which is all that CVode works for */
31: /* compute the shifted Jacobian (1/gm)*I + Jrest */
32: TSComputeIJacobian(ts,ts->ptime,yy,yydot,1/gm,J,P,PETSC_FALSE);
33: VecResetArray(yy);
34: MatScale(P,gm); /* turn into I-gm*Jrest, J is not used by Sundials */
35: *jcurPtr = TRUE;
36: TSSundialsGetPC(ts,&pc);
37: PCSetOperators(pc,J,P);
38: return(0);
39: }
41: /*
42: TSPSolve_Sundials - routine that we provide to Sundials that applies the preconditioner.
43: */
44: PetscErrorCode TSPSolve_Sundials(realtype tn,N_Vector y,N_Vector fy,N_Vector r,N_Vector z,
45: realtype _gamma,realtype delta,int lr,void *P_data,N_Vector vtemp)
46: {
47: TS ts = (TS) P_data;
48: TS_Sundials *cvode = (TS_Sundials*)ts->data;
49: PC pc;
50: Vec rr = cvode->w1,zz = cvode->w2;
52: PetscScalar *r_data,*z_data;
55: /* Make the PETSc work vectors rr and zz point to the arrays in the SUNDIALS vectors r and z respectively*/
56: r_data = (PetscScalar*) N_VGetArrayPointer(r);
57: z_data = (PetscScalar*) N_VGetArrayPointer(z);
58: VecPlaceArray(rr,r_data);
59: VecPlaceArray(zz,z_data);
61: /* Solve the Px=r and put the result in zz */
62: TSSundialsGetPC(ts,&pc);
63: PCApply(pc,rr,zz);
64: VecResetArray(rr);
65: VecResetArray(zz);
66: return(0);
67: }
69: /*
70: TSFunction_Sundials - routine that we provide to Sundials that applies the right hand side.
71: */
72: int TSFunction_Sundials(realtype t,N_Vector y,N_Vector ydot,void *ctx)
73: {
74: TS ts = (TS) ctx;
75: DM dm;
76: DMTS tsdm;
77: TSIFunction ifunction;
78: MPI_Comm comm;
79: TS_Sundials *cvode = (TS_Sundials*)ts->data;
80: Vec yy = cvode->w1,yyd = cvode->w2,yydot = cvode->ydot;
81: PetscScalar *y_data,*ydot_data;
85: PetscObjectGetComm((PetscObject)ts,&comm);
86: /* Make the PETSc work vectors yy and yyd point to the arrays in the SUNDIALS vectors y and ydot respectively*/
87: y_data = (PetscScalar*) N_VGetArrayPointer(y);
88: ydot_data = (PetscScalar*) N_VGetArrayPointer(ydot);
89: VecPlaceArray(yy,y_data);CHKERRABORT(comm,ierr);
90: VecPlaceArray(yyd,ydot_data);CHKERRABORT(comm,ierr);
92: /* Now compute the right hand side function, via IFunction unless only the more efficient RHSFunction is set */
93: TSGetDM(ts,&dm);
94: DMGetDMTS(dm,&tsdm);
95: DMTSGetIFunction(dm,&ifunction,NULL);
96: if (!ifunction) {
97: TSComputeRHSFunction(ts,t,yy,yyd);
98: } else { /* If rhsfunction is also set, this computes both parts and shifts them to the right */
99: VecZeroEntries(yydot);
100: TSComputeIFunction(ts,t,yy,yydot,yyd,PETSC_FALSE);CHKERRABORT(comm,ierr);
101: VecScale(yyd,-1.);
102: }
103: VecResetArray(yy);CHKERRABORT(comm,ierr);
104: VecResetArray(yyd);CHKERRABORT(comm,ierr);
105: return(0);
106: }
108: /*
109: TSStep_Sundials - Calls Sundials to integrate the ODE.
110: */
111: PetscErrorCode TSStep_Sundials(TS ts)
112: {
113: TS_Sundials *cvode = (TS_Sundials*)ts->data;
115: PetscInt flag;
116: long int nits,lits,nsteps;
117: realtype t,tout;
118: PetscScalar *y_data;
119: void *mem;
122: mem = cvode->mem;
123: tout = ts->max_time;
124: VecGetArray(ts->vec_sol,&y_data);
125: N_VSetArrayPointer((realtype*)y_data,cvode->y);
126: VecRestoreArray(ts->vec_sol,NULL);
128: /* We would like to TSPreStage() and TSPostStage()
129: * before each stage solve but CVode does not appear to support this. */
130: if (cvode->monitorstep)
131: flag = CVode(mem,tout,cvode->y,&t,CV_ONE_STEP);
132: else
133: flag = CVode(mem,tout,cvode->y,&t,CV_NORMAL);
135: if (flag) { /* display error message */
136: switch (flag) {
137: case CV_ILL_INPUT:
138: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ILL_INPUT");
139: break;
140: case CV_TOO_CLOSE:
141: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_CLOSE");
142: break;
143: case CV_TOO_MUCH_WORK: {
144: PetscReal tcur;
145: CVodeGetNumSteps(mem,&nsteps);
146: CVodeGetCurrentTime(mem,&tcur);
147: SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_WORK. At t=%g, nsteps %D exceeds maxstep %D. Increase '-ts_max_steps <>' or modify TSSetMaxSteps()",(double)tcur,nsteps,ts->max_steps);
148: } break;
149: case CV_TOO_MUCH_ACC:
150: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_TOO_MUCH_ACC");
151: break;
152: case CV_ERR_FAILURE:
153: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_ERR_FAILURE");
154: break;
155: case CV_CONV_FAILURE:
156: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_CONV_FAILURE");
157: break;
158: case CV_LINIT_FAIL:
159: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LINIT_FAIL");
160: break;
161: case CV_LSETUP_FAIL:
162: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSETUP_FAIL");
163: break;
164: case CV_LSOLVE_FAIL:
165: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_LSOLVE_FAIL");
166: break;
167: case CV_RHSFUNC_FAIL:
168: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RHSFUNC_FAIL");
169: break;
170: case CV_FIRST_RHSFUNC_ERR:
171: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_FIRST_RHSFUNC_ERR");
172: break;
173: case CV_REPTD_RHSFUNC_ERR:
174: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_REPTD_RHSFUNC_ERR");
175: break;
176: case CV_UNREC_RHSFUNC_ERR:
177: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_UNREC_RHSFUNC_ERR");
178: break;
179: case CV_RTFUNC_FAIL:
180: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, CV_RTFUNC_FAIL");
181: break;
182: default:
183: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVode() fails, flag %d",flag);
184: }
185: }
187: /* log inner nonlinear and linear iterations */
188: CVodeGetNumNonlinSolvIters(mem,&nits);
189: CVSpilsGetNumLinIters(mem,&lits);
190: ts->snes_its += nits; ts->ksp_its = lits;
192: /* copy the solution from cvode->y to cvode->update and sol */
193: VecPlaceArray(cvode->w1,y_data);
194: VecCopy(cvode->w1,cvode->update);
195: VecResetArray(cvode->w1);
196: VecCopy(cvode->update,ts->vec_sol);
198: ts->time_step = t - ts->ptime;
199: ts->ptime = t;
201: CVodeGetNumSteps(mem,&nsteps);
202: if (!cvode->monitorstep) ts->steps += nsteps - 1; /* TSStep() increments the step counter by one */
203: return(0);
204: }
206: static PetscErrorCode TSInterpolate_Sundials(TS ts,PetscReal t,Vec X)
207: {
208: TS_Sundials *cvode = (TS_Sundials*)ts->data;
209: N_Vector y;
211: PetscScalar *x_data;
212: PetscInt glosize,locsize;
215: /* get the vector size */
216: VecGetSize(X,&glosize);
217: VecGetLocalSize(X,&locsize);
219: /* allocate the memory for N_Vec y */
220: y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
221: if (!y) SETERRQ(PETSC_COMM_SELF,1,"Interpolated y is not allocated");
223: VecGetArray(X,&x_data);
224: N_VSetArrayPointer((realtype*)x_data,y);
225: CVodeGetDky(cvode->mem,t,0,y);
226: VecRestoreArray(X,&x_data);
227: return(0);
228: }
230: PetscErrorCode TSReset_Sundials(TS ts)
231: {
232: TS_Sundials *cvode = (TS_Sundials*)ts->data;
236: VecDestroy(&cvode->update);
237: VecDestroy(&cvode->ydot);
238: VecDestroy(&cvode->w1);
239: VecDestroy(&cvode->w2);
240: if (cvode->mem) CVodeFree(&cvode->mem);
241: return(0);
242: }
244: PetscErrorCode TSDestroy_Sundials(TS ts)
245: {
246: TS_Sundials *cvode = (TS_Sundials*)ts->data;
250: TSReset_Sundials(ts);
251: MPI_Comm_free(&(cvode->comm_sundials));
252: PetscFree(ts->data);
253: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",NULL);
254: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",NULL);
255: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",NULL);
256: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",NULL);
257: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",NULL);
258: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",NULL);
259: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",NULL);
260: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",NULL);
261: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",NULL);
262: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",NULL);
263: return(0);
264: }
266: PetscErrorCode TSSetUp_Sundials(TS ts)
267: {
268: TS_Sundials *cvode = (TS_Sundials*)ts->data;
270: PetscInt glosize,locsize,i,flag;
271: PetscScalar *y_data,*parray;
272: void *mem;
273: PC pc;
274: PCType pctype;
275: PetscBool pcnone;
278: if (ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for exact final time option 'MATCHSTEP' when using Sundials");
280: /* get the vector size */
281: VecGetSize(ts->vec_sol,&glosize);
282: VecGetLocalSize(ts->vec_sol,&locsize);
284: /* allocate the memory for N_Vec y */
285: cvode->y = N_VNew_Parallel(cvode->comm_sundials,locsize,glosize);
286: if (!cvode->y) SETERRQ(PETSC_COMM_SELF,1,"cvode->y is not allocated");
288: /* initialize N_Vec y: copy ts->vec_sol to cvode->y */
289: VecGetArray(ts->vec_sol,&parray);
290: y_data = (PetscScalar*) N_VGetArrayPointer(cvode->y);
291: for (i = 0; i < locsize; i++) y_data[i] = parray[i];
292: VecRestoreArray(ts->vec_sol,NULL);
294: VecDuplicate(ts->vec_sol,&cvode->update);
295: VecDuplicate(ts->vec_sol,&cvode->ydot);
296: PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->update);
297: PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->ydot);
299: /*
300: Create work vectors for the TSPSolve_Sundials() routine. Note these are
301: allocated with zero space arrays because the actual array space is provided
302: by Sundials and set using VecPlaceArray().
303: */
304: VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,0,&cvode->w1);
305: VecCreateMPIWithArray(PetscObjectComm((PetscObject)ts),1,locsize,PETSC_DECIDE,0,&cvode->w2);
306: PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w1);
307: PetscLogObjectParent((PetscObject)ts,(PetscObject)cvode->w2);
309: /* Call CVodeCreate to create the solver memory and the use of a Newton iteration */
310: mem = CVodeCreate(cvode->cvode_type, CV_NEWTON);
311: if (!mem) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"CVodeCreate() fails");
312: cvode->mem = mem;
314: /* Set the pointer to user-defined data */
315: flag = CVodeSetUserData(mem, ts);
316: if (flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetUserData() fails");
318: /* Sundials may choose to use a smaller initial step, but will never use a larger step. */
319: flag = CVodeSetInitStep(mem,(realtype)ts->time_step);
320: if (flag) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetInitStep() failed");
321: if (cvode->mindt > 0) {
322: flag = CVodeSetMinStep(mem,(realtype)cvode->mindt);
323: if (flag) {
324: if (flag == CV_MEM_NULL) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, cvode_mem pointer is NULL");
325: else if (flag == CV_ILL_INPUT) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed, hmin is nonpositive or it exceeds the maximum allowable step size");
326: else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMinStep() failed");
327: }
328: }
329: if (cvode->maxdt > 0) {
330: flag = CVodeSetMaxStep(mem,(realtype)cvode->maxdt);
331: if (flag) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_LIB,"CVodeSetMaxStep() failed");
332: }
334: /* Call CVodeInit to initialize the integrator memory and specify the
335: * user's right hand side function in u'=f(t,u), the inital time T0, and
336: * the initial dependent variable vector cvode->y */
337: flag = CVodeInit(mem,TSFunction_Sundials,ts->ptime,cvode->y);
338: if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeInit() fails, flag %d",flag);
340: /* specifies scalar relative and absolute tolerances */
341: flag = CVodeSStolerances(mem,cvode->reltol,cvode->abstol);
342: if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSStolerances() fails, flag %d",flag);
344: /* Specify max num of steps to be taken by cvode in its attempt to reach the next output time */
345: flag = CVodeSetMaxNumSteps(mem,ts->max_steps);
346: if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVodeSetMaxNumSteps() fails, flag %d",flag);
348: /* call CVSpgmr to use GMRES as the linear solver. */
349: /* setup the ode integrator with the given preconditioner */
350: TSSundialsGetPC(ts,&pc);
351: PCGetType(pc,&pctype);
352: PetscObjectTypeCompare((PetscObject)pc,PCNONE,&pcnone);
353: if (pcnone) {
354: flag = CVSpgmr(mem,PREC_NONE,0);
355: if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
356: } else {
357: flag = CVSpgmr(mem,PREC_LEFT,cvode->maxl);
358: if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmr() fails, flag %d",flag);
360: /* Set preconditioner and solve routines Precond and PSolve,
361: and the pointer to the user-defined block data */
362: flag = CVSpilsSetPreconditioner(mem,TSPrecond_Sundials,TSPSolve_Sundials);
363: if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpilsSetPreconditioner() fails, flag %d", flag);
364: }
366: flag = CVSpilsSetGSType(mem, MODIFIED_GS);
367: if (flag) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CVSpgmrSetGSType() fails, flag %d",flag);
368: return(0);
369: }
371: /* type of CVODE linear multistep method */
372: const char *const TSSundialsLmmTypes[] = {"","ADAMS","BDF","TSSundialsLmmType","SUNDIALS_",0};
373: /* type of G-S orthogonalization used by CVODE linear solver */
374: const char *const TSSundialsGramSchmidtTypes[] = {"","MODIFIED","CLASSICAL","TSSundialsGramSchmidtType","SUNDIALS_",0};
376: PetscErrorCode TSSetFromOptions_Sundials(PetscOptionItems *PetscOptionsObject,TS ts)
377: {
378: TS_Sundials *cvode = (TS_Sundials*)ts->data;
380: int indx;
381: PetscBool flag;
382: PC pc;
385: PetscOptionsHead(PetscOptionsObject,"SUNDIALS ODE solver options");
386: PetscOptionsEList("-ts_sundials_type","Scheme","TSSundialsSetType",TSSundialsLmmTypes,3,TSSundialsLmmTypes[cvode->cvode_type],&indx,&flag);
387: if (flag) {
388: TSSundialsSetType(ts,(TSSundialsLmmType)indx);
389: }
390: PetscOptionsEList("-ts_sundials_gramschmidt_type","Type of orthogonalization","TSSundialsSetGramSchmidtType",TSSundialsGramSchmidtTypes,3,TSSundialsGramSchmidtTypes[cvode->gtype],&indx,&flag);
391: if (flag) {
392: TSSundialsSetGramSchmidtType(ts,(TSSundialsGramSchmidtType)indx);
393: }
394: PetscOptionsReal("-ts_sundials_atol","Absolute tolerance for convergence","TSSundialsSetTolerance",cvode->abstol,&cvode->abstol,NULL);
395: PetscOptionsReal("-ts_sundials_rtol","Relative tolerance for convergence","TSSundialsSetTolerance",cvode->reltol,&cvode->reltol,NULL);
396: PetscOptionsReal("-ts_sundials_mindt","Minimum step size","TSSundialsSetMinTimeStep",cvode->mindt,&cvode->mindt,NULL);
397: PetscOptionsReal("-ts_sundials_maxdt","Maximum step size","TSSundialsSetMaxTimeStep",cvode->maxdt,&cvode->maxdt,NULL);
398: PetscOptionsReal("-ts_sundials_linear_tolerance","Convergence tolerance for linear solve","TSSundialsSetLinearTolerance",cvode->linear_tol,&cvode->linear_tol,NULL);
399: PetscOptionsInt("-ts_sundials_maxl","Max dimension of the Krylov subspace","TSSundialsSetMaxl",cvode->maxl,&cvode->maxl,NULL);
400: PetscOptionsBool("-ts_sundials_monitor_steps","Monitor SUNDIALS internal steps","TSSundialsMonitorInternalSteps",cvode->monitorstep,&cvode->monitorstep,NULL);
401: PetscOptionsTail();
402: TSSundialsGetPC(ts,&pc);
403: PCSetFromOptions(pc);
404: return(0);
405: }
407: PetscErrorCode TSView_Sundials(TS ts,PetscViewer viewer)
408: {
409: TS_Sundials *cvode = (TS_Sundials*)ts->data;
411: char *type;
412: char atype[] = "Adams";
413: char btype[] = "BDF: backward differentiation formula";
414: PetscBool iascii,isstring;
415: long int nsteps,its,nfevals,nlinsetups,nfails,itmp;
416: PetscInt qlast,qcur;
417: PetscReal hinused,hlast,hcur,tcur,tolsfac;
418: PC pc;
421: if (cvode->cvode_type == SUNDIALS_ADAMS) type = atype;
422: else type = btype;
424: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
425: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
426: if (iascii) {
427: PetscViewerASCIIPrintf(viewer,"Sundials integrater does not use SNES!\n");
428: PetscViewerASCIIPrintf(viewer,"Sundials integrater type %s\n",type);
429: PetscViewerASCIIPrintf(viewer,"Sundials abs tol %g rel tol %g\n",cvode->abstol,cvode->reltol);
430: PetscViewerASCIIPrintf(viewer,"Sundials linear solver tolerance factor %g\n",cvode->linear_tol);
431: PetscViewerASCIIPrintf(viewer,"Sundials max dimension of Krylov subspace %D\n",cvode->maxl);
432: if (cvode->gtype == SUNDIALS_MODIFIED_GS) {
433: PetscViewerASCIIPrintf(viewer,"Sundials using modified Gram-Schmidt for orthogonalization in GMRES\n");
434: } else {
435: PetscViewerASCIIPrintf(viewer,"Sundials using unmodified (classical) Gram-Schmidt for orthogonalization in GMRES\n");
436: }
437: if (cvode->mindt > 0) {PetscViewerASCIIPrintf(viewer,"Sundials minimum time step %g\n",cvode->mindt);}
438: if (cvode->maxdt > 0) {PetscViewerASCIIPrintf(viewer,"Sundials maximum time step %g\n",cvode->maxdt);}
440: /* Outputs from CVODE, CVSPILS */
441: CVodeGetTolScaleFactor(cvode->mem,&tolsfac);
442: PetscViewerASCIIPrintf(viewer,"Sundials suggested factor for tolerance scaling %g\n",tolsfac);
443: CVodeGetIntegratorStats(cvode->mem,&nsteps,&nfevals,
444: &nlinsetups,&nfails,&qlast,&qcur,
445: &hinused,&hlast,&hcur,&tcur);
446: PetscViewerASCIIPrintf(viewer,"Sundials cumulative number of internal steps %D\n",nsteps);
447: PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to rhs function %D\n",nfevals);
448: PetscViewerASCIIPrintf(viewer,"Sundials no. of calls to linear solver setup function %D\n",nlinsetups);
449: PetscViewerASCIIPrintf(viewer,"Sundials no. of error test failures %D\n",nfails);
451: CVodeGetNonlinSolvStats(cvode->mem,&its,&nfails);
452: PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear solver iterations %D\n",its);
453: PetscViewerASCIIPrintf(viewer,"Sundials no. of nonlinear convergence failure %D\n",nfails);
455: CVSpilsGetNumLinIters(cvode->mem, &its); /* its = no. of calls to TSPrecond_Sundials() */
456: PetscViewerASCIIPrintf(viewer,"Sundials no. of linear iterations %D\n",its);
457: CVSpilsGetNumConvFails(cvode->mem,&itmp);
458: PetscViewerASCIIPrintf(viewer,"Sundials no. of linear convergence failures %D\n",itmp);
460: TSSundialsGetPC(ts,&pc);
461: PCView(pc,viewer);
462: CVSpilsGetNumPrecEvals(cvode->mem,&itmp);
463: PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner evaluations %D\n",itmp);
464: CVSpilsGetNumPrecSolves(cvode->mem,&itmp);
465: PetscViewerASCIIPrintf(viewer,"Sundials no. of preconditioner solves %D\n",itmp);
467: CVSpilsGetNumJtimesEvals(cvode->mem,&itmp);
468: PetscViewerASCIIPrintf(viewer,"Sundials no. of Jacobian-vector product evaluations %D\n",itmp);
469: CVSpilsGetNumRhsEvals(cvode->mem,&itmp);
470: PetscViewerASCIIPrintf(viewer,"Sundials no. of rhs calls for finite diff. Jacobian-vector evals %D\n",itmp);
471: } else if (isstring) {
472: PetscViewerStringSPrintf(viewer,"Sundials type %s",type);
473: }
474: return(0);
475: }
478: /* --------------------------------------------------------------------------*/
479: PetscErrorCode TSSundialsSetType_Sundials(TS ts,TSSundialsLmmType type)
480: {
481: TS_Sundials *cvode = (TS_Sundials*)ts->data;
484: cvode->cvode_type = type;
485: return(0);
486: }
488: PetscErrorCode TSSundialsSetMaxl_Sundials(TS ts,PetscInt maxl)
489: {
490: TS_Sundials *cvode = (TS_Sundials*)ts->data;
493: cvode->maxl = maxl;
494: return(0);
495: }
497: PetscErrorCode TSSundialsSetLinearTolerance_Sundials(TS ts,double tol)
498: {
499: TS_Sundials *cvode = (TS_Sundials*)ts->data;
502: cvode->linear_tol = tol;
503: return(0);
504: }
506: PetscErrorCode TSSundialsSetGramSchmidtType_Sundials(TS ts,TSSundialsGramSchmidtType type)
507: {
508: TS_Sundials *cvode = (TS_Sundials*)ts->data;
511: cvode->gtype = type;
512: return(0);
513: }
515: PetscErrorCode TSSundialsSetTolerance_Sundials(TS ts,double aabs,double rel)
516: {
517: TS_Sundials *cvode = (TS_Sundials*)ts->data;
520: if (aabs != PETSC_DECIDE) cvode->abstol = aabs;
521: if (rel != PETSC_DECIDE) cvode->reltol = rel;
522: return(0);
523: }
525: PetscErrorCode TSSundialsSetMinTimeStep_Sundials(TS ts,PetscReal mindt)
526: {
527: TS_Sundials *cvode = (TS_Sundials*)ts->data;
530: cvode->mindt = mindt;
531: return(0);
532: }
534: PetscErrorCode TSSundialsSetMaxTimeStep_Sundials(TS ts,PetscReal maxdt)
535: {
536: TS_Sundials *cvode = (TS_Sundials*)ts->data;
539: cvode->maxdt = maxdt;
540: return(0);
541: }
542: PetscErrorCode TSSundialsGetPC_Sundials(TS ts,PC *pc)
543: {
544: SNES snes;
545: KSP ksp;
549: TSGetSNES(ts,&snes);
550: SNESGetKSP(snes,&ksp);
551: KSPGetPC(ksp,pc);
552: return(0);
553: }
555: PetscErrorCode TSSundialsGetIterations_Sundials(TS ts,int *nonlin,int *lin)
556: {
558: if (nonlin) *nonlin = ts->snes_its;
559: if (lin) *lin = ts->ksp_its;
560: return(0);
561: }
563: PetscErrorCode TSSundialsMonitorInternalSteps_Sundials(TS ts,PetscBool s)
564: {
565: TS_Sundials *cvode = (TS_Sundials*)ts->data;
568: cvode->monitorstep = s;
569: return(0);
570: }
571: /* -------------------------------------------------------------------------------------------*/
573: /*@C
574: TSSundialsGetIterations - Gets the number of nonlinear and linear iterations used so far by Sundials.
576: Not Collective
578: Input parameters:
579: . ts - the time-step context
581: Output Parameters:
582: + nonlin - number of nonlinear iterations
583: - lin - number of linear iterations
585: Level: advanced
587: Notes:
588: These return the number since the creation of the TS object
590: .keywords: non-linear iterations, linear iterations
592: .seealso: TSSundialsSetType(), TSSundialsSetMaxl(),
593: TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
594: TSSundialsGetIterations(), TSSundialsSetType(),
595: TSSundialsSetLinearTolerance(), TSSundialsGetPC(), TSSetExactFinalTime()
597: @*/
598: PetscErrorCode TSSundialsGetIterations(TS ts,int *nonlin,int *lin)
599: {
603: PetscUseMethod(ts,"TSSundialsGetIterations_C",(TS,int*,int*),(ts,nonlin,lin));
604: return(0);
605: }
607: /*@
608: TSSundialsSetType - Sets the method that Sundials will use for integration.
610: Logically Collective on TS
612: Input parameters:
613: + ts - the time-step context
614: - type - one of SUNDIALS_ADAMS or SUNDIALS_BDF
616: Level: intermediate
618: .keywords: Adams, backward differentiation formula
620: .seealso: TSSundialsGetIterations(), TSSundialsSetMaxl(),
621: TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
622: TSSundialsGetIterations(), TSSundialsSetType(),
623: TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
624: TSSetExactFinalTime()
625: @*/
626: PetscErrorCode TSSundialsSetType(TS ts,TSSundialsLmmType type)
627: {
631: PetscTryMethod(ts,"TSSundialsSetType_C",(TS,TSSundialsLmmType),(ts,type));
632: return(0);
633: }
635: /*@
636: TSSundialsSetMaxl - Sets the dimension of the Krylov space used by
637: GMRES in the linear solver in SUNDIALS. SUNDIALS DOES NOT use restarted GMRES so
638: this is the maximum number of GMRES steps that will be used.
640: Logically Collective on TS
642: Input parameters:
643: + ts - the time-step context
644: - maxl - number of direction vectors (the dimension of Krylov subspace).
646: Level: advanced
648: .keywords: GMRES
650: .seealso: TSSundialsGetIterations(), TSSundialsSetType(),
651: TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
652: TSSundialsGetIterations(), TSSundialsSetType(),
653: TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
654: TSSetExactFinalTime()
656: @*/
657: PetscErrorCode TSSundialsSetMaxl(TS ts,PetscInt maxl)
658: {
663: PetscTryMethod(ts,"TSSundialsSetMaxl_C",(TS,PetscInt),(ts,maxl));
664: return(0);
665: }
667: /*@
668: TSSundialsSetLinearTolerance - Sets the tolerance used to solve the linear
669: system by SUNDIALS.
671: Logically Collective on TS
673: Input parameters:
674: + ts - the time-step context
675: - tol - the factor by which the tolerance on the nonlinear solver is
676: multiplied to get the tolerance on the linear solver, .05 by default.
678: Level: advanced
680: .keywords: GMRES, linear convergence tolerance, SUNDIALS
682: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
683: TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
684: TSSundialsGetIterations(), TSSundialsSetType(),
685: TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
686: TSSetExactFinalTime()
688: @*/
689: PetscErrorCode TSSundialsSetLinearTolerance(TS ts,double tol)
690: {
695: PetscTryMethod(ts,"TSSundialsSetLinearTolerance_C",(TS,double),(ts,tol));
696: return(0);
697: }
699: /*@
700: TSSundialsSetGramSchmidtType - Sets type of orthogonalization used
701: in GMRES method by SUNDIALS linear solver.
703: Logically Collective on TS
705: Input parameters:
706: + ts - the time-step context
707: - type - either SUNDIALS_MODIFIED_GS or SUNDIALS_CLASSICAL_GS
709: Level: advanced
711: .keywords: Sundials, orthogonalization
713: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
714: TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(),
715: TSSundialsGetIterations(), TSSundialsSetType(),
716: TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
717: TSSetExactFinalTime()
719: @*/
720: PetscErrorCode TSSundialsSetGramSchmidtType(TS ts,TSSundialsGramSchmidtType type)
721: {
725: PetscTryMethod(ts,"TSSundialsSetGramSchmidtType_C",(TS,TSSundialsGramSchmidtType),(ts,type));
726: return(0);
727: }
729: /*@
730: TSSundialsSetTolerance - Sets the absolute and relative tolerance used by
731: Sundials for error control.
733: Logically Collective on TS
735: Input parameters:
736: + ts - the time-step context
737: . aabs - the absolute tolerance
738: - rel - the relative tolerance
740: See the Cvode/Sundials users manual for exact details on these parameters. Essentially
741: these regulate the size of the error for a SINGLE timestep.
743: Level: intermediate
745: .keywords: Sundials, tolerance
747: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetGMRESMaxl(),
748: TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(),
749: TSSundialsGetIterations(), TSSundialsSetType(),
750: TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC(),
751: TSSetExactFinalTime()
753: @*/
754: PetscErrorCode TSSundialsSetTolerance(TS ts,double aabs,double rel)
755: {
759: PetscTryMethod(ts,"TSSundialsSetTolerance_C",(TS,double,double),(ts,aabs,rel));
760: return(0);
761: }
763: /*@
764: TSSundialsGetPC - Extract the PC context from a time-step context for Sundials.
766: Input Parameter:
767: . ts - the time-step context
769: Output Parameter:
770: . pc - the preconditioner context
772: Level: advanced
774: .seealso: TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
775: TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
776: TSSundialsGetIterations(), TSSundialsSetType(),
777: TSSundialsSetLinearTolerance(), TSSundialsSetTolerance()
778: @*/
779: PetscErrorCode TSSundialsGetPC(TS ts,PC *pc)
780: {
784: PetscUseMethod(ts,"TSSundialsGetPC_C",(TS,PC*),(ts,pc));
785: return(0);
786: }
788: /*@
789: TSSundialsSetMinTimeStep - Smallest time step to be chosen by the adaptive controller.
791: Input Parameter:
792: + ts - the time-step context
793: - mindt - lowest time step if positive, negative to deactivate
795: Note:
796: Sundials will error if it is not possible to keep the estimated truncation error below
797: the tolerance set with TSSundialsSetTolerance() without going below this step size.
799: Level: beginner
801: .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
802: @*/
803: PetscErrorCode TSSundialsSetMinTimeStep(TS ts,PetscReal mindt)
804: {
808: PetscTryMethod(ts,"TSSundialsSetMinTimeStep_C",(TS,PetscReal),(ts,mindt));
809: return(0);
810: }
812: /*@
813: TSSundialsSetMaxTimeStep - Largest time step to be chosen by the adaptive controller.
815: Input Parameter:
816: + ts - the time-step context
817: - maxdt - lowest time step if positive, negative to deactivate
819: Level: beginner
821: .seealso: TSSundialsSetType(), TSSundialsSetTolerance(),
822: @*/
823: PetscErrorCode TSSundialsSetMaxTimeStep(TS ts,PetscReal maxdt)
824: {
828: PetscTryMethod(ts,"TSSundialsSetMaxTimeStep_C",(TS,PetscReal),(ts,maxdt));
829: return(0);
830: }
832: /*@
833: TSSundialsMonitorInternalSteps - Monitor Sundials internal steps (Defaults to false).
835: Input Parameter:
836: + ts - the time-step context
837: - ft - PETSC_TRUE if monitor, else PETSC_FALSE
839: Level: beginner
841: .seealso:TSSundialsGetIterations(), TSSundialsSetType(), TSSundialsSetMaxl(),
842: TSSundialsSetLinearTolerance(), TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(),
843: TSSundialsGetIterations(), TSSundialsSetType(),
844: TSSundialsSetLinearTolerance(), TSSundialsSetTolerance(), TSSundialsGetPC()
845: @*/
846: PetscErrorCode TSSundialsMonitorInternalSteps(TS ts,PetscBool ft)
847: {
851: PetscTryMethod(ts,"TSSundialsMonitorInternalSteps_C",(TS,PetscBool),(ts,ft));
852: return(0);
853: }
854: /* -------------------------------------------------------------------------------------------*/
855: /*MC
856: TSSUNDIALS - ODE solver using the LLNL CVODE/SUNDIALS package (now called SUNDIALS)
858: Options Database:
859: + -ts_sundials_type <bdf,adams> -
860: . -ts_sundials_gramschmidt_type <modified, classical> - type of orthogonalization inside GMRES
861: . -ts_sundials_atol <tol> - Absolute tolerance for convergence
862: . -ts_sundials_rtol <tol> - Relative tolerance for convergence
863: . -ts_sundials_linear_tolerance <tol> -
864: . -ts_sundials_maxl <maxl> - Max dimension of the Krylov subspace
865: - -ts_sundials_monitor_steps - Monitor SUNDIALS internal steps
868: Notes:
869: This uses its own nonlinear solver and Krylov method so PETSc SNES and KSP options do not apply,
870: only PETSc PC options.
872: Level: beginner
874: .seealso: TSCreate(), TS, TSSetType(), TSSundialsSetType(), TSSundialsSetMaxl(), TSSundialsSetLinearTolerance(),
875: TSSundialsSetGramSchmidtType(), TSSundialsSetTolerance(), TSSundialsGetPC(), TSSundialsGetIterations(), TSSetExactFinalTime()
877: M*/
878: PETSC_EXTERN PetscErrorCode TSCreate_Sundials(TS ts)
879: {
880: TS_Sundials *cvode;
882: PC pc;
885: ts->ops->reset = TSReset_Sundials;
886: ts->ops->destroy = TSDestroy_Sundials;
887: ts->ops->view = TSView_Sundials;
888: ts->ops->setup = TSSetUp_Sundials;
889: ts->ops->step = TSStep_Sundials;
890: ts->ops->interpolate = TSInterpolate_Sundials;
891: ts->ops->setfromoptions = TSSetFromOptions_Sundials;
892: ts->default_adapt_type = TSADAPTNONE;
894: PetscNewLog(ts,&cvode);
896: ts->usessnes = PETSC_TRUE;
898: ts->data = (void*)cvode;
899: cvode->cvode_type = SUNDIALS_BDF;
900: cvode->gtype = SUNDIALS_CLASSICAL_GS;
901: cvode->maxl = 5;
902: cvode->linear_tol = .05;
903: cvode->monitorstep = PETSC_TRUE;
905: MPI_Comm_dup(PetscObjectComm((PetscObject)ts),&(cvode->comm_sundials));
907: cvode->mindt = -1.;
908: cvode->maxdt = -1.;
910: /* set tolerance for Sundials */
911: cvode->reltol = 1e-6;
912: cvode->abstol = 1e-6;
914: /* set PCNONE as default pctype */
915: TSSundialsGetPC_Sundials(ts,&pc);
916: PCSetType(pc,PCNONE);
918: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetType_C",TSSundialsSetType_Sundials);
919: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxl_C",TSSundialsSetMaxl_Sundials);
920: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetLinearTolerance_C",TSSundialsSetLinearTolerance_Sundials);
921: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetGramSchmidtType_C",TSSundialsSetGramSchmidtType_Sundials);
922: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetTolerance_C",TSSundialsSetTolerance_Sundials);
923: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMinTimeStep_C",TSSundialsSetMinTimeStep_Sundials);
924: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsSetMaxTimeStep_C",TSSundialsSetMaxTimeStep_Sundials);
925: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetPC_C",TSSundialsGetPC_Sundials);
926: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsGetIterations_C",TSSundialsGetIterations_Sundials);
927: PetscObjectComposeFunction((PetscObject)ts,"TSSundialsMonitorInternalSteps_C",TSSundialsMonitorInternalSteps_Sundials);
928: return(0);
929: }