for the a given face before the compute step. Since face data
Type | Intent | Optional | 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. |
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.
pure function tem_isRecvFace(facePos, faces ) result( isRecvFace )
! --------------------------------------------------------------------------
!> 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 :: isRecvFace
! --------------------------------------------------------------------------
integer :: leftPrp, rightPrp
! --------------------------------------------------------------------------
! get the left and right property of the faces.
leftPrp = faces%faceList%leftPrp%val(facePos)
rightPrp = faces%faceList%rightPrp%val(facePos)
isrecvface = -1 ! default value (error)
! So, check for all conditions under which face info has to be received.
! Fluid-fluid face (purely local)
! -> nothing to be sent
if ( (leftPrp == tem_fluidFace_prp) &
& .and. (rightPrp == tem_fluidFace_prp) ) then
isRecvFace = 0
! Fluid-remote face (the halo is not refined)
! -> left element is on my rank, so I will compute this face. Therefore,
! I will also receive the remote information. In this case I will
! receive info from the right adjacent element.
elseif ( (leftPrp == tem_fluidFace_prp) &
& .and. (rightPrp == tem_remoteFace_prp) ) then
isRecvFace = tem_right
! Remote-fluid face (the halo is not refined)
! -> left element is NOT on my rank, so another rank will compute this face.
! So, I have to send info for the left adjacent element, but I do not
! not receive information.
elseif( (leftPrp == tem_remoteFace_prp) &
& .and. (rightPrp == tem_fluidFace_prp) ) then
isRecvFace = 0
! 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 recv any info about any element adjacent to this face.
elseif( leftPrp.eq.tem_remoteFace_prp .and. &
& rightPrp.eq.tem_notExist_prp ) then
isRecvFace = 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 recv any info about any element adjacent to this face.
elseif( leftPrp.eq.tem_notExist_prp .and. &
& rightPrp.eq.tem_remoteFace_prp ) then
isRecvFace = 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
isRecvFace = 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
isRecvFace = 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
isRecvFace = 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
isRecvFace = 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
isRecvFace = 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
isRecvFace = 0
! Remote-fromCoarser face
! -> My rank holds the right elements. Therefore, my element is on a coarser
! level.
! However, the remote elements are refined, but not located on a single
! rank (otherwise it would have the notExist property instead of the
! remoteFace property). Therefore, my rank will send data to the remote
! ranks (as they are on a finer level), but we do not receive any info.
elseif( leftPrp.eq.tem_remoteFace_prp .and. &
& rightPrp.eq.tem_fromCoarserFace_prp ) then
isRecvFace = 0
! FormCoarser-Remote face
! -> My rank holds the left element. Therefore, my element is on a coarser
! level. However, the remote elements are refined, but not located on a
! single rank (otherwise it would have the notExist property instead of
! the remoteFace property). Therefore, my rank will send data to the
! remote ranks (as they are on a finer level), but we do not receive any
! info.
elseif( leftPrp.eq.tem_fromCoarserFace_prp .and. &
& rightPrp.eq.tem_remoteFace_prp ) then
isRecvFace = 0
! Remote-Remote face
! -> My rank holds none of the elements located next to this face. So we do
! not receive any information for this face.
elseif( leftPrp.eq.tem_remoteFace_prp .and. &
& rightPrp.eq.tem_remoteFace_prp ) then
isRecvFace = 0
! notExist - remote face
elseif( leftPrp.eq.tem_notExist_prp .and. &
& rightPrp.eq.tem_remoteFace_prp) then
isRecvFace = 0
! remote - notExist face
elseif( leftPrp.eq.tem_remoteFace_prp .and. &
& rightPrp.eq.tem_notExist_prp ) then
isRecvFace = 0
! fluid-remote+fromCoaser face
elseif( leftPrp.eq.tem_fluidFace_prp .and. &
& rightPrp.eq.ior( tem_remoteFace_prp, tem_fromCoarserFace_prp )) then
isRecvFace = tem_right
! remote+fromCoarser-fluid face
elseif( leftPrp.eq.ior(tem_remoteFace_prp,tem_fromCoarserFace_prp) .and. &
& rightPrp.eq.tem_fluidFace_prp ) then
isRecvFace = tem_left
! notExist-remote+fromCoaser face (nothing to be received, this rank does
! not have any element shared by this face.)
elseif( leftPrp.eq.tem_notExist_prp .and. &
& rightPrp.eq.ior(tem_remoteFace_prp, tem_fromCoarserFace_prp)) then
isRecvFace = 0
! remote+fromCoarser-remote face (nothing to be received, everything is
! remote)
elseif( leftPrp.eq.ior(tem_remoteFace_prp, tem_fromCoarserFace_prp) .and. &
& rightPrp.eq.tem_notExist_prp ) then
isRecvFace = 0
! remote+fromCoarser-remote+fromCoaser face (nothing to be received,
! everything is remote)
elseif( leftPrp.eq.ior(tem_remoteFace_prp, tem_fromCoarserFace_prp) .and. &
& rightPrp.eq.ior(tem_remoteFace_prp, tem_fromCoarserFace_prp)) then
isRecvFace = 0
! fromFiner - not exist face (nothing to be done)
elseif( leftPrp.eq.tem_fromFinerFace_prp .and. &
& rightPrp.eq.tem_notExist_prp ) then
isRecvFace = 0
! notExist - fromFiner face (nothing to be done)
elseif( leftPrp.eq.tem_notExist_prp .and. &
& rightPrp.eq.tem_fromFinerFace_prp ) then
isRecvFace = 0
! fromFiner+remote - fluid face (we do not have to recv info here, but we
! send it)
elseif( leftPrp.eq.ior(tem_fromFinerFace_prp,tem_remoteFace_prp) .and. &
& rightPrp.eq.tem_fluidFace_prp ) then
isRecvFace = 0
! fluid - remote+fromFiner face (we do not recv info here, but we send it)
elseif( leftPrp.eq.tem_fluidFace_prp .and. &
& rightPrp.eq.ior(tem_fromFinerFace_prp,tem_remoteFace_prp)) then
isRecvFace = 0
! fromCoarser - notExist face (nothing to receive)
elseif( leftPrp.eq.tem_fromCoarserFace_prp .and. &
& rightPrp.eq.tem_notExist_prp ) then
isRecvFace = 0
! notExist - fromCoarser face (nothing to receive)
elseif( leftPrp.eq.tem_notExist_prp .and. &
& rightPrp.eq.tem_fromCoarserFace_prp ) then
isRecvFace = 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
isRecvFace = 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
isRecvFace = 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
isRecvFace = 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
isRecvFace = 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
isRecvFace = 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
isRecvFace = 0
! fromFiner+remote - fromCoarser face (we do not have to recv info here,
! but we send it)
elseif( leftPrp.eq.ior(tem_fromFinerFace_prp,tem_remoteFace_prp) .and. &
& rightPrp.eq.tem_fromCoarserFace_prp ) then
isRecvFace = 0
! fromCoarser - remote+fromFiner face (we do not recv info here, but we
! send it)
elseif( leftPrp.eq.tem_fromCoarserFace_prp .and. &
& rightPrp.eq.ior(tem_fromFinerFace_prp,tem_remoteFace_prp)) then
isRecvFace = 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
isRecvFace = 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
isRecvFace = 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
isRecvFace = 0
! fromFiner - remote face (nothing to be done)
elseif( leftPrp.eq.tem_fromFinerFace_prp .and. &
& rightPrp.eq.tem_remoteFace_prp ) then
isRecvFace = 0
! remote - fromFiner face (nothing to be done)
elseif( leftPrp.eq.tem_remoteFace_prp .and. &
& rightPrp.eq.tem_fromFinerFace_prp ) then
isRecvFace = 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
isRecvFace = 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
isRecvFace = 0
! fluid - boundary face (nothing to be done)
elseif(leftPrp.eq.tem_fluidFace_prp .and. rightPrp.eq.tem_bndFace_prp) then
isRecvFace = 0
! boundary - fluid face (nothing to be done)
elseif(leftPrp.eq.tem_bndFace_prp .and. rightPrp.eq.tem_fluidFace_prp) then
isRecvFace = 0
! remote - boundary face (nothing to be done)
elseif(leftPrp.eq.tem_remoteFace_prp .and. rightPrp.eq.tem_bndFace_prp) then
isRecvFace = 0
! boundary - remote face (nothing to be done)
elseif(leftPrp.eq.tem_bndFace_prp .and. rightPrp.eq.tem_remoteFace_prp) then
isRecvFace = 0
! fromFiner - boundary face (nothing to be done)
elseif( leftPrp.eq.tem_fromFinerFace_prp .and. &
& rightPrp.eq.tem_bndFace_prp ) then
isRecvFace = 0
! boundary - fromFiner face (nothing to be done)
elseif( leftPrp.eq.tem_bndFace_prp .and. &
& rightPrp.eq.tem_fromFinerFace_prp ) then
isRecvFace = 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
isRecvFace = 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
isRecvFace = 0
! fromCoarser - boundary face (nothing to be done)
elseif( leftPrp.eq.tem_fromCoarserFace_prp .and. &
& rightPrp.eq.tem_bndFace_prp ) then
isRecvFace = 0
! boundary - fromCoarser face (nothing to be done)
elseif( leftPrp.eq.tem_bndFace_prp .and. &
& rightPrp.eq.tem_fromCoarserFace_prp ) then
isRecvFace = 0
end if
end function tem_isRecvFace