28 public :: comm_data_transfer
32 interface comm_data_transfer
35 end interface comm_data_transfer
40 end interface comm_var
42 public :: comm_stat_sum
43 public :: comm_stat_sum_eachlayer
44 public :: comm_stat_avg
45 public :: comm_stat_max
46 public :: comm_stat_min
48 interface comm_stat_sum
51 end interface comm_stat_sum
53 interface comm_stat_sum_eachlayer
56 end interface comm_stat_sum_eachlayer
58 interface comm_stat_avg
61 end interface comm_stat_avg
63 interface comm_stat_max
66 end interface comm_stat_max
68 interface comm_stat_min
71 end interface comm_stat_min
132 private :: comm_list_generate
133 private :: comm_sortdest
134 private :: comm_sortdest_pl
135 private :: comm_sortdest_singular
137 private :: comm_debugtest
143 logical,
private :: comm_apply_barrier = .false.
144 integer,
private :: comm_varmax = 50
145 logical,
private :: debug = .false.
146 logical,
private :: testonly = .false.
149 integer,
private,
parameter :: rellist_vindex = 6
150 integer,
private,
parameter :: i_recv_grid = 1
151 integer,
private,
parameter :: i_recv_rgn = 2
152 integer,
private,
parameter :: i_recv_prc = 3
153 integer,
private,
parameter :: i_send_grid = 4
154 integer,
private,
parameter :: i_send_rgn = 5
155 integer,
private,
parameter :: i_send_prc = 6
157 integer,
private,
allocatable :: rellist(:,:)
158 integer,
private :: rellist_nmax
161 integer,
private,
parameter :: recv_nlim = 20
162 integer,
private,
parameter :: send_nlim = 20
164 integer,
private :: copy_nmax_r2r = 0
165 integer,
private :: recv_nmax_r2r = 0
166 integer,
private :: send_nmax_r2r = 0
168 integer,
private :: copy_nmax_p2r = 0
169 integer,
private :: recv_nmax_p2r = 0
170 integer,
private :: send_nmax_p2r = 0
172 integer,
private :: copy_nmax_r2p = 0
173 integer,
private :: recv_nmax_r2p = 0
174 integer,
private :: send_nmax_r2p = 0
176 integer,
private :: singular_nmax = 0
178 integer,
private,
parameter :: info_vindex = 3
179 integer,
private,
parameter :: i_size = 1
180 integer,
private,
parameter :: i_prc_from = 2
181 integer,
private,
parameter :: i_prc_to = 3
184 integer,
private,
parameter :: list_vindex = 4
185 integer,
private,
parameter :: i_grid_from = 1
186 integer,
private,
parameter :: i_l_from = 2
187 integer,
private,
parameter :: i_grid_to = 3
188 integer,
private,
parameter :: i_l_to = 4
191 integer,
private :: req_count
206 namelist / param_comm_icoa / &
207 comm_apply_barrier, &
217 if(
io_l )
write(
io_fid_log,*)
'+++ Module[comm]/Category[common share]'
221 if(
io_l )
write(
io_fid_log,*)
'*** PARAM_COMM_ICOA is not specified. use default.'
222 elseif( ierr > 0 )
then
223 write(*,*)
'xxx Not appropriate names in namelist PARAM_COMM_ICOA. STOP.'
230 elseif(
rp ==
sp )
then
233 write(*,*)
'xxx precision is not supportd'
243 if(
io_l )
write(
io_fid_log,*)
'====== communication information ======'
245 call comm_list_generate
248 call comm_sortdest_pl
249 call comm_sortdest_singular
251 allocate(
req_list( recv_nmax_r2r + send_nmax_r2r &
252 + recv_nmax_p2r + send_nmax_p2r &
253 + recv_nmax_r2p + send_nmax_r2p ) )
255 if( testonly )
call comm_debugtest
262 subroutine comm_list_generate
290 integer :: prc, prc_rmt
291 integer :: rgnid, rgnid_rmt
292 integer :: i, j, i_rmt, j_rmt
321 rellist(i_recv_grid,cnt) =
suf(i,j)
322 rellist(i_recv_rgn, cnt) = rgnid
323 rellist(i_recv_prc, cnt) = prc
324 rellist(i_send_grid,cnt) =
suf(i_rmt,j_rmt)
325 rellist(i_send_rgn, cnt) = rgnid_rmt
326 rellist(i_send_prc, cnt) = prc_rmt
343 rellist(i_recv_grid,cnt) =
suf(i,j)
344 rellist(i_recv_rgn, cnt) = rgnid
345 rellist(i_recv_prc, cnt) = prc
346 rellist(i_send_grid,cnt) =
suf(i_rmt,j_rmt)
347 rellist(i_send_rgn, cnt) = rgnid_rmt
348 rellist(i_send_prc, cnt) = prc_rmt
367 rellist(i_recv_grid,cnt) =
suf(i,j)
368 rellist(i_recv_rgn, cnt) = rgnid
369 rellist(i_recv_prc, cnt) = prc
370 rellist(i_send_grid,cnt) =
suf(i_rmt,j_rmt)
371 rellist(i_send_rgn, cnt) = rgnid_rmt
372 rellist(i_send_prc, cnt) = prc_rmt
389 rellist(i_recv_grid,cnt) =
suf(i,j)
390 rellist(i_recv_rgn, cnt) = rgnid
391 rellist(i_recv_prc, cnt) = prc
392 rellist(i_send_grid,cnt) =
suf(i_rmt,j_rmt)
393 rellist(i_send_rgn, cnt) = rgnid_rmt
394 rellist(i_send_prc, cnt) = prc_rmt
413 rellist(i_recv_grid,cnt) =
suf(i,j)
414 rellist(i_recv_rgn, cnt) = rgnid
415 rellist(i_recv_prc, cnt) = prc
416 rellist(i_send_grid,cnt) =
suf(i_rmt,j_rmt)
417 rellist(i_send_rgn, cnt) = rgnid_rmt
418 rellist(i_send_prc, cnt) = prc_rmt
435 rellist(i_recv_grid,cnt) =
suf(i,j)
436 rellist(i_recv_rgn, cnt) = rgnid
437 rellist(i_recv_prc, cnt) = prc
438 rellist(i_send_grid,cnt) =
suf(i_rmt,j_rmt)
439 rellist(i_send_rgn, cnt) = rgnid_rmt
440 rellist(i_send_prc, cnt) = prc_rmt
459 rellist(i_recv_grid,cnt) =
suf(i,j)
460 rellist(i_recv_rgn, cnt) = rgnid
461 rellist(i_recv_prc, cnt) = prc
462 rellist(i_send_grid,cnt) =
suf(i_rmt,j_rmt)
463 rellist(i_send_rgn, cnt) = rgnid_rmt
464 rellist(i_send_prc, cnt) = prc_rmt
481 rellist(i_recv_grid,cnt) =
suf(i,j)
482 rellist(i_recv_rgn, cnt) = rgnid
483 rellist(i_recv_prc, cnt) = prc
484 rellist(i_send_grid,cnt) =
suf(i_rmt,j_rmt)
485 rellist(i_send_rgn, cnt) = rgnid_rmt
486 rellist(i_send_prc, cnt) = prc_rmt
513 rellist(i_recv_grid,cnt) =
suf(i,j)
514 rellist(i_recv_rgn, cnt) = rgnid
515 rellist(i_recv_prc, cnt) = prc
516 rellist(i_send_grid,cnt) =
suf(i_rmt,j_rmt)
517 rellist(i_send_rgn, cnt) = rgnid_rmt
518 rellist(i_send_prc, cnt) = prc_rmt
535 rellist(i_recv_grid,cnt) =
suf(i,j)
536 rellist(i_recv_rgn, cnt) = rgnid
537 rellist(i_recv_prc, cnt) = prc
538 rellist(i_send_grid,cnt) =
suf(i_rmt,j_rmt)
539 rellist(i_send_rgn, cnt) = rgnid_rmt
540 rellist(i_send_prc, cnt) = prc_rmt
563 rellist(i_recv_grid,cnt) =
suf(i,j)
564 rellist(i_recv_rgn, cnt) = rgnid
565 rellist(i_recv_prc, cnt) = prc
566 rellist(i_send_grid,cnt) =
suf(i_rmt,j_rmt)
567 rellist(i_send_rgn, cnt) = rgnid_rmt
568 rellist(i_send_prc, cnt) = prc_rmt
584 rellist(i_recv_grid,cnt) =
suf(i,j)
585 rellist(i_recv_rgn, cnt) = rgnid
586 rellist(i_recv_prc, cnt) = prc
587 rellist(i_send_grid,cnt) =
suf(i_rmt,j_rmt)
588 rellist(i_send_rgn, cnt) = rgnid_rmt
589 rellist(i_send_prc, cnt) = prc_rmt
607 rellist(i_recv_grid,cnt) =
suf(i,j)
608 rellist(i_recv_rgn, cnt) = rgnid
609 rellist(i_recv_prc, cnt) = prc
610 rellist(i_send_grid,cnt) =
suf(i_rmt,j_rmt)
611 rellist(i_send_rgn, cnt) = rgnid_rmt
612 rellist(i_send_prc, cnt) = prc_rmt
635 rellist(i_recv_grid,cnt) =
suf(i,j)
636 rellist(i_recv_rgn, cnt) = rgnid
637 rellist(i_recv_prc, cnt) = prc
638 rellist(i_send_grid,cnt) =
suf(i_rmt,j_rmt)
639 rellist(i_send_rgn, cnt) = rgnid_rmt
640 rellist(i_send_prc, cnt) = prc_rmt
651 if(
io_l )
write(
io_fid_log,
'(7(A10))')
'Count',
'|recv_grid',
'| recv_rgn',
'| recv_prc', &
652 '|send_grid',
'| send_rgn',
'| send_prc'
653 do cnt = 1, rellist_nmax
659 end subroutine comm_list_generate
663 subroutine comm_sortdest
674 integer :: sendbuf1(1)
675 integer :: recvbuf1(1)
677 integer,
allocatable :: sendbuf_info(:)
678 integer,
allocatable :: recvbuf_info(:)
679 integer,
allocatable :: sendbuf_list(:,:,:)
680 integer,
allocatable :: recvbuf_list(:,:,:)
681 integer,
allocatable :: req_list_r2r(:)
683 integer :: recv_nglobal_r2r
684 integer :: send_size_nglobal
686 integer :: cnt, irank, ipos
687 integer :: totalsize, rank, tag
689 integer :: i_from, j_from, r_from, g_from, l_from, p_from
690 integer :: i_to, j_to, r_to, g_to, l_to, p_to
707 allocate(
recv_list_r2r(list_vindex,rellist_nmax,recv_nlim) )
708 allocate(
send_list_r2r(list_vindex,rellist_nmax,send_nlim) )
714 do cnt = 1, rellist_nmax
716 if ( rellist(i_recv_prc,cnt) == rellist(i_send_prc,cnt) )
then
727 do n = 1, recv_nmax_r2r
728 if (
recv_info_r2r(i_prc_from,n) == rellist(i_send_prc,cnt) )
then
734 if ( irank < 0 )
then
735 recv_nmax_r2r = recv_nmax_r2r + 1
736 irank = recv_nmax_r2r
745 recv_list_r2r(i_grid_from,ipos,irank) = rellist(i_send_grid,cnt)
747 recv_list_r2r(i_grid_to ,ipos,irank) = rellist(i_recv_grid,cnt)
761 sendbuf1(1) = recv_nmax_r2r
763 call mpi_allreduce( sendbuf1(1), &
771 recv_nglobal_r2r = recvbuf1(1)
774 allocate( sendbuf_info(info_vindex*recv_nglobal_r2r) )
775 allocate( recvbuf_info(info_vindex*recv_nglobal_r2r*
prc_nprocs) )
778 do irank = 1, recv_nmax_r2r
779 n = (irank-1) * info_vindex
786 totalsize = info_vindex * recv_nglobal_r2r
788 if ( totalsize > 0 )
then
789 call mpi_allgather( sendbuf_info(1), &
799 send_size_nglobal = 0
803 n = (p-1) * info_vindex
805 if ( recvbuf_info(n+i_prc_from) ==
prc_myrank )
then
806 send_nmax_r2r = send_nmax_r2r + 1
807 irank = send_nmax_r2r
814 send_size_nglobal = max( send_size_nglobal, recvbuf_info(n+i_size) )
818 if(
io_l )
write(
io_fid_log,*)
'*** Recv_nmax_r2r(global) = ', recv_nglobal_r2r
819 if(
io_l )
write(
io_fid_log,*)
'*** Recv_nmax_r2r(local) = ', recv_nmax_r2r
820 if(
io_l )
write(
io_fid_log,*)
'*** Send_nmax_r2r(local) = ', send_nmax_r2r
821 if(
io_l )
write(
io_fid_log,*)
'*** Send_size_r2r(global) = ', send_size_nglobal
823 if(
io_l )
write(
io_fid_log,
'(A)')
'|---------------------------------------'
826 do irank = 1, recv_nmax_r2r
829 do irank = 1, send_nmax_r2r
834 allocate( req_list_r2r(recv_nmax_r2r+send_nmax_r2r) )
836 allocate( sendbuf_list(list_vindex,send_size_nglobal,recv_nmax_r2r) )
837 allocate( recvbuf_list(list_vindex,send_size_nglobal,send_nmax_r2r) )
838 sendbuf_list(:,:,:) = -1
842 do irank = 1, send_nmax_r2r
843 req_count = req_count + 1
848 call mpi_irecv( recvbuf_list(1,1,irank), &
854 req_list_r2r(req_count), &
858 do irank = 1, recv_nmax_r2r
863 req_count = req_count + 1
868 call mpi_isend( sendbuf_list(1,1,irank), &
874 req_list_r2r(req_count), &
879 if ( recv_nmax_r2r+send_nmax_r2r > 0 )
then
880 call mpi_waitall( recv_nmax_r2r+send_nmax_r2r, &
882 mpi_statuses_ignore, &
886 do irank = 1, send_nmax_r2r
897 '|ifrom',
'|jfrom',
'|rfrom',
'|gfrom',
'|lfrom',
'|pfrom', &
898 '| ito',
'| jto',
'| rto',
'| gto',
'| lto',
'| pto'
913 if(
io_l )
write(
io_fid_log,
'(13(I6))') ipos, i_from, j_from, r_from, g_from, l_from, p_from, &
914 i_to , j_to , r_to , g_to , l_to , p_to
919 do irank = 1, recv_nmax_r2r
921 '|ifrom',
'|jfrom',
'|rfrom',
'|gfrom',
'|lfrom',
'|pfrom', &
922 '| ito',
'| jto',
'| rto',
'| gto',
'| lto',
'| pto'
937 if(
io_l )
write(
io_fid_log,
'(13(I6))') ipos, i_from, j_from, r_from, g_from, l_from, p_from, &
938 i_to , j_to , r_to , g_to , l_to , p_to
944 do irank = 1, send_nmax_r2r
946 '|ifrom',
'|jfrom',
'|rfrom',
'|gfrom',
'|lfrom',
'|pfrom', &
947 '| ito',
'| jto',
'| rto',
'| gto',
'| lto',
'| pto'
962 if(
io_l )
write(
io_fid_log,
'(13(I6))') ipos, i_from, j_from, r_from, g_from, l_from, p_from, &
963 i_to , j_to , r_to , g_to , l_to , p_to
974 end subroutine comm_sortdest
978 subroutine comm_sortdest_pl
1000 integer :: prc, prc_rmt
1001 integer :: rgnid, rgnid_rmt
1002 logical :: check_vert_pl
1005 integer :: irank, ipos
1007 integer :: send_size_nglobal_pl
1009 integer :: l, l_pl, n, v, vv
1010 integer :: i_from, j_from, r_from, g_from, l_from, p_from
1011 integer :: i_to, j_to, r_to, g_to, l_to, p_to
1027 allocate(
recv_list_p2r(list_vindex,send_size_nglobal_pl,recv_nlim) )
1028 allocate(
send_list_p2r(list_vindex,send_size_nglobal_pl,send_nlim) )
1044 allocate(
recv_list_r2p(list_vindex,send_size_nglobal_pl,recv_nlim) )
1045 allocate(
send_list_r2p(list_vindex,send_size_nglobal_pl,send_nlim) )
1059 if ( rgnid_rmt ==
i_npl )
then
1065 elseif( rgnid_rmt ==
i_spl )
then
1073 if ( check_vert_pl )
then
1084 if ( prc == prc_rmt )
then
1096 do n = 1, recv_nmax_p2r
1103 if ( irank < 0 )
then
1104 recv_nmax_p2r = recv_nmax_p2r + 1
1105 irank = recv_nmax_p2r
1121 do n = 1, send_nmax_r2p
1128 if ( irank < 0 )
then
1129 send_nmax_r2p = send_nmax_r2p + 1
1130 irank = send_nmax_r2p
1163 if ( rgnid ==
i_npl )
then
1168 elseif( rgnid ==
i_spl )
then
1178 if ( prc == prc_rmt )
then
1190 do n = 1, recv_nmax_r2p
1197 if ( irank < 0 )
then
1198 recv_nmax_r2p = recv_nmax_r2p + 1
1199 irank = recv_nmax_r2p
1215 do n = 1, send_nmax_p2r
1222 if ( irank < 0 )
then
1223 send_nmax_p2r = send_nmax_p2r + 1
1224 irank = send_nmax_p2r
1258 if(
io_l )
write(
io_fid_log,*)
'*** Recv_nmax_p2r(local) = ', recv_nmax_p2r
1259 if(
io_l )
write(
io_fid_log,*)
'*** Send_nmax_p2r(local) = ', send_nmax_p2r
1261 if(
io_l )
write(
io_fid_log,
'(A)')
'|---------------------------------------'
1264 do irank = 1, recv_nmax_p2r
1267 do irank = 1, send_nmax_p2r
1276 '|rfrom',
'|gfrom',
'|lfrom',
'|pfrom', &
1277 '| ito',
'| jto',
'| rto',
'| gto',
'| lto',
'| pto'
1290 if(
io_l )
write(
io_fid_log,
'(11(I6))') ipos, r_from, g_from, l_from, p_from, &
1291 i_to , j_to , r_to , g_to , l_to , p_to
1296 do irank = 1, recv_nmax_p2r
1298 '|rfrom',
'|gfrom',
'|lfrom',
'|pfrom', &
1299 '| ito',
'| jto',
'| rto',
'| gto',
'| lto',
'| pto'
1312 if(
io_l )
write(
io_fid_log,
'(11(I6))') ipos, r_from, g_from, l_from, p_from, &
1313 i_to , j_to , r_to , g_to , l_to , p_to
1319 do irank = 1, send_nmax_p2r
1321 '|rfrom',
'|gfrom',
'|lfrom',
'|pfrom', &
1322 '| ito',
'| jto',
'| rto',
'| gto',
'| lto',
'| pto'
1335 if(
io_l )
write(
io_fid_log,
'(11(I6))') ipos, r_from, g_from, l_from, p_from, &
1336 i_to , j_to , r_to , g_to , l_to , p_to
1342 if(
io_l )
write(
io_fid_log,*)
'*** Recv_nmax_r2p(local) = ', recv_nmax_r2p
1343 if(
io_l )
write(
io_fid_log,*)
'*** Send_nmax_r2p(local) = ', send_nmax_r2p
1345 if(
io_l )
write(
io_fid_log,
'(A)')
'|---------------------------------------'
1348 do irank = 1, recv_nmax_r2p
1351 do irank = 1, send_nmax_r2p
1360 '|ifrom',
'|jfrom',
'|rfrom',
'|gfrom',
'|lfrom',
'|pfrom', &
1361 '| rto',
'| gto',
'| lto',
'| pto'
1374 if(
io_l )
write(
io_fid_log,
'(11(I6))') ipos, i_from, j_from, r_from, g_from, l_from, p_from, &
1375 r_to , g_to , l_to , p_to
1380 do irank = 1, recv_nmax_r2p
1382 '|ifrom',
'|jfrom',
'|rfrom',
'|gfrom',
'|lfrom',
'|pfrom', &
1383 '| rto',
'| gto',
'| lto',
'| pto'
1396 if(
io_l )
write(
io_fid_log,
'(11(I6))') ipos, i_from, j_from, r_from, g_from, l_from, p_from, &
1397 r_to , g_to , l_to , p_to
1403 do irank = 1, send_nmax_r2p
1405 '|ifrom',
'|jfrom',
'|rfrom',
'|gfrom',
'|lfrom',
'|pfrom', &
1406 '| rto',
'| gto',
'| lto',
'| pto'
1419 if(
io_l )
write(
io_fid_log,
'(11(I6))') ipos, i_from, j_from, r_from, g_from, l_from, p_from, &
1420 r_to , g_to , l_to , p_to
1425 if(
io_l )
write(
io_fid_log,*)
'*** Send_size_p2r,r2p = ', send_size_nglobal_pl
1437 end subroutine comm_sortdest_pl
1441 subroutine comm_sortdest_singular
1458 integer :: i, j, i_rmt, j_rmt
1461 integer :: i_from, j_from, r_from, g_from, l_from, p_from
1462 integer :: i_to, j_to, r_to, g_to, l_to, p_to
1529 if(
io_l )
write(
io_fid_log,
'(A)')
'|---------------------------------------'
1538 '|ifrom',
'|jfrom',
'|rfrom',
'|gfrom',
'|lfrom',
'|pfrom', &
1539 '| ito',
'| jto',
'| rto',
'| gto',
'| lto',
'| pto'
1554 if(
io_l )
write(
io_fid_log,
'(13(I6))') ipos, i_from, j_from, r_from, g_from, l_from, p_from, &
1555 i_to , j_to , r_to , g_to , l_to , p_to
1560 end subroutine comm_sortdest_singular
1572 real(SP),
intent(inout) :: var (:,:,:,:)
1573 real(SP),
intent(inout) :: var_pl(:,:,:,:)
1575 integer :: shp(4), kmax, vmax
1576 integer :: totalsize, rank, tag
1577 integer :: irank, ipos, imax
1578 integer :: ij_from, l_from, ij_to, l_to
1580 integer :: k, v, ikv
1587 if ( comm_apply_barrier )
then
1601 if ( kmax * vmax >
adm_kall * comm_varmax )
then
1602 write(*,*)
'xxx [COMM_data_transfer] kmax * vmax exceeds ADM_kall * COMM_varmax, stop!'
1603 write(*,*)
'xxx kmax * vmax = ', kmax * vmax
1604 write(*,*)
'xxx ADM_kall * COMM_varmax = ',
adm_kall * comm_varmax
1628 do irank = 1, recv_nmax_r2r
1629 req_count = req_count + 1
1645 do irank = 1, recv_nmax_p2r
1646 req_count = req_count + 1
1662 do irank = 1, recv_nmax_r2p
1663 req_count = req_count + 1
1683 do irank = 1, send_nmax_r2r
1699 ikv = (v-1) * imax * kmax &
1710 req_count = req_count + 1
1711 totalsize = imax * kmax * vmax
1727 do irank = 1, send_nmax_p2r
1743 ikv = (v-1) * imax * kmax &
1754 req_count = req_count + 1
1755 totalsize = imax * kmax * vmax
1771 do irank = 1, send_nmax_r2p
1787 ikv = (v-1) * imax * kmax &
1798 req_count = req_count + 1
1799 totalsize = imax * kmax * vmax
1826 do irank = 1, copy_nmax_r2r
1842 var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
1851 do irank = 1, copy_nmax_p2r
1867 var(ij_to,k,l_to,v) = var_pl(ij_from,k,l_from,v)
1876 do irank = 1, copy_nmax_r2p
1892 var_pl(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
1903 if ( req_count > 0 )
then
1904 call mpi_waitall( req_count, &
1906 mpi_statuses_ignore, &
1923 do irank = 1, recv_nmax_r2r
1937 ikv = (v-1) * imax * kmax &
1950 do irank = 1, recv_nmax_p2r
1964 ikv = (v-1) * imax * kmax &
1977 do irank = 1, recv_nmax_r2p
1991 ikv = (v-1) * imax * kmax &
2007 do irank = 1, singular_nmax
2023 var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2051 real(DP),
intent(inout) :: var (:,:,:,:)
2052 real(DP),
intent(inout) :: var_pl(:,:,:,:)
2054 integer :: shp(4), kmax, vmax
2055 integer :: totalsize, rank, tag
2056 integer :: irank, ipos, imax
2057 integer :: ij_from, l_from, ij_to, l_to
2059 integer :: k, v, ikv
2066 if ( comm_apply_barrier )
then
2080 if ( kmax * vmax >
adm_kall * comm_varmax )
then
2081 write(*,*)
'xxx [COMM_data_transfer] kmax * vmax exceeds ADM_kall * COMM_varmax, stop!'
2082 write(*,*)
'xxx kmax * vmax = ', kmax * vmax
2083 write(*,*)
'xxx ADM_kall * COMM_varmax = ',
adm_kall * comm_varmax
2107 do irank = 1, recv_nmax_r2r
2108 req_count = req_count + 1
2115 mpi_double_precision, &
2124 do irank = 1, recv_nmax_p2r
2125 req_count = req_count + 1
2132 mpi_double_precision, &
2141 do irank = 1, recv_nmax_r2p
2142 req_count = req_count + 1
2149 mpi_double_precision, &
2162 do irank = 1, send_nmax_r2r
2178 ikv = (v-1) * imax * kmax &
2189 req_count = req_count + 1
2190 totalsize = imax * kmax * vmax
2197 mpi_double_precision, &
2206 do irank = 1, send_nmax_p2r
2222 ikv = (v-1) * imax * kmax &
2233 req_count = req_count + 1
2234 totalsize = imax * kmax * vmax
2241 mpi_double_precision, &
2250 do irank = 1, send_nmax_r2p
2266 ikv = (v-1) * imax * kmax &
2277 req_count = req_count + 1
2278 totalsize = imax * kmax * vmax
2285 mpi_double_precision, &
2305 do irank = 1, copy_nmax_r2r
2321 var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2330 do irank = 1, copy_nmax_p2r
2346 var(ij_to,k,l_to,v) = var_pl(ij_from,k,l_from,v)
2355 do irank = 1, copy_nmax_r2p
2371 var_pl(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2382 if ( req_count > 0 )
then
2383 call mpi_waitall( req_count, &
2385 mpi_statuses_ignore, &
2402 do irank = 1, recv_nmax_r2r
2416 ikv = (v-1) * imax * kmax &
2429 do irank = 1, recv_nmax_p2r
2443 ikv = (v-1) * imax * kmax &
2456 do irank = 1, recv_nmax_r2p
2470 ikv = (v-1) * imax * kmax &
2486 do irank = 1, singular_nmax
2502 var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2529 real(
dp),
intent(inout) :: var (:,:,:,:)
2531 integer :: shp(4), kmax, vmax
2532 integer :: totalsize, rank, tag
2533 integer :: irank, ipos,
imax
2534 integer :: ij_from, l_from, ij_to, l_to
2536 integer :: k, v, ikv
2540 if ( comm_apply_barrier )
then
2552 if ( kmax * vmax >
adm_kall * comm_varmax )
then
2553 write(*,*)
'xxx [COMM_data_transfer_nopl] kmax * vmax exceeds ADM_kall * COMM_varmax, stop!'
2554 write(*,*)
'xxx kmax * vmax = ', kmax * vmax
2555 write(*,*)
'xxx ADM_kall * COMM_varmax = ',
adm_kall * comm_varmax
2570 do irank = 1, recv_nmax_r2r
2571 req_count = req_count + 1
2578 mpi_double_precision, &
2587 do irank = 1, send_nmax_r2r
2598 ikv = (v-1) *
imax * kmax &
2608 req_count = req_count + 1
2609 totalsize =
imax * kmax * vmax
2615 mpi_double_precision, &
2626 do irank = 1, copy_nmax_r2r
2637 var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2645 if ( req_count > 0 )
then
2646 call mpi_waitall( req_count, &
2648 mpi_statuses_ignore, &
2655 do irank = 1, recv_nmax_r2r
2664 ikv = (v-1) *
imax * kmax &
2676 do irank = 1, singular_nmax
2687 var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2715 integer,
intent(in) :: kmax
2716 integer,
intent(in) :: vmax
2717 real(SP),
intent(inout) :: var (ADM_gall, kmax,ADM_lall, vmax)
2718 real(SP),
intent(inout) :: var_pl(ADM_gall_pl,kmax,ADM_lall_pl,vmax)
2720 real(SP) :: sendbuf_h2p_SP(kmax*vmax,PRC_RGN_total_pl)
2721 real(SP) :: recvbuf_h2p_SP(kmax*vmax,PRC_RGN_total_pl)
2723 integer :: totalsize, rank, tag
2724 integer :: irank, ipos
2725 integer :: ij_from, l_from, ij_to, l_to
2726 integer :: r_from, r_to
2732 if ( comm_apply_barrier )
then
2745 do irank = 1, send_nmax_p2r
2751 req_count = req_count + 1
2752 totalsize = kmax * vmax
2756 call mpi_irecv( recvbuf_h2p_sp(1,
i_npl), &
2767 req_count = req_count + 1
2768 totalsize = kmax * vmax
2772 call mpi_irecv( recvbuf_h2p_sp(1,
i_spl), &
2785 do irank = 1, recv_nmax_p2r
2794 kk = (v-1) * kmax + k
2795 sendbuf_h2p_sp(kk,
i_npl) = var(ij_from,k,l_from,v)
2799 req_count = req_count + 1
2800 totalsize = kmax * vmax
2804 call mpi_isend( sendbuf_h2p_sp(1,
i_npl), &
2817 kk = (v-1) * kmax + k
2818 sendbuf_h2p_sp(kk,
i_spl) = var(ij_from,k,l_from,v)
2822 req_count = req_count + 1
2823 totalsize = kmax * vmax
2827 call mpi_isend( sendbuf_h2p_sp(1,
i_spl), &
2840 do irank = 1, copy_nmax_p2r
2853 var_pl(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2861 if ( req_count > 0 )
then
2862 call mpi_waitall( req_count, &
2864 mpi_statuses_ignore, &
2869 do irank = 1, send_nmax_p2r
2880 kk = (v-1) * kmax + k
2881 var_pl(ij_to,k,l_to,v) = recvbuf_h2p_sp(kk,
i_npl)
2889 kk = (v-1) * kmax + k
2890 var_pl(ij_to,k,l_to,v) = recvbuf_h2p_sp(kk,
i_spl)
2923 integer,
intent(in) :: kmax
2924 integer,
intent(in) :: vmax
2925 real(DP),
intent(inout) :: var (ADM_gall, kmax,ADM_lall, vmax)
2926 real(DP),
intent(inout) :: var_pl(ADM_gall_pl,kmax,ADM_lall_pl,vmax)
2928 real(DP) :: sendbuf_h2p_DP(kmax*vmax,PRC_RGN_total_pl)
2929 real(DP) :: recvbuf_h2p_DP(kmax*vmax,PRC_RGN_total_pl)
2931 integer :: totalsize, rank, tag
2932 integer :: irank, ipos
2933 integer :: ij_from, l_from, ij_to, l_to
2934 integer :: r_from, r_to
2940 if ( comm_apply_barrier )
then
2953 do irank = 1, send_nmax_p2r
2959 req_count = req_count + 1
2960 totalsize = kmax * vmax
2964 call mpi_irecv( recvbuf_h2p_dp(1,
i_npl), &
2966 mpi_double_precision, &
2975 req_count = req_count + 1
2976 totalsize = kmax * vmax
2980 call mpi_irecv( recvbuf_h2p_dp(1,
i_spl), &
2982 mpi_double_precision, &
2993 do irank = 1, recv_nmax_p2r
3002 kk = (v-1) * kmax + k
3003 sendbuf_h2p_dp(kk,
i_npl) = var(ij_from,k,l_from,v)
3007 req_count = req_count + 1
3008 totalsize = kmax * vmax
3012 call mpi_isend( sendbuf_h2p_dp(1,
i_npl), &
3014 mpi_double_precision, &
3025 kk = (v-1) * kmax + k
3026 sendbuf_h2p_dp(kk,
i_spl) = var(ij_from,k,l_from,v)
3030 req_count = req_count + 1
3031 totalsize = kmax * vmax
3035 call mpi_isend( sendbuf_h2p_dp(1,
i_spl), &
3037 mpi_double_precision, &
3048 do irank = 1, copy_nmax_p2r
3061 var_pl(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
3069 if ( req_count > 0 )
then
3070 call mpi_waitall( req_count, &
3072 mpi_statuses_ignore, &
3077 do irank = 1, send_nmax_p2r
3088 kk = (v-1) * kmax + k
3089 var_pl(ij_to,k,l_to,v) = recvbuf_h2p_dp(kk,
i_npl)
3097 kk = (v-1) * kmax + k
3098 var_pl(ij_to,k,l_to,v) = recvbuf_h2p_dp(kk,
i_spl)
3117 function suf(i,j)
result(suffix)
3129 subroutine comm_debugtest
3146 integer :: i, j, k, l, ij, rgnid, prc
3152 var(:,:,:,:) = -999.0_rp
3153 var_pl(:,:,:,:) = -999.0_rp
3164 var(ij,k,l,1) = real(prc, kind=
rp)
3165 var(ij,k,l,2) = real(rgnid,kind=
rp)
3166 var(ij,k,l,3) = real(i, kind=
rp)
3167 var(ij,k,l,4) = real(j, kind=
rp)
3174 var(1,k,l,:) = -1.0_rp
3185 var_pl(ij,k,l,1) = real(-prc, kind=
rp)
3186 var_pl(ij,k,l,2) = real(-rgnid,kind=
rp)
3187 var_pl(ij,k,l,3) = real(-ij, kind=
rp)
3188 var_pl(ij,k,l,4) = real(-ij, kind=
rp)
3208 '(',int(var(ij,k,l,1)),
',',int(var(ij,k,l,2)),
')'
3227 '(',int(var_pl(ij,k,l,1)),
',',int(var_pl(ij,k,l,2)),
')'
3248 '(',int(var(ij,k,l,3)),
',',int(var(ij,k,l,4)),
')'
3267 '(',int(var_pl(ij,k,l,3)),
',',int(var_pl(ij,k,l,4)),
')'
3278 call comm_data_transfer( var(:,:,:,:), var_pl(:,:,:,:) )
3295 '(',int(var(ij,k,l,1)),
',',int(var(ij,k,l,2)),
')'
3314 '(',int(var_pl(ij,k,l,1)),
',',int(var_pl(ij,k,l,2)),
')'
3335 '(',int(var(ij,k,l,3)),
',',int(var(ij,k,l,4)),
')'
3354 '(',int(var_pl(ij,k,l,3)),
',',int(var_pl(ij,k,l,4)),
')'
3372 var(ij,k,l,1) = real(prc, kind=
rp)
3373 var(ij,k,l,2) = real(rgnid,kind=
rp)
3374 var(ij,k,l,3) = real(i, kind=
rp)
3375 var(ij,k,l,4) = real(j, kind=
rp)
3385 var(ij,k,l,1) = real(prc, kind=
rp)
3386 var(ij,k,l,2) = real(rgnid,kind=
rp)
3387 var(ij,k,l,3) = real(i, kind=
rp)
3388 var(ij,k,l,4) = real(j, kind=
rp)
3397 call comm_var( var(:,:,:,:), var_pl(:,:,:,:),
adm_kall, 4 )
3414 '(',int(var(ij,k,l,1)),
',',int(var(ij,k,l,2)),
')'
3433 '(',int(var_pl(ij,k,l,1)),
',',int(var_pl(ij,k,l,2)),
')'
3454 '(',int(var(ij,k,l,3)),
',',int(var(ij,k,l,4)),
')'
3473 '(',int(var_pl(ij,k,l,3)),
',',int(var_pl(ij,k,l,4)),
')'
3482 end subroutine comm_debugtest
3491 real(SP),
intent(in) :: localsum
3492 real(SP),
intent(out) :: globalsum
3494 real(SP) :: sendbuf(1)
3495 real(SP) :: recvbuf(PRC_nprocs)
3501 sendbuf(1) = localsum
3503 call mpi_allgather( sendbuf, &
3512 globalsum = sum( recvbuf(:) )
3514 globalsum = localsum
3527 real(DP),
intent(in) :: localsum
3528 real(DP),
intent(out) :: globalsum
3530 real(DP) :: sendbuf(1)
3531 real(DP) :: recvbuf(PRC_nprocs)
3537 sendbuf(1) = localsum
3539 call mpi_allgather( sendbuf, &
3541 mpi_double_precision, &
3544 mpi_double_precision, &
3548 globalsum = sum( recvbuf(:) )
3550 globalsum = localsum
3563 integer,
intent(in) :: kall
3564 real(SP),
intent(in) :: localsum (kall)
3565 real(SP),
intent(out) :: globalsum(kall)
3567 real(SP) :: sendbuf(kall)
3568 integer :: displs (PRC_nprocs)
3569 integer :: counts (PRC_nprocs)
3570 real(SP) :: recvbuf(kall,PRC_nprocs)
3576 do p = 1, prc_nprocs
3577 displs(p) = (p-1) * kall
3582 sendbuf(:) = localsum(:)
3584 call mpi_allgatherv( sendbuf, &
3595 globalsum(k) = sum( recvbuf(k,:) )
3599 globalsum(k) = localsum(k)
3613 integer,
intent(in) :: kall
3614 real(DP),
intent(in) :: localsum (kall)
3615 real(DP),
intent(out) :: globalsum(kall)
3617 real(DP) :: sendbuf(kall)
3618 integer :: displs (PRC_nprocs)
3619 integer :: counts (PRC_nprocs)
3620 real(DP) :: recvbuf(kall,PRC_nprocs)
3626 do p = 1, prc_nprocs
3627 displs(p) = (p-1) * kall
3632 sendbuf(:) = localsum(:)
3634 call mpi_allgatherv( sendbuf, &
3636 mpi_double_precision, &
3640 mpi_double_precision, &
3645 globalsum(k) = sum( recvbuf(k,:) )
3649 globalsum(k) = localsum(k)
3663 real(SP),
intent(in) :: localavg
3664 real(SP),
intent(out) :: globalavg
3666 real(SP) :: sendbuf(1)
3667 real(SP) :: recvbuf(PRC_nprocs)
3673 sendbuf(1) = localavg
3675 call mpi_allgather( sendbuf, &
3684 globalavg = sum( recvbuf(:) ) / real(prc_nprocs,kind=sp)
3686 globalavg = localavg
3699 real(DP),
intent(in) :: localavg
3700 real(DP),
intent(out) :: globalavg
3702 real(DP) :: sendbuf(1)
3703 real(DP) :: recvbuf(PRC_nprocs)
3709 sendbuf(1) = localavg
3711 call mpi_allgather( sendbuf, &
3713 mpi_double_precision, &
3716 mpi_double_precision, &
3720 globalavg = sum( recvbuf(:) ) / real(prc_nprocs,kind=dp)
3722 globalavg = localavg
3735 real(SP),
intent(in) :: localmax
3736 real(SP),
intent(out) :: globalmax
3738 real(SP) :: sendbuf(1)
3739 real(SP) :: recvbuf(PRC_nprocs)
3744 sendbuf(1) = localmax
3746 call mpi_allgather( sendbuf, &
3755 globalmax = maxval( recvbuf(:) )
3767 real(DP),
intent(in) :: localmax
3768 real(DP),
intent(out) :: globalmax
3770 real(DP) :: sendbuf(1)
3771 real(DP) :: recvbuf(PRC_nprocs)
3776 sendbuf(1) = localmax
3778 call mpi_allgather( sendbuf, &
3780 mpi_double_precision, &
3783 mpi_double_precision, &
3787 globalmax = maxval( recvbuf(:) )
3799 real(SP),
intent(in) :: localmin
3800 real(SP),
intent(out) :: globalmin
3802 real(SP) :: sendbuf(1)
3803 real(SP) :: recvbuf(PRC_nprocs)
3808 sendbuf(1) = localmin
3810 call mpi_allgather( sendbuf, &
3819 globalmin = minval( recvbuf(:) )
3831 real(DP),
intent(in) :: localmin
3832 real(DP),
intent(out) :: globalmin
3834 real(DP) :: sendbuf(1)
3835 real(DP) :: recvbuf(PRC_nprocs)
3840 sendbuf(1) = localmin
3842 call mpi_allgather( sendbuf, &
3844 mpi_double_precision, &
3847 mpi_double_precision, &
3851 globalmin = minval( recvbuf(:) )