request_remoteHalos Subroutine

private subroutine request_remoteHalos(levelDesc, proc, tree, iLevel, stencil, pathFirst, pathLast)

Inverse Communication: Communicate, which elements each process needs from me.

In this routine, we send the treeIDs of the halo elements to the processes, where they are located. Later on, we fill these halos locally with information from these processes (sourceProcs). In this routine however, we now SEND information to these sourceProcs, so do not get confused here.

Arguments

Type IntentOptional Attributes Name
type(tem_levelDesc_type), intent(inout) :: levelDesc(tree%global%minlevel:)

the level descriptor to be filled

type(tem_comm_env_type), intent(in) :: proc

Process description to use.

type(treelmesh_type), intent(in) :: tree

the global tree

integer, intent(in) :: iLevel

current level

type(tem_stencilHeader_type), intent(in) :: stencil

stencil definition

type(tem_path_type), intent(in) :: pathFirst(:)

first and last treeID path in every process

type(tem_path_type), intent(in) :: pathLast(:)

first and last treeID path in every process


Calls

proc~~request_remotehalos~~CallsGraph proc~request_remotehalos request_remoteHalos interface~append~11 append proc~request_remotehalos->interface~append~11 interface~destroy~25 destroy proc~request_remotehalos->interface~destroy~25 interface~init~24 init proc~request_remotehalos->interface~init~24 interface~positionofval~5 positionofval proc~request_remotehalos->interface~positionofval~5 mpi_irecv mpi_irecv proc~request_remotehalos->mpi_irecv mpi_isend mpi_isend proc~request_remotehalos->mpi_isend mpi_waitall mpi_waitall proc~request_remotehalos->mpi_waitall proc~create_allparentneighbors create_allParentNeighbors proc~request_remotehalos->proc~create_allparentneighbors proc~identify_halo identify_halo proc~request_remotehalos->proc~identify_halo proc~identify_stencilneigh identify_stencilNeigh proc~request_remotehalos->proc~identify_stencilneigh proc~tem_abort tem_abort proc~request_remotehalos->proc~tem_abort proc~tem_horizontalspacer tem_horizontalSpacer proc~request_remotehalos->proc~tem_horizontalspacer proc~append_ga_dynlong append_ga_dynlong interface~append~11->proc~append_ga_dynlong proc~append_ga_dynlong_vec append_ga_dynlong_vec interface~append~11->proc~append_ga_dynlong_vec proc~destroy_ga2d_real destroy_ga2d_real interface~destroy~25->proc~destroy_ga2d_real proc~init_ga2d_real init_ga2d_real interface~init~24->proc~init_ga2d_real proc~posofval_label posofval_label interface~positionofval~5->proc~posofval_label proc~create_allparentneighbors->interface~append~11 proc~create_allparentneighbors->proc~identify_stencilneigh interface~tem_parentof tem_ParentOf proc~create_allparentneighbors->interface~tem_parentof proc~identify_elements identify_elements proc~create_allparentneighbors->proc~identify_elements proc~identify_halo->interface~positionofval~5 proc~identify_local_element identify_local_element proc~identify_halo->proc~identify_local_element proc~tem_levelof tem_LevelOf proc~identify_halo->proc~tem_levelof proc~identify_stencilneigh->proc~identify_elements mpi_abort mpi_abort proc~tem_abort->mpi_abort

Called by

proc~~request_remotehalos~~CalledByGraph proc~request_remotehalos request_remoteHalos proc~communicate_elements communicate_elements proc~communicate_elements->proc~request_remotehalos proc~tem_find_allelements tem_find_allElements proc~tem_find_allelements->proc~communicate_elements proc~tem_create_leveldesc tem_create_levelDesc proc~tem_create_leveldesc->proc~tem_find_allelements proc~tem_dimbydim_construction tem_dimByDim_construction proc~tem_dimbydim_construction->proc~tem_create_leveldesc

