particle_set_ReadHDF5 Subroutine

private impure subroutine particle_set_ReadHDF5(this, iter, time)

Uses

  • proc~~particle_set_readhdf5~~UsesGraph proc~particle_set_readhdf5 particle_set%particle_set_ReadHDF5 module~leaputils leapUtils proc~particle_set_readhdf5->module~leaputils module~leapkinds leapKinds module~leaputils->module~leapkinds iso_fortran_env iso_fortran_env module~leapkinds->iso_fortran_env

Reads particle data from file in parallel using HDF5.

Type Bound

particle_set

Arguments

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

Set of resolved particles

integer, intent(out) :: iter

Iteration at write

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

Time at write


Calls

proc~~particle_set_readhdf5~~CallsGraph proc~particle_set_readhdf5 particle_set%particle_set_ReadHDF5 c c proc~particle_set_readhdf5->c fc fc proc~particle_set_readhdf5->fc fh fh proc~particle_set_readhdf5->fh none~readattributes hdf5_obj%ReadAttributes proc~particle_set_readhdf5->none~readattributes p p proc~particle_set_readhdf5->p particles particles proc~particle_set_readhdf5->particles proc~hdf5_obj_final hdf5_obj%hdf5_obj_Final proc~particle_set_readhdf5->proc~hdf5_obj_final proc~hdf5_obj_init hdf5_obj%hdf5_obj_Init proc~particle_set_readhdf5->proc~hdf5_obj_init proc~lagrangian_set_resize lagrangian_set%lagrangian_set_Resize proc~particle_set_readhdf5->proc~lagrangian_set_resize proc~stringtool_obj_removeextension stringtool_obj%stringtool_obj_RemoveExtension proc~particle_set_readhdf5->proc~stringtool_obj_removeextension tc tc proc~particle_set_readhdf5->tc th th proc~particle_set_readhdf5->th v v proc~particle_set_readhdf5->v w w proc~particle_set_readhdf5->w 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 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~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

