tem_isSendFace Function

private pure function tem_isSendFace(facePos, faces) result(isSendFace)

for the a given face before the compute step.

Arguments

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

The position of the face in faces to be checked.

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

The description of the faces.

Return Value integer

Returns 0 if no information has to be received from remote ranks before a compute step can be done. In case that information has to be received it will return a non-zero positive value. Furhtermore, it determines whether information for the left face value (left adjacent element) or the right face value (right adjacent element) has to be received. In the first case it will return tem_left and in the second case it will return tem_right. If something went wrong, e.g. an unknown combination occurs, it will return -1.


Called by

proc~~tem_issendface~~CalledByGraph proc~tem_issendface tem_isSendFace proc~tem_build_facesendbuffers tem_build_faceSendBuffers proc~tem_build_facesendbuffers->proc~tem_issendface proc~tem_build_facebuffers tem_build_faceBuffers proc~tem_build_facebuffers->proc~tem_build_facesendbuffers proc~tem_build_face_info tem_build_face_info proc~tem_build_face_info->proc~tem_build_facebuffers

Source Code

  pure function tem_isSendFace(facePos, faces) result( isSendFace )
    ! --------------------------------------------------------------------------
    !> The position of the face in faces to be checked.
    integer , intent(in) :: facePos
    !> The description of the faces.
    type(tem_face_descriptor_type),intent(in)  :: faces
    !> Returns 0 if no information has to be received from remote ranks
    !! before a compute step can be done. In case that information has to be
    !! received it will return a non-zero positive value. Furhtermore,
    !! it determines whether information for the left face value
    !! (left adjacent element)
    !! or the right face value (right adjacent element) has to be received.
    !! In the first case it will return
    !! [[tem_faceData_module:tem_left]] and in the second case
    !! it will return
    !! [[tem_faceData_module:tem_right]].
    !! If something went wrong, e.g. an unknown combination occurs,
    !! it will return -1.
    integer :: isSendFace
    ! --------------------------------------------------------------------------
    integer :: leftPrp, rightPrp
    ! --------------------------------------------------------------------------

    leftPrp = faces%faceList%leftPrp%val(facePos)
    rightPrp = faces%faceList%rightPrp%val(facePos)

    issendface = -1 ! default value (error)

    ! So, check for all conditions under which face info has to be send.

    ! Fluid-fluid face (purely local)
    ! -> nothing to send
    if( leftPrp.eq.tem_fluidFace_prp .and. rightPrp.eq.tem_fluidFace_prp ) then
      isSendFace = 0

    ! Fluid-remote face (the halo is not refined)
    ! -> left element is on my rank, so I will compute it. Therefore, I will
    !    receive information about the right element, but I do not have to send
    !    anything.
    elseif( leftPrp.eq.tem_fluidFace_prp .and.                                 &
      &     rightPrp.eq.tem_remoteFace_prp ) then
      isSendFace = 0

    ! Remote-fluid face (the halo is not refined)
    ! -> left element is on a remote rank, so I will not compute it. Therfore
    !    I will send the element to the remote rank.
    elseif( leftPrp.eq.tem_remoteFace_prp .and.                                &
      &     rightPrp.eq.tem_fluidFace_prp ) then
      isSendFace = tem_right

    ! Remote-notExisting face (the remote is not refined)
    ! -> remote is not located on my rank and the other element does not
    !    exist. This face does not have to be computed on any rank, so my rank
    !    is not sending any info about any element adjacent to this face.
    elseif( leftPrp.eq.tem_remoteFace_prp .and.                                &
      &     rightPrp.eq.tem_notExist_prp ) then
      isSendFace = 0

    ! notExisting-remote face (the remote is not refined)
    ! -> remote is not located on my rank and the other element does not
    !    exist. This face does not have to be computed on any rank, so my rank
    !    is not send any info about any element adjacent to this face.
    elseif( leftPrp.eq.tem_notExist_prp .and.                                  &
      &     rightPrp.eq.tem_remoteFace_prp ) then
      isSendFace = 0

    ! Fluid-fromFiner face (purely local)
    ! -> everything is local, so there is no need to send or recv information
    !    from another rank.
    elseif( leftPrp.eq.tem_fluidFace_prp .and.                                 &
      &     rightPrp.eq.tem_fromFinerFace_prp ) then
      isSendFace = 0

    ! FromFiner-fluid face (purely local)
    ! -> everything is local, so there is no need to send or recv information
    !    from another rank.
    elseif( leftPrp.eq.tem_fromFinerFace_prp .and.                             &
      &     rightPrp.eq.tem_fluidFace_prp ) then
      isSendFace = 0

    ! FromFiner-FromFiner face (purely local)
    ! -> everything is local. Somehow two ghost elements meet each other, but
    !    we do not have to recv or send any data for these faces.
    elseif( leftPrp.eq.tem_fromFinerFace_prp .and.                             &
      &     rightPrp.eq.tem_fromFinerFace_prp ) then
      isSendFace = 0

    ! FromCoarser-fluid face (purely local)
    ! -> everything is local, so no data has to be send or received.
    elseif( leftPrp.eq.tem_fromCoarserFace_prp .and.                           &
      &     rightPrp.eq.tem_fluidFace_prp ) then
      isSendFace = 0

    ! Fluid-fromCoarser face (purely local)
    ! -> everything is local, so no data has to be send or received.
    elseif( leftPrp.eq.tem_fluidFace_prp .and.                                 &
      &     rightPrp.eq.tem_fromCoarserFace_prp ) then
      isSendFace = 0

    ! FromCoarser-fromCoarser face (purely local)
    ! -> everything is local, so no data has to be send or received.
    elseif( leftPrp.eq.tem_fromCoarserFace_prp .and.                           &
      &     rightPrp.eq.tem_fromCoarserFace_prp ) then
      isSendFace = 0

    ! Remote-fromCoarser face
    ! -> The remote fluid elements are located on a finer level
    elseif( leftPrp.eq.tem_remoteFace_prp .and.                                &
      &     rightPrp.eq.tem_fromCoarserFace_prp ) then
      isSendFace = tem_right

    ! FromCoarser-remote face
    ! -> The remote fluid elements are located on a finer level
    elseif( leftPrp.eq.tem_fromCoarserFace_prp .and.                           &
      &     rightPrp.eq.tem_remoteFace_prp ) then
      isSendFace = tem_left

    ! Remote-remote face
    ! -> The halo elements meet each other. This can happen, but we do not have
    !    to send or receive information for this face, because we do not
    !    compute anything here.
    elseif( leftPrp.eq.tem_remoteFace_prp .and.                                &
      &     rightPrp.eq.tem_remoteFace_prp ) then
      isSendFace = 0

    ! notExist - remote face (nothing to be send, we do not have info about
    ! adjacent elements)
    elseif( leftPrp.eq.tem_notExist_prp .and.                                  &
      &     rightPrp.eq.tem_remoteFace_prp ) then
      isSendFace = 0

    ! remote - notExist face (nothing to be send, we do not have any info about
    ! adjacent elements)
    elseif( leftPrp.eq.tem_remoteFace_prp .and.                                &
      &     rightPrp.eq.tem_notExist_prp ) then
      isSendFace = 0

    ! fluid-remote+fromCoaser face (here, we receive data in some situations,
    ! but we never send info in this case)
    elseif( leftPrp.eq.tem_fluidFace_prp .and.                                 &
      &     rightPrp.eq.ior(tem_remoteFace_prp, tem_fromCoarserFace_prp)) then
      isSendFace = 0

    ! remote+fromCoarser-fluid face (here, we receive data in some situations,
    ! but we never send info for this faces)
    elseif( leftPrp.eq.ior(tem_remoteFace_prp,tem_fromCoarserFace_prp) .and.   &
      &     rightPrp.eq.tem_fluidFace_prp ) then
      isSendFace = 0

    ! notExist-remote+fromCoaser face (nothing to be send here)
    elseif( leftPrp.eq.tem_notExist_prp .and.                                  &
      &     rightPrp.eq.ior(tem_remoteFace_prp, tem_fromCoarserFace_prp)) then
      isSendFace = 0

    ! remote+fromCoarser-remote face (everything is remote, so nothing to be
    ! sent)
    elseif( leftPrp.eq.ior(tem_remoteFace_prp, tem_fromCoarserFace_prp) .and.  &
      &     rightPrp.eq.tem_notExist_prp ) then
      isSendFace = 0

    ! remote+fromCoarser-remote+fromCoaser face (everything is remote, so
    ! nothing to be sent)
    elseif( leftPrp.eq.ior(tem_remoteFace_prp, tem_fromCoarserFace_prp) .and.  &
      &     rightPrp.eq.ior(tem_remoteFace_prp, tem_fromCoarserFace_prp)) then
      isSendFace = 0

    ! fromFiner - not exist face (nothing to be send here).
    elseif( leftPrp.eq.tem_fromFinerFace_prp .and.                             &
      &     rightPrp.eq.tem_notExist_prp ) then
      isSendFace = 0

    ! notExist - fromFiner face (nothing to be send here)
    elseif( leftPrp.eq.tem_notExist_prp .and.                                  &
      &     rightPrp.eq.tem_fromFinerFace_prp ) then
      isSendFace = 0

    ! fromFiner+remote - fluid face (will be send on the finer level)
    elseif( leftPrp.eq.ior(tem_fromFinerFace_prp,tem_remoteFace_prp) .and.     &
      &     rightPrp.eq.tem_fluidFace_prp ) then
      isSendFace = 0

    ! fluid - remote+fromFiner face (will be send on the finer level)
    elseif( leftPrp.eq.tem_fluidFace_prp .and.                                 &
      &     rightPrp.eq.ior(tem_fromFinerFace_prp,tem_remoteFace_prp)) then
      isSendFace = 0

    ! fromCoarser - notExist face (nothing to send)
    elseif( leftPrp.eq.tem_fromCoarserFace_prp .and.                           &
      &     rightPrp.eq.tem_notExist_prp ) then
      isSendFace = 0

    ! notExist - fromCoarser face (nothing to send)
    elseif( leftPrp.eq.tem_notExist_prp .and.                                  &
      &     rightPrp.eq.tem_fromCoarserFace_prp ) then
      isSendFace = 0

    ! remote+fromFiner-remote face (everything is remote)
    elseif( leftPrp.eq.ior(tem_fromFinerFace_prp, tem_remoteFace_prp) .and.    &
      &     rightPrp.eq.tem_remoteFace_prp ) then
      isSendFace = 0

    ! remote-remote+fromFiner face (everything is remote)
    elseif( leftPrp.eq. tem_remoteFace_prp .and.                               &
      &     rightPrp.eq.ior(tem_remoteFace_prp, tem_fromFinerFace_prp)) then
      isSendFace = 0

    ! notExist-remote+fromFiner face (nothing is on my rank)
    elseif( leftPrp.eq. tem_notExist_prp .and.                                 &
      &     rightPrp.eq.ior(tem_remoteFace_prp, tem_fromFinerFace_prp)) then
      isSendFace = 0

    ! remote+fromFiner-notExist face (nothing is on my rank)
    elseif( leftPrp.eq.ior(tem_remoteFace_prp, tem_fromFinerFace_prp) .and.    &
      &     rightPrp.eq.tem_notExist_prp ) then
      isSendFace = 0

    ! fromFiner-fromCoarser face (everything is local, this is an intermediate
    ! face level (i.e. leveljumps of difference larger than 1 occure here) )
    elseif( leftPrp.eq.tem_fromFinerFace_prp .and.                             &
      &     rightPrp.eq.tem_fromCoarserFace_prp ) then
      isSendFace = 0

    ! fromCoarser-fromFiner face (everything is local, this is an intermediate
    ! face level (i.e. leveljumps of difference larger than 1 occure here) )
    elseif( leftPrp.eq.tem_fromCoarserFace_prp .and.                           &
      &     rightPrp.eq.tem_fromFinerFace_prp ) then
      isSendFace = 0

    ! fromFiner+remote - fromCoarser face (send on a finer level)
    elseif( leftPrp.eq.ior(tem_fromFinerFace_prp,tem_remoteFace_prp) .and.     &
      &     rightPrp.eq.tem_fromCoarserFace_prp ) then
      isSendFace = 0

    ! fromCoarser - remote+fromFiner face (send on a finer level)
    elseif( leftPrp.eq.tem_fromCoarserFace_prp .and.                           &
      &     rightPrp.eq.ior(tem_fromFinerFace_prp,tem_remoteFace_prp)) then
      isSendFace = 0

    ! remote+fromFiner - remote+fromFiner face (everything is remote, nothing
    ! to be done)
    elseif( leftPrp.eq.ior(tem_fromFinerFace_prp,tem_remoteFace_prp) .and.     &
      &     rightPrp.eq.ior(tem_fromFinerFace_prp,tem_remoteFace_prp)) then
      isSendFace = 0

    ! remote+fromFiner - fromFiner face (everything is from finer, nothing to
    ! be done)
    elseif( leftPrp.eq.ior(tem_fromFinerFace_prp,tem_remoteFace_prp) .and.     &
      &     rightPrp.eq.tem_fromFinerFace_prp ) then
      isSendFace = 0

    ! fromFiner - remote+fromFiner face (everything is from finer, nothing to
    ! be done)
    elseif( leftPrp.eq.tem_fromFinerFace_prp .and.                             &
      &     rightPrp.eq.ior(tem_fromFinerFace_prp,tem_remoteFace_prp)) then
      isSendFace = 0

    ! fromFiner - remote face (nothing to be done)
    elseif( leftPrp.eq.tem_fromFinerFace_prp .and.                             &
      &     rightPrp.eq.tem_remoteFace_prp ) then
      isSendFace = 0

    ! remote - fromFiner face (nothing to be done)
    elseif( leftPrp.eq.tem_remoteFace_prp .and.                                &
      &     rightPrp.eq.tem_fromFinerFace_prp ) then
      isSendFace = 0

    ! remote+fromCoarser - fromCoarser face (elements are on coarser
    ! level, here and remote so nothing to be done)
    elseif( leftPrp.eq.ior(tem_fromCoarserFace_prp,tem_remoteFace_prp) .and.   &
      &     rightPrp.eq.tem_fromCoarserFace_prp ) then
      isSendFace = 0

    ! fromCoarser - remote+fromCoarser face (elements are on coarser
    ! level, here and remote so nothing to be done)
    elseif( leftPrp.eq.tem_fromCoarserFace_prp .and.                           &
      &     rightPrp.eq.ior(tem_fromCoarserFace_prp,tem_remoteFace_prp)) then
      isSendFace = 0

    ! fluid - boundary face (nothing to be done)
    elseif(leftPrp.eq.tem_fluidFace_prp .and. rightPrp.eq.tem_bndFace_prp) then
      isSendFace = 0

    ! boundary - fluid face (nothing to be done)
    elseif(leftPrp.eq.tem_bndFace_prp .and. rightPrp.eq.tem_fluidFace_prp) then
      isSendFace = 0

    ! remote - boundary face (nothing to be done)
    elseif(leftPrp.eq.tem_remoteFace_prp .and. rightPrp.eq.tem_bndFace_prp) then
      isSendFace = 0

    ! boundary - remote face (nothing to be done)
    elseif(leftPrp.eq.tem_bndFace_prp .and. rightPrp.eq.tem_remoteFace_prp) then
      isSendFace = 0

    ! fromFiner - boundary face (nothing to be done)
    elseif( leftPrp.eq.tem_fromFinerFace_prp .and.                             &
      &     rightPrp.eq.tem_bndFace_prp ) then
      isSendFace = 0

    ! boundary - fromFiner face (nothing to be done)
    elseif( leftPrp.eq.tem_bndFace_prp .and.                                   &
      &     rightPrp.eq.tem_fromFinerFace_prp ) then
      isSendFace = 0

    ! fromFiner+remote - boundary face (nothing to be done)
    elseif( leftPrp.eq.ior(tem_fromFinerFace_prp,tem_remoteFace_prp) .and.     &
      &     rightPrp.eq.tem_bndFace_prp ) then
      isSendFace = 0

    ! boundary - fromFiner+remote face (nothing to be done)
    elseif( leftPrp.eq.tem_bndFace_prp .and.                                   &
      &     rightPrp.eq.ior(tem_fromFinerFace_prp,tem_remoteFace_prp)) then
      isSendFace = 0

    ! fromCoarser - boundary face (nothing to be done)
    elseif( leftPrp.eq.tem_fromCoarserFace_prp .and.                           &
      &     rightPrp.eq.tem_bndFace_prp ) then
      isSendFace = 0

    ! boundary - fromCoarser face (nothing to be done)
    elseif( leftPrp.eq.tem_bndFace_prp .and.                                   &
      &     rightPrp.eq.tem_fromCoarserFace_prp ) then
      isSendFace = 0

    end if

  end function tem_isSendFace