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