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