Routine to build communication buffer using elemRanks. This routine can be used only if all elements need to be communicated but they need process-wise seperation. Uses nScalars to get position in the value array to communicate. For send buffer: elemRanks contains target ranks to send data to For recv buffer: elemRanks contains source ranks to recv data from
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(tem_communication_type), | intent(out) | :: | commBuffer |
send or recv communication buffer to be created |
||
integer, | intent(in) | :: | nScalars |
Number of scalars per element |
||
integer, | intent(in) | :: | nElems |
Total number of elements or points to communicate |
||
integer, | intent(in) | :: | elemRanks(nElems) |
Target or source rank for each element or point |
subroutine tem_comm_createBuffer( commBuffer, nScalars, nElems, elemRanks )
! -------------------------------------------------------------------------!
!> send or recv communication buffer to be created
type(tem_communication_type), intent(out) :: commBuffer
!> Number of scalars per element
integer, intent(in) :: nScalars
!> Total number of elements or points to communicate
integer, intent(in) :: nElems
!> Target or source rank for each element or point
integer, intent(in) :: elemRanks(nElems)
! -------------------------------------------------------------------------!
integer :: iElem, iProc, iVar, counter, pntPos
type(dyn_intArray_type) :: partnerProc
integer, allocatable :: pos(:)
! -------------------------------------------------------------------------!
! Create dynamic array of rank ids to communicate to initialize commBuffer
do iElem = 1, nElems
call append( me = partnerProc, &
& val = elemRanks(iElem) )
end do
! Initialize commBuffer
call tem_comm_init( me = commBuffer, &
& nProcs = partnerProc%nVals )
! Store rank id to communicate
do iProc = 1, partnerProc%nVals
commBuffer%proc(iProc) = partnerProc%val(iProc)
end do
! Create map from commBuffer to elem array
do iElem = 1, nElems
iProc = PositionOfVal( me = partnerProc, &
& val = elemRanks(iElem) )
call append( me = commBuffer%elemPos(iProc), &
& val = iElem )
end do
call destroy(partnerProc)
commBuffer%nElemsProc(:) = commBuffer%elemPos(:)%nVals
! Get position in state array to send or recv data
allocate( pos( nScalars * maxVal( commBuffer%nElemsProc(:) ) ) )
! Assign comm buffer positions
do iProc = 1, commBuffer%nProcs
counter = 0
! loop of nElems per proc
do iElem = 1, commBuffer%nElemsProc(iProc)
! position of this proc point in the point array
pntPos = commBuffer%elemPos(iProc)%val(iElem)
do iVar = 1, nScalars
counter = counter + 1
! position in evalVal array which has size: nElems*nScalars
pos(counter) = (pntPos-1)*nScalars + iVar
end do !iVar
end do !iElem
! copy position array to me%pos, allocate me%val array
call tem_commbuf_real_fillPos( me = commBuffer%buf_real(iProc), &
& pos = pos, &
& nVals = counter )
end do !iProc
deallocate(pos)
end subroutine tem_comm_createBuffer