build_levelElements Subroutine

private subroutine build_levelElements(levelDesc, tree, proc, stencil, pathFirst, pathLast)

Assemble the fluid list and identify neighbor relations start with building up the ghost and halo element collection as well

Arguments

Type IntentOptional Attributes Name
type(tem_levelDesc_type), intent(inout) :: levelDesc(tree%global%minLevel:)

the level descriptor to be filled

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

the global tree

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

Process description to use.

type(tem_stencilHeader_type) :: stencil

array of all stencils used in the simulation

type(tem_path_type), intent(in) :: pathFirst(:)

first and last treeID path in every process

type(tem_path_type), intent(in) :: pathLast(:)

first and last treeID path in every process


Calls

proc~~build_levelelements~~CallsGraph proc~build_levelelements build_levelElements interface~tem_logging_isactive tem_logging_isActive proc~build_levelelements->interface~tem_logging_isactive proc~identify_additionalneigh identify_additionalNeigh proc~build_levelelements->proc~identify_additionalneigh proc~identify_elements identify_elements proc~build_levelelements->proc~identify_elements proc~tem_elemlist_dump tem_elemList_dump proc~build_levelelements->proc~tem_elemlist_dump proc~tem_horizontalspacer tem_horizontalSpacer proc~build_levelelements->proc~tem_horizontalspacer proc~tem_logging_isactive_for tem_logging_isActive_for interface~tem_logging_isactive->proc~tem_logging_isactive_for proc~tem_logging_isactive_primary tem_logging_isActive_primary interface~tem_logging_isactive->proc~tem_logging_isactive_primary proc~identify_additionalneigh->proc~identify_elements proc~identify_additionalneigh->proc~tem_horizontalspacer interface~positionofval~5 positionofval proc~identify_additionalneigh->interface~positionofval~5 proc~identify_elements->proc~identify_elements interface~append~11 append proc~identify_elements->interface~append~11 interface~init~24 init proc~identify_elements->interface~init~24 proc~create_allparentneighbors create_allParentNeighbors proc~identify_elements->proc~create_allparentneighbors proc~identify_stencilneigh identify_stencilNeigh proc~identify_elements->proc~identify_stencilneigh proc~single_process_element single_process_element proc~identify_elements->proc~single_process_element proc~tem_directchildren tem_directChildren proc~identify_elements->proc~tem_directchildren proc~tem_find_depproc tem_find_depProc proc~identify_elements->proc~tem_find_depproc proc~tem_levelof tem_LevelOf proc~identify_elements->proc~tem_levelof proc~tem_pathof tem_PathOf proc~identify_elements->proc~tem_pathof proc~tem_tidinfo tem_tIDinfo proc~identify_elements->proc~tem_tidinfo proc~tem_element_dump tem_element_dump proc~tem_elemlist_dump->proc~tem_element_dump

Called by

