tem_build_fromFinerList Subroutine

private subroutine tem_build_fromFinerList(minLevel, maxLevel, nEligibleChildren, faces)

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

Arguments

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

The min refinement level of the mesh.

integer, intent(in) :: maxLevel

The max refinement level of the mesh.

integer, intent(in) :: nEligibleChildren

The number of eligible children for the vertical face dependency

type(tem_face_type), intent(inout) :: faces(minLevel:maxLevel)

The created face descriptor.


Calls

proc~~tem_build_fromfinerlist~~CallsGraph proc~tem_build_fromfinerlist tem_build_fromFinerList proc~tem_isfromfinerface tem_isFromFinerFace proc~tem_build_fromfinerlist->proc~tem_isfromfinerface proc~tem_abort tem_abort proc~tem_isfromfinerface->proc~tem_abort mpi_abort mpi_abort proc~tem_abort->mpi_abort

Called by

proc~~tem_build_fromfinerlist~~CalledByGraph proc~tem_build_fromfinerlist tem_build_fromFinerList proc~tem_build_facelists tem_build_faceLists proc~tem_build_facelists->proc~tem_build_fromfinerlist proc~tem_build_face_info tem_build_face_info proc~tem_build_face_info->proc~tem_build_facelists

Source Code

  subroutine tem_build_fromFinerList( minLevel, maxLevel, nEligibleChildren, &
                                    & faces )
    ! --------------------------------------------------------------------------
    !> The min refinement level of the mesh.
    integer, intent(in) :: minLevel
    !> The max refinement level of the mesh.
    integer, intent(in) :: maxLevel
    !> The number of eligible children for the vertical face dependency
    integer, intent(in) :: nEligibleChildren
    !> The created face descriptor.
    type(tem_face_type),intent(inout) :: faces(minLevel:maxLevel)
    ! --------------------------------------------------------------------------
    integer :: nFromFinerFaces(3,2)
    integer :: iDir, iFace, iFromFinerFace(2), iDep, iLevel
    integer :: childFace, childFaceOp, childPos, childPosOp, isFromFiner
    ! --------------------------------------------------------------------------

    ! Loop over all the levels
    levelLoop: do iLevel = minLevel, maxLevel

      ! First, we count all the element which will be covered by this rank.
      nFromFinerFaces = 0
      do iDir = 1,3
        do iFace = 1, faces(iLevel)%faces(iDir)%faceList%faceId%nVals

          ! decide if the element will be computed on this level and this rank.
          isFromFiner = tem_isFromFinerFace(iFace, faces(iLevel)%faces(iDir))
          if( isFromFiner .eq. tem_left) then
            ! If left element is refined, its right face has to be interpolated.
            nFromFinerFaces(iDir,tem_right) =                          &
              &                     nFromFinerFaces(iDir,tem_right) + 1
          elseif( isFromFiner .eq. tem_right) then
            ! If right element is refined, its left face has to be interpolated.
            nFromFinerFaces(iDir,tem_left) =                           &
              &                      nFromFinerFaces(iDir,tem_left) + 1
          end if

        end do
      end do

      ! Allocate memory for the from finer faces
      do iDir = 1,3
        do iFace = 1,2
          allocate( faces(iLevel)%faces(iDir)%fromFinerFace(iFace)%            &
            &                          elemPos( nFromFinerFaces(iDir, iFace) ))
          allocate( faces(iLevel)%faces(iDir)%fromFinerFace(iFace)%            &
            &      childPos( nEligibleChildren, nFromFinerFaces(iDir, iFace) ))
          allocate( faces(iLevel)%faces(iDir)%fromFinerFace(iFace)%            &
            &                        elemPosOp( nFromFinerFaces(iDir, iFace) ))
          allocate( faces(iLevel)%faces(iDir)%fromFinerFace(iFace)%            &
            &    childPosOp( nEligibleChildren, nFromFinerFaces(iDir, iFace) ))
          allocate( faces(iLevel)%faces(iDir)%fromFinerFace(iFace)%            &
            &                          facePos( nFromFinerFaces(iDir, iFace) ))
        end do
      end do

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

          ! Decide if the face has to be interpolated.
          isFromFiner = tem_isFromFinerFace(iFace, faces(iLevel)%faces(iDir))
          if( isFromFiner .eq. tem_left) then
            ! The left element is refined -> add the left element to the right
            ! faces that have to be interpolated.
            faces(iLevel)%faces(iDir)%fromFinerFace(tem_right)%        &
              &       elemPos(iFromFinerFace(tem_right)) =             &
              &       faces(iLevel)%faces(iDir)%faceList%leftElemPos%val(iFace)
            faces(iLevel)%faces(iDir)%fromFinerFace(tem_right)%        &
              &      elemPosOp(iFromFinerFace(tem_right)) =            &
              &      faces(iLevel)%faces(iDir)%faceList%rightElemPos%val(iFace)
            faces(iLevel)%faces(iDir)%fromFinerFace(tem_right)%        &
              &              facePos(iFromFinerFace(tem_right)) = iFace
            ! Now, we set the dependencies, i.e. we identify the positions of
            ! the child elements and set them in the from finer description.
            do iDep = 1,nEligibleChildren
              childFace = faces(iLevel)%faces(iDir)%faceDep%                   &
                &                                      childFacePos(iDep,iFace)
              childPos = faces(iLevel+1)%faces(iDir)%faceList%                 &
                &                                    leftElemPos%val(childFace)
              faces(iLevel)%faces(iDir)%fromFinerFace(tem_right)%      &
                &   childPos(iDep,iFromFinerFace(tem_right)) = childPos
              childFaceOp = faces(iLevel)%faces(iDir)%faceDep%                 &
                &                                      childFacePosOp(iDep,iFace)
              childPosOp = faces(iLevel+1)%faces(iDir)%faceList%               &
                &                                   rightElemPos%val(childFaceOp)
              faces(iLevel)%faces(iDir)%fromFinerFace(tem_right)%      &
                &            childPosOp(iDep,iFromFinerFace(tem_right))&
                & = childPosOp
            end do

            ! Count the number of right faces
            iFromFinerFace(tem_right) =                                &
              &                           iFromFinerFace(tem_right) + 1

          elseif( isFromFiner .eq. tem_right) then
            ! The right element is refined -> add the right element to the left
            ! faces that have to be interpolated.
            faces(iLevel)%faces(iDir)%fromFinerFace(tem_left)%         &
              &      elemPos(iFromFinerFace(tem_left)) =               &
              &      faces(iLevel)%faces(iDir)%faceList%rightElemPos%val(iFace)
            faces(iLevel)%faces(iDir)%fromFinerFace(tem_left)%         &
              &       elemPosOp(iFromFinerFace(tem_left)) =            &
              &       faces(iLevel)%faces(iDir)%faceList%leftElemPos%val(iFace)
            faces(iLevel)%faces(iDir)%fromFinerFace(tem_left)%         &
              &               facePos(iFromFinerFace(tem_left)) = iFace
            ! Now, we set the dependencies.
            ! Now, we set the dependencies, i.e. we identify the positions of
            ! the child elements and set them in the from finer description.
            do iDep = 1,nEligibleChildren
              childFace = faces(iLevel)%faces(iDir)%faceDep%                   &
                &                                      childFacePos(iDep,iFace)
              childPos = faces(iLevel+1)%faces(iDir)%faceList%                 &
                &                                   rightElemPos%val(childFace)
              faces(iLevel)%faces(iDir)%fromFinerFace(tem_left)%       &
                &    childPos(iDep,iFromFinerFace(tem_left)) = childPos
              childFaceOp = faces(iLevel)%faces(iDir)%faceDep%                 &
                &                                      childFacePosOp(iDep,iFace)
              childPosOp = faces(iLevel+1)%faces(iDir)%faceList%               &
                &                                    leftElemPos%val(childFaceOp)
              faces(iLevel)%faces(iDir)%fromFinerFace(tem_left)%       &
                & childPosOp(iDep,iFromFinerFace(tem_left)) = childPosOp
            end do

            ! Count the number of left faces
            iFromFinerFace(tem_left) =                                 &
              &                            iFromFinerFace(tem_left) + 1

          end if
        end do
      end do
    end do levelLoop

  end subroutine tem_build_fromFinerList