Source Code

    impure subroutine particle_set_ReadHDF5(this,iter,time)
      !> Reads particle data from file in parallel using HDF5.
      use leapUtils, only: stringtool_obj
      implicit none
      class(particle_set), intent(inout) :: this                               !! Set of resolved particles
      integer,             intent(out)   :: iter                               !! Iteration at write
      real(wp),            intent(out)   :: time                               !! Time at write
      ! Work variables
      type(stringtool_obj):: stringtool
      character(len=:),      &
            allocatable   :: filename
      type(hdf5_obj)      :: hdf5
      integer,        allocatable :: buffi(:)
      integer(leapI8),allocatable :: buffi8(:)
      real(wp),       allocatable :: buffr(:)
      integer              :: n

      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)
        this%count_=int(this%count/mpi%nproc)
        if (mpi%rank%mine.le.mod(this%count,mpi%nproc)) this%count_=this%count_+1
      end associate

      ! Resize the data array to contain read Lagrangian objects
      call this%Resize(this%count_)

      ! Allocate buffers
      allocate(buffi8(this%count_),buffi(this%count_),buffr(this%count_))

      ! Read the base Lagrangian object structure
      if (this%count_.eq.0) then
        call hdf5%Read("/","id", buffi8)
        call hdf5%Read("/","x",  buffr )
        call hdf5%Read("/","y",  buffr )
        call hdf5%Read("/","z",  buffr )
        call hdf5%Read("/","i",  buffi )
        call hdf5%Read("/","j",  buffi )
        call hdf5%Read("/","k",  buffi )
        call hdf5%Read("/","vx", buffr )
        call hdf5%Read("/","vy", buffr )
        call hdf5%Read("/","vz", buffr )
        call hdf5%Read("/","wx", buffr )
        call hdf5%Read("/","wy", buffr )
        call hdf5%Read("/","wz", buffr )
        call hdf5%Read("/","d",  buffr )
        call hdf5%Read("/","s",  buffi )
        call hdf5%Read("/","rho",buffr )
        call hdf5%Read("/","Fhx",buffr )
        call hdf5%Read("/","Fhy",buffr )
        call hdf5%Read("/","Fhz",buffr )
        call hdf5%Read("/","Thx",buffr )
        call hdf5%Read("/","Thy",buffr )
        call hdf5%Read("/","Thz",buffr )
        call hdf5%Read("/","Fcx",buffr )
        call hdf5%Read("/","Fcy",buffr )
        call hdf5%Read("/","Fcz",buffr )
        call hdf5%Read("/","Tcx",buffr )
        call hdf5%Read("/","Tcy",buffr )
        call hdf5%Read("/","Tcz",buffr )
      else
        select type (particles=>this%p)
        type is (particle_obj)
          call hdf5%Read("/","id", buffi8); particles(:)%id   = buffi8
          call hdf5%Read("/","x",  buffr ); particles(:)%p(1) = buffr
          call hdf5%Read("/","y",  buffr ); particles(:)%p(2) = buffr
          call hdf5%Read("/","z",  buffr ); particles(:)%p(3) = buffr
          call hdf5%Read("/","i",  buffi ); particles(:)%c(1) = buffi
          call hdf5%Read("/","j",  buffi ); particles(:)%c(2) = buffi
          call hdf5%Read("/","k",  buffi ); particles(:)%c(3) = buffi
          call hdf5%Read("/","vx", buffr ); particles(:)%v(1) = buffr
          call hdf5%Read("/","vy", buffr ); particles(:)%v(2) = buffr
          call hdf5%Read("/","vz", buffr ); particles(:)%v(3) = buffr
          call hdf5%Read("/","wx", buffr ); particles(:)%w(1) = buffr
          call hdf5%Read("/","wy", buffr ); particles(:)%w(2) = buffr
          call hdf5%Read("/","wz", buffr ); particles(:)%w(3) = buffr
          call hdf5%Read("/","d",  buffr ); particles(:)%d    = buffr
          call hdf5%Read("/","s",  buffi ); particles(:)%s    = buffi
          call hdf5%Read("/","rho",buffr ); particles(:)%rho  = buffr
          call hdf5%Read("/","Fhx",buffr ); particles(:)%Fh(1)= buffr
          call hdf5%Read("/","Fhy",buffr ); particles(:)%Fh(2)= buffr
          call hdf5%Read("/","Fhz",buffr ); particles(:)%Fh(3)= buffr
          call hdf5%Read("/","Thx",buffr ); particles(:)%Th(1)= buffr
          call hdf5%Read("/","Thy",buffr ); particles(:)%Th(2)= buffr
          call hdf5%Read("/","Thz",buffr ); particles(:)%Th(3)= buffr
          call hdf5%Read("/","Fcx",buffr ); particles(:)%Fc(1)= buffr
          call hdf5%Read("/","Fcy",buffr ); particles(:)%Fc(2)= buffr
          call hdf5%Read("/","Fcz",buffr ); particles(:)%Fc(3)= buffr
          call hdf5%Read("/","Tcx",buffr ); particles(:)%Tc(1)= buffr
          call hdf5%Read("/","Tcy",buffr ); particles(:)%Tc(2)= buffr
          call hdf5%Read("/","Tcz",buffr ); particles(:)%Tc(3)= buffr
        end select
      end if

      deallocate(buffi8,buffi,buffr)

      ! Close file
      call hdf5%Close()

      ! Finalize HDF5
      call hdf5%Finalize()

      ! Set unread variables to zero
      select type (particles=>this%p)
      class is (particle_obj)
        do n=1,this%count_
          particles(n)%pold  = 0.0_wp
          particles(n)%vold  = 0.0_wp
          particles(n)%wold  = 0.0_wp
          particles(n)%Fhold = 0.0_wp
          particles(n)%Thold = 0.0_wp
          particles(n)%Fcold = 0.0_wp
          particles(n)%Tcold = 0.0_wp
        end do
      end select

      return
    end subroutine particle_set_ReadHDF5