marker_set_ReadHDF5 Subroutine

private impure subroutine marker_set_ReadHDF5(this, iter, time)

Uses

  • proc~~marker_set_readhdf5~~UsesGraph proc~marker_set_readhdf5 marker_set%marker_set_ReadHDF5 module~leaputils leapUtils proc~marker_set_readhdf5->module~leaputils module~leapkinds leapKinds module~leaputils->module~leapkinds iso_fortran_env iso_fortran_env module~leapkinds->iso_fortran_env

Reads marker data from file in parallel using HDF5.

Type Bound

marker_set

Arguments

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

A collection of tessellation elements

integer, intent(out) :: iter

Iteration at write

real(kind=wp), intent(out) :: time

Time at write


Calls

proc~~marker_set_readhdf5~~CallsGraph proc~marker_set_readhdf5 marker_set%marker_set_ReadHDF5 c c proc~marker_set_readhdf5->c f f proc~marker_set_readhdf5->f markers markers proc~marker_set_readhdf5->markers n n proc~marker_set_readhdf5->n none~readattributes hdf5_obj%ReadAttributes proc~marker_set_readhdf5->none~readattributes p p proc~marker_set_readhdf5->p proc~hdf5_obj_final hdf5_obj%hdf5_obj_Final proc~marker_set_readhdf5->proc~hdf5_obj_final proc~hdf5_obj_init hdf5_obj%hdf5_obj_Init proc~marker_set_readhdf5->proc~hdf5_obj_init proc~lagrangian_set_communicate lagrangian_set%lagrangian_set_Communicate proc~marker_set_readhdf5->proc~lagrangian_set_communicate proc~lagrangian_set_localize lagrangian_set%lagrangian_set_Localize proc~marker_set_readhdf5->proc~lagrangian_set_localize proc~lagrangian_set_resize lagrangian_set%lagrangian_set_Resize proc~marker_set_readhdf5->proc~lagrangian_set_resize proc~stringtool_obj_removeextension stringtool_obj%stringtool_obj_RemoveExtension proc~marker_set_readhdf5->proc~stringtool_obj_removeextension v v proc~marker_set_readhdf5->v proc~hdf5_obj_readattributes0d hdf5_obj%hdf5_obj_ReadAttributes0D none~readattributes->proc~hdf5_obj_readattributes0d proc~hdf5_obj_readattributes1d hdf5_obj%hdf5_obj_ReadAttributes1D none~readattributes->proc~hdf5_obj_readattributes1d 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 proc~lagrangian_set_communicate->proc~lagrangian_set_resize mpi_gather mpi_gather proc~lagrangian_set_communicate->mpi_gather mpi_recv mpi_recv proc~lagrangian_set_communicate->mpi_recv mpi_send mpi_send proc~lagrangian_set_communicate->mpi_send proc~lagrangian_set_recycle lagrangian_set%lagrangian_set_Recycle proc~lagrangian_set_communicate->proc~lagrangian_set_recycle proc~block_obj_getowningcell block_obj%block_obj_GetOwningCell proc~lagrangian_set_localize->proc~block_obj_getowningcell 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_readattributes1d->h5aclose_f proc~hdf5_obj_readattributes1d->h5aget_space_f proc~hdf5_obj_readattributes1d->h5aget_type_f proc~hdf5_obj_readattributes1d->h5aopen_f proc~hdf5_obj_readattributes1d->h5aread_f proc~hdf5_obj_readattributes1d->h5sclose_f proc~hdf5_obj_readattributes1d->h5sget_simple_extent_dims_f proc~hdf5_obj_readattributes1d->proc~hdf5_obj_closegroup proc~hdf5_obj_readattributes1d->proc~hdf5_obj_getgroupobject proc~hdf5_obj_readattributes1d->proc~hdf5_obj_opengroup proc~lagrangian_set_recycle->proc~lagrangian_set_resize 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~~marker_set_readhdf5~~CalledByGraph proc~marker_set_readhdf5 marker_set%marker_set_ReadHDF5 proc~respart_set_readhdf5 ResPart_set%ResPart_set_ReadHDF5 proc~respart_set_readhdf5->proc~marker_set_readhdf5 proc~solid_set_readhdf5 solid_set%solid_set_ReadHDF5 proc~solid_set_readhdf5->proc~marker_set_readhdf5

