lagrangian_set_UpdateGhostObjectsDir Subroutine

private impure subroutine lagrangian_set_UpdateGhostObjectsDir(this, dist, idir)

Updates ghost objects in the idir direction. Copies objects that lie "dist"-away from the block's boundaries in idir-direction to neighboring MPI-ranks. Copied objects get a negative ID to designate them as ghost objects

Arguments

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

A set of Lagrangian objects

real(kind=wp), intent(in) :: dist

Distance from boundaries

integer, intent(in) :: idir

Direction of communication


Calls

proc~~lagrangian_set_updateghostobjectsdir~~CallsGraph proc~lagrangian_set_updateghostobjectsdir lagrangian_set_UpdateGhostObjectsDir mpi_irecv mpi_irecv proc~lagrangian_set_updateghostobjectsdir->mpi_irecv mpi_isend mpi_isend proc~lagrangian_set_updateghostobjectsdir->mpi_isend mpi_wait mpi_wait proc~lagrangian_set_updateghostobjectsdir->mpi_wait mpi_waitall mpi_waitall proc~lagrangian_set_updateghostobjectsdir->mpi_waitall proc~lagrangian_set_resize lagrangian_set%lagrangian_set_Resize proc~lagrangian_set_updateghostobjectsdir->proc~lagrangian_set_resize

Called by

