tem_commbuf_int_gatherindexed Subroutine

private subroutine tem_commbuf_int_gatherindexed(me, pos, nvals)

gather the indexed mpi datatype, which describes how the data in the state vector relates to the entries in the buffer. in contrast to the simple indexed type above, we try to minimize the number of blocks here, and gather contiguous blocks of memory together.

Arguments

Type IntentOptional Attributes Name
type(tem_intbuffer_type), intent(inout) :: me
integer, intent(in) :: pos(nvals)
integer, intent(in) :: nvals

Calls

proc~~tem_commbuf_int_gatherindexed~~CallsGraph proc~tem_commbuf_int_gatherindexed tem_commbuf_int_gatherindexed interface~append~9 append proc~tem_commbuf_int_gatherindexed->interface~append~9 mpi_type_commit mpi_type_commit proc~tem_commbuf_int_gatherindexed->mpi_type_commit proc~check_mpi_error check_mpi_error proc~tem_commbuf_int_gatherindexed->proc~check_mpi_error interface~init~22 init proc~tem_commbuf_int_gatherindexed->interface~init~22 mpi_type_indexed mpi_type_indexed proc~tem_commbuf_int_gatherindexed->mpi_type_indexed interface~destroy~9 destroy proc~tem_commbuf_int_gatherindexed->interface~destroy~9 proc~append_da_label append_da_label interface~append~9->proc~append_da_label proc~append_da_veclabel append_da_veclabel interface~append~9->proc~append_da_veclabel mpi_error_string mpi_error_string proc~check_mpi_error->mpi_error_string proc~tem_abort tem_abort proc~check_mpi_error->proc~tem_abort proc~init_ga2d_real init_ga2d_real interface~init~22->proc~init_ga2d_real proc~destroy_da_label destroy_da_label interface~destroy~9->proc~destroy_da_label interface~sortedposofval~5 sortedposofval proc~append_da_label->interface~sortedposofval~5 interface~expand~9 expand proc~append_da_label->interface~expand~9 mpi_abort mpi_abort proc~tem_abort->mpi_abort proc~append_da_veclabel->interface~expand~9 proc~sortposofval_label sortposofval_label interface~sortedposofval~5->proc~sortposofval_label proc~expand_da_label expand_da_label interface~expand~9->proc~expand_da_label

Contents


Source Code

  subroutine tem_commbuf_int_gatherindexed( me, pos, nvals )
    ! -------------------------------------------------------------------- !
    type(tem_intbuffer_type), intent(inout) :: me
    integer, intent(in) :: nvals
    integer, intent(in) :: pos(nvals)
    ! -------------------------------------------------------------------- !
    type(grw_intarray_type) :: blocklength
    type(grw_intarray_type) :: displ
    integer :: ival, counter
    integer :: ierror
    ! -------------------------------------------------------------------- !

    me%nvals = nvals

    ! initialize growing arrays, a kb should be fine to start with...
    call init(blocklength, 256)
    call init(displ, 256)

    if (nvals > 0) then

      ! start with the displacement of the first entry in the list
      call append(displ, pos(1)-1)
      counter = 1

      do ival=2,nvals
        if (pos(ival) == pos(ival-1)+1) then
          ! contiguous memory location following the previous one, increase the
          ! the blocklength.
          counter = counter + 1
        else
          ! new block encountered, record the block found so far
          call append(blocklength, counter)

          ! start new block
          call append(displ, pos(ival)-1)
          counter = 1
        end if
      end do

      ! finish the last block, by recording its found length:
      call append(blocklength, counter)

    end if

    ! call mpi_type_indexed(count, array_of_blocklengths, &
    !   &                   array_of_displacements, oldtype, newtype, ierror)
    call mpi_type_indexed( displ%nvals, blocklength%val, displ%val, &
      &                    mpi_integer, me%memindexed, ierror     )
    call check_mpi_error(ierror,'type indexed in tem_commbuf_int_gatherindexed')
    call mpi_type_commit(me%memindexed, ierror)
    call check_mpi_error(ierror,'commit memindexed in tem_commbuf_int_gatherindexed')

    call destroy(displ)
    call destroy(blocklength)

  end subroutine tem_commbuf_int_gatherindexed