eulerian_set_ReadSingleHDF5 Subroutine

private impure subroutine eulerian_set_ReadSingleHDF5(this, hdf5, ind)

Reads one Eulerian object based on name using LEAP's HDF5.

Type Bound

eulerian_set

Arguments

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

An Eulerian Set

type(hdf5_obj), intent(inout) :: hdf5
integer, intent(in) :: ind

Index of Eulerian object


Calls

proc~~eulerian_set_readsinglehdf5~~CallsGraph proc~eulerian_set_readsinglehdf5 eulerian_set%eulerian_set_ReadSingleHDF5 cell cell proc~eulerian_set_readsinglehdf5->cell proc~eulerian_obj_updateghostcells eulerian_obj_base%eulerian_obj_UpdateGhostCells proc~eulerian_set_readsinglehdf5->proc~eulerian_obj_updateghostcells proc~eulerian_obj_updateghostcells_x eulerian_obj_base%eulerian_obj_UpdateGhostCells_x proc~eulerian_obj_updateghostcells->proc~eulerian_obj_updateghostcells_x proc~eulerian_obj_updateghostcells_y eulerian_obj_base%eulerian_obj_UpdateGhostCells_y proc~eulerian_obj_updateghostcells->proc~eulerian_obj_updateghostcells_y proc~eulerian_obj_updateghostcells_z eulerian_obj_base%eulerian_obj_UpdateGhostCells_z proc~eulerian_obj_updateghostcells->proc~eulerian_obj_updateghostcells_z proc~eulerian_obj_updateghostcells_x->cell mpi_irecv mpi_irecv proc~eulerian_obj_updateghostcells_x->mpi_irecv mpi_isend mpi_isend proc~eulerian_obj_updateghostcells_x->mpi_isend mpi_waitall mpi_waitall proc~eulerian_obj_updateghostcells_x->mpi_waitall proc~eulerian_obj_updateghostcells_y->cell proc~eulerian_obj_updateghostcells_y->mpi_irecv proc~eulerian_obj_updateghostcells_y->mpi_isend proc~eulerian_obj_updateghostcells_y->mpi_waitall proc~eulerian_obj_updateghostcells_z->cell proc~eulerian_obj_updateghostcells_z->mpi_irecv proc~eulerian_obj_updateghostcells_z->mpi_isend proc~eulerian_obj_updateghostcells_z->mpi_waitall

Called by

proc~~eulerian_set_readsinglehdf5~~CalledByGraph proc~eulerian_set_readsinglehdf5 eulerian_set%eulerian_set_ReadSingleHDF5 none~readsingle eulerian_set%ReadSingle none~readsingle->proc~eulerian_set_readsinglehdf5 proc~eulerian_set_readh5hut eulerian_set%eulerian_set_ReadH5HUT proc~eulerian_set_readh5hut->none~readsingle proc~eulerian_set_readhdf5 eulerian_set%eulerian_set_ReadHDF5 proc~eulerian_set_readhdf5->none~readsingle

Source Code

    impure subroutine eulerian_set_ReadSingleHDF5(this,hdf5,ind)
      !> Reads one Eulerian object based on name using LEAP's HDF5.
      implicit none
      class(eulerian_set), intent(inout) :: this                               !! An Eulerian Set
      type(hdf5_obj),      intent(inout) :: hdf5
      integer,             intent(in)    :: ind                                !! Index of Eulerian object
      ! Work variable
      real(wp),allocatable :: buff_r(:,:,:)                                    !! Buffer to write
      integer, allocatable :: buff_i(:,:,:)                                    !! Buffer to write

      associate (lo => this%block%lo, hi => this%block%hi)
        select type (my_field => this%field(ind)%p)
        type is (eulerian_obj_r)
          ! Allocate buffer
          allocate(buff_r(lo(1):hi(1), lo(2):hi(2), lo(3):hi(3)))

          ! Read field
          call hdf5%Read('/',trim(adjustl(my_field%name)),buff_r,this%block%lo,this%block%hi)
          my_field%cell(lo(1):hi(1), lo(2):hi(2), lo(3):hi(3))=buff_r

          ! Clean
          deallocate(buff_r)

        type is (eulerian_obj_i)
          ! Allocate buffer
          allocate(buff_i(lo(1):hi(1), lo(2):hi(2), lo(3):hi(3)))

          ! Read field
          call hdf5%Read('/',trim(adjustl(my_field%name)),buff_i,lo,hi)
          my_field%cell(lo(1):hi(1), lo(2):hi(2), lo(3):hi(3))=buff_i

          ! Clean
          deallocate(buff_i)
        end select

        ! Update ghost cells
        call this%field(ind)%p%UpdateGhostCells
      end associate

      return
    end subroutine eulerian_set_ReadSingleHDF5