proc~~build_levelelements~~CalledByGraph proc~build_levelelements build_levelElements proc~tem_find_allelements tem_find_allElements proc~tem_find_allelements->proc~build_levelelements proc~tem_create_leveldesc tem_create_levelDesc proc~tem_create_leveldesc->proc~tem_find_allelements 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 build_levelElements( levelDesc, tree, proc, stencil,       &
    &                             pathFirst, pathLast )
    ! -------------------------------------------------------------------- !
    !> the global tree
    type(treelmesh_type),        intent(in) :: tree
    !> the level descriptor to be filled
    type(tem_levelDesc_type), intent(inout) :: levelDesc( tree%global%minLevel:)
    !> Process description to use.
    type(tem_comm_env_type), intent(in)     :: proc
    !> array of all stencils used in the simulation
    type(tem_stencilHeader_type)            :: stencil
    !> first and last treeID path in every process
    type(tem_path_type), intent(in)         :: pathFirst(:), pathLast(:)
    ! -------------------------------------------------------------------- !
    integer(kind=long_k) :: neighID   ! neighboring neighID
    integer :: elemPos, minLevel, maxLevel
    integer :: iElem, iNeighElem, iLevel, iStencil
    ! position of where to read the stencil neighbor neighID in the element
    integer :: neighPos
    ! -------------------------------------------------------------------- !

    minLevel = tree%global%minLevel
    maxLevel = tree%global%maxLevel

    write(logUnit(3),*) 'Building level-wise fluid list ...'

    call tem_horizontalSpacer( fUnit = dbgUnit(1), before = 1 )
    write(dbgUnit(3),*) 'Inside routine: build_levelElements'

    ! step 1. Identify the neighbors of each elements including ghost from finer
    ! and ghost from coarser
    ! Now iterate over all the neighbor lists and add the required elements if
    ! necessary
    do iLevel = minLevel, maxLevel
      write(dbgUnit(7),*) 'Level: ', iLevel
      do iElem = 1, levelDesc( iLevel )%elem%nElems( eT_fluid )
        do iStencil = 1, levelDesc( iLevel )%elem%stencil%val( iElem )%nVals
          do iNeighElem = 1, levelDesc( iLevel )%elem                          &
            &                         %stencil%val( iElem )%val( iStencil )%QQN
            ! Same as in identify_additionalNeigh
            neighPos = levelDesc( iLevel )%elem%stencil%val( iElem )           &
              &                           %val( iStencil )%tIDpos( iNeighElem )
            neighID = levelDesc(iLevel)%elem%neighID%val(iElem)%val( neighPos )

! write(dbgUnit(7),"(6(A,I0))") 'iElem: ', iElem, &
!   &                 ', treeID: ', levelDesc( iLevel )%elem%tID%val( iElem ), &
!   &                 ', iStencil: ', iStencil, &
!   &                 ', iNeighElem: ', iNeighElem, &
!   &                 ', neighPos: ', neighPos, &
!   &                 ', neighID: ', neighID

            if( neighID > 0_long_k ) then
              ! identify the process in which requested neighID is located
              !   - if remote process add to elem list as halo element.
              !   - if different level add to elem list as ghostFromFiner or
              !     ghostFromCoarser.
              ! Also return position (elemPos) of neighID in elem%tID list
              call identify_elements( treeID    = neighID,      &
                &                     tree      = tree,         &
                &                     pathFirst = pathFirst,    &
                &                     pathLast  = pathLast,     &
                &                     levelDesc = levelDesc,    &
                &                     elemPos   = elemPos,      &
                &                     proc      = proc,         &
                &                     nesting   = 0,            &
                &                     stencil   = stencil )
            else ! neighID < 0, i.e. it is a bcID
              elemPos = int( neighID )
            end if ! neighPos > 0
            levelDesc( iLevel )%elem%stencil%val( iElem )%val( iStencil )      &
              &                               %totalPos( iNeighElem ) = elemPos
          end do ! iNeighElem
        end do ! iStencil
      end do ! iElem
    end do ! iLevel

    ! JQ: halo and ghost (with empty stencil) are added to elem list.
    !     shall we calculate possible neighbors for them?

    ! Find neighbors of neighbors for elements in the require list.
    ! Same procedure as above but identify neighbors of elements in require
    ! list along 1st stencil directions.
    ! What exactly is the require list for?
    ! - Used ONLY for boundary stencil with higher order neighbors i.e
    !   only when require nVals > 0
    call identify_additionalNeigh( tree       = tree,             &
      &                            proc       = proc,             &
      &                            levelDesc  = levelDesc,        &
      &                            pathFirst  = pathFirst,        &
      &                            pathLast   = pathLast,         &
      &                            stencil    = stencil )

    ! dump elemList into debug unit
    if (tem_logging_isActive(main_debug%logger, 5)) then
      do iLevel = minLevel, maxLevel
          call tem_elemList_dump( me      = levelDesc( iLevel )%elem, &
            &                     nUnit   = dbgUnit(5),               &
            &                     stencil = .true.,                   &
            &                     string  = 'after build level elm'   )
      end do ! iLevel
    end if

    write(dbgUnit(1),*) 'Leave routine: build_levelElements'
    call tem_horizontalSpacer( fUnit = dbgUnit(1), after = 1 )

  end subroutine build_levelElements