PROGRESS  master
prg_parallel_mod.F90
Go to the documentation of this file.
1 
3 !
4 !
5 
7 
8  use bml
9 
10 #ifdef DO_MPI
11  use mpi
12 #endif
13 
14 #ifdef DO_MPI
15 #ifdef SINGLE
16 #define REAL_MPI_TYPE MPI_FLOAT
17 #else
18 #define REAL_MPI_TYPE MPI_DOUBLE
19 #endif
20 #endif
21 
22  implicit none
23 
24  private !Everything is private by default
25 
26  integer, parameter :: dp = kind(1.0d0)
27 
28  integer :: myrank, nranks
29  integer :: ierr, reqcount
30  integer, allocatable :: requestlist(:), rused(:)
31 
32  public :: rankreducedata_t
33  public :: getnranks
34  public :: getmyrank
35  public :: printrank
36  public :: prg_initparallel
37  public :: prg_shutdownparallel
38  public :: prg_barrierparallel
39  public :: sendreceiveparallel
40  public :: isendparallel
41  public :: sendparallel
42  public :: prg_iprg_recvparallel
43  public :: prg_recvparallel
44  public :: sumintparallel
45  public :: sumrealparallel
46  public :: maxintparallel
47  public :: maxrealparallel
48  public :: minintparallel
49  public :: minrealparallel
50  public :: prg_minrealreduce
51  public :: prg_maxrealreduce
52  public :: prg_maxintreduce2
53  public :: prg_sumintreduce2
54  public :: prg_sumintreducen
55  public :: prg_sumrealreduce
56  public :: prg_sumrealreduce2
57  public :: prg_sumrealreduce3
58  public :: prg_sumrealreducen
59  public :: minrankrealparallel
60  public :: maxrankrealparallel
61  public :: prg_bcastparallel
62  public :: prg_bcastrealparallel
63  public :: prg_bcastintparallel
64  public :: allgatherrealparallel
65  public :: allgatherintparallel
66  public :: allgathervrealparallel
67  public :: allgathervintparallel
70  public :: prg_allgatherparallel
71  public :: prg_wait
72 
75 
77  real(dp) :: val
78 
80  integer :: rank
81 
82  end type rankreducedata_t
83 
84 contains
85 
86  !
87  ! Return total number of ranks
88  !
89  function getnranks() result(nR)
90  integer :: nr
91 
92  nr = nranks
93  return
94 
95  end function getnranks
96 
97  !
98  ! Return the local rank
99  !
100  function getmyrank() result(mR)
101  integer :: mr
102 
103  mr = myrank
104 
105  return
106 
107  end function getmyrank
108 
109  !
110  ! Check if this is rank 0 for printing
111  !
112  function printrank() result(pR)
113  integer :: pr
114 
115  pr = 0
116 
117  if (myrank .eq. 0) then
118  pr = 1
119  endif
120 
121  return
122 
123  end function printrank
124 
125  !
126  ! Initialize MPI
127  !
128  subroutine prg_initparallel()
129 
130 #ifdef DO_MPI
131  call mpi_init(ierr)
132  call mpi_comm_rank(mpi_comm_world, myrank, ierr)
133  call mpi_comm_size(mpi_comm_world, nranks, ierr)
134 
135  call bml_initf(mpi_comm_world)
136 
137  ! if (printRank() .eq. 1) then
138  write(*,*) "MPI started in progress, rank ", myrank, " out of ", nranks, " ranks"
139  ! endif
140 
141 #else
142  myrank = 0
143  nranks = 1
144 #endif
145 
146  allocate(requestlist(nranks))
147  allocate(rused(nranks))
148  rused = 0
149 
150  end subroutine prg_initparallel
151 
152  !
153  ! Shutdown MPI
154  !
155  subroutine prg_shutdownparallel()
156 
157  if(allocated(requestlist))deallocate(requestlist)
158  if(allocated(rused))deallocate(rused)
159 
160 #ifdef DO_MPI
161  call bml_shutdownf()
162 
163  call mpi_finalize(ierr)
164 #endif
165 
166  end subroutine prg_shutdownparallel
167 
168  !
169  ! Save request
170  !
171  function saverequest(irequest)
172 
173  integer :: saverequest
174  integer,intent(in) :: irequest
175  integer :: i
176 
177 #ifdef DO_MPI
178  do i = 1, nranks
179  if (rused(i) == 0) then
180  requestlist(i) = irequest
181  rused(i) = 1
182 
183  saverequest = i
184  return
185  endif
186  enddo
187 #endif
188 
189  saverequest = -1
190  return
191 
192  end function saverequest
193 
194  !
195  ! Barrier - all ranks synchronized
196  !
197  subroutine prg_barrierparallel()
198 
199 #ifdef DO_MPI
200  call mpi_barrier(mpi_comm_world, ierr)
201 #endif
202 
203  end subroutine prg_barrierparallel
204 
205  !
206  ! Coordinated send and receive
207  !
208  subroutine sendreceiveparallel(sendBuf, sendLen, dest, recvBuf, recvLen, source, nreceived)
209 
210  real(dp), intent(in) :: sendbuf(*)
211  real(dp), intent(out) :: recvbuf(*)
212  integer, intent(in) :: sendlen, recvlen, dest, source
213  integer, intent(out) :: nreceived
214 #ifdef DO_MPI
215  integer :: mstatus(mpi_status_size)
216 
217  call mpi_sendrecv(sendbuf, sendlen, real_mpi_type, dest, 0, &
218  recvbuf, recvlen, real_mpi_type, source, 0, &
219  mpi_comm_world, mstatus, ierr)
220  call mpi_get_count(mstatus, real_mpi_type, nreceived, ierr)
221 #else
222  recvbuf(1:sendlen) = sendbuf(1:sendlen)
223  nreceived = sendlen
224 #endif
225 
226  end subroutine sendreceiveparallel
227 
228  !
229  ! Non-blocking send
230  !
231  subroutine isendparallel(sendBuf, sendLen, dest)
232 
233  real(dp),intent(in) :: sendbuf(*)
234  integer,intent(in) :: sendlen, dest
235  integer :: request
236 
237 #ifdef DO_MPI
238  call mpi_isend(sendbuf, sendlen, real_mpi_type, dest, &
239  0, mpi_comm_world, request, ierr)
240 #endif
241 
242  end subroutine isendparallel
243 
244  !
245  ! Blocking send
246  !
247  subroutine sendparallel(sendBuf, sendLen, dest)
248 
249  real(dp),intent(in) :: sendbuf(*)
250  integer,intent(in) :: sendlen, dest
251 
252 #ifdef DO_MPI
253  call mpi_send(sendbuf, sendlen, real_mpi_type, dest, &
254  0, mpi_comm_world, ierr)
255 #endif
256 
257  end subroutine sendparallel
258 
259  !
260  ! Non-blocking receive from any other rank
261  !
262  subroutine prg_iprg_recvparallel(recvBuf, recvLen, rind)
263 
264  real(dp) :: recvbuf(*)
265  integer,intent(in) :: recvlen
266  integer :: rind
267 #ifdef DO_MPI
268  integer :: request
269 
270  call mpi_irecv(recvbuf, recvlen, real_mpi_type, mpi_any_source, &
271  0, mpi_comm_world, request, ierr)
272  rind = saverequest(request)
273 #endif
274 
275  end subroutine prg_iprg_recvparallel
276 
277  !
278  ! Blocking receive from any other rank
279  !
280  subroutine prg_recvparallel(recvBuf, recvLen)
281 
282  real(dp) :: recvbuf(*)
283  integer,intent(in) :: recvlen
284 #ifdef DO_MPI
285  integer :: mstatus(mpi_status_size)
286 
287  call mpi_recv(recvbuf, recvlen, real_mpi_type, mpi_any_source, &
288  0, mpi_comm_world, mstatus, ierr)
289 #endif
290 
291  end subroutine prg_recvparallel
292 
293  !
294  ! Integer sum reduce
295  !
296  subroutine sumintparallel(sendBuf, recvBuf, icount)
297 
298  integer,intent(in) :: sendbuf(*)
299  integer :: recvbuf(*)
300  integer,intent(in) :: icount
301  integer :: i
302 
303 #ifdef DO_MPI
304  call mpi_allreduce(sendbuf, recvbuf, icount, mpi_int, &
305  mpi_sum, mpi_comm_world, ierr)
306 #else
307  do i = 1, icount
308  recvbuf(i) = sendbuf(i)
309  enddo
310 #endif
311 
312  end subroutine sumintparallel
313 
314  !
315  ! Real sum reduce
316  !
317  subroutine sumrealparallel(sendBuf, recvBuf, icount)
318 
319  real(dp),intent(in) :: sendbuf(*)
320  real(dp),intent(out) :: recvbuf(*)
321  integer,intent(in) :: icount
322  integer :: i
323 
324 #ifdef DO_MPI
325  call mpi_allreduce(sendbuf, recvbuf, icount, real_mpi_type, &
326  mpi_sum, mpi_comm_world, ierr)
327 #else
328  do i = 1, icount
329  recvbuf(i) = sendbuf(i)
330  enddo
331 #endif
332 
333  end subroutine sumrealparallel
334 
335  !
336  ! Integer max reduce
337  !
338  subroutine maxintparallel(sendBuf, recvBuf, icount)
339 
340  integer,intent(in) :: sendbuf(*)
341  integer,intent(out) :: recvbuf(*)
342  integer,intent(in) :: icount
343  integer :: i
344 
345 #ifdef DO_MPI
346  call mpi_allreduce(sendbuf, recvbuf, icount, mpi_int, &
347  mpi_max, mpi_comm_world, ierr)
348 #else
349  do i = 1, icount
350  recvbuf(i) = sendbuf(i)
351  enddo
352 #endif
353 
354  end subroutine maxintparallel
355 
356  !
357  ! Real max reduce
358  !
359  subroutine maxrealparallel(sendBuf, recvBuf, icount)
360 
361  real(dp),intent(in) :: sendbuf(*)
362  real(dp),intent(out) :: recvbuf(*)
363  integer,intent(in) :: icount
364  integer :: i
365 
366 #ifdef DO_MPI
367  call mpi_allreduce(sendbuf, recvbuf, icount, real_mpi_type, &
368  mpi_max, mpi_comm_world, ierr)
369 #else
370  do i = 1, icount
371  recvbuf(i) = sendbuf(i)
372  enddo
373 #endif
374 
375  end subroutine maxrealparallel
376 
377  !
378  ! Integer min reduce
379  !
380  subroutine minintparallel(sendBuf, recvBuf, icount)
381 
382  integer,intent(in) :: sendbuf(*)
383  integer,intent(out) :: recvbuf(*)
384  integer,intent(in) :: icount
385  integer :: i
386 
387 #ifdef DO_MPI
388  call mpi_allreduce(sendbuf, recvbuf, icount, mpi_int, &
389  mpi_min, mpi_comm_world, ierr)
390 #else
391  do i = 1, icount
392  recvbuf(i) = sendbuf(i)
393  enddo
394 #endif
395 
396  end subroutine minintparallel
397 
398  !
399  ! Real min reduce
400  !
401  subroutine minrealparallel(sendBuf, recvBuf, icount)
402 
403  real(dp),intent(in) :: sendbuf(*)
404  real(dp),intent(out) :: recvbuf(*)
405  integer,intent(in) :: icount
406  integer :: i
407 
408 #ifdef DO_MPI
409  call mpi_allreduce(sendbuf, recvbuf, icount, real_mpi_type, &
410  mpi_min, mpi_comm_world, ierr)
411 #else
412  do i = 1, icount
413  recvbuf(i) = sendbuf(i)
414  enddo
415 #endif
416 
417  end subroutine minrealparallel
418 
419  !
420  ! Real min reduce for 1 value
421  !
422  subroutine prg_minrealreduce(rvalue)
423 
424  real(dp),intent(inout) :: rvalue
425  real(dp) :: slocal(1),sglobal(1)
426 
427  slocal(1) = rvalue
428 
429  call minrealparallel(slocal, sglobal, 1);
430 
431  rvalue = sglobal(1)
432 
433  end subroutine prg_minrealreduce
434 
435  !
436  ! Real max reduce for 1 value
437  !
438  subroutine prg_maxrealreduce(rvalue)
439 
440  real(dp),intent(inout) :: rvalue
441  real(dp) :: slocal(1), sglobal(1)
442 
443  slocal(1) = rvalue
444 
445  call maxrealparallel(slocal, sglobal, 1);
446 
447  rvalue = sglobal(1)
448 
449  end subroutine prg_maxrealreduce
450 
451  !
452  ! Integer max reduce for 2 values
453  !
454  subroutine prg_maxintreduce2(value1, value2)
455 
456  integer,intent(inout) :: value1, value2
457  integer :: slocal(2), sglobal(2)
458 
459  slocal(1) = value1
460  slocal(2) = value2
461 
462  call maxintparallel(slocal, sglobal, 2);
463 
464  value1 = sglobal(1)
465  value2 = sglobal(2)
466 
467  end subroutine prg_maxintreduce2
468 
469  !
470  ! Integer sum reduce for 2 values
471  !
472  subroutine prg_sumintreduce2(value1, value2)
473 
474  integer,intent(inout) :: value1, value2
475  integer :: slocal(2), sglobal(2)
476 
477  slocal(1) = value1
478  slocal(2) = value2
479 
480  call sumintparallel(slocal, sglobal, 2);
481 
482  value1 = sglobal(1)
483  value2 = sglobal(2)
484 
485  end subroutine prg_sumintreduce2
486 
487  !
488  ! Real sum reduce for 1 value
489  !
490  subroutine prg_sumrealreduce(value1)
491 
492  real(dp),intent(inout) :: value1
493  real(dp):: slocal(1), sglobal(1)
494 
495  slocal(1) = value1
496 
497  call sumrealparallel(slocal, sglobal, 1);
498 
499  value1 = sglobal(1)
500 
501  end subroutine prg_sumrealreduce
502 
503  !
504  ! Real sum reduce for 2 values
505  !
506  subroutine prg_sumrealreduce2(value1, value2)
507 
508  real(dp),intent(inout) :: value1, value2
509  real(dp):: slocal(2), sglobal(2)
510 
511  slocal(1) = value1
512  slocal(2) = value2
513 
514  call sumrealparallel(slocal, sglobal, 2);
515 
516  value1 = sglobal(1)
517  value2 = sglobal(2)
518 
519  end subroutine prg_sumrealreduce2
520 
521  !
522  ! Real sum reduce for 3 values
523  !
524  subroutine prg_sumrealreduce3(value1, value2, value3)
525 
526  real(dp),intent(inout) :: value1, value2, value3
527  real(dp):: slocal(3), sglobal(3)
528 
529  slocal(1) = value1
530  slocal(2) = value2
531  slocal(3) = value3
532 
533  call sumrealparallel(slocal, sglobal, 3);
534 
535  value1 = sglobal(1)
536  value2 = sglobal(2)
537  value3 = sglobal(3)
538 
539  end subroutine prg_sumrealreduce3
540 
541  !
542  ! Real sum reduce for N values
543  !
544  subroutine prg_sumrealreducen(valueVec, N)
545 
546  real(dp),intent(inout) :: valuevec(n)
547  integer, intent(in) :: n
548 
549  real(dp), allocatable :: sglobal(:)
550 
551  allocate(sglobal(n))
552 
553  call sumrealparallel(valuevec, sglobal, n);
554 
555  valuevec = sglobal
556 
557  deallocate(sglobal)
558 
559  end subroutine prg_sumrealreducen
560 
561 
562  !
563  ! Real sum reduce for Int values
564  !
565  subroutine prg_sumintreducen(valueVec, N)
566 
567  integer, intent(inout) :: valuevec(n)
568  integer, intent(in) :: n
569 
570  integer, allocatable :: sglobal(:)
571 
572  allocate(sglobal(n))
573 
574  call sumintparallel(valuevec, sglobal, n);
575 
576  valuevec = sglobal
577 
578  deallocate(sglobal)
579 
580  end subroutine prg_sumintreducen
581 
582  !
583  ! Wrapper for MPI_Allreduce real min with rank.
584  !
585  subroutine minrankrealparallel(sendBuf, recvBuf, icount)
586 
587  type(rankreducedata_t), intent(in) :: sendbuf(*)
588  type(rankreducedata_t), intent(out) :: recvbuf(*)
589  integer, intent(in) :: icount
590 
591 #ifdef DO_MPI
592  call mpi_allreduce(sendbuf, recvbuf, icount, mpi_double_int, &
593  mpi_minloc, mpi_comm_world, ierr)
594 #else
595  integer :: i
596 
597  do i = 1, icount
598  recvbuf(i)%val = sendbuf(i)%val
599  recvbuf(i)%rank = sendbuf(i)%rank
600  enddo
601 #endif
602 
603  end subroutine minrankrealparallel
604 
605  !
606  ! Wrapper for MPI_Allreduce real max with rank.
607  !
608  subroutine maxrankrealparallel(sendBuf, recvBuf, icount)
609 
610  type(rankreducedata_t), intent(in) :: sendbuf(*)
611  type(rankreducedata_t), intent(out) :: recvbuf(*)
612  integer, intent(in) :: icount
613 
614 #ifdef DO_MPI
615  call mpi_allreduce(sendbuf, recvbuf, icount, mpi_double_int, &
616  mpi_maxloc, mpi_comm_world, ierr)
617 #else
618  integer :: i
619 
620  do i = 1, icount
621  recvbuf(i)%val = sendbuf(i)%val
622  recvbuf(i)%rank = sendbuf(i)%rank
623  enddo
624 #endif
625 
626  end subroutine maxrankrealparallel
627 
628  !
629  ! Wrapper for MPI broadcast
630  !
631  subroutine prg_bcastparallel(buf, blen, root)
632 
633  character, intent(inout) :: buf(*)
634  integer, intent(in) :: blen, root
635 
636 #ifdef DO_MPI
637  call mpi_bcast(buf, blen, mpi_byte, root, mpi_comm_world, ierr)
638 #endif
639 
640  end subroutine prg_bcastparallel
641 
642  !
643  ! Wrapper for MPI broadcast
644  !
645  subroutine prg_bcastrealparallel(buf, blen, root)
646 
647  real(dp), intent(inout) :: buf(*)
648  integer, intent(in) :: blen, root
649 
650 #ifdef DO_MPI
651  call mpi_bcast(buf, blen, real_mpi_type, root, mpi_comm_world, ierr)
652 #endif
653 
654  end subroutine prg_bcastrealparallel
655 
656  !
657  ! Wrapper for MPI broadcast
658  !
659  subroutine prg_bcastintparallel(buf, blen, root)
660 
661  integer, intent(inout) :: buf(*)
662  integer, intent(in) :: blen, root
663 
664 #ifdef DO_MPI
665  call mpi_bcast(buf, blen, mpi_int, root, mpi_comm_world, ierr)
666 #endif
667 
668  end subroutine prg_bcastintparallel
669 
670  !
671  ! Wrapper for real MPI_AllGather
672  !
673  subroutine allgatherrealparallel(sendBuf, sendLen, recvBuf, recvLen)
674 
675  real(dp), intent(in) :: sendbuf(*)
676  real(dp), intent(out) :: recvbuf(*)
677  integer, intent(in) :: sendlen, recvlen
678 
679 #ifdef DO_MPI
680  call mpi_allgather(sendbuf, sendlen, real_mpi_type, recvbuf, recvlen, &
681  real_mpi_type, mpi_comm_world, ierr)
682 #endif
683 
684  end subroutine allgatherrealparallel
685 
686  !
687  ! Wrapper for integer MPI_AllGather
688  !
689  subroutine allgatherintparallel(sendBuf, sendLen, recvBuf, recvLen)
690 
691  integer, intent(in) :: sendbuf(*)
692  integer, intent(out) :: recvbuf(*)
693  integer, intent(in) :: sendlen, recvlen
694 
695 #ifdef DO_MPI
696  call mpi_allgather(sendbuf, sendlen, mpi_int, recvbuf, recvlen, &
697  mpi_int, mpi_comm_world, ierr)
698 #endif
699 
700  end subroutine allgatherintparallel
701 
702  !
703  ! Wrapper for real MPI_AllGatherV
704  !
705  subroutine allgathervrealparallel(sendBuf, sendLen, recvBuf, recvLen, recvDispl)
706 
707  real(dp), intent(in) :: sendbuf(*)
708  real(dp), intent(out) :: recvbuf(*)
709  integer, intent(in) :: sendlen
710  integer, intent(in) :: recvlen(*)
711  integer, intent(in) :: recvdispl(*)
712 
713 #ifdef DO_MPI
714  call mpi_allgatherv(sendbuf, sendlen, real_mpi_type, recvbuf, &
715  recvlen, recvdispl, real_mpi_type, &
716  mpi_comm_world, ierr)
717 #endif
718 
719 
720  end subroutine allgathervrealparallel
721 
722  !
723  ! Wrapper for integer MPI_AllGatherV
724  !
725  subroutine allgathervintparallel(sendBuf, sendLen, recvBuf, recvLen, recvDispl)
726 
727  integer, intent(in) :: sendbuf(*)
728  integer, intent(out) :: recvbuf(*)
729  integer, intent(in) :: sendlen
730  integer, intent(in) :: recvlen(*)
731  integer, intent(in) :: recvdispl(*)
732 
733 #ifdef DO_MPI
734  call mpi_allgatherv(sendbuf, sendlen, mpi_int, recvbuf, recvlen, &
735  recvdispl, mpi_int, mpi_comm_world, ierr)
736 #endif
737 
738  end subroutine allgathervintparallel
739 
740  !
741  ! Wrapper for real MPI_AllReduce
742  !
743  subroutine prg_allsumrealreduceparallel(buf, buflen)
744 
745  real(dp), intent(inout) :: buf(*)
746  integer, intent(in) :: buflen
747 
748 #ifdef DO_MPI
749  call mpi_allreduce(mpi_in_place, buf, buflen, real_mpi_type, &
750  mpi_sum, mpi_comm_world, ierr)
751 #endif
752 
753  end subroutine prg_allsumrealreduceparallel
754 
755  !
756  ! Wrapper for integer MPI_AllReduce
757  !
758  subroutine prg_allsumintreduceparallel(buf, buflen)
759 
760  integer, intent(inout) :: buf(*)
761  integer, intent(in) :: buflen
762 
763 #ifdef DO_MPI
764  call mpi_allreduce(mpi_in_place, buf, buflen, mpi_int, &
765  mpi_sum, mpi_comm_world, ierr)
766 #endif
767 
768  end subroutine prg_allsumintreduceparallel
769 
770  !
771  ! Wrapper for bml_allGatherVParallel
772  !
773  subroutine prg_allgatherparallel(a)
774 
775  type (bml_matrix_t), intent(inout) :: a
776 
777 #ifdef DO_MPI
778  call bml_allgathervparallel(a)
779 #endif
780 
781  end subroutine prg_allgatherparallel
782 
783 
784  !
785  ! Wrapper for MPI_prg_wait
786  !
787  subroutine prg_wait()
788 
789  integer :: status(3)
790  integer :: request, ierr
791 
792 #ifdef DO_MPI
793  ! call MPI_WAIT(request, status, ierr)
794 #endif
795 
796  end subroutine prg_wait
797 
798 
799 
800 end module prg_parallel_mod
The parallel module.
subroutine, public prg_iprg_recvparallel(recvBuf, recvLen, rind)
subroutine, public minintparallel(sendBuf, recvBuf, icount)
subroutine, public maxrankrealparallel(sendBuf, recvBuf, icount)
subroutine, public prg_barrierparallel()
integer, dimension(:), allocatable requestlist
subroutine, public prg_maxintreduce2(value1, value2)
integer function saverequest(irequest)
subroutine, public prg_sumrealreduce3(value1, value2, value3)
subroutine, public minrankrealparallel(sendBuf, recvBuf, icount)
subroutine, public sumrealparallel(sendBuf, recvBuf, icount)
subroutine, public prg_sumrealreduce(value1)
subroutine, public prg_sumrealreducen(valueVec, N)
subroutine, public allgatherrealparallel(sendBuf, sendLen, recvBuf, recvLen)
subroutine, public allgatherintparallel(sendBuf, sendLen, recvBuf, recvLen)
subroutine, public prg_maxrealreduce(rvalue)
subroutine, public sendparallel(sendBuf, sendLen, dest)
subroutine, public prg_shutdownparallel()
subroutine, public allgathervrealparallel(sendBuf, sendLen, recvBuf, recvLen, recvDispl)
subroutine, public isendparallel(sendBuf, sendLen, dest)
subroutine, public prg_allgatherparallel(a)
subroutine, public prg_bcastparallel(buf, blen, root)
integer function, public printrank()
subroutine, public maxrealparallel(sendBuf, recvBuf, icount)
subroutine, public prg_bcastintparallel(buf, blen, root)
subroutine, public prg_sumintreduce2(value1, value2)
subroutine, public sendreceiveparallel(sendBuf, sendLen, dest, recvBuf, recvLen, source, nreceived)
integer, parameter dp
subroutine, public prg_bcastrealparallel(buf, blen, root)
subroutine, public allgathervintparallel(sendBuf, sendLen, recvBuf, recvLen, recvDispl)
subroutine, public prg_initparallel()
subroutine, public maxintparallel(sendBuf, recvBuf, icount)
subroutine, public prg_allsumrealreduceparallel(buf, buflen)
subroutine, public prg_allsumintreduceparallel(buf, buflen)
subroutine, public prg_sumintreducen(valueVec, N)
integer function, public getnranks()
subroutine, public minrealparallel(sendBuf, recvBuf, icount)
subroutine, public prg_minrealreduce(rvalue)
subroutine, public prg_recvparallel(recvBuf, recvLen)
integer, dimension(:), allocatable rused
integer function, public getmyrank()
subroutine, public prg_wait()
subroutine, public sumintparallel(sendBuf, recvBuf, icount)
subroutine, public prg_sumrealreduce2(value1, value2)
Data structure for rection over MPI ranks.