FIXME : Add documentation
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(redist_type), | intent(inout) | :: | r | |||
complex, | intent(in), | dimension (r%to_low(1):, r%to_low(2):) | :: | from_here | ||
complex, | intent(inout), | dimension (r%from_low(1):, r%from_low(2):, r%from_low(3):, r%from_low(4):) | :: | to_here |
subroutine c_redist_42_inv (r, from_here, to_here)
use job_manage, only: time_message
use mp, only: iproc, nproc, send, receive, get_mp_times
implicit none
type (redist_type), intent (in out) :: r
complex, dimension (r%to_low(1):, &
r%to_low(2):), intent (in) :: from_here
complex, dimension (r%from_low(1):, &
r%from_low(2):, &
r%from_low(3):, &
r%from_low(4):), intent (in out) :: to_here
integer :: i, idp, ipto, ipfrom, iadp
real :: mp_total, mp_total_after
if (.not. using_measure_scatter) then
call time_message(.false.,time_redist,' Redistribute')
call get_mp_times(total_time = mp_total)
end if
! redistribute from local processor to local processor
do i = 1, r%to(iproc)%nn
to_here(r%from(iproc)%k(i), &
r%from(iproc)%l(i), &
r%from(iproc)%m(i), &
r%from(iproc)%n(i)) &
= from_here(r%to(iproc)%k(i), &
r%to(iproc)%l(i))
end do
! 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%to(ipto)%nn > 0) then
do i = 1, r%to(ipto)%nn
r%complex_buff(i) = from_here(r%to(ipto)%k(i), &
r%to(ipto)%l(i))
end do
call send (r%complex_buff(1:r%to(ipto)%nn), ipto, idp)
end if
! receive from idpth preceding processor
if (r%from(ipfrom)%nn > 0) then
call receive (r%complex_buff(1:r%from(ipfrom)%nn), ipfrom, idp)
do i = 1, r%from(ipfrom)%nn
to_here(r%from(ipfrom)%k(i), &
r%from(ipfrom)%l(i), &
r%from(ipfrom)%m(i), &
r%from(ipfrom)%n(i)) &
= r%complex_buff(i)
end do
end if
else
! receive from idpth preceding processor
if (r%from(ipfrom)%nn > 0) then
call receive (r%complex_buff(1:r%from(ipfrom)%nn), ipfrom, idp)
do i = 1, r%from(ipfrom)%nn
to_here(r%from(ipfrom)%k(i), &
r%from(ipfrom)%l(i), &
r%from(ipfrom)%m(i), &
r%from(ipfrom)%n(i)) &
= r%complex_buff(i)
end do
end if
! send to idpth next processor
if (r%to(ipto)%nn > 0) then
do i = 1, r%to(ipto)%nn
r%complex_buff(i) = from_here(r%to(ipto)%k(i), &
r%to(ipto)%l(i))
end do
call send (r%complex_buff(1:r%to(ipto)%nn), ipto, idp)
end if
end if
end do
if (.not. using_measure_scatter) then
call time_message(.false.,time_redist,' Redistribute')
call get_mp_times(total_time = mp_total_after)
time_redist_mpi = time_redist_mpi + (mp_total_after - mp_total)
end if
end subroutine c_redist_42_inv