hdf5_obj_WriteAttributes0D Subroutine

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

Uses

  • proc~~hdf5_obj_writeattributes0d~~UsesGraph proc~hdf5_obj_writeattributes0d hdf5_obj%hdf5_obj_WriteAttributes0D iso_c_binding iso_c_binding proc~hdf5_obj_writeattributes0d->iso_c_binding

Writes a scalar attribute.

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(in), target :: val

Attribute label


Calls

proc~~hdf5_obj_writeattributes0d~~CallsGraph proc~hdf5_obj_writeattributes0d hdf5_obj%hdf5_obj_WriteAttributes0D h5aclose_f h5aclose_f proc~hdf5_obj_writeattributes0d->h5aclose_f h5acreate_f h5acreate_f proc~hdf5_obj_writeattributes0d->h5acreate_f h5awrite_f h5awrite_f proc~hdf5_obj_writeattributes0d->h5awrite_f h5sclose_f h5sclose_f proc~hdf5_obj_writeattributes0d->h5sclose_f h5screate_simple_f h5screate_simple_f proc~hdf5_obj_writeattributes0d->h5screate_simple_f proc~hdf5_obj_closegroup hdf5_obj%hdf5_obj_CloseGroup proc~hdf5_obj_writeattributes0d->proc~hdf5_obj_closegroup proc~hdf5_obj_getgroupobject hdf5_obj%hdf5_obj_GetGroupObject proc~hdf5_obj_writeattributes0d->proc~hdf5_obj_getgroupobject proc~hdf5_obj_opengroup hdf5_obj%hdf5_obj_OpenGroup proc~hdf5_obj_writeattributes0d->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_writeattributes0d~~CalledByGraph proc~hdf5_obj_writeattributes0d hdf5_obj%hdf5_obj_WriteAttributes0D none~writeattributes hdf5_obj%WriteAttributes none~writeattributes->proc~hdf5_obj_writeattributes0d proc~bc_set_write bc_set%bc_set_Write proc~bc_set_write->none~writeattributes proc~block_obj_write block_obj%block_obj_Write proc~block_obj_write->none~writeattributes proc~eulerian_set_writehdf5 eulerian_set%eulerian_set_WriteHDF5 proc~eulerian_set_writehdf5->none~writeattributes none~writesingle eulerian_set%WriteSingle proc~eulerian_set_writehdf5->none~writesingle proc~h5hut_obj_newtimestep h5hut_obj%h5hut_obj_NewTimeStep proc~h5hut_obj_newtimestep->none~writeattributes proc~h5hut_obj_writeattributes0d h5hut_obj%h5hut_obj_WriteAttributes0D proc~h5hut_obj_writeattributes0d->none~writeattributes proc~h5hut_obj_writeattributes1d h5hut_obj%h5hut_obj_WriteAttributes1D proc~h5hut_obj_writeattributes1d->none~writeattributes proc~h5hut_obj_writegrid h5hut_obj%h5hut_obj_WriteGrid proc~h5hut_obj_writegrid->none~writeattributes proc~marker_set_writehdf5 marker_set%marker_set_WriteHDF5 proc~marker_set_writehdf5->none~writeattributes proc~particle_set_writehdf5 particle_set%particle_set_WriteHDF5 proc~particle_set_writehdf5->none~writeattributes proc~respart_set_writehdf5 ResPart_set%ResPart_set_WriteHDF5 proc~respart_set_writehdf5->none~writeattributes proc~respart_set_writehdf5->proc~marker_set_writehdf5 proc~solid_set_writehdf5 solid_set%solid_set_WriteHDF5 proc~solid_set_writehdf5->none~writeattributes proc~solid_set_writehdf5->proc~marker_set_writehdf5 none~writeattributes~2 h5hut_obj%WriteAttributes none~writeattributes~2->proc~h5hut_obj_writeattributes0d none~writeattributes~2->proc~h5hut_obj_writeattributes1d proc~eulerian_set_writeh5hut eulerian_set%eulerian_set_WriteH5HUT proc~eulerian_set_writeh5hut->proc~h5hut_obj_newtimestep proc~eulerian_set_writeh5hut->none~writesingle proc~eulerian_set_writesingleh5hut eulerian_set%eulerian_set_WriteSingleH5HUT proc~eulerian_set_writesingleh5hut->proc~h5hut_obj_writegrid proc~marker_set_writeh5hut marker_set%marker_set_WriteH5HUT proc~marker_set_writeh5hut->proc~h5hut_obj_newtimestep proc~particle_set_writeh5hut particle_set%particle_set_WriteH5HUT proc~particle_set_writeh5hut->proc~h5hut_obj_newtimestep proc~respart_set_writeh5hut ResPart_set%ResPart_set_WriteH5HUT proc~respart_set_writeh5hut->proc~h5hut_obj_newtimestep proc~respart_set_writeh5hut->proc~marker_set_writeh5hut none~writesingle->proc~eulerian_set_writesingleh5hut proc~solid_set_writeh5hut solid_set%solid_set_WriteH5HUT proc~solid_set_writeh5hut->none~writeattributes~2 proc~eulerian_set_writesilo eulerian_set%eulerian_set_WriteSILO proc~eulerian_set_writesilo->none~writesingle

Source Code

    impure subroutine hdf5_obj_WriteAttributes0D(this,groupname,label,val)
      !> Writes a scalar attribute.
      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(in),   &
                               target :: val                                   !! Attribute label
      ! Work variables
      integer         :: ierr
      integer(HID_T)  :: obj
      integer, target :: log2int
      type(c_ptr)     :: fptr
      integer(HID_T)  :: sid
      integer(HID_T)  :: dtype
      integer(HID_T)  :: attrid
      integer(HSIZE_T):: dims(1)
      integer         :: ndim

      ! Open group
      call this%OpenGroup(groupname)

      select type(val)
      type is (real(leapDP))
        dtype = H5T_NATIVE_DOUBLE
      type is (real(leapSP))
        dtype = H5T_NATIVE_REAL
      type is (integer(leapI4))
        dtype = H5T_NATIVE_INTEGER
      type is (integer(leapI8))
        dtype = H5T_STD_I64LE
      type is (logical)
        dtype = H5T_NATIVE_INTEGER
      end select

      ! Create data space
      dims(1) = 1
      ndim    = 1
      call H5Screate_simple_f(ndim, dims, sid, ierr)

      ! Create Attribute ID
      if (trim(adjustl(groupname)).eq.'/') then
        call H5Acreate_f(this%fid, trim(adjustl(label)), dtype, sid, attrid, ierr)
      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 H5Acreate_f(obj, trim(adjustl(label)), dtype, sid, attrid, ierr)
      end if

      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)
        ! Special treatment of FORTRAN logicals
        log2int = 0
        if (val.eqv..true.) log2int = 1
        fptr = c_loc(log2int)
      end select

      call H5Awrite_f(attrid,dtype,fptr,ierr)

      ! 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_WriteAttributes0D