1 #include <petsc/private/viewerhdf5impl.h> 2 #include <petsclayouthdf5.h> /*I "petsclayoutdf5.h" I*/ 3 #include <petscis.h> /*I "petscis.h" I*/ 4 5 #if defined(PETSC_HAVE_HDF5) 6 7 struct _n_HDF5ReadCtx { 8 const char *name; 9 hid_t file, group, dataset, dataspace; 10 int lenInd, bsInd, complexInd, rdim; 11 hsize_t *dims; 12 PetscBool complexVal, dim2; 13 14 // Needed for compression 15 PetscInt runs; 16 PetscInt *cind; 17 }; 18 typedef struct _n_HDF5ReadCtx *HDF5ReadCtx; 19 20 PetscErrorCode PetscViewerHDF5CheckTimestepping_Internal(PetscViewer viewer, const char name[]) 21 { 22 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 23 PetscBool timestepping = PETSC_FALSE; 24 25 PetscFunctionBegin; 26 PetscCall(PetscViewerHDF5ReadAttribute(viewer, name, "timestepping", PETSC_BOOL, &hdf5->defTimestepping, ×tepping)); 27 if (timestepping != hdf5->timestepping) { 28 char *group; 29 30 PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group)); 31 SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Dataset %s/%s stored with timesteps? %s Timestepping pushed? %s", group, name, PetscBools[timestepping], PetscBools[hdf5->timestepping]); 32 } 33 PetscFunctionReturn(PETSC_SUCCESS); 34 } 35 36 static PetscErrorCode PetscViewerHDF5ReadInitialize_Private(PetscViewer viewer, const char name[], HDF5ReadCtx *ctx) 37 { 38 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 39 HDF5ReadCtx h = NULL; 40 41 PetscFunctionBegin; 42 PetscCall(PetscViewerHDF5CheckTimestepping_Internal(viewer, name)); 43 PetscCall(PetscNew(&h)); 44 h->name = name; 45 PetscCall(PetscViewerHDF5OpenGroup(viewer, NULL, &h->file, &h->group)); 46 PetscCallHDF5Return(h->dataset, H5Dopen2, (h->group, name, H5P_DEFAULT)); 47 PetscCallHDF5Return(h->dataspace, H5Dget_space, (h->dataset)); 48 PetscCall(PetscViewerHDF5ReadAttribute(viewer, name, "complex", PETSC_BOOL, &h->complexVal, &h->complexVal)); 49 if (!hdf5->horizontal) { 50 /* MATLAB stores column vectors horizontally */ 51 PetscCall(PetscViewerHDF5HasAttribute(viewer, name, "MATLAB_class", &hdf5->horizontal)); 52 } 53 h->runs = 0; 54 h->cind = NULL; 55 *ctx = h; 56 PetscFunctionReturn(PETSC_SUCCESS); 57 } 58 59 static PetscErrorCode PetscViewerHDF5ReadFinalize_Private(PetscViewer viewer, HDF5ReadCtx *ctx) 60 { 61 HDF5ReadCtx h; 62 63 PetscFunctionBegin; 64 h = *ctx; 65 PetscCallHDF5(H5Gclose, (h->group)); 66 PetscCallHDF5(H5Sclose, (h->dataspace)); 67 PetscCallHDF5(H5Dclose, (h->dataset)); 68 PetscCall(PetscFree((*ctx)->dims)); 69 PetscCall(PetscFree((*ctx)->cind)); 70 PetscCall(PetscFree(*ctx)); 71 PetscFunctionReturn(PETSC_SUCCESS); 72 } 73 74 // Need forward declaration because we have a cyclic call chain 75 static PetscErrorCode PetscViewerHDF5Load_Internal(PetscViewer, const char[], PetscBool, PetscLayout, hid_t, void **); 76 77 static PetscErrorCode PetscViewerHDF5ReadSizes_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscBool uncompress, PetscBool setup, PetscLayout *map_) 78 { 79 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 80 PetscInt bs, N; 81 PetscLayout map; 82 PetscBool compressed; 83 84 PetscFunctionBegin; 85 if (!*map_) PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)viewer), map_)); 86 map = *map_; 87 88 PetscCall(PetscViewerHDF5HasAttribute(viewer, ctx->name, "compressed", &compressed)); 89 if (compressed && uncompress) { 90 hid_t inttype; 91 PetscLayout cmap; 92 PetscInt *lcind; 93 PetscMPIInt *counts, *displs; 94 const PetscInt *range; 95 PetscInt N = 0; 96 PetscMPIInt size; 97 MPI_Comm comm; 98 99 #if defined(PETSC_USE_64BIT_INDICES) 100 inttype = H5T_NATIVE_LLONG; 101 #else 102 inttype = H5T_NATIVE_INT; 103 #endif 104 PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm)); 105 PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)viewer), &cmap)); 106 cmap->bs = 3; 107 PetscCall(PetscViewerHDF5Load_Internal(viewer, ctx->name, PETSC_FALSE, cmap, inttype, (void **)&lcind)); 108 PetscCheck(!(cmap->n % 3), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Compressed IS must have an even number of entries, not %" PetscInt_FMT, cmap->n); 109 for (PetscInt i = 0; i < cmap->n / 3; ++i) N += lcind[i * 3 + 0]; 110 PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &N, 1, MPIU_INT, MPIU_SUM, comm)); 111 ctx->runs = cmap->N / 3; 112 PetscCall(PetscMalloc1(cmap->N, &ctx->cind)); 113 PetscCallMPI(MPI_Comm_size(comm, &size)); 114 PetscCall(PetscLayoutGetRanges(cmap, &range)); 115 PetscCall(PetscMalloc2(size, &counts, size, &displs)); 116 for (PetscInt r = 0; r < size; ++r) { 117 PetscCall(PetscMPIIntCast(range[r + 1] - range[r], &counts[r])); 118 PetscCall(PetscMPIIntCast(range[r], &displs[r])); 119 } 120 PetscCallMPI(MPI_Allgatherv(lcind, cmap->n, MPIU_INT, ctx->cind, counts, displs, MPIU_INT, comm)); 121 PetscCall(PetscFree2(counts, displs)); 122 PetscCall(PetscFree(lcind)); 123 PetscCall(PetscLayoutDestroy(&cmap)); 124 125 ctx->dim2 = PETSC_FALSE; 126 ctx->rdim = 1; 127 ctx->lenInd = 0; 128 PetscCall(PetscMalloc1(ctx->rdim, &ctx->dims)); 129 ctx->dims[0] = N; 130 bs = 1; 131 goto layout; 132 } 133 134 /* Get actual number of dimensions in dataset */ 135 PetscCallHDF5Return(ctx->rdim, H5Sget_simple_extent_dims, (ctx->dataspace, NULL, NULL)); 136 PetscCall(PetscMalloc1(ctx->rdim, &ctx->dims)); 137 PetscCallHDF5Return(ctx->rdim, H5Sget_simple_extent_dims, (ctx->dataspace, ctx->dims, NULL)); 138 139 /* 140 Dimensions are in this order: 141 [0] timesteps (optional) 142 [lenInd] entries (numbers or blocks) 143 ... 144 [bsInd] entries of blocks (optional) 145 [bsInd+1] real & imaginary part (optional) 146 = rdim-1 147 */ 148 149 /* Get entries dimension index */ 150 ctx->lenInd = 0; 151 if (hdf5->timestepping) ++ctx->lenInd; 152 153 /* Get block dimension index */ 154 if (ctx->complexVal) { 155 ctx->bsInd = ctx->rdim - 2; 156 ctx->complexInd = ctx->rdim - 1; 157 } else { 158 ctx->bsInd = ctx->rdim - 1; 159 ctx->complexInd = -1; 160 } 161 PetscCheck(ctx->lenInd <= ctx->bsInd, PetscObjectComm((PetscObject)viewer), PETSC_ERR_PLIB, "Calculated block dimension index = %d < %d = length dimension index.", ctx->bsInd, ctx->lenInd); 162 PetscCheck(ctx->bsInd <= ctx->rdim - 1, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Calculated block dimension index = %d > %d = total number of dimensions - 1.", ctx->bsInd, ctx->rdim - 1); 163 PetscCheck(!ctx->complexVal || ctx->dims[ctx->complexInd] == 2, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Complex numbers must have exactly 2 parts (%" PRIuHSIZE ")", ctx->dims[ctx->complexInd]); 164 165 if (hdf5->horizontal) { 166 /* support horizontal 1D arrays (MATLAB vectors) - swap meaning of blocks and entries */ 167 int t = ctx->lenInd; 168 ctx->lenInd = ctx->bsInd; 169 ctx->bsInd = t; 170 } 171 172 /* Get block size */ 173 ctx->dim2 = PETSC_FALSE; 174 if (ctx->lenInd == ctx->bsInd) { 175 bs = 1; /* support vectors stored as 1D array */ 176 } else { 177 bs = (PetscInt)ctx->dims[ctx->bsInd]; 178 if (bs == 1) ctx->dim2 = PETSC_TRUE; /* vector with blocksize of 1, still stored as 2D array */ 179 } 180 181 layout: 182 /* Get global size */ 183 PetscCall(PetscIntCast(bs * ctx->dims[ctx->lenInd], &N)); 184 185 /* Set global size, blocksize and type if not yet set */ 186 if (map->bs < 0) { 187 PetscCall(PetscLayoutSetBlockSize(map, bs)); 188 } else PetscCheck(map->bs == bs, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Block size of array in file is %" PetscInt_FMT ", not %" PetscInt_FMT " as expected", bs, map->bs); 189 if (map->N < 0) { 190 PetscCall(PetscLayoutSetSize(map, N)); 191 } else PetscCheck(map->N == N, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Global size of array %s in file is %" PetscInt_FMT ", not %" PetscInt_FMT " as expected", ctx->name, N, map->N); 192 if (setup) PetscCall(PetscLayoutSetUp(map)); 193 PetscFunctionReturn(PETSC_SUCCESS); 194 } 195 196 static PetscErrorCode PetscViewerHDF5ReadSelectHyperslab_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout map, hid_t *memspace) 197 { 198 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 199 hsize_t *count, *offset; 200 PetscInt bs, n, low; 201 int i; 202 203 PetscFunctionBegin; 204 /* Compute local size and ownership range */ 205 PetscCall(PetscLayoutSetUp(map)); 206 PetscCall(PetscLayoutGetBlockSize(map, &bs)); 207 PetscCall(PetscLayoutGetLocalSize(map, &n)); 208 PetscCall(PetscLayoutGetRange(map, &low, NULL)); 209 210 /* Each process defines a dataset and reads it from the hyperslab in the file */ 211 PetscCall(PetscMalloc2(ctx->rdim, &count, ctx->rdim, &offset)); 212 for (i = 0; i < ctx->rdim; i++) { 213 /* By default, select all entries with no offset */ 214 offset[i] = 0; 215 count[i] = ctx->dims[i]; 216 } 217 if (hdf5->timestepping) { 218 count[0] = 1; 219 offset[0] = hdf5->timestep; 220 } 221 { 222 PetscCall(PetscHDF5IntCast(n / bs, &count[ctx->lenInd])); 223 PetscCall(PetscHDF5IntCast(low / bs, &offset[ctx->lenInd])); 224 } 225 PetscCallHDF5Return(*memspace, H5Screate_simple, (ctx->rdim, count, NULL)); 226 PetscCallHDF5(H5Sselect_hyperslab, (ctx->dataspace, H5S_SELECT_SET, offset, NULL, count, NULL)); 227 PetscCall(PetscFree2(count, offset)); 228 PetscFunctionReturn(PETSC_SUCCESS); 229 } 230 231 static PetscErrorCode PetscViewerHDF5ReadArray_Private(PetscViewer viewer, HDF5ReadCtx h, hid_t datatype, hid_t memspace, void *arr) 232 { 233 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 234 235 PetscFunctionBegin; 236 PetscCallHDF5(H5Dread, (h->dataset, datatype, memspace, h->dataspace, hdf5->dxpl_id, arr)); 237 PetscFunctionReturn(PETSC_SUCCESS); 238 } 239 240 static PetscErrorCode PetscViewerHDF5Load_Internal(PetscViewer viewer, const char name[], PetscBool uncompress, PetscLayout map, hid_t datatype, void **newarr) 241 { 242 PetscBool has; 243 char *group; 244 HDF5ReadCtx h = NULL; 245 hid_t memspace = 0; 246 size_t unitsize; 247 void *arr; 248 249 PetscFunctionBegin; 250 PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group)); 251 PetscCall(PetscViewerHDF5HasDataset(viewer, name, &has)); 252 PetscCheck(has, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", name, group); 253 PetscCall(PetscViewerHDF5ReadInitialize_Private(viewer, name, &h)); 254 #if defined(PETSC_USE_COMPLEX) 255 if (!h->complexVal) { 256 H5T_class_t clazz = H5Tget_class(datatype); 257 PetscCheck(clazz != H5T_FLOAT, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Dataset %s/%s is marked as real but PETSc is configured for complex scalars. The conversion is not yet implemented. Configure with --with-scalar-type=real to read this dataset", group ? group : "", name); 258 } 259 #else 260 PetscCheck(!h->complexVal, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Dataset %s/%s is marked as complex but PETSc is configured for real scalars. Configure with --with-scalar-type=complex to read this dataset", group, name); 261 #endif 262 263 PetscCall(PetscViewerHDF5ReadSizes_Private(viewer, h, uncompress, PETSC_TRUE, &map)); 264 PetscCall(PetscViewerHDF5ReadSelectHyperslab_Private(viewer, h, map, &memspace)); 265 266 if (h->runs && uncompress) { 267 PetscInt *ind; 268 269 PetscCall(PetscInfo(viewer, "Read compressed object with name %s of size %" PetscInt_FMT ":%" PetscInt_FMT "\n", name, map->n, map->N)); 270 // Each process stores the whole compression, so skip any leading parts 271 PetscCall(PetscMalloc1(map->n, &ind)); 272 for (PetscInt i = 0, off = 0; i < h->runs; ++i) { 273 for (PetscInt j = 0, inc = 0; j < h->cind[i * 3 + 0]; ++j, ++off, inc += h->cind[i * 3 + 1]) { 274 if (off >= map->rend) { 275 i = h->runs; 276 break; 277 } 278 if (off >= map->rstart) ind[off - map->rstart] = h->cind[i * 3 + 2] + inc; 279 } 280 } 281 *newarr = ind; 282 goto cleanup; 283 } 284 285 unitsize = H5Tget_size(datatype); 286 if (h->complexVal) unitsize *= 2; 287 /* unitsize is size_t i.e. always unsigned, so the negative check is pointless? */ 288 PetscCheck(unitsize > 0 && unitsize <= PetscMax(sizeof(PetscInt), sizeof(PetscScalar)), PETSC_COMM_SELF, PETSC_ERR_LIB, "Sanity check failed: HDF5 function H5Tget_size(datatype) returned suspicious value %zu", unitsize); 289 PetscCall(PetscMalloc(map->n * unitsize, &arr)); 290 291 PetscCall(PetscViewerHDF5ReadArray_Private(viewer, h, datatype, memspace, arr)); 292 *newarr = arr; 293 294 cleanup: 295 PetscCallHDF5(H5Sclose, (memspace)); 296 PetscCall(PetscViewerHDF5ReadFinalize_Private(viewer, &h)); 297 PetscCall(PetscFree(group)); 298 PetscFunctionReturn(PETSC_SUCCESS); 299 } 300 301 /*@C 302 PetscViewerHDF5Load - Read a raw array from the `PETSCVIEWERHDF5` dataset in parallel 303 304 Collective; No Fortran Support 305 306 Input Parameters: 307 + viewer - The `PETSCVIEWERHDF5` viewer 308 . name - The dataset name 309 - datatype - The HDF5 datatype of the items in the dataset 310 311 Input/Output Parameter: 312 . map - The layout which specifies array partitioning, on output the 313 set up layout (with global size and blocksize according to dataset) 314 315 Output Parameter: 316 . newarr - The partitioned array, a memory image of the given dataset 317 318 Level: developer 319 320 Notes: 321 This is intended mainly for internal use; users should use higher level routines such as `ISLoad()`, `VecLoad()`, `DMLoad()`. 322 323 The array is partitioned according to the given `PetscLayout` which is converted to an HDF5 hyperslab. 324 325 This name is relative to the current group returned by `PetscViewerHDF5OpenGroup()`. 326 327 .seealso: `PetscViewer`, `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5ReadSizes()`, 328 `VecLoad()`, `ISLoad()`, `PetscLayout` 329 @*/ 330 PetscErrorCode PetscViewerHDF5Load(PetscViewer viewer, const char name[], PetscLayout map, hid_t datatype, void **newarr) 331 { 332 PetscFunctionBegin; 333 PetscCall(PetscViewerHDF5Load_Internal(viewer, name, PETSC_TRUE, map, datatype, newarr)); 334 PetscFunctionReturn(PETSC_SUCCESS); 335 } 336 337 /*@C 338 PetscViewerHDF5ReadSizes - Read block size and global size of a `Vec` or `IS` stored in an HDF5 file. 339 340 Input Parameters: 341 + viewer - The `PETSCVIEWERHDF5` viewer 342 - name - The dataset name 343 344 Output Parameters: 345 + bs - block size 346 - N - global size 347 348 Level: advanced 349 350 Notes: 351 The dataset is stored as an HDF5 dataspace with 1-4 dimensions in the order 352 1) # timesteps (optional), 2) # blocks, 3) # elements per block (optional), 4) real and imaginary part (only for complex). 353 354 The dataset can be stored as a 2D dataspace even if its blocksize is 1; see `PetscViewerHDF5SetBaseDimension2()`. 355 356 .seealso: `PetscViewer`, `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `VecLoad()`, `ISLoad()`, `VecGetSize()`, `ISGetSize()`, `PetscViewerHDF5SetBaseDimension2()` 357 @*/ 358 PetscErrorCode PetscViewerHDF5ReadSizes(PetscViewer viewer, const char name[], PetscInt *bs, PetscInt *N) 359 { 360 HDF5ReadCtx h = NULL; 361 PetscLayout map = NULL; 362 363 PetscFunctionBegin; 364 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 365 PetscCall(PetscViewerHDF5ReadInitialize_Private(viewer, name, &h)); 366 PetscCall(PetscViewerHDF5ReadSizes_Private(viewer, h, PETSC_TRUE, PETSC_FALSE, &map)); 367 PetscCall(PetscViewerHDF5ReadFinalize_Private(viewer, &h)); 368 if (bs) *bs = map->bs; 369 if (N) *N = map->N; 370 PetscCall(PetscLayoutDestroy(&map)); 371 PetscFunctionReturn(PETSC_SUCCESS); 372 } 373 374 #endif /* defined(PETSC_HAVE_HDF5) */ 375