Source Code

    impure subroutine marker_set_ReadHDF5(this,iter,time)
      !> Reads marker data from file in parallel using HDF5.
      use leapUtils, only: stringtool_obj
      implicit none
      class(marker_set), intent(inout) :: this                                 !! A collection of tessellation elements
      integer,           intent(out)   :: iter                                 !! Iteration at write
      real(wp),          intent(out)   :: time                                 !! Time at write
      ! Work variables
      integer,          allocatable :: buffi(:)                                !! Integer buffer
      integer(leapI8),  allocatable :: buffi8(:)                               !! Integer (kind=8) buffer
      real(wp),         allocatable :: buffr(:)                                !! Real buffer
      integer,          allocatable :: size_per_rank(:)                        !! Number of markers per rank
      integer                       :: max_chunk                               !! Max number of chunks on any MPI rank
      integer                       :: chunk_size                              !! Number of markers in chunk
      integer                       :: goff                                    !! Global offset
      integer                       :: loff                                    !! Local offset
      integer                       :: n                                       !! Iterator

      type(stringtool_obj):: stringtool
      character(len=:),      &
            allocatable   :: filename
      type(hdf5_obj)      :: hdf5

      filename = 'HDF5/'//stringtool%RemoveExtension(this%read_file)

      ! Initialize HDF5
      call hdf5%Initialize(this%parallel)

      ! Open file
      call hdf5%Open(filename,'R')

      ! Read attributes
      call hdf5%ReadAttributes("/","Time", time)
      call hdf5%ReadAttributes("/","Iter", iter)
      call hdf5%ReadAttributes("/","Count",this%count)

      ! Distribute equally the data on each mpi rank
      associate(mpi=>this%parallel)
        allocate(size_per_rank(mpi%nproc))
        do n=1,mpi%nproc
          size_per_rank(n) = int(this%count/mpi%nproc)
          if (n.le.mod(this%count,mpi%nproc)) size_per_rank(n)=size_per_rank(n)+1
        end do
      end associate

      ! Allocate buffers
      call this%parallel%Max(ceiling(size_per_rank(this%parallel%rank%mine)/real(this%R_chunk_size,wp)),max_chunk)

      loff  = 0
      chunk_size = 0
      do n=1,max_chunk
        ! Local offset
        loff    = loff + chunk_size

        ! Determine next chunk size
        if (loff.lt.size_per_rank(this%parallel%rank%mine)) then
          chunk_size = min(size_per_rank(this%parallel%rank%mine) - loff,this%R_chunk_size)

          ! Resize the data array to contain read Lagrangian objects
          call this%Resize(this%count_+chunk_size)
        else
          chunk_size = 0
        end if

        ! Global offset
        goff = sum(size_per_rank(1:this%parallel%rank%mine-1)) + loff + 1

        allocate(buffi8(chunk_size),buffi(chunk_size),buffr(chunk_size))

        ! Read the base Lagrangian object structure
        if (chunk_size.eq.0) then
          call hdf5%Read("/","id",buffi8,offset=goff)
          call hdf5%Read("/","x", buffr, offset=goff)
          call hdf5%Read("/","y", buffr, offset=goff)
          call hdf5%Read("/","z", buffr, offset=goff)
          call hdf5%Read("/","i", buffi, offset=goff)
          call hdf5%Read("/","j", buffi, offset=goff)
          call hdf5%Read("/","k", buffi, offset=goff)
          call hdf5%Read("/","s", buffi, offset=goff)
          call hdf5%Read("/","SA",buffr, offset=goff)
          call hdf5%Read("/","nx",buffr, offset=goff)
          call hdf5%Read("/","ny",buffr, offset=goff)
          call hdf5%Read("/","nz",buffr, offset=goff)
          call hdf5%Read("/","vx",buffr, offset=goff)
          call hdf5%Read("/","vy",buffr, offset=goff)
          call hdf5%Read("/","vz",buffr, offset=goff)
          call hdf5%Read("/","fx",buffr, offset=goff)
          call hdf5%Read("/","fy",buffr, offset=goff)
          call hdf5%Read("/","fz",buffr, offset=goff)
        else
          select type (markers=>this%p)
          type is (marker_obj)
            call hdf5%Read("/","id",buffi8,offset=goff); markers(this%count_-chunk_size+1:this%count_)%id   = buffi8
            call hdf5%Read("/","x", buffr, offset=goff); markers(this%count_-chunk_size+1:this%count_)%p(1) = buffr
            call hdf5%Read("/","y", buffr, offset=goff); markers(this%count_-chunk_size+1:this%count_)%p(2) = buffr
            call hdf5%Read("/","z", buffr, offset=goff); markers(this%count_-chunk_size+1:this%count_)%p(3) = buffr
            call hdf5%Read("/","i", buffi, offset=goff); markers(this%count_-chunk_size+1:this%count_)%c(1) = buffi
            call hdf5%Read("/","j", buffi, offset=goff); markers(this%count_-chunk_size+1:this%count_)%c(2) = buffi
            call hdf5%Read("/","k", buffi, offset=goff); markers(this%count_-chunk_size+1:this%count_)%c(3) = buffi
            call hdf5%Read("/","s", buffi, offset=goff); markers(this%count_-chunk_size+1:this%count_)%s    = buffi
            call hdf5%Read("/","SA",buffr, offset=goff); markers(this%count_-chunk_size+1:this%count_)%SA   = buffr
            call hdf5%Read("/","nx",buffr, offset=goff); markers(this%count_-chunk_size+1:this%count_)%n(1) = buffr
            call hdf5%Read("/","ny",buffr, offset=goff); markers(this%count_-chunk_size+1:this%count_)%n(2) = buffr
            call hdf5%Read("/","nz",buffr, offset=goff); markers(this%count_-chunk_size+1:this%count_)%n(3) = buffr
            call hdf5%Read("/","vx",buffr, offset=goff); markers(this%count_-chunk_size+1:this%count_)%v(1) = buffr
            call hdf5%Read("/","vy",buffr, offset=goff); markers(this%count_-chunk_size+1:this%count_)%v(2) = buffr
            call hdf5%Read("/","vz",buffr, offset=goff); markers(this%count_-chunk_size+1:this%count_)%v(3) = buffr
            call hdf5%Read("/","fx",buffr, offset=goff); markers(this%count_-chunk_size+1:this%count_)%f(1) = buffr
            call hdf5%Read("/","fy",buffr, offset=goff); markers(this%count_-chunk_size+1:this%count_)%f(2) = buffr
            call hdf5%Read("/","fz",buffr, offset=goff); markers(this%count_-chunk_size+1:this%count_)%f(3) = buffr
          end select
        end if

        ! Free buffers
        deallocate(buffi8,buffi,buffr)

        ! Send markers to their mpiranks
        call this%Communicate()
      end do

      ! Free remaining arrays, close file and finalize HDF5
      deallocate(size_per_rank)
      call hdf5%Close()
      call hdf5%Finalize()

      ! Set unread variables to zero
      select type (markers=>this%p)
      class is (marker_obj)
        do n=1,this%count_
          markers(n)%pold  = 0.0_wp
          markers(n)%vold  = 0.0_wp
        end do
      end select

      ! Localize markers on the grid
      call this%Localize()

      return
    end subroutine marker_set_ReadHDF5