tem_build_computeList Subroutine

private subroutine tem_build_computeList(faces, nEligibleChildren)

Extracts all the faces which will be computed by current rank.

Arguments

Type IntentOptional Attributes Name
type(tem_face_type), intent(inout) :: faces

The created face descriptor.

integer, intent(in) :: nEligibleChildren

The number of eligible children for the vertical face dependency


Calls

proc~~tem_build_computelist~~CallsGraph proc~tem_build_computelist tem_build_computeList proc~tem_iscomputeface tem_isComputeFace proc~tem_build_computelist->proc~tem_iscomputeface proc~tem_abort tem_abort proc~tem_iscomputeface->proc~tem_abort mpi_abort mpi_abort proc~tem_abort->mpi_abort

Called by

proc~~tem_build_computelist~~CalledByGraph proc~tem_build_computelist tem_build_computeList proc~tem_build_facelists tem_build_faceLists proc~tem_build_facelists->proc~tem_build_computelist proc~tem_build_face_info tem_build_face_info proc~tem_build_face_info->proc~tem_build_facelists

Source Code

  subroutine tem_build_computeList( faces, nEligibleChildren )
    ! --------------------------------------------------------------------------
    !> The created face descriptor.
    type(tem_face_type),intent(inout) :: faces
    !> The number of eligible children for the vertical face dependency
    integer, intent(in) :: nEligibleChildren
    ! --------------------------------------------------------------------------
    integer :: nComputeFaces(3)
    integer :: iDir, iFace, iComputeFace
    ! --------------------------------------------------------------------------


    ! First, we count all the element which will be covered by this rank.
    nComputeFaces = 0
    do iDir = 1,3

      do iFace = 1, faces%faces(iDir)%faceList%faceId%nVals

        ! decide if the element will be computed on this level and this rank.
        if ( tem_isComputeFace(iFace,             &
          &                    faces%faces(iDir), &
          &                    nEligibleChildren) ) then
          nComputeFaces(iDir) = nComputeFaces(iDir) + 1
        end if

      end do

    end do

    ! Allocate memory for the compute faces
    do iDir = 1,3
      allocate( faces%faces(iDir)%computeFace%leftPos( nComputeFaces(iDir) ) )
      allocate( faces%faces(iDir)%computeFace%rightPos( nComputeFaces(iDir) ) )
      allocate( faces%faces(iDir)%computeFace%facePos( nComputeFaces(iDir) ) )
    end do

    ! Now, we add all the compute elements to the iterator
    do iDir = 1,3
      iComputeFace = 1
      do iFace = 1, faces%faces(iDir)%faceList%faceId%nVals

        ! Decide if the element will be computed on this level and this rank.
        ! If yes, we add it to the list of compute faces.
        if ( tem_isComputeFace(iFace,             &
          &                    faces%faces(iDir), &
          &                    nEligibleChildren) ) then
          faces%faces(iDir)%computeFace%leftPos(iComputeFace) =                &
            &                 faces%faces(iDir)%faceList%leftElemPos%val(iFace)
          faces%faces(iDir)%computeFace%rightPos(iComputeFace) =               &
            &                faces%faces(iDir)%faceList%rightElemPos%val(iFace)
          faces%faces(iDir)%computeFace%facePos(iComputeFace) = iFace
          iComputeFace = iComputeFace + 1
        end if

      end do
    end do
  end subroutine tem_build_computeList