Actual source code: gr2.c
1: #define PETSCDM_DLL
3: /*
4: Plots vectors obtained with DACreate2d()
5: */
7: #include private/daimpl.h
8: #include private/vecimpl.h
10: #if defined(PETSC_HAVE_PNETCDF)
12: #include "pnetcdf.h"
14: #endif
17: /*
18: The data that is passed into the graphics callback
19: */
20: typedef struct {
21: PetscInt m,n,step,k;
22: PetscReal min,max,scale;
23: PetscScalar *xy,*v;
24: PetscTruth showgrid;
25: } ZoomCtx;
27: /*
28: This does the drawing for one particular field
29: in one particular set of coordinates. It is a callback
30: called from PetscDrawZoom()
31: */
34: PetscErrorCode VecView_MPI_Draw_DA2d_Zoom(PetscDraw draw,void *ctx)
35: {
36: ZoomCtx *zctx = (ZoomCtx*)ctx;
38: PetscInt m,n,i,j,k,step,id,c1,c2,c3,c4;
39: PetscReal s,min,x1,x2,x3,x4,y_1,y2,y3,y4;
40: PetscScalar *v,*xy;
43: m = zctx->m;
44: n = zctx->n;
45: step = zctx->step;
46: k = zctx->k;
47: v = zctx->v;
48: xy = zctx->xy;
49: s = zctx->scale;
50: min = zctx->min;
51:
52: /* PetscDraw the contour plot patch */
53: for (j=0; j<n-1; j++) {
54: for (i=0; i<m-1; i++) {
55: #if !defined(PETSC_USE_COMPLEX)
56: 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));
57: id = i+j*m+1; x2 = xy[2*id];y2 = y_1; c2 = (int)(PETSC_DRAW_BASIC_COLORS+s*(v[k+step*id]-min));
58: 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));
59: id = i+j*m+m; x4 = x1; y4 = y3; c4 = (int)(PETSC_DRAW_BASIC_COLORS+s*(v[k+step*id]-min));
60: #else
61: 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));
62: 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));
63: 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));
64: id = i+j*m+m; x4 = x1; y4 = y3; c4 = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscRealPart(v[k+step*id])-min));
65: #endif
66: PetscDrawTriangle(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3);
67: PetscDrawTriangle(draw,x1,y_1,x3,y3,x4,y4,c1,c3,c4);
68: if (zctx->showgrid) {
69: PetscDrawLine(draw,x1,y_1,x2,y2,PETSC_DRAW_BLACK);
70: PetscDrawLine(draw,x2,y2,x3,y3,PETSC_DRAW_BLACK);
71: PetscDrawLine(draw,x3,y3,x4,y4,PETSC_DRAW_BLACK);
72: PetscDrawLine(draw,x4,y4,x1,y_1,PETSC_DRAW_BLACK);
73: }
74: }
75: }
76: return(0);
77: }
81: PetscErrorCode VecView_MPI_Draw_DA2d(Vec xin,PetscViewer viewer)
82: {
83: DA da,dac,dag;
84: PetscErrorCode ierr;
85: PetscMPIInt rank;
86: PetscInt igstart,N,s,M,istart,isize,jgstart,w;
87: const PetscInt *lx,*ly;
88: PetscReal coors[4],ymin,ymax,xmin,xmax;
89: PetscDraw draw,popup;
90: PetscTruth isnull,useports = PETSC_FALSE;
91: MPI_Comm comm;
92: Vec xlocal,xcoor,xcoorl;
93: DAPeriodicType periodic;
94: DAStencilType st;
95: ZoomCtx zctx;
96: PetscDrawViewPorts *ports;
97: PetscViewerFormat format;
100: zctx.showgrid = PETSC_FALSE;
101: PetscViewerDrawGetDraw(viewer,0,&draw);
102: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
104: PetscObjectQuery((PetscObject)xin,"DA",(PetscObject*)&da);
105: if (!da) SETERRQ(PETSC_ERR_ARG_WRONG,"Vector not generated from a DA");
107: PetscObjectGetComm((PetscObject)xin,&comm);
108: MPI_Comm_rank(comm,&rank);
110: DAGetInfo(da,0,&M,&N,0,&zctx.m,&zctx.n,0,&w,&s,&periodic,&st);
111: DAGetOwnershipRanges(da,&lx,&ly,PETSC_NULL);
113: /*
114: Obtain a sequential vector that is going to contain the local values plus ONE layer of
115: ghosted values to draw the graphics from. We also need its corresponding DA (dac) that will
116: update the local values pluse ONE layer of ghost values.
117: */
118: PetscObjectQuery((PetscObject)da,"GraphicsGhosted",(PetscObject*)&xlocal);
119: if (!xlocal) {
120: if (periodic != DA_NONPERIODIC || s != 1 || st != DA_STENCIL_BOX) {
121: /*
122: if original da is not of stencil width one, or periodic or not a box stencil then
123: create a special DA to handle one level of ghost points for graphics
124: */
125: DACreate2d(comm,DA_NONPERIODIC,DA_STENCIL_BOX,M,N,zctx.m,zctx.n,w,1,lx,ly,&dac);
126: PetscInfo(da,"Creating auxilary DA for managing graphics ghost points\n");
127: } else {
128: /* otherwise we can use the da we already have */
129: dac = da;
130: }
131: /* create local vector for holding ghosted values used in graphics */
132: DACreateLocalVector(dac,&xlocal);
133: if (dac != da) {
134: /* don't keep any public reference of this DA, is is only available through xlocal */
135: DADestroy(dac);
136: } else {
137: /* remove association between xlocal and da, because below we compose in the opposite
138: direction and if we left this connect we'd get a loop, so the objects could
139: never be destroyed */
140: PetscObjectCompose((PetscObject)xlocal,"DA",0);
141: }
142: PetscObjectCompose((PetscObject)da,"GraphicsGhosted",(PetscObject)xlocal);
143: PetscObjectDereference((PetscObject)xlocal);
144: } else {
145: if (periodic == DA_NONPERIODIC && s == 1 && st == DA_STENCIL_BOX) {
146: dac = da;
147: } else {
148: PetscObjectQuery((PetscObject)xlocal,"DA",(PetscObject*)&dac);
149: }
150: }
152: /*
153: Get local (ghosted) values of vector
154: */
155: DAGlobalToLocalBegin(dac,xin,INSERT_VALUES,xlocal);
156: DAGlobalToLocalEnd(dac,xin,INSERT_VALUES,xlocal);
157: VecGetArray(xlocal,&zctx.v);
159: /* get coordinates of nodes */
160: DAGetCoordinates(da,&xcoor);
161: if (!xcoor) {
162: DASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,0.0);
163: DAGetCoordinates(da,&xcoor);
164: }
166: /*
167: Determine the min and max coordinates in plot
168: */
169: VecStrideMin(xcoor,0,PETSC_NULL,&xmin);
170: VecStrideMax(xcoor,0,PETSC_NULL,&xmax);
171: VecStrideMin(xcoor,1,PETSC_NULL,&ymin);
172: VecStrideMax(xcoor,1,PETSC_NULL,&ymax);
173: coors[0] = xmin - .05*(xmax- xmin); coors[2] = xmax + .05*(xmax - xmin);
174: coors[1] = ymin - .05*(ymax- ymin); coors[3] = ymax + .05*(ymax - ymin);
175: PetscInfo4(da,"Preparing DA 2d contour plot coordinates %G %G %G %G\n",coors[0],coors[1],coors[2],coors[3]);
177: /*
178: get local ghosted version of coordinates
179: */
180: PetscObjectQuery((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject*)&xcoorl);
181: if (!xcoorl) {
182: /* create DA to get local version of graphics */
183: DACreate2d(comm,DA_NONPERIODIC,DA_STENCIL_BOX,M,N,zctx.m,zctx.n,2,1,lx,ly,&dag);
184: PetscInfo(dag,"Creating auxilary DA for managing graphics coordinates ghost points\n");
185: DACreateLocalVector(dag,&xcoorl);
186: PetscObjectCompose((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject)xcoorl);
187: DADestroy(dag);/* dereference dag */
188: PetscObjectDereference((PetscObject)xcoorl);
189: } else {
190: PetscObjectQuery((PetscObject)xcoorl,"DA",(PetscObject*)&dag);
191: }
192: DAGlobalToLocalBegin(dag,xcoor,INSERT_VALUES,xcoorl);
193: DAGlobalToLocalEnd(dag,xcoor,INSERT_VALUES,xcoorl);
194: VecGetArray(xcoorl,&zctx.xy);
195:
196: /*
197: Get information about size of area each processor must do graphics for
198: */
199: DAGetInfo(dac,0,&M,&N,0,0,0,0,&zctx.step,0,&periodic,0);
200: DAGetGhostCorners(dac,&igstart,&jgstart,0,&zctx.m,&zctx.n,0);
201: DAGetCorners(dac,&istart,0,0,&isize,0,0);
203: PetscOptionsGetTruth(PETSC_NULL,"-draw_contour_grid",&zctx.showgrid,PETSC_NULL);
205: PetscViewerGetFormat(viewer,&format);
206: PetscOptionsGetTruth(PETSC_NULL,"-draw_ports",&useports,PETSC_NULL);
207: if (useports || format == PETSC_VIEWER_DRAW_PORTS){
208: PetscDrawSynchronizedClear(draw);
209: PetscDrawViewPortsCreate(draw,zctx.step,&ports);
210: }
211: /*
212: Loop over each field; drawing each in a different window
213: */
214: for (zctx.k=0; zctx.k<zctx.step; zctx.k++) {
215: if (useports) {
216: PetscDrawViewPortsSet(ports,zctx.k);
217: } else {
218: PetscViewerDrawGetDraw(viewer,zctx.k,&draw);
219: PetscDrawSynchronizedClear(draw);
220: }
222: /*
223: Determine the min and max color in plot
224: */
225: VecStrideMin(xin,zctx.k,PETSC_NULL,&zctx.min);
226: VecStrideMax(xin,zctx.k,PETSC_NULL,&zctx.max);
227: if (zctx.min == zctx.max) {
228: zctx.min -= 1.e-12;
229: zctx.max += 1.e-12;
230: }
232: if (!rank) {
233: char *title;
235: DAGetFieldName(da,zctx.k,&title);
236: if (title) {
237: PetscDrawSetTitle(draw,title);
238: }
239: }
240: PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);
241: PetscInfo2(da,"DA 2d contour plot min %G max %G\n",zctx.min,zctx.max);
243: PetscDrawGetPopup(draw,&popup);
244: if (popup) {PetscDrawScalePopup(popup,zctx.min,zctx.max);}
246: zctx.scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/(zctx.max - zctx.min);
248: PetscDrawZoom(draw,VecView_MPI_Draw_DA2d_Zoom,&zctx);
249: }
250: if (useports){
251: PetscDrawViewPortsDestroy(ports);
252: }
254: VecRestoreArray(xcoorl,&zctx.xy);
255: VecRestoreArray(xlocal,&zctx.v);
256: VecDestroy(xcoor);
257: return(0);
258: }
261: #if defined(PETSC_HAVE_HDF5)
264: PetscErrorCode VecView_MPI_HDF5_DA(Vec xin,PetscViewer viewer)
265: {
267: DA da;
268: hsize_t dim,dims[5];
269: hsize_t count[5];
270: hsize_t offset[5];
271: PetscInt cnt = 0;
272: PetscScalar *x;
273: const char *vecname;
274: hid_t filespace; /* file dataspace identifier */
275: hid_t plist_id; /* property list identifier */
276: hid_t dset_id; /* dataset identifier */
277: hid_t memspace; /* memory dataspace identifier */
278: hid_t file_id;
279: herr_t status;
282: PetscViewerHDF5GetFileId(viewer, &file_id);
283: PetscObjectQuery((PetscObject)xin,"DA",(PetscObject*)&da);
284: if (!da) SETERRQ(PETSC_ERR_ARG_WRONG,"Vector not generated from a DA");
286: /* Create the dataspace for the dataset */
287: dim = PetscHDF5IntCast(da->dim + ((da->w == 1) ? 0 : 1));
288: if (da->dim == 3) dims[cnt++] = PetscHDF5IntCast(da->P);
289: if (da->dim > 1) dims[cnt++] = PetscHDF5IntCast(da->N);
290: dims[cnt++] = PetscHDF5IntCast(da->M);
291: if (da->w > 1) dims[cnt++] = PetscHDF5IntCast(da->w);
292: #if defined(PETSC_USE_COMPLEX)
293: dim++;
294: dims[cnt++] = 2;
295: #endif
296: filespace = H5Screate_simple(dim, dims, NULL);
297: if (filespace == -1) SETERRQ(PETSC_ERR_LIB,"Cannot H5Screate_simple()");
299: /* Create the dataset with default properties and close filespace */
300: PetscObjectGetName((PetscObject)xin,&vecname);
301: #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
302: dset_id = H5Dcreate2(file_id, vecname, H5T_NATIVE_DOUBLE, filespace, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
303: #else
304: dset_id = H5Dcreate(file_id, vecname, H5T_NATIVE_DOUBLE, filespace, H5P_DEFAULT);
305: #endif
306: if (dset_id == -1) SETERRQ(PETSC_ERR_LIB,"Cannot H5Dcreate2()");
307: status = H5Sclose(filespace);CHKERRQ(status);
309: /* Each process defines a dataset and writes it to the hyperslab in the file */
310: cnt = 0;
311: if (da->dim == 3) offset[cnt++] = PetscHDF5IntCast(da->zs);
312: if (da->dim > 1) offset[cnt++] = PetscHDF5IntCast(da->ys);
313: offset[cnt++] = PetscHDF5IntCast(da->xs/da->w);
314: if (da->w > 1) offset[cnt++] = 0;
315: #if defined(PETSC_USE_COMPLEX)
316: offset[cnt++] = 0;
317: #endif
318: cnt = 0;
319: if (da->dim == 3) count[cnt++] = PetscHDF5IntCast(da->ze - da->zs);
320: if (da->dim > 1) count[cnt++] = PetscHDF5IntCast(da->ye - da->ys);
321: count[cnt++] = PetscHDF5IntCast((da->xe - da->xs)/da->w);
322: if (da->w > 1) count[cnt++] = PetscHDF5IntCast(da->w);
323: #if defined(PETSC_USE_COMPLEX)
324: count[cnt++] = 2;
325: #endif
326: memspace = H5Screate_simple(dim, count, NULL);
327: if (memspace == -1) SETERRQ(PETSC_ERR_LIB,"Cannot H5Screate_simple()");
330: filespace = H5Dget_space(dset_id);
331: if (filespace == -1) SETERRQ(PETSC_ERR_LIB,"Cannot H5Dget_space()");
332: status = H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);CHKERRQ(status);
334: /* Create property list for collective dataset write */
335: plist_id = H5Pcreate(H5P_DATASET_XFER);
336: if (plist_id == -1) SETERRQ(PETSC_ERR_LIB,"Cannot H5Pcreate()");
337: #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
338: status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);CHKERRQ(status);
339: #endif
340: /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */
342: VecGetArray(xin, &x);
343: status = H5Dwrite(dset_id, H5T_NATIVE_DOUBLE, memspace, filespace, plist_id, x);CHKERRQ(status);
344: status = H5Fflush(file_id, H5F_SCOPE_GLOBAL);CHKERRQ(status);
345: VecRestoreArray(xin, &x);
347: /* Close/release resources */
348: status = H5Pclose(plist_id);CHKERRQ(status);
349: status = H5Sclose(filespace);CHKERRQ(status)
350: status = H5Sclose(memspace);CHKERRQ(status);
351: status = H5Dclose(dset_id);CHKERRQ(status);
352: PetscInfo1(xin,"Wrote Vec object with name %s\n",vecname);
353: return(0);
354: }
355: #endif
357: #if defined(PETSC_HAVE_PNETCDF)
360: PetscErrorCode VecView_MPI_Netcdf_DA(Vec xin,PetscViewer viewer)
361: {
363: PetscInt ncid,xstart,xdim_num=1;
364: PetscInt dim,m,n,p,dof,swidth,M,N,P;
365: PetscInt xin_dim,xin_id,xin_n,xin_N,xyz_dim,xyz_id,xyz_n,xyz_N;
366: const PetscInt *lx,*ly,*lz;
367: PetscScalar *xarray;
368: DA da,dac;
369: Vec natural,xyz;
370: DAStencilType stencil;
371: DAPeriodicType periodic;
372: MPI_Comm comm;
375: PetscObjectGetComm((PetscObject)xin,&comm);
376: PetscObjectQuery((PetscObject)xin,"DA",(PetscObject*)&da);
377: if (!da) SETERRQ(PETSC_ERR_ARG_WRONG,"Vector not generated from a DA");
378: DAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&dof,&swidth,&periodic,&stencil);
380: /* create the appropriate DA to map the coordinates to natural ordering */
381: DAGetOwnershipRanges(da,&lx,&ly,&lz);
382: if (dim == 1) {
383: DACreate1d(comm,DA_NONPERIODIC,m,dim,0,lx,&dac);
384: } else if (dim == 2) {
385: DACreate2d(comm,DA_NONPERIODIC,DA_STENCIL_BOX,m,n,M,N,dim,0,lx,ly,&dac);
386: } else if (dim == 3) {
387: DACreate3d(comm,DA_NONPERIODIC,DA_STENCIL_BOX,m,n,p,M,N,P,dim,0,lx,ly,lz,&dac);
388: } else {
389: SETERRQ1(PETSC_ERR_ARG_CORRUPT,"Dimension is not 1 2 or 3: %D\n",dim);
390: }
391: DACreateNaturalVector(dac,&xyz);
392: PetscObjectSetOptionsPrefix((PetscObject)xyz,"coor_");
393: DAGlobalToNaturalBegin(dac,da->coordinates,INSERT_VALUES,xyz);
394: DAGlobalToNaturalEnd(dac,da->coordinates,INSERT_VALUES,xyz);
395: /* Create the DA vector in natural ordering */
396: DACreateNaturalVector(da,&natural);
397: DAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural);
398: DAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural);
399: /* Write the netCDF dataset */
400: PetscViewerNetcdfGetID(viewer,&ncid);
401: if (ncid < 0) SETERRQ(PETSC_ERR_ORDER,"First call PetscViewerNetcdfOpen to create NetCDF dataset");
402: /* define dimensions */
403: VecGetSize(xin,&xin_N);
404: VecGetLocalSize(xin,&xin_n);
405: ncmpi_def_dim(ncid,"PETSc_DA_Vector_Global_Size",xin_N,&xin_dim);
406: VecGetSize(xyz,&xyz_N);
407: VecGetLocalSize(xyz,&xyz_n);
408: ncmpi_def_dim(ncid,"PETSc_DA_Coordinate_Vector_Global_Size",xyz_N,&xyz_dim);
409: /* define variables */
410: ncmpi_def_var(ncid,"PETSc_DA_Vector",NC_DOUBLE,xdim_num,&xin_dim,&xin_id);
411: ncmpi_def_var(ncid,"PETSc_DA_Coordinate_Vector",NC_DOUBLE,xdim_num,&xyz_dim,&xyz_id);
412: /* leave define mode */
413: ncmpi_enddef(ncid);
414: /* store the vector */
415: VecGetArray(xin,&xarray);
416: VecGetOwnershipRange(xin,&xstart,PETSC_NULL);
417: ncmpi_put_vara_double_all(ncid,xin_id,(const MPI_Offset*)&xstart,(const MPI_Offset*)&xin_n,xarray);
418: VecRestoreArray(xin,&xarray);
419: /* store the coordinate vector */
420: VecGetArray(xyz,&xarray);
421: VecGetOwnershipRange(xyz,&xstart,PETSC_NULL);
422: ncmpi_put_vara_double_all(ncid,xyz_id,(const MPI_Offset*)&xstart,(const MPI_Offset*)&xyz_n,xarray);
423: VecRestoreArray(xyz,&xarray);
424: /* destroy the vectors and da */
425: VecDestroy(natural);
426: VecDestroy(xyz);
427: DADestroy(dac);
428: return(0);
429: }
430: #endif
432: EXTERN PetscErrorCode VecView_MPI_Draw_DA1d(Vec,PetscViewer);
434: #if defined(PETSC_HAVE_MPIIO)
437: static PetscErrorCode DAArrayMPIIO(DA da,PetscViewer viewer,Vec xin,PetscTruth write)
438: {
440: MPI_File mfdes;
441: PetscMPIInt gsizes[4],lsizes[4],lstarts[4],asiz,dof;
442: MPI_Datatype view;
443: PetscScalar *array;
444: MPI_Offset off;
445: MPI_Aint ub,ul;
446: PetscInt type,rows,vecrows,tr[2];
449: VecGetSize(xin,&vecrows);
450: if (!write) {
451: /* Read vector header. */
452: PetscViewerBinaryRead(viewer,tr,2,PETSC_INT);
453: type = tr[0];
454: rows = tr[1];
455: if (type != VEC_FILE_COOKIE) {
456: SETERRQ(PETSC_ERR_ARG_WRONG,"Not vector next in file");
457: }
458: if (rows != vecrows) SETERRQ(PETSC_ERR_ARG_SIZ,"Vector in file not same size as DA vector");
459: } else {
460: tr[0] = VEC_FILE_COOKIE;
461: tr[1] = vecrows;
462: PetscViewerBinaryWrite(viewer,tr,2,PETSC_INT,PETSC_TRUE);
463: }
465: dof = PetscMPIIntCast(da->w);
466: gsizes[0] = dof; gsizes[1] = PetscMPIIntCast(da->M); gsizes[2] = PetscMPIIntCast(da->N); gsizes[3] = PetscMPIIntCast(da->P);
467: lsizes[0] = dof;lsizes[1] = PetscMPIIntCast((da->xe-da->xs)/dof); lsizes[2] = PetscMPIIntCast(da->ye-da->ys); lsizes[3] = PetscMPIIntCast(da->ze-da->zs);
468: lstarts[0] = 0; lstarts[1] = PetscMPIIntCast(da->xs/dof); lstarts[2] = PetscMPIIntCast(da->ys); lstarts[3] = PetscMPIIntCast(da->zs);
469: MPI_Type_create_subarray(da->dim+1,gsizes,lsizes,lstarts,MPI_ORDER_FORTRAN,MPIU_SCALAR,&view);
470: MPI_Type_commit(&view);
471:
472: PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);
473: PetscViewerBinaryGetMPIIOOffset(viewer,&off);
474: MPI_File_set_view(mfdes,off,MPIU_SCALAR,view,(char *)"native",MPI_INFO_NULL);
475: VecGetArray(xin,&array);
476: asiz = lsizes[1]*(lsizes[2] > 0 ? lsizes[2] : 1)*(lsizes[3] > 0 ? lsizes[3] : 1)*dof;
477: if (write) {
478: MPIU_File_write_all(mfdes,array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);
479: } else {
480: MPIU_File_read_all(mfdes,array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);
481: }
482: MPI_Type_get_extent(view,&ul,&ub);
483: PetscViewerBinaryAddMPIIOOffset(viewer,ub);
484: VecRestoreArray(xin,&array);
485: MPI_Type_free(&view);
486: return(0);
487: }
488: #endif
493: PetscErrorCode VecView_MPI_DA(Vec xin,PetscViewer viewer)
494: {
495: DA da;
497: PetscInt dim;
498: Vec natural;
499: PetscTruth isdraw;
500: #if defined(PETSC_HAVE_HDF5)
501: PetscTruth ishdf5;
502: #endif
503: #if defined(PETSC_HAVE_PNETCDF)
504: PetscTruth isnetcdf;
505: #endif
506: const char *prefix;
509: PetscObjectQuery((PetscObject)xin,"DA",(PetscObject*)&da);
510: if (!da) SETERRQ(PETSC_ERR_ARG_WRONG,"Vector not generated from a DA");
511: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
512: #if defined(PETSC_HAVE_HDF5)
513: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_HDF5,&ishdf5);
514: #endif
515: #if defined(PETSC_HAVE_PNETCDF)
516: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_NETCDF,&isnetcdf);
517: #endif
518: if (isdraw) {
519: DAGetInfo(da,&dim,0,0,0,0,0,0,0,0,0,0);
520: if (dim == 1) {
521: VecView_MPI_Draw_DA1d(xin,viewer);
522: } else if (dim == 2) {
523: VecView_MPI_Draw_DA2d(xin,viewer);
524: } else {
525: SETERRQ1(PETSC_ERR_SUP,"Cannot graphically view vector associated with this dimensional DA %D",dim);
526: }
527: #if defined(PETSC_HAVE_HDF5)
528: } else if (ishdf5) {
529: VecView_MPI_HDF5_DA(xin,viewer);
530: #endif
531: #if defined(PETSC_HAVE_PNETCDF)
532: } else if (isnetcdf) {
533: VecView_MPI_Netcdf_DA(xin,viewer);
534: #endif
535: } else {
536: #if defined(PETSC_HAVE_MPIIO)
537: PetscTruth isbinary,isMPIIO;
539: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
540: if (isbinary) {
541: PetscViewerBinaryGetMPIIO(viewer,&isMPIIO);
542: if (isMPIIO) {
543: DAArrayMPIIO(da,viewer,xin,PETSC_TRUE);
544: return(0);
545: }
546: }
547: #endif
548:
549: /* call viewer on natural ordering */
550: PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);
551: DACreateNaturalVector(da,&natural);
552: PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);
553: DAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural);
554: DAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural);
555: PetscObjectName((PetscObject)xin);
556: PetscObjectSetName((PetscObject)natural,((PetscObject)xin)->name);
557: VecView(natural,viewer);
558: VecDestroy(natural);
559: }
560: return(0);
561: }
564: #if defined(PETSC_HAVE_HDF5)
568: PetscErrorCode VecLoadIntoVector_HDF5_DA(PetscViewer viewer,Vec xin)
569: {
570: DA da;
572: hsize_t dim,dims[5];
573: hsize_t count[5];
574: hsize_t offset[5];
575: PetscInt cnt = 0;
576: PetscScalar *x;
577: const char *vecname;
578: hid_t filespace; /* file dataspace identifier */
579: hid_t plist_id; /* property list identifier */
580: hid_t dset_id; /* dataset identifier */
581: hid_t memspace; /* memory dataspace identifier */
582: hid_t file_id;
583: herr_t status;
586: PetscViewerHDF5GetFileId(viewer, &file_id);
587: PetscObjectQuery((PetscObject)xin,"DA",(PetscObject*)&da);
588: if (!da) SETERRQ(PETSC_ERR_ARG_WRONG,"Vector not generated from a DA");
590: /* Create the dataspace for the dataset */
591: dim = PetscHDF5IntCast(da->dim + ((da->w == 1) ? 0 : 1));
592: dims[cnt++] = PetscHDF5IntCast(da->P);
593: dims[cnt++] = PetscHDF5IntCast(da->N);
594: dims[cnt++] = PetscHDF5IntCast(da->M);
595: if (da->w > 1) PetscHDF5IntCast(dims[cnt++] = da->w);
596: #if defined(PETSC_USE_COMPLEX)
597: dim++;
598: dims[cnt++] = 2;
599: #endif
601: /* Create the dataset with default properties and close filespace */
602: PetscObjectGetName((PetscObject)xin,&vecname);
603: #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
604: dset_id = H5Dopen2(file_id, vecname, H5P_DEFAULT);
605: #else
606: dset_id = H5Dopen(file_id, vecname);
607: #endif
608: if (dset_id == -1) SETERRQ1(PETSC_ERR_LIB,"Cannot H5Dopen2() with Vec named %s",vecname);
609: filespace = H5Dget_space(dset_id);
611: /* Each process defines a dataset and reads it from the hyperslab in the file */
612: cnt = 0;
613: offset[cnt++] = PetscHDF5IntCast(da->zs);
614: offset[cnt++] = PetscHDF5IntCast(da->ys);
615: offset[cnt++] = PetscHDF5IntCast(da->xs/da->w);
616: if (da->w > 1) offset[cnt++] = 0;
617: #if defined(PETSC_USE_COMPLEX)
618: offset[cnt++] = 0;
619: #endif
620: cnt = 0;
621: count[cnt++] = PetscHDF5IntCast(da->ze - da->zs);
622: count[cnt++] = PetscHDF5IntCast(da->ye - da->ys);
623: count[cnt++] = PetscHDF5IntCast((da->xe - da->xs)/da->w);
624: if (da->w > 1) count[cnt++] = PetscHDF5IntCast(da->w);
625: #if defined(PETSC_USE_COMPLEX)
626: count[cnt++] = 2;
627: #endif
628: memspace = H5Screate_simple(dim, count, NULL);
629: if (memspace == -1) SETERRQ(PETSC_ERR_LIB,"Cannot H5Screate_simple()");
631: status = H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);CHKERRQ(status);
633: /* Create property list for collective dataset write */
634: plist_id = H5Pcreate(H5P_DATASET_XFER);
635: if (plist_id == -1) SETERRQ(PETSC_ERR_LIB,"Cannot H5Pcreate()");
636: #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
637: status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);CHKERRQ(status);
638: #endif
639: /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */
641: VecGetArray(xin, &x);
642: status = H5Dread(dset_id, H5T_NATIVE_DOUBLE, memspace, filespace, plist_id, x);CHKERRQ(status);
643: VecRestoreArray(xin, &x);
645: /* Close/release resources */
646: status = H5Pclose(plist_id);CHKERRQ(status);
647: status = H5Sclose(filespace);CHKERRQ(status);
648: status = H5Sclose(memspace);CHKERRQ(status);
649: status = H5Dclose(dset_id);CHKERRQ(status);
650: return(0);
651: }
653: #endif
659: PetscErrorCode VecLoadIntoVector_Binary_DA(PetscViewer viewer,Vec xin)
660: {
661: DA da;
663: Vec natural;
664: const char *prefix;
665: PetscInt bs;
666: PetscTruth flag;
667: #if defined(PETSC_HAVE_HDF5)
668: PetscTruth ishdf5;
669: #endif
670: #if defined(PETSC_HAVE_MPIIO)
671: PetscTruth isMPIIO;
672: #endif
675: PetscObjectQuery((PetscObject)xin,"DA",(PetscObject*)&da);
676: if (!da) SETERRQ(PETSC_ERR_ARG_WRONG,"Vector not generated from a DA");
678: #if defined(PETSC_HAVE_HDF5)
679: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_HDF5,&ishdf5);
680: if (ishdf5) {
681: VecLoadIntoVector_HDF5_DA(viewer,xin);
682: return(0);
683: }
684: #endif
686: #if defined(PETSC_HAVE_MPIIO)
687: PetscViewerBinaryGetMPIIO(viewer,&isMPIIO);
688: if (isMPIIO) {
689: DAArrayMPIIO(da,viewer,xin,PETSC_FALSE);
690: return(0);
691: }
692: #endif
694: PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);
695: DACreateNaturalVector(da,&natural);
696: PetscObjectSetName((PetscObject)natural,((PetscObject)xin)->name);
697: PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);
698: VecLoadIntoVector(viewer,natural);
699: DANaturalToGlobalBegin(da,natural,INSERT_VALUES,xin);
700: DANaturalToGlobalEnd(da,natural,INSERT_VALUES,xin);
701: VecDestroy(natural);
702: PetscInfo(xin,"Loading vector from natural ordering into DA\n");
703: PetscOptionsGetInt(((PetscObject)xin)->prefix,"-vecload_block_size",&bs,&flag);
704: if (flag && bs != da->w) {
705: PetscInfo2(xin,"Block size in file %D not equal to DA's dof %D\n",bs,da->w);
706: }
707: return(0);
708: }