Actual source code: gr2.c
petsc-3.3-p0 2012-06-05
2: /*
3: Plots vectors obtained with DMDACreate2d()
4: */
6: #include <petsc-private/daimpl.h> /*I "petscdmda.h" I*/
7: #include <petsc-private/vecimpl.h>
9: /*
10: The data that is passed into the graphics callback
11: */
12: typedef struct {
13: PetscInt m,n,step,k;
14: PetscReal min,max,scale;
15: PetscScalar *xy,*v;
16: PetscBool showgrid;
17: } ZoomCtx;
19: /*
20: This does the drawing for one particular field
21: in one particular set of coordinates. It is a callback
22: called from PetscDrawZoom()
23: */
26: PetscErrorCode VecView_MPI_Draw_DA2d_Zoom(PetscDraw draw,void *ctx)
27: {
28: ZoomCtx *zctx = (ZoomCtx*)ctx;
30: PetscInt m,n,i,j,k,step,id,c1,c2,c3,c4;
31: PetscReal s,min,x1,x2,x3,x4,y_1,y2,y3,y4;
32: PetscScalar *v,*xy;
35: m = zctx->m;
36: n = zctx->n;
37: step = zctx->step;
38: k = zctx->k;
39: v = zctx->v;
40: xy = zctx->xy;
41: s = zctx->scale;
42: min = zctx->min;
43:
44: /* PetscDraw the contour plot patch */
45: for (j=0; j<n-1; j++) {
46: for (i=0; i<m-1; i++) {
47: #if !defined(PETSC_USE_COMPLEX)
48: id = i+j*m; x1 = xy[2*id];y_1 = xy[2*id+1];c1 = (int)(PETSC_DRAW_BASIC_COLORS+s*(v[k+step*id]-min));
49: id = i+j*m+1; x2 = xy[2*id];y2 = y_1; c2 = (int)(PETSC_DRAW_BASIC_COLORS+s*(v[k+step*id]-min));
50: id = i+j*m+1+m;x3 = x2; y3 = xy[2*id+1];c3 = (int)(PETSC_DRAW_BASIC_COLORS+s*(v[k+step*id]-min));
51: id = i+j*m+m; x4 = x1; y4 = y3; c4 = (int)(PETSC_DRAW_BASIC_COLORS+s*(v[k+step*id]-min));
52: #else
53: id = i+j*m; x1 = PetscRealPart(xy[2*id]);y_1 = PetscRealPart(xy[2*id+1]);c1 = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscRealPart(v[k+step*id])-min));
54: id = i+j*m+1; x2 = PetscRealPart(xy[2*id]);y2 = y_1; c2 = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscRealPart(v[k+step*id])-min));
55: id = i+j*m+1+m;x3 = x2; y3 = PetscRealPart(xy[2*id+1]);c3 = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscRealPart(v[k+step*id])-min));
56: id = i+j*m+m; x4 = x1; y4 = y3; c4 = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscRealPart(v[k+step*id])-min));
57: #endif
58: PetscDrawTriangle(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3);
59: PetscDrawTriangle(draw,x1,y_1,x3,y3,x4,y4,c1,c3,c4);
60: if (zctx->showgrid) {
61: PetscDrawLine(draw,x1,y_1,x2,y2,PETSC_DRAW_BLACK);
62: PetscDrawLine(draw,x2,y2,x3,y3,PETSC_DRAW_BLACK);
63: PetscDrawLine(draw,x3,y3,x4,y4,PETSC_DRAW_BLACK);
64: PetscDrawLine(draw,x4,y4,x1,y_1,PETSC_DRAW_BLACK);
65: }
66: }
67: }
68: return(0);
69: }
73: PetscErrorCode VecView_MPI_Draw_DA2d(Vec xin,PetscViewer viewer)
74: {
75: DM da,dac,dag;
76: PetscErrorCode ierr;
77: PetscMPIInt rank;
78: PetscInt N,s,M,w;
79: const PetscInt *lx,*ly;
80: PetscReal coors[4],ymin,ymax,xmin,xmax;
81: PetscDraw draw,popup;
82: PetscBool isnull,useports = PETSC_FALSE;
83: MPI_Comm comm;
84: Vec xlocal,xcoor,xcoorl;
85: DMDABoundaryType bx,by;
86: DMDAStencilType st;
87: ZoomCtx zctx;
88: PetscDrawViewPorts *ports = PETSC_NULL;
89: PetscViewerFormat format;
90: PetscInt *displayfields;
91: PetscInt ndisplayfields,i,nbounds;
92: PetscBool flg;
93: const PetscReal *bounds;
96: zctx.showgrid = PETSC_FALSE;
97: PetscViewerDrawGetDraw(viewer,0,&draw);
98: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
99: PetscViewerDrawGetBounds(viewer,&nbounds,&bounds);
101: PetscObjectQuery((PetscObject)xin,"DM",(PetscObject*)&da);
102: if (!da) SETERRQ(((PetscObject)xin)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
104: PetscObjectGetComm((PetscObject)xin,&comm);
105: MPI_Comm_rank(comm,&rank);
107: DMDAGetInfo(da,0,&M,&N,0,&zctx.m,&zctx.n,0,&w,&s,&bx,&by,0,&st);
108: DMDAGetOwnershipRanges(da,&lx,&ly,PETSC_NULL);
110: /*
111: Obtain a sequential vector that is going to contain the local values plus ONE layer of
112: ghosted values to draw the graphics from. We also need its corresponding DMDA (dac) that will
113: update the local values pluse ONE layer of ghost values.
114: */
115: PetscObjectQuery((PetscObject)da,"GraphicsGhosted",(PetscObject*)&xlocal);
116: if (!xlocal) {
117: if (bx != DMDA_BOUNDARY_NONE || by != DMDA_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) {
118: /*
119: if original da is not of stencil width one, or periodic or not a box stencil then
120: create a special DMDA to handle one level of ghost points for graphics
121: */
122: DMDACreate2d(comm,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,zctx.m,zctx.n,w,1,lx,ly,&dac);
123: PetscInfo(da,"Creating auxilary DMDA for managing graphics ghost points\n");
124: } else {
125: /* otherwise we can use the da we already have */
126: dac = da;
127: }
128: /* create local vector for holding ghosted values used in graphics */
129: DMCreateLocalVector(dac,&xlocal);
130: if (dac != da) {
131: /* don't keep any public reference of this DMDA, is is only available through xlocal */
132: PetscObjectDereference((PetscObject)dac);
133: } else {
134: /* remove association between xlocal and da, because below we compose in the opposite
135: direction and if we left this connect we'd get a loop, so the objects could
136: never be destroyed */
137: PetscObjectRemoveReference((PetscObject)xlocal,"DM");
138: }
139: PetscObjectCompose((PetscObject)da,"GraphicsGhosted",(PetscObject)xlocal);
140: PetscObjectDereference((PetscObject)xlocal);
141: } else {
142: if (bx != DMDA_BOUNDARY_NONE || by != DMDA_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) {
143: PetscObjectQuery((PetscObject)xlocal,"DM",(PetscObject*)&dac);
144: } else {
145: dac = da;
146: }
147: }
149: /*
150: Get local (ghosted) values of vector
151: */
152: DMGlobalToLocalBegin(dac,xin,INSERT_VALUES,xlocal);
153: DMGlobalToLocalEnd(dac,xin,INSERT_VALUES,xlocal);
154: VecGetArray(xlocal,&zctx.v);
156: /* get coordinates of nodes */
157: DMDAGetCoordinates(da,&xcoor);
158: if (!xcoor) {
159: DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,0.0);
160: DMDAGetCoordinates(da,&xcoor);
161: }
163: /*
164: Determine the min and max coordinates in plot
165: */
166: VecStrideMin(xcoor,0,PETSC_NULL,&xmin);
167: VecStrideMax(xcoor,0,PETSC_NULL,&xmax);
168: VecStrideMin(xcoor,1,PETSC_NULL,&ymin);
169: VecStrideMax(xcoor,1,PETSC_NULL,&ymax);
170: coors[0] = xmin - .05*(xmax- xmin); coors[2] = xmax + .05*(xmax - xmin);
171: coors[1] = ymin - .05*(ymax- ymin); coors[3] = ymax + .05*(ymax - ymin);
172: PetscInfo4(da,"Preparing DMDA 2d contour plot coordinates %G %G %G %G\n",coors[0],coors[1],coors[2],coors[3]);
174: /*
175: get local ghosted version of coordinates
176: */
177: PetscObjectQuery((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject*)&xcoorl);
178: if (!xcoorl) {
179: /* create DMDA to get local version of graphics */
180: DMDACreate2d(comm,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,zctx.m,zctx.n,2,1,lx,ly,&dag);
181: PetscInfo(dag,"Creating auxilary DMDA for managing graphics coordinates ghost points\n");
182: DMCreateLocalVector(dag,&xcoorl);
183: PetscObjectCompose((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject)xcoorl);
184: PetscObjectDereference((PetscObject)dag);
185: PetscObjectDereference((PetscObject)xcoorl);
186: } else {
187: PetscObjectQuery((PetscObject)xcoorl,"DM",(PetscObject*)&dag);
188: }
189: DMGlobalToLocalBegin(dag,xcoor,INSERT_VALUES,xcoorl);
190: DMGlobalToLocalEnd(dag,xcoor,INSERT_VALUES,xcoorl);
191: VecGetArray(xcoorl,&zctx.xy);
192:
193: /*
194: Get information about size of area each processor must do graphics for
195: */
196: DMDAGetInfo(dac,0,&M,&N,0,0,0,0,&zctx.step,0,&bx,&by,0,0);
197: DMDAGetGhostCorners(dac,0,0,0,&zctx.m,&zctx.n,0);
199: PetscOptionsGetBool(PETSC_NULL,"-draw_contour_grid",&zctx.showgrid,PETSC_NULL);
200: PetscMalloc(zctx.step*sizeof(PetscInt),&displayfields);
201: for (i=0; i<zctx.step; i++) displayfields[i] = i;
202: ndisplayfields = zctx.step;
203: PetscOptionsGetIntArray(PETSC_NULL,"-draw_fields",displayfields,&ndisplayfields,&flg);
204: if (!flg) ndisplayfields = zctx.step;
206: PetscViewerGetFormat(viewer,&format);
207: PetscOptionsGetBool(PETSC_NULL,"-draw_ports",&useports,PETSC_NULL);
208: if (useports || format == PETSC_VIEWER_DRAW_PORTS){
209: PetscDrawSynchronizedClear(draw);
210: PetscDrawViewPortsCreate(draw,ndisplayfields,&ports);
211: }
213: /*
214: Loop over each field; drawing each in a different window
215: */
216: for (i=0; i<ndisplayfields; i++) {
217: zctx.k = displayfields[i];
218: if (useports) {
219: PetscDrawViewPortsSet(ports,i);
220: } else {
221: PetscViewerDrawGetDraw(viewer,i,&draw);
222: PetscDrawSynchronizedClear(draw);
223: }
225: /*
226: Determine the min and max color in plot
227: */
228: VecStrideMin(xin,zctx.k,PETSC_NULL,&zctx.min);
229: VecStrideMax(xin,zctx.k,PETSC_NULL,&zctx.max);
230: if (zctx.k < nbounds) {
231: zctx.min = PetscMin(zctx.min,bounds[2*zctx.k]);
232: zctx.max = PetscMax(zctx.max,bounds[2*zctx.k+1]);
233: }
234: if (zctx.min == zctx.max) {
235: zctx.min -= 1.e-12;
236: zctx.max += 1.e-12;
237: }
239: if (!rank) {
240: const char *title;
242: DMDAGetFieldName(da,zctx.k,&title);
243: if (title) {
244: PetscDrawSetTitle(draw,title);
245: }
246: }
247: PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);
248: PetscInfo2(da,"DMDA 2d contour plot min %G max %G\n",zctx.min,zctx.max);
250: PetscDrawGetPopup(draw,&popup);
251: if (popup) {PetscDrawScalePopup(popup,zctx.min,zctx.max);}
253: zctx.scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/(zctx.max - zctx.min);
255: PetscDrawZoom(draw,VecView_MPI_Draw_DA2d_Zoom,&zctx);
256: }
257: PetscFree(displayfields);
258: PetscDrawViewPortsDestroy(ports);
260: VecRestoreArray(xcoorl,&zctx.xy);
261: VecRestoreArray(xlocal,&zctx.v);
262: return(0);
263: }
266: #if defined(PETSC_HAVE_HDF5)
269: PetscErrorCode VecView_MPI_HDF5_DA(Vec xin,PetscViewer viewer)
270: {
271: DM dm;
272: DM_DA *da;
273: hsize_t dim,dims[5];
274: hsize_t count[5];
275: hsize_t offset[5];
276: PetscInt cnt = 0;
277: PetscScalar *x;
278: const char *vecname;
279: hid_t filespace; /* file dataspace identifier */
280: hid_t plist_id; /* property list identifier */
281: hid_t dset_id; /* dataset identifier */
282: hid_t memspace; /* memory dataspace identifier */
283: hid_t scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
284: hid_t file_id;
285: hid_t group;
286: herr_t status;
287: PetscInt timestep;
291: PetscViewerHDF5OpenGroup(viewer, &file_id, &group);
292: //PetscViewerHDF5GetTimestep(viewer, ×tep);
294: PetscObjectQuery((PetscObject)xin,"DM",(PetscObject*)&dm);
295: if (!dm) SETERRQ(((PetscObject)xin)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
296: da = (DM_DA*)dm->data;
298: /* Create the dataspace for the dataset */
299: dim = PetscHDF5IntCast(da->dim + ((da->w == 1) ? 0 : 1));
300: if (da->dim == 3) dims[cnt++] = PetscHDF5IntCast(da->P);
301: if (da->dim > 1) dims[cnt++] = PetscHDF5IntCast(da->N);
302: dims[cnt++] = PetscHDF5IntCast(da->M);
303: if (da->w > 1) dims[cnt++] = PetscHDF5IntCast(da->w);
304: #if defined(PETSC_USE_COMPLEX)
305: dim++;
306: dims[cnt++] = 2;
307: #endif
308: filespace = H5Screate_simple(dim, dims, NULL);
309: if (filespace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Screate_simple()");
311: #if defined(PETSC_USE_REAL_SINGLE)
312: scalartype = H5T_NATIVE_FLOAT;
313: #elif defined(PETSC_USE_REAL___FLOAT128)
314: #error "HDF5 output with 128 bit floats not supported."
315: #else
316: scalartype = H5T_NATIVE_DOUBLE;
317: #endif
319: /* Create the dataset with default properties and close filespace */
320: PetscObjectGetName((PetscObject)xin,&vecname);
321: if (!H5Lexists(group, vecname, H5P_DEFAULT)) {
322: #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
323: dset_id = H5Dcreate2(group, vecname, scalartype, filespace, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
324: #else
325: dset_id = H5Dcreate(group, vecname, scalartype, filespace, H5P_DEFAULT);
326: #endif
327: if (dset_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Dcreate2()");
328: } else {
329: dset_id = H5Dopen2(group, vecname, H5P_DEFAULT);
330: status = H5Dset_extent(dset_id, dims);CHKERRQ(status);
331: }
332: status = H5Sclose(filespace);CHKERRQ(status);
334: /* Each process defines a dataset and writes it to the hyperslab in the file */
335: cnt = 0;
336: if (da->dim == 3) offset[cnt++] = PetscHDF5IntCast(da->zs);
337: if (da->dim > 1) offset[cnt++] = PetscHDF5IntCast(da->ys);
338: offset[cnt++] = PetscHDF5IntCast(da->xs/da->w);
339: if (da->w > 1) offset[cnt++] = 0;
340: #if defined(PETSC_USE_COMPLEX)
341: offset[cnt++] = 0;
342: #endif
343: cnt = 0;
344: if (da->dim == 3) count[cnt++] = PetscHDF5IntCast(da->ze - da->zs);
345: if (da->dim > 1) count[cnt++] = PetscHDF5IntCast(da->ye - da->ys);
346: count[cnt++] = PetscHDF5IntCast((da->xe - da->xs)/da->w);
347: if (da->w > 1) count[cnt++] = PetscHDF5IntCast(da->w);
348: #if defined(PETSC_USE_COMPLEX)
349: count[cnt++] = 2;
350: #endif
351: memspace = H5Screate_simple(dim, count, NULL);
352: if (memspace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Screate_simple()");
355: filespace = H5Dget_space(dset_id);
356: if (filespace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Dget_space()");
357: status = H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);CHKERRQ(status);
359: /* Create property list for collective dataset write */
360: plist_id = H5Pcreate(H5P_DATASET_XFER);
361: if (plist_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Pcreate()");
362: #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
363: status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);CHKERRQ(status);
364: #endif
365: /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */
367: VecGetArray(xin, &x);
368: status = H5Dwrite(dset_id, scalartype, memspace, filespace, plist_id, x);CHKERRQ(status);
369: status = H5Fflush(file_id, H5F_SCOPE_GLOBAL);CHKERRQ(status);
370: VecRestoreArray(xin, &x);
372: /* Close/release resources */
373: if (group != file_id) {
374: status = H5Gclose(group);CHKERRQ(status);
375: }
376: status = H5Pclose(plist_id);CHKERRQ(status);
377: status = H5Sclose(filespace);CHKERRQ(status);
378: status = H5Sclose(memspace);CHKERRQ(status);
379: status = H5Dclose(dset_id);CHKERRQ(status);
380: PetscInfo1(xin,"Wrote Vec object with name %s\n",vecname);
381: return(0);
382: }
383: #endif
385: extern PetscErrorCode VecView_MPI_Draw_DA1d(Vec,PetscViewer);
387: #if defined(PETSC_HAVE_MPIIO)
390: static PetscErrorCode DMDAArrayMPIIO(DM da,PetscViewer viewer,Vec xin,PetscBool write)
391: {
392: PetscErrorCode ierr;
393: MPI_File mfdes;
394: PetscMPIInt gsizes[4],lsizes[4],lstarts[4],asiz,dof;
395: MPI_Datatype view;
396: const PetscScalar *array;
397: MPI_Offset off;
398: MPI_Aint ub,ul;
399: PetscInt type,rows,vecrows,tr[2];
400: DM_DA *dd = (DM_DA*)da->data;
403: VecGetSize(xin,&vecrows);
404: if (!write) {
405: /* Read vector header. */
406: PetscViewerBinaryRead(viewer,tr,2,PETSC_INT);
407: type = tr[0];
408: rows = tr[1];
409: if (type != VEC_FILE_CLASSID) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Not vector next in file");
410: if (rows != vecrows) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Vector in file not same size as DMDA vector");
411: } else {
412: tr[0] = VEC_FILE_CLASSID;
413: tr[1] = vecrows;
414: PetscViewerBinaryWrite(viewer,tr,2,PETSC_INT,PETSC_TRUE);
415: }
417: dof = PetscMPIIntCast(dd->w);
418: gsizes[0] = dof; gsizes[1] = PetscMPIIntCast(dd->M); gsizes[2] = PetscMPIIntCast(dd->N); gsizes[3] = PetscMPIIntCast(dd->P);
419: lsizes[0] = dof;lsizes[1] = PetscMPIIntCast((dd->xe-dd->xs)/dof); lsizes[2] = PetscMPIIntCast(dd->ye-dd->ys); lsizes[3] = PetscMPIIntCast(dd->ze-dd->zs);
420: lstarts[0] = 0; lstarts[1] = PetscMPIIntCast(dd->xs/dof); lstarts[2] = PetscMPIIntCast(dd->ys); lstarts[3] = PetscMPIIntCast(dd->zs);
421: MPI_Type_create_subarray(dd->dim+1,gsizes,lsizes,lstarts,MPI_ORDER_FORTRAN,MPIU_SCALAR,&view);
422: MPI_Type_commit(&view);
423:
424: PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);
425: PetscViewerBinaryGetMPIIOOffset(viewer,&off);
426: MPI_File_set_view(mfdes,off,MPIU_SCALAR,view,(char *)"native",MPI_INFO_NULL);
427: VecGetArrayRead(xin,&array);
428: asiz = lsizes[1]*(lsizes[2] > 0 ? lsizes[2] : 1)*(lsizes[3] > 0 ? lsizes[3] : 1)*dof;
429: if (write) {
430: MPIU_File_write_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);
431: } else {
432: MPIU_File_read_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);
433: }
434: MPI_Type_get_extent(view,&ul,&ub);
435: PetscViewerBinaryAddMPIIOOffset(viewer,ub);
436: VecRestoreArrayRead(xin,&array);
437: MPI_Type_free(&view);
438: return(0);
439: }
440: #endif
442: EXTERN_C_BEGIN
445: PetscErrorCode VecView_MPI_DA(Vec xin,PetscViewer viewer)
446: {
447: DM da;
449: PetscInt dim;
450: Vec natural;
451: PetscBool isdraw,isvtk;
452: #if defined(PETSC_HAVE_HDF5)
453: PetscBool ishdf5;
454: #endif
455: const char *prefix,*name;
458: PetscObjectQuery((PetscObject)xin,"DM",(PetscObject*)&da);
459: if (!da) SETERRQ(((PetscObject)xin)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
460: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
461: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk);
462: #if defined(PETSC_HAVE_HDF5)
463: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);
464: #endif
465: if (isdraw) {
466: DMDAGetInfo(da,&dim,0,0,0,0,0,0,0,0,0,0,0,0);
467: if (dim == 1) {
468: VecView_MPI_Draw_DA1d(xin,viewer);
469: } else if (dim == 2) {
470: VecView_MPI_Draw_DA2d(xin,viewer);
471: } else {
472: SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Cannot graphically view vector associated with this dimensional DMDA %D",dim);
473: }
474: } else if (isvtk) { /* Duplicate the Vec and hold a reference to the DM */
475: Vec Y;
476: PetscObjectReference((PetscObject)da);
477: VecDuplicate(xin,&Y);
478: if (((PetscObject)xin)->name) {
479: /* If xin was named, copy the name over to Y. The duplicate names are safe because nobody else will ever see Y. */
480: PetscObjectSetName((PetscObject)Y,((PetscObject)xin)->name);
481: }
482: VecCopy(xin,Y);
483: PetscViewerVTKAddField(viewer,(PetscObject)da,DMDAVTKWriteAll,PETSC_VTK_POINT_FIELD,(PetscObject)Y);
484: #if defined(PETSC_HAVE_HDF5)
485: } else if (ishdf5) {
486: VecView_MPI_HDF5_DA(xin,viewer);
487: #endif
488: } else {
489: #if defined(PETSC_HAVE_MPIIO)
490: PetscBool isbinary,isMPIIO;
492: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
493: if (isbinary) {
494: PetscViewerBinaryGetMPIIO(viewer,&isMPIIO);
495: if (isMPIIO) {
496: DMDAArrayMPIIO(da,viewer,xin,PETSC_TRUE);
497: return(0);
498: }
499: }
500: #endif
501:
502: /* call viewer on natural ordering */
503: PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);
504: DMDACreateNaturalVector(da,&natural);
505: PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);
506: DMDAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural);
507: DMDAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural);
508: PetscObjectGetName((PetscObject)xin,&name);
509: PetscObjectSetName((PetscObject)natural,name);
510: VecView(natural,viewer);
511: VecDestroy(&natural);
512: }
513: return(0);
514: }
515: EXTERN_C_END
517: #if defined(PETSC_HAVE_HDF5)
520: PetscErrorCode VecLoad_HDF5_DA(Vec xin, PetscViewer viewer)
521: {
522: DM da;
524: hsize_t dim;
525: hsize_t count[5];
526: hsize_t offset[5];
527: PetscInt cnt = 0;
528: PetscScalar *x;
529: const char *vecname;
530: hid_t filespace; /* file dataspace identifier */
531: hid_t plist_id; /* property list identifier */
532: hid_t dset_id; /* dataset identifier */
533: hid_t memspace; /* memory dataspace identifier */
534: hid_t file_id;
535: herr_t status;
536: DM_DA *dd;
539: PetscViewerHDF5GetFileId(viewer, &file_id);
540: PetscObjectQuery((PetscObject)xin,"DM",(PetscObject*)&da);
541: dd = (DM_DA*)da->data;
543: /* Create the dataspace for the dataset */
544: dim = PetscHDF5IntCast(dd->dim + ((dd->w == 1) ? 0 : 1));
545: #if defined(PETSC_USE_COMPLEX)
546: dim++;
547: #endif
549: /* Create the dataset with default properties and close filespace */
550: PetscObjectGetName((PetscObject)xin,&vecname);
551: #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
552: dset_id = H5Dopen2(file_id, vecname, H5P_DEFAULT);
553: #else
554: dset_id = H5Dopen(file_id, vecname);
555: #endif
556: if (dset_id == -1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Dopen2() with Vec named %s",vecname);
557: filespace = H5Dget_space(dset_id);
559: /* Each process defines a dataset and reads it from the hyperslab in the file */
560: cnt = 0;
561: if (dd->dim == 3) offset[cnt++] = PetscHDF5IntCast(dd->zs);
562: if (dd->dim > 1) offset[cnt++] = PetscHDF5IntCast(dd->ys);
563: offset[cnt++] = PetscHDF5IntCast(dd->xs/dd->w);
564: if (dd->w > 1) offset[cnt++] = 0;
565: #if defined(PETSC_USE_COMPLEX)
566: offset[cnt++] = 0;
567: #endif
568: cnt = 0;
569: if (dd->dim == 3) count[cnt++] = PetscHDF5IntCast(dd->ze - dd->zs);
570: if (dd->dim > 1) count[cnt++] = PetscHDF5IntCast(dd->ye - dd->ys);
571: count[cnt++] = PetscHDF5IntCast((dd->xe - dd->xs)/dd->w);
572: if (dd->w > 1) count[cnt++] = PetscHDF5IntCast(dd->w);
573: #if defined(PETSC_USE_COMPLEX)
574: count[cnt++] = 2;
575: #endif
576: memspace = H5Screate_simple(dim, count, NULL);
577: if (memspace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Screate_simple()");
579: status = H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);CHKERRQ(status);
581: /* Create property list for collective dataset write */
582: plist_id = H5Pcreate(H5P_DATASET_XFER);
583: if (plist_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Pcreate()");
584: #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
585: status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);CHKERRQ(status);
586: #endif
587: /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */
589: VecGetArray(xin, &x);
590: status = H5Dread(dset_id, H5T_NATIVE_DOUBLE, memspace, filespace, plist_id, x);CHKERRQ(status);
591: VecRestoreArray(xin, &x);
593: /* Close/release resources */
594: status = H5Pclose(plist_id);CHKERRQ(status);
595: status = H5Sclose(filespace);CHKERRQ(status);
596: status = H5Sclose(memspace);CHKERRQ(status);
597: status = H5Dclose(dset_id);CHKERRQ(status);
598: return(0);
599: }
600: #endif
604: PetscErrorCode VecLoad_Binary_DA(Vec xin, PetscViewer viewer)
605: {
606: DM da;
608: Vec natural;
609: const char *prefix;
610: PetscInt bs;
611: PetscBool flag;
612: DM_DA *dd;
613: #if defined(PETSC_HAVE_MPIIO)
614: PetscBool isMPIIO;
615: #endif
618: PetscObjectQuery((PetscObject)xin,"DM",(PetscObject*)&da);
619: dd = (DM_DA*)da->data;
620: #if defined(PETSC_HAVE_MPIIO)
621: PetscViewerBinaryGetMPIIO(viewer,&isMPIIO);
622: if (isMPIIO) {
623: DMDAArrayMPIIO(da,viewer,xin,PETSC_FALSE);
624: return(0);
625: }
626: #endif
628: PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);
629: DMDACreateNaturalVector(da,&natural);
630: PetscObjectSetName((PetscObject)natural,((PetscObject)xin)->name);
631: PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);
632: VecLoad_Binary(natural,viewer);
633: DMDANaturalToGlobalBegin(da,natural,INSERT_VALUES,xin);
634: DMDANaturalToGlobalEnd(da,natural,INSERT_VALUES,xin);
635: VecDestroy(&natural);
636: PetscInfo(xin,"Loading vector from natural ordering into DMDA\n");
637: PetscOptionsGetInt(((PetscObject)xin)->prefix,"-vecload_block_size",&bs,&flag);
638: if (flag && bs != dd->w) {
639: PetscInfo2(xin,"Block size in file %D not equal to DMDA's dof %D\n",bs,dd->w);
640: }
641: return(0);
642: }
644: EXTERN_C_BEGIN
647: PetscErrorCode VecLoad_Default_DA(Vec xin, PetscViewer viewer)
648: {
650: DM da;
651: PetscBool isbinary;
652: #if defined(PETSC_HAVE_HDF5)
653: PetscBool ishdf5;
654: #endif
657: PetscObjectQuery((PetscObject)xin,"DM",(PetscObject*)&da);
658: if (!da) SETERRQ(((PetscObject)xin)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA");
660: #if defined(PETSC_HAVE_HDF5)
661: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);
662: #endif
663: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
665: if (isbinary) {
666: VecLoad_Binary_DA(xin,viewer);
667: #if defined(PETSC_HAVE_HDF5)
668: } else if (ishdf5) {
669: VecLoad_HDF5_DA(xin,viewer);
670: #endif
671: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for vector loading", ((PetscObject)viewer)->type_name);
672: return(0);
673: }
674: EXTERN_C_END