Building neighor array
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | iStencil |
Index of your neighbor list. |
||
type(tem_levelDesc_type), | intent(inout) | :: | levelDesc(tree%global%minLevel:) |
Level descriptor for each level of your mesh (starting from min level). |
||
type(treelmesh_type), | intent(in) | :: | tree |
Tree representation of your mesh. |
||
type(tem_stencilHeader_type) | :: | computeStencil |
The stencil you build the horizontal dependencies for. |
subroutine tem_debug_HorizontalDependencies( iStencil, levelDesc, tree, &
& computeStencil )
! ---------------------------------------------------------------------------
!> Tree representation of your mesh.
type(treelmesh_type), intent(in) :: tree
!> Level descriptor for each level of your mesh (starting from min level).
type(tem_levelDesc_type), intent(inout) :: levelDesc(tree%global%minLevel:)
!> Index of your neighbor list.
integer, intent(in) :: iStencil
!> The stencil you build the horizontal dependencies for.
type(tem_stencilHeader_type) :: computeStencil
! ---------------------------------------------------------------------------
integer :: iElem, iLevel, nElems, nEntries, iVal, elemPos, tIDpos
integer :: neighPos, stencilSize, stencilPos
integer(kind=long_k) :: tID, bufID, refID
! ---------------------------------------------------------------------------
call tem_horizontalSpacer( fUnit = dbgUnit(1), before = 1 )
write(dbgUnit(1),"(A )") 'Get into tem_debug_horizontalDependencies routine'
write(dbgUnit(1),"(A,I3)") 'Check horizontal dependencies stencil:',iStencil
write(logUnit(1),*)' Checking, if the horizontal dependencies were set'// &
& ' correctly'
do iLevel = tree%global%minLevel,tree%global%maxLevel
nEntries = size( levelDesc( iLevel )%neigh(iStencil)%nghElems, 1)
nElems = size( levelDesc( iLevel )%neigh(iStencil)%nghElems, 2)
write(dbgUnit(1),"(3(A,I0))") ' iLevel: ', iLevel, &
& ' nEntries: ', nEntries, ', nElems: ', nElems
if( nEntries > 0 ) then
write(dbgUnit(1),"(A11, A10)") 'fluid, ', 'neighbors:'
do iElem = 1, nElems
! First get the treeID
if( computeStencil%useAll ) then
tIDPos = iElem
tID = levelDesc( iLevel )%total( tIDpos )
else
if( computeStencil%elemLvl( iLevel )%nVals <= 0 ) exit
tIDPos = computeStencil%elemLvl( iLevel )%val( iElem )
tID = tree%treeID( tIDpos )
end if
! Get the position of the treeID in the element list
elemPos = PositionOfVal(levelDesc( iLevel )%elem%tID, tID )
if( elemPos <= levelDesc( iLevel )%elem%nElems(eT_fluid) )then
! find the right stencil
stencilPos = tem_stencil_getHeaderPos( &
& me = levelDesc( iLevel )%elem%stencil%val( elemPos ), &
& val = iStencil )
else
! set the stencil position to be the fluid stencil
stencilPos = 1
end if
!write(*,*)'in debug_HorDeps, level ', iLevel,' elemPos: ', elemPos
!write(*,*)'stencilPos: ', stencilPos, ' iStencil: ', iStencil
!write(*,*)'nVals: ', levelDesc( iLevel )%elem%stencil%val( elemPos )%nVals
!write(*,*)'val%pos', levelDesc( iLevel )%elem%stencil%val( elemPos )%val(:)%headerPos
stencilSize = size( levelDesc( iLevel )%elem%stencil%val( elemPos ) &
& %val( stencilPos )%tIDpos )
!write(*,*)'Stencil entries: ', nEntries, ' stencilSize: ', stencilSize
if (nEntries /= stencilSize) then
write(dbgUnit(1),*)'stencil entries do not match for tID ',tID, &
& ': ', nEntries, ' vs. ', stencilSize
end if
write( dbgUnit(3), '(I10,A)', advance='no' ) tID, ':'
do iVal = 1, nEntries
! required neighbor ID
tIDpos = levelDesc( iLevel )%elem%stencil%val( elemPos )%val( &
& stencilPos )%tIDpos( iVal )
if( tIDpos > 0 ) then
refID = levelDesc( iLevel )%elem%neighID%val( elemPos )%val( &
& tIDpos )
else
refID = 0
end if
! actually obtained neighbor ID
neighPos = levelDesc( iLevel )%neigh(iStencil)% &
& nghElems( iVal, iElem )
if( neighPos > 0 ) then
bufID = levelDesc( iLevel )%total( neighPos )
else
bufID = neighPos
endif
write(dbgUnit(3),'(I10)',advance='no') bufID
if( refID /= bufID ) then
write(dbgUnit(3),'(a,i0,a)', advance='no') ' (',refID,')'
if( bufID > 0_long_k ) then
write(dbgUnit(1),*)' Non matching stencil entries!!'
write(dbgUnit(1),*)' original treeID: ', tID
write(dbgUnit(1),*)' required stencil treeID: ', refID
write(dbgUnit(1),*) ' encountered treeID: ', bufID
write(logUnit(1),*)' Found non-matching treeID in the neighborhood'
write(logUnit(1),*)' original treeID: ', tID
write(logUnit(1),*)' required stencil treeID: ', refID
write(logUnit(1),*)' encountered treeID: ', bufID
call tem_abort()
end if
end if
end do ! iVal = 1, nEntries
write(dbgUnit(3),*)''
end do ! iElem
write(dbgUnit(3),*)''
end if
end do ! iLevel
write(dbgUnit(1),*) 'Leave tem_debug_HorizontalDependencies'
call tem_horizontalSpacer( fUnit = dbgUnit(1), after = 1 )
end subroutine tem_debug_HorizontalDependencies