1*235f7792SMatthew G. Knepley #include <petscis.h> /*I "petscis.h" I*/ 2*235f7792SMatthew G. Knepley #include <petsc-private/isimpl.h> 3*235f7792SMatthew G. Knepley #include <petscviewerhdf5.h> 4*235f7792SMatthew G. Knepley 5*235f7792SMatthew G. Knepley #if defined(PETSC_HAVE_HDF5) 6*235f7792SMatthew G. Knepley #undef __FUNCT__ 7*235f7792SMatthew G. Knepley #define __FUNCT__ "ISLoad_HDF5" 8*235f7792SMatthew G. Knepley /* 9*235f7792SMatthew G. Knepley This should handle properly the cases where PetscInt is 32 or 64 and hsize_t is 32 or 64. These means properly casting with 10*235f7792SMatthew G. Knepley checks back and forth between the two types of variables. 11*235f7792SMatthew G. Knepley */ 12*235f7792SMatthew G. Knepley PetscErrorCode ISLoad_HDF5(IS is, PetscViewer viewer) 13*235f7792SMatthew G. Knepley { 14*235f7792SMatthew G. Knepley hid_t inttype; /* int type (H5T_NATIVE_INT or H5T_NATIVE_LLONG) */ 15*235f7792SMatthew G. Knepley hid_t file_id, group, dset_id, filespace, memspace, plist_id; 16*235f7792SMatthew G. Knepley hsize_t rdim, dim; 17*235f7792SMatthew G. Knepley hsize_t dims[3], count[3], offset[3]; 18*235f7792SMatthew G. Knepley herr_t status; 19*235f7792SMatthew G. Knepley PetscInt n, N, bs = 1, bsInd, lenInd, low, timestep; 20*235f7792SMatthew G. Knepley const PetscInt *ind; 21*235f7792SMatthew G. Knepley const char *isname; 22*235f7792SMatthew G. Knepley PetscErrorCode ierr; 23*235f7792SMatthew G. Knepley 24*235f7792SMatthew G. Knepley PetscFunctionBegin; 25*235f7792SMatthew G. Knepley ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr); 26*235f7792SMatthew G. Knepley ierr = PetscViewerHDF5GetTimestep(viewer, ×tep);CHKERRQ(ierr); 27*235f7792SMatthew G. Knepley ierr = ISGetBlockSize(is, &bs);CHKERRQ(ierr); 28*235f7792SMatthew G. Knepley /* Create the dataset with default properties and close filespace */ 29*235f7792SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) is, &isname);CHKERRQ(ierr); 30*235f7792SMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 31*235f7792SMatthew G. Knepley dset_id = H5Dopen2(group, isname, H5P_DEFAULT); 32*235f7792SMatthew G. Knepley #else 33*235f7792SMatthew G. Knepley dset_id = H5Dopen(group, isname); 34*235f7792SMatthew G. Knepley #endif 35*235f7792SMatthew G. Knepley if (dset_id == -1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Could not H5Dopen() with IS named %s", isname); 36*235f7792SMatthew G. Knepley /* Retrieve the dataspace for the dataset */ 37*235f7792SMatthew G. Knepley filespace = H5Dget_space(dset_id); 38*235f7792SMatthew G. Knepley if (filespace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Could not H5Dget_space()"); 39*235f7792SMatthew G. Knepley dim = 0; 40*235f7792SMatthew G. Knepley if (timestep >= 0) ++dim; 41*235f7792SMatthew G. Knepley ++dim; 42*235f7792SMatthew G. Knepley if (bs >= 1) ++dim; 43*235f7792SMatthew G. Knepley rdim = H5Sget_simple_extent_dims(filespace, dims, NULL); 44*235f7792SMatthew G. Knepley bsInd = rdim-1; 45*235f7792SMatthew G. Knepley lenInd = timestep >= 0 ? 1 : 0; 46*235f7792SMatthew G. Knepley if (rdim != dim) { 47*235f7792SMatthew G. Knepley if (rdim == dim+1 && bs == -1) bs = dims[bsInd]; 48*235f7792SMatthew G. Knepley else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Dimension of array in file %d not %d as expected",rdim,dim); 49*235f7792SMatthew G. Knepley } else if (bs >= 1 && bs != (PetscInt) dims[bsInd]) { 50*235f7792SMatthew G. Knepley ierr = ISSetBlockSize(is, dims[bsInd]);CHKERRQ(ierr); 51*235f7792SMatthew G. Knepley if (ierr) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Block size %d specified for IS does not match blocksize in file %d",bs,dims[bsInd]); 52*235f7792SMatthew G. Knepley bs = dims[bsInd]; 53*235f7792SMatthew G. Knepley } 54*235f7792SMatthew G. Knepley 55*235f7792SMatthew G. Knepley /* Set Vec sizes,blocksize,and type if not already set */ 56*235f7792SMatthew G. Knepley ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr); 57*235f7792SMatthew G. Knepley ierr = ISGetSize(is, &N);CHKERRQ(ierr); 58*235f7792SMatthew G. Knepley if (n < 0 && N < 0) {ierr = PetscLayoutSetSize(is->map, dims[lenInd]*bs);CHKERRQ(ierr);} 59*235f7792SMatthew G. Knepley ierr = PetscLayoutSetUp(is->map);CHKERRQ(ierr); 60*235f7792SMatthew G. Knepley /* If sizes and type already set,check if the vector global size is correct */ 61*235f7792SMatthew G. Knepley ierr = ISGetSize(is, &N);CHKERRQ(ierr); 62*235f7792SMatthew G. Knepley if (N/bs != (PetscInt) dims[lenInd]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "IS in file different length (%d) then input vector (%d)", (PetscInt) dims[lenInd], N/bs); 63*235f7792SMatthew G. Knepley 64*235f7792SMatthew G. Knepley /* Each process defines a dataset and reads it from the hyperslab in the file */ 65*235f7792SMatthew G. Knepley ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr); 66*235f7792SMatthew G. Knepley dim = 0; 67*235f7792SMatthew G. Knepley if (timestep >= 0) { 68*235f7792SMatthew G. Knepley count[dim] = 1; 69*235f7792SMatthew G. Knepley ++dim; 70*235f7792SMatthew G. Knepley } 71*235f7792SMatthew G. Knepley ierr = PetscHDF5IntCast(n/bs,count + dim);CHKERRQ(ierr); 72*235f7792SMatthew G. Knepley ++dim; 73*235f7792SMatthew G. Knepley if (bs >= 1) { 74*235f7792SMatthew G. Knepley count[dim] = bs; 75*235f7792SMatthew G. Knepley ++dim; 76*235f7792SMatthew G. Knepley } 77*235f7792SMatthew G. Knepley memspace = H5Screate_simple(dim, count, NULL); 78*235f7792SMatthew G. Knepley if (memspace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Could not H5Screate_simple()"); 79*235f7792SMatthew G. Knepley 80*235f7792SMatthew G. Knepley /* Select hyperslab in the file */ 81*235f7792SMatthew G. Knepley ierr = PetscLayoutGetRange(is->map, &low, NULL);CHKERRQ(ierr); 82*235f7792SMatthew G. Knepley dim = 0; 83*235f7792SMatthew G. Knepley if (timestep >= 0) { 84*235f7792SMatthew G. Knepley offset[dim] = timestep; 85*235f7792SMatthew G. Knepley ++dim; 86*235f7792SMatthew G. Knepley } 87*235f7792SMatthew G. Knepley ierr = PetscHDF5IntCast(low/bs,offset + dim);CHKERRQ(ierr); 88*235f7792SMatthew G. Knepley ++dim; 89*235f7792SMatthew G. Knepley if (bs >= 1) { 90*235f7792SMatthew G. Knepley offset[dim] = 0; 91*235f7792SMatthew G. Knepley ++dim; 92*235f7792SMatthew G. Knepley } 93*235f7792SMatthew G. Knepley status = H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);CHKERRQ(status); 94*235f7792SMatthew G. Knepley 95*235f7792SMatthew G. Knepley /* Create property list for collective dataset read */ 96*235f7792SMatthew G. Knepley plist_id = H5Pcreate(H5P_DATASET_XFER); 97*235f7792SMatthew G. Knepley if (plist_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Could not H5Pcreate()"); 98*235f7792SMatthew G. Knepley #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 99*235f7792SMatthew G. Knepley status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);CHKERRQ(status); 100*235f7792SMatthew G. Knepley #endif 101*235f7792SMatthew G. Knepley /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 102*235f7792SMatthew G. Knepley 103*235f7792SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 104*235f7792SMatthew G. Knepley inttype = H5T_NATIVE_LLONG; 105*235f7792SMatthew G. Knepley #else 106*235f7792SMatthew G. Knepley inttype = H5T_NATIVE_INT; 107*235f7792SMatthew G. Knepley #endif 108*235f7792SMatthew G. Knepley ierr = PetscMalloc1(n,&ind);CHKERRQ(ierr); 109*235f7792SMatthew G. Knepley status = H5Dread(dset_id, inttype, memspace, filespace, plist_id, (void *) ind);CHKERRQ(status); 110*235f7792SMatthew G. Knepley ierr = ISGeneralSetIndices(is, n, ind, PETSC_OWN_POINTER);CHKERRQ(ierr); 111*235f7792SMatthew G. Knepley 112*235f7792SMatthew G. Knepley /* Close/release resources */ 113*235f7792SMatthew G. Knepley if (group != file_id) {status = H5Gclose(group);CHKERRQ(status);} 114*235f7792SMatthew G. Knepley status = H5Pclose(plist_id);CHKERRQ(status); 115*235f7792SMatthew G. Knepley status = H5Sclose(filespace);CHKERRQ(status); 116*235f7792SMatthew G. Knepley status = H5Sclose(memspace);CHKERRQ(status); 117*235f7792SMatthew G. Knepley status = H5Dclose(dset_id);CHKERRQ(status); 118*235f7792SMatthew G. Knepley PetscFunctionReturn(0); 119*235f7792SMatthew G. Knepley } 120*235f7792SMatthew G. Knepley #endif 121*235f7792SMatthew G. Knepley 122*235f7792SMatthew G. Knepley #undef __FUNCT__ 123*235f7792SMatthew G. Knepley #define __FUNCT__ "ISLoad_Default" 124*235f7792SMatthew G. Knepley PetscErrorCode ISLoad_Default(IS is, PetscViewer viewer) 125*235f7792SMatthew G. Knepley { 126*235f7792SMatthew G. Knepley PetscBool isbinary; 127*235f7792SMatthew G. Knepley #if defined(PETSC_HAVE_HDF5) 128*235f7792SMatthew G. Knepley PetscBool ishdf5; 129*235f7792SMatthew G. Knepley #endif 130*235f7792SMatthew G. Knepley PetscErrorCode ierr; 131*235f7792SMatthew G. Knepley 132*235f7792SMatthew G. Knepley PetscFunctionBegin; 133*235f7792SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 134*235f7792SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr); 135*235f7792SMatthew G. Knepley if (isbinary) { 136*235f7792SMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject) is), PETSC_ERR_SUP, "This should be implemented"); 137*235f7792SMatthew G. Knepley #if defined(PETSC_HAVE_HDF5) 138*235f7792SMatthew G. Knepley } else if (ishdf5) { 139*235f7792SMatthew G. Knepley if (!((PetscObject) is)->name) { 140*235f7792SMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Since HDF5 format gives ASCII name for each object in file; must use ISLoad() after setting name of Vec with PetscObjectSetName()"); 141*235f7792SMatthew G. Knepley } 142*235f7792SMatthew G. Knepley ierr = ISLoad_HDF5(is, viewer);CHKERRQ(ierr); 143*235f7792SMatthew G. Knepley #endif 144*235f7792SMatthew G. Knepley } 145*235f7792SMatthew G. Knepley PetscFunctionReturn(0); 146*235f7792SMatthew G. Knepley } 147