identify_elements Subroutine

private recursive subroutine identify_elements(treeID, tree, pathFirst, pathLast, levelDesc, elemPos, proc, Stencil, nesting, skip_add_additionalGhost)

Check, on which partition a given element is located add required elements to corresponding lists: if remote, add to halo if ghost, add to resp. ghost list

Arguments

Type IntentOptional Attributes Name
integer(kind=long_k), intent(in) :: treeID

treeID to identify

type(treelmesh_type), intent(in) :: tree

tree information

type(tem_path_type), intent(in) :: pathFirst(:)

first treeID path in every process

type(tem_path_type), intent(in) :: pathLast(:)

last treeID path in every process

type(tem_levelDesc_type), intent(inout) :: levelDesc(tree%global%minLevel:)

the level descriptor to be filled

integer, intent(out) :: elemPos

nTreeID element position in the levelDesc % elem list

type(tem_comm_env_type), intent(in) :: proc

Process description to use.

type(tem_stencilHeader_type), intent(in) :: Stencil

current stencil definition

integer, intent(in) :: nesting

nesting level

logical, intent(in), optional :: skip_add_additionalGhost

logical, optional, if true no ghosts are added


Source Code

  recursive subroutine identify_elements( treeID, tree, pathFirst, pathLast,  &
    &                                     levelDesc, elemPos, proc,           &
    &                                     Stencil, nesting,                   &
    &                                     skip_add_additionalGhost )
    ! -------------------------------------------------------------------- !
    !> treeID to identify
    integer(kind=long_k), intent(in) :: treeID
    !> tree information
    type(treelmesh_type), intent(in) :: tree
    !> first treeID path in every process
    type(tem_path_type), intent(in) :: pathFirst(:)
    !> last treeID path in every process
    type(tem_path_type), intent(in) :: pathLast(:)
    !> the level descriptor to be filled
    type(tem_levelDesc_type), intent(inout) :: levelDesc(tree%global%minLevel:)
    !> Process description to use.
    type(tem_comm_env_type), intent(in) :: proc
    !> nTreeID element position in the levelDesc % elem list
    integer, intent(out) :: elemPos
    !> current stencil definition
    type(tem_stencilHeader_type), intent(in) :: Stencil
    !> nesting level
    integer, intent(in) :: nesting
    !> logical, optional, if true no ghosts are added
    logical, intent(in), optional :: skip_add_additionalGhost
    ! -------------------------------------------------------------------- !
    integer(kind=long_k) :: children(8) ! child elements
    integer :: nDepProcs ! processes that this neigh depend on
    integer :: depProc
    integer :: iChild, neighLevel
    type(tem_path_type) :: elemPath
    type(tem_stencilElement_type) :: emptyStencil(1)
    integer :: childPos   ! position of child element
    integer :: nNesting ! new nesting level
    ! Position of the current element in the hash
    integer :: hashpos, elemNesting
    logical :: cacheHit
    logical :: updated ! was the element updated during identify_local_element
    logical :: l_skip_add_additionalGhost
    ! -------------------------------------------------------------------- !
    
    if (present(skip_add_additionalGhost)) then 
      l_skip_add_additionalGhost = skip_add_additionalGhost
    else 
      l_skip_add_additionalGhost = .false.
    end if 
    
    ! Increase the nesting
    nNesting = nesting + 1
    cacheHit = .false.
    elemNesting = -999

    if (treeID > 0) then ! it is a element, otherwise it is a bcID
      neighLevel = tem_LevelOf(TreeID)
      elemPath  = tem_PathOf(TreeID)

      hashpos = int(mod(TreeID, nHashes))

      hashmatch: if (hash(hashpos) == TreeID) then

        if ( nesting == levelDesc(neighlevel)         &
          &               %elem                       &
          &               %haloNesting                &
          &               %val(hash_elemPos(hashPos)) ) then
          cacheHit = .true.
        end if
        if ( levelDesc(neighlevel)         &
          &    %elem                       &
          &    %needsUpdate                &
          &    %val(hash_elemPos(hashPos)) ) then
          cacheHit = .false.
          ! Set the needs update to false as it was now updated
          levelDesc(neighLevel)            &
            &  %elem                       &
            &  %needsUpdate                &
            &  %val(hash_elemPos(hashPos)) = .false.
        end if

      end if hashmatch

      cachemiss: if (.not. cacheHit) then
        ! Cache miss
        ! Probably did not hit this element yet, put it in the hash now,
        ! and identify it.
        ! (If this treeID is already in the hash, we definitely have
        !  already processed this treeID, and need nothing to do anymore).
        ! This leaves elements to be identified multiple times only in
        ! the case of hash collisions, which should be pretty few, for
        ! a sufficiently large hash.
        ! -----------------------------------------------
        call tem_find_depProc( depProc   = depProc,   &
          &                    nDepProcs = nDepProcs, &
          &                    tree      = tree,      &
          &                    elemPath  = elemPath,  &
          &                    PathFirst = PathFirst, &
          &                    PathLast  = PathLast   )

        ! How many processes possess a part of the requested treeID
        if (nDepProcs == 1) then
          updated = .false.
          ! Might be a local or halo of same level or a ghostFromCoarser
          ! or a ghostFromFiner. If halo and ghost it is distributed
          ! in single process.
          ! elemPos here is position of TreeID in levelDesc elem list
          call single_process_element( targetID       = TreeID,               &
            &                          levelDesc      = levelDesc,            &
            &                          tree           = tree,                 &
            &                          proc           = proc,                 &
            &                          iProc          = depProc,              &
            &                          minLevel       = tree%global%minlevel, &
            &                          elemPos        = elemPos,              &
            &                          updated        = updated,              &
            &                          stencil        = Stencil,              &
            &                          nesting        = nesting,               &
            &              skip_add_additionalGhost = l_skip_add_additionalGhost)

          if (elemPos > 0) then
          
            elemNesting = min( nesting,              &
              &                levelDesc(neighLevel) &
              &                  %elem               &
              &                  %haloNesting        &
              &                  %val(elemPos)       )

            if ( nesting < levelDesc(neighLevel)%elem%haloNesting  &
              &                                      %val(elemPos) ) then
              ! element needs updating as the current nesting is smaller than the
              ! one in the element property
              levelDesc(neighLevel)%elem%needsUpdate%val(elemPos) = .true.
              levelDesc(neighLevel)%elem%haloNesting%val(elemPos) &
                &  = elemNesting
              updated = .true.
            end if
            
          endif

          ! it is ghost from coarser element in current process proc%rank
          if (updated) then
            if ( levelDesc( neighLevel )%elem%eType%val( elemPos )   &
              &                               == eT_ghostFromCoarser ) then
              if (nesting < nestingLimit) then
                ! Create the direct neighbors of the ghostFromCoarser

                ! identify all the compute neighbors of the current element
                call identify_stencilNeigh( iElem          = elemPos,        &
                  &                         iLevel         = neighLevel,     &
                  &                         tree           = tree,           &
                  &                         iStencil       = 1,              &
                  &                         pathFirst      = pathFirst,      &
                  &                         pathLast       = pathLast,       &
                  &                         levelDesc      = levelDesc,      &
                  &                         proc           = proc,           &
                  &                         stencil        = Stencil,        &
                  &                         nesting        = elemNesting + 1 )

              end if ! nesting < 1?

              ! Create all elements required up to the actual existing fluid
              ! element these include the neighbors of the parents. In a level
              ! jump >1, these intermediate levels have to provide valid
              ! quantities over two of their computation updates to account for
              ! the recursive algorithm.
              call create_allParentNeighbors( &
                &    targetID  = levelDesc(neighLevel)%elem%tID%val(elemPos), &
                &    level     = neighLevel,     &
                &    tree      = tree,           &
                &    stencil   = Stencil, &
                &    levelDesc = levelDesc,      &
                &    pathFirst = pathFirst,      &
                &    pathLast  = pathLast,       &
                &    proc      = proc )

            end if ! ghostFromCoarser?
          end if ! updated?

        else if ( nDepProcs > 1 ) then
          ! more than one depending processes
          ! create ghost and find children (can only have finer children)
          call init( me        = emptyStencil(1),   &
            &        QQN       = Stencil%QQN,       &
            &        headerPos = 1 )
          ! Add to the level-wise ghost list
          ! append and store position of element for return
          call append( me              = levelDesc( neighLevel )%elem, &
            &          tID             = TreeID,                       &
            &          property        = 0_long_k,                     &
            &          eType           = eT_distributedGhostFromFiner, &
            &          stencilElements = emptyStencil,                 &
            &          pos             = elemPos                       )

          !... and store the position in the ghost list !
          ! Now find all children of this ghost
          ! childPos is return but not used anymore
          children = tem_directChildren( TreeID )
          do iChild = 1,8
            call identify_elements( TreeID    = children( iChild ), &
              &                     tree      = tree,               &
              &                     pathFirst = pathFirst,          &
              &                     pathLast  = pathLast,           &
              &                     levelDesc = levelDesc,          &
              &                     proc      = proc,               &
              &                     elemPos   = childPos,           &
              &                     stencil   = Stencil,            &
              &                     nesting   = nNesting,            &
              &  skip_add_additionalGhost = l_skip_add_additionalGhost)
          end do
        else ! nDepProcs < 1
          ! cell not existing. stop.
          write(logUnit(6),*) '---------- WARNING!!! -------------------------'
          write(logUnit(6),*) 'cell', TreeID ,' not existing on proc ',       &
            &        proc%rank, '. Number of dependencies', nDepProcs
          call tem_tIDinfo( me = TreeID, tree = tree, nUnit = logUnit(6) )
          write(logUnit(6),*) 'WARNING!!! This should never occur and '        &
            &                 // 'points to a bug or a buggy mesh.'
        end if

        if (elemPos > 0) then
          ! Set the encountered element position to the hash
          hash(hashpos) = TreeID
          hash_elemPos( hashpos ) = elemPos
        end if
        ! end cache miss
        ! -----------------------------------------------
      else cachemiss
        ! -----------------------------------------------
        ! Cache hit, i.e. element already in levelDesc%elem, just update nesting
        elemPos = hash_elemPos( hashpos )
        levelDesc( neighLevel )%elem%haloNesting%val( elemPos )            &
          &  = min( levelDesc( neighLevel )%elem%haloNesting%val(elemPos), &
          &         nesting                                                )
        ! -----------------------------------------------
      end if cachemiss
    else ! treeID <= 0, i.e. it is a bcID
      elemPos = int( TreeID )
    end if ! TreeID > 0

  end subroutine identify_elements