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.
Type | Intent | Optional | 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.
|
||
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. |
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