boxBoxOverlap Function

private function boxBoxOverlap(center_a, dim_a, center_b, dim_b, norm_b) result(overlap)

This function checks for intersection of a axis aligned box and a parallelepiped.

Arguments

Type IntentOptional Attributes Name
real(kind=rk), intent(in) :: center_a(3)
real(kind=rk), intent(in) :: dim_a(3)
real(kind=rk), intent(in) :: center_b(3)
real(kind=rk), intent(in) :: dim_b(3)

Halflength of the parallelepiped in each direction.

Todo

HK: Dim is a bad naming here, halflength would be better.

real(kind=rk), intent(in) :: norm_b(3,3)

Return Value logical


Called by

proc~~boxboxoverlap~~CalledByGraph proc~boxboxoverlap boxBoxOverlap proc~tem_boxcubeoverlap tem_boxCubeOverlap proc~tem_boxcubeoverlap->proc~boxboxoverlap proc~tem_cano_initsubtree tem_cano_initSubTree proc~tem_cano_initsubtree->proc~tem_boxcubeoverlap proc~tem_shape_subtreefromgeominters tem_shape_subTreeFromGeomInters proc~tem_shape_subtreefromgeominters->proc~tem_cano_initsubtree proc~tem_shape2subtree tem_shape2subTree proc~tem_shape2subtree->proc~tem_shape_subtreefromgeominters

Source Code

  function boxBoxOverlap(center_a, dim_a, center_b, dim_b, norm_b) &
    &      result(overlap)

    real(kind=rk), intent(in) :: center_a(3) !< Center of the axis aligned box
    real(kind=rk), intent(in) :: dim_a(3) !< Length of the AAB in each direction
    real(kind=rk), intent(in) :: center_b(3) !< Center of the parallelepiped

    !> Halflength of the parallelepiped in each direction.
    !!@todo HK: Dim is a bad naming here, halflength would be better.
    real(kind=rk), intent(in) :: dim_b(3)

    real(kind=rk), intent(in) :: norm_b(3,3) !< normal vector of parallelepiped
    logical :: overlap


    ! Local Variables

    real(kind=rk)             :: vec_ab(3) ! Vector connecting center_a and
                                           ! center_b
!HK!    real(kind=rk)             :: coeffs(3,3)  ! Matrix for Coefficients needed
!HK!                                              ! for the algorithm.
!HK!
!HK!    real(kind=rk)             :: norm_a(3,3) ! edge directions of the axis
    !                                          aligned box

    ! Initialize variable
    overlap = .true.

    vec_ab = center_b - center_a

!HK!    norm_a(:,1) = (/1., 0., 0./)
!HK!    norm_a(:,2) = (/0., 1., 0./)
!HK!    norm_a(:,3) = (/0., 0., 1./)

    ! coeffs Matrix is structured as follows:
    ! --> second dimension
    ! xx xy xz
    ! yx yy yz
    ! zx zy zz


    ! Case 1:
    ! Calculation of so far needed Coefficients:
!HK!    coeffs(1,1) = dot_product(norm_a(:,1),norm_b(:,1))
!HK!    coeffs(1,2) = dot_product(norm_a(:,1),norm_b(:,2))
!HK!    coeffs(1,3) = dot_product(norm_a(:,1),norm_b(:,3))

    if ( abs(vec_ab(1)) &
      &  .gt. dim_a(1) + abs(dim_b(1)*norm_b(1,1)) &
      &                + abs(dim_b(2)*norm_b(1,2)) &
      &                + abs(dim_b(3)*norm_b(1,3))) then
      overlap = .false.
      return
    end if


    ! Case 2:
    ! Calculation of so far needed Coefficients:
!HK!    coeffs(2,1) = dot_product(norm_a(:,2),norm_b(:,1))
!HK!    coeffs(2,2) = dot_product(norm_a(:,2),norm_b(:,2))
!HK!    coeffs(2,3) = dot_product(norm_a(:,2),norm_b(:,3))

    if ( abs(vec_ab(2)) &
      &  .gt. dim_a(2) + abs(dim_b(1)*norm_b(2,1)) &
      &                + abs(dim_b(2)*norm_b(2,2)) &
      &                + abs(dim_b(3)*norm_b(2,3))) then
      overlap = .false.
      return
    end if

    ! Case 3:
    ! Calculation of so far needed Coefficients:
