1235f7792SMatthew G. Knepley #include <petscis.h> /*I "petscis.h" I*/ 2235f7792SMatthew G. Knepley #include <petsc-private/isimpl.h> 3235f7792SMatthew G. Knepley #include <petscviewerhdf5.h> 4235f7792SMatthew G. Knepley 5235f7792SMatthew G. Knepley #if defined(PETSC_HAVE_HDF5) 6235f7792SMatthew G. Knepley #undef __FUNCT__ 7235f7792SMatthew G. Knepley #define __FUNCT__ "ISLoad_HDF5" 8235f7792SMatthew G. Knepley /* 9235f7792SMatthew 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 10235f7792SMatthew G. Knepley checks back and forth between the two types of variables. 11235f7792SMatthew G. Knepley */ 12235f7792SMatthew G. Knepley PetscErrorCode ISLoad_HDF5(IS is, PetscViewer viewer) 13235f7792SMatthew G. Knepley { 14235f7792SMatthew G. Knepley hid_t inttype; /* int type (H5T_NATIVE_INT or H5T_NATIVE_LLONG) */ 15235f7792SMatthew G. Knepley hid_t file_id, group, dset_id, filespace, memspace, plist_id; 16*794a961bSBarry Smith int rdim, dim; 17235f7792SMatthew G. Knepley hsize_t dims[3], count[3], offset[3]; 1833334011SMatthew G. Knepley PetscInt n, N, bs, bsInd, lenInd, low, timestep; 19235f7792SMatthew G. Knepley const PetscInt *ind; 20235f7792SMatthew G. Knepley const char *isname; 21235f7792SMatthew G. Knepley PetscErrorCode ierr; 22235f7792SMatthew G. Knepley 23235f7792SMatthew G. Knepley PetscFunctionBegin; 24235f7792SMatthew G. Knepley ierr = PetscViewerHDF5OpenGroup(viewer, &file_id, &group);CHKERRQ(ierr); 25235f7792SMatthew G. Knepley ierr = PetscViewerHDF5GetTimestep(viewer, ×tep);CHKERRQ(ierr); 26235f7792SMatthew G. Knepley ierr = ISGetBlockSize(is, &bs);CHKERRQ(ierr); 27235f7792SMatthew G. Knepley /* Create the dataset with default properties and close filespace */ 28235f7792SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) is, &isname);CHKERRQ(ierr); 29235f7792SMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 30729ab6d9SBarry Smith PetscStackCallHDF5Return(dset_id,H5Dopen2,(group, isname, H5P_DEFAULT)); 31235f7792SMatthew G. Knepley #else 32729ab6d9SBarry Smith PetscStackCallHDF5Return(dset_id,H5Dopen,(group, isname)); 33235f7792SMatthew G. Knepley #endif 34235f7792SMatthew G. Knepley /* Retrieve the dataspace for the dataset */ 35729ab6d9SBarry Smith PetscStackCallHDF5Return(filespace,H5Dget_space,(dset_id)); 36235f7792SMatthew G. Knepley dim = 0; 37235f7792SMatthew G. Knepley if (timestep >= 0) ++dim; 38235f7792SMatthew G. Knepley ++dim; 3933334011SMatthew G. Knepley ++dim; 40729ab6d9SBarry Smith PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(filespace, dims, NULL)); 41235f7792SMatthew G. Knepley bsInd = rdim-1; 42235f7792SMatthew G. Knepley lenInd = timestep >= 0 ? 1 : 0; 4333334011SMatthew G. Knepley if (rdim != dim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Dimension of array in file %d not %d as expected",rdim,dim); 4433334011SMatthew G. Knepley else if (bs != (PetscInt) dims[bsInd]) { 4533334011SMatthew G. Knepley ierr = ISSetBlockSize(is, dims[bsInd]); 46235f7792SMatthew 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]); 47235f7792SMatthew G. Knepley bs = dims[bsInd]; 48235f7792SMatthew G. Knepley } 49235f7792SMatthew G. Knepley 50235f7792SMatthew G. Knepley /* Set Vec sizes,blocksize,and type if not already set */ 51235f7792SMatthew G. Knepley ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr); 52235f7792SMatthew G. Knepley ierr = ISGetSize(is, &N);CHKERRQ(ierr); 53235f7792SMatthew G. Knepley if (n < 0 && N < 0) {ierr = PetscLayoutSetSize(is->map, dims[lenInd]*bs);CHKERRQ(ierr);} 54235f7792SMatthew G. Knepley ierr = PetscLayoutSetUp(is->map);CHKERRQ(ierr); 55235f7792SMatthew G. Knepley /* If sizes and type already set,check if the vector global size is correct */ 56235f7792SMatthew G. Knepley ierr = ISGetSize(is, &N);CHKERRQ(ierr); 57235f7792SMatthew 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); 58235f7792SMatthew G. Knepley 59235f7792SMatthew G. Knepley /* Each process defines a dataset and reads it from the hyperslab in the file */ 60235f7792SMatthew G. Knepley ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr); 61235f7792SMatthew G. Knepley dim = 0; 62235f7792SMatthew G. Knepley if (timestep >= 0) { 63235f7792SMatthew G. Knepley count[dim] = 1; 64235f7792SMatthew G. Knepley ++dim; 65235f7792SMatthew G. Knepley } 66235f7792SMatthew G. Knepley ierr = PetscHDF5IntCast(n/bs,count + dim);CHKERRQ(ierr); 67235f7792SMatthew G. Knepley ++dim; 68235f7792SMatthew G. Knepley if (bs >= 1) { 69235f7792SMatthew G. Knepley count[dim] = bs; 70235f7792SMatthew G. Knepley ++dim; 71235f7792SMatthew G. Knepley } 72729ab6d9SBarry Smith PetscStackCallHDF5Return(memspace,H5Screate_simple,(dim, count, NULL)); 73235f7792SMatthew G. Knepley 74235f7792SMatthew G. Knepley /* Select hyperslab in the file */ 75235f7792SMatthew G. Knepley ierr = PetscLayoutGetRange(is->map, &low, NULL);CHKERRQ(ierr); 76235f7792SMatthew G. Knepley dim = 0; 77235f7792SMatthew G. Knepley if (timestep >= 0) { 78235f7792SMatthew G. Knepley offset[dim] = timestep; 79235f7792SMatthew G. Knepley ++dim; 80235f7792SMatthew G. Knepley } 81235f7792SMatthew G. Knepley ierr = PetscHDF5IntCast(low/bs,offset + dim);CHKERRQ(ierr); 82235f7792SMatthew G. Knepley ++dim; 83235f7792SMatthew G. Knepley if (bs >= 1) { 84235f7792SMatthew G. Knepley offset[dim] = 0; 85235f7792SMatthew G. Knepley ++dim; 86235f7792SMatthew G. Knepley } 87729ab6d9SBarry Smith PetscStackCallHDF5(H5Sselect_hyperslab,(filespace, H5S_SELECT_SET, offset, NULL, count, NULL)); 88235f7792SMatthew G. Knepley 89235f7792SMatthew G. Knepley /* Create property list for collective dataset read */ 90729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_DATASET_XFER)); 91235f7792SMatthew G. Knepley #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 92729ab6d9SBarry Smith PetscStackCallHDF5(H5Pset_dxpl_mpio,(plist_id, H5FD_MPIO_COLLECTIVE)); 93235f7792SMatthew G. Knepley #endif 94235f7792SMatthew G. Knepley /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 95235f7792SMatthew G. Knepley 96235f7792SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 97235f7792SMatthew G. Knepley inttype = H5T_NATIVE_LLONG; 98235f7792SMatthew G. Knepley #else 99235f7792SMatthew G. Knepley inttype = H5T_NATIVE_INT; 100235f7792SMatthew G. Knepley #endif 101235f7792SMatthew G. Knepley ierr = PetscMalloc1(n,&ind);CHKERRQ(ierr); 102729ab6d9SBarry Smith PetscStackCallHDF5(H5Dread,(dset_id, inttype, memspace, filespace, plist_id, (void *) ind)); 103235f7792SMatthew G. Knepley ierr = ISGeneralSetIndices(is, n, ind, PETSC_OWN_POINTER);CHKERRQ(ierr); 104235f7792SMatthew G. Knepley 105235f7792SMatthew G. Knepley /* Close/release resources */ 106729ab6d9SBarry Smith if (group != file_id) PetscStackCallHDF5(H5Gclose,(group)); 107729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 108729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(filespace)); 109729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(memspace)); 110729ab6d9SBarry Smith PetscStackCallHDF5(H5Dclose,(dset_id)); 111235f7792SMatthew G. Knepley PetscFunctionReturn(0); 112235f7792SMatthew G. Knepley } 113235f7792SMatthew G. Knepley #endif 114235f7792SMatthew G. Knepley 115235f7792SMatthew G. Knepley #undef __FUNCT__ 116235f7792SMatthew G. Knepley #define __FUNCT__ "ISLoad_Default" 117235f7792SMatthew G. Knepley PetscErrorCode ISLoad_Default(IS is, PetscViewer viewer) 118235f7792SMatthew G. Knepley { 119235f7792SMatthew G. Knepley PetscBool isbinary; 120235f7792SMatthew G. Knepley #if defined(PETSC_HAVE_HDF5) 121235f7792SMatthew G. Knepley PetscBool ishdf5; 122235f7792SMatthew G. Knepley #endif 123235f7792SMatthew G. Knepley PetscErrorCode ierr; 124235f7792SMatthew G. Knepley 125235f7792SMatthew G. Knepley PetscFunctionBegin; 126235f7792SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 12782be971eSMatthew G. Knepley #if defined(PETSC_HAVE_HDF5) 128235f7792SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);CHKERRQ(ierr); 12982be971eSMatthew G. Knepley #endif 130235f7792SMatthew G. Knepley if (isbinary) { 131235f7792SMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject) is), PETSC_ERR_SUP, "This should be implemented"); 132235f7792SMatthew G. Knepley #if defined(PETSC_HAVE_HDF5) 133235f7792SMatthew G. Knepley } else if (ishdf5) { 134235f7792SMatthew G. Knepley if (!((PetscObject) is)->name) { 135235f7792SMatthew 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()"); 136235f7792SMatthew G. Knepley } 137235f7792SMatthew G. Knepley ierr = ISLoad_HDF5(is, viewer);CHKERRQ(ierr); 138235f7792SMatthew G. Knepley #endif 139235f7792SMatthew G. Knepley } 140235f7792SMatthew G. Knepley PetscFunctionReturn(0); 141235f7792SMatthew G. Knepley } 142