m_Transfer.F90 24 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818
  1. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  2. ! Math and Computer Science Division, Argonne National Laboratory !
  3. !-----------------------------------------------------------------------
  4. ! CVS m_Transfer.F90,v 1.13 2008-01-28 06:48:03 jacob Exp
  5. ! CVS MCT_2_8_0
  6. !BOP -------------------------------------------------------------------
  7. !
  8. ! !MODULE: m_Transfer - Routines for the MxN transfer of Attribute Vectors
  9. !
  10. ! !DESCRIPTION:
  11. ! This module provides routines for doing MxN transfer of data in an
  12. ! Attribute Vector between two components on separate sets of MPI processes.
  13. ! Uses the Router datatype.
  14. !
  15. ! !SEE ALSO:
  16. ! m_Rearranger
  17. ! !INTERFACE:
  18. module m_Transfer
  19. ! !USES:
  20. use m_MCTWorld, only : MCTWorld
  21. use m_MCTWorld, only : ThisMCTWorld
  22. use m_AttrVect, only : AttrVect
  23. use m_AttrVect, only : nIAttr,nRAttr
  24. use m_AttrVect, only : Permute, Unpermute
  25. use m_AttrVect, only : AttrVect_init => init
  26. use m_AttrVect, only : AttrVect_copy => copy
  27. use m_AttrVect, only : AttrVect_clean => clean
  28. use m_AttrVect, only : lsize
  29. use m_Router, only : Router
  30. use m_mpif90
  31. use m_die
  32. use m_stdio
  33. implicit none
  34. private ! except
  35. ! !PUBLIC MEMBER FUNCTIONS:
  36. public :: isend
  37. public :: send
  38. public :: waitsend
  39. public :: irecv
  40. public :: recv
  41. public :: waitrecv
  42. interface isend ; module procedure isend_ ; end interface
  43. interface send ; module procedure send_ ; end interface
  44. interface waitsend ; module procedure waitsend_ ; end interface
  45. interface irecv ; module procedure irecv_ ; end interface
  46. interface recv ; module procedure recv_ ; end interface
  47. interface waitrecv ; module procedure waitrecv_ ; end interface
  48. ! !DEFINED PARAMETERS:
  49. integer,parameter :: DefaultTag = 600
  50. ! !REVISION HISTORY:
  51. ! 08Nov02 - R. Jacob <jacob@mcs.anl.gov> - make new module by combining
  52. ! MCT_Send, MCT_Recv and MCT_Recvsum
  53. ! 11Nov02 - R. Jacob <jacob@mcs.anl.gov> - Remove MCT_Recvsum and use
  54. ! optional argument in recv_ to do the same thing.
  55. ! 23Jul03 - R. Jacob <jacob@mcs.anl.gov> - Move buffers for data and
  56. ! MPI_Reqest and MPI_Status arrays to Router. Use them.
  57. ! 24Jul03 - R. Jacob <jacob@mcs.anl.gov> - Split send_ into isend_ and
  58. ! waitsend_. Redefine send_.
  59. ! 22Jan08 - R. Jacob <jacob@mcs.anl.gov> - Handle unordered GSMaps
  60. !EOP ___________________________________________________________________
  61. character(len=*),parameter :: myname='MCT::m_Transfer'
  62. contains
  63. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  64. ! Math and Computer Science Division, Argonne National Laboratory !
  65. !BOP -------------------------------------------------------------------
  66. !
  67. ! !IROUTINE: isend_ - Distributed non-blocking send of an Attribute Vector
  68. !
  69. ! !DESCRIPTION:
  70. ! Send the the data in the {\tt AttrVect} {\tt aV} to the
  71. ! component specified in the {\tt Router} {\tt Rout}. An error will
  72. ! result if the size of the attribute vector does not match the size
  73. ! parameter stored in the {\tt Router}.
  74. !
  75. ! Requires a corresponding {\tt recv\_} or {\tt irecv\_} to be called on the other component.
  76. !
  77. ! The optional argument {\tt Tag} can be used to set the tag value used in
  78. ! the data transfer. DefaultTag will be used otherwise. {\tt Tag} must be
  79. ! the same in the matching {\tt recv\_} or {\tt irecv\_}.
  80. !
  81. ! {\bf N.B.:} The {\tt AttrVect} argument in the corresponding
  82. ! {\tt recv\_} call is assumed to have exactly the same attributes
  83. ! in exactly the same order as {\tt aV}.
  84. !
  85. ! !INTERFACE:
  86. subroutine isend_(aVin, Rout, Tag)
  87. !
  88. ! !USES:
  89. !
  90. implicit none
  91. ! !INPUT PARAMETERS:
  92. !
  93. Type(AttrVect),target,intent(in) :: aVin
  94. Type(Router), intent(inout) :: Rout
  95. integer,optional, intent(in) :: Tag
  96. ! !REVISION HISTORY:
  97. ! 07Feb01 - R. Jacob <jacob@mcs.anl.gov> - initial prototype
  98. ! 08Feb01 - R. Jacob <jacob@mcs.anl.gov> - First working code
  99. ! 18May01 - R. Jacob <jacob@mcs.anl.gov> - use MP_Type to determine type in mpi_send
  100. ! 07Jun01 - R. Jacob <jacob@mcs.anl.gov> - remove logic to check "direction" of Router.
  101. ! remove references to ThisMCTWorld%mylrank
  102. ! 03Aug01 - E. Ong <eong@mcs.anl.gov> - Explicitly specify the starting address in mpi_send.
  103. ! 15Feb02 - R. Jacob <jacob@mcs.anl.gov> - Use MCT_comm
  104. ! 26Mar02 - E. Ong <eong@mcs.anl.gov> - Apply faster copy order
  105. ! 26Sep02 - R. Jacob <jacob@mcs.anl.gov> - Check Av against Router lAvsize
  106. ! 05Nov02 - R. Jacob <jacob@mcs.anl.gov> - Remove iList, rList arguments.
  107. ! 08Nov02 - R. Jacob <jacob@mcs.anl.gov> - MCT_Send is now send_ in m_Transfer
  108. ! 11Nov02 - R. Jacob <jacob@mcs.anl.gov> - Use DefaultTag and add optional Tag argument
  109. ! 25Jul03 - R. Jacob <jacob@mcs.anl.gov> - Split into isend_ and waitsend_
  110. ! 22Jan08 - R. Jacob <jacob@mcs.anl.gov> - Handle unordered GSMaps by permuting before send.
  111. ! remove special case for sending one segment directly from Av which probably
  112. ! wasn't safe.
  113. !EOP ___________________________________________________________________
  114. character(len=*),parameter :: myname_=myname//'::isend_'
  115. integer :: numi,numr,i,j,k,ier
  116. integer :: mycomp,othercomp
  117. integer :: AttrIndex,VectIndex,seg_start,seg_end
  118. integer :: proc,nseg,mytag
  119. integer :: mp_Type_rp1
  120. logical :: unordered
  121. type(AttrVect),pointer :: Av
  122. type(AttrVect),target :: Avtmp
  123. !--------------------------------------------------------
  124. ! Return if no one to send to
  125. if(Rout%nprocs .eq. 0 ) RETURN
  126. ! set up Av to send from
  127. unordered = associated(Rout%permarr)
  128. if (unordered) then
  129. call AttrVect_init(Avtmp,Avin,lsize(Avin))
  130. call AttrVect_copy(Avin,aVtmp)
  131. call Permute(aVtmp,Rout%permarr)
  132. Av => Avtmp
  133. else
  134. Av => Avin
  135. endif
  136. !check Av size against Router
  137. !
  138. if(lsize(aV) /= Rout%lAvsize) then
  139. write(stderr,'(2a)') myname_, &
  140. ' MCTERROR: AV size not appropriate for this Router...exiting'
  141. call die(myname_)
  142. endif
  143. ! get ids of components involved in this communication
  144. mycomp=Rout%comp1id
  145. othercomp=Rout%comp2id
  146. ! find total number of real and integer vectors
  147. ! for now, assume we are sending all of them
  148. Rout%numiatt = nIAttr(aV)
  149. Rout%numratt = nRAttr(aV)
  150. numi = Rout%numiatt
  151. numr = Rout%numratt
  152. !!!!!!!!!!!!!! IF SENDING INTEGER DATA
  153. if(numi .ge. 1) then
  154. ! allocate buffers to hold all outgoing data
  155. do proc=1,Rout%nprocs
  156. allocate(Rout%ip1(proc)%pi(Rout%locsize(proc)*numi),stat=ier)
  157. if(ier/=0) call die(myname_,'allocate(Rout%ip1%pi)',ier)
  158. enddo
  159. endif
  160. !!!!!!!!!!!!!! IF SENDING REAL DATA
  161. if(numr .ge. 1) then
  162. ! allocate buffers to hold all outgoing data
  163. do proc=1,Rout%nprocs
  164. allocate(Rout%rp1(proc)%pr(Rout%locsize(proc)*numr),stat=ier)
  165. if(ier/=0) call die(myname_,'allocate(Rout%rp1%pr)',ier)
  166. enddo
  167. mp_Type_rp1=MP_Type(Rout%rp1(1)%pr(1))
  168. endif
  169. ! Load data going to each processor
  170. do proc = 1,Rout%nprocs
  171. j=1
  172. k=1
  173. ! load the correct pieces of the integer and real vectors
  174. ! if Rout%num_segs(proc)=1, then this will do one loop
  175. do nseg = 1,Rout%num_segs(proc)
  176. seg_start = Rout%seg_starts(proc,nseg)
  177. seg_end = seg_start + Rout%seg_lengths(proc,nseg)-1
  178. do VectIndex = seg_start,seg_end
  179. do AttrIndex = 1,numi
  180. Rout%ip1(proc)%pi(j) = aV%iAttr(AttrIndex,VectIndex)
  181. j=j+1
  182. enddo
  183. do AttrIndex = 1,numr
  184. Rout%rp1(proc)%pr(k) = aV%rAttr(AttrIndex,VectIndex)
  185. k=k+1
  186. enddo
  187. enddo
  188. enddo
  189. ! Send the integer data
  190. if(numi .ge. 1) then
  191. ! set tag
  192. mytag = DefaultTag
  193. if(present(Tag)) mytag=Tag
  194. call MPI_ISEND(Rout%ip1(proc)%pi(1), &
  195. Rout%locsize(proc)*numi,MP_INTEGER,Rout%pe_list(proc), &
  196. mytag,ThisMCTWorld%MCT_comm,Rout%ireqs(proc),ier)
  197. if(ier /= 0) call MP_perr_die(myname_,'MPI_ISEND(ints)',ier)
  198. endif
  199. ! Send the real data
  200. if(numr .ge. 1) then
  201. ! set tag
  202. mytag = DefaultTag + 1
  203. if(present(Tag)) mytag=Tag +1
  204. call MPI_ISEND(Rout%rp1(proc)%pr(1), &
  205. Rout%locsize(proc)*numr,mp_Type_rp1,Rout%pe_list(proc), &
  206. mytag,ThisMCTWorld%MCT_comm,Rout%rreqs(proc),ier)
  207. if(ier /= 0) call MP_perr_die(myname_,'MPI_ISEND(reals)',ier)
  208. endif
  209. enddo
  210. if (unordered) then
  211. call AttrVect_clean(aVtmp)
  212. nullify(aV)
  213. else
  214. nullify(aV)
  215. endif
  216. end subroutine isend_
  217. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  218. ! Math and Computer Science Division, Argonne National Laboratory !
  219. !BOP -------------------------------------------------------------------
  220. !
  221. ! !IROUTINE: waitsend_ - Wait for a distributed non-blocking send to complete
  222. !
  223. ! !DESCRIPTION:
  224. ! Wait for the data being sent with the {\tt Router} {\tt Rout} to complete.
  225. !
  226. ! !INTERFACE:
  227. subroutine waitsend_(Rout)
  228. !
  229. ! !USES:
  230. !
  231. implicit none
  232. ! !INPUT PARAMETERS:
  233. !
  234. Type(Router), intent(inout) :: Rout
  235. ! !REVISION HISTORY:
  236. ! 24Jul03 - R. Jacob <jacob@mcs.anl.gov> - First working version is
  237. ! the wait part of original send_
  238. !EOP ___________________________________________________________________
  239. character(len=*),parameter :: myname_=myname//'::waitsend_'
  240. integer :: proc,ier
  241. ! Return if nothing to wait for
  242. if(Rout%nprocs .eq. 0 ) RETURN
  243. ! wait for all sends to complete
  244. if(Rout%numiatt .ge. 1) then
  245. call MPI_WAITALL(Rout%nprocs,Rout%ireqs,Rout%istatus,ier)
  246. if(ier /= 0) call MP_perr_die(myname_,'MPI_WAITALL(ints)',ier)
  247. do proc=1,Rout%nprocs
  248. deallocate(Rout%ip1(proc)%pi,stat=ier)
  249. if(ier/=0) call die(myname_,'deallocate(ip1%pi)',ier)
  250. enddo
  251. endif
  252. if(Rout%numratt .ge. 1) then
  253. call MPI_WAITALL(Rout%nprocs,Rout%rreqs,Rout%rstatus,ier)
  254. if(ier /= 0) call MP_perr_die(myname_,'MPI_WAITALL(reals)',ier)
  255. do proc=1,Rout%nprocs
  256. deallocate(Rout%rp1(proc)%pr,stat=ier)
  257. if(ier/=0) call die(myname_,'deallocate(rp1%pi)',ier)
  258. enddo
  259. endif
  260. end subroutine waitsend_
  261. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  262. ! Math and Computer Science Division, Argonne National Laboratory !
  263. !BOP -------------------------------------------------------------------
  264. !
  265. ! !IROUTINE: send_ - Distributed blocking send of an Attribute Vector
  266. !
  267. ! !DESCRIPTION:
  268. ! Send the the data in the {\tt AttrVect} {\tt aV} to the
  269. ! component specified in the {\tt Router} {\tt Rout}. An error will
  270. ! result if the size of the attribute vector does not match the size
  271. ! parameter stored in the {\tt Router}.
  272. !
  273. ! Requires a corresponding {\tt recv\_} or {\tt irecv\_} to be called on the other
  274. ! component.
  275. !
  276. ! The optional argument {\tt Tag} can be used to set the tag value used in
  277. ! the data transfer. DefaultTag will be used otherwise. {\tt Tag} must be
  278. ! the same in the matching {\tt recv\_} or {\tt irecv\_}.
  279. !
  280. ! {\bf N.B.:} The {\tt AttrVect} argument in the corresponding
  281. ! {\tt recv} call is assumed to have exactly the same attributes
  282. ! in exactly the same order as {\tt aV}.
  283. !
  284. ! !INTERFACE:
  285. subroutine send_(aV, Rout, Tag)
  286. !
  287. ! !USES:
  288. !
  289. implicit none
  290. ! !INPUT PARAMETERS:
  291. !
  292. Type(AttrVect), intent(in) :: aV
  293. Type(Router), intent(inout) :: Rout
  294. integer,optional, intent(in) :: Tag
  295. ! !REVISION HISTORY:
  296. ! 24Jul03 - R. Jacob <jacob@mcs.anl.gov> - New version uses isend and waitsend
  297. !EOP ___________________________________________________________________
  298. character(len=*),parameter :: myname_=myname//'::send_'
  299. call isend_(aV,Rout,Tag)
  300. call waitsend_(Rout)
  301. end subroutine send_
  302. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  303. ! Math and Computer Science Division, Argonne National Laboratory !
  304. !BOP -------------------------------------------------------------------
  305. !
  306. ! !IROUTINE: irecv_ - Distributed receive of an Attribute Vector
  307. !
  308. ! !DESCRIPTION:
  309. ! Recieve into the {\tt AttrVect} {\tt aV} the data coming from the
  310. ! component specified in the {\tt Router} {\tt Rout}. An error will
  311. ! result if the size of the attribute vector does not match the size
  312. ! parameter stored in the {\tt Router}.
  313. !
  314. ! Requires a corresponding {\tt send\_} or {\tt isend\_} to be called
  315. ! on the other component.
  316. !
  317. ! The optional argument {\tt Tag} can be used to set the tag value used in
  318. ! the data transfer. DefaultTag will be used otherwise. {\tt Tag} must be
  319. ! the same in the matching {\tt send\_} or {\tt isend\_}.
  320. !
  321. ! If data for a grid point is coming from more than one process, {\tt recv\_}
  322. ! will overwrite the duplicate values leaving the last received value
  323. ! in the output aV. If the optional argument {\tt Sum} is invoked, the output
  324. ! will contain the sum of any duplicate values received for the same grid point.
  325. !
  326. ! Will return as soon as MPI\_IRECV's are posted. Call {\tt waitrecv\_} to
  327. ! complete the receive operation.
  328. !
  329. ! {\bf N.B.:} The {\tt AttrVect} argument in the corresponding
  330. ! {\tt send\_} call is assumed to have exactly the same attributes
  331. ! in exactly the same order as {\tt aV}.
  332. !
  333. ! !INTERFACE:
  334. subroutine irecv_(aV, Rout, Tag, Sum)
  335. !
  336. ! !USES:
  337. !
  338. implicit none
  339. ! !INPUT/OUTPUT PARAMETERS:
  340. !
  341. Type(AttrVect), intent(inout) :: aV
  342. ! !INPUT PARAMETERS:
  343. !
  344. Type(Router), intent(inout) :: Rout
  345. integer,optional, intent(in) :: Tag
  346. logical,optional, intent(in) :: Sum
  347. ! !REVISION HISTORY:
  348. ! 07Feb01 - R. Jacob <jacob@mcs.anl.gov> - initial prototype
  349. ! 07Jun01 - R. Jacob <jacob@mcs.anl.gov> - remove logic to
  350. ! check "direction" of Router. remove references
  351. ! to ThisMCTWorld%mylrank
  352. ! 03Aug01 - E.T. Ong <eong@mcs.anl.gov> - explicity specify starting
  353. ! address in MPI_RECV
  354. ! 27Nov01 - E.T. Ong <eong@mcs.anl.gov> - deallocated to prevent
  355. ! memory leaks
  356. ! 15Feb02 - R. Jacob <jacob@mcs.anl.gov> - Use MCT_comm
  357. ! 26Mar02 - E. Ong <eong@mcs.anl.gov> - Apply faster copy order.
  358. ! 26Sep02 - R. Jacob <jacob@mcs.anl.gov> - Check Av against Router lAvsize
  359. ! 08Nov02 - R. Jacob <jacob@mcs.anl.gov> - MCT_Recv is now recv_ in m_Transfer
  360. ! 11Nov02 - R. Jacob <jacob@mcs.anl.gov> - Add optional Sum argument to
  361. ! tell recv_ to sum data for the same point received from multiple
  362. ! processors. Replaces recvsum_ which had replaced MCT_Recvsum.
  363. ! Use DefaultTag and add optional Tag argument
  364. ! 25Jul03 - R. Jacob <jacob@mcs.anl.gov> - break into irecv_ and waitrecv_
  365. !EOP ___________________________________________________________________
  366. character(len=*),parameter :: myname_=myname//'::irecv_'
  367. integer :: numi,numr,i,j,k,ier
  368. integer :: mycomp,othercomp
  369. integer :: seg_start,seg_end
  370. integer :: proc,numprocs,nseg,mytag
  371. integer :: mp_Type_rp1
  372. logical :: DoSum
  373. !--------------------------------------------------------
  374. ! Return if no one to receive from
  375. if(Rout%nprocs .eq. 0 ) RETURN
  376. !check Av size against Router
  377. !
  378. if(lsize(aV) /= Rout%lAvsize) then
  379. write(stderr,'(2a)') myname_, &
  380. ' MCTERROR: AV size not appropriate for this Router...exiting'
  381. call die(myname_)
  382. endif
  383. DoSum = .false.
  384. if(present(Sum)) DoSum=Sum
  385. mycomp=Rout%comp1id
  386. othercomp=Rout%comp2id
  387. ! find total number of real and integer vectors
  388. ! for now, assume we are receiving all of them
  389. Rout%numiatt = nIAttr(aV)
  390. Rout%numratt = nRAttr(aV)
  391. numi = Rout%numiatt
  392. numr = Rout%numratt
  393. !!!!!!!!!!!!!! IF RECEVING INTEGER DATA
  394. if(numi .ge. 1) then
  395. ! allocate buffers to hold all incoming data
  396. do proc=1,Rout%nprocs
  397. allocate(Rout%ip1(proc)%pi(Rout%locsize(proc)*numi),stat=ier)
  398. if(ier/=0) call die(myname_,'allocate(Rout%ip1%pi)',ier)
  399. enddo
  400. endif
  401. !!!!!!!!!!!!!! IF RECEIVING REAL DATA
  402. if(numr .ge. 1) then
  403. ! allocate buffers to hold all incoming data
  404. do proc=1,Rout%nprocs
  405. allocate(Rout%rp1(proc)%pr(Rout%locsize(proc)*numr),stat=ier)
  406. if(ier/=0) call die(myname_,'allocate(Rout%rp1%pr)',ier)
  407. enddo
  408. mp_Type_rp1=MP_Type(Rout%rp1(1)%pr(1))
  409. endif
  410. ! Post all MPI_IRECV
  411. do proc=1,Rout%nprocs
  412. ! receive the integer data
  413. if(numi .ge. 1) then
  414. ! set tag
  415. mytag = DefaultTag
  416. if(present(Tag)) mytag=Tag
  417. if( Rout%num_segs(proc) > 1 .or. DoSum ) then
  418. call MPI_IRECV(Rout%ip1(proc)%pi(1), &
  419. Rout%locsize(proc)*numi,MP_INTEGER,Rout%pe_list(proc), &
  420. mytag,ThisMCTWorld%MCT_comm,Rout%ireqs(proc),ier)
  421. else
  422. call MPI_IRECV(aV%iAttr(1,Rout%seg_starts(proc,1)), &
  423. Rout%locsize(proc)*numi,MP_INTEGER,Rout%pe_list(proc), &
  424. mytag,ThisMCTWorld%MCT_comm,Rout%ireqs(proc),ier)
  425. endif
  426. if(ier /= 0) call MP_perr_die(myname_,'MPI_IRECV(ints)',ier)
  427. endif
  428. ! receive the real data
  429. if(numr .ge. 1) then
  430. ! corresponding tag logic must be in send_
  431. mytag = DefaultTag + 1
  432. if(present(Tag)) mytag=Tag +1
  433. if( Rout%num_segs(proc) > 1 .or. DoSum ) then
  434. call MPI_IRECV(Rout%rp1(proc)%pr(1), &
  435. Rout%locsize(proc)*numr,mp_Type_rp1,Rout%pe_list(proc), &
  436. mytag,ThisMCTWorld%MCT_comm,Rout%rreqs(proc),ier)
  437. else
  438. call MPI_IRECV(aV%rAttr(1,Rout%seg_starts(proc,1)), &
  439. Rout%locsize(proc)*numr,mp_Type_rp1,Rout%pe_list(proc), &
  440. mytag,ThisMCTWorld%MCT_comm,Rout%rreqs(proc),ier)
  441. endif
  442. if(ier /= 0) call MP_perr_die(myname_,'MPI_IRECV(reals)',ier)
  443. endif
  444. enddo
  445. end subroutine irecv_
  446. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  447. ! Math and Computer Science Division, Argonne National Laboratory !
  448. !BOP -------------------------------------------------------------------
  449. !
  450. ! !IROUTINE: waitrecv_ - Wait for a distributed non-blocking recv to complete
  451. !
  452. ! !DESCRIPTION:
  453. ! Wait for the data being received with the {\tt Router} {\tt Rout} to complete.
  454. ! When done, copy the data into the {\tt AttrVect} {\tt aV}.
  455. !
  456. ! !INTERFACE:
  457. subroutine waitrecv_(aV, Rout, Sum)
  458. !
  459. ! !USES:
  460. !
  461. implicit none
  462. ! !INPUT/OUTPUT PARAMETERS:
  463. !
  464. Type(AttrVect), intent(inout) :: aV
  465. Type(Router), intent(inout) :: Rout
  466. ! !INPUT PARAMETERS:
  467. !
  468. logical,optional, intent(in) :: Sum
  469. ! !REVISION HISTORY:
  470. ! 25Jul03 - R. Jacob <jacob@mcs.anl.gov> - First working version is the wait
  471. ! and copy parts from old recv_.
  472. ! 25Jan08 - R. Jacob <jacob@mcs.anl.gov> - Handle unordered GSMaps by
  473. ! applying permutation to received array.
  474. !EOP ___________________________________________________________________
  475. character(len=*),parameter :: myname_=myname//'::waitrecv_'
  476. integer :: proc,ier,j,k,nseg
  477. integer :: AttrIndex,VectIndex,seg_start,seg_end
  478. logical :: DoSum
  479. logical :: unordered
  480. ! Return if nothing to wait for
  481. if(Rout%nprocs .eq. 0 ) RETURN
  482. !check Av size against Router
  483. !
  484. if(lsize(aV) /= Rout%lAvsize) then
  485. write(stderr,'(2a)') myname_, &
  486. ' MCTERROR: AV size not appropriate for this Router...exiting'
  487. call die(myname_)
  488. endif
  489. unordered = associated(Rout%permarr)
  490. DoSum = .false.
  491. if(present(Sum)) DoSum=Sum
  492. ! wait for all recieves to complete
  493. if(Rout%numiatt .ge. 1) then
  494. call MPI_WAITALL(Rout%nprocs,Rout%ireqs,Rout%istatus,ier)
  495. if(ier /= 0) call MP_perr_die(myname_,'MPI_WAITALL(ints)',ier)
  496. endif
  497. if(Rout%numratt .ge. 1) then
  498. call MPI_WAITALL(Rout%nprocs,Rout%rreqs,Rout%rstatus,ier)
  499. if(ier /= 0) call MP_perr_die(myname_,'MPI_WAITALL(reals)',ier)
  500. endif
  501. ! Load data which came from each processor
  502. do proc=1,Rout%nprocs
  503. if( (Rout%num_segs(proc) > 1) .or. DoSum ) then
  504. j=1
  505. k=1
  506. if(DoSum) then
  507. ! sum the correct pieces of the integer and real vectors
  508. do nseg = 1,Rout%num_segs(proc)
  509. seg_start = Rout%seg_starts(proc,nseg)
  510. seg_end = seg_start + Rout%seg_lengths(proc,nseg)-1
  511. do VectIndex = seg_start,seg_end
  512. do AttrIndex = 1,Rout%numiatt
  513. aV%iAttr(AttrIndex,VectIndex)= &
  514. aV%iAttr(AttrIndex,VectIndex)+Rout%ip1(proc)%pi(j)
  515. j=j+1
  516. enddo
  517. do AttrIndex = 1,Rout%numratt
  518. aV%rAttr(AttrIndex,VectIndex)= &
  519. aV%rAttr(AttrIndex,VectIndex)+Rout%rp1(proc)%pr(k)
  520. k=k+1
  521. enddo
  522. enddo
  523. enddo
  524. else
  525. ! load the correct pieces of the integer and real vectors
  526. do nseg = 1,Rout%num_segs(proc)
  527. seg_start = Rout%seg_starts(proc,nseg)
  528. seg_end = seg_start + Rout%seg_lengths(proc,nseg)-1
  529. do VectIndex = seg_start,seg_end
  530. do AttrIndex = 1,Rout%numiatt
  531. aV%iAttr(AttrIndex,VectIndex)=Rout%ip1(proc)%pi(j)
  532. j=j+1
  533. enddo
  534. do AttrIndex = 1,Rout%numratt
  535. aV%rAttr(AttrIndex,VectIndex)=Rout%rp1(proc)%pr(k)
  536. k=k+1
  537. enddo
  538. enddo
  539. enddo
  540. endif
  541. endif
  542. enddo
  543. !........................WAITANY METHOD................................
  544. !
  545. !....NOTE: Make status argument a 1-dimensional array
  546. ! ! Load data which came from each processor
  547. ! do numprocs = 1,Rout%nprocs
  548. ! ! Load the integer data
  549. ! if(Rout%numiatt .ge. 1) then
  550. ! call MPI_WAITANY(Rout%nprocs,Rout%ireqs,proc,Rout%istatus,ier)
  551. ! if(ier /= 0) call MP_perr_die(myname_,'MPI_WAITANY(ints)',ier)
  552. ! j=1
  553. ! ! load the correct pieces of the integer vectors
  554. ! do nseg = 1,Rout%num_segs(proc)
  555. ! seg_start = Rout%seg_starts(proc,nseg)
  556. ! seg_end = seg_start + Rout%seg_lengths(proc,nseg)-1
  557. ! do VectIndex = seg_start,seg_end
  558. ! do AttrIndex = 1,Rout%numiatt
  559. ! aV%iAttr(AttrIndex,VectIndex)=Rout%ip1(proc)%pi(j)
  560. ! j=j+1
  561. ! enddo
  562. ! enddo
  563. ! enddo
  564. ! endif
  565. ! ! Load the real data
  566. ! if(numr .ge. 1) then
  567. ! call MPI_WAITANY(Rout%nprocs,Rout%rreqs,proc,Rout%rstatus,ier)
  568. ! if(ier /= 0) call MP_perr_die(myname_,'MPI_WAITANY(reals)',ier)
  569. ! k=1
  570. ! ! load the correct pieces of the real vectors
  571. ! do nseg = 1,Rout%num_segs(proc)
  572. ! seg_start = Rout%seg_starts(proc,nseg)
  573. ! seg_end = seg_start + Rout%seg_lengths(proc,nseg)-1
  574. ! do VectIndex = seg_start,seg_end
  575. ! do AttrIndex = 1,numr
  576. ! aV%rAttr(AttrIndex,VectIndex)=Rout%rp1(proc)%pr(k)
  577. ! k=k+1
  578. ! enddo
  579. ! enddo
  580. ! enddo
  581. ! endif
  582. ! enddo
  583. !........................................................................
  584. ! Deallocate all structures
  585. if(Rout%numiatt .ge. 1) then
  586. ! Deallocate the receive buffers
  587. do proc=1,Rout%nprocs
  588. deallocate(Rout%ip1(proc)%pi,stat=ier)
  589. if(ier/=0) call die(myname_,'deallocate(Rout%ip1%pi)',ier)
  590. enddo
  591. endif
  592. if(Rout%numratt .ge. 1) then
  593. ! Deallocate the receive buffers
  594. do proc=1,Rout%nprocs
  595. deallocate(Rout%rp1(proc)%pr,stat=ier)
  596. if(ier/=0) call die(myname_,'deallocate(Rout%rp1%pr)',ier)
  597. enddo
  598. endif
  599. if (unordered) call Unpermute(aV,Rout%permarr)
  600. end subroutine waitrecv_
  601. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  602. ! Math and Computer Science Division, Argonne National Laboratory !
  603. !BOP -------------------------------------------------------------------
  604. !
  605. ! !IROUTINE: recv_ - Distributed receive of an Attribute Vector
  606. !
  607. ! !DESCRIPTION:
  608. ! Recieve into the {\tt AttrVect} {\tt aV} the data coming from the
  609. ! component specified in the {\tt Router} {\tt Rout}. An error will
  610. ! result if the size of the attribute vector does not match the size
  611. ! parameter stored in the {\tt Router}.
  612. !
  613. ! Requires a corresponding {\tt send\_} or {\tt isend\_}to be called
  614. ! on the other component.
  615. !
  616. ! The optional argument {\tt Tag} can be used to set the tag value used in
  617. ! the data transfer. DefaultTag will be used otherwise. {\tt Tag} must be
  618. ! the same in the matching {\tt send\_}
  619. !
  620. ! If data for a grid point is coming from more than one process, {\tt recv\_}
  621. ! will overwrite the duplicate values leaving the last received value
  622. ! in the output aV. If the optional argument {\tt Sum} is invoked, the output
  623. ! will contain the sum of any duplicate values received for the same grid point.
  624. !
  625. ! Will not return until all data has been received.
  626. !
  627. ! {\bf N.B.:} The {\tt AttrVect} argument in the corresponding
  628. ! {\tt send\_} call is assumed to have exactly the same attributes
  629. ! in exactly the same order as {\tt aV}.
  630. !
  631. ! !INTERFACE:
  632. subroutine recv_(aV, Rout, Tag, Sum)
  633. !
  634. ! !USES:
  635. !
  636. implicit none
  637. ! !INPUT/OUTPUT PARAMETERS:
  638. !
  639. Type(AttrVect), intent(inout) :: aV
  640. ! !INPUT PARAMETERS:
  641. !
  642. Type(Router), intent(inout) :: Rout
  643. integer,optional, intent(in) :: Tag
  644. logical,optional, intent(in) :: Sum
  645. ! !REVISION HISTORY:
  646. ! 25Jul03 - R. Jacob <jacob@mcs.anl.gov> - Rewrite using irecv and waitrecv
  647. !EOP ___________________________________________________________________
  648. character(len=*),parameter :: myname_=myname//'::recv_'
  649. call irecv_(aV,Rout,Tag,Sum)
  650. call waitrecv_(aV,Rout,Sum)
  651. end subroutine recv_
  652. end module m_Transfer