proc~~lagrangian_set_updateghostobjectsdir~~CalledByGraph proc~lagrangian_set_updateghostobjectsdir lagrangian_set_UpdateGhostObjectsDir proc~lagrangian_set_updateghostobjects lagrangian_set%lagrangian_set_UpdateGhostObjects proc~lagrangian_set_updateghostobjects->proc~lagrangian_set_updateghostobjectsdir proc~collision_obj_updateghostobjects collision_obj%collision_obj_UpdateGhostObjects proc~collision_obj_updateghostobjects->proc~lagrangian_set_updateghostobjects proc~cdifs_obj_updatecollisions cdifs_obj_UpdateCollisions proc~cdifs_obj_updatecollisions->proc~collision_obj_updateghostobjects proc~grans_obj_advancesolution_computecollisionforces grans_obj_AdvanceSolution_ComputeCollisionForces proc~grans_obj_advancesolution_computecollisionforces->proc~collision_obj_updateghostobjects 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 lagrangian_set_UpdateGhostObjectsDir(this,dist,idir)
      !> Updates ghost objects in the idir direction.
      ! Copies objects that lie "dist"-away from the
      ! block's boundaries in idir-direction to
      ! neighboring MPI-ranks. Copied objects get a
      ! negative ID to designate them as ghost objects
      implicit none
      class(lagrangian_set), intent(inout) :: this                             !! A set of Lagrangian objects
      real(wp),              intent(in)    :: dist                             !! Distance from boundaries
      integer,               intent(in)    :: idir                             !! Direction of communication
      ! Work variables
      integer :: nsendL,nsendR                                                 !! Number of ghost objects to send to left/right rank
      integer :: nrecvL,nrecvR                                                 !! Number of ghost objects to receive from left/right rank
      class(lagrangian_obj),allocatable :: bufL(:),bufR(:)                     !! Buffers containing ghost objects to be sent
      type(MPI_STATUS)  :: statuses(4)                                         !! For MPI non-blocking send/receive
      type(MPI_REQUEST) :: requests(4)                                         !! For MPI non-blocking send/receive
      integer :: n                                                             !! Iterator
      integer :: nobjects_old
      integer :: ierr

      ! Figure out how many objects need to be sent
      nsendL=0
      nsendR=0
      do n=1,this%count_
        ! Determine wehther this object is close enough to the boundary
        ! to become a ghost object for a neighboring block
        if ((this%p(n)%p(idir).lt.this%block%xlo(idir)+dist) .and. this%parallel%rank%L(idir)-1.ne.MPI_PROC_NULL) nsendL=nsendL+1
        if ((this%p(n)%p(idir).gt.this%block%xhi(idir)-dist) .and. this%parallel%rank%R(idir)-1.ne.MPI_PROC_NULL) nsendR=nsendR+1
      end do

      ! Allocate buffers to send objects
      if (nsendL.gt.0) allocate(bufL(nsendL),source=this%sample)
      if (nsendR.gt.0) allocate(bufR(nsendR),source=this%sample)

      ! Copy objects in buffers and make ID negative
      ! to signify that this is a ghost object
      nsendL=0
      nsendR=0
      associate(periodicity => this%block%periods, pmin => this%block%pmin, pmax => this%block%pmax)
        do n=1,this%count_

          ! Send object left, if it's close enough to the boundary
          if ((this%p(n)%p(idir).lt.this%block%xlo(idir)+dist) .and. this%parallel%rank%L(idir)-1.ne.MPI_PROC_NULL) then
            nsendL=nsendL+1

            ! Store object in buffer
            bufL(nsendL)=this%p(n)

            ! Make id sign negative to mark this as ghost object
            bufL(nsendL)%id=-abs(bufL(nsendL)%id)
            if (periodicity(idir) .and.  &
                this%parallel%rank%dir(idir) .eq. 1)  then
                bufL(nsendL)%p(idir) = bufL(nsendL)%p(idir) + (pmax(idir)-pmin(idir))
            end if
          end if

          ! Send object right, if it's close enough to the boundary
          if ((this%p(n)%p(idir).gt.this%block%xhi(idir)-dist) .and. this%parallel%rank%R(idir)-1.ne.MPI_PROC_NULL) then
            nsendR=nsendR+1

            ! Store object in buffer
            bufR(nsendR)=this%p(n)

            ! Make id sign negative to mark this as ghost object
            bufR(nsendR)%id=-abs(bufR(nsendR)%id)

            if (periodicity(idir) .and.  &
                this%parallel%rank%dir(idir) .eq. this%parallel%np(idir))  then
                bufR(nsendR)%p(idir) = bufR(nsendR)%p(idir) - (pmax(idir)-pmin(idir))
            end if
          end if
        end do
      end associate

      ! Communicate sizes with non-blocking MPI directives
      nrecvR = 0
      nrecvL = 0
      associate(mpi=>this%parallel)
        ! Post receives from Left and right ranks
        call MPI_IRECV( nrecvL, 1, mpi%INTEGER, &
                        mpi%rank%L(idir)-1, 0, mpi%comm%g, requests(1), ierr)
        call MPI_IRECV( nrecvR, 1, mpi%INTEGER, &
                        mpi%rank%R(idir)-1, 0, mpi%comm%g, requests(2), ierr)
        ! Send sizes to left and right ranks
        call MPI_ISEND( nsendR, 1, mpi%INTEGER, &
                        mpi%rank%R(idir)-1, 0, mpi%comm%g, requests(3), ierr)
        call MPI_ISEND( nsendL, 1, mpi%INTEGER, &
                        mpi%rank%L(idir)-1, 0, mpi%comm%g, requests(4), ierr)
        ! Synchronize
        call MPI_WAITALL( 4, requests, statuses, ierr )
      end associate

      ! Store old size and resize the set
      nobjects_old = this%count_
      call this%Resize(this%count_+nrecvL+nrecvR)

      ! Communicate sizes with non-blocking MPI directives
      associate(mpi=>this%parallel)
        ! Post receives from Left and right ranks
        if (nrecvL.gt.0) then
          call MPI_IRECV( this%p(nobjects_old+1)%id,       nrecvL, this%MPI_TYPE, &
                        mpi%rank%L(idir)-1, 0, mpi%comm%g, requests(1), ierr)
        end if
        if (nrecvR.gt.0) then
          call MPI_IRECV( this%p(nobjects_old+nrecvL+1)%id,nrecvR, this%MPI_TYPE, &
                        mpi%rank%R(idir)-1, 0, mpi%comm%g, requests(2), ierr)
        end if
        ! Send to left and right ranks
        if (nsendR.gt.0) then
          call MPI_ISEND( bufR(1)%id,                      nsendR, this%MPI_TYPE, &
                        mpi%rank%R(idir)-1, 0, mpi%comm%g, requests(3), ierr)
        end if
        if (nsendL.gt.0) then
          call MPI_ISEND( bufL(1)%id,                      nsendL, this%MPI_TYPE, &
                        mpi%rank%L(idir)-1, 0, mpi%comm%g, requests(4), ierr)
        end if
        ! Synchronize
        if (nrecvL.gt.0) call MPI_WAIT(requests(1),statuses(1),ierr)
        if (nrecvR.gt.0) call MPI_WAIT(requests(2),statuses(2),ierr)
        if (nsendR.gt.0) call MPI_WAIT(requests(3),statuses(3),ierr)
        if (nsendL.gt.0) call MPI_WAIT(requests(4),statuses(4),ierr)
      end associate

      return
    end subroutine lagrangian_set_UpdateGhostObjectsDir