block_obj_Read Subroutine

private impure subroutine block_obj_Read(this, name)

Reads block data using HDF5.

Type Bound

block_obj

Arguments

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

A block object

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

Name of file to write


Calls

proc~~block_obj_read~~CallsGraph proc~block_obj_read block_obj%block_obj_Read none~readattributes hdf5_obj%ReadAttributes proc~block_obj_read->none~readattributes proc~block_obj_setconveniencepointers block_obj%block_obj_SetConveniencePointers proc~block_obj_read->proc~block_obj_setconveniencepointers proc~block_obj_setupmpitypes block_obj%block_obj_SetupMPITypes proc~block_obj_read->proc~block_obj_setupmpitypes proc~block_obj_updategridghostcells block_obj%block_obj_UpdateGridGhostCells proc~block_obj_read->proc~block_obj_updategridghostcells proc~block_obj_updatemidpoints block_obj%block_obj_UpdateMidPoints proc~block_obj_read->proc~block_obj_updatemidpoints proc~block_obj_updatespacing block_obj%block_obj_UpdateSpacing proc~block_obj_read->proc~block_obj_updatespacing proc~hdf5_obj_final hdf5_obj%hdf5_obj_Final proc~block_obj_read->proc~hdf5_obj_final proc~hdf5_obj_init hdf5_obj%hdf5_obj_Init proc~block_obj_read->proc~hdf5_obj_init proc~hdf5_obj_readcoord hdf5_obj%hdf5_obj_ReadCoord proc~block_obj_read->proc~hdf5_obj_readcoord 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 mpi_type_commit mpi_type_commit proc~block_obj_setupmpitypes->mpi_type_commit mpi_type_free mpi_type_free proc~block_obj_setupmpitypes->mpi_type_free mpi_type_vector mpi_type_vector proc~block_obj_setupmpitypes->mpi_type_vector mpi_irecv mpi_irecv proc~block_obj_updategridghostcells->mpi_irecv mpi_isend mpi_isend proc~block_obj_updategridghostcells->mpi_isend mpi_wait mpi_wait proc~block_obj_updategridghostcells->mpi_wait mpi_waitall mpi_waitall proc~block_obj_updategridghostcells->mpi_waitall proc~block_obj_updateextents block_obj%block_obj_UpdateExtents proc~block_obj_updategridghostcells->proc~block_obj_updateextents 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_readcoord->h5dclose_f h5dget_space_f h5dget_space_f proc~hdf5_obj_readcoord->h5dget_space_f h5dget_type_f h5dget_type_f proc~hdf5_obj_readcoord->h5dget_type_f h5dopen_f h5dopen_f proc~hdf5_obj_readcoord->h5dopen_f h5dread_f h5dread_f proc~hdf5_obj_readcoord->h5dread_f h5sclose_f h5sclose_f proc~hdf5_obj_readcoord->h5sclose_f h5screate_simple_f h5screate_simple_f proc~hdf5_obj_readcoord->h5screate_simple_f h5sselect_hyperslab_f h5sselect_hyperslab_f proc~hdf5_obj_readcoord->h5sselect_hyperslab_f h5tclose_f h5tclose_f proc~hdf5_obj_readcoord->h5tclose_f proc~hdf5_obj_closegroup hdf5_obj%hdf5_obj_CloseGroup proc~hdf5_obj_readcoord->proc~hdf5_obj_closegroup proc~hdf5_obj_getgroupobject hdf5_obj%hdf5_obj_GetGroupObject proc~hdf5_obj_readcoord->proc~hdf5_obj_getgroupobject proc~hdf5_obj_opengroup hdf5_obj%hdf5_obj_OpenGroup proc~hdf5_obj_readcoord->proc~hdf5_obj_opengroup proc~axis_obj_init axis_obj%axis_obj_Init proc~block_obj_updateextents->proc~axis_obj_init 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_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~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_Read(this,name)
      !> Reads block data using HDF5.
      implicit none
      class(block_obj), intent(inout) :: 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         :: dir
      character(len=2):: label

      ! 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)),'R')

      associate (lo => this%lo,  hi => this%hi, periods => this%periods)
        ! Get periodicity and number of grid points
        call hdf5%ReadAttributes('/','Periodicity',periods)
        call hdf5%ReadAttributes('/','Ng',         Ng     )

        ! Set bounds
        lo = 1
        hi = Ng

        ! Initialize axes
        do dir=1,3
          call this%axis(dir)%Initialize(this%lo(dir),this%hi(dir),this%ngc)
        end do

        ! Associate pointers
        call this%SetConveniencePointers

        ! Read coordinates from file
        do dir=1,3
          ! Form label for this axis: x1, x2, x3
          write(label,fmt='(a,i1)') 'x',dir

          ! Read coordinates
          call hdf5%ReadCoord('/',label,this%axis(dir)%x(lo(dir):hi(dir)+1))
        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

      ! Update values within ghostcells
      call this%UpdateGridGhostCells()

      ! Update mid points (xm)
      call this%UpdateMidPoints()

      ! Update spacing (dxm)
      call this%UpdateSpacing()

      ! Create MPI type for ghostcell communication
      call this%SetupMPITypes()

      return
    end subroutine block_obj_Read