Find all the virtual children of the sourceID down to the targetLevel and add to the level-wise ghostFromCoarser list in the level descriptor
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer(kind=long_k), | intent(in) | :: | sourceID |
source treeID (existing founded ID in tree%treeID list or children ID from recursion) |
||
integer(kind=long_k), | intent(in) | :: | sourceProperty |
property of source element |
||
integer, | intent(in) | :: | foundPos |
position of this sourceID in elem%tID list |
||
type(tem_path_type), | intent(in) | :: | elemPath |
element path |
||
integer, | intent(in) | :: | targetLevel |
level upto which virtual children must be created |
||
type(tem_levelDesc_type), | intent(inout) | :: | levelDesc(minLevel:) |
the level descriptor to be filled |
||
integer, | intent(in) | :: | minLevel |
minimum level in the tree |
||
type(treelmesh_type), | intent(in) | :: | tree |
tree information |
||
type(tem_stencilHeader_type), | intent(in) | :: | stencil |
current stencil definition |
||
integer, | intent(in) | :: | nesting |
nesting level |
||
logical, | intent(out) | :: | updated |
was the element updated in this call? |
recursive subroutine add_all_virtual_children( & & sourceID, sourceProperty, foundPos, elemPath, & & targetLevel, levelDesc, minLevel, tree, & & stencil, nesting, updated ) ! -------------------------------------------------------------------- ! !> source treeID (existing founded ID in tree%treeID list or children ID !! from recursion) integer(kind=long_k), intent(in) :: sourceID !> property of source element integer(kind=long_k), intent(in) :: sourceProperty !> element path type(tem_path_type), intent(in) :: elemPath !> nesting level integer, intent(in) :: nesting !> position of this sourceID in elem%tID list integer, intent(in) :: foundPos !> level upto which virtual children must be created integer,intent(in) :: targetLevel !> minimum level in the tree integer, intent(in) :: minLevel !> the level descriptor to be filled type(tem_levelDesc_type), intent(inout) :: levelDesc(minLevel:) !> tree information type(treelmesh_type), intent(in) :: tree !> current stencil definition type(tem_stencilHeader_type), intent(in) :: stencil !> was the element updated in this call? logical, intent(out) :: updated ! -------------------------------------------------------------------- ! integer :: targetPos, iChild, iDir, nVals, sourceLevel ! position of the existing (source) tID in the elem list integer :: sourcePos logical :: wasAdded logical :: childUpdated integer(kind=long_k) :: cTreeID ! current treeID to identify integer(kind=long_k) :: ichildID ! child treeID integer(kind=long_k) :: curNeighborID integer(kind=long_k), allocatable :: tNeighID(:) type(tem_stencilElement_type) :: tStencil(1) integer :: iChildCoord(4), curLevel, offset(4), xc(4), childCoord(4) integer :: addedPos ! -------------------------------------------------------------------- ! !Position of the coarser source element sourceLevel = tem_LevelOf(sourceID) sourcePos = PositionOfVal( me = levelDesc( sourceLevel )%elem%tID, & & val = sourceID ) allocate( tNeighID(stencil%QQN) ) offset(4) = 0 ! By default, set that no element was updated updated = .false. childUpdated = .false. ! create virual children until target level is reached if ( sourceLevel < targetLevel ) then curLevel = sourceLevel + 1 call init( me = tStencil(1), QQN = stencil%QQN, headerPos = 1 ) ! Add to the level-wise ghost list cTreeID = elemPath%node( targetLevel - sourceLevel ) call append( me = levelDesc( curLevel )%elem, & & tID = cTreeID, & & property = sourceProperty, & & eType = eT_ghostFromCoarser, & & sourceProc = tree%global%myPart+1, & & stencilElements = tStencil, & & pos = targetPos, & & haloNesting = nesting, & & wasAdded = wasAdded ) if (wasAdded) then write(dbgUnit(7),"(2(A,I0))") 'Added as a GhostFromCoarser: ', & & cTreeID, ', to level: ', curLevel updated = .true. iChild = tem_childNumber(cTreeID) iChildCoord(:) = tem_coordOfId( int(iChild, long_k) ) ! inherit the boundary infos from the parent tNeighID = 0_long_k ! Get neighIds of (source) coarse element do iDir = 1, 3 call tem_find_BCs_fromCoarser( dir = iDir, & & childCoord = iChildCoord, & & sourceLevel = sourceLevel, & & sourcePos = sourcePos, & & neighID = tNeighID, & & minLevel = minLevel, & & levelDesc = levelDesc, & & computeStencil = Stencil ) end do ! loop over all directions to determine neighIDs for ghost child ! (targetPos) do iDir = 1, stencil%QQN ! compute virtual neighbor child from coarser neighbor in iDir ! direction if( tNeighID(iDir) > 0_long_k ) then ! find the corresponding children of the neighbor defined for my ! parent. Find the child in the corresponding direction ! get the corresponding neighbor in the offset direction, which ! is a child of the coarser neighbor element. Use +2 offset to get ! only positive numbers ! (child coordinates always range between {0..1}) offset(1:3) = stencil%cxDir(1:3, iDir ) xc(1:3) = mod( iChildCoord(1:3) + 2 + offset(1:3), 2 ) xc(4) = 1 ! get the child number of the corresponding neighbor within the ! coarser element. (is {1...8}) ichildID = tem_IdOfCoord( xc ) ! calculate the child treeID from the coarser neighbor and the ! childID curNeighborID = tNeighID(iDir)*8_long_k + ichildID else if (tNeighID(iDir) == 0_long_k) then ! virtual neighbor child exist in current parent. ! If the neighbor is still 0, it must be a direct sibling, ! so we can directly compute its tID. ! find the neighbor according to stencil offset from me childCoord = tem_coordOfId( cTreeID ) offset(1:3) = stencil%cxDir(1:3, iDir ) curNeighborID = tem_IdOfCoord( childCoord + offset ) else ! tNeighID(iDir) < 0_long_k ! inherit the boundary ID curNeighborID = tNeighID(iDir) end if ! append the neighbor ID ... call append( me = levelDesc( curLevel )%elem%neighID & & %val( targetPos ), & & val = curNeighborID, & & pos = addedPos ) ! ... and store this position in the stencil levelDesc( curLevel )%elem%stencil%val( targetPos ) & & %val(1)%tIDpos( iDir ) & & = addedPos ! to append the neighbors of ghosts to required list !call append( me = levelDesc( curLevel )%require, & ! & val = curNeighborID, & ! & pos = addedPos, & ! & wasAdded = wasAdded ) end do !iDir QQN ! Set prp_hasBnd if any of neighbors is boundary nVals = levelDesc( curLevel )%elem%neighID%val( targetPos )%nVals if ( minval( levelDesc(curLevel)%elem%neighID%val(targetPos) & & %val(1:nVals) ) < 0 ) then ! Found boundary levelDesc( curLevel )%elem%property%val( targetPos ) & & = ibset( levelDesc( curLevel )%elem%property%val(targetPos), & & prp_hasBnd ) else ! Unset the property bit levelDesc( curLevel )%elem%property%val( targetPos ) & & = ibclr( levelDesc( curLevel )%elem%property%val(targetPos), & & prp_hasBnd ) end if else ! Existing element encountered. updated = .false. end if ! Overwrite the eventually existing nesting with the smallest value. ! The smallest nesting determines if further neighbors have to be ! retrieved in communicate_elements if ( nesting < levelDesc(curLevel)%elem%haloNesting & & %val(targetPos) ) then ! needs update updated = .true. levelDesc( curLevel )%elem%needsUpdate%val( targetPos ) = .true. levelDesc( curLevel )%elem%haloNesting%val( targetPos ) & & = min( nesting, & & levelDesc( curLevel )%elem%haloNesting%val(targetPos) ) end if ! In any case we have to recurse down to the target level ! lower levels might not yet exist. call add_all_virtual_children( sourceID = cTreeID, & & foundPos = foundPos, & & elemPath = elemPath, & & nesting = nesting, & & targetLevel = targetLevel, & & levelDesc = levelDesc, & & minLevel = minLevel, & & updated = childUpdated, & & tree = tree, & & sourceProperty = sourceProperty, & & Stencil = Stencil ) end if ! sourceLevel < targetLevel updated = ( updated .or. childUpdated ) end subroutine add_all_virtual_children