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