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