Source Code

  subroutine request_remoteHalos( levelDesc, proc, tree, iLevel, stencil,&
    &                             pathFirst, pathLast )
    ! ---------------------------------------------------------------------------
    !> the global tree
    type(treelmesh_type), intent(in) :: tree
    !> the level descriptor to be filled
    type(tem_levelDesc_type), intent(inout) :: levelDesc(tree%global%minlevel: )
    !> Process description to use.
    type(tem_comm_env_type), intent(in) :: proc
    !> current level
    integer, intent(in) :: iLevel
    !> stencil definition
    type(tem_stencilHeader_type), intent(in) :: stencil
    !> first and last treeID path in every process
    type(tem_path_type), intent(in) :: pathFirst(:), pathLast(:)
    ! ---------------------------------------------------------------------------
    integer :: iProc, iErr, iElem, elemPos, procPos
    integer :: haloLevel
    integer, allocatable :: rq_handle(:)
    integer, allocatable :: status(:,: )
    integer :: nCommunications, nesting
    type( grw_longArray_type ), allocatable ::  treeIDs_fromTarget(:)
    type( grw_intArray_type ),  allocatable :: nestings_fromTarget(:)
    type( grw_longArray_type ), allocatable ::  treeIDs_toSource(:)
    type( grw_intArray_type ),  allocatable :: nestings_toSource(:)
    logical :: updated
    integer, parameter :: message_flag_long = 24
    integer, parameter :: message_flag_int  = 25
    ! ---------------------------------------------------------------------------

    call tem_horizontalSpacer( fUnit = dbgUnit(1), before = 1 )
    write(dbgUnit(1),*) "Get into routine: request_remoteHalos"
    write(dbgUnit(1),*) 'Requesting remote halos on level: ', iLevel

    ! two communications: treeID and nesting
    nCommunications = 2 * (  levelDesc( iLevel )%sendbuffer%nProcs &
      &                    + levelDesc( iLevel )%recvbuffer%nProcs )
    allocate( rq_handle( nCommunications ) )
    allocate( status( mpi_status_size, nCommunications ) )
    rq_handle(:) = MPI_REQUEST_NULL

    ! Warning: Inverse Communication ! (send to source, recv from target)
    ! ---------------------------------------------------------------------
    ! I receive from target what elements are needed by them
    ! SendBuffer contains my elements required by remote targets
    allocate(  treeIDs_fromTarget( levelDesc( iLevel )%sendbuffer%nProcs ))
    allocate( nestings_fromTarget( levelDesc( iLevel )%sendbuffer%nProcs ))

    do iProc = 1, levelDesc( iLevel )%sendbuffer%nProcs

      ! Allocate the buffers
      call init( me     = treeIDs_fromTarget( iProc ),                      &
        &        length = levelDesc( iLevel )%sendbuffer%nElemsProc( iProc ))
      call init( me     = nestings_fromTarget( iProc ),                     &
        &        length = levelDesc( iLevel )%sendbuffer%nElemsProc( iProc ))

      ! Receive the element tree IDs
      call mpi_irecv( treeIDs_fromTarget( iProc )%val,                         &
        &             treeIDs_fromTarget( iProc )%ContainerSize,               &
        &             mpi_integer8,                                            &
        &             levelDesc( iLevel )%sendbuffer%proc(iProc),              &
        &             message_flag_long,                                       &
        &             proc%comm,                                               &
        &             rq_handle( iProc),                                       &
        &             iErr )
      ! Receive the element nestings
      call mpi_irecv( nestings_fromTarget( iProc )%val,                        &
        &             nestings_fromTarget( iProc )%ContainerSize,              &
        &             mpi_integer,                                             &
        &             levelDesc( iLevel )%sendbuffer%proc(iProc),              &
        &             message_flag_int,                                        &
        &             proc%comm,                                               &
        &             rq_handle( iProc+nCommunications/2 ),                    &
        &             iErr  )
      ! Update the number of elements inside the growing array of the recv
      ! buffer
       treeIDs_fromTarget( iProc )%nVals = treeIDs_fromTarget( iProc )%ContainerSize
      nestings_fromTarget( iProc )%nVals = treeIDs_fromTarget( iProc )%ContainerSize
    end do ! iProc
    ! ---------------------------------------------------------------------

    ! ---------------------------------------------------------------------
    ! I send to source what elements I need from them
    allocate(  treeIDs_toSource( levelDesc( iLevel )%recvbuffer%nProcs ) )
    allocate( nestings_toSource( levelDesc( iLevel )%recvbuffer%nProcs ) )

    do iProc = 1, levelDesc( iLevel )%recvbuffer%nProcs

      ! Get the position of the process in the dynamic halos list (might be
      ! unordered)
      procPos = PositionOfVal( me  = levelDesc(iLevel)%haloList%PartnerProc, &
        &                      val = levelDesc(iLevel)%recvbuffer            &
        &                                             %proc(iProc) + 1       )

      call init(  me     = treeIDs_toSource( iProc ),                        &
        &         length = levelDesc( iLevel )%recvbuffer%nElemsProc( iProc ))
      call init(  me     = nestings_toSource( iProc ),                       &
        &         length = levelDesc( iLevel )%recvbuffer%nElemsProc( iProc ))

      ! Collect the halo treeIDs into send buffers for all the processes to
      ! request from
      do iElem = 1, levelDesc( iLevel )%recvbuffer%nElemsProc( iProc )
        elemPos = levelDesc( iLevel )%haloList%halos%val(procPos)%val(iElem)
        call append( me  = treeIDs_toSource( iProc ),                  &
          &          val = levelDesc( iLevel )%elem%tID%val( elemPos ) )
        call append( me  = nestings_toSource( iProc ),                         &
          &          val = levelDesc( iLevel )%elem%haloNesting%val( elemPos ) )
      enddo

      ! Send treeIDs
      call mpi_isend(                                                           &
        &        treeIDs_toSource(iProc)%val,                                   &
        &        treeIDs_toSource(iProc)%nVals,                                 &
        &        mpi_integer8,                                                  &
        &        levelDesc( iLevel )%recvbuffer%proc( iProc ),                  &
        &        message_flag_long,                                             &
        &        proc%comm,                                                     &
        &        rq_handle( iProc + levelDesc( iLevel )%sendbuffer%nProcs),     &
        &        ierr  )
      ! Send nesting
      call mpi_isend(                                                           &
        &        nestings_toSource(iProc)%val,                                  &
        &        nestings_toSource(iProc)%nVals,                                &
        &        mpi_integer,                                                   &
        &        levelDesc( iLevel )%recvbuffer%proc( iProc ),                  &
        &        message_flag_int,                                              &
        &        proc%comm,                                                     &
        &        rq_handle( iProc + levelDesc( iLevel )%sendbuffer%nProcs       &
        &                         + nCommunications/2),                         &
        &        ierr  )

    end do ! iProc

    call mpi_waitall( nCommunications, rq_handle, status, ierr)
    write(logUnit(5),*) '  Received requested halo elements successfully.'
    deallocate(  treeIDs_toSource )
    deallocate( nestings_toSource )
    deallocate(         rq_handle )
    deallocate(         status    )
    ! Requested halo elements were received.
    ! ---------------------------------------------------------------------------

    ! ---------------------------------------------------------------------------
    ! Now identify the requested halos.
    if( allocated( levelDesc( iLevel )%sendbufferFromCoarser%elemPos))         &
      &         deallocate( levelDesc( iLevel )%sendbufferFromCoarser%elemPos )
    allocate( levelDesc( iLevel )%sendbufferFromCoarser%elemPos(               &
      &                                levelDesc( iLevel )%sendbuffer%nProcs ))
    if( allocated( levelDesc( iLevel )%sendbufferFromFiner%elemPos ))          &
      &           deallocate( levelDesc( iLevel )%sendbufferFromFiner%elemPos )
    allocate( levelDesc( iLevel )%sendbufferFromFiner%elemPos(                 &
      &                                levelDesc( iLevel )%sendbuffer%nProcs ))

    ! Add elements of received buffers to elem
    do iProc = 1, levelDesc( iLevel )%sendbuffer%nProcs
      ! Allocate the buffer for the element position indices
      call init( me = levelDesc( iLevel )%sendbuffer%elemPos(iProc) )
      call init( me = levelDesc( iLevel )%sendbufferFromCoarser%elemPos(iProc) )
      call init( me = levelDesc( iLevel )%sendbufferFromFiner%elemPos(iProc) )

      do iElem = 1, treeIDs_fromTarget( iProc )%nVals
        nesting = nestings_fromTarget( iProc )%val( iElem )
        ! identify requested halo treeID in local process
        call identify_halo( haloTreeID = treeIDs_fromTarget( iProc )%val(iElem), &
          &                 elemPos    = elemPos,                          &
          &                 halolevel  = haloLevel,                        &
          &                 levelDesc  = levelDesc,                        &
          &                 nesting    = nesting,                          &
          &                 updated    = updated,                          &
          &                 tree       = tree,                             &
          &                 minLevel   = tree%global%minLevel,             &
          &                 stencil    = stencil )
        if( elemPos > 0 ) then
          ! if requested halo is ghostFromCoarser then find stencil neighbors of
          ! this halo element
          if ( (nesting < nestingLimit)                           &
            &  .and. (levelDesc( iLevel )%elem%eType%val(elemPos) &
            &         == eT_ghostFromCoarser)                     ) then
            ! identify all the compute neighbors of the current element
            call identify_stencilNeigh( iElem          = elemPos,    &
              &                         iLevel         = iLevel,     &
              &                         tree           = tree,       &
              &                         iStencil       = 1,          &
              &                         pathFirst      = pathFirst,  &
              &                         pathLast       = pathLast,   &
              &                         levelDesc      = levelDesc,  &
              &                         proc           = proc,       &
              &                         stencil        = stencil,    &
              &                         nesting        = nesting + 1 )
          end if

          ! if requested halo element haloNesting < found halo element (elemPos)
          ! haloNesting
          if ( nestings_fromTarget(iProc)%val(iElem)                  &
            &  < levelDesc( haloLevel )%elem%haloNesting%val(elemPos) ) then
            levelDesc(haloLevel)%elem%needsUpdate%val(elemPos) = .true.
            levelDesc(haloLevel)                                             &
              &  %elem                                                       &
              &  %haloNesting                                                &
              &  %val(elemPos) = min( nestings_fromTarget(iProc)%val(iElem), &
              &                       levelDesc(haloLevel)                   &
              &                         %elem                                &
              &                         %haloNesting                         &
              &                         %val(elemPos)                        )
          end if

          ! only add, if the element was added locally
          select case( levelDesc(iLevel)%elem%eType%val(elemPos) )
          ! Depending on the type of the element, add to the
          ! regular buffer, bufferFromCoarser, bufferFromFiner
          case( eT_fluid )
            call append( me  = levelDesc( iLevel )%sendbuffer%elemPos( iProc ),&
              &          val = elemPos )
          case( eT_ghostFromCoarser )
            call append( me  = levelDesc( iLevel )%sendbufferFromCoarser       &
              &                                             %elemPos( iProc ), &
              &          val = elemPos )

            ! for ghostFromCoarser determine neighbors of coarser element
            if( levelDesc( haloLevel )%elem%haloNesting%val( elemPos )         &
              &                                        < nestingLimit ) then
              call create_allParentNeighbors(                                  &
                &      targetID       = levelDesc(iLevel)%elem%tID%val( elemPos ),&
                &      level          = iLevel,                                &
                &      tree           = tree,                                  &
                &      stencil        = stencil,                        &
                &      levelDesc      = levelDesc,                             &
                &      pathFirst      = pathFirst,                             &
                &      pathLast       = pathLast,                              &
                &      proc           = proc )
            end if

          case( eT_ghostFromFiner )
            call append( me  = levelDesc( iLevel )%sendbufferFromFiner         &
              &                %elemPos( iProc ),                              &
              &          val = elemPos )
          case( eT_distributedGhostFromFiner)
             write(logUnit(1),*)' Found distributed ghost From Finer in '//    &
               &                'request remote Halos'
             write(logUnit(1),*)' This case should not occur!'
             call tem_abort()
          end select
        end if ! elemPos > 0
      end do ! ielem recvbuffer

      levelDesc( iLevel )%sendbuffer%nElemsProc( iProc )                       &
        &   = levelDesc( iLevel )%sendbuffer%elemPos( iProc )%nVals
      levelDesc( iLevel )%sendbufferFromCoarser%nElemsProc( iProc )            &
        &   = levelDesc( iLevel )%sendbufferFromCoarser%elemPos( iProc )%nVals
      levelDesc( iLevel )%sendbufferFromFiner%nElemsProc( iProc )              &
        &   = levelDesc( iLevel )%sendbufferFromFiner%elemPos( iProc )%nVals

      ! destroy temp variables
      call destroy( me =  treeIDs_fromTarget( iProc ) )
      call destroy( me = nestings_fromTarget( iProc ) )
    end do ! iProc
    deallocate(  treeIDs_fromTarget )
    deallocate( nestings_fromTarget )
    ! Now each Process knows, which elements to send to others
    ! ---------------------------------------------------------------------
    write(logUnit(5),*) 'Finished requesting remote halos'

    write(dbgUnit(1),*) "Leave  routine: request_remoteHalos"
    call tem_horizontalSpacer( fUnit = dbgUnit(1), after = 1 )

  end subroutine request_remoteHalos