tem_get_faces Subroutine

private subroutine tem_get_faces(levelDesc, direction, faces)

direction. It attaches also two-sided properties to the face.

Arguments

Type IntentOptional Attributes Name
type(tem_levelDesc_type), intent(in) :: levelDesc

Level descriptor of the level of the mesh you want to collect the faces for.

integer, intent(in) :: direction

The spatial direction to collect the faces for: 1 -> x direction 2 -> y direction 3 -> z direction

type(tem_face_descriptor_type), intent(inout) :: faces

Description of the faces on the current level.


Calls

proc~~tem_get_faces~~CallsGraph proc~tem_get_faces tem_get_faces proc~tem_init_facelist tem_init_faceList proc~tem_get_faces->proc~tem_init_facelist proc~tem_addface tem_addFace proc~tem_get_faces->proc~tem_addface proc~tem_get_elemprp tem_get_elemPrp proc~tem_get_faces->proc~tem_get_elemprp proc~tem_get_faceneigh tem_get_faceNeigh proc~tem_get_faces->proc~tem_get_faceneigh interface~init~22 init proc~tem_init_facelist->interface~init~22 interface~append~23 append proc~tem_addface->interface~append~23 proc~tem_abort tem_abort proc~tem_get_elemprp->proc~tem_abort proc~tem_treeidintotal tem_treeIDinTotal proc~tem_get_faceneigh->proc~tem_treeidintotal proc~tem_idofcoord tem_IdOfCoord proc~tem_get_faceneigh->proc~tem_idofcoord proc~tem_get_faceneigh->proc~tem_abort proc~tem_coordofid tem_CoordOfId proc~tem_get_faceneigh->proc~tem_coordofid proc~tem_etypeofid tem_eTypeOfId proc~tem_treeidintotal->proc~tem_etypeofid tem_positioninsorted tem_positioninsorted proc~tem_treeidintotal->tem_positioninsorted proc~init_ga2d_real init_ga2d_real interface~init~22->proc~init_ga2d_real proc~append_arrayga2d_real append_arrayga2d_real interface~append~23->proc~append_arrayga2d_real proc~append_singlega2d_real append_singlega2d_real interface~append~23->proc~append_singlega2d_real mpi_abort mpi_abort proc~tem_abort->mpi_abort proc~tem_levelof tem_LevelOf proc~tem_coordofid->proc~tem_levelof interface~positionofval~5 positionofval proc~tem_etypeofid->interface~positionofval~5 interface~expand~22 expand proc~append_arrayga2d_real->interface~expand~22 proc~append_singlega2d_real->interface~expand~22

Called by

proc~~tem_get_faces~~CalledByGraph proc~tem_get_faces tem_get_faces proc~tem_collect_faces tem_collect_faces proc~tem_collect_faces->proc~tem_get_faces proc~tem_build_face_info tem_build_face_info proc~tem_build_face_info->proc~tem_collect_faces

Contents

Source Code


Source Code

  subroutine tem_get_faces(levelDesc, direction, faces)
    ! --------------------------------------------------------------------------
    !> Level descriptor of the level of the mesh you want to collect the
    !! faces for.
    type(tem_levelDesc_type), intent(in) :: levelDesc
    !> The spatial direction to collect the faces for:
    !! 1 -> x direction
    !! 2 -> y direction
    !! 3 -> z direction
    integer, intent(in) :: direction
    !> Description of the faces on the current level.
    type(tem_face_descriptor_type), intent(inout) :: faces
    ! --------------------------------------------------------------------------
    integer :: iElem
    integer :: neighIndex
    integer :: neighPos
    integer :: elemPrp, neighPrp
    integer(kind=long_k) :: elemId, neighId
    ! --------------------------------------------------------------------------

    ! Initialize the facelist (with a dynamic array for the faceIDs)
    call tem_init_faceList( faces%faceList )

    ! Loop over all the elements (includes fluid, ghost and halos) have a look
    ! to the left and to the right face.
    ! Try to add both faces of the element to the list of all faces.
    elemLoop: do iElem = 1, levelDesc%nElems

      ! Get the current element property
      elemPrp = tem_get_elemPrp(levelDesc, iElem)
      elemId = levelDesc%total(iElem)

      ! Look for the left face of the current element.
      ! We use the compute stencil to determine the position of the neighboring
      ! left element.
      neighIndex = tem_left
      call tem_get_faceNeigh( levelDesc, iElem, direction, neighIndex, &
        &                     neighId, neighPos                        )
      neighPrp = tem_get_elemPrp(levelDesc, neighPos)

      ! For the left face, the current element will be on its right side.
      call tem_addFace( leftElemId   = neighId,  &
        &               leftElemPos  = neighPos, &
        &               leftElemPrp  = neighPrp, &
        &               rightElemId  = elemId,   &
        &               rightElemPos = iElem,    &
        &               rightElemPrp = elemPrp,  &
        &               faces        = faces     )

      ! Look for the right face of the current element.
      ! Get the right neighbor of the element.
      neighIndex = tem_right
      call tem_get_faceNeigh( levelDesc, iElem, direction, neighIndex, &
        &                     neighId, neighPos                        )
      neighPrp = tem_get_elemPrp(levelDesc, neighPos)

      ! For the right face, the current element will be on its left side.
      call tem_addFace( leftElemId   = elemId,   &
        &               leftElemPos  = iElem,    &
        &               leftElemPrp  = elemPrp,  &
        &               rightElemId  = neighId,  &
        &               rightElemPos = neighPos, &
        &               rightElemPrp = neighPrp, &
        &               faces        = faces     )

    end do elemLoop

  end subroutine tem_get_faces