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