tem_build_faceSendBuffers Subroutine

private subroutine tem_build_faceSendBuffers(levelDesc, faces, buf)

current rank has to send information for (before the other ranks can make the compute step).

Arguments

Type IntentOptional Attributes Name
type(tem_levelDesc_type), intent(in) :: levelDesc

Dimension-by-dimension level descriptor for the current level and direction.

type(tem_face_descriptor_type), intent(inout) :: faces

The communication pattern you want use for the buffer. The created face descriptor.

type(tem_communication_type), intent(out) :: buf(2)

The created send buffer. Size is two, due to left and right limes of the face values. To access it use tem_left and tem_right.


Calls

proc~~tem_build_facesendbuffers~~CallsGraph proc~tem_build_facesendbuffers tem_build_faceSendBuffers interface~positionofval~5 positionofval proc~tem_build_facesendbuffers->interface~positionofval~5 interface~append~23 append proc~tem_build_facesendbuffers->interface~append~23 proc~tem_abort tem_abort proc~tem_build_facesendbuffers->proc~tem_abort proc~tem_issendface tem_isSendFace proc~tem_build_facesendbuffers->proc~tem_issendface interface~init~22 init proc~tem_build_facesendbuffers->interface~init~22 proc~posofval_label posofval_label interface~positionofval~5->proc~posofval_label proc~append_arrayga2d_real append_arrayga2d_real interface~append~23->proc~append_arrayga2d_real proc~append_singlega2d_real append_singlega2d_real interface~append~23->proc~append_singlega2d_real mpi_abort mpi_abort proc~tem_abort->mpi_abort proc~init_ga2d_real init_ga2d_real interface~init~22->proc~init_ga2d_real interface~sortedposofval~5 sortedposofval proc~posofval_label->interface~sortedposofval~5 interface~expand~22 expand proc~append_arrayga2d_real->interface~expand~22 proc~append_singlega2d_real->interface~expand~22 proc~sortposofval_label sortposofval_label interface~sortedposofval~5->proc~sortposofval_label proc~expand_ga2d_real expand_ga2d_real interface~expand~22->proc~expand_ga2d_real

Called by

proc~~tem_build_facesendbuffers~~CalledByGraph proc~tem_build_facesendbuffers tem_build_faceSendBuffers proc~tem_build_facebuffers tem_build_faceBuffers proc~tem_build_facebuffers->proc~tem_build_facesendbuffers proc~tem_build_face_info tem_build_face_info proc~tem_build_face_info->proc~tem_build_facebuffers

Contents


