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