m_GlobalSegMapComms.F90 18 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555
  1. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  2. ! Math and Computer Science Division, Argonne National Laboratory !
  3. !-----------------------------------------------------------------------
  4. ! CVS m_GlobalSegMapComms.F90,v 1.6 2004-04-23 23:58:47 jacob Exp
  5. ! CVS MCT_2_8_0
  6. !BOP -------------------------------------------------------------------
  7. !
  8. ! !MODULE: m_GlobalSegMapComms - GlobalSegMap Communications Support
  9. !
  10. ! !DESCRIPTION:
  11. !
  12. ! This module provides communications support for the {\tt GlobalSegMap}
  13. ! datatype. Both blocking and non-blocking point-to-point communications
  14. ! are provided for send (analogues to {\tt MPI\_SEND()/MPI\_ISEND()})
  15. ! A receive and broadcast method is also supplied.
  16. !
  17. ! !INTERFACE:
  18. module m_GlobalSegMapComms
  19. implicit none
  20. private ! except
  21. ! !PUBLIC MEMBER FUNCTIONS:
  22. public :: send
  23. public :: recv
  24. public :: isend
  25. public :: bcast
  26. interface bcast ; module procedure bcast_ ; end interface
  27. interface send ; module procedure send_ ; end interface
  28. interface recv ; module procedure recv_ ; end interface
  29. interface isend ; module procedure isend_ ; end interface
  30. ! !REVISION HISTORY:
  31. ! 11Aug03 - J.W. Larson <larson@mcs.anl.gov> - initial version
  32. !
  33. !EOP ___________________________________________________________________
  34. character(len=*),parameter :: myname='MCT::m_GlobalSegMapComms'
  35. contains
  36. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  37. ! Math and Computer Science Division, Argonne National Laboratory !
  38. !BOP -------------------------------------------------------------------
  39. !
  40. ! !IROUTINE: send_ - Point-to-point blocking Send of a GlobalSegMap
  41. !
  42. ! !DESCRIPTION:
  43. ! This routine performs a blocking send of a {\tt GlobalSegMap} (the
  44. ! input argument {\tt outgoingGSMap}) to the root processor on component
  45. ! {\tt comp\_id}. The input {\tt INTEGER} argument {\tt TagBase}
  46. ! is used to generate tags for the messages associated with this operation;
  47. ! there are six messages involved, so the user should avoid using tag
  48. ! values {\tt TagBase} and {\tt TagBase + 5}. All six messages are blocking.
  49. ! The success (failure) of this operation is reported in the zero
  50. ! (non-zero) value of the optional {\tt INTEGER} output variable {\tt status}.
  51. !
  52. ! !INTERFACE:
  53. subroutine send_(outgoingGSMap, comp_id, TagBase, status)
  54. !
  55. ! !USES:
  56. !
  57. use m_mpif90
  58. use m_die, only : MP_perr_die,die
  59. use m_stdio
  60. use m_GlobalSegMap, only : GlobalSegMap
  61. use m_GlobalSegMap, only : GlobalSegMap_ngseg => ngseg
  62. use m_GlobalSegMap, only : GlobalSegMap_comp_id => comp_ID
  63. use m_GlobalSegMap, only : GlobalSegMap_gsize => gsize
  64. use m_MCTWorld, only : ComponentToWorldRank
  65. use m_MCTWorld, only : ThisMCTWorld
  66. implicit none
  67. ! !INPUT PARAMETERS:
  68. type(GlobalSegMap), intent(IN) :: outgoingGSMap
  69. integer, intent(IN) :: comp_id
  70. integer, intent(IN) :: TagBase
  71. ! !OUTPUT PARAMETERS:
  72. integer, optional, intent(OUT) :: status
  73. ! !REVISION HISTORY:
  74. ! 13Aug03 - J.W. Larson <larson@mcs.anl.gov> - API and initial version.
  75. ! 26Aug03 - R. Jacob <jacob@mcs.anl.gov> - use same method as isend_
  76. ! 05Mar04 - R. Jacob <jacob@mcs.anl.gov> - match new isend_ method.
  77. !EOP ___________________________________________________________________
  78. character(len=*),parameter :: myname_=myname//'::send_'
  79. integer :: ierr
  80. integer :: destID
  81. integer :: nsegs
  82. if(present(status)) status = 0 ! the success value
  83. destID = ComponentToWorldRank(0, comp_id, ThisMCTWorld)
  84. ! Next, send the buffer size to destID so it can prepare a
  85. ! receive buffer of the correct size.
  86. nsegs = GlobalSegMap_ngseg(outgoingGSMap)
  87. call MPI_SEND(outgoingGSMap%comp_id, 1, MP_Type(outgoingGSMap%comp_id), destID, &
  88. TagBase, ThisMCTWorld%MCT_comm, ierr)
  89. if(ierr /= 0) then
  90. call MP_perr_die(myname_, 'Send compid failed',ierr)
  91. endif
  92. call MPI_SEND(outgoingGSMap%ngseg, 1, MP_Type(outgoingGSMap%ngseg), destID, &
  93. TagBase+1, ThisMCTWorld%MCT_comm, ierr)
  94. if(ierr /= 0) then
  95. call MP_perr_die(myname_, 'Send ngseg failed',ierr)
  96. endif
  97. call MPI_SEND(outgoingGSMap%gsize, 1, MP_Type(outgoingGSMap%gsize), destID, &
  98. TagBase+2, ThisMCTWorld%MCT_comm, ierr)
  99. if(ierr /= 0) then
  100. call MP_perr_die(myname_, 'Send gsize failed',ierr)
  101. endif
  102. ! Send segment information data (3 messages)
  103. call MPI_SEND(outgoingGSMap%start, nsegs, &
  104. MP_Type(outgoingGSMap%start(1)), &
  105. destID, TagBase+3, ThisMCTWorld%MCT_comm, ierr)
  106. if(ierr /= 0) then
  107. call MP_perr_die(myname_, 'Send outgoingGSMap%start failed',ierr)
  108. endif
  109. call MPI_SEND(outgoingGSMap%length, nsegs, &
  110. MP_Type(outgoingGSMap%length(1)), &
  111. destID, TagBase+4, ThisMCTWorld%MCT_comm, ierr)
  112. if(ierr /= 0) then
  113. call MP_perr_die(myname_, 'Send outgoingGSMap%length failed',ierr)
  114. endif
  115. call MPI_SEND(outgoingGSMap%pe_loc, nsegs, &
  116. MP_Type(outgoingGSMap%pe_loc(1)), &
  117. destID, TagBase+5, ThisMCTWorld%MCT_comm, ierr)
  118. if(ierr /= 0) then
  119. call MP_perr_die(myname_, 'Send outgoingGSMap%pe_loc failed',ierr)
  120. endif
  121. end subroutine send_
  122. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  123. ! Math and Computer Science Division, Argonne National Laboratory !
  124. !BOP -------------------------------------------------------------------
  125. !
  126. ! !IROUTINE: isend_ - Point-to-point Non-blocking Send of a GlobalSegMap
  127. !
  128. ! !DESCRIPTION:
  129. ! This routine performs a non-blocking send of a {\tt GlobalSegMap} (the
  130. ! input argument {\tt outgoingGSMap}) to the root processor on component
  131. ! {\tt comp\_id} The input {\tt INTEGER} argument {\tt TagBase}
  132. ! is used to generate tags for the messages associated with this operation;
  133. ! there are six messages involved, so the user should avoid using tag
  134. ! values {\tt TagBase} and {\tt TagBase + 5}. All six messages are non-
  135. ! blocking, and the request handles for them are returned in the output
  136. ! {\tt INTEGER} array {\tt reqHandle}, which can be checked for completion
  137. ! using any of MPI's wait functions. The success (failure) of
  138. ! this operation is reported in the zero (non-zero) value of the optional
  139. ! {\tt INTEGER} output variable {\tt status}.
  140. !
  141. ! {\bf N.B.}: Data is sent directly out of {\tt outgoingGSMap} so it
  142. ! must not be deleted until the send has completed.
  143. !
  144. ! {\bf N.B.}: The array {\tt reqHandle} represents allocated memory that
  145. ! must be deallocated when it is no longer needed. Failure to do so will
  146. ! create a memory leak.
  147. !
  148. ! !INTERFACE:
  149. subroutine isend_(outgoingGSMap, comp_id, TagBase, reqHandle, status)
  150. !
  151. ! !USES:
  152. !
  153. use m_mpif90
  154. use m_die, only : MP_perr_die,die
  155. use m_stdio
  156. use m_GlobalSegMap, only : GlobalSegMap
  157. use m_GlobalSegMap, only : GlobalSegMap_ngseg => ngseg
  158. use m_MCTWorld, only : ComponentToWorldRank
  159. use m_MCTWorld, only : ThisMCTWorld
  160. implicit none
  161. ! !INPUT PARAMETERS:
  162. type(GlobalSegMap), intent(IN) :: outgoingGSMap
  163. integer, intent(IN) :: comp_id
  164. integer, intent(IN) :: TagBase
  165. ! !OUTPUT PARAMETERS:
  166. integer, dimension(:), pointer :: reqHandle
  167. integer, optional, intent(OUT) :: status
  168. ! !REVISION HISTORY:
  169. ! 13Aug03 - J.W. Larson <larson@mcs.anl.gov> - API and initial version.
  170. ! 05Mar04 - R. Jacob <jacob@mcs.anl.gov> - Send everything directly out
  171. ! of input GSMap. Don't use a SendBuffer.
  172. !
  173. !EOP ___________________________________________________________________
  174. character(len=*),parameter :: myname_=myname//'::isend_'
  175. integer :: ierr,destID,nsegs
  176. if(present(status)) status = 0 ! the success value
  177. destID = ComponentToWorldRank(0, comp_id, ThisMCTWorld)
  178. allocate(reqHandle(6), stat=ierr)
  179. if(ierr /= 0) then
  180. write(stderr,'(2a,i8)') myname_, &
  181. 'FATAL--allocation of send buffer failed with ierr=',ierr
  182. call die(myname_)
  183. endif
  184. ! Next, send the buffer size to destID so it can prepare a
  185. ! receive buffer of the correct size (3 messages).
  186. nsegs = GlobalSegMap_ngseg(outgoingGSMap)
  187. call MPI_ISEND(outgoingGSMap%comp_id, 1, MP_Type(outgoingGSMap%comp_id), destID, &
  188. TagBase, ThisMCTWorld%MCT_comm, reqHandle(1), ierr)
  189. if(ierr /= 0) then
  190. call MP_perr_die(myname_, 'Send compid failed',ierr)
  191. endif
  192. call MPI_ISEND(outgoingGSMap%ngseg, 1, MP_Type(outgoingGSMap%ngseg), destID, &
  193. TagBase+1, ThisMCTWorld%MCT_comm, reqHandle(2), ierr)
  194. if(ierr /= 0) then
  195. call MP_perr_die(myname_, 'Send ngseg failed',ierr)
  196. endif
  197. call MPI_ISEND(outgoingGSMap%gsize, 1, MP_Type(outgoingGSMap%gsize), destID, &
  198. TagBase+2, ThisMCTWorld%MCT_comm, reqHandle(3), ierr)
  199. if(ierr /= 0) then
  200. call MP_perr_die(myname_, 'Send gsize failed',ierr)
  201. endif
  202. ! Send segment information data (3 messages)
  203. call MPI_ISEND(outgoingGSMap%start, nsegs, &
  204. MP_Type(outgoingGSMap%start(1)), &
  205. destID, TagBase+3, ThisMCTWorld%MCT_comm, reqHandle(4), ierr)
  206. if(ierr /= 0) then
  207. call MP_perr_die(myname_, 'Send outgoingGSMap%start failed',ierr)
  208. endif
  209. call MPI_ISEND(outgoingGSMap%length, nsegs, &
  210. MP_Type(outgoingGSMap%length(1)), &
  211. destID, TagBase+4, ThisMCTWorld%MCT_comm, reqHandle(5), ierr)
  212. if(ierr /= 0) then
  213. call MP_perr_die(myname_, 'Send outgoingGSMap%length failed',ierr)
  214. endif
  215. call MPI_ISEND(outgoingGSMap%pe_loc, nsegs, &
  216. MP_Type(outgoingGSMap%pe_loc(1)), &
  217. destID, TagBase+5, ThisMCTWorld%MCT_comm, reqHandle(6), ierr)
  218. if(ierr /= 0) then
  219. call MP_perr_die(myname_, 'Send outgoingGSMap%pe_loc failed',ierr)
  220. endif
  221. end subroutine isend_
  222. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  223. ! Math and Computer Science Division, Argonne National Laboratory !
  224. !BOP -------------------------------------------------------------------
  225. !
  226. ! !IROUTINE: recv_ - Point-to-point blocking Receive of a GlobalSegMap
  227. !
  228. ! !DESCRIPTION:
  229. ! This routine performs a blocking receive of a {\tt GlobalSegMap} (the
  230. ! input argument {\tt outgoingGSMap}) from the root processor on component
  231. ! {\tt comp\_id}. The input {\tt INTEGER} argument {\tt TagBase}
  232. ! is used to generate tags for the messages associated with this operation;
  233. ! there are six messages involved, so the user should avoid using tag
  234. ! values {\tt TagBase} and {\tt TagBase + 5}. The success (failure) of this
  235. ! operation is reported in the zero (non-zero) value of the optional {\tt INTEGER}
  236. ! output variable {\tt status}.
  237. !
  238. ! !INTERFACE:
  239. subroutine recv_(incomingGSMap, comp_id, TagBase, status)
  240. !
  241. ! !USES:
  242. !
  243. use m_mpif90
  244. use m_die, only : MP_perr_die, die
  245. use m_stdio
  246. use m_GlobalSegMap, only : GlobalSegMap
  247. use m_GlobalSegMap, only : GlobalSegMap_init => init
  248. use m_MCTWorld, only : ComponentToWorldRank
  249. use m_MCTWorld, only : ThisMCTWorld
  250. implicit none
  251. ! !INPUT PARAMETERS:
  252. integer, intent(IN) :: comp_id
  253. integer, intent(IN) :: TagBase
  254. ! !OUTPUT PARAMETERS:
  255. type(GlobalSegMap), intent(OUT) :: incomingGSMap
  256. integer, optional, intent(OUT) :: status
  257. ! !REVISION HISTORY:
  258. ! 13Aug03 - J.W. Larson <larson@mcs.anl.gov> - API and initial version.
  259. ! 25Aug03 - R.Jacob <larson@mcs.anl.gov> - rename to recv_.
  260. !EOP ___________________________________________________________________
  261. character(len=*),parameter :: myname_=myname//'::recv_'
  262. integer :: ierr,sourceID
  263. integer :: MPstatus(MP_STATUS_SIZE)
  264. integer :: RecvBuffer(3)
  265. if(present(status)) status = 0 ! the success value
  266. sourceID = ComponentToWorldRank(0, comp_id, ThisMCTWorld)
  267. ! Receive the GlobalSegMap's basic constants: component id,
  268. ! grid size, and number of segments. The number of segments
  269. ! is needed to construct the arrays into which segment
  270. ! information will be received. Thus, this receive blocks.
  271. call MPI_RECV(RecvBuffer(1), 1, MP_Type(RecvBuffer(1)), sourceID, &
  272. TagBase, ThisMCTWorld%MCT_comm, MPstatus, ierr)
  273. if(ierr /= 0) then
  274. call MP_perr_die(myname_, 'Receive of compid failed',ierr)
  275. endif
  276. call MPI_RECV(RecvBuffer(2), 1, MP_Type(RecvBuffer(2)), sourceID, &
  277. TagBase+1, ThisMCTWorld%MCT_comm, MPstatus, ierr)
  278. if(ierr /= 0) then
  279. call MP_perr_die(myname_, 'Receive of ngseg failed',ierr)
  280. endif
  281. call MPI_RECV(RecvBuffer(3), 1, MP_Type(RecvBuffer(3)), sourceID, &
  282. TagBase+2, ThisMCTWorld%MCT_comm, MPstatus, ierr)
  283. if(ierr /= 0) then
  284. call MP_perr_die(myname_, 'Receive of gsize failed',ierr)
  285. endif
  286. ! Create Empty GlobaSegMap into which segment information
  287. ! will be received
  288. call GlobalSegMap_init(incomingGSMap, RecvBuffer(1), RecvBuffer(2), &
  289. RecvBuffer(3))
  290. ! Receive segment information data (3 messages)
  291. call MPI_RECV(incomingGSMap%start, RecvBuffer(2), &
  292. MP_Type(incomingGSMap%start(1)), &
  293. sourceID, TagBase+3, ThisMCTWorld%MCT_comm, MPstatus, ierr)
  294. if(ierr /= 0) then
  295. call MP_perr_die(myname_, 'Recv incomingGSMap%start failed',ierr)
  296. endif
  297. call MPI_RECV(incomingGSMap%length, RecvBuffer(2), &
  298. MP_Type(incomingGSMap%length(1)), &
  299. sourceID, TagBase+4, ThisMCTWorld%MCT_comm, MPstatus, ierr)
  300. if(ierr /= 0) then
  301. call MP_perr_die(myname_, 'Recv incomingGSMap%length failed',ierr)
  302. endif
  303. call MPI_RECV(incomingGSMap%pe_loc, RecvBuffer(2), &
  304. MP_Type(incomingGSMap%pe_loc(1)), &
  305. sourceID, TagBase+5, ThisMCTWorld%MCT_comm, MPstatus, ierr)
  306. if(ierr /= 0) then
  307. call MP_perr_die(myname_, 'Recv incomingGSMap%pe_loc failed',ierr)
  308. endif
  309. end subroutine recv_
  310. !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  311. ! Math and Computer Science Division, Argonne National Laboratory !
  312. !BOP -------------------------------------------------------------------
  313. !
  314. ! !IROUTINE: bcast_ - broadcast a GlobalSegMap object
  315. !
  316. ! !DESCRIPTION:
  317. !
  318. ! The routine {\tt bcast\_()} takes the input/output {\em GlobalSegMap}
  319. ! argument {\tt GSMap} (on input valid only on the {\tt root} process,
  320. ! on output valid on all processes) and broadcasts it to all processes
  321. ! on the communicator associated with the F90 handle {\tt comm}. The
  322. ! success (failure) of this operation is returned as a zero (non-zero)
  323. ! value of the optional output {\tt INTEGER} argument {\tt status}.
  324. !
  325. ! !INTERFACE:
  326. subroutine bcast_(GSMap, root, comm, status)
  327. !
  328. ! !USES:
  329. !
  330. use m_mpif90
  331. use m_die, only : MP_perr_die,die
  332. use m_stdio
  333. use m_GlobalSegMap, only : GlobalSegMap
  334. implicit none
  335. ! !INPUT PARAMETERS:
  336. integer, intent(in) :: root
  337. integer, intent(in) :: comm
  338. ! !INPUT/OUTPUT PARAMETERS:
  339. type(GlobalSegMap), intent(inout) :: GSMap ! Output GlobalSegMap
  340. ! !OUTPUT PARAMETERS:
  341. integer, optional, intent(out) :: status ! global vector size
  342. ! !REVISION HISTORY:
  343. ! 17Oct01 - J.W. Larson <larson@mcs.anl.gov> - Initial version.
  344. ! 11Aug03 - J.W. Larson <larson@mcs.anl.gov> - Relocated from original
  345. ! location in m_GlobalSegMap.
  346. !EOP ___________________________________________________________________
  347. character(len=*),parameter :: myname_=myname//'::bcast_'
  348. integer :: myID, ierr, n
  349. integer, dimension(:), allocatable :: IntBuffer
  350. ! Step One: which process am I?
  351. call MP_COMM_RANK(comm, myID, ierr)
  352. if(ierr /= 0) call MP_perr_die(myname_,'MP_comm_rank()',ierr)
  353. ! Step Two: Broadcast the scalar bits of the GlobalSegMap from
  354. ! the root.
  355. allocate(IntBuffer(3), stat=ierr) ! allocate buffer space (all PEs)
  356. if(ierr /= 0) then
  357. if(.not. present(status)) then
  358. call die(myname_,'allocate(IntBuffer)',ierr)
  359. else
  360. write(stderr,*) myname_,':: error during allocate(IntBuffer)'
  361. status = 2
  362. return
  363. endif
  364. endif
  365. if(myID == root) then ! pack the buffer
  366. IntBuffer(1) = GSMap%comp_id
  367. IntBuffer(2) = GSMap%ngseg
  368. IntBuffer(3) = GSMap%gsize
  369. endif
  370. call MPI_BCAST(IntBuffer, 3, MP_type(IntBuffer(1)), root, comm, ierr)
  371. if(ierr /= 0) call MP_perr_die(myname_,'MPI_BCAST(IntBuffer)',ierr)
  372. if(myID /= root) then ! unpack from buffer to GSMap
  373. GSMap%comp_id = IntBuffer(1)
  374. GSMap%ngseg = IntBuffer(2)
  375. GSMap%gsize = IntBuffer(3)
  376. endif
  377. deallocate(IntBuffer, stat=ierr) ! deallocate buffer space
  378. if(ierr /= 0) then
  379. if(.not. present(status)) then
  380. call die(myname_,'deallocate(IntBuffer)',ierr)
  381. else
  382. write(stderr,*) myname_,':: error during deallocate(IntBuffer)'
  383. status = 4
  384. return
  385. endif
  386. endif
  387. ! Step Three: Broadcast the vector bits of GSMap from the root.
  388. ! Pack them into one big array to save latency costs associated
  389. ! with multiple broadcasts.
  390. allocate(IntBuffer(3*GSMap%ngseg), stat=ierr) ! allocate buffer space (all PEs)
  391. if(ierr /= 0) then
  392. if(.not. present(status)) then
  393. call die(myname_,'second allocate(IntBuffer)',ierr)
  394. else
  395. write(stderr,*) myname_,':: error during second allocate(IntBuffer)'
  396. status = 5
  397. return
  398. endif
  399. endif
  400. if(myID == root) then ! pack outgoing broadcast buffer
  401. do n=1,GSMap%ngseg
  402. IntBuffer(n) = GSMap%start(n)
  403. IntBuffer(GSMap%ngseg+n) = GSMap%length(n)
  404. IntBuffer(2*GSMap%ngseg+n) = GSMap%pe_loc(n)
  405. end do
  406. endif
  407. call MPI_BCAST(IntBuffer, 3*GSMap%ngseg, MP_Type(IntBuffer(1)), root, comm, ierr)
  408. if(ierr /= 0) call MP_perr_die(myname_,'Error in second MPI_BCAST(IntBuffer)',ierr)
  409. if(myID /= root) then ! Allocate GSMap%start, GSMap%length,...and fill them
  410. allocate(GSMap%start(GSMap%ngseg), GSMap%length(GSMap%ngseg), &
  411. GSMap%pe_loc(GSMap%ngseg), stat=ierr)
  412. if(ierr /= 0) then
  413. if(.not. present(status)) then
  414. call die(myname_,'off-root allocate(GSMap%start...)',ierr)
  415. else
  416. write(stderr,*) myname_,':: error during off-root allocate(GSMap%start...)'
  417. status = 7
  418. return
  419. endif
  420. endif
  421. do n=1,GSMap%ngseg ! unpack the buffer into the GlobalSegMap
  422. GSMap%start(n) = IntBuffer(n)
  423. GSMap%length(n) = IntBuffer(GSMap%ngseg+n)
  424. GSMap%pe_loc(n) = IntBuffer(2*GSMap%ngseg+n)
  425. end do
  426. endif
  427. ! Clean up buffer space:
  428. deallocate(IntBuffer, stat=ierr)
  429. if(ierr /= 0) then
  430. if(.not. present(status)) then
  431. call die(myname_,'second deallocate(IntBuffer)',ierr)
  432. else
  433. write(stderr,*) myname_,':: error during second deallocate(IntBuffer)'
  434. status = 8
  435. return
  436. endif
  437. endif
  438. end subroutine bcast_
  439. end module m_GlobalSegMapComms