c_redist_36_mpi_copy Subroutine

private subroutine c_redist_36_mpi_copy(r, from_here, to_here, ntgrid)

Uses

FIXME : Add documentation

Arguments

Type IntentOptional Attributes Name
type(redist_type), intent(inout) :: r
complex, intent(in), dimension (r%from_low(1):, r%from_low(2):, r%from_low(3):) :: from_here
complex, intent(inout), dimension (r%to_low(1):, r%to_low(2):, r%to_low(3):, r%to_low(4):, r%to_low(5):, r%to_low(6):) :: to_here
integer, intent(in) :: ntgrid

Contents

Source Code


Source Code

  subroutine c_redist_36_mpi_copy(r, from_here, to_here, ntgrid)
    use mp, only: iproc, nproc, send, receive
    implicit none
    integer, intent(in) :: ntgrid

    type (redist_type), intent (in out) :: r

    complex, dimension (r%from_low(1):, &
                        r%from_low(2):, &
                        r%from_low(3):), intent (in) :: from_here

    complex, dimension (r%to_low(1):, &
                        r%to_low(2):, &
                        r%to_low(3):, &
                        r%to_low(4):, &
                        r%to_low(5):, &
                        r%to_low(6):), intent (in out) :: to_here

    integer :: i, idp, ipto, ipfrom, iadp
    !AJ Added as the r%complex_buf is not designed for sending contigous blocks, rather to pack individual elements
    !AJ into the buffer.  As each element in my redist object is actually (2*ntgrid+1)*2 (sigma) then I need a bigger
    !AJ send and receive buffer.
    complex, dimension(:,:,:), allocatable :: temp_buff

    allocate(temp_buff(-ntgrid:ntgrid, 1:2, size(r%complex_buff)))

    ! redistribute to idpth next processor from idpth preceding processor
    ! or redistribute from idpth preceding processor to idpth next processor
    ! to avoid deadlocks
    do idp = 1, nproc-1
       ipto = mod(iproc + idp, nproc)
       ipfrom = mod(iproc + nproc - idp, nproc)
       iadp = min(idp, nproc - idp)
       ! avoid deadlock AND ensure mostly parallel resolution
       if (mod(iproc/iadp,2) == 0) then
          ! send to idpth next processor
          if (r%from(ipto)%nn > 0) then
             do i = 1, r%from(ipto)%nn 
                temp_buff(:,:,i) = from_here(:,:,r%from(ipto)%m(i))
             end do
             call send (temp_buff(:,:,1:r%from(ipto)%nn), ipto, idp)
          end if

          !AJ FIX TO DO CONTIGIOUS BLOCKS
          ! receive from idpth preceding processor
          if (r%to(ipfrom)%nn > 0) then
             call receive (temp_buff(:,:,1:r%to(ipfrom)%nn), ipfrom, idp)
             do i = 1, r%to(ipfrom)%nn
                to_here(:,:,r%to(ipfrom)%m(i), &
                            r%to(ipfrom)%n(i), &
                            r%to(ipfrom)%o(i), &
                            r%to(ipfrom)%p(i)) &
                            = temp_buff(:,:,i)
             end do
          end if
       else
          ! receive from idpth preceding processor
          if (r%to(ipfrom)%nn > 0) then
             call receive (temp_buff(:,:,1:r%to(ipfrom)%nn), ipfrom, idp)
             do i = 1, r%to(ipfrom)%nn
                to_here(:,:,r%to(ipfrom)%m(i), &
                            r%to(ipfrom)%n(i), &
                            r%to(ipfrom)%o(i), &
                            r%to(ipfrom)%p(i)) &
                            = temp_buff(:,:,i)
             end do
          end if
          ! send to idpth next processor
          if (r%from(ipto)%nn > 0) then
             do i = 1, r%from(ipto)%nn
                temp_buff(:,:,i) = from_here(:,:,r%from(ipto)%m(i))
             end do
             call send (temp_buff(:,:,1:r%from(ipto)%nn), ipto, idp)
          end if
       end if
    end do

    deallocate(temp_buff)

  end subroutine c_redist_36_mpi_copy