1 2 /* 3 Plots vectors obtained with DMDACreate2d() 4 */ 5 6 #include <petsc-private/daimpl.h> /*I "petscdmda.h" I*/ 7 #include <petsc-private/vecimpl.h> 8 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; 18 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 */ 24 #undef __FUNCT__ 25 #define __FUNCT__ "VecView_MPI_Draw_DA2d_Zoom" 26 PetscErrorCode VecView_MPI_Draw_DA2d_Zoom(PetscDraw draw,void *ctx) 27 { 28 ZoomCtx *zctx = (ZoomCtx*)ctx; 29 PetscErrorCode ierr; 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; 33 34 PetscFunctionBegin; 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 ierr = PetscDrawTriangle(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3);CHKERRQ(ierr); 59 ierr = PetscDrawTriangle(draw,x1,y_1,x3,y3,x4,y4,c1,c3,c4);CHKERRQ(ierr); 60 if (zctx->showgrid) { 61 ierr = PetscDrawLine(draw,x1,y_1,x2,y2,PETSC_DRAW_BLACK);CHKERRQ(ierr); 62 ierr = PetscDrawLine(draw,x2,y2,x3,y3,PETSC_DRAW_BLACK);CHKERRQ(ierr); 63 ierr = PetscDrawLine(draw,x3,y3,x4,y4,PETSC_DRAW_BLACK);CHKERRQ(ierr); 64 ierr = PetscDrawLine(draw,x4,y4,x1,y_1,PETSC_DRAW_BLACK);CHKERRQ(ierr); 65 } 66 } 67 } 68 PetscFunctionReturn(0); 69 } 70 71 #undef __FUNCT__ 72 #define __FUNCT__ "VecView_MPI_Draw_DA2d" 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; 94 95 PetscFunctionBegin; 96 zctx.showgrid = PETSC_FALSE; 97 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 98 ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 99 ierr = PetscViewerDrawGetBounds(viewer,&nbounds,&bounds);CHKERRQ(ierr); 100 101 ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 102 if (!da) SETERRQ(((PetscObject)xin)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 103 104 ierr = PetscObjectGetComm((PetscObject)xin,&comm);CHKERRQ(ierr); 105 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 106 107 ierr = DMDAGetInfo(da,0,&M,&N,0,&zctx.m,&zctx.n,0,&w,&s,&bx,&by,0,&st);CHKERRQ(ierr); 108 ierr = DMDAGetOwnershipRanges(da,&lx,&ly,PETSC_NULL);CHKERRQ(ierr); 109 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 ierr = PetscObjectQuery((PetscObject)da,"GraphicsGhosted",(PetscObject*)&xlocal);CHKERRQ(ierr); 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 ierr = DMDACreate2d(comm,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,zctx.m,zctx.n,w,1,lx,ly,&dac);CHKERRQ(ierr); 123 ierr = PetscInfo(da,"Creating auxilary DMDA for managing graphics ghost points\n");CHKERRQ(ierr); 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 ierr = DMCreateLocalVector(dac,&xlocal);CHKERRQ(ierr); 130 if (dac != da) { 131 /* don't keep any public reference of this DMDA, is is only available through xlocal */ 132 ierr = PetscObjectDereference((PetscObject)dac);CHKERRQ(ierr); 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 ierr = PetscObjectRemoveReference((PetscObject)xlocal,"__PETSc_dm");CHKERRQ(ierr); 138 } 139 ierr = PetscObjectCompose((PetscObject)da,"GraphicsGhosted",(PetscObject)xlocal);CHKERRQ(ierr); 140 ierr = PetscObjectDereference((PetscObject)xlocal);CHKERRQ(ierr); 141 } else { 142 if (bx != DMDA_BOUNDARY_NONE || by != DMDA_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) { 143 ierr = VecGetDM(xlocal, &dac);CHKERRQ(ierr); 144 } else { 145 dac = da; 146 } 147 } 148 149 /* 150 Get local (ghosted) values of vector 151 */ 152 ierr = DMGlobalToLocalBegin(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr); 153 ierr = DMGlobalToLocalEnd(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr); 154 ierr = VecGetArray(xlocal,&zctx.v);CHKERRQ(ierr); 155 156 /* get coordinates of nodes */ 157 ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr); 158 if (!xcoor) { 159 ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,0.0);CHKERRQ(ierr); 160 ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr); 161 } 162 163 /* 164 Determine the min and max coordinates in plot 165 */ 166 ierr = VecStrideMin(xcoor,0,PETSC_NULL,&xmin);CHKERRQ(ierr); 167 ierr = VecStrideMax(xcoor,0,PETSC_NULL,&xmax);CHKERRQ(ierr); 168 ierr = VecStrideMin(xcoor,1,PETSC_NULL,&ymin);CHKERRQ(ierr); 169 ierr = VecStrideMax(xcoor,1,PETSC_NULL,&ymax);CHKERRQ(ierr); 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 ierr = PetscInfo4(da,"Preparing DMDA 2d contour plot coordinates %G %G %G %G\n",coors[0],coors[1],coors[2],coors[3]);CHKERRQ(ierr); 173 174 /* 175 get local ghosted version of coordinates 176 */ 177 ierr = PetscObjectQuery((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject*)&xcoorl);CHKERRQ(ierr); 178 if (!xcoorl) { 179 /* create DMDA to get local version of graphics */ 180 ierr = DMDACreate2d(comm,DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,zctx.m,zctx.n,2,1,lx,ly,&dag);CHKERRQ(ierr); 181 ierr = PetscInfo(dag,"Creating auxilary DMDA for managing graphics coordinates ghost points\n");CHKERRQ(ierr); 182 ierr = DMCreateLocalVector(dag,&xcoorl);CHKERRQ(ierr); 183 ierr = PetscObjectCompose((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject)xcoorl);CHKERRQ(ierr); 184 ierr = PetscObjectDereference((PetscObject)dag);CHKERRQ(ierr); 185 ierr = PetscObjectDereference((PetscObject)xcoorl);CHKERRQ(ierr); 186 } else { 187 ierr = VecGetDM(xcoorl,&dag);CHKERRQ(ierr); 188 } 189 ierr = DMGlobalToLocalBegin(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr); 190 ierr = DMGlobalToLocalEnd(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr); 191 ierr = VecGetArray(xcoorl,&zctx.xy);CHKERRQ(ierr); 192 193 /* 194 Get information about size of area each processor must do graphics for 195 */ 196 ierr = DMDAGetInfo(dac,0,&M,&N,0,0,0,0,&zctx.step,0,&bx,&by,0,0);CHKERRQ(ierr); 197 ierr = DMDAGetGhostCorners(dac,0,0,0,&zctx.m,&zctx.n,0);CHKERRQ(ierr); 198 199 ierr = PetscOptionsGetBool(PETSC_NULL,"-draw_contour_grid",&zctx.showgrid,PETSC_NULL);CHKERRQ(ierr); 200 ierr = PetscMalloc(zctx.step*sizeof(PetscInt),&displayfields);CHKERRQ(ierr); 201 for (i=0; i<zctx.step; i++) displayfields[i] = i; 202 ndisplayfields = zctx.step; 203 ierr = PetscOptionsGetIntArray(PETSC_NULL,"-draw_fields",displayfields,&ndisplayfields,&flg);CHKERRQ(ierr); 204 if (!flg) ndisplayfields = zctx.step; 205 206 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 207 ierr = PetscOptionsGetBool(PETSC_NULL,"-draw_ports",&useports,PETSC_NULL);CHKERRQ(ierr); 208 if (useports || format == PETSC_VIEWER_DRAW_PORTS){ 209 ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr); 210 ierr = PetscDrawViewPortsCreate(draw,ndisplayfields,&ports);CHKERRQ(ierr); 211 } 212 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 ierr = PetscDrawViewPortsSet(ports,i);CHKERRQ(ierr); 220 } else { 221 ierr = PetscViewerDrawGetDraw(viewer,i,&draw);CHKERRQ(ierr); 222 } 223 224 /* 225 Determine the min and max color in plot 226 */ 227 ierr = VecStrideMin(xin,zctx.k,PETSC_NULL,&zctx.min);CHKERRQ(ierr); 228 ierr = VecStrideMax(xin,zctx.k,PETSC_NULL,&zctx.max);CHKERRQ(ierr); 229 if (zctx.k < nbounds) { 230 zctx.min = PetscMin(zctx.min,bounds[2*zctx.k]); 231 zctx.max = PetscMax(zctx.max,bounds[2*zctx.k+1]); 232 } 233 if (zctx.min == zctx.max) { 234 zctx.min -= 1.e-12; 235 zctx.max += 1.e-12; 236 } 237 238 if (!rank) { 239 const char *title; 240 241 ierr = DMDAGetFieldName(da,zctx.k,&title);CHKERRQ(ierr); 242 if (title) { 243 ierr = PetscDrawSetTitle(draw,title);CHKERRQ(ierr); 244 } 245 } 246 ierr = PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);CHKERRQ(ierr); 247 ierr = PetscInfo2(da,"DMDA 2d contour plot min %G max %G\n",zctx.min,zctx.max);CHKERRQ(ierr); 248 249 ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 250 if (popup) {ierr = PetscDrawScalePopup(popup,zctx.min,zctx.max);CHKERRQ(ierr);} 251 252 zctx.scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/(zctx.max - zctx.min); 253 254 ierr = PetscDrawZoom(draw,VecView_MPI_Draw_DA2d_Zoom,&zctx);CHKERRQ(ierr); 255 } 256 ierr = PetscFree(displayfields);CHKERRQ(ierr); 257 ierr = PetscDrawViewPortsDestroy(ports);CHKERRQ(ierr); 258 259 ierr = VecRestoreArray(xcoorl,&zctx.xy);CHKERRQ(ierr); 260 ierr = VecRestoreArray(xlocal,&zctx.v);CHKERRQ(ierr); 261 PetscFunctionReturn(0); 262 } 263 264 265 #if defined(PETSC_HAVE_HDF5) 266 #undef __FUNCT__ 267 #define __FUNCT__ "VecView_MPI_HDF5_DA" 268 PetscErrorCode VecView_MPI_HDF5_DA(Vec xin,PetscViewer viewer) 269 { 270 DM dm; 271 DM_DA *da; 272 hid_t filespace; /* file dataspace identifier */ 273 hid_t chunkspace; /* chunk dataset property 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 hid_t group; 279 hid_t scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */ 280 herr_t status; 281 hsize_t i, dim; 282 hsize_t maxDims[6], dims[6], chunkDims[6], count[6], offset[6]; 283 PetscInt timestep; 284 PetscScalar *x; 285 const char *vecname; 286 PetscErrorCode ierr; 287 288 PetscFunctionBegin; 289 ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr); 290 ierr = PetscViewerHDF5GetTimestep(viewer, ×tep);CHKERRQ(ierr); 291 292 ierr = VecGetDM(xin,&dm);CHKERRQ(ierr); 293 if (!dm) SETERRQ(((PetscObject)xin)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 294 da = (DM_DA*)dm->data; 295 296 /* Create the dataspace for the dataset. 297 * 298 * dims - holds the current dimensions of the dataset 299 * 300 * maxDims - holds the maximum dimensions of the dataset (unlimited 301 * for the number of time steps with the current dimensions for the 302 * other dimensions; so only additional time steps can be added). 303 * 304 * chunkDims - holds the size of a single time step (required to 305 * permit extending dataset). 306 */ 307 dim = 0; 308 if (timestep >= 0) { 309 dims[dim] = timestep+1; 310 maxDims[dim] = H5S_UNLIMITED; 311 chunkDims[dim] = 1; 312 ++dim; 313 } 314 if (da->dim == 3) { 315 dims[dim] = PetscHDF5IntCast(da->P); 316 maxDims[dim] = dims[dim]; 317 chunkDims[dim] = dims[dim]; 318 ++dim; 319 } 320 if (da->dim > 1) { 321 dims[dim] = PetscHDF5IntCast(da->N); 322 maxDims[dim] = dims[dim]; 323 chunkDims[dim] = dims[dim]; 324 ++dim; 325 } 326 dims[dim] = PetscHDF5IntCast(da->M); 327 maxDims[dim] = dims[dim]; 328 chunkDims[dim] = dims[dim]; 329 ++dim; 330 if (da->w > 1) { 331 dims[dim] = PetscHDF5IntCast(da->w); 332 maxDims[dim] = dims[dim]; 333 chunkDims[dim] = dims[dim]; 334 ++dim; 335 } 336 #if defined(PETSC_USE_COMPLEX) 337 dims[dim] = 2; 338 maxDims[dim] = dims[dim]; 339 chunkDims[dim] = dims[dim]; 340 ++dim; 341 #endif 342 for (i=0; i < dim; ++i) 343 filespace = H5Screate_simple(dim, dims, maxDims); 344 if (filespace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Screate_simple()"); 345 346 #if defined(PETSC_USE_REAL_SINGLE) 347 scalartype = H5T_NATIVE_FLOAT; 348 #elif defined(PETSC_USE_REAL___FLOAT128) 349 #error "HDF5 output with 128 bit floats not supported." 350 #else 351 scalartype = H5T_NATIVE_DOUBLE; 352 #endif 353 354 /* Create the dataset with default properties and close filespace */ 355 ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr); 356 if (!H5Lexists(group, vecname, H5P_DEFAULT)) { 357 /* Create chunk */ 358 chunkspace = H5Pcreate(H5P_DATASET_CREATE); 359 if (chunkspace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Pcreate()"); 360 status = H5Pset_chunk(chunkspace, dim, chunkDims); CHKERRQ(status); 361 362 #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 363 dset_id = H5Dcreate2(group, vecname, scalartype, filespace, H5P_DEFAULT, chunkspace, H5P_DEFAULT); 364 #else 365 dset_id = H5Dcreate(group, vecname, scalartype, filespace, H5P_DEFAULT); 366 #endif 367 if (dset_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Dcreate2()"); 368 } else { 369 dset_id = H5Dopen2(group, vecname, H5P_DEFAULT); 370 status = H5Dset_extent(dset_id, dims);CHKERRQ(status); 371 } 372 status = H5Sclose(filespace);CHKERRQ(status); 373 374 /* Each process defines a dataset and writes it to the hyperslab in the file */ 375 dim = 0; 376 if (timestep >= 0) { 377 offset[dim] = timestep; 378 ++dim; 379 } 380 if (da->dim == 3) offset[dim++] = PetscHDF5IntCast(da->zs); 381 if (da->dim > 1) offset[dim++] = PetscHDF5IntCast(da->ys); 382 offset[dim++] = PetscHDF5IntCast(da->xs/da->w); 383 if (da->w > 1) offset[dim++] = 0; 384 #if defined(PETSC_USE_COMPLEX) 385 offset[dim++] = 0; 386 #endif 387 dim = 0; 388 if (timestep >= 0) { 389 count[dim] = 1; 390 ++dim; 391 } 392 if (da->dim == 3) count[dim++] = PetscHDF5IntCast(da->ze - da->zs); 393 if (da->dim > 1) count[dim++] = PetscHDF5IntCast(da->ye - da->ys); 394 count[dim++] = PetscHDF5IntCast((da->xe - da->xs)/da->w); 395 if (da->w > 1) count[dim++] = PetscHDF5IntCast(da->w); 396 #if defined(PETSC_USE_COMPLEX) 397 count[dim++] = 2; 398 #endif 399 memspace = H5Screate_simple(dim, count, NULL); 400 if (memspace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Screate_simple()"); 401 402 filespace = H5Dget_space(dset_id); 403 if (filespace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Dget_space()"); 404 status = H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);CHKERRQ(status); 405 406 /* Create property list for collective dataset write */ 407 plist_id = H5Pcreate(H5P_DATASET_XFER); 408 if (plist_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Pcreate()"); 409 #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 410 status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);CHKERRQ(status); 411 #endif 412 /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 413 414 ierr = VecGetArray(xin, &x);CHKERRQ(ierr); 415 status = H5Dwrite(dset_id, scalartype, memspace, filespace, plist_id, x);CHKERRQ(status); 416 status = H5Fflush(file_id, H5F_SCOPE_GLOBAL);CHKERRQ(status); 417 ierr = VecRestoreArray(xin, &x);CHKERRQ(ierr); 418 419 /* Close/release resources */ 420 if (group != file_id) { 421 status = H5Gclose(group);CHKERRQ(status); 422 } 423 status = H5Pclose(plist_id);CHKERRQ(status); 424 status = H5Sclose(filespace);CHKERRQ(status); 425 status = H5Sclose(memspace);CHKERRQ(status); 426 status = H5Dclose(dset_id);CHKERRQ(status); 427 ierr = PetscInfo1(xin,"Wrote Vec object with name %s\n",vecname);CHKERRQ(ierr); 428 PetscFunctionReturn(0); 429 } 430 #endif 431 432 extern PetscErrorCode VecView_MPI_Draw_DA1d(Vec,PetscViewer); 433 434 #if defined(PETSC_HAVE_MPIIO) 435 #undef __FUNCT__ 436 #define __FUNCT__ "DMDAArrayMPIIO" 437 static PetscErrorCode DMDAArrayMPIIO(DM da,PetscViewer viewer,Vec xin,PetscBool write) 438 { 439 PetscErrorCode ierr; 440 MPI_File mfdes; 441 PetscMPIInt gsizes[4],lsizes[4],lstarts[4],asiz,dof; 442 MPI_Datatype view; 443 const PetscScalar *array; 444 MPI_Offset off; 445 MPI_Aint ub,ul; 446 PetscInt type,rows,vecrows,tr[2]; 447 DM_DA *dd = (DM_DA*)da->data; 448 449 PetscFunctionBegin; 450 ierr = VecGetSize(xin,&vecrows);CHKERRQ(ierr); 451 if (!write) { 452 /* Read vector header. */ 453 ierr = PetscViewerBinaryRead(viewer,tr,2,PETSC_INT);CHKERRQ(ierr); 454 type = tr[0]; 455 rows = tr[1]; 456 if (type != VEC_FILE_CLASSID) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Not vector next in file"); 457 if (rows != vecrows) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_SIZ,"Vector in file not same size as DMDA vector"); 458 } else { 459 tr[0] = VEC_FILE_CLASSID; 460 tr[1] = vecrows; 461 ierr = PetscViewerBinaryWrite(viewer,tr,2,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 462 } 463 464 dof = PetscMPIIntCast(dd->w); 465 gsizes[0] = dof; gsizes[1] = PetscMPIIntCast(dd->M); gsizes[2] = PetscMPIIntCast(dd->N); gsizes[3] = PetscMPIIntCast(dd->P); 466 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); 467 lstarts[0] = 0; lstarts[1] = PetscMPIIntCast(dd->xs/dof); lstarts[2] = PetscMPIIntCast(dd->ys); lstarts[3] = PetscMPIIntCast(dd->zs); 468 ierr = MPI_Type_create_subarray(dd->dim+1,gsizes,lsizes,lstarts,MPI_ORDER_FORTRAN,MPIU_SCALAR,&view);CHKERRQ(ierr); 469 ierr = MPI_Type_commit(&view);CHKERRQ(ierr); 470 471 ierr = PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);CHKERRQ(ierr); 472 ierr = PetscViewerBinaryGetMPIIOOffset(viewer,&off);CHKERRQ(ierr); 473 ierr = MPI_File_set_view(mfdes,off,MPIU_SCALAR,view,(char *)"native",MPI_INFO_NULL);CHKERRQ(ierr); 474 ierr = VecGetArrayRead(xin,&array);CHKERRQ(ierr); 475 asiz = lsizes[1]*(lsizes[2] > 0 ? lsizes[2] : 1)*(lsizes[3] > 0 ? lsizes[3] : 1)*dof; 476 if (write) { 477 ierr = MPIU_File_write_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr); 478 } else { 479 ierr = MPIU_File_read_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr); 480 } 481 ierr = MPI_Type_get_extent(view,&ul,&ub);CHKERRQ(ierr); 482 ierr = PetscViewerBinaryAddMPIIOOffset(viewer,ub);CHKERRQ(ierr); 483 ierr = VecRestoreArrayRead(xin,&array);CHKERRQ(ierr); 484 ierr = MPI_Type_free(&view);CHKERRQ(ierr); 485 PetscFunctionReturn(0); 486 } 487 #endif 488 489 EXTERN_C_BEGIN 490 #undef __FUNCT__ 491 #define __FUNCT__ "VecView_MPI_DA" 492 PetscErrorCode VecView_MPI_DA(Vec xin,PetscViewer viewer) 493 { 494 DM da; 495 PetscErrorCode ierr; 496 PetscInt dim; 497 Vec natural; 498 PetscBool isdraw,isvtk; 499 #if defined(PETSC_HAVE_HDF5) 500 PetscBool ishdf5; 501 #endif 502 const char *prefix,*name; 503 504 PetscFunctionBegin; 505 ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 506 if (!da) SETERRQ(((PetscObject)xin)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 507 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 508 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk);CHKERRQ(ierr); 509 #if defined(PETSC_HAVE_HDF5) 510 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr); 511 #endif 512 if (isdraw) { 513 ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 514 if (dim == 1) { 515 ierr = VecView_MPI_Draw_DA1d(xin,viewer);CHKERRQ(ierr); 516 } else if (dim == 2) { 517 ierr = VecView_MPI_Draw_DA2d(xin,viewer);CHKERRQ(ierr); 518 } else { 519 SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Cannot graphically view vector associated with this dimensional DMDA %D",dim); 520 } 521 } else if (isvtk) { /* Duplicate the Vec and hold a reference to the DM */ 522 Vec Y; 523 ierr = PetscObjectReference((PetscObject)da);CHKERRQ(ierr); 524 ierr = VecDuplicate(xin,&Y);CHKERRQ(ierr); 525 if (((PetscObject)xin)->name) { 526 /* If xin was named, copy the name over to Y. The duplicate names are safe because nobody else will ever see Y. */ 527 ierr = PetscObjectSetName((PetscObject)Y,((PetscObject)xin)->name);CHKERRQ(ierr); 528 } 529 ierr = VecCopy(xin,Y);CHKERRQ(ierr); 530 ierr = PetscViewerVTKAddField(viewer,(PetscObject)da,DMDAVTKWriteAll,PETSC_VTK_POINT_FIELD,(PetscObject)Y);CHKERRQ(ierr); 531 #if defined(PETSC_HAVE_HDF5) 532 } else if (ishdf5) { 533 ierr = VecView_MPI_HDF5_DA(xin,viewer);CHKERRQ(ierr); 534 #endif 535 } else { 536 #if defined(PETSC_HAVE_MPIIO) 537 PetscBool isbinary,isMPIIO; 538 539 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 540 if (isbinary) { 541 ierr = PetscViewerBinaryGetMPIIO(viewer,&isMPIIO);CHKERRQ(ierr); 542 if (isMPIIO) { 543 ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_TRUE);CHKERRQ(ierr); 544 PetscFunctionReturn(0); 545 } 546 } 547 #endif 548 549 /* call viewer on natural ordering */ 550 ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr); 551 ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr); 552 ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr); 553 ierr = DMDAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr); 554 ierr = DMDAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr); 555 ierr = PetscObjectGetName((PetscObject)xin,&name);CHKERRQ(ierr); 556 ierr = PetscObjectSetName((PetscObject)natural,name);CHKERRQ(ierr); 557 ierr = VecView(natural,viewer);CHKERRQ(ierr); 558 ierr = VecDestroy(&natural);CHKERRQ(ierr); 559 } 560 PetscFunctionReturn(0); 561 } 562 EXTERN_C_END 563 564 #if defined(PETSC_HAVE_HDF5) 565 #undef __FUNCT__ 566 #define __FUNCT__ "VecLoad_HDF5_DA" 567 PetscErrorCode VecLoad_HDF5_DA(Vec xin, PetscViewer viewer) 568 { 569 DM da; 570 PetscErrorCode ierr; 571 hsize_t dim; 572 hsize_t count[5]; 573 hsize_t offset[5]; 574 PetscInt cnt = 0; 575 PetscScalar *x; 576 const char *vecname; 577 hid_t filespace; /* file dataspace identifier */ 578 hid_t plist_id; /* property list identifier */ 579 hid_t dset_id; /* dataset identifier */ 580 hid_t memspace; /* memory dataspace identifier */ 581 hid_t file_id; 582 herr_t status; 583 DM_DA *dd; 584 585 PetscFunctionBegin; 586 ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 587 ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 588 dd = (DM_DA*)da->data; 589 590 /* Create the dataspace for the dataset */ 591 dim = PetscHDF5IntCast(dd->dim + ((dd->w == 1) ? 0 : 1)); 592 #if defined(PETSC_USE_COMPLEX) 593 dim++; 594 #endif 595 596 /* Create the dataset with default properties and close filespace */ 597 ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr); 598 #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 599 dset_id = H5Dopen2(file_id, vecname, H5P_DEFAULT); 600 #else 601 dset_id = H5Dopen(file_id, vecname); 602 #endif 603 if (dset_id == -1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Dopen2() with Vec named %s",vecname); 604 filespace = H5Dget_space(dset_id); 605 606 /* Each process defines a dataset and reads it from the hyperslab in the file */ 607 cnt = 0; 608 if (dd->dim == 3) offset[cnt++] = PetscHDF5IntCast(dd->zs); 609 if (dd->dim > 1) offset[cnt++] = PetscHDF5IntCast(dd->ys); 610 offset[cnt++] = PetscHDF5IntCast(dd->xs/dd->w); 611 if (dd->w > 1) offset[cnt++] = 0; 612 #if defined(PETSC_USE_COMPLEX) 613 offset[cnt++] = 0; 614 #endif 615 cnt = 0; 616 if (dd->dim == 3) count[cnt++] = PetscHDF5IntCast(dd->ze - dd->zs); 617 if (dd->dim > 1) count[cnt++] = PetscHDF5IntCast(dd->ye - dd->ys); 618 count[cnt++] = PetscHDF5IntCast((dd->xe - dd->xs)/dd->w); 619 if (dd->w > 1) count[cnt++] = PetscHDF5IntCast(dd->w); 620 #if defined(PETSC_USE_COMPLEX) 621 count[cnt++] = 2; 622 #endif 623 memspace = H5Screate_simple(dim, count, NULL); 624 if (memspace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Screate_simple()"); 625 626 status = H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);CHKERRQ(status); 627 628 /* Create property list for collective dataset write */ 629 plist_id = H5Pcreate(H5P_DATASET_XFER); 630 if (plist_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Pcreate()"); 631 #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 632 status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);CHKERRQ(status); 633 #endif 634 /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 635 636 ierr = VecGetArray(xin, &x);CHKERRQ(ierr); 637 status = H5Dread(dset_id, H5T_NATIVE_DOUBLE, memspace, filespace, plist_id, x);CHKERRQ(status); 638 ierr = VecRestoreArray(xin, &x);CHKERRQ(ierr); 639 640 /* Close/release resources */ 641 status = H5Pclose(plist_id);CHKERRQ(status); 642 status = H5Sclose(filespace);CHKERRQ(status); 643 status = H5Sclose(memspace);CHKERRQ(status); 644 status = H5Dclose(dset_id);CHKERRQ(status); 645 PetscFunctionReturn(0); 646 } 647 #endif 648 649 #undef __FUNCT__ 650 #define __FUNCT__ "VecLoad_Binary_DA" 651 PetscErrorCode VecLoad_Binary_DA(Vec xin, PetscViewer viewer) 652 { 653 DM da; 654 PetscErrorCode ierr; 655 Vec natural; 656 const char *prefix; 657 PetscInt bs; 658 PetscBool flag; 659 DM_DA *dd; 660 #if defined(PETSC_HAVE_MPIIO) 661 PetscBool isMPIIO; 662 #endif 663 664 PetscFunctionBegin; 665 ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 666 dd = (DM_DA*)da->data; 667 #if defined(PETSC_HAVE_MPIIO) 668 ierr = PetscViewerBinaryGetMPIIO(viewer,&isMPIIO);CHKERRQ(ierr); 669 if (isMPIIO) { 670 ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_FALSE);CHKERRQ(ierr); 671 PetscFunctionReturn(0); 672 } 673 #endif 674 675 ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr); 676 ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr); 677 ierr = PetscObjectSetName((PetscObject)natural,((PetscObject)xin)->name);CHKERRQ(ierr); 678 ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr); 679 ierr = VecLoad_Binary(natural,viewer);CHKERRQ(ierr); 680 ierr = DMDANaturalToGlobalBegin(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr); 681 ierr = DMDANaturalToGlobalEnd(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr); 682 ierr = VecDestroy(&natural);CHKERRQ(ierr); 683 ierr = PetscInfo(xin,"Loading vector from natural ordering into DMDA\n");CHKERRQ(ierr); 684 ierr = PetscOptionsGetInt(((PetscObject)xin)->prefix,"-vecload_block_size",&bs,&flag);CHKERRQ(ierr); 685 if (flag && bs != dd->w) { 686 ierr = PetscInfo2(xin,"Block size in file %D not equal to DMDA's dof %D\n",bs,dd->w);CHKERRQ(ierr); 687 } 688 PetscFunctionReturn(0); 689 } 690 691 EXTERN_C_BEGIN 692 #undef __FUNCT__ 693 #define __FUNCT__ "VecLoad_Default_DA" 694 PetscErrorCode VecLoad_Default_DA(Vec xin, PetscViewer viewer) 695 { 696 PetscErrorCode ierr; 697 DM da; 698 PetscBool isbinary; 699 #if defined(PETSC_HAVE_HDF5) 700 PetscBool ishdf5; 701 #endif 702 703 PetscFunctionBegin; 704 ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 705 if (!da) SETERRQ(((PetscObject)xin)->comm,PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 706 707 #if defined(PETSC_HAVE_HDF5) 708 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr); 709 #endif 710 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 711 712 if (isbinary) { 713 ierr = VecLoad_Binary_DA(xin,viewer);CHKERRQ(ierr); 714 #if defined(PETSC_HAVE_HDF5) 715 } else if (ishdf5) { 716 ierr = VecLoad_HDF5_DA(xin,viewer);CHKERRQ(ierr); 717 #endif 718 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for vector loading", ((PetscObject)viewer)->type_name); 719 PetscFunctionReturn(0); 720 } 721 EXTERN_C_END 722