subroutine MAPL_LocStreamCreateXform ( Xform, LocStreamOut, LocStreamIn, NAME, MASK_OUT, &
UseFCollect, RC )
type(MAPL_LocStreamXform), intent(OUT) :: Xform
type(MAPL_LocStream), intent(IN ) :: LocStreamOut
type(MAPL_LocStream), intent(IN ) :: LocStreamIn
character(len=*), intent(IN ) :: NAME
logical, optional, intent(IN ) :: MASK_OUT(:)
logical, optional, intent(IN ) :: UseFCollect
integer, optional, intent(OUT) :: RC
! Local variables
integer :: STATUS
integer :: N, M, MM
logical :: DONE(LocStreamOut%Ptr%NT_local)
logical, pointer :: ISDONE(:)
logical :: dn(1)
type (ESMF_VM) :: vm
integer :: NDES, hash
integer, pointer :: GLOBAL_IdByPe(:) =>null() ! All Location Ids in PE order
integer :: I, First, Last
logical, allocatable :: IsNeeded(:)
integer, allocatable :: PELens(:), Begs(:), Ends(:)
integer :: NumSenders
integer :: myId
integer :: MyLen(1)
#if defined(ONE_SIDED_COMM)
integer :: SizeOfReal
#endif
! Both streams must be subsets of same parent.
! The parent stream is usually an exchange grid.
!-----------------------------------------------
_ASSERT(trim(LocStreamOut%PTR%ROOTNAME)==trim(LocStreamIn%PTR%ROOTNAME),'needs informative message')
allocate(XFORM%Ptr, STAT=STATUS)
_VERIFY(STATUS)
Xform%Ptr%InputStream = LocStreamIn
Xform%Ptr%OutputStream = LocStreamOut
Xform%Ptr%Name = NAME
Xform%Ptr%do_not_use_fcollect = .false. ! defaults to FCOLLECT for now
#ifdef DO_NOT_USE_FCOLLECT
Xform%Ptr%do_not_use_fcollect = .true.
#endif
if (present(UseFCollect)) then
Xform%Ptr%do_not_use_fcollect = .not. UseFCollect
end if
! We have to fill all output locations where mask is true
!--------------------------------------------------------
if(present(MASK_OUT)) then
DONE = .not. MASK_OUT
else
DONE = .false.
end if
Xform%Ptr%count = count(.not.DONE)
ALLOCATE(Xform%Ptr%IndexOut(Xform%Ptr%count), stat=STATUS)
_VERIFY(STATUS)
ALLOCATE(Xform%Ptr%IndexIn (Xform%Ptr%count), stat=STATUS)
_VERIFY(STATUS)
MM=1
Hash = MAPL_HashCreate(8*1024)
do M = 1, LocStreamIn%Ptr%NT_local
n = MAPL_HashIncrement(Hash,LocStreamIn%Ptr%Local_Id(M))
_ASSERT(N==M,'needs informative message')
enddo
do N = 1, LocStreamOut%Ptr%NT_local
if(DONE(N)) cycle
M = MAPL_HashIncrement(Hash,LocStreamOut%Ptr%Local_Id(N))
if(m<=LocStreamIn%Ptr%NT_local) then
Xform%Ptr%IndexOut(MM) = N
Xform%Ptr%IndexIn (MM) = M
DONE (N) = .TRUE.
MM=MM+1
endif
end do
call MAPL_HashDestroy(Hash)
Xform%Ptr%LastLocal = MM-1
! Otherwise, assume nothing and do a full collect.
!-------------------------------------------------
call ESMF_VMGetCurrent ( vm, rc=status )
_VERIFY(STATUS)
call ESMF_VMGet ( vm, petCount=nDEs, &
mpiCommunicator=Xform%Ptr%Comm, rc=status )
_VERIFY(STATUS)
allocate(IsDone(NDES))
dn(1) = all(done)
call MAPL_CommsAllGather(vm, dn, 1, &
isdone, 1, rc=status)
_VERIFY(STATUS)
Xform%Ptr%Local = all(isdone)
deallocate(IsDone)
NEED_COMM: if(.not.Xform%Ptr%Local) then
if (Xform%Ptr%do_not_use_fcollect) then
allocate(PELens(NDES),Begs(NDES),Ends(NDES),IsNeeded(Ndes))
#if defined(ONE_SIDED_COMM)
allocate(Xform%Ptr%Buff(LocStreamIn%Ptr%NT_LOCAL))
allocate(Xform%Ptr%Len(NDES), stat=status)
_VERIFY(STATUS)
CALL MPI_TYPE_GET_EXTENT(MPI_REAL, lb, SizeOfReal, status)
_VERIFY(STATUS)
call mpi_Win_Create(Xform%Ptr%Buff,LocStreamIn%Ptr%NT_LOCAL*SizeOfReal, &
SizeOfReal,MPI_INFO_NULL,Xform%Ptr%Comm,Xform%Ptr%Window,status)
_VERIFY(STATUS)
#endif
MyLen(1) = LocStreamIn%Ptr%NT_LOCAL
call MAPL_CommsAllGather(vm, MyLen, 1, &
PELens, 1, rc=status)
_VERIFY(STATUS)
Begs(1) = 1
Ends(1) = PELens(1)
do i=2,NDES
Begs(i) = Ends(i-1) + 1
Ends(i) = Ends(i-1) + PELens(i)
end do
_ASSERT(Ends(NDES) == LocStreamIn%Ptr%NT_GLOBAL,'needs informative message')
endif
allocate(GLOBAL_IdByPe(LocStreamIn%Ptr%NT_GLOBAL), STAT=STATUS)
_VERIFY(STATUS)
! Collect all tile ides in the input stream's pe order
!-----------------------------------------------------
call ESMFL_FCOLLECT(LocStreamIn%Ptr%TILEGRID, GLOBAL_IdByPe, &
LocStreamIn%Ptr%LOCAL_ID, RC=STATUS)
_VERIFY(STATUS)
! Make a Hash of the tile locations by input order
!-------------------------------------------------
Hash = MAPL_HashCreate(8*1024)
do M = 1, LocStreamIn%Ptr%NT_global
n = MAPL_HashIncrement(Hash,Global_IdByPe(M))
_ASSERT(N==M,'needs informative message')
enddo
if(Xform%Ptr%do_not_use_fcollect) then
! Find out which processors have output tiles we need
!----------------------------------------------------
IsNeeded = .false.
do N = 1, LocStreamOut%Ptr%NT_local
if(.not.DONE(N)) then
M = MAPL_HashIncrement(Hash,LocStreamOut%Ptr%Local_Id(N))
do i=1,ndes
if(M>=Begs(i) .and. M<=Ends(i)) then
IsNeeded(i) = .true.
exit
end if
enddo
end if
end do
! Allocate my senders and their size in the Xform.
! Note that fullinput has all tiles from all of those
! pes that have tiles we need.
!-----------------------------------------------------
NumSenders = count(IsNeeded)
allocate(Xform%Ptr%senders(NumSenders), stat=status)
_VERIFY(STATUS)
allocate(Xform%Ptr% len(NumSenders), stat=status)
_VERIFY(STATUS)
First = 1
Last = 0
M = 0
do I=1,NDES
if(Isneeded(i)) then
m = m + 1
Last = Last + PELens(i)
Global_IdByPe(First:Last) &
= Global_IdByPe(Begs(i):Ends(i))
First = First + PELens(i)
Xform%Ptr%senders(m) = i-1
Xform%Ptr% len(m) = PELens(i)
end if
end do
deallocate(PELens,Begs,Ends,IsNeeded)
Xform%Ptr%InputLen = Last
call ESMF_VmGet(VM, localPet=MYID, rc=status)
_VERIFY(STATUS)
Xform%Ptr%myId = myid
#if defined(TWO_SIDED_COMM)
block
integer, allocatable :: allSenders(:,:)
integer :: NumReceivers
integer :: k
#if 0
allocate(allSenders(ndes,ndes), stat=status)
_VERIFY(STATUS)
allSenders(:,myId+1) = -1
if (m>0) allSenders(1:M,myId+1) = Xform%Ptr%senders
do I=1,NDES
call MAPL_CommsBcast(vm, DATA=allSenders(:,I), N=ndes, ROOT=I-1, RC=status)
_VERIFY(STATUS)
end do
NumReceivers = count(allSenders == myId)
allocate(Xform%Ptr%receivers(NumReceivers), stat=status)
_VERIFY(STATUS)
M = 0
do I=1,NDES
if(myId == I-1) cycle ! skip myself
do K=1,NDES
if(allSenders(K,I) < 0) exit !senders are packed. we have reached end ...
if(allSenders(K,I) == myId) then
M = M+1
Xform%Ptr%receivers(m) = i-1
exit
end if
end do
end do
_ASSERT(NumReceivers==M,'needs informative message')
deallocate(allSenders)
#else
allocate(allSenders(ndes,1), stat=status)
_VERIFY(STATUS)
block
integer :: lNumReceivers
do I=1,NDES
lNumReceivers = 0
if (m>0) lNumReceivers = count(Xform%Ptr%senders == I-1)
call MPI_GATHER( lNumReceivers, 1, MPI_INTEGER, &
allSenders(:,1), 1, MPI_INTEGER, &
I-1, Xform%Ptr%Comm, status )
_VERIFY(status)
enddo
end block
call ESMF_VMBarrier(vm, rc=status)
_VERIFY(STATUS)
NumReceivers = 0
do I=1,NDES
NumReceivers = NumReceivers + allSenders(I,1)
end do
allocate(Xform%Ptr%receivers(NumReceivers), stat=status)
_VERIFY(STATUS)
M = 0
do I=1,NDES
if(myId == I-1) cycle ! skip myself
do K=1,allSenders(I,1)
M = M+1
Xform%Ptr%receivers(M) = I-1
end do
end do
_ASSERT(NumReceivers==M,'needs informative message')
deallocate(allSenders)
#endif
end block
#endif
! Put the tiles we being brought over into a hash table
!------------------------------------------------------
call MAPL_HashDestroy(Hash)
Hash = MAPL_HashCreate(8*1024)
do M = 1, Xform%Ptr%InputLen
n = MAPL_HashIncrement(Hash,Global_IdByPe(M))
_ASSERT(N==M,'needs informative message')
enddo
! Pick out the ones we need fromthose brought over
!-------------------------------------------------
endif
do N = 1, LocStreamOut%Ptr%NT_local
if(.not.DONE(N)) then
M = MAPL_HashIncrement(Hash,LocStreamOut%Ptr%Local_Id(N))
Xform%Ptr%IndexOut(MM) = N
Xform%Ptr%IndexIn (MM) = M
DONE (N) = .TRUE.
MM=MM+1
end if
end do
call MAPL_HashDestroy(Hash)
deallocate(GLOBAL_IdByPe)
end if NEED_COMM
! Make sure that did it
!----------------------
_ASSERT(all(DONE),'needs informative message')
_RETURN(ESMF_SUCCESS)
end subroutine MAPL_LocStreamCreateXform