h5hut_obj_WriteGrid Subroutine

private impure subroutine h5hut_obj_WriteGrid(this, name, xlo, dx)

Writes grid attributes.

Type Bound

h5hut_obj

Arguments

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

A H5hut object

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

Variable name

real(kind=wp), intent(in) :: xlo(3)

Coordinate of the low bound

real(kind=wp), intent(in) :: dx(3)

Grid spacing


Calls

proc~~h5hut_obj_writegrid~~CallsGraph proc~h5hut_obj_writegrid h5hut_obj%h5hut_obj_WriteGrid none~writeattributes hdf5_obj%WriteAttributes proc~h5hut_obj_writegrid->none~writeattributes proc~hdf5_obj_creategroup hdf5_obj%hdf5_obj_CreateGroup proc~h5hut_obj_writegrid->proc~hdf5_obj_creategroup proc~hdf5_obj_writeattributes0d hdf5_obj%hdf5_obj_WriteAttributes0D none~writeattributes->proc~hdf5_obj_writeattributes0d proc~hdf5_obj_writeattributes1d hdf5_obj%hdf5_obj_WriteAttributes1D none~writeattributes->proc~hdf5_obj_writeattributes1d h5gcreate_f h5gcreate_f proc~hdf5_obj_creategroup->h5gcreate_f proc~hashtbl_obj_hashstring hashtbl_obj%hashtbl_obj_HashString proc~hdf5_obj_creategroup->proc~hashtbl_obj_hashstring proc~hashtbl_obj_put hashtbl_obj%hashtbl_obj_Put proc~hdf5_obj_creategroup->proc~hashtbl_obj_put proc~hdf5_obj_closegroup hdf5_obj%hdf5_obj_CloseGroup proc~hdf5_obj_creategroup->proc~hdf5_obj_closegroup proc~hdf5_obj_fixgroupname hdf5_obj%hdf5_obj_FixGroupName proc~hdf5_obj_creategroup->proc~hdf5_obj_fixgroupname proc~sllist_obj_put sllist_obj%sllist_obj_Put proc~hashtbl_obj_put->proc~sllist_obj_put proc~hdf5_obj_closegroup->proc~hashtbl_obj_hashstring proc~hdf5_obj_closegroup->proc~hdf5_obj_fixgroupname h5gclose_f h5gclose_f proc~hdf5_obj_closegroup->h5gclose_f proc~hashtbl_obj_remove hashtbl_obj%hashtbl_obj_Remove proc~hdf5_obj_closegroup->proc~hashtbl_obj_remove proc~hdf5_obj_getgroupobject hdf5_obj%hdf5_obj_GetGroupObject proc~hdf5_obj_closegroup->proc~hdf5_obj_getgroupobject proc~hdf5_obj_writeattributes0d->proc~hdf5_obj_closegroup 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_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_writeattributes1d->proc~hdf5_obj_closegroup proc~hdf5_obj_writeattributes1d->h5aclose_f proc~hdf5_obj_writeattributes1d->h5acreate_f proc~hdf5_obj_writeattributes1d->h5awrite_f proc~hdf5_obj_writeattributes1d->h5sclose_f proc~hdf5_obj_writeattributes1d->h5screate_simple_f proc~hdf5_obj_writeattributes1d->proc~hdf5_obj_getgroupobject proc~hdf5_obj_writeattributes1d->proc~hdf5_obj_opengroup proc~sllist_obj_remove sllist_obj%sllist_obj_Remove proc~hashtbl_obj_remove->proc~sllist_obj_remove proc~hdf5_obj_getgroupobject->proc~hashtbl_obj_hashstring proc~hdf5_obj_getgroupobject->proc~hdf5_obj_fixgroupname none~get~3 hashtbl_obj%Get proc~hdf5_obj_getgroupobject->none~get~3 proc~hdf5_obj_opengroup->proc~hashtbl_obj_hashstring proc~hdf5_obj_opengroup->proc~hashtbl_obj_put proc~hdf5_obj_opengroup->proc~hdf5_obj_fixgroupname proc~hdf5_obj_opengroup->proc~hdf5_obj_getgroupobject h5oopen_f h5oopen_f proc~hdf5_obj_opengroup->h5oopen_f proc~sllist_obj_put->proc~sllist_obj_put 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 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_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~~h5hut_obj_writegrid~~CalledByGraph proc~h5hut_obj_writegrid h5hut_obj%h5hut_obj_WriteGrid proc~eulerian_set_writesingleh5hut eulerian_set%eulerian_set_WriteSingleH5HUT proc~eulerian_set_writesingleh5hut->proc~h5hut_obj_writegrid none~writesingle eulerian_set%WriteSingle none~writesingle->proc~eulerian_set_writesingleh5hut proc~eulerian_set_writeh5hut eulerian_set%eulerian_set_WriteH5HUT proc~eulerian_set_writeh5hut->none~writesingle proc~eulerian_set_writehdf5 eulerian_set%eulerian_set_WriteHDF5 proc~eulerian_set_writehdf5->none~writesingle proc~eulerian_set_writesilo eulerian_set%eulerian_set_WriteSILO proc~eulerian_set_writesilo->none~writesingle

Source Code

    impure subroutine h5hut_obj_WriteGrid(this,name,xlo,dx)
      !> Writes grid attributes.
      implicit none
      class(h5hut_obj), intent(inout) :: this                                  !! A H5hut object
      character(len=*), intent(in) :: name                                     !! Variable name
      real(wp),         intent(in) :: xlo(3)                                   !! Coordinate of the low bound
      real(wp),         intent(in) :: dx(3)                                    !! Grid spacing
      ! Work variables
      real(leapDP)     :: xlo_(3),dx_(3)

      xlo_=real(xlo, leapDP)
      dx_ =real(dx,  leapDP)

      ! Create paraent block group if needed
      if (.not.this%block_group_exists)  then
        call this%hdf5%CreateGroup(trim(adjustl(this%step_name))//'Block/')
        this%block_group_exists = .true.
      end if

      ! Set block attributes
      call this%hdf5%WriteAttributes(trim(adjustl(this%step_name))//'Block/'//trim(adjustl(name))//'/','__Origin__',xlo_)
      call this%hdf5%WriteAttributes(trim(adjustl(this%step_name))//'Block/'//trim(adjustl(name))//'/','__Spacing__',dx_)

      return
    end subroutine h5hut_obj_WriteGrid