This routine creates subTree from geometry intersection
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(tem_shape_type), | intent(in) | :: | me |
shape 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 |
||
integer, | intent(inout) | :: | countPoints |
How many points there will be |
||
type(tem_grwPoints_type), | intent(inout) | :: | grwPnts |
growing array to store tracking points |
||
logical, | intent(in) | :: | storePnts |
to Store points in grwPnts |
||
type(dyn_intarray_type), | intent(inout) | :: | map2global |
growing array for the map2global |
subroutine tem_shape_subTreeFromGeomInters( me, inTree, countElems, &
& countPoints, grwPnts, &
& storePnts, map2global )
! --------------------------------------------------------------------------
!> shape objects on which to work
type(tem_shape_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 )
!> How many points there will be
integer, intent( inout ) :: countPoints
!> growing array for the map2global
type(dyn_intArray_type), intent(inout) :: map2global
!> growing array to store tracking points
type(tem_grwPoints_type), intent(inout) :: grwPnts
!> to Store points in grwPnts
logical, intent(in) :: storePnts
! --------------------------------------------------------------------------
real(kind=rk) :: tStart, tEnd
integer :: iElem, dPos, tLevel, iObj
logical :: wasAdded, intersects, addToSubTree
logical :: foundAny, foundGlobally
type(tem_cube_type) :: cube
integer(kind=long_k) :: treeID
integer :: iError
! --------------------------------------------------------------------------
write(logUnit(5),"(A)") ' Extracting subTree from geometrical shape' &
& // ' - tree intersection'
tStart = mpi_wtime()
select case (trim(me%kind))
case ('canoND')
! treat canoND seperately to efficiently handle cano kind = points
call tem_cano_initSubTree( me = me%canoND(:), &
& inTree = inTree, &
& countElems = countElems, &
& map2global = map2global, &
& shapeInverted = me%inverted )
foundany = (sum(countElems) > 0)
call MPI_Allreduce( foundAny, foundGlobally, 1, &
& MPI_LOGICAL, MPI_LOR, inTree%global%comm, iError )
if (storePnts .and. (.not. foundGlobally)) then
write(logUnit(5),"(A)") ' WARNING: did not find any element for canonical shape'
write(logUnit(5),"(A)") ' looking at potential neighboring elements'
! No elements were found for the canoND objects, check whether
! any neighboring element is close, to extrapolate nearby points.
call tem_cano_checkNeigh( me = me%canoND(:), &
& inTree = inTree, &
& countElems = countElems, &
& map2global = map2global )
write(logUnit(5),*) ' --> after considering neighbors found ', sum(countElems)
write(logUnit(5),*) ' elements on this process (', inTree%global%myPart, ')'
end if
case ('triangle', 'stl', 'sphere', 'ellipsoid', 'cylinder')
do iElem = 1, inTree%nElems
treeID = inTree%treeID(iElem)
call tem_convertTreeIDtoCube(cube, inTree, treeID)
intersects = .false.
! Treat each object of same geometry kind
select case (trim(me%kind))
case ('triangle')
do iObj = 1, size(me%triangle)
intersects = intersects &
& .or. tem_triangleCubeOverlap(me%triangle(iObj), cube)
end do
case ('stl')
intersects = tem_stlCubeOverlap(me%stl_data, cube)
case ('sphere')
do iObj = 1, size(me%sphere)
intersects = intersects &
& .or. tem_sphereCubeOverlap(me%sphere(iObj), cube)
end do
case ('ellipsoid')
do iObj = 1, size(me%ellipsoid)
intersects = intersects &
& .or. tem_ellipsoidCubeOverlap(me%ellipsoid(iObj), cube)
end do
case ('cylinder')
do iObj = 1, size(me%cylinder)
intersects = intersects &
& .or. tem_cylinderCubeOverlap(me%cylinder(iObj), cube)
end do
end select
addToSubTree = .false.
if (.not. me%inverted .and. intersects) then
! Shape intersects with current element and not inverted
addToSubTree = .true.
else if (me%inverted .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
case default
call tem_abort('In tem_shape_subTreeFromGeomInters: Unknown shape kind')
end select
countPoints = 0
call init(me=grwPnts)
if ( storePnts .and. (map2global%nVals > 0) ) then
select case (trim(me%kind))
case ('canoND')
call tem_cano_storePntsInSubTree( me = me%canoND(:), &
& inTree = inTree, &
& countPoints = countPoints, &
& grwPnts = grwPnts, &
& map2global = map2global )
case default
call tem_abort('Use get_points supported only for canoND')
end select
end if
tEnd = mpi_wtime()
write(logunit(4),"(A,E12.6)") 'Done. This process cost: ', tEnd-tStart
end subroutine tem_shape_subTreeFromGeomInters