tem_find_BCs_fromFiner Subroutine

private subroutine tem_find_BCs_fromFiner(childPos, sourceLevel, targetLevel, targetPos, levelDesc, minLevel, stencil)

Inherit the neighborhood from the sourceELem to the targetElem Note that targetElem is inout, as it might have already values assigned.

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: childPos(8)

position of all childs in the levelDesc elem tID list

integer, intent(in) :: sourceLevel

level of child

integer, intent(in) :: targetLevel

level of parent

integer, intent(in) :: targetPos

added position of parent in the levelDesc elem tID list

type(tem_levelDesc_type) :: levelDesc(minLevel:)

the level descriptor to be filled

integer, intent(in) :: minLevel

minimum level in the tree

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

current stencil definition


Calls

proc~~tem_find_bcs_fromfiner~~CallsGraph proc~tem_find_bcs_fromfiner tem_find_BCs_fromFiner interface~append~11 append proc~tem_find_bcs_fromfiner->interface~append~11 interface~init~24 init proc~tem_find_bcs_fromfiner->interface~init~24 proc~stenciltochild stencilToChild proc~tem_find_bcs_fromfiner->proc~stenciltochild proc~update_childneighborid update_childNeighborID proc~tem_find_bcs_fromfiner->proc~update_childneighborid proc~append_ga_dynlong append_ga_dynlong interface~append~11->proc~append_ga_dynlong proc~append_ga_dynlong_vec append_ga_dynlong_vec interface~append~11->proc~append_ga_dynlong_vec proc~init_ga2d_real init_ga2d_real interface~init~24->proc~init_ga2d_real interface~tem_parentof tem_ParentOf proc~update_childneighborid->interface~tem_parentof proc~tem_idofcoord tem_IdOfCoord proc~update_childneighborid->proc~tem_idofcoord proc~tem_directparent tem_directParent interface~tem_parentof->proc~tem_directparent proc~tem_parentatlevel tem_ParentAtLevel interface~tem_parentof->proc~tem_parentatlevel interface~expand~9 expand proc~append_ga_dynlong->interface~expand~9 proc~append_ga_dynlong_vec->interface~expand~9 proc~expand_ga_dynlong expand_ga_dynlong interface~expand~9->proc~expand_ga_dynlong proc~tem_levelof tem_LevelOf proc~tem_parentatlevel->proc~tem_levelof

Called by

