xref: /petsc/src/vec/is/utils/hdf5/hdf5io.c (revision 9ad2ceda86f92eef2ad3194ad31c6a1b35bb29cb)
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 {
810b424faSBarry Smith   hid_t     file, group, dataset, dataspace;
910b424faSBarry Smith   int       lenInd, bsInd, complexInd, rdim;
1010b424faSBarry Smith   hsize_t  *dims;
1110b424faSBarry Smith   PetscBool complexVal, dim2;
1210b424faSBarry Smith };
1310b424faSBarry Smith typedef struct _n_HDF5ReadCtx *HDF5ReadCtx;
1410b424faSBarry Smith 
1510b424faSBarry Smith PetscErrorCode PetscViewerHDF5CheckTimestepping_Internal(PetscViewer viewer, const char name[])
1610b424faSBarry Smith {
1710b424faSBarry Smith   PetscViewer_HDF5 *hdf5         = (PetscViewer_HDF5 *)viewer->data;
1810b424faSBarry Smith   PetscBool         timestepping = PETSC_FALSE;
1910b424faSBarry Smith 
2010b424faSBarry Smith   PetscFunctionBegin;
2110b424faSBarry Smith   PetscCall(PetscViewerHDF5ReadAttribute(viewer, name, "timestepping", PETSC_BOOL, &hdf5->defTimestepping, &timestepping));
2210b424faSBarry Smith   if (timestepping != hdf5->timestepping) {
2310b424faSBarry Smith     char *group;
2410b424faSBarry Smith 
2510b424faSBarry Smith     PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group));
2610b424faSBarry 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]);
2710b424faSBarry Smith   }
2810b424faSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
2910b424faSBarry Smith }
3010b424faSBarry Smith 
3110b424faSBarry Smith static PetscErrorCode PetscViewerHDF5ReadInitialize_Private(PetscViewer viewer, const char name[], HDF5ReadCtx *ctx)
3210b424faSBarry Smith {
3310b424faSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
3410b424faSBarry Smith   HDF5ReadCtx       h    = NULL;
3510b424faSBarry Smith 
3610b424faSBarry Smith   PetscFunctionBegin;
3710b424faSBarry Smith   PetscCall(PetscViewerHDF5CheckTimestepping_Internal(viewer, name));
3810b424faSBarry Smith   PetscCall(PetscNew(&h));
3910b424faSBarry Smith   PetscCall(PetscViewerHDF5OpenGroup(viewer, NULL, &h->file, &h->group));
4010b424faSBarry Smith   PetscCallHDF5Return(h->dataset, H5Dopen2, (h->group, name, H5P_DEFAULT));
4110b424faSBarry Smith   PetscCallHDF5Return(h->dataspace, H5Dget_space, (h->dataset));
4210b424faSBarry Smith   PetscCall(PetscViewerHDF5ReadAttribute(viewer, name, "complex", PETSC_BOOL, &h->complexVal, &h->complexVal));
4310b424faSBarry Smith   if (!hdf5->horizontal) {
4410b424faSBarry Smith     /* MATLAB stores column vectors horizontally */
4510b424faSBarry Smith     PetscCall(PetscViewerHDF5HasAttribute(viewer, name, "MATLAB_class", &hdf5->horizontal));
4610b424faSBarry Smith   }
4710b424faSBarry Smith   *ctx = h;
4810b424faSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
4910b424faSBarry Smith }
5010b424faSBarry Smith 
5110b424faSBarry Smith static PetscErrorCode PetscViewerHDF5ReadFinalize_Private(PetscViewer viewer, HDF5ReadCtx *ctx)
5210b424faSBarry Smith {
5310b424faSBarry Smith   HDF5ReadCtx h;
5410b424faSBarry Smith 
5510b424faSBarry Smith   PetscFunctionBegin;
5610b424faSBarry Smith   h = *ctx;
5710b424faSBarry Smith   PetscCallHDF5(H5Gclose, (h->group));
5810b424faSBarry Smith   PetscCallHDF5(H5Sclose, (h->dataspace));
5910b424faSBarry Smith   PetscCallHDF5(H5Dclose, (h->dataset));
6010b424faSBarry Smith   PetscCall(PetscFree((*ctx)->dims));
6110b424faSBarry Smith   PetscCall(PetscFree(*ctx));
6210b424faSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
6310b424faSBarry Smith }
6410b424faSBarry Smith 
6510b424faSBarry Smith static PetscErrorCode PetscViewerHDF5ReadSizes_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscBool setup, PetscLayout *map_)
6610b424faSBarry Smith {
6710b424faSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
686497c311SBarry Smith   PetscInt          bs, N;
6910b424faSBarry Smith   PetscLayout       map;
7010b424faSBarry Smith 
7110b424faSBarry Smith   PetscFunctionBegin;
724ad8454bSPierre Jolivet   if (!*map_) PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)viewer), map_));
7310b424faSBarry Smith   map = *map_;
7410b424faSBarry Smith 
7510b424faSBarry Smith   /* Get actual number of dimensions in dataset */
7610b424faSBarry Smith   PetscCallHDF5Return(ctx->rdim, H5Sget_simple_extent_dims, (ctx->dataspace, NULL, NULL));
7710b424faSBarry Smith   PetscCall(PetscMalloc1(ctx->rdim, &ctx->dims));
7810b424faSBarry Smith   PetscCallHDF5Return(ctx->rdim, H5Sget_simple_extent_dims, (ctx->dataspace, ctx->dims, NULL));
7910b424faSBarry Smith 
8010b424faSBarry Smith   /*
8110b424faSBarry Smith      Dimensions are in this order:
8210b424faSBarry Smith      [0]        timesteps (optional)
8310b424faSBarry Smith      [lenInd]   entries (numbers or blocks)
8410b424faSBarry Smith      ...
8510b424faSBarry Smith      [bsInd]    entries of blocks (optional)
8610b424faSBarry Smith      [bsInd+1]  real & imaginary part (optional)
8710b424faSBarry Smith       = rdim-1
8810b424faSBarry Smith    */
8910b424faSBarry Smith 
9010b424faSBarry Smith   /* Get entries dimension index */
9110b424faSBarry Smith   ctx->lenInd = 0;
9210b424faSBarry Smith   if (hdf5->timestepping) ++ctx->lenInd;
9310b424faSBarry Smith 
9410b424faSBarry Smith   /* Get block dimension index */
9510b424faSBarry Smith   if (ctx->complexVal) {
9610b424faSBarry Smith     ctx->bsInd      = ctx->rdim - 2;
9710b424faSBarry Smith     ctx->complexInd = ctx->rdim - 1;
9810b424faSBarry Smith   } else {
9910b424faSBarry Smith     ctx->bsInd      = ctx->rdim - 1;
10010b424faSBarry Smith     ctx->complexInd = -1;
10110b424faSBarry Smith   }
10210b424faSBarry 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);
10310b424faSBarry 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);
10410b424faSBarry 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]);
10510b424faSBarry Smith 
10610b424faSBarry Smith   if (hdf5->horizontal) {
10710b424faSBarry Smith     /* support horizontal 1D arrays (MATLAB vectors) - swap meaning of blocks and entries */
108*9ad2cedaSPierre Jolivet     int t       = ctx->lenInd;
10910b424faSBarry Smith     ctx->lenInd = ctx->bsInd;
11010b424faSBarry Smith     ctx->bsInd  = t;
11110b424faSBarry Smith   }
11210b424faSBarry Smith 
11310b424faSBarry Smith   /* Get block size */
11410b424faSBarry Smith   ctx->dim2 = PETSC_FALSE;
11510b424faSBarry Smith   if (ctx->lenInd == ctx->bsInd) {
11610b424faSBarry Smith     bs = 1; /* support vectors stored as 1D array */
11710b424faSBarry Smith   } else {
11810b424faSBarry Smith     bs = (PetscInt)ctx->dims[ctx->bsInd];
11910b424faSBarry Smith     if (bs == 1) ctx->dim2 = PETSC_TRUE; /* vector with blocksize of 1, still stored as 2D array */
12010b424faSBarry Smith   }
12110b424faSBarry Smith 
12210b424faSBarry Smith   /* Get global size */
1236497c311SBarry Smith   PetscCall(PetscIntCast(bs * ctx->dims[ctx->lenInd], &N));
12410b424faSBarry Smith 
12510b424faSBarry Smith   /* Set global size, blocksize and type if not yet set */
12610b424faSBarry Smith   if (map->bs < 0) {
12710b424faSBarry Smith     PetscCall(PetscLayoutSetBlockSize(map, bs));
12810b424faSBarry 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);
12910b424faSBarry Smith   if (map->N < 0) {
13010b424faSBarry Smith     PetscCall(PetscLayoutSetSize(map, N));
13110b424faSBarry Smith   } else PetscCheck(map->N == N, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Global size of array in file is %" PetscInt_FMT ", not %" PetscInt_FMT " as expected", N, map->N);
13210b424faSBarry Smith   if (setup) PetscCall(PetscLayoutSetUp(map));
13310b424faSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
13410b424faSBarry Smith }
13510b424faSBarry Smith 
13610b424faSBarry Smith static PetscErrorCode PetscViewerHDF5ReadSelectHyperslab_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout map, hid_t *memspace)
13710b424faSBarry Smith {
13810b424faSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
13910b424faSBarry Smith   hsize_t          *count, *offset;
14010b424faSBarry Smith   PetscInt          bs, n, low;
14110b424faSBarry Smith   int               i;
14210b424faSBarry Smith 
14310b424faSBarry Smith   PetscFunctionBegin;
14410b424faSBarry Smith   /* Compute local size and ownership range */
14510b424faSBarry Smith   PetscCall(PetscLayoutSetUp(map));
14610b424faSBarry Smith   PetscCall(PetscLayoutGetBlockSize(map, &bs));
14710b424faSBarry Smith   PetscCall(PetscLayoutGetLocalSize(map, &n));
14810b424faSBarry Smith   PetscCall(PetscLayoutGetRange(map, &low, NULL));
14910b424faSBarry Smith 
15010b424faSBarry Smith   /* Each process defines a dataset and reads it from the hyperslab in the file */
15110b424faSBarry Smith   PetscCall(PetscMalloc2(ctx->rdim, &count, ctx->rdim, &offset));
15210b424faSBarry Smith   for (i = 0; i < ctx->rdim; i++) {
15310b424faSBarry Smith     /* By default, select all entries with no offset */
15410b424faSBarry Smith     offset[i] = 0;
15510b424faSBarry Smith     count[i]  = ctx->dims[i];
15610b424faSBarry Smith   }
15710b424faSBarry Smith   if (hdf5->timestepping) {
15810b424faSBarry Smith     count[0]  = 1;
15910b424faSBarry Smith     offset[0] = hdf5->timestep;
16010b424faSBarry Smith   }
16110b424faSBarry Smith   {
16210b424faSBarry Smith     PetscCall(PetscHDF5IntCast(n / bs, &count[ctx->lenInd]));
16310b424faSBarry Smith     PetscCall(PetscHDF5IntCast(low / bs, &offset[ctx->lenInd]));
16410b424faSBarry Smith   }
16510b424faSBarry Smith   PetscCallHDF5Return(*memspace, H5Screate_simple, (ctx->rdim, count, NULL));
16610b424faSBarry Smith   PetscCallHDF5(H5Sselect_hyperslab, (ctx->dataspace, H5S_SELECT_SET, offset, NULL, count, NULL));
16710b424faSBarry Smith   PetscCall(PetscFree2(count, offset));
16810b424faSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
16910b424faSBarry Smith }
17010b424faSBarry Smith 
17110b424faSBarry Smith static PetscErrorCode PetscViewerHDF5ReadArray_Private(PetscViewer viewer, HDF5ReadCtx h, hid_t datatype, hid_t memspace, void *arr)
17210b424faSBarry Smith {
17310b424faSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
17410b424faSBarry Smith 
17510b424faSBarry Smith   PetscFunctionBegin;
17610b424faSBarry Smith   PetscCallHDF5(H5Dread, (h->dataset, datatype, memspace, h->dataspace, hdf5->dxpl_id, arr));
17710b424faSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
17810b424faSBarry Smith }
17910b424faSBarry Smith 
18010b424faSBarry Smith /*@C
181a24573d6SBarry Smith   PetscViewerHDF5Load - Read a raw array from the `PETSCVIEWERHDF5` dataset in parallel
18210b424faSBarry Smith 
18310b424faSBarry Smith   Collective; No Fortran Support
18410b424faSBarry Smith 
18510b424faSBarry Smith   Input Parameters:
18610b424faSBarry Smith + viewer   - The `PETSCVIEWERHDF5` viewer
18710b424faSBarry Smith . name     - The dataset name
18810b424faSBarry Smith - datatype - The HDF5 datatype of the items in the dataset
18910b424faSBarry Smith 
19010b424faSBarry Smith   Input/Output Parameter:
19110b424faSBarry Smith . map - The layout which specifies array partitioning, on output the
19210b424faSBarry Smith              set up layout (with global size and blocksize according to dataset)
19310b424faSBarry Smith 
19410b424faSBarry Smith   Output Parameter:
19510b424faSBarry Smith . newarr - The partitioned array, a memory image of the given dataset
19610b424faSBarry Smith 
19710b424faSBarry Smith   Level: developer
19810b424faSBarry Smith 
19910b424faSBarry Smith   Notes:
20010b424faSBarry Smith   This is intended mainly for internal use; users should use higher level routines such as `ISLoad()`, `VecLoad()`, `DMLoad()`.
20110b424faSBarry Smith 
20210b424faSBarry Smith   The array is partitioned according to the given `PetscLayout` which is converted to an HDF5 hyperslab.
20310b424faSBarry Smith 
20410b424faSBarry Smith   This name is relative to the current group returned by `PetscViewerHDF5OpenGroup()`.
20510b424faSBarry Smith 
20610b424faSBarry Smith .seealso: `PetscViewer`, `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5ReadSizes()`,
207a24573d6SBarry Smith           `VecLoad()`, `ISLoad()`, `PetscLayout`
20810b424faSBarry Smith @*/
209cc4c1da9SBarry Smith PetscErrorCode PetscViewerHDF5Load(PetscViewer viewer, const char name[], PetscLayout map, hid_t datatype, void **newarr)
21010b424faSBarry Smith {
21110b424faSBarry Smith   PetscBool   has;
21210b424faSBarry Smith   char       *group;
21310b424faSBarry Smith   HDF5ReadCtx h        = NULL;
21410b424faSBarry Smith   hid_t       memspace = 0;
21510b424faSBarry Smith   size_t      unitsize;
21610b424faSBarry Smith   void       *arr;
21710b424faSBarry Smith 
21810b424faSBarry Smith   PetscFunctionBegin;
21910b424faSBarry Smith   PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group));
22010b424faSBarry Smith   PetscCall(PetscViewerHDF5HasDataset(viewer, name, &has));
22110b424faSBarry Smith   PetscCheck(has, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", name, group);
22210b424faSBarry Smith   PetscCall(PetscViewerHDF5ReadInitialize_Private(viewer, name, &h));
22310b424faSBarry Smith   #if defined(PETSC_USE_COMPLEX)
22410b424faSBarry Smith   if (!h->complexVal) {
22510b424faSBarry Smith     H5T_class_t clazz = H5Tget_class(datatype);
22610b424faSBarry Smith     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);
22710b424faSBarry Smith   }
22810b424faSBarry Smith   #else
22910b424faSBarry Smith   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);
23010b424faSBarry Smith   #endif
23110b424faSBarry Smith 
23210b424faSBarry Smith   PetscCall(PetscViewerHDF5ReadSizes_Private(viewer, h, PETSC_TRUE, &map));
23310b424faSBarry Smith   PetscCall(PetscViewerHDF5ReadSelectHyperslab_Private(viewer, h, map, &memspace));
23410b424faSBarry Smith 
23510b424faSBarry Smith   unitsize = H5Tget_size(datatype);
23610b424faSBarry Smith   if (h->complexVal) unitsize *= 2;
23710b424faSBarry Smith   /* unitsize is size_t i.e. always unsigned, so the negative check is pointless? */
23810b424faSBarry Smith   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);
23910b424faSBarry Smith   PetscCall(PetscMalloc(map->n * unitsize, &arr));
24010b424faSBarry Smith 
24110b424faSBarry Smith   PetscCall(PetscViewerHDF5ReadArray_Private(viewer, h, datatype, memspace, arr));
24210b424faSBarry Smith   PetscCallHDF5(H5Sclose, (memspace));
24310b424faSBarry Smith   PetscCall(PetscViewerHDF5ReadFinalize_Private(viewer, &h));
24410b424faSBarry Smith   PetscCall(PetscFree(group));
24510b424faSBarry Smith   *newarr = arr;
24610b424faSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
24710b424faSBarry Smith }
24810b424faSBarry Smith 
24910b424faSBarry Smith /*@C
25010b424faSBarry Smith   PetscViewerHDF5ReadSizes - Read block size and global size of a `Vec` or `IS` stored in an HDF5 file.
25110b424faSBarry Smith 
25210b424faSBarry Smith   Input Parameters:
25310b424faSBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
25410b424faSBarry Smith - name   - The dataset name
25510b424faSBarry Smith 
25610b424faSBarry Smith   Output Parameters:
25710b424faSBarry Smith + bs - block size
25810b424faSBarry Smith - N  - global size
25910b424faSBarry Smith 
26010b424faSBarry Smith   Level: advanced
26110b424faSBarry Smith 
26210b424faSBarry Smith   Notes:
26310b424faSBarry Smith   The dataset is stored as an HDF5 dataspace with 1-4 dimensions in the order
26410b424faSBarry Smith   1) # timesteps (optional), 2) # blocks, 3) # elements per block (optional), 4) real and imaginary part (only for complex).
26510b424faSBarry Smith 
26610b424faSBarry Smith   The dataset can be stored as a 2D dataspace even if its blocksize is 1; see `PetscViewerHDF5SetBaseDimension2()`.
26710b424faSBarry Smith 
26810b424faSBarry Smith .seealso: `PetscViewer`, `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `VecLoad()`, `ISLoad()`, `VecGetSize()`, `ISGetSize()`, `PetscViewerHDF5SetBaseDimension2()`
26910b424faSBarry Smith @*/
27010b424faSBarry Smith PetscErrorCode PetscViewerHDF5ReadSizes(PetscViewer viewer, const char name[], PetscInt *bs, PetscInt *N)
27110b424faSBarry Smith {
27210b424faSBarry Smith   HDF5ReadCtx h   = NULL;
27310b424faSBarry Smith   PetscLayout map = NULL;
27410b424faSBarry Smith 
27510b424faSBarry Smith   PetscFunctionBegin;
27610b424faSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
27710b424faSBarry Smith   PetscCall(PetscViewerHDF5ReadInitialize_Private(viewer, name, &h));
27810b424faSBarry Smith   PetscCall(PetscViewerHDF5ReadSizes_Private(viewer, h, PETSC_FALSE, &map));
27910b424faSBarry Smith   PetscCall(PetscViewerHDF5ReadFinalize_Private(viewer, &h));
28010b424faSBarry Smith   if (bs) *bs = map->bs;
28110b424faSBarry Smith   if (N) *N = map->N;
28210b424faSBarry Smith   PetscCall(PetscLayoutDestroy(&map));
28310b424faSBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
28410b424faSBarry Smith }
28510b424faSBarry Smith 
28610b424faSBarry Smith #endif /* defined(PETSC_HAVE_HDF5) */
287