tem_get_faceNeigh Subroutine

private subroutine tem_get_faceNeigh(levelDesc, elemPos, dir, leftOrRight, neighId, neighPos)

Function to get the face neighbor of a certain element in the level descriptor. Even non-existing face neighbors can be handled by this routine.

Arguments

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

The level descriptor the element is located in.

integer, intent(in) :: elemPos

The element position in the level descriptor.

integer, intent(in) :: dir

The spatial direction of the face.

  1. -> x direction
  2. -> y direction
  3. -> z direction
integer, intent(in) :: leftOrRight

Find left or right face neighbor of the element. Use \ref tem_face_module::tem_left or \ref tem_face_module::tem_right.

integer(kind=long_k), intent(out) :: neighId

The treeid of the face neighbor element.

integer, intent(out) :: neighPos

Position of the neighbor in the total list, or 0 if element does not exist.


Calls

proc~~tem_get_faceneigh~~CallsGraph proc~tem_get_faceneigh tem_get_faceNeigh proc~tem_coordofid tem_CoordOfId proc~tem_get_faceneigh->proc~tem_coordofid proc~tem_treeidintotal tem_treeIDinTotal proc~tem_get_faceneigh->proc~tem_treeidintotal proc~tem_abort tem_abort proc~tem_get_faceneigh->proc~tem_abort proc~tem_idofcoord tem_IdOfCoord proc~tem_get_faceneigh->proc~tem_idofcoord proc~tem_levelof tem_LevelOf proc~tem_coordofid->proc~tem_levelof proc~tem_etypeofid tem_eTypeOfId proc~tem_treeidintotal->proc~tem_etypeofid tem_positioninsorted tem_positioninsorted proc~tem_treeidintotal->tem_positioninsorted mpi_abort mpi_abort proc~tem_abort->mpi_abort interface~positionofval~5 positionofval proc~tem_etypeofid->interface~positionofval~5 proc~posofval_label posofval_label interface~positionofval~5->proc~posofval_label

Called by

proc~~tem_get_faceneigh~~CalledByGraph proc~tem_get_faceneigh tem_get_faceNeigh proc~tem_get_faces tem_get_faces proc~tem_get_faces->proc~tem_get_faceneigh proc~tem_extend_commfromfinerprp tem_extend_commFromFinerPrp proc~tem_extend_commfromfinerprp->proc~tem_get_faceneigh proc~tem_extend_commfromcoarserprp tem_extend_commFromCoarserPrp proc~tem_extend_commfromcoarserprp->proc~tem_get_faceneigh proc~tem_collect_faces tem_collect_faces proc~tem_collect_faces->proc~tem_get_faces proc~tem_extend_remoteprp tem_extend_remotePrp proc~tem_extend_remoteprp->proc~tem_extend_commfromfinerprp proc~tem_extend_remoteprp->proc~tem_extend_commfromcoarserprp proc~tem_build_face_info tem_build_face_info proc~tem_build_face_info->proc~tem_collect_faces proc~tem_build_face_info->proc~tem_extend_remoteprp

Contents

Source Code


Source Code

  subroutine tem_get_faceNeigh( levelDesc, elemPos, dir, leftOrRight, neighId, &
    &                           neighPos )
    ! --------------------------------------------------------------------------
    !> The level descriptor the element is located in.
    type(tem_levelDesc_type), intent(in) :: levelDesc
    !> The element position in the level descriptor.
    integer, intent(in) :: elemPos
    !> The spatial direction of the face.
    !!
    !! 1. -> x direction
    !! 2. -> y direction
    !! 3. -> z direction
    integer, intent(in) :: dir
    !> Find left or right face neighbor of the element.
    !! Use
    !! \ref tem_face_module::tem_left
    !! or
    !! \ref tem_face_module::tem_right.
    integer, intent(in) :: leftOrRight
    !> The treeid of the face neighbor element.
    integer(kind=long_k), intent(out) :: neighId
    !> Position of the neighbor in the total list, or 0 if element does not
    !! exist.
    integer, intent(out) :: neighPos
    ! --------------------------------------------------------------------------
    integer(kind=long_k) :: elemId
    integer :: elemCoord(4)
    integer :: dirOffset
    ! --------------------------------------------------------------------------

    ! We check the neighbor position of the element.
    neighPos = levelDesc%neigh(1)%nghElems( leftOrRight, elemPos )
    elemId = levelDesc%total(elemPos)

    ! Can we find the neighbor in the level descriptor. If not, we have
    ! to build the tree id our own.
    if (neighPos>0 .and. neighPos<=levelDesc%nElems) then
      neighId = levelDesc%total(neighPos)
    else
      ! So, the neighbor element does not exist in level descriptor.
      ! Let's create the tree id (of the neighbor on this level) by our own.
      elemCoord = tem_coordOfId(elemId)
      ! Initialize dirOffset just to silence compiler warning
      dirOffset = 0
      select case(leftOrRight)
      case(tem_left)
        dirOffset = -1
      case(tem_right)
        dirOffset = 1
      case default
        write(*,*) 'ERROR in tem_get_faceNeigh: unknown face neighbor '//      &
          &        'direction stopping ... '
        write(*,*) 'Given direction was: ', dir
        call tem_abort()
      end select
      elemCoord(dir) = elemCoord(dir) + dirOffset
      neighId = tem_idOfCoord(elemCoord)

      ! Check if it is boundary element or not
      if (neighPos<0) then
        ! In case of a boundary element, we keep the negative boundary
        ! unmodified.
      else
        ! If no, we try to locate the neighbor in the total list of the level
        ! descriptor
        neighPos = tem_treeIDinTotal( neighId, levelDesc)
      end if
    endif

  end subroutine tem_get_faceNeigh