tem_faceDep_verticalDown Subroutine

private subroutine tem_faceDep_verticalDown(coarseFaces, fineFaces, dir, nEligibleChildren)

Arguments

Type IntentOptional Attributes Name
type(tem_face_descriptor_type), intent(inout) :: coarseFaces

Face description on the coarse level. The dependencies to the finer level will be appended to this face descriptor.

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

Face description on the finer 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


Calls

proc~~tem_facedep_verticaldown~~CallsGraph proc~tem_facedep_verticaldown tem_faceDep_verticalDown proc~tem_adddep_down tem_addDep_down proc~tem_facedep_verticaldown->proc~tem_adddep_down proc~tem_getface_prp tem_getFace_prp proc~tem_facedep_verticaldown->proc~tem_getface_prp proc~tem_reqdep_down tem_reqDep_down proc~tem_facedep_verticaldown->proc~tem_reqdep_down interface~positionofval~5 positionofval proc~tem_adddep_down->interface~positionofval~5 proc~tem_abort tem_abort proc~tem_adddep_down->proc~tem_abort proc~tem_coordofid tem_CoordOfId proc~tem_adddep_down->proc~tem_coordofid proc~tem_directchildren tem_directChildren proc~tem_adddep_down->proc~tem_directchildren proc~tem_eligiblechildren tem_eligibleChildren proc~tem_adddep_down->proc~tem_eligiblechildren proc~tem_idofcoord tem_IdOfCoord proc~tem_adddep_down->proc~tem_idofcoord proc~posofval_label posofval_label interface~positionofval~5->proc~posofval_label mpi_abort mpi_abort proc~tem_abort->mpi_abort proc~tem_levelof tem_LevelOf proc~tem_coordofid->proc~tem_levelof interface~sortedposofval~5 sortedposofval proc~posofval_label->interface~sortedposofval~5

Called by

proc~~tem_facedep_verticaldown~~CalledByGraph proc~tem_facedep_verticaldown tem_faceDep_verticalDown proc~tem_facedep_vertical tem_faceDep_vertical proc~tem_facedep_vertical->proc~tem_facedep_verticaldown proc~tem_build_face_info tem_build_face_info proc~tem_build_face_info->proc~tem_facedep_vertical

Source Code

  subroutine tem_faceDep_verticalDown( coarseFaces, fineFaces, dir, &
                                     & nEligibleChildren )
    ! --------------------------------------------------------------------------
    !> Face description on the coarse level. The dependencies to the finer
    !! level will be appended to this face descriptor.
    type(tem_face_descriptor_type), intent(inout) :: coarseFaces
    !> Face description on the finer 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 :: iFace
    integer :: leftPrp, rightPrp
    ! --------------------------------------------------------------------------

    ! Run over all the faces and add the downward depndency (depending on the
    ! properties of the face).
    faceLoop: do iFace = 1, coarseFaces%faceList%faceId%nVals

      ! get its properties from left and right
      call tem_getFace_prp( coarseFaces, iFace, leftPrp, rightPrp )

      ! check if this face requires downwards dependencies. If yes, we just add
      ! this dependency to the coarse face descriptor.
      if( tem_reqDep_down(leftPrp, rightPrp) ) then
        call tem_addDep_down( iFace, coarseFaces, fineFaces, &
                            & dir, nEligibleChildren )
      end if

    end do faceLoop

  end subroutine tem_faceDep_verticalDown