call this routine from your geometry initialization routine in the solver create all the necessary level-wise objects, such as element lists, dependencies
1.) build all dependencies for halos and ghost which are needed for interpolation/reconstruction (including MPI communications) 2.) build the pointers for each element to its neighbors/stencil elements. All this information is stored in tem_levelDesc_type
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(treelmesh_type), | intent(inout) | :: | tree |
the global tree |
||
type(tem_levelDesc_type), | intent(inout) | :: | levelDesc(tree%global%minLevel:) |
the level descriptor to be filled |
||
integer, | intent(out), | allocatable | :: | levelPointer(:) |
Pointer from treeIDlist entry to level-wise fluid part of total list |
|
type(tem_stencilHeader_type) | :: | computeStencil(:) |
array of all stencils used in the simulation |
|||
type(tem_commPattern_type), | intent(in) | :: | commPattern |
the communication pattern used |
||
logical, | intent(in), | optional | :: | cleanup |
cleanup arrays afterwards? |
|
integer, | intent(in), | optional | :: | reqNesting |
nesting level |
|
type(tem_comm_env_type), | intent(in) | :: | proc |
Process description to use. |
subroutine tem_find_allElements( tree, levelDesc, levelPointer, &
& computeStencil, commPattern, cleanup, &
& reqNesting, proc )
! -------------------------------------------------------------------- !
!> the global tree
type(treelmesh_type), intent(inout) :: tree
!> the level descriptor to be filled
type(tem_levelDesc_type), intent(inout) :: levelDesc(tree%global%minLevel:)
!> Pointer from treeIDlist entry to level-wise fluid part of total list
integer, allocatable, intent(out) :: levelPointer(:)
!> array of all stencils used in the simulation
type(tem_stencilHeader_type) :: computeStencil(:)
!> the communication pattern used
type(tem_commPattern_type), intent(in) :: commPattern
!> nesting level
integer, intent(in), optional :: reqNesting
!> cleanup arrays afterwards?
logical, intent(in), optional :: cleanup
!> Process description to use.
type(tem_comm_env_type), intent(in) :: proc
! -------------------------------------------------------------------- !
! Geometry and Tree related variables
integer :: iLevel
type( tem_path_type ), allocatable :: pathFirst(:), pathLast(:)
logical :: doAdditional
logical :: clean_constructionArrays
! -------------------------------------------------------------------- !
call tem_horizontalSpacer( fUnit = dbgUnit(1), before = 1 )
write(dbgUnit(1),*) 'Inside routine: tem_find_allElements'
if( present( cleanup )) then
clean_constructionArrays = cleanup
else
clean_constructionArrays = .false.
end if
if( present( reqNesting )) then
nestingLimit = reqNesting
else
nestingLimit = 1
end if
write(logUnit(3),*) 'Building up the element and total list.'
! Build the pathlist for improved performance when finding local element
! positions since the tree is the same for every scheme the allocation only
! has to be performed once
if ( .not. allocated(tree%pathList))then
allocate( tree%pathList(tree%nElems) )
end if
! Fill the pathList for each element in the treeID list
tree%pathList = tem_PathOf( tree%treeID )
! first and last path in every process
allocate(pathFirst( tree%global%nParts ))
allocate(pathLast( tree%global%nParts ))
pathFirst = tem_PathOf( tree%Part_First )
pathLast = tem_PathOf( tree%Part_Last )
! Step 2: build levelDesc element list including identification of neighbor
! elements
call build_levelElements( levelDesc = levelDesc, &
& tree = tree, &
& proc = proc, &
& stencil = computeStencil(1), &
& pathFirst = pathFirst, &
& pathLast = pathLast )
! Step 3: assign totalPnt to elem%tID in sorted fashion and prepare
! haloPrc list
do iLevel = tree%global%minLevel, tree%global%maxLevel
call identify_lists( levelDesc(iLevel) )
end do
! Step 4: Communicate halo elements
! Here, first local halos are send to their corresponding process to
! check for existence of halo elements.
! Each process checks whether requested element exist and returns
! list of available elements.
! With this new information halo list is redefined.
call communicate_elements( tree, proc, levelDesc, commPattern, &
& pathFirst, pathLast, computeStencil )
! Step 5: do additional communication if there are elements in require list
! which are neighbors of higher order boundaries.
call check_additionalComm( levelDesc, proc, doAdditional, &
& tree%global%minlevel )
! If doAdditional then identify neighbors of higher order boundary
! neighbor elements.
! After this identification, new halo elements might have to be added so
! we need to communicate all halo elements again.
if( doAdditional ) then
! passing only first stencil as this is the required compute stencil
call identify_additionalNeigh( tree, proc, levelDesc, &
& pathFirst, pathLast, computeStencil(1) )
call communicate_elements( tree, proc, levelDesc, commPattern, &
& pathFirst, pathLast, computeStencil )
end if
! Step 6: assemble levelDesc total(treeID) list, property list and
! pntTID list in sorted fashion (fluids+ghosts+halos)
! which are pre-assembled in element type
call assemble_lists( levelDesc, &
& tree%global%minLevel, tree%global%maxLevel, &
& tree )
! Step 7:
call tem_build_levelPointer( levelPointer, tree, levelDesc )
call update_elemPosToTotalPos( levelDesc, levelPointer, tree, &
& computeStencil )
! Warning: Truncation introduces a memory peak because of copy
! operations! Better not use...
!call truncate_lists( levelDesc, tree%global%minLevel )
if( clean_constructionArrays ) then
call tem_cleanupDependencyArrays( levelDesc )
endif
write(logUnit(3),*) 'Finished building up element and total list. '
deallocate(pathFirst)
deallocate(pathLast)
deallocate(hash)
deallocate(hash_elempos)
write(dbgUnit(1),*) 'Leave routine: tem_find_allElements'
call tem_horizontalSpacer( fUnit = dbgUnit(1), after = 1 )
end subroutine tem_find_allElements