110b424faSBarry Smith #include <petsc/private/viewerhdf5impl.h> 210b424faSBarry Smith #include <petsclayouthdf5.h> /*I "petsclayoutdf5.h" I*/ 310b424faSBarry Smith #include <petscis.h> /*I "petscis.h" I*/ 410b424faSBarry Smith 510b424faSBarry Smith #if defined(PETSC_HAVE_HDF5) 610b424faSBarry Smith 710b424faSBarry Smith struct _n_HDF5ReadCtx { 8*21c42226SMatthew G. Knepley const char *name; 910b424faSBarry Smith hid_t file, group, dataset, dataspace; 1010b424faSBarry Smith int lenInd, bsInd, complexInd, rdim; 1110b424faSBarry Smith hsize_t *dims; 1210b424faSBarry Smith PetscBool complexVal, dim2; 13*21c42226SMatthew G. Knepley 14*21c42226SMatthew G. Knepley // Needed for compression 15*21c42226SMatthew G. Knepley PetscInt runs; 16*21c42226SMatthew G. Knepley PetscInt *cind; 1710b424faSBarry Smith }; 1810b424faSBarry Smith typedef struct _n_HDF5ReadCtx *HDF5ReadCtx; 1910b424faSBarry Smith 2010b424faSBarry Smith PetscErrorCode PetscViewerHDF5CheckTimestepping_Internal(PetscViewer viewer, const char name[]) 2110b424faSBarry Smith { 2210b424faSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 2310b424faSBarry Smith PetscBool timestepping = PETSC_FALSE; 2410b424faSBarry Smith 2510b424faSBarry Smith PetscFunctionBegin; 2610b424faSBarry Smith PetscCall(PetscViewerHDF5ReadAttribute(viewer, name, "timestepping", PETSC_BOOL, &hdf5->defTimestepping, ×tepping)); 2710b424faSBarry Smith if (timestepping != hdf5->timestepping) { 2810b424faSBarry Smith char *group; 2910b424faSBarry Smith 3010b424faSBarry Smith PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group)); 3110b424faSBarry Smith 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]); 3210b424faSBarry Smith } 3310b424faSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 3410b424faSBarry Smith } 3510b424faSBarry Smith 3610b424faSBarry Smith static PetscErrorCode PetscViewerHDF5ReadInitialize_Private(PetscViewer viewer, const char name[], HDF5ReadCtx *ctx) 3710b424faSBarry Smith { 3810b424faSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 3910b424faSBarry Smith HDF5ReadCtx h = NULL; 4010b424faSBarry Smith 4110b424faSBarry Smith PetscFunctionBegin; 4210b424faSBarry Smith PetscCall(PetscViewerHDF5CheckTimestepping_Internal(viewer, name)); 4310b424faSBarry Smith PetscCall(PetscNew(&h)); 44*21c42226SMatthew G. Knepley h->name = name; 4510b424faSBarry Smith PetscCall(PetscViewerHDF5OpenGroup(viewer, NULL, &h->file, &h->group)); 4610b424faSBarry Smith PetscCallHDF5Return(h->dataset, H5Dopen2, (h->group, name, H5P_DEFAULT)); 4710b424faSBarry Smith PetscCallHDF5Return(h->dataspace, H5Dget_space, (h->dataset)); 4810b424faSBarry Smith PetscCall(PetscViewerHDF5ReadAttribute(viewer, name, "complex", PETSC_BOOL, &h->complexVal, &h->complexVal)); 4910b424faSBarry Smith if (!hdf5->horizontal) { 5010b424faSBarry Smith /* MATLAB stores column vectors horizontally */ 5110b424faSBarry Smith PetscCall(PetscViewerHDF5HasAttribute(viewer, name, "MATLAB_class", &hdf5->horizontal)); 5210b424faSBarry Smith } 53*21c42226SMatthew G. Knepley h->runs = 0; 54*21c42226SMatthew G. Knepley h->cind = NULL; 5510b424faSBarry Smith *ctx = h; 5610b424faSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 5710b424faSBarry Smith } 5810b424faSBarry Smith 5910b424faSBarry Smith static PetscErrorCode PetscViewerHDF5ReadFinalize_Private(PetscViewer viewer, HDF5ReadCtx *ctx) 6010b424faSBarry Smith { 6110b424faSBarry Smith HDF5ReadCtx h; 6210b424faSBarry Smith 6310b424faSBarry Smith PetscFunctionBegin; 6410b424faSBarry Smith h = *ctx; 6510b424faSBarry Smith PetscCallHDF5(H5Gclose, (h->group)); 6610b424faSBarry Smith PetscCallHDF5(H5Sclose, (h->dataspace)); 6710b424faSBarry Smith PetscCallHDF5(H5Dclose, (h->dataset)); 6810b424faSBarry Smith PetscCall(PetscFree((*ctx)->dims)); 69*21c42226SMatthew G. Knepley PetscCall(PetscFree((*ctx)->cind)); 7010b424faSBarry Smith PetscCall(PetscFree(*ctx)); 7110b424faSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 7210b424faSBarry Smith } 7310b424faSBarry Smith 74*21c42226SMatthew G. Knepley // Need forward declaration because we have a cyclic call chain 75*21c42226SMatthew G. Knepley static PetscErrorCode PetscViewerHDF5Load_Internal(PetscViewer, const char[], PetscBool, PetscLayout, hid_t, void **); 76*21c42226SMatthew G. Knepley 77*21c42226SMatthew G. Knepley static PetscErrorCode PetscViewerHDF5ReadSizes_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscBool uncompress, PetscBool setup, PetscLayout *map_) 7810b424faSBarry Smith { 7910b424faSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 806497c311SBarry Smith PetscInt bs, N; 8110b424faSBarry Smith PetscLayout map; 82*21c42226SMatthew G. Knepley PetscBool compressed; 8310b424faSBarry Smith 8410b424faSBarry Smith PetscFunctionBegin; 854ad8454bSPierre Jolivet if (!*map_) PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)viewer), map_)); 8610b424faSBarry Smith map = *map_; 8710b424faSBarry Smith 88*21c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5HasAttribute(viewer, ctx->name, "compressed", &compressed)); 89*21c42226SMatthew G. Knepley if (compressed && uncompress) { 90*21c42226SMatthew G. Knepley hid_t inttype; 91*21c42226SMatthew G. Knepley PetscLayout cmap; 92*21c42226SMatthew G. Knepley PetscInt *lcind; 93*21c42226SMatthew G. Knepley PetscMPIInt *counts, *displs; 94*21c42226SMatthew G. Knepley const PetscInt *range; 95*21c42226SMatthew G. Knepley PetscInt N = 0; 96*21c42226SMatthew G. Knepley PetscMPIInt size; 97*21c42226SMatthew G. Knepley MPI_Comm comm; 98*21c42226SMatthew G. Knepley 99*21c42226SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 100*21c42226SMatthew G. Knepley inttype = H5T_NATIVE_LLONG; 101*21c42226SMatthew G. Knepley #else 102*21c42226SMatthew G. Knepley inttype = H5T_NATIVE_INT; 103*21c42226SMatthew G. Knepley #endif 104*21c42226SMatthew G. Knepley PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm)); 105*21c42226SMatthew G. Knepley PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)viewer), &cmap)); 106*21c42226SMatthew G. Knepley cmap->bs = 3; 107*21c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5Load_Internal(viewer, ctx->name, PETSC_FALSE, cmap, inttype, (void **)&lcind)); 108*21c42226SMatthew G. Knepley 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*21c42226SMatthew G. Knepley for (PetscInt i = 0; i < cmap->n / 3; ++i) N += lcind[i * 3 + 0]; 110*21c42226SMatthew G. Knepley PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &N, 1, MPIU_INT, MPIU_SUM, comm)); 111*21c42226SMatthew G. Knepley ctx->runs = cmap->N / 3; 112*21c42226SMatthew G. Knepley PetscCall(PetscMalloc1(cmap->N, &ctx->cind)); 113*21c42226SMatthew G. Knepley PetscCallMPI(MPI_Comm_size(comm, &size)); 114*21c42226SMatthew G. Knepley PetscCall(PetscLayoutGetRanges(cmap, &range)); 115*21c42226SMatthew G. Knepley PetscCall(PetscMalloc2(size, &counts, size, &displs)); 116*21c42226SMatthew G. Knepley for (PetscInt r = 0; r < size; ++r) { 117*21c42226SMatthew G. Knepley PetscCall(PetscMPIIntCast(range[r + 1] - range[r], &counts[r])); 118*21c42226SMatthew G. Knepley PetscCall(PetscMPIIntCast(range[r], &displs[r])); 119*21c42226SMatthew G. Knepley } 120*21c42226SMatthew G. Knepley PetscCallMPI(MPI_Allgatherv(lcind, cmap->n, MPIU_INT, ctx->cind, counts, displs, MPIU_INT, comm)); 121*21c42226SMatthew G. Knepley PetscCall(PetscFree2(counts, displs)); 122*21c42226SMatthew G. Knepley PetscCall(PetscFree(lcind)); 123*21c42226SMatthew G. Knepley PetscCall(PetscLayoutDestroy(&cmap)); 124*21c42226SMatthew G. Knepley 125*21c42226SMatthew G. Knepley ctx->dim2 = PETSC_FALSE; 126*21c42226SMatthew G. Knepley ctx->rdim = 1; 127*21c42226SMatthew G. Knepley ctx->lenInd = 0; 128*21c42226SMatthew G. Knepley PetscCall(PetscMalloc1(ctx->rdim, &ctx->dims)); 129*21c42226SMatthew G. Knepley ctx->dims[0] = N; 130*21c42226SMatthew G. Knepley bs = 1; 131*21c42226SMatthew G. Knepley goto layout; 132*21c42226SMatthew G. Knepley } 133*21c42226SMatthew G. Knepley 13410b424faSBarry Smith /* Get actual number of dimensions in dataset */ 13510b424faSBarry Smith PetscCallHDF5Return(ctx->rdim, H5Sget_simple_extent_dims, (ctx->dataspace, NULL, NULL)); 13610b424faSBarry Smith PetscCall(PetscMalloc1(ctx->rdim, &ctx->dims)); 13710b424faSBarry Smith PetscCallHDF5Return(ctx->rdim, H5Sget_simple_extent_dims, (ctx->dataspace, ctx->dims, NULL)); 13810b424faSBarry Smith 13910b424faSBarry Smith /* 14010b424faSBarry Smith Dimensions are in this order: 14110b424faSBarry Smith [0] timesteps (optional) 14210b424faSBarry Smith [lenInd] entries (numbers or blocks) 14310b424faSBarry Smith ... 14410b424faSBarry Smith [bsInd] entries of blocks (optional) 14510b424faSBarry Smith [bsInd+1] real & imaginary part (optional) 14610b424faSBarry Smith = rdim-1 14710b424faSBarry Smith */ 14810b424faSBarry Smith 14910b424faSBarry Smith /* Get entries dimension index */ 15010b424faSBarry Smith ctx->lenInd = 0; 15110b424faSBarry Smith if (hdf5->timestepping) ++ctx->lenInd; 15210b424faSBarry Smith 15310b424faSBarry Smith /* Get block dimension index */ 15410b424faSBarry Smith if (ctx->complexVal) { 15510b424faSBarry Smith ctx->bsInd = ctx->rdim - 2; 15610b424faSBarry Smith ctx->complexInd = ctx->rdim - 1; 15710b424faSBarry Smith } else { 15810b424faSBarry Smith ctx->bsInd = ctx->rdim - 1; 15910b424faSBarry Smith ctx->complexInd = -1; 16010b424faSBarry Smith } 16110b424faSBarry Smith PetscCheck(ctx->lenInd <= ctx->bsInd, PetscObjectComm((PetscObject)viewer), PETSC_ERR_PLIB, "Calculated block dimension index = %d < %d = length dimension index.", ctx->bsInd, ctx->lenInd); 16210b424faSBarry Smith 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); 16310b424faSBarry Smith 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]); 16410b424faSBarry Smith 16510b424faSBarry Smith if (hdf5->horizontal) { 16610b424faSBarry Smith /* support horizontal 1D arrays (MATLAB vectors) - swap meaning of blocks and entries */ 1679ad2cedaSPierre Jolivet int t = ctx->lenInd; 16810b424faSBarry Smith ctx->lenInd = ctx->bsInd; 16910b424faSBarry Smith ctx->bsInd = t; 17010b424faSBarry Smith } 17110b424faSBarry Smith 17210b424faSBarry Smith /* Get block size */ 17310b424faSBarry Smith ctx->dim2 = PETSC_FALSE; 17410b424faSBarry Smith if (ctx->lenInd == ctx->bsInd) { 17510b424faSBarry Smith bs = 1; /* support vectors stored as 1D array */ 17610b424faSBarry Smith } else { 17710b424faSBarry Smith bs = (PetscInt)ctx->dims[ctx->bsInd]; 17810b424faSBarry Smith if (bs == 1) ctx->dim2 = PETSC_TRUE; /* vector with blocksize of 1, still stored as 2D array */ 17910b424faSBarry Smith } 18010b424faSBarry Smith 181*21c42226SMatthew G. Knepley layout: 18210b424faSBarry Smith /* Get global size */ 1836497c311SBarry Smith PetscCall(PetscIntCast(bs * ctx->dims[ctx->lenInd], &N)); 18410b424faSBarry Smith 18510b424faSBarry Smith /* Set global size, blocksize and type if not yet set */ 18610b424faSBarry Smith if (map->bs < 0) { 18710b424faSBarry Smith PetscCall(PetscLayoutSetBlockSize(map, bs)); 18810b424faSBarry Smith } 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); 18910b424faSBarry Smith if (map->N < 0) { 19010b424faSBarry Smith PetscCall(PetscLayoutSetSize(map, N)); 191*21c42226SMatthew G. Knepley } 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); 19210b424faSBarry Smith if (setup) PetscCall(PetscLayoutSetUp(map)); 19310b424faSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 19410b424faSBarry Smith } 19510b424faSBarry Smith 19610b424faSBarry Smith static PetscErrorCode PetscViewerHDF5ReadSelectHyperslab_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout map, hid_t *memspace) 19710b424faSBarry Smith { 19810b424faSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 19910b424faSBarry Smith hsize_t *count, *offset; 20010b424faSBarry Smith PetscInt bs, n, low; 20110b424faSBarry Smith int i; 20210b424faSBarry Smith 20310b424faSBarry Smith PetscFunctionBegin; 20410b424faSBarry Smith /* Compute local size and ownership range */ 20510b424faSBarry Smith PetscCall(PetscLayoutSetUp(map)); 20610b424faSBarry Smith PetscCall(PetscLayoutGetBlockSize(map, &bs)); 20710b424faSBarry Smith PetscCall(PetscLayoutGetLocalSize(map, &n)); 20810b424faSBarry Smith PetscCall(PetscLayoutGetRange(map, &low, NULL)); 20910b424faSBarry Smith 21010b424faSBarry Smith /* Each process defines a dataset and reads it from the hyperslab in the file */ 21110b424faSBarry Smith PetscCall(PetscMalloc2(ctx->rdim, &count, ctx->rdim, &offset)); 21210b424faSBarry Smith for (i = 0; i < ctx->rdim; i++) { 21310b424faSBarry Smith /* By default, select all entries with no offset */ 21410b424faSBarry Smith offset[i] = 0; 21510b424faSBarry Smith count[i] = ctx->dims[i]; 21610b424faSBarry Smith } 21710b424faSBarry Smith if (hdf5->timestepping) { 21810b424faSBarry Smith count[0] = 1; 21910b424faSBarry Smith offset[0] = hdf5->timestep; 22010b424faSBarry Smith } 22110b424faSBarry Smith { 22210b424faSBarry Smith PetscCall(PetscHDF5IntCast(n / bs, &count[ctx->lenInd])); 22310b424faSBarry Smith PetscCall(PetscHDF5IntCast(low / bs, &offset[ctx->lenInd])); 22410b424faSBarry Smith } 22510b424faSBarry Smith PetscCallHDF5Return(*memspace, H5Screate_simple, (ctx->rdim, count, NULL)); 22610b424faSBarry Smith PetscCallHDF5(H5Sselect_hyperslab, (ctx->dataspace, H5S_SELECT_SET, offset, NULL, count, NULL)); 22710b424faSBarry Smith PetscCall(PetscFree2(count, offset)); 22810b424faSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 22910b424faSBarry Smith } 23010b424faSBarry Smith 23110b424faSBarry Smith static PetscErrorCode PetscViewerHDF5ReadArray_Private(PetscViewer viewer, HDF5ReadCtx h, hid_t datatype, hid_t memspace, void *arr) 23210b424faSBarry Smith { 23310b424faSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 23410b424faSBarry Smith 23510b424faSBarry Smith PetscFunctionBegin; 23610b424faSBarry Smith PetscCallHDF5(H5Dread, (h->dataset, datatype, memspace, h->dataspace, hdf5->dxpl_id, arr)); 23710b424faSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 23810b424faSBarry Smith } 23910b424faSBarry Smith 240*21c42226SMatthew G. Knepley static PetscErrorCode PetscViewerHDF5Load_Internal(PetscViewer viewer, const char name[], PetscBool uncompress, PetscLayout map, hid_t datatype, void **newarr) 241*21c42226SMatthew G. Knepley { 242*21c42226SMatthew G. Knepley PetscBool has; 243*21c42226SMatthew G. Knepley char *group; 244*21c42226SMatthew G. Knepley HDF5ReadCtx h = NULL; 245*21c42226SMatthew G. Knepley hid_t memspace = 0; 246*21c42226SMatthew G. Knepley size_t unitsize; 247*21c42226SMatthew G. Knepley void *arr; 248*21c42226SMatthew G. Knepley 249*21c42226SMatthew G. Knepley PetscFunctionBegin; 250*21c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group)); 251*21c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5HasDataset(viewer, name, &has)); 252*21c42226SMatthew G. Knepley PetscCheck(has, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", name, group); 253*21c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5ReadInitialize_Private(viewer, name, &h)); 254*21c42226SMatthew G. Knepley #if defined(PETSC_USE_COMPLEX) 255*21c42226SMatthew G. Knepley if (!h->complexVal) { 256*21c42226SMatthew G. Knepley H5T_class_t clazz = H5Tget_class(datatype); 257*21c42226SMatthew G. Knepley 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*21c42226SMatthew G. Knepley } 259*21c42226SMatthew G. Knepley #else 260*21c42226SMatthew G. Knepley 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*21c42226SMatthew G. Knepley #endif 262*21c42226SMatthew G. Knepley 263*21c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5ReadSizes_Private(viewer, h, uncompress, PETSC_TRUE, &map)); 264*21c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5ReadSelectHyperslab_Private(viewer, h, map, &memspace)); 265*21c42226SMatthew G. Knepley 266*21c42226SMatthew G. Knepley if (h->runs && uncompress) { 267*21c42226SMatthew G. Knepley PetscInt *ind; 268*21c42226SMatthew G. Knepley 269*21c42226SMatthew G. Knepley PetscCall(PetscInfo(viewer, "Read compressed object with name %s of size %" PetscInt_FMT ":%" PetscInt_FMT "\n", name, map->n, map->N)); 270*21c42226SMatthew G. Knepley // Each process stores the whole compression, so skip any leading parts 271*21c42226SMatthew G. Knepley PetscCall(PetscMalloc1(map->n, &ind)); 272*21c42226SMatthew G. Knepley for (PetscInt i = 0, off = 0; i < h->runs; ++i) { 273*21c42226SMatthew G. Knepley for (PetscInt j = 0, inc = 0; j < h->cind[i * 3 + 0]; ++j, ++off, inc += h->cind[i * 3 + 1]) { 274*21c42226SMatthew G. Knepley if (off >= map->rend) { 275*21c42226SMatthew G. Knepley i = h->runs; 276*21c42226SMatthew G. Knepley break; 277*21c42226SMatthew G. Knepley } 278*21c42226SMatthew G. Knepley if (off >= map->rstart) ind[off - map->rstart] = h->cind[i * 3 + 2] + inc; 279*21c42226SMatthew G. Knepley } 280*21c42226SMatthew G. Knepley } 281*21c42226SMatthew G. Knepley *newarr = ind; 282*21c42226SMatthew G. Knepley goto cleanup; 283*21c42226SMatthew G. Knepley } 284*21c42226SMatthew G. Knepley 285*21c42226SMatthew G. Knepley unitsize = H5Tget_size(datatype); 286*21c42226SMatthew G. Knepley if (h->complexVal) unitsize *= 2; 287*21c42226SMatthew G. Knepley /* unitsize is size_t i.e. always unsigned, so the negative check is pointless? */ 288*21c42226SMatthew G. Knepley 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*21c42226SMatthew G. Knepley PetscCall(PetscMalloc(map->n * unitsize, &arr)); 290*21c42226SMatthew G. Knepley 291*21c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5ReadArray_Private(viewer, h, datatype, memspace, arr)); 292*21c42226SMatthew G. Knepley *newarr = arr; 293*21c42226SMatthew G. Knepley 294*21c42226SMatthew G. Knepley cleanup: 295*21c42226SMatthew G. Knepley PetscCallHDF5(H5Sclose, (memspace)); 296*21c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5ReadFinalize_Private(viewer, &h)); 297*21c42226SMatthew G. Knepley PetscCall(PetscFree(group)); 298*21c42226SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 299*21c42226SMatthew G. Knepley } 300*21c42226SMatthew G. Knepley 30110b424faSBarry Smith /*@C 302a24573d6SBarry Smith PetscViewerHDF5Load - Read a raw array from the `PETSCVIEWERHDF5` dataset in parallel 30310b424faSBarry Smith 30410b424faSBarry Smith Collective; No Fortran Support 30510b424faSBarry Smith 30610b424faSBarry Smith Input Parameters: 30710b424faSBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 30810b424faSBarry Smith . name - The dataset name 30910b424faSBarry Smith - datatype - The HDF5 datatype of the items in the dataset 31010b424faSBarry Smith 31110b424faSBarry Smith Input/Output Parameter: 31210b424faSBarry Smith . map - The layout which specifies array partitioning, on output the 31310b424faSBarry Smith set up layout (with global size and blocksize according to dataset) 31410b424faSBarry Smith 31510b424faSBarry Smith Output Parameter: 31610b424faSBarry Smith . newarr - The partitioned array, a memory image of the given dataset 31710b424faSBarry Smith 31810b424faSBarry Smith Level: developer 31910b424faSBarry Smith 32010b424faSBarry Smith Notes: 32110b424faSBarry Smith This is intended mainly for internal use; users should use higher level routines such as `ISLoad()`, `VecLoad()`, `DMLoad()`. 32210b424faSBarry Smith 32310b424faSBarry Smith The array is partitioned according to the given `PetscLayout` which is converted to an HDF5 hyperslab. 32410b424faSBarry Smith 32510b424faSBarry Smith This name is relative to the current group returned by `PetscViewerHDF5OpenGroup()`. 32610b424faSBarry Smith 32710b424faSBarry Smith .seealso: `PetscViewer`, `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5ReadSizes()`, 328a24573d6SBarry Smith `VecLoad()`, `ISLoad()`, `PetscLayout` 32910b424faSBarry Smith @*/ 330cc4c1da9SBarry Smith PetscErrorCode PetscViewerHDF5Load(PetscViewer viewer, const char name[], PetscLayout map, hid_t datatype, void **newarr) 33110b424faSBarry Smith { 33210b424faSBarry Smith PetscFunctionBegin; 333*21c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5Load_Internal(viewer, name, PETSC_TRUE, map, datatype, newarr)); 33410b424faSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 33510b424faSBarry Smith } 33610b424faSBarry Smith 33710b424faSBarry Smith /*@C 33810b424faSBarry Smith PetscViewerHDF5ReadSizes - Read block size and global size of a `Vec` or `IS` stored in an HDF5 file. 33910b424faSBarry Smith 34010b424faSBarry Smith Input Parameters: 34110b424faSBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 34210b424faSBarry Smith - name - The dataset name 34310b424faSBarry Smith 34410b424faSBarry Smith Output Parameters: 34510b424faSBarry Smith + bs - block size 34610b424faSBarry Smith - N - global size 34710b424faSBarry Smith 34810b424faSBarry Smith Level: advanced 34910b424faSBarry Smith 35010b424faSBarry Smith Notes: 35110b424faSBarry Smith The dataset is stored as an HDF5 dataspace with 1-4 dimensions in the order 35210b424faSBarry Smith 1) # timesteps (optional), 2) # blocks, 3) # elements per block (optional), 4) real and imaginary part (only for complex). 35310b424faSBarry Smith 35410b424faSBarry Smith The dataset can be stored as a 2D dataspace even if its blocksize is 1; see `PetscViewerHDF5SetBaseDimension2()`. 35510b424faSBarry Smith 35610b424faSBarry Smith .seealso: `PetscViewer`, `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `VecLoad()`, `ISLoad()`, `VecGetSize()`, `ISGetSize()`, `PetscViewerHDF5SetBaseDimension2()` 35710b424faSBarry Smith @*/ 35810b424faSBarry Smith PetscErrorCode PetscViewerHDF5ReadSizes(PetscViewer viewer, const char name[], PetscInt *bs, PetscInt *N) 35910b424faSBarry Smith { 36010b424faSBarry Smith HDF5ReadCtx h = NULL; 36110b424faSBarry Smith PetscLayout map = NULL; 36210b424faSBarry Smith 36310b424faSBarry Smith PetscFunctionBegin; 36410b424faSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 36510b424faSBarry Smith PetscCall(PetscViewerHDF5ReadInitialize_Private(viewer, name, &h)); 366*21c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5ReadSizes_Private(viewer, h, PETSC_TRUE, PETSC_FALSE, &map)); 36710b424faSBarry Smith PetscCall(PetscViewerHDF5ReadFinalize_Private(viewer, &h)); 36810b424faSBarry Smith if (bs) *bs = map->bs; 36910b424faSBarry Smith if (N) *N = map->N; 37010b424faSBarry Smith PetscCall(PetscLayoutDestroy(&map)); 37110b424faSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 37210b424faSBarry Smith } 37310b424faSBarry Smith 37410b424faSBarry Smith #endif /* defined(PETSC_HAVE_HDF5) */ 375