proc~~tem_find_bcs_fromfiner~~CalledByGraph proc~tem_find_bcs_fromfiner tem_find_BCs_fromFiner proc~add_ghostfromfiner add_ghostFromFiner proc~add_ghostfromfiner->proc~tem_find_bcs_fromfiner proc~add_ghostfromfiner->proc~add_ghostfromfiner proc~identify_local_element identify_local_element proc~identify_local_element->proc~add_ghostfromfiner 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_fromFiner( childPos, sourceLevel, targetLevel, &
    &                                targetPos, levelDesc, minLevel,     &
    &                                stencil                             )
    ! -------------------------------------------------------------------- !
    !> position of all childs in the levelDesc elem tID list
    integer, intent(in) :: childPos(8)
    !> level of child
    integer, intent(in) :: sourceLevel
    !> level of parent
    integer, intent(in) :: targetLevel
    !> added position of parent in the levelDesc elem tID list
    integer, intent(in) :: targetPos
    !> minimum level in the tree
    integer, intent(in) :: minLevel
    !> the level descriptor to be filled
    type(tem_levelDesc_type ) :: levelDesc(minLevel:)
    !> current stencil definition
    type(tem_stencilHeader_type), intent(in) :: stencil
    ! -------------------------------------------------------------------- !
    integer :: dir
    ! Tangential direction iterators
    integer :: iDirX, iDirY, iDir
    integer :: dirX !< first tangential direction
    integer :: dirY !< second tangential direction
    integer :: iStencilElem  !< stencil element iterator
    integer :: childCoord(4)
    ! the type of the current stencil link (side, edge, corner)
    integer :: linktype
    integer :: iDirNormal
    integer :: iStencil, addedPos
    type(tem_stencilElement_type) :: tStencil
    integer(kind=long_k) :: tNeighID
    ! -------------------------------------------------------------------- !

    if ( .not. allocated( levelDesc( targetLevel )%elem%stencil%              &
      &                                              val(targetPos)%val )) then
      call init( me     = levelDesc( targetLevel )%elem%stencil%val(targetPos),&
        &        length = 1 )
    end if

    call init( me = tStencil, QQN = stencil%QQN, headerPos = 1 )
    call append( me  = levelDesc( targetLevel )%elem%stencil%val(targetPos),   &
      &          val = tStencil )
    if( levelDesc( targetLevel )%elem%neighID%val(targetPos)%nVals > 0 ) then
      write(dbgUnit(2),*) 'Warning: tem_find_BCs_fromFiner has already' &
        &                       //' encountered'
      write(dbgUnit(2),*) '         existing neighbor IDs for the stencil'
    end if
    iStencil = levelDesc( targetLevel )%elem%stencil%val( targetPos )%nVals

    childCoord(4) = 1
    do iStencilElem = 1, stencil%QQN
      ! Reset current direction
      ! only 0 can not be used for inheritance of BCs later because boundaries
      ! are stored as negative treeIDs. Transformation of 0 into a HUGE
      ! negative value enables easy inheritance.
      tNeighID = -huge( 0_long_k )
      ! Count the non-zero entries to decide the kind of the current
      ! link direction in the stencil. (only direct neighbors)
      ! Distinguishing three different kinds of neighbors:
      ! SIDE, EDGE and CORNER
      linktype = sum( abs(stencil%cxDir( :, iStencilElem )))
      select case(linktype)
      case(1)
        ! 1 non-zero
        ! SIDE, needs four children
        ! identify the normal direction of the side
        dir = maxloc( abs( stencil%cxDir(:, iStencilElem )), 1)
        ! Convert the stencil direction to child direction
        iDirNormal = stencilToChild( stencil%cxDir(dir, iStencilElem) )
        childCoord( dir ) = iDirNormal
        ! identify the tangential directions
        dirX = mod( dir,   3) + 1
        dirY = mod( dir+1, 3) + 1
        do iDirY = 0, 1
          childCoord( dirY ) = iDirY
          do iDirX = 0, 1
            childCoord( dirX ) = iDirX
            call update_childNeighborID( neighID      = tNeighID,            &
              &                          childCoord   = childCoord,          &
              &                          childPos     = childPos,            &
              &                          iStencil     = iStencil,            &
              &                          iStencilElem = iStencilElem,        &
              &                          elem  = levelDesc(sourceLevel)%elem )
          end do
        end do

      case(2)
        ! 2 non-zeroes
        ! EDGE, needs two children
        ! Determine the 0-direction (dimension with the zero entry)
        dir = minloc( abs( stencil%cxDir(:,iStencilElem )), 1)
        ! Get the edge directions
        dirX = mod( dir,   3) + 1
        dirY = mod( dir+1, 3) + 1
        childCoord( dirX ) = stencilToChild(stencil%cxDir(dirX,iStencilElem))
        childCoord( dirY ) = stencilToChild(stencil%cxDir(dirY,iStencilElem))
        do iDir = 0, 1
          childCoord( dir ) = iDir
          call update_childNeighborID( neighID      = tNeighID,            &
            &                          childCoord   = childCoord,          &
            &                          childPos     = childPos,            &
            &                          iStencil     = iStencil,            &
            &                          iStencilElem = iStencilElem,        &
            &                          elem  = levelDesc(sourceLevel)%elem )
        enddo

      case(3)
        ! No zero at all, all three directions have a offset
        ! CORNER, just a single child is connected to this link
        childCoord( 1:3 ) = stencilToChild( stencil%cxDir(:,iStencilElem) )
        call update_childNeighborID( neighID      = tNeighID,            &
          &                          childCoord   = childCoord,          &
          &                          childPos     = childPos,            &
          &                          iStencil     = iStencil,            &
          &                          iStencilElem = iStencilElem,        &
          &                          elem  = levelDesc(sourceLevel)%elem )
      end select

      ! Append the neighID of virtual parent
      call append( me  = levelDesc(targetLevel)%elem%neighID%val(targetPos), &
        &          val = tNeighID,                                           &
        &          pos = addedPos                                            )
      levelDesc(targetLevel)     &
        &  %elem                 &
        &  %stencil              &
        &  %val(targetPos)       &
        &  %val(iStencil)        &
        &  %tIDpos(iStencilElem) = addedPos

      ! to append the neighbors of ghosts to required list
      !call append( me       = levelDesc( targetLevel )%require, &
      !  &          val      = tNeighID,                         &
      !  &          pos      = addedPos )

    end do

  end subroutine tem_find_BCs_fromFiner