descriptor.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | coarseFacePos |
The position of the coarse face in coarseFaces you want to add the child dependencies for. |
||
type(tem_face_descriptor_type), | intent(inout) | :: | coarseFaces |
The description of the faces on the coarse level. The dependency will be added to this face descriptor. |
||
type(tem_face_descriptor_type), | intent(in) | :: | fineFaces |
The descriptor of the faces on the fine level. |
||
integer, | intent(in) | :: | dir |
The spatial direction of the face to add the downward dependencies for. 1 --> x direction 2 --> y direction 3 --> z direction |
||
integer, | intent(in) | :: | nEligibleChildren |
The number of eligible children for the vertical face dependency |
subroutine tem_addDep_down( coarseFacePos, coarseFaces, fineFaces, dir, &
& nEligibleChildren )
! --------------------------------------------------------------------------
!> The position of the coarse face in coarseFaces you want to add
!! the child dependencies for.
integer, intent(in) :: coarseFacePos
!> The description of the faces on the coarse level. The dependency will be
!! added to this face descriptor.
type(tem_face_descriptor_type), intent(inout) :: coarseFaces
!> The descriptor of the faces on the fine level.
type(tem_face_descriptor_type), intent(in) :: fineFaces
!> The spatial direction of the face to add the downward dependencies for.
!! 1 --> x direction
!! 2 --> y direction
!! 3 --> z direction
integer, intent(in) :: dir
!> The number of eligible children for the vertical face dependency
integer, intent(in) :: nEligibleChildren
! --------------------------------------------------------------------------
integer(kind=long_k) :: faceId, rightElemId, childFaceId, childFaceIdOp, &
& leftOf_childFaceIdOp
integer :: leftOf_coord(4)
integer(kind=long_k) :: children(8), childrenOp(8)
integer, allocatable :: eligibleChildren(:), eligibleChildrenOp(:)
integer :: eligChildPar, eligChildParOp, iChild, childFacePos, childFacePosOp
! --------------------------------------------------------------------------
! Left element id of the face.
faceId = coarseFaces%faceList%faceId%val(coarseFacePos)
rightElemId = coarseFaces%faceList%rightElemId%val(coarseFacePos)
! Identify the all the eligible child elements.
select case(nEligibleChildren)
case(4)
select case(dir)
case(1) ! x - direction
eligChildPar = q__E
eligChildParOp = q__W
case(2) ! y - direction
eligChildPar = q__N
eligChildParOp = q__S
case(3) ! z - direction
eligChildPar = q__T
eligChildParOp = q__B
case default
write(*,*) 'ERROR in tem_addDep_down: unknown spatial direction, '// &
& 'stopping ...'
call tem_abort()
end select
call tem_eligibleChildren(eligibleChildren, eligChildPar)
call tem_eligibleChildren(eligibleChildrenOp, eligChildParOp)
children = tem_directChildren(faceId)
childrenOp = tem_directChildren(rightElemId)
case(8)
allocate(eligibleChildren(nEligibleChildren))
allocate(eligibleChildrenOp(nEligibleChildren))
eligibleChildren = (/ 1, 2, 3, 4, 5, 6, 7, 8 /)
eligibleChildrenOp = (/ 1, 2, 3, 4, 5, 6, 7, 8 /)
children = tem_directChildren(faceId)
childrenOp = tem_directChildren(rightElemId)
case default
write(logUnit(1),*) 'ERROR in tem_addDep_down: Unknown number of'
write(logUnit(1),*) ' eligible children, stopping ...'
call tem_abort()
end select
! Loop over the relevant, eligible children.
childLoop: do iChild = 1,nEligibleChildren
! The tree id of the current eligible child
childFaceId = children(eligibleChildren(iChild))
childFaceIdOp = childrenOp(eligibleChildrenOp(iChild))
! Get the TreeID which is left of childFaceIdOp for face identification
leftOf_coord = tem_CoordOfId( childFaceIdOp )
leftOf_coord(dir) = leftOf_coord(dir)-1
leftOf_childFaceIdOp = tem_IdOfCoord(leftOf_coord)
! Get its position in the fine face descriptor.
childFacePos = PositionOfVal( me = fineFaces%faceList%faceId, &
& val = childFaceId )
! Get the op position in the fine face descriptor.
childFacePosOp = PositionOfVal( me = fineFaces%faceList%faceId, &
& val = leftOf_childFaceIdOp )
! If the child element exists in the fine face descriptor, add a
! dependency from coarse element to its child element.
if(childFacePos<=0) then
coarseFaces%faceDep%childFaceId(iChild, coarseFacePos) = -1_long_k
coarseFaces%faceDep%childFacePos(iChild, coarseFacePos) = -1
else
coarseFaces%faceDep%childFaceId(iChild, coarseFacePos) = childFaceId
coarseFaces%faceDep%childFacePos(iChild, coarseFacePos) = childFacePos
end if
if(childFacePosOp<=0) then
coarseFaces%faceDep%childFaceIdOp(iChild, coarseFacePos) = -1_long_k
coarseFaces%faceDep%childFacePosOp(iChild, coarseFacePos) = -1
else
coarseFaces%faceDep%childFaceIdOp(iChild, coarseFacePos) = leftOf_childFaceIdOp
coarseFaces%faceDep%childFacePosOp(iChild, coarseFacePos) = childFacePosOp
end if
end do childLoop
end subroutine tem_addDep_down