tem_find_BCs_fromCoarser Subroutine

private subroutine tem_find_BCs_fromCoarser(dir, childCoord, sourceLevel, sourcePos, neighID, computeStencil, levelDesc, minLevel)

Inherit the neighborhood from the sourceELem to the targetElem

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: dir
integer, intent(in) :: childCoord(4)

coordinate of virtual(ghost) child

integer, intent(in) :: sourceLevel

coarse element level

integer, intent(in) :: sourcePos

position of coarser element in original treeID list

integer(kind=long_k), intent(inout) :: neighID(:)

neighbor treeIDs of child element

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

current stencil definition

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

the level descriptor to be filled

integer, intent(in) :: minLevel

minimum level in the tree


Calls

proc~~tem_find_bcs_fromcoarser~~CallsGraph proc~tem_find_bcs_fromcoarser tem_find_BCs_fromCoarser proc~childtostencil childToStencil proc~tem_find_bcs_fromcoarser->proc~childtostencil

Called by

proc~~tem_find_bcs_fromcoarser~~CalledByGraph proc~tem_find_bcs_fromcoarser tem_find_BCs_fromCoarser proc~add_all_virtual_children add_all_virtual_children proc~add_all_virtual_children->proc~tem_find_bcs_fromcoarser proc~add_all_virtual_children->proc~add_all_virtual_children proc~identify_local_element identify_local_element proc~identify_local_element->proc~add_all_virtual_children proc~identify_halo identify_halo proc~identify_halo->proc~identify_local_element proc~single_process_element single_process_element proc~single_process_element->proc~identify_local_element proc~identify_elements identify_elements proc~identify_elements->proc~single_process_element proc~request_remotehalos request_remoteHalos proc~request_remotehalos->proc~identify_halo

Source Code

  subroutine tem_find_BCs_fromCoarser( dir, childCoord, sourceLevel,       &
    &                                  sourcePos, neighID, computeStencil, &
    &                                  levelDesc, minLevel                 )
    ! -------------------------------------------------------------------- !
    !> coarse element level
    integer, intent(in) :: sourceLevel
    !> position of coarser element in original treeID list
    integer, intent(in) :: sourcePos
    !>
    integer, intent(in) :: dir
    !> minimum level in the tree
    integer, intent(in) :: minLevel
    !> neighbor treeIDs of child element
    integer(kind=long_k), intent(inout) :: neighID(:)
    !> coordinate of virtual(ghost) child
    integer, intent(in) :: childCoord(4)
    !> current stencil definition
    type(tem_stencilHeader_type), intent(in) :: computeStencil
    !> the level descriptor to be filled
    type(tem_levelDesc_type), intent(in) :: levelDesc(minLevel:)
    ! -------------------------------------------------------------------- !
    ! Tangential direction iterators
    integer :: iDirX, iDirY
    integer :: dirX ! first tangential direction
    integer :: dirY ! second tangential direction
    integer :: curDir ! child face normal direction, always {-1,1}
    integer :: myLink(4) ! child link direction
    integer :: parentLink(3) ! each entry must be either {-1,0,1}
    integer :: iChildStencil ! child stencil element iterator
    integer :: iStencilElem  ! stencil element iterator
    integer :: iStencil
    integer :: posInNeighID
    ! -------------------------------------------------------------------- !

    ! curDir is -1 if childCoord(dir) is 0
    ! curDir is 1 if childCoord(dir) is 1
    ! it is needed to find in which direction of child's neighIDs are to be
    ! found
    curDir = childToStencil( childCoord( dir ))
    myLink( dir ) = curDir
    myLink(4) = 0
    parentLink(dir) = curDir

    ! identify the tangential directions
    dirX = mod( dir, 3) + 1
    dirY = mod( dir+1, 3) + 1
    iStencil = 1

    ! Iterate over all nine child links of the surface with the normal pointing
    ! into dir
    do iDirY = -1, 1
      myLink( dirY ) = iDirY
      ! select the valid directions to look for the neighbors in the second
      ! tangential direction
      if( childCoord( dirY ) == 0) then
        parentLink( dirY ) = min( iDirY, 0 )
      else
        parentLink( dirY ) = max( iDirY, 0 )
      end if
      do iDirX = -1, 1
        myLink( dirX ) = iDirX
        if( childCoord( dirX ) == 0) then
          parentLink( dirX ) = min( iDirX, 0 )
        else
          parentLink( dirX ) = max( iDirX, 0 )
        end if
        ! matching direction of child to compute stencil
        iChildStencil = 0
        do iStencilElem = 1, computeStencil%QQN
          if (      computeStencil%cxDir(1, iStencilElem) == myLink(1)  &
            & .and. computeStencil%cxDir(2, iStencilElem) == myLink(2)  &
            & .and. computeStencil%cxDir(3, iStencilElem) == myLink(3)  &
            &                                                           ) then
            ! Found matching stencil entry for current child direction
            iChildStencil = iStencilElem
          end if
        end do

        if ( iChildStencil > 0 ) then
          ! matching direction of parent to compute stencil
          do iStencilElem = 1, computeStencil%QQN
            if (      computeStencil%cxDir(1, iStencilElem) == parentLink(1) &
              & .and. computeStencil%cxDir(2, iStencilElem) == parentLink(2) &
              & .and. computeStencil%cxDir(3, iStencilElem) == parentLink(3) &
              &                                                              ) then
              ! Found matching stencil entry for current child direction
              ! Set the parent's neighbor here. Later we replace the
              ! parent's neighbor by the current level's neighbor
              posInNeighID = levelDesc(sourceLevel) &
                &              %elem                &
                &              %stencil             &
                &              %val(sourcePos)      &
                &              %val(iStencil)       &
                &              %tIDpos(iStencilElem)
              neighID( iChildStencil )        &
                &  = levelDesc( sourceLevel ) &
                &      %elem                  &
                &      %neighID               &
                &      %val(sourcePos)        &
                &      %val(posInNeighID)
            end if
          end do ! iStencilElem
        end if ! iChildStencil > 0
      end do ! iDirX
    end do ! iDirY

  end subroutine tem_find_BCs_fromCoarser