hdf5_obj_Read3D Subroutine

private impure subroutine hdf5_obj_Read3D(this, groupname, name, array, lo, hi)

Uses

  • proc~~hdf5_obj_read3d~~UsesGraph proc~hdf5_obj_read3d hdf5_obj%hdf5_obj_Read3D iso_c_binding iso_c_binding proc~hdf5_obj_read3d->iso_c_binding

Reads a 3D dataset located under groupname and given by name.

Type Bound

hdf5_obj

Arguments

Type IntentOptional 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) :: lo(3)

Low bounds

integer, intent(in) :: hi(3)

High bounds


Calls

proc~~hdf5_obj_read3d~~CallsGraph proc~hdf5_obj_read3d hdf5_obj%hdf5_obj_Read3D h5dclose_f h5dclose_f proc~hdf5_obj_read3d->h5dclose_f h5dget_space_f h5dget_space_f proc~hdf5_obj_read3d->h5dget_space_f h5dget_type_f h5dget_type_f proc~hdf5_obj_read3d->h5dget_type_f h5dopen_f h5dopen_f proc~hdf5_obj_read3d->h5dopen_f h5dread_f h5dread_f proc~hdf5_obj_read3d->h5dread_f h5pclose_f h5pclose_f proc~hdf5_obj_read3d->h5pclose_f h5pcreate_f h5pcreate_f proc~hdf5_obj_read3d->h5pcreate_f h5pset_dxpl_mpio_f h5pset_dxpl_mpio_f proc~hdf5_obj_read3d->h5pset_dxpl_mpio_f h5sclose_f h5sclose_f proc~hdf5_obj_read3d->h5sclose_f h5screate_simple_f h5screate_simple_f proc~hdf5_obj_read3d->h5screate_simple_f h5sget_simple_extent_ndims_f h5sget_simple_extent_ndims_f proc~hdf5_obj_read3d->h5sget_simple_extent_ndims_f h5sselect_hyperslab_f h5sselect_hyperslab_f proc~hdf5_obj_read3d->h5sselect_hyperslab_f h5tclose_f h5tclose_f proc~hdf5_obj_read3d->h5tclose_f proc~hdf5_obj_closegroup hdf5_obj%hdf5_obj_CloseGroup proc~hdf5_obj_read3d->proc~hdf5_obj_closegroup proc~hdf5_obj_getgroupobject hdf5_obj%hdf5_obj_GetGroupObject proc~hdf5_obj_read3d->proc~hdf5_obj_getgroupobject proc~hdf5_obj_opengroup hdf5_obj%hdf5_obj_OpenGroup proc~hdf5_obj_read3d->proc~hdf5_obj_opengroup proc~hdf5_obj_closegroup->proc~hdf5_obj_getgroupobject h5gclose_f h5gclose_f proc~hdf5_obj_closegroup->h5gclose_f proc~hashtbl_obj_hashstring hashtbl_obj%hashtbl_obj_HashString proc~hdf5_obj_closegroup->proc~hashtbl_obj_hashstring proc~hashtbl_obj_remove hashtbl_obj%hashtbl_obj_Remove proc~hdf5_obj_closegroup->proc~hashtbl_obj_remove proc~hdf5_obj_fixgroupname hdf5_obj%hdf5_obj_FixGroupName proc~hdf5_obj_closegroup->proc~hdf5_obj_fixgroupname none~get~3 hashtbl_obj%Get proc~hdf5_obj_getgroupobject->none~get~3 proc~hdf5_obj_getgroupobject->proc~hashtbl_obj_hashstring proc~hdf5_obj_getgroupobject->proc~hdf5_obj_fixgroupname proc~hdf5_obj_opengroup->proc~hdf5_obj_getgroupobject h5oopen_f h5oopen_f proc~hdf5_obj_opengroup->h5oopen_f proc~hdf5_obj_opengroup->proc~hashtbl_obj_hashstring proc~hashtbl_obj_put hashtbl_obj%hashtbl_obj_Put proc~hdf5_obj_opengroup->proc~hashtbl_obj_put proc~hdf5_obj_opengroup->proc~hdf5_obj_fixgroupname proc~hashtbl_obj_get_int4 hashtbl_obj%hashtbl_obj_Get_int4 none~get~3->proc~hashtbl_obj_get_int4 proc~hashtbl_obj_get_int8 hashtbl_obj%hashtbl_obj_Get_int8 none~get~3->proc~hashtbl_obj_get_int8 proc~hashtbl_obj_get_real_dp hashtbl_obj%hashtbl_obj_Get_real_dp none~get~3->proc~hashtbl_obj_get_real_dp proc~hashtbl_obj_get_real_sp hashtbl_obj%hashtbl_obj_Get_real_sp none~get~3->proc~hashtbl_obj_get_real_sp proc~sllist_obj_put sllist_obj%sllist_obj_Put proc~hashtbl_obj_put->proc~sllist_obj_put proc~sllist_obj_remove sllist_obj%sllist_obj_Remove proc~hashtbl_obj_remove->proc~sllist_obj_remove none~get~2 sllist_obj%Get proc~hashtbl_obj_get_int4->none~get~2 proc~hashtbl_obj_get_int8->none~get~2 proc~hashtbl_obj_get_real_dp->none~get~2 proc~hashtbl_obj_get_real_sp->none~get~2 proc~sllist_obj_put->proc~sllist_obj_put proc~sllist_obj_get_int4 sllist_obj%sllist_obj_Get_int4 none~get~2->proc~sllist_obj_get_int4 proc~sllist_obj_get_int8 sllist_obj%sllist_obj_Get_int8 none~get~2->proc~sllist_obj_get_int8 proc~sllist_obj_get_real_dp sllist_obj%sllist_obj_Get_real_dp none~get~2->proc~sllist_obj_get_real_dp proc~sllist_obj_get_real_sp sllist_obj%sllist_obj_Get_real_sp none~get~2->proc~sllist_obj_get_real_sp proc~sllist_obj_get_int4->proc~sllist_obj_get_int4 proc~sllist_obj_get_int8->proc~sllist_obj_get_int8 proc~sllist_obj_get_real_dp->proc~sllist_obj_get_real_dp proc~sllist_obj_get_real_sp->proc~sllist_obj_get_real_sp

