1 2 /* 3 Plots vectors obtained with DMDACreate2d() 4 */ 5 6 #include <petsc-private/dmdaimpl.h> /*I "petscdmda.h" I*/ 7 #include <petsc-private/vecimpl.h> 8 #include <petscdraw.h> 9 #include <petscviewerhdf5.h> 10 11 /* 12 The data that is passed into the graphics callback 13 */ 14 typedef struct { 15 PetscInt m,n,step,k; 16 PetscReal min,max,scale; 17 const PetscScalar *xy,*v; 18 PetscBool showgrid; 19 const char *name0,*name1; 20 } ZoomCtx; 21 22 /* 23 This does the drawing for one particular field 24 in one particular set of coordinates. It is a callback 25 called from PetscDrawZoom() 26 */ 27 #undef __FUNCT__ 28 #define __FUNCT__ "VecView_MPI_Draw_DA2d_Zoom" 29 PetscErrorCode VecView_MPI_Draw_DA2d_Zoom(PetscDraw draw,void *ctx) 30 { 31 ZoomCtx *zctx = (ZoomCtx*)ctx; 32 PetscErrorCode ierr; 33 PetscInt m,n,i,j,k,step,id,c1,c2,c3,c4; 34 PetscReal s,min,max,x1,x2,x3,x4,y_1,y2,y3,y4,xmin = PETSC_MAX_REAL,xmax = PETSC_MIN_REAL,ymin = PETSC_MAX_REAL,ymax = PETSC_MIN_REAL; 35 PetscReal xminf,xmaxf,yminf,ymaxf,w; 36 const PetscScalar *v,*xy; 37 char value[16]; 38 size_t len; 39 40 PetscFunctionBegin; 41 m = zctx->m; 42 n = zctx->n; 43 step = zctx->step; 44 k = zctx->k; 45 v = zctx->v; 46 xy = zctx->xy; 47 s = zctx->scale; 48 min = zctx->min; 49 max = zctx->max; 50 51 /* PetscDraw the contour plot patch */ 52 for (j=0; j<n-1; j++) { 53 for (i=0; i<m-1; i++) { 54 id = i+j*m; 55 x1 = PetscRealPart(xy[2*id]); 56 y_1 = PetscRealPart(xy[2*id+1]); 57 c1 = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscClipInterval(PetscRealPart(v[k+step*id]),min,max)-min)); 58 xmin = PetscMin(xmin,x1); 59 ymin = PetscMin(ymin,y_1); 60 xmax = PetscMax(xmax,x1); 61 ymax = PetscMax(ymax,y_1); 62 63 id = i+j*m+1; 64 x2 = PetscRealPart(xy[2*id]); 65 y2 = PetscRealPart(xy[2*id+1]); 66 c2 = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscClipInterval(PetscRealPart(v[k+step*id]),min,max)-min)); 67 xmin = PetscMin(xmin,x2); 68 ymin = PetscMin(ymin,y2); 69 xmax = PetscMax(xmax,x2); 70 ymax = PetscMax(ymax,y2); 71 72 id = i+j*m+1+m; 73 x3 = PetscRealPart(xy[2*id]); 74 y3 = PetscRealPart(xy[2*id+1]); 75 c3 = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscClipInterval(PetscRealPart(v[k+step*id]),min,max)-min)); 76 xmin = PetscMin(xmin,x3); 77 ymin = PetscMin(ymin,y3); 78 xmax = PetscMax(xmax,x3); 79 ymax = PetscMax(ymax,y3); 80 81 id = i+j*m+m; 82 x4 = PetscRealPart(xy[2*id]); 83 y4 = PetscRealPart(xy[2*id+1]); 84 c4 = (int)(PETSC_DRAW_BASIC_COLORS+s*(PetscClipInterval(PetscRealPart(v[k+step*id]),min,max)-min)); 85 xmin = PetscMin(xmin,x4); 86 ymin = PetscMin(ymin,y4); 87 xmax = PetscMax(xmax,x4); 88 ymax = PetscMax(ymax,y4); 89 90 ierr = PetscDrawTriangle(draw,x1,y_1,x2,y2,x3,y3,c1,c2,c3);CHKERRQ(ierr); 91 ierr = PetscDrawTriangle(draw,x1,y_1,x3,y3,x4,y4,c1,c3,c4);CHKERRQ(ierr); 92 if (zctx->showgrid) { 93 ierr = PetscDrawLine(draw,x1,y_1,x2,y2,PETSC_DRAW_BLACK);CHKERRQ(ierr); 94 ierr = PetscDrawLine(draw,x2,y2,x3,y3,PETSC_DRAW_BLACK);CHKERRQ(ierr); 95 ierr = PetscDrawLine(draw,x3,y3,x4,y4,PETSC_DRAW_BLACK);CHKERRQ(ierr); 96 ierr = PetscDrawLine(draw,x4,y4,x1,y_1,PETSC_DRAW_BLACK);CHKERRQ(ierr); 97 } 98 } 99 } 100 if (zctx->name0) { 101 PetscReal xl,yl,xr,yr,x,y; 102 ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 103 x = xl + .3*(xr - xl); 104 xl = xl + .01*(xr - xl); 105 y = yr - .3*(yr - yl); 106 yl = yl + .01*(yr - yl); 107 ierr = PetscDrawString(draw,x,yl,PETSC_DRAW_BLACK,zctx->name0);CHKERRQ(ierr); 108 ierr = PetscDrawStringVertical(draw,xl,y,PETSC_DRAW_BLACK,zctx->name1);CHKERRQ(ierr); 109 } 110 /* 111 Ideally we would use the PetscDrawAxis object to manage displaying the coordinate limits 112 but that may require some refactoring. 113 */ 114 ierr = MPI_Allreduce(&xmin,&xminf,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)draw));CHKERRQ(ierr); 115 ierr = MPI_Allreduce(&xmax,&xmaxf,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)draw));CHKERRQ(ierr); 116 ierr = MPI_Allreduce(&ymin,&yminf,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)draw));CHKERRQ(ierr); 117 ierr = MPI_Allreduce(&ymax,&ymaxf,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)draw));CHKERRQ(ierr); 118 ierr = PetscSNPrintf(value,16,"%f",xminf);CHKERRQ(ierr); 119 ierr = PetscDrawString(draw,xminf,yminf - .05*(ymaxf - yminf),PETSC_DRAW_BLACK,value);CHKERRQ(ierr); 120 ierr = PetscSNPrintf(value,16,"%f",xmaxf);CHKERRQ(ierr); 121 ierr = PetscStrlen(value,&len);CHKERRQ(ierr); 122 ierr = PetscDrawStringGetSize(draw,&w,NULL);CHKERRQ(ierr); 123 ierr = PetscDrawString(draw,xmaxf - len*w,yminf - .05*(ymaxf - yminf),PETSC_DRAW_BLACK,value);CHKERRQ(ierr); 124 ierr = PetscSNPrintf(value,16,"%f",yminf);CHKERRQ(ierr); 125 ierr = PetscDrawString(draw,xminf - .05*(xmaxf - xminf),yminf,PETSC_DRAW_BLACK,value);CHKERRQ(ierr); 126 ierr = PetscSNPrintf(value,16,"%f",ymaxf);CHKERRQ(ierr); 127 ierr = PetscDrawString(draw,xminf - .05*(xmaxf - xminf),ymaxf,PETSC_DRAW_BLACK,value);CHKERRQ(ierr); 128 PetscFunctionReturn(0); 129 } 130 131 #undef __FUNCT__ 132 #define __FUNCT__ "VecView_MPI_Draw_DA2d" 133 PetscErrorCode VecView_MPI_Draw_DA2d(Vec xin,PetscViewer viewer) 134 { 135 DM da,dac,dag; 136 PetscErrorCode ierr; 137 PetscMPIInt rank; 138 PetscInt N,s,M,w,ncoors = 4; 139 const PetscInt *lx,*ly; 140 PetscReal coors[4],ymin,ymax,xmin,xmax; 141 PetscDraw draw,popup; 142 PetscBool isnull,useports = PETSC_FALSE; 143 MPI_Comm comm; 144 Vec xlocal,xcoor,xcoorl; 145 DMBoundaryType bx,by; 146 DMDAStencilType st; 147 ZoomCtx zctx; 148 PetscDrawViewPorts *ports = NULL; 149 PetscViewerFormat format; 150 PetscInt *displayfields; 151 PetscInt ndisplayfields,i,nbounds; 152 const PetscReal *bounds; 153 154 PetscFunctionBegin; 155 zctx.showgrid = PETSC_FALSE; 156 157 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 158 ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); 159 ierr = PetscViewerDrawGetBounds(viewer,&nbounds,&bounds);CHKERRQ(ierr); 160 161 ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 162 if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 163 164 ierr = PetscObjectGetComm((PetscObject)xin,&comm);CHKERRQ(ierr); 165 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 166 167 ierr = DMDAGetInfo(da,0,&M,&N,0,&zctx.m,&zctx.n,0,&w,&s,&bx,&by,0,&st);CHKERRQ(ierr); 168 ierr = DMDAGetOwnershipRanges(da,&lx,&ly,NULL);CHKERRQ(ierr); 169 170 /* 171 Obtain a sequential vector that is going to contain the local values plus ONE layer of 172 ghosted values to draw the graphics from. We also need its corresponding DMDA (dac) that will 173 update the local values pluse ONE layer of ghost values. 174 */ 175 ierr = PetscObjectQuery((PetscObject)da,"GraphicsGhosted",(PetscObject*)&xlocal);CHKERRQ(ierr); 176 if (!xlocal) { 177 if (bx != DM_BOUNDARY_NONE || by != DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) { 178 /* 179 if original da is not of stencil width one, or periodic or not a box stencil then 180 create a special DMDA to handle one level of ghost points for graphics 181 */ 182 ierr = DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,zctx.m,zctx.n,w,1,lx,ly,&dac);CHKERRQ(ierr); 183 ierr = PetscInfo(da,"Creating auxilary DMDA for managing graphics ghost points\n");CHKERRQ(ierr); 184 } else { 185 /* otherwise we can use the da we already have */ 186 dac = da; 187 } 188 /* create local vector for holding ghosted values used in graphics */ 189 ierr = DMCreateLocalVector(dac,&xlocal);CHKERRQ(ierr); 190 if (dac != da) { 191 /* don't keep any public reference of this DMDA, is is only available through xlocal */ 192 ierr = PetscObjectDereference((PetscObject)dac);CHKERRQ(ierr); 193 } else { 194 /* remove association between xlocal and da, because below we compose in the opposite 195 direction and if we left this connect we'd get a loop, so the objects could 196 never be destroyed */ 197 ierr = PetscObjectRemoveReference((PetscObject)xlocal,"__PETSc_dm");CHKERRQ(ierr); 198 } 199 ierr = PetscObjectCompose((PetscObject)da,"GraphicsGhosted",(PetscObject)xlocal);CHKERRQ(ierr); 200 ierr = PetscObjectDereference((PetscObject)xlocal);CHKERRQ(ierr); 201 } else { 202 if (bx != DM_BOUNDARY_NONE || by != DM_BOUNDARY_NONE || s != 1 || st != DMDA_STENCIL_BOX) { 203 ierr = VecGetDM(xlocal, &dac);CHKERRQ(ierr); 204 } else { 205 dac = da; 206 } 207 } 208 209 /* 210 Get local (ghosted) values of vector 211 */ 212 ierr = DMGlobalToLocalBegin(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr); 213 ierr = DMGlobalToLocalEnd(dac,xin,INSERT_VALUES,xlocal);CHKERRQ(ierr); 214 ierr = VecGetArrayRead(xlocal,&zctx.v);CHKERRQ(ierr); 215 216 /* get coordinates of nodes */ 217 ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr); 218 if (!xcoor) { 219 ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,0.0);CHKERRQ(ierr); 220 ierr = DMGetCoordinates(da,&xcoor);CHKERRQ(ierr); 221 } 222 223 /* 224 Determine the min and max coordinates in plot 225 */ 226 ierr = VecStrideMin(xcoor,0,NULL,&xmin);CHKERRQ(ierr); 227 ierr = VecStrideMax(xcoor,0,NULL,&xmax);CHKERRQ(ierr); 228 ierr = VecStrideMin(xcoor,1,NULL,&ymin);CHKERRQ(ierr); 229 ierr = VecStrideMax(xcoor,1,NULL,&ymax);CHKERRQ(ierr); 230 coors[0] = xmin - .05*(xmax- xmin); coors[2] = xmax + .05*(xmax - xmin); 231 coors[1] = ymin - .05*(ymax- ymin); coors[3] = ymax + .05*(ymax - ymin); 232 ierr = PetscInfo4(da,"Preparing DMDA 2d contour plot coordinates %g %g %g %g\n",(double)coors[0],(double)coors[1],(double)coors[2],(double)coors[3]);CHKERRQ(ierr); 233 234 ierr = PetscOptionsGetRealArray(NULL,"-draw_coordinates",coors,&ncoors,NULL);CHKERRQ(ierr); 235 236 /* 237 get local ghosted version of coordinates 238 */ 239 ierr = PetscObjectQuery((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject*)&xcoorl);CHKERRQ(ierr); 240 if (!xcoorl) { 241 /* create DMDA to get local version of graphics */ 242 ierr = DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,zctx.m,zctx.n,2,1,lx,ly,&dag);CHKERRQ(ierr); 243 ierr = PetscInfo(dag,"Creating auxilary DMDA for managing graphics coordinates ghost points\n");CHKERRQ(ierr); 244 ierr = DMCreateLocalVector(dag,&xcoorl);CHKERRQ(ierr); 245 ierr = PetscObjectCompose((PetscObject)da,"GraphicsCoordinateGhosted",(PetscObject)xcoorl);CHKERRQ(ierr); 246 ierr = PetscObjectDereference((PetscObject)dag);CHKERRQ(ierr); 247 ierr = PetscObjectDereference((PetscObject)xcoorl);CHKERRQ(ierr); 248 } else { 249 ierr = VecGetDM(xcoorl,&dag);CHKERRQ(ierr); 250 } 251 ierr = DMGlobalToLocalBegin(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr); 252 ierr = DMGlobalToLocalEnd(dag,xcoor,INSERT_VALUES,xcoorl);CHKERRQ(ierr); 253 ierr = VecGetArrayRead(xcoorl,&zctx.xy);CHKERRQ(ierr); 254 255 /* 256 Get information about size of area each processor must do graphics for 257 */ 258 ierr = DMDAGetInfo(dac,0,&M,&N,0,0,0,0,&zctx.step,0,&bx,&by,0,0);CHKERRQ(ierr); 259 ierr = DMDAGetGhostCorners(dac,0,0,0,&zctx.m,&zctx.n,0);CHKERRQ(ierr); 260 261 ierr = PetscOptionsGetBool(NULL,"-draw_contour_grid",&zctx.showgrid,NULL);CHKERRQ(ierr); 262 263 ierr = DMDASelectFields(da,&ndisplayfields,&displayfields);CHKERRQ(ierr); 264 265 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 266 ierr = PetscOptionsGetBool(NULL,"-draw_ports",&useports,NULL);CHKERRQ(ierr); 267 if (useports || format == PETSC_VIEWER_DRAW_PORTS) { 268 ierr = PetscDrawSynchronizedClear(draw);CHKERRQ(ierr); 269 ierr = PetscDrawViewPortsCreate(draw,ndisplayfields,&ports);CHKERRQ(ierr); 270 zctx.name0 = 0; 271 zctx.name1 = 0; 272 } else { 273 ierr = DMDAGetCoordinateName(da,0,&zctx.name0);CHKERRQ(ierr); 274 ierr = DMDAGetCoordinateName(da,1,&zctx.name1);CHKERRQ(ierr); 275 } 276 277 /* 278 Loop over each field; drawing each in a different window 279 */ 280 for (i=0; i<ndisplayfields; i++) { 281 zctx.k = displayfields[i]; 282 if (useports) { 283 ierr = PetscDrawViewPortsSet(ports,i);CHKERRQ(ierr); 284 } else { 285 ierr = PetscViewerDrawGetDraw(viewer,i,&draw);CHKERRQ(ierr); 286 } 287 288 /* 289 Determine the min and max color in plot 290 */ 291 ierr = VecStrideMin(xin,zctx.k,NULL,&zctx.min);CHKERRQ(ierr); 292 ierr = VecStrideMax(xin,zctx.k,NULL,&zctx.max);CHKERRQ(ierr); 293 if (zctx.k < nbounds) { 294 zctx.min = bounds[2*zctx.k]; 295 zctx.max = bounds[2*zctx.k+1]; 296 } 297 if (zctx.min == zctx.max) { 298 zctx.min -= 1.e-12; 299 zctx.max += 1.e-12; 300 } 301 302 if (!rank) { 303 const char *title; 304 305 ierr = DMDAGetFieldName(da,zctx.k,&title);CHKERRQ(ierr); 306 if (title) { 307 ierr = PetscDrawSetTitle(draw,title);CHKERRQ(ierr); 308 } 309 } 310 ierr = PetscDrawSetCoordinates(draw,coors[0],coors[1],coors[2],coors[3]);CHKERRQ(ierr); 311 ierr = PetscInfo2(da,"DMDA 2d contour plot min %g max %g\n",(double)zctx.min,(double)zctx.max);CHKERRQ(ierr); 312 313 ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 314 if (popup) {ierr = PetscDrawScalePopup(popup,zctx.min,zctx.max);CHKERRQ(ierr);} 315 316 zctx.scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/(zctx.max - zctx.min); 317 318 ierr = PetscDrawZoom(draw,VecView_MPI_Draw_DA2d_Zoom,&zctx);CHKERRQ(ierr); 319 } 320 ierr = PetscFree(displayfields);CHKERRQ(ierr); 321 ierr = PetscDrawViewPortsDestroy(ports);CHKERRQ(ierr); 322 323 ierr = VecRestoreArrayRead(xcoorl,&zctx.xy);CHKERRQ(ierr); 324 ierr = VecRestoreArrayRead(xlocal,&zctx.v);CHKERRQ(ierr); 325 PetscFunctionReturn(0); 326 } 327 328 #if defined(PETSC_HAVE_HDF5) 329 #undef __FUNCT__ 330 #define __FUNCT__ "VecGetHDF5ChunkSize" 331 static PetscErrorCode VecGetHDF5ChunkSize(DM_DA *da, Vec xin, PetscInt dimension, PetscInt timestep, hsize_t *chunkDims) 332 { 333 PetscMPIInt comm_size; 334 PetscErrorCode ierr; 335 hsize_t chunk_size, target_size, dim; 336 hsize_t vec_size = sizeof(PetscScalar)*da->M*da->N*da->P*da->w; 337 hsize_t avg_local_vec_size,KiB = 1024,MiB = KiB*KiB,GiB = MiB*KiB,min_size = MiB; 338 hsize_t max_chunks = 64*KiB; /* HDF5 internal limitation */ 339 hsize_t max_chunk_size = 4*GiB; /* HDF5 internal limitation */ 340 PetscInt zslices=da->p, yslices=da->n, xslices=da->m; 341 342 PetscFunctionBegin; 343 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)xin), &comm_size);CHKERRQ(ierr); 344 avg_local_vec_size = (hsize_t) ceil(vec_size*1.0/comm_size); /* we will attempt to use this as the chunk size */ 345 346 target_size = (hsize_t) ceil(PetscMin(vec_size,PetscMin(max_chunk_size,PetscMax(avg_local_vec_size,PetscMax(ceil(vec_size*1.0/max_chunks),min_size))))); 347 /* following line uses sizeof(PetscReal) instead of sizeof(PetscScalar) because the last dimension of chunkDims[] captures the 2* when complex numbers are being used */ 348 chunk_size = (hsize_t) PetscMax(1,chunkDims[0])*PetscMax(1,chunkDims[1])*PetscMax(1,chunkDims[2])*PetscMax(1,chunkDims[3])*PetscMax(1,chunkDims[4])*PetscMax(1,chunkDims[5])*sizeof(PetscReal); 349 350 /* 351 if size/rank > max_chunk_size, we need radical measures: even going down to 352 avg_local_vec_size is not enough, so we simply use chunk size of 4 GiB no matter 353 what, composed in the most efficient way possible. 354 N.B. this minimises the number of chunks, which may or may not be the optimal 355 solution. In a BG, for example, the optimal solution is probably to make # chunks = # 356 IO nodes involved, but this author has no access to a BG to figure out how to 357 reliably find the right number. And even then it may or may not be enough. 358 */ 359 if (avg_local_vec_size > max_chunk_size) { 360 /* check if we can just split local z-axis: is that enough? */ 361 zslices = (PetscInt)ceil(vec_size*1.0/(da->p*max_chunk_size))*zslices; 362 if (zslices > da->P) { 363 /* lattice is too large in xy-directions, splitting z only is not enough */ 364 zslices = da->P; 365 yslices= (PetscInt)ceil(vec_size*1.0/(zslices*da->n*max_chunk_size))*yslices; 366 if (yslices > da->N) { 367 /* lattice is too large in x-direction, splitting along z, y is not enough */ 368 yslices = da->N; 369 xslices= (PetscInt)ceil(vec_size*1.0/(zslices*yslices*da->m*max_chunk_size))*xslices; 370 } 371 } 372 dim = 0; 373 if (timestep >= 0) { 374 ++dim; 375 } 376 /* prefer to split z-axis, even down to planar slices */ 377 if (dimension == 3) { 378 chunkDims[dim++] = (hsize_t) da->P/zslices; 379 chunkDims[dim++] = (hsize_t) da->N/yslices; 380 chunkDims[dim++] = (hsize_t) da->M/xslices; 381 } else { 382 /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */ 383 chunkDims[dim++] = (hsize_t) da->N/yslices; 384 chunkDims[dim++] = (hsize_t) da->M/xslices; 385 } 386 chunk_size = (hsize_t) PetscMax(1,chunkDims[0])*PetscMax(1,chunkDims[1])*PetscMax(1,chunkDims[2])*PetscMax(1,chunkDims[3])*PetscMax(1,chunkDims[4])*PetscMax(1,chunkDims[5])*sizeof(double); 387 } else { 388 if (target_size < chunk_size) { 389 /* only change the defaults if target_size < chunk_size */ 390 dim = 0; 391 if (timestep >= 0) { 392 ++dim; 393 } 394 /* prefer to split z-axis, even down to planar slices */ 395 if (dimension == 3) { 396 /* try splitting the z-axis to core-size bits, i.e. divide chunk size by # comm_size in z-direction */ 397 if (target_size >= chunk_size/da->p) { 398 /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */ 399 chunkDims[dim] = (hsize_t) ceil(da->P*1.0/da->p); 400 } else { 401 /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be 402 radical and let everyone write all they've got */ 403 chunkDims[dim++] = (hsize_t) ceil(da->P*1.0/da->p); 404 chunkDims[dim++] = (hsize_t) ceil(da->N*1.0/da->n); 405 chunkDims[dim++] = (hsize_t) ceil(da->M*1.0/da->m); 406 } 407 } else { 408 /* This is a 2D world exceeding 4GiB in size; yes, I've seen them, even used myself */ 409 if (target_size >= chunk_size/da->n) { 410 /* just make chunks the size of <local_z>x<whole_world_y>x<whole_world_x>x<dof> */ 411 chunkDims[dim] = (hsize_t) ceil(da->N*1.0/da->n); 412 } else { 413 /* oops, just splitting the z-axis is NOT ENOUGH, need to split more; let's be 414 radical and let everyone write all they've got */ 415 chunkDims[dim++] = (hsize_t) ceil(da->N*1.0/da->n); 416 chunkDims[dim++] = (hsize_t) ceil(da->M*1.0/da->m); 417 } 418 419 } 420 chunk_size = (hsize_t) PetscMax(1,chunkDims[0])*PetscMax(1,chunkDims[1])*PetscMax(1,chunkDims[2])*PetscMax(1,chunkDims[3])*PetscMax(1,chunkDims[4])*PetscMax(1,chunkDims[5])*sizeof(double); 421 } else { 422 /* precomputed chunks are fine, we don't need to do anything */ 423 } 424 } 425 PetscFunctionReturn(0); 426 } 427 #endif 428 429 #if defined(PETSC_HAVE_HDF5) 430 #undef __FUNCT__ 431 #define __FUNCT__ "VecView_MPI_HDF5_DA" 432 PetscErrorCode VecView_MPI_HDF5_DA(Vec xin,PetscViewer viewer) 433 { 434 DM dm; 435 DM_DA *da; 436 hid_t filespace; /* file dataspace identifier */ 437 hid_t chunkspace; /* chunk dataset property identifier */ 438 hid_t plist_id; /* property list identifier */ 439 hid_t dset_id; /* dataset identifier */ 440 hid_t memspace; /* memory dataspace identifier */ 441 hid_t file_id; 442 hid_t group; 443 hid_t scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */ 444 hsize_t dim; 445 hsize_t maxDims[6]={0}, dims[6]={0}, chunkDims[6]={0}, count[6]={0}, offset[6]={0}; /* we depend on these being sane later on */ 446 PetscInt timestep, dimension; 447 const PetscScalar *x; 448 const char *vecname; 449 PetscErrorCode ierr; 450 451 PetscFunctionBegin; 452 ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr); 453 ierr = PetscViewerHDF5GetTimestep(viewer, ×tep);CHKERRQ(ierr); 454 455 ierr = VecGetDM(xin,&dm);CHKERRQ(ierr); 456 if (!dm) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 457 da = (DM_DA*)dm->data; 458 ierr = DMGetDimension(dm, &dimension);CHKERRQ(ierr); 459 460 /* Create the dataspace for the dataset. 461 * 462 * dims - holds the current dimensions of the dataset 463 * 464 * maxDims - holds the maximum dimensions of the dataset (unlimited 465 * for the number of time steps with the current dimensions for the 466 * other dimensions; so only additional time steps can be added). 467 * 468 * chunkDims - holds the size of a single time step (required to 469 * permit extending dataset). 470 */ 471 dim = 0; 472 if (timestep >= 0) { 473 dims[dim] = timestep+1; 474 maxDims[dim] = H5S_UNLIMITED; 475 chunkDims[dim] = 1; 476 ++dim; 477 } 478 if (dimension == 3) { 479 ierr = PetscHDF5IntCast(da->P,dims+dim);CHKERRQ(ierr); 480 maxDims[dim] = dims[dim]; 481 chunkDims[dim] = dims[dim]; 482 ++dim; 483 } 484 if (dimension > 1) { 485 ierr = PetscHDF5IntCast(da->N,dims+dim);CHKERRQ(ierr); 486 maxDims[dim] = dims[dim]; 487 chunkDims[dim] = dims[dim]; 488 ++dim; 489 } 490 ierr = PetscHDF5IntCast(da->M,dims+dim);CHKERRQ(ierr); 491 maxDims[dim] = dims[dim]; 492 chunkDims[dim] = dims[dim]; 493 ++dim; 494 if (da->w > 1) { 495 ierr = PetscHDF5IntCast(da->w,dims+dim);CHKERRQ(ierr); 496 maxDims[dim] = dims[dim]; 497 chunkDims[dim] = dims[dim]; 498 ++dim; 499 } 500 #if defined(PETSC_USE_COMPLEX) 501 dims[dim] = 2; 502 maxDims[dim] = dims[dim]; 503 chunkDims[dim] = dims[dim]; 504 ++dim; 505 #endif 506 507 ierr = VecGetHDF5ChunkSize(da, xin, dimension, timestep, chunkDims); CHKERRQ(ierr); 508 509 PetscStackCallHDF5Return(filespace,H5Screate_simple,(dim, dims, maxDims)); 510 511 #if defined(PETSC_USE_REAL_SINGLE) 512 scalartype = H5T_NATIVE_FLOAT; 513 #elif defined(PETSC_USE_REAL___FLOAT128) 514 #error "HDF5 output with 128 bit floats not supported." 515 #else 516 scalartype = H5T_NATIVE_DOUBLE; 517 #endif 518 519 /* Create the dataset with default properties and close filespace */ 520 ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr); 521 if (!H5Lexists(group, vecname, H5P_DEFAULT)) { 522 /* Create chunk */ 523 PetscStackCallHDF5Return(chunkspace,H5Pcreate,(H5P_DATASET_CREATE)); 524 PetscStackCallHDF5(H5Pset_chunk,(chunkspace, dim, chunkDims)); 525 526 #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 527 PetscStackCallHDF5Return(dset_id,H5Dcreate2,(group, vecname, scalartype, filespace, H5P_DEFAULT, chunkspace, H5P_DEFAULT)); 528 #else 529 PetscStackCallHDF5Return(dset_id,H5Dcreate,(group, vecname, scalartype, filespace, H5P_DEFAULT)); 530 #endif 531 } else { 532 PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, vecname, H5P_DEFAULT)); 533 PetscStackCallHDF5(H5Dset_extent,(dset_id, dims)); 534 } 535 PetscStackCallHDF5(H5Sclose,(filespace)); 536 537 /* Each process defines a dataset and writes it to the hyperslab in the file */ 538 dim = 0; 539 if (timestep >= 0) { 540 offset[dim] = timestep; 541 ++dim; 542 } 543 if (dimension == 3) {ierr = PetscHDF5IntCast(da->zs,offset + dim++);CHKERRQ(ierr);} 544 if (dimension > 1) {ierr = PetscHDF5IntCast(da->ys,offset + dim++);CHKERRQ(ierr);} 545 ierr = PetscHDF5IntCast(da->xs/da->w,offset + dim++);CHKERRQ(ierr); 546 if (da->w > 1) offset[dim++] = 0; 547 #if defined(PETSC_USE_COMPLEX) 548 offset[dim++] = 0; 549 #endif 550 dim = 0; 551 if (timestep >= 0) { 552 count[dim] = 1; 553 ++dim; 554 } 555 if (dimension == 3) {ierr = PetscHDF5IntCast(da->ze - da->zs,count + dim++);CHKERRQ(ierr);} 556 if (dimension > 1) {ierr = PetscHDF5IntCast(da->ye - da->ys,count + dim++);CHKERRQ(ierr);} 557 ierr = PetscHDF5IntCast((da->xe - da->xs)/da->w,count + dim++);CHKERRQ(ierr); 558 if (da->w > 1) {ierr = PetscHDF5IntCast(da->w,count + dim++);CHKERRQ(ierr);} 559 #if defined(PETSC_USE_COMPLEX) 560 count[dim++] = 2; 561 #endif 562 PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL)); 563 PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id)); 564 PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL)); 565 566 /* Create property list for collective dataset write */ 567 PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_DATASET_XFER)); 568 #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 569 PetscStackCallHDF5(H5Pset_dxpl_mpio,(plist_id, H5FD_MPIO_COLLECTIVE)); 570 #endif 571 /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 572 573 ierr = VecGetArrayRead(xin, &x);CHKERRQ(ierr); 574 PetscStackCallHDF5(H5Dwrite,(dset_id, scalartype, memspace, filespace, plist_id, x)); 575 PetscStackCallHDF5(H5Fflush,(file_id, H5F_SCOPE_GLOBAL)); 576 ierr = VecRestoreArrayRead(xin, &x);CHKERRQ(ierr); 577 578 /* Close/release resources */ 579 if (group != file_id) { 580 PetscStackCallHDF5(H5Gclose,(group)); 581 } 582 PetscStackCallHDF5(H5Pclose,(plist_id)); 583 PetscStackCallHDF5(H5Sclose,(filespace)); 584 PetscStackCallHDF5(H5Sclose,(memspace)); 585 PetscStackCallHDF5(H5Dclose,(dset_id)); 586 ierr = PetscInfo1(xin,"Wrote Vec object with name %s\n",vecname);CHKERRQ(ierr); 587 PetscFunctionReturn(0); 588 } 589 #endif 590 591 extern PetscErrorCode VecView_MPI_Draw_DA1d(Vec,PetscViewer); 592 593 #if defined(PETSC_HAVE_MPIIO) 594 #undef __FUNCT__ 595 #define __FUNCT__ "DMDAArrayMPIIO" 596 static PetscErrorCode DMDAArrayMPIIO(DM da,PetscViewer viewer,Vec xin,PetscBool write) 597 { 598 PetscErrorCode ierr; 599 MPI_File mfdes; 600 PetscMPIInt gsizes[4],lsizes[4],lstarts[4],asiz,dof; 601 MPI_Datatype view; 602 const PetscScalar *array; 603 MPI_Offset off; 604 MPI_Aint ub,ul; 605 PetscInt type,rows,vecrows,tr[2]; 606 DM_DA *dd = (DM_DA*)da->data; 607 PetscBool skipheader; 608 PetscInt num = 2; 609 610 PetscFunctionBegin; 611 ierr = VecGetSize(xin,&vecrows);CHKERRQ(ierr); 612 ierr = PetscViewerBinaryGetSkipHeader(viewer,&skipheader);CHKERRQ(ierr); 613 if (!write) { 614 /* Read vector header. */ 615 if (!skipheader) { 616 ierr = PetscViewerBinaryRead(viewer,tr,&num,PETSC_INT);CHKERRQ(ierr); 617 type = tr[0]; 618 rows = tr[1]; 619 if (type != VEC_FILE_CLASSID) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Not vector next in file"); 620 if (rows != vecrows) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_SIZ,"Vector in file not same size as DMDA vector"); 621 } 622 } else { 623 tr[0] = VEC_FILE_CLASSID; 624 tr[1] = vecrows; 625 if (!skipheader) { 626 ierr = PetscViewerBinaryWrite(viewer,tr,2,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 627 } 628 } 629 630 ierr = PetscMPIIntCast(dd->w,&dof);CHKERRQ(ierr); 631 gsizes[0] = dof; 632 ierr = PetscMPIIntCast(dd->M,gsizes+1);CHKERRQ(ierr); 633 ierr = PetscMPIIntCast(dd->N,gsizes+2);CHKERRQ(ierr); 634 ierr = PetscMPIIntCast(dd->P,gsizes+3);CHKERRQ(ierr); 635 lsizes[0] = dof; 636 ierr = PetscMPIIntCast((dd->xe-dd->xs)/dof,lsizes+1);CHKERRQ(ierr); 637 ierr = PetscMPIIntCast(dd->ye-dd->ys,lsizes+2);CHKERRQ(ierr); 638 ierr = PetscMPIIntCast(dd->ze-dd->zs,lsizes+3);CHKERRQ(ierr); 639 lstarts[0] = 0; 640 ierr = PetscMPIIntCast(dd->xs/dof,lstarts+1);CHKERRQ(ierr); 641 ierr = PetscMPIIntCast(dd->ys,lstarts+2);CHKERRQ(ierr); 642 ierr = PetscMPIIntCast(dd->zs,lstarts+3);CHKERRQ(ierr); 643 ierr = MPI_Type_create_subarray(da->dim+1,gsizes,lsizes,lstarts,MPI_ORDER_FORTRAN,MPIU_SCALAR,&view);CHKERRQ(ierr); 644 ierr = MPI_Type_commit(&view);CHKERRQ(ierr); 645 646 ierr = PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);CHKERRQ(ierr); 647 ierr = PetscViewerBinaryGetMPIIOOffset(viewer,&off);CHKERRQ(ierr); 648 ierr = MPI_File_set_view(mfdes,off,MPIU_SCALAR,view,(char*)"native",MPI_INFO_NULL);CHKERRQ(ierr); 649 ierr = VecGetArrayRead(xin,&array);CHKERRQ(ierr); 650 asiz = lsizes[1]*(lsizes[2] > 0 ? lsizes[2] : 1)*(lsizes[3] > 0 ? lsizes[3] : 1)*dof; 651 if (write) { 652 ierr = MPIU_File_write_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr); 653 } else { 654 ierr = MPIU_File_read_all(mfdes,(PetscScalar*)array,asiz,MPIU_SCALAR,MPI_STATUS_IGNORE);CHKERRQ(ierr); 655 } 656 ierr = MPI_Type_get_extent(view,&ul,&ub);CHKERRQ(ierr); 657 ierr = PetscViewerBinaryAddMPIIOOffset(viewer,ub);CHKERRQ(ierr); 658 ierr = VecRestoreArrayRead(xin,&array);CHKERRQ(ierr); 659 ierr = MPI_Type_free(&view);CHKERRQ(ierr); 660 PetscFunctionReturn(0); 661 } 662 #endif 663 664 #undef __FUNCT__ 665 #define __FUNCT__ "VecView_MPI_DA" 666 PetscErrorCode VecView_MPI_DA(Vec xin,PetscViewer viewer) 667 { 668 DM da; 669 PetscErrorCode ierr; 670 PetscInt dim; 671 Vec natural; 672 PetscBool isdraw,isvtk; 673 #if defined(PETSC_HAVE_HDF5) 674 PetscBool ishdf5; 675 #endif 676 const char *prefix,*name; 677 PetscViewerFormat format; 678 679 PetscFunctionBegin; 680 ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 681 if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 682 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 683 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk);CHKERRQ(ierr); 684 #if defined(PETSC_HAVE_HDF5) 685 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr); 686 #endif 687 if (isdraw) { 688 ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 689 if (dim == 1) { 690 ierr = VecView_MPI_Draw_DA1d(xin,viewer);CHKERRQ(ierr); 691 } else if (dim == 2) { 692 ierr = VecView_MPI_Draw_DA2d(xin,viewer);CHKERRQ(ierr); 693 } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot graphically view vector associated with this dimensional DMDA %D",dim); 694 } else if (isvtk) { /* Duplicate the Vec and hold a reference to the DM */ 695 Vec Y; 696 ierr = PetscObjectReference((PetscObject)da);CHKERRQ(ierr); 697 ierr = VecDuplicate(xin,&Y);CHKERRQ(ierr); 698 if (((PetscObject)xin)->name) { 699 /* If xin was named, copy the name over to Y. The duplicate names are safe because nobody else will ever see Y. */ 700 ierr = PetscObjectSetName((PetscObject)Y,((PetscObject)xin)->name);CHKERRQ(ierr); 701 } 702 ierr = VecCopy(xin,Y);CHKERRQ(ierr); 703 ierr = PetscViewerVTKAddField(viewer,(PetscObject)da,DMDAVTKWriteAll,PETSC_VTK_POINT_FIELD,(PetscObject)Y);CHKERRQ(ierr); 704 #if defined(PETSC_HAVE_HDF5) 705 } else if (ishdf5) { 706 ierr = VecView_MPI_HDF5_DA(xin,viewer);CHKERRQ(ierr); 707 #endif 708 } else { 709 #if defined(PETSC_HAVE_MPIIO) 710 PetscBool isbinary,isMPIIO; 711 712 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 713 if (isbinary) { 714 ierr = PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO);CHKERRQ(ierr); 715 if (isMPIIO) { 716 ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_TRUE);CHKERRQ(ierr); 717 PetscFunctionReturn(0); 718 } 719 } 720 #endif 721 722 /* call viewer on natural ordering */ 723 ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr); 724 ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr); 725 ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr); 726 ierr = DMDAGlobalToNaturalBegin(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr); 727 ierr = DMDAGlobalToNaturalEnd(da,xin,INSERT_VALUES,natural);CHKERRQ(ierr); 728 ierr = PetscObjectGetName((PetscObject)xin,&name);CHKERRQ(ierr); 729 ierr = PetscObjectSetName((PetscObject)natural,name);CHKERRQ(ierr); 730 731 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 732 if (format == PETSC_VIEWER_BINARY_MATLAB) { 733 /* temporarily remove viewer format so it won't trigger in the VecView() */ 734 ierr = PetscViewerSetFormat(viewer,PETSC_VIEWER_DEFAULT);CHKERRQ(ierr); 735 } 736 737 ierr = VecView(natural,viewer);CHKERRQ(ierr); 738 739 if (format == PETSC_VIEWER_BINARY_MATLAB) { 740 MPI_Comm comm; 741 FILE *info; 742 const char *fieldname; 743 char fieldbuf[256]; 744 PetscInt dim,ni,nj,nk,pi,pj,pk,dof,n; 745 746 /* set the viewer format back into the viewer */ 747 ierr = PetscViewerSetFormat(viewer,format);CHKERRQ(ierr); 748 ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr); 749 ierr = PetscViewerBinaryGetInfoPointer(viewer,&info);CHKERRQ(ierr); 750 ierr = DMDAGetInfo(da,&dim,&ni,&nj,&nk,&pi,&pj,&pk,&dof,0,0,0,0,0);CHKERRQ(ierr); 751 ierr = PetscFPrintf(comm,info,"#--- begin code written by PetscViewerBinary for MATLAB format ---#\n");CHKERRQ(ierr); 752 ierr = PetscFPrintf(comm,info,"#$$ tmp = PetscBinaryRead(fd); \n");CHKERRQ(ierr); 753 if (dim == 1) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d);\n",dof,ni);CHKERRQ(ierr); } 754 if (dim == 2) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d);\n",dof,ni,nj);CHKERRQ(ierr); } 755 if (dim == 3) { ierr = PetscFPrintf(comm,info,"#$$ tmp = reshape(tmp,%d,%d,%d,%d);\n",dof,ni,nj,nk);CHKERRQ(ierr); } 756 757 for (n=0; n<dof; n++) { 758 ierr = DMDAGetFieldName(da,n,&fieldname);CHKERRQ(ierr); 759 if (!fieldname || !fieldname[0]) { 760 ierr = PetscSNPrintf(fieldbuf,sizeof fieldbuf,"field%D",n);CHKERRQ(ierr); 761 fieldname = fieldbuf; 762 } 763 if (dim == 1) { ierr = PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:))';\n",name,fieldname,n+1);CHKERRQ(ierr); } 764 if (dim == 2) { ierr = PetscFPrintf(comm,info,"#$$ Set.%s.%s = squeeze(tmp(%d,:,:))';\n",name,fieldname,n+1);CHKERRQ(ierr); } 765 if (dim == 3) { ierr = PetscFPrintf(comm,info,"#$$ Set.%s.%s = permute(squeeze(tmp(%d,:,:,:)),[2 1 3]);\n",name,fieldname,n+1);CHKERRQ(ierr);} 766 } 767 ierr = PetscFPrintf(comm,info,"#$$ clear tmp; \n");CHKERRQ(ierr); 768 ierr = PetscFPrintf(comm,info,"#--- end code written by PetscViewerBinary for MATLAB format ---#\n\n");CHKERRQ(ierr); 769 } 770 771 ierr = VecDestroy(&natural);CHKERRQ(ierr); 772 } 773 PetscFunctionReturn(0); 774 } 775 776 #if defined(PETSC_HAVE_HDF5) 777 #undef __FUNCT__ 778 #define __FUNCT__ "VecLoad_HDF5_DA" 779 PetscErrorCode VecLoad_HDF5_DA(Vec xin, PetscViewer viewer) 780 { 781 DM da; 782 PetscErrorCode ierr; 783 hsize_t dim; 784 hsize_t count[5]; 785 hsize_t offset[5]; 786 PetscInt cnt = 0, dimension; 787 PetscScalar *x; 788 const char *vecname; 789 hid_t filespace; /* file dataspace identifier */ 790 hid_t plist_id; /* property list identifier */ 791 hid_t dset_id; /* dataset identifier */ 792 hid_t memspace; /* memory dataspace identifier */ 793 hid_t file_id,group; 794 hid_t scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */ 795 DM_DA *dd; 796 797 PetscFunctionBegin; 798 #if defined(PETSC_USE_REAL_SINGLE) 799 scalartype = H5T_NATIVE_FLOAT; 800 #elif defined(PETSC_USE_REAL___FLOAT128) 801 #error "HDF5 output with 128 bit floats not supported." 802 #else 803 scalartype = H5T_NATIVE_DOUBLE; 804 #endif 805 806 ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr); 807 ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 808 dd = (DM_DA*)da->data; 809 ierr = DMGetDimension(da, &dimension);CHKERRQ(ierr); 810 811 /* Create the dataspace for the dataset */ 812 ierr = PetscHDF5IntCast(dimension + ((dd->w == 1) ? 0 : 1),&dim);CHKERRQ(ierr); 813 #if defined(PETSC_USE_COMPLEX) 814 dim++; 815 #endif 816 817 /* Create the dataset with default properties and close filespace */ 818 ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr); 819 #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 820 PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, vecname, H5P_DEFAULT)); 821 #else 822 PetscStackCallHDF5Return(dset_id,H5Dopen,(group, vecname)); 823 #endif 824 PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id)); 825 826 /* Each process defines a dataset and reads it from the hyperslab in the file */ 827 cnt = 0; 828 if (dimension == 3) {ierr = PetscHDF5IntCast(dd->zs,offset + cnt++);CHKERRQ(ierr);} 829 if (dimension > 1) {ierr = PetscHDF5IntCast(dd->ys,offset + cnt++);CHKERRQ(ierr);} 830 ierr = PetscHDF5IntCast(dd->xs/dd->w,offset + cnt++);CHKERRQ(ierr); 831 if (dd->w > 1) offset[cnt++] = 0; 832 #if defined(PETSC_USE_COMPLEX) 833 offset[cnt++] = 0; 834 #endif 835 cnt = 0; 836 if (dimension == 3) {ierr = PetscHDF5IntCast(dd->ze - dd->zs,count + cnt++);CHKERRQ(ierr);} 837 if (dimension > 1) {ierr = PetscHDF5IntCast(dd->ye - dd->ys,count + cnt++);CHKERRQ(ierr);} 838 ierr = PetscHDF5IntCast((dd->xe - dd->xs)/dd->w,count + cnt++);CHKERRQ(ierr); 839 if (dd->w > 1) {ierr = PetscHDF5IntCast(dd->w,count + cnt++);CHKERRQ(ierr);} 840 #if defined(PETSC_USE_COMPLEX) 841 count[cnt++] = 2; 842 #endif 843 PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL)); 844 PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL)); 845 846 /* Create property list for collective dataset write */ 847 PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_DATASET_XFER)); 848 #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 849 PetscStackCallHDF5(H5Pset_dxpl_mpio,(plist_id, H5FD_MPIO_COLLECTIVE)); 850 #endif 851 /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 852 853 ierr = VecGetArray(xin, &x);CHKERRQ(ierr); 854 PetscStackCallHDF5(H5Dread,(dset_id, scalartype, memspace, filespace, plist_id, x)); 855 ierr = VecRestoreArray(xin, &x);CHKERRQ(ierr); 856 857 /* Close/release resources */ 858 if (group != file_id) { 859 PetscStackCallHDF5(H5Gclose,(group)); 860 } 861 PetscStackCallHDF5(H5Pclose,(plist_id)); 862 PetscStackCallHDF5(H5Sclose,(filespace)); 863 PetscStackCallHDF5(H5Sclose,(memspace)); 864 PetscStackCallHDF5(H5Dclose,(dset_id)); 865 PetscFunctionReturn(0); 866 } 867 #endif 868 869 #undef __FUNCT__ 870 #define __FUNCT__ "VecLoad_Binary_DA" 871 PetscErrorCode VecLoad_Binary_DA(Vec xin, PetscViewer viewer) 872 { 873 DM da; 874 PetscErrorCode ierr; 875 Vec natural; 876 const char *prefix; 877 PetscInt bs; 878 PetscBool flag; 879 DM_DA *dd; 880 #if defined(PETSC_HAVE_MPIIO) 881 PetscBool isMPIIO; 882 #endif 883 884 PetscFunctionBegin; 885 ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 886 dd = (DM_DA*)da->data; 887 #if defined(PETSC_HAVE_MPIIO) 888 ierr = PetscViewerBinaryGetUseMPIIO(viewer,&isMPIIO);CHKERRQ(ierr); 889 if (isMPIIO) { 890 ierr = DMDAArrayMPIIO(da,viewer,xin,PETSC_FALSE);CHKERRQ(ierr); 891 PetscFunctionReturn(0); 892 } 893 #endif 894 895 ierr = PetscObjectGetOptionsPrefix((PetscObject)xin,&prefix);CHKERRQ(ierr); 896 ierr = DMDACreateNaturalVector(da,&natural);CHKERRQ(ierr); 897 ierr = PetscObjectSetName((PetscObject)natural,((PetscObject)xin)->name);CHKERRQ(ierr); 898 ierr = PetscObjectSetOptionsPrefix((PetscObject)natural,prefix);CHKERRQ(ierr); 899 ierr = VecLoad(natural,viewer);CHKERRQ(ierr); 900 ierr = DMDANaturalToGlobalBegin(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr); 901 ierr = DMDANaturalToGlobalEnd(da,natural,INSERT_VALUES,xin);CHKERRQ(ierr); 902 ierr = VecDestroy(&natural);CHKERRQ(ierr); 903 ierr = PetscInfo(xin,"Loading vector from natural ordering into DMDA\n");CHKERRQ(ierr); 904 ierr = PetscOptionsGetInt(((PetscObject)xin)->prefix,"-vecload_block_size",&bs,&flag);CHKERRQ(ierr); 905 if (flag && bs != dd->w) { 906 ierr = PetscInfo2(xin,"Block size in file %D not equal to DMDA's dof %D\n",bs,dd->w);CHKERRQ(ierr); 907 } 908 PetscFunctionReturn(0); 909 } 910 911 #undef __FUNCT__ 912 #define __FUNCT__ "VecLoad_Default_DA" 913 PetscErrorCode VecLoad_Default_DA(Vec xin, PetscViewer viewer) 914 { 915 PetscErrorCode ierr; 916 DM da; 917 PetscBool isbinary; 918 #if defined(PETSC_HAVE_HDF5) 919 PetscBool ishdf5; 920 #endif 921 922 PetscFunctionBegin; 923 ierr = VecGetDM(xin,&da);CHKERRQ(ierr); 924 if (!da) SETERRQ(PetscObjectComm((PetscObject)xin),PETSC_ERR_ARG_WRONG,"Vector not generated from a DMDA"); 925 926 #if defined(PETSC_HAVE_HDF5) 927 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr); 928 #endif 929 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 930 931 if (isbinary) { 932 ierr = VecLoad_Binary_DA(xin,viewer);CHKERRQ(ierr); 933 #if defined(PETSC_HAVE_HDF5) 934 } else if (ishdf5) { 935 ierr = VecLoad_HDF5_DA(xin,viewer);CHKERRQ(ierr); 936 #endif 937 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for vector loading", ((PetscObject)viewer)->type_name); 938 PetscFunctionReturn(0); 939 } 940