Actual source code: gr2.c
1: #define PETSCDM_DLL
3: /*
4: Plots vectors obtained with DACreate2d()
5: */
7: #include ../src/dm/da/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;
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: PetscViewerDrawGetDraw(viewer,0,&draw);
101: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
103: PetscObjectQuery((PetscObject)xin,"DA",(PetscObject*)&da);
104: if (!da) SETERRQ(PETSC_ERR_ARG_WRONG,"Vector not generated from a DA");
106: PetscObjectGetComm((PetscObject)xin,&comm);
107: MPI_Comm_rank(comm,&rank);
109: DAGetInfo(da,0,&M,&N,0,&zctx.m,&zctx.n,0,&w,&s,&periodic,&st);
110: DAGetOwnershipRanges(da,&lx,&ly,PETSC_NULL);
112: /*
113: Obtain a sequential vector that is going to contain the local values plus ONE layer of
114: ghosted values to draw the graphics from. We also need its corresponding DA (dac) that will
115: update the local values pluse ONE layer of ghost values.
116: */
117: PetscObjectQuery((PetscObject)da,"GraphicsGhosted",(PetscObject*)&xlocal);
118: if (!xlocal) {
119: if (periodic != DA_NONPERIODIC || s != 1 || st != DA_STENCIL_BOX) {
120: /*
121: if original da is not of stencil width one, or periodic or not a box stencil then
122: create a special DA to handle one level of ghost points for graphics
123: */
124: DACreate2d(comm,DA_NONPERIODIC,DA_STENCIL_BOX,M,N,zctx.m,zctx.n,w,1,lx,ly,&dac);
125: PetscInfo(da,"Creating auxilary DA for managing graphics ghost points\n");
126: } else {
127: /* otherwise we can use the da we already have */
128: dac = da;
129: }
130: /* create local vector for holding ghosted values used in graphics */
131: DACreateLocalVector(dac,&xlocal);
132: if (dac != da) {
133: /* don't keep any public reference of this DA, is is only available through xlocal */
134: DADestroy(dac);
135: } else {
136: /* remove association between xlocal and da, because below we compose in the opposite
137: direction and if we left this connect we'd get a loop, so the objects could
138: never be destroyed */
139: PetscObjectCompose((PetscObject)xlocal,"DA",0);
140: }
141: PetscObjectCompose((PetscObject)da,"GraphicsGhosted",(PetscObject)xlocal);
142: PetscObjectDereference((PetscObject)xlocal);
143: } else {
144: if (periodic == DA_NONPERIODIC && s == 1 && st == DA_STENCIL_BOX) {
145: dac = da;
146: } else {
147: PetscObjectQuery((PetscObject)xlocal,"DA",(PetscObject*)&dac);
148: }
149: }
151: /*
152: Get local (ghosted) values of vector
153: */
154: DAGlobalToLocalBegin(dac,xin,INSERT_VALUES,xlocal);
155: DAGlobalToLocalEnd(dac,xin,INSERT_VALUES,xlocal);
156: VecGetArray(xlocal,&zctx.v);
158: /* get coordinates of nodes */
159: DAGetCoordinates(da,&xcoor);
160: if (!xcoor) {
161: DASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,0.0);
162: DAGetCoordinates(da,&xcoor);
163: }
165: /*
166: Determine the min and max coordinates in plot
167: */
168: VecStrideMin(xcoor,0,PETSC_NULL,&xmin);
169: VecStrideMax(xcoor,0,PETSC_NULL,&xmax);
170: VecStrideMin(xcoor,1,PETSC_NULL,&ymin);
171: VecStrideMax(xcoor,1,PETSC_NULL,&ymax);
172: coors[0] = xmin - .05*(xmax- xmin); coors[2] = xmax + .05*(xmax - xmin);
173: coors[1] = ymin - .05*(ymax- ymin); coors[3] = ymax + .05*(ymax - ymin);
174: PetscInfo4(da,"Preparing DA 2d contour plot coordinates %G %G %G %G\n",coors[0],coors[1],coors[2],coors[3]);
176: /*
177: get local ghosted version of coordinates
178: */
179: PetscObjectQuery((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject*)&xcoorl);
180: if (!xcoorl) {
181: /* create DA to get local version of graphics */
182: DACreate2d(comm,DA_NONPERIODIC,DA_STENCIL_BOX,M,N,zctx.m,zctx.n,2,1,lx,ly,&dag);
183: PetscInfo(dag,"Creating auxilary DA for managing graphics coordinates ghost points\n");
184: DACreateLocalVector(dag,&xcoorl);
185: PetscObjectCompose((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject)xcoorl);
186: DADestroy(dag);/* dereference dag */
187: PetscObjectDereference((PetscObject)xcoorl);
188: } else {
189: PetscObjectQuery((PetscObject)xcoorl,"DA",(PetscObject*)&dag);
190: }
191: DAGlobalToLocalBegin(dag,xcoor,INSERT_VALUES,xcoorl);
192: DAGlobalToLocalEnd(dag,xcoor,INSERT_VALUES,xcoorl);
193: VecGetArray(xcoorl,&zctx.xy);
194:
195: /*
196: Get information about size of area each processor must do graphics for
197: */
198: DAGetInfo(dac,0,&M,&N,0,0,0,0,&zctx.step,0,&periodic,0);
199: DAGetGhostCorners(dac,&igstart,&jgstart,0,&zctx.m,&zctx.n,0);
200: DAGetCorners(dac,&istart,0,0,&isize,0,0);
202: PetscOptionsHasName(PETSC_NULL,"-draw_contour_grid",&zctx.showgrid);
204: PetscViewerGetFormat(viewer,&format);
205: PetscOptionsHasName(PETSC_NULL,"-draw_ports",&useports);
206: if (useports || format == PETSC_VIEWER_DRAW_PORTS){
207: PetscDrawSynchronizedClear(draw);
208: PetscDrawViewPortsCreate(draw,zctx.step,&ports);
209: }
210: /*
211: Loop over each field; drawing each in a different window
212: */
213: for (zctx.k=0; zctx.k<zctx.step; zctx.k++) {
214: if (useports) {
215: PetscDrawViewPortsSet(ports,zctx.k);
216: } else {
217: PetscViewerDrawGetDraw(viewer,zctx.k,&draw);
218: PetscDrawSynchronizedClear(draw);
219: }
221: /*
222: Determine the min and max color in plot
223: */
224: VecStrideMin(xin,zctx.k,PETSC_NULL,&zctx.min);
225: VecStrideMax(xin,zctx.k,PETSC_NULL,&zctx.max);
226: if (zctx.min == zctx.max) {
227: zctx.min -= 1.e-12;
228: zctx.max += 1.e-12;
229: }
231: if (!rank) {
232: char *title;
234: DAGetFieldName(da,zctx.k,&title);
235: if (title) {
236: PetscDrawSetTitle(draw,title);
237: }
238: }
239: PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);
240: PetscInfo2(da,"DA 2d contour plot min %G max %G\n",zctx.min,zctx.max);
242: PetscDrawGetPopup(draw,&popup);
243: if (popup) {PetscDrawScalePopup(popup,zctx.min,zctx.max);}
245: zctx.scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/(zctx.max - zctx.min);
247: PetscDrawZoom(draw,VecView_MPI_Draw_DA2d_Zoom,&zctx);
248: }
249: if (useports){
250: PetscDrawViewPortsDestroy(ports);
251: }
253: VecRestoreArray(xcoorl,&zctx.xy);
254: VecRestoreArray(xlocal,&zctx.v);
255: VecDestroy(xcoor);
256: return(0);
257: }
260: #if defined(PETSC_HAVE_HDF5)
263: PetscErrorCode VecView_MPI_HDF5_DA(Vec xin,PetscViewer viewer)
264: {
266: DA da;
267: hsize_t dim,dims[4];
268: hsize_t count[4];
269: hsize_t offset[4];
270: PetscInt cnt = 0;
271: PetscScalar *x;
272: const char *vecname;
273: hid_t filespace; /* file dataspace identifier */
274: hid_t plist_id; /* property list identifier */
275: hid_t dset_id; /* dataset identifier */
276: hid_t memspace; /* memory dataspace identifier */
277: hid_t file_id;
278: herr_t status;
281: PetscViewerHDF5GetFileId(viewer, &file_id);
282: PetscObjectQuery((PetscObject)xin,"DA",(PetscObject*)&da);
283: if (!da) SETERRQ(PETSC_ERR_ARG_WRONG,"Vector not generated from a DA");
285: /* Create the dataspace for the dataset */
286: dim = da->dim + ((da->w == 1) ? 0 : 1);
287: if (da->w > 1) dims[cnt++] = da->w;
288: dims[cnt++] = da->M;
289: dims[cnt++] = da->N;
290: dims[cnt] = da->P;
291: filespace = H5Screate_simple(dim, dims, NULL);
292: if (filespace == -1) SETERRQ(PETSC_ERR_LIB,"Cannot H5Screate_simple()");
294: /* Create the dataset with default properties and close filespace */
295: PetscObjectGetName((PetscObject)xin,&vecname);
296: #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
297: dset_id = H5Dcreate2(file_id, vecname, H5T_NATIVE_DOUBLE, filespace, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
298: #else
299: dset_id = H5Dcreate(file_id, vecname, H5T_NATIVE_DOUBLE, filespace, H5P_DEFAULT);
300: #endif
301: if (dset_id == -1) SETERRQ(PETSC_ERR_LIB,"Cannot H5Dcreate2()");
302: status = H5Sclose(filespace);CHKERRQ(status);
304: /* Each process defines a dataset and writes it to the hyperslab in the file */
305: cnt = 0;
306: if (da->w > 1) offset[cnt++] = 0;
307: offset[cnt++] = da->xs/da->w;
308: offset[cnt++] = da->ys;
309: offset[cnt++] = da->zs;
310: cnt = 0;
311: if (da->w > 1) count[cnt++] = da->w;
312: count[cnt++] = (da->xe - da->xs)/da->w;
313: count[cnt++] = da->ye - da->ys;
314: count[cnt++] = da->ze - da->zs;
315: memspace = H5Screate_simple(dim, count, NULL);
316: if (memspace == -1) SETERRQ(PETSC_ERR_LIB,"Cannot H5Screate_simple()");
319: filespace = H5Dget_space(dset_id);
320: if (filespace == -1) SETERRQ(PETSC_ERR_LIB,"Cannot H5Dget_space()");
321: status = H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);CHKERRQ(status);
323: /* Create property list for collective dataset write */
324: plist_id = H5Pcreate(H5P_DATASET_XFER);
325: if (plist_id == -1) SETERRQ(PETSC_ERR_LIB,"Cannot H5Pcreate()");
326: #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
327: status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);CHKERRQ(status);
328: #endif
329: /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */
331: VecGetArray(xin, &x);
332: status = H5Dwrite(dset_id, H5T_NATIVE_DOUBLE, memspace, filespace, plist_id, x);CHKERRQ(status);
333: status = H5Fflush(file_id, H5F_SCOPE_GLOBAL);CHKERRQ(status);
334: VecRestoreArray(xin, &x);
336: /* Close/release resources */
337: status = H5Pclose(plist_id);CHKERRQ(status);
338: status = H5Sclose(filespace);CHKERRQ(status)
339: status = H5Sclose(memspace);CHKERRQ(status);
340: status = H5Dclose(dset_id);CHKERRQ(status);
341: PetscInfo1(xin,"Wrote Vec object with name %s\n",vecname);
342: return(0);
343: }
344: #endif
346: #if defined(PETSC_HAVE_PNETCDF)
349: PetscErrorCode VecView_MPI_Netcdf_DA(Vec xin,PetscViewer viewer)
350: {
352: PetscInt ncid,xstart,xdim_num=1;
353: PetscInt dim,m,n,p,dof,swidth,M,N,P;
354: PetscInt xin_dim,xin_id,xin_n,xin_N,xyz_dim,xyz_id,xyz_n,xyz_N;
355: const PetscInt *lx,*ly,*lz;
356: PetscScalar *xarray;
357: DA da,dac;
358: Vec natural,xyz;
359: DAStencilType stencil;
360: DAPeriodicType periodic;
361: MPI_Comm comm;
364: PetscObjectGetComm((PetscObject)xin,&comm);
365: PetscObjectQuery((PetscObject)xin,"DA",(PetscObject*)&da);
366: if (!da) SETERRQ(PETSC_ERR_ARG_WRONG,"Vector not generated from a DA");
367: DAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&dof,&swidth,&periodic,&stencil);
369: /* create the appropriate DA to map the coordinates to natural ordering */
370: DAGetOwnershipRanges(da,&lx,&ly,&lz);
371: if (dim == 1) {
372: DACreate1d(comm,DA_NONPERIODIC,m,dim,0,lx,&dac);
373: } else if (dim == 2) {
374: DACreate2d(comm,DA_NONPERIODIC,DA_STENCIL_BOX,m,n,M,N,dim,0,lx,ly,&dac);
375: } else if (dim == 3) {
376: DACreate3d(comm,DA_NONPERIODIC,DA_STENCIL_BOX,m,n,p,M,N,P,dim,0,lx,ly,lz,&dac);
377: } else {
378: SETERRQ1(PETSC_ERR_ARG_CORRUPT,"Dimension is not 1 2 or 3: %D\n",dim);
379: }
380: DACreateNaturalVector(dac,&xyz);
381: PetscObjectSetOptionsPrefix((PetscObject)xyz,"coor_");
382: DAGlobalToNaturalBegin(dac,da->coordinates,INSERT_VALUES,xyz);
383: DAGlobalToNaturalEnd(dac,da->coordinates,INSERT_VALUES,xyz);
384: /* Create the DA vector in natural ordering */
385: DACreateNaturalVector(da,&natural);
386: DAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural);
387: DAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural);
388: /* Write the netCDF dataset */
389: PetscViewerNetcdfGetID(viewer,&ncid);
390: if (ncid < 0) SETERRQ(PETSC_ERR_ORDER,"First call PetscViewerNetcdfOpen to create NetCDF dataset");
391: /* define dimensions */
392: VecGetSize(xin,&xin_N);
393: VecGetLocalSize(xin,&xin_n);
394: ncmpi_def_dim(ncid,"PETSc_DA_Vector_Global_Size",xin_N,&xin_dim);
395: VecGetSize(xyz,&xyz_N);
396: VecGetLocalSize(xyz,&xyz_n);
397: ncmpi_def_dim(ncid,"PETSc_DA_Coordinate_Vector_Global_Size",xyz_N,&xyz_dim);
398: /* define variables */
399: ncmpi_def_var(ncid,"PETSc_DA_Vector",NC_DOUBLE,xdim_num,&xin_dim,&xin_id);
400: ncmpi_def_var(ncid,"PETSc_DA_Coordinate_Vector",NC_DOUBLE,xdim_num,&xyz_dim,&xyz_id);
401: /* leave define mode */
402: ncmpi_enddef(ncid);
403: /* store the vector */
404: VecGetArray(xin,&xarray);
405: VecGetOwnershipRange(xin,&xstart,PETSC_NULL);
406: ncmpi_put_vara_double_all(ncid,xin_id,(const MPI_Offset*)&xstart,(const MPI_Offset*)&xin_n,xarray);
407: VecRestoreArray(xin,&xarray);
408: /* store the coordinate vector */
409: VecGetArray(xyz,&xarray);
410: VecGetOwnershipRange(xyz,&xstart,PETSC_NULL);
411: ncmpi_put_vara_double_all(ncid,xyz_id,(const MPI_Offset*)&xstart,(const MPI_Offset*)&xyz_n,xarray);
412: VecRestoreArray(xyz,&xarray);
413: /* destroy the vectors and da */
414: VecDestroy(natural);
415: VecDestroy(xyz);
416: DADestroy(dac);
417: return(0);
418: }
419: #endif
421: EXTERN PetscErrorCode VecView_MPI_Draw_DA1d(Vec,PetscViewer);
423: #if defined(PETSC_HAVE_MPIIO)
426: static PetscErrorCode DAArrayMPIIO(DA da,PetscViewer viewer,Vec xin,PetscTruth write)
427: {
429: MPI_File mfdes;
430: PetscMPIInt gsizes[4],lsizes[4],lstarts[4],asiz,dof;
431: MPI_Datatype view;
432: PetscScalar *array;
433: MPI_Offset off;
434: MPI_Aint ub,ul;
435: PetscInt type,rows,vecrows,tr[2];
438: VecGetSize(xin,&vecrows);
439: if (!write) {
440: /* Read vector header. */
441: PetscViewerBinaryRead(viewer,tr,2,PETSC_INT);
442: type = tr[0];
443: rows = tr[1];
444: if (type != VEC_FILE_COOKIE) {
445: SETERRQ(PETSC_ERR_ARG_WRONG,"Not vector next in file");
446: }
447: if (rows != vecrows) SETERRQ(PETSC_ERR_ARG_SIZ,"Vector in file not same size as DA vector");
448: } else {
449: tr[0] = VEC_FILE_COOKIE;
450: tr[1] = vecrows;
451: PetscViewerBinaryWrite(viewer,tr,2,PETSC_INT,PETSC_TRUE);
452: }
454: dof = PetscMPIIntCast(da->w);
455: gsizes[0] = dof; gsizes[1] = PetscMPIIntCast(da->M); gsizes[2] = PetscMPIIntCast(da->N); gsizes[3] = PetscMPIIntCast(da->P);
456: 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);
457: lstarts[0] = 0; lstarts[1] = PetscMPIIntCast(da->xs/dof); lstarts[2] = PetscMPIIntCast(da->ys); lstarts[3] = PetscMPIIntCast(da->zs);
458: MPI_Type_create_subarray(da->dim+1,gsizes,lsizes,lstarts,MPI_ORDER_FORTRAN,MPIU_SCALAR,&view);
459: MPI_Type_commit(&view);
460:
461: PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);
462: PetscViewerBinaryGetMPIIOOffset(viewer,&off);
463: MPI_File_set_view(mfdes,off,MPIU_SCALAR,view,(char *)"native",MPI_INFO_NULL);
464: VecGetArray(xin,&array);
465: asiz = lsizes[1]*(lsizes[2] > 0 ? lsizes[2] : 1)*(lsizes[3] > 0 ? lsizes[3] : 1)*dof;
466: if (write) {
467: MPIU_File_write_all(mfdes,array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);
468: } else {
469: MPIU_File_read_all(mfdes,array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);
470: }
471: MPI_Type_get_extent(view,&ul,&ub);
472: PetscViewerBinaryAddMPIIOOffset(viewer,ub);
473: VecRestoreArray(xin,&array);
474: MPI_Type_free(&view);
475: return(0);
476: }
477: #endif
482: PetscErrorCode VecView_MPI_DA(Vec xin,PetscViewer viewer)
483: {
484: DA da;
486: PetscInt dim;
487: Vec natural;
488: PetscTruth isdraw;
489: #if defined(PETSC_HAVE_HDF5)
490: PetscTruth ishdf5;
491: #endif
492: #if defined(PETSC_HAVE_PNETCDF)
493: PetscTruth isnetcdf;
494: #endif
495: const char *prefix;
498: PetscObjectQuery((PetscObject)xin,"DA",(PetscObject*)&da);
499: if (!da) SETERRQ(PETSC_ERR_ARG_WRONG,"Vector not generated from a DA");
500: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
501: #if defined(PETSC_HAVE_HDF5)
502: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_HDF5,&ishdf5);
503: #endif
504: #if defined(PETSC_HAVE_PNETCDF)
505: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_NETCDF,&isnetcdf);
506: #endif
507: if (isdraw) {
508: DAGetInfo(da,&dim,0,0,0,0,0,0,0,0,0,0);
509: if (dim == 1) {
510: VecView_MPI_Draw_DA1d(xin,viewer);
511: } else if (dim == 2) {
512: VecView_MPI_Draw_DA2d(xin,viewer);
513: } else {
514: SETERRQ1(PETSC_ERR_SUP,"Cannot graphically view vector associated with this dimensional DA %D",dim);
515: }
516: #if defined(PETSC_HAVE_HDF5)
517: } else if (ishdf5) {
518: VecView_MPI_HDF5_DA(xin,viewer);
519: #endif
520: #if defined(PETSC_HAVE_PNETCDF)
521: } else if (isnetcdf) {
522: VecView_MPI_Netcdf_DA(xin,viewer);
523: #endif
524: } else {
525: #if defined(PETSC_HAVE_MPIIO)
526: PetscTruth isbinary,isMPIIO;
528: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
529: if (isbinary) {
530: PetscViewerBinaryGetMPIIO(viewer,&isMPIIO);
531: if (isMPIIO) {
532: DAArrayMPIIO(da,viewer,xin,PETSC_TRUE);
533: return(0);
534: }
535: }
536: #endif
537:
538: /* call viewer on natural ordering */
539: PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);
540: DACreateNaturalVector(da,&natural);
541: PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);
542: DAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural);
543: DAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural);
544: PetscObjectName((PetscObject)xin);
545: PetscObjectSetName((PetscObject)natural,((PetscObject)xin)->name);
546: VecView(natural,viewer);
547: VecDestroy(natural);
548: }
549: return(0);
550: }
553: #if defined(PETSC_HAVE_HDF5)
557: PetscErrorCode VecLoadIntoVector_HDF5_DA(PetscViewer viewer,Vec xin)
558: {
559: DA da;
561: hsize_t dim,dims[4];
562: hsize_t count[4];
563: hsize_t offset[4];
564: PetscInt cnt = 0;
565: PetscScalar *x;
566: const char *vecname;
567: hid_t filespace; /* file dataspace identifier */
568: hid_t plist_id; /* property list identifier */
569: hid_t dset_id; /* dataset identifier */
570: hid_t memspace; /* memory dataspace identifier */
571: hid_t file_id;
572: herr_t status;
575: PetscViewerHDF5GetFileId(viewer, &file_id);
576: PetscObjectQuery((PetscObject)xin,"DA",(PetscObject*)&da);
577: if (!da) SETERRQ(PETSC_ERR_ARG_WRONG,"Vector not generated from a DA");
579: /* Create the dataspace for the dataset */
580: dim = da->dim + ((da->w == 1) ? 0 : 1);
581: if (da->w > 1) dims[cnt++] = da->w;
582: dims[cnt++] = da->M;
583: dims[cnt++] = da->N;
584: dims[cnt] = da->P;
586: /* Create the dataset with default properties and close filespace */
587: PetscObjectGetName((PetscObject)xin,&vecname);
588: #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
589: dset_id = H5Dopen2(file_id, vecname, H5P_DEFAULT);
590: #else
591: dset_id = H5Dopen(file_id, vecname);
592: #endif
593: if (dset_id == -1) SETERRQ1(PETSC_ERR_LIB,"Cannot H5Dopen2() with Vec named %s",vecname);
594: filespace = H5Dget_space(dset_id);
596: /* Each process defines a dataset and reads it from the hyperslab in the file */
597: cnt = 0;
598: if (da->w > 1) offset[cnt++] = 0;
599: offset[cnt++] = da->xs/da->w;
600: offset[cnt++] = da->ys;
601: offset[cnt++] = da->zs;
602: cnt = 0;
603: if (da->w > 1) count[cnt++] = da->w;
604: count[cnt++] = (da->xe - da->xs)/da->w;
605: count[cnt++] = da->ye - da->ys;
606: count[cnt++] = da->ze - da->zs;
607: memspace = H5Screate_simple(dim, count, NULL);
608: if (memspace == -1) SETERRQ(PETSC_ERR_LIB,"Cannot H5Screate_simple()");
610: status = H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);CHKERRQ(status);
612: /* Create property list for collective dataset write */
613: plist_id = H5Pcreate(H5P_DATASET_XFER);
614: if (plist_id == -1) SETERRQ(PETSC_ERR_LIB,"Cannot H5Pcreate()");
615: #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
616: status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);CHKERRQ(status);
617: #endif
618: /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */
620: VecGetArray(xin, &x);
621: status = H5Dread(dset_id, H5T_NATIVE_DOUBLE, memspace, filespace, plist_id, x);CHKERRQ(status);
622: VecRestoreArray(xin, &x);
624: /* Close/release resources */
625: status = H5Pclose(plist_id);CHKERRQ(status);
626: status = H5Sclose(filespace);CHKERRQ(status);
627: status = H5Sclose(memspace);CHKERRQ(status);
628: status = H5Dclose(dset_id);CHKERRQ(status);
629: return(0);
630: }
632: #endif
638: PetscErrorCode VecLoadIntoVector_Binary_DA(PetscViewer viewer,Vec xin)
639: {
640: DA da;
642: Vec natural;
643: const char *prefix;
644: PetscInt bs;
645: PetscTruth flag;
646: #if defined(PETSC_HAVE_HDF5)
647: PetscTruth ishdf5;
648: #endif
649: #if defined(PETSC_HAVE_MPIIO)
650: PetscTruth isMPIIO;
651: #endif
654: PetscObjectQuery((PetscObject)xin,"DA",(PetscObject*)&da);
655: if (!da) SETERRQ(PETSC_ERR_ARG_WRONG,"Vector not generated from a DA");
657: #if defined(PETSC_HAVE_HDF5)
658: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_HDF5,&ishdf5);
659: if (ishdf5) {
660: VecLoadIntoVector_HDF5_DA(viewer,xin);
661: return(0);
662: }
663: #endif
665: #if defined(PETSC_HAVE_MPIIO)
666: PetscViewerBinaryGetMPIIO(viewer,&isMPIIO);
667: if (isMPIIO) {
668: DAArrayMPIIO(da,viewer,xin,PETSC_FALSE);
669: return(0);
670: }
671: #endif
673: PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);
674: DACreateNaturalVector(da,&natural);
675: PetscObjectSetName((PetscObject)natural,((PetscObject)xin)->name);
676: PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);
677: VecLoadIntoVector(viewer,natural);
678: DANaturalToGlobalBegin(da,natural,INSERT_VALUES,xin);
679: DANaturalToGlobalEnd(da,natural,INSERT_VALUES,xin);
680: VecDestroy(natural);
681: PetscInfo(xin,"Loading vector from natural ordering into DA\n");
682: PetscOptionsGetInt(((PetscObject)xin)->prefix,"-vecload_block_size",&bs,&flag);
683: if (flag && bs != da->w) {
684: PetscInfo2(xin,"Block size in file %D not equal to DA's dof %D\n",bs,da->w);
685: }
686: return(0);
687: }