ResPart_set_Regroup Subroutine

private impure subroutine ResPart_set_Regroup(this)

Regroup markers with their respective centroids on the same MPI block

Type Bound

ResPart_set

Arguments

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

Set of resolved particles


Calls

proc~~respart_set_regroup~~CallsGraph proc~respart_set_regroup ResPart_set%ResPart_set_Regroup markers markers proc~respart_set_regroup->markers proc~block_obj_locate block_obj%block_obj_Locate proc~respart_set_regroup->proc~block_obj_locate proc~lagrangian_set_communicate lagrangian_set%lagrangian_set_Communicate proc~respart_set_regroup->proc~lagrangian_set_communicate proc~lagrangian_set_updatecount lagrangian_set%lagrangian_set_UpdateCount proc~respart_set_regroup->proc~lagrangian_set_updatecount mpi_cart_rank mpi_cart_rank proc~block_obj_locate->mpi_cart_rank mpi_gather mpi_gather proc~lagrangian_set_communicate->mpi_gather mpi_recv mpi_recv proc~lagrangian_set_communicate->mpi_recv mpi_send mpi_send proc~lagrangian_set_communicate->mpi_send proc~lagrangian_set_recycle lagrangian_set%lagrangian_set_Recycle proc~lagrangian_set_communicate->proc~lagrangian_set_recycle proc~lagrangian_set_resize lagrangian_set%lagrangian_set_Resize proc~lagrangian_set_communicate->proc~lagrangian_set_resize mpi_allgather mpi_allgather proc~lagrangian_set_updatecount->mpi_allgather proc~lagrangian_set_recycle->proc~lagrangian_set_resize

Called by

proc~~respart_set_regroup~~CalledByGraph proc~respart_set_regroup ResPart_set%ResPart_set_Regroup proc~cdifs_obj_advancesolutionrp cdifs_obj_AdvanceSolutionRP proc~cdifs_obj_advancesolutionrp->proc~respart_set_regroup proc~respart_set_gethydroforces ResPart_set%ResPart_set_GetHydroForces proc~cdifs_obj_advancesolutionrp->proc~respart_set_gethydroforces proc~grans_obj_advancesolution grans_obj_AdvanceSolution proc~grans_obj_advancesolution->proc~respart_set_regroup proc~respart_set_gethydroforces->proc~respart_set_regroup 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 ResPart_set_Regroup(this)
      !> Regroup markers with their respective centroids on
      ! the same MPI block
      implicit none
      class(ResPart_set), intent(inout):: this                                 !! Set of resolved particles
      ! Work variables
      integer, allocatable :: buf1d(:)
      integer :: n
      integer(kind=8):: id

      ! Update this%count
      call this%UpdateCount()

      ! Allocate global arrays of centroids MPI ranks
      if (.not.allocated(this%ranks)) allocate(this%ranks(this%count))
      if (size(this%ranks).ne. this%count) then
        deallocate(this%ranks)
        allocate(this%ranks(this%count))
      end if

      allocate(buf1D(this%count))

      ! Set to zero
      this%ranks(:)=0
      ! Update Gloabl ranks
      do n=1,this%count_
        id=this%p(n)%id
        this%ranks(id)=this%block%Locate(this%p(n)%p)
      end do

      call this%parallel%Sum(this%ranks,buf1D); this%ranks=buf1D

      select type (markers=>this%ib%p)
      type is (marker_obj)
        do n=1,this%ib%count_
          ! %s stores global id of centroid, %t stores owning MPI rank
          markers(n)%t=this%ranks(markers(n)%s)
        end do
      end select
      ! Communicate with specified particle locator
      call this%ib%Communicate(ResPart_set_GetOwnerRankByRP)

      deallocate(buf1D)

      return
    end subroutine ResPart_set_Regroup