Update the neighor arrays depending on what is given in the element stencil
The array levelDesc( iLevel )%neigh( iStenci )%nghElems( 1:QQN, 1:nElems ) is being filled up here
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | iStencil |
Index of your neighbor list. |
||
type(tem_levelDesc_type), | intent(inout) | :: | levelDesc |
Level descriptor for each level of your mesh (starting from min level). |
||
type(tem_stencilHeader_type), | intent(in) | :: | computeStencil |
The stencil you build the horizontal dependencies for. |
||
integer, | intent(in) | :: | list(:) |
stencil elemLvl points to sorted original treeID list |
||
integer, | intent(in) | :: | nElems |
number of elements |
||
type(treelmesh_type), | intent(in) | :: | tree |
tree information |
subroutine tem_build_treeHorizontalDep( iStencil, levelDesc, computeStencil, &
& list, nElems, tree )
! ---------------------------------------------------------------------------
!> Level descriptor for each level of your mesh (starting from min level).
type(tem_levelDesc_type), intent(inout) :: levelDesc
!> tree information
type(treelmesh_type), intent(in) :: tree
!> Index of your neighbor list.
integer, intent(in) :: iStencil
!> The stencil you build the horizontal dependencies for.
type(tem_stencilHeader_type), intent(in) :: computeStencil
!> number of elements
integer, intent(in) :: nElems
!> stencil elemLvl points to sorted original treeID list
integer, intent(in) :: list(:)
! ---------------------------------------------------------------------------
integer :: iElem, iStencilElem
integer :: elemPos, stencilPos, levelPos, totalPos
integer(kind=long_k) :: tID
logical :: missingNeigh
real(kind=rk) :: bary(3) ! bary center of element
! ---------------------------------------------------------------------------
! run over the given element list to update
do iElem = 1, nElems
missingNeigh = .false.
tID = tree%treeID( list( iElem ))
! Get the position of the treeID in the element list
elemPos = PositionOfVal( me = levelDesc%elem%tID, val = tID )
if( elemPos > 0 ) then
totalPos = tem_treeIDinTotal( &
& tID = tID, &
& levelDesc = levelDesc, &
& eType = levelDesc%elem%eType%val( elemPos ))
stencilPos = tem_stencil_getHeaderPos( &
& me = levelDesc%elem%stencil%val( elemPos ), &
& val = iStencil )
if( stencilPos > 0 ) then
do iStencilElem = 1, levelDesc%elem%stencil%val(elemPos) &
& %val(stencilPos)%QQN
! Totalpos has before already been updated, so just read from
! totalPos
levelPos = levelDesc%elem%stencil%val(elemPos) &
& %val(stencilPos) &
& %totalPos(iStencilElem)
if( levelPos <= 0 .and. computeStencil%requireNeighNeigh ) then
missingNeigh = .true.
end if
levelDesc%neigh( iStencil )%nghElems( iStencilElem, iElem ) &
& = levelPos
end do
end if
if( missingNeigh ) then
bary = tem_BaryOfId(tree, tID)
write(logUnit(10),"(A,I0,A,I0,A,3ES12.5)") &
& 'missing neighbor iElem: ', iElem,&
& ', tID: ', tID, ', bary: ', bary
! Set the property of the current element to solid
levelDesc%property( totalPos ) &
& = ibset( levelDesc%property( totalPos ), prp_solid )
end if
else ! elemPos <= 0, i.e. tID can not be found in levelDesc%elem%tID
write(dbgUnit(1),*) 'Warning!!!'
write(dbgUnit(1),*) tID, ' not found during build_treehorizontalDep'
write(dbgUnit(1),*) 'That should not occur.'
end if
end do
end subroutine tem_build_treeHorizontalDep