PetscErrorCode ISView_General_HDF5(IS is, PetscViewer viewer) { hid_t filespace; /* file dataspace identifier */ hid_t chunkspace; /* chunk dataset property identifier */ hid_t plist_id; /* property list identifier */ hid_t dset_id; /* dataset identifier */ hid_t memspace; /* memory dataspace identifier */ hid_t inttype; /* int type (H5T_NATIVE_INT or H5T_NATIVE_LLONG) */ hid_t file_id, group; herr_t status; hsize_t dim, maxDims[3], dims[3], chunkDims[3], count[3],offset[3]; PetscInt bs, N, n, timestep, low; const PetscInt *ind; const char *isname; PetscErrorCode ierr; PetscFunctionBegin; ierr = ISGetBlockSize(is,&bs);CHKERRQ(ierr); ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr); ierr = PetscViewerHDF5GetTimestep(viewer, ×tep);CHKERRQ(ierr); /* Create the dataspace for the dataset. * * dims - holds the current dimensions of the dataset * * maxDims - holds the maximum dimensions of the dataset (unlimited * for the number of time steps with the current dimensions for the * other dimensions; so only additional time steps can be added). * * chunkDims - holds the size of a single time step (required to * permit extending dataset). */ dim = 0; if (timestep >= 0) { dims[dim] = timestep+1; maxDims[dim] = H5S_UNLIMITED; chunkDims[dim] = 1; ++dim; } ierr = ISGetSize(is, &N);CHKERRQ(ierr); ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr); ierr = PetscHDF5IntCast(N/bs,dims + dim);CHKERRQ(ierr); maxDims[dim] = dims[dim]; chunkDims[dim] = dims[dim]; ++dim; if (bs >= 1) { dims[dim] = bs; maxDims[dim] = dims[dim]; chunkDims[dim] = dims[dim]; ++dim; } filespace = H5Screate_simple(dim, dims, maxDims); if (filespace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Cannot H5Screate_simple()"); #if defined(PETSC_USE_64BIT_INDICES) inttype = H5T_NATIVE_LLONG; #else inttype = H5T_NATIVE_INT; #endif /* Create the dataset with default properties and close filespace */ ierr = PetscObjectGetName((PetscObject) is, &isname);CHKERRQ(ierr); if (!H5Lexists(group, isname, H5P_DEFAULT)) { /* Create chunk */ chunkspace = H5Pcreate(H5P_DATASET_CREATE); if (chunkspace == -1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Cannot H5Pcreate()"); status = H5Pset_chunk(chunkspace, dim, chunkDims);CHKERRQ(status); #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) dset_id = H5Dcreate2(group, isname, inttype, filespace, H5P_DEFAULT, chunkspace, H5P_DEFAULT); #else dset_id = H5Dcreate(group, isname, inttype, filespace, H5P_DEFAULT); #endif if (dset_id == -1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Cannot H5Dcreate2()"); status = H5Pclose(chunkspace);CHKERRQ(status); } else { dset_id = H5Dopen2(group, isname, H5P_DEFAULT); status = H5Dset_extent(dset_id, dims);CHKERRQ(status); } status = H5Sclose(filespace);CHKERRQ(status); /* Each process defines a dataset and writes it to the hyperslab in the file */ dim = 0; if (timestep >= 0) { count[dim] = 1; ++dim; } ierr = PetscHDF5IntCast(n/bs,count + dim);CHKERRQ(ierr); ++dim; if (bs >= 1) { count[dim] = bs; ++dim; } if (n > 0) { memspace = H5Screate_simple(dim, count, NULL); if (memspace == -1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Cannot H5Screate_simple()"); } else { /* Can't create dataspace with zero for any dimension, so create null dataspace. */ memspace = H5Screate(H5S_NULL); if (memspace == -1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Cannot H5Screate()"); } /* Select hyperslab in the file */ ierr = PetscLayoutGetRange(is->map, &low, NULL);CHKERRQ(ierr); dim = 0; if (timestep >= 0) { offset[dim] = timestep; ++dim; } ierr = PetscHDF5IntCast(low/bs,offset + dim);CHKERRQ(ierr); ++dim; if (bs >= 1) { offset[dim] = 0; ++dim; } if (n > 0) { filespace = H5Dget_space(dset_id); if (filespace == -1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Cannot H5Dget_space()"); status = H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);CHKERRQ(status); } else { /* Create null filespace to match null memspace. */ filespace = H5Screate(H5S_NULL); if (filespace == -1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Cannot H5Screate(H5S_NULL)"); } /* Create property list for collective dataset write */ plist_id = H5Pcreate(H5P_DATASET_XFER); if (plist_id == -1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Cannot H5Pcreate()"); #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);CHKERRQ(status); #endif /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ ierr = ISGetIndices(is, &ind);CHKERRQ(ierr); status = H5Dwrite(dset_id, inttype, memspace, filespace, plist_id, ind);CHKERRQ(status); status = H5Fflush(file_id, H5F_SCOPE_GLOBAL);CHKERRQ(status); ierr = ISGetIndices(is, &ind);CHKERRQ(ierr); /* Close/release resources */ if (group != file_id) {status = H5Gclose(group);CHKERRQ(status);} status = H5Pclose(plist_id);CHKERRQ(status); status = H5Sclose(filespace);CHKERRQ(status); status = H5Sclose(memspace);CHKERRQ(status); status = H5Dclose(dset_id);CHKERRQ(status); ierr = PetscInfo1(is, "Wrote IS object with name %s\n", isname);CHKERRQ(ierr); PetscFunctionReturn(0); }
/* This should handle properly the cases where PetscInt is 32 or 64 and hsize_t is 32 or 64. These means properly casting with checks back and forth between the two types of variables. */ PetscErrorCode VecLoad_HDF5(Vec xin, PetscViewer viewer) { hid_t file_id, group, dset_id, filespace, memspace, plist_id; hsize_t rdim, dim; hsize_t dims[4], count[4], offset[4]; herr_t status; PetscInt n, N, bs = 1, bsInd, lenInd, low, timestep; PetscScalar *x; const char *vecname; PetscErrorCode ierr; PetscFunctionBegin; ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr); ierr = PetscViewerHDF5GetTimestep(viewer, ×tep);CHKERRQ(ierr); ierr = VecGetBlockSize(xin,&bs);CHKERRQ(ierr); /* Create the dataset with default properties and close filespace */ ierr = PetscObjectGetName((PetscObject)xin,&vecname);CHKERRQ(ierr); #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) dset_id = H5Dopen2(group, vecname, H5P_DEFAULT); #else dset_id = H5Dopen(group, vecname); #endif if (dset_id == -1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Could not H5Dopen() with Vec named %s",vecname); /* Retrieve the dataspace for the dataset */ filespace = H5Dget_space(dset_id); if (filespace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Could not H5Dget_space()"); dim = 0; if (timestep >= 0) ++dim; ++dim; if (bs >= 1) ++dim; #if defined(PETSC_USE_COMPLEX) ++dim; #endif rdim = H5Sget_simple_extent_dims(filespace, dims, NULL); #if defined(PETSC_USE_COMPLEX) bsInd = rdim-2; #else bsInd = rdim-1; #endif lenInd = timestep >= 0 ? 1 : 0; if (rdim != dim) { if (rdim == dim+1 && bs == -1) bs = dims[bsInd]; else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Dimension of array in file %d not %d as expected",rdim,dim); } else if (bs >= 1 && bs != (PetscInt) dims[bsInd]) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Block size %d specified for vector does not match blocksize in file %d",bs,dims[bsInd]); /* Set Vec sizes,blocksize,and type if not already set */ if ((xin)->map->n < 0 && (xin)->map->N < 0) { ierr = VecSetSizes(xin, PETSC_DECIDE, dims[lenInd]*bs);CHKERRQ(ierr); } /* If sizes and type already set,check if the vector global size is correct */ ierr = VecGetSize(xin, &N);CHKERRQ(ierr); if (N/bs != (PetscInt) dims[lenInd]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Vector in file different length (%d) then input vector (%d)", (PetscInt) dims[lenInd], N/bs); /* Each process defines a dataset and reads it from the hyperslab in the file */ ierr = VecGetLocalSize(xin, &n);CHKERRQ(ierr); dim = 0; if (timestep >= 0) { count[dim] = 1; ++dim; } ierr = PetscHDF5IntCast(n/bs,count + dim);CHKERRQ(ierr); ++dim; if (bs >= 1) { count[dim] = bs; ++dim; } #if defined(PETSC_USE_COMPLEX) count[dim] = 2; ++dim; #endif memspace = H5Screate_simple(dim, count, NULL); if (memspace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Could not H5Screate_simple()"); /* Select hyperslab in the file */ ierr = VecGetOwnershipRange(xin, &low, NULL);CHKERRQ(ierr); dim = 0; if (timestep >= 0) { offset[dim] = timestep; ++dim; } ierr = PetscHDF5IntCast(low/bs,offset + dim);CHKERRQ(ierr); ++dim; if (bs >= 1) { offset[dim] = 0; ++dim; } #if defined(PETSC_USE_COMPLEX) offset[dim] = 0; ++dim; #endif status = H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);CHKERRQ(status); /* Create property list for collective dataset read */ plist_id = H5Pcreate(H5P_DATASET_XFER); if (plist_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Could not H5Pcreate()"); #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);CHKERRQ(status); #endif /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ ierr = VecGetArray(xin, &x);CHKERRQ(ierr); status = H5Dread(dset_id, H5T_NATIVE_DOUBLE, memspace, filespace, plist_id, x);CHKERRQ(status); ierr = VecRestoreArray(xin, &x);CHKERRQ(ierr); /* Close/release resources */ if (group != file_id) { status = H5Gclose(group);CHKERRQ(status); } status = H5Pclose(plist_id);CHKERRQ(status); status = H5Sclose(filespace);CHKERRQ(status); status = H5Sclose(memspace);CHKERRQ(status); status = H5Dclose(dset_id);CHKERRQ(status); ierr = VecAssemblyBegin(xin);CHKERRQ(ierr); ierr = VecAssemblyEnd(xin);CHKERRQ(ierr); PetscFunctionReturn(0); }