hdf5_obj_ReadAttributes0D Subroutine

private impure subroutine hdf5_obj_ReadAttributes0D(this, groupname, label, val)

Uses

  • proc~~hdf5_obj_readattributes0d~~UsesGraph proc~hdf5_obj_readattributes0d hdf5_obj%hdf5_obj_ReadAttributes0D iso_c_binding iso_c_binding proc~hdf5_obj_readattributes0d->iso_c_binding

Read a scalar attribute under a given group.

Type Bound

hdf5_obj

Arguments

Type IntentOptional Attributes Name
class(hdf5_obj), intent(inout) :: this

A HDF5 object

character(len=*), intent(in) :: groupname

Groupname

character(len=*), intent(in) :: label

Attribute label

class(*), intent(out), target :: val

Attribute value


Calls

proc~~hdf5_obj_readattributes0d~~CallsGraph proc~hdf5_obj_readattributes0d hdf5_obj%hdf5_obj_ReadAttributes0D h5aclose_f h5aclose_f proc~hdf5_obj_readattributes0d->h5aclose_f h5aget_space_f h5aget_space_f proc~hdf5_obj_readattributes0d->h5aget_space_f h5aget_type_f h5aget_type_f proc~hdf5_obj_readattributes0d->h5aget_type_f h5aopen_f h5aopen_f proc~hdf5_obj_readattributes0d->h5aopen_f h5aread_f h5aread_f proc~hdf5_obj_readattributes0d->h5aread_f h5sclose_f h5sclose_f proc~hdf5_obj_readattributes0d->h5sclose_f h5sget_simple_extent_dims_f h5sget_simple_extent_dims_f proc~hdf5_obj_readattributes0d->h5sget_simple_extent_dims_f proc~hdf5_obj_closegroup hdf5_obj%hdf5_obj_CloseGroup proc~hdf5_obj_readattributes0d->proc~hdf5_obj_closegroup proc~hdf5_obj_getgroupobject hdf5_obj%hdf5_obj_GetGroupObject proc~hdf5_obj_readattributes0d->proc~hdf5_obj_getgroupobject proc~hdf5_obj_opengroup hdf5_obj%hdf5_obj_OpenGroup proc~hdf5_obj_readattributes0d->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

Called by

proc~~hdf5_obj_readattributes0d~~CalledByGraph proc~hdf5_obj_readattributes0d hdf5_obj%hdf5_obj_ReadAttributes0D none~readattributes hdf5_obj%ReadAttributes none~readattributes->proc~hdf5_obj_readattributes0d proc~bc_set_read bc_set%bc_set_Read proc~bc_set_read->none~readattributes proc~block_obj_read block_obj%block_obj_Read proc~block_obj_read->none~readattributes proc~eulerian_set_readhdf5 eulerian_set%eulerian_set_ReadHDF5 proc~eulerian_set_readhdf5->none~readattributes proc~h5hut_obj_jumptostep h5hut_obj%h5hut_obj_JumpToStep proc~h5hut_obj_jumptostep->none~readattributes proc~h5hut_obj_lasttimestep h5hut_obj%h5hut_obj_LastTimeStep proc~h5hut_obj_lasttimestep->none~readattributes proc~h5hut_obj_readattributes0d h5hut_obj%h5hut_obj_ReadAttributes0D proc~h5hut_obj_readattributes0d->none~readattributes proc~h5hut_obj_readattributes1d h5hut_obj%h5hut_obj_ReadAttributes1D proc~h5hut_obj_readattributes1d->none~readattributes proc~marker_set_readhdf5 marker_set%marker_set_ReadHDF5 proc~marker_set_readhdf5->none~readattributes proc~particle_set_readhdf5 particle_set%particle_set_ReadHDF5 proc~particle_set_readhdf5->none~readattributes proc~respart_set_readhdf5 ResPart_set%ResPart_set_ReadHDF5 proc~respart_set_readhdf5->none~readattributes proc~respart_set_readhdf5->proc~marker_set_readhdf5 proc~solid_set_readhdf5 solid_set%solid_set_ReadHDF5 proc~solid_set_readhdf5->none~readattributes proc~solid_set_readhdf5->proc~marker_set_readhdf5 none~readattributes~2 h5hut_obj%ReadAttributes none~readattributes~2->proc~h5hut_obj_readattributes0d none~readattributes~2->proc~h5hut_obj_readattributes1d proc~eulerian_set_readh5hut eulerian_set%eulerian_set_ReadH5HUT proc~eulerian_set_readh5hut->proc~h5hut_obj_jumptostep proc~eulerian_set_readh5hut->proc~h5hut_obj_lasttimestep proc~marker_set_readh5hut marker_set%marker_set_ReadH5HUT proc~marker_set_readh5hut->proc~h5hut_obj_jumptostep proc~marker_set_readh5hut->proc~h5hut_obj_lasttimestep proc~particle_set_readh5hut particle_set%particle_set_ReadH5HUT proc~particle_set_readh5hut->proc~h5hut_obj_jumptostep proc~particle_set_readh5hut->proc~h5hut_obj_lasttimestep proc~respart_set_readh5hut ResPart_set%ResPart_set_ReadH5HUT proc~respart_set_readh5hut->proc~h5hut_obj_jumptostep proc~respart_set_readh5hut->proc~h5hut_obj_lasttimestep proc~solid_set_readh5hut solid_set%solid_set_ReadH5HUT proc~solid_set_readh5hut->none~readattributes~2