Source Code

    impure subroutine hdf5_obj_Read3D(this,groupname,name,array,lo,hi)
      !> Reads a 3D dataset located under groupname and given by name.
      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)    :: lo(3)                                 !! Low bounds
      integer,          intent(in)    :: hi(3)                                 !! High bounds
      ! 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(3)
      integer(HSIZE_T) :: offset(3)
      integer          :: ndim
      integer          :: npoint
      integer          :: maxhi(3)
      integer          :: minlo(3)
      integer(HID_T)   :: obj
      integer          :: ierr
      type(c_ptr)      :: fptr = c_null_ptr

      ! Get number of grid points in each direction
      call this%parallel%max(hi,maxhi)
      call this%parallel%min(lo,minlo)

      ! Check for invalid arrays
      npoint = sum(hi-lo+1)
      if (npoint.lt.0) then
        dims = int(0,kind=HSIZE_T)
      else
        dims = int(hi-lo+1,kind=HSIZE_T)
      end if

      ! Open group, if not already open
      call this%OpenGroup(groupname)

      ! Find parent 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 3D
      call H5Sget_simple_extent_ndims_f(sidg,ndim,ierr)
      if (ndim.ne.3) &
        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)

      ! Create the local data space
      call H5Screate_simple_f(ndim,dims,sidc,ierr)

      ! Select a subset of the global data space
      offset = int(lo - minlo,kind=HSIZE_T)
      call H5Sselect_hyperslab_f(sidg,H5S_SELECT_SET_F,offset,dims,ierr)

      ! Select all elements from the local data space
      offset = int(0,kind=HSIZE_T)
      call H5Sselect_hyperslab_f(sidc,H5S_SELECT_SET_F,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_Read3D