m_FcComms.F90 20 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685
  1. !BOP -------------------------------------------------------------------
  2. !
  3. ! !MODULE: m_FcComms - MPI collective communication operators
  4. ! with explict flow control
  5. !
  6. ! !DESCRIPTION:
  7. !
  8. ! This module includes implementations of MPI collective operators that
  9. ! have proven problematic on certain systems when run at scale. By
  10. ! introducing additonal flow control, these problems (exhausting internal
  11. ! system resources) can be avoided. These routines were ported from
  12. ! the Community Atmosphere Model's spmd_utils.F90.
  13. !
  14. ! !INTERFACE:
  15. !
  16. ! Workaround for performance issue with rsend on cray systems with
  17. ! gemini interconnect
  18. !
  19. #ifdef _NO_MPI_RSEND
  20. #define MPI_RSEND MPI_SEND
  21. #define mpi_rsend mpi_send
  22. #define MPI_IRSEND MPI_ISEND
  23. #define mpi_irsend mpi_isend
  24. #endif
  25. module m_FcComms
  26. implicit none
  27. private ! except
  28. public :: fc_gather_int ! flow control version of mpi_gather for integer vectors
  29. public :: fc_gather_fp ! flow control version of mpi_gather for FP vectors
  30. public :: fc_gatherv_int ! flow control version of mpi_gatherv for integer vectors
  31. public :: fc_gatherv_fp ! flow control version of mpi_gatherv for integer vectors
  32. public :: get_fcblocksize ! get current value of max_gather_block_size
  33. public :: set_fcblocksize ! set current value of max_gather_block_size
  34. ! !REVISION HISTORY:
  35. ! 30Jan09 - P.H. Worley <worleyph@ornl.gov> - imported routines
  36. ! from CAM's spmd_utils to create this module.
  37. integer, public :: max_gather_block_size = 64
  38. character(len=*),parameter :: myname='MCT(MPEU)::m_FcComms'
  39. contains
  40. !BOP -------------------------------------------------------------------
  41. !
  42. ! !IROUTINE: fc_gather_int - Gather an array of type integer
  43. !
  44. ! !DESCRIPTION:
  45. ! This routine gathers a {\em distributed} array of type {\em integer}
  46. ! to the {\tt root} process. Explicit handshaking messages are used
  47. ! to control the number of processes communicating with the root
  48. ! at any one time.
  49. !
  50. ! If flow_cntl optional parameter
  51. ! < 0 : use MPI_Gather
  52. ! >= 0: use point-to-point with handshaking messages and
  53. ! preposting receive requests up to
  54. ! min(max(1,flow_cntl),max_gather_block_size)
  55. ! ahead if optional flow_cntl parameter is present.
  56. ! Otherwise, max_gather_block_size is used in its place.
  57. ! Default value is max_gather_block_size.
  58. ! !INTERFACE:
  59. !
  60. subroutine fc_gather_int (sendbuf, sendcnt, sendtype, &
  61. recvbuf, recvcnt, recvtype, &
  62. root, comm, flow_cntl )
  63. !
  64. ! !USES:
  65. !
  66. use m_die
  67. use m_mpif90
  68. !
  69. ! !INPUT PARAMETERS:
  70. !
  71. integer, intent(in) :: sendbuf(*)
  72. integer, intent(in) :: sendcnt
  73. integer, intent(in) :: sendtype
  74. integer, intent(in) :: recvcnt
  75. integer, intent(in) :: recvtype
  76. integer, intent(in) :: root
  77. integer, intent(in) :: comm
  78. integer, optional, intent(in) :: flow_cntl
  79. ! !OUTPUT PARAMETERS:
  80. !
  81. integer, intent(out) :: recvbuf(*)
  82. ! !REVISION HISTORY:
  83. ! 30Jan09 - P.H. Worley - imported from spmd_utils.F90
  84. !EOP ___________________________________________________________________
  85. character(len=*),parameter :: myname_=myname//'::fc_gather_int'
  86. integer :: signal
  87. logical fc_gather ! use explicit flow control?
  88. integer gather_block_size ! number of preposted receive requests
  89. integer :: mytid, mysize, mtag, p, i, count, displs
  90. integer :: preposts, head, tail
  91. integer :: rcvid(max_gather_block_size)
  92. integer :: status(MP_STATUS_SIZE)
  93. integer :: ier ! MPI error code
  94. signal = 1
  95. if ( present(flow_cntl) ) then
  96. if (flow_cntl >= 0) then
  97. gather_block_size = min(max(1,flow_cntl),max_gather_block_size)
  98. fc_gather = .true.
  99. else
  100. fc_gather = .false.
  101. endif
  102. else
  103. gather_block_size = max(1,max_gather_block_size)
  104. fc_gather = .true.
  105. endif
  106. if (fc_gather) then
  107. call mpi_comm_rank (comm, mytid, ier)
  108. call mpi_comm_size (comm, mysize, ier)
  109. mtag = 0
  110. if (root .eq. mytid) then
  111. ! prepost gather_block_size irecvs, and start receiving data
  112. preposts = min(mysize-1, gather_block_size)
  113. head = 0
  114. count = 0
  115. do p=0, mysize-1
  116. if (p .ne. root) then
  117. if (recvcnt > 0) then
  118. count = count + 1
  119. if (count > preposts) then
  120. tail = mod(head,preposts) + 1
  121. call mpi_wait (rcvid(tail), status, ier)
  122. end if
  123. head = mod(head,preposts) + 1
  124. displs = p*recvcnt
  125. call mpi_irecv ( recvbuf(displs+1), recvcnt, &
  126. recvtype, p, mtag, comm, rcvid(head), &
  127. ier )
  128. call mpi_send ( signal, 1, recvtype, p, mtag, comm, ier )
  129. end if
  130. end if
  131. end do
  132. ! copy local data
  133. displs = mytid*recvcnt
  134. do i=1,sendcnt
  135. recvbuf(displs+i) = sendbuf(i)
  136. enddo
  137. ! wait for final data
  138. do i=1,min(count,preposts)
  139. call mpi_wait (rcvid(i), status, ier)
  140. enddo
  141. else
  142. if (sendcnt > 0) then
  143. call mpi_recv ( signal, 1, sendtype, root, mtag, comm, &
  144. status, ier )
  145. call mpi_rsend ( sendbuf, sendcnt, sendtype, root, mtag, &
  146. comm, ier )
  147. end if
  148. endif
  149. if (ier /= 0) then
  150. call MP_perr_die(myname_,':: (point-to-point implementation)',ier)
  151. end if
  152. else
  153. call mpi_gather (sendbuf, sendcnt, sendtype, &
  154. recvbuf, recvcnt, recvtype, &
  155. root, comm, ier)
  156. if (ier /= 0) then
  157. call MP_perr_die(myname_,':: MPI_GATHER',ier)
  158. end if
  159. endif
  160. return
  161. end subroutine fc_gather_int
  162. !BOP -------------------------------------------------------------------
  163. !
  164. ! !IROUTINE: fc_gather_fp - Gather an array of type FP
  165. !
  166. ! !DESCRIPTION:
  167. ! This routine gathers a {\em distributed} array of type {\em FP} to
  168. ! the {\tt root} process. Explicit handshaking messages are used
  169. ! to control the number of processes communicating with the root
  170. ! at any one time.
  171. !
  172. ! If flow_cntl optional parameter
  173. ! < 0 : use MPI_Gather
  174. ! >= 0: use point-to-point with handshaking messages and
  175. ! preposting receive requests up to
  176. ! min(max(1,flow_cntl),max_gather_block_size)
  177. ! ahead if optional flow_cntl parameter is present.
  178. ! Otherwise, max_gather_block_size is used in its place.
  179. ! Default value is max_gather_block_size.
  180. ! !INTERFACE:
  181. !
  182. subroutine fc_gather_fp (sendbuf, sendcnt, sendtype, &
  183. recvbuf, recvcnt, recvtype, &
  184. root, comm, flow_cntl )
  185. !
  186. ! !USES:
  187. !
  188. use m_realkinds, only : FP
  189. use m_die
  190. use m_mpif90
  191. !
  192. ! !INPUT PARAMETERS:
  193. !
  194. real (FP), intent(in) :: sendbuf(*)
  195. integer, intent(in) :: sendcnt
  196. integer, intent(in) :: sendtype
  197. integer, intent(in) :: recvcnt
  198. integer, intent(in) :: recvtype
  199. integer, intent(in) :: root
  200. integer, intent(in) :: comm
  201. integer, optional, intent(in) :: flow_cntl
  202. ! !OUTPUT PARAMETERS:
  203. !
  204. real (FP), intent(out) :: recvbuf(*)
  205. ! !REVISION HISTORY:
  206. ! 30Jan09 - P.H. Worley - imported from spmd_utils.F90
  207. !EOP ___________________________________________________________________
  208. character(len=*),parameter :: myname_=myname//'::fc_gather_fp'
  209. real (FP) :: signal
  210. logical fc_gather ! use explicit flow control?
  211. integer gather_block_size ! number of preposted receive requests
  212. integer :: mytid, mysize, mtag, p, i, count, displs
  213. integer :: preposts, head, tail
  214. integer :: rcvid(max_gather_block_size)
  215. integer :: status(MP_STATUS_SIZE)
  216. integer :: ier ! MPI error code
  217. signal = 1.0
  218. if ( present(flow_cntl) ) then
  219. if (flow_cntl >= 0) then
  220. gather_block_size = min(max(1,flow_cntl),max_gather_block_size)
  221. fc_gather = .true.
  222. else
  223. fc_gather = .false.
  224. endif
  225. else
  226. gather_block_size = max(1,max_gather_block_size)
  227. fc_gather = .true.
  228. endif
  229. if (fc_gather) then
  230. call mpi_comm_rank (comm, mytid, ier)
  231. call mpi_comm_size (comm, mysize, ier)
  232. mtag = 0
  233. if (root .eq. mytid) then
  234. ! prepost gather_block_size irecvs, and start receiving data
  235. preposts = min(mysize-1, gather_block_size)
  236. head = 0
  237. count = 0
  238. do p=0, mysize-1
  239. if (p .ne. root) then
  240. if (recvcnt > 0) then
  241. count = count + 1
  242. if (count > preposts) then
  243. tail = mod(head,preposts) + 1
  244. call mpi_wait (rcvid(tail), status, ier)
  245. end if
  246. head = mod(head,preposts) + 1
  247. displs = p*recvcnt
  248. call mpi_irecv ( recvbuf(displs+1), recvcnt, &
  249. recvtype, p, mtag, comm, rcvid(head), &
  250. ier )
  251. call mpi_send ( signal, 1, recvtype, p, mtag, comm, ier )
  252. end if
  253. end if
  254. end do
  255. ! copy local data
  256. displs = mytid*recvcnt
  257. do i=1,sendcnt
  258. recvbuf(displs+i) = sendbuf(i)
  259. enddo
  260. ! wait for final data
  261. do i=1,min(count,preposts)
  262. call mpi_wait (rcvid(i), status, ier)
  263. enddo
  264. else
  265. if (sendcnt > 0) then
  266. call mpi_recv ( signal, 1, sendtype, root, mtag, comm, &
  267. status, ier )
  268. call mpi_rsend ( sendbuf, sendcnt, sendtype, root, mtag, &
  269. comm, ier )
  270. end if
  271. endif
  272. if (ier /= 0) then
  273. call MP_perr_die(myname_,':: (point-to-point implementation)',ier)
  274. end if
  275. else
  276. call mpi_gather (sendbuf, sendcnt, sendtype, &
  277. recvbuf, recvcnt, recvtype, &
  278. root, comm, ier)
  279. if (ier /= 0) then
  280. call MP_perr_die(myname_,':: MPI_GATHER',ier)
  281. end if
  282. endif
  283. return
  284. end subroutine fc_gather_fp
  285. !BOP -------------------------------------------------------------------
  286. !
  287. ! !IROUTINE: fc_gatherv_int - Gather an array of type integer
  288. !
  289. ! !DESCRIPTION:
  290. ! This routine gathers a {\em distributed} array of type {\em integer}
  291. ! to the {\tt root} process. Explicit handshaking messages are used
  292. ! to control the number of processes communicating with the root
  293. ! at any one time.
  294. !
  295. ! If flow_cntl optional parameter
  296. ! < 0 : use MPI_Gatherv
  297. ! >= 0: use point-to-point with handshaking messages and
  298. ! preposting receive requests up to
  299. ! min(max(1,flow_cntl),max_gather_block_size)
  300. ! ahead if optional flow_cntl parameter is present.
  301. ! Otherwise, max_gather_block_size is used in its place.
  302. ! Default value is max_gather_block_size.
  303. ! !INTERFACE:
  304. !
  305. subroutine fc_gatherv_int (sendbuf, sendcnt, sendtype, &
  306. recvbuf, recvcnts, displs, recvtype, &
  307. root, comm, flow_cntl )
  308. !
  309. ! !USES:
  310. !
  311. use m_die
  312. use m_mpif90
  313. !
  314. ! !INPUT PARAMETERS:
  315. !
  316. integer, intent(in) :: sendbuf(*)
  317. integer, intent(in) :: sendcnt
  318. integer, intent(in) :: sendtype
  319. integer, intent(in) :: recvcnts(*)
  320. integer, intent(in) :: displs(*)
  321. integer, intent(in) :: recvtype
  322. integer, intent(in) :: root
  323. integer, intent(in) :: comm
  324. integer, optional, intent(in) :: flow_cntl
  325. ! !OUTPUT PARAMETERS:
  326. !
  327. integer, intent(out) :: recvbuf(*)
  328. ! !REVISION HISTORY:
  329. ! 30Jan09 - P.H. Worley - imported from spmd_utils.F90
  330. !EOP ___________________________________________________________________
  331. character(len=*),parameter :: myname_=myname//'::fc_gatherv_int'
  332. integer :: signal
  333. logical fc_gather ! use explicit flow control?
  334. integer gather_block_size ! number of preposted receive requests
  335. integer :: mytid, mysize, mtag, p, q, i, count
  336. integer :: preposts, head, tail
  337. integer :: rcvid(max_gather_block_size)
  338. integer :: status(MP_STATUS_SIZE)
  339. integer :: ier ! MPI error code
  340. signal = 1
  341. if ( present(flow_cntl) ) then
  342. if (flow_cntl >= 0) then
  343. gather_block_size = min(max(1,flow_cntl),max_gather_block_size)
  344. fc_gather = .true.
  345. else
  346. fc_gather = .false.
  347. endif
  348. else
  349. gather_block_size = max(1,max_gather_block_size)
  350. fc_gather = .true.
  351. endif
  352. if (fc_gather) then
  353. call mpi_comm_rank (comm, mytid, ier)
  354. call mpi_comm_size (comm, mysize, ier)
  355. mtag = 0
  356. if (root .eq. mytid) then
  357. ! prepost gather_block_size irecvs, and start receiving data
  358. preposts = min(mysize-1, gather_block_size)
  359. head = 0
  360. count = 0
  361. do p=0, mysize-1
  362. if (p .ne. root) then
  363. q = p+1
  364. if (recvcnts(q) > 0) then
  365. count = count + 1
  366. if (count > preposts) then
  367. tail = mod(head,preposts) + 1
  368. call mpi_wait (rcvid(tail), status, ier)
  369. end if
  370. head = mod(head,preposts) + 1
  371. call mpi_irecv ( recvbuf(displs(q)+1), recvcnts(q), &
  372. recvtype, p, mtag, comm, rcvid(head), &
  373. ier )
  374. call mpi_send ( signal, 1, recvtype, p, mtag, comm, ier )
  375. end if
  376. end if
  377. end do
  378. ! copy local data
  379. q = mytid+1
  380. do i=1,sendcnt
  381. recvbuf(displs(q)+i) = sendbuf(i)
  382. enddo
  383. ! wait for final data
  384. do i=1,min(count,preposts)
  385. call mpi_wait (rcvid(i), status, ier)
  386. enddo
  387. else
  388. if (sendcnt > 0) then
  389. call mpi_recv ( signal, 1, sendtype, root, mtag, comm, &
  390. status, ier )
  391. call mpi_rsend ( sendbuf, sendcnt, sendtype, root, mtag, &
  392. comm, ier )
  393. end if
  394. endif
  395. if (ier /= 0) then
  396. call MP_perr_die(myname_,':: (point-to-point implementation)',ier)
  397. end if
  398. else
  399. call mpi_gatherv (sendbuf, sendcnt, sendtype, &
  400. recvbuf, recvcnts, displs, recvtype, &
  401. root, comm, ier)
  402. if (ier /= 0) then
  403. call MP_perr_die(myname_,':: MPI_GATHERV',ier)
  404. end if
  405. endif
  406. return
  407. end subroutine fc_gatherv_int
  408. !BOP -------------------------------------------------------------------
  409. !
  410. ! !IROUTINE: fc_gatherv_fp - Gather an array of type FP
  411. !
  412. ! !DESCRIPTION:
  413. ! This routine gathers a {\em distributed} array of type {\em FP} to
  414. ! the {\tt root} process. Explicit handshaking messages are used
  415. ! to control the number of processes communicating with the root
  416. ! at any one time.
  417. !
  418. ! If flow_cntl optional parameter
  419. ! < 0 : use MPI_Gatherv
  420. ! >= 0: use point-to-point with handshaking messages and
  421. ! preposting receive requests up to
  422. ! min(max(1,flow_cntl),max_gather_block_size)
  423. ! ahead if optional flow_cntl parameter is present.
  424. ! Otherwise, max_gather_block_size is used in its place.
  425. ! Default value is max_gather_block_size.
  426. ! !INTERFACE:
  427. !
  428. subroutine fc_gatherv_fp (sendbuf, sendcnt, sendtype, &
  429. recvbuf, recvcnts, displs, recvtype, &
  430. root, comm, flow_cntl )
  431. !
  432. ! !USES:
  433. !
  434. use m_realkinds, only : FP
  435. use m_die
  436. use m_mpif90
  437. !
  438. ! !INPUT PARAMETERS:
  439. !
  440. real (FP), intent(in) :: sendbuf(*)
  441. integer, intent(in) :: sendcnt
  442. integer, intent(in) :: sendtype
  443. integer, intent(in) :: recvcnts(*)
  444. integer, intent(in) :: displs(*)
  445. integer, intent(in) :: recvtype
  446. integer, intent(in) :: root
  447. integer, intent(in) :: comm
  448. integer, optional, intent(in) :: flow_cntl
  449. ! !OUTPUT PARAMETERS:
  450. !
  451. real (FP), intent(out) :: recvbuf(*)
  452. ! !REVISION HISTORY:
  453. ! 30Jan09 - P.H. Worley - imported from spmd_utils.F90
  454. !EOP ___________________________________________________________________
  455. character(len=*),parameter :: myname_=myname//'::fc_gatherv_fp'
  456. real (FP) :: signal
  457. logical fc_gather ! use explicit flow control?
  458. integer gather_block_size ! number of preposted receive requests
  459. integer :: mytid, mysize, mtag, p, q, i, count
  460. integer :: preposts, head, tail
  461. integer :: rcvid(max_gather_block_size)
  462. integer :: status(MP_STATUS_SIZE)
  463. integer :: ier ! MPI error code
  464. signal = 1.0
  465. if ( present(flow_cntl) ) then
  466. if (flow_cntl >= 0) then
  467. gather_block_size = min(max(1,flow_cntl),max_gather_block_size)
  468. fc_gather = .true.
  469. else
  470. fc_gather = .false.
  471. endif
  472. else
  473. gather_block_size = max(1,max_gather_block_size)
  474. fc_gather = .true.
  475. endif
  476. if (fc_gather) then
  477. call mpi_comm_rank (comm, mytid, ier)
  478. call mpi_comm_size (comm, mysize, ier)
  479. mtag = 0
  480. if (root .eq. mytid) then
  481. ! prepost gather_block_size irecvs, and start receiving data
  482. preposts = min(mysize-1, gather_block_size)
  483. head = 0
  484. count = 0
  485. do p=0, mysize-1
  486. if (p .ne. root) then
  487. q = p+1
  488. if (recvcnts(q) > 0) then
  489. count = count + 1
  490. if (count > preposts) then
  491. tail = mod(head,preposts) + 1
  492. call mpi_wait (rcvid(tail), status, ier)
  493. end if
  494. head = mod(head,preposts) + 1
  495. call mpi_irecv ( recvbuf(displs(q)+1), recvcnts(q), &
  496. recvtype, p, mtag, comm, rcvid(head), &
  497. ier )
  498. call mpi_send ( signal, 1, recvtype, p, mtag, comm, ier )
  499. end if
  500. end if
  501. end do
  502. ! copy local data
  503. q = mytid+1
  504. do i=1,sendcnt
  505. recvbuf(displs(q)+i) = sendbuf(i)
  506. enddo
  507. ! wait for final data
  508. do i=1,min(count,preposts)
  509. call mpi_wait (rcvid(i), status, ier)
  510. enddo
  511. else
  512. if (sendcnt > 0) then
  513. call mpi_recv ( signal, 1, sendtype, root, mtag, comm, &
  514. status, ier )
  515. call mpi_rsend ( sendbuf, sendcnt, sendtype, root, mtag, &
  516. comm, ier )
  517. end if
  518. endif
  519. if (ier /= 0) then
  520. call MP_perr_die(myname_,':: (point-to-point implementation)',ier)
  521. end if
  522. else
  523. call mpi_gatherv (sendbuf, sendcnt, sendtype, &
  524. recvbuf, recvcnts, displs, recvtype, &
  525. root, comm, ier)
  526. if (ier /= 0) then
  527. call MP_perr_die(myname_,':: MPI_GATHERV',ier)
  528. end if
  529. endif
  530. return
  531. end subroutine fc_gatherv_fp
  532. !BOP -------------------------------------------------------------------
  533. !
  534. ! !IROUTINE: get_fcblocksize - return max_gather_block_size
  535. !
  536. ! !DESCRIPTION:
  537. ! This function returns the current value of max_gather_block_size
  538. !
  539. ! !INTERFACE:
  540. function get_fcblocksize()
  541. ! !USES:
  542. !
  543. ! No external modules are used by this function.
  544. implicit none
  545. ! !INPUT PARAMETERS:
  546. !
  547. ! !OUTPUT PARAMETERS:
  548. !
  549. integer :: get_fcblocksize
  550. ! !REVISION HISTORY:
  551. ! 03Mar09 - R. Jacob (jacob@mcs.anl.gov) -- intial version
  552. !EOP ___________________________________________________________________
  553. character(len=*),parameter :: myname_=myname//'::get_fcblocksize'
  554. get_fcblocksize = max_gather_block_size
  555. end function get_fcblocksize
  556. !BOP -------------------------------------------------------------------
  557. !
  558. ! !IROUTINE: set_fcblocksize - set max_gather_block_size
  559. !
  560. ! !DESCRIPTION:
  561. ! This function sets the current value of max_gather_block_size
  562. !
  563. ! !INTERFACE:
  564. subroutine set_fcblocksize(gather_block_size)
  565. ! !USES:
  566. !
  567. ! No external modules are used by this function.
  568. implicit none
  569. ! !INPUT PARAMETERS:
  570. !
  571. integer :: gather_block_size
  572. ! !OUTPUT PARAMETERS:
  573. !
  574. ! !REVISION HISTORY:
  575. ! 03Mar09 - R. Jacob (jacob@mcs.anl.gov) -- intial version
  576. !EOP ___________________________________________________________________
  577. character(len=*),parameter :: myname_=myname//':: set_fcblocksize'
  578. max_gather_block_size = gather_block_size
  579. end subroutine set_fcblocksize
  580. end module m_FcComms