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