ResPart_set_ReadH5HUT Subroutine

private impure subroutine ResPart_set_ReadH5HUT(this, iter, time, step)

Reads ResPart data from file using H5HUT.

Type Bound

ResPart_set

Arguments

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

Set of resolved particles

integer, intent(out) :: iter

Iteration at write

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

Time at write

integer, intent(in), optional :: step

User supplied step to open


Calls

proc~~respart_set_readh5hut~~CallsGraph proc~respart_set_readh5hut ResPart_set%ResPart_set_ReadH5HUT c c proc~respart_set_readh5hut->c fc fc proc~respart_set_readh5hut->fc fh fh proc~respart_set_readh5hut->fh p p proc~respart_set_readh5hut->p particles particles proc~respart_set_readh5hut->particles proc~h5hut_obj_final h5hut_obj%h5hut_obj_Final proc~respart_set_readh5hut->proc~h5hut_obj_final proc~h5hut_obj_getnpoints h5hut_obj%h5hut_obj_GetNPoints proc~respart_set_readh5hut->proc~h5hut_obj_getnpoints proc~h5hut_obj_init h5hut_obj%h5hut_obj_Init proc~respart_set_readh5hut->proc~h5hut_obj_init proc~h5hut_obj_jumptostep h5hut_obj%h5hut_obj_JumpToStep proc~respart_set_readh5hut->proc~h5hut_obj_jumptostep proc~h5hut_obj_lasttimestep h5hut_obj%h5hut_obj_LastTimeStep proc~respart_set_readh5hut->proc~h5hut_obj_lasttimestep proc~lagrangian_set_communicate lagrangian_set%lagrangian_set_Communicate proc~respart_set_readh5hut->proc~lagrangian_set_communicate proc~lagrangian_set_localize lagrangian_set%lagrangian_set_Localize proc~respart_set_readh5hut->proc~lagrangian_set_localize proc~lagrangian_set_resize lagrangian_set%lagrangian_set_Resize proc~respart_set_readh5hut->proc~lagrangian_set_resize proc~respart_set_updatelookup ResPart_set%ResPart_set_UpdateLookup proc~respart_set_readh5hut->proc~respart_set_updatelookup tc tc proc~respart_set_readh5hut->tc th th proc~respart_set_readh5hut->th v v proc~respart_set_readh5hut->v w w proc~respart_set_readh5hut->w proc~hdf5_obj_final hdf5_obj%hdf5_obj_Final proc~h5hut_obj_final->proc~hdf5_obj_final proc~hdf5_obj_getnpoints hdf5_obj%hdf5_obj_GetNPoints proc~h5hut_obj_getnpoints->proc~hdf5_obj_getnpoints proc~hdf5_obj_readdatasetnames hdf5_obj%hdf5_obj_ReadDatasetNames proc~h5hut_obj_getnpoints->proc~hdf5_obj_readdatasetnames proc~hdf5_obj_init hdf5_obj%hdf5_obj_Init proc~h5hut_obj_init->proc~hdf5_obj_init none~readattributes hdf5_obj%ReadAttributes proc~h5hut_obj_jumptostep->none~readattributes proc~h5hut_obj_setstep h5hut_obj%h5hut_obj_SetStep proc~h5hut_obj_jumptostep->proc~h5hut_obj_setstep proc~h5hut_obj_lasttimestep->none~readattributes proc~h5hut_obj_getnsteps h5hut_obj%h5hut_obj_GetNSteps proc~h5hut_obj_lasttimestep->proc~h5hut_obj_getnsteps proc~h5hut_obj_lasttimestep->proc~h5hut_obj_setstep 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 proc~respart_set_updatelookup->particles 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 proc~hdf5_obj_readgroupnames hdf5_obj%hdf5_obj_ReadGroupNames proc~h5hut_obj_getnsteps->proc~hdf5_obj_readgroupnames proc~hdf5_obj_creategroup hdf5_obj%hdf5_obj_CreateGroup proc~h5hut_obj_setstep->proc~hdf5_obj_creategroup 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 h5dclose_f h5dclose_f proc~hdf5_obj_getnpoints->h5dclose_f h5dget_space_f h5dget_space_f proc~hdf5_obj_getnpoints->h5dget_space_f h5dopen_f h5dopen_f proc~hdf5_obj_getnpoints->h5dopen_f h5sclose_f h5sclose_f proc~hdf5_obj_getnpoints->h5sclose_f h5sget_simple_extent_npoints_f h5sget_simple_extent_npoints_f proc~hdf5_obj_getnpoints->h5sget_simple_extent_npoints_f proc~hdf5_obj_closegroup hdf5_obj%hdf5_obj_CloseGroup proc~hdf5_obj_getnpoints->proc~hdf5_obj_closegroup proc~hdf5_obj_getgroupobject hdf5_obj%hdf5_obj_GetGroupObject proc~hdf5_obj_getnpoints->proc~hdf5_obj_getgroupobject proc~hdf5_obj_opengroup hdf5_obj%hdf5_obj_OpenGroup proc~hdf5_obj_getnpoints->proc~hdf5_obj_opengroup 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 h5gn_members_f h5gn_members_f proc~hdf5_obj_readdatasetnames->h5gn_members_f h5iget_name_f h5iget_name_f proc~hdf5_obj_readdatasetnames->h5iget_name_f h5oclose_f h5oclose_f proc~hdf5_obj_readdatasetnames->h5oclose_f h5oget_info_f h5oget_info_f proc~hdf5_obj_readdatasetnames->h5oget_info_f h5oopen_by_idx_f h5oopen_by_idx_f proc~hdf5_obj_readdatasetnames->h5oopen_by_idx_f proc~hdf5_obj_fixgroupname hdf5_obj%hdf5_obj_FixGroupName proc~hdf5_obj_readdatasetnames->proc~hdf5_obj_fixgroupname proc~lagrangian_set_recycle->proc~lagrangian_set_resize proc~hdf5_obj_closegroup->proc~hdf5_obj_fixgroupname 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_creategroup->proc~hdf5_obj_closegroup proc~hdf5_obj_creategroup->proc~hdf5_obj_fixgroupname h5gcreate_f h5gcreate_f proc~hdf5_obj_creategroup->h5gcreate_f 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_getgroupobject->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_opengroup->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~hdf5_obj_opengroup->proc~hashtbl_obj_put proc~hdf5_obj_readattributes0d->h5sclose_f proc~hdf5_obj_readattributes0d->proc~hdf5_obj_closegroup proc~hdf5_obj_readattributes0d->proc~hdf5_obj_getgroupobject proc~hdf5_obj_readattributes0d->proc~hdf5_obj_opengroup 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 h5sget_simple_extent_dims_f h5sget_simple_extent_dims_f proc~hdf5_obj_readattributes0d->h5sget_simple_extent_dims_f proc~hdf5_obj_readattributes1d->h5sclose_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_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->h5sget_simple_extent_dims_f proc~hdf5_obj_readgroupnames->h5gn_members_f proc~hdf5_obj_readgroupnames->h5iget_name_f proc~hdf5_obj_readgroupnames->h5oclose_f proc~hdf5_obj_readgroupnames->h5oget_info_f proc~hdf5_obj_readgroupnames->h5oopen_by_idx_f proc~hdf5_obj_readgroupnames->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 ResPart_set_ReadH5HUT(this,iter,time,step)
      !> Reads ResPart data from file using H5HUT.
      implicit none
      class(ResPart_set), intent(inout) :: this                                !! Set of resolved particles
      integer,            intent(out)   :: iter                                !! Iteration at write
      real(wp),           intent(out)   :: time                                !! Time at write
      integer,            intent(in),    &
                               optional :: step                                !! User supplied step to open
      ! Work variables
      type(h5hut_obj) :: h5
      integer,        allocatable :: buffi(:)
      integer(leapI8),allocatable :: buffi8(:)
      real(wp),       allocatable :: buffr(:)
      integer         :: n

      ! Open the file
      call h5%initialize(trim(adjustl(this%read_file)),"R",this%parallel)

      if (.not.present(step)) then
        ! Jump to last step
        call h5%LastTimeStep(iter,time)
      else
        ! Jump to desired step
        call h5%JumpToStep(step,iter,time)
      end if

      ! Get total number of point objects
      call h5%getnpoints(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 h5%Read("id",  buffi8)
        call h5%Read("x",   buffr )
        call h5%Read("y",   buffr )
        call h5%Read("z",   buffr )
        call h5%Read("i",   buffi )
        call h5%Read("j",   buffi )
        call h5%Read("k",   buffi )
        call h5%Read("u",   buffr )
        call h5%Read("v",   buffr )
        call h5%Read("w",   buffr )
        call h5%Read("wx",  buffr )
        call h5%Read("wy",  buffr )
        call h5%Read("wz",  buffr )
        call h5%Read("d",   buffr )
        call h5%Read("s",   buffi )
        call h5%Read("rho", buffr )
        call h5%Read("Fhx", buffr )
        call h5%Read("Fhy", buffr )
        call h5%Read("Fhz", buffr )
        call h5%Read("Thx", buffr )
        call h5%Read("Thy", buffr )
        call h5%Read("Thz", buffr )
        call h5%Read("Fcx", buffr )
        call h5%Read("Fcy", buffr )
        call h5%Read("Fcz", buffr )
        call h5%Read("Tcx", buffr )
        call h5%Read("Tcy", buffr )
        call h5%Read("Tcz", buffr )
      else
        select type (particles=>this%p)
        type is (ResPart_obj)
          call h5%Read("id",  buffi8); particles(:)%id   = buffi8
          call h5%Read("x",   buffr ); particles(:)%p(1) = buffr
          call h5%Read("y",   buffr ); particles(:)%p(2) = buffr
          call h5%Read("z",   buffr ); particles(:)%p(3) = buffr
          call h5%Read("i",   buffi ); particles(:)%c(1) = buffi
          call h5%Read("j",   buffi ); particles(:)%c(2) = buffi
          call h5%Read("k",   buffi ); particles(:)%c(3) = buffi
          call h5%Read("u",   buffr ); particles(:)%v(1) = buffr
          call h5%Read("v",   buffr ); particles(:)%v(2) = buffr
          call h5%Read("w",   buffr ); particles(:)%v(3) = buffr
          call h5%Read("wx",  buffr ); particles(:)%w(1) = buffr
          call h5%Read("wy",  buffr ); particles(:)%w(2) = buffr
          call h5%Read("wz",  buffr ); particles(:)%w(3) = buffr
          call h5%Read("d",   buffr ); particles(:)%d    = buffr
          call h5%Read("s",   buffi ); particles(:)%s    = buffi
          call h5%Read("rho", buffr ); particles(:)%rho  = buffr
          call h5%Read("Fhx", buffr ); particles(:)%Fh(1)= buffr
          call h5%Read("Fhy", buffr ); particles(:)%Fh(2)= buffr
          call h5%Read("Fhz", buffr ); particles(:)%Fh(3)= buffr
          call h5%Read("Thx", buffr ); particles(:)%Th(1)= buffr
          call h5%Read("Thy", buffr ); particles(:)%Th(2)= buffr
          call h5%Read("Thz", buffr ); particles(:)%Th(3)= buffr
          call h5%Read("Fcx", buffr ); particles(:)%Fc(1)= buffr
          call h5%Read("Fcy", buffr ); particles(:)%Fc(2)= buffr
          call h5%Read("Fcz", buffr ); particles(:)%Fc(3)= buffr
          call h5%Read("Tcx", buffr ); particles(:)%Tc(1)= buffr
          call h5%Read("Tcy", buffr ); particles(:)%Tc(2)= buffr
          call h5%Read("Tcz", buffr ); particles(:)%Tc(3)= buffr
        end select
      end if

      ! Free data
      deallocate(buffi8,buffi,buffr)

      call h5%Finalize

      ! Set unread variables to zero
      select type (particles=>this%p)
      class is (ResPart_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

      ! Send Lagrangian centroids to their mpiranks
      call this%Communicate

      ! Localize markers on the grid
      call this%Localize

      ! Read the associated IB data
      if (.not.present(step)) then
        call this%ib%Read(iter,time)
      else
        call this%ib%Read(iter,time,step)
      end if

      ! Update lookup table
      call this%UpdateLookup

      return
    end subroutine ResPart_set_ReadH5HUT