Source Code

    impure subroutine hdf5_obj_ReadAttributes0D(this,groupname,label,val)
      !> Read a scalar attribute under a given group.
      use iso_c_binding
      implicit none
      class(hdf5_obj),  intent(inout) :: this                                  !! A HDF5 object
      character(len=*), intent(in)    :: groupname                             !! Groupname
      character(len=*), intent(in)    :: label                                 !! Attribute label
      class(*),         intent(out),   &
                               target :: val                                   !! Attribute value
      ! Work variables
      integer         :: ierr
      integer(HID_T)  :: obj
      integer, target :: int2log
      type(c_ptr)     :: fptr
      integer(HID_T)  :: dtype
      integer(HID_T)  :: sid
      integer(HID_T)  :: attrid
      integer(HSIZE_T):: dims(1)
      integer(HSIZE_T):: maxdims(1)

      ! Open group
      call this%OpenGroup(groupname)

      ! Open the attribute
      if (trim(adjustl(groupname)).eq.'/') then
        call H5Aopen_f(this%fid, trim(adjustl(label)), attrid, ierr, aapl_id=H5P_DEFAULT_F)
      else
        ! Get this group's object
        obj = this%GetGroupObject(groupname)
        if (obj.eq.-1_HID_T) &
          call this%parallel%Stop("Unable to find group "//trim(adjustl(groupname))//" in hdf5 file "//this%filename)

        call H5Aopen_f(obj, trim(adjustl(label)), attrid, ierr, aapl_id=H5P_DEFAULT_F)
      end if

      ! Get datatype and spaceid
      call H5Aget_type_f(attrid,dtype,ierr)
      call H5Aget_space_f(attrid,sid,ierr)

      ! Check dimensions
      call H5Sget_simple_extent_dims_f(sid, dims, maxdims, ierr)
      if (dims(1).ne.1) &
          call this%parallel%Stop("Size mismatch when reading attribute"//trim(adjustl(label))//" under "//trim(adjustl(groupname))//" in hdf5 file "//this%filename)

      ! Read attribute
      select type(val)
      type is (real(leapDP))
        fptr = c_loc(val)
      type is (real(leapSP))
        fptr = c_loc(val)
      type is (integer(leapI4))
        fptr = c_loc(val)
      type is (integer(leapI8))
        fptr = c_loc(val)
      type is (logical)
        fptr = c_loc(int2log)
      end select
      call H5Aread_f(attrid, dtype,fptr, ierr)

      select type(val)
      type is (logical)
        if (int2log.eq.1) then
          val = .true.
        else
          val = .false.
        end if
      end select

      ! Deallocate pointer
      fptr = c_null_ptr

      ! Remove attribute ID
      call H5Aclose_f(attrid, ierr)

      ! Close data space
      call H5Sclose_f(sid, ierr)

      ! Close group
      call this%CloseGroup(groupname)

      return
    end subroutine hdf5_obj_ReadAttributes0D