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: &
983  prc_have_pl, &
984  prc_rgn_local, &
985  prc_rgn_vert_pl, &
987  prc_rgn_lp2r, &
988  prc_rgn_r2lp, &
989  prc_rgn_l2r, &
990  prc_rgn_r2p_pl, &
991  i_l, &
992  i_prc, &
993  i_rgnid, &
994  i_n, &
995  i_s, &
996  i_npl, &
997  i_spl
998  implicit none
999 
1000  integer :: prc, prc_rmt
1001  integer :: rgnid, rgnid_rmt
1002  logical :: check_vert_pl
1003  integer :: pl_to
1004 
1005  integer :: irank, ipos
1006 
1007  integer :: send_size_nglobal_pl
1008 
1009  integer :: l, l_pl, n, v, vv
1010  integer :: i_from, j_from, r_from, g_from, l_from, p_from
1011  integer :: i_to, j_to, r_to, g_to, l_to, p_to
1012  !---------------------------------------------------------------------------
1013 
1014  send_size_nglobal_pl = prc_rgn_ndiamond
1015 
1016  allocate( copy_info_p2r(info_vindex) )
1017  allocate( recv_info_p2r(info_vindex,recv_nlim) )
1018  allocate( send_info_p2r(info_vindex,send_nlim) )
1019  copy_info_p2r(:) = -1
1020  recv_info_p2r(:,:) = -1
1021  send_info_p2r(:,:) = -1
1022  copy_info_p2r(i_size) = 0
1023  recv_info_p2r(i_size,:) = 0
1024  send_info_p2r(i_size,:) = 0
1025 
1026  allocate( copy_list_p2r(list_vindex,send_size_nglobal_pl) )
1027  allocate( recv_list_p2r(list_vindex,send_size_nglobal_pl,recv_nlim) )
1028  allocate( send_list_p2r(list_vindex,send_size_nglobal_pl,send_nlim) )
1029  copy_list_p2r(:,:) = -1
1030  recv_list_p2r(:,:,:) = -1
1031  send_list_p2r(:,:,:) = -1
1032 
1033  allocate( copy_info_r2p(info_vindex) )
1034  allocate( recv_info_r2p(info_vindex,recv_nlim) )
1035  allocate( send_info_r2p(info_vindex,send_nlim) )
1036  copy_info_r2p(:) = -1
1037  recv_info_r2p(:,:) = -1
1038  send_info_r2p(:,:) = -1
1039  copy_info_r2p(i_size) = 0
1040  recv_info_r2p(i_size,:) = 0
1041  send_info_r2p(i_size,:) = 0
1042 
1043  allocate( copy_list_r2p(list_vindex,send_size_nglobal_pl) )
1044  allocate( recv_list_r2p(list_vindex,send_size_nglobal_pl,recv_nlim) )
1045  allocate( send_list_r2p(list_vindex,send_size_nglobal_pl,send_nlim) )
1046  copy_list_r2p(:,:) = -1
1047  recv_list_r2p(:,:,:) = -1
1048  send_list_r2p(:,:,:) = -1
1049 
1050  ! Search in regular region
1051  do l = 1, prc_rgn_local
1052  rgnid = prc_rgn_l2r(l)
1053  prc = prc_myrank
1054 
1055  do l_pl = i_npl, i_spl
1056  rgnid_rmt = l_pl
1057  prc_rmt = prc_rgn_r2p_pl(rgnid_rmt)
1058 
1059  if ( rgnid_rmt == i_npl ) then
1060  check_vert_pl = prc_rgn_vert_pl(i_npl,rgnid)
1061  i_from = adm_gmin
1062  j_from = adm_gmax
1063  i_to = adm_gmin
1064  j_to = adm_gmax + 1
1065  elseif( rgnid_rmt == i_spl ) then
1066  check_vert_pl = prc_rgn_vert_pl(i_spl,rgnid)
1067  i_from = adm_gmax
1068  j_from = adm_gmin
1069  i_to = adm_gmax + 1
1070  j_to = adm_gmin
1071  endif
1072 
1073  if ( check_vert_pl ) then
1074 
1075  ! search destination in the pole halo
1076  do v = 1, adm_vlink
1077  if ( rgnid == prc_rgn_vert_tab_pl(i_rgnid,rgnid_rmt,v) ) then
1078  pl_to = v
1079  if( pl_to < adm_gmin_pl ) pl_to = adm_gmax_pl
1080  exit
1081  endif
1082  enddo
1083 
1084  if ( prc == prc_rmt ) then ! no communication
1085  ! copy region inner -> pole halo
1086  copy_info_r2p(i_size) = copy_info_r2p(i_size) + 1
1087  ipos = copy_info_r2p(i_size)
1088 
1089  copy_list_r2p(i_grid_from,ipos) = suf(i_from,j_from)
1090  copy_list_r2p(i_l_from ,ipos) = l
1091  copy_list_r2p(i_grid_to ,ipos) = pl_to
1092  copy_list_r2p(i_l_to ,ipos) = l_pl
1093  else ! node-to-node communication
1094  ! receive pole center -> region halo
1095  irank = -1
1096  do n = 1, recv_nmax_p2r
1097  if ( recv_info_p2r(i_prc_from,n) == prc_rmt ) then
1098  irank = n
1099  exit
1100  endif
1101  enddo
1102 
1103  if ( irank < 0 ) then ! register new rank id
1104  recv_nmax_p2r = recv_nmax_p2r + 1
1105  irank = recv_nmax_p2r
1106 
1107  recv_info_p2r(i_prc_from,irank) = prc_rmt
1108  recv_info_p2r(i_prc_to ,irank) = prc
1109  endif
1110 
1111  recv_info_p2r(i_size,irank) = recv_info_p2r(i_size,irank) + 1
1112  ipos = recv_info_p2r(i_size,irank)
1113 
1114  recv_list_p2r(i_grid_from,ipos,irank) = adm_gslf_pl
1115  recv_list_p2r(i_l_from ,ipos,irank) = l_pl
1116  recv_list_p2r(i_grid_to ,ipos,irank) = suf(i_to,j_to)
1117  recv_list_p2r(i_l_to ,ipos,irank) = l
1118 
1119  ! send region inner -> pole halo
1120  irank = -1
1121  do n = 1, send_nmax_r2p
1122  if ( send_info_r2p(i_prc_to,n) == prc_rmt ) then
1123  irank = n
1124  exit
1125  endif
1126  enddo
1127 
1128  if ( irank < 0 ) then ! register new rank id
1129  send_nmax_r2p = send_nmax_r2p + 1
1130  irank = send_nmax_r2p
1131 
1132  send_info_r2p(i_prc_from,irank) = prc
1133  send_info_r2p(i_prc_to ,irank) = prc_rmt
1134  endif
1135 
1136  send_info_r2p(i_size,irank) = send_info_r2p(i_size,irank) + 1
1137  ipos = send_info_r2p(i_size,irank)
1138 
1139  send_list_r2p(i_grid_from,ipos,irank) = suf(i_from,j_from)
1140  send_list_r2p(i_l_from ,ipos,irank) = l
1141  send_list_r2p(i_grid_to ,ipos,irank) = pl_to
1142  send_list_r2p(i_l_to ,ipos,irank) = l_pl
1143  endif
1144  endif
1145 
1146  enddo ! loop l_pl
1147  enddo ! loop l
1148 
1149  ! Search in pole
1150  if ( prc_have_pl ) then
1151 
1152  do l_pl = i_npl, i_spl
1153  rgnid = l_pl
1154  prc = prc_myrank
1155 
1156  do v = adm_vlink+1, 2, -1
1157  vv = v
1158  if( v == adm_vlink+1 ) vv = 1
1159 
1160  rgnid_rmt = prc_rgn_vert_tab_pl(i_rgnid,rgnid,vv)
1161  prc_rmt = prc_rgn_r2lp(i_prc,rgnid_rmt)
1162 
1163  if ( rgnid == i_npl ) then
1164  i_from = adm_gmin
1165  j_from = adm_gmax
1166  i_to = adm_gmin
1167  j_to = adm_gmax + 1
1168  elseif( rgnid == i_spl ) then
1169  i_from = adm_gmax
1170  j_from = adm_gmin
1171  i_to = adm_gmax + 1
1172  j_to = adm_gmin
1173  endif
1174 
1175  pl_to = vv
1176  if( pl_to < adm_gmin_pl ) pl_to = adm_gmax_pl
1177 
1178  if ( prc == prc_rmt ) then ! no communication
1179  ! copy pole center -> region halo
1180  copy_info_p2r(i_size) = copy_info_p2r(i_size) + 1
1181  ipos = copy_info_p2r(i_size)
1182 
1183  copy_list_p2r(i_grid_from,ipos) = adm_gslf_pl
1184  copy_list_p2r(i_l_from ,ipos) = l_pl
1185  copy_list_p2r(i_grid_to ,ipos) = suf(i_to,j_to)
1186  copy_list_p2r(i_l_to ,ipos) = prc_rgn_r2lp(i_l,rgnid_rmt)
1187  else ! node-to-node communication
1188  ! send region inner -> pole halo
1189  irank = -1
1190  do n = 1, recv_nmax_r2p
1191  if ( recv_info_r2p(i_prc_from,n) == prc_rmt ) then
1192  irank = n
1193  exit
1194  endif
1195  enddo
1196 
1197  if ( irank < 0 ) then ! register new rank id
1198  recv_nmax_r2p = recv_nmax_r2p + 1
1199  irank = recv_nmax_r2p
1200 
1201  recv_info_r2p(i_prc_from,irank) = prc_rmt
1202  recv_info_r2p(i_prc_to ,irank) = prc
1203  endif
1204 
1205  recv_info_r2p(i_size,irank) = recv_info_r2p(i_size,irank) + 1
1206  ipos = recv_info_r2p(i_size,irank)
1207 
1208  recv_list_r2p(i_grid_from,ipos,irank) = suf(i_from,j_from)
1209  recv_list_r2p(i_l_from ,ipos,irank) = prc_rgn_r2lp(i_l,rgnid_rmt)
1210  recv_list_r2p(i_grid_to ,ipos,irank) = pl_to
1211  recv_list_r2p(i_l_to ,ipos,irank) = l_pl
1212 
1213  ! receive pole center -> region halo
1214  irank = -1
1215  do n = 1, send_nmax_p2r
1216  if ( send_info_p2r(i_prc_to,n) == prc_rmt ) then
1217  irank = n
1218  exit
1219  endif
1220  enddo
1221 
1222  if ( irank < 0 ) then ! register new rank id
1223  send_nmax_p2r = send_nmax_p2r + 1
1224  irank = send_nmax_p2r
1225 
1226  send_info_p2r(i_prc_from,irank) = prc
1227  send_info_p2r(i_prc_to ,irank) = prc_rmt
1228  endif
1229 
1230  send_info_p2r(i_size,irank) = send_info_p2r(i_size,irank) + 1
1231  ipos = send_info_p2r(i_size,irank)
1232 
1233  send_list_p2r(i_grid_from,ipos,irank) = adm_gslf_pl
1234  send_list_p2r(i_l_from ,ipos,irank) = l_pl
1235  send_list_p2r(i_grid_to ,ipos,irank) = suf(i_to,j_to)
1236  send_list_p2r(i_l_to ,ipos,irank) = prc_rgn_r2lp(i_l,rgnid_rmt)
1237  endif
1238 
1239  enddo ! loop pole halo
1240 
1241  enddo ! loop l_pl
1242 
1243  if ( copy_info_p2r(i_size) > 0 ) then
1244  copy_nmax_p2r = 1
1245  copy_info_p2r(i_prc_from) = prc_myrank
1246  copy_info_p2r(i_prc_to ) = prc_myrank
1247  endif
1248 
1249  if ( copy_info_r2p(i_size) > 0 ) then
1250  copy_nmax_r2p = 1
1251  copy_info_r2p(i_prc_from) = prc_myrank
1252  copy_info_r2p(i_prc_to ) = prc_myrank
1253  endif
1254 
1255  endif
1256 
1257  if( io_l ) write(io_fid_log,*)
1258  if( io_l ) write(io_fid_log,*) '*** Recv_nmax_p2r(local) = ', recv_nmax_p2r
1259  if( io_l ) write(io_fid_log,*) '*** Send_nmax_p2r(local) = ', send_nmax_p2r
1260  if( io_l ) write(io_fid_log,*)
1261  if( io_l ) write(io_fid_log,'(A)') '|---------------------------------------'
1262  if( io_l ) write(io_fid_log,'(A)') '| size prc_from prc_to'
1263  if( io_l ) write(io_fid_log,'(A10,3(I10))') '| Copy_p2r', copy_info_p2r(:)
1264  do irank = 1, recv_nmax_p2r
1265  if( io_l ) write(io_fid_log,'(A10,3(I10))') '| Recv_p2r', recv_info_p2r(:,irank)
1266  enddo
1267  do irank = 1, send_nmax_p2r
1268  if( io_l ) write(io_fid_log,'(A10,3(I10))') '| Send_p2r', send_info_p2r(:,irank)
1269  enddo
1270 
1271  if ( debug ) then
1272  if( io_l ) write(io_fid_log,*)
1273  if( io_l ) write(io_fid_log,*) '--- Copy_list_p2r'
1274  if( io_l ) write(io_fid_log,*)
1275  if( io_l ) write(io_fid_log,'(11(A6))') 'number', &
1276  '|rfrom','|gfrom','|lfrom','|pfrom', &
1277  '| ito','| jto','| rto','| gto','| lto','| pto'
1278  do ipos = 1, copy_info_p2r(i_size)
1279  g_from = copy_list_p2r(i_grid_from,ipos)
1280  l_from = copy_list_p2r(i_l_from,ipos)
1281  p_from = copy_info_p2r(i_prc_from)
1282  r_from = l_from
1283  g_to = copy_list_p2r(i_grid_to,ipos)
1284  l_to = copy_list_p2r(i_l_to,ipos)
1285  p_to = copy_info_p2r(i_prc_to)
1286  i_to = mod(g_to-1,adm_gall_1d) + 1
1287  j_to = (g_to-i_to) / adm_gall_1d + 1
1288  r_to = prc_rgn_lp2r(l_to,p_to)
1289 
1290  if( io_l ) write(io_fid_log,'(11(I6))') ipos, r_from, g_from, l_from, p_from, &
1291  i_to , j_to , r_to , g_to , l_to , p_to
1292  enddo
1293 
1294  if( io_l ) write(io_fid_log,*)
1295  if( io_l ) write(io_fid_log,*) '--- Recv_list_p2r'
1296  do irank = 1, recv_nmax_p2r
1297  if( io_l ) write(io_fid_log,'(11(A6))') 'number', &
1298  '|rfrom','|gfrom','|lfrom','|pfrom', &
1299  '| ito','| jto','| rto','| gto','| lto','| pto'
1300  do ipos = 1, recv_info_p2r(i_size,irank)
1301  g_from = recv_list_p2r(i_grid_from,ipos,irank)
1302  l_from = recv_list_p2r(i_l_from,ipos,irank)
1303  p_from = recv_info_p2r(i_prc_from,irank)
1304  r_from = l_from
1305  g_to = recv_list_p2r(i_grid_to,ipos,irank)
1306  l_to = recv_list_p2r(i_l_to,ipos,irank)
1307  p_to = recv_info_p2r(i_prc_to,irank)
1308  i_to = mod(g_to-1,adm_gall_1d) + 1
1309  j_to = (g_to-i_to) / adm_gall_1d + 1
1310  r_to = prc_rgn_lp2r(l_to,p_to)
1311 
1312  if( io_l ) write(io_fid_log,'(11(I6))') ipos, r_from, g_from, l_from, p_from, &
1313  i_to , j_to , r_to , g_to , l_to , p_to
1314  enddo
1315  enddo
1316 
1317  if( io_l ) write(io_fid_log,*)
1318  if( io_l ) write(io_fid_log,*) '--- Send_list_p2r'
1319  do irank = 1, send_nmax_p2r
1320  if( io_l ) write(io_fid_log,'(11(A6))') 'number', &
1321  '|rfrom','|gfrom','|lfrom','|pfrom', &
1322  '| ito','| jto','| rto','| gto','| lto','| pto'
1323  do ipos = 1, send_info_p2r(i_size,irank)
1324  g_from = send_list_p2r(i_grid_from,ipos,irank)
1325  l_from = send_list_p2r(i_l_from,ipos,irank)
1326  p_from = send_info_p2r(i_prc_from,irank)
1327  r_from = l_from
1328  g_to = send_list_p2r(i_grid_to,ipos,irank)
1329  l_to = send_list_p2r(i_l_to,ipos,irank)
1330  p_to = send_info_p2r(i_prc_to,irank)
1331  i_to = mod(g_to-1,adm_gall_1d) + 1
1332  j_to = (g_to-i_to) / adm_gall_1d + 1
1333  r_to = prc_rgn_lp2r(l_to,p_to)
1334 
1335  if( io_l ) write(io_fid_log,'(11(I6))') ipos, r_from, g_from, l_from, p_from, &
1336  i_to , j_to , r_to , g_to , l_to , p_to
1337  enddo
1338  enddo
1339  endif
1340 
1341  if( io_l ) write(io_fid_log,*)
1342  if( io_l ) write(io_fid_log,*) '*** Recv_nmax_r2p(local) = ', recv_nmax_r2p
1343  if( io_l ) write(io_fid_log,*) '*** Send_nmax_r2p(local) = ', send_nmax_r2p
1344  if( io_l ) write(io_fid_log,*)
1345  if( io_l ) write(io_fid_log,'(A)') '|---------------------------------------'
1346  if( io_l ) write(io_fid_log,'(A)') '| size prc_from prc_to'
1347  if( io_l ) write(io_fid_log,'(A10,3(I10))') '| Copy_r2p', copy_info_r2p(:)
1348  do irank = 1, recv_nmax_r2p
1349  if( io_l ) write(io_fid_log,'(A10,3(I10))') '| Recv_r2p', recv_info_r2p(:,irank)
1350  enddo
1351  do irank = 1, send_nmax_r2p
1352  if( io_l ) write(io_fid_log,'(A10,3(I10))') '| Send_r2p', send_info_r2p(:,irank)
1353  enddo
1354 
1355  if ( debug ) then
1356  if( io_l ) write(io_fid_log,*)
1357  if( io_l ) write(io_fid_log,*) '--- Copy_list_r2p'
1358  if( io_l ) write(io_fid_log,*)
1359  if( io_l ) write(io_fid_log,'(11(A6))') 'number', &
1360  '|ifrom','|jfrom','|rfrom','|gfrom','|lfrom','|pfrom', &
1361  '| rto','| gto','| lto','| pto'
1362  do ipos = 1, copy_info_r2p(i_size)
1363  g_from = copy_list_r2p(i_grid_from,ipos)
1364  l_from = copy_list_r2p(i_l_from,ipos)
1365  p_from = copy_info_r2p(i_prc_from)
1366  i_from = mod(g_from-1,adm_gall_1d) + 1
1367  j_from = (g_from-i_from) / adm_gall_1d + 1
1368  r_from = prc_rgn_lp2r(l_from,p_from)
1369  g_to = copy_list_r2p(i_grid_to,ipos)
1370  l_to = copy_list_r2p(i_l_to,ipos)
1371  p_to = copy_info_r2p(i_prc_to)
1372  r_to = l_to
1373 
1374  if( io_l ) write(io_fid_log,'(11(I6))') ipos, i_from, j_from, r_from, g_from, l_from, p_from, &
1375  r_to , g_to , l_to , p_to
1376  enddo
1377 
1378  if( io_l ) write(io_fid_log,*)
1379  if( io_l ) write(io_fid_log,*) '--- Recv_list_r2p'
1380  do irank = 1, recv_nmax_r2p
1381  if( io_l ) write(io_fid_log,'(11(A6))') 'number', &
1382  '|ifrom','|jfrom','|rfrom','|gfrom','|lfrom','|pfrom', &
1383  '| rto','| gto','| lto','| pto'
1384  do ipos = 1, recv_info_r2p(i_size,irank)
1385  g_from = recv_list_r2p(i_grid_from,ipos,irank)
1386  l_from = recv_list_r2p(i_l_from,ipos,irank)
1387  p_from = recv_info_r2p(i_prc_from,irank)
1388  i_from = mod(g_from-1,adm_gall_1d) + 1
1389  j_from = (g_from-i_from) / adm_gall_1d + 1
1390  r_from = prc_rgn_lp2r(l_from,p_from)
1391  g_to = recv_list_r2p(i_grid_to,ipos,irank)
1392  l_to = recv_list_r2p(i_l_to,ipos,irank)
1393  p_to = recv_info_r2p(i_prc_to,irank)
1394  r_to = l_to
1395 
1396  if( io_l ) write(io_fid_log,'(11(I6))') ipos, i_from, j_from, r_from, g_from, l_from, p_from, &
1397  r_to , g_to , l_to , p_to
1398  enddo
1399  enddo
1400 
1401  if( io_l ) write(io_fid_log,*)
1402  if( io_l ) write(io_fid_log,*) '--- Send_list_r2p'
1403  do irank = 1, send_nmax_r2p
1404  if( io_l ) write(io_fid_log,'(11(A6))') 'number', &
1405  '|ifrom','|jfrom','|rfrom','|gfrom','|lfrom','|pfrom', &
1406  '| rto','| gto','| lto','| pto'
1407  do ipos = 1, send_info_r2p(i_size,irank)
1408  g_from = send_list_r2p(i_grid_from,ipos,irank)
1409  l_from = send_list_r2p(i_l_from,ipos,irank)
1410  p_from = send_info_r2p(i_prc_from,irank)
1411  i_from = mod(g_from-1,adm_gall_1d) + 1
1412  j_from = (g_from-i_from) / adm_gall_1d + 1
1413  r_from = prc_rgn_lp2r(l_from,p_from)
1414  g_to = send_list_r2p(i_grid_to,ipos,irank)
1415  l_to = send_list_r2p(i_l_to,ipos,irank)
1416  p_to = send_info_r2p(i_prc_to,irank)
1417  r_to = l_to
1418 
1419  if( io_l ) write(io_fid_log,'(11(I6))') ipos, i_from, j_from, r_from, g_from, l_from, p_from, &
1420  r_to , g_to , l_to , p_to
1421  enddo
1422  enddo
1423  endif
1424  if( io_l ) write(io_fid_log,*)
1425  if( io_l ) write(io_fid_log,*) '*** Send_size_p2r,r2p = ', send_size_nglobal_pl
1426 
1427  allocate( sendbuf_r2p_sp(send_size_nglobal_pl*adm_kall*comm_varmax,send_nmax_r2p) )
1428  allocate( recvbuf_r2p_sp(send_size_nglobal_pl*adm_kall*comm_varmax,recv_nmax_r2p) )
1429  allocate( sendbuf_p2r_sp(send_size_nglobal_pl*adm_kall*comm_varmax,send_nmax_p2r) )
1430  allocate( recvbuf_p2r_sp(send_size_nglobal_pl*adm_kall*comm_varmax,recv_nmax_p2r) )
1431  allocate( sendbuf_r2p_dp(send_size_nglobal_pl*adm_kall*comm_varmax,send_nmax_r2p) )
1432  allocate( recvbuf_r2p_dp(send_size_nglobal_pl*adm_kall*comm_varmax,recv_nmax_r2p) )
1433  allocate( sendbuf_p2r_dp(send_size_nglobal_pl*adm_kall*comm_varmax,send_nmax_p2r) )
1434  allocate( recvbuf_p2r_dp(send_size_nglobal_pl*adm_kall*comm_varmax,recv_nmax_p2r) )
1435 
1436  return
1437  end subroutine comm_sortdest_pl
1438 
1439  !-----------------------------------------------------------------------------
1441  subroutine comm_sortdest_singular
1442  use scale_prc, only: &
1443  prc_myrank
1444  use scale_prc_icoa, only: &
1445  prc_rgn_local, &
1446  prc_rgn_vert_num, &
1447  prc_rgn_vert_pl, &
1448  prc_rgn_lp2r, &
1449  prc_rgn_l2r, &
1450  i_w, &
1451  i_n, &
1452  i_s, &
1453  i_npl, &
1454  i_spl
1455  implicit none
1456 
1457  integer :: rgnid
1458  integer :: i, j, i_rmt, j_rmt
1459 
1460  integer :: l, ipos
1461  integer :: i_from, j_from, r_from, g_from, l_from, p_from
1462  integer :: i_to, j_to, r_to, g_to, l_to, p_to
1463  !---------------------------------------------------------------------------
1464 
1465  allocate( singular_info(info_vindex) )
1466  singular_info(:) = -1
1467  singular_info(i_size) = 0
1468 
1469  allocate( singular_list(list_vindex,4*adm_lall) )
1470  singular_list(:,:) = -1
1471 
1472  do l = 1, prc_rgn_local
1473  rgnid = prc_rgn_l2r(l)
1474 
1475  if ( prc_rgn_vert_num(i_w,rgnid) == 3 ) then
1476  singular_info(i_size) = singular_info(i_size) + 1
1477  ipos = singular_info(i_size)
1478 
1479  i = adm_gmin
1480  j = adm_gmin - 1
1481  i_rmt = adm_gmin - 1
1482  j_rmt = adm_gmin - 1
1483 
1484  singular_list(i_grid_from,ipos) = suf(i,j)
1485  singular_list(i_l_from ,ipos) = l
1486  singular_list(i_grid_to ,ipos) = suf(i_rmt,j_rmt)
1487  singular_list(i_l_to ,ipos) = l
1488  endif
1489 
1490  if ( prc_rgn_vert_num(i_n,rgnid) /= 4 .OR. prc_rgn_vert_pl(i_npl,rgnid) ) then
1491  singular_info(i_size) = singular_info(i_size) + 1
1492  ipos = singular_info(i_size)
1493 
1494  i = adm_gmin
1495  j = adm_gmax + 1
1496  i_rmt = adm_gmin - 1
1497  j_rmt = adm_gmax + 1
1498 
1499  singular_list(i_grid_from,ipos) = suf(i,j)
1500  singular_list(i_l_from ,ipos) = l
1501  singular_list(i_grid_to ,ipos) = suf(i_rmt,j_rmt)
1502  singular_list(i_l_to ,ipos) = l
1503  endif
1504 
1505  if ( prc_rgn_vert_num(i_s,rgnid) /= 4 .OR. prc_rgn_vert_pl(i_spl,rgnid) ) then
1506  singular_info(i_size) = singular_info(i_size) + 1
1507  ipos = singular_info(i_size)
1508 
1509  i = adm_gmax + 1
1510  j = adm_gmin
1511  i_rmt = adm_gmax + 1
1512  j_rmt = adm_gmin - 1
1513 
1514  singular_list(i_grid_from,ipos) = suf(i,j)
1515  singular_list(i_l_from ,ipos) = l
1516  singular_list(i_grid_to ,ipos) = suf(i_rmt,j_rmt)
1517  singular_list(i_l_to ,ipos) = l
1518  endif
1519 
1520  enddo ! loop l
1521 
1522  if ( singular_info(i_size) > 0 ) then
1523  singular_nmax = 1
1524  singular_info(i_prc_from) = prc_myrank
1525  singular_info(i_prc_to ) = prc_myrank
1526  endif
1527 
1528  if( io_l ) write(io_fid_log,*)
1529  if( io_l ) write(io_fid_log,'(A)') '|---------------------------------------'
1530  if( io_l ) write(io_fid_log,'(A)') '| size prc_from prc_to'
1531  if( io_l ) write(io_fid_log,'(A10,3(I10))') '| Singular', singular_info(:)
1532 
1533  if ( debug ) then
1534  if( io_l ) write(io_fid_log,*)
1535  if( io_l ) write(io_fid_log,*) '--- Singular_list'
1536  if( io_l ) write(io_fid_log,*)
1537  if( io_l ) write(io_fid_log,'(13(A6))') 'number', &
1538  '|ifrom','|jfrom','|rfrom','|gfrom','|lfrom','|pfrom', &
1539  '| ito','| jto','| rto','| gto','| lto','| pto'
1540  do ipos = 1, singular_info(i_size)
1541  g_from = singular_list(i_grid_from,ipos)
1542  l_from = singular_list(i_l_from,ipos)
1543  p_from = singular_info(i_prc_from)
1544  i_from = mod(g_from-1,adm_gall_1d) + 1
1545  j_from = (g_from-i_from) / adm_gall_1d + 1
1546  r_from = prc_rgn_lp2r(l_from,p_from)
1547  g_to = singular_list(i_grid_to,ipos)
1548  l_to = singular_list(i_l_to,ipos)
1549  p_to = singular_info(i_prc_to)
1550  i_to = mod(g_to-1,adm_gall_1d) + 1
1551  j_to = (g_to-i_to) / adm_gall_1d + 1
1552  r_to = prc_rgn_lp2r(l_to,p_to)
1553 
1554  if( io_l ) write(io_fid_log,'(13(I6))') ipos, i_from, j_from, r_from, g_from, l_from, p_from, &
1555  i_to , j_to , r_to , g_to , l_to , p_to
1556  enddo
1557  endif
1558 
1559  return
1560  end subroutine comm_sortdest_singular
1561 
1562  !-----------------------------------------------------------------------------
1564  subroutine comm_data_transfer_sp( &
1565  var, &
1566  var_pl )
1567  use scale_prc, only: &
1569  prc_abort
1570  implicit none
1571 
1572  real(SP), intent(inout) :: var (:,:,:,:)
1573  real(SP), intent(inout) :: var_pl(:,:,:,:)
1574 
1575  integer :: shp(4), kmax, vmax
1576  integer :: totalsize, rank, tag
1577  integer :: irank, ipos, imax
1578  integer :: ij_from, l_from, ij_to, l_to
1579 
1580  integer :: k, v, ikv
1581  integer :: ierr
1582  !---------------------------------------------------------------------------
1583  !$acc data &
1584  !$acc pcopy(var) &
1585  !$acc pcopy(var_pl)
1586 
1587  if ( comm_apply_barrier ) then
1588  call prof_rapstart('COMM_barrier',2)
1589  call mpi_barrier(prc_local_comm_world,ierr)
1590  call prof_rapend ('COMM_barrier',2)
1591  endif
1592 
1593  !$acc wait
1594 
1595  call prof_rapstart('COMM_data_transfer',2)
1596 
1597  shp = shape(var)
1598  kmax = shp(2)
1599  vmax = shp(4)
1600 
1601  if ( kmax * vmax > adm_kall * comm_varmax ) then
1602  write(*,*) 'xxx [COMM_data_transfer] kmax * vmax exceeds ADM_kall * COMM_varmax, stop!'
1603  write(*,*) 'xxx kmax * vmax = ', kmax * vmax
1604  write(*,*) 'xxx ADM_kall * COMM_varmax = ', adm_kall * comm_varmax
1605  call prc_abort
1606  endif
1607 
1608  !---< start communication >---
1609  ! Theres no p2r & r2p communication without calling COMM_sortdest_pl.
1610  ! receive pole => region
1611  ! receive region => pole
1612  ! receive region => region
1613  ! pack and send pole => region
1614  ! pack and send region => pole
1615  ! pack and send region => region
1616  ! copy pole => region
1617  ! copy region => pole
1618  ! copy region => region
1619  ! wait all
1620  ! unpack pole => region
1621  ! unpack region => pole
1622  ! unpack region => region
1623  ! copy region halo => region halo (singular point)
1624 
1625  req_count = 0
1626 
1627  !--- receive r2r
1628  do irank = 1, recv_nmax_r2r
1629  req_count = req_count + 1
1630  totalsize = recv_info_r2r(i_size ,irank) * kmax * vmax
1631  rank = recv_info_r2r(i_prc_from,irank)
1632  tag = recv_info_r2r(i_prc_from,irank)
1633 
1634  call mpi_irecv( recvbuf_r2r_sp(1,irank), & ! [OUT]
1635  totalsize, & ! [IN]
1636  mpi_real, & ! [IN]
1637  rank, & ! [IN]
1638  tag, & ! [IN]
1639  prc_local_comm_world, & ! [IN]
1640  req_list(req_count), & ! [OUT]
1641  ierr ) ! [OUT]
1642  enddo
1643 
1644  !--- receive p2r
1645  do irank = 1, recv_nmax_p2r
1646  req_count = req_count + 1
1647  totalsize = recv_info_p2r(i_size ,irank) * kmax * vmax
1648  rank = recv_info_p2r(i_prc_from,irank)
1649  tag = recv_info_p2r(i_prc_from,irank) + 100000
1650 
1651  call mpi_irecv( recvbuf_p2r_sp(1,irank), & ! [OUT]
1652  totalsize, & ! [IN]
1653  mpi_real, & ! [IN]
1654  rank, & ! [IN]
1655  tag, & ! [IN]
1656  prc_local_comm_world, & ! [IN]
1657  req_list(req_count), & ! [OUT]
1658  ierr ) ! [OUT]
1659  enddo
1660 
1661  !--- receive r2p
1662  do irank = 1, recv_nmax_r2p
1663  req_count = req_count + 1
1664  totalsize = recv_info_r2p(i_size ,irank) * kmax * vmax
1665  rank = recv_info_r2p(i_prc_from,irank)
1666  tag = recv_info_r2p(i_prc_from,irank) + 200000
1667 
1668  call mpi_irecv( recvbuf_r2p_sp(1,irank), & ! [OUT]
1669  totalsize, & ! [IN]
1670  mpi_real, & ! [IN]
1671  rank, & ! [IN]
1672  tag, & ! [IN]
1673  prc_local_comm_world, & ! [IN]
1674  req_list(req_count), & ! [OUT]
1675  ierr ) ! [OUT]
1676  enddo
1677 
1678  !$acc data &
1679  !$acc pcopyin(Send_list_r2r,Send_list_p2r,Send_list_r2p)
1680  !$acc wait
1681 
1682  !--- pack and send r2r
1683  do irank = 1, send_nmax_r2r
1684  imax = send_info_r2r(i_size,irank)
1685 
1686  !$acc kernels pcopy(sendbuf_r2r_SP)
1687  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
1688  !$omp shared(irank,imax,kmax,vmax,Send_list_r2r,sendbuf_r2r_SP,var) &
1689  !$omp collapse(2)
1690  !$acc loop independent
1691  do v = 1, vmax
1692  !$acc loop independent
1693  do k = 1, kmax
1694  !$acc loop independent
1695  do ipos = 1, imax
1696  ij_from = send_list_r2r(i_grid_from,ipos,irank)
1697  l_from = send_list_r2r(i_l_from ,ipos,irank)
1698 
1699  ikv = (v-1) * imax * kmax &
1700  + (k-1) * imax &
1701  + ipos
1702 
1703  sendbuf_r2r_sp(ikv,irank) = var(ij_from,k,l_from,v)
1704  enddo
1705  enddo
1706  enddo
1707  !$omp end parallel do
1708  !$acc end kernels
1709 
1710  req_count = req_count + 1
1711  totalsize = imax * kmax * vmax
1712  rank = send_info_r2r(i_prc_to ,irank)
1713  tag = send_info_r2r(i_prc_from,irank)
1714 
1715  !$acc wait
1716  call mpi_isend( sendbuf_r2r_sp(1,irank), & ! [IN]
1717  totalsize, & ! [IN]
1718  mpi_real, & ! [IN]
1719  rank, & ! [IN]
1720  tag, & ! [IN]
1721  prc_local_comm_world, & ! [IN]
1722  req_list(req_count), & ! [OUT]
1723  ierr ) ! [OUT]
1724  enddo
1725 
1726  !--- pack and send p2r
1727  do irank = 1, send_nmax_p2r
1728  imax = send_info_p2r(i_size,irank)
1729 
1730  !$acc kernels pcopy(sendbuf_p2r_SP)
1731  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
1732  !$omp shared(irank,imax,kmax,vmax,Send_list_p2r,sendbuf_p2r_SP,var_pl) &
1733  !$omp collapse(2)
1734  !$acc loop independent
1735  do v = 1, vmax
1736  !$acc loop independent
1737  do k = 1, kmax
1738  !$acc loop independent
1739  do ipos = 1, imax
1740  ij_from = send_list_p2r(i_grid_from,ipos,irank)
1741  l_from = send_list_p2r(i_l_from ,ipos,irank)
1742 
1743  ikv = (v-1) * imax * kmax &
1744  + (k-1) * imax &
1745  + ipos
1746 
1747  sendbuf_p2r_sp(ikv,irank) = var_pl(ij_from,k,l_from,v)
1748  enddo
1749  enddo
1750  enddo
1751  !$omp end parallel do
1752  !$acc end kernels
1753 
1754  req_count = req_count + 1
1755  totalsize = imax * kmax * vmax
1756  rank = send_info_p2r(i_prc_to ,irank)
1757  tag = send_info_p2r(i_prc_from,irank) + 100000
1758 
1759  !$acc wait
1760  call mpi_isend( sendbuf_p2r_sp(1,irank), & ! [IN]
1761  totalsize, & ! [IN]
1762  mpi_real, & ! [IN]
1763  rank, & ! [IN]
1764  tag, & ! [IN]
1765  prc_local_comm_world, & ! [IN]
1766  req_list(req_count), & ! [OUT]
1767  ierr ) ! [OUT]
1768  enddo
1769 
1770  !--- pack and send r2p
1771  do irank = 1, send_nmax_r2p
1772  imax = send_info_r2p(i_size,irank)
1773 
1774  !$acc kernels pcopy(sendbuf_r2p_SP)
1775  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
1776  !$omp shared(irank,imax,kmax,vmax,Send_list_r2p,sendbuf_r2p_SP,var) &
1777  !$omp collapse(2)
1778  !$acc loop independent
1779  do v = 1, vmax
1780  !$acc loop independent
1781  do k = 1, kmax
1782  !$acc loop independent
1783  do ipos = 1, imax
1784  ij_from = send_list_r2p(i_grid_from,ipos,irank)
1785  l_from = send_list_r2p(i_l_from ,ipos,irank)
1786 
1787  ikv = (v-1) * imax * kmax &
1788  + (k-1) * imax &
1789  + ipos
1790 
1791  sendbuf_r2p_sp(ikv,irank) = var(ij_from,k,l_from,v)
1792  enddo
1793  enddo
1794  enddo
1795  !$omp end parallel do
1796  !$acc end kernels
1797 
1798  req_count = req_count + 1
1799  totalsize = imax * kmax * vmax
1800  rank = send_info_r2p(i_prc_to ,irank)
1801  tag = send_info_r2p(i_prc_from,irank) + 200000
1802 
1803  !$acc wait
1804  call mpi_isend( sendbuf_r2p_sp(1,irank), & ! [IN]
1805  totalsize, & ! [IN]
1806  mpi_real, & ! [IN]
1807  rank, & ! [IN]
1808  tag, & ! [IN]
1809  prc_local_comm_world, & ! [IN]
1810  req_list(req_count), & ! [OUT]
1811  ierr ) ! [OUT]
1812  enddo
1813 
1814  !$acc end data
1815  !$acc data &
1816  !$acc pcopyin(Copy_list_r2r,Copy_list_p2r,Copy_list_r2p)
1817  !$acc wait
1818 
1819  !$omp parallel default(none),private(ipos,k,v,imax,irank,ij_from,l_from,ij_to,l_to) &
1820  !$omp shared(kmax,vmax,var,var_pl, &
1821  !$omp Copy_nmax_p2r,Copy_info_p2r,Copy_list_p2r, &
1822  !$omp Copy_nmax_r2p,Copy_info_r2p,Copy_list_r2p, &
1823  !$omp Copy_nmax_r2r,Copy_info_r2r,Copy_list_r2r )
1824 
1825  !--- copy r2r
1826  do irank = 1, copy_nmax_r2r
1827  imax = copy_info_r2r(i_size)
1828 
1829  !$acc kernels
1830  !$omp do collapse(2)
1831  !$acc loop independent
1832  do v = 1, vmax
1833  !$acc loop independent
1834  do k = 1, kmax
1835  !$acc loop independent
1836  do ipos = 1, imax
1837  ij_from = copy_list_r2r(i_grid_from,ipos)
1838  l_from = copy_list_r2r(i_l_from ,ipos)
1839  ij_to = copy_list_r2r(i_grid_to ,ipos)
1840  l_to = copy_list_r2r(i_l_to ,ipos)
1841 
1842  var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
1843  enddo
1844  enddo
1845  enddo
1846  !$omp end do nowait
1847  !$acc end kernels
1848  enddo
1849 
1850  !--- copy p2r
1851  do irank = 1, copy_nmax_p2r
1852  imax = copy_info_p2r(i_size)
1853 
1854  !$acc kernels
1855  !$omp do collapse(2)
1856  !$acc loop independent
1857  do v = 1, vmax
1858  !$acc loop independent
1859  do k = 1, kmax
1860  !$acc loop independent
1861  do ipos = 1, imax
1862  ij_from = copy_list_p2r(i_grid_from,ipos)
1863  l_from = copy_list_p2r(i_l_from ,ipos)
1864  ij_to = copy_list_p2r(i_grid_to ,ipos)
1865  l_to = copy_list_p2r(i_l_to ,ipos)
1866 
1867  var(ij_to,k,l_to,v) = var_pl(ij_from,k,l_from,v)
1868  enddo
1869  enddo
1870  enddo
1871  !$omp end do nowait
1872  !$acc end kernels
1873  enddo
1874 
1875  !--- copy r2p
1876  do irank = 1, copy_nmax_r2p
1877  imax = copy_info_r2p(i_size)
1878 
1879  !$acc kernels
1880  !$omp do collapse(2)
1881  !$acc loop independent
1882  do v = 1, vmax
1883  !$acc loop independent
1884  do k = 1, kmax
1885  !$acc loop independent
1886  do ipos = 1, imax
1887  ij_from = copy_list_r2p(i_grid_from,ipos)
1888  l_from = copy_list_r2p(i_l_from ,ipos)
1889  ij_to = copy_list_r2p(i_grid_to ,ipos)
1890  l_to = copy_list_r2p(i_l_to ,ipos)
1891 
1892  var_pl(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
1893  enddo
1894  enddo
1895  enddo
1896  !$omp end do nowait
1897  !$acc end kernels
1898  enddo
1899 
1900  !$omp end parallel
1901 
1902  !--- wait all
1903  if ( req_count > 0 ) then
1904  call mpi_waitall( req_count, & ! [IN]
1905  req_list(:), & ! [IN]
1906  mpi_statuses_ignore, & ! [OUT]
1907  ierr ) ! [OUT]
1908  endif
1909 
1910  !$acc end data
1911  !$acc data &
1912  !$acc pcopyin(Recv_list_r2r,Recv_list_p2r,Recv_list_r2p)
1913  !$acc wait
1914 
1915  !$omp parallel default(none),private(ipos,k,v,imax,irank,ikv,ij_from,l_from,ij_to,l_to) &
1916  !$omp shared(kmax,vmax,var,var_pl, &
1917  !$omp Recv_nmax_p2r,Recv_info_p2r,Recv_list_p2r,recvbuf_p2r_SP, &
1918  !$omp Recv_nmax_r2p,Recv_info_r2p,Recv_list_r2p,recvbuf_r2p_SP, &
1919  !$omp Recv_nmax_r2r,Recv_info_r2r,Recv_list_r2r,recvbuf_r2r_SP, &
1920  !$omp Singular_nmax,Singular_info,Singular_list )
1921 
1922  !--- unpack r2r
1923  do irank = 1, recv_nmax_r2r
1924  imax = recv_info_r2r(i_size,irank)
1925 
1926  !$acc kernels pcopyin(recvbuf_r2r_SP)
1927  !$omp do collapse(2)
1928  !$acc loop independent
1929  do v = 1, vmax
1930  !$acc loop independent
1931  do k = 1, kmax
1932  !$acc loop independent
1933  do ipos = 1, imax
1934  ij_to = recv_list_r2r(i_grid_to,ipos,irank)
1935  l_to = recv_list_r2r(i_l_to ,ipos,irank)
1936 
1937  ikv = (v-1) * imax * kmax &
1938  + (k-1) * imax &
1939  + ipos
1940 
1941  var(ij_to,k,l_to,v) = recvbuf_r2r_sp(ikv,irank)
1942  enddo
1943  enddo
1944  enddo
1945  !$omp end do nowait
1946  !$acc end kernels
1947  enddo
1948 
1949  !--- unpack p2r
1950  do irank = 1, recv_nmax_p2r
1951  imax = recv_info_p2r(i_size,irank)
1952 
1953  !$acc kernels pcopyin(recvbuf_p2r_SP)
1954  !$omp do collapse(2)
1955  !$acc loop independent
1956  do v = 1, vmax
1957  !$acc loop independent
1958  do k = 1, kmax
1959  !$acc loop independent
1960  do ipos = 1, imax
1961  ij_to = recv_list_p2r(i_grid_to,ipos,irank)
1962  l_to = recv_list_p2r(i_l_to ,ipos,irank)
1963 
1964  ikv = (v-1) * imax * kmax &
1965  + (k-1) * imax &
1966  + ipos
1967 
1968  var(ij_to,k,l_to,v) = recvbuf_p2r_sp(ikv,irank)
1969  enddo
1970  enddo
1971  enddo
1972  !$omp end do nowait
1973  !$acc end kernels
1974  enddo
1975 
1976  !--- unpack r2p
1977  do irank = 1, recv_nmax_r2p
1978  imax = recv_info_r2p(i_size,irank)
1979 
1980  !$acc kernels pcopyin(recvbuf_r2p_SP)
1981  !$omp do collapse(2)
1982  !$acc loop independent
1983  do v = 1, vmax
1984  !$acc loop independent
1985  do k = 1, kmax
1986  !$acc loop independent
1987  do ipos = 1, imax
1988  ij_to = recv_list_r2p(i_grid_to,ipos,irank)
1989  l_to = recv_list_r2p(i_l_to ,ipos,irank)
1990 
1991  ikv = (v-1) * imax * kmax &
1992  + (k-1) * imax &
1993  + ipos
1994 
1995  var_pl(ij_to,k,l_to,v) = recvbuf_r2p_sp(ikv,irank)
1996  enddo
1997  enddo
1998  enddo
1999  !$omp end do nowait
2000  !$acc end kernels
2001  enddo
2002 
2003  !$acc end data
2004  !$acc wait
2005 
2006  !--- singular point (halo to halo)
2007  do irank = 1, singular_nmax
2008  imax = singular_info(i_size)
2009 
2010  !$acc kernels pcopyin(Singular_list)
2011  !$omp do collapse(2)
2012  !$acc loop independent
2013  do v = 1, vmax
2014  !$acc loop independent
2015  do k = 1, kmax
2016  !$acc loop independent
2017  do ipos = 1, imax
2018  ij_from = singular_list(i_grid_from,ipos)
2019  l_from = singular_list(i_l_from ,ipos)
2020  ij_to = singular_list(i_grid_to ,ipos)
2021  l_to = singular_list(i_l_to ,ipos)
2022 
2023  var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2024  enddo
2025  enddo
2026  enddo
2027  !$omp end do nowait
2028  !$acc end kernels
2029  enddo
2030 
2031  !$omp end parallel
2032 
2033  call prof_rapend('COMM_data_transfer',2)
2034 
2035  !$acc end data
2036  !$acc wait
2037 
2038  return
2039  end subroutine comm_data_transfer_sp
2040 
2041  !-----------------------------------------------------------------------------
2043  subroutine comm_data_transfer_dp( &
2044  var, &
2045  var_pl )
2046  use scale_prc, only: &
2048  prc_abort
2049  implicit none
2050 
2051  real(DP), intent(inout) :: var (:,:,:,:)
2052  real(DP), intent(inout) :: var_pl(:,:,:,:)
2053 
2054  integer :: shp(4), kmax, vmax
2055  integer :: totalsize, rank, tag
2056  integer :: irank, ipos, imax
2057  integer :: ij_from, l_from, ij_to, l_to
2058 
2059  integer :: k, v, ikv
2060  integer :: ierr
2061  !---------------------------------------------------------------------------
2062  !$acc data &
2063  !$acc pcopy(var) &
2064  !$acc pcopy(var_pl)
2065 
2066  if ( comm_apply_barrier ) then
2067  call prof_rapstart('COMM_barrier',2)
2068  call mpi_barrier(prc_local_comm_world,ierr)
2069  call prof_rapend ('COMM_barrier',2)
2070  endif
2071 
2072  !$acc wait
2073 
2074  call prof_rapstart('COMM_data_transfer',2)
2075 
2076  shp = shape(var)
2077  kmax = shp(2)
2078  vmax = shp(4)
2079 
2080  if ( kmax * vmax > adm_kall * comm_varmax ) then
2081  write(*,*) 'xxx [COMM_data_transfer] kmax * vmax exceeds ADM_kall * COMM_varmax, stop!'
2082  write(*,*) 'xxx kmax * vmax = ', kmax * vmax
2083  write(*,*) 'xxx ADM_kall * COMM_varmax = ', adm_kall * comm_varmax
2084  call prc_abort
2085  endif
2086 
2087  !---< start communication >---
2088  ! Theres no p2r & r2p communication without calling COMM_sortdest_pl.
2089  ! receive pole => region
2090  ! receive region => pole
2091  ! receive region => region
2092  ! pack and send pole => region
2093  ! pack and send region => pole
2094  ! pack and send region => region
2095  ! copy pole => region
2096  ! copy region => pole
2097  ! copy region => region
2098  ! wait all
2099  ! unpack pole => region
2100  ! unpack region => pole
2101  ! unpack region => region
2102  ! copy region halo => region halo (singular point)
2103 
2104  req_count = 0
2105 
2106  !--- receive r2r
2107  do irank = 1, recv_nmax_r2r
2108  req_count = req_count + 1
2109  totalsize = recv_info_r2r(i_size ,irank) * kmax * vmax
2110  rank = recv_info_r2r(i_prc_from,irank)
2111  tag = recv_info_r2r(i_prc_from,irank)
2112 
2113  call mpi_irecv( recvbuf_r2r_dp(1,irank), & ! [OUT]
2114  totalsize, & ! [IN]
2115  mpi_double_precision, & ! [IN]
2116  rank, & ! [IN]
2117  tag, & ! [IN]
2118  prc_local_comm_world, & ! [IN]
2119  req_list(req_count), & ! [OUT]
2120  ierr ) ! [OUT]
2121  enddo
2122 
2123  !--- receive p2r
2124  do irank = 1, recv_nmax_p2r
2125  req_count = req_count + 1
2126  totalsize = recv_info_p2r(i_size ,irank) * kmax * vmax
2127  rank = recv_info_p2r(i_prc_from,irank)
2128  tag = recv_info_p2r(i_prc_from,irank) + 100000
2129 
2130  call mpi_irecv( recvbuf_p2r_dp(1,irank), & ! [OUT]
2131  totalsize, & ! [IN]
2132  mpi_double_precision, & ! [IN]
2133  rank, & ! [IN]
2134  tag, & ! [IN]
2135  prc_local_comm_world, & ! [IN]
2136  req_list(req_count), & ! [OUT]
2137  ierr ) ! [OUT]
2138  enddo
2139 
2140  !--- receive r2p
2141  do irank = 1, recv_nmax_r2p
2142  req_count = req_count + 1
2143  totalsize = recv_info_r2p(i_size ,irank) * kmax * vmax
2144  rank = recv_info_r2p(i_prc_from,irank)
2145  tag = recv_info_r2p(i_prc_from,irank) + 200000
2146 
2147  call mpi_irecv( recvbuf_r2p_dp(1,irank), & ! [OUT]
2148  totalsize, & ! [IN]
2149  mpi_double_precision, & ! [IN]
2150  rank, & ! [IN]
2151  tag, & ! [IN]
2152  prc_local_comm_world, & ! [IN]
2153  req_list(req_count), & ! [OUT]
2154  ierr ) ! [OUT]
2155  enddo
2156 
2157  !$acc data &
2158  !$acc pcopyin(Send_list_r2r,Send_list_p2r,Send_list_r2p)
2159  !$acc wait
2160 
2161  !--- pack and send r2r
2162  do irank = 1, send_nmax_r2r
2163  imax = send_info_r2r(i_size,irank)
2164 
2165  !$acc kernels pcopy(sendbuf_r2r_DP)
2166  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
2167  !$omp shared(irank,imax,kmax,vmax,Send_list_r2r,sendbuf_r2r_DP,var) &
2168  !$omp collapse(2)
2169  !$acc loop independent
2170  do v = 1, vmax
2171  !$acc loop independent
2172  do k = 1, kmax
2173  !$acc loop independent
2174  do ipos = 1, imax
2175  ij_from = send_list_r2r(i_grid_from,ipos,irank)
2176  l_from = send_list_r2r(i_l_from ,ipos,irank)
2177 
2178  ikv = (v-1) * imax * kmax &
2179  + (k-1) * imax &
2180  + ipos
2181 
2182  sendbuf_r2r_dp(ikv,irank) = var(ij_from,k,l_from,v)
2183  enddo
2184  enddo
2185  enddo
2186  !$omp end parallel do
2187  !$acc end kernels
2188 
2189  req_count = req_count + 1
2190  totalsize = imax * kmax * vmax
2191  rank = send_info_r2r(i_prc_to ,irank)
2192  tag = send_info_r2r(i_prc_from,irank)
2193 
2194  !$acc wait
2195  call mpi_isend( sendbuf_r2r_dp(1,irank), & ! [IN]
2196  totalsize, & ! [IN]
2197  mpi_double_precision, & ! [IN]
2198  rank, & ! [IN]
2199  tag, & ! [IN]
2200  prc_local_comm_world, & ! [IN]
2201  req_list(req_count), & ! [OUT]
2202  ierr ) ! [OUT]
2203  enddo
2204 
2205  !--- pack and send p2r
2206  do irank = 1, send_nmax_p2r
2207  imax = send_info_p2r(i_size,irank)
2208 
2209  !$acc kernels pcopy(sendbuf_p2r_DP)
2210  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
2211  !$omp shared(irank,imax,kmax,vmax,Send_list_p2r,sendbuf_p2r_DP,var_pl) &
2212  !$omp collapse(2)
2213  !$acc loop independent
2214  do v = 1, vmax
2215  !$acc loop independent
2216  do k = 1, kmax
2217  !$acc loop independent
2218  do ipos = 1, imax
2219  ij_from = send_list_p2r(i_grid_from,ipos,irank)
2220  l_from = send_list_p2r(i_l_from ,ipos,irank)
2221 
2222  ikv = (v-1) * imax * kmax &
2223  + (k-1) * imax &
2224  + ipos
2225 
2226  sendbuf_p2r_dp(ikv,irank) = var_pl(ij_from,k,l_from,v)
2227  enddo
2228  enddo
2229  enddo
2230  !$omp end parallel do
2231  !$acc end kernels
2232 
2233  req_count = req_count + 1
2234  totalsize = imax * kmax * vmax
2235  rank = send_info_p2r(i_prc_to ,irank)
2236  tag = send_info_p2r(i_prc_from,irank) + 100000
2237 
2238  !$acc wait
2239  call mpi_isend( sendbuf_p2r_dp(1,irank), & ! [IN]
2240  totalsize, & ! [IN]
2241  mpi_double_precision, & ! [IN]
2242  rank, & ! [IN]
2243  tag, & ! [IN]
2244  prc_local_comm_world, & ! [IN]
2245  req_list(req_count), & ! [OUT]
2246  ierr ) ! [OUT]
2247  enddo
2248 
2249  !--- pack and send r2p
2250  do irank = 1, send_nmax_r2p
2251  imax = send_info_r2p(i_size,irank)
2252 
2253  !$acc kernels pcopy(sendbuf_r2p_DP)
2254  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
2255  !$omp shared(irank,imax,kmax,vmax,Send_list_r2p,sendbuf_r2p_DP,var) &
2256  !$omp collapse(2)
2257  !$acc loop independent
2258  do v = 1, vmax
2259  !$acc loop independent
2260  do k = 1, kmax
2261  !$acc loop independent
2262  do ipos = 1, imax
2263  ij_from = send_list_r2p(i_grid_from,ipos,irank)
2264  l_from = send_list_r2p(i_l_from ,ipos,irank)
2265 
2266  ikv = (v-1) * imax * kmax &
2267  + (k-1) * imax &
2268  + ipos
2269 
2270  sendbuf_r2p_dp(ikv,irank) = var(ij_from,k,l_from,v)
2271  enddo
2272  enddo
2273  enddo
2274  !$omp end parallel do
2275  !$acc end kernels
2276 
2277  req_count = req_count + 1
2278  totalsize = imax * kmax * vmax
2279  rank = send_info_r2p(i_prc_to ,irank)
2280  tag = send_info_r2p(i_prc_from,irank) + 200000
2281 
2282  !$acc wait
2283  call mpi_isend( sendbuf_r2p_dp(1,irank), & ! [IN]
2284  totalsize, & ! [IN]
2285  mpi_double_precision, & ! [IN]
2286  rank, & ! [IN]
2287  tag, & ! [IN]
2288  prc_local_comm_world, & ! [IN]
2289  req_list(req_count), & ! [OUT]
2290  ierr ) ! [OUT]
2291  enddo
2292 
2293  !$acc end data
2294  !$acc data &
2295  !$acc pcopyin(Copy_list_r2r,Copy_list_p2r,Copy_list_r2p)
2296  !$acc wait
2297 
2298  !$omp parallel default(none),private(ipos,k,v,imax,irank,ij_from,l_from,ij_to,l_to) &
2299  !$omp shared(kmax,vmax,var,var_pl, &
2300  !$omp Copy_nmax_p2r,Copy_info_p2r,Copy_list_p2r, &
2301  !$omp Copy_nmax_r2p,Copy_info_r2p,Copy_list_r2p, &
2302  !$omp Copy_nmax_r2r,Copy_info_r2r,Copy_list_r2r )
2303 
2304  !--- copy r2r
2305  do irank = 1, copy_nmax_r2r
2306  imax = copy_info_r2r(i_size)
2307 
2308  !$acc kernels
2309  !$omp do collapse(2)
2310  !$acc loop independent
2311  do v = 1, vmax
2312  !$acc loop independent
2313  do k = 1, kmax
2314  !$acc loop independent
2315  do ipos = 1, imax
2316  ij_from = copy_list_r2r(i_grid_from,ipos)
2317  l_from = copy_list_r2r(i_l_from ,ipos)
2318  ij_to = copy_list_r2r(i_grid_to ,ipos)
2319  l_to = copy_list_r2r(i_l_to ,ipos)
2320 
2321  var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2322  enddo
2323  enddo
2324  enddo
2325  !$omp end do nowait
2326  !$acc end kernels
2327  enddo
2328 
2329  !--- copy p2r
2330  do irank = 1, copy_nmax_p2r
2331  imax = copy_info_p2r(i_size)
2332 
2333  !$acc kernels
2334  !$omp do collapse(2)
2335  !$acc loop independent
2336  do v = 1, vmax
2337  !$acc loop independent
2338  do k = 1, kmax
2339  !$acc loop independent
2340  do ipos = 1, imax
2341  ij_from = copy_list_p2r(i_grid_from,ipos)
2342  l_from = copy_list_p2r(i_l_from ,ipos)
2343  ij_to = copy_list_p2r(i_grid_to ,ipos)
2344  l_to = copy_list_p2r(i_l_to ,ipos)
2345 
2346  var(ij_to,k,l_to,v) = var_pl(ij_from,k,l_from,v)
2347  enddo
2348  enddo
2349  enddo
2350  !$omp end do nowait
2351  !$acc end kernels
2352  enddo
2353 
2354  !--- copy r2p
2355  do irank = 1, copy_nmax_r2p
2356  imax = copy_info_r2p(i_size)
2357 
2358  !$acc kernels
2359  !$omp do collapse(2)
2360  !$acc loop independent
2361  do v = 1, vmax
2362  !$acc loop independent
2363  do k = 1, kmax
2364  !$acc loop independent
2365  do ipos = 1, imax
2366  ij_from = copy_list_r2p(i_grid_from,ipos)
2367  l_from = copy_list_r2p(i_l_from ,ipos)
2368  ij_to = copy_list_r2p(i_grid_to ,ipos)
2369  l_to = copy_list_r2p(i_l_to ,ipos)
2370 
2371  var_pl(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2372  enddo
2373  enddo
2374  enddo
2375  !$omp end do nowait
2376  !$acc end kernels
2377  enddo
2378 
2379  !$omp end parallel
2380 
2381  !--- wait all
2382  if ( req_count > 0 ) then
2383  call mpi_waitall( req_count, & ! [IN]
2384  req_list(:), & ! [IN]
2385  mpi_statuses_ignore, & ! [OUT]
2386  ierr ) ! [OUT]
2387  endif
2388 
2389  !$acc end data
2390  !$acc data &
2391  !$acc pcopyin(Recv_list_r2r,Recv_list_p2r,Recv_list_r2p)
2392  !$acc wait
2393 
2394  !$omp parallel default(none),private(ipos,k,v,imax,irank,ikv,ij_from,l_from,ij_to,l_to) &
2395  !$omp shared(kmax,vmax,var,var_pl, &
2396  !$omp Recv_nmax_p2r,Recv_info_p2r,Recv_list_p2r,recvbuf_p2r_DP, &
2397  !$omp Recv_nmax_r2p,Recv_info_r2p,Recv_list_r2p,recvbuf_r2p_DP, &
2398  !$omp Recv_nmax_r2r,Recv_info_r2r,Recv_list_r2r,recvbuf_r2r_DP, &
2399  !$omp Singular_nmax,Singular_info,Singular_list )
2400 
2401  !--- unpack r2r
2402  do irank = 1, recv_nmax_r2r
2403  imax = recv_info_r2r(i_size,irank)
2404 
2405  !$acc kernels pcopyin(recvbuf_r2r_DP)
2406  !$omp do collapse(2)
2407  !$acc loop independent
2408  do v = 1, vmax
2409  !$acc loop independent
2410  do k = 1, kmax
2411  !$acc loop independent
2412  do ipos = 1, imax
2413  ij_to = recv_list_r2r(i_grid_to,ipos,irank)
2414  l_to = recv_list_r2r(i_l_to ,ipos,irank)
2415 
2416  ikv = (v-1) * imax * kmax &
2417  + (k-1) * imax &
2418  + ipos
2419 
2420  var(ij_to,k,l_to,v) = recvbuf_r2r_dp(ikv,irank)
2421  enddo
2422  enddo
2423  enddo
2424  !$omp end do nowait
2425  !$acc end kernels
2426  enddo
2427 
2428  !--- unpack p2r
2429  do irank = 1, recv_nmax_p2r
2430  imax = recv_info_p2r(i_size,irank)
2431 
2432  !$acc kernels pcopyin(recvbuf_p2r_DP)
2433  !$omp do collapse(2)
2434  !$acc loop independent
2435  do v = 1, vmax
2436  !$acc loop independent
2437  do k = 1, kmax
2438  !$acc loop independent
2439  do ipos = 1, imax
2440  ij_to = recv_list_p2r(i_grid_to,ipos,irank)
2441  l_to = recv_list_p2r(i_l_to ,ipos,irank)
2442 
2443  ikv = (v-1) * imax * kmax &
2444  + (k-1) * imax &
2445  + ipos
2446 
2447  var(ij_to,k,l_to,v) = recvbuf_p2r_dp(ikv,irank)
2448  enddo
2449  enddo
2450  enddo
2451  !$omp end do nowait
2452  !$acc end kernels
2453  enddo
2454 
2455  !--- unpack r2p
2456  do irank = 1, recv_nmax_r2p
2457  imax = recv_info_r2p(i_size,irank)
2458 
2459  !$acc kernels pcopyin(recvbuf_r2p_DP)
2460  !$omp do collapse(2)
2461  !$acc loop independent
2462  do v = 1, vmax
2463  !$acc loop independent
2464  do k = 1, kmax
2465  !$acc loop independent
2466  do ipos = 1, imax
2467  ij_to = recv_list_r2p(i_grid_to,ipos,irank)
2468  l_to = recv_list_r2p(i_l_to ,ipos,irank)
2469 
2470  ikv = (v-1) * imax * kmax &
2471  + (k-1) * imax &
2472  + ipos
2473 
2474  var_pl(ij_to,k,l_to,v) = recvbuf_r2p_dp(ikv,irank)
2475  enddo
2476  enddo
2477  enddo
2478  !$omp end do nowait
2479  !$acc end kernels
2480  enddo
2481 
2482  !$acc end data
2483  !$acc wait
2484 
2485  !--- singular point (halo to halo)
2486  do irank = 1, singular_nmax
2487  imax = singular_info(i_size)
2488 
2489  !$acc kernels pcopyin(Singular_list)
2490  !$omp do collapse(2)
2491  !$acc loop independent
2492  do v = 1, vmax
2493  !$acc loop independent
2494  do k = 1, kmax
2495  !$acc loop independent
2496  do ipos = 1, imax
2497  ij_from = singular_list(i_grid_from,ipos)
2498  l_from = singular_list(i_l_from ,ipos)
2499  ij_to = singular_list(i_grid_to ,ipos)
2500  l_to = singular_list(i_l_to ,ipos)
2501 
2502  var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2503  enddo
2504  enddo
2505  enddo
2506  !$omp end do nowait
2507  !$acc end kernels
2508  enddo
2509 
2510  !$omp end parallel
2511 
2512  call prof_rapend('COMM_data_transfer',2)
2513 
2514  !$acc end data
2515  !$acc wait
2516 
2517  return
2518  end subroutine comm_data_transfer_dp
2519 
2520  !-----------------------------------------------------------------------------
2522  subroutine comm_data_transfer_nopl( &
2523  var )
2524  use scale_prc, only: &
2526  prc_abort
2527  implicit none
2528 
2529  real(dp), intent(inout) :: var (:,:,:,:)
2530 
2531  integer :: shp(4), kmax, vmax
2532  integer :: totalsize, rank, tag
2533  integer :: irank, ipos, imax
2534  integer :: ij_from, l_from, ij_to, l_to
2535 
2536  integer :: k, v, ikv
2537  integer :: ierr
2538  !---------------------------------------------------------------------------
2539 
2540  if ( comm_apply_barrier ) then
2541  call prof_rapstart('COMM_barrier',2)
2542  call mpi_barrier(prc_local_comm_world,ierr)
2543  call prof_rapend ('COMM_barrier',2)
2544  endif
2545 
2546  call prof_rapstart('COMM_data_transfer',2)
2547 
2548  shp = shape(var)
2549  kmax = shp(2)
2550  vmax = shp(4)
2551 
2552  if ( kmax * vmax > adm_kall * comm_varmax ) then
2553  write(*,*) 'xxx [COMM_data_transfer_nopl] kmax * vmax exceeds ADM_kall * COMM_varmax, stop!'
2554  write(*,*) 'xxx kmax * vmax = ', kmax * vmax
2555  write(*,*) 'xxx ADM_kall * COMM_varmax = ', adm_kall * comm_varmax
2556  call prc_abort
2557  endif
2558 
2559  !---< start communication >---
2560  ! receive region => region
2561  ! pack and send region => region
2562  ! copy region => region
2563  ! wait all
2564  ! unpack region => region
2565  ! copy region halo => region halo (singular point)
2566 
2567  req_count = 0
2568 
2569  !--- receive r2r
2570  do irank = 1, recv_nmax_r2r
2571  req_count = req_count + 1
2572  totalsize = recv_info_r2r(i_size ,irank) * kmax * vmax
2573  rank = recv_info_r2r(i_prc_from,irank)
2574  tag = recv_info_r2r(i_prc_from,irank)
2575 
2576  call mpi_irecv( recvbuf_r2r_dp(1,irank), & ! [OUT]
2577  totalsize, & ! [IN]
2578  mpi_double_precision, & ! [IN]
2579  rank, & ! [IN]
2580  tag, & ! [IN]
2581  prc_local_comm_world, & ! [IN]
2582  req_list(req_count), & ! [OUT]
2583  ierr ) ! [OUT]
2584  enddo
2585 
2586  !--- pack and send r2r
2587  do irank = 1, send_nmax_r2r
2588  imax = send_info_r2r(i_size,irank)
2589 
2590  !$omp parallel do default(none),private(ipos,k,v,ikv,ij_from,l_from) &
2591  !$omp shared(irank,imax,kmax,vmax,Send_list_r2r,sendbuf_r2r_DP,var)
2592  do v = 1, vmax
2593  do k = 1, kmax
2594  do ipos = 1, imax
2595  ij_from = send_list_r2r(i_grid_from,ipos,irank)
2596  l_from = send_list_r2r(i_l_from ,ipos,irank)
2597 
2598  ikv = (v-1) * imax * kmax &
2599  + (k-1) * imax &
2600  + ipos
2601 
2602  sendbuf_r2r_dp(ikv,irank) = var(ij_from,k,l_from,v)
2603  enddo
2604  enddo
2605  enddo
2606  !$omp end parallel do
2607 
2608  req_count = req_count + 1
2609  totalsize = imax * kmax * vmax
2610  rank = send_info_r2r(i_prc_to ,irank)
2611  tag = send_info_r2r(i_prc_from,irank)
2612 
2613  call mpi_isend( sendbuf_r2r_dp(1,irank), & ! [IN]
2614  totalsize, & ! [IN]
2615  mpi_double_precision, & ! [IN]
2616  rank, & ! [IN]
2617  tag, & ! [IN]
2618  prc_local_comm_world, & ! [IN]
2619  req_list(req_count), & ! [OUT]
2620  ierr ) ! [OUT]
2621  enddo
2622 
2623  !--- copy r2r
2624  !$omp parallel do default(none),private(ipos,k,v,imax,irank,ij_from,l_from,ij_to,l_to) &
2625  !$omp shared(kmax,vmax,Copy_nmax_r2r,Copy_info_r2r,Copy_list_r2r,var)
2626  do irank = 1, copy_nmax_r2r
2627  imax = copy_info_r2r(i_size)
2628 
2629  do v = 1, vmax
2630  do k = 1, kmax
2631  do ipos = 1, imax
2632  ij_from = copy_list_r2r(i_grid_from,ipos)
2633  l_from = copy_list_r2r(i_l_from ,ipos)
2634  ij_to = copy_list_r2r(i_grid_to ,ipos)
2635  l_to = copy_list_r2r(i_l_to ,ipos)
2636 
2637  var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2638  enddo
2639  enddo
2640  enddo
2641  enddo
2642  !$omp end parallel do
2643 
2644  !--- wait all
2645  if ( req_count > 0 ) then
2646  call mpi_waitall( req_count, & ! [IN]
2647  req_list(:), & ! [IN]
2648  mpi_statuses_ignore, & ! [OUT]
2649  ierr ) ! [OUT]
2650  endif
2651 
2652  !--- unpack r2r
2653  !$omp parallel do default(none),private(ipos,k,v,imax,irank,ij_to,l_to,ikv) &
2654  !$omp shared(kmax,vmax,Recv_nmax_r2r,Recv_info_r2r,Recv_list_r2r,recvbuf_r2r_DP,var)
2655  do irank = 1, recv_nmax_r2r
2656  imax = recv_info_r2r(i_size,irank)
2657 
2658  do v = 1, vmax
2659  do k = 1, kmax
2660  do ipos = 1, imax
2661  ij_to = recv_list_r2r(i_grid_to,ipos,irank)
2662  l_to = recv_list_r2r(i_l_to ,ipos,irank)
2663 
2664  ikv = (v-1) * imax * kmax &
2665  + (k-1) * imax &
2666  + ipos
2667 
2668  var(ij_to,k,l_to,v) = recvbuf_r2r_dp(ikv,irank)
2669  enddo
2670  enddo
2671  enddo
2672  enddo
2673  !$omp end parallel do
2674 
2675  !--- singular point (halo to halo)
2676  do irank = 1, singular_nmax
2677  imax = singular_info(i_size)
2678 
2679  do v = 1, vmax
2680  do k = 1, kmax
2681  do ipos = 1, imax
2682  ij_from = singular_list(i_grid_from,ipos)
2683  l_from = singular_list(i_l_from ,ipos)
2684  ij_to = singular_list(i_grid_to ,ipos)
2685  l_to = singular_list(i_l_to ,ipos)
2686 
2687  var(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2688  enddo
2689  enddo
2690  enddo
2691  enddo
2692 
2693  call prof_rapend('COMM_data_transfer',2)
2694 
2695  return
2696  end subroutine comm_data_transfer_nopl
2697 
2698  !-----------------------------------------------------------------------------
2700  subroutine comm_var_sp( &
2701  var, &
2702  var_pl, &
2703  kmax, &
2704  vmax )
2705  use scale_prc, only: &
2707  use scale_prc_icoa, only: &
2708  prc_rgn_total_pl, &
2709  prc_rgn_rgn4pl, &
2710  prc_rgn_lp2r, &
2711  i_npl, &
2712  i_spl
2713  implicit none
2714 
2715  integer, intent(in) :: kmax
2716  integer, intent(in) :: vmax
2717  real(SP), intent(inout) :: var (ADM_gall, kmax,ADM_lall, vmax)
2718  real(SP), intent(inout) :: var_pl(ADM_gall_pl,kmax,ADM_lall_pl,vmax)
2719 
2720  real(SP) :: sendbuf_h2p_SP(kmax*vmax,PRC_RGN_total_pl)
2721  real(SP) :: recvbuf_h2p_SP(kmax*vmax,PRC_RGN_total_pl)
2722 
2723  integer :: totalsize, rank, tag
2724  integer :: irank, ipos
2725  integer :: ij_from, l_from, ij_to, l_to
2726  integer :: r_from, r_to
2727 
2728  integer :: k, v, kk
2729  integer :: ierr
2730  !---------------------------------------------------------------------------
2731 
2732  if ( comm_apply_barrier ) then
2733  call prof_rapstart('COMM_barrier',2)
2734  call mpi_barrier(prc_local_comm_world,ierr)
2735  call prof_rapend ('COMM_barrier',2)
2736  endif
2737 
2738  call prof_rapstart('COMM_var',2)
2739 
2740  if ( comm_pl ) then
2741 
2742  req_count = 0
2743 
2744  !--- receive p2r-reverse
2745  do irank = 1, send_nmax_p2r
2746  do ipos = 1, send_info_p2r(i_size,irank)
2747  l_from = send_list_p2r(i_l_to ,ipos,irank)
2748  r_from = prc_rgn_lp2r(l_from,send_info_p2r(i_prc_to,irank))
2749 
2750  if ( r_from == prc_rgn_rgn4pl(i_npl) ) then
2751  req_count = req_count + 1
2752  totalsize = kmax * vmax
2753  rank = send_info_p2r(i_prc_to ,irank)
2754  tag = send_info_p2r(i_prc_from,irank) + 100000
2755 
2756  call mpi_irecv( recvbuf_h2p_sp(1,i_npl), & ! [OUT]
2757  totalsize, & ! [IN]
2758  mpi_real, & ! [IN]
2759  rank, & ! [IN]
2760  tag, & ! [IN]
2761  prc_local_comm_world, & ! [IN]
2762  req_list(req_count), & ! [OUT]
2763  ierr ) ! [OUT]
2764  endif
2765 
2766  if ( r_from == prc_rgn_rgn4pl(i_spl) ) then
2767  req_count = req_count + 1
2768  totalsize = kmax * vmax
2769  rank = send_info_p2r(i_prc_to ,irank)
2770  tag = send_info_p2r(i_prc_from,irank) + 200000
2771 
2772  call mpi_irecv( recvbuf_h2p_sp(1,i_spl), & ! [OUT]
2773  totalsize, & ! [IN]
2774  mpi_real, & ! [IN]
2775  rank, & ! [IN]
2776  tag, & ! [IN]
2777  prc_local_comm_world, & ! [IN]
2778  req_list(req_count), & ! [OUT]
2779  ierr ) ! [OUT]
2780  endif
2781  enddo
2782  enddo
2783 
2784  !--- pack and send p2r-reverse
2785  do irank = 1, recv_nmax_p2r
2786  do ipos = 1, recv_info_p2r(i_size,irank)
2787  ij_from = recv_list_p2r(i_grid_to,ipos,irank)
2788  l_from = recv_list_p2r(i_l_to ,ipos,irank)
2789  r_from = prc_rgn_lp2r(l_from,recv_info_p2r(i_prc_to,irank))
2790 
2791  if ( r_from == prc_rgn_rgn4pl(i_npl) ) then
2792  do v = 1, vmax
2793  do k = 1, kmax
2794  kk = (v-1) * kmax + k
2795  sendbuf_h2p_sp(kk,i_npl) = var(ij_from,k,l_from,v)
2796  enddo
2797  enddo
2798 
2799  req_count = req_count + 1
2800  totalsize = kmax * vmax
2801  rank = recv_info_p2r(i_prc_from,irank)
2802  tag = recv_info_p2r(i_prc_from,irank) + 100000
2803 
2804  call mpi_isend( sendbuf_h2p_sp(1,i_npl), & ! [IN]
2805  totalsize, & ! [IN]
2806  mpi_real, & ! [IN]
2807  rank, & ! [IN]
2808  tag, & ! [IN]
2809  prc_local_comm_world, & ! [IN]
2810  req_list(req_count), & ! [OUT]
2811  ierr ) ! [OUT]
2812  endif
2813 
2814  if ( r_from == prc_rgn_rgn4pl(i_spl) ) then
2815  do v = 1, vmax
2816  do k = 1, kmax
2817  kk = (v-1) * kmax + k
2818  sendbuf_h2p_sp(kk,i_spl) = var(ij_from,k,l_from,v)
2819  enddo
2820  enddo
2821 
2822  req_count = req_count + 1
2823  totalsize = kmax * vmax
2824  rank = recv_info_p2r(i_prc_from,irank)
2825  tag = recv_info_p2r(i_prc_from,irank) + 200000
2826 
2827  call mpi_isend( sendbuf_h2p_sp(1,i_spl), & ! [IN]
2828  totalsize, & ! [IN]
2829  mpi_real, & ! [IN]
2830  rank, & ! [IN]
2831  tag, & ! [IN]
2832  prc_local_comm_world, & ! [IN]
2833  req_list(req_count), & ! [OUT]
2834  ierr ) ! [OUT]
2835  endif
2836  enddo
2837  enddo
2838 
2839  !--- copy p2r-reverse
2840  do irank = 1, copy_nmax_p2r
2841  do ipos = 1, copy_info_p2r(i_size)
2842  ij_from = copy_list_p2r(i_grid_to ,ipos)
2843  l_from = copy_list_p2r(i_l_to ,ipos)
2844  r_from = prc_rgn_lp2r(l_from,copy_info_p2r(i_prc_to))
2845  ij_to = copy_list_p2r(i_grid_from,ipos)
2846  l_to = copy_list_p2r(i_l_from ,ipos)
2847  r_to = prc_rgn_lp2r(l_to,copy_info_p2r(i_prc_from))
2848 
2849  if ( r_from == prc_rgn_rgn4pl(i_npl) &
2850  .OR. r_from == prc_rgn_rgn4pl(i_spl) ) then
2851  do v = 1, vmax
2852  do k = 1, kmax
2853  var_pl(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
2854  enddo
2855  enddo
2856  endif
2857  enddo
2858  enddo
2859 
2860  !--- wait all
2861  if ( req_count > 0 ) then
2862  call mpi_waitall( req_count, & ! [IN]
2863  req_list(:), & ! [IN]
2864  mpi_statuses_ignore, & ! [OUT]
2865  ierr ) ! [OUT]
2866  endif
2867 
2868  !--- unpack p2r-reverse
2869  do irank = 1, send_nmax_p2r
2870  do ipos = 1, send_info_p2r(i_size,irank)
2871  l_from = send_list_p2r(i_l_to ,ipos,irank)
2872  r_from = prc_rgn_lp2r(l_from,send_info_p2r(i_prc_to,irank))
2873 
2874  ij_to = send_list_p2r(i_grid_from,ipos,irank)
2875  l_to = send_list_p2r(i_l_from ,ipos,irank)
2876 
2877  if ( r_from == prc_rgn_rgn4pl(i_npl) ) then
2878  do v = 1, vmax
2879  do k = 1, kmax
2880  kk = (v-1) * kmax + k
2881  var_pl(ij_to,k,l_to,v) = recvbuf_h2p_sp(kk,i_npl)
2882  enddo
2883  enddo
2884  endif
2885 
2886  if ( r_from == prc_rgn_rgn4pl(i_spl) ) then
2887  do v = 1, vmax
2888  do k = 1, kmax
2889  kk = (v-1) * kmax + k
2890  var_pl(ij_to,k,l_to,v) = recvbuf_h2p_sp(kk,i_spl)
2891  enddo
2892  enddo
2893  endif
2894  enddo
2895  enddo
2896 
2897  endif
2898 
2899  call comm_data_transfer_sp(var,var_pl)
2900 
2901  call prof_rapend('COMM_var',2)
2902 
2903  return
2904  end subroutine comm_var_sp
2905 
2906  !-----------------------------------------------------------------------------
2908  subroutine comm_var_dp( &
2909  var, &
2910  var_pl, &
2911  kmax, &
2912  vmax )
2913  use scale_prc, only: &
2915  use scale_prc_icoa, only: &
2916  prc_rgn_total_pl, &
2917  prc_rgn_rgn4pl, &
2918  prc_rgn_lp2r, &
2919  i_npl, &
2920  i_spl
2921  implicit none
2922 
2923  integer, intent(in) :: kmax
2924  integer, intent(in) :: vmax
2925  real(DP), intent(inout) :: var (ADM_gall, kmax,ADM_lall, vmax)
2926  real(DP), intent(inout) :: var_pl(ADM_gall_pl,kmax,ADM_lall_pl,vmax)
2927 
2928  real(DP) :: sendbuf_h2p_DP(kmax*vmax,PRC_RGN_total_pl)
2929  real(DP) :: recvbuf_h2p_DP(kmax*vmax,PRC_RGN_total_pl)
2930 
2931  integer :: totalsize, rank, tag
2932  integer :: irank, ipos
2933  integer :: ij_from, l_from, ij_to, l_to
2934  integer :: r_from, r_to
2935 
2936  integer :: k, v, kk
2937  integer :: ierr
2938  !---------------------------------------------------------------------------
2939 
2940  if ( comm_apply_barrier ) then
2941  call prof_rapstart('COMM_barrier',2)
2942  call mpi_barrier(prc_local_comm_world,ierr)
2943  call prof_rapend ('COMM_barrier',2)
2944  endif
2945 
2946  call prof_rapstart('COMM_var',2)
2947 
2948  if ( comm_pl ) then
2949 
2950  req_count = 0
2951 
2952  !--- receive p2r-reverse
2953  do irank = 1, send_nmax_p2r
2954  do ipos = 1, send_info_p2r(i_size,irank)
2955  l_from = send_list_p2r(i_l_to ,ipos,irank)
2956  r_from = prc_rgn_lp2r(l_from,send_info_p2r(i_prc_to,irank))
2957 
2958  if ( r_from == prc_rgn_rgn4pl(i_npl) ) then
2959  req_count = req_count + 1
2960  totalsize = kmax * vmax
2961  rank = send_info_p2r(i_prc_to ,irank)
2962  tag = send_info_p2r(i_prc_from,irank) + 100000
2963 
2964  call mpi_irecv( recvbuf_h2p_dp(1,i_npl), & ! [OUT]
2965  totalsize, & ! [IN]
2966  mpi_double_precision, & ! [IN]
2967  rank, & ! [IN]
2968  tag, & ! [IN]
2969  prc_local_comm_world, & ! [IN]
2970  req_list(req_count), & ! [OUT]
2971  ierr ) ! [OUT]
2972  endif
2973 
2974  if ( r_from == prc_rgn_rgn4pl(i_spl) ) then
2975  req_count = req_count + 1
2976  totalsize = kmax * vmax
2977  rank = send_info_p2r(i_prc_to ,irank)
2978  tag = send_info_p2r(i_prc_from,irank) + 200000
2979 
2980  call mpi_irecv( recvbuf_h2p_dp(1,i_spl), & ! [OUT]
2981  totalsize, & ! [IN]
2982  mpi_double_precision, & ! [IN]
2983  rank, & ! [IN]
2984  tag, & ! [IN]
2985  prc_local_comm_world, & ! [IN]
2986  req_list(req_count), & ! [OUT]
2987  ierr ) ! [OUT]
2988  endif
2989  enddo
2990  enddo
2991 
2992  !--- pack and send p2r-reverse
2993  do irank = 1, recv_nmax_p2r
2994  do ipos = 1, recv_info_p2r(i_size,irank)
2995  ij_from = recv_list_p2r(i_grid_to,ipos,irank)
2996  l_from = recv_list_p2r(i_l_to ,ipos,irank)
2997  r_from = prc_rgn_lp2r(l_from,recv_info_p2r(i_prc_to,irank))
2998 
2999  if ( r_from == prc_rgn_rgn4pl(i_npl) ) then
3000  do v = 1, vmax
3001  do k = 1, kmax
3002  kk = (v-1) * kmax + k
3003  sendbuf_h2p_dp(kk,i_npl) = var(ij_from,k,l_from,v)
3004  enddo
3005  enddo
3006 
3007  req_count = req_count + 1
3008  totalsize = kmax * vmax
3009  rank = recv_info_p2r(i_prc_from,irank)
3010  tag = recv_info_p2r(i_prc_from,irank) + 100000
3011 
3012  call mpi_isend( sendbuf_h2p_dp(1,i_npl), & ! [IN]
3013  totalsize, & ! [IN]
3014  mpi_double_precision, & ! [IN]
3015  rank, & ! [IN]
3016  tag, & ! [IN]
3017  prc_local_comm_world, & ! [IN]
3018  req_list(req_count), & ! [OUT]
3019  ierr ) ! [OUT]
3020  endif
3021 
3022  if ( r_from == prc_rgn_rgn4pl(i_spl) ) then
3023  do v = 1, vmax
3024  do k = 1, kmax
3025  kk = (v-1) * kmax + k
3026  sendbuf_h2p_dp(kk,i_spl) = var(ij_from,k,l_from,v)
3027  enddo
3028  enddo
3029 
3030  req_count = req_count + 1
3031  totalsize = kmax * vmax
3032  rank = recv_info_p2r(i_prc_from,irank)
3033  tag = recv_info_p2r(i_prc_from,irank) + 200000
3034 
3035  call mpi_isend( sendbuf_h2p_dp(1,i_spl), & ! [IN]
3036  totalsize, & ! [IN]
3037  mpi_double_precision, & ! [IN]
3038  rank, & ! [IN]
3039  tag, & ! [IN]
3040  prc_local_comm_world, & ! [IN]
3041  req_list(req_count), & ! [OUT]
3042  ierr ) ! [OUT]
3043  endif
3044  enddo
3045  enddo
3046 
3047  !--- copy p2r-reverse
3048  do irank = 1, copy_nmax_p2r
3049  do ipos = 1, copy_info_p2r(i_size)
3050  ij_from = copy_list_p2r(i_grid_to ,ipos)
3051  l_from = copy_list_p2r(i_l_to ,ipos)
3052  r_from = prc_rgn_lp2r(l_from,copy_info_p2r(i_prc_to))
3053  ij_to = copy_list_p2r(i_grid_from,ipos)
3054  l_to = copy_list_p2r(i_l_from ,ipos)
3055  r_to = prc_rgn_lp2r(l_to,copy_info_p2r(i_prc_from))
3056 
3057  if ( r_from == prc_rgn_rgn4pl(i_npl) &
3058  .OR. r_from == prc_rgn_rgn4pl(i_spl) ) then
3059  do v = 1, vmax
3060  do k = 1, kmax
3061  var_pl(ij_to,k,l_to,v) = var(ij_from,k,l_from,v)
3062  enddo
3063  enddo
3064  endif
3065  enddo
3066  enddo
3067 
3068  !--- wait all
3069  if ( req_count > 0 ) then
3070  call mpi_waitall( req_count, & ! [IN]
3071  req_list(:), & ! [IN]
3072  mpi_statuses_ignore, & ! [OUT]
3073  ierr ) ! [OUT]
3074  endif
3075 
3076  !--- unpack p2r-reverse
3077  do irank = 1, send_nmax_p2r
3078  do ipos = 1, send_info_p2r(i_size,irank)
3079  l_from = send_list_p2r(i_l_to ,ipos,irank)
3080  r_from = prc_rgn_lp2r(l_from,send_info_p2r(i_prc_to,irank))
3081 
3082  ij_to = send_list_p2r(i_grid_from,ipos,irank)
3083  l_to = send_list_p2r(i_l_from ,ipos,irank)
3084 
3085  if ( r_from == prc_rgn_rgn4pl(i_npl) ) then
3086  do v = 1, vmax
3087  do k = 1, kmax
3088  kk = (v-1) * kmax + k
3089  var_pl(ij_to,k,l_to,v) = recvbuf_h2p_dp(kk,i_npl)
3090  enddo
3091  enddo
3092  endif
3093 
3094  if ( r_from == prc_rgn_rgn4pl(i_spl) ) then
3095  do v = 1, vmax
3096  do k = 1, kmax
3097  kk = (v-1) * kmax + k
3098  var_pl(ij_to,k,l_to,v) = recvbuf_h2p_dp(kk,i_spl)
3099  enddo
3100  enddo
3101  endif
3102  enddo
3103  enddo
3104 
3105  endif
3106 
3107  call comm_data_transfer_dp(var,var_pl)
3108 
3109  call prof_rapend('COMM_var',2)
3110 
3111  return
3112  end subroutine comm_var_dp
3113 
3114  !-----------------------------------------------------------------------------
3117  function suf(i,j) result(suffix)
3118  implicit none
3119 
3120  integer :: suffix
3121  integer :: i, j
3122  !---------------------------------------------------------------------------
3123 
3124  suffix = adm_gall_1d * (j-1) + i
3125 
3126  end function suf
3127 
3128  !-----------------------------------------------------------------------------
3129  subroutine comm_debugtest
3130  use scale_prc, only: &
3131  prc_myrank, &
3133  use scale_prc_icoa, only: &
3134  prc_rgn_l2r, &
3135  prc_rgn_vert_pl, &
3136  prc_rgn_have_sgp, &
3137  i_n, &
3138  i_s, &
3139  i_npl, &
3140  i_spl
3141  implicit none
3142 
3143  real(rp) :: var (adm_gall ,adm_kall,adm_lall ,4)
3144  real(rp) :: var_pl(adm_gall_pl,adm_kall,adm_lall_pl,4)
3145 
3146  integer :: i, j, k, l, ij, rgnid, prc
3147  !---------------------------------------------------------------------------
3148 
3149  if( io_l ) write(io_fid_log,*)
3150  if( io_l ) write(io_fid_log,*) '+++ TEST start'
3151 
3152  var(:,:,:,:) = -999.0_rp
3153  var_pl(:,:,:,:) = -999.0_rp
3154 
3155  do l = 1, adm_lall
3156  rgnid = prc_rgn_l2r(l)
3157  prc = prc_myrank
3158 
3159  do k = adm_kmin, adm_kmax
3160  do j = adm_gmin, adm_gmax
3161  do i = adm_gmin, adm_gmax
3162  ij = adm_gall_1d * (j-1) + i
3163 
3164  var(ij,k,l,1) = real(prc, kind=rp)
3165  var(ij,k,l,2) = real(rgnid,kind=rp)
3166  var(ij,k,l,3) = real(i, kind=rp)
3167  var(ij,k,l,4) = real(j, kind=rp)
3168  enddo
3169  enddo
3170  enddo
3171 
3172  if ( prc_rgn_have_sgp(l) ) then
3173  do k = adm_kmin, adm_kmax
3174  var(1,k,l,:) = -1.0_rp
3175  enddo
3176  endif
3177  enddo
3178 
3179  do l = 1, adm_lall_pl
3180  rgnid = l
3181  prc = prc_myrank
3182 
3183  do k = adm_kmin, adm_kmax
3184  do ij = 1, adm_gall_pl
3185  var_pl(ij,k,l,1) = real(-prc, kind=rp)
3186  var_pl(ij,k,l,2) = real(-rgnid,kind=rp)
3187  var_pl(ij,k,l,3) = real(-ij, kind=rp)
3188  var_pl(ij,k,l,4) = real(-ij, kind=rp)
3189  enddo
3190  enddo
3191  enddo
3192 
3193  if( io_l ) write(io_fid_log,*) '##### (prc,rgnid) #####'
3194  do l = 1, adm_lall
3195  do k = adm_kmin, adm_kmin
3196  if( io_l ) write(io_fid_log,*)
3197  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3198  do i = 1, adm_gall_1d
3199  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3200  enddo
3201  if( io_l ) write(io_fid_log,*)
3202 
3203  do j = adm_gall_1d, 1, -1
3204  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3205  do i = 1, adm_gall_1d
3206  ij = adm_gall_1d * (j-1) + i
3207  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3208  '(',int(var(ij,k,l,1)),',',int(var(ij,k,l,2)),')'
3209  enddo
3210  if( io_l ) write(io_fid_log,*)
3211  enddo
3212  enddo
3213  enddo
3214 
3215  do l = 1, adm_lall_pl
3216  do k = adm_kmin, adm_kmin
3217  if( io_l ) write(io_fid_log,*)
3218  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3219  do i = 1, adm_gall_pl
3220  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3221  enddo
3222  if( io_l ) write(io_fid_log,*)
3223 
3224  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3225  do ij = 1, adm_gall_pl
3226  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3227  '(',int(var_pl(ij,k,l,1)),',',int(var_pl(ij,k,l,2)),')'
3228  enddo
3229  if( io_l ) write(io_fid_log,*)
3230  enddo
3231  enddo
3232 
3233  if( io_l ) write(io_fid_log,*) '##### (i,j) #####'
3234  do l = 1, adm_lall
3235  do k = adm_kmin, adm_kmin
3236  if( io_l ) write(io_fid_log,*)
3237  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3238  do i = 1, adm_gall_1d
3239  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3240  enddo
3241  if( io_l ) write(io_fid_log,*)
3242 
3243  do j = adm_gall_1d, 1, -1
3244  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3245  do i = 1, adm_gall_1d
3246  ij = adm_gall_1d * (j-1) + i
3247  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3248  '(',int(var(ij,k,l,3)),',',int(var(ij,k,l,4)),')'
3249  enddo
3250  if( io_l ) write(io_fid_log,*)
3251  enddo
3252  enddo
3253  enddo
3254 
3255  do l = 1, adm_lall_pl
3256  do k = adm_kmin, adm_kmin
3257  if( io_l ) write(io_fid_log,*)
3258  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3259  do i = 1, adm_gall_pl
3260  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3261  enddo
3262  if( io_l ) write(io_fid_log,*)
3263 
3264  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3265  do ij = 1, adm_gall_pl
3266  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3267  '(',int(var_pl(ij,k,l,3)),',',int(var_pl(ij,k,l,4)),')'
3268  enddo
3269  if( io_l ) write(io_fid_log,*)
3270  enddo
3271  enddo
3272 
3273 
3274 
3275  if( io_l ) write(io_fid_log,*)
3276  if( io_l ) write(io_fid_log,*) '+++ Communication start'
3277 
3278  call comm_data_transfer( var(:,:,:,:), var_pl(:,:,:,:) )
3279 
3280  if( io_l ) write(io_fid_log,*) '##### (prc,rgnid) #####'
3281  do l = 1, adm_lall
3282  do k = adm_kmin, adm_kmin
3283  if( io_l ) write(io_fid_log,*)
3284  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3285  do i = 1, adm_gall_1d
3286  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3287  enddo
3288  if( io_l ) write(io_fid_log,*)
3289 
3290  do j = adm_gall_1d, 1, -1
3291  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3292  do i = 1, adm_gall_1d
3293  ij = adm_gall_1d * (j-1) + i
3294  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3295  '(',int(var(ij,k,l,1)),',',int(var(ij,k,l,2)),')'
3296  enddo
3297  if( io_l ) write(io_fid_log,*)
3298  enddo
3299  enddo
3300  enddo
3301 
3302  do l = 1, adm_lall_pl
3303  do k = adm_kmin, adm_kmin
3304  if( io_l ) write(io_fid_log,*)
3305  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3306  do i = 1, adm_gall_pl
3307  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3308  enddo
3309  if( io_l ) write(io_fid_log,*)
3310 
3311  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3312  do ij = 1, adm_gall_pl
3313  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3314  '(',int(var_pl(ij,k,l,1)),',',int(var_pl(ij,k,l,2)),')'
3315  enddo
3316  if( io_l ) write(io_fid_log,*)
3317  enddo
3318  enddo
3319 
3320  if( io_l ) write(io_fid_log,*) '##### (i,j) #####'
3321  do l = 1, adm_lall
3322  do k = adm_kmin, adm_kmin
3323  if( io_l ) write(io_fid_log,*)
3324  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3325  do i = 1, adm_gall_1d
3326  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3327  enddo
3328  if( io_l ) write(io_fid_log,*)
3329 
3330  do j = adm_gall_1d, 1, -1
3331  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3332  do i = 1, adm_gall_1d
3333  ij = adm_gall_1d * (j-1) + i
3334  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3335  '(',int(var(ij,k,l,3)),',',int(var(ij,k,l,4)),')'
3336  enddo
3337  if( io_l ) write(io_fid_log,*)
3338  enddo
3339  enddo
3340  enddo
3341 
3342  do l = 1, adm_lall_pl
3343  do k = adm_kmin, adm_kmin
3344  if( io_l ) write(io_fid_log,*)
3345  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3346  do i = 1, adm_gall_pl
3347  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3348  enddo
3349  if( io_l ) write(io_fid_log,*)
3350 
3351  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3352  do ij = 1, adm_gall_pl
3353  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3354  '(',int(var_pl(ij,k,l,3)),',',int(var_pl(ij,k,l,4)),')'
3355  enddo
3356  if( io_l ) write(io_fid_log,*)
3357  enddo
3358  enddo
3359 
3360 
3361 
3362  do l = 1, adm_lall
3363  rgnid = prc_rgn_l2r(l)
3364  prc = prc_myrank
3365 
3366  if ( prc_rgn_vert_pl(i_npl,rgnid) ) then
3367  do k = adm_kmin, adm_kmax
3368  j = adm_gmax+1
3369  i = adm_gmin
3370  ij = adm_gall_1d * (j-1) + i
3371 
3372  var(ij,k,l,1) = real(prc, kind=rp)
3373  var(ij,k,l,2) = real(rgnid,kind=rp)
3374  var(ij,k,l,3) = real(i, kind=rp)
3375  var(ij,k,l,4) = real(j, kind=rp)
3376  enddo
3377  endif
3378 
3379  if ( prc_rgn_vert_pl(i_spl,rgnid) ) then
3380  do k = adm_kmin, adm_kmax
3381  j = adm_gmin
3382  i = adm_gmax+1
3383  ij = adm_gall_1d * (j-1) + i
3384 
3385  var(ij,k,l,1) = real(prc, kind=rp)
3386  var(ij,k,l,2) = real(rgnid,kind=rp)
3387  var(ij,k,l,3) = real(i, kind=rp)
3388  var(ij,k,l,4) = real(j, kind=rp)
3389  enddo
3390  endif
3391 
3392  enddo
3393 
3394  if( io_l ) write(io_fid_log,*)
3395  if( io_l ) write(io_fid_log,*) '+++ pole fill start'
3396 
3397  call comm_var( var(:,:,:,:), var_pl(:,:,:,:), adm_kall, 4 )
3398 
3399  if( io_l ) write(io_fid_log,*) '##### (prc,rgnid) #####'
3400  do l = 1, adm_lall
3401  do k = adm_kmin, adm_kmin
3402  if( io_l ) write(io_fid_log,*)
3403  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3404  do i = 1, adm_gall_1d
3405  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3406  enddo
3407  if( io_l ) write(io_fid_log,*)
3408 
3409  do j = adm_gall_1d, 1, -1
3410  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3411  do i = 1, adm_gall_1d
3412  ij = adm_gall_1d * (j-1) + i
3413  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3414  '(',int(var(ij,k,l,1)),',',int(var(ij,k,l,2)),')'
3415  enddo
3416  if( io_l ) write(io_fid_log,*)
3417  enddo
3418  enddo
3419  enddo
3420 
3421  do l = 1, adm_lall_pl
3422  do k = adm_kmin, adm_kmin
3423  if( io_l ) write(io_fid_log,*)
3424  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3425  do i = 1, adm_gall_pl
3426  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3427  enddo
3428  if( io_l ) write(io_fid_log,*)
3429 
3430  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3431  do ij = 1, adm_gall_pl
3432  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3433  '(',int(var_pl(ij,k,l,1)),',',int(var_pl(ij,k,l,2)),')'
3434  enddo
3435  if( io_l ) write(io_fid_log,*)
3436  enddo
3437  enddo
3438 
3439  if( io_l ) write(io_fid_log,*) '##### (i,j) #####'
3440  do l = 1, adm_lall
3441  do k = adm_kmin, adm_kmin
3442  if( io_l ) write(io_fid_log,*)
3443  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3444  do i = 1, adm_gall_1d
3445  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3446  enddo
3447  if( io_l ) write(io_fid_log,*)
3448 
3449  do j = adm_gall_1d, 1, -1
3450  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3451  do i = 1, adm_gall_1d
3452  ij = adm_gall_1d * (j-1) + i
3453  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3454  '(',int(var(ij,k,l,3)),',',int(var(ij,k,l,4)),')'
3455  enddo
3456  if( io_l ) write(io_fid_log,*)
3457  enddo
3458  enddo
3459  enddo
3460 
3461  do l = 1, adm_lall_pl
3462  do k = adm_kmin, adm_kmin
3463  if( io_l ) write(io_fid_log,*)
3464  if( io_l ) write(io_fid_log,'(A9)',advance='no') ' |'
3465  do i = 1, adm_gall_pl
3466  if( io_l ) write(io_fid_log,'(I9)',advance='no') i
3467  enddo
3468  if( io_l ) write(io_fid_log,*)
3469 
3470  if( io_l ) write(io_fid_log,'(I8,A1)',advance='no') j, '|'
3471  do ij = 1, adm_gall_pl
3472  if( io_l ) write(io_fid_log,'(A1,I3,A1,I3,A1)',advance='no') &
3473  '(',int(var_pl(ij,k,l,3)),',',int(var_pl(ij,k,l,4)),')'
3474  enddo
3475  if( io_l ) write(io_fid_log,*)
3476  enddo
3477  enddo
3478 
3479  call prc_mpifinish
3480 
3481  stop
3482  end subroutine comm_debugtest
3483 
3484  !-----------------------------------------------------------------------------
3485  subroutine comm_stat_sum_sp( localsum, globalsum )
3486  use scale_prc, only: &
3488  prc_nprocs
3489  implicit none
3490 
3491  real(SP), intent(in) :: localsum
3492  real(SP), intent(out) :: globalsum
3493 
3494  real(SP) :: sendbuf(1)
3495  real(SP) :: recvbuf(PRC_nprocs)
3496 
3497  integer :: ierr
3498  !---------------------------------------------------------------------------
3499 
3500  if ( comm_pl ) then
3501  sendbuf(1) = localsum
3502 
3503  call mpi_allgather( sendbuf, &
3504  1, &
3505  mpi_real, &
3506  recvbuf, &
3507  1, &
3508  mpi_real, &
3510  ierr )
3511 
3512  globalsum = sum( recvbuf(:) )
3513  else
3514  globalsum = localsum
3515  endif
3516 
3517  return
3518  end subroutine comm_stat_sum_sp
3519 
3520  !-----------------------------------------------------------------------------
3521  subroutine comm_stat_sum_dp( localsum, globalsum )
3522  use scale_prc, only: &
3524  prc_nprocs
3525  implicit none
3526 
3527  real(DP), intent(in) :: localsum
3528  real(DP), intent(out) :: globalsum
3529 
3530  real(DP) :: sendbuf(1)
3531  real(DP) :: recvbuf(PRC_nprocs)
3532 
3533  integer :: ierr
3534  !---------------------------------------------------------------------------
3535 
3536  if ( comm_pl ) then
3537  sendbuf(1) = localsum
3538 
3539  call mpi_allgather( sendbuf, &
3540  1, &
3541  mpi_double_precision, &
3542  recvbuf, &
3543  1, &
3544  mpi_double_precision, &
3546  ierr )
3547 
3548  globalsum = sum( recvbuf(:) )
3549  else
3550  globalsum = localsum
3551  endif
3552 
3553  return
3554  end subroutine comm_stat_sum_dp
3555 
3556  !-----------------------------------------------------------------------------
3557  subroutine comm_stat_sum_eachlayer_sp( kall, localsum, globalsum )
3558  use scale_prc, only: &
3560  prc_nprocs
3561  implicit none
3562 
3563  integer, intent(in) :: kall
3564  real(SP), intent(in) :: localsum (kall)
3565  real(SP), intent(out) :: globalsum(kall)
3566 
3567  real(SP) :: sendbuf(kall)
3568  integer :: displs (PRC_nprocs)
3569  integer :: counts (PRC_nprocs)
3570  real(SP) :: recvbuf(kall,PRC_nprocs)
3571 
3572  integer :: ierr
3573  integer :: k, p
3574  !---------------------------------------------------------------------------
3575 
3576  do p = 1, prc_nprocs
3577  displs(p) = (p-1) * kall
3578  counts(p) = kall
3579  enddo
3580 
3581  if ( comm_pl ) then
3582  sendbuf(:) = localsum(:)
3583 
3584  call mpi_allgatherv( sendbuf, &
3585  kall, &
3586  mpi_real, &
3587  recvbuf, &
3588  counts, &
3589  displs, &
3590  mpi_real, &
3592  ierr )
3593 
3594  do k = 1, kall
3595  globalsum(k) = sum( recvbuf(k,:) )
3596  enddo
3597  else
3598  do k = 1, kall
3599  globalsum(k) = localsum(k)
3600  enddo
3601  endif
3602 
3603  return
3604  end subroutine comm_stat_sum_eachlayer_sp
3605 
3606  !-----------------------------------------------------------------------------
3607  subroutine comm_stat_sum_eachlayer_dp( kall, localsum, globalsum )
3608  use scale_prc, only: &
3610  prc_nprocs
3611  implicit none
3612 
3613  integer, intent(in) :: kall
3614  real(DP), intent(in) :: localsum (kall)
3615  real(DP), intent(out) :: globalsum(kall)
3616 
3617  real(DP) :: sendbuf(kall)
3618  integer :: displs (PRC_nprocs)
3619  integer :: counts (PRC_nprocs)
3620  real(DP) :: recvbuf(kall,PRC_nprocs)
3621 
3622  integer :: ierr
3623  integer :: k, p
3624  !---------------------------------------------------------------------------
3625 
3626  do p = 1, prc_nprocs
3627  displs(p) = (p-1) * kall
3628  counts(p) = kall
3629  enddo
3630 
3631  if ( comm_pl ) then
3632  sendbuf(:) = localsum(:)
3633 
3634  call mpi_allgatherv( sendbuf, &
3635  kall, &
3636  mpi_double_precision, &
3637  recvbuf, &
3638  counts, &
3639  displs, &
3640  mpi_double_precision, &
3642  ierr )
3643 
3644  do k = 1, kall
3645  globalsum(k) = sum( recvbuf(k,:) )
3646  enddo
3647  else
3648  do k = 1, kall
3649  globalsum(k) = localsum(k)
3650  enddo
3651  endif
3652 
3653  return
3654  end subroutine comm_stat_sum_eachlayer_dp
3655 
3656  !-----------------------------------------------------------------------------
3657  subroutine comm_stat_avg_sp( localavg, globalavg )
3658  use scale_prc, only: &
3660  prc_nprocs
3661  implicit none
3662 
3663  real(SP), intent(in) :: localavg
3664  real(SP), intent(out) :: globalavg
3665 
3666  real(SP) :: sendbuf(1)
3667  real(SP) :: recvbuf(PRC_nprocs)
3668 
3669  integer :: ierr
3670  !---------------------------------------------------------------------------
3671 
3672  if ( comm_pl ) then
3673  sendbuf(1) = localavg
3674 
3675  call mpi_allgather( sendbuf, &
3676  1, &
3677  mpi_real, &
3678  recvbuf, &
3679  1, &
3680  mpi_real, &
3682  ierr )
3683 
3684  globalavg = sum( recvbuf(:) ) / real(prc_nprocs,kind=sp)
3685  else
3686  globalavg = localavg
3687  endif
3688 
3689  return
3690  end subroutine comm_stat_avg_sp
3691 
3692  !-----------------------------------------------------------------------------
3693  subroutine comm_stat_avg_dp( localavg, globalavg )
3694  use scale_prc, only: &
3696  prc_nprocs
3697  implicit none
3698 
3699  real(DP), intent(in) :: localavg
3700  real(DP), intent(out) :: globalavg
3701 
3702  real(DP) :: sendbuf(1)
3703  real(DP) :: recvbuf(PRC_nprocs)
3704 
3705  integer :: ierr
3706  !---------------------------------------------------------------------------
3707 
3708  if ( comm_pl ) then
3709  sendbuf(1) = localavg
3710 
3711  call mpi_allgather( sendbuf, &
3712  1, &
3713  mpi_double_precision, &
3714  recvbuf, &
3715  1, &
3716  mpi_double_precision, &
3718  ierr )
3719 
3720  globalavg = sum( recvbuf(:) ) / real(prc_nprocs,kind=dp)
3721  else
3722  globalavg = localavg
3723  endif
3724 
3725  return
3726  end subroutine comm_stat_avg_dp
3727 
3728  !-----------------------------------------------------------------------------
3729  subroutine comm_stat_max_sp( localmax, globalmax )
3730  use scale_prc, only: &
3732  prc_nprocs
3733  implicit none
3734 
3735  real(SP), intent(in) :: localmax
3736  real(SP), intent(out) :: globalmax
3737 
3738  real(SP) :: sendbuf(1)
3739  real(SP) :: recvbuf(PRC_nprocs)
3740 
3741  integer :: ierr
3742  !---------------------------------------------------------------------------
3743 
3744  sendbuf(1) = localmax
3745 
3746  call mpi_allgather( sendbuf, &
3747  1, &
3748  mpi_real, &
3749  recvbuf, &
3750  1, &
3751  mpi_real, &
3753  ierr )
3754 
3755  globalmax = maxval( recvbuf(:) )
3756 
3757  return
3758  end subroutine comm_stat_max_sp
3759 
3760  !-----------------------------------------------------------------------------
3761  subroutine comm_stat_max_dp( localmax, globalmax )
3762  use scale_prc, only: &
3764  prc_nprocs
3765  implicit none
3766 
3767  real(DP), intent(in) :: localmax
3768  real(DP), intent(out) :: globalmax
3769 
3770  real(DP) :: sendbuf(1)
3771  real(DP) :: recvbuf(PRC_nprocs)
3772 
3773  integer :: ierr
3774  !---------------------------------------------------------------------------
3775 
3776  sendbuf(1) = localmax
3777 
3778  call mpi_allgather( sendbuf, &
3779  1, &
3780  mpi_double_precision, &
3781  recvbuf, &
3782  1, &
3783  mpi_double_precision, &
3785  ierr )
3786 
3787  globalmax = maxval( recvbuf(:) )
3788 
3789  return
3790  end subroutine comm_stat_max_dp
3791 
3792  !-----------------------------------------------------------------------------
3793  subroutine comm_stat_min_sp( localmin, globalmin )
3794  use scale_prc, only: &
3796  prc_nprocs
3797  implicit none
3798 
3799  real(SP), intent(in) :: localmin
3800  real(SP), intent(out) :: globalmin
3801 
3802  real(SP) :: sendbuf(1)
3803  real(SP) :: recvbuf(PRC_nprocs)
3804 
3805  integer :: ierr
3806  !---------------------------------------------------------------------------
3807 
3808  sendbuf(1) = localmin
3809 
3810  call mpi_allgather( sendbuf, &
3811  1, &
3812  mpi_real, &
3813  recvbuf, &
3814  1, &
3815  mpi_real, &
3817  ierr )
3818 
3819  globalmin = minval( recvbuf(:) )
3820 
3821  return
3822  end subroutine comm_stat_min_sp
3823 
3824  !-----------------------------------------------------------------------------
3825  subroutine comm_stat_min_dp( localmin, globalmin )
3826  use scale_prc, only: &
3828  prc_nprocs
3829  implicit none
3830 
3831  real(DP), intent(in) :: localmin
3832  real(DP), intent(out) :: globalmin
3833 
3834  real(DP) :: sendbuf(1)
3835  real(DP) :: recvbuf(PRC_nprocs)
3836 
3837  integer :: ierr
3838  !---------------------------------------------------------------------------
3839 
3840  sendbuf(1) = localmin
3841 
3842  call mpi_allgather( sendbuf, &
3843  1, &
3844  mpi_double_precision, &
3845  recvbuf, &
3846  1, &
3847  mpi_double_precision, &
3849  ierr )
3850 
3851  globalmin = minval( recvbuf(:) )
3852 
3853  return
3854  end subroutine comm_stat_min_dp
3855 
3856 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:350
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:1567
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:89
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:3118
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_prc_icoa::prc_rgn_ndiamond
integer, public prc_rgn_ndiamond
number of diamonds
Definition: scale_prc_icoA.F90:67
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:174
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:91
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:2524
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:3658
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:3762
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:3558
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:2705
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:3486
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:58
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:2913
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:3826
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:3522
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:90
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:3794
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:3730
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:63
scale_prof::prof_rapend
subroutine, public prof_rapend(rapname_base, level, disable_barrier)
Save raptime.
Definition: scale_prof.F90:246
scale_io::io_nml
logical, public io_nml
output log or not? (for namelist, this process)
Definition: scale_io.F90:64
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:3608
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:366
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:3694
scale_io::io_fid_conf
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:57
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:2046
scale_comm_icoa::send_info_r2p
integer, dimension(:,:), allocatable, public send_info_r2p
Definition: scale_comm_icoA.F90:92