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