!HK!    coeffs(3,1) = dot_product(norm_a(:,3),norm_b(:,1))
!HK!    coeffs(3,2) = dot_product(norm_a(:,3),norm_b(:,2))
!HK!    coeffs(3,3) = dot_product(norm_a(:,3),norm_b(:,3))

    if ( abs(vec_ab(3)) &
      &  .gt. dim_a(3) + abs(dim_b(1)*norm_b(3,1)) &
      &                + abs(dim_b(2)*norm_b(3,2)) &
      &                + abs(dim_b(3)*norm_b(3,3))) then
      overlap = .false.
      return
    end if


    ! Case 4:
    if ( abs(dot_product(vec_ab, norm_b(:,1))) &
      &  .gt. dim_b(1) + abs(dim_a(1)*norm_b(1,1)) &
      &                + abs(dim_a(2)*norm_b(2,1)) &
      &                + abs(dim_a(3)*norm_b(3,1))) then
      overlap = .false.
      return
    end if

    ! Case 5:
    if ( abs(dot_product(vec_ab, norm_b(:,2))) &
      &  .gt. dim_b(2) + abs(dim_a(1)*norm_b(1,2)) &
      &                + abs(dim_a(2)*norm_b(2,2)) &
      &                + abs(dim_a(3)*norm_b(3,2))) then
      overlap = .false.
      return
    end if

    ! Case 6:
    if ( abs(dot_product(vec_ab, norm_b(:,3))) &
      &  .gt. dim_b(3) + abs(dim_a(1)*norm_b(1,3)) &
      &                + abs(dim_a(2)*norm_b(2,3)) &
      &                + abs(dim_a(3)*norm_b(3,3))) then
      overlap = .false.
      return
    end if


    ! Case 7:
    if ( abs(vec_ab(3))*norm_b(2,1) - vec_ab(2)*norm_b(3,1) &
      &  .gt. abs(dim_a(2)*norm_b(3,1)) + abs(dim_a(3)*norm_b(2,1)) &
      &       + abs(dim_b(2)*norm_b(1,3)) + abs(dim_b(3)*norm_b(1,2)) ) then
      overlap = .false.
      return
    end if

    ! Case 8:
    if ( abs(vec_ab(3))*norm_b(2,2) - vec_ab(2)*norm_b(3,2) &
      &  .gt. abs(dim_a(2)*norm_b(3,2)) + abs(dim_a(3)*norm_b(2,2)) &
      &       + abs(dim_b(1)*norm_b(1,3)) + abs(dim_b(3)*norm_b(1,1)) ) then
      overlap = .false.
      return
    end if

    ! Case 9:
    if ( abs(vec_ab(3))*norm_b(2,3) - vec_ab(2)*norm_b(3,3) &
      &  .gt. abs(dim_a(2)*norm_b(3,3)) + abs(dim_a(3)*norm_b(2,3)) &
      &       + abs(dim_b(1)*norm_b(1,2)) + abs(dim_b(2)*norm_b(1,1)) ) then
      overlap = .false.
      return
    end if

    ! Case 10:
    if ( abs(vec_ab(1))*norm_b(3,1) - vec_ab(3)*norm_b(1,1) &
      &  .gt. abs(dim_a(1)*norm_b(3,1)) + abs(dim_a(3)*norm_b(1,1)) &
      &       + abs(dim_b(2)*norm_b(2,3)) + abs(dim_b(3)*norm_b(2,2)) ) then
      overlap = .false.
      return
    end if

    ! Case 11:
    if ( abs(vec_ab(1))*norm_b(3,2) - vec_ab(3)*norm_b(1,2) &
      &  .gt. abs(dim_a(1)*norm_b(3,2)) + abs(dim_a(3)*norm_b(1,2)) &
      &       + abs(dim_b(1)*norm_b(2,3)) + abs(dim_b(3)*norm_b(2,1)) ) then
      overlap = .false.
      return
    end if

    ! Case 12:
    if ( abs(vec_ab(1))*norm_b(3,3) - vec_ab(3)*norm_b(1,3) &
      &  .gt. abs(dim_a(1)*norm_b(3,3)) + abs(dim_a(3)*norm_b(1,3)) &
      &       + abs(dim_b(1)*norm_b(2,2)) + abs(dim_b(2)*norm_b(2,1)) ) then
      overlap = .false.
      return
    end if

    ! Case 13:
    if ( abs(vec_ab(2))*norm_b(1,1) - vec_ab(1)*norm_b(2,1) &
      &  .gt. abs(dim_a(1)*norm_b(2,1)) + abs(dim_a(2)*norm_b(1,1)) &
      &       + abs(dim_b(2)*norm_b(3,3)) + abs(dim_b(3)*norm_b(3,2)) ) then
      overlap = .false.
      return
    end if

    ! Case 14:
    if ( abs(vec_ab(2))*norm_b(1,2) - vec_ab(1)*norm_b(2,2) &
      &  .gt. abs(dim_a(1)*norm_b(2,2)) + abs(dim_a(2)*norm_b(1,2)) &
      &       + abs(dim_b(1)*norm_b(3,3)) + abs(dim_b(3)*norm_b(3,1)) ) then
      overlap = .false.
      return
    end if

    ! Case 15:
    if ( abs(vec_ab(2))*norm_b(1,3) - vec_ab(1)*norm_b(2,3) &
      &  .gt. abs(dim_a(1)*norm_b(2,3)) + abs(dim_a(2)*norm_b(1,3)) &
      &       + abs(dim_b(1)*norm_b(3,2)) + abs(dim_b(2)*norm_b(3,1)) ) then
      overlap = .false.
      return
    end if

  end function boxBoxOverlap