Source Code

  subroutine tem_build_faceSendBuffers( levelDesc, faces, buf )
    ! --------------------------------------------------------------------------
    !> Dimension-by-dimension level descriptor for the current level and
    !! direction.
    type(tem_levelDesc_type), intent(in) :: levelDesc
    !> The communication pattern you want use for the buffer.
    ! type(tem_commpattern_type), intent(in) :: commPattern
    !> The created face descriptor.
    type(tem_face_descriptor_type),intent(inout)  :: faces
    !> The created send buffer. Size is two, due to left and right limes of
    !! the face values. To access it use tem_left and tem_right.
    type(tem_communication_type), intent(out) :: buf(2)
    ! --------------------------------------------------------------------------
    integer :: iFace, iProc, iSide, elemPos, faceSide, tIdPos, rank, rankPos,  &
      &        faceIndex, elemAddPos
    logical :: wasAdded
    integer(kind=long_k) :: elemId, trgElemId
    ! The element treeids I will send data for. We have to differentiate the
    ! left and right face of a face here.
    type(dyn_longArray_type) :: elemIds(2)
    ! The elements positions I will send data for. We have to differentiate the
    ! left and right face of a face here.
    type(grw_intArray_type) :: elemPositions(2)
    ! The source ranks.
    type(dyn_intArray_type) :: trgRank(2)
    ! The position of the source ranks.
    type(grw_intArray_type) :: trgRankPos(2)
    ! --------------------------------------------------------------------------

    ! Init our growing arrays.
    call init( me = elemIds(1) , length = 16 )
    call init( me = elemIds(2) , length = 16 )
    call init( me = elemPositions(1) , length = 16)
    call init( me = elemPositions(2) , length = 16)
    call init( me = trgRankPos(1) , length = 16)
    call init( me = trgRankPos(2) , length = 16)
    call init( me = trgRank(1) , length = 16 )
    call init( me = trgRank(2) , length = 16 )

    ! Iterate over all the faces and count the communicate faces.
    do iFace = 1, faces%faceList%faceId%nVals

      ! Check whether we will send info of one of the adjacent elements
      ! of this face.
      faceSide = tem_isSendFace(iFace, faces)
      if ( faceSide /= 0 ) then

        ! Check for which adjacent element (left or right) the current rank
        ! will send information about and count it.
        if (faceSide == tem_left) then
          ! The position and id of the elment I will send data for
          elemPos = faces%faceList%leftElemPos%val(iFace)
          elemId = faces%faceList%faceId%val(iFace)
          ! I will send information about my left element, so I will send
          ! the right face info of this element.
          faceSide = tem_right
          ! Get the element id in the opposite direction of the face. This
          ! will be used to determine the target rank.
          trgElemId = faces%faceList%rightElemId%val(iFace)
        else
          ! The position and id of the elment I will send data for
          elemPos = faces%faceList%rightElemPos%val(iFace)
          elemId = faces%faceList%rightElemId%val(iFace)
          ! I will send info about the right element, so I send info about
          ! it's left face.
          faceSide = tem_left
          ! Get the element id in the opposite direction of the face. This
          ! will be used to determine the target rank.
          trgElemId = faces%faceList%faceId%val(iFace)
        end if

        ! Get the source rank of the target element id. This will give us
        ! the target rank of the face.
        ! Therefore we have a look at the level
        ! descriptor of the current spatial direction and level.
        tIdPos = PositionOfVal( me = levelDesc%elem%tId, val = trgElemId )
        if (tIdPos <= 0) then
          write(*,*) 'ERROR in tem_build_faceSendBuffers: not able to '//  &
            &        'identify target rank for a certain face, stopping ...'
          call tem_abort()
        end if
        rank = levelDesc%elem%sourceProc%val(tIdPos)

        ! Store element position and soruce proc in a dynamic array.
        call append( me       = elemIds(faceSide), &
          &          val      = elemId,            &
          &          pos      = elemAddPos,        &
          &          wasAdded = wasAdded           )
        if(wasAdded) then
          call append( me = elemPositions(faceSide), val = elemPos )
          call append( me = trgRank(faceSide), val = rank, pos = rankPos )
          call append( me = trgRankPos(faceSide), val = rankPos )
        end if

      end if
    end do

    ! Now, we build the buffer for the faces (left and right values).
    do iSide = 1, 2
      ! Init the datastructures in the buffer
      buf(iSide)%nProcs = trgRank(iSide)%nVals
      allocate(buf(iSide)%proc( trgRank(iSide)%nVals ))
      allocate(buf(iSide)%nElemsProc( trgRank(iSide)%nVals ))
      allocate(buf(iSide)%elemPos( trgRank(iSide)%nVals ))
      do iProc = 1, buf(iSide)%nProcs
        call init( me = buf(iSide)%elemPos(iProc), length = 16 )
      end do

      ! Set the positions of the elements I will receive.
      do iFace = 1, elemPositions(iSide)%nVals
        faceIndex = elemIds(iSide)%sorted(iFace)
        elemPos = elemPositions(iSide)%val(faceIndex)
        rankPos = trgRankPos(iSide)%val(faceIndex)
        call append( me = buf(iSide)%elemPos(rankPos), val = elemPos )
      end do

      ! Count the number of elements we recv from each proc.
      do iProc = 1, buf(iSide)%nProcs
        buf(iSide)%proc(iProc) = trgRank(iSide)%val(iProc) - 1
        buf(iSide)%nElemsProc(iProc) = buf(iSide)%elemPos(iProc)%nVals
      end do
    end do

  end subroutine tem_build_faceSendBuffers