tem_build_treeHorizontalDep Subroutine

private subroutine tem_build_treeHorizontalDep(iStencil, levelDesc, computeStencil, list, nElems, tree)

Update the neighor arrays depending on what is given in the element stencil

The array levelDesc( iLevel )%neigh( iStenci )%nghElems( 1:QQN, 1:nElems ) is being filled up here

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: iStencil

Index of your neighbor list.

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

Level descriptor for each level of your mesh (starting from min level).

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

The stencil you build the horizontal dependencies for.

integer, intent(in) :: list(:)

stencil elemLvl points to sorted original treeID list

integer, intent(in) :: nElems

number of elements

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

tree information


Calls

proc~~tem_build_treehorizontaldep~~CallsGraph proc~tem_build_treehorizontaldep tem_build_treeHorizontalDep interface~positionofval~5 positionofval proc~tem_build_treehorizontaldep->interface~positionofval~5 proc~tem_treeidintotal tem_treeIDinTotal proc~tem_build_treehorizontaldep->proc~tem_treeidintotal proc~tem_baryofid tem_BaryOfId proc~tem_build_treehorizontaldep->proc~tem_baryofid proc~tem_stencil_getheaderpos tem_stencil_getHeaderPos proc~tem_build_treehorizontaldep->proc~tem_stencil_getheaderpos proc~posofval_label posofval_label interface~positionofval~5->proc~posofval_label proc~tem_etypeofid tem_eTypeOfId proc~tem_treeidintotal->proc~tem_etypeofid tem_positioninsorted tem_positioninsorted proc~tem_treeidintotal->tem_positioninsorted proc~tem_elemsizelevel tem_ElemSizeLevel proc~tem_baryofid->proc~tem_elemsizelevel proc~tem_coordofid tem_CoordOfId proc~tem_baryofid->proc~tem_coordofid interface~sortedposofval~5 sortedposofval proc~posofval_label->interface~sortedposofval~5 proc~tem_etypeofid->interface~positionofval~5 proc~tem_levelof tem_LevelOf proc~tem_coordofid->proc~tem_levelof proc~sortposofval_label sortposofval_label interface~sortedposofval~5->proc~sortposofval_label

Called by

proc~~tem_build_treehorizontaldep~~CalledByGraph proc~tem_build_treehorizontaldep tem_build_treeHorizontalDep proc~tem_build_horizontaldependencies tem_build_horizontalDependencies proc~tem_build_horizontaldependencies->proc~tem_build_treehorizontaldep proc~tem_create_leveldesc tem_create_levelDesc proc~tem_create_leveldesc->proc~tem_build_horizontaldependencies proc~tem_dimbydim_construction tem_dimByDim_construction proc~tem_dimbydim_construction->proc~tem_create_leveldesc proc~tem_build_face_info tem_build_face_info proc~tem_build_face_info->proc~tem_dimbydim_construction

Contents


Source Code

  subroutine tem_build_treeHorizontalDep( iStencil, levelDesc, computeStencil, &
    &                                     list, nElems, tree )
    ! ---------------------------------------------------------------------------
    !> Level descriptor for each level of your mesh (starting from min level).
    type(tem_levelDesc_type), intent(inout) :: levelDesc
    !> tree information
    type(treelmesh_type), intent(in) :: tree
    !> Index of your neighbor list.
    integer, intent(in) :: iStencil
    !> The stencil you build the horizontal dependencies for.
    type(tem_stencilHeader_type), intent(in) :: computeStencil
    !> number of elements
    integer, intent(in) :: nElems
    !> stencil elemLvl points to sorted original treeID list
    integer, intent(in) :: list(:)
    ! ---------------------------------------------------------------------------
    integer :: iElem, iStencilElem
    integer :: elemPos, stencilPos, levelPos, totalPos
    integer(kind=long_k) :: tID
    logical :: missingNeigh
    real(kind=rk) :: bary(3) ! bary center of element
    ! ---------------------------------------------------------------------------

    ! run over the given element list to update
    do iElem = 1, nElems
      missingNeigh = .false.
      tID = tree%treeID( list( iElem ))
      ! Get the position of the treeID in the element list
      elemPos = PositionOfVal( me = levelDesc%elem%tID, val = tID )
      if( elemPos > 0 ) then
        totalPos = tem_treeIDinTotal(                                          &
          &                    tID       = tID,                                &
          &                    levelDesc = levelDesc,                          &
          &                    eType     = levelDesc%elem%eType%val( elemPos ))
        stencilPos = tem_stencil_getHeaderPos(                                 &
          &                       me  = levelDesc%elem%stencil%val( elemPos ), &
          &                       val = iStencil )
        if( stencilPos > 0 ) then
          do iStencilElem = 1, levelDesc%elem%stencil%val(elemPos)             &
            &                                        %val(stencilPos)%QQN
            ! Totalpos has before already been updated, so just read from
            ! totalPos
            levelPos = levelDesc%elem%stencil%val(elemPos)                     &
              &                              %val(stencilPos)                  &
              &                              %totalPos(iStencilElem)
            if( levelPos <= 0 .and. computeStencil%requireNeighNeigh ) then
              missingNeigh = .true.
            end if
            levelDesc%neigh( iStencil )%nghElems( iStencilElem, iElem )        &
              &                                                     = levelPos
          end do
        end if
        if( missingNeigh ) then
          bary = tem_BaryOfId(tree, tID)
          write(logUnit(10),"(A,I0,A,I0,A,3ES12.5)") &
            &                   'missing neighbor iElem: ', iElem,&
            &                   ', tID: ', tID, ', bary: ', bary
          ! Set the property of the current element to solid
          levelDesc%property( totalPos )                                       &
            &  = ibset( levelDesc%property( totalPos ), prp_solid )
        end if
      else ! elemPos <= 0, i.e. tID can not be found in levelDesc%elem%tID
        write(dbgUnit(1),*) 'Warning!!!'
        write(dbgUnit(1),*) tID, ' not found during build_treehorizontalDep'
        write(dbgUnit(1),*) 'That should not occur.'
      end if
    end do

  end subroutine tem_build_treeHorizontalDep