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