collision_obj_ComputeCollisionsRPvPP Subroutine

private pure subroutine collision_obj_ComputeCollisionsRPvPP(this)

Computes collisions between Point and Resolved Particles.

Type Bound

collision_obj

Arguments

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

Collision utility


Calls

proc~~collision_obj_computecollisionsrpvpp~~CallsGraph proc~collision_obj_computecollisionsrpvpp collision_obj%collision_obj_ComputeCollisionsRPvPP none~get~2 sllist_obj%Get proc~collision_obj_computecollisionsrpvpp->none~get~2 proc~cross_product~3 cross_product proc~collision_obj_computecollisionsrpvpp->proc~cross_product~3 proc~dem_col DEM_col proc~collision_obj_computecollisionsrpvpp->proc~dem_col type1 type1 proc~collision_obj_computecollisionsrpvpp->type1 type2 type2 proc~collision_obj_computecollisionsrpvpp->type2 proc~sllist_obj_get_int4 sllist_obj%sllist_obj_Get_int4 none~get~2->proc~sllist_obj_get_int4 proc~sllist_obj_get_int8 sllist_obj%sllist_obj_Get_int8 none~get~2->proc~sllist_obj_get_int8 proc~sllist_obj_get_real_dp sllist_obj%sllist_obj_Get_real_dp none~get~2->proc~sllist_obj_get_real_dp proc~sllist_obj_get_real_sp sllist_obj%sllist_obj_Get_real_sp none~get~2->proc~sllist_obj_get_real_sp proc~dem_col->proc~cross_product~3 proc~sllist_obj_get_int4->proc~sllist_obj_get_int4 proc~sllist_obj_get_int8->proc~sllist_obj_get_int8 proc~sllist_obj_get_real_dp->proc~sllist_obj_get_real_dp proc~sllist_obj_get_real_sp->proc~sllist_obj_get_real_sp

Called by

proc~~collision_obj_computecollisionsrpvpp~~CalledByGraph proc~collision_obj_computecollisionsrpvpp collision_obj%collision_obj_ComputeCollisionsRPvPP proc~collision_obj_computecollisions collision_obj%collision_obj_ComputeCollisions proc~collision_obj_computecollisions->proc~collision_obj_computecollisionsrpvpp proc~cdifs_obj_updatecollisions cdifs_obj_UpdateCollisions proc~cdifs_obj_updatecollisions->proc~collision_obj_computecollisions proc~grans_obj_advancesolution_computecollisionforces grans_obj_AdvanceSolution_ComputeCollisionForces proc~grans_obj_advancesolution_computecollisionforces->proc~collision_obj_computecollisions 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

    pure subroutine collision_obj_ComputeCollisionsRPvPP(this)
      !> Computes collisions between Point and Resolved Particles.
      implicit none
      class(collision_obj), intent(inout) :: this                              !! Collision utility
      ! Work variable
      integer :: ic,jc,kc
      integer :: i,j,k
      integer :: key1,key2
      integer :: id1,id2
      real(wp):: r1(3), r2(3)
      real(wp):: d1   , d2
      real(wp):: m1   , m2
      real(wp):: v1(3), v2(3)
      real(wp):: w1(3), w2(3)
      real(wp):: n12(3)
      real(wp):: Fcol(3)
      real(wp):: Tcol(3)
      real(wp), parameter :: Pi=4.0_wp*atan(1.0_wp)


      associate (lo =>this%cblock%lo, hi=>this%cblock%hi)
        ! Loop over bins
        select type(type1=>this%RP%p)
        type is(ResPart_obj)
          select type (type2=>this%PP%p)
          type is (particle_obj)

            ! Loop over collision grid
            do kc=lo(3),hi(3)
              do jc=lo(2),hi(2)
                do ic=lo(1),hi(1)

                  ! Loop over objects in this cell
                  do key1 = 1,this%RPobjincell(ic,jc,kc)
                    ! Get particle ID
                    call this%RPneighbors(ic,jc,kc)%get(key=key1,val=id1)

                    r1 = type1(id1)%p
                    d1 = type1(id1)%d
                    m1 = type1(id1)%rho*Pi/6.0_wp*type1(id1)%d**3
                    v1 = type1(id1)%v
                    w1 = type1(id1)%w

                    ! Check collisions with neighbors
                    do k = kc-1,kc+1
                      do j = jc-1,jc+1
                        do i = ic-1,ic+1
                          do key2 = 1,this%PPobjincell(i,j,k)
                            call this%PPneighbors(i,j,k)%get(key=key2,val=id2)

                            r2 = type2(id2)%p
                            d2 = type2(id2)%d
                            m2 = type2(id2)%rho*Pi/6.0_wp*type2(id2)%d**3
                            v2 = type2(id2)%v
                            w2 = type2(id2)%w

                            ! Force
                            Fcol = DEM_col(r1,r2,d1,d2,m1,m2,v1,v2,w1,w2,this%tcol,this%edry,this%muc)
                            ! Torque
                            n12 = (r2-r1)/(norm2(r2-r1)+epsilon(1.0_wp))
                            Tcol = cross_product (0.5_wp*d1*n12,Fcol)

                            type1(id1)%Fc = type1(id1)%Fc + Fcol
                            type1(id1)%Tc = type1(id1)%Tc + Tcol
                          end do
                        end do
                      end do
                    end do
                  end do

                  ! Loop over objects in this cell
                  do key1 = 1,this%PPobjincell(ic,jc,kc)
                    ! Get particle ID
                    call this%PPneighbors(ic,jc,kc)%get(key=key1,val=id1)

                    r1 = type2(id1)%p
                    d1 = type2(id1)%d
                    m1 = type2(id1)%rho*Pi/6.0_wp*type2(id1)%d**3
                    v1 = type2(id1)%v
                    w1 = type2(id1)%w

                    ! Check collisions with neighbors
                    do k = kc-1,kc+1
                      do j = jc-1,jc+1
                        do i = ic-1,ic+1
                          do key2 = 1,this%RPobjincell(i,j,k)
                            call this%RPneighbors(i,j,k)%get(key=key2,val=id2)

                            r2 = type1(id2)%p
                            d2 = type1(id2)%d
                            m2 = type1(id2)%rho*Pi/6.0_wp*type1(id2)%d**3
                            v2 = type1(id2)%v
                            w2 = type1(id2)%w

                            ! Force
                            Fcol = DEM_col(r1,r2,d1,d2,m1,m2,v1,v2,w1,w2,this%tcol,this%edry,this%muc)
                            ! Torque
                            n12 = (r2-r1)/(norm2(r2-r1)+epsilon(1.0_wp))
                            Tcol = cross_product (0.5_wp*d1*n12,Fcol)

                            type2(id1)%Fc = type2(id1)%Fc + Fcol
                            type2(id1)%Tc = type2(id1)%Tc + Tcol
                          end do
                        end do
                      end do
                    end do

                  end do
                end do
              end do
            end do

          end select
        end select
      end associate

      return
    end subroutine collision_obj_ComputeCollisionsRPvPP