Create a subTree from a given inTree and an array of shapes.
All elements to the inTree%treeID are collected into the subTree%map2global array and a sub-communicator with participating processes is created. Additionally how many and which elements exist on my local process and are requested from the shapes is identified.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(treelmesh_type), | intent(in) | :: | inTree | Global mesh from which the elements are identified and then stored to |
||
type(tem_subTree_type), | intent(out) | :: | subTree | new mesh |
||
type(tem_shape_type), | intent(in) | :: | inShape(:) | shape objects on which to work |
||
type(tem_levelDesc_type), | intent(in), | optional | :: | levelDesc(:) | optional level descriptor needed for local shape |
|
type(tem_BC_prop_type), | intent(in) | :: | bc_prop | bc property which is used to identify elements belong to certain BCs |
||
logical, | intent(in), | optional | :: | storePnts | To store space points in subTree |
|
type(tem_stencilHeader_type), | intent(in), | optional | :: | stencil | stencil used to find bcID on certain links |
|
character(len=*), | intent(in), | optional | :: | prefix | prefix for the subTree |
subroutine tem_create_subTree_of( inTree, subTree, inShape, levelDesc, &
& bc_prop, storePnts, stencil, prefix )
! ---------------------------------------------------------------------------
!> Global mesh from which the elements are identified and then stored to
type(treelmesh_type), intent(in) :: inTree
!> new mesh
type(tem_subTree_type), intent(out) :: subTree
!> shape objects on which to work
type(tem_shape_type ),intent(in) :: inShape(:)
!> optional level descriptor needed for local shape
type(tem_levelDesc_type ), optional, intent(in) :: levelDesc(:)
!> bc property which is used to identify elements belong to certain BCs
type( tem_bc_prop_type ), intent(in) :: bc_prop
!> To store space points in subTree
logical, optional, intent(in) :: storePnts
!> stencil used to find bcID on certain links
type( tem_stencilHeader_type ), optional, intent(in) :: stencil
!> prefix for the subTree
character(len=*), optional, intent(in) :: prefix
! ---------------------------------------------------------------------------
! local growing array for the map to the global tree
type(dyn_intArray_type) :: local_map2global
type(tem_grwPoints_type) :: local_grwPnts
integer, allocatable :: map2globalTemp(:)
integer(kind=long_k), allocatable :: tmp_treeID(:)
integer :: iShape, iElem
! mpi variables for parallel writing of the log-wise trees
integer :: color, iError
! if the local rank has part of the tracking object
logical :: participateInMesh, globalParticipateInMesh
! if any element be found for this iShape
logical :: foundAny
integer :: local_countElems( globalMaxLevels )
integer :: local_countPnts
logical :: local_storePnts
! --------------------------------------------------------------------------
if (present(storePnts)) then
local_storePnts = storePnts
else
local_storePnts = .false.
end if
subTree%useGlobalMesh = .false.
subTree%useLocalMesh = .false.
subTree%created_new_comm = .false.
! Prepare output structure
! init here the sub-tree mesh
subTree%global = inTree%global
! initialize the number of elements on this partition with 0
subTree%nElems = 0
subTree%global%nElems = 0
if ( present( prefix )) subTree%global%dirname = trim(prefix)
! Reset counter for actual elements for this shape object
! @todo: indeed only sum(local_countElems) is used later.
! thus no need to calculate countElems in following routines.
local_countElems = 0
local_countPnts = 0
call init( local_map2global, length = 128)
do iShape = 1, size( inShape )
select case( inShape(iShape)%shapeID )
! Use elements belonging to a canonical n-dimensional shape
case( tem_geometrical_shape )
call tem_log(5, 'iShape '//trim(tem_toStr(iShape)) &
& //' is geometrical shape')
call tem_shape_subTreeFromGeomInters( me = inShape(iShape), &
& inTree = inTree, &
& countElems = local_countElems, &
& countPoints = local_countPnts, &
& grwPnts = local_grwPnts, &
& storePnts = local_storePnts, &
& map2global = local_map2global )
! Use all elements for this shape
case( tem_global_shape )
call tem_log(5, 'iShape '//trim(tem_toStr(iShape))//' is a global shape')
! Global shape means, all elements are part of the shape.
! so nothing to do, just continue with the next shape
local_countElems = 1
subTree%useGlobalMesh = .true.
exit
! Only use elements with a certain property
case( tem_property_shape )
call tem_log(5, 'iShape '//trim(tem_toStr(iShape))//' is a property shape')
call tem_shape_initPropElements( inShape(iShape)%propBits, inTree, &
& local_countElems, local_map2global )
! Only use elements belong to certain boundaries
case( tem_boundary_shape )
call tem_log(5, 'iShape '//trim(tem_toStr(iShape))//' is a boundary shape')
call tem_shape_findElemByBCLabels( bcLabels = inShape(iShape)%bcLabels,&
& bc_prop = bc_prop, &
& foundAny = foundAny, &
& map2global = local_map2Global, &
& bcIDs = subTree%bc_ID, &
& stencil = stencil )
if ( foundAny ) local_countElems = 1
! elements on the local partition
case( tem_local_shape )
call tem_log(5, 'iShape '//trim(tem_toStr(iShape))//' is a local shape')
if(present(levelDesc)) then
subTree%useLocalMesh = .true.
call tem_shape_initLocal( levelDesc, tmp_treeID )
else
write(logUnit(1),*)'ERROR:levelDesc is not passed to '// &
& 'tem_creat_subTree_of'
write(logUnit(1),*)' levelDesc is required to create '// &
& 'subTree for tem_local_shape'
call tem_abort()
end if
local_countElems = 1
! update the communicator
subTree%global%comm = mpi_comm_self
exit
case( tem_level_shape )
call tem_log(5, 'iShape '//trim(tem_toStr(iShape))//' is a level shape')
call tem_shape_initByLevels( inTree, inShape(iShape)%minLevel, &
& inShape(iShape)%maxLevel, &
& local_countElems, local_map2global )
end select ! shapeID
end do ! iShape
if( .not. subTree%useGlobalMesh ) then
subTree%global%predefined = ''
end if
! Exchange information about the mesh
! mpi_exscan: track( iLog )%tree%offset
! Now create a new communicator with ranks that own parts of the sub-mesh
! This is done with the color
participateInMesh = ( sum( local_countElems ) > 0 .or. local_countPnts > 0 )
! abort if no elements are found for the given shape
call MPI_ALLREDUCE( participateInMesh, globalParticipateInMesh, 1, &
& MPI_LOGICAL, MPI_LOR, inTree%global%comm, iError )
if (.not. globalParticipateInMesh) then
write(logUnit(1),*) 'Error: No elements found for '// &
& 'defined shape '//trim(prefix)//' in the mesh'
call tem_abort()
end if
if( participateInMesh ) then
color = 1
else
! Not participating, set color to undefined
color = MPI_UNDEFINED
end if
! Input rank as a key to keep sorting as in global communicator
! reorder process ranks into tracking communicator
if( .not. subTree%useLocalMesh ) then
call mpi_comm_split( inTree%global%comm, color, inTree%global%myPart, &
& subTree%global%comm, iError )
if( iError .ne. mpi_success ) then
write(logUnit(1),*)' There was an error splitting the communicator'
write(logUnit(1),*)' for the tracking subset. Aborting...'
call tem_abort()
end if
subTree%created_new_comm = .true.
end if
if( participateInMesh ) then
call mpi_comm_size( subTree%global%comm, subTree%global%nParts, iError )
call mpi_comm_rank( subTree%global%comm, subTree%global%myPart, iError )
if( .not. subTree%useGlobalmesh) then
if( subTree%useLocalMesh ) then
call tem_subTree_from( me = subTree, &
& treeID = tmp_treeID )
! update the property list of the newly defined tree according to the
! the local shape (simply copy the bits)
call tem_copyPropertyBits( levelDesc = levelDesc, subTree = subTree )
else ! NOT useLocalMesh
! Copy the created dynamic local_tIDlist to a long integer array
allocate( map2globalTemp( local_map2global%nVals ))
! Assign list of links to the temporary list, first sorting
do iElem = 1, local_map2global%nVals
map2globalTemp( iElem ) = local_map2global%val( &
& local_map2global%sorted( iElem ))
end do
call tem_subTree_from( me = subTree, &
& map2global = map2globalTemp, &
& grwPnts = local_grwPnts )
call destroy( me = local_grwPnts )
! copy property bits from the global tree
call tem_copyPropertyBits( inTree = inTree, subTree = subTree )
deallocate( map2globalTemp )
end if ! useLocalMesh with level desc
! Calculate the effective extents of the tracking sub-tree bounding box
call tem_setEffBoundingBox( subTree, inTree )
else ! useGlobalMesh
! assign the global properties to the subTree
subTree%nElems = inTree%nElems
subTree%global%nElems = inTree%global%nElems
end if ! not global mesh?
if (subTree%global%myPart == 0) then
if (subTree%global%nElems == 1 .and. allocated(subTree%map2global) ) then
write(*,"(A,3F12.6)") 'Barycenter of tracking element: ', &
& tem_baryOfId(inTree, inTree%treeID(subTree%map2global(1)))
end if
end if !subTree root
else ! NOT participate in mesh
subTree%nElems = 0
subTree%global%nElems = 0
allocate( subTree%map2global(0) )
subTree%nPoints = 0
subTree%glob_nPoints = 0_long_k
allocate( subTree%points(0,0) )
end if ! Participate in mesh?
write(logUnit(5),"(A,I0)") 'Found corresponding nElems: ', &
& subTree%global%nElems
call destroy( me = local_map2global )
call tem_horizontalSpacer( fUnit= logUnit(5) )
end subroutine tem_create_subTree_of