22 #include "inc_openmp.h" 65 public :: comm_horizontal_max
66 public :: comm_horizontal_min
72 module procedure comm_vars_2d
73 module procedure comm_vars_3d
74 end interface comm_vars
77 module procedure comm_vars8_2d
78 module procedure comm_vars8_3d
79 end interface comm_vars8
82 module procedure comm_wait_2d
83 module procedure comm_wait_3d
84 end interface comm_wait
86 interface comm_horizontal_max
88 module procedure comm_horizontal_max_3d
89 end interface comm_horizontal_max
91 interface comm_horizontal_min
93 module procedure comm_horizontal_min_3d
94 end interface comm_horizontal_min
99 end interface comm_gather
111 end interface comm_bcast
129 integer,
private :: comm_nreq_max
130 integer,
private :: comm_vsize_max
131 integer,
private :: comm_vsize_max_pc
133 logical,
private :: comm_isallperiodic
135 integer,
private :: comm_size2d_ns4
136 integer,
private :: comm_size2d_ns8
137 integer,
private :: comm_size2d_we
138 integer,
private :: comm_size2d_4c
140 integer,
private :: comm_vars_id = 0
142 logical,
private :: comm_use_mpi_pc = .true.
144 real(RP),
private,
allocatable :: recvpack_w2p(:,:)
145 real(RP),
private,
allocatable :: recvpack_e2p(:,:)
146 real(RP),
private,
allocatable :: sendpack_p2w(:,:)
147 real(RP),
private,
allocatable :: sendpack_p2e(:,:)
149 logical,
private,
allocatable :: use_packbuf(:)
152 integer,
private,
allocatable :: req_cnt (:)
153 integer,
private,
allocatable :: req_list(:,:)
154 integer,
private,
allocatable :: preq_cnt (:)
155 integer,
private,
allocatable :: preq_list(:,:)
156 integer,
private,
allocatable :: pseqid(:)
171 namelist / param_comm / &
176 integer :: nreq_ns, nreq_we, nreq_4c
181 if( io_l )
write(io_fid_log,*)
182 if( io_l )
write(io_fid_log,*)
'++++++ Module[COMM] / Categ[ATMOS-RM COMM] / Origin[SCALElib]' 184 comm_vsize_max = max( 10 +
qa*2, 25 )
185 comm_vsize_max_pc = 50 +
qa*2
191 if( io_l )
write(io_fid_log,*)
'*** Not found namelist. Default used.' 192 elseif( ierr > 0 )
then 193 write(*,*)
'xxx Not appropriate names in namelist PARAM_COMM. Check!' 196 if( io_nml )
write(io_fid_nml,nml=param_comm)
202 if ( comm_use_mpi_pc )
then 203 comm_nreq_max = 2 * nreq_ns + 2 * nreq_we + 4 * nreq_4c + 1
205 comm_nreq_max = 2 * nreq_ns + 2 * nreq_we + 4 * nreq_4c
209 comm_size2d_ns8 =
imax 211 comm_size2d_4c =
ihalo 213 allocate( recvpack_w2p(comm_size2d_we*
ka,comm_vsize_max) )
214 allocate( recvpack_e2p(comm_size2d_we*
ka,comm_vsize_max) )
215 allocate( sendpack_p2w(comm_size2d_we*
ka,comm_vsize_max) )
216 allocate( sendpack_p2e(comm_size2d_we*
ka,comm_vsize_max) )
218 allocate( use_packbuf(comm_vsize_max) )
219 use_packbuf(:) = .false.
222 allocate( req_cnt( comm_vsize_max) )
223 allocate( req_list(comm_nreq_max,comm_vsize_max) )
225 req_list(:,:) = mpi_request_null
227 if ( comm_use_mpi_pc )
then 228 allocate( preq_cnt( comm_vsize_max_pc) )
229 allocate( preq_list(comm_nreq_max+1,comm_vsize_max_pc) )
231 preq_list(:,:) = mpi_request_null
233 allocate( pseqid(comm_vsize_max_pc) )
237 comm_isallperiodic = .true.
239 comm_isallperiodic = .false.
242 if (
rp == kind(0.d0) )
then 244 elseif(
rp == kind(0.0) )
then 247 write(*,*)
'xxx precision is not supportd' 253 if( io_l )
write(io_fid_log,*)
254 if( io_l )
write(io_fid_log,*)
'*** Maximum number of vars for one communication: ', &
256 if( io_l )
write(io_fid_log,*)
'*** Data size of var (3D,including halo) [byte] : ', &
258 if( io_l )
write(io_fid_log,*)
'*** Data size of halo [byte] : ', &
260 if( io_l )
write(io_fid_log,*)
'*** Ratio of halo against the whole 3D grid : ', &
261 real(2*IA*JHALO+2*JMAX*IHALO) /
real(
ia*
ja)
262 if( io_l )
write(io_fid_log,*)
'*** All side is periodic? : ', comm_isallperiodic
275 character(len=*),
intent(in) :: varname
276 real(RP),
intent(inout) :: var(:,:,:)
277 integer,
intent(inout) :: vid
280 if ( vid > comm_vsize_max )
then 281 write(*,*)
'xxx vid exceeds max', vid, comm_vsize_max
285 if ( comm_use_mpi_pc )
then 287 comm_vars_id = comm_vars_id + 1
288 if ( comm_vars_id > comm_vsize_max_pc )
then 289 write(*,*)
'xxx number of variable for MPI PC exceeds max', comm_vars_id, comm_vsize_max_pc
297 vid = comm_vars_id + comm_vsize_max
299 if(
io_l )
write(
io_fid_log,
'(1x,A,I3.3,2A)')
'*** [Pers.COMM] Initialize variable : ID = ', vid, &
300 ', name = ', trim(varname)
315 character(len=*),
intent(in) :: varname
316 real(RP),
intent(inout) :: var(:,:,:)
317 integer,
intent(inout) :: vid
320 if ( vid > comm_vsize_max )
then 321 write(*,*)
'xxx vid exceeds max', vid, comm_vsize_max
325 if ( comm_use_mpi_pc )
then 327 comm_vars_id = comm_vars_id + 1
328 if ( comm_vars_id > comm_vsize_max_pc )
then 329 write(*,*)
'xxx number of variable for MPI PC exceeds max', comm_vars_id, comm_vsize_max_pc
334 call vars8_init_mpi_pc(var, comm_vars_id, vid)
337 vid = comm_vars_id + comm_vsize_max
339 if(
io_l )
write(
io_fid_log,
'(1x,A,I3.3,2A)')
'*** [Pers.COMM] Initialize variable : ID = ', vid, &
340 ', name = ', trim(varname)
348 subroutine comm_vars_3d(var, vid)
351 real(RP),
intent(inout) :: var(:,:,:)
352 integer,
intent(in) :: vid
355 if ( vid > comm_vsize_max )
then 361 call vars_3d_mpi(var, vid)
366 end subroutine comm_vars_3d
369 subroutine comm_vars8_3d(var, vid)
372 real(RP),
intent(inout) :: var(:,:,:)
373 integer,
intent(in) :: vid
376 if ( vid > comm_vsize_max )
then 387 end subroutine comm_vars8_3d
390 subroutine comm_wait_3d(var, vid, FILL_BND)
393 real(RP),
intent(inout) :: var(:,:,:)
394 integer,
intent(in) :: vid
395 logical,
intent(in),
optional :: fill_bnd
401 if (
present(fill_bnd) ) fill_bnd_ = fill_bnd
403 if ( vid > comm_vsize_max )
then 405 call wait_3d_mpi_pc(var, vid-comm_vsize_max)
414 if ( .NOT. comm_isallperiodic )
then 415 if ( fill_bnd_ )
then 416 call copy_boundary_3d(var)
421 end subroutine comm_wait_3d
424 subroutine comm_vars_2d(var, vid)
426 real(RP),
intent(inout) :: var(:,:)
427 integer,
intent(in) :: vid
435 end subroutine comm_vars_2d
438 subroutine comm_vars8_2d(var, vid)
441 real(RP),
intent(inout) :: var(:,:)
442 integer,
intent(in) :: vid
450 end subroutine comm_vars8_2d
453 subroutine comm_wait_2d(var, vid, FILL_BND)
456 real(RP),
intent(inout) :: var(:,:)
457 integer,
intent(in) :: vid
458 logical,
intent(in),
optional :: fill_bnd
464 if (
present(fill_bnd) ) fill_bnd_ = fill_bnd
467 call wait_2d_mpi(var, vid)
470 if( .NOT. comm_isallperiodic )
then 471 if ( fill_bnd_ )
then 472 call copy_boundary_2d(var)
477 end subroutine comm_wait_2d
486 real(RP),
intent(out) :: varmean(
ka)
487 real(RP),
intent(in) :: var (
ka,
ia,
ja)
489 real(RP) :: statval (
ka)
490 real(RP) :: statcnt (
ka)
491 real(RP) :: allstatval(
ka)
492 real(RP) :: allstatcnt(
ka)
505 statval(k) = statval(k) + var(k,i,j)
506 statcnt(k) = statcnt(k) + 1.0_rp
515 call mpi_allreduce( statval(1), &
523 call mpi_allreduce( statcnt(1), &
534 zerosw = 0.5_rp - sign(0.5_rp, allstatcnt(k) - 1.e-12_rp )
535 varmean(k) = allstatval(k) / ( allstatcnt(k) + zerosw ) * ( 1.0_rp - zerosw )
547 real(RP),
intent(out) :: varmax
548 real(RP),
intent(in) :: var(IA,JA)
551 real(RP) :: allstatval
561 call mpi_allreduce( statval, &
578 subroutine comm_horizontal_max_3d( varmax, var )
583 real(RP),
intent(out) :: varmax(KA)
584 real(RP),
intent(in) :: var (KA,IA,JA)
586 real(RP) :: statval (KA)
587 real(RP) :: allstatval(KA)
593 statval(:) = -1.e19_rp
595 statval(k) = maxval(var(k,
is:
ie,
js:
je))
601 call mpi_allreduce( statval(1), &
612 varmax(k) = allstatval(k)
618 end subroutine comm_horizontal_max_3d
625 real(RP),
intent(out) :: varmin
626 real(RP),
intent(in) :: var(IA,JA)
629 real(RP) :: allstatval
639 call mpi_allreduce( statval, &
656 subroutine comm_horizontal_min_3d( varmin, var )
661 real(RP),
intent(out) :: varmin(KA)
662 real(RP),
intent(in) :: var (KA,IA,JA)
664 real(RP) :: statval (KA)
665 real(RP) :: allstatval(KA)
671 statval(:) = -1.e19_rp
673 statval(k) = minval(var(k,
is:
ie,
js:
je))
679 call mpi_allreduce( statval(1), &
690 varmin(k) = allstatval(k)
696 end subroutine comm_horizontal_min_3d
705 real(RP),
intent(out) :: recv(:,:)
706 real(RP),
intent(in) :: send(:,:)
707 integer,
intent(in) :: gIA
708 integer,
intent(in) :: gJA
710 integer :: sendcounts, recvcounts
714 sendcounts = gia * gja
715 recvcounts = gia * gja
717 call mpi_gather( send(:,:), &
737 real(RP),
intent(out) :: recv(:,:,:)
738 real(RP),
intent(in) :: send(:,:,:)
739 integer,
intent(in) :: gIA
740 integer,
intent(in) :: gJA
741 integer,
intent(in) :: gKA
743 integer :: sendcounts, recvcounts
747 sendcounts = gia * gja * gka
748 recvcounts = gia * gja * gka
750 call mpi_gather( send(:,:,:), &
770 real(RP),
intent(inout) :: var
780 call mpi_bcast( var, &
799 real(RP),
intent(inout) :: var(:)
800 integer,
intent(in) :: gIA
810 call mpi_bcast( var(:), &
829 real(RP),
intent(inout) :: var(:,:)
830 integer,
intent(in) :: gIA
831 integer,
intent(in) :: gJA
841 call mpi_bcast( var(:,:), &
860 real(RP),
intent(inout) :: var(:,:,:)
861 integer,
intent(in) :: gIA
862 integer,
intent(in) :: gJA
863 integer,
intent(in) :: gKA
871 counts = gia * gja * gka
873 call mpi_bcast( var(:,:,:), &
892 real(RP),
intent(inout) :: var(:,:,:,:)
893 integer,
intent(in) :: gIA
894 integer,
intent(in) :: gJA
895 integer,
intent(in) :: gKA
896 integer,
intent(in) :: gTime
904 counts = gia * gja * gka * gtime
905 if ( gia>0 .AND. gja>0 .AND. gka>0 .AND. gtime>0 .AND. &
907 write(*,*)
'xxx counts overflow' 911 call mpi_bcast( var(:,:,:,:), &
930 integer,
intent(inout) :: var
940 call mpi_bcast( var, &
959 integer,
intent(inout) :: var(:)
960 integer,
intent(in) :: gIA
970 call mpi_bcast( var(:), &
989 integer,
intent(inout) :: var(:,:)
990 integer,
intent(in) :: gIA
991 integer,
intent(in) :: gJA
1001 call mpi_bcast( var(:,:), &
1020 logical,
intent(inout) :: var
1030 call mpi_bcast( var, &
1048 real(RP),
intent(inout) :: var(:,:,:)
1049 integer,
intent(in) :: vid
1050 integer,
intent(in) :: seqid
1052 integer :: ireq, tag, ierr
1066 mpi_proc_null, tag+comm_nreq_max+1,
comm_world, &
1067 preq_list(comm_nreq_max+1,vid), ierr )
1079 call mpi_recv_init( recvpack_e2p(:,seqid), comm_size2d_we*kd,
comm_datatype, &
1083 call mpi_recv_init( recvpack_w2p(:,seqid), comm_size2d_we*kd,
comm_datatype, &
1089 call mpi_send_init( sendpack_p2w(:,seqid), comm_size2d_we*kd,
comm_datatype, &
1093 call mpi_send_init( sendpack_p2e(:,seqid), comm_size2d_we*kd,
comm_datatype, &
1105 preq_cnt(vid) = ireq - 1
1110 call mpi_testall( preq_cnt(vid), preq_list(1:preq_cnt(vid),vid), &
1111 flag, mpi_statuses_ignore, ierr )
1117 subroutine vars8_init_mpi_pc(var, vid, seqid)
1120 real(RP),
intent(inout) :: var(:,:,:)
1121 integer,
intent(in) :: vid
1122 integer,
intent(in) :: seqid
1124 integer :: ireq, tag, tagc
1139 mpi_proc_null, tag+comm_nreq_max+1,
comm_world, &
1140 preq_list(comm_nreq_max+1,vid), ierr )
1143 if ( comm_isallperiodic )
then 1149 call mpi_recv_init( var(1,
ie+1,j), comm_size2d_4c*kd,
comm_datatype, &
1165 call mpi_recv_init( var(1,
ie+1,j), comm_size2d_4c*kd,
comm_datatype, &
1196 call mpi_recv_init( recvpack_e2p(:,seqid), comm_size2d_we*kd,
comm_datatype, &
1201 call mpi_recv_init( recvpack_w2p(:,seqid), comm_size2d_we*kd,
comm_datatype, &
1208 call mpi_send_init( sendpack_p2w(:,seqid), comm_size2d_we*kd,
comm_datatype, &
1213 call mpi_send_init( sendpack_p2e(:,seqid), comm_size2d_we*kd,
comm_datatype, &
1272 call mpi_recv_init( var(1,
ie+1,j), comm_size2d_4c*kd,
comm_datatype, &
1280 call mpi_recv_init( var(1,
ie+1,j), comm_size2d_4c*kd,
comm_datatype, &
1288 call mpi_recv_init( var(1,
ie+1,j), comm_size2d_4c*kd,
comm_datatype, &
1324 call mpi_recv_init( var(1,
ie+1,j), comm_size2d_4c*kd,
comm_datatype, &
1332 call mpi_recv_init( var(1,
ie+1,j), comm_size2d_4c*kd,
comm_datatype, &
1340 call mpi_recv_init( var(1,
ie+1,j), comm_size2d_4c*kd,
comm_datatype, &
1395 call mpi_recv_init( recvpack_e2p(:,seqid), comm_size2d_we*kd,
comm_datatype, &
1402 call mpi_recv_init( recvpack_w2p(:,seqid), comm_size2d_we*kd,
comm_datatype, &
1411 call mpi_send_init( sendpack_p2w(:,seqid), comm_size2d_we*kd,
comm_datatype, &
1418 call mpi_send_init( sendpack_p2e(:,seqid), comm_size2d_we*kd,
comm_datatype, &
1480 call mpi_send_init( var(1,
ie+1,j), comm_size2d_4c*kd,
comm_datatype, &
1532 call mpi_send_init( var(1,
ie+1,j), comm_size2d_4c*kd,
comm_datatype, &
1549 preq_cnt(vid) = ireq - 1
1554 call mpi_testall( preq_cnt(vid), preq_list(1:preq_cnt(vid),vid), &
1555 flag, mpi_statuses_ignore, ierr )
1559 end subroutine vars8_init_mpi_pc
1561 subroutine vars_3d_mpi(var, vid)
1566 real(RP),
intent(inout) :: var(:,:,:)
1567 integer,
intent(in) :: vid
1570 integer :: ireq, tag
1582 if ( use_packbuf(vid) )
then 1583 write(*,*)
'packing buffer is already used', vid
1586 use_packbuf(vid) = .true.
1589 if ( comm_isallperiodic )
then 1601 call mpi_irecv( recvpack_e2p(:,vid), comm_size2d_we*kd,
comm_datatype, &
1605 call mpi_irecv( recvpack_w2p(:,vid), comm_size2d_we*kd,
comm_datatype, &
1609 call pack_3d(var, vid)
1613 call mpi_isend( sendpack_p2w(:,vid), comm_size2d_we*kd,
comm_datatype, &
1617 call mpi_isend( sendpack_p2e(:,vid), comm_size2d_we*kd,
comm_datatype, &
1646 call mpi_irecv( recvpack_e2p(:,vid), comm_size2d_we*kd,
comm_datatype, &
1652 call mpi_irecv( recvpack_w2p(:,vid), comm_size2d_we*kd,
comm_datatype, &
1657 call pack_3d(var, vid)
1662 call mpi_isend( sendpack_p2w(:,vid), comm_size2d_we*kd,
comm_datatype, &
1668 call mpi_isend( sendpack_p2e(:,vid), comm_size2d_we*kd,
comm_datatype, &
1687 req_cnt(vid) = ireq - 1
1690 end subroutine vars_3d_mpi
1697 real(RP),
intent(inout) :: var(:,:,:)
1698 integer,
intent(in) :: vid
1700 integer :: ireq, tag, tagc
1714 if ( use_packbuf(vid) )
then 1715 write(*,*)
'packing buffer is already used', vid
1718 use_packbuf(vid) = .true.
1721 if ( comm_isallperiodic )
then 1774 call mpi_irecv( recvpack_e2p(:,vid), comm_size2d_we*kd,
comm_datatype, &
1779 call mpi_irecv( recvpack_w2p(:,vid), comm_size2d_we*kd,
comm_datatype, &
1783 call pack_3d(var, vid)
1789 call mpi_isend( sendpack_p2w(:,vid), comm_size2d_we*kd,
comm_datatype, &
1794 call mpi_isend( sendpack_p2e(:,vid), comm_size2d_we*kd,
comm_datatype, &
1976 call mpi_irecv( recvpack_e2p(:,vid), comm_size2d_we*kd,
comm_datatype, &
1983 call mpi_irecv( recvpack_w2p(:,vid), comm_size2d_we*kd,
comm_datatype, &
1988 call pack_3d(var, vid)
1994 call mpi_isend( sendpack_p2w(:,vid), comm_size2d_we*kd,
comm_datatype, &
2001 call mpi_isend( sendpack_p2e(:,vid), comm_size2d_we*kd,
comm_datatype, &
2037 call mpi_isend( var(1,1,j), comm_size2d_4c*kd,
comm_datatype, &
2089 call mpi_isend( var(1,1,j), comm_size2d_4c*kd,
comm_datatype, &
2132 req_cnt(vid) = ireq - 1
2142 real(RP),
intent(inout) :: var(:,:)
2143 integer,
intent(in) :: vid
2145 integer :: ireq, tag
2153 if ( use_packbuf(vid) )
then 2154 write(*,*)
'packing buffer is already used', vid
2157 use_packbuf(vid) = .true.
2160 if ( comm_isallperiodic )
then 2164 call mpi_irecv( var(:,
js-
jhalo:
js-1), comm_size2d_ns4, &
2170 call mpi_irecv( var(:,
je+1:
je+
jhalo), comm_size2d_ns4, &
2176 call mpi_irecv( recvpack_e2p(:,vid), comm_size2d_we, &
2182 call mpi_irecv( recvpack_w2p(:,vid), comm_size2d_we, &
2187 call pack_2d(var, vid)
2190 call mpi_isend( sendpack_p2w(:,vid), comm_size2d_we, &
2196 call mpi_isend( sendpack_p2e(:,vid), comm_size2d_we, &
2202 call mpi_isend( var(:,
je-
jhalo+1:
je), comm_size2d_ns4, &
2208 call mpi_isend( var(:,
js:
js+
jhalo-1), comm_size2d_ns4, &
2218 call mpi_irecv( var(:,
js-
jhalo:
js-1), comm_size2d_ns4, &
2226 call mpi_irecv( var(:,
je+1:
je+
jhalo), comm_size2d_ns4, &
2234 call mpi_irecv( recvpack_e2p(:,vid), comm_size2d_we, &
2242 call mpi_irecv( recvpack_w2p(:,vid), comm_size2d_we, &
2248 call pack_2d(var, vid)
2252 call mpi_isend( sendpack_p2w(:,vid), comm_size2d_we, &
2260 call mpi_isend( sendpack_p2e(:,vid), comm_size2d_we, &
2268 call mpi_isend( var(:,
je-
jhalo+1:
je), comm_size2d_ns4, &
2276 call mpi_isend( var(:,
js:
js+
jhalo-1), comm_size2d_ns4, &
2284 req_cnt(vid) = ireq - 1
2294 real(RP),
intent(inout) :: var(:,:)
2295 integer,
intent(in) :: vid
2297 integer :: ireq, tag, tagc
2307 if ( use_packbuf(vid) )
then 2308 write(*,*)
'packing buffer is already used', vid
2311 use_packbuf(vid) = .true.
2314 if ( comm_isallperiodic )
then 2320 call mpi_irecv( var(
ie+1,j), comm_size2d_4c, &
2329 call mpi_irecv( var(
is-
ihalo,j), comm_size2d_4c, &
2338 call mpi_irecv( var(
ie+1,j), comm_size2d_4c, &
2347 call mpi_irecv( var(
is-
ihalo,j), comm_size2d_4c, &
2356 call mpi_irecv( var(
is,j), comm_size2d_ns8, &
2365 call mpi_irecv( var(
is,j), comm_size2d_ns8, &
2372 call mpi_irecv( recvpack_e2p(:,vid), comm_size2d_we, &
2377 call mpi_irecv( recvpack_w2p(:,vid), comm_size2d_we, &
2382 call pack_2d(var, vid)
2385 call mpi_isend( sendpack_p2w(:,vid), comm_size2d_we, &
2391 call mpi_isend( sendpack_p2e(:,vid), comm_size2d_we, &
2399 call mpi_isend( var(
is,j), comm_size2d_ns8, &
2409 call mpi_isend( var(
is,j), comm_size2d_ns8, &
2419 call mpi_isend( var(
is,j), comm_size2d_4c, &
2429 call mpi_isend( var(
ie-
ihalo+1,j), comm_size2d_4c, &
2439 call mpi_isend( var(
is,j), comm_size2d_4c, &
2449 call mpi_isend( var(
ie-
ihalo+1,j), comm_size2d_4c, &
2462 call mpi_irecv( var(
ie+1,j), comm_size2d_4c, &
2471 call mpi_irecv( var(
ie+1,j), comm_size2d_4c, &
2480 call mpi_irecv( var(
ie+1,j), comm_size2d_4c, &
2492 call mpi_irecv( var(
is-
ihalo,j), comm_size2d_4c, &
2501 call mpi_irecv( var(
is-
ihalo,j), comm_size2d_4c, &
2510 call mpi_irecv( var(
is-
ihalo,j), comm_size2d_4c, &
2522 call mpi_irecv( var(
ie+1,j), comm_size2d_4c, &
2531 call mpi_irecv( var(
ie+1,j), comm_size2d_4c, &
2540 call mpi_irecv( var(
ie+1,j), comm_size2d_4c, &
2552 call mpi_irecv( var(
is-
ihalo,j), comm_size2d_4c, &
2561 call mpi_irecv( var(
is-
ihalo,j), comm_size2d_4c, &
2570 call mpi_irecv( var(
is-
ihalo,j), comm_size2d_4c, &
2582 call mpi_irecv( var(
is,j), comm_size2d_ns8, &
2594 call mpi_irecv( var(
is,j), comm_size2d_ns8, &
2604 call mpi_irecv( recvpack_e2p(:,vid), comm_size2d_we, &
2612 call mpi_irecv( recvpack_w2p(:,vid), comm_size2d_we, &
2618 call pack_2d(var, vid)
2622 call mpi_isend( sendpack_p2w(:,vid), comm_size2d_we, &
2630 call mpi_isend( sendpack_p2e(:,vid), comm_size2d_we, &
2640 call mpi_isend( var(
is,j), comm_size2d_ns8, &
2652 call mpi_isend( var(
is,j), comm_size2d_ns8, &
2664 call mpi_isend( var(
is,j), comm_size2d_4c, &
2673 call mpi_isend( var(
is,j), comm_size2d_4c, &
2682 call mpi_isend( var(
is,j), comm_size2d_4c, &
2694 call mpi_isend( var(
ie-
ihalo+1,j), comm_size2d_4c, &
2703 call mpi_isend( var(
ie-
ihalo+1,j), comm_size2d_4c, &
2712 call mpi_isend( var(
ie-
ihalo+1,j), comm_size2d_4c, &
2724 call mpi_isend( var(
is,j), comm_size2d_4c, &
2733 call mpi_isend( var(
is,j), comm_size2d_4c, &
2742 call mpi_isend( var(
is,j), comm_size2d_4c, &
2754 call mpi_isend( var(
ie-
ihalo+1,j), comm_size2d_4c, &
2763 call mpi_isend( var(
ie-
ihalo+1,j), comm_size2d_4c, &
2772 call mpi_isend( var(
ie-
ihalo+1,j), comm_size2d_4c, &
2782 req_cnt(vid) = ireq - 1
2792 real(RP),
intent(inout) :: var(:,:,:)
2793 integer,
intent(in) :: vid
2798 if ( use_packbuf(pseqid(vid)) )
then 2799 write(*,*)
'packing buffer is already used', vid, pseqid(vid)
2802 use_packbuf(pseqid(vid)) = .true.
2805 call pack_3d(var, pseqid(vid))
2807 call mpi_startall(preq_cnt(vid), preq_list(1:preq_cnt(vid),vid), ierr)
2814 real(RP),
intent(inout) :: var(:,:,:)
2815 integer,
intent(in) :: vid
2821 call mpi_waitall( req_cnt(vid), &
2822 req_list(1:req_cnt(vid),vid), &
2823 mpi_statuses_ignore, &
2825 call unpack_3d(var, vid)
2828 use_packbuf(vid) = .false.
2834 subroutine wait_2d_mpi(var, vid)
2837 real(RP),
intent(inout) :: var(:,:)
2838 integer,
intent(in) :: vid
2844 call mpi_waitall( req_cnt(vid), &
2845 req_list(1:req_cnt(vid),vid), &
2846 mpi_statuses_ignore, &
2848 call unpack_2d(var, vid)
2851 use_packbuf(vid) = .false.
2855 end subroutine wait_2d_mpi
2857 subroutine wait_3d_mpi_pc(var, vid)
2860 real(RP),
intent(inout) :: var(:,:,:)
2861 integer,
intent(in) :: vid
2866 call mpi_waitall( preq_cnt(vid), &
2867 preq_list(1:preq_cnt(vid),vid), &
2868 mpi_statuses_ignore, &
2870 call unpack_3d(var, pseqid(vid))
2873 use_packbuf(pseqid(vid)) = .false.
2877 end subroutine wait_3d_mpi_pc
2879 subroutine pack_3d(var, vid)
2882 real(RP),
intent(in) :: var(:,:,:)
2883 integer,
intent(in) :: vid
2886 integer :: k, i, j, n
2892 if ( comm_isallperiodic )
then 2903 sendpack_p2w(n,vid) = var(k,i,j)
2916 sendpack_p2e(n,vid) = var(k,i,j)
2934 sendpack_p2w(n,vid) = var(k,i,j)
2950 sendpack_p2e(n,vid) = var(k,i,j)
2961 end subroutine pack_3d
2963 subroutine pack_2d(var, vid)
2966 real(RP),
intent(in) :: var(:,:)
2967 integer,
intent(in) :: vid
2973 if ( comm_isallperiodic )
then 2983 sendpack_p2w(n,vid) = var(i,j)
2994 sendpack_p2e(n,vid) = var(i,j)
3009 sendpack_p2w(n,vid) = var(i,j)
3022 sendpack_p2e(n,vid) = var(i,j)
3032 end subroutine pack_2d
3034 subroutine unpack_3d(var, vid)
3037 real(RP),
intent(inout) :: var(:,:,:)
3038 integer,
intent(in) :: vid
3041 integer :: i, j, k, n
3048 if ( comm_isallperiodic )
then 3058 var(k,i,j) = recvpack_e2p(n,vid)
3071 var(k,i,j) = recvpack_w2p(n,vid)
3088 var(k,i,j) = recvpack_e2p(n,vid)
3104 var(k,i,j) = recvpack_w2p(n,vid)
3115 end subroutine unpack_3d
3117 subroutine unpack_2d(var, vid)
3120 real(RP),
intent(inout) :: var(:,:)
3121 integer,
intent(in) :: vid
3128 if( comm_isallperiodic )
then 3136 var(i,j) = recvpack_e2p(n,vid)
3146 var(i,j) = recvpack_w2p(n,vid)
3160 var(i,j) = recvpack_e2p(n,vid)
3173 var(i,j) = recvpack_w2p(n,vid)
3183 end subroutine unpack_2d
3185 subroutine copy_boundary_3d(var)
3188 real(RP),
intent(inout) :: var(:,:,:)
3197 var(:,i,j) = var(:,i,
je)
3206 var(:,i,j) = var(:,i,
js)
3215 var(:,i,j) = var(:,
ie,j)
3224 var(:,i,j) = var(:,
is,j)
3234 var(:,i,j) = var(:,
is,
je)
3240 var(:,i,j) = var(:,i,
je)
3246 var(:,i,j) = var(:,
is,j)
3256 var(:,i,j) = var(:,
is,
js)
3262 var(:,i,j) = var(:,i,
js)
3268 var(:,i,j) = var(:,
is,j)
3278 var(:,i,j) = var(:,
ie,
je)
3284 var(:,i,j) = var(:,i,
je)
3290 var(:,i,j) = var(:,
ie,j)
3300 var(:,i,j) = var(:,
ie,
js)
3306 var(:,i,j) = var(:,i,
js)
3312 var(:,i,j) = var(:,
ie,j)
3318 end subroutine copy_boundary_3d
3320 subroutine copy_boundary_2d(var)
3323 real(RP),
intent(inout) :: var(:,:)
3332 var(i,j) = var(i,
je)
3342 var(i,j) = var(i,
js)
3351 var(i,j) = var(
ie,j)
3360 var(i,j) = var(
is,j)
3370 var(i,j) = var(
is,
je)
3377 var(i,j) = var(i,
je)
3384 var(i,j) = var(
is,j)
3394 var(i,j) = var(
is,
js)
3401 var(i,j) = var(i,
js)
3408 var(i,j) = var(
is,j)
3418 var(i,j) = var(
ie,
je)
3425 var(i,j) = var(i,
je)
3432 var(i,j) = var(
ie,j)
3442 var(i,j) = var(
ie,
js)
3449 var(i,j) = var(i,
js)
3456 var(i,j) = var(
ie,j)
3462 end subroutine copy_boundary_2d
3468 namelist / param_comm / &
3469 comm_vsize_max_pc, &
3472 integer :: i, j, ierr
3475 deallocate( recvpack_w2p )
3476 deallocate( recvpack_e2p )
3477 deallocate( sendpack_p2w )
3478 deallocate( sendpack_p2e )
3480 deallocate( use_packbuf )
3483 deallocate( req_cnt )
3484 deallocate( req_list )
3486 if ( comm_use_mpi_pc )
then 3487 do j=1, comm_vsize_max_pc
3488 do i=1, comm_nreq_max+1
3489 if (preq_list(i,j) .NE. mpi_request_null) &
3490 call mpi_request_free(preq_list(i,j), ierr)
3493 deallocate( preq_cnt )
3494 deallocate( preq_list )
3495 deallocate( pseqid )
integer, public imax
of computational cells: x, local
integer, public is
start point of inner domain: x, local
integer, public comm_datatype
datatype of variable
integer, parameter, public prc_s
[node direction] south
subroutine comm_bcast_scr(var)
Broadcast data for whole process value in scalar field.
integer, public je
end point of inner domain: y, local
logical, public prc_has_n
integer, parameter, public prc_w
[node direction] west
real(rp), public const_huge
huge number
integer, public prc_local_comm_world
local communicator
subroutine, public prc_mpistop
Abort MPI.
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
logical, public io_l
output log or not? (this process)
subroutine wait_3d_mpi(var, vid)
logical, public prc_has_e
integer, public ke
end point of inner domain: z, local
subroutine comm_bcast_int_1d(var, gIA)
Broadcast data for whole process value in 1D field (integer)
subroutine comm_bcast_2d(var, gIA, gJA)
Broadcast data for whole process value in 2D field.
integer, dimension(8), public prc_next
node ID of 8 neighbour process
logical, public prc_has_s
real(rp), public const_undef
integer, parameter, public prc_n
[node direction] north
integer, public ia
of whole cells: x, local, with HALO
subroutine, public comm_horizontal_mean(varmean, var)
calculate horizontal mean (global total with communication)
integer, public ka
of whole cells: z, local, with HALO
subroutine comm_gather_2d(recv, send, gIA, gJA)
Get data from whole process value in 2D field.
integer, public comm_world
communication world ID
integer, public jhalo
of halo cells: y
subroutine, public comm_cleanup
integer, parameter, public prc_nw
[node direction] northwest
integer, public js
start point of inner domain: y, local
logical, public comm_fill_bnd
switch whether fill boundary data
integer, parameter, public prc_e
[node direction] east
subroutine vars_3d_mpi_pc(var, vid)
subroutine comm_horizontal_min_2d(varmin, var)
Get minimum value in horizontal area.
subroutine comm_bcast_logical_scr(var)
Broadcast data for whole process value in scalar (logical)
subroutine vars_2d_mpi(var, vid)
subroutine vars_init_mpi_pc(var, vid, seqid)
integer, parameter, public prc_masterrank
master process in each communicator
integer, public ks
start point of inner domain: z, local
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
subroutine comm_bcast_4d(var, gIA, gJA, gKA, gTime)
Broadcast data for whole process value in 4D field.
integer, public ie
end point of inner domain: x, local
subroutine comm_horizontal_max_2d(varmax, var)
Get maximum value in horizontal area.
subroutine comm_bcast_int_2d(var, gIA, gJA)
Broadcast data for whole process value in 2D field (integer)
integer, parameter, public prc_sw
[node direction] southwest
subroutine comm_bcast_int_scr(var)
Broadcast data for whole process value in scalar (integer)
subroutine comm_bcast_3d(var, gIA, gJA, gKA)
Broadcast data for whole process value in 3D field.
integer, parameter, public prc_ne
[node direction] northeast
integer, public io_fid_conf
Config file ID.
subroutine, public comm_setup
Setup.
integer, public io_fid_log
Log file ID.
subroutine comm_gather_3d(recv, send, gIA, gJA, gKA)
Get data from whole process value in 3D field.
integer, parameter, public prc_se
[node direction] southeast
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
logical, public prc_has_w
integer, parameter, public rp
subroutine vars8_2d_mpi(var, vid)
subroutine vars8_3d_mpi(var, vid)
integer, public jmax
of computational cells: y, local
subroutine comm_bcast_1d(var, gIA)
Broadcast data for whole process value in 1D field.
integer, public ihalo
of halo cells: x
integer, public ja
of whole cells: y, local, with HALO
subroutine, public comm_vars_init(varname, var, vid)
Register variables.