collision_obj_UpdateNeighborList Subroutine

private impure subroutine collision_obj_UpdateNeighborList(this)

Updates neighbor lists in preparation for collisions.

Type Bound

collision_obj

Arguments

Type IntentOptional Attributes Name
class(collision_obj), intent(inout) :: this

Collision utility


Calls

proc~~collision_obj_updateneighborlist~~CallsGraph proc~collision_obj_updateneighborlist collision_obj%collision_obj_UpdateNeighborList proc~block_obj_getowningcellwgc block_obj%block_obj_GetOwningCellWGC proc~collision_obj_updateneighborlist->proc~block_obj_getowningcellwgc proc~sllist_obj_put sllist_obj%sllist_obj_Put proc~collision_obj_updateneighborlist->proc~sllist_obj_put proc~sllist_obj_put->proc~sllist_obj_put

Called by

proc~~collision_obj_updateneighborlist~~CalledByGraph proc~collision_obj_updateneighborlist collision_obj%collision_obj_UpdateNeighborList proc~cdifs_obj_updatecollisions cdifs_obj_UpdateCollisions proc~cdifs_obj_updatecollisions->proc~collision_obj_updateneighborlist proc~grans_obj_advancesolution_computecollisionforces grans_obj_AdvanceSolution_ComputeCollisionForces proc~grans_obj_advancesolution_computecollisionforces->proc~collision_obj_updateneighborlist proc~cdifs_obj_advancesolutionrp cdifs_obj_AdvanceSolutionRP proc~cdifs_obj_advancesolutionrp->proc~cdifs_obj_updatecollisions proc~grans_obj_advancesolution grans_obj_AdvanceSolution proc~grans_obj_advancesolution->proc~grans_obj_advancesolution_computecollisionforces interface~grans_obj_advancesolution grans_obj%grans_obj_AdvanceSolution interface~grans_obj_advancesolution->proc~grans_obj_advancesolution proc~cdifs_obj_advancesolution cdifs_obj_AdvanceSolution proc~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolutionrp interface~cdifs_obj_advancesolution cdifs_obj%cdifs_obj_AdvanceSolution interface~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolution

Source Code

    impure subroutine collision_obj_UpdateNeighborList(this)
      !> Updates neighbor lists in preparation for collisions.
      implicit none
      class(collision_obj), intent(inout) :: this                              !! Collision utility
      ! Work variables
      integer :: coor(3)
      integer :: n
      integer :: count

      ! Resolved Particles
      if (associated(this%RP)) then
        if (.not.allocated(this%RPobjincell) &
            .or. (.not.allocated(this%RPneighbors))) then
            call this%parallel%Stop("Using collisions utility without proper initialization")
        end if

        this%RPobjincell = 0
        do n=1,this%RP%count_
          ! Get coordinate on collision block
          coor = this%cblock%GetOwningCellWGC(this%RP%p(n)%p)
          ! Increment counter
          count = this%RPobjincell(coor(1),coor(2),coor(3)) + 1
          this%RPobjincell(coor(1),coor(2),coor(3)) = count
          ! Add to neighbor list
          call this%RPneighbors(coor(1),coor(2),coor(3))%put(key=count,val=n)
        end do
      end if

      ! Point Particles
      if (associated(this%PP)) then
        if (.not.allocated(this%PPobjincell) &
            .or. (.not.allocated(this%PPneighbors))) then
            call this%parallel%Stop("Using collisions utility without proper initialization")
        end if

        this%PPobjincell = 0
        do n=1,this%PP%count_
          ! Get coordinate on collision block
          coor = this%cblock%GetOwningCellWGC(this%PP%p(n)%p)
          ! Increment counter
          count = this%PPobjincell(coor(1),coor(2),coor(3)) + 1
          this%PPobjincell(coor(1),coor(2),coor(3)) = count
          ! Add to neighbor list
          call this%PPneighbors(coor(1),coor(2),coor(3))%put(key=count,val=n)
        end do
      end if

      ! Immersed Boundary Points
      if (associated(this%IB)) then
        if (.not.allocated(this%IBobjincell) &
            .or. (.not.allocated(this%IBneighbors))) then
            call this%parallel%Stop("Using collisions utility without proper initialization")
        end if

        this%IBobjincell = 0
        do n=1,this%IB%count_
          ! Get coordinate on collision block
          coor = this%cblock%GetOwningCellWGC(this%IB%p(n)%p)
          ! Increment counter
          count = this%IBobjincell(coor(1),coor(2),coor(3)) + 1
          this%IBobjincell(coor(1),coor(2),coor(3)) = count
          ! Add to neighbor list
          call this%IBneighbors(coor(1),coor(2),coor(3))%put(key=count,val=n)
        end do
      end if

      return
    end subroutine collision_obj_UpdateNeighborList