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