tem_build_listHorizontalDep Subroutine

private subroutine tem_build_listHorizontalDep(iStencil, levelDesc, posInSortElem, nElems, iIndex)

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).

integer, intent(in) :: posInSortElem(:)

Positions in sorted elem%tID list

integer, intent(in) :: nElems

number of elements

integer, intent(inout) :: iIndex

Calls

proc~~tem_build_listhorizontaldep~~CallsGraph proc~tem_build_listhorizontaldep tem_build_listHorizontalDep proc~tem_treeidintotal tem_treeIDinTotal proc~tem_build_listhorizontaldep->proc~tem_treeidintotal proc~tem_etypeofid tem_eTypeOfId proc~tem_treeidintotal->proc~tem_etypeofid tem_positioninsorted tem_positioninsorted proc~tem_treeidintotal->tem_positioninsorted interface~positionofval~5 positionofval proc~tem_etypeofid->interface~positionofval~5 proc~posofval_label posofval_label interface~positionofval~5->proc~posofval_label

Called by

proc~~tem_build_listhorizontaldep~~CalledByGraph proc~tem_build_listhorizontaldep tem_build_listHorizontalDep proc~tem_build_horizontaldependencies tem_build_horizontalDependencies proc~tem_build_horizontaldependencies->proc~tem_build_listhorizontaldep 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

Source Code

  subroutine tem_build_listHorizontalDep( iStencil, levelDesc, &
    &                                     posInSortElem, nElems, iIndex )
    ! ---------------------------------------------------------------------------
    !> Level descriptor for each level of your mesh (starting from min level).
    type(tem_levelDesc_type), intent(inout) :: levelDesc
    !> Index of your neighbor list.
    integer, intent(in) :: iStencil
    !>
    integer, intent(inout) :: iIndex
    !> number of elements
    integer, intent(in) :: nElems
    !> Positions in sorted elem%tID list
    integer, intent(in) :: posInSortElem(:)
    ! ---------------------------------------------------------------------------
    integer :: iElem, iStencilElem
    integer :: posInElem, levelPos, neighPos, neighVal
    integer(kind=long_k) :: nTreeID ! neighbor treeID value to treat
    logical :: invalid
    ! ---------------------------------------------------------------------------
    ! run over the given element list to update
    do iElem = 1, nElems
      posInElem = posInSortElem( iElem )
!write(*,*)'in build_listHorDep nVals: ', levelDesc%elem%stencil%val( posInElem )%nVals
!write(*,*)'in build_listHorDep iStencil: ', iStencil
      ! Check if there are enough stencils given in the element
      if( levelDesc%elem%stencil%val( posInElem )%nVals >= iStencil ) then
        ! Increase the index element counter, if a stencil was encountered
        iIndex = iIndex + 1
        ! Run over each entry of the stencil
        do iStencilElem = 1, levelDesc%elem%stencil%val( posInElem )%            &
          &                                                val( iStencil )%QQN
          neighPos = levelDesc%elem%stencil%val( posInElem )%val( iStencil )     &
            &                                           %tIDpos( iStencilElem )
          invalid = .false.
          if( neighPos > 0 ) then
            nTreeID = levelDesc%elem%neighID%val( posInElem )%val( neighPos )
            if( nTreeID < 0_long_k ) then
              ! If the entry in the stencil is < 0, it must be a boundary ID
              ! -> simply set
              neighVal = int(nTreeID)
            else
              ! If the entry is > 0, test if there is also an entry for the
              ! element position directly (%elem array)
              if( allocated( levelDesc%elem%stencil%val( posInElem )             &
                &                            %val( iStencil )%totalPos) ) then
                ! If there is an entry, test, if it is valid
                levelPos = levelDesc%elem%stencil%val( posInElem )               &
                  &                    %val( iStencil )%totalPos( iStencilElem )
                if( levelPos > 0 ) then
                  ! If there is a valid entry, get the position of this element
                  ! in the total list
                  neighVal = levelPos
                else
                  ! If the entry in %elem was not valid, use the tID to find
                  invalid = .true.
                end if ! levelPos > 0
              else
                ! if there was no %elem entry given at all, directly use the
                ! tID to find the element position in the total list
                invalid = .true.
              end if ! allocated totalPos
            end if ! nTreeID > 0

            ! If above anything went wrong, ...
            if( invalid ) then
              ! ... try to find it with the slowest method: just look up the
              ! tID in the total list
              neighVal = tem_treeIDinTotal( nTreeID, levelDesc )
            end if
          else ! neighPos < 0
            neighVal = neighPos
          end if ! neighPos > 0
          levelDesc%neigh(iStencil)%nghElems( iStencilElem, iIndex ) = neighVal
        end do ! stencil element
      end if ! nStencils > iStencil
    end do ! iElem

  end subroutine tem_build_listHorizontalDep