Reads a 1D dataset located under groupname and given by name. If no offset is provided, uses default file view. Otherwise, sets file view manually.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(hdf5_obj), | intent(inout) | :: | this |
A HDF5 object |
||
| character(len=*), | intent(in) | :: | groupname |
Parent group containing the dataset |
||
| character(len=*), | intent(in) | :: | name |
Dataset name |
||
| class(*), | intent(out), | target | :: | array(:) |
Data array |
|
| integer, | intent(in), | optional | :: | offset |
Indicates number of elements to skip before reading |
impure subroutine hdf5_obj_Read1D(this,groupname,name,array,offset) !> Reads a 1D dataset located under groupname and given by name. ! If no offset is provided, uses default file view. ! Otherwise, sets file view manually. use iso_c_binding implicit none class(hdf5_obj), intent(inout) :: this !! A HDF5 object character(len=*), intent(in) :: groupname !! Parent group containing the dataset character(len=*), intent(in) :: name !! Dataset name class(*), intent(out), & target :: array(:) !! Data array integer, intent(in), & optional :: offset !! Indicates number of elements to skip before reading ! Work variables integer(HID_T) :: did integer(HID_T) :: sidc integer(HID_T) :: sidg integer(HID_T) :: plist integer(HID_T) :: dtype integer(HSIZE_T) :: dims(1) integer(HSIZE_T) :: my_offset(1) integer :: ndim integer, & allocatable :: npoint_proc(:),buff(:) integer(HID_T) :: obj integer :: ierr type(c_ptr) :: fptr = c_null_ptr ! Open group, if not already open call this%OpenGroup(groupname) ! Get group object obj = this%GetGroupObject(groupname) ! Open dataset call H5Dopen_f(obj, trim(adjustl(name)), did, ierr) if (ierr.ne.0) & call this%parallel%Stop("Unable to open dataset "//trim(adjustl(name))//" under "//trim(adjustl(groupname))//" in HDF5 file "//this%filename) ! Get the global data space from file call H5Dget_space_f(did,sidg,ierr) ! Check data is 1D call H5Sget_simple_extent_ndims_f(sidg,ndim,ierr) if (ndim.ne.1) & call this%parallel%Stop("Dimension mismatch for dataset "//trim(adjustl(name))//" under "//trim(adjustl(groupname))//" in HDF5 file "//this%filename) ! Get datatype from file call H5Dget_type_f(did,dtype,ierr) ! Get the global data space from file call H5Dget_space_f(did,sidg,ierr) ! Create the local data space dims(1) = int(size(array),kind=HSIZE_T) call H5Screate_simple_f(ndim,dims,sidc,ierr) ! Select a subset of the global data space if (present(offset)) then my_offset(1) = int(offset,HID_T) else ! Get number of grid points in each direction allocate(npoint_proc(this%parallel%nproc), source = 0) allocate(buff(this%parallel%nproc)) npoint_proc(this%parallel%rank%mine) = size(array) call this%parallel%sum(npoint_proc, buff); npoint_proc=buff deallocate(buff) my_offset(1) = int(sum(npoint_proc(1:this%parallel%rank%mine-1)),kind=HSIZE_T) end if call H5Sselect_hyperslab_f(sidg,H5S_SELECT_SET_F,my_offset,dims,ierr) ! Select all elements from the local data space my_offset = int(0,kind=HSIZE_T) call H5Sselect_hyperslab_f(sidc,H5S_SELECT_SET_F,my_offset,dims,ierr) ! Create property list call H5Pcreate_f(H5P_DATASET_XFER_F, plist, ierr ) call H5Pset_dxpl_mpio_f(plist,H5FD_MPIO_COLLECTIVE_F,ierr ) if (size(array).ne.0) then select type(array) type is (real(leapDP)) fptr = c_loc(array) type is (real(leapSP)) fptr = c_loc(array) type is (integer(leapI4)) fptr = c_loc(array) type is (integer(leapI8)) fptr = c_loc(array) end select end if ! Read dataset call H5Dread_f(did,dtype,fptr,ierr,mem_space_id=sidc,file_space_id=sidg,xfer_prp=plist) if (ierr.ne.0) & call this%parallel%Stop("Unable to read dataset "//trim(adjustl(name))//" under "//trim(adjustl(groupname))//" in HDF5 file "//this%filename) ! Deallocate pointer fptr = c_null_ptr ! Close dataset call H5Dclose_f(did,ierr) ! Close data spaces call H5Sclose_f(sidc,ierr) call H5Sclose_f(sidg,ierr) ! Close property list call H5Pclose_f(plist,ierr) ! Close datatype call H5Tclose_f(dtype,ierr) ! Close group call this%CloseGroup(groupname) return end subroutine hdf5_obj_Read1D