single_process_element Subroutine

private subroutine single_process_element(targetID, levelDesc, tree, proc, iProc, minLevel, elemPos, stencil, nesting, updated, skip_add_additionalGhost)

Determine the location (which process) of a requested element, which was identified to be located on one single process (can be local or remote) If it is located on a remote process: add to halo list local process: identify if ghost or fluid

Arguments

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

neighboring treeID

type(tem_levelDesc_type), intent(inout) :: levelDesc(minLevel:)

the level descriptor to be filled

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

tree information

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

Process description to use.

integer, intent(in) :: iProc

Process on which targetID is located

integer, intent(in) :: minLevel

minimum level fluid element in the tree

integer, intent(out) :: elemPos

targetID element position in the levelDesc % elem list

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?

logical, intent(in), optional :: skip_add_additionalGhost

logical, optional, if true no ghosts are added


Calls

proc~~single_process_element~~CallsGraph proc~single_process_element single_process_element interface~append~29 append proc~single_process_element->interface~append~29 interface~init~20 init proc~single_process_element->interface~init~20 proc~identify_local_element identify_local_element proc~single_process_element->proc~identify_local_element proc~tem_abort tem_abort proc~single_process_element->proc~tem_abort proc~tem_levelof tem_LevelOf proc~single_process_element->proc~tem_levelof proc~append_ga_dynlong append_ga_dynlong interface~append~29->proc~append_ga_dynlong proc~append_ga_dynlong_vec append_ga_dynlong_vec interface~append~29->proc~append_ga_dynlong_vec proc~init_ga2d_real init_ga2d_real interface~init~20->proc~init_ga2d_real proc~identify_local_element->interface~append~29 proc~identify_local_element->proc~tem_levelof interface~positionofval~5 positionofval proc~identify_local_element->interface~positionofval~5 proc~add_all_virtual_children add_all_virtual_children proc~identify_local_element->proc~add_all_virtual_children proc~add_ghostfromfiner add_ghostFromFiner proc~identify_local_element->proc~add_ghostfromfiner proc~tem_pathof tem_PathOf proc~identify_local_element->proc~tem_pathof proc~tem_posofpath tem_PosOfPath proc~identify_local_element->proc~tem_posofpath proc~tem_tidinfo tem_tIDinfo proc~identify_local_element->proc~tem_tidinfo mpi_abort mpi_abort proc~tem_abort->mpi_abort proc~posofval_label posofval_label interface~positionofval~5->proc~posofval_label proc~add_all_virtual_children->interface~append~29 proc~add_all_virtual_children->interface~init~20 proc~add_all_virtual_children->proc~tem_levelof proc~add_all_virtual_children->interface~positionofval~5 proc~add_all_virtual_children->proc~add_all_virtual_children proc~tem_childnumber tem_childNumber proc~add_all_virtual_children->proc~tem_childnumber proc~tem_coordofid tem_CoordOfId proc~add_all_virtual_children->proc~tem_coordofid proc~tem_find_bcs_fromcoarser tem_find_BCs_fromCoarser proc~add_all_virtual_children->proc~tem_find_bcs_fromcoarser proc~tem_idofcoord tem_IdOfCoord proc~add_all_virtual_children->proc~tem_idofcoord proc~add_ghostfromfiner->interface~append~29 proc~add_ghostfromfiner->proc~tem_levelof proc~add_ghostfromfiner->interface~positionofval~5 proc~add_ghostfromfiner->proc~add_ghostfromfiner proc~add_ghostfromfiner->proc~tem_pathof proc~add_ghostfromfiner->proc~tem_posofpath proc~tem_directchildren tem_directChildren proc~add_ghostfromfiner->proc~tem_directchildren proc~tem_find_bcs_fromfiner tem_find_BCs_fromFiner proc~add_ghostfromfiner->proc~tem_find_bcs_fromfiner interface~expand~25 expand proc~append_ga_dynlong->interface~expand~25 proc~append_ga_dynlong_vec->interface~expand~25 proc~tem_baryofid tem_BaryOfId proc~tem_tidinfo->proc~tem_baryofid proc~tem_tidinfo->proc~tem_coordofid proc~tem_elemsize tem_ElemSize proc~tem_tidinfo->proc~tem_elemsize proc~expand_ga_dynlong expand_ga_dynlong interface~expand~25->proc~expand_ga_dynlong interface~sortedposofval~5 sortedposofval proc~posofval_label->interface~sortedposofval~5 proc~tem_baryofid->proc~tem_coordofid proc~tem_elemsizelevel tem_ElemSizeLevel proc~tem_baryofid->proc~tem_elemsizelevel proc~tem_coordofid->proc~tem_levelof proc~tem_elemsize->proc~tem_levelof proc~tem_elemsize->proc~tem_elemsizelevel proc~childtostencil childToStencil proc~tem_find_bcs_fromcoarser->proc~childtostencil proc~tem_find_bcs_fromfiner->interface~append~29 proc~tem_find_bcs_fromfiner->interface~init~20 proc~stenciltochild stencilToChild proc~tem_find_bcs_fromfiner->proc~stenciltochild proc~update_childneighborid update_childNeighborID proc~tem_find_bcs_fromfiner->proc~update_childneighborid

Called by

