block_obj_Write Subroutine

private impure subroutine block_obj_Write(this, name)

Writes block data using HDF5.

Type Bound

block_obj

Arguments

Type IntentOptional Attributes Name
class(block_obj), intent(in) :: this

A block object

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

Name of file to write


Calls

proc~~block_obj_write~~CallsGraph proc~block_obj_write block_obj%block_obj_Write none~writeattributes hdf5_obj%WriteAttributes proc~block_obj_write->none~writeattributes proc~hdf5_obj_final hdf5_obj%hdf5_obj_Final proc~block_obj_write->proc~hdf5_obj_final proc~hdf5_obj_init hdf5_obj%hdf5_obj_Init proc~block_obj_write->proc~hdf5_obj_init proc~hdf5_obj_writecoord hdf5_obj%hdf5_obj_WriteCoord proc~block_obj_write->proc~hdf5_obj_writecoord 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 h5garbage_collect_f h5garbage_collect_f proc~hdf5_obj_final->h5garbage_collect_f proc~hashtbl_obj_final hashtbl_obj%hashtbl_obj_Final proc~hdf5_obj_final->proc~hashtbl_obj_final h5open_f h5open_f proc~hdf5_obj_init->h5open_f proc~hashtbl_obj_init hashtbl_obj%hashtbl_obj_Init proc~hdf5_obj_init->proc~hashtbl_obj_init h5dclose_f h5dclose_f proc~hdf5_obj_writecoord->h5dclose_f h5dcreate_f h5dcreate_f proc~hdf5_obj_writecoord->h5dcreate_f h5dwrite_f h5dwrite_f proc~hdf5_obj_writecoord->h5dwrite_f h5pclose_f h5pclose_f proc~hdf5_obj_writecoord->h5pclose_f h5pcreate_f h5pcreate_f proc~hdf5_obj_writecoord->h5pcreate_f h5pset_dxpl_mpio_f h5pset_dxpl_mpio_f proc~hdf5_obj_writecoord->h5pset_dxpl_mpio_f h5sclose_f h5sclose_f proc~hdf5_obj_writecoord->h5sclose_f h5screate_simple_f h5screate_simple_f proc~hdf5_obj_writecoord->h5screate_simple_f h5sselect_hyperslab_f h5sselect_hyperslab_f proc~hdf5_obj_writecoord->h5sselect_hyperslab_f proc~hdf5_obj_closegroup hdf5_obj%hdf5_obj_CloseGroup proc~hdf5_obj_writecoord->proc~hdf5_obj_closegroup proc~hdf5_obj_getgroupobject hdf5_obj%hdf5_obj_GetGroupObject proc~hdf5_obj_writecoord->proc~hdf5_obj_getgroupobject proc~hdf5_obj_opengroup hdf5_obj%hdf5_obj_OpenGroup proc~hdf5_obj_writecoord->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~hdf5_obj_writeattributes0d->h5sclose_f proc~hdf5_obj_writeattributes0d->h5screate_simple_f proc~hdf5_obj_writeattributes0d->proc~hdf5_obj_closegroup proc~hdf5_obj_writeattributes0d->proc~hdf5_obj_getgroupobject proc~hdf5_obj_writeattributes0d->proc~hdf5_obj_opengroup 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 proc~hdf5_obj_writeattributes1d->h5sclose_f proc~hdf5_obj_writeattributes1d->h5screate_simple_f proc~hdf5_obj_writeattributes1d->proc~hdf5_obj_closegroup proc~hdf5_obj_writeattributes1d->proc~hdf5_obj_getgroupobject proc~hdf5_obj_writeattributes1d->proc~hdf5_obj_opengroup proc~hdf5_obj_writeattributes1d->h5aclose_f proc~hdf5_obj_writeattributes1d->h5acreate_f proc~hdf5_obj_writeattributes1d->h5awrite_f 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 block_obj_Write(this,name)
      !> Writes block data using HDF5.
      implicit none
      class(block_obj), intent(in) :: this                                     !! A block object
      character(len=*), intent(in) :: name                                     !! Name of file to write
      ! Work variables
      type(hdf5_obj), &
              pointer :: hdf5
      type(hdf5_obj), &
              target  :: hdf5t
      integer         :: Ng(3)
      integer         :: maxhi(3)
      integer         :: minlo(3)
      character(len=2):: label
      integer         :: nsize
      integer         :: dir
      integer         :: Nblock(3)
      real(wp),       &
          allocatable :: xg(:)
      real(wp),     &
          allocatable :: buff(:)

      ! Use internal pointer, if associated.
      ! Otherwise, use a new object.
      if (associated(this%hdf5)) then
        hdf5 => this%hdf5
      else
        call hdf5t%Initialize(this%parallel)
        hdf5 => hdf5t
      end if

      ! Open file
      call hdf5%Open(trim(adjustl(name)),'W')

      associate (lo => this%lo, hi => this%hi, x => this%x, y => this%y, &
         z =>this%z, periods => this%periods, rank => this%parallel%rank,&
         np => this%parallel%np )

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

        ! Total number of grid points
        Ng = maxhi-minlo+1

        ! Write header info
        call hdf5%WriteAttributes('/','Periodicity',periods)
        call hdf5%WriteAttributes('/','Ng',         Ng     )

        do dir=1,3
          ! Build global coordinate array
          allocate(xg  (minlo(dir):maxhi(dir)+1))
          allocate(buff(Ng(dir)+1))

          xg = 0.0_wp
          xg(lo(dir):hi(dir)) = this%axis(dir)%x(lo(dir):hi(dir))

          if (rank%dir(dir).eq.np(dir)) xg(hi(dir)+1) = this%axis(dir)%x(hi(dir)+1)

          Nblock      = np
          Nblock(dir) = 1
          nsize       = Nblock(1)*Nblock(2)*Nblock(3)

          call this%parallel%Sum(xg,buff); xg=buff/nsize

          ! Form label for this axis: x1, x2, x3
          write(label,fmt='(a,i1)') 'x',dir

          ! Write coordinates along this axis
          call hdf5%WriteCoord('/', label, xg)

          deallocate(xg)
          deallocate(buff)
        end do

      end associate

      ! Close file
      call hdf5%Close()

      ! Free data
      if (associated(this%hdf5)) then
        hdf5 => null()
      else
        call hdf5t%Finalize()
        hdf5 => null()
      end if

      return
    end subroutine block_obj_Write