SCALE-RM
scale_comm_cartesC.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 #include "scalelib.h"
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use mpi
17  use scale_precision
18  use scale_io
19  use scale_prof
21  use scale_tracer
22 
23  use scale_prc, only: &
24  prc_abort
25  use scale_prc_cartesc, only: &
26  prc_next, &
27  prc_w, &
28  prc_n, &
29  prc_e, &
30  prc_s, &
31  prc_nw, &
32  prc_ne, &
33  prc_sw, &
34  prc_se, &
35  prc_has_w, &
36  prc_has_n, &
37  prc_has_e, &
38  prc_has_s
39  !-----------------------------------------------------------------------------
40  implicit none
41  private
42  !-----------------------------------------------------------------------------
43  !
44  !++ Public procedure
45  !
46  public :: comm_setup
47  public :: comm_vars_init
48  public :: comm_vars8_init
49  public :: comm_vars
50  public :: comm_vars8
51  public :: comm_wait
52  public :: comm_gather
53  public :: comm_bcast
54  public :: comm_cleanup
55 
56  interface comm_vars
57  module procedure comm_vars_2d
58  module procedure comm_vars_3d
59  end interface comm_vars
60 
61  interface comm_vars8
62  module procedure comm_vars8_2d
63  module procedure comm_vars8_3d
64  end interface comm_vars8
65 
66  interface comm_wait
67  module procedure comm_wait_2d
68  module procedure comm_wait_3d
69  end interface comm_wait
70 
71  interface comm_gather
72  module procedure comm_gather_2d
73  module procedure comm_gather_3d
74  end interface comm_gather
75 
76  interface comm_bcast
77  module procedure comm_bcast_scr
78  module procedure comm_bcast_1d
79  module procedure comm_bcast_2d
80  module procedure comm_bcast_3d
81  module procedure comm_bcast_4d
82  module procedure comm_bcast_int_scr
83  module procedure comm_bcast_int_1d
84  module procedure comm_bcast_int_2d
85  module procedure comm_bcast_logical_scr
86  module procedure comm_bcast_character
87  end interface comm_bcast
88 
89  !-----------------------------------------------------------------------------
90  !
91  !++ Public parameters & variables
92  !
93  integer, public :: comm_datatype
94  integer, public :: comm_world
95 
96  !-----------------------------------------------------------------------------
97  !
98  !++ Private procedure
99  !
100  !-----------------------------------------------------------------------------
101  !
102  !++ Private parameters & variables
103  !
104  integer, private :: comm_nreq_max
105  integer, private :: comm_vsize_max
106  integer, private :: comm_vsize_max_pc
107 
108  logical, private :: comm_isallperiodic
109 
110  integer, private :: comm_size2d_ns4
111  integer, private :: comm_size2d_ns8
112  integer, private :: comm_size2d_we
113  integer, private :: comm_size2d_4c
114 
115  integer, private :: comm_vars_id = 0
116 
117  logical, private :: comm_use_mpi_pc = .true.
118 
119  real(RP), private, allocatable :: recvpack_w2p(:,:)
120  real(RP), private, allocatable :: recvpack_e2p(:,:)
121  real(RP), private, allocatable :: sendpack_p2w(:,:)
122  real(RP), private, allocatable :: sendpack_p2e(:,:)
123 #ifdef DEBUG
124  logical, private, allocatable :: use_packbuf(:)
125 #endif
126 
127  integer, private, allocatable :: req_cnt (:)
128  integer, private, allocatable :: req_list(:,:)
129  integer, private, allocatable :: preq_cnt (:)
130  integer, private, allocatable :: preq_list(:,:)
131  integer, private, allocatable :: pseqid(:)
132 
133  !-----------------------------------------------------------------------------
134 contains
135 
136  !-----------------------------------------------------------------------------
138  !-----------------------------------------------------------------------------
139  subroutine comm_setup
140  use scale_prc, only: &
142  implicit none
143 
144  namelist / param_comm_cartesc / &
145  comm_vsize_max, &
146  comm_vsize_max_pc, &
147  comm_use_mpi_pc
148 
149  integer :: nreq_NS, nreq_WE, nreq_4C
150 
151  logical, save :: initialized = .false.
152 
153  integer :: ierr
154  !---------------------------------------------------------------------------
155 
156  if ( initialized ) return
157 
158  log_newline
159  log_info("COMM_setup",*) 'Setup'
160 
161  comm_vsize_max = max( 10 + qa*2, 25 )
162  comm_vsize_max_pc = 50 + qa*2
163 
164  !--- read namelist
165  rewind(io_fid_conf)
166  read(io_fid_conf,nml=param_comm_cartesc,iostat=ierr)
167  if( ierr < 0 ) then !--- missing
168  log_info("COMM_setup",*) 'Not found namelist. Default used.'
169  elseif( ierr > 0 ) then !--- fatal error
170  log_error("COMM_setup",*) 'Not appropriate names in namelist PARAM_COMM_CARTESC. Check!'
171  call prc_abort
172  endif
173  log_nml(param_comm_cartesc)
174 
175  nreq_ns = 2 * jhalo !--- send x JHALO, recv x JHALO
176  nreq_we = 2 !--- send x 1 , recv x 1
177  nreq_4c = 2 * jhalo !--- send x JHALO, recv x JHALO
178 
179  if ( comm_use_mpi_pc ) then
180  comm_nreq_max = 2 * nreq_ns + 2 * nreq_we + 4 * nreq_4c + 1
181  else
182  comm_nreq_max = 2 * nreq_ns + 2 * nreq_we + 4 * nreq_4c
183  end if
184 
185  comm_size2d_ns4 = ia * jhalo
186  comm_size2d_ns8 = imax
187  comm_size2d_we = jmax * ihalo
188  comm_size2d_4c = ihalo
189 
190  allocate( recvpack_w2p(comm_size2d_we*ka,comm_vsize_max) )
191  allocate( recvpack_e2p(comm_size2d_we*ka,comm_vsize_max) )
192  allocate( sendpack_p2w(comm_size2d_we*ka,comm_vsize_max) )
193  allocate( sendpack_p2e(comm_size2d_we*ka,comm_vsize_max) )
194 #ifdef DEBUG
195  allocate( use_packbuf(comm_vsize_max) )
196  use_packbuf(:) = .false.
197 #endif
198 
199  allocate( req_cnt( comm_vsize_max) )
200  allocate( req_list(comm_nreq_max,comm_vsize_max) )
201  req_cnt(:) = -1
202  req_list(:,:) = mpi_request_null
203 
204  if ( comm_use_mpi_pc ) then
205  allocate( preq_cnt( comm_vsize_max_pc) )
206  allocate( preq_list(comm_nreq_max+1,comm_vsize_max_pc) )
207  preq_cnt(:) = -1
208  preq_list(:,:) = mpi_request_null
209 
210  allocate( pseqid(comm_vsize_max_pc) )
211  end if
212 
213  if ( prc_has_n .AND. prc_has_s .AND. prc_has_w .AND. prc_has_e ) then
214  comm_isallperiodic = .true.
215  else
216  comm_isallperiodic = .false.
217  endif
218 
219  if ( rp == kind(0.d0) ) then
220  comm_datatype = mpi_double_precision
221  elseif( rp == kind(0.0) ) then
222  comm_datatype = mpi_real
223  else
224  log_error("COMM_setup",*) 'precision is not supportd'
225  call prc_abort
226  endif
227 
229 
230  log_newline
231  log_info("COMM_setup",*) 'Communication information'
232  log_info_cont(*) 'Maximum number of vars for one communication: ', comm_vsize_max
233  log_info_cont(*) 'Data size of var (3D,including halo) [byte] : ', rp*ka*ia*ja
234  log_info_cont(*) 'Data size of halo [byte] : ', rp*ka*(2*ia*jhalo+2*jmax*ihalo)
235  log_info_cont(*) 'Ratio of halo against the whole 3D grid : ', real(2*IA*JHALO+2*JMAX*IHALO) / real(ia*ja)
236  log_info_cont(*) 'All side is periodic? : ', comm_isallperiodic
237 
238  initialized = .true.
239 
240  return
241  end subroutine comm_setup
242 
243  !-----------------------------------------------------------------------------
245  subroutine comm_vars_init( &
246  varname, &
247  var, &
248  vid )
249  implicit none
250 
251  character(len=*), intent(in) :: varname
252  real(RP), intent(inout) :: var(:,:,:)
253  integer, intent(inout) :: vid
254  !---------------------------------------------------------------------------
255 
256  if ( vid > comm_vsize_max ) then
257  log_error("COMM_vars_init",*) 'vid exceeds max', vid, comm_vsize_max
258  call prc_abort
259  end if
260 
261  if ( comm_use_mpi_pc ) then
262 
263  comm_vars_id = comm_vars_id + 1
264  if ( comm_vars_id > comm_vsize_max_pc ) then
265  log_error("COMM_vars_init",*) 'number of variable for MPI PC exceeds max', comm_vars_id, comm_vsize_max_pc
266  call prc_abort
267  end if
268 
269  call prof_rapstart('COMM_init_pers', 2)
270  call vars_init_mpi_pc(var, comm_vars_id, vid)
271  call prof_rapend ('COMM_init_pers', 2)
272 
273  vid = comm_vars_id + comm_vsize_max
274 
275  log_info("COMM_vars_init",'(1x,A,I3.3,2A)') 'Initialize variable : ID = ', vid, &
276  ', name = ', trim(varname)
277 
278  end if
279 
280  return
281  end subroutine comm_vars_init
282 
283  !-----------------------------------------------------------------------------
285  subroutine comm_vars8_init( &
286  varname, &
287  var, &
288  vid )
289  implicit none
290 
291  character(len=*), intent(in) :: varname
292  real(RP), intent(inout) :: var(:,:,:)
293  integer, intent(inout) :: vid
294  !---------------------------------------------------------------------------
295 
296  if ( vid > comm_vsize_max ) then
297  log_error("COMM_vars8_init",*) 'vid exceeds max', vid, comm_vsize_max
298  call prc_abort
299  end if
300 
301  if ( comm_use_mpi_pc ) then
302 
303  comm_vars_id = comm_vars_id + 1
304  if ( comm_vars_id > comm_vsize_max_pc ) then
305  log_error("COMM_vars8_init",*) 'number of variable for MPI PC exceeds max', comm_vars_id, comm_vsize_max_pc
306  call prc_abort
307  end if
308 
309  call prof_rapstart('COMM_init_pers', 2)
310  call vars8_init_mpi_pc(var, comm_vars_id, vid)
311  call prof_rapend ('COMM_init_pers', 2)
312 
313  vid = comm_vars_id + comm_vsize_max
314 
315  log_info("COMM_vars8_init",'(1x,A,I3.3,2A)') 'Initialize variable : ID = ', vid, &
316  ', name = ', trim(varname)
317 
318  end if
319 
320  return
321  end subroutine comm_vars8_init
322 
323  !-----------------------------------------------------------------------------
324  subroutine comm_vars_3d(var, vid)
325  implicit none
326 
327  real(RP), intent(inout) :: var(:,:,:)
328  integer, intent(in) :: vid
329  !---------------------------------------------------------------------------
330 
331  if ( vid > comm_vsize_max ) then
332  call prof_rapstart('COMM_vars_pers', 2)
333  call vars_3d_mpi_pc(var, vid-comm_vsize_max)
334  call prof_rapend ('COMM_vars_pers', 2)
335  else
336  call prof_rapstart('COMM_vars', 2)
337  call vars_3d_mpi(var, vid)
338  call prof_rapend ('COMM_vars', 2)
339  end if
340 
341  return
342  end subroutine comm_vars_3d
343 
344  !-----------------------------------------------------------------------------
345  subroutine comm_vars8_3d(var, vid)
346  implicit none
347 
348  real(RP), intent(inout) :: var(:,:,:)
349  integer, intent(in) :: vid
350  !---------------------------------------------------------------------------
351 
352  if ( vid > comm_vsize_max ) then
353  call prof_rapstart('COMM_vars_pers', 2)
354  call vars_3d_mpi_pc(var, vid-comm_vsize_max)
355  call prof_rapend ('COMM_vars_pers', 2)
356  else
357  call prof_rapstart('COMM_vars', 2)
358  call vars8_3d_mpi(var, vid)
359  call prof_rapend ('COMM_vars', 2)
360  end if
361 
362  return
363  end subroutine comm_vars8_3d
364 
365  !-----------------------------------------------------------------------------
366  subroutine comm_wait_3d(var, vid, FILL_BND)
367  implicit none
368 
369  real(RP), intent(inout) :: var(:,:,:)
370  integer, intent(in) :: vid
371  logical, intent(in), optional :: FILL_BND
372 
373  logical :: FILL_BND_
374  !---------------------------------------------------------------------------
375 
376  fill_bnd_ = .true.
377  if ( present(fill_bnd) ) fill_bnd_ = fill_bnd
378 
379  if ( vid > comm_vsize_max ) then
380  call prof_rapstart('COMM_wait_pers', 2)
381  call wait_3d_mpi_pc(var, vid-comm_vsize_max)
382  call prof_rapend ('COMM_wait_pers', 2)
383  else
384  call prof_rapstart('COMM_wait', 2)
385  call wait_3d_mpi(var, vid)
386  call prof_rapend ('COMM_wait', 2)
387  end if
388 
389  ! copy inner data to boundary
390  if ( .NOT. comm_isallperiodic ) then
391  if ( fill_bnd_ ) then
392  call copy_boundary_3d(var)
393  end if
394  end if
395 
396  return
397  end subroutine comm_wait_3d
398 
399  !-----------------------------------------------------------------------------
400  subroutine comm_vars_2d(var, vid)
401  implicit none
402  real(RP), intent(inout) :: var(:,:)
403  integer, intent(in) :: vid
404  !---------------------------------------------------------------------------
405 
406  call prof_rapstart('COMM_vars', 2)
407  call vars_2d_mpi(var, vid)
408  call prof_rapend ('COMM_vars', 2)
409 
410  return
411  end subroutine comm_vars_2d
412 
413  !-----------------------------------------------------------------------------
414  subroutine comm_vars8_2d(var, vid)
415  implicit none
416 
417  real(RP), intent(inout) :: var(:,:)
418  integer, intent(in) :: vid
419  !---------------------------------------------------------------------------
420 
421  call prof_rapstart('COMM_vars', 2)
422  call vars8_2d_mpi(var, vid)
423  call prof_rapend ('COMM_vars', 2)
424 
425  return
426  end subroutine comm_vars8_2d
427 
428  !-----------------------------------------------------------------------------
429  subroutine comm_wait_2d(var, vid, FILL_BND)
430  implicit none
431 
432  real(RP), intent(inout) :: var(:,:)
433  integer, intent(in) :: vid
434  logical, intent(in), optional :: FILL_BND
435 
436  logical :: FILL_BND_
437  !---------------------------------------------------------------------------
438 
439  fill_bnd_ = .true.
440  if ( present(fill_bnd) ) fill_bnd_ = fill_bnd
441 
442  call prof_rapstart('COMM_wait', 2)
443  call wait_2d_mpi(var, vid)
444  call prof_rapend ('COMM_wait', 2)
445 
446  if( .NOT. comm_isallperiodic ) then
447  if ( fill_bnd_ ) then
448  call copy_boundary_2d(var)
449  end if
450  end if
451 
452  return
453  end subroutine comm_wait_2d
454 
455  !-----------------------------------------------------------------------------
457  subroutine comm_horizontal_mean_2d( varmean, var )
458  use scale_const, only: &
460  implicit none
461 
462  real(RP), intent(out) :: varmean
463  real(RP), intent(in) :: var (ia,ja)
464 
465  real(RP) :: statval
466  real(RP) :: statcnt
467  real(RP) :: allstatval
468  real(RP) :: allstatcnt
469  real(RP) :: zerosw
470 
471  integer :: ierr
472  integer :: i, j
473  !---------------------------------------------------------------------------
474 
475  statval = 0.0_rp
476  statcnt = 0.0_rp
477  do j = js, je
478  do i = is, ie
479  if ( abs(var(i,j)) < abs(const_undef) ) then
480  statval = statval + var(i,j)
481  statcnt = statcnt + 1.0_rp
482  endif
483  enddo
484  enddo
485 
486  ! [NOTE] always communicate globally
487  call prof_rapstart('COMM_Allreduce', 2)
488  ! All reduce
489  call mpi_allreduce( statval, &
490  allstatval, &
491  1, &
492  comm_datatype, &
493  mpi_sum, &
494  comm_world, &
495  ierr )
496  ! All reduce
497  call mpi_allreduce( statcnt, &
498  allstatcnt, &
499  1, &
500  comm_datatype, &
501  mpi_sum, &
502  comm_world, &
503  ierr )
504 
505  call prof_rapend ('COMM_Allreduce', 2)
506 
507  zerosw = 0.5_rp - sign(0.5_rp, allstatcnt - 1.e-12_rp )
508  varmean = allstatval / ( allstatcnt + zerosw ) * ( 1.0_rp - zerosw )
509  !LOG_INFO("COMM_horizontal_mean_2D",*) varmean, allstatval, allstatcnt
510 
511  return
512  end subroutine comm_horizontal_mean_2d
513 
514  !-----------------------------------------------------------------------------
516  subroutine comm_horizontal_mean_3d( varmean, var )
517  use scale_const, only: &
519  implicit none
520 
521  real(RP), intent(out) :: varmean(ka)
522  real(RP), intent(in) :: var (ka,ia,ja)
523 
524  real(RP) :: statval (ka)
525  real(RP) :: statcnt (ka)
526  real(RP) :: allstatval(ka)
527  real(RP) :: allstatcnt(ka)
528  real(RP) :: zerosw
529 
530  integer :: ierr
531  integer :: k, i, j
532  !---------------------------------------------------------------------------
533 
534  statval(:) = 0.0_rp
535  statcnt(:) = 0.0_rp
536  do j = js, je
537  do i = is, ie
538  do k = 1, ka
539  if ( abs(var(k,i,j)) < abs(const_undef) ) then
540  statval(k) = statval(k) + var(k,i,j)
541  statcnt(k) = statcnt(k) + 1.0_rp
542  endif
543  enddo
544  enddo
545  enddo
546 
547  ! [NOTE] always communicate globally
548  call prof_rapstart('COMM_Allreduce', 2)
549  ! All reduce
550  call mpi_allreduce( statval(1), &
551  allstatval(1), &
552  ka, &
553  comm_datatype, &
554  mpi_sum, &
555  comm_world, &
556  ierr )
557  ! All reduce
558  call mpi_allreduce( statcnt(1), &
559  allstatcnt(1), &
560  ka, &
561  comm_datatype, &
562  mpi_sum, &
563  comm_world, &
564  ierr )
565 
566  call prof_rapend ('COMM_Allreduce', 2)
567 
568  do k = 1, ka
569  zerosw = 0.5_rp - sign(0.5_rp, allstatcnt(k) - 1.e-12_rp )
570  varmean(k) = allstatval(k) / ( allstatcnt(k) + zerosw ) * ( 1.0_rp - zerosw )
571  !LOG_INFO("COMM_horizontal_mean_3D",*) k, varmean(k), allstatval(k), allstatcnt(k)
572  enddo
573 
574  return
575  end subroutine comm_horizontal_mean_3d
576 
577  !-----------------------------------------------------------------------------
579  subroutine comm_gather_2d( recv, send, gIA, gJA )
580  use scale_prc, only: &
582  implicit none
583 
584  real(RP), intent(out) :: recv(:,:,:)
585  real(RP), intent(in) :: send(:,:)
586  integer, intent(in) :: gIA
587  integer, intent(in) :: gJA
588 
589  integer :: sendcounts, recvcounts
590  integer :: ierr
591  !---------------------------------------------------------------------------
592 
593  sendcounts = gia * gja
594  recvcounts = gia * gja
595 
596  call mpi_gather( send(:,:), &
597  sendcounts, &
598  comm_datatype, &
599  recv(:,:,:), &
600  recvcounts, &
601  comm_datatype, &
602  prc_masterrank, &
603  comm_world, &
604  ierr )
605 
606  return
607  end subroutine comm_gather_2d
608 
609  !-----------------------------------------------------------------------------
611  subroutine comm_gather_3d( recv, send, gIA, gJA, gKA )
612  use scale_prc, only: &
614  implicit none
615 
616  real(RP), intent(out) :: recv(:,:,:,:)
617  real(RP), intent(in) :: send(:,:,:)
618  integer, intent(in) :: gIA
619  integer, intent(in) :: gJA
620  integer, intent(in) :: gKA
621 
622  integer :: sendcounts, recvcounts
623  integer :: ierr
624  !---------------------------------------------------------------------------
625 
626  sendcounts = gia * gja * gka
627  recvcounts = gia * gja * gka
628 
629  call mpi_gather( send(:,:,:), &
630  sendcounts, &
631  comm_datatype, &
632  recv(:,:,:,:), &
633  recvcounts, &
634  comm_datatype, &
635  prc_masterrank, &
636  comm_world, &
637  ierr )
638 
639  return
640  end subroutine comm_gather_3d
641 
642  !-----------------------------------------------------------------------------
644  subroutine comm_bcast_scr( var )
645  use scale_prc, only: &
647  implicit none
648 
649  real(RP), intent(inout) :: var
650 
651  integer :: counts
652  integer :: ierr
653  !---------------------------------------------------------------------------
654 
655  call prof_rapstart('COMM_Bcast', 2)
656 
657  counts = 1
658 
659  call mpi_bcast( var, &
660  counts, &
661  comm_datatype, &
662  prc_masterrank, &
663  comm_world, &
664  ierr )
665 
666  call prof_rapend('COMM_Bcast', 2)
667 
668  return
669  end subroutine comm_bcast_scr
670 
671  !-----------------------------------------------------------------------------
673  subroutine comm_bcast_1d( var, gIA )
674  use scale_prc, only: &
676  implicit none
677 
678  real(RP), intent(inout) :: var(:)
679  integer, intent(in) :: gIA
680 
681  integer :: counts
682  integer :: ierr
683  !---------------------------------------------------------------------------
684 
685  call prof_rapstart('COMM_Bcast', 2)
686 
687  counts = gia
688 
689  call mpi_bcast( var(:), &
690  counts, &
691  comm_datatype, &
692  prc_masterrank, &
693  comm_world, &
694  ierr )
695 
696  call prof_rapend('COMM_Bcast', 2)
697 
698  return
699  end subroutine comm_bcast_1d
700 
701  !-----------------------------------------------------------------------------
703  subroutine comm_bcast_2d( var, gIA, gJA )
704  use scale_prc, only: &
706  implicit none
707 
708  real(RP), intent(inout) :: var(:,:)
709  integer, intent(in) :: gIA
710  integer, intent(in) :: gJA
711 
712  integer :: counts
713  integer :: ierr
714  !---------------------------------------------------------------------------
715 
716  call prof_rapstart('COMM_Bcast', 2)
717 
718  counts = gia * gja
719 
720  call mpi_bcast( var(:,:), &
721  counts, &
722  comm_datatype, &
723  prc_masterrank, &
724  comm_world, &
725  ierr )
726 
727  call prof_rapend('COMM_Bcast', 2)
728 
729  return
730  end subroutine comm_bcast_2d
731 
732  !-----------------------------------------------------------------------------
734  subroutine comm_bcast_3d( var, gIA, gJA, gKA )
735  use scale_prc, only: &
737  implicit none
738 
739  real(RP), intent(inout) :: var(:,:,:)
740  integer, intent(in) :: gIA
741  integer, intent(in) :: gJA
742  integer, intent(in) :: gKA
743 
744  integer :: counts
745  integer :: ierr
746  !---------------------------------------------------------------------------
747 
748  call prof_rapstart('COMM_Bcast', 2)
749 
750  counts = gia * gja * gka
751 
752  call mpi_bcast( var(:,:,:), &
753  counts, &
754  comm_datatype, &
755  prc_masterrank, &
756  comm_world, &
757  ierr )
758 
759  call prof_rapend('COMM_Bcast', 2)
760 
761  return
762  end subroutine comm_bcast_3d
763 
764  !-----------------------------------------------------------------------------
766  subroutine comm_bcast_4d( var, gIA, gJA, gKA, gTime )
767  use scale_prc, only: &
769  implicit none
770 
771  real(RP), intent(inout) :: var(:,:,:,:)
772  integer, intent(in) :: gIA
773  integer, intent(in) :: gJA
774  integer, intent(in) :: gKA
775  integer, intent(in) :: gTime
776 
777  integer :: counts
778  integer :: ierr
779  !---------------------------------------------------------------------------
780 
781  call prof_rapstart('COMM_Bcast', 2)
782 
783  counts = gia * gja * gka * gtime
784  if ( gia>0 .AND. gja>0 .AND. gka>0 .AND. gtime>0 .AND. &
785  counts < 0 ) then
786  log_error("COMM_bcast_4D",*) 'counts overflow'
787  call prc_abort
788  end if
789 
790  call mpi_bcast( var(:,:,:,:), &
791  counts, &
792  comm_datatype, &
793  prc_masterrank, &
794  comm_world, &
795  ierr )
796 
797  call prof_rapend('COMM_Bcast', 2)
798 
799  return
800  end subroutine comm_bcast_4d
801 
802  !-----------------------------------------------------------------------------
804  subroutine comm_bcast_int_scr( var )
805  use scale_prc, only: &
807  implicit none
808 
809  integer, intent(inout) :: var
810 
811  integer :: counts
812  integer :: ierr
813  !---------------------------------------------------------------------------
814 
815  call prof_rapstart('COMM_Bcast', 2)
816 
817  counts = 1
818 
819  call mpi_bcast( var, &
820  counts, &
821  mpi_integer, &
822  prc_masterrank, &
823  comm_world, &
824  ierr )
825 
826  call prof_rapend('COMM_Bcast', 2)
827 
828  return
829  end subroutine comm_bcast_int_scr
830 
831  !-----------------------------------------------------------------------------
833  subroutine comm_bcast_int_1d( var, gIA )
834  use scale_prc, only: &
836  implicit none
837 
838  integer, intent(inout) :: var(:)
839  integer, intent(in) :: gIA
840 
841  integer :: counts
842  integer :: ierr
843  !---------------------------------------------------------------------------
844 
845  call prof_rapstart('COMM_Bcast', 2)
846 
847  counts = gia
848 
849  call mpi_bcast( var(:), &
850  counts, &
851  mpi_integer, &
852  prc_masterrank, &
853  comm_world, &
854  ierr )
855 
856  call prof_rapend('COMM_Bcast', 2)
857 
858  return
859  end subroutine comm_bcast_int_1d
860 
861  !-----------------------------------------------------------------------------
863  subroutine comm_bcast_int_2d( var, gIA, gJA )
864  use scale_prc, only: &
866  implicit none
867 
868  integer, intent(inout) :: var(:,:)
869  integer, intent(in) :: gIA
870  integer, intent(in) :: gJA
871 
872  integer :: counts
873  integer :: ierr
874  !---------------------------------------------------------------------------
875 
876  call prof_rapstart('COMM_Bcast', 2)
877 
878  counts = gia * gja
879 
880  call mpi_bcast( var(:,:), &
881  counts, &
882  mpi_integer, &
883  prc_masterrank, &
884  comm_world, &
885  ierr )
886 
887  call prof_rapend('COMM_Bcast', 2)
888 
889  return
890  end subroutine comm_bcast_int_2d
891 
892  !-----------------------------------------------------------------------------
894  subroutine comm_bcast_logical_scr( var )
895  use scale_prc, only: &
897  implicit none
898 
899  logical, intent(inout) :: var
900 
901  integer :: counts
902  integer :: ierr
903  !---------------------------------------------------------------------------
904 
905  call prof_rapstart('COMM_Bcast', 2)
906 
907  counts = 1
908 
909  call mpi_bcast( var, &
910  counts, &
911  mpi_logical, &
912  prc_masterrank, &
913  comm_world, &
914  ierr )
915 
916  call prof_rapend('COMM_Bcast', 2)
917 
918  return
919  end subroutine comm_bcast_logical_scr
920 
921  !-----------------------------------------------------------------------------
923  subroutine comm_bcast_character( var )
924  use scale_prc, only: &
926  implicit none
927 
928  character(len=*), intent(inout) :: var
929 
930  integer :: counts
931  integer :: ierr
932  !---------------------------------------------------------------------------
933 
934  call prof_rapstart('COMM_Bcast', 2)
935 
936  counts = len(var)
937 
938  call mpi_bcast( var, &
939  counts, &
940  mpi_character, &
941  prc_masterrank, &
942  comm_world, &
943  ierr )
944 
945  call prof_rapend('COMM_Bcast', 2)
946 
947  return
948  end subroutine comm_bcast_character
949 
950 !-------------------------------------------------------------------------------
951 ! private routines
952 !-------------------------------------------------------------------------------
953  subroutine vars_init_mpi_pc(var, vid, seqid)
954  implicit none
955 
956  real(RP), intent(inout) :: var(:,:,:)
957  integer, intent(in) :: vid
958  integer, intent(in) :: seqid
959 
960  integer :: ireq, tag, ierr
961  logical :: flag
962 
963  integer :: kd
964  integer :: i
965 
966  tag = vid * 100
967  ireq = 1
968 
969  kd = size(var, 1)
970 
971  ! register whole array to inner table of MPI and/or lower library
972  ! otherwise a lot of sub small segments would be registered
973  call mpi_send_init( var(:,:,:), size(var), comm_datatype, &
974  mpi_proc_null, tag+comm_nreq_max+1, comm_world, &
975  preq_list(comm_nreq_max+1,vid), ierr )
976 
977  !--- From 4-Direction HALO communicate
978  ! From S
979  call mpi_recv_init( var(:,:,js-jhalo:js-1), comm_size2d_ns4*kd, comm_datatype, &
980  prc_next(prc_s), tag+1, comm_world, preq_list(ireq,vid), ierr )
981  ireq = ireq + 1
982  ! From N
983  call mpi_recv_init( var(:,:,je+1:je+jhalo), comm_size2d_ns4*kd, comm_datatype, &
984  prc_next(prc_n), tag+2, comm_world, preq_list(ireq,vid), ierr )
985  ireq = ireq + 1
986  ! From E
987  call mpi_recv_init( recvpack_e2p(:,seqid), comm_size2d_we*kd, comm_datatype, &
988  prc_next(prc_e), tag+3, comm_world, preq_list(ireq,vid), ierr )
989  ireq = ireq + 1
990  ! From W
991  call mpi_recv_init( recvpack_w2p(:,seqid), comm_size2d_we*kd, comm_datatype, &
992  prc_next(prc_w), tag+4, comm_world, preq_list(ireq,vid), ierr )
993  ireq = ireq + 1
994 
995  !--- To 4-Direction HALO communicate
996  ! To W HALO
997  call mpi_send_init( sendpack_p2w(:,seqid), comm_size2d_we*kd, comm_datatype, &
998  prc_next(prc_w), tag+3, comm_world, preq_list(ireq,vid), ierr )
999  ireq = ireq + 1
1000  ! To E HALO
1001  call mpi_send_init( sendpack_p2e(:,seqid), comm_size2d_we*kd, comm_datatype, &
1002  prc_next(prc_e), tag+4, comm_world, preq_list(ireq,vid), ierr )
1003  ireq = ireq + 1
1004  ! To N HALO
1005  call mpi_send_init( var(:,:,je-jhalo+1:je), comm_size2d_ns4*kd, comm_datatype, &
1006  prc_next(prc_n), tag+1, comm_world, preq_list(ireq,vid), ierr )
1007  ireq = ireq + 1
1008  ! To S HALO
1009  call mpi_send_init( var(:,:,js:js+jhalo-1), comm_size2d_ns4*kd, comm_datatype, &
1010  prc_next(prc_s), tag+2, comm_world, preq_list(ireq,vid), ierr )
1011  ireq = ireq + 1
1012 
1013  preq_cnt(vid) = ireq - 1
1014  pseqid(vid) = seqid
1015 
1016  ! to finish initial processes of MPIa
1017  do i = 1, 32
1018  call mpi_testall( preq_cnt(vid), preq_list(1:preq_cnt(vid),vid), &
1019  flag, mpi_statuses_ignore, ierr )
1020  enddo
1021 
1022  return
1023  end subroutine vars_init_mpi_pc
1024 
1025  subroutine vars8_init_mpi_pc(var, vid, seqid)
1026  implicit none
1027 
1028  real(RP), intent(inout) :: var(:,:,:)
1029  integer, intent(in) :: vid
1030  integer, intent(in) :: seqid
1031 
1032  integer :: ireq, tag, tagc
1033  integer :: ierr
1034  integer :: kd
1035  logical :: flag
1036 
1037  integer :: i, j
1038 
1039  kd = size(var, 1)
1040 
1041  tag = vid * 100
1042  ireq = 1
1043 
1044  ! register whole array to inner table of MPI and/or lower library
1045  ! otherwise a lot of sub small segments would be registered
1046  call mpi_send_init( var(:,:,:), size(var), comm_datatype, &
1047  mpi_proc_null, tag+comm_nreq_max+1, comm_world, &
1048  preq_list(comm_nreq_max+1,vid), ierr )
1049 
1050 
1051  if ( comm_isallperiodic ) then ! periodic condition
1052 
1053  !--- From 8-Direction HALO communicate
1054  ! From SE
1055  tagc = 0
1056  do j = js-jhalo, js-1
1057  call mpi_recv_init( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
1058  prc_next(prc_se), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1059  ireq = ireq + 1
1060  tagc = tagc + 1
1061  enddo
1062  ! From SW
1063  tagc = 10
1064  do j = js-jhalo, js-1
1065  call mpi_recv_init( var(1,is-ihalo,j), comm_size2d_4c*kd, comm_datatype, &
1066  prc_next(prc_sw), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1067  ireq = ireq + 1
1068  tagc = tagc + 1
1069  enddo
1070  ! From NE
1071  tagc = 20
1072  do j = je+1, je+jhalo
1073  call mpi_recv_init( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
1074  prc_next(prc_ne), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1075  ireq = ireq + 1
1076  tagc = tagc + 1
1077  enddo
1078  ! From NW
1079  tagc = 30
1080  do j = je+1, je+jhalo
1081  call mpi_recv_init( var(1,is-ihalo,j), comm_size2d_4c*kd, comm_datatype, &
1082  prc_next(prc_nw), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1083  ireq = ireq + 1
1084  tagc = tagc + 1
1085  enddo
1086  ! From S
1087  tagc = 40
1088  do j = js-jhalo, js-1
1089  call mpi_recv_init( var(1,is,j), comm_size2d_ns8*kd, comm_datatype, &
1090  prc_next(prc_s), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1091  ireq = ireq + 1
1092  tagc = tagc + 1
1093  enddo
1094  ! From N
1095  tagc = 50
1096  do j = je+1, je+jhalo
1097  call mpi_recv_init( var(1,is,j), comm_size2d_ns8*kd, comm_datatype, &
1098  prc_next(prc_n), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1099  ireq = ireq + 1
1100  tagc = tagc + 1
1101  enddo
1102  ! From E
1103  tagc = 60
1104  call mpi_recv_init( recvpack_e2p(:,seqid), comm_size2d_we*kd, comm_datatype, &
1105  prc_next(prc_e), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1106  ireq = ireq + 1
1107  ! From W
1108  tagc = 70
1109  call mpi_recv_init( recvpack_w2p(:,seqid), comm_size2d_we*kd, comm_datatype, &
1110  prc_next(prc_w), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1111  ireq = ireq + 1
1112 
1113  !--- To 8-Direction HALO communicate
1114  ! To W HALO
1115  tagc = 60
1116  call mpi_send_init( sendpack_p2w(:,seqid), comm_size2d_we*kd, comm_datatype, &
1117  prc_next(prc_w), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1118  ireq = ireq + 1
1119  ! To E HALO
1120  tagc = 70
1121  call mpi_send_init( sendpack_p2e(:,seqid), comm_size2d_we*kd, comm_datatype, &
1122  prc_next(prc_e), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1123  ireq = ireq + 1
1124  ! To N HALO
1125  tagc = 40
1126  do j = je-jhalo+1, je
1127  call mpi_send_init( var(1,is,j), comm_size2d_ns8*kd, comm_datatype, &
1128  prc_next(prc_n), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1129  ireq = ireq + 1
1130  tagc = tagc + 1
1131  enddo
1132  ! To S HALO
1133  tagc = 50
1134  do j = js, js+jhalo-1
1135  call mpi_send_init( var(1,is,j), comm_size2d_ns8*kd, comm_datatype, &
1136  prc_next(prc_s), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1137  ireq = ireq + 1
1138  tagc = tagc + 1
1139  enddo
1140  ! To NW HALO
1141  tagc = 0
1142  do j = je-jhalo+1, je
1143  call mpi_send_init( var(1,is,j), comm_size2d_4c*kd, comm_datatype, &
1144  prc_next(prc_nw), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1145  ireq = ireq + 1
1146  tagc = tagc + 1
1147  enddo
1148  ! To NE HALO
1149  tagc = 10
1150  do j = je-jhalo+1, je
1151  call mpi_send_init( var(1,ie-ihalo+1,j), comm_size2d_4c*kd, comm_datatype, &
1152  prc_next(prc_ne), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1153  ireq = ireq + 1
1154  tagc = tagc + 1
1155  enddo
1156  ! To SW HALO
1157  tagc = 20
1158  do j = js, js+jhalo-1
1159  call mpi_send_init( var(1,is,j), comm_size2d_4c*kd, comm_datatype, &
1160  prc_next(prc_sw), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1161  ireq = ireq + 1
1162  tagc = tagc + 1
1163  enddo
1164  ! To SE HALO
1165  tagc = 30
1166  do j = js, js+jhalo-1
1167  call mpi_send_init( var(1,ie-ihalo+1,j), comm_size2d_4c*kd, comm_datatype, &
1168  prc_next(prc_se), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1169  ireq = ireq + 1
1170  tagc = tagc + 1
1171  enddo
1172 
1173  else ! non-periodic condition
1174 
1175  !--- From 8-Direction HALO communicate
1176  ! From SE
1177  if ( prc_has_s .AND. prc_has_e ) then
1178  tagc = 0
1179  do j = js-jhalo, js-1
1180  call mpi_recv_init( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
1181  prc_next(prc_se), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1182  ireq = ireq + 1
1183  tagc = tagc + 1
1184  enddo
1185  else if ( prc_has_s ) then
1186  tagc = 0
1187  do j = js-jhalo, js-1
1188  call mpi_recv_init( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
1189  prc_next(prc_s), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1190  ireq = ireq + 1
1191  tagc = tagc + 1
1192  enddo
1193  else if ( prc_has_e ) then
1194  tagc = 0
1195  do j = js-jhalo, js-1
1196  call mpi_recv_init( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
1197  prc_next(prc_e), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1198  ireq = ireq + 1
1199  tagc = tagc + 1
1200  enddo
1201  endif
1202  ! From SW
1203  if ( prc_has_s .AND. prc_has_w ) then
1204  tagc = 10
1205  do j = js-jhalo, js-1
1206  call mpi_recv_init( var(1,is-ihalo,j), comm_size2d_4c*kd, comm_datatype, &
1207  prc_next(prc_sw), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1208  ireq = ireq + 1
1209  tagc = tagc + 1
1210  enddo
1211  else if ( prc_has_s ) then
1212  tagc = 10
1213  do j = js-jhalo, js-1
1214  call mpi_recv_init( var(1,is-ihalo,j), comm_size2d_4c*kd, comm_datatype, &
1215  prc_next(prc_s), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1216  ireq = ireq + 1
1217  tagc = tagc + 1
1218  enddo
1219  else if ( prc_has_w ) then
1220  tagc = 10
1221  do j = js-jhalo, js-1
1222  call mpi_recv_init( var(1,is-ihalo,j), comm_size2d_4c*kd, comm_datatype, &
1223  prc_next(prc_w), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1224  ireq = ireq + 1
1225  tagc = tagc + 1
1226  enddo
1227  endif
1228  ! From NE
1229  if ( prc_has_n .AND. prc_has_e ) then
1230  tagc = 20
1231  do j = je+1, je+jhalo
1232  call mpi_recv_init( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
1233  prc_next(prc_ne), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1234  ireq = ireq + 1
1235  tagc = tagc + 1
1236  enddo
1237  else if ( prc_has_n ) then
1238  tagc = 20
1239  do j = je+1, je+jhalo
1240  call mpi_recv_init( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
1241  prc_next(prc_n), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1242  ireq = ireq + 1
1243  tagc = tagc + 1
1244  enddo
1245  else if ( prc_has_e ) then
1246  tagc = 20
1247  do j = je+1, je+jhalo
1248  call mpi_recv_init( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
1249  prc_next(prc_e), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1250  ireq = ireq + 1
1251  tagc = tagc + 1
1252  enddo
1253  endif
1254  ! From NW
1255  if ( prc_has_n .AND. prc_has_w ) then
1256  tagc = 30
1257  do j = je+1, je+jhalo
1258  call mpi_recv_init( var(1,is-ihalo,j), comm_size2d_4c*kd, comm_datatype, &
1259  prc_next(prc_nw), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1260  ireq = ireq + 1
1261  tagc = tagc + 1
1262  enddo
1263  else if ( prc_has_n ) then
1264  tagc = 30
1265  do j = je+1, je+jhalo
1266  call mpi_recv_init( var(1,is-ihalo,j), comm_size2d_4c*kd, comm_datatype, &
1267  prc_next(prc_n), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1268  ireq = ireq + 1
1269  tagc = tagc + 1
1270  enddo
1271  else if ( prc_has_w ) then
1272  tagc = 30
1273  do j = je+1, je+jhalo
1274  call mpi_recv_init( var(1,is-ihalo,j), comm_size2d_4c*kd, comm_datatype, &
1275  prc_next(prc_w), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1276  ireq = ireq + 1
1277  tagc = tagc + 1
1278  enddo
1279  endif
1280  ! From S
1281  if ( prc_has_s ) then
1282  tagc = 40
1283  do j = js-jhalo, js-1
1284  call mpi_recv_init( var(1,is,j), comm_size2d_ns8*kd, comm_datatype, &
1285  prc_next(prc_s), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1286  ireq = ireq + 1
1287  tagc = tagc + 1
1288  enddo
1289  endif
1290  ! From N
1291  if ( prc_has_n ) then
1292  tagc = 50
1293  do j = je+1, je+jhalo
1294  call mpi_recv_init( var(1,is,j), comm_size2d_ns8*kd, comm_datatype, &
1295  prc_next(prc_n), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1296  ireq = ireq + 1
1297  tagc = tagc + 1
1298  enddo
1299  endif
1300  ! From E
1301  if ( prc_has_e ) then
1302  tagc = 60
1303  call mpi_recv_init( recvpack_e2p(:,seqid), comm_size2d_we*kd, comm_datatype, &
1304  prc_next(prc_e), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1305  ireq = ireq + 1
1306  endif
1307  ! From W
1308  if ( prc_has_w ) then
1309  tagc = 70
1310  call mpi_recv_init( recvpack_w2p(:,seqid), comm_size2d_we*kd, comm_datatype, &
1311  prc_next(prc_w), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1312  ireq = ireq + 1
1313  endif
1314 
1315  !--- To 8-Direction HALO communicate
1316  ! To W HALO
1317  if ( prc_has_w ) then
1318  tagc = 60
1319  call mpi_send_init( sendpack_p2w(:,seqid), comm_size2d_we*kd, comm_datatype, &
1320  prc_next(prc_w), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1321  ireq = ireq + 1
1322  endif
1323  ! To E HALO
1324  if ( prc_has_e ) then
1325  tagc = 70
1326  call mpi_send_init( sendpack_p2e(:,seqid), comm_size2d_we*kd, comm_datatype, &
1327  prc_next(prc_e), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1328  ireq = ireq + 1
1329  endif
1330  ! To N HALO
1331  if ( prc_has_n ) then
1332  tagc = 40
1333  do j = je-jhalo+1, je
1334  call mpi_send_init( var(1,is,j), comm_size2d_ns8*kd, comm_datatype, &
1335  prc_next(prc_n), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1336  ireq = ireq + 1
1337  tagc = tagc + 1
1338  enddo
1339  endif
1340  ! To S HALO
1341  if ( prc_has_s ) then
1342  tagc = 50
1343  do j = js, js+jhalo-1
1344  call mpi_send_init( var(1,is,j), comm_size2d_ns8*kd, comm_datatype, &
1345  prc_next(prc_s), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1346  ireq = ireq + 1
1347  tagc = tagc + 1
1348  enddo
1349  endif
1350  ! To NW HALO
1351  if ( prc_has_n .AND. prc_has_w ) then
1352  tagc = 0
1353  do j = je-jhalo+1, je
1354  call mpi_send_init( var(1,is,j), comm_size2d_4c*kd, comm_datatype, &
1355  prc_next(prc_nw), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1356  ireq = ireq + 1
1357  tagc = tagc + 1
1358  enddo
1359  else if ( prc_has_n ) then
1360  tagc = 10
1361  do j = je-jhalo+1, je
1362  call mpi_send_init( var(1,is-ihalo,j), comm_size2d_4c*kd, comm_datatype, &
1363  prc_next(prc_n), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1364  ireq = ireq + 1
1365  tagc = tagc + 1
1366  enddo
1367  else if ( prc_has_w ) then
1368  tagc = 20
1369  do j = je+1, je+jhalo
1370  call mpi_send_init( var(1,is,j), comm_size2d_4c*kd, comm_datatype, &
1371  prc_next(prc_w), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1372  ireq = ireq + 1
1373  tagc = tagc + 1
1374  enddo
1375  endif
1376  ! To NE HALO
1377  if ( prc_has_n .AND. prc_has_e ) then
1378  tagc = 10
1379  do j = je-jhalo+1, je
1380  call mpi_send_init( var(1,ie-ihalo+1,j), comm_size2d_4c*kd, comm_datatype, &
1381  prc_next(prc_ne), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1382  ireq = ireq + 1
1383  tagc = tagc + 1
1384  enddo
1385  else if ( prc_has_n ) then
1386  tagc = 0
1387  do j = je-jhalo+1, je
1388  call mpi_send_init( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
1389  prc_next(prc_n), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1390  ireq = ireq + 1
1391  tagc = tagc + 1
1392  enddo
1393  else if ( prc_has_e ) then
1394  tagc = 30
1395  do j = je+1, je+jhalo
1396  call mpi_send_init( var(1,ie-ihalo+1,j), comm_size2d_4c*kd, comm_datatype, &
1397  prc_next(prc_e), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1398  ireq = ireq + 1
1399  tagc = tagc + 1
1400  enddo
1401  endif
1402  ! To SW HALO
1403  if ( prc_has_s .AND. prc_has_w ) then
1404  tagc = 20
1405  do j = js, js+jhalo-1
1406  call mpi_send_init( var(1,is,j), comm_size2d_4c*kd, comm_datatype, &
1407  prc_next(prc_sw), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1408  ireq = ireq + 1
1409  tagc = tagc + 1
1410  enddo
1411  else if ( prc_has_s ) then
1412  tagc = 30
1413  do j = js, js+jhalo-1
1414  call mpi_send_init( var(1,is-ihalo,j), comm_size2d_4c*kd, comm_datatype, &
1415  prc_next(prc_s), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1416  ireq = ireq + 1
1417  tagc = tagc + 1
1418  enddo
1419  else if ( prc_has_w ) then
1420  tagc = 0
1421  do j = js-jhalo, js-1
1422  call mpi_send_init( var(1,is,j), comm_size2d_4c*kd, comm_datatype, &
1423  prc_next(prc_w), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1424  ireq = ireq + 1
1425  tagc = tagc + 1
1426  enddo
1427  endif
1428  ! To SE HALO
1429  if ( prc_has_s .AND. prc_has_e ) then
1430  tagc = 30
1431  do j = js, js+jhalo-1
1432  call mpi_send_init( var(1,ie-ihalo+1,j), comm_size2d_4c*kd, comm_datatype, &
1433  prc_next(prc_se), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1434  ireq = ireq + 1
1435  tagc = tagc + 1
1436  enddo
1437  else if ( prc_has_s ) then
1438  tagc = 20
1439  do j = js, js+jhalo-1
1440  call mpi_send_init( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
1441  prc_next(prc_s), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1442  ireq = ireq + 1
1443  tagc = tagc + 1
1444  enddo
1445  else if ( prc_has_e ) then
1446  tagc = 10
1447  do j = js-jhalo, js-1
1448  call mpi_send_init( var(1,ie-ihalo+1,j), comm_size2d_4c*kd, comm_datatype, &
1449  prc_next(prc_e), tag+tagc, comm_world, preq_list(ireq,vid), ierr )
1450  ireq = ireq + 1
1451  tagc = tagc + 1
1452  enddo
1453  endif
1454 
1455  endif
1456 
1457  preq_cnt(vid) = ireq - 1
1458  pseqid(vid) = seqid
1459 
1460  ! to finish initial processes of MPIa
1461  do i = 1, 32
1462  call mpi_testall( preq_cnt(vid), preq_list(1:preq_cnt(vid),vid), &
1463  flag, mpi_statuses_ignore, ierr )
1464  enddo
1465 
1466  return
1467  end subroutine vars8_init_mpi_pc
1468 
1469  subroutine vars_3d_mpi(var, vid)
1470  use scale_prc, only: &
1471  prc_abort
1472  implicit none
1473 
1474  real(RP), intent(inout) :: var(:,:,:)
1475  integer, intent(in) :: vid
1476 
1477 
1478  integer :: ireq, tag
1479 
1480  integer :: kd
1481  integer :: ierr
1482  !---------------------------------------------------------------------------
1483 
1484  tag = vid * 100
1485  ireq = 1
1486 
1487  kd = size(var, 1)
1488 
1489 #ifdef DEBUG
1490  if ( use_packbuf(vid) ) then
1491  log_error("vars_3D_mpi",*) 'packing buffer is already used', vid
1492  call prc_abort
1493  end if
1494  use_packbuf(vid) = .true.
1495 #endif
1496 
1497  if ( comm_isallperiodic ) then ! periodic condition
1498 
1499  !--- From 4-Direction HALO communicate
1500  ! From S
1501  call mpi_irecv( var(:,:,js-jhalo:js-1), comm_size2d_ns4*kd, comm_datatype, &
1502  prc_next(prc_s), tag+1, comm_world, req_list(ireq,vid), ierr )
1503  ireq = ireq + 1
1504  ! From N
1505  call mpi_irecv( var(:,:,je+1:je+jhalo), comm_size2d_ns4*kd, comm_datatype, &
1506  prc_next(prc_n), tag+2, comm_world, req_list(ireq,vid), ierr )
1507  ireq = ireq + 1
1508  ! From E
1509  call mpi_irecv( recvpack_e2p(:,vid), comm_size2d_we*kd, comm_datatype, &
1510  prc_next(prc_e), tag+3, comm_world, req_list(ireq,vid), ierr )
1511  ireq = ireq + 1
1512  ! From W
1513  call mpi_irecv( recvpack_w2p(:,vid), comm_size2d_we*kd, comm_datatype, &
1514  prc_next(prc_w), tag+4, comm_world, req_list(ireq,vid), ierr )
1515  ireq = ireq + 1
1516 
1517  call pack_3d(var, vid)
1518 
1519  !--- To 4-Direction HALO communicate
1520  ! To W HALO
1521  call mpi_isend( sendpack_p2w(:,vid), comm_size2d_we*kd, comm_datatype, &
1522  prc_next(prc_w), tag+3, comm_world, req_list(ireq,vid), ierr )
1523  ireq = ireq + 1
1524  ! To E HALO
1525  call mpi_isend( sendpack_p2e(:,vid), comm_size2d_we*kd, comm_datatype, &
1526  prc_next(prc_e), tag+4, comm_world, req_list(ireq,vid), ierr )
1527  ireq = ireq + 1
1528  ! To N HALO
1529  call mpi_isend( var(:,:,je-jhalo+1:je), comm_size2d_ns4*kd, comm_datatype, &
1530  prc_next(prc_n), tag+1, comm_world, req_list(ireq,vid), ierr )
1531  ireq = ireq + 1
1532  ! To S HALO
1533  call mpi_isend( var(:,:,js:js+jhalo-1), comm_size2d_ns4*kd, comm_datatype, &
1534  prc_next(prc_s), tag+2, comm_world, req_list(ireq,vid), ierr )
1535  ireq = ireq + 1
1536 
1537  else ! non-periodic condition
1538 
1539  !--- From 4-Direction HALO communicate
1540  ! From S
1541  if ( prc_has_s ) then
1542  call mpi_irecv( var(:,:,js-jhalo:js-1), comm_size2d_ns4*kd, comm_datatype, &
1543  prc_next(prc_s), tag+1, comm_world, req_list(ireq,vid), ierr )
1544  ireq = ireq + 1
1545  endif
1546  ! From N
1547  if ( prc_has_n ) then
1548  call mpi_irecv( var(:,:,je+1:je+jhalo), comm_size2d_ns4*kd, comm_datatype, &
1549  prc_next(prc_n), tag+2, comm_world, req_list(ireq,vid), ierr )
1550  ireq = ireq + 1
1551  endif
1552  ! From E
1553  if ( prc_has_e ) then
1554  call mpi_irecv( recvpack_e2p(:,vid), comm_size2d_we*kd, comm_datatype, &
1555  prc_next(prc_e), tag+3, comm_world, req_list(ireq,vid), ierr )
1556  ireq = ireq + 1
1557  endif
1558  ! From W
1559  if ( prc_has_w ) then
1560  call mpi_irecv( recvpack_w2p(:,vid), comm_size2d_we*kd, comm_datatype, &
1561  prc_next(prc_w), tag+4, comm_world, req_list(ireq,vid), ierr )
1562  ireq = ireq + 1
1563  endif
1564 
1565  call pack_3d(var, vid)
1566 
1567  !--- To 4-Direction HALO communicate
1568  ! To W HALO
1569  if ( prc_has_w ) then
1570  call mpi_isend( sendpack_p2w(:,vid), comm_size2d_we*kd, comm_datatype, &
1571  prc_next(prc_w), tag+3, comm_world, req_list(ireq,vid), ierr )
1572  ireq = ireq + 1
1573  endif
1574  ! To E HALO
1575  if ( prc_has_e ) then
1576  call mpi_isend( sendpack_p2e(:,vid), comm_size2d_we*kd, comm_datatype, &
1577  prc_next(prc_e), tag+4, comm_world, req_list(ireq,vid), ierr )
1578  ireq = ireq + 1
1579  endif
1580  ! To N HALO
1581  if ( prc_has_n ) then
1582  call mpi_isend( var(:,:,je-jhalo+1:je), comm_size2d_ns4*kd, comm_datatype, &
1583  prc_next(prc_n), tag+1, comm_world, req_list(ireq,vid), ierr )
1584  ireq = ireq + 1
1585  endif
1586  ! To S HALO
1587  if ( prc_has_s ) then
1588  call mpi_isend( var(:,:,js:js+jhalo-1), comm_size2d_ns4*kd, comm_datatype, &
1589  prc_next(prc_s), tag+2, comm_world, req_list(ireq,vid), ierr )
1590  ireq = ireq + 1
1591  endif
1592 
1593  endif
1594 
1595  req_cnt(vid) = ireq - 1
1596 
1597  return
1598  end subroutine vars_3d_mpi
1599 
1600  subroutine vars8_3d_mpi(var, vid)
1601  use scale_prc, only: &
1602  prc_abort
1603  implicit none
1604 
1605  real(RP), intent(inout) :: var(:,:,:)
1606  integer, intent(in) :: vid
1607 
1608  integer :: ireq, tag, tagc
1609 
1610  integer :: kd
1611 
1612  integer :: ierr
1613  integer :: j
1614  !---------------------------------------------------------------------------
1615 
1616  tag = vid * 100
1617  ireq = 1
1618 
1619  kd = size(var, 1)
1620 
1621 #ifdef DEBUG
1622  if ( use_packbuf(vid) ) then
1623  log_error("vars8_3D_mpi",*) 'packing buffer is already used', vid
1624  call prc_abort
1625  end if
1626  use_packbuf(vid) = .true.
1627 #endif
1628 
1629  if ( comm_isallperiodic ) then ! periodic condition
1630 
1631  !--- From 8-Direction HALO communicate
1632  ! From SE
1633  tagc = 0
1634  do j = js-jhalo, js-1
1635  call mpi_irecv( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
1636  prc_next(prc_se), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1637  ireq = ireq + 1
1638  tagc = tagc + 1
1639  enddo
1640  ! From SW
1641  tagc = 10
1642  do j = js-jhalo, js-1
1643  call mpi_irecv( var(1,is-ihalo,j), comm_size2d_4c*kd, comm_datatype, &
1644  prc_next(prc_sw), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1645  ireq = ireq + 1
1646  tagc = tagc + 1
1647  enddo
1648  ! From NE
1649  tagc = 20
1650  do j = je+1, je+jhalo
1651  call mpi_irecv( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
1652  prc_next(prc_ne), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1653  ireq = ireq + 1
1654  tagc = tagc + 1
1655  enddo
1656  ! From NW
1657  tagc = 30
1658  do j = je+1, je+jhalo
1659  call mpi_irecv( var(1,is-ihalo,j), comm_size2d_4c*kd, comm_datatype, &
1660  prc_next(prc_nw), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1661  ireq = ireq + 1
1662  tagc = tagc + 1
1663  enddo
1664  ! From S
1665  tagc = 40
1666  do j = js-jhalo, js-1
1667  call mpi_irecv( var(1,is,j), comm_size2d_ns8*kd, comm_datatype, &
1668  prc_next(prc_s), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1669  ireq = ireq + 1
1670  tagc = tagc + 1
1671  enddo
1672  ! From N
1673  tagc = 50
1674  do j = je+1, je+jhalo
1675  call mpi_irecv( var(1,is,j), comm_size2d_ns8*kd, comm_datatype, &
1676  prc_next(prc_n), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1677  ireq = ireq + 1
1678  tagc = tagc + 1
1679  enddo
1680  ! From E
1681  tagc = 60
1682  call mpi_irecv( recvpack_e2p(:,vid), comm_size2d_we*kd, comm_datatype, &
1683  prc_next(prc_e), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1684  ireq = ireq + 1
1685  ! From W
1686  tagc = 70
1687  call mpi_irecv( recvpack_w2p(:,vid), comm_size2d_we*kd, comm_datatype, &
1688  prc_next(prc_w), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1689  ireq = ireq + 1
1690 
1691  call pack_3d(var, vid)
1692 
1693 
1694  !--- To 8-Direction HALO communicate
1695  ! To W HALO
1696  tagc = 60
1697  call mpi_isend( sendpack_p2w(:,vid), comm_size2d_we*kd, comm_datatype, &
1698  prc_next(prc_w), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1699  ireq = ireq + 1
1700  ! To E HALO
1701  tagc = 70
1702  call mpi_isend( sendpack_p2e(:,vid), comm_size2d_we*kd, comm_datatype, &
1703  prc_next(prc_e), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1704  ireq = ireq + 1
1705  ! To N HALO
1706  tagc = 40
1707  do j = je-jhalo+1, je
1708  call mpi_isend( var(1,is,j), comm_size2d_ns8*kd, comm_datatype, &
1709  prc_next(prc_n), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1710  ireq = ireq + 1
1711  tagc = tagc + 1
1712  enddo
1713  ! To S HALO
1714  tagc = 50
1715  do j = js, js+jhalo-1
1716  call mpi_isend( var(1,is,j), comm_size2d_ns8*kd, comm_datatype, &
1717  prc_next(prc_s), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1718  ireq = ireq + 1
1719  tagc = tagc + 1
1720  enddo
1721  ! To NW HALO
1722  tagc = 0
1723  do j = je-jhalo+1, je
1724  call mpi_isend( var(1,is,j), comm_size2d_4c*kd, comm_datatype, &
1725  prc_next(prc_nw), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1726  ireq = ireq + 1
1727  tagc = tagc + 1
1728  enddo
1729  ! To NE HALO
1730  tagc = 10
1731  do j = je-jhalo+1, je
1732  call mpi_isend( var(1,ie-ihalo+1,j), comm_size2d_4c*kd, comm_datatype, &
1733  prc_next(prc_ne), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1734  ireq = ireq + 1
1735  tagc = tagc + 1
1736  enddo
1737  ! To SW HALO
1738  tagc = 20
1739  do j = js, js+jhalo-1
1740  call mpi_isend( var(1,is,j), comm_size2d_4c*kd, comm_datatype, &
1741  prc_next(prc_sw), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1742  ireq = ireq + 1
1743  tagc = tagc + 1
1744  enddo
1745  ! To SE HALO
1746  tagc = 30
1747  do j = js, js+jhalo-1
1748  call mpi_isend( var(1,ie-ihalo+1,j), comm_size2d_4c*kd, comm_datatype, &
1749  prc_next(prc_se), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1750  ireq = ireq + 1
1751  tagc = tagc + 1
1752  enddo
1753 
1754  else ! non-periodic condition
1755 
1756  !--- From 8-Direction HALO communicate
1757  ! From SE
1758  if ( prc_has_s .AND. prc_has_e ) then
1759  tagc = 0
1760  do j = js-jhalo, js-1
1761  call mpi_irecv( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
1762  prc_next(prc_se), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1763  ireq = ireq + 1
1764  tagc = tagc + 1
1765  enddo
1766  else if ( prc_has_s ) then
1767  tagc = 0
1768  do j = js-jhalo, js-1
1769  call mpi_irecv( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
1770  prc_next(prc_s), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1771  ireq = ireq + 1
1772  tagc = tagc + 1
1773  enddo
1774  else if ( prc_has_e ) then
1775  tagc = 0
1776  do j = js-jhalo, js-1
1777  call mpi_irecv( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
1778  prc_next(prc_e), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1779  ireq = ireq + 1
1780  tagc = tagc + 1
1781  enddo
1782  endif
1783  ! From SW
1784  if ( prc_has_s .AND. prc_has_w ) then
1785  tagc = 10
1786  do j = js-jhalo, js-1
1787  call mpi_irecv( var(1,is-ihalo,j), comm_size2d_4c*kd, comm_datatype, &
1788  prc_next(prc_sw), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1789  ireq = ireq + 1
1790  tagc = tagc + 1
1791  enddo
1792  else if ( prc_has_s ) then
1793  tagc = 10
1794  do j = js-jhalo, js-1
1795  call mpi_irecv( var(1,is-ihalo,j), comm_size2d_4c*kd, comm_datatype, &
1796  prc_next(prc_s), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1797  ireq = ireq + 1
1798  tagc = tagc + 1
1799  enddo
1800  else if ( prc_has_w ) then
1801  tagc = 10
1802  do j = js-jhalo, js-1
1803  call mpi_irecv( var(1,is-ihalo,j), comm_size2d_4c*kd, comm_datatype, &
1804  prc_next(prc_w), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1805  ireq = ireq + 1
1806  tagc = tagc + 1
1807  enddo
1808  endif
1809  ! From NE
1810  if ( prc_has_n .AND. prc_has_e ) then
1811  tagc = 20
1812  do j = je+1, je+jhalo
1813  call mpi_irecv( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
1814  prc_next(prc_ne), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1815  ireq = ireq + 1
1816  tagc = tagc + 1
1817  enddo
1818  else if ( prc_has_n ) then
1819  tagc = 20
1820  do j = je+1, je+jhalo
1821  call mpi_irecv( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
1822  prc_next(prc_n), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1823  ireq = ireq + 1
1824  tagc = tagc + 1
1825  enddo
1826  else if ( prc_has_e ) then
1827  tagc = 20
1828  do j = je+1, je+jhalo
1829  call mpi_irecv( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
1830  prc_next(prc_e), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1831  ireq = ireq + 1
1832  tagc = tagc + 1
1833  enddo
1834  endif
1835  ! From NW
1836  if ( prc_has_n .AND. prc_has_w ) then
1837  tagc = 30
1838  do j = je+1, je+jhalo
1839  call mpi_irecv( var(1,is-ihalo,j), comm_size2d_4c*kd, comm_datatype, &
1840  prc_next(prc_nw), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1841  ireq = ireq + 1
1842  tagc = tagc + 1
1843  enddo
1844  else if ( prc_has_n ) then
1845  tagc = 30
1846  do j = je+1, je+jhalo
1847  call mpi_irecv( var(1,is-ihalo,j), comm_size2d_4c*kd, comm_datatype, &
1848  prc_next(prc_n), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1849  ireq = ireq + 1
1850  tagc = tagc + 1
1851  enddo
1852  else if ( prc_has_w ) then
1853  tagc = 30
1854  do j = je+1, je+jhalo
1855  call mpi_irecv( var(1,is-ihalo,j), comm_size2d_4c*kd, comm_datatype, &
1856  prc_next(prc_w), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1857  ireq = ireq + 1
1858  tagc = tagc + 1
1859  enddo
1860  endif
1861  ! From S
1862  if ( prc_has_s ) then
1863  tagc = 40
1864  do j = js-jhalo, js-1
1865  call mpi_irecv( var(1,is,j), comm_size2d_ns8*kd, comm_datatype, &
1866  prc_next(prc_s), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1867  ireq = ireq + 1
1868  tagc = tagc + 1
1869  enddo
1870  endif
1871  ! From N
1872  if ( prc_has_n ) then
1873  tagc = 50
1874  do j = je+1, je+jhalo
1875  call mpi_irecv( var(1,is,j), comm_size2d_ns8*kd, comm_datatype, &
1876  prc_next(prc_n), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1877  ireq = ireq + 1
1878  tagc = tagc + 1
1879  enddo
1880  endif
1881  ! From E
1882  if ( prc_has_e ) then
1883  tagc = 60
1884  call mpi_irecv( recvpack_e2p(:,vid), comm_size2d_we*kd, comm_datatype, &
1885  prc_next(prc_e), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1886  ireq = ireq + 1
1887  endif
1888  ! From W
1889  if ( prc_has_w ) then
1890  tagc = 70
1891  call mpi_irecv( recvpack_w2p(:,vid), comm_size2d_we*kd, comm_datatype, &
1892  prc_next(prc_w), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1893  ireq = ireq + 1
1894  endif
1895 
1896  call pack_3d(var, vid)
1897 
1898  !--- To 8-Direction HALO communicate
1899  ! To W HALO
1900  if ( prc_has_w ) then
1901  tagc = 60
1902  call mpi_isend( sendpack_p2w(:,vid), comm_size2d_we*kd, comm_datatype, &
1903  prc_next(prc_w), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1904  ireq = ireq + 1
1905  endif
1906  ! To E HALO
1907  if ( prc_has_e ) then
1908  tagc = 70
1909  call mpi_isend( sendpack_p2e(:,vid), comm_size2d_we*kd, comm_datatype, &
1910  prc_next(prc_e), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1911  ireq = ireq + 1
1912  endif
1913  ! To N HALO
1914  if ( prc_has_n ) then
1915  tagc = 40
1916  do j = je-jhalo+1, je
1917  call mpi_isend( var(1,is,j), comm_size2d_ns8*kd, comm_datatype, &
1918  prc_next(prc_n), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1919  ireq = ireq + 1
1920  tagc = tagc + 1
1921  enddo
1922  endif
1923  ! To S HALO
1924  if ( prc_has_s ) then
1925  tagc = 50
1926  do j = js, js+jhalo-1
1927  call mpi_isend( var(1,is,j), comm_size2d_ns8*kd, comm_datatype, &
1928  prc_next(prc_s), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1929  ireq = ireq + 1
1930  tagc = tagc + 1
1931  enddo
1932  endif
1933  ! To NW HALO
1934  if ( prc_has_n .AND. prc_has_w ) then
1935  tagc = 0
1936  do j = je-jhalo+1, je
1937  call mpi_isend( var(1,is,j), comm_size2d_4c*kd, comm_datatype, &
1938  prc_next(prc_nw), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1939  ireq = ireq + 1
1940  tagc = tagc + 1
1941  enddo
1942  else if ( prc_has_n ) then
1943  tagc = 10
1944  do j = je-jhalo+1, je
1945  call mpi_isend( var(1,1,j), comm_size2d_4c*kd, comm_datatype, &
1946  prc_next(prc_n), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1947  ireq = ireq + 1
1948  tagc = tagc + 1
1949  enddo
1950  else if ( prc_has_w ) then
1951  tagc = 20
1952  do j = je+1, je+jhalo
1953  call mpi_isend( var(1,is,j), comm_size2d_4c*kd, comm_datatype, &
1954  prc_next(prc_w), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1955  ireq = ireq + 1
1956  tagc = tagc + 1
1957  enddo
1958  endif
1959  ! To NE HALO
1960  if ( prc_has_n .AND. prc_has_e ) then
1961  tagc = 10
1962  do j = je-jhalo+1, je
1963  call mpi_isend( var(1,ie-ihalo+1,j), comm_size2d_4c*kd, comm_datatype, &
1964  prc_next(prc_ne), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1965  ireq = ireq + 1
1966  tagc = tagc + 1
1967  enddo
1968  else if ( prc_has_n ) then
1969  tagc = 0
1970  do j = je-jhalo+1, je
1971  call mpi_isend( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
1972  prc_next(prc_n), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1973  ireq = ireq + 1
1974  tagc = tagc + 1
1975  enddo
1976  else if ( prc_has_e ) then
1977  tagc = 30
1978  do j = je+1, je+jhalo
1979  call mpi_isend( var(1,ie-ihalo+1,j), comm_size2d_4c*kd, comm_datatype, &
1980  prc_next(prc_e), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1981  ireq = ireq + 1
1982  tagc = tagc + 1
1983  enddo
1984  endif
1985  ! To SW HALO
1986  if ( prc_has_s .AND. prc_has_w ) then
1987  tagc = 20
1988  do j = js, js+jhalo-1
1989  call mpi_isend( var(1,is,j), comm_size2d_4c*kd, comm_datatype, &
1990  prc_next(prc_sw), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1991  ireq = ireq + 1
1992  tagc = tagc + 1
1993  enddo
1994  else if ( prc_has_s ) then
1995  tagc = 30
1996  do j = js, js+jhalo-1
1997  call mpi_isend( var(1,1,j), comm_size2d_4c*kd, comm_datatype, &
1998  prc_next(prc_s), tag+tagc, comm_world, req_list(ireq,vid), ierr )
1999  ireq = ireq + 1
2000  tagc = tagc + 1
2001  enddo
2002  else if ( prc_has_w ) then
2003  tagc = 0
2004  do j = js-jhalo, js-1
2005  call mpi_isend( var(1,is,j), comm_size2d_4c*kd, comm_datatype, &
2006  prc_next(prc_w), tag+tagc, comm_world, req_list(ireq,vid), ierr )
2007  ireq = ireq + 1
2008  tagc = tagc + 1
2009  enddo
2010  endif
2011  ! To SE HALO
2012  if ( prc_has_s .AND. prc_has_e ) then
2013  tagc = 30
2014  do j = js, js+jhalo-1
2015  call mpi_isend( var(1,ie-ihalo+1,j), comm_size2d_4c*kd, comm_datatype, &
2016  prc_next(prc_se), tag+tagc, comm_world, req_list(ireq,vid), ierr )
2017  ireq = ireq + 1
2018  tagc = tagc + 1
2019  enddo
2020  else if ( prc_has_s ) then
2021  tagc = 20
2022  do j = js, js+jhalo-1
2023  call mpi_isend( var(1,ie+1,j), comm_size2d_4c*kd, comm_datatype, &
2024  prc_next(prc_s), tag+tagc, comm_world, req_list(ireq,vid), ierr )
2025  ireq = ireq + 1
2026  tagc = tagc + 1
2027  enddo
2028  else if ( prc_has_e ) then
2029  tagc = 10
2030  do j = js-jhalo, js-1
2031  call mpi_isend( var(1,ie-ihalo+1,j), comm_size2d_4c*kd, comm_datatype, &
2032  prc_next(prc_e), tag+tagc, comm_world, req_list(ireq,vid), ierr )
2033  ireq = ireq + 1
2034  tagc = tagc + 1
2035  enddo
2036  endif
2037 
2038  endif
2039 
2040  req_cnt(vid) = ireq - 1
2041 
2042  return
2043  end subroutine vars8_3d_mpi
2044 
2045  subroutine vars_2d_mpi(var, vid)
2046  use scale_prc, only: &
2047  prc_abort
2048  implicit none
2049 
2050  real(RP), intent(inout) :: var(:,:)
2051  integer, intent(in) :: vid
2052 
2053  integer :: ireq, tag
2054  integer :: ierr
2055  !---------------------------------------------------------------------------
2056 
2057  tag = vid * 100
2058  ireq = 1
2059 
2060 #ifdef DEBUG
2061  if ( use_packbuf(vid) ) then
2062  log_error("vars_2D_mpi",*) 'packing buffer is already used', vid
2063  call prc_abort
2064  end if
2065  use_packbuf(vid) = .true.
2066 #endif
2067 
2068  if ( comm_isallperiodic ) then
2069  !--- periodic condition
2070  !--- From 4-Direction HALO communicate
2071  ! From S
2072  call mpi_irecv( var(:,js-jhalo:js-1), comm_size2d_ns4, &
2073  comm_datatype, prc_next(prc_s), tag+1, &
2074  comm_world, req_list(ireq,vid), ierr )
2075  ireq = ireq + 1
2076 
2077  ! From N
2078  call mpi_irecv( var(:,je+1:je+jhalo), comm_size2d_ns4, &
2079  comm_datatype, prc_next(prc_n), tag+2, &
2080  comm_world, req_list(ireq,vid), ierr )
2081  ireq = ireq + 1
2082 
2083  ! From E
2084  call mpi_irecv( recvpack_e2p(:,vid), comm_size2d_we, &
2085  comm_datatype, prc_next(prc_e), tag+3, &
2086  comm_world, req_list(ireq,vid), ierr )
2087  ireq = ireq + 1
2088 
2089  ! From W
2090  call mpi_irecv( recvpack_w2p(:,vid), comm_size2d_we, &
2091  comm_datatype, prc_next(prc_w), tag+4, &
2092  comm_world, req_list(ireq,vid), ierr )
2093  ireq = ireq + 1
2094 
2095  call pack_2d(var, vid)
2096 
2097  ! To W HALO communicate
2098  call mpi_isend( sendpack_p2w(:,vid), comm_size2d_we, &
2099  comm_datatype, prc_next(prc_w), tag+3, &
2100  comm_world, req_list(ireq,vid), ierr )
2101  ireq = ireq + 1
2102 
2103  ! To E HALO communicate
2104  call mpi_isend( sendpack_p2e(:,vid), comm_size2d_we, &
2105  comm_datatype, prc_next(prc_e), tag+4, &
2106  comm_world, req_list(ireq,vid), ierr )
2107  ireq = ireq + 1
2108 
2109  ! To N HALO communicate
2110  call mpi_isend( var(:,je-jhalo+1:je), comm_size2d_ns4, &
2111  comm_datatype, prc_next(prc_n), tag+1, &
2112  comm_world, req_list(ireq,vid), ierr )
2113  ireq = ireq + 1
2114 
2115  ! To S HALO communicate
2116  call mpi_isend( var(:,js:js+jhalo-1), comm_size2d_ns4, &
2117  comm_datatype, prc_next(prc_s), tag+2, &
2118  comm_world, req_list(ireq,vid), ierr )
2119  ireq = ireq + 1
2120 
2121  else
2122  !--- non-periodic condition
2123  !--- From 4-Direction HALO communicate
2124  ! From S
2125  if ( prc_has_s ) then
2126  call mpi_irecv( var(:,js-jhalo:js-1), comm_size2d_ns4, &
2127  comm_datatype, prc_next(prc_s), tag+1, &
2128  comm_world, req_list(ireq,vid), ierr )
2129  ireq = ireq + 1
2130  endif
2131 
2132  ! From N
2133  if ( prc_has_n ) then
2134  call mpi_irecv( var(:,je+1:je+jhalo), comm_size2d_ns4, &
2135  comm_datatype, prc_next(prc_n), tag+2, &
2136  comm_world, req_list(ireq,vid), ierr )
2137  ireq = ireq + 1
2138  endif
2139 
2140  ! From E
2141  if ( prc_has_e ) then
2142  call mpi_irecv( recvpack_e2p(:,vid), comm_size2d_we, &
2143  comm_datatype, prc_next(prc_e), tag+3, &
2144  comm_world, req_list(ireq,vid), ierr )
2145  ireq = ireq + 1
2146  endif
2147 
2148  ! From W
2149  if ( prc_has_w ) then
2150  call mpi_irecv( recvpack_w2p(:,vid), comm_size2d_we, &
2151  comm_datatype, prc_next(prc_w), tag+4, &
2152  comm_world, req_list(ireq,vid), ierr )
2153  ireq = ireq + 1
2154  endif
2155 
2156  call pack_2d(var, vid)
2157 
2158  ! To W HALO communicate
2159  if ( prc_has_w ) then
2160  call mpi_isend( sendpack_p2w(:,vid), comm_size2d_we, &
2161  comm_datatype, prc_next(prc_w), tag+3, &
2162  comm_world, req_list(ireq,vid), ierr )
2163  ireq = ireq + 1
2164  endif
2165 
2166  ! To E HALO communicate
2167  if ( prc_has_e ) then
2168  call mpi_isend( sendpack_p2e(:,vid), comm_size2d_we, &
2169  comm_datatype, prc_next(prc_e), tag+4, &
2170  comm_world, req_list(ireq,vid), ierr )
2171  ireq = ireq + 1
2172  endif
2173 
2174  ! To N HALO communicate
2175  if ( prc_has_n ) then
2176  call mpi_isend( var(:,je-jhalo+1:je), comm_size2d_ns4, &
2177  comm_datatype, prc_next(prc_n), tag+1, &
2178  comm_world, req_list(ireq,vid), ierr )
2179  ireq = ireq + 1
2180  endif
2181 
2182  ! To S HALO communicate
2183  if ( prc_has_s ) then
2184  call mpi_isend( var(:,js:js+jhalo-1), comm_size2d_ns4, &
2185  comm_datatype, prc_next(prc_s), tag+2, &
2186  comm_world, req_list(ireq,vid), ierr )
2187  ireq = ireq + 1
2188  endif
2189 
2190  endif
2191 
2192  req_cnt(vid) = ireq - 1
2193 
2194  return
2195  end subroutine vars_2d_mpi
2196 
2197  subroutine vars8_2d_mpi(var, vid)
2198  use scale_prc, only: &
2199  prc_abort
2200  implicit none
2201 
2202  real(RP), intent(inout) :: var(:,:)
2203  integer, intent(in) :: vid
2204 
2205  integer :: ireq, tag, tagc
2206 
2207  integer :: ierr
2208  integer :: j
2209  !---------------------------------------------------------------------------
2210 
2211  tag = vid * 100
2212  ireq = 1
2213 
2214 #ifdef DEBUG
2215  if ( use_packbuf(vid) ) then
2216  log_error("vars8_2D_mpi",*) 'packing buffer is already used', vid
2217  call prc_abort
2218  end if
2219  use_packbuf(vid) = .true.
2220 #endif
2221 
2222  if ( comm_isallperiodic ) then
2223  !--- periodic condition
2224  !--- From 8-Direction HALO communicate
2225  ! From SE
2226  tagc = 0
2227  do j = js-jhalo, js-1
2228  call mpi_irecv( var(ie+1,j), comm_size2d_4c, &
2229  comm_datatype, prc_next(prc_se), tag+tagc, &
2230  comm_world, req_list(ireq,vid), ierr )
2231  ireq = ireq + 1
2232  tagc = tagc + 1
2233  enddo
2234  ! From SW
2235  tagc = 10
2236  do j = js-jhalo, js-1
2237  call mpi_irecv( var(is-ihalo,j), comm_size2d_4c, &
2238  comm_datatype, prc_next(prc_sw), tag+tagc, &
2239  comm_world, req_list(ireq,vid), ierr )
2240  ireq = ireq + 1
2241  tagc = tagc + 1
2242  enddo
2243  ! From NE
2244  tagc = 20
2245  do j = je+1, je+jhalo
2246  call mpi_irecv( var(ie+1,j), comm_size2d_4c, &
2247  comm_datatype, prc_next(prc_ne), tag+tagc, &
2248  comm_world, req_list(ireq,vid), ierr )
2249  ireq = ireq + 1
2250  tagc = tagc + 1
2251  enddo
2252  ! From NW
2253  tagc = 30
2254  do j = je+1, je+jhalo
2255  call mpi_irecv( var(is-ihalo,j), comm_size2d_4c, &
2256  comm_datatype, prc_next(prc_nw), tag+tagc, &
2257  comm_world, req_list(ireq,vid), ierr )
2258  ireq = ireq + 1
2259  tagc = tagc + 1
2260  enddo
2261  ! From S
2262  tagc = 40
2263  do j = js-jhalo, js-1
2264  call mpi_irecv( var(is,j), comm_size2d_ns8, &
2265  comm_datatype, prc_next(prc_s), tag+tagc, &
2266  comm_world, req_list(ireq,vid), ierr )
2267  ireq = ireq + 1
2268  tagc = tagc + 1
2269  enddo
2270  ! From N
2271  tagc = 50
2272  do j = je+1, je+jhalo
2273  call mpi_irecv( var(is,j), comm_size2d_ns8, &
2274  comm_datatype, prc_next(prc_n), tag+tagc, &
2275  comm_world, req_list(ireq,vid), ierr )
2276  ireq = ireq + 1
2277  tagc = tagc + 1
2278  enddo
2279  ! From E
2280  call mpi_irecv( recvpack_e2p(:,vid), comm_size2d_we, &
2281  comm_datatype, prc_next(prc_e), tag+60, &
2282  comm_world, req_list(ireq,vid), ierr )
2283  ireq = ireq + 1
2284  ! From W
2285  call mpi_irecv( recvpack_w2p(:,vid), comm_size2d_we, &
2286  comm_datatype, prc_next(prc_w), tag+70, &
2287  comm_world, req_list(ireq,vid), ierr )
2288  ireq = ireq + 1
2289 
2290  call pack_2d(var, vid)
2291 
2292  ! To W HALO communicate
2293  call mpi_isend( sendpack_p2w(:,vid), comm_size2d_we, &
2294  comm_datatype, prc_next(prc_w), tag+60, &
2295  comm_world, req_list(ireq,vid), ierr )
2296  ireq = ireq + 1
2297 
2298  ! To E HALO communicate
2299  call mpi_isend( sendpack_p2e(:,vid), comm_size2d_we, &
2300  comm_datatype, prc_next(prc_e), tag+70, &
2301  comm_world, req_list(ireq,vid), ierr )
2302  ireq = ireq + 1
2303 
2304  ! To N HALO communicate
2305  tagc = 40
2306  do j = je-jhalo+1, je
2307  call mpi_isend( var(is,j), comm_size2d_ns8, &
2308  comm_datatype, prc_next(prc_n), tag+tagc, &
2309  comm_world, req_list(ireq,vid), ierr )
2310  ireq = ireq + 1
2311  tagc = tagc + 1
2312  enddo
2313 
2314  ! To S HALO communicate
2315  tagc = 50
2316  do j = js, js+jhalo-1
2317  call mpi_isend( var(is,j), comm_size2d_ns8, &
2318  comm_datatype, prc_next(prc_s), tag+tagc, &
2319  comm_world, req_list(ireq,vid), ierr )
2320  ireq = ireq + 1
2321  tagc = tagc + 1
2322  enddo
2323 
2324  ! To NW HALO communicate
2325  tagc = 0
2326  do j = je-jhalo+1, je
2327  call mpi_isend( var(is,j), comm_size2d_4c, &
2328  comm_datatype, prc_next(prc_nw), tag+tagc, &
2329  comm_world, req_list(ireq,vid), ierr )
2330  ireq = ireq + 1
2331  tagc = tagc + 1
2332  enddo
2333 
2334  ! To NE HALO communicate
2335  tagc = 10
2336  do j = je-jhalo+1, je
2337  call mpi_isend( var(ie-ihalo+1,j), comm_size2d_4c, &
2338  comm_datatype, prc_next(prc_ne), tag+tagc, &
2339  comm_world, req_list(ireq,vid), ierr )
2340  ireq = ireq + 1
2341  tagc = tagc + 1
2342  enddo
2343 
2344  ! To SW HALO communicate
2345  tagc = 20
2346  do j = js, js+jhalo-1
2347  call mpi_isend( var(is,j), comm_size2d_4c, &
2348  comm_datatype, prc_next(prc_sw), tag+tagc, &
2349  comm_world, req_list(ireq,vid), ierr )
2350  ireq = ireq + 1
2351  tagc = tagc + 1
2352  enddo
2353 
2354  ! To SE HALO communicate
2355  tagc = 30
2356  do j = js, js+jhalo-1
2357  call mpi_isend( var(ie-ihalo+1,j), comm_size2d_4c, &
2358  comm_datatype, prc_next(prc_se), tag+tagc, &
2359  comm_world, req_list(ireq,vid), ierr )
2360  ireq = ireq + 1
2361  tagc = tagc + 1
2362  enddo
2363  else
2364  !--- non-periodic condition
2365  !--- From 8-Direction HALO communicate
2366  ! From SE
2367  if ( prc_has_s .AND. prc_has_e ) then
2368  tagc = 0
2369  do j = js-jhalo, js-1
2370  call mpi_irecv( var(ie+1,j), comm_size2d_4c, &
2371  comm_datatype, prc_next(prc_se), tag+tagc, &
2372  comm_world, req_list(ireq,vid), ierr )
2373  ireq = ireq + 1
2374  tagc = tagc + 1
2375  enddo
2376  else if ( prc_has_s ) then
2377  tagc = 0
2378  do j = js-jhalo, js-1
2379  call mpi_irecv( var(ie+1,j), comm_size2d_4c, &
2380  comm_datatype, prc_next(prc_s), tag+tagc, &
2381  comm_world, req_list(ireq,vid), ierr )
2382  ireq = ireq + 1
2383  tagc = tagc + 1
2384  enddo
2385  else if ( prc_has_e ) then
2386  tagc = 0
2387  do j = js-jhalo, js-1
2388  call mpi_irecv( var(ie+1,j), comm_size2d_4c, &
2389  comm_datatype, prc_next(prc_e), tag+tagc, &
2390  comm_world, req_list(ireq,vid), ierr )
2391  ireq = ireq + 1
2392  tagc = tagc + 1
2393  enddo
2394  endif
2395 
2396  ! From SW
2397  if ( prc_has_s .AND. prc_has_w ) then
2398  tagc = 10
2399  do j = js-jhalo, js-1
2400  call mpi_irecv( var(is-ihalo,j), comm_size2d_4c, &
2401  comm_datatype, prc_next(prc_sw), tag+tagc, &
2402  comm_world, req_list(ireq,vid), ierr )
2403  ireq = ireq + 1
2404  tagc = tagc + 1
2405  enddo
2406  else if ( prc_has_s ) then
2407  tagc = 10
2408  do j = js-jhalo, js-1
2409  call mpi_irecv( var(is-ihalo,j), comm_size2d_4c, &
2410  comm_datatype, prc_next(prc_s), tag+tagc, &
2411  comm_world, req_list(ireq,vid), ierr )
2412  ireq = ireq + 1
2413  tagc = tagc + 1
2414  enddo
2415  else if ( prc_has_w ) then
2416  tagc = 10
2417  do j = js-jhalo, js-1
2418  call mpi_irecv( var(is-ihalo,j), comm_size2d_4c, &
2419  comm_datatype, prc_next(prc_w), tag+tagc, &
2420  comm_world, req_list(ireq,vid), ierr )
2421  ireq = ireq + 1
2422  tagc = tagc + 1
2423  enddo
2424  endif
2425 
2426  ! From NE
2427  if ( prc_has_n .AND. prc_has_e ) then
2428  tagc = 20
2429  do j = je+1, je+jhalo
2430  call mpi_irecv( var(ie+1,j), comm_size2d_4c, &
2431  comm_datatype, prc_next(prc_ne), tag+tagc, &
2432  comm_world, req_list(ireq,vid), ierr )
2433  ireq = ireq + 1
2434  tagc = tagc + 1
2435  enddo
2436  else if ( prc_has_n ) then
2437  tagc = 20
2438  do j = je+1, je+jhalo
2439  call mpi_irecv( var(ie+1,j), comm_size2d_4c, &
2440  comm_datatype, prc_next(prc_n), tag+tagc, &
2441  comm_world, req_list(ireq,vid), ierr )
2442  ireq = ireq + 1
2443  tagc = tagc + 1
2444  enddo
2445  else if ( prc_has_e ) then
2446  tagc = 20
2447  do j = je+1, je+jhalo
2448  call mpi_irecv( var(ie+1,j), comm_size2d_4c, &
2449  comm_datatype, prc_next(prc_e), tag+tagc, &
2450  comm_world, req_list(ireq,vid), ierr )
2451  ireq = ireq + 1
2452  tagc = tagc + 1
2453  enddo
2454  endif
2455 
2456  ! From NW
2457  if ( prc_has_n .AND. prc_has_w ) then
2458  tagc = 30
2459  do j = je+1, je+jhalo
2460  call mpi_irecv( var(is-ihalo,j), comm_size2d_4c, &
2461  comm_datatype, prc_next(prc_nw), tag+tagc, &
2462  comm_world, req_list(ireq,vid), ierr )
2463  ireq = ireq + 1
2464  tagc = tagc + 1
2465  enddo
2466  else if ( prc_has_n ) then
2467  tagc = 30
2468  do j = je+1, je+jhalo
2469  call mpi_irecv( var(is-ihalo,j), comm_size2d_4c, &
2470  comm_datatype, prc_next(prc_n), tag+tagc, &
2471  comm_world, req_list(ireq,vid), ierr )
2472  ireq = ireq + 1
2473  tagc = tagc + 1
2474  enddo
2475  else if ( prc_has_w ) then
2476  tagc = 30
2477  do j = je+1, je+jhalo
2478  call mpi_irecv( var(is-ihalo,j), comm_size2d_4c, &
2479  comm_datatype, prc_next(prc_w), tag+tagc, &
2480  comm_world, req_list(ireq,vid), ierr )
2481  ireq = ireq + 1
2482  tagc = tagc + 1
2483  enddo
2484  endif
2485 
2486  ! From S
2487  if ( prc_has_s ) then
2488  tagc = 40
2489  do j = js-jhalo, js-1
2490  call mpi_irecv( var(is,j), comm_size2d_ns8, &
2491  comm_datatype, prc_next(prc_s), tag+tagc, &
2492  comm_world, req_list(ireq,vid), ierr )
2493  ireq = ireq + 1
2494  tagc = tagc + 1
2495  enddo
2496  endif
2497 
2498  ! From N
2499  if ( prc_has_n ) then
2500  tagc = 50
2501  do j = je+1, je+jhalo
2502  call mpi_irecv( var(is,j), comm_size2d_ns8, &
2503  comm_datatype, prc_next(prc_n), tag+tagc, &
2504  comm_world, req_list(ireq,vid), ierr )
2505  ireq = ireq + 1
2506  tagc = tagc + 1
2507  enddo
2508  endif
2509 
2510  ! From E
2511  if ( prc_has_e ) then
2512  call mpi_irecv( recvpack_e2p(:,vid), comm_size2d_we, &
2513  comm_datatype, prc_next(prc_e), tag+60, &
2514  comm_world, req_list(ireq,vid), ierr )
2515  ireq = ireq + 1
2516  endif
2517 
2518  ! From W
2519  if ( prc_has_w ) then
2520  call mpi_irecv( recvpack_w2p(:,vid), comm_size2d_we, &
2521  comm_datatype, prc_next(prc_w), tag+70, &
2522  comm_world, req_list(ireq,vid), ierr )
2523  ireq = ireq + 1
2524  endif
2525 
2526  call pack_2d(var, vid)
2527 
2528  ! To W HALO communicate
2529  if ( prc_has_w ) then
2530  call mpi_isend( sendpack_p2w(:,vid), comm_size2d_we, &
2531  comm_datatype, prc_next(prc_w), tag+60, &
2532  comm_world, req_list(ireq,vid), ierr )
2533  ireq = ireq + 1
2534  endif
2535 
2536  ! To E HALO communicate
2537  if ( prc_has_e ) then
2538  call mpi_isend( sendpack_p2e(:,vid), comm_size2d_we, &
2539  comm_datatype, prc_next(prc_e), tag+70, &
2540  comm_world, req_list(ireq,vid), ierr )
2541  ireq = ireq + 1
2542  endif
2543 
2544  ! To N HALO communicate
2545  if ( prc_has_n ) then
2546  tagc = 40
2547  do j = je-jhalo+1, je
2548  call mpi_isend( var(is,j), comm_size2d_ns8, &
2549  comm_datatype, prc_next(prc_n), tag+tagc, &
2550  comm_world, req_list(ireq,vid), ierr )
2551  ireq = ireq + 1
2552  tagc = tagc + 1
2553  enddo
2554  endif
2555 
2556  ! To S HALO communicate
2557  if ( prc_has_s ) then
2558  tagc = 50
2559  do j = js, js+jhalo-1
2560  call mpi_isend( var(is,j), comm_size2d_ns8, &
2561  comm_datatype, prc_next(prc_s), tag+tagc, &
2562  comm_world, req_list(ireq,vid), ierr )
2563  ireq = ireq + 1
2564  tagc = tagc + 1
2565  enddo
2566  endif
2567 
2568  ! To NW HALO communicate
2569  if ( prc_has_n .AND. prc_has_w ) then
2570  tagc = 0
2571  do j = je-jhalo+1, je
2572  call mpi_isend( var(is,j), comm_size2d_4c, &
2573  comm_datatype, prc_next(prc_nw), tag+tagc, &
2574  comm_world, req_list(ireq,vid), ierr )
2575  ireq = ireq + 1
2576  tagc = tagc + 1
2577  enddo
2578  else if ( prc_has_n ) then
2579  tagc = 10
2580  do j = je-jhalo+1, je
2581  call mpi_isend( var(is,j), comm_size2d_4c, &
2582  comm_datatype, prc_next(prc_n), tag+tagc, &
2583  comm_world, req_list(ireq,vid), ierr )
2584  ireq = ireq + 1
2585  tagc = tagc + 1
2586  enddo
2587  else if ( prc_has_w ) then
2588  tagc = 20
2589  do j = je-jhalo+1, je
2590  call mpi_isend( var(is,j), comm_size2d_4c, &
2591  comm_datatype, prc_next(prc_w), tag+tagc, &
2592  comm_world, req_list(ireq,vid), ierr )
2593  ireq = ireq + 1
2594  tagc = tagc + 1
2595  enddo
2596  endif
2597 
2598  ! To NE HALO communicate
2599  if ( prc_has_n .AND. prc_has_e ) then
2600  tagc = 10
2601  do j = je-jhalo+1, je
2602  call mpi_isend( var(ie-ihalo+1,j), comm_size2d_4c, &
2603  comm_datatype, prc_next(prc_ne), tag+tagc, &
2604  comm_world, req_list(ireq,vid), ierr )
2605  ireq = ireq + 1
2606  tagc = tagc + 1
2607  enddo
2608  else if ( prc_has_n ) then
2609  tagc = 0
2610  do j = je-jhalo+1, je
2611  call mpi_isend( var(ie-ihalo+1,j), comm_size2d_4c, &
2612  comm_datatype, prc_next(prc_n), tag+tagc, &
2613  comm_world, req_list(ireq,vid), ierr )
2614  ireq = ireq + 1
2615  tagc = tagc + 1
2616  enddo
2617  else if ( prc_has_e ) then
2618  tagc = 30
2619  do j = je-jhalo+1, je
2620  call mpi_isend( var(ie-ihalo+1,j), comm_size2d_4c, &
2621  comm_datatype, prc_next(prc_e), tag+tagc, &
2622  comm_world, req_list(ireq,vid), ierr )
2623  ireq = ireq + 1
2624  tagc = tagc + 1
2625  enddo
2626  endif
2627 
2628  ! To SW HALO communicate
2629  if ( prc_has_s .AND. prc_has_w ) then
2630  tagc = 20
2631  do j = js, js+jhalo-1
2632  call mpi_isend( var(is,j), comm_size2d_4c, &
2633  comm_datatype, prc_next(prc_sw), tag+tagc, &
2634  comm_world, req_list(ireq,vid), ierr )
2635  ireq = ireq + 1
2636  tagc = tagc + 1
2637  enddo
2638  else if ( prc_has_s ) then
2639  tagc = 30
2640  do j = js, js+jhalo-1
2641  call mpi_isend( var(is,j), comm_size2d_4c, &
2642  comm_datatype, prc_next(prc_s), tag+tagc, &
2643  comm_world, req_list(ireq,vid), ierr )
2644  ireq = ireq + 1
2645  tagc = tagc + 1
2646  enddo
2647  else if ( prc_has_w ) then
2648  tagc = 0
2649  do j = js, js+jhalo-1
2650  call mpi_isend( var(is,j), comm_size2d_4c, &
2651  comm_datatype, prc_next(prc_w), tag+tagc, &
2652  comm_world, req_list(ireq,vid), ierr )
2653  ireq = ireq + 1
2654  tagc = tagc + 1
2655  enddo
2656  endif
2657 
2658  ! To SE HALO communicate
2659  if ( prc_has_s .AND. prc_has_e ) then
2660  tagc = 30
2661  do j = js, js+jhalo-1
2662  call mpi_isend( var(ie-ihalo+1,j), comm_size2d_4c, &
2663  comm_datatype, prc_next(prc_se), tag+tagc, &
2664  comm_world, req_list(ireq,vid), ierr )
2665  ireq = ireq + 1
2666  tagc = tagc + 1
2667  enddo
2668  else if ( prc_has_s ) then
2669  tagc = 20
2670  do j = js, js+jhalo-1
2671  call mpi_isend( var(ie-ihalo+1,j), comm_size2d_4c, &
2672  comm_datatype, prc_next(prc_s), tag+tagc, &
2673  comm_world, req_list(ireq,vid), ierr )
2674  ireq = ireq + 1
2675  tagc = tagc + 1
2676  enddo
2677  else if ( prc_has_e ) then
2678  tagc = 10
2679  do j = js, js+jhalo-1
2680  call mpi_isend( var(ie-ihalo+1,j), comm_size2d_4c, &
2681  comm_datatype, prc_next(prc_e), tag+tagc, &
2682  comm_world, req_list(ireq,vid), ierr )
2683  ireq = ireq + 1
2684  tagc = tagc + 1
2685  enddo
2686  endif
2687 
2688  endif
2689 
2690  req_cnt(vid) = ireq - 1
2691 
2692  return
2693  end subroutine vars8_2d_mpi
2694 
2695  subroutine vars_3d_mpi_pc(var, vid)
2696  use scale_prc, only: &
2697  prc_abort
2698  implicit none
2699 
2700  real(RP), intent(inout) :: var(:,:,:)
2701  integer, intent(in) :: vid
2702  integer :: ierr
2703  !---------------------------------------------------------------------------
2704 
2705 #ifdef DEBUG
2706  if ( use_packbuf(pseqid(vid)) ) then
2707  log_error("vars_3D_mpi_pc",*) 'packing buffer is already used', vid, pseqid(vid)
2708  call prc_abort
2709  end if
2710  use_packbuf(pseqid(vid)) = .true.
2711 #endif
2712 
2713  call pack_3d(var, pseqid(vid))
2714 
2715  call mpi_startall(preq_cnt(vid), preq_list(1:preq_cnt(vid),vid), ierr)
2716 
2717  return
2718  end subroutine vars_3d_mpi_pc
2719 
2720  subroutine wait_3d_mpi(var, vid)
2721  implicit none
2722  real(RP), intent(inout) :: var(:,:,:)
2723  integer, intent(in) :: vid
2724 
2725  integer :: ierr
2726  !---------------------------------------------------------------------------
2727 
2728  !--- wait packets
2729  call mpi_waitall( req_cnt(vid), &
2730  req_list(1:req_cnt(vid),vid), &
2731  mpi_statuses_ignore, &
2732  ierr )
2733  call unpack_3d(var, vid)
2734 
2735 #ifdef DEBUG
2736  use_packbuf(vid) = .false.
2737 #endif
2738 
2739  return
2740  end subroutine wait_3d_mpi
2741 
2742  subroutine wait_2d_mpi(var, vid)
2743  implicit none
2744 
2745  real(RP), intent(inout) :: var(:,:)
2746  integer, intent(in) :: vid
2747 
2748  integer :: ierr
2749  !---------------------------------------------------------------------------
2750 
2751  !--- wait packets
2752  call mpi_waitall( req_cnt(vid), &
2753  req_list(1:req_cnt(vid),vid), &
2754  mpi_statuses_ignore, &
2755  ierr )
2756  call unpack_2d(var, vid)
2757 
2758 #ifdef DEBUG
2759  use_packbuf(vid) = .false.
2760 #endif
2761 
2762  return
2763  end subroutine wait_2d_mpi
2764 
2765  subroutine wait_3d_mpi_pc(var, vid)
2766  implicit none
2767 
2768  real(RP), intent(inout) :: var(:,:,:)
2769  integer, intent(in) :: vid
2770 
2771  integer :: ierr
2772 
2773  !--- wait packets
2774  call mpi_waitall( preq_cnt(vid), &
2775  preq_list(1:preq_cnt(vid),vid), &
2776  mpi_statuses_ignore, &
2777  ierr )
2778  call unpack_3d(var, pseqid(vid))
2779 
2780 #ifdef DEBUG
2781  use_packbuf(pseqid(vid)) = .false.
2782 #endif
2783 
2784  return
2785  end subroutine wait_3d_mpi_pc
2786 
2787  subroutine pack_3d(var, vid)
2788  implicit none
2789 
2790  real(RP), intent(in) :: var(:,:,:)
2791  integer, intent(in) :: vid
2792 
2793  integer :: kd
2794  integer :: k, i, j, n
2795 
2796  kd = size(var, 1)
2797 
2798  call prof_rapstart('COMM_pack', 3)
2799 
2800  if ( comm_isallperiodic ) then ! periodic condition
2801 
2802  !--- packing packets to West
2803 !OCL NORECURRENCE(sendpack_P2W)
2804  !$omp parallel do private(i,j,k,n) OMP_SCHEDULE_ collapse(2)
2805  do j = js, je
2806  do i = is, is+ihalo-1
2807  do k = 1, kd
2808  n = (j-js) * kd * ihalo &
2809  + (i-is) * kd &
2810  + k
2811  sendpack_p2w(n,vid) = var(k,i,j)
2812  enddo
2813  enddo
2814  enddo
2815  !--- packing packets to East
2816 !OCL NORECURRENCE(sendpack_P2E)
2817  !$omp parallel do private(i,j,k,n) OMP_SCHEDULE_ collapse(2)
2818  do j = js, je
2819  do i = ie-ihalo+1, ie
2820  do k = 1, kd
2821  n = (j-js) * kd * ihalo &
2822  + (i-ie+ihalo-1) * kd &
2823  + k
2824  sendpack_p2e(n,vid) = var(k,i,j)
2825  enddo
2826  enddo
2827  enddo
2828 
2829  else
2830 
2831  if ( prc_has_w ) then
2832  !--- packing packets to West
2833 !OCL NORECURRENCE(sendpack_P2W)
2834  !$omp parallel do default(none) private(i,j,k,n) OMP_SCHEDULE_ collapse(2) &
2835  !$omp shared(JS,JE,IS,IHALO,kd,var,sendpack_P2W,vid)
2836  do j = js, je
2837  do i = is, is+ihalo-1
2838  do k = 1, kd
2839  n = (j-js) * kd * ihalo &
2840  + (i-is) * kd &
2841  + k
2842  sendpack_p2w(n,vid) = var(k,i,j)
2843  enddo
2844  enddo
2845  enddo
2846  endif
2847  if ( prc_has_e ) then
2848  !--- packing packets to East
2849 !OCL NORECURRENCE(sendpack_P2E)
2850  !$omp parallel do default(none) private(i,j,k,n) OMP_SCHEDULE_ collapse(2) &
2851  !$omp shared(JS,JE,IE,IHALO,kd,var,sendpack_P2E,vid)
2852  do j = js, je
2853  do i = ie-ihalo+1, ie
2854  do k = 1, kd
2855  n = (j-js) * kd * ihalo &
2856  + (i-ie+ihalo-1) * kd &
2857  + k
2858  sendpack_p2e(n,vid) = var(k,i,j)
2859  enddo
2860  enddo
2861  enddo
2862  endif
2863 
2864  end if
2865 
2866  call prof_rapend('COMM_pack', 3)
2867 
2868  return
2869  end subroutine pack_3d
2870 
2871  subroutine pack_2d(var, vid)
2872  implicit none
2873 
2874  real(RP), intent(in) :: var(:,:)
2875  integer, intent(in) :: vid
2876 
2877  integer :: i, j, n
2878 
2879  call prof_rapstart('COMM_pack', 3)
2880 
2881  if ( comm_isallperiodic ) then ! periodic condition
2882 
2883  !--- To 4-Direction HALO communicate
2884  !--- packing packets to West
2885 !OCL NORECURRENCE(sendpack_P2W)
2886  !$omp parallel do private(i,j,n) OMP_SCHEDULE_ collapse(2)
2887  do j = js, je
2888  do i = is, is+ihalo-1
2889  n = (j-js) * ihalo &
2890  + (i-is) + 1
2891  sendpack_p2w(n,vid) = var(i,j)
2892  enddo
2893  enddo
2894 
2895  !--- packing packets to East
2896  !$omp parallel do private(i,j,n) OMP_SCHEDULE_ collapse(2)
2897 !OCL NORECURRENCE(sendpack_P2E)
2898  do j = js, je
2899  do i = ie-ihalo+1, ie
2900  n = (j-js) * ihalo &
2901  + (i-ie+ihalo-1) + 1
2902  sendpack_p2e(n,vid) = var(i,j)
2903  enddo
2904  enddo
2905 
2906  else
2907 
2908  !--- To 4-Direction HALO communicate
2909  !--- packing packets to West
2910  if ( prc_has_w ) then
2911 !OCL NORECURRENCE(sendpack_P2W)
2912  !$omp parallel do private(i,j,n) OMP_SCHEDULE_ collapse(2)
2913  do j = js, je
2914  do i = is, is+ihalo-1
2915  n = (j-js) * ihalo &
2916  + (i-is) + 1
2917  sendpack_p2w(n,vid) = var(i,j)
2918  enddo
2919  enddo
2920  endif
2921 
2922  !--- packing packets to East
2923  if ( prc_has_e ) then
2924 !OCL NORECURRENCE(sendpack_P2E)
2925  !$omp parallel do private(i,j,n) OMP_SCHEDULE_ collapse(2)
2926  do j = js, je
2927  do i = ie-ihalo+1, ie
2928  n = (j-js) * ihalo &
2929  + (i-ie+ihalo-1) + 1
2930  sendpack_p2e(n,vid) = var(i,j)
2931  enddo
2932  enddo
2933  endif
2934 
2935  end if
2936 
2937  call prof_rapend('COMM_pack', 3)
2938 
2939  return
2940  end subroutine pack_2d
2941 
2942  subroutine unpack_3d(var, vid)
2943  implicit none
2944 
2945  real(RP), intent(inout) :: var(:,:,:)
2946  integer, intent(in) :: vid
2947 
2948  integer :: kd
2949  integer :: i, j, k, n
2950  !---------------------------------------------------------------------------
2951 
2952  kd = size(var, 1)
2953 
2954  call prof_rapstart('COMM_unpack', 3)
2955 
2956  if ( comm_isallperiodic ) then ! periodic condition
2957 
2958  !--- unpacking packets from East
2959  !$omp parallel do private(i,j,k,n) OMP_SCHEDULE_ collapse(2)
2960  do j = js, je
2961  do i = ie+1, ie+ihalo
2962  do k = 1, kd
2963  n = (j-js) * kd * ihalo &
2964  + (i-ie-1) * kd &
2965  + k
2966  var(k,i,j) = recvpack_e2p(n,vid)
2967  enddo
2968  enddo
2969  enddo
2970 
2971  !--- unpacking packets from West
2972  !$omp parallel do private(i,j,k,n) OMP_SCHEDULE_ collapse(2)
2973  do j = js, je
2974  do i = is-ihalo, is-1
2975  do k = 1, kd
2976  n = (j-js) * kd * ihalo &
2977  + (i-is+ihalo) * kd &
2978  + k
2979  var(k,i,j) = recvpack_w2p(n,vid)
2980  enddo
2981  enddo
2982  enddo
2983 
2984  else ! non-periodic condition
2985 
2986  if ( prc_has_e ) then
2987  !--- unpacking packets from East
2988  !$omp parallel do default(none) private(i,j,k,n) OMP_SCHEDULE_ collapse(2) &
2989  !$omp shared(JS,JE,IE,IHALO,kd,var,recvpack_E2P,vid)
2990  do j = js, je
2991  do i = ie+1, ie+ihalo
2992  do k = 1, kd
2993  n = (j-js) * kd * ihalo &
2994  + (i-ie-1) * kd &
2995  + k
2996  var(k,i,j) = recvpack_e2p(n,vid)
2997  enddo
2998  enddo
2999  enddo
3000  endif
3001 
3002  if ( prc_has_w ) then
3003  !--- unpacking packets from West
3004  !$omp parallel do default(none) private(i,j,k,n) OMP_SCHEDULE_ collapse(2) &
3005  !$omp shared(JS,JE,IS,IHALO,kd,var,recvpack_W2P,vid)
3006  do j = js, je
3007  do i = is-ihalo, is-1
3008  do k = 1, kd
3009  n = (j-js) * kd * ihalo &
3010  + (i-is+ihalo) * kd &
3011  + k
3012  var(k,i,j) = recvpack_w2p(n,vid)
3013  enddo
3014  enddo
3015  enddo
3016  endif
3017 
3018  end if
3019 
3020  call prof_rapend('COMM_unpack', 3)
3021 
3022  return
3023  end subroutine unpack_3d
3024 
3025  subroutine unpack_2d(var, vid)
3026  implicit none
3027 
3028  real(RP), intent(inout) :: var(:,:)
3029  integer, intent(in) :: vid
3030 
3031  integer :: i, j, n
3032  !---------------------------------------------------------------------------
3033 
3034  call prof_rapstart('COMM_unpack', 3)
3035 
3036  if( comm_isallperiodic ) then
3037  !--- periodic condition
3038  !--- unpacking packets from East
3039  !$omp parallel do private(i,j,n) OMP_SCHEDULE_ collapse(2)
3040  do j = js, je
3041  do i = ie+1, ie+ihalo
3042  n = (j-js) * ihalo &
3043  + (i-ie-1) + 1
3044  var(i,j) = recvpack_e2p(n,vid)
3045  enddo
3046  enddo
3047 
3048  !--- unpacking packets from West
3049  !$omp parallel do private(i,j,n) OMP_SCHEDULE_ collapse(2)
3050  do j = js, je
3051  do i = is-ihalo, is-1
3052  n = (j-js) * ihalo &
3053  + (i-is+ihalo) + 1
3054  var(i,j) = recvpack_w2p(n,vid)
3055  enddo
3056  enddo
3057 
3058  else
3059  !--- non-periodic condition
3060 
3061  !--- unpacking packets from East / copy inner data to HALO(East)
3062  if( prc_has_e ) then
3063  !$omp parallel do private(i,j,n) OMP_SCHEDULE_ collapse(2)
3064  do j = js, je
3065  do i = ie+1, ie+ihalo
3066  n = (j-js) * ihalo &
3067  + (i-ie-1) + 1
3068  var(i,j) = recvpack_e2p(n,vid)
3069  enddo
3070  enddo
3071  end if
3072 
3073 
3074  !--- unpacking packets from West / copy inner data to HALO(West)
3075  if( prc_has_w ) then
3076  !$omp parallel do private(i,j,n) OMP_SCHEDULE_ collapse(2)
3077  do j = js, je
3078  do i = is-ihalo, is-1
3079  n = (j-js) * ihalo &
3080  + (i-is+ihalo) + 1
3081  var(i,j) = recvpack_w2p(n,vid)
3082  enddo
3083  enddo
3084  end if
3085 
3086  end if
3087 
3088  call prof_rapend('COMM_unpack', 3)
3089 
3090  return
3091  end subroutine unpack_2d
3092 
3093  subroutine copy_boundary_3d(var)
3094  implicit none
3095 
3096  real(RP), intent(inout) :: var(:,:,:)
3097 
3098  integer :: i, j
3099  !---------------------------------------------------------------------------
3100 
3101  !--- copy inner data to HALO(North)
3102  if ( .NOT. prc_has_n ) then
3103  do j = je+1, je+jhalo
3104  do i = is, ie
3105  var(:,i,j) = var(:,i,je)
3106  enddo
3107  enddo
3108  endif
3109 
3110  !--- copy inner data to HALO(South)
3111  if ( .NOT. prc_has_s ) then
3112  do j = js-jhalo, js-1
3113  do i = is, ie
3114  var(:,i,j) = var(:,i,js)
3115  enddo
3116  enddo
3117  endif
3118 
3119  !--- copy inner data to HALO(East)
3120  if ( .NOT. prc_has_e ) then
3121  do j = js, je
3122  do i = ie+1, ie+ihalo
3123  var(:,i,j) = var(:,ie,j)
3124  enddo
3125  enddo
3126  end if
3127 
3128  !--- copy inner data to HALO(West)
3129  if ( .NOT. prc_has_w ) then
3130  do j = js, je
3131  do i = is-ihalo, is-1
3132  var(:,i,j) = var(:,is,j)
3133  enddo
3134  enddo
3135  end if
3136 
3137  !--- copy inner data to HALO(NorthWest)
3138  if ( .NOT. prc_has_n .AND. &
3139  .NOT. prc_has_w ) then
3140  do j = je+1, je+jhalo
3141  do i = is-ihalo, is-1
3142  var(:,i,j) = var(:,is,je)
3143  enddo
3144  enddo
3145  elseif( .NOT. prc_has_n ) then
3146  do j = je+1, je+jhalo
3147  do i = is-ihalo, is-1
3148  var(:,i,j) = var(:,i,je)
3149  enddo
3150  enddo
3151  elseif( .NOT. prc_has_w ) then
3152  do j = je+1, je+jhalo
3153  do i = is-ihalo, is-1
3154  var(:,i,j) = var(:,is,j)
3155  enddo
3156  enddo
3157  endif
3158 
3159  !--- copy inner data to HALO(SouthWest)
3160  if ( .NOT. prc_has_s .AND. &
3161  .NOT. prc_has_w ) then
3162  do j = js-ihalo, js-1
3163  do i = is-ihalo, is-1
3164  var(:,i,j) = var(:,is,js)
3165  enddo
3166  enddo
3167  elseif( .NOT. prc_has_s ) then
3168  do j = js-ihalo, js-1
3169  do i = is-ihalo, is-1
3170  var(:,i,j) = var(:,i,js)
3171  enddo
3172  enddo
3173  elseif( .NOT. prc_has_w ) then
3174  do j = js-ihalo, js-1
3175  do i = is-ihalo, is-1
3176  var(:,i,j) = var(:,is,j)
3177  enddo
3178  enddo
3179  endif
3180 
3181  !--- copy inner data to HALO(NorthEast)
3182  if ( .NOT. prc_has_n .AND. &
3183  .NOT. prc_has_e ) then
3184  do j = je+1, je+jhalo
3185  do i = ie+1, ie+ihalo
3186  var(:,i,j) = var(:,ie,je)
3187  enddo
3188  enddo
3189  elseif( .NOT. prc_has_n ) then
3190  do j = je+1, je+jhalo
3191  do i = ie+1, ie+ihalo
3192  var(:,i,j) = var(:,i,je)
3193  enddo
3194  enddo
3195  elseif( .NOT. prc_has_e ) then
3196  do j = je+1, je+jhalo
3197  do i = ie+1, ie+ihalo
3198  var(:,i,j) = var(:,ie,j)
3199  enddo
3200  enddo
3201  endif
3202 
3203  !--- copy inner data to HALO(SouthEast)
3204  if ( .NOT. prc_has_s .AND. &
3205  .NOT. prc_has_e ) then
3206  do j = js-ihalo, js-1
3207  do i = ie+1, ie+ihalo
3208  var(:,i,j) = var(:,ie,js)
3209  enddo
3210  enddo
3211  elseif( .NOT. prc_has_s ) then
3212  do j = js-ihalo, js-1
3213  do i = ie+1, ie+ihalo
3214  var(:,i,j) = var(:,i,js)
3215  enddo
3216  enddo
3217  elseif( .NOT. prc_has_e ) then
3218  do j = js-ihalo, js-1
3219  do i = ie+1, ie+ihalo
3220  var(:,i,j) = var(:,ie,j)
3221  enddo
3222  enddo
3223  endif
3224 
3225  return
3226  end subroutine copy_boundary_3d
3227 
3228  subroutine copy_boundary_2d(var)
3229  implicit none
3230 
3231  real(RP), intent(inout) :: var(:,:)
3232 
3233  integer :: i, j
3234  !---------------------------------------------------------------------------
3235  !--- copy inner data to HALO(North)
3236  if( .NOT. prc_has_n ) then
3237  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
3238  do j = je+1, je+jhalo
3239  do i = is, ie
3240  var(i,j) = var(i,je)
3241  enddo
3242  enddo
3243  endif
3244 
3245  !--- copy inner data to HALO(South)
3246  if( .NOT. prc_has_s ) then
3247  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
3248  do j = js-jhalo, js-1
3249  do i = is, ie
3250  var(i,j) = var(i,js)
3251  enddo
3252  enddo
3253  endif
3254 
3255  if( .NOT. prc_has_e ) then
3256  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
3257  do j = js, je
3258  do i = ie+1, ie+ihalo
3259  var(i,j) = var(ie,j)
3260  enddo
3261  enddo
3262  endif
3263 
3264  if( .NOT. prc_has_w ) then
3265  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
3266  do j = js, je
3267  do i = is-ihalo, is-1
3268  var(i,j) = var(is,j)
3269  enddo
3270  enddo
3271  endif
3272 
3273  !--- copy inner data to HALO(NorthWest)
3274  if( .NOT. prc_has_n .AND. .NOT. prc_has_w ) then
3275  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
3276  do j = je+1, je+jhalo
3277  do i = is-ihalo, is-1
3278  var(i,j) = var(is,je)
3279  enddo
3280  enddo
3281  elseif( .NOT. prc_has_n ) then
3282  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
3283  do j = je+1, je+jhalo
3284  do i = is-ihalo, is-1
3285  var(i,j) = var(i,je)
3286  enddo
3287  enddo
3288  elseif( .NOT. prc_has_w ) then
3289  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
3290  do j = je+1, je+jhalo
3291  do i = is-ihalo, is-1
3292  var(i,j) = var(is,j)
3293  enddo
3294  enddo
3295  endif
3296 
3297  !--- copy inner data to HALO(SouthWest)
3298  if( .NOT. prc_has_s .AND. .NOT. prc_has_w ) then
3299  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
3300  do j = js-ihalo, js-1
3301  do i = is-ihalo, is-1
3302  var(i,j) = var(is,js)
3303  enddo
3304  enddo
3305  elseif( .NOT. prc_has_s ) then
3306  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
3307  do j = js-ihalo, js-1
3308  do i = is-ihalo, is-1
3309  var(i,j) = var(i,js)
3310  enddo
3311  enddo
3312  elseif( .NOT. prc_has_w ) then
3313  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
3314  do j = js-ihalo, js-1
3315  do i = is-ihalo, is-1
3316  var(i,j) = var(is,j)
3317  enddo
3318  enddo
3319  endif
3320 
3321  !--- copy inner data to HALO(NorthEast)
3322  if( .NOT. prc_has_n .AND. .NOT. prc_has_e ) then
3323  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
3324  do j = je+1, je+jhalo
3325  do i = ie+1, ie+ihalo
3326  var(i,j) = var(ie,je)
3327  enddo
3328  enddo
3329  elseif( .NOT. prc_has_n ) then
3330  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
3331  do j = je+1, je+jhalo
3332  do i = ie+1, ie+ihalo
3333  var(i,j) = var(i,je)
3334  enddo
3335  enddo
3336  elseif( .NOT. prc_has_e ) then
3337  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
3338  do j = je+1, je+jhalo
3339  do i = ie+1, ie+ihalo
3340  var(i,j) = var(ie,j)
3341  enddo
3342  enddo
3343  endif
3344 
3345  !--- copy inner data to HALO(SouthEast)
3346  if( .NOT. prc_has_s .AND. .NOT. prc_has_e ) then
3347  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
3348  do j = js-ihalo, js-1
3349  do i = ie+1, ie+ihalo
3350  var(i,j) = var(ie,js)
3351  enddo
3352  enddo
3353  elseif( .NOT. prc_has_s ) then
3354  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
3355  do j = js-ihalo, js-1
3356  do i = ie+1, ie+ihalo
3357  var(i,j) = var(i,js)
3358  enddo
3359  enddo
3360  elseif( .NOT. prc_has_e ) then
3361  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
3362  do j = js-ihalo, js-1
3363  do i = ie+1, ie+ihalo
3364  var(i,j) = var(ie,j)
3365  enddo
3366  enddo
3367  endif
3368 
3369  return
3370  end subroutine copy_boundary_2d
3371 
3372  subroutine comm_cleanup
3373  use mpi
3374  implicit none
3375 
3376  integer :: i, j, ierr
3377  !---------------------------------------------------------------------------
3378 
3379  deallocate( recvpack_w2p )
3380  deallocate( recvpack_e2p )
3381  deallocate( sendpack_p2w )
3382  deallocate( sendpack_p2e )
3383 #ifdef DEBUG
3384  deallocate( use_packbuf )
3385 #endif
3386 
3387  deallocate( req_cnt )
3388  deallocate( req_list )
3389 
3390  if ( comm_use_mpi_pc ) then
3391  do j=1, comm_vsize_max_pc
3392  do i=1, comm_nreq_max+1
3393  if (preq_list(i,j) .NE. mpi_request_null) &
3394  call mpi_request_free(preq_list(i,j), ierr)
3395  enddo
3396  enddo
3397  deallocate( preq_cnt )
3398  deallocate( preq_list )
3399  deallocate( pseqid )
3400  end if
3401  end subroutine comm_cleanup
3402 
3403 end module scale_comm_cartesc
subroutine vars_3d_mpi_pc(var, vid)
integer, public jmax
of computational cells: y, local
subroutine vars_init_mpi_pc(var, vid, seqid)
integer, parameter, public prc_sw
[node direction] southwest
integer, public comm_world
communication world ID
subroutine vars8_3d_mpi(var, vid)
integer, dimension(8), public prc_next
node ID of 8 neighbour process
subroutine comm_bcast_int_2d(var, gIA, gJA)
Broadcast data for whole process value in 2D field (integer)
integer, public imax
of computational cells: x, local
integer, parameter, public prc_n
[node direction] north
subroutine wait_3d_mpi(var, vid)
integer, public ia
of whole cells: x, local, with HALO
integer, public comm_datatype
datatype of variable
subroutine comm_bcast_1d(var, gIA)
Broadcast data for whole process value in 1D field.
subroutine comm_bcast_3d(var, gIA, gJA, gKA)
Broadcast data for whole process value in 3D field.
integer, public qa
integer, public ja
of whole cells: y, local, with HALO
subroutine vars8_2d_mpi(var, vid)
module process / cartesC
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
subroutine comm_gather_2d(recv, send, gIA, gJA)
Get data from whole process value in 2D field.
logical, public prc_has_s
integer, parameter, public prc_se
[node direction] southeast
logical, public prc_has_n
logical, public prc_has_e
real(rp), public const_undef
Definition: scale_const.F90:41
module COMMUNICATION
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
subroutine comm_bcast_2d(var, gIA, gJA)
Broadcast data for whole process value in 2D field.
module TRACER
subroutine comm_bcast_4d(var, gIA, gJA, gKA, gTime)
Broadcast data for whole process value in 4D field.
module atmosphere / grid / cartesC index
subroutine comm_horizontal_mean_3d(varmean, var)
calculate horizontal mean (global total with communication) 3D
module PROCESS
Definition: scale_prc.F90:11
integer, public je
end point of inner domain: y, local
subroutine comm_bcast_character(var)
Broadcast data for whole process value in character.
integer, parameter, public prc_w
[node direction] west
integer, parameter, public prc_masterrank
master process in each communicator
Definition: scale_prc.F90:65
subroutine, public comm_vars_init(varname, var, vid)
Register variables.
subroutine, public comm_vars8_init(varname, var, vid)
Register variables.
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module CONSTANT
Definition: scale_const.F90:11
integer, public js
start point of inner domain: y, local
subroutine vars_2d_mpi(var, vid)
subroutine comm_bcast_int_1d(var, gIA)
Broadcast data for whole process value in 1D field (integer)
integer, parameter, public prc_nw
[node direction] northwest
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:157
subroutine, public comm_cleanup
integer, parameter, public prc_e
[node direction] east
module profiler
Definition: scale_prof.F90:11
subroutine comm_bcast_int_scr(var)
Broadcast data for whole process value in scalar (integer)
module PRECISION
integer, public ka
of whole cells: z, local, with HALO
integer, parameter, public prc_s
[node direction] south
subroutine comm_gather_3d(recv, send, gIA, gJA, gKA)
Get data from whole process value in 3D field.
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:87
subroutine, public comm_setup
Setup.
module STDIO
Definition: scale_io.F90:10
subroutine comm_bcast_scr(var)
Broadcast data for whole process value in scalar field.
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:210
integer, parameter, public rp
subroutine comm_bcast_logical_scr(var)
Broadcast data for whole process value in scalar (logical)
integer, parameter, public prc_ne
[node direction] northeast
logical, public prc_has_w