collision_obj_ComputeCollisionsRPvIB Subroutine

private pure subroutine collision_obj_ComputeCollisionsRPvIB(this)

Computes collisions between Resolved Particles and Immersed Boundaries.

Type Bound

collision_obj

Arguments

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

Collision utility


Calls

proc~~collision_obj_computecollisionsrpvib~~CallsGraph proc~collision_obj_computecollisionsrpvib collision_obj%collision_obj_ComputeCollisionsRPvIB center center proc~collision_obj_computecollisionsrpvib->center none~get~2 sllist_obj%Get proc~collision_obj_computecollisionsrpvib->none~get~2 proc~cross_product~3 cross_product proc~collision_obj_computecollisionsrpvib->proc~cross_product~3 proc~dem_col DEM_col proc~collision_obj_computecollisionsrpvib->proc~dem_col triangle triangle proc~collision_obj_computecollisionsrpvib->triangle 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_computecollisionsrpvib~~CalledByGraph proc~collision_obj_computecollisionsrpvib collision_obj%collision_obj_ComputeCollisionsRPvIB proc~collision_obj_computecollisions collision_obj%collision_obj_ComputeCollisions proc~collision_obj_computecollisions->proc~collision_obj_computecollisionsrpvib 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_ComputeCollisionsRPvIB(this)
      !> Computes collisions between Resolved Particles and Immersed Boundaries.
      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):: Fc(3)
      real(wp):: Tcol(3)
      integer :: ncol
      real(wp), parameter :: Pi=4.0_wp*atan(1.0_wp)

      associate (lo =>this%cblock%lo, hi=>this%cblock%hi)
        select type(center=>this%RP%p)
        type is (ResPart_obj)
          select type(triangle=>this%IB%p)
          type is (marker_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 RP 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)

                    ! Ignore Ghost Particles
                    if (center(id1)%id.le.0) cycle

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

                    ! Check collisions with neighbor IB triangles
                    Fcol = 0.0_wp
                    Tcol = 0.0_wp
                    ncol = 0
                    do k = kc-1,kc+1
                      do j = jc-1,jc+1
                        do i = ic-1,ic+1
                          do key2 = 1,this%IBobjincell(i,j,k)
                            call this%IBneighbors(i,j,k)%get(key=key2,val=id2)

                            ! IB info (static IB for now)
                            ! r2 is the projection of r1 on the plane going through
                            ! triangle(id2)%p and oriented along tirange(id2)%n
                            r2 = r1 - dot_product((r1-triangle(id2)%p),triangle(id2)%n)*triangle(id2)%n
                            d2 = 0.0_wp
                            m2 = m1*1e6_wp
                            v2 = 0.0_wp
                            w2 = 0.0_wp

                            ! Force due to this IB triangle
                            Fc = DEM_col(r1,r2,d1,d2,m1,m2,v1,v2,w1,w2,this%tcol,this%edry,this%muc)

                            if (norm2(Fc).ne.0.0_wp) then
                              ! Force
                              Fcol = Fcol + Fc
                              ! Torque
                              n12 = (r2-r1)/(norm2(r2-r1)+epsilon(1.0_wp))
                              Tcol = Tcol + cross_product (0.5_wp*d1*n12,Fc)

                              ncol = ncol + 1
                            end if
                          end do
                        end do
                      end do
                    end do

                    ncol = max(ncol,1)

                    Fcol = Fcol/real(ncol,wp)
                    Tcol = Tcol/real(ncol,wp)

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

                end do
              end do
            end do
          end select
        end select
      end associate

      return
    end subroutine collision_obj_ComputeCollisionsRPvIB