110b424faSBarry Smith #include <petsc/private/viewerhdf5impl.h> 2*2e1d0745SJose E. Roman #include <petsclayouthdf5.h> /*I "petsclayouthdf5.h" I*/ 310b424faSBarry Smith 410b424faSBarry Smith struct _n_HDF5ReadCtx { 521c42226SMatthew G. Knepley const char *name; 610b424faSBarry Smith hid_t file, group, dataset, dataspace; 710b424faSBarry Smith int lenInd, bsInd, complexInd, rdim; 810b424faSBarry Smith hsize_t *dims; 910b424faSBarry Smith PetscBool complexVal, dim2; 1021c42226SMatthew G. Knepley 1121c42226SMatthew G. Knepley // Needed for compression 1221c42226SMatthew G. Knepley PetscInt runs; 1321c42226SMatthew G. Knepley PetscInt *cind; 1410b424faSBarry Smith }; 1510b424faSBarry Smith typedef struct _n_HDF5ReadCtx *HDF5ReadCtx; 1610b424faSBarry Smith 1710b424faSBarry Smith PetscErrorCode PetscViewerHDF5CheckTimestepping_Internal(PetscViewer viewer, const char name[]) 1810b424faSBarry Smith { 1910b424faSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 2010b424faSBarry Smith PetscBool timestepping = PETSC_FALSE; 2110b424faSBarry Smith 2210b424faSBarry Smith PetscFunctionBegin; 2310b424faSBarry Smith PetscCall(PetscViewerHDF5ReadAttribute(viewer, name, "timestepping", PETSC_BOOL, &hdf5->defTimestepping, ×tepping)); 2410b424faSBarry Smith if (timestepping != hdf5->timestepping) { 2510b424faSBarry Smith char *group; 2610b424faSBarry Smith 2710b424faSBarry Smith PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group)); 2810b424faSBarry 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]); 2910b424faSBarry Smith } 3010b424faSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 3110b424faSBarry Smith } 3210b424faSBarry Smith 3310b424faSBarry Smith static PetscErrorCode PetscViewerHDF5ReadInitialize_Private(PetscViewer viewer, const char name[], HDF5ReadCtx *ctx) 3410b424faSBarry Smith { 3510b424faSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 3610b424faSBarry Smith HDF5ReadCtx h = NULL; 3710b424faSBarry Smith 3810b424faSBarry Smith PetscFunctionBegin; 3910b424faSBarry Smith PetscCall(PetscViewerHDF5CheckTimestepping_Internal(viewer, name)); 4010b424faSBarry Smith PetscCall(PetscNew(&h)); 4121c42226SMatthew G. Knepley h->name = name; 4210b424faSBarry Smith PetscCall(PetscViewerHDF5OpenGroup(viewer, NULL, &h->file, &h->group)); 4310b424faSBarry Smith PetscCallHDF5Return(h->dataset, H5Dopen2, (h->group, name, H5P_DEFAULT)); 4410b424faSBarry Smith PetscCallHDF5Return(h->dataspace, H5Dget_space, (h->dataset)); 4510b424faSBarry Smith PetscCall(PetscViewerHDF5ReadAttribute(viewer, name, "complex", PETSC_BOOL, &h->complexVal, &h->complexVal)); 4610b424faSBarry Smith if (!hdf5->horizontal) { 4710b424faSBarry Smith /* MATLAB stores column vectors horizontally */ 4810b424faSBarry Smith PetscCall(PetscViewerHDF5HasAttribute(viewer, name, "MATLAB_class", &hdf5->horizontal)); 4910b424faSBarry Smith } 5021c42226SMatthew G. Knepley h->runs = 0; 5121c42226SMatthew G. Knepley h->cind = NULL; 5210b424faSBarry Smith *ctx = h; 5310b424faSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 5410b424faSBarry Smith } 5510b424faSBarry Smith 5610b424faSBarry Smith static PetscErrorCode PetscViewerHDF5ReadFinalize_Private(PetscViewer viewer, HDF5ReadCtx *ctx) 5710b424faSBarry Smith { 5810b424faSBarry Smith HDF5ReadCtx h; 5910b424faSBarry Smith 6010b424faSBarry Smith PetscFunctionBegin; 6110b424faSBarry Smith h = *ctx; 6210b424faSBarry Smith PetscCallHDF5(H5Gclose, (h->group)); 6310b424faSBarry Smith PetscCallHDF5(H5Sclose, (h->dataspace)); 6410b424faSBarry Smith PetscCallHDF5(H5Dclose, (h->dataset)); 6510b424faSBarry Smith PetscCall(PetscFree((*ctx)->dims)); 6621c42226SMatthew G. Knepley PetscCall(PetscFree((*ctx)->cind)); 6710b424faSBarry Smith PetscCall(PetscFree(*ctx)); 6810b424faSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 6910b424faSBarry Smith } 7010b424faSBarry Smith 7121c42226SMatthew G. Knepley // Need forward declaration because we have a cyclic call chain 7221c42226SMatthew G. Knepley static PetscErrorCode PetscViewerHDF5Load_Internal(PetscViewer, const char[], PetscBool, PetscLayout, hid_t, void **); 7321c42226SMatthew G. Knepley 7421c42226SMatthew G. Knepley static PetscErrorCode PetscViewerHDF5ReadSizes_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscBool uncompress, PetscBool setup, PetscLayout *map_) 7510b424faSBarry Smith { 7610b424faSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 776497c311SBarry Smith PetscInt bs, N; 7810b424faSBarry Smith PetscLayout map; 7921c42226SMatthew G. Knepley PetscBool compressed; 8010b424faSBarry Smith 8110b424faSBarry Smith PetscFunctionBegin; 824ad8454bSPierre Jolivet if (!*map_) PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)viewer), map_)); 8310b424faSBarry Smith map = *map_; 8410b424faSBarry Smith 8521c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5HasAttribute(viewer, ctx->name, "compressed", &compressed)); 8621c42226SMatthew G. Knepley if (compressed && uncompress) { 8721c42226SMatthew G. Knepley hid_t inttype; 8821c42226SMatthew G. Knepley PetscLayout cmap; 89fef1ebd0SPierre Jolivet PetscInt *lcind, N = 0; 90fef1ebd0SPierre Jolivet PetscMPIInt *counts, *displs, size, n; 9121c42226SMatthew G. Knepley const PetscInt *range; 9221c42226SMatthew G. Knepley MPI_Comm comm; 9321c42226SMatthew G. Knepley 9421c42226SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 9521c42226SMatthew G. Knepley inttype = H5T_NATIVE_LLONG; 9621c42226SMatthew G. Knepley #else 9721c42226SMatthew G. Knepley inttype = H5T_NATIVE_INT; 9821c42226SMatthew G. Knepley #endif 9921c42226SMatthew G. Knepley PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm)); 10021c42226SMatthew G. Knepley PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)viewer), &cmap)); 10121c42226SMatthew G. Knepley cmap->bs = 3; 10221c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5Load_Internal(viewer, ctx->name, PETSC_FALSE, cmap, inttype, (void **)&lcind)); 10321c42226SMatthew 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); 10421c42226SMatthew G. Knepley for (PetscInt i = 0; i < cmap->n / 3; ++i) N += lcind[i * 3 + 0]; 10521c42226SMatthew G. Knepley PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &N, 1, MPIU_INT, MPIU_SUM, comm)); 10621c42226SMatthew G. Knepley ctx->runs = cmap->N / 3; 10721c42226SMatthew G. Knepley PetscCall(PetscMalloc1(cmap->N, &ctx->cind)); 10821c42226SMatthew G. Knepley PetscCallMPI(MPI_Comm_size(comm, &size)); 10921c42226SMatthew G. Knepley PetscCall(PetscLayoutGetRanges(cmap, &range)); 11021c42226SMatthew G. Knepley PetscCall(PetscMalloc2(size, &counts, size, &displs)); 11121c42226SMatthew G. Knepley for (PetscInt r = 0; r < size; ++r) { 11221c42226SMatthew G. Knepley PetscCall(PetscMPIIntCast(range[r + 1] - range[r], &counts[r])); 11321c42226SMatthew G. Knepley PetscCall(PetscMPIIntCast(range[r], &displs[r])); 11421c42226SMatthew G. Knepley } 115fef1ebd0SPierre Jolivet PetscCall(PetscMPIIntCast(cmap->n, &n)); 116fef1ebd0SPierre Jolivet PetscCallMPI(MPI_Allgatherv(lcind, n, MPIU_INT, ctx->cind, counts, displs, MPIU_INT, comm)); 11721c42226SMatthew G. Knepley PetscCall(PetscFree2(counts, displs)); 11821c42226SMatthew G. Knepley PetscCall(PetscFree(lcind)); 11921c42226SMatthew G. Knepley PetscCall(PetscLayoutDestroy(&cmap)); 12021c42226SMatthew G. Knepley 12121c42226SMatthew G. Knepley ctx->dim2 = PETSC_FALSE; 12221c42226SMatthew G. Knepley ctx->rdim = 1; 12321c42226SMatthew G. Knepley ctx->lenInd = 0; 12421c42226SMatthew G. Knepley PetscCall(PetscMalloc1(ctx->rdim, &ctx->dims)); 12521c42226SMatthew G. Knepley ctx->dims[0] = N; 12621c42226SMatthew G. Knepley bs = 1; 12721c42226SMatthew G. Knepley goto layout; 12821c42226SMatthew G. Knepley } 12921c42226SMatthew G. Knepley 13010b424faSBarry Smith /* Get actual number of dimensions in dataset */ 13110b424faSBarry Smith PetscCallHDF5Return(ctx->rdim, H5Sget_simple_extent_dims, (ctx->dataspace, NULL, NULL)); 13210b424faSBarry Smith PetscCall(PetscMalloc1(ctx->rdim, &ctx->dims)); 13310b424faSBarry Smith PetscCallHDF5Return(ctx->rdim, H5Sget_simple_extent_dims, (ctx->dataspace, ctx->dims, NULL)); 13410b424faSBarry Smith 13510b424faSBarry Smith /* 13610b424faSBarry Smith Dimensions are in this order: 13710b424faSBarry Smith [0] timesteps (optional) 13810b424faSBarry Smith [lenInd] entries (numbers or blocks) 13910b424faSBarry Smith ... 14010b424faSBarry Smith [bsInd] entries of blocks (optional) 14110b424faSBarry Smith [bsInd+1] real & imaginary part (optional) 14210b424faSBarry Smith = rdim-1 14310b424faSBarry Smith */ 14410b424faSBarry Smith 14510b424faSBarry Smith /* Get entries dimension index */ 14610b424faSBarry Smith ctx->lenInd = 0; 14710b424faSBarry Smith if (hdf5->timestepping) ++ctx->lenInd; 14810b424faSBarry Smith 14910b424faSBarry Smith /* Get block dimension index */ 15010b424faSBarry Smith if (ctx->complexVal) { 15110b424faSBarry Smith ctx->bsInd = ctx->rdim - 2; 15210b424faSBarry Smith ctx->complexInd = ctx->rdim - 1; 15310b424faSBarry Smith } else { 15410b424faSBarry Smith ctx->bsInd = ctx->rdim - 1; 15510b424faSBarry Smith ctx->complexInd = -1; 15610b424faSBarry Smith } 15710b424faSBarry 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); 15810b424faSBarry 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); 15910b424faSBarry 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]); 16010b424faSBarry Smith 16110b424faSBarry Smith if (hdf5->horizontal) { 16210b424faSBarry Smith /* support horizontal 1D arrays (MATLAB vectors) - swap meaning of blocks and entries */ 1639ad2cedaSPierre Jolivet int t = ctx->lenInd; 16410b424faSBarry Smith ctx->lenInd = ctx->bsInd; 16510b424faSBarry Smith ctx->bsInd = t; 16610b424faSBarry Smith } 16710b424faSBarry Smith 16810b424faSBarry Smith /* Get block size */ 16910b424faSBarry Smith ctx->dim2 = PETSC_FALSE; 17010b424faSBarry Smith if (ctx->lenInd == ctx->bsInd) { 17110b424faSBarry Smith bs = 1; /* support vectors stored as 1D array */ 17210b424faSBarry Smith } else { 17310b424faSBarry Smith bs = (PetscInt)ctx->dims[ctx->bsInd]; 17410b424faSBarry Smith if (bs == 1) ctx->dim2 = PETSC_TRUE; /* vector with blocksize of 1, still stored as 2D array */ 17510b424faSBarry Smith } 17610b424faSBarry Smith 17721c42226SMatthew G. Knepley layout: 17810b424faSBarry Smith /* Get global size */ 1796497c311SBarry Smith PetscCall(PetscIntCast(bs * ctx->dims[ctx->lenInd], &N)); 18010b424faSBarry Smith 18110b424faSBarry Smith /* Set global size, blocksize and type if not yet set */ 18210b424faSBarry Smith if (map->bs < 0) { 18310b424faSBarry Smith PetscCall(PetscLayoutSetBlockSize(map, bs)); 18410b424faSBarry 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); 18510b424faSBarry Smith if (map->N < 0) { 18610b424faSBarry Smith PetscCall(PetscLayoutSetSize(map, N)); 18721c42226SMatthew 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); 18810b424faSBarry Smith if (setup) PetscCall(PetscLayoutSetUp(map)); 18910b424faSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 19010b424faSBarry Smith } 19110b424faSBarry Smith 19210b424faSBarry Smith static PetscErrorCode PetscViewerHDF5ReadSelectHyperslab_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout map, hid_t *memspace) 19310b424faSBarry Smith { 19410b424faSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 19510b424faSBarry Smith hsize_t *count, *offset; 19610b424faSBarry Smith PetscInt bs, n, low; 19710b424faSBarry Smith int i; 19810b424faSBarry Smith 19910b424faSBarry Smith PetscFunctionBegin; 20010b424faSBarry Smith /* Compute local size and ownership range */ 20110b424faSBarry Smith PetscCall(PetscLayoutSetUp(map)); 20210b424faSBarry Smith PetscCall(PetscLayoutGetBlockSize(map, &bs)); 20310b424faSBarry Smith PetscCall(PetscLayoutGetLocalSize(map, &n)); 20410b424faSBarry Smith PetscCall(PetscLayoutGetRange(map, &low, NULL)); 20510b424faSBarry Smith 20610b424faSBarry Smith /* Each process defines a dataset and reads it from the hyperslab in the file */ 20710b424faSBarry Smith PetscCall(PetscMalloc2(ctx->rdim, &count, ctx->rdim, &offset)); 20810b424faSBarry Smith for (i = 0; i < ctx->rdim; i++) { 20910b424faSBarry Smith /* By default, select all entries with no offset */ 21010b424faSBarry Smith offset[i] = 0; 21110b424faSBarry Smith count[i] = ctx->dims[i]; 21210b424faSBarry Smith } 21310b424faSBarry Smith if (hdf5->timestepping) { 21410b424faSBarry Smith count[0] = 1; 21510b424faSBarry Smith offset[0] = hdf5->timestep; 21610b424faSBarry Smith } 21710b424faSBarry Smith { 21810b424faSBarry Smith PetscCall(PetscHDF5IntCast(n / bs, &count[ctx->lenInd])); 21910b424faSBarry Smith PetscCall(PetscHDF5IntCast(low / bs, &offset[ctx->lenInd])); 22010b424faSBarry Smith } 22110b424faSBarry Smith PetscCallHDF5Return(*memspace, H5Screate_simple, (ctx->rdim, count, NULL)); 22210b424faSBarry Smith PetscCallHDF5(H5Sselect_hyperslab, (ctx->dataspace, H5S_SELECT_SET, offset, NULL, count, NULL)); 22310b424faSBarry Smith PetscCall(PetscFree2(count, offset)); 22410b424faSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 22510b424faSBarry Smith } 22610b424faSBarry Smith 22710b424faSBarry Smith static PetscErrorCode PetscViewerHDF5ReadArray_Private(PetscViewer viewer, HDF5ReadCtx h, hid_t datatype, hid_t memspace, void *arr) 22810b424faSBarry Smith { 22910b424faSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 23010b424faSBarry Smith 23110b424faSBarry Smith PetscFunctionBegin; 23210b424faSBarry Smith PetscCallHDF5(H5Dread, (h->dataset, datatype, memspace, h->dataspace, hdf5->dxpl_id, arr)); 23310b424faSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 23410b424faSBarry Smith } 23510b424faSBarry Smith 23621c42226SMatthew G. Knepley static PetscErrorCode PetscViewerHDF5Load_Internal(PetscViewer viewer, const char name[], PetscBool uncompress, PetscLayout map, hid_t datatype, void **newarr) 23721c42226SMatthew G. Knepley { 23821c42226SMatthew G. Knepley PetscBool has; 23921c42226SMatthew G. Knepley char *group; 24021c42226SMatthew G. Knepley HDF5ReadCtx h = NULL; 24121c42226SMatthew G. Knepley hid_t memspace = 0; 24221c42226SMatthew G. Knepley size_t unitsize; 24321c42226SMatthew G. Knepley void *arr; 24421c42226SMatthew G. Knepley 24521c42226SMatthew G. Knepley PetscFunctionBegin; 24621c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group)); 24721c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5HasDataset(viewer, name, &has)); 24821c42226SMatthew G. Knepley PetscCheck(has, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", name, group); 24921c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5ReadInitialize_Private(viewer, name, &h)); 25021c42226SMatthew G. Knepley #if defined(PETSC_USE_COMPLEX) 25121c42226SMatthew G. Knepley if (!h->complexVal) { 25221c42226SMatthew G. Knepley H5T_class_t clazz = H5Tget_class(datatype); 25321c42226SMatthew 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); 25421c42226SMatthew G. Knepley } 25521c42226SMatthew G. Knepley #else 25621c42226SMatthew 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); 25721c42226SMatthew G. Knepley #endif 25821c42226SMatthew G. Knepley 25921c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5ReadSizes_Private(viewer, h, uncompress, PETSC_TRUE, &map)); 26021c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5ReadSelectHyperslab_Private(viewer, h, map, &memspace)); 26121c42226SMatthew G. Knepley 26221c42226SMatthew G. Knepley if (h->runs && uncompress) { 26321c42226SMatthew G. Knepley PetscInt *ind; 26421c42226SMatthew G. Knepley 26521c42226SMatthew G. Knepley PetscCall(PetscInfo(viewer, "Read compressed object with name %s of size %" PetscInt_FMT ":%" PetscInt_FMT "\n", name, map->n, map->N)); 26621c42226SMatthew G. Knepley // Each process stores the whole compression, so skip any leading parts 26721c42226SMatthew G. Knepley PetscCall(PetscMalloc1(map->n, &ind)); 26821c42226SMatthew G. Knepley for (PetscInt i = 0, off = 0; i < h->runs; ++i) { 26921c42226SMatthew G. Knepley for (PetscInt j = 0, inc = 0; j < h->cind[i * 3 + 0]; ++j, ++off, inc += h->cind[i * 3 + 1]) { 27021c42226SMatthew G. Knepley if (off >= map->rend) { 27121c42226SMatthew G. Knepley i = h->runs; 27221c42226SMatthew G. Knepley break; 27321c42226SMatthew G. Knepley } 27421c42226SMatthew G. Knepley if (off >= map->rstart) ind[off - map->rstart] = h->cind[i * 3 + 2] + inc; 27521c42226SMatthew G. Knepley } 27621c42226SMatthew G. Knepley } 27721c42226SMatthew G. Knepley *newarr = ind; 27821c42226SMatthew G. Knepley goto cleanup; 27921c42226SMatthew G. Knepley } 28021c42226SMatthew G. Knepley 28121c42226SMatthew G. Knepley unitsize = H5Tget_size(datatype); 28221c42226SMatthew G. Knepley if (h->complexVal) unitsize *= 2; 28321c42226SMatthew G. Knepley /* unitsize is size_t i.e. always unsigned, so the negative check is pointless? */ 28421c42226SMatthew 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); 28521c42226SMatthew G. Knepley PetscCall(PetscMalloc(map->n * unitsize, &arr)); 28621c42226SMatthew G. Knepley 28721c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5ReadArray_Private(viewer, h, datatype, memspace, arr)); 28821c42226SMatthew G. Knepley *newarr = arr; 28921c42226SMatthew G. Knepley 29021c42226SMatthew G. Knepley cleanup: 29121c42226SMatthew G. Knepley PetscCallHDF5(H5Sclose, (memspace)); 29221c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5ReadFinalize_Private(viewer, &h)); 29321c42226SMatthew G. Knepley PetscCall(PetscFree(group)); 29421c42226SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 29521c42226SMatthew G. Knepley } 29621c42226SMatthew G. Knepley 29710b424faSBarry Smith /*@C 298a24573d6SBarry Smith PetscViewerHDF5Load - Read a raw array from the `PETSCVIEWERHDF5` dataset in parallel 29910b424faSBarry Smith 30010b424faSBarry Smith Collective; No Fortran Support 30110b424faSBarry Smith 30210b424faSBarry Smith Input Parameters: 30310b424faSBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 30410b424faSBarry Smith . name - The dataset name 30510b424faSBarry Smith - datatype - The HDF5 datatype of the items in the dataset 30610b424faSBarry Smith 30710b424faSBarry Smith Input/Output Parameter: 30810b424faSBarry Smith . map - The layout which specifies array partitioning, on output the 30910b424faSBarry Smith set up layout (with global size and blocksize according to dataset) 31010b424faSBarry Smith 31110b424faSBarry Smith Output Parameter: 31210b424faSBarry Smith . newarr - The partitioned array, a memory image of the given dataset 31310b424faSBarry Smith 31410b424faSBarry Smith Level: developer 31510b424faSBarry Smith 31610b424faSBarry Smith Notes: 31710b424faSBarry Smith This is intended mainly for internal use; users should use higher level routines such as `ISLoad()`, `VecLoad()`, `DMLoad()`. 31810b424faSBarry Smith 31910b424faSBarry Smith The array is partitioned according to the given `PetscLayout` which is converted to an HDF5 hyperslab. 32010b424faSBarry Smith 32110b424faSBarry Smith This name is relative to the current group returned by `PetscViewerHDF5OpenGroup()`. 32210b424faSBarry Smith 32310b424faSBarry Smith .seealso: `PetscViewer`, `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5ReadSizes()`, 324a24573d6SBarry Smith `VecLoad()`, `ISLoad()`, `PetscLayout` 32510b424faSBarry Smith @*/ 326cc4c1da9SBarry Smith PetscErrorCode PetscViewerHDF5Load(PetscViewer viewer, const char name[], PetscLayout map, hid_t datatype, void **newarr) 32710b424faSBarry Smith { 32810b424faSBarry Smith PetscFunctionBegin; 32921c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5Load_Internal(viewer, name, PETSC_TRUE, map, datatype, newarr)); 33010b424faSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 33110b424faSBarry Smith } 33210b424faSBarry Smith 333*2e1d0745SJose E. Roman /*@ 33410b424faSBarry Smith PetscViewerHDF5ReadSizes - Read block size and global size of a `Vec` or `IS` stored in an HDF5 file. 33510b424faSBarry Smith 33610b424faSBarry Smith Input Parameters: 33710b424faSBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 33810b424faSBarry Smith - name - The dataset name 33910b424faSBarry Smith 34010b424faSBarry Smith Output Parameters: 34110b424faSBarry Smith + bs - block size 34210b424faSBarry Smith - N - global size 34310b424faSBarry Smith 34410b424faSBarry Smith Level: advanced 34510b424faSBarry Smith 34610b424faSBarry Smith Notes: 34710b424faSBarry Smith The dataset is stored as an HDF5 dataspace with 1-4 dimensions in the order 34810b424faSBarry Smith 1) # timesteps (optional), 2) # blocks, 3) # elements per block (optional), 4) real and imaginary part (only for complex). 34910b424faSBarry Smith 35010b424faSBarry Smith The dataset can be stored as a 2D dataspace even if its blocksize is 1; see `PetscViewerHDF5SetBaseDimension2()`. 35110b424faSBarry Smith 35210b424faSBarry Smith .seealso: `PetscViewer`, `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `VecLoad()`, `ISLoad()`, `VecGetSize()`, `ISGetSize()`, `PetscViewerHDF5SetBaseDimension2()` 35310b424faSBarry Smith @*/ 35410b424faSBarry Smith PetscErrorCode PetscViewerHDF5ReadSizes(PetscViewer viewer, const char name[], PetscInt *bs, PetscInt *N) 35510b424faSBarry Smith { 35610b424faSBarry Smith HDF5ReadCtx h = NULL; 35710b424faSBarry Smith PetscLayout map = NULL; 35810b424faSBarry Smith 35910b424faSBarry Smith PetscFunctionBegin; 36010b424faSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 36110b424faSBarry Smith PetscCall(PetscViewerHDF5ReadInitialize_Private(viewer, name, &h)); 36221c42226SMatthew G. Knepley PetscCall(PetscViewerHDF5ReadSizes_Private(viewer, h, PETSC_TRUE, PETSC_FALSE, &map)); 36310b424faSBarry Smith PetscCall(PetscViewerHDF5ReadFinalize_Private(viewer, &h)); 36410b424faSBarry Smith if (bs) *bs = map->bs; 36510b424faSBarry Smith if (N) *N = map->N; 36610b424faSBarry Smith PetscCall(PetscLayoutDestroy(&map)); 36710b424faSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 36810b424faSBarry Smith } 369