Communicate with all existing process the number of requested halo elements
After this routine, each process knows how many processes there are to communicate with and how many elements have to be transferred
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(tem_levelDesc_type), | intent(inout) | :: | me(minlevel:maxLevel) |
the level descriptor to be filled |
||
type(tem_comm_env_type), | intent(in) | :: | proc |
Process description to use. |
||
integer, | intent(in) | :: | minLevel |
level range |
||
integer, | intent(in) | :: | maxLevel |
level range |
subroutine communicate_nElemsToTransfer( me, proc, minLevel, maxLevel )
! ---------------------------------------------------------------------------
!> level range
integer, intent(in) :: minLevel, maxLevel
!> Process description to use.
type(tem_comm_env_type), intent(in) :: proc
!> the level descriptor to be filled
type(tem_levelDesc_type), intent(inout) :: me(minlevel:maxLevel)
! ---------------------------------------------------------------------------
! the number of halos which will be sent to this mpi process for each level
! and for each process:
! first dimension is the level, second one is the mpi process
integer, allocatable :: nHalosHis(:,:)
integer, allocatable :: nHalosMine(:,:)
integer, allocatable :: sources(:)
integer, allocatable :: sourceHalos(:)
integer :: iProc, iLevel, ierr
integer :: nLevels, sourceProc, nProcs
real(kind=rk) :: tStart, tEnd
! ---------------------------------------------------------------------------
! ---------------------------------------------------------------------------
! Count number of halos on each level for each other process
! these values will be sent in the following
nLevels = maxLevel - minLevel + 1
tStart = mpi_wtime()
if (use_sparse_alltoall) then
write(logunit(4),*) 'Starting sparse(all-to-all) exchange of halo counts'
! Using levelwise sparse all-to-all here for now out of simplicity.
! If there is a need, we could first count the number of levels to
! exchange with each process, and then have only this in a single,
! sparse all-to-all.
! However, this is a little involved and probably requires us to
! employ a dynamic array for the processes...
! In most cases the levelwise sparse alltoall as done now is probably
! good enough.
do iLevel = minLevel, maxLevel
write(logunit(10),*) ' Level: ', iLevel
nProcs = me(iLevel)%haloList%PartnerProc%nVals
allocate( nHalosMine(nProcs,1) )
nHalosMine(:,1) = me(iLevel)%haloList%halos%val(:nProcs)%nVals
call tem_comm_init( me = me(iLevel)%recvbuffer, &
& nProcs = nProcs )
me(iLevel)%recvBuffer%nElemsProc = nHalosMine(:,1)
me(iLevel)%recvBuffer%proc &
& = me(iLevel)%halolist%partnerProc%val(:nProcs) - 1
call tem_sparse_alltoall_int( targets = me(iLevel)%recvBuffer &
& %Proc, &
& send_buffer = nHalosMine(:,1), &
& sources = sources, &
& recv_buffer = sourceHalos, &
& tag = iLevel, &
& comm = proc%comm )
call tem_comm_init( me = me( iLevel )%sendbuffer, &
& nProcs = size(sources) )
me(iLevel)%sendBuffer%nElemsProc = sourceHalos
me(iLevel)%sendBuffer%proc = sources
deallocate(sourceHalos)
deallocate(sources)
deallocate(nHalosMine)
end do
else
write(logunit(4),*) 'Starting all-to-all exchange of halo counts'
! allocate the send buffer.
! Contains the number of elements for each level, which are requested
! from the specific process
allocate( nHalosMine( minLevel:maxLevel,proc%comm_size ))
! allocate the recv buffer.
! Will contain the number of elements for each level, which other
! processes request from me. This one is filled in the mpi_alltoall call
allocate( nHalosHis( minLevel:maxLevel, proc%comm_size ))
! Reset the number of requested and requesting halos
nHalosMine(:,:) = 0
nHalosHis(:,:) = 0
do iLevel = minLevel, maxLevel
nProcs = me(iLevel)%haloList%PartnerProc%nVals
do iProc = 1, nProcs
sourceProc = me(iLevel)%haloList%partnerProc%val(iProc)
nHalosMine(iLevel, sourceProc) = &
& me(iLevel)%haloList%halos%val(iProc)%nVals
end do ! iProc
call tem_comm_init( me(iLevel)%recvbuffer, nProcs )
! count valid procs and nElemsProc
call tem_comm_count( me(iLevel)%recvBuffer, proc%comm_size, &
& nHalosMine(iLevel,:) )
end do ! iLevel
! number of halos which are sent to a process
! might be different from number of halos
! received from the same process for that we use mpi_alltoall
! Exchange the number of requested treeIDs
call mpi_alltoall( nHalosMine, nLevels, mpi_integer, &
& nHalosHis, nLevels, mpi_integer, &
& proc%comm, ierr )
do iLevel = minLevel, maxLevel
! now check which processs really ask for information and allocate
! buffer
nProcs = count( nHalosHis( iLevel, : ) > 0 )
call tem_comm_init( me( iLevel )%sendbuffer, nProcs )
! count valid procs and nElemsProc
call tem_comm_count( me(iLevel)%sendBuffer, proc%comm_size, &
& nHalosHis(iLevel,:))
end do
deallocate( nHalosMine )
deallocate( nHalosHis )
end if
tEnd = mpi_wtime()
write(logunit(4),"(A)") 'Finished all-to-all exchanging number of elements'
write(logunit(4),"(A,E12.6)") 'All to all cost: ', tEnd-tStart
do iLevel = minLevel, maxLevel
! Copy all headers (information about target and source procs)
! to bufferFromCoarser and bufferFromFiner
me( iLevel )%sendbufferFromCoarser = me( iLevel )%sendbuffer
me( iLevel )%sendbufferFromFiner = me( iLevel )%sendbuffer
me( iLevel )%recvbufferFromCoarser = me( iLevel )%recvbuffer
me( iLevel )%recvbufferFromFiner = me( iLevel )%recvbuffer
end do ! iLevel
write(logUnit(5),*) ' Finished counting valid procs and nElemsProcs.'
end subroutine communicate_nElemsToTransfer