eulerian_set_WriteSingleH5HUT Subroutine

private impure subroutine eulerian_set_WriteSingleH5HUT(this, h5, ind)

Write a single Eulerian object to file using H5HUT.

Type Bound

eulerian_set

Arguments

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

An Eulerian Set

type(h5hut_obj), intent(inout) :: h5

H5hut structure

integer, intent(in) :: ind

Index of Eulerian object


Calls

proc~~eulerian_set_writesingleh5hut~~CallsGraph proc~eulerian_set_writesingleh5hut eulerian_set%eulerian_set_WriteSingleH5HUT cell cell proc~eulerian_set_writesingleh5hut->cell proc~h5hut_obj_writegrid h5hut_obj%h5hut_obj_WriteGrid proc~eulerian_set_writesingleh5hut->proc~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~~eulerian_set_writesingleh5hut~~CalledByGraph proc~eulerian_set_writesingleh5hut eulerian_set%eulerian_set_WriteSingleH5HUT 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 eulerian_set_WriteSingleH5HUT(this,h5,ind)
      !> Write a single Eulerian object to file using H5HUT.
      implicit none
      class(eulerian_set), intent(inout) :: this                               !! An Eulerian Set
      type(h5hut_obj),     intent(inout) :: h5                                 !! H5hut structure
      integer,             intent(in)    :: ind                                !! Index of Eulerian object
      ! Work variables
      real(wp),allocatable :: buff_r(:,:,:)                                    !! Buffer to write
      integer, allocatable :: buff_i(:,:,:)                                    !! Buffer to write
      real(wp)        :: shift(3)                                              !! Shift to account for staggering
      real(wp)        :: origin(3)                                             !! Absolute origin


      ! Get grid origin and displacement
      origin=this%block%pmin

      ! Write grid attributes
      shift=0.5_wp*this%block%dx
      select case(this%field(ind)%p%staggering)
      case (0)
        ! cell centered
        ! Do nothing
      case (1)
        ! X1-face centered
        shift(1)=0.0_wp
      case (2)
        ! X2-face centered
        shift(2)=0.0_wp
      case (3)
        ! X3-face centered
        shift(3)=0.0_wp
      end select

      associate (lo => this%block%lo, hi => this%block%hi, ngc=> this%block%ngc)
        select type (my_field=>this%field(ind)%p)
        type is (eulerian_obj_r)
          ! Allocate buffer
          allocate(buff_r(lo(1):hi(1), lo(2):hi(2), lo(3):hi(3)))
          ! Write data
          buff_r = my_field%cell(lo(1):hi(1), lo(2):hi(2), lo(3):hi(3))
          call h5%Write(my_field%name,buff_r,lo,hi)
          deallocate(buff_r)
        type is (eulerian_obj_i)
          ! Allocate buffer
          allocate(buff_i(lo(1):hi(1), lo(2):hi(2), lo(3):hi(3)))
          ! Write data
          buff_i = my_field%cell(lo(1):hi(1), lo(2):hi(2), lo(3):hi(3))
          call h5%Write(my_field%name,buff_i,lo,hi)
          deallocate(buff_i)
        end select
      end associate

      ! Write grid attributes
      call h5%WriteGrid(this%field(ind)%p%name,origin+shift,this%block%dx)

      return
    end subroutine eulerian_set_WriteSingleH5HUT