Create subtree from the intersection of canonical shapes and inTree
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(tem_canonicalND_type), | intent(in) | :: | me(:) |
canonicalND objects on which to work |
||
type(treelmesh_type), | intent(in) | :: | inTree |
Global tree |
||
integer, | intent(inout) | :: | countElems(globalMaxLevels) |
How many elements there will be for each level in the track |
||
type(dyn_intarray_type), | intent(inout) | :: | map2global |
growing array for the map2global |
||
logical, | intent(in) | :: | shapeInverted |
If true then elements not intersected are added to subTree |
subroutine tem_cano_initSubTree( me, inTree, countElems, map2global, &
& shapeInverted )
! --------------------------------------------------------------------------
!> canonicalND objects on which to work
type(tem_canonicalND_type ),intent(in) :: me(:)
!> Global tree
type(treelmesh_type), intent(in) :: inTree
!> How many elements there will be for each level in the track
integer, intent( inout ) :: countElems( globalMaxLevels )
!> growing array for the map2global
type(dyn_intArray_type), intent(inout) :: map2global
!> If true then elements not intersected are added to subTree
logical, intent(in) :: shapeInverted
! --------------------------------------------------------------------------
integer :: iCano, tLevel, elemPos, iElem, dPos, maxLevel
integer(kind=long_k) :: treeID
logical :: wasAdded, intersects, addToSubTree
type(tem_cube_type) :: cube
! --------------------------------------------------------------------------
write(logUnit(4),"(A)") 'Initializing canonical shapes'
maxLevel = inTree%global%maxLevel
! Create subTree by intersecting canoND objects with the inTree
do iCano = 1, size(me)
! for a point there is no need to check against all elements
! instead just convert the point into treeID and check for its
! exitence in given tree
if (trim(me(iCano)%kind) == 'point') then
treeID = tem_IdOfCoord( tem_CoordOfReal(inTree, &
& me(iCano)%point%coord(1:3), &
& maxLevel) )
! get position of the treeID in inTree
elemPos = tem_PosOfId( treeID, inTree%treeID)
found_elem: if (elemPos > 0) then
! append the position in inTree to the map (note that already existing
! ones are omitted)
call append( me = map2global, &
& pos = dpos, &
& val = elemPos, &
& wasAdded = wasAdded )
! Count up if it was added
if (wasAdded) then
tLevel = tem_levelOf( treeID )
countElems( tLevel ) = countElems( tLevel ) + 1
end if ! wasAdded
end if found_elem
else
do iElem = 1, inTree%nElems
treeID = inTree%treeID(iElem)
call tem_convertTreeIDtoCube(cube, inTree, treeID)
intersects = .false.
select case ( trim(me(iCano)%kind) )
case ('line')
intersects = tem_lineCubeOverlap(me(iCano)%line, cube)
case ('plane')
intersects = tem_planeCubeOverlap(me(iCano)%plane, cube)
case ('box')
intersects = tem_boxCubeOverlap(me(iCano)%box, cube)
end select
addToSubTree = .false.
if (.not. shapeInverted .and. intersects) then
! Shape intersects with current element and not inverted
addToSubTree = .true.
else if (shapeInverted .and. .not. intersects) then
! shape not intersected and is inverted shape so add this to subTree
addToSubTree = .true.
end if
if (addToSubTree) then
! append iElem in inTree to the map (note that already existing
! ones are omitted)
call append( me = map2global, &
& pos = dpos, &
& val = iElem , &
& wasAdded = wasAdded )
! Count up if it was added
if( wasAdded ) then
tLevel = tem_levelOf( treeID )
countElems( tLevel ) = countElems( tLevel ) + 1
end if ! wasAdded
end if !intersects
end do !iElem
end if
end do !iCano
end subroutine tem_cano_initSubTree