tem_extend_commFromFinerPrp Subroutine

private subroutine tem_extend_commFromFinerPrp(levelDesc, direction, faces)

the neighboring halo element is from finer on the remote partition.

Arguments

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

Level descriptor of the level of the mesh you want to collect the faces for.

integer, intent(in) :: direction

The face direction.

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

Description of the faces on the current level.


Calls

proc~~tem_extend_commfromfinerprp~~CallsGraph proc~tem_extend_commfromfinerprp tem_extend_commFromFinerPrp proc~tem_get_faceneigh tem_get_faceNeigh proc~tem_extend_commfromfinerprp->proc~tem_get_faceneigh proc~tem_appendface_prp tem_appendFace_prp proc~tem_extend_commfromfinerprp->proc~tem_appendface_prp proc~tem_treeidintotal tem_treeIDinTotal proc~tem_get_faceneigh->proc~tem_treeidintotal proc~tem_idofcoord tem_IdOfCoord proc~tem_get_faceneigh->proc~tem_idofcoord proc~tem_coordofid tem_CoordOfId proc~tem_get_faceneigh->proc~tem_coordofid proc~tem_abort tem_abort proc~tem_get_faceneigh->proc~tem_abort interface~positionofval~5 positionofval proc~tem_appendface_prp->interface~positionofval~5 proc~tem_melt_faceprp tem_melt_facePrp proc~tem_appendface_prp->proc~tem_melt_faceprp proc~tem_etypeofid tem_eTypeOfId proc~tem_treeidintotal->proc~tem_etypeofid tem_positioninsorted tem_positioninsorted proc~tem_treeidintotal->tem_positioninsorted proc~posofval_label posofval_label interface~positionofval~5->proc~posofval_label proc~tem_levelof tem_LevelOf proc~tem_coordofid->proc~tem_levelof mpi_abort mpi_abort proc~tem_abort->mpi_abort proc~tem_etypeofid->interface~positionofval~5 interface~sortedposofval~5 sortedposofval proc~posofval_label->interface~sortedposofval~5

Called by

proc~~tem_extend_commfromfinerprp~~CalledByGraph proc~tem_extend_commfromfinerprp tem_extend_commFromFinerPrp proc~tem_extend_remoteprp tem_extend_remotePrp proc~tem_extend_remoteprp->proc~tem_extend_commfromfinerprp proc~tem_build_face_info tem_build_face_info proc~tem_build_face_info->proc~tem_extend_remoteprp

Contents


Source Code

  subroutine tem_extend_commFromFinerPrp( levelDesc, direction, faces )
    ! --------------------------------------------------------------------------
    !> Level descriptor of the level of the mesh you want to collect the
    !! faces for.
    type(tem_levelDesc_type), intent(in) :: levelDesc
    !> The face direction.
    integer, intent(in) :: direction
    !> Description of the faces on the current level.
    type(tem_face_descriptor_type), intent(inout) :: faces
    ! --------------------------------------------------------------------------
    integer :: iProc, iElem
    integer :: neighIndex
    integer(kind=long_k) :: elemId, neighId
    integer :: elemPos, neighPos
    ! --------------------------------------------------------------------------

    ! We run over all the elements in the from finer buffer and attach
    ! the from finer property to the face the element is aligned to.
    processLoop: do iProc = 1, levelDesc%recvbufferFromFiner%nProcs

      elementLoop: do iElem = 1, levelDesc%recvbufferFromFiner%nElemsProc(iProc)

        ! The element position in the total list of the level descriptor.
        elemPos = levelDesc%recvbufferFromFiner%elemPos(iProc)%val(iElem)
        elemId = levelDesc%total(elemPos)

        ! Attach elem property for its right face (checks also if face exists).
        neighIndex = tem_right
        call tem_get_faceNeigh( levelDesc, elemPos, direction, neighIndex, &
          &                     neighId, neighPos                          )
        call tem_appendFace_prp( elemId, tem_fromFinerFace_prp, faces, &
          &                      tem_left                      )

        ! Attach elem property for its left face (checks also if face exists).
        neighIndex = tem_left
        call tem_get_faceNeigh( levelDesc, elemPos, direction, neighIndex, &
          &                     neighId, neighPos                          )
        call tem_appendFace_prp( neighId, tem_fromFinerFace_prp, faces, &
          &                      tem_right                      )

    end do elementLoop
  end do processLoop

  end subroutine tem_extend_commFromFinerPrp