SCALE-RM
scale_comm_icoA.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
9 !-------------------------------------------------------------------------------
11  !-----------------------------------------------------------------------------
12  !
13  !++ Used modules
14  !
15  use mpi
16  use scale_precision
17  use scale_io
18  use scale_prof
20  !-----------------------------------------------------------------------------
21  implicit none
22  private
23  !-----------------------------------------------------------------------------
24  !
25  !++ Public procedures
26  !
27  public :: comm_setup
28  public :: comm_data_transfer
29  public :: comm_data_transfer_nopl
30  public :: comm_var
31 
32  interface comm_data_transfer
33  module procedure comm_data_transfer_sp
34  module procedure comm_data_transfer_dp
35  end interface comm_data_transfer
36 
37  interface comm_var
38  module procedure comm_var_sp
39  module procedure comm_var_dp
40  end interface comm_var
41 
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
47 
48  interface comm_stat_sum
49  module procedure comm_stat_sum_sp
50  module procedure comm_stat_sum_dp
51  end interface comm_stat_sum
52 
53  interface comm_stat_sum_eachlayer
54  module procedure comm_stat_sum_eachlayer_sp
55  module procedure comm_stat_sum_eachlayer_dp
56  end interface comm_stat_sum_eachlayer
57 
58  interface comm_stat_avg
59  module procedure comm_stat_avg_sp
60  module procedure comm_stat_avg_dp
61  end interface comm_stat_avg
62 
63  interface comm_stat_max
64  module procedure comm_stat_max_sp
65  module procedure comm_stat_max_dp
66  end interface comm_stat_max
67 
68  interface comm_stat_min
69  module procedure comm_stat_min_sp
70  module procedure comm_stat_min_dp
71  end interface comm_stat_min
72 
73  !-----------------------------------------------------------------------------
74  !
75  !++ Public parameters & variables
76  !
77  logical, public :: comm_pl = .true.
78 
79  integer, public :: comm_datatype
80 
81  ! node-node communication/copy info
82  integer, public, allocatable :: copy_info_r2r(:)
83  integer, public, allocatable :: recv_info_r2r(:,:)
84  integer, public, allocatable :: send_info_r2r(:,:)
85 
86  integer, public, allocatable :: copy_info_p2r(:)
87  integer, public, allocatable :: recv_info_p2r(:,:)
88  integer, public, allocatable :: send_info_p2r(:,:)
89 
90  integer, public, allocatable :: copy_info_r2p(:)
91  integer, public, allocatable :: recv_info_r2p(:,:)
92  integer, public, allocatable :: send_info_r2p(:,:)
93 
94  integer, public, allocatable :: singular_info(:)
95 
96  ! node-node communication/copy list
97  integer, public, allocatable :: copy_list_r2r(:,:)
98  integer, public, allocatable :: recv_list_r2r(:,:,:)
99  integer, public, allocatable :: send_list_r2r(:,:,:)
100 
101  integer, public, allocatable :: copy_list_p2r(:,:)
102  integer, public, allocatable :: recv_list_p2r(:,:,:)
103  integer, public, allocatable :: send_list_p2r(:,:,:)
104 
105  integer, public, allocatable :: copy_list_r2p(:,:)
106  integer, public, allocatable :: recv_list_r2p(:,:,:)
107  integer, public, allocatable :: send_list_r2p(:,:,:)
108 
109  integer, public, allocatable :: singular_list(:,:)
110 
111  ! working buffer
112  integer, public, allocatable :: req_list(:)
113 
114  real(SP), public, allocatable :: sendbuf_r2r_sp(:,:)
115  real(SP), public, allocatable :: recvbuf_r2r_sp(:,:)
116  real(SP), public, allocatable :: sendbuf_p2r_sp(:,:)
117  real(SP), public, allocatable :: recvbuf_p2r_sp(:,:)
118  real(SP), public, allocatable :: sendbuf_r2p_sp(:,:)
119  real(SP), public, allocatable :: recvbuf_r2p_sp(:,:)
120 
121  real(DP), public, allocatable :: sendbuf_r2r_dp(:,:)
122  real(DP), public, allocatable :: recvbuf_r2r_dp(:,:)
123  real(DP), public, allocatable :: sendbuf_p2r_dp(:,:)
124  real(DP), public, allocatable :: recvbuf_p2r_dp(:,:)
125  real(DP), public, allocatable :: sendbuf_r2p_dp(:,:)
126  real(DP), public, allocatable :: recvbuf_r2p_dp(:,:)
127 
128  !-----------------------------------------------------------------------------
129  !
130  !++ Private procedures
131  !
132  private :: comm_list_generate
133  private :: comm_sortdest
134  private :: comm_sortdest_pl
135  private :: comm_sortdest_singular
136 
137  private :: comm_debugtest
138 
139  !-----------------------------------------------------------------------------
140  !
141  !++ Private parameters & variables
142  !
143  logical, private :: comm_apply_barrier = .false.
144  integer, private :: comm_varmax = 50
145  logical, private :: debug = .false.
146  logical, private :: testonly = .false.
147 
148  ! send/recv relationship
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
156 
157  integer, private, allocatable :: rellist(:,:)
158  integer, private :: rellist_nmax
159 
160  ! node-node communication/copy info
161  integer, private, parameter :: recv_nlim = 20
162  integer, private, parameter :: send_nlim = 20
163 
164  integer, private :: copy_nmax_r2r = 0
165  integer, private :: recv_nmax_r2r = 0
166  integer, private :: send_nmax_r2r = 0
167 
168  integer, private :: copy_nmax_p2r = 0
169  integer, private :: recv_nmax_p2r = 0
170  integer, private :: send_nmax_p2r = 0
171 
172  integer, private :: copy_nmax_r2p = 0
173  integer, private :: recv_nmax_r2p = 0
174  integer, private :: send_nmax_r2p = 0
175 
176  integer, private :: singular_nmax = 0
177 
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
182 
183  ! node-node communication/copy list
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
189 
190  ! working buffer
191  integer, private :: req_count
192 
193  !-----------------------------------------------------------------------------
194 contains
195  !-----------------------------------------------------------------------------
197  subroutine comm_setup
198  use scale_prc, only: &
199  prc_abort
200  use scale_prc_icoa, only: &
201  prc_rgn_r2p_pl, &
202  i_npl, &
203  i_spl
204  implicit none
205 
206  namelist / param_comm_icoa / &
207  comm_apply_barrier, &
208  comm_varmax, &
209  debug, &
210  testonly
211 
212  integer :: ierr
213  !---------------------------------------------------------------------------
214 
215  !--- read parameters
216  if( io_l ) write(io_fid_log,*)
217  if( io_l ) write(io_fid_log,*) '+++ Module[comm]/Category[common share]'
218  rewind(io_fid_conf)
219  read(io_fid_conf,nml=param_comm_icoa,iostat=ierr)
220  if ( ierr < 0 ) then
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.'
224  call prc_abort
225  endif
226  if( io_nml ) write(io_fid_log,nml=param_comm_icoa)
227 
228  if ( rp == dp ) then
229  comm_datatype = mpi_double_precision
230  elseif( rp == sp ) then
231  comm_datatype = mpi_real
232  else
233  write(*,*) 'xxx precision is not supportd'
234  call prc_abort
235  endif
236 
237  if ( prc_rgn_r2p_pl(i_npl) < 0 &
238  .AND. prc_rgn_r2p_pl(i_spl) < 0 ) then
239  comm_pl = .false.
240  endif
241 
242  if( io_l ) write(io_fid_log,*)
243  if( io_l ) write(io_fid_log,*) '====== communication information ======'
244 
245  call comm_list_generate
246 
247  call comm_sortdest
248  call comm_sortdest_pl
249  call comm_sortdest_singular
250 
251  allocate( req_list( recv_nmax_r2r + send_nmax_r2r &
252  + recv_nmax_p2r + send_nmax_p2r &
253  + recv_nmax_r2p + send_nmax_r2p ) )
254 
255  if( testonly ) call comm_debugtest
256 
257  return
258  end subroutine comm_setup
259 
260  !-----------------------------------------------------------------------------
262  subroutine comm_list_generate
263  use scale_prc, only: &
264  prc_myrank
265  use scale_prc_icoa, only: &
266  prc_rgn_local, &
267  prc_rgn_r2lp, &
268  prc_rgn_l2r, &
272  i_prc, &
273  i_rgnid, &
274  i_dir, &
275  i_sw, &
276  i_nw, &
277  i_ne, &
278  i_se, &
279  i_w, &
280  i_n, &
281  i_e, &
282  i_s
283  implicit none
284 
285  integer :: ginnar
286 
287  integer :: prc, prc_rmt
288  integer :: rgnid, rgnid_rmt
289  integer :: i, j, i_rmt, j_rmt
290 
291  integer :: n, l, cnt
292  !---------------------------------------------------------------------------
293 
294  ginnar = adm_gmax - adm_gmin + 1
295 
296  allocate( rellist(rellist_vindex,adm_gall*adm_lall) )
297 
298  cnt = 0
299  do l = 1, prc_rgn_local
300  rgnid = prc_rgn_l2r(l)
301  prc = prc_myrank
302 
303  !---< South West >---
304 
305  ! NE -> SW halo
306  if ( prc_rgn_edge_tab(i_dir,i_sw,rgnid) == i_ne ) then
307  rgnid_rmt = prc_rgn_edge_tab(i_rgnid,i_sw,rgnid)
308  prc_rmt = prc_rgn_r2lp(i_prc,rgnid_rmt)
309 
310  do n = 1, ginnar
311  cnt = cnt + 1
312 
313  i = adm_gmin - 1 + n
314  j = adm_gmin - 1
315  i_rmt = adm_gmin - 1 + n
316  j_rmt = adm_gmax
317 
318  rellist(i_recv_grid,cnt) = suf(i,j)
319  rellist(i_recv_rgn, cnt) = rgnid
320  rellist(i_recv_prc, cnt) = prc
321  rellist(i_send_grid,cnt) = suf(i_rmt,j_rmt)
322  rellist(i_send_rgn, cnt) = rgnid_rmt
323  rellist(i_send_prc, cnt) = prc_rmt
324  enddo
325  endif
326 
327  ! SE -> SW halo (Southern Hemisphere, Edge of diamond)
328  if ( prc_rgn_edge_tab(i_dir,i_sw,rgnid) == i_se ) then
329  rgnid_rmt = prc_rgn_edge_tab(i_rgnid,i_sw,rgnid)
330  prc_rmt = prc_rgn_r2lp(i_prc,rgnid_rmt)
331 
332  do n = 1, ginnar
333  cnt = cnt + 1
334 
335  i = adm_gmin - 1 + n
336  j = adm_gmin - 1
337  i_rmt = adm_gmax
338  j_rmt = adm_gmax + 1 - n ! reverse order
339 
340  rellist(i_recv_grid,cnt) = suf(i,j)
341  rellist(i_recv_rgn, cnt) = rgnid
342  rellist(i_recv_prc, cnt) = prc
343  rellist(i_send_grid,cnt) = suf(i_rmt,j_rmt)
344  rellist(i_send_rgn, cnt) = rgnid_rmt
345  rellist(i_send_prc, cnt) = prc_rmt
346  enddo
347  endif
348 
349  !---< North West >---
350 
351  ! SE -> NW
352  if ( prc_rgn_edge_tab(i_dir,i_nw,rgnid) == i_se ) then
353  rgnid_rmt = prc_rgn_edge_tab(i_rgnid,i_nw,rgnid)
354  prc_rmt = prc_rgn_r2lp(i_prc,rgnid_rmt)
355 
356  do n = 1, ginnar
357  cnt = cnt + 1
358 
359  i = adm_gmin - 1
360  j = adm_gmin - 1 + n
361  i_rmt = adm_gmax
362  j_rmt = adm_gmin - 1 + n
363 
364  rellist(i_recv_grid,cnt) = suf(i,j)
365  rellist(i_recv_rgn, cnt) = rgnid
366  rellist(i_recv_prc, cnt) = prc
367  rellist(i_send_grid,cnt) = suf(i_rmt,j_rmt)
368  rellist(i_send_rgn, cnt) = rgnid_rmt
369  rellist(i_send_prc, cnt) = prc_rmt
370  enddo
371  endif
372 
373  ! NE -> NW (Northern Hemisphere, Edge of diamond)
374  if ( prc_rgn_edge_tab(i_dir,i_nw,rgnid) == i_ne ) then
375  rgnid_rmt = prc_rgn_edge_tab(i_rgnid,i_nw,rgnid)
376  prc_rmt = prc_rgn_r2lp(i_prc,rgnid_rmt)
377 
378  do n = 1, ginnar
379  cnt = cnt + 1
380 
381  i = adm_gmin - 1
382  j = adm_gmin - 1 + n
383  i_rmt = adm_gmax + 1 - n ! reverse order
384  j_rmt = adm_gmax
385 
386  rellist(i_recv_grid,cnt) = suf(i,j)
387  rellist(i_recv_rgn, cnt) = rgnid
388  rellist(i_recv_prc, cnt) = prc
389  rellist(i_send_grid,cnt) = suf(i_rmt,j_rmt)
390  rellist(i_send_rgn, cnt) = rgnid_rmt
391  rellist(i_send_prc, cnt) = prc_rmt
392  enddo
393  endif
394 
395  !---< North East >---
396 
397  ! SW -> NE
398  if ( prc_rgn_edge_tab(i_dir,i_ne,rgnid) == i_sw ) then
399  rgnid_rmt = prc_rgn_edge_tab(i_rgnid,i_ne,rgnid)
400  prc_rmt = prc_rgn_r2lp(i_prc,rgnid_rmt)
401 
402  do n = 1, ginnar
403  cnt = cnt + 1
404 
405  i = adm_gmin - 1 + n
406  j = adm_gmax + 1
407  i_rmt = adm_gmin - 1 + n
408  j_rmt = adm_gmin
409 
410  rellist(i_recv_grid,cnt) = suf(i,j)
411  rellist(i_recv_rgn, cnt) = rgnid
412  rellist(i_recv_prc, cnt) = prc
413  rellist(i_send_grid,cnt) = suf(i_rmt,j_rmt)
414  rellist(i_send_rgn, cnt) = rgnid_rmt
415  rellist(i_send_prc, cnt) = prc_rmt
416  enddo
417  endif
418 
419  ! NW -> NE (Northern Hemisphere, Edge of diamond)
420  if ( prc_rgn_edge_tab(i_dir,i_ne,rgnid) == i_nw ) then
421  rgnid_rmt = prc_rgn_edge_tab(i_rgnid,i_ne,rgnid)
422  prc_rmt = prc_rgn_r2lp(i_prc,rgnid_rmt)
423 
424  do n = 1, ginnar
425  cnt = cnt + 1
426 
427  i = adm_gmin + n ! shift 1 grid
428  j = adm_gmax + 1
429  i_rmt = adm_gmin
430  j_rmt = adm_gmax + 1 - n ! reverse order
431 
432  rellist(i_recv_grid,cnt) = suf(i,j)
433  rellist(i_recv_rgn, cnt) = rgnid
434  rellist(i_recv_prc, cnt) = prc
435  rellist(i_send_grid,cnt) = suf(i_rmt,j_rmt)
436  rellist(i_send_rgn, cnt) = rgnid_rmt
437  rellist(i_send_prc, cnt) = prc_rmt
438  enddo
439  endif
440 
441  !---< South East >---
442 
443  ! NW -> SE
444  if ( prc_rgn_edge_tab(i_dir,i_se,rgnid) == i_nw ) then
445  rgnid_rmt = prc_rgn_edge_tab(i_rgnid,i_se,rgnid)
446  prc_rmt = prc_rgn_r2lp(i_prc,rgnid_rmt)
447 
448  do n = 1, ginnar
449  cnt = cnt + 1
450 
451  i = adm_gmax + 1
452  j = adm_gmin - 1 + n
453  i_rmt = adm_gmin
454  j_rmt = adm_gmin - 1 + n
455 
456  rellist(i_recv_grid,cnt) = suf(i,j)
457  rellist(i_recv_rgn, cnt) = rgnid
458  rellist(i_recv_prc, cnt) = prc
459  rellist(i_send_grid,cnt) = suf(i_rmt,j_rmt)
460  rellist(i_send_rgn, cnt) = rgnid_rmt
461  rellist(i_send_prc, cnt) = prc_rmt
462  enddo
463  endif
464 
465  ! SW -> SE (Southern Hemisphere, Edge of diamond)
466  if ( prc_rgn_edge_tab(i_dir,i_se,rgnid) == i_sw ) then
467  rgnid_rmt = prc_rgn_edge_tab(i_rgnid,i_se,rgnid)
468  prc_rmt = prc_rgn_r2lp(i_prc,rgnid_rmt)
469 
470  do n = 1, ginnar
471  cnt = cnt + 1
472 
473  i = adm_gmax + 1
474  j = adm_gmin + n ! shift 1 grid
475  i_rmt = adm_gmax + 1 - n ! reverse order
476  j_rmt = adm_gmin
477 
478  rellist(i_recv_grid,cnt) = suf(i,j)
479  rellist(i_recv_rgn, cnt) = rgnid
480  rellist(i_recv_prc, cnt) = prc
481  rellist(i_send_grid,cnt) = suf(i_rmt,j_rmt)
482  rellist(i_send_rgn, cnt) = rgnid_rmt
483  rellist(i_send_prc, cnt) = prc_rmt
484  enddo
485  endif
486 
487  !---< Vertex : link to the next next region >---
488 
489  ! West Vertex
490  if ( prc_rgn_vert_num(i_w,rgnid) == 4 ) then
491  rgnid_rmt = prc_rgn_vert_tab(i_rgnid,i_w,rgnid,3)
492  prc_rmt = prc_rgn_r2lp(i_prc,rgnid_rmt)
493 
494  cnt = cnt + 1
495 
496  i = adm_gmin - 1
497  j = adm_gmin - 1
498 
499  if ( prc_rgn_vert_tab(i_dir,i_w,rgnid,3) == i_n ) then
500  i_rmt = adm_gmin
501  j_rmt = adm_gmax
502  elseif( prc_rgn_vert_tab(i_dir,i_w,rgnid,3) == i_e ) then
503  i_rmt = adm_gmax
504  j_rmt = adm_gmax
505  elseif( prc_rgn_vert_tab(i_dir,i_w,rgnid,3) == i_s ) then
506  i_rmt = adm_gmax
507  j_rmt = adm_gmin
508  endif
509 
510  rellist(i_recv_grid,cnt) = suf(i,j)
511  rellist(i_recv_rgn, cnt) = rgnid
512  rellist(i_recv_prc, cnt) = prc
513  rellist(i_send_grid,cnt) = suf(i_rmt,j_rmt)
514  rellist(i_send_rgn, cnt) = rgnid_rmt
515  rellist(i_send_prc, cnt) = prc_rmt
516  endif
517 
518  ! North Vertex
519  if ( prc_rgn_vert_num(i_n,rgnid) == 4 ) then
520  rgnid_rmt = prc_rgn_vert_tab(i_rgnid,i_n,rgnid,3)
521  prc_rmt = prc_rgn_r2lp(i_prc,rgnid_rmt)
522 
523  ! known as north pole point
524  if ( prc_rgn_vert_tab(i_dir,i_n,rgnid,3) == i_w ) then
525  cnt = cnt + 1
526 
527  i = adm_gmin
528  j = adm_gmax + 1
529  i_rmt = adm_gmin
530  j_rmt = adm_gmin
531 
532  rellist(i_recv_grid,cnt) = suf(i,j)
533  rellist(i_recv_rgn, cnt) = rgnid
534  rellist(i_recv_prc, cnt) = prc
535  rellist(i_send_grid,cnt) = suf(i_rmt,j_rmt)
536  rellist(i_send_rgn, cnt) = rgnid_rmt
537  rellist(i_send_prc, cnt) = prc_rmt
538  endif
539 
540  ! unused vertex point
541  cnt = cnt + 1
542 
543  i = adm_gmin - 1
544  j = adm_gmax + 1
545 
546  if ( prc_rgn_vert_tab(i_dir,i_n,rgnid,3) == i_w ) then
547  i_rmt = adm_gmin
548  j_rmt = adm_gmin + 1
549  elseif( prc_rgn_vert_tab(i_dir,i_n,rgnid,3) == i_n ) then
550  i_rmt = adm_gmin
551  j_rmt = adm_gmax
552  elseif( prc_rgn_vert_tab(i_dir,i_n,rgnid,3) == i_e ) then
553  i_rmt = adm_gmax
554  j_rmt = adm_gmax
555  elseif( prc_rgn_vert_tab(i_dir,i_n,rgnid,3) == i_s ) then
556  i_rmt = adm_gmax
557  j_rmt = adm_gmin
558  endif
559 
560  rellist(i_recv_grid,cnt) = suf(i,j)
561  rellist(i_recv_rgn, cnt) = rgnid
562  rellist(i_recv_prc, cnt) = prc
563  rellist(i_send_grid,cnt) = suf(i_rmt,j_rmt)
564  rellist(i_send_rgn, cnt) = rgnid_rmt
565  rellist(i_send_prc, cnt) = prc_rmt
566  endif
567 
568  ! East Vertex
569  if ( prc_rgn_vert_num(i_e,rgnid) == 4 ) then
570  if ( prc_rgn_vert_tab(i_dir,i_e,rgnid,3) == i_w ) then
571  rgnid_rmt = prc_rgn_vert_tab(i_rgnid,i_e,rgnid,3)
572  prc_rmt = prc_rgn_r2lp(i_prc,rgnid_rmt)
573 
574  cnt = cnt + 1
575 
576  i = adm_gmax + 1
577  j = adm_gmax + 1
578  i_rmt = adm_gmin
579  j_rmt = adm_gmin
580 
581  rellist(i_recv_grid,cnt) = suf(i,j)
582  rellist(i_recv_rgn, cnt) = rgnid
583  rellist(i_recv_prc, cnt) = prc
584  rellist(i_send_grid,cnt) = suf(i_rmt,j_rmt)
585  rellist(i_send_rgn, cnt) = rgnid_rmt
586  rellist(i_send_prc, cnt) = prc_rmt
587  endif
588  endif
589 
590  ! South Vertex
591  if ( prc_rgn_vert_num(i_s,rgnid) == 4 ) then
592  rgnid_rmt = prc_rgn_vert_tab(i_rgnid,i_s,rgnid,3)
593  prc_rmt = prc_rgn_r2lp(i_prc,rgnid_rmt)
594 
595  ! known as south pole point
596  if ( prc_rgn_vert_tab(i_dir,i_s,rgnid,3) == i_w ) then
597  cnt = cnt + 1
598 
599  i = adm_gmax + 1
600  j = adm_gmin
601  i_rmt = adm_gmin
602  j_rmt = adm_gmin
603 
604  rellist(i_recv_grid,cnt) = suf(i,j)
605  rellist(i_recv_rgn, cnt) = rgnid
606  rellist(i_recv_prc, cnt) = prc
607  rellist(i_send_grid,cnt) = suf(i_rmt,j_rmt)
608  rellist(i_send_rgn, cnt) = rgnid_rmt
609  rellist(i_send_prc, cnt) = prc_rmt
610  endif
611 
612  ! unused vertex point
613  cnt = cnt + 1
614 
615  i = adm_gmax + 1
616  j = adm_gmin - 1
617 
618  if ( prc_rgn_vert_tab(i_dir,i_s,rgnid,3) == i_w ) then
619  i_rmt = adm_gmin + 1
620  j_rmt = adm_gmin
621  elseif( prc_rgn_vert_tab(i_dir,i_s,rgnid,3) == i_n ) then
622  i_rmt = adm_gmin
623  j_rmt = adm_gmax
624  elseif( prc_rgn_vert_tab(i_dir,i_s,rgnid,3) == i_e ) then
625  i_rmt = adm_gmax
626  j_rmt = adm_gmax
627  elseif( prc_rgn_vert_tab(i_dir,i_s,rgnid,3) == i_s ) then
628  i_rmt = adm_gmax
629  j_rmt = adm_gmin
630  endif
631 
632  rellist(i_recv_grid,cnt) = suf(i,j)
633  rellist(i_recv_rgn, cnt) = rgnid
634  rellist(i_recv_prc, cnt) = prc
635  rellist(i_send_grid,cnt) = suf(i_rmt,j_rmt)
636  rellist(i_send_rgn, cnt) = rgnid_rmt
637  rellist(i_send_prc, cnt) = prc_rmt
638  endif
639  enddo ! loop l
640 
641  rellist_nmax = cnt
642 
643  if( io_l ) write(io_fid_log,*)
644  if( io_l ) write(io_fid_log,*) '*** rellist_nmax:', rellist_nmax
645 
646  if ( debug ) then
647  if( io_l ) write(io_fid_log,*) '--- Relation Table'
648  if( io_l ) write(io_fid_log,'(7(A10))') 'Count', '|recv_grid', '| recv_rgn', '| recv_prc', &
649  '|send_grid', '| send_rgn', '| send_prc'
650  do cnt = 1, rellist_nmax
651  if( io_l ) write(io_fid_log,'(7(I10))') cnt, rellist(:,cnt)
652  enddo
653  endif
654 
655  return
656  end subroutine comm_list_generate
657 
658  !-----------------------------------------------------------------------------
660  subroutine comm_sortdest
661  use scale_prc, only: &
663  prc_myrank, &
664  prc_nprocs
665  use scale_prc_icoa, only: &
666  prc_rgn_lp2r, &
667  prc_rgn_r2lp, &
668  i_l
669  implicit none
670 
671  integer :: sendbuf1(1)
672  integer :: recvbuf1(1)
673 
674  integer, allocatable :: sendbuf_info(:)
675  integer, allocatable :: recvbuf_info(:)
676  integer, allocatable :: sendbuf_list(:,:,:)
677  integer, allocatable :: recvbuf_list(:,:,:)
678  integer, allocatable :: REQ_list_r2r(:)
679 
680  integer :: Recv_nglobal_r2r
681  integer :: Send_size_nglobal
682 
683  integer :: cnt, irank, ipos
684  integer :: totalsize, rank, tag
685 
686  integer :: i_from, j_from, r_from, g_from, l_from, p_from
687  integer :: i_to, j_to, r_to, g_to, l_to, p_to
688 
689  integer :: ierr
690  integer :: n, p
691  !---------------------------------------------------------------------------
692 
693  allocate( copy_info_r2r(info_vindex) )
694  allocate( recv_info_r2r(info_vindex,recv_nlim) )
695  allocate( send_info_r2r(info_vindex,send_nlim) )
696  copy_info_r2r(:) = -1
697  recv_info_r2r(:,:) = -1
698  send_info_r2r(:,:) = -1
699  copy_info_r2r(i_size) = 0
700  recv_info_r2r(i_size,:) = 0
701  send_info_r2r(i_size,:) = 0
702 
703  allocate( copy_list_r2r(list_vindex,rellist_nmax) )
704  allocate( recv_list_r2r(list_vindex,rellist_nmax,recv_nlim) )
705  allocate( send_list_r2r(list_vindex,rellist_nmax,send_nlim) )
706  copy_list_r2r(:,:) = -1
707  recv_list_r2r(:,:,:) = -1
708  send_list_r2r(:,:,:) = -1
709 
710  ! sorting list according to destination
711  do cnt = 1, rellist_nmax
712 
713  if ( rellist(i_recv_prc,cnt) == rellist(i_send_prc,cnt) ) then ! no communication
714  copy_info_r2r(i_size) = copy_info_r2r(i_size) + 1
715  ipos = copy_info_r2r(i_size)
716 
717  copy_list_r2r(i_grid_from,ipos) = rellist(i_send_grid,cnt)
718  copy_list_r2r(i_l_from ,ipos) = prc_rgn_r2lp(i_l,rellist(i_send_rgn,cnt))
719  copy_list_r2r(i_grid_to ,ipos) = rellist(i_recv_grid,cnt)
720  copy_list_r2r(i_l_to ,ipos) = prc_rgn_r2lp(i_l,rellist(i_recv_rgn,cnt))
721  else ! node-to-node communication
722  !--- search existing rank id (identify by prc_from)
723  irank = -1
724  do n = 1, recv_nmax_r2r
725  if ( recv_info_r2r(i_prc_from,n) == rellist(i_send_prc,cnt) ) then
726  irank = n
727  exit
728  endif
729  enddo
730 
731  if ( irank < 0 ) then ! register new rank id
732  recv_nmax_r2r = recv_nmax_r2r + 1
733  irank = recv_nmax_r2r
734 
735  recv_info_r2r(i_prc_from,irank) = rellist(i_send_prc,cnt)
736  recv_info_r2r(i_prc_to ,irank) = rellist(i_recv_prc,cnt)
737  endif
738 
739  recv_info_r2r(i_size,irank) = recv_info_r2r(i_size,irank) + 1
740  ipos = recv_info_r2r(i_size,irank)
741 
742  recv_list_r2r(i_grid_from,ipos,irank) = rellist(i_send_grid,cnt)
743  recv_list_r2r(i_l_from ,ipos,irank) = prc_rgn_r2lp(i_l,rellist(i_send_rgn,cnt))
744  recv_list_r2r(i_grid_to ,ipos,irank) = rellist(i_recv_grid,cnt)
745  recv_list_r2r(i_l_to ,ipos,irank) = prc_rgn_r2lp(i_l,rellist(i_recv_rgn,cnt))
746  endif
747  enddo
748 
749  if ( copy_info_r2r(i_size) > 0 ) then
750  copy_nmax_r2r = 1
751  copy_info_r2r(i_prc_from) = prc_myrank
752  copy_info_r2r(i_prc_to ) = prc_myrank
753  endif
754 
755 
756 
757  ! get maximum number of rank to communication
758  sendbuf1(1) = recv_nmax_r2r
759 
760  call mpi_allreduce( sendbuf1(1), & ! [IN]
761  recvbuf1(1), & ! [OUT]
762  1, & ! [IN]
763  mpi_integer, & ! [IN]
764  mpi_max, & ! [IN]
765  prc_local_comm_world, & ! [IN]
766  ierr ) ! [OUT]
767 
768  recv_nglobal_r2r = recvbuf1(1)
769 
770  ! Distribute receive request from each rank to all
771  allocate( sendbuf_info(info_vindex*recv_nglobal_r2r) )
772  allocate( recvbuf_info(info_vindex*recv_nglobal_r2r*prc_nprocs) )
773  sendbuf_info(:) = -1
774 
775  do irank = 1, recv_nmax_r2r
776  n = (irank-1) * info_vindex
777 
778  sendbuf_info(n+i_size ) = recv_info_r2r(i_size ,irank)
779  sendbuf_info(n+i_prc_from) = recv_info_r2r(i_prc_from,irank)
780  sendbuf_info(n+i_prc_to ) = recv_info_r2r(i_prc_to ,irank)
781  enddo
782 
783  totalsize = info_vindex * recv_nglobal_r2r
784 
785  if ( totalsize > 0 ) then
786  call mpi_allgather( sendbuf_info(1), & ! [IN]
787  totalsize, & ! [IN]
788  mpi_integer, & ! [IN]
789  recvbuf_info(1), & ! [OUT]
790  totalsize, & ! [IN]
791  mpi_integer, & ! [IN]
792  prc_local_comm_world, & ! [IN]
793  ierr ) ! [OUT]
794  endif
795 
796  send_size_nglobal = 0
797 
798  ! Accept receive request to my rank
799  do p = 1, recv_nglobal_r2r*prc_nprocs
800  n = (p-1) * info_vindex
801 
802  if ( recvbuf_info(n+i_prc_from) == prc_myrank ) then
803  send_nmax_r2r = send_nmax_r2r + 1
804  irank = send_nmax_r2r
805 
806  send_info_r2r(i_size ,irank) = recvbuf_info(n+i_size )
807  send_info_r2r(i_prc_from,irank) = recvbuf_info(n+i_prc_from)
808  send_info_r2r(i_prc_to ,irank) = recvbuf_info(n+i_prc_to )
809  endif
810 
811  send_size_nglobal = max( send_size_nglobal, recvbuf_info(n+i_size) )
812  enddo
813 
814  if( io_l ) write(io_fid_log,*)
815  if( io_l ) write(io_fid_log,*) '*** Recv_nmax_r2r(global) = ', recv_nglobal_r2r
816  if( io_l ) write(io_fid_log,*) '*** Recv_nmax_r2r(local) = ', recv_nmax_r2r
817  if( io_l ) write(io_fid_log,*) '*** Send_nmax_r2r(local) = ', send_nmax_r2r
818  if( io_l ) write(io_fid_log,*) '*** Send_size_r2r(global) = ', send_size_nglobal
819  if( io_l ) write(io_fid_log,*)
820  if( io_l ) write(io_fid_log,'(A)') '|---------------------------------------'
821  if( io_l ) write(io_fid_log,'(A)') '| size prc_from prc_to'
822  if( io_l ) write(io_fid_log,'(A10,3(I10))') '| Copy_r2r', copy_info_r2r(:)
823  do irank = 1, recv_nmax_r2r
824  if( io_l ) write(io_fid_log,'(A10,3(I10))') '| Recv_r2r', recv_info_r2r(:,irank)
825  enddo
826  do irank = 1, send_nmax_r2r
827  if( io_l ) write(io_fid_log,'(A10,3(I10))') '| Send_r2r', send_info_r2r(:,irank)
828  enddo
829 
830  ! Communicate detailed information in each pair
831  allocate( req_list_r2r(recv_nmax_r2r+send_nmax_r2r) )
832 
833  allocate( sendbuf_list(list_vindex,send_size_nglobal,recv_nmax_r2r) )
834  allocate( recvbuf_list(list_vindex,send_size_nglobal,send_nmax_r2r) )
835  sendbuf_list(:,:,:) = -1
836 
837  req_count = 0
838 
839  do irank = 1, send_nmax_r2r
840  req_count = req_count + 1
841  totalsize = send_info_r2r(i_size,irank) * list_vindex
842  rank = send_info_r2r(i_prc_to ,irank)
843  tag = send_info_r2r(i_prc_from,irank)
844 
845  call mpi_irecv( recvbuf_list(1,1,irank), & ! [OUT]
846  totalsize, & ! [IN]
847  mpi_integer, & ! [IN]
848  rank, & ! [IN]
849  tag, & ! [IN]
850  prc_local_comm_world, & ! [IN]
851  req_list_r2r(req_count), & ! [OUT]
852  ierr ) ! [OUT]
853  enddo
854 
855  do irank = 1, recv_nmax_r2r
856  do ipos = 1, recv_info_r2r(i_size,irank)
857  sendbuf_list(:,ipos,irank) = recv_list_r2r(:,ipos,irank)
858  enddo
859 
860  req_count = req_count + 1
861  totalsize = recv_info_r2r(i_size,irank) * list_vindex
862  rank = recv_info_r2r(i_prc_from,irank)
863  tag = recv_info_r2r(i_prc_from,irank)
864 
865  call mpi_isend( sendbuf_list(1,1,irank), & ! [IN]
866  totalsize, & ! [IN]
867  mpi_integer, & ! [IN]
868  rank, & ! [IN]
869  tag, & ! [IN]
870  prc_local_comm_world, & ! [IN]
871  req_list_r2r(req_count), & ! [OUT]
872  ierr ) ! [OUT]
873  enddo
874 
875  !--- wait packets
876  if ( recv_nmax_r2r+send_nmax_r2r > 0 ) then
877  call mpi_waitall( recv_nmax_r2r+send_nmax_r2r, & ! [IN]
878  req_list_r2r(:), & ! [IN]
879  mpi_statuses_ignore, & ! [OUT]
880  ierr ) ! [OUT]
881  endif
882 
883  do irank = 1, send_nmax_r2r
884  do ipos = 1, send_info_r2r(i_size,irank)
885  send_list_r2r(:,ipos,irank) = recvbuf_list(:,ipos,irank)
886  enddo
887  enddo
888 
889  if ( debug ) then
890  if( io_l ) write(io_fid_log,*)
891  if( io_l ) write(io_fid_log,*) '--- Copy_list_r2r'
892  if( io_l ) write(io_fid_log,*)
893  if( io_l ) write(io_fid_log,'(13(A6))') 'number', &
894  '|ifrom','|jfrom','|rfrom','|gfrom','|lfrom','|pfrom', &
895  '| ito','| jto','| rto','| gto','| lto','| pto'
896  do ipos = 1, copy_info_r2r(i_size)
897  g_from = copy_list_r2r(i_grid_from,ipos)
898  l_from = copy_list_r2r(i_l_from,ipos)
899  p_from = copy_info_r2r(i_prc_from)
900  i_from = mod(g_from-1,adm_gall_1d) + 1
901  j_from = (g_from-i_from) / adm_gall_1d + 1
902  r_from = prc_rgn_lp2r(l_from,p_from)
903  g_to = copy_list_r2r(i_grid_to,ipos)
904  l_to = copy_list_r2r(i_l_to,ipos)
905  p_to = copy_info_r2r(i_prc_to)
906  i_to = mod(g_to-1,adm_gall_1d) + 1
907  j_to = (g_to-i_to) / adm_gall_1d + 1
908  r_to = prc_rgn_lp2r(l_to,p_to)
909 
910  if( io_l ) write(io_fid_log,'(13(I6))') ipos, i_from, j_from, r_from, g_from, l_from, p_from, &
911  i_to , j_to , r_to , g_to , l_to , p_to
912  enddo
913 
914  if( io_l ) write(io_fid_log,*)
915  if( io_l ) write(io_fid_log,*) '--- Recv_list_r2r'
916  do irank = 1, recv_nmax_r2r
917  if( io_l ) write(io_fid_log,'(13(A6))') 'number', &
918  '|ifrom','|jfrom','|rfrom','|gfrom','|lfrom','|pfrom', &
919  '| ito','| jto','| rto','| gto','| lto','| pto'
920  do ipos = 1, recv_info_r2r(i_size,irank)
921  g_from = recv_list_r2r(i_grid_from,ipos,irank)
922  l_from = recv_list_r2r(i_l_from,ipos,irank)
923  p_from = recv_info_r2r(i_prc_from,irank)
924  i_from = mod(g_from-1,adm_gall_1d) + 1
925  j_from = (g_from-i_from) / adm_gall_1d + 1
926  r_from = prc_rgn_lp2r(l_from,p_from)
927  g_to = recv_list_r2r(i_grid_to,ipos,irank)
928  l_to = recv_list_r2r(i_l_to,ipos,irank)
929  p_to = recv_info_r2r(i_prc_to,irank)
930  i_to = mod(g_to-1,adm_gall_1d) + 1
931  j_to = (g_to-i_to) / adm_gall_1d + 1
932  r_to = prc_rgn_lp2r(l_to,p_to)
933 
934  if( io_l ) write(io_fid_log,'(13(I6))') ipos, i_from, j_from, r_from, g_from, l_from, p_from, &
935  i_to , j_to , r_to , g_to , l_to , p_to
936  enddo
937  enddo
938 
939  if( io_l ) write(io_fid_log,*)
940  if( io_l ) write(io_fid_log,*) '--- Send_list_r2r'
941  do irank = 1, send_nmax_r2r
942  if( io_l ) write(io_fid_log,'(13(A6))') 'number', &
943  '|ifrom','|jfrom','|rfrom','|gfrom','|lfrom','|pfrom', &
944  '| ito','| jto','| rto','| gto','| lto','| pto'
945  do ipos = 1, send_info_r2r(i_size,irank)
946  g_from = send_list_r2r(i_grid_from,ipos,irank)
947  l_from = send_list_r2r(i_l_from,ipos,irank)
948  p_from = send_info_r2r(i_prc_from,irank)
949  i_from = mod(g_from-1,adm_gall_1d) + 1
950  j_from = (g_from-i_from) / adm_gall_1d + 1
951  r_from = prc_rgn_lp2r(l_from,p_from)
952  g_to = send_list_r2r(i_grid_to,ipos,irank)
953  l_to = send_list_r2r(i_l_to,ipos,irank)
954  p_to = send_info_r2r(i_prc_to,irank)
955  i_to = mod(g_to-1,adm_gall_1d) + 1
956  j_to = (g_to-i_to) / adm_gall_1d + 1
957  r_to = prc_rgn_lp2r(l_to,p_to)
958 
959  if( io_l ) write(io_fid_log,'(13(I6))') ipos, i_from, j_from, r_from, g_from, l_from, p_from, &
960  i_to , j_to , r_to , g_to , l_to , p_to
961  enddo
962  enddo
963  endif
964 
965  allocate( sendbuf_r2r_sp(send_size_nglobal*adm_kall*comm_varmax,send_nmax_r2r) )
966  allocate( recvbuf_r2r_sp(send_size_nglobal*adm_kall*comm_varmax,recv_nmax_r2r) )
967  allocate( sendbuf_r2r_dp(send_size_nglobal*adm_kall*comm_varmax,send_nmax_r2r) )
968  allocate( recvbuf_r2r_dp(send_size_nglobal*adm_kall*comm_varmax,recv_nmax_r2r) )
969 
970  return
971  end subroutine comm_sortdest
972 
973  !-----------------------------------------------------------------------------
975  subroutine comm_sortdest_pl
976  use scale_prc, only: &
977  prc_myrank
978  use scale_prc_icoa, only: &
979  prc_have_pl, &
980  prc_rgn_local, &
983  prc_rgn_lp2r, &
984  prc_rgn_r2lp, &
985  prc_rgn_l2r, &
986  prc_rgn_r2p_pl, &
987  i_l, &
988  i_prc, &
989  i_rgnid, &
990  i_n, &
991  i_s, &
992  i_npl, &
993  i_spl
994  implicit none
995 
996  integer :: prc, prc_rmt
997  integer :: rgnid, rgnid_rmt
998  integer :: check_vert_num
999  integer :: pl_to
1000 
1001  integer :: irank, ipos
1002 
1003  integer, parameter :: Send_size_nglobal_pl = 10
1004 
1005  integer :: l, l_pl, n, v, vv
1006  integer :: i_from, j_from, r_from, g_from, l_from, p_from
1007  integer :: i_to, j_to, r_to, g_to, l_to, p_to
1008  !---------------------------------------------------------------------------
1009 
1010  allocate( copy_info_p2r(info_vindex) )
1011  allocate( recv_info_p2r(info_vindex,recv_nlim) )
1012  allocate( send_info_p2r(info_vindex,send_nlim) )
1013  copy_info_p2r(:) = -1
1014  recv_info_p2r(:,:) = -1
1015  send_info_p2r(:,:) = -1
1016  copy_info_p2r(i_size) = 0
1017  recv_info_p2r(i_size,:) = 0
1018  send_info_p2r(i_size,:) = 0
1019 
1020  allocate( copy_list_p2r(list_vindex,send_size_nglobal_pl) )
1021  allocate( recv_list_p2r(list_vindex,send_size_nglobal_pl,recv_nlim) )
1022  allocate( send_list_p2r(list_vindex,send_size_nglobal_pl,send_nlim) )
1023  copy_list_p2r(:,:) = -1
1024  recv_list_p2r(:,:,:) = -1
1025  send_list_p2r(:,:,:) = -1
1026 
1027  allocate( copy_info_r2p(info_vindex) )
1028  allocate( recv_info_r2p(info_vindex,recv_nlim) )
1029  allocate( send_info_r2p(info_vindex,send_nlim) )
1030  copy_info_r2p(:) = -1
1031  recv_info_r2p(:,:) = -1
1032  send_info_r2p(:,:) = -1
1033  copy_info_r2p(i_size) = 0
1034  recv_info_r2p(i_size,:) = 0
1035  send_info_r2p(i_size,:) = 0
1036 
1037  allocate( copy_list_r2p(list_vindex,send_size_nglobal_pl) )
1038  allocate( recv_list_r2p(list_vindex,send_size_nglobal_pl,recv_nlim) )
1039  allocate( send_list_r2p(list_vindex,send_size_nglobal_pl,send_nlim) )
1040  copy_list_r2p(:,:) = -1
1041  recv_list_r2p(:,:,:) = -1
1042  send_list_r2p(:,:,:) = -1
1043 
1044  ! Search in regular region
1045  do l = 1, prc_rgn_local
1046  rgnid = prc_rgn_l2r(l)
1047  prc = prc_myrank
1048 
1049  do l_pl = i_npl, i_spl
1050  rgnid_rmt = l_pl
1051  prc_rmt = prc_rgn_r2p_pl(rgnid_rmt)
1052 
1053  if ( rgnid_rmt == i_npl ) then
1054  check_vert_num = prc_rgn_vert_num(i_n,rgnid)
1055  i_from = adm_gmin
1056  j_from = adm_gmax
1057  i_to = adm_gmin
1058  j_to = adm_gmax + 1
1059  elseif( rgnid_rmt == i_spl ) then
1060  check_vert_num = prc_rgn_vert_num(i_s,rgnid)
1061  i_from = adm_gmax
1062  j_from = adm_gmin
1063  i_to = adm_gmax + 1
1064  j_to = adm_gmin
1065  endif
1066 
1067  if ( check_vert_num == adm_vlink ) then
1068 
1069  ! search destination in the pole halo
1070  do v = 1, adm_vlink
1071  if ( rgnid == prc_rgn_vert_tab_pl(i_rgnid,rgnid_rmt,v) ) then
1072  pl_to = v
1073  if( pl_to < adm_gmin_pl ) pl_to = adm_gmax_pl
1074  exit
1075  endif
1076  enddo
1077 
1078  if ( prc == prc_rmt ) then ! no communication
1079  ! copy region inner -> pole halo
1080  copy_info_r2p(i_size) = copy_info_r2p(i_size) + 1
1081  ipos = copy_info_r2p(i_size)
1082 
1083  copy_list_r2p(i_grid_from,ipos) = suf(i_from,j_from)
1084  copy_list_r2p(i_l_from ,ipos) = l
1085  copy_list_r2p(i_grid_to ,ipos) = pl_to
1086  copy_list_r2p(i_l_to ,ipos) = l_pl
1087  else ! node-to-node communication
1088  ! receive pole center -> region halo
1089  irank = -1
1090  do n = 1, recv_nmax_p2r
1091  if ( recv_info_p2r(i_prc_from,n) == prc_rmt ) then
1092  irank = n
1093  exit
1094  endif
1095  enddo
1096 
1097  if ( irank < 0 ) then ! register new rank id
1098  recv_nmax_p2r = recv_nmax_p2r + 1
1099  irank = recv_nmax_p2r
1100 
1101  recv_info_p2r(i_prc_from,irank) = prc_rmt
1102  recv_info_p2r(i_prc_to ,irank) = prc
1103  endif
1104 
1105  recv_info_p2r(i_size,irank) = recv_info_p2r(i_size,irank) + 1
1106  ipos = recv_info_p2r(i_size,irank)
1107 
1108  recv_list_p2r(i_grid_from,ipos,irank) = adm_gslf_pl
1109  recv_list_p2r(i_l_from ,ipos,irank) = l_pl
1110  recv_list_p2r(i_grid_to ,ipos,irank) = suf(i_to,j_to)
1111  recv_list_p2r(i_l_to ,ipos,irank) = l
1112 
1113  ! send region inner -> pole halo
1114  irank = -1
1115  do n = 1, send_nmax_r2p
1116  if ( send_info_r2p(i_prc_to,n) == prc_rmt ) then
1117  irank = n
1118  exit
1119  endif
1120  enddo
1121 
1122  if ( irank < 0 ) then ! register new rank id
1123  send_nmax_r2p = send_nmax_r2p + 1
1124  irank = send_nmax_r2p
1125 
1126  send_info_r2p(i_prc_from,irank) = prc
1127  send_info_r2p(i_prc_to ,irank) = prc_rmt
1128  endif
1129 
1130  send_info_r2p(i_size,irank) = send_info_r2p(i_size,irank) + 1
1131  ipos = send_info_r2p(i_size,irank)
1132 
1133  send_list_r2p(i_grid_from,ipos,irank) = suf(i_from,j_from)
1134  send_list_r2p(i_l_from ,ipos,irank) = l
1135  send_list_r2p(i_grid_to ,ipos,irank) = pl_to
1136  send_list_r2p(i_l_to ,ipos,irank) = l_pl
1137  endif
1138  endif
1139 
1140  enddo ! loop l_pl
1141  enddo ! loop l
1142 
1143  ! Search in pole
1144  if ( prc_have_pl ) then
1145 
1146  do l_pl = i_npl, i_spl
1147  rgnid = l_pl
1148  prc = prc_myrank
1149 
1150  do v = adm_vlink+1, 2, -1
1151  vv = v
1152  if( v == adm_vlink+1 ) vv = 1
1153 
1154  rgnid_rmt = prc_rgn_vert_tab_pl(i_rgnid,rgnid,vv)
1155  prc_rmt = prc_rgn_r2lp(i_prc,rgnid_rmt)
1156 
1157  if ( rgnid == i_npl ) then
1158  i_from = adm_gmin
1159  j_from = adm_gmax
1160  i_to = adm_gmin
1161  j_to = adm_gmax + 1
1162  elseif( rgnid == i_spl ) then
1163  i_from = adm_gmax
1164  j_from = adm_gmin
1165  i_to = adm_gmax + 1
1166  j_to = adm_gmin
1167  endif
1168 
1169  pl_to = vv
1170  if( pl_to < adm_gmin_pl ) pl_to = adm_gmax_pl
1171 
1172  if ( prc == prc_rmt ) then ! no communication
1173  ! copy pole center -> region halo
1174  copy_info_p2r(i_size) = copy_info_p2r(i_size) + 1
1175  ipos = copy_info_p2r(i_size)
1176 
1177  copy_list_p2r(i_grid_from,ipos) = adm_gslf_pl
1178  copy_list_p2r(i_l_from ,ipos) = l_pl
1179  copy_list_p2r(i_grid_to ,ipos) = suf(i_to,j_to)
1180  copy_list_p2r(i_l_to ,ipos) = prc_rgn_r2lp(i_l,rgnid_rmt)
1181  else ! node-to-node communication
1182  ! send region inner -> pole halo
1183  irank = -1
1184  do n = 1, recv_nmax_r2p
1185  if ( recv_info_r2p(i_prc_from,n) == prc_rmt ) then
1186  irank = n
1187  exit
1188  endif
1189  enddo
1190 
1191  if ( irank < 0 ) then ! register new rank id
1192  recv_nmax_r2p = recv_nmax_r2p + 1
1193  irank = recv_nmax_r2p
1194 
1195  recv_info_r2p(i_prc_from,irank) = prc_rmt
1196  recv_info_r2p(i_prc_to ,irank) = prc
1197  endif
1198 
1199  recv_info_r2p(i_size,irank) = recv_info_r2p(i_size,irank) + 1
1200  ipos = recv_info_r2p(i_size,irank)
1201 
1202  recv_list_r2p(i_grid_from,ipos,irank) = suf(i_from,j_from)
1203  recv_list_r2p(i_l_from ,ipos,irank) = prc_rgn_r2lp(i_l,rgnid_rmt)
1204  recv_list_r2p(i_grid_to ,ipos,irank) = pl_to
1205  recv_list_r2p(i_l_to ,ipos,irank) = l_pl
1206 
1207  ! receive pole center -> region halo
1208  irank = -1
1209  do n = 1, send_nmax_p2r
1210  if ( send_info_p2r(i_prc_to,n) == prc_rmt ) then
1211  irank = n
1212  exit
1213  endif
1214  enddo
1215 
1216  if ( irank < 0 ) then ! register new rank id
1217  send_nmax_p2r = send_nmax_p2r + 1
1218  irank = send_nmax_p2r
1219 
1220  send_info_p2r(i_prc_from,irank) = prc
1221  send_info_p2r(i_prc_to ,irank) = prc_rmt
1222  endif
1223 
1224  send_info_p2r(i_size,irank) = send_info_p2r(i_size,irank) + 1
1225  ipos = send_info_p2r(i_size,irank)
1226 
1227  send_list_p2r(i_grid_from,ipos,irank) = adm_gslf_pl
1228  send_list_p2r(i_l_from ,ipos,irank) = l_pl
1229  send_list_p2r(i_grid_to ,ipos,irank) = suf(i_to,j_to)
1230  send_list_p2r(i_l_to ,ipos,irank) = prc_rgn_r2lp(i_l,rgnid_rmt)
1231  endif
1232 
1233  enddo ! loop pole halo
1234 
1235  enddo ! loop l_pl
1236 
1237  if ( copy_info_p2r(i_size) > 0 ) then
1238  copy_nmax_p2r = 1
1239  copy_info_p2r(i_prc_from) = prc_myrank
1240  copy_info_p2r(i_prc_to ) = prc_myrank
1241  endif
1242 
1243  if ( copy_info_r2p(i_size) > 0 ) then
1244  copy_nmax_r2p = 1
1245  copy_info_r2p(i_prc_from) = prc_myrank
1246  copy_info_r2p(i_prc_to ) = prc_myrank
1247  endif
1248 
1249  endif
1250 
1251  if( io_l ) write(io_fid_log,*)
1252  if( io_l ) write(io_fid_log,*) '*** Recv_nmax_p2r(local) = ', recv_nmax_p2r
1253  if( io_l ) write(io_fid_log,*) '*** Send_nmax_p2r(local) = ', send_nmax_p2r
1254  if( io_l ) write(io_fid_log,*)
1255  if( io_l ) write(io_fid_log,'(A)') '|---------------------------------------'
1256  if( io_l ) write(io_fid_log,'(A)') '| size prc_from prc_to'
1257  if( io_l ) write(io_fid_log,'(A10,3(I10))') '| Copy_p2r', copy_info_p2r(:)
1258  do irank = 1, recv_nmax_p2r
1259  if( io_l ) write(io_fid_log,'(A10,3(I10))') '| Recv_p2r', recv_info_p2r(:,irank)
1260  enddo
1261  do irank = 1, send_nmax_p2r
1262  if( io_l ) write(io_fid_log,'(A10,3(I10))') '| Send_p2r', send_info_p2r(:,irank)
1263  enddo
1264 
1265  if ( debug ) then
1266  if( io_l ) write(io_fid_log,*)
1267  if( io_l ) write(io_fid_log,*) '--- Copy_list_p2r'
1268  if( io_l ) write(io_fid_log,*)
1269  if( io_l ) write(io_fid_log,'(11(A6))') 'number', &
1270  '|rfrom','|gfrom','|lfrom','|pfrom', &
1271  '| ito','| jto','| rto','| gto','| lto','| pto'
1272  do ipos = 1, copy_info_p2r(i_size)
1273  g_from = copy_list_p2r(i_grid_from,ipos)
1274  l_from = copy_list_p2r(i_l_from,ipos)
1275  p_from = copy_info_p2r(i_prc_from)
1276  r_from = l_from
1277  g_to = copy_list_p2r(i_grid_to,ipos)
1278  l_to = copy_list_p2r(i_l_to,ipos)
1279  p_to = copy_info_p2r(i_prc_to)
1280  i_to = mod(g_to-1,adm_gall_1d) + 1
1281  j_to = (g_to-i_to) / adm_gall_1d + 1
1282  r_to = prc_rgn_lp2r(l_to,p_to)
1283 
1284  if( io_l ) write(io_fid_log,'(11(I6))') ipos, r_from, g_from, l_from, p_from, &
1285  i_to , j_to , r_to , g_to , l_to , p_to
1286  enddo
1287 
1288  if( io_l ) write(io_fid_log,*)
1289  if( io_l ) write(io_fid_log,*) '--- Recv_list_p2r'
1290  do irank = 1, recv_nmax_p2r
1291  if( io_l ) write(io_fid_log,'(11(A6))') 'number', &
1292  '|rfrom','|gfrom','|lfrom','|pfrom', &
1293  '| ito','| jto','| rto','| gto','| lto','| pto'
1294  do ipos = 1, recv_info_p2r(i_size,irank)
1295  g_from = recv_list_p2r(i_grid_from,ipos,irank)
1296  l_from = recv_list_p2r(i_l_from,ipos,irank)
1297  p_from = recv_info_p2r(i_prc_from,irank)
1298  r_from = l_from
1299  g_to = recv_list_p2r(i_grid_to,ipos,irank)
1300  l_to = recv_list_p2r(i_l_to,ipos,irank)
1301  p_to = recv_info_p2r(i_prc_to,irank)
1302  i_to = mod(g_to-1,adm_gall_1d) + 1
1303  j_to = (g_to-i_to) / adm_gall_1d + 1
1304  r_to = prc_rgn_lp2r(l_to,p_to)
1305 
1306  if( io_l ) write(io_fid_log,'(11(I6))') ipos, r_from, g_from, l_from, p_from, &
1307  i_to , j_to , r_to , g_to , l_to , p_to
1308  enddo
1309  enddo
1310 
1311  if( io_l ) write(io_fid_log,*)
1312  if( io_l ) write(io_fid_log,*) '--- Send_list_p2r'
1313  do irank = 1, send_nmax_p2r
1314  if( io_l ) write(io_fid_log,'(11(A6))') 'number', &
1315  '|rfrom','|gfrom','|lfrom','|pfrom', &
1316  '| ito','| jto','| rto','| gto','| lto','| pto'
1317  do ipos = 1, send_info_p2r(i_size,irank)
1318  g_from = send_list_p2r(i_grid_from,ipos,irank)
1319  l_from = send_list_p2r(i_l_from,ipos,irank)
1320  p_from = send_info_p2r(i_prc_from,irank)
1321  r_from = l_from
1322  g_to = send_list_p2r(i_grid_to,ipos,irank)
1323  l_to = send_list_p2r(i_l_to,ipos,irank)
1324  p_to = send_info_p2r(i_prc_to,irank)
1325  i_to = mod(g_to-1,adm_gall_1d) + 1
1326  j_to = (g_to-i_to) / adm_gall_1d + 1
1327  r_to = prc_rgn_lp2r(l_to,p_to)
1328 
1329  if( io_l ) write(io_fid_log,'(11(I6))') ipos, r_from, g_from, l_from, p_from, &
1330  i_to , j_to , r_to , g_to , l_to , p_to
1331  enddo
1332  enddo
1333  endif
1334 
1335  if( io_l ) write(io_fid_log,*)
1336  if( io_l ) write(io_fid_log,*) '*** Recv_nmax_r2p(local) = ', recv_nmax_r2p
1337  if( io_l ) write(io_fid_log,*) '*** Send_nmax_r2p(local) = ', send_nmax_r2p
1338  if( io_l ) write(io_fid_log,*)
1339  if( io_l ) write(io_fid_log,'(A)') '|---------------------------------------'
1340  if( io_l ) write(io_fid_log,'(A)') '| size prc_from prc_to'
1341  if( io_l ) write(io_fid_log,'(A10,3(I10))') '| Copy_r2p', copy_info_r2p(:)
1342  do irank = 1, recv_nmax_r2p
1343  if( io_l ) write(io_fid_log,'(A10,3(I10))') '| Recv_r2p', recv_info_r2p(:,irank)
1344  enddo
1345  do irank = 1, send_nmax_r2p
1346  if( io_l ) write(io_fid_log,'(A10,3(I10))') '| Send_r2p', send_info_r2p(:,irank)
1347  enddo
1348 
1349  if ( debug ) then
1350  if( io_l ) write(io_fid_log,*)
1351  if( io_l ) write(io_fid_log,*) '--- Copy_list_r2p'
1352  if( io_l ) write(io_fid_log,*)
1353  if( io_l ) write(io_fid_log,'(11(A6))') 'number', &
1354  '|ifrom','|jfrom','|rfrom','|gfrom','|lfrom','|pfrom', &
1355  '| rto','| gto','| lto','| pto'
1356  do ipos = 1, copy_info_r2p(i_size)
1357  g_from = copy_list_r2p(i_grid_from,ipos)
1358  l_from = copy_list_r2p(i_l_from,ipos)
1359  p_from = copy_info_r2p(i_prc_from)
1360  i_from = mod(g_from-1,adm_gall_1d) + 1
1361  j_from = (g_from-i_from) / adm_gall_1d + 1
1362  r_from = prc_rgn_lp2r(l_from,p_from)
1363  g_to = copy_list_r2p(i_grid_to,ipos)
1364  l_to = copy_list_r2p(i_l_to,ipos)
1365  p_to = copy_info_r2p(i_prc_to)
1366  r_to = l_to
1367 
1368  if( io_l ) write(io_fid_log,'(11(I6))') ipos, i_from, j_from, r_from, g_from, l_from, p_from, &
1369  r_to , g_to , l_to , p_to
1370  enddo
1371 
1372  if( io_l ) write(io_fid_log,*)
1373  if( io_l ) write(io_fid_log,*) '--- Recv_list_r2p'
1374  do irank = 1, recv_nmax_r2p
1375  if( io_l ) write(io_fid_log,'(11(A6))') 'number', &
1376  '|ifrom','|jfrom','|rfrom','|gfrom','|lfrom','|pfrom', &
1377  '| rto','| gto','| lto','| pto'
1378  do ipos = 1, recv_info_r2p(i_size,irank)
1379  g_from = recv_list_r2p(i_grid_from,ipos,irank)
1380  l_from = recv_list_r2p(i_l_from,ipos,irank)
1381  p_from = recv_info_r2p(i_prc_from,irank)
1382  i_from = mod(g_from-1,adm_gall_1d) + 1
1383  j_from = (g_from-i_from) / adm_gall_1d + 1
1384  r_from = prc_rgn_lp2r(l_from,p_from)
1385  g_to = recv_list_r2p(i_grid_to,ipos,irank)
1386  l_to = recv_list_r2p(i_l_to,ipos,irank)
1387  p_to = recv_info_r2p(i_prc_to,irank)
1388  r_to = l_to
1389 
1390  if( io_l ) write(io_fid_log,'(11(I6))') ipos, i_from, j_from, r_from, g_from, l_from, p_from, &
1391  r_to , g_to , l_to , p_to
1392  enddo
1393  enddo
1394 
1395  if( io_l ) write(io_fid_log,*)
1396  if( io_l ) write(io_fid_log,*) '--- Send_list_r2p'
1397  do irank = 1, send_nmax_r2p
1398  if( io_l ) write(io_fid_log,'(11(A6))') 'number', &
1399  '|ifrom','|jfrom','|rfrom','|gfrom','|lfrom','|pfrom', &
1400  '| rto','| gto','| lto','| pto'
1401  do ipos = 1, send_info_r2p(i_size,irank)
1402  g_from = send_list_r2p(i_grid_from,ipos,irank)
1403  l_from = send_list_r2p(i_l_from,ipos,irank)
1404  p_from = send_info_r2p(i_prc_from,irank)
1405  i_from = mod(g_from-1,adm_gall_1d) + 1
1406  j_from = (g_from-i_from) / adm_gall_1d + 1
1407  r_from = prc_rgn_lp2r(l_from,p_from)
1408  g_to = send_list_r2p(i_grid_to,ipos,irank)
1409  l_to = send_list_r2p(i_l_to,ipos,irank)
1410  p_to = send_info_r2p(i_prc_to,irank)
1411  r_to = l_to
1412 
1413  if( io_l ) write(io_fid_log,'(11(I6))') ipos, i_from, j_from, r_from, g_from, l_from, p_from, &
1414  r_to , g_to , l_to , p_to
1415  enddo
1416  enddo
1417  endif
1418  if( io_l ) write(io_fid_log,*)
1419  if( io_l ) write(io_fid_log,*) '*** Send_size_p2r,r2p = ', send_size_nglobal_pl
1420 
1421  allocate( sendbuf_r2p_sp(send_size_nglobal_pl*adm_kall*comm_varmax,send_nmax_r2p) )
1422  allocate( recvbuf_r2p_sp(send_size_nglobal_pl*adm_kall*comm_varmax,recv_nmax_r2p) )
1423  allocate( sendbuf_p2r_sp(send_size_nglobal_pl*adm_kall*comm_varmax,send_nmax_p2r) )
1424  allocate( recvbuf_p2r_sp(send_size_nglobal_pl*adm_kall*comm_varmax,recv_nmax_p2r) )
1425  allocate( sendbuf_r2p_dp(send_size_nglobal_pl*adm_kall*comm_varmax,send_nmax_r2p) )
1426  allocate( recvbuf_r2p_dp(send_size_nglobal_pl*adm_kall*comm_varmax,recv_nmax_r2p) )
1427  allocate( sendbuf_p2r_dp(send_size_nglobal_pl*adm_kall*comm_varmax,send_nmax_p2r) )
1428  allocate( recvbuf_p2r_dp(send_size_nglobal_pl*adm_kall*comm_varmax,recv_nmax_p2r) )
1429 
1430  return
1431  end subroutine comm_sortdest_pl
1432 
1433  !-----------------------------------------------------------------------------
1435  subroutine comm_sortdest_singular
1436  use scale_prc, only: &
1437  prc_myrank
1438  use scale_prc_icoa, only: &
1439  prc_rgn_local, &
1440  prc_rgn_vert_num, &
1441  prc_rgn_lp2r, &
1442  prc_rgn_l2r, &
1443  i_w, &
1444  i_n, &
1445  i_s
1446  implicit none
1447 
1448  integer :: rgnid
1449  integer :: i, j, i_rmt, j_rmt
1450 
1451  integer :: l, ipos
1452  integer :: i_from, j_from, r_from, g_from, l_from, p_from
1453  integer :: i_to, j_to, r_to, g_to, l_to, p_to
1454  !---------------------------------------------------------------------------
1455 
1456  allocate( singular_info(info_vindex) )
1457  singular_info(:) = -1
1458  singular_info(i_size) = 0
1459 
1460  allocate( singular_list(list_vindex,4*adm_lall) )
1461  singular_list(:,:) = -1
1462 
1463  do l = 1, prc_rgn_local
1464  rgnid = prc_rgn_l2r(l)
1465 
1466  if ( prc_rgn_vert_num(i_w,rgnid) == 3 ) then
1467  singular_info(i_size) = singular_info(i_size) + 1
1468  ipos = singular_info(i_size)
1469 
1470  i = adm_gmin
1471  j = adm_gmin - 1
1472  i_rmt = adm_gmin - 1
1473  j_rmt = adm_gmin - 1
1474 
1475  singular_list(i_grid_from,ipos) = suf(i,j)
1476  singular_list(i_l_from ,ipos) = l
1477  singular_list(i_grid_to ,ipos) = suf(i_rmt,j_rmt)
1478  singular_list(i_l_to ,ipos) = l
1479  endif
1480 
1481  if ( prc_rgn_vert_num(i_n,rgnid) /= 4 ) then
1482  singular_info(i_size) = singular_info(i_size) + 1
1483  ipos = singular_info(i_size)
1484 
1485  i = adm_gmin
1486  j = adm_gmax + 1
1487  i_rmt = adm_gmin - 1
1488  j_rmt = adm_gmax + 1
1489 
1490  singular_list(i_grid_from,ipos) = suf(i,j)
1491  singular_list(i_l_from ,ipos) = l
1492  singular_list(i_grid_to ,ipos) = suf(i_rmt,j_rmt)
1493  singular_list(i_l_to ,ipos) = l
1494  endif
1495 
1496  if ( prc_rgn_vert_num(i_s,rgnid) /= 4 ) then
1497  singular_info(i_size) = singular_info(i_size) + 1
1498  ipos = singular_info(i_size)
1499 
1500  i = adm_gmax + 1
1501  j = adm_gmin
1502  i_rmt = adm_gmax + 1
1503  j_rmt = adm_gmin - 1
1504 
1505  singular_list(i_grid_from,ipos) = suf(i,j)
1506  singular_list(i_l_from ,ipos) = l
1507  singular_list(i_grid_to ,ipos) = suf(i_rmt,j_rmt)
1508  singular_list(i_l_to ,ipos) = l
1509  endif
1510 
1511  enddo ! loop l
1512 
1513  if ( singular_info(i_size) > 0 ) then
1514  singular_nmax = 1
1515  singular_info(i_prc_from) = prc_myrank
1516  singular_info(i_prc_to ) = prc_myrank
1517  endif
1518 
1519  if( io_l ) write(io_fid_log,*)
1520  if( io_l ) write(io_fid_log,'(A)') '|---------------------------------------'
1521  if( io_l ) write(io_fid_log,'(A)') '| size prc_from prc_to'
1522  if( io_l ) write(io_fid_log,'(A10,3(I10))') '| Singular', singular_info(:)
1523 
1524  if ( debug ) then
1525  if( io_l ) write(io_fid_log,*)
1526  if( io_l ) write(io_fid_log,*) '--- Singular_list'
1527  if( io_l ) write(io_fid_log,*)
1528  if( io_l ) write(io_fid_log,'(13(A6))') 'number', &
1529  '|ifrom','|jfrom','|rfrom','|gfrom','|lfrom','|pfrom', &
1530  '| ito','| jto','| rto','| gto','| lto','| pto'
1531  do ipos = 1, singular_info(i_size)
1532  g_from = singular_list(i_grid_from,ipos)
1533  l_from = singular_list(i_l_from,ipos)
1534  p_from = singular_info(i_prc_from)
1535  i_from = mod(g_from-1,adm_gall_1d) + 1
1536  j_from = (g_from-i_from) / adm_gall_1d + 1
1537  r_from = prc_rgn_lp2r(l_from,p_from)
1538  g_to = singular_list(i_grid_to,ipos)
1539  l_to = singular_list(i_l_to,ipos)
1540  p_to = singular_info(i_prc_to)
1541  i_to = mod(g_to-1,adm_gall_1d) + 1
1542  j_to = (g_to-i_to) / adm_gall_1d + 1
1543  r_to = prc_rgn_lp2r(l_to,p_to)
1544 
1545  if( io_l ) write(io_fid_log,'(13(I6))') ipos, i_from, j_from, r_from, g_from, l_from, p_from, &
1546  i_to , j_to , r_to , g_to , l_to , p_to
1547  enddo
1548  endif
1549 
1550  return
1551  end subroutine comm_sortdest_singular
1552 
1553  !-----------------------------------------------------------------------------
1555  subroutine comm_data_transfer_sp( &
1556  var, &
1557  var_pl )
1558  use scale_prc, only: &
1560  prc_abort
1561  implicit none
1562 
1563  real(SP), intent(inout) :: var (:,:,:,:)
1564  real(SP), intent(inout) :: var_pl(:,:,:,:)
1565 
1566  integer :: shp(4), kmax, vmax
1567  integer :: totalsize, rank, tag
1568  integer :: irank, ipos, imax
1569  integer :: ij_from, l_from, ij_to, l_to
1570 
1571  integer :: k, v, ikv
1572  integer :: ierr
1573  !---------------------------------------------------------------------------
1574  !$acc data &
1575  !$acc pcopy(var) &
1576  !$acc pcopy(var_pl)
1577 
1578  if ( comm_apply_barrier ) then
1579  call prof_rapstart('COMM_barrier',2)
1580  call mpi_barrier(prc_local_comm_world,ierr)
1581  call prof_rapend ('COMM_barrier',2)
1582  endif
1583 
1584  !$acc wait
1585 
1586  call prof_rapstart('COMM_data_transfer',2)
1587 
1588  shp = shape(var)
1589  kmax = shp(2)
1590  vmax = shp(4)
1591 
1592  if ( kmax * vmax > adm_kall * comm_varmax ) then
1593  write(*,*) 'xxx [COMM_data_transfer] kmax * vmax exceeds ADM_kall * COMM_varmax, stop!'
1594  write(*,*) 'xxx kmax * vmax = ', kmax * vmax
1595  write(*,*) 'xxx ADM_kall * COMM_varmax = ', adm_kall * comm_varmax
1596  call prc_abort
1597  endif
1598 
1599  !---< start communication >---
1600  ! Theres no p2r & r2p communication without calling COMM_sortdest_pl.
1601  ! receive pole => region
1602  ! receive region => pole
1603  ! receive region => region
1604  ! pack and send pole => region
1605  ! pack and send region => pole
1606  ! pack and send region => region
1607  ! copy pole => region
1608  ! copy region => pole
1609  ! copy region => region
1610  ! wait all
1611  ! unpack pole => region
1612  ! unpack region => pole
1613  ! unpack region => region
1614  ! copy region halo => region halo (singular point)
1615 
1616  req_count = 0
1617 
1618  !--- receive r2r
1619  do irank = 1, recv_nmax_r2r
1620  req_count = req_count + 1
1621  totalsize = recv_info_r2r(i_size ,irank) * kmax * vmax
1622  rank = recv_info_r2r(i_prc_from,irank)
1623  tag = recv_info_r2r(i_prc_from,irank)
1624 
1625  call mpi_irecv( recvbuf_r2r_sp(1,irank), & ! [OUT]
1626  totalsize, & ! [IN]
1627  mpi_real, & ! [IN]
1628  rank, & ! [IN]
1629  tag, & ! [IN]
1630  prc_local_comm_world, & ! [IN]
1631  req_list(req_count), & ! [OUT]
1632  ierr ) ! [OUT]
1633  enddo
1634 
1635  !--- receive p2r
1636  do irank = 1, recv_nmax_p2r
1637  req_count = req_count + 1
1638  totalsize = recv_info_p2r(i_size ,irank) * kmax * vmax
1639  rank = recv_info_p2r(i_prc_from,irank)
1640  tag = recv_info_p2r(i_prc_from,irank) + 1000000
1641 
1642  call mpi_irecv( recvbuf_p2r_sp(1,irank), & ! [OUT]
1643  totalsize, & ! [IN]
1644  mpi_real, & ! [IN]
1645  rank, & ! [IN]
1646  tag, & ! [IN]
1647  prc_local_comm_world, & ! [IN]
1648  req_list(req_count), & ! [OUT]
1649  ierr ) ! [OUT]
1650  enddo
1651 
1652  !--- receive r2p
1653  do irank = 1, recv_nmax_r2p
1654  req_count = req_count + 1
1655  totalsize = recv_info_r2p(i_size ,irank) * kmax * vmax
1656  rank = recv_info_r2p(i_prc_from,irank)
1657  tag = recv_info_r2p(i_prc_from,irank) + 2000000
1658 
1659  call mpi_irecv( recvbuf_r2p_sp(1,irank), & ! [OUT]
1660  totalsize, & ! [IN]
1661  mpi_real, & ! [IN]
1662  rank, & ! [IN]
1663  tag, & ! [IN]
1664  prc_local_comm_world, & ! [IN]
1665  req_list(req_count), & ! [OUT]
1666  ierr ) ! [OUT]
1667  enddo
1668 
1669  !$acc data &
1670  !$acc pcopyin(Send_list_r2r,Send_list_p2r,Send_list_r2p)
1671  !$acc wait
1672 
1673  !--- pack and send r2r
1674  do irank = 1, send_nmax_r2r
1675  imax = send_info_r2r(i_size,irank)
1676 
1677  !$acc kernels pcopy(sendbuf_r2r_SP)
1678  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
1679  !$omp shared(irank,imax,kmax,vmax,Send_list_r2r,sendbuf_r2r_SP,var) &
1680  !$omp collapse(2)
1681  !$acc loop independent
1682  do v = 1, vmax
1683  !$acc loop independent
1684  do k = 1, kmax
1685  !$acc loop independent
1686  do ipos = 1, imax
1687  ij_from = send_list_r2r(i_grid_from,ipos,irank)
1688  l_from = send_list_r2r(i_l_from ,ipos,irank)
1689 
1690  ikv = (v-1) * imax * kmax &
1691  + (k-1) * imax &
1692  + ipos
1693 
1694  sendbuf_r2r_sp(ikv,irank) = var(ij_from,k,l_from,v)
1695  enddo
1696  enddo
1697  enddo
1698  !$omp end parallel do
1699  !$acc end kernels
1700 
1701  req_count = req_count + 1
1702  totalsize = imax * kmax * vmax
1703  rank = send_info_r2r(i_prc_to ,irank)
1704  tag = send_info_r2r(i_prc_from,irank)
1705 
1706  !$acc wait
1707  call mpi_isend( sendbuf_r2r_sp(1,irank), & ! [IN]
1708  totalsize, & ! [IN]
1709  mpi_real, & ! [IN]
1710  rank, & ! [IN]
1711  tag, & ! [IN]
1712  prc_local_comm_world, & ! [IN]
1713  req_list(req_count), & ! [OUT]
1714  ierr ) ! [OUT]
1715  enddo
1716 
1717  !--- pack and send p2r
1718  do irank = 1, send_nmax_p2r
1719  imax = send_info_p2r(i_size,irank)
1720 
1721  !$acc kernels pcopy(sendbuf_p2r_SP)
1722  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
1723  !$omp shared(irank,imax,kmax,vmax,Send_list_p2r,sendbuf_p2r_SP,var_pl) &
1724  !$omp collapse(2)
1725  !$acc loop independent
1726  do v = 1, vmax
1727  !$acc loop independent
1728  do k = 1, kmax
1729  !$acc loop independent
1730  do ipos = 1, imax
1731  ij_from = send_list_p2r(i_grid_from,ipos,irank)
1732  l_from = send_list_p2r(i_l_from ,ipos,irank)
1733 
1734  ikv = (v-1) * imax * kmax &
1735  + (k-1) * imax &
1736  + ipos
1737 
1738  sendbuf_p2r_sp(ikv,irank) = var_pl(ij_from,k,l_from,v)
1739  enddo
1740  enddo
1741  enddo
1742  !$omp end parallel do
1743  !$acc end kernels
1744 
1745  req_count = req_count + 1
1746  totalsize = imax * kmax * vmax
1747  rank = send_info_p2r(i_prc_to ,irank)
1748  tag = send_info_p2r(i_prc_from,irank) + 1000000
1749 
1750  !$acc wait
1751  call mpi_isend( sendbuf_p2r_sp(1,irank), & ! [IN]
1752  totalsize, & ! [IN]
1753  mpi_real, & ! [IN]
1754  rank, & ! [IN]
1755  tag, & ! [IN]
1756  prc_local_comm_world, & ! [IN]
1757  req_list(req_count), & ! [OUT]
1758  ierr ) ! [OUT]
1759  enddo
1760 
1761  !--- pack and send r2p
1762  do irank = 1, send_nmax_r2p
1763  imax = send_info_r2p(i_size,irank)
1764 
1765  !$acc kernels pcopy(sendbuf_r2p_SP)
1766  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
1767  !$omp shared(irank,imax,kmax,vmax,Send_list_r2p,sendbuf_r2p_SP,var) &
1768  !$omp collapse(2)
1769  !$acc loop independent
1770  do v = 1, vmax
1771  !$acc loop independent
1772  do k = 1, kmax
1773  !$acc loop independent
1774  do ipos = 1, imax
1775  ij_from = send_list_r2p(i_grid_from,ipos,irank)
1776  l_from = send_list_r2p(i_l_from ,ipos,irank)
1777 
1778  ikv = (v-1) * imax * kmax &
1779  + (k-1) * imax &
1780  + ipos
1781 
1782  sendbuf_r2p_sp(ikv,irank) = var(ij_from,k,l_from,v)
1783  enddo
1784  enddo
1785  enddo
1786  !$omp end parallel do
1787  !$acc end kernels
1788 
1789  req_count = req_count + 1
1790  totalsize = imax * kmax * vmax
1791  rank = send_info_r2p(i_prc_to ,irank)
1792  tag = send_info_r2p(i_prc_from,irank) + 2000000
1793 
1794  !$acc wait
1795  call mpi_isend( sendbuf_r2p_sp(1,irank), & ! [IN]
1796  totalsize, & ! [IN]
1797  mpi_real, & ! [IN]
1798  rank, & ! [IN]
1799  tag, & ! [IN]
1800  prc_local_comm_world, & ! [IN]
1801  req_list(req_count), & ! [OUT]
1802  ierr ) ! [OUT]
1803  enddo
1804 
1805  !$acc end data
1806  !$acc data &
1807  !$acc pcopyin(Copy_list_r2r,Copy_list_p2r,Copy_list_r2p)
1808  !$acc wait
1809 
1810  !$omp parallel default(none),private(ipos,k,v,imax,irank,ij_from,l_from,ij_to,l_to) &
1811  !$omp shared(kmax,vmax,var,var_pl, &
1812  !$omp Copy_nmax_p2r,Copy_info_p2r,Copy_list_p2r, &
1813  !$omp Copy_nmax_r2p,Copy_info_r2p,Copy_list_r2p, &
1814  !$omp Copy_nmax_r2r,Copy_info_r2r,Copy_list_r2r )
1815 
1816  !--- copy r2r
1817  do irank = 1, copy_nmax_r2r
1818  imax = copy_info_r2r(i_size)
1819 
1820  !$acc kernels
1821  !$omp do collapse(2)
1822  !$acc loop independent
1823  do v = 1, vmax
1824  !$acc loop independent
1825  do k = 1, kmax
1826  !$acc loop independent
1827  do ipos = 1, imax
1828  ij_from = copy_list_r2r(i_grid_from,ipos)
1829  l_from = copy_list_r2r(i_l_from ,ipos)
1830  ij_to = copy_list_r2r(i_grid_to ,ipos)
1831  l_to = copy_list_r2r(i_l_to ,ipos)
1832 
1833  var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
1834  enddo
1835  enddo
1836  enddo
1837  !$omp end do nowait
1838  !$acc end kernels
1839  enddo
1840 
1841  !--- copy p2r
1842  do irank = 1, copy_nmax_p2r
1843  imax = copy_info_p2r(i_size)
1844 
1845  !$acc kernels
1846  !$omp do collapse(2)
1847  !$acc loop independent
1848  do v = 1, vmax
1849  !$acc loop independent
1850  do k = 1, kmax
1851  !$acc loop independent
1852  do ipos = 1, imax
1853  ij_from = copy_list_p2r(i_grid_from,ipos)
1854  l_from = copy_list_p2r(i_l_from ,ipos)
1855  ij_to = copy_list_p2r(i_grid_to ,ipos)
1856  l_to = copy_list_p2r(i_l_to ,ipos)
1857 
1858  var(ij_to,k,l_to,v) = var_pl(ij_from,k,l_from,v)
1859  enddo
1860  enddo
1861  enddo
1862  !$omp end do nowait
1863  !$acc end kernels
1864  enddo
1865 
1866  !--- copy r2p
1867  do irank = 1, copy_nmax_r2p
1868  imax = copy_info_r2p(i_size)
1869 
1870  !$acc kernels
1871  !$omp do collapse(2)
1872  !$acc loop independent
1873  do v = 1, vmax
1874  !$acc loop independent
1875  do k = 1, kmax
1876  !$acc loop independent
1877  do ipos = 1, imax
1878  ij_from = copy_list_r2p(i_grid_from,ipos)
1879  l_from = copy_list_r2p(i_l_from ,ipos)
1880  ij_to = copy_list_r2p(i_grid_to ,ipos)
1881  l_to = copy_list_r2p(i_l_to ,ipos)
1882 
1883  var_pl(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
1884  enddo
1885  enddo
1886  enddo
1887  !$omp end do nowait
1888  !$acc end kernels
1889  enddo
1890 
1891  !$omp end parallel
1892 
1893  !--- wait all
1894  if ( req_count > 0 ) then
1895  call mpi_waitall( req_count, & ! [IN]
1896  req_list(:), & ! [IN]
1897  mpi_statuses_ignore, & ! [OUT]
1898  ierr ) ! [OUT]
1899  endif
1900 
1901  !$acc end data
1902  !$acc data &
1903  !$acc pcopyin(Recv_list_r2r,Recv_list_p2r,Recv_list_r2p)
1904  !$acc wait
1905 
1906  !$omp parallel default(none),private(ipos,k,v,imax,irank,ikv,ij_from,l_from,ij_to,l_to) &
1907  !$omp shared(kmax,vmax,var,var_pl, &
1908  !$omp Recv_nmax_p2r,Recv_info_p2r,Recv_list_p2r,recvbuf_p2r_SP, &
1909  !$omp Recv_nmax_r2p,Recv_info_r2p,Recv_list_r2p,recvbuf_r2p_SP, &
1910  !$omp Recv_nmax_r2r,Recv_info_r2r,Recv_list_r2r,recvbuf_r2r_SP, &
1911  !$omp Singular_nmax,Singular_info,Singular_list )
1912 
1913  !--- unpack r2r
1914  do irank = 1, recv_nmax_r2r
1915  imax = recv_info_r2r(i_size,irank)
1916 
1917  !$acc kernels pcopyin(recvbuf_r2r_SP)
1918  !$omp do collapse(2)
1919  !$acc loop independent
1920  do v = 1, vmax
1921  !$acc loop independent
1922  do k = 1, kmax
1923  !$acc loop independent
1924  do ipos = 1, imax
1925  ij_to = recv_list_r2r(i_grid_to,ipos,irank)
1926  l_to = recv_list_r2r(i_l_to ,ipos,irank)
1927 
1928  ikv = (v-1) * imax * kmax &
1929  + (k-1) * imax &
1930  + ipos
1931 
1932  var(ij_to,k,l_to,v) = recvbuf_r2r_sp(ikv,irank)
1933  enddo
1934  enddo
1935  enddo
1936  !$omp end do nowait
1937  !$acc end kernels
1938  enddo
1939 
1940  !--- unpack p2r
1941  do irank = 1, recv_nmax_p2r
1942  imax = recv_info_p2r(i_size,irank)
1943 
1944  !$acc kernels pcopyin(recvbuf_p2r_SP)
1945  !$omp do collapse(2)
1946  !$acc loop independent
1947  do v = 1, vmax
1948  !$acc loop independent
1949  do k = 1, kmax
1950  !$acc loop independent
1951  do ipos = 1, imax
1952  ij_to = recv_list_p2r(i_grid_to,ipos,irank)
1953  l_to = recv_list_p2r(i_l_to ,ipos,irank)
1954 
1955  ikv = (v-1) * imax * kmax &
1956  + (k-1) * imax &
1957  + ipos
1958 
1959  var(ij_to,k,l_to,v) = recvbuf_p2r_sp(ikv,irank)
1960  enddo
1961  enddo
1962  enddo
1963  !$omp end do nowait
1964  !$acc end kernels
1965  enddo
1966 
1967  !--- unpack r2p
1968  do irank = 1, recv_nmax_r2p
1969  imax = recv_info_r2p(i_size,irank)
1970 
1971  !$acc kernels pcopyin(recvbuf_r2p_SP)
1972  !$omp do collapse(2)
1973  !$acc loop independent
1974  do v = 1, vmax
1975  !$acc loop independent
1976  do k = 1, kmax
1977  !$acc loop independent
1978  do ipos = 1, imax
1979  ij_to = recv_list_r2p(i_grid_to,ipos,irank)
1980  l_to = recv_list_r2p(i_l_to ,ipos,irank)
1981 
1982  ikv = (v-1) * imax * kmax &
1983  + (k-1) * imax &
1984  + ipos
1985 
1986  var_pl(ij_to,k,l_to,v) = recvbuf_r2p_sp(ikv,irank)
1987  enddo
1988  enddo
1989  enddo
1990  !$omp end do nowait
1991  !$acc end kernels
1992  enddo
1993 
1994  !$acc end data
1995  !$acc wait
1996 
1997  !--- singular point (halo to halo)
1998  do irank = 1, singular_nmax
1999  imax = singular_info(i_size)
2000 
2001  !$acc kernels pcopyin(Singular_list)
2002  !$omp do collapse(2)
2003  !$acc loop independent
2004  do v = 1, vmax
2005  !$acc loop independent
2006  do k = 1, kmax
2007  !$acc loop independent
2008  do ipos = 1, imax
2009  ij_from = singular_list(i_grid_from,ipos)
2010  l_from = singular_list(i_l_from ,ipos)
2011  ij_to = singular_list(i_grid_to ,ipos)
2012  l_to = singular_list(i_l_to ,ipos)
2013 
2014  var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2015  enddo
2016  enddo
2017  enddo
2018  !$omp end do nowait
2019  !$acc end kernels
2020  enddo
2021 
2022  !$omp end parallel
2023 
2024  call prof_rapend('COMM_data_transfer',2)
2025 
2026  !$acc end data
2027  !$acc wait
2028 
2029  return
2030  end subroutine comm_data_transfer_sp
2031 
2032  !-----------------------------------------------------------------------------
2034  subroutine comm_data_transfer_dp( &
2035  var, &
2036  var_pl )
2037  use scale_prc, only: &
2039  prc_abort
2040  implicit none
2041 
2042  real(DP), intent(inout) :: var (:,:,:,:)
2043  real(DP), intent(inout) :: var_pl(:,:,:,:)
2044 
2045  integer :: shp(4), kmax, vmax
2046  integer :: totalsize, rank, tag
2047  integer :: irank, ipos, imax
2048  integer :: ij_from, l_from, ij_to, l_to
2049 
2050  integer :: k, v, ikv
2051  integer :: ierr
2052  !---------------------------------------------------------------------------
2053  !$acc data &
2054  !$acc pcopy(var) &
2055  !$acc pcopy(var_pl)
2056 
2057  if ( comm_apply_barrier ) then
2058  call prof_rapstart('COMM_barrier',2)
2059  call mpi_barrier(prc_local_comm_world,ierr)
2060  call prof_rapend ('COMM_barrier',2)
2061  endif
2062 
2063  !$acc wait
2064 
2065  call prof_rapstart('COMM_data_transfer',2)
2066 
2067  shp = shape(var)
2068  kmax = shp(2)
2069  vmax = shp(4)
2070 
2071  if ( kmax * vmax > adm_kall * comm_varmax ) then
2072  write(*,*) 'xxx [COMM_data_transfer] kmax * vmax exceeds ADM_kall * COMM_varmax, stop!'
2073  write(*,*) 'xxx kmax * vmax = ', kmax * vmax
2074  write(*,*) 'xxx ADM_kall * COMM_varmax = ', adm_kall * comm_varmax
2075  call prc_abort
2076  endif
2077 
2078  !---< start communication >---
2079  ! Theres no p2r & r2p communication without calling COMM_sortdest_pl.
2080  ! receive pole => region
2081  ! receive region => pole
2082  ! receive region => region
2083  ! pack and send pole => region
2084  ! pack and send region => pole
2085  ! pack and send region => region
2086  ! copy pole => region
2087  ! copy region => pole
2088  ! copy region => region
2089  ! wait all
2090  ! unpack pole => region
2091  ! unpack region => pole
2092  ! unpack region => region
2093  ! copy region halo => region halo (singular point)
2094 
2095  req_count = 0
2096 
2097  !--- receive r2r
2098  do irank = 1, recv_nmax_r2r
2099  req_count = req_count + 1
2100  totalsize = recv_info_r2r(i_size ,irank) * kmax * vmax
2101  rank = recv_info_r2r(i_prc_from,irank)
2102  tag = recv_info_r2r(i_prc_from,irank)
2103 
2104  call mpi_irecv( recvbuf_r2r_dp(1,irank), & ! [OUT]
2105  totalsize, & ! [IN]
2106  mpi_double_precision, & ! [IN]
2107  rank, & ! [IN]
2108  tag, & ! [IN]
2109  prc_local_comm_world, & ! [IN]
2110  req_list(req_count), & ! [OUT]
2111  ierr ) ! [OUT]
2112  enddo
2113 
2114  !--- receive p2r
2115  do irank = 1, recv_nmax_p2r
2116  req_count = req_count + 1
2117  totalsize = recv_info_p2r(i_size ,irank) * kmax * vmax
2118  rank = recv_info_p2r(i_prc_from,irank)
2119  tag = recv_info_p2r(i_prc_from,irank) + 1000000
2120 
2121  call mpi_irecv( recvbuf_p2r_dp(1,irank), & ! [OUT]
2122  totalsize, & ! [IN]
2123  mpi_double_precision, & ! [IN]
2124  rank, & ! [IN]
2125  tag, & ! [IN]
2126  prc_local_comm_world, & ! [IN]
2127  req_list(req_count), & ! [OUT]
2128  ierr ) ! [OUT]
2129  enddo
2130 
2131  !--- receive r2p
2132  do irank = 1, recv_nmax_r2p
2133  req_count = req_count + 1
2134  totalsize = recv_info_r2p(i_size ,irank) * kmax * vmax
2135  rank = recv_info_r2p(i_prc_from,irank)
2136  tag = recv_info_r2p(i_prc_from,irank) + 2000000
2137 
2138  call mpi_irecv( recvbuf_r2p_dp(1,irank), & ! [OUT]
2139  totalsize, & ! [IN]
2140  mpi_double_precision, & ! [IN]
2141  rank, & ! [IN]
2142  tag, & ! [IN]
2143  prc_local_comm_world, & ! [IN]
2144  req_list(req_count), & ! [OUT]
2145  ierr ) ! [OUT]
2146  enddo
2147 
2148  !$acc data &
2149  !$acc pcopyin(Send_list_r2r,Send_list_p2r,Send_list_r2p)
2150  !$acc wait
2151 
2152  !--- pack and send r2r
2153  do irank = 1, send_nmax_r2r
2154  imax = send_info_r2r(i_size,irank)
2155 
2156  !$acc kernels pcopy(sendbuf_r2r_DP)
2157  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
2158  !$omp shared(irank,imax,kmax,vmax,Send_list_r2r,sendbuf_r2r_DP,var) &
2159  !$omp collapse(2)
2160  !$acc loop independent
2161  do v = 1, vmax
2162  !$acc loop independent
2163  do k = 1, kmax
2164  !$acc loop independent
2165  do ipos = 1, imax
2166  ij_from = send_list_r2r(i_grid_from,ipos,irank)
2167  l_from = send_list_r2r(i_l_from ,ipos,irank)
2168 
2169  ikv = (v-1) * imax * kmax &
2170  + (k-1) * imax &
2171  + ipos
2172 
2173  sendbuf_r2r_dp(ikv,irank) = var(ij_from,k,l_from,v)
2174  enddo
2175  enddo
2176  enddo
2177  !$omp end parallel do
2178  !$acc end kernels
2179 
2180  req_count = req_count + 1
2181  totalsize = imax * kmax * vmax
2182  rank = send_info_r2r(i_prc_to ,irank)
2183  tag = send_info_r2r(i_prc_from,irank)
2184 
2185  !$acc wait
2186  call mpi_isend( sendbuf_r2r_dp(1,irank), & ! [IN]
2187  totalsize, & ! [IN]
2188  mpi_double_precision, & ! [IN]
2189  rank, & ! [IN]
2190  tag, & ! [IN]
2191  prc_local_comm_world, & ! [IN]
2192  req_list(req_count), & ! [OUT]
2193  ierr ) ! [OUT]
2194  enddo
2195 
2196  !--- pack and send p2r
2197  do irank = 1, send_nmax_p2r
2198  imax = send_info_p2r(i_size,irank)
2199 
2200  !$acc kernels pcopy(sendbuf_p2r_DP)
2201  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
2202  !$omp shared(irank,imax,kmax,vmax,Send_list_p2r,sendbuf_p2r_DP,var_pl) &
2203  !$omp collapse(2)
2204  !$acc loop independent
2205  do v = 1, vmax
2206  !$acc loop independent
2207  do k = 1, kmax
2208  !$acc loop independent
2209  do ipos = 1, imax
2210  ij_from = send_list_p2r(i_grid_from,ipos,irank)
2211  l_from = send_list_p2r(i_l_from ,ipos,irank)
2212 
2213  ikv = (v-1) * imax * kmax &
2214  + (k-1) * imax &
2215  + ipos
2216 
2217  sendbuf_p2r_dp(ikv,irank) = var_pl(ij_from,k,l_from,v)
2218  enddo
2219  enddo
2220  enddo
2221  !$omp end parallel do
2222  !$acc end kernels
2223 
2224  req_count = req_count + 1
2225  totalsize = imax * kmax * vmax
2226  rank = send_info_p2r(i_prc_to ,irank)
2227  tag = send_info_p2r(i_prc_from,irank) + 1000000
2228 
2229  !$acc wait
2230  call mpi_isend( sendbuf_p2r_dp(1,irank), & ! [IN]
2231  totalsize, & ! [IN]
2232  mpi_double_precision, & ! [IN]
2233  rank, & ! [IN]
2234  tag, & ! [IN]
2235  prc_local_comm_world, & ! [IN]
2236  req_list(req_count), & ! [OUT]
2237  ierr ) ! [OUT]
2238  enddo
2239 
2240  !--- pack and send r2p
2241  do irank = 1, send_nmax_r2p
2242  imax = send_info_r2p(i_size,irank)
2243 
2244  !$acc kernels pcopy(sendbuf_r2p_DP)
2245  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
2246  !$omp shared(irank,imax,kmax,vmax,Send_list_r2p,sendbuf_r2p_DP,var) &
2247  !$omp collapse(2)
2248  !$acc loop independent
2249  do v = 1, vmax
2250  !$acc loop independent
2251  do k = 1, kmax
2252  !$acc loop independent
2253  do ipos = 1, imax
2254  ij_from = send_list_r2p(i_grid_from,ipos,irank)
2255  l_from = send_list_r2p(i_l_from ,ipos,irank)
2256 
2257  ikv = (v-1) * imax * kmax &
2258  + (k-1) * imax &
2259  + ipos
2260 
2261  sendbuf_r2p_dp(ikv,irank) = var(ij_from,k,l_from,v)
2262  enddo
2263  enddo
2264  enddo
2265  !$omp end parallel do
2266  !$acc end kernels
2267 
2268  req_count = req_count + 1
2269  totalsize = imax * kmax * vmax
2270  rank = send_info_r2p(i_prc_to ,irank)
2271  tag = send_info_r2p(i_prc_from,irank) + 2000000
2272 
2273  !$acc wait
2274  call mpi_isend( sendbuf_r2p_dp(1,irank), & ! [IN]
2275  totalsize, & ! [IN]
2276  mpi_double_precision, & ! [IN]
2277  rank, & ! [IN]
2278  tag, & ! [IN]
2279  prc_local_comm_world, & ! [IN]
2280  req_list(req_count), & ! [OUT]
2281  ierr ) ! [OUT]
2282  enddo
2283 
2284  !$acc end data
2285  !$acc data &
2286  !$acc pcopyin(Copy_list_r2r,Copy_list_p2r,Copy_list_r2p)
2287  !$acc wait
2288 
2289  !$omp parallel default(none),private(ipos,k,v,imax,irank,ij_from,l_from,ij_to,l_to) &
2290  !$omp shared(kmax,vmax,var,var_pl, &
2291  !$omp Copy_nmax_p2r,Copy_info_p2r,Copy_list_p2r, &
2292  !$omp Copy_nmax_r2p,Copy_info_r2p,Copy_list_r2p, &
2293  !$omp Copy_nmax_r2r,Copy_info_r2r,Copy_list_r2r )
2294 
2295  !--- copy r2r
2296  do irank = 1, copy_nmax_r2r
2297  imax = copy_info_r2r(i_size)
2298 
2299  !$acc kernels
2300  !$omp do collapse(2)
2301  !$acc loop independent
2302  do v = 1, vmax
2303  !$acc loop independent
2304  do k = 1, kmax
2305  !$acc loop independent
2306  do ipos = 1, imax
2307  ij_from = copy_list_r2r(i_grid_from,ipos)
2308  l_from = copy_list_r2r(i_l_from ,ipos)
2309  ij_to = copy_list_r2r(i_grid_to ,ipos)
2310  l_to = copy_list_r2r(i_l_to ,ipos)
2311 
2312  var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2313  enddo
2314  enddo
2315  enddo
2316  !$omp end do nowait
2317  !$acc end kernels
2318  enddo
2319 
2320  !--- copy p2r
2321  do irank = 1, copy_nmax_p2r
2322  imax = copy_info_p2r(i_size)
2323 
2324  !$acc kernels
2325  !$omp do collapse(2)
2326  !$acc loop independent
2327  do v = 1, vmax
2328  !$acc loop independent
2329  do k = 1, kmax
2330  !$acc loop independent
2331  do ipos = 1, imax
2332  ij_from = copy_list_p2r(i_grid_from,ipos)
2333  l_from = copy_list_p2r(i_l_from ,ipos)
2334  ij_to = copy_list_p2r(i_grid_to ,ipos)
2335  l_to = copy_list_p2r(i_l_to ,ipos)
2336 
2337  var(ij_to,k,l_to,v) = var_pl(ij_from,k,l_from,v)
2338  enddo
2339  enddo
2340  enddo
2341  !$omp end do nowait
2342  !$acc end kernels
2343  enddo
2344 
2345  !--- copy r2p
2346  do irank = 1, copy_nmax_r2p
2347  imax = copy_info_r2p(i_size)
2348 
2349  !$acc kernels
2350  !$omp do collapse(2)
2351  !$acc loop independent
2352  do v = 1, vmax
2353  !$acc loop independent
2354  do k = 1, kmax
2355  !$acc loop independent
2356  do ipos = 1, imax
2357  ij_from = copy_list_r2p(i_grid_from,ipos)
2358  l_from = copy_list_r2p(i_l_from ,ipos)
2359  ij_to = copy_list_r2p(i_grid_to ,ipos)
2360  l_to = copy_list_r2p(i_l_to ,ipos)
2361 
2362  var_pl(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2363  enddo
2364  enddo
2365  enddo
2366  !$omp end do nowait
2367  !$acc end kernels
2368  enddo
2369 
2370  !$omp end parallel
2371 
2372  !--- wait all
2373  if ( req_count > 0 ) then
2374  call mpi_waitall( req_count, & ! [IN]
2375  req_list(:), & ! [IN]
2376  mpi_statuses_ignore, & ! [OUT]
2377  ierr ) ! [OUT]
2378  endif
2379 
2380  !$acc end data
2381  !$acc data &
2382  !$acc pcopyin(Recv_list_r2r,Recv_list_p2r,Recv_list_r2p)
2383  !$acc wait
2384 
2385  !$omp parallel default(none),private(ipos,k,v,imax,irank,ikv,ij_from,l_from,ij_to,l_to) &
2386  !$omp shared(kmax,vmax,var,var_pl, &
2387  !$omp Recv_nmax_p2r,Recv_info_p2r,Recv_list_p2r,recvbuf_p2r_DP, &
2388  !$omp Recv_nmax_r2p,Recv_info_r2p,Recv_list_r2p,recvbuf_r2p_DP, &
2389  !$omp Recv_nmax_r2r,Recv_info_r2r,Recv_list_r2r,recvbuf_r2r_DP, &
2390  !$omp Singular_nmax,Singular_info,Singular_list )
2391 
2392  !--- unpack r2r
2393  do irank = 1, recv_nmax_r2r
2394  imax = recv_info_r2r(i_size,irank)
2395 
2396  !$acc kernels pcopyin(recvbuf_r2r_DP)
2397  !$omp do collapse(2)
2398  !$acc loop independent
2399  do v = 1, vmax
2400  !$acc loop independent
2401  do k = 1, kmax
2402  !$acc loop independent
2403  do ipos = 1, imax
2404  ij_to = recv_list_r2r(i_grid_to,ipos,irank)
2405  l_to = recv_list_r2r(i_l_to ,ipos,irank)
2406 
2407  ikv = (v-1) * imax * kmax &
2408  + (k-1) * imax &
2409  + ipos
2410 
2411  var(ij_to,k,l_to,v) = recvbuf_r2r_dp(ikv,irank)
2412  enddo
2413  enddo
2414  enddo
2415  !$omp end do nowait
2416  !$acc end kernels
2417  enddo
2418 
2419  !--- unpack p2r
2420  do irank = 1, recv_nmax_p2r
2421  imax = recv_info_p2r(i_size,irank)
2422 
2423  !$acc kernels pcopyin(recvbuf_p2r_DP)
2424  !$omp do collapse(2)
2425  !$acc loop independent
2426  do v = 1, vmax
2427  !$acc loop independent
2428  do k = 1, kmax
2429  !$acc loop independent
2430  do ipos = 1, imax
2431  ij_to = recv_list_p2r(i_grid_to,ipos,irank)
2432  l_to = recv_list_p2r(i_l_to ,ipos,irank)
2433 
2434  ikv = (v-1) * imax * kmax &
2435  + (k-1) * imax &
2436  + ipos
2437 
2438  var(ij_to,k,l_to,v) = recvbuf_p2r_dp(ikv,irank)
2439  enddo
2440  enddo
2441  enddo
2442  !$omp end do nowait
2443  !$acc end kernels
2444  enddo
2445 
2446  !--- unpack r2p
2447  do irank = 1, recv_nmax_r2p
2448  imax = recv_info_r2p(i_size,irank)
2449 
2450  !$acc kernels pcopyin(recvbuf_r2p_DP)
2451  !$omp do collapse(2)
2452  !$acc loop independent
2453  do v = 1, vmax
2454  !$acc loop independent
2455  do k = 1, kmax
2456  !$acc loop independent
2457  do ipos = 1, imax
2458  ij_to = recv_list_r2p(i_grid_to,ipos,irank)
2459  l_to = recv_list_r2p(i_l_to ,ipos,irank)
2460 
2461  ikv = (v-1) * imax * kmax &
2462  + (k-1) * imax &
2463  + ipos
2464 
2465  var_pl(ij_to,k,l_to,v) = recvbuf_r2p_dp(ikv,irank)
2466  enddo
2467  enddo
2468  enddo
2469  !$omp end do nowait
2470  !$acc end kernels
2471  enddo
2472 
2473  !$acc end data
2474  !$acc wait
2475 
2476  !--- singular point (halo to halo)
2477  do irank = 1, singular_nmax
2478  imax = singular_info(i_size)
2479 
2480  !$acc kernels pcopyin(Singular_list)
2481  !$omp do collapse(2)
2482  !$acc loop independent
2483  do v = 1, vmax
2484  !$acc loop independent
2485  do k = 1, kmax
2486  !$acc loop independent
2487  do ipos = 1, imax
2488  ij_from = singular_list(i_grid_from,ipos)
2489  l_from = singular_list(i_l_from ,ipos)
2490  ij_to = singular_list(i_grid_to ,ipos)
2491  l_to = singular_list(i_l_to ,ipos)
2492 
2493  var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2494  enddo
2495  enddo
2496  enddo
2497  !$omp end do nowait
2498  !$acc end kernels
2499  enddo
2500 
2501  !$omp end parallel
2502 
2503  call prof_rapend('COMM_data_transfer',2)
2504 
2505  !$acc end data
2506  !$acc wait
2507 
2508  return
2509  end subroutine comm_data_transfer_dp
2510 
2511  !-----------------------------------------------------------------------------
2513  subroutine comm_data_transfer_nopl( &
2514  var )
2515  use scale_prc, only: &
2517  prc_abort
2518  implicit none
2519 
2520  real(DP), intent(inout) :: var (:,:,:,:)
2521 
2522  integer :: shp(4), kmax, vmax
2523  integer :: totalsize, rank, tag
2524  integer :: irank, ipos, imax
2525  integer :: ij_from, l_from, ij_to, l_to
2526 
2527  integer :: k, v, ikv
2528  integer :: ierr
2529  !---------------------------------------------------------------------------
2530 
2531  if ( comm_apply_barrier ) then
2532  call prof_rapstart('COMM_barrier',2)
2533  call mpi_barrier(prc_local_comm_world,ierr)
2534  call prof_rapend ('COMM_barrier',2)
2535  endif
2536 
2537  call prof_rapstart('COMM_data_transfer',2)
2538 
2539  shp = shape(var)
2540  kmax = shp(2)
2541  vmax = shp(4)
2542 
2543  if ( kmax * vmax > adm_kall * comm_varmax ) then
2544  write(*,*) 'xxx [COMM_data_transfer_nopl] kmax * vmax exceeds ADM_kall * COMM_varmax, stop!'
2545  write(*,*) 'xxx kmax * vmax = ', kmax * vmax
2546  write(*,*) 'xxx ADM_kall * COMM_varmax = ', adm_kall * comm_varmax
2547  call prc_abort
2548  endif
2549 
2550  !---< start communication >---
2551  ! receive region => region
2552  ! pack and send region => region
2553  ! copy region => region
2554  ! wait all
2555  ! unpack region => region
2556  ! copy region halo => region halo (singular point)
2557 
2558  req_count = 0
2559 
2560  !--- receive r2r
2561  do irank = 1, recv_nmax_r2r
2562  req_count = req_count + 1
2563  totalsize = recv_info_r2r(i_size ,irank) * kmax * vmax
2564  rank = recv_info_r2r(i_prc_from,irank)
2565  tag = recv_info_r2r(i_prc_from,irank)
2566 
2567  call mpi_irecv( recvbuf_r2r_dp(1,irank), & ! [OUT]
2568  totalsize, & ! [IN]
2569  mpi_double_precision, & ! [IN]
2570  rank, & ! [IN]
2571  tag, & ! [IN]
2572  prc_local_comm_world, & ! [IN]
2573  req_list(req_count), & ! [OUT]
2574  ierr ) ! [OUT]
2575  enddo
2576 
2577  !--- pack and send r2r
2578  do irank = 1, send_nmax_r2r
2579  imax = send_info_r2r(i_size,irank)
2580 
2581  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
2582  !$omp shared(irank,imax,kmax,vmax,Send_list_r2r,sendbuf_r2r_DP,var)
2583  do v = 1, vmax
2584  do k = 1, kmax
2585  do ipos = 1, imax
2586  ij_from = send_list_r2r(i_grid_from,ipos,irank)
2587  l_from = send_list_r2r(i_l_from ,ipos,irank)
2588 
2589  ikv = (v-1) * imax * kmax &
2590  + (k-1) * imax &
2591  + ipos
2592 
2593  sendbuf_r2r_dp(ikv,irank) = var(ij_from,k,l_from,v)
2594  enddo
2595  enddo
2596  enddo
2597  !$omp end parallel do
2598 
2599  req_count = req_count + 1
2600  totalsize = imax * kmax * vmax
2601  rank = send_info_r2r(i_prc_to ,irank)
2602  tag = send_info_r2r(i_prc_from,irank)
2603 
2604  call mpi_isend( sendbuf_r2r_dp(1,irank), & ! [IN]
2605  totalsize, & ! [IN]
2606  mpi_double_precision, & ! [IN]
2607  rank, & ! [IN]
2608  tag, & ! [IN]
2609  prc_local_comm_world, & ! [IN]
2610  req_list(req_count), & ! [OUT]
2611  ierr ) ! [OUT]
2612  enddo
2613 
2614  !--- copy r2r
2615  !$omp parallel do default(none),private(ipos,k,v,imax,irank,ij_from,l_from,ij_to,l_to) &
2616  !$omp shared(kmax,vmax,Copy_nmax_r2r,Copy_info_r2r,Copy_list_r2r,var)
2617  do irank = 1, copy_nmax_r2r
2618  imax = copy_info_r2r(i_size)
2619 
2620  do v = 1, vmax
2621  do k = 1, kmax
2622  do ipos = 1, imax
2623  ij_from = copy_list_r2r(i_grid_from,ipos)
2624  l_from = copy_list_r2r(i_l_from ,ipos)
2625  ij_to = copy_list_r2r(i_grid_to ,ipos)
2626  l_to = copy_list_r2r(i_l_to ,ipos)
2627 
2628  var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2629  enddo
2630  enddo
2631  enddo
2632  enddo
2633  !$omp end parallel do
2634 
2635  !--- wait all
2636  if ( req_count > 0 ) then
2637  call mpi_waitall( req_count, & ! [IN]
2638  req_list(:), & ! [IN]
2639  mpi_statuses_ignore, & ! [OUT]
2640  ierr ) ! [OUT]
2641  endif
2642 
2643  !--- unpack r2r
2644  !$omp parallel do default(none),private(ipos,k,v,imax,irank,ij_to,l_to,ikv) &
2645  !$omp shared(kmax,vmax,Recv_nmax_r2r,Recv_info_r2r,Recv_list_r2r,recvbuf_r2r_DP,var)
2646  do irank = 1, recv_nmax_r2r
2647  imax = recv_info_r2r(i_size,irank)
2648 
2649  do v = 1, vmax
2650  do k = 1, kmax
2651  do ipos = 1, imax
2652  ij_to = recv_list_r2r(i_grid_to,ipos,irank)
2653  l_to = recv_list_r2r(i_l_to ,ipos,irank)
2654 
2655  ikv = (v-1) * imax * kmax &
2656  + (k-1) * imax &
2657  + ipos
2658 
2659  var(ij_to,k,l_to,v) = recvbuf_r2r_dp(ikv,irank)
2660  enddo
2661  enddo
2662  enddo
2663  enddo
2664  !$omp end parallel do
2665 
2666  !--- singular point (halo to halo)
2667  do irank = 1, singular_nmax
2668  imax = singular_info(i_size)
2669 
2670  do v = 1, vmax
2671  do k = 1, kmax
2672  do ipos = 1, imax
2673  ij_from = singular_list(i_grid_from,ipos)
2674  l_from = singular_list(i_l_from ,ipos)
2675  ij_to = singular_list(i_grid_to ,ipos)
2676  l_to = singular_list(i_l_to ,ipos)
2677 
2678  var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2679  enddo
2680  enddo
2681  enddo
2682  enddo
2683 
2684  call prof_rapend('COMM_data_transfer',2)
2685 
2686  return
2687  end subroutine comm_data_transfer_nopl
2688 
2689  !-----------------------------------------------------------------------------
2691  subroutine comm_var_sp( &
2692  var, &
2693  var_pl, &
2694  kmax, &
2695  vmax )
2696  use scale_prc, only: &
2698  use scale_prc_icoa, only: &
2699  prc_rgn_total_pl, &
2700  prc_rgn_rgn4pl, &
2701  prc_rgn_lp2r, &
2702  i_npl, &
2703  i_spl
2704  implicit none
2705 
2706  integer, intent(in) :: kmax
2707  integer, intent(in) :: vmax
2708  real(SP), intent(inout) :: var (adm_gall, kmax,adm_lall, vmax)
2709  real(SP), intent(inout) :: var_pl(adm_gall_pl,kmax,adm_lall_pl,vmax)
2710 
2711  real(SP) :: sendbuf_h2p_SP(kmax*vmax,prc_rgn_total_pl)
2712  real(SP) :: recvbuf_h2p_SP(kmax*vmax,prc_rgn_total_pl)
2713 
2714  integer :: totalsize, rank, tag
2715  integer :: irank, ipos
2716  integer :: ij_from, l_from, ij_to, l_to
2717  integer :: r_from, r_to
2718 
2719  integer :: k, v, kk
2720  integer :: ierr
2721  !---------------------------------------------------------------------------
2722 
2723  if ( comm_apply_barrier ) then
2724  call prof_rapstart('COMM_barrier',2)
2725  call mpi_barrier(prc_local_comm_world,ierr)
2726  call prof_rapend ('COMM_barrier',2)
2727  endif
2728 
2729  call prof_rapstart('COMM_var',2)
2730 
2731  if ( comm_pl ) then
2732 
2733  req_count = 0
2734 
2735  !--- receive p2r-reverse
2736  do irank = 1, send_nmax_p2r
2737  do ipos = 1, send_info_p2r(i_size,irank)
2738  l_from = send_list_p2r(i_l_to ,ipos,irank)
2739  r_from = prc_rgn_lp2r(l_from,send_info_p2r(i_prc_to,irank))
2740 
2741  if ( r_from == prc_rgn_rgn4pl(i_npl) ) then
2742  req_count = req_count + 1
2743  totalsize = kmax * vmax
2744  rank = send_info_p2r(i_prc_to ,irank)
2745  tag = send_info_p2r(i_prc_from,irank) + 1000000
2746 
2747  call mpi_irecv( recvbuf_h2p_sp(1,i_npl), & ! [OUT]
2748  totalsize, & ! [IN]
2749  mpi_real, & ! [IN]
2750  rank, & ! [IN]
2751  tag, & ! [IN]
2752  prc_local_comm_world, & ! [IN]
2753  req_list(req_count), & ! [OUT]
2754  ierr ) ! [OUT]
2755  endif
2756 
2757  if ( r_from == prc_rgn_rgn4pl(i_spl) ) then
2758  req_count = req_count + 1
2759  totalsize = kmax * vmax
2760  rank = send_info_p2r(i_prc_to ,irank)
2761  tag = send_info_p2r(i_prc_from,irank) + 2000000
2762 
2763  call mpi_irecv( recvbuf_h2p_sp(1,i_spl), & ! [OUT]
2764  totalsize, & ! [IN]
2765  mpi_real, & ! [IN]
2766  rank, & ! [IN]
2767  tag, & ! [IN]
2768  prc_local_comm_world, & ! [IN]
2769  req_list(req_count), & ! [OUT]
2770  ierr ) ! [OUT]
2771  endif
2772  enddo
2773  enddo
2774 
2775  !--- pack and send p2r-reverse
2776  do irank = 1, recv_nmax_p2r
2777  do ipos = 1, recv_info_p2r(i_size,irank)
2778  ij_from = recv_list_p2r(i_grid_to,ipos,irank)
2779  l_from = recv_list_p2r(i_l_to ,ipos,irank)
2780  r_from = prc_rgn_lp2r(l_from,recv_info_p2r(i_prc_to,irank))
2781 
2782  if ( r_from == prc_rgn_rgn4pl(i_npl) ) then
2783  do v = 1, vmax
2784  do k = 1, kmax
2785  kk = (v-1) * kmax + k
2786  sendbuf_h2p_sp(kk,i_npl) = var(ij_from,k,l_from,v)
2787  enddo
2788  enddo
2789 
2790  req_count = req_count + 1
2791  totalsize = kmax * vmax
2792  rank = recv_info_p2r(i_prc_from,irank)
2793  tag = recv_info_p2r(i_prc_from,irank) + 1000000
2794 
2795  call mpi_isend( sendbuf_h2p_sp(1,i_npl), & ! [IN]
2796  totalsize, & ! [IN]
2797  mpi_real, & ! [IN]
2798  rank, & ! [IN]
2799  tag, & ! [IN]
2800  prc_local_comm_world, & ! [IN]
2801  req_list(req_count), & ! [OUT]
2802  ierr ) ! [OUT]
2803  endif
2804 
2805  if ( r_from == prc_rgn_rgn4pl(i_spl) ) then
2806  do v = 1, vmax
2807  do k = 1, kmax
2808  kk = (v-1) * kmax + k
2809  sendbuf_h2p_sp(kk,i_spl) = var(ij_from,k,l_from,v)
2810  enddo
2811  enddo
2812 
2813  req_count = req_count + 1
2814  totalsize = kmax * vmax
2815  rank = recv_info_p2r(i_prc_from,irank)
2816  tag = recv_info_p2r(i_prc_from,irank) + 2000000
2817 
2818  call mpi_isend( sendbuf_h2p_sp(1,i_spl), & ! [IN]
2819  totalsize, & ! [IN]
2820  mpi_real, & ! [IN]
2821  rank, & ! [IN]
2822  tag, & ! [IN]
2823  prc_local_comm_world, & ! [IN]
2824  req_list(req_count), & ! [OUT]
2825  ierr ) ! [OUT]
2826  endif
2827  enddo
2828  enddo
2829 
2830  !--- copy p2r-reverse
2831  do irank = 1, copy_nmax_p2r
2832  do ipos = 1, copy_info_p2r(i_size)
2833  ij_from = copy_list_p2r(i_grid_to ,ipos)
2834  l_from = copy_list_p2r(i_l_to ,ipos)
2835  r_from = prc_rgn_lp2r(l_from,copy_info_p2r(i_prc_to))
2836  ij_to = copy_list_p2r(i_grid_from,ipos)
2837  l_to = copy_list_p2r(i_l_from ,ipos)
2838  r_to = prc_rgn_lp2r(l_to,copy_info_p2r(i_prc_from))
2839 
2840  if ( r_from == prc_rgn_rgn4pl(i_npl) &
2841  .OR. r_from == prc_rgn_rgn4pl(i_spl) ) then
2842  do v = 1, vmax
2843  do k = 1, kmax
2844  var_pl(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2845  enddo
2846  enddo
2847  endif
2848  enddo
2849  enddo
2850 
2851  !--- wait all
2852  if ( req_count > 0 ) then
2853  call mpi_waitall( req_count, & ! [IN]
2854  req_list(:), & ! [IN]
2855  mpi_statuses_ignore, & ! [OUT]
2856  ierr ) ! [OUT]
2857  endif
2858 
2859  !--- unpack p2r-reverse
2860  do irank = 1, send_nmax_p2r
2861  do ipos = 1, send_info_p2r(i_size,irank)
2862  l_from = send_list_p2r(i_l_to ,ipos,irank)
2863  r_from = prc_rgn_lp2r(l_from,send_info_p2r(i_prc_to,irank))
2864 
2865  ij_to = send_list_p2r(i_grid_from,ipos,irank)
2866  l_to = send_list_p2r(i_l_from ,ipos,irank)
2867 
2868  if ( r_from == prc_rgn_rgn4pl(i_npl) ) then
2869  do v = 1, vmax
2870  do k = 1, kmax
2871  kk = (v-1) * kmax + k
2872  var_pl(ij_to,k,l_to,v) = recvbuf_h2p_sp(kk,i_npl)
2873  enddo
2874  enddo
2875  endif
2876 
2877  if ( r_from == prc_rgn_rgn4pl(i_spl) ) then
2878  do v = 1, vmax
2879  do k = 1, kmax
2880  kk = (v-1) * kmax + k
2881  var_pl(ij_to,k,l_to,v) = recvbuf_h2p_sp(kk,i_spl)
2882  enddo
2883  enddo
2884  endif
2885  enddo
2886  enddo
2887 
2888  endif
2889 
2890  call comm_data_transfer_sp(var,var_pl)
2891 
2892  call prof_rapend('COMM_var',2)
2893 
2894  return
2895  end subroutine comm_var_sp
2896 
2897  !-----------------------------------------------------------------------------
2899  subroutine comm_var_dp( &
2900  var, &
2901  var_pl, &
2902  kmax, &
2903  vmax )
2904  use scale_prc, only: &
2906  use scale_prc_icoa, only: &
2907  prc_rgn_total_pl, &
2908  prc_rgn_rgn4pl, &
2909  prc_rgn_lp2r, &
2910  i_npl, &
2911  i_spl
2912  implicit none
2913 
2914  integer, intent(in) :: kmax
2915  integer, intent(in) :: vmax
2916  real(DP), intent(inout) :: var (adm_gall, kmax,adm_lall, vmax)
2917  real(DP), intent(inout) :: var_pl(adm_gall_pl,kmax,adm_lall_pl,vmax)
2918 
2919  real(DP) :: sendbuf_h2p_DP(kmax*vmax,prc_rgn_total_pl)
2920  real(DP) :: recvbuf_h2p_DP(kmax*vmax,prc_rgn_total_pl)
2921 
2922  integer :: totalsize, rank, tag
2923  integer :: irank, ipos
2924  integer :: ij_from, l_from, ij_to, l_to
2925  integer :: r_from, r_to
2926 
2927  integer :: k, v, kk
2928  integer :: ierr
2929  !---------------------------------------------------------------------------
2930 
2931  if ( comm_apply_barrier ) then
2932  call prof_rapstart('COMM_barrier',2)
2933  call mpi_barrier(prc_local_comm_world,ierr)
2934  call prof_rapend ('COMM_barrier',2)
2935  endif
2936 
2937  call prof_rapstart('COMM_var',2)
2938 
2939  if ( comm_pl ) then
2940 
2941  req_count = 0
2942 
2943  !--- receive p2r-reverse
2944  do irank = 1, send_nmax_p2r
2945  do ipos = 1, send_info_p2r(i_size,irank)
2946  l_from = send_list_p2r(i_l_to ,ipos,irank)
2947  r_from = prc_rgn_lp2r(l_from,send_info_p2r(i_prc_to,irank))
2948 
2949  if ( r_from == prc_rgn_rgn4pl(i_npl) ) then
2950  req_count = req_count + 1
2951  totalsize = kmax * vmax
2952  rank = send_info_p2r(i_prc_to ,irank)
2953  tag = send_info_p2r(i_prc_from,irank) + 1000000
2954 
2955  call mpi_irecv( recvbuf_h2p_dp(1,i_npl), & ! [OUT]
2956  totalsize, & ! [IN]
2957  mpi_double_precision, & ! [IN]
2958  rank, & ! [IN]
2959  tag, & ! [IN]
2960  prc_local_comm_world, & ! [IN]
2961  req_list(req_count), & ! [OUT]
2962  ierr ) ! [OUT]
2963  endif
2964 
2965  if ( r_from == prc_rgn_rgn4pl(i_spl) ) then
2966  req_count = req_count + 1
2967  totalsize = kmax * vmax
2968  rank = send_info_p2r(i_prc_to ,irank)
2969  tag = send_info_p2r(i_prc_from,irank) + 2000000
2970 
2971  call mpi_irecv( recvbuf_h2p_dp(1,i_spl), & ! [OUT]
2972  totalsize, & ! [IN]
2973  mpi_double_precision, & ! [IN]
2974  rank, & ! [IN]
2975  tag, & ! [IN]
2976  prc_local_comm_world, & ! [IN]
2977  req_list(req_count), & ! [OUT]
2978  ierr ) ! [OUT]
2979  endif
2980  enddo
2981  enddo
2982 
2983  !--- pack and send p2r-reverse
2984  do irank = 1, recv_nmax_p2r
2985  do ipos = 1, recv_info_p2r(i_size,irank)
2986  ij_from = recv_list_p2r(i_grid_to,ipos,irank)
2987  l_from = recv_list_p2r(i_l_to ,ipos,irank)
2988  r_from = prc_rgn_lp2r(l_from,recv_info_p2r(i_prc_to,irank))
2989 
2990  if ( r_from == prc_rgn_rgn4pl(i_npl) ) then
2991  do v = 1, vmax
2992  do k = 1, kmax
2993  kk = (v-1) * kmax + k
2994  sendbuf_h2p_dp(kk,i_npl) = var(ij_from,k,l_from,v)
2995  enddo
2996  enddo
2997 
2998  req_count = req_count + 1
2999  totalsize = kmax * vmax
3000  rank = recv_info_p2r(i_prc_from,irank)
3001  tag = recv_info_p2r(i_prc_from,irank) + 1000000
3002 
3003  call mpi_isend( sendbuf_h2p_dp(1,i_npl), & ! [IN]
3004  totalsize, & ! [IN]
3005  mpi_double_precision, & ! [IN]
3006  rank, & ! [IN]
3007  tag, & ! [IN]
3008  prc_local_comm_world, & ! [IN]
3009  req_list(req_count), & ! [OUT]
3010  ierr ) ! [OUT]
3011  endif
3012 
3013  if ( r_from == prc_rgn_rgn4pl(i_spl) ) then
3014  do v = 1, vmax
3015  do k = 1, kmax
3016  kk = (v-1) * kmax + k
3017  sendbuf_h2p_dp(kk,i_spl) = var(ij_from,k,l_from,v)
3018  enddo
3019  enddo
3020 
3021  req_count = req_count + 1
3022  totalsize = kmax * vmax
3023  rank = recv_info_p2r(i_prc_from,irank)
3024  tag = recv_info_p2r(i_prc_from,irank) + 2000000
3025 
3026  call mpi_isend( sendbuf_h2p_dp(1,i_spl), & ! [IN]
3027  totalsize, & ! [IN]
3028  mpi_double_precision, & ! [IN]
3029  rank, & ! [IN]
3030  tag, & ! [IN]
3031  prc_local_comm_world, & ! [IN]
3032  req_list(req_count), & ! [OUT]
3033  ierr ) ! [OUT]
3034  endif
3035  enddo
3036  enddo
3037 
3038  !--- copy p2r-reverse
3039  do irank = 1, copy_nmax_p2r
3040  do ipos = 1, copy_info_p2r(i_size)
3041  ij_from = copy_list_p2r(i_grid_to ,ipos)
3042  l_from = copy_list_p2r(i_l_to ,ipos)
3043  r_from = prc_rgn_lp2r(l_from,copy_info_p2r(i_prc_to))
3044  ij_to = copy_list_p2r(i_grid_from,ipos)
3045  l_to = copy_list_p2r(i_l_from ,ipos)
3046  r_to = prc_rgn_lp2r(l_to,copy_info_p2r(i_prc_from))
3047 
3048  if ( r_from == prc_rgn_rgn4pl(i_npl) &
3049  .OR. r_from == prc_rgn_rgn4pl(i_spl) ) then
3050  do v = 1, vmax
3051  do k = 1, kmax
3052  var_pl(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
3053  enddo
3054  enddo
3055  endif
3056  enddo
3057  enddo
3058 
3059  !--- wait all
3060  if ( req_count > 0 ) then
3061  call mpi_waitall( req_count, & ! [IN]
3062  req_list(:), & ! [IN]
3063  mpi_statuses_ignore, & ! [OUT]
3064  ierr ) ! [OUT]
3065  endif
3066 
3067  !--- unpack p2r-reverse
3068  do irank = 1, send_nmax_p2r
3069  do ipos = 1, send_info_p2r(i_size,irank)
3070  l_from = send_list_p2r(i_l_to ,ipos,irank)
3071  r_from = prc_rgn_lp2r(l_from,send_info_p2r(i_prc_to,irank))
3072 
3073  ij_to = send_list_p2r(i_grid_from,ipos,irank)
3074  l_to = send_list_p2r(i_l_from ,ipos,irank)
3075 
3076  if ( r_from == prc_rgn_rgn4pl(i_npl) ) then
3077  do v = 1, vmax
3078  do k = 1, kmax
3079  kk = (v-1) * kmax + k
3080  var_pl(ij_to,k,l_to,v) = recvbuf_h2p_dp(kk,i_npl)
3081  enddo
3082  enddo
3083  endif
3084 
3085  if ( r_from == prc_rgn_rgn4pl(i_spl) ) then
3086  do v = 1, vmax
3087  do k = 1, kmax
3088  kk = (v-1) * kmax + k
3089  var_pl(ij_to,k,l_to,v) = recvbuf_h2p_dp(kk,i_spl)
3090  enddo
3091  enddo
3092  endif
3093  enddo
3094  enddo
3095 
3096  endif
3097 
3098  call comm_data_transfer_dp(var,var_pl)
3099 
3100  call prof_rapend('COMM_var',2)
3101 
3102  return
3103  end subroutine comm_var_dp
3104 
3105  !-----------------------------------------------------------------------------
3108  function suf(i,j) result(suffix)
3109  implicit none
3110 
3111  integer :: suffix
3112  integer :: i, j
3113  !---------------------------------------------------------------------------
3114 
3115  suffix = adm_gall_1d * (j-1) + i
3116 
3117  end function suf
3118 
3119  !-----------------------------------------------------------------------------
3120  subroutine comm_debugtest
3121  use scale_prc, only: &
3122  prc_myrank, &
3124  use scale_prc_icoa, only: &
3125  prc_rgn_l2r, &
3126  prc_rgn_vert_num, &
3127  prc_rgn_have_sgp, &
3128  i_n, &
3129  i_s
3130  implicit none
3131 
3132  real(RP) :: var (adm_gall ,adm_kall,adm_lall ,4)
3133  real(RP) :: var_pl(adm_gall_pl,adm_kall,adm_lall_pl,4)
3134 
3135  integer :: i, j, k, l, ij, rgnid, prc
3136  !---------------------------------------------------------------------------
3137 
3138  if( io_l ) write(io_fid_log,*)
3139  if( io_l ) write(io_fid_log,*) '+++ TEST start'
3140 
3141  var(:,:,:,:) = -999.0_rp
3142  var_pl(:,:,:,:) = -999.0_rp
3143 
3144  do l = 1, adm_lall
3145  rgnid = prc_rgn_l2r(l)
3146  prc = prc_myrank
3147 
3148  do k = adm_kmin, adm_kmax
3149  do j = adm_gmin, adm_gmax
3150  do i = adm_gmin, adm_gmax
3151  ij = adm_gall_1d * (j-1) + i
3152 
3153  var(ij,k,l,1) = real(prc, kind=rp)
3154  var(ij,k,l,2) = real(rgnid,kind=rp)
3155  var(ij,k,l,3) = real(i, kind=rp)
3156  var(ij,k,l,4) = real(j, kind=rp)
3157  enddo
3158  enddo
3159  enddo
3160 
3161  if ( prc_rgn_have_sgp(l) ) then
3162  do k = adm_kmin, adm_kmax
3163  var(1,k,l,:) = -1.0_rp
3164  enddo
3165  endif
3166  enddo
3167 
3168  do l = 1, adm_lall_pl
3169  rgnid = l
3170  prc = prc_myrank
3171 
3172  do k = adm_kmin, adm_kmax
3173  do ij = 1, adm_gall_pl
3174  var_pl(ij,k,l,1) = real(-prc, kind=rp)
3175  var_pl(ij,k,l,2) = real(-rgnid,kind=rp)
3176  var_pl(ij,k,l,3) = real(-ij, kind=rp)
3177  var_pl(ij,k,l,4) = real(-ij, kind=rp)
3178  enddo
3179  enddo
3180  enddo
3181 
3182  if( io_l ) write(io_fid_log,*) '##### (prc,rgnid) #####'
3183  do l = 1, adm_lall
3184  do k = adm_kmin, adm_kmin
3185  if( io_l ) write(io_fid_log,*)
3186  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3187  do i = 1, adm_gall_1d
3188  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3189  enddo
3190  if( io_l ) write(io_fid_log,*)
3191 
3192  do j = adm_gall_1d, 1, -1
3193  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3194  do i = 1, adm_gall_1d
3195  ij = adm_gall_1d * (j-1) + i
3196  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3197  '(',int(var(ij,k,l,1)),',',int(var(ij,k,l,2)),')'
3198  enddo
3199  if( io_l ) write(io_fid_log,*)
3200  enddo
3201  enddo
3202  enddo
3203 
3204  do l = 1, adm_lall_pl
3205  do k = adm_kmin, adm_kmin
3206  if( io_l ) write(io_fid_log,*)
3207  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3208  do i = 1, adm_gall_pl
3209  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3210  enddo
3211  if( io_l ) write(io_fid_log,*)
3212 
3213  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3214  do ij = 1, adm_gall_pl
3215  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3216  '(',int(var_pl(ij,k,l,1)),',',int(var_pl(ij,k,l,2)),')'
3217  enddo
3218  if( io_l ) write(io_fid_log,*)
3219  enddo
3220  enddo
3221 
3222  if( io_l ) write(io_fid_log,*) '##### (i,j) #####'
3223  do l = 1, adm_lall
3224  do k = adm_kmin, adm_kmin
3225  if( io_l ) write(io_fid_log,*)
3226  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3227  do i = 1, adm_gall_1d
3228  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3229  enddo
3230  if( io_l ) write(io_fid_log,*)
3231 
3232  do j = adm_gall_1d, 1, -1
3233  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3234  do i = 1, adm_gall_1d
3235  ij = adm_gall_1d * (j-1) + i
3236  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3237  '(',int(var(ij,k,l,3)),',',int(var(ij,k,l,4)),')'
3238  enddo
3239  if( io_l ) write(io_fid_log,*)
3240  enddo
3241  enddo
3242  enddo
3243 
3244  do l = 1, adm_lall_pl
3245  do k = adm_kmin, adm_kmin
3246  if( io_l ) write(io_fid_log,*)
3247  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3248  do i = 1, adm_gall_pl
3249  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3250  enddo
3251  if( io_l ) write(io_fid_log,*)
3252 
3253  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3254  do ij = 1, adm_gall_pl
3255  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3256  '(',int(var_pl(ij,k,l,3)),',',int(var_pl(ij,k,l,4)),')'
3257  enddo
3258  if( io_l ) write(io_fid_log,*)
3259  enddo
3260  enddo
3261 
3262 
3263 
3264  if( io_l ) write(io_fid_log,*)
3265  if( io_l ) write(io_fid_log,*) '+++ Communication start'
3266 
3267  call comm_data_transfer( var(:,:,:,:), var_pl(:,:,:,:) )
3268 
3269  if( io_l ) write(io_fid_log,*) '##### (prc,rgnid) #####'
3270  do l = 1, adm_lall
3271  do k = adm_kmin, adm_kmin
3272  if( io_l ) write(io_fid_log,*)
3273  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3274  do i = 1, adm_gall_1d
3275  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3276  enddo
3277  if( io_l ) write(io_fid_log,*)
3278 
3279  do j = adm_gall_1d, 1, -1
3280  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3281  do i = 1, adm_gall_1d
3282  ij = adm_gall_1d * (j-1) + i
3283  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3284  '(',int(var(ij,k,l,1)),',',int(var(ij,k,l,2)),')'
3285  enddo
3286  if( io_l ) write(io_fid_log,*)
3287  enddo
3288  enddo
3289  enddo
3290 
3291  do l = 1, adm_lall_pl
3292  do k = adm_kmin, adm_kmin
3293  if( io_l ) write(io_fid_log,*)
3294  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3295  do i = 1, adm_gall_pl
3296  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3297  enddo
3298  if( io_l ) write(io_fid_log,*)
3299 
3300  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3301  do ij = 1, adm_gall_pl
3302  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3303  '(',int(var_pl(ij,k,l,1)),',',int(var_pl(ij,k,l,2)),')'
3304  enddo
3305  if( io_l ) write(io_fid_log,*)
3306  enddo
3307  enddo
3308 
3309  if( io_l ) write(io_fid_log,*) '##### (i,j) #####'
3310  do l = 1, adm_lall
3311  do k = adm_kmin, adm_kmin
3312  if( io_l ) write(io_fid_log,*)
3313  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3314  do i = 1, adm_gall_1d
3315  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3316  enddo
3317  if( io_l ) write(io_fid_log,*)
3318 
3319  do j = adm_gall_1d, 1, -1
3320  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3321  do i = 1, adm_gall_1d
3322  ij = adm_gall_1d * (j-1) + i
3323  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3324  '(',int(var(ij,k,l,3)),',',int(var(ij,k,l,4)),')'
3325  enddo
3326  if( io_l ) write(io_fid_log,*)
3327  enddo
3328  enddo
3329  enddo
3330 
3331  do l = 1, adm_lall_pl
3332  do k = adm_kmin, adm_kmin
3333  if( io_l ) write(io_fid_log,*)
3334  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3335  do i = 1, adm_gall_pl
3336  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3337  enddo
3338  if( io_l ) write(io_fid_log,*)
3339 
3340  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3341  do ij = 1, adm_gall_pl
3342  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3343  '(',int(var_pl(ij,k,l,3)),',',int(var_pl(ij,k,l,4)),')'
3344  enddo
3345  if( io_l ) write(io_fid_log,*)
3346  enddo
3347  enddo
3348 
3349 
3350 
3351  do l = 1, adm_lall
3352  rgnid = prc_rgn_l2r(l)
3353  prc = prc_myrank
3354 
3355  if ( prc_rgn_vert_num(i_n,rgnid) == adm_vlink ) then
3356  do k = adm_kmin, adm_kmax
3357  j = adm_gmax+1
3358  i = adm_gmin
3359  ij = adm_gall_1d * (j-1) + i
3360 
3361  var(ij,k,l,1) = real(prc, kind=rp)
3362  var(ij,k,l,2) = real(rgnid,kind=rp)
3363  var(ij,k,l,3) = real(i, kind=rp)
3364  var(ij,k,l,4) = real(j, kind=rp)
3365  enddo
3366  endif
3367 
3368  if ( prc_rgn_vert_num(i_s,rgnid) == adm_vlink ) then
3369  do k = adm_kmin, adm_kmax
3370  j = adm_gmin
3371  i = adm_gmax+1
3372  ij = adm_gall_1d * (j-1) + i
3373 
3374  var(ij,k,l,1) = real(prc, kind=rp)
3375  var(ij,k,l,2) = real(rgnid,kind=rp)
3376  var(ij,k,l,3) = real(i, kind=rp)
3377  var(ij,k,l,4) = real(j, kind=rp)
3378  enddo
3379  endif
3380 
3381  enddo
3382 
3383  if( io_l ) write(io_fid_log,*)
3384  if( io_l ) write(io_fid_log,*) '+++ pole fill start'
3385 
3386  call comm_var( var(:,:,:,:), var_pl(:,:,:,:), adm_kall, 4 )
3387 
3388  if( io_l ) write(io_fid_log,*) '##### (prc,rgnid) #####'
3389  do l = 1, adm_lall
3390  do k = adm_kmin, adm_kmin
3391  if( io_l ) write(io_fid_log,*)
3392  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3393  do i = 1, adm_gall_1d
3394  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3395  enddo
3396  if( io_l ) write(io_fid_log,*)
3397 
3398  do j = adm_gall_1d, 1, -1
3399  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3400  do i = 1, adm_gall_1d
3401  ij = adm_gall_1d * (j-1) + i
3402  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3403  '(',int(var(ij,k,l,1)),',',int(var(ij,k,l,2)),')'
3404  enddo
3405  if( io_l ) write(io_fid_log,*)
3406  enddo
3407  enddo
3408  enddo
3409 
3410  do l = 1, adm_lall_pl
3411  do k = adm_kmin, adm_kmin
3412  if( io_l ) write(io_fid_log,*)
3413  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3414  do i = 1, adm_gall_pl
3415  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3416  enddo
3417  if( io_l ) write(io_fid_log,*)
3418 
3419  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3420  do ij = 1, adm_gall_pl
3421  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3422  '(',int(var_pl(ij,k,l,1)),',',int(var_pl(ij,k,l,2)),')'
3423  enddo
3424  if( io_l ) write(io_fid_log,*)
3425  enddo
3426  enddo
3427 
3428  if( io_l ) write(io_fid_log,*) '##### (i,j) #####'
3429  do l = 1, adm_lall
3430  do k = adm_kmin, adm_kmin
3431  if( io_l ) write(io_fid_log,*)
3432  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3433  do i = 1, adm_gall_1d
3434  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3435  enddo
3436  if( io_l ) write(io_fid_log,*)
3437 
3438  do j = adm_gall_1d, 1, -1
3439  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3440  do i = 1, adm_gall_1d
3441  ij = adm_gall_1d * (j-1) + i
3442  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3443  '(',int(var(ij,k,l,3)),',',int(var(ij,k,l,4)),')'
3444  enddo
3445  if( io_l ) write(io_fid_log,*)
3446  enddo
3447  enddo
3448  enddo
3449 
3450  do l = 1, adm_lall_pl
3451  do k = adm_kmin, adm_kmin
3452  if( io_l ) write(io_fid_log,*)
3453  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3454  do i = 1, adm_gall_pl
3455  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3456  enddo
3457  if( io_l ) write(io_fid_log,*)
3458 
3459  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3460  do ij = 1, adm_gall_pl
3461  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3462  '(',int(var_pl(ij,k,l,3)),',',int(var_pl(ij,k,l,4)),')'
3463  enddo
3464  if( io_l ) write(io_fid_log,*)
3465  enddo
3466  enddo
3467 
3468  call prc_mpifinish
3469 
3470  stop
3471  end subroutine comm_debugtest
3472 
3473  !-----------------------------------------------------------------------------
3474  subroutine comm_stat_sum_sp( localsum, globalsum )
3475  use scale_prc, only: &
3477  prc_nprocs
3478  implicit none
3479 
3480  real(SP), intent(in) :: localsum
3481  real(SP), intent(out) :: globalsum
3482 
3483  real(SP) :: sendbuf(1)
3484  real(SP) :: recvbuf(prc_nprocs)
3485 
3486  integer :: ierr
3487  !---------------------------------------------------------------------------
3488 
3489  if ( comm_pl ) then
3490  sendbuf(1) = localsum
3491 
3492  call mpi_allgather( sendbuf, &
3493  1, &
3494  mpi_real, &
3495  recvbuf, &
3496  1, &
3497  mpi_real, &
3499  ierr )
3500 
3501  globalsum = sum( recvbuf(:) )
3502  else
3503  globalsum = localsum
3504  endif
3505 
3506  return
3507  end subroutine comm_stat_sum_sp
3508 
3509  !-----------------------------------------------------------------------------
3510  subroutine comm_stat_sum_dp( localsum, globalsum )
3511  use scale_prc, only: &
3513  prc_nprocs
3514  implicit none
3515 
3516  real(DP), intent(in) :: localsum
3517  real(DP), intent(out) :: globalsum
3518 
3519  real(DP) :: sendbuf(1)
3520  real(DP) :: recvbuf(prc_nprocs)
3521 
3522  integer :: ierr
3523  !---------------------------------------------------------------------------
3524 
3525  if ( comm_pl ) then
3526  sendbuf(1) = localsum
3527 
3528  call mpi_allgather( sendbuf, &
3529  1, &
3530  mpi_double_precision, &
3531  recvbuf, &
3532  1, &
3533  mpi_double_precision, &
3535  ierr )
3536 
3537  globalsum = sum( recvbuf(:) )
3538  else
3539  globalsum = localsum
3540  endif
3541 
3542  return
3543  end subroutine comm_stat_sum_dp
3544 
3545  !-----------------------------------------------------------------------------
3546  subroutine comm_stat_sum_eachlayer_sp( kall, localsum, globalsum )
3547  use scale_prc, only: &
3549  prc_nprocs
3550  implicit none
3551 
3552  integer, intent(in) :: kall
3553  real(SP), intent(in) :: localsum (kall)
3554  real(SP), intent(out) :: globalsum(kall)
3555 
3556  real(SP) :: sendbuf(kall)
3557  integer :: displs (prc_nprocs)
3558  integer :: counts (prc_nprocs)
3559  real(SP) :: recvbuf(kall,prc_nprocs)
3560 
3561  integer :: ierr
3562  integer :: k, p
3563  !---------------------------------------------------------------------------
3564 
3565  do p = 1, prc_nprocs
3566  displs(p) = (p-1) * kall
3567  counts(p) = kall
3568  enddo
3569 
3570  if ( comm_pl ) then
3571  sendbuf(:) = localsum(:)
3572 
3573  call mpi_allgatherv( sendbuf, &
3574  kall, &
3575  mpi_real, &
3576  recvbuf, &
3577  counts, &
3578  displs, &
3579  mpi_real, &
3581  ierr )
3582 
3583  do k = 1, kall
3584  globalsum(k) = sum( recvbuf(k,:) )
3585  enddo
3586  else
3587  do k = 1, kall
3588  globalsum(k) = localsum(k)
3589  enddo
3590  endif
3591 
3592  return
3593  end subroutine comm_stat_sum_eachlayer_sp
3594 
3595  !-----------------------------------------------------------------------------
3596  subroutine comm_stat_sum_eachlayer_dp( kall, localsum, globalsum )
3597  use scale_prc, only: &
3599  prc_nprocs
3600  implicit none
3601 
3602  integer, intent(in) :: kall
3603  real(DP), intent(in) :: localsum (kall)
3604  real(DP), intent(out) :: globalsum(kall)
3605 
3606  real(DP) :: sendbuf(kall)
3607  integer :: displs (prc_nprocs)
3608  integer :: counts (prc_nprocs)
3609  real(DP) :: recvbuf(kall,prc_nprocs)
3610 
3611  integer :: ierr
3612  integer :: k, p
3613  !---------------------------------------------------------------------------
3614 
3615  do p = 1, prc_nprocs
3616  displs(p) = (p-1) * kall
3617  counts(p) = kall
3618  enddo
3619 
3620  if ( comm_pl ) then
3621  sendbuf(:) = localsum(:)
3622 
3623  call mpi_allgatherv( sendbuf, &
3624  kall, &
3625  mpi_double_precision, &
3626  recvbuf, &
3627  counts, &
3628  displs, &
3629  mpi_double_precision, &
3631  ierr )
3632 
3633  do k = 1, kall
3634  globalsum(k) = sum( recvbuf(k,:) )
3635  enddo
3636  else
3637  do k = 1, kall
3638  globalsum(k) = localsum(k)
3639  enddo
3640  endif
3641 
3642  return
3643  end subroutine comm_stat_sum_eachlayer_dp
3644 
3645  !-----------------------------------------------------------------------------
3646  subroutine comm_stat_avg_sp( localavg, globalavg )
3647  use scale_prc, only: &
3649  prc_nprocs
3650  implicit none
3651 
3652  real(SP), intent(in) :: localavg
3653  real(SP), intent(out) :: globalavg
3654 
3655  real(SP) :: sendbuf(1)
3656  real(SP) :: recvbuf(prc_nprocs)
3657 
3658  integer :: ierr
3659  !---------------------------------------------------------------------------
3660 
3661  if ( comm_pl ) then
3662  sendbuf(1) = localavg
3663 
3664  call mpi_allgather( sendbuf, &
3665  1, &
3666  mpi_real, &
3667  recvbuf, &
3668  1, &
3669  mpi_real, &
3671  ierr )
3672 
3673  globalavg = sum( recvbuf(:) ) / real(prc_nprocs,kind=sp)
3674  else
3675  globalavg = localavg
3676  endif
3677 
3678  return
3679  end subroutine comm_stat_avg_sp
3680 
3681  !-----------------------------------------------------------------------------
3682  subroutine comm_stat_avg_dp( localavg, globalavg )
3683  use scale_prc, only: &
3685  prc_nprocs
3686  implicit none
3687 
3688  real(DP), intent(in) :: localavg
3689  real(DP), intent(out) :: globalavg
3690 
3691  real(DP) :: sendbuf(1)
3692  real(DP) :: recvbuf(prc_nprocs)
3693 
3694  integer :: ierr
3695  !---------------------------------------------------------------------------
3696 
3697  if ( comm_pl ) then
3698  sendbuf(1) = localavg
3699 
3700  call mpi_allgather( sendbuf, &
3701  1, &
3702  mpi_double_precision, &
3703  recvbuf, &
3704  1, &
3705  mpi_double_precision, &
3707  ierr )
3708 
3709  globalavg = sum( recvbuf(:) ) / real(prc_nprocs,kind=dp)
3710  else
3711  globalavg = localavg
3712  endif
3713 
3714  return
3715  end subroutine comm_stat_avg_dp
3716 
3717  !-----------------------------------------------------------------------------
3718  subroutine comm_stat_max_sp( localmax, globalmax )
3719  use scale_prc, only: &
3721  prc_nprocs
3722  implicit none
3723 
3724  real(SP), intent(in) :: localmax
3725  real(SP), intent(out) :: globalmax
3726 
3727  real(SP) :: sendbuf(1)
3728  real(SP) :: recvbuf(prc_nprocs)
3729 
3730  integer :: ierr
3731  !---------------------------------------------------------------------------
3732 
3733  sendbuf(1) = localmax
3734 
3735  call mpi_allgather( sendbuf, &
3736  1, &
3737  mpi_real, &
3738  recvbuf, &
3739  1, &
3740  mpi_real, &
3742  ierr )
3743 
3744  globalmax = maxval( recvbuf(:) )
3745 
3746  return
3747  end subroutine comm_stat_max_sp
3748 
3749  !-----------------------------------------------------------------------------
3750  subroutine comm_stat_max_dp( localmax, globalmax )
3751  use scale_prc, only: &
3753  prc_nprocs
3754  implicit none
3755 
3756  real(DP), intent(in) :: localmax
3757  real(DP), intent(out) :: globalmax
3758 
3759  real(DP) :: sendbuf(1)
3760  real(DP) :: recvbuf(prc_nprocs)
3761 
3762  integer :: ierr
3763  !---------------------------------------------------------------------------
3764 
3765  sendbuf(1) = localmax
3766 
3767  call mpi_allgather( sendbuf, &
3768  1, &
3769  mpi_double_precision, &
3770  recvbuf, &
3771  1, &
3772  mpi_double_precision, &
3774  ierr )
3775 
3776  globalmax = maxval( recvbuf(:) )
3777 
3778  return
3779  end subroutine comm_stat_max_dp
3780 
3781  !-----------------------------------------------------------------------------
3782  subroutine comm_stat_min_sp( localmin, globalmin )
3783  use scale_prc, only: &
3785  prc_nprocs
3786  implicit none
3787 
3788  real(SP), intent(in) :: localmin
3789  real(SP), intent(out) :: globalmin
3790 
3791  real(SP) :: sendbuf(1)
3792  real(SP) :: recvbuf(prc_nprocs)
3793 
3794  integer :: ierr
3795  !---------------------------------------------------------------------------
3796 
3797  sendbuf(1) = localmin
3798 
3799  call mpi_allgather( sendbuf, &
3800  1, &
3801  mpi_real, &
3802  recvbuf, &
3803  1, &
3804  mpi_real, &
3806  ierr )
3807 
3808  globalmin = minval( recvbuf(:) )
3809 
3810  return
3811  end subroutine comm_stat_min_sp
3812 
3813  !-----------------------------------------------------------------------------
3814  subroutine comm_stat_min_dp( localmin, globalmin )
3815  use scale_prc, only: &
3817  prc_nprocs
3818  implicit none
3819 
3820  real(DP), intent(in) :: localmin
3821  real(DP), intent(out) :: globalmin
3822 
3823  real(DP) :: sendbuf(1)
3824  real(DP) :: recvbuf(prc_nprocs)
3825 
3826  integer :: ierr
3827  !---------------------------------------------------------------------------
3828 
3829  sendbuf(1) = localmin
3830 
3831  call mpi_allgather( sendbuf, &
3832  1, &
3833  mpi_double_precision, &
3834  recvbuf, &
3835  1, &
3836  mpi_double_precision, &
3838  ierr )
3839 
3840  globalmin = minval( recvbuf(:) )
3841 
3842  return
3843  end subroutine comm_stat_min_dp
3844 
3845 end module scale_comm_icoa
real(dp), dimension(:,:), allocatable, public recvbuf_r2r_dp
integer, parameter, public i_prc
process
integer, dimension(:,:,:), allocatable, public send_list_p2r
subroutine comm_var_sp(var, var_pl, kmax, vmax)
Data transfer with region halo => pole center.
integer, dimension(:,:,:), allocatable, public recv_list_r2r
subroutine comm_stat_sum_dp(localsum, globalsum)
subroutine, public comm_setup
Setup.
integer, dimension(:,:,:), allocatable, public prc_rgn_edge_tab
region link information (for 4 edges)
integer, parameter, public i_e
east
integer, parameter, public i_l
local region
integer, parameter, public adm_gmin_pl
integer, dimension(:,:,:), allocatable, public recv_list_p2r
integer, parameter, public i_se
south east
integer, parameter, public i_spl
south pole
real(sp), dimension(:,:), allocatable, public sendbuf_r2p_sp
integer, parameter, public adm_gslf_pl
module COMMUNICATION
integer, dimension(:,:,:,:), allocatable, public prc_rgn_vert_tab
region link information (for 4 vertexes)
integer, dimension(:), allocatable, public copy_info_p2r
integer, dimension(:,:), allocatable, public recv_info_p2r
integer, dimension(:), allocatable, public prc_rgn_l2r
l,prc_me => rgn
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
integer, dimension(:,:), allocatable, public send_info_r2p
integer, parameter, public adm_lall_pl
subroutine comm_stat_sum_sp(localsum, globalsum)
real(dp), dimension(:,:), allocatable, public recvbuf_r2p_dp
integer, dimension(:), allocatable, public singular_info
integer, dimension(:,:), allocatable, public recv_info_r2p
integer, public prc_nprocs
myrank in local communicator
Definition: scale_prc.F90:88
integer, parameter, public i_dir
direction
integer, dimension(:,:), allocatable, public singular_list
subroutine comm_data_transfer_dp(var, var_pl)
Data transfer kernel.
subroutine comm_stat_avg_sp(localavg, globalavg)
integer, parameter, public i_sw
south west
logical, public io_l
output log or not? (this process)
Definition: scale_io.F90:61
integer, dimension(:,:), allocatable, public send_info_p2r
integer, parameter, public i_w
west
real(sp), dimension(:,:), allocatable, public sendbuf_r2r_sp
integer, parameter, public i_nw
north west
integer, dimension(:), allocatable, public copy_info_r2p
integer, dimension(:,:), allocatable, public recv_info_r2r
real(dp), dimension(:,:), allocatable, public recvbuf_p2r_dp
integer, dimension(:,:), allocatable, public prc_rgn_vert_num
number of region around the vertex (4 vertexes)
integer, dimension(:,:), allocatable, public copy_list_p2r
integer, dimension(:,:,:), allocatable, public recv_list_r2p
module PROCESS
Definition: scale_prc.F90:11
integer, dimension(prc_rgn_total_pl), public prc_rgn_r2p_pl
process ID which have the pole regions
subroutine comm_stat_max_dp(localmax, globalmax)
integer, public comm_datatype
integer, dimension(:), allocatable, public req_list
integer, dimension(prc_rgn_total_pl), public prc_rgn_rgn4pl
region, having pole data in the halo
real(dp), dimension(:,:), allocatable, public sendbuf_r2r_dp
subroutine comm_stat_sum_eachlayer_sp(kall, localsum, globalsum)
integer, dimension(:), allocatable, public copy_info_r2r
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:89
integer, parameter, public i_ne
north east
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
real(sp), dimension(:,:), allocatable, public recvbuf_p2r_sp
integer, dimension(:,:), allocatable, public send_info_r2r
real(dp), dimension(:,:), allocatable, public sendbuf_p2r_dp
real(sp), dimension(:,:), allocatable, public recvbuf_r2p_sp
subroutine comm_var_dp(var, var_pl, kmax, vmax)
Data transfer with region halo => pole center.
subroutine comm_stat_max_sp(localmax, globalmax)
subroutine, public comm_data_transfer_nopl(var)
Data transfer kernel.
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:157
integer, public prc_rgn_local
number of regular region (local)
integer, parameter, public i_npl
north pole
module profiler
Definition: scale_prof.F90:11
integer, parameter, public i_rgnid
region id
integer, dimension(:,:,:), allocatable, public prc_rgn_vert_tab_pl
region link information (for 4 vertexes)
integer function suf(i, j)
suffix calculation
module atmosphere / grid / icosahedralA / index
integer, parameter, public prc_rgn_total_pl
number of pole region
integer, dimension(:,:), allocatable, public prc_rgn_lp2r
l,prc => rgn
subroutine comm_stat_avg_dp(localavg, globalavg)
subroutine, public prc_mpifinish
Stop MPI peacefully.
Definition: scale_prc.F90:354
real(dp), dimension(:,:), allocatable, public sendbuf_r2p_dp
module PRECISION
integer, dimension(:,:,:), allocatable, public send_list_r2r
subroutine comm_data_transfer_sp(var, var_pl)
Data transfer kernel.
logical, public comm_pl
integer, parameter, public i_s
south
subroutine comm_stat_min_sp(localmin, globalmin)
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:87
integer, parameter, public sp
logical, public io_nml
output log or not? (for namelist, this process)
Definition: scale_io.F90:62
module STDIO
Definition: scale_io.F90:10
subroutine comm_stat_sum_eachlayer_dp(kall, localsum, globalsum)
integer, parameter, public i_n
north
logical, public prc_have_pl
this ID manages pole region?
subroutine comm_stat_min_dp(localmin, globalmin)
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:210
integer, public io_fid_log
Log file ID.
Definition: scale_io.F90:56
integer, dimension(:,:,:), allocatable, public send_list_r2p
integer, parameter, public rp
real(sp), dimension(:,:), allocatable, public sendbuf_p2r_sp
module process / icoA
logical, dimension(:), allocatable, public prc_rgn_have_sgp
region have singlar point?
integer, parameter, public dp
integer, dimension(:,:), allocatable, public prc_rgn_r2lp
rgn => l,prc
real(sp), dimension(:,:), allocatable, public recvbuf_r2r_sp
integer, dimension(:,:), allocatable, public copy_list_r2r
integer, dimension(:,:), allocatable, public copy_list_r2p