1 #include <petscis.h> /*I "petscis.h" I*/ 2 #include <petsc/private/isimpl.h> 3 #include <petscviewerhdf5.h> 4 5 #if defined(PETSC_HAVE_HDF5) 6 #undef __FUNCT__ 7 #define __FUNCT__ "ISLoad_HDF5" 8 /* 9 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 checks back and forth between the two types of variables. 11 */ 12 PetscErrorCode ISLoad_HDF5(IS is, PetscViewer viewer) 13 { 14 hid_t inttype; /* int type (H5T_NATIVE_INT or H5T_NATIVE_LLONG) */ 15 hid_t file_id, group, dset_id, filespace, memspace, plist_id; 16 int rdim, dim; 17 hsize_t dims[3], count[3], offset[3]; 18 PetscInt n, N, bs, bsInd, lenInd, low, timestep; 19 const PetscInt *ind; 20 const char *isname; 21 PetscErrorCode ierr; 22 23 PetscFunctionBegin; 24 if (!((PetscObject)is)->name) SETERRQ(PetscObjectComm((PetscObject)is), PETSC_ERR_SUP, "Since HDF5 format gives ASCII name for each object in file; must use ISLoad() after setting name of Vec with PetscObjectSetName()"); 25 ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr); 26 ierr = PetscViewerHDF5GetTimestep(viewer, ×tep);CHKERRQ(ierr); 27 ierr = ISGetBlockSize(is, &bs);CHKERRQ(ierr); 28 /* Create the dataset with default properties and close filespace */ 29 ierr = PetscObjectGetName((PetscObject) is, &isname);CHKERRQ(ierr); 30 #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 31 PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, isname, H5P_DEFAULT)); 32 #else 33 PetscStackCallHDF5Return(dset_id,H5Dopen,(group, isname)); 34 #endif 35 /* Retrieve the dataspace for the dataset */ 36 PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id)); 37 dim = 0; 38 if (timestep >= 0) ++dim; 39 ++dim; 40 ++dim; 41 PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(filespace, dims, NULL)); 42 bsInd = rdim-1; 43 lenInd = timestep >= 0 ? 1 : 0; 44 if (rdim != dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Dimension of array in file %d not %d as expected",rdim,dim); 45 else if (bs != (PetscInt) dims[bsInd]) { 46 ierr = ISSetBlockSize(is, dims[bsInd]); 47 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]); 48 bs = dims[bsInd]; 49 } 50 51 /* Set Vec sizes,blocksize,and type if not already set */ 52 ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr); 53 ierr = ISGetSize(is, &N);CHKERRQ(ierr); 54 if (n < 0 && N < 0) {ierr = PetscLayoutSetSize(is->map, dims[lenInd]*bs);CHKERRQ(ierr);} 55 ierr = PetscLayoutSetUp(is->map);CHKERRQ(ierr); 56 /* If sizes and type already set,check if the vector global size is correct */ 57 ierr = ISGetSize(is, &N);CHKERRQ(ierr); 58 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); 59 60 /* Each process defines a dataset and reads it from the hyperslab in the file */ 61 ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr); 62 dim = 0; 63 if (timestep >= 0) { 64 count[dim] = 1; 65 ++dim; 66 } 67 ierr = PetscHDF5IntCast(n/bs,count + dim);CHKERRQ(ierr); 68 ++dim; 69 if (bs >= 1) { 70 count[dim] = bs; 71 ++dim; 72 } 73 PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL)); 74 75 /* Select hyperslab in the file */ 76 ierr = PetscLayoutGetRange(is->map, &low, NULL);CHKERRQ(ierr); 77 dim = 0; 78 if (timestep >= 0) { 79 offset[dim] = timestep; 80 ++dim; 81 } 82 ierr = PetscHDF5IntCast(low/bs,offset + dim);CHKERRQ(ierr); 83 ++dim; 84 if (bs >= 1) { 85 offset[dim] = 0; 86 ++dim; 87 } 88 PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL)); 89 90 /* Create property list for collective dataset read */ 91 PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_DATASET_XFER)); 92 #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 93 PetscStackCallHDF5(H5Pset_dxpl_mpio,(plist_id, H5FD_MPIO_COLLECTIVE)); 94 #endif 95 /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 96 97 #if defined(PETSC_USE_64BIT_INDICES) 98 inttype = H5T_NATIVE_LLONG; 99 #else 100 inttype = H5T_NATIVE_INT; 101 #endif 102 ierr = PetscMalloc1(n,&ind);CHKERRQ(ierr); 103 PetscStackCallHDF5(H5Dread,(dset_id, inttype, memspace, filespace, plist_id, (void *) ind)); 104 ierr = ISGeneralSetIndices(is, n, ind, PETSC_OWN_POINTER);CHKERRQ(ierr); 105 106 /* Close/release resources */ 107 if (group != file_id) PetscStackCallHDF5(H5Gclose,(group)); 108 PetscStackCallHDF5(H5Pclose,(plist_id)); 109 PetscStackCallHDF5(H5Sclose,(filespace)); 110 PetscStackCallHDF5(H5Sclose,(memspace)); 111 PetscStackCallHDF5(H5Dclose,(dset_id)); 112 PetscFunctionReturn(0); 113 } 114 #endif 115 116 #undef __FUNCT__ 117 #define __FUNCT__ "ISLoad_Binary" 118 PetscErrorCode ISLoad_Binary(IS is, PetscViewer viewer) 119 { 120 PetscErrorCode ierr; 121 PetscBool isgeneral,skipheader; 122 int fd; 123 PetscInt tr[2],N,ln,*idx; 124 MPI_Request request; 125 MPI_Status status; 126 MPI_Comm comm; 127 PetscMPIInt rank,size,tag; 128 129 PetscFunctionBegin; 130 ierr = PetscObjectGetComm((PetscObject)is,&comm);CHKERRQ(ierr); 131 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 132 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 133 ierr = PetscObjectGetNewTag((PetscObject)viewer,&tag);CHKERRQ(ierr); 134 /* force binary viewer to load .info file if it has not yet done so */ 135 ierr = PetscViewerSetUp(viewer);CHKERRQ(ierr); 136 ierr = PetscObjectTypeCompare((PetscObject)is,ISGENERAL,&isgeneral);CHKERRQ(ierr); 137 if (!isgeneral) SETERRQ(comm,PETSC_ERR_ARG_INCOMP,"IS must be of type ISGENERAL to load into it"); 138 ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 139 ierr = PetscViewerBinaryGetSkipHeader(viewer,&skipheader);CHKERRQ(ierr); 140 if (skipheader) SETERRQ(comm,PETSC_ERR_USER, "Currently no support for binary files without headers"); 141 142 ierr = PetscViewerBinaryRead(viewer,tr,2,NULL,PETSC_INT);CHKERRQ(ierr); 143 if (tr[0] != IS_FILE_CLASSID) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"Not an IS next in file"); 144 145 /* Has IS already had its layout defined */ 146 ierr = PetscLayoutGetSize(is->map,&N);CHKERRQ(ierr); 147 if (N > -1 && N != tr[1]) SETERRQ2(comm,PETSC_ERR_ARG_SIZ,"Size of IS in file %D does not match size of IS provided",tr[1],N); 148 if (N == -1) { 149 N = tr[1]; 150 ierr = PetscLayoutSetSize(is->map,N);CHKERRQ(ierr); 151 ierr = PetscLayoutSetUp(is->map);CHKERRQ(ierr); 152 } 153 ierr = PetscLayoutGetLocalSize(is->map,&ln);CHKERRQ(ierr); 154 ierr = PetscMalloc1(ln,&idx);CHKERRQ(ierr); 155 if (!rank) { 156 ierr = PetscBinaryRead(fd,idx,ln,PETSC_INT);CHKERRQ(ierr); 157 158 if (size > 1) { 159 PetscInt *range,n,i,*idxwork; 160 161 /* read in other chuncks and send to other processors */ 162 /* determine maximum chunck owned by other */ 163 range = is->map->range; 164 n = 1; 165 for (i=1; i<size; i++) n = PetscMax(n,range[i+1] - range[i]); 166 167 ierr = PetscMalloc1(n,&idxwork);CHKERRQ(ierr); 168 for (i=1; i<size; i++) { 169 n = range[i+1] - range[i]; 170 ierr = PetscBinaryRead(fd,idxwork,n,PETSC_INT);CHKERRQ(ierr); 171 ierr = MPI_Isend(idxwork,n,MPIU_INT,i,tag,comm,&request);CHKERRQ(ierr); 172 ierr = MPI_Wait(&request,&status);CHKERRQ(ierr); 173 } 174 ierr = PetscFree(idxwork);CHKERRQ(ierr); 175 } 176 } else { 177 ierr = MPI_Recv(idx,ln,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr); 178 } 179 ierr = ISGeneralSetIndices(is,ln,idx,PETSC_OWN_POINTER);CHKERRQ(ierr); 180 PetscFunctionReturn(0); 181 } 182 183 #undef __FUNCT__ 184 #define __FUNCT__ "ISLoad_Default" 185 PetscErrorCode ISLoad_Default(IS is, PetscViewer viewer) 186 { 187 PetscBool isbinary; 188 #if defined(PETSC_HAVE_HDF5) 189 PetscBool ishdf5; 190 #endif 191 PetscErrorCode ierr; 192 193 PetscFunctionBegin; 194 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 195 #if defined(PETSC_HAVE_HDF5) 196 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr); 197 #endif 198 if (isbinary) { 199 ierr = ISLoad_Binary(is, viewer);CHKERRQ(ierr); 200 #if defined(PETSC_HAVE_HDF5) 201 } else if (ishdf5) { 202 ierr = ISLoad_HDF5(is, viewer);CHKERRQ(ierr); 203 #endif 204 } 205 PetscFunctionReturn(0); 206 } 207