proc~~single_process_element~~CalledByGraph proc~single_process_element single_process_element proc~identify_elements identify_elements proc~identify_elements->proc~single_process_element proc~identify_elements->proc~identify_elements proc~create_allparentneighbors create_allParentNeighbors proc~identify_elements->proc~create_allparentneighbors proc~identify_stencilneigh identify_stencilNeigh proc~identify_elements->proc~identify_stencilneigh proc~build_levelelements build_levelElements proc~build_levelelements->proc~identify_elements proc~identify_additionalneigh identify_additionalNeigh proc~build_levelelements->proc~identify_additionalneigh proc~create_allparentneighbors->proc~identify_elements proc~create_allparentneighbors->proc~identify_stencilneigh proc~identify_additionalneigh->proc~identify_elements proc~identify_stencilneigh->proc~identify_elements proc~request_remotehalos request_remoteHalos proc~request_remotehalos->proc~create_allparentneighbors proc~request_remotehalos->proc~identify_stencilneigh proc~tem_find_allelements tem_find_allElements proc~tem_find_allelements->proc~build_levelelements proc~tem_find_allelements->proc~identify_additionalneigh proc~communicate_elements communicate_elements proc~communicate_elements->proc~request_remotehalos proc~tem_create_leveldesc tem_create_levelDesc proc~tem_create_leveldesc->proc~tem_find_allelements

Source Code

  subroutine single_process_element( targetID, levelDesc, tree, proc, iProc, &
    &                                minLevel, elemPos, stencil, nesting,    &
    &                                updated, skip_add_additionalGhost )
    ! -------------------------------------------------------------------- !
    !> neighboring treeID
    integer(kind=long_k), intent(in)         :: targetID
    !> minimum level fluid element in the tree
    integer, intent(in)                      :: minLevel
    !> the level descriptor to be filled
    type(tem_levelDesc_type), intent(inout)  :: levelDesc(minLevel:)
    !> Process description to use.
    type(tem_comm_env_type), intent(in)      :: proc
    !> Process on which targetID is located
    integer, intent(in)                      :: iProc
    !> tree information
    type(treelmesh_type), intent(in)         :: tree
    !> current stencil definition
    type(tem_stencilHeader_type), intent(in) :: stencil
    !> targetID element position in the levelDesc % elem list
    integer, intent(out)                     :: elemPos
    !> nesting level
    integer, intent(in)                      :: nesting
    !> was the element updated in this call?
    logical, intent(out)                     :: updated
    !> logical, optional, if true no ghosts are added
    logical, intent(in), optional :: skip_add_additionalGhost
    ! -------------------------------------------------------------------- !
    type(tem_stencilElement_type) :: emptyStencil(1)
    integer :: targetLevel
    logical :: wasAdded
    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 

    targetLevel = tem_LevelOf(targetID) ! Has to be same as tLevel!?
    if ( (targetLevel < minLevel)                 &
      &  .or. (targetLevel > uBound(levelDesc,1)) ) then
      write(logUnit(1),*) ' ERROR: level which is not included in the fluid'// &
        &                 ' tree was demanded in singleProcNeigh'
      write(logUnit(1),"(2(A,I0))") 'treeID: ', targetID, ', level: ', &
        &                           targetLevel
      call tem_abort()
    end if

    ! Set the element updated flag as a default to false
    updated = .false.

    ! If it is a remote cell on only one process -> regular halo
    if (iProc /= proc%rank + 1) then
      ! REMOTE
      call init( me = emptyStencil(1), QQN=stencil%QQN )

      ! append this targetID as halo element to levelDesc elem list
      call append( me              = levelDesc(targetLevel)%elem, &
        &          tID             = targetID,                    &
        &          eType           = eT_halo,                     &
        &          property        = 0_long_k,                    &
        &          sourceProc      = iProc,                       &
        &          haloNesting     = nesting,                     &
        &          stencilElements = emptyStencil,                &
        &          pos             = elemPos,                     &
        &          wasAdded        = wasAdded                     )

      if (.not. wasAdded) then
        ! If this element was already there, make sure we use the minimal
        ! nesting level requested for this element.
        updated = ( nesting < levelDesc(targetLevel) &
          &                     %elem                &
          &                     %haloNesting         &
          &                     %val(elemPos)        )
        ! If the nesting has been updated (decreased, we need to revisit this
        ! element in search for its neighbors).
        levelDesc(targetLevel)%elem%needsUpdate%val(elemPos) = updated

        levelDesc(targetLevel)%elem%haloNesting%val(elemPos)              &
          &  = min( levelDesc(targetLevel)%elem%haloNesting%val(elemPos), &
          &         nesting                                               )
      else
        ! New halo element added
        updated = .true.
        write(dbgUnit(7),"(A,I0)") 'Added as a Halo: ', targetID, &
          &                        'to level: ', targetLevel
      end if ! wasAdded

    else ! iProc == proc%rank + 1
      ! LOCAL

      ! Either a local ghost or fluid cell
      call identify_local_element(                                 &
        &    targetID                 = targetID,                  &
        &    levelDesc                = levelDesc,                 &
        &    tree                     = tree,                      &
        &    elemPos                  = elemPos,                   &
        &    minLevel                 = minLevel,                  &
        &    nesting                  = nesting,                   &
        &    updated                  = updated,                   &
        &    stencil                  = stencil,                   &
        &    skip_add_additionalGhost = l_skip_add_additionalGhost )
    end if ! iProc /= proc%rank + 1

  end subroutine single_process_element