SCALE-RM
scale_comm_cartesC_nest.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 #include "scalelib.h"
13  !-----------------------------------------------------------------------------
14  !
15  !++ used modules
16  !
17  use mpi
18  use scale_precision
19  use scale_io
20  use scale_prof
21  use scale_debug
23  use scale_index
24  use scale_tracer
25  !-----------------------------------------------------------------------------
26  implicit none
27  private
28  !-----------------------------------------------------------------------------
29  !
30  !++ Public procedure
31  !
32  public :: comm_cartesc_nest_setup
38  public :: comm_cartesc_nest_test
40 
41  !-----------------------------------------------------------------------------
42  !
43  !++ Public parameters & variables
44  !
45  integer, public :: intercomm_parent ! inter-communicator to parent
46  integer, public :: intercomm_daughter ! inter-communicator to daughter
47 
48  integer, public :: comm_cartesc_nest_filiation(10)
49  integer, public :: handling_num
50  integer, public :: comm_cartesc_nest_tile_num_x
51  integer, public :: comm_cartesc_nest_tile_num_y
52  integer, public, allocatable :: comm_cartesc_nest_tile_id(:)
53 
54  integer, public :: parent_kmax(2)
55  integer, public :: parent_imax(2)
56  integer, public :: parent_jmax(2)
57  integer, public :: parent_ka(2)
58  integer, public :: parent_ia(2)
59  integer, public :: parent_ja(2)
60  integer, public :: parent_okmax(2)
61  integer, public :: parent_lkmax(2)
62  real(DP), public :: parent_dtsec(2)
63  integer, public :: parent_nstep(2)
64 
65  integer, public :: daughter_kmax(2)
66  integer, public :: daughter_imax(2)
67  integer, public :: daughter_jmax(2)
68  integer, public :: daughter_ka(2)
69  integer, public :: daughter_ia(2)
70  integer, public :: daughter_ja(2)
71  integer, public :: daughter_okmax(2)
72  integer, public :: daughter_lkmax(2)
73  real(DP), public :: daughter_dtsec(2)
74  integer, public :: daughter_nstep(2)
75 
76  integer, public :: prnt_ks(2)
77  integer, public :: prnt_ke(2)
78  integer, public :: prnt_is(2)
79  integer, public :: prnt_ie(2)
80  integer, public :: prnt_js(2)
81  integer, public :: prnt_je(2)
82 
83  integer, public :: datr_ks(2)
84  integer, public :: datr_ke(2)
85  integer, public :: datr_is(2)
86  integer, public :: datr_ie(2)
87  integer, public :: datr_js(2)
88  integer, public :: datr_je(2)
89 
90  integer, public :: tileal_ka(2)
91  integer, public :: tileal_ia(2)
92  integer, public :: tileal_ja(2)
93 
94  integer, public :: comm_cartesc_nest_bnd_qa = 1
95  integer, public :: comm_cartesc_nest_interp_level = 5
97 
98  logical, public :: use_nesting = .false.
99  logical, public :: offline = .false.
100  logical, public :: online_iam_parent = .false.
101  logical, public :: online_iam_daughter = .false.
102  integer, public :: online_domain_num = 1
103  logical, public :: online_use_velz = .false.
104  logical, public :: online_no_rotate = .false.
105  logical, public :: online_boundary_use_qhyd = .false.
106  logical, public :: online_boundary_diagqhyd = .false.
107 
108  !-----------------------------------------------------------------------------
109  !
110  !++ Private procedure
111  !
112  private :: comm_cartesc_nest_parentsize
113  private :: comm_cartesc_nest_catalogue
114  private :: comm_cartesc_nest_ping
115  private :: comm_cartesc_nest_setup_nestdown
116  private :: comm_cartesc_nest_importgrid_nestdown
117  private :: comm_cartesc_nest_intercomm_nestdown
118  private :: comm_cartesc_nest_issuer_of_receive
119  private :: comm_cartesc_nest_issuer_of_wait
120 
121  interface comm_cartesc_nest_intercomm_nestdown
123  end interface comm_cartesc_nest_intercomm_nestdown
124 
125  interface comm_cartesc_nest_issuer_of_receive
127  end interface comm_cartesc_nest_issuer_of_receive
128 
129  interface comm_cartesc_nest_issuer_of_wait
130  module procedure comm_cartesc_nest_issuer_of_wait_3d
131  end interface comm_cartesc_nest_issuer_of_wait
132 
133  !-----------------------------------------------------------------------------
134  !
135  !++ Private parameters & variables
136  !
137  real(RP), private, allocatable :: latlon_catalog(:,:,:)
138  real(RP), private :: latlon_local (4,2)
139 
140  integer, private :: parent_prc_num_x(2)
141  integer, private :: parent_prc_num_y(2)
142  integer, private :: parent_prc_nprocs(2)
143 
144  integer, private :: daughter_prc_num_x(2)
145  integer, private :: daughter_prc_num_y(2)
146  integer, private :: daughter_prc_nprocs(2)
147 
148  integer, private :: comm_cartesc_nest_tile_all
149  integer, private :: comm_cartesc_nest_tile_allmax_p
150  integer, private :: comm_cartesc_nest_tile_allmax_d
151  integer, private, allocatable :: comm_cartesc_nest_tile_list_p(:,:)
152  integer, private, allocatable :: comm_cartesc_nest_tile_list_d(:,:)
153  integer, private, allocatable :: comm_cartesc_nest_tile_list_yp(:)
154  integer, private :: num_yp
155 
156  character(len=H_LONG), private :: offline_parent_basename
157  integer, private :: offline_parent_prc_num_x
158  integer, private :: offline_parent_prc_num_y
159  integer, private :: offline_parent_kmax
160  integer, private :: offline_parent_imax
161  integer, private :: offline_parent_jmax
162  integer, private :: offline_parent_lkmax
163  integer, private :: offline_parent_okmax
164  integer(8), private :: online_wait_limit
165  logical, private :: online_daughter_use_velz
166  logical, private :: online_daughter_no_rotate
167  logical, private :: online_aggressive_comm
168 
169  integer, parameter :: i_lon = 1
170  integer, parameter :: i_lat = 2
171 
172  integer, parameter :: i_min = 1
173  integer, parameter :: i_max = 2
174  integer, parameter :: i_bndqa = 20
175 
176  integer, parameter :: i_sclr = 1
177  integer, parameter :: i_zstg = 2
178  integer, parameter :: i_xstg = 3
179  integer, parameter :: i_ystg = 4
180 
181  integer, parameter :: itp_ng = 4
182  integer, private :: itp_nh = 4
183  integer, private :: itp_nv = 2
184 
185  integer, parameter :: tag_lon = 1
186  integer, parameter :: tag_lat = 2
187  integer, parameter :: tag_lonuy = 3
188  integer, parameter :: tag_latuy = 4
189  integer, parameter :: tag_lonxv = 5
190  integer, parameter :: tag_latxv = 6
191  integer, parameter :: tag_cz = 7
192  integer, parameter :: tag_fz = 8
193 
194  integer, parameter :: tag_dens = 1
195  integer, parameter :: tag_momz = 2
196  integer, parameter :: tag_momx = 3
197  integer, parameter :: tag_momy = 4
198  integer, parameter :: tag_rhot = 5
199  integer, parameter :: tag_qx = 6
200 
201  integer, parameter :: order_tag_comm = 100000
202  integer, parameter :: order_tag_var = 1000
203  ! intercomm tag id: IC | VAR | YP
204  ! (total: 6columns) X X X X X X
205 
206  integer, private, parameter :: interp_search_divnum = 10
207 
208  integer, private :: intercomm_id(2)
209 
210  integer, private, parameter :: max_isu = 100 ! maximum number of receive/wait issue
211  integer, private, parameter :: max_isuf = 20 ! maximum number of receive/wait issue (z-stag)
212  integer, private :: max_rq = 1000 ! maximum number of req: tentative approach
213  integer, private :: rq_ctl_p ! for control request id (counting)
214  integer, private :: rq_ctl_d ! for control request id (counting)
215  integer, private :: rq_tot_p ! for control request id (total number)
216  integer, private :: rq_tot_d ! for control request id (total number)
217  integer, private, allocatable :: ireq_p(:) ! buffer of request-id for parent
218  integer, private, allocatable :: ireq_d(:) ! buffer of request-id for daughter
219  integer, private, allocatable :: call_order(:) ! calling order from parent
220 
221  real(RP), private, allocatable :: buffer_2d (:,:) ! buffer of communicator: 2D (with HALO)
222  real(RP), private, allocatable :: buffer_3d (:,:,:) ! buffer of communicator: 3D (with HALO)
223  real(RP), private, allocatable :: buffer_3df (:,:,:) ! buffer of communicator: 3D-Kface (with HALO)
224  real(RP), private, allocatable :: recvbuf_3d (:,:,:,:) ! buffer of receiver: 3D (with HALO)
225  real(RP), private, allocatable :: recvbuf_3df(:,:,:,:) ! buffer of receiver: 3D-Kface (with HALO)
226 
227  real(RP), private, allocatable :: buffer_ref_lon (:,:) ! buffer of communicator: LON
228  real(RP), private, allocatable :: buffer_ref_lonuy(:,:) ! buffer of communicator: LONUY
229  real(RP), private, allocatable :: buffer_ref_lonxv(:,:) ! buffer of communicator: LONXV
230  real(RP), private, allocatable :: buffer_ref_lat (:,:) ! buffer of communicator: LAT
231  real(RP), private, allocatable :: buffer_ref_latuy(:,:) ! buffer of communicator: LATUY
232  real(RP), private, allocatable :: buffer_ref_latxv(:,:) ! buffer of communicator: LATXV
233  real(RP), private, allocatable :: buffer_ref_cz (:,:,:) ! buffer of communicator: CZ
234  real(RP), private, allocatable :: buffer_ref_fz (:,:,:) ! buffer of communicator: FZ
235 
236 
237  real(RP), private, allocatable :: buffer_ref_3d (:,:,:) ! buffer of communicator: 3D data (with HALO)
238  real(RP), private, allocatable :: buffer_ref_3df(:,:,:) ! buffer of communicator: 3D at z-Face (with HALO)
239 
240  real(RP), private, allocatable :: org_dens(:,:,:) ! buffer of communicator: DENS
241  real(RP), private, allocatable :: org_momz(:,:,:) ! buffer of communicator: MOMZ
242  real(RP), private, allocatable :: org_momx(:,:,:) ! buffer of communicator: MOMX
243  real(RP), private, allocatable :: org_momy(:,:,:) ! buffer of communicator: MOMY
244  real(RP), private, allocatable :: org_u_ll(:,:,:) ! buffer of communicator: U_ll
245  real(RP), private, allocatable :: org_v_ll(:,:,:) ! buffer of communicator: V_ll
246  real(RP), private, allocatable :: org_rhot(:,:,:) ! buffer of communicator: RHOT
247  real(RP), private, allocatable :: org_qtrc(:,:,:,:) ! buffer of communicator: QTRC
248 
249  integer, private, allocatable :: igrd (:,:,:,:) ! interpolation target grids in x-axis
250  integer, private, allocatable :: jgrd (:,:,:,:) ! interpolation target grids in y-axis
251  real(RP), private, allocatable :: hfact(:,:,:,:) ! interpolation factor for horizontal direction
252  integer, private, allocatable :: kgrd (:,:,:,:,:,:) ! interpolation target grids in z-axis
253  real(RP), private, allocatable :: vfact(:,:,:,:,:,:) ! interpolation factor for vertical direction
254 
255  integer(8), private :: nwait_p, nwait_d, nrecv, nsend
256 
257  character(len=H_SHORT) :: mp_type
258 
259  !-----------------------------------------------------------------------------
260 contains
261  !-----------------------------------------------------------------------------
263  subroutine comm_cartesc_nest_setup ( &
264  QA_MP, &
265  MP_TYPE_in, &
266  inter_parent, &
267  inter_child )
268  use scale_file, only: &
269  file_open, &
270  file_get_attribute, &
271  file_get_shape
272  use scale_const, only: &
273  d2r => const_d2r
274  use scale_prc, only: &
275  prc_abort, &
278  use scale_interp, only: &
279  interp_setup, &
281  use scale_comm_cartesc, only: &
282  comm_bcast
283  use scale_atmos_grid_cartesc_real, only: &
294  use scale_atmos_hydrometeor, only: &
296  implicit none
297 
298  integer, intent(in) :: QA_MP
299  character(len=*), intent(in) :: MP_TYPE_in
300  integer, intent(in), optional :: inter_parent
301  integer, intent(in), optional :: inter_child
302 
303 
304  character(len=H_LONG) :: LATLON_CATALOGUE_FNAME = 'latlon_domain_catalogue.txt'
305 
306  integer :: ONLINE_SPECIFIED_MAXRQ = 0
307  integer :: i
308  integer :: fid, ierr
309  integer :: parent_id
310 
311  logical :: flag_parent = .false.
312  logical :: flag_child = .false.
313 
314  integer :: imaxg(1), jmaxg(1)
315  integer :: pnum_x(1), pnum_y(1)
316  integer :: dims(1)
317  logical :: error, existed
318 
319  namelist / param_comm_cartesc_nest / &
320  latlon_catalogue_fname, &
321  offline_parent_basename, &
322  offline_parent_prc_num_x, &
323  offline_parent_prc_num_y, &
327  online_use_velz, &
330  online_aggressive_comm, &
331  online_wait_limit, &
332  online_specified_maxrq, &
335 
336  !---------------------------------------------------------------------------
337 
338  log_newline
339  log_info("COMM_CARTESC_NEST_setup",*) 'Setup'
340 
341  if ( present(inter_parent) ) then
342  if( inter_parent /= mpi_comm_null ) flag_child = .true. ! exist parent, so work as a child
343  endif
344  if ( present(inter_child) ) then
345  if( inter_child /= mpi_comm_null ) flag_parent = .true. ! exist child, so work as a parent
346  endif
347 
348  offline_parent_basename = ""
349 
350  nwait_p = 0
351  nwait_d = 0
352  nrecv = 0
353  nsend = 0
354 
355  handling_num = 0
357  online_wait_limit = 999999999
358  online_aggressive_comm = .true.
359 
360  !--- read namelist
361  rewind(io_fid_conf)
362  read(io_fid_conf,nml=param_comm_cartesc_nest,iostat=ierr)
363  if( ierr < 0 ) then !--- missing
364  log_info("COMM_CARTESC_NEST_setup",*) 'Not found namelist. Default used.'
365  elseif( ierr > 0 ) then !--- fatal error
366  log_error("COMM_CARTESC_NEST_setup",*) 'Not appropriate names in namelist PARAM_COMM_CARTESC_NEST. Check!'
367  call prc_abort
368  endif
369  log_nml(param_comm_cartesc_nest)
370 
372 
373  if ( offline_parent_basename /= "" ) then
374 
375  offline = .true.
376  use_nesting = .true.
377 
378  if ( prc_ismaster ) then
379  call file_open( offline_parent_basename, & ! (in)
380  fid, & ! (out)
381  aggregate = .false. ) ! (in)
382 
383  call file_get_attribute( fid, "global", "scale_atmos_grid_cartesC_index_imaxg", &
384  imaxg(:), existed=existed )
385  if ( existed ) then
386  call file_get_attribute( fid, "global", "scale_cartesC_prc_num_x", &
387  pnum_x(:), existed=existed )
388  end if
389  if ( existed ) then
390  offline_parent_imax = imaxg(1) / pnum_x(1)
391  else
392  ! for old file
393  call file_get_shape( fid, "CX", dims(:) )
394  offline_parent_imax = dims(1)-ihalo*2
395  end if
396 
397  call file_get_attribute( fid, "global", "scale_atmos_grid_cartesC_index_jmaxg", &
398  jmaxg(:), existed=existed )
399  if ( existed ) then
400  call file_get_attribute( fid, "global", "scale_cartesC_prc_num_y", &
401  pnum_y(:), existed=existed )
402  end if
403  if ( existed ) then
404  offline_parent_jmax = jmaxg(1) / pnum_y(1)
405  else
406  ! for old file
407  call file_get_shape( fid, "CY", dims(:) )
408  offline_parent_jmax = dims(1)-jhalo*2
409  end if
410 
411  call file_get_attribute( fid, "global", "scale_atmos_grid_cartesC_index_kmax", &
412  dims(:), existed=existed )
413  if ( existed ) then
414  offline_parent_kmax = dims(1)
415  else
416  call file_get_shape( fid, "z", dims(:), error=error )
417  if ( error ) then
418  offline_parent_kmax = 0
419  else
420  offline_parent_kmax = dims(1)
421  endif
422  end if
423 
424  call file_get_attribute( fid, "global", "scale_ocean_grid_cartesC_index_kmax", &
425  dims(:), existed=existed )
426  if ( existed ) then
427  offline_parent_okmax = dims(1)
428  else
429  ! for old file
430  call file_get_shape( fid, "oz", dims(:), error=error )
431  if ( error ) then
432  offline_parent_okmax = 0
433  else
434  offline_parent_okmax = dims(1)
435  endif
436  end if
437 
438  call file_get_attribute( fid, "global", "scale_land_grid_cartesC_index_kmax", &
439  dims(:), existed=existed )
440  if ( existed ) then
441  offline_parent_lkmax = dims(1)
442  else
443  ! for old file
444  call file_get_shape( fid, "lz", dims(:), error=error )
445  if ( error ) then
446  offline_parent_lkmax = 0
447  else
448  offline_parent_lkmax = dims(1)
449  endif
450  end if
451 
452  endif
453  call comm_bcast( offline_parent_imax )
454  call comm_bcast( offline_parent_jmax )
455  call comm_bcast( offline_parent_kmax )
456  call comm_bcast( offline_parent_okmax )
457  call comm_bcast( offline_parent_lkmax )
458  endif
459 
460  if ( online_iam_daughter .or. online_iam_parent ) then
461 
462  if ( offline ) then
463  log_error("COMM_CARTESC_NEST_setup",*) 'OFFLINE and ONLINE cannot be use at the same time'
464  call prc_abort
465  endif
466 
467  use_nesting = .true.
468  endif
469 
471 
473  itp_nv = 2
474 
476  if( online_specified_maxrq > max_rq ) max_rq = online_specified_maxrq
477 
478  allocate( ireq_p(max_rq) )
479  allocate( ireq_d(max_rq) )
480  allocate( call_order(max_rq) )
481  ireq_p(:) = mpi_request_null
482  ireq_d(:) = mpi_request_null
483 
484  if ( use_nesting ) then
485 
486  if ( offline .OR. online_iam_daughter ) then
487  latlon_local(i_min,i_lon) = minval(atmos_grid_cartesc_real_lonuv(:,:)) / d2r
488  latlon_local(i_max,i_lon) = maxval(atmos_grid_cartesc_real_lonuv(:,:)) / d2r
489  latlon_local(i_min,i_lat) = minval(atmos_grid_cartesc_real_latuv(:,:)) / d2r
490  latlon_local(i_max,i_lat) = maxval(atmos_grid_cartesc_real_latuv(:,:)) / d2r
491  endif
492 
493  if ( offline ) then
494 
495  handling_num = 1
496  parent_prc_num_x(handling_num) = offline_parent_prc_num_x
497  parent_prc_num_y(handling_num) = offline_parent_prc_num_y
498  parent_kmax(handling_num) = offline_parent_kmax
499  parent_imax(handling_num) = offline_parent_imax
500  parent_jmax(handling_num) = offline_parent_jmax
501  parent_okmax(handling_num) = offline_parent_okmax
502  parent_lkmax(handling_num) = offline_parent_lkmax
503 
504  parent_prc_nprocs(handling_num) = parent_prc_num_x(handling_num) * parent_prc_num_y(handling_num)
505  allocate( latlon_catalog(parent_prc_nprocs(handling_num),2,2) )
506 
507  !--- read latlon catalogue
508  fid = io_get_available_fid()
509  open( fid, &
510  file = trim(latlon_catalogue_fname), &
511  form = 'formatted', &
512  status = 'old', &
513  iostat = ierr )
514 
515  if ( ierr /= 0 ) then
516  log_error("COMM_CARTESC_NEST_setup",*) 'cannot open latlon-catalogue file!'
517  call prc_abort
518  endif
519 
520  do i = 1, parent_prc_nprocs(handling_num)
521  read(fid,'(i8,4f32.24)',iostat=ierr) parent_id, &
522  latlon_catalog(i,i_min,i_lon), latlon_catalog(i,i_max,i_lon), & ! LON: MIN, MAX
523  latlon_catalog(i,i_min,i_lat), latlon_catalog(i,i_max,i_lat) ! LAT: MIN, MAX
524  if ( i /= parent_id ) then
525  log_error("COMM_CARTESC_NEST_setup",*) 'internal error: parent mpi id'
526  call prc_abort
527  endif
528  if ( ierr /= 0 ) exit
529  enddo
530  close(fid)
531 
533 
534  else ! ONLINE RELATIONSHIP
535 ! if ( present(flag_parent) .AND. present(flag_child) ) then
536 ! LOG_INFO("COMM_CARTESC_NEST_setup",'(1x,A)') &
537 ! '*** Setup Online Nesting Inter-Domain Communicator (IDC)'
538 ! else
539 ! LOG_ERROR("COMM_CARTESC_NEST_setup",*) 'Internal Error:'
540 ! LOG_ERROR_CONT(*) 'The flag_parent and flag_child are needed.'
541 ! LOG_WARN("COMM_CARTESC_NEST_setup",*) ' domain: ', ONLINE_DOMAIN_NUM
542 ! call PRC_abort
543 ! endif
544 
545  if( online_boundary_use_qhyd ) then
546  mp_type = mp_type_in
548  elseif ( atmos_hydrometeor_dry ) then
549  mp_type = "DRY"
551  else
552  mp_type = "QV"
554  endif
555 
556  log_info("COMM_CARTESC_NEST_setup",*) "flag_parent", flag_parent, "flag_child", flag_child
557  log_info("COMM_CARTESC_NEST_setup",*) "ONLINE_IAM_PARENT", online_iam_parent, "ONLINE_IAM_DAUGHTER", online_iam_daughter
558 
559  if( flag_parent ) then ! must do first before daughter processes
560  !-------------------------------------------------
561  if ( .NOT. online_iam_parent ) then
562  log_error("COMM_CARTESC_NEST_setup",*) '[NEST_setup] Parent Flag from launcher is not consistent with namelist!'
563  log_error_cont(*) 'PARENT - domain : ', online_domain_num
564  call prc_abort
565  endif
566 
567  handling_num = 1 !HANDLING_NUM + 1
568  intercomm_id(handling_num) = online_domain_num
569  comm_cartesc_nest_filiation(intercomm_id(handling_num)) = 1
570 
571  intercomm_daughter = inter_child
572  log_info("COMM_CARTESC_NEST_setup",'(1x,A,I2,A)') 'Online Nesting - PARENT [INTERCOMM_ID:', &
573  intercomm_id(handling_num), ' ]'
574  log_info("COMM_CARTESC_NEST_setup",*) 'Online Nesting - INTERCOMM :', intercomm_daughter
575 
576  call comm_cartesc_nest_ping( handling_num )
577 
578  call comm_cartesc_nest_parentsize( handling_num )
579 
580  call comm_cartesc_nest_catalogue( handling_num )
581  call mpi_barrier(intercomm_daughter, ierr)
582 
592 
593  log_info("COMM_CARTESC_NEST_setup",'(1x,A)' ) 'Informations of Parent Domain [me]'
594  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_nprocs :', parent_prc_nprocs(handling_num)
595  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_NUM_X :', parent_prc_num_x(handling_num)
596  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_NUM_Y :', parent_prc_num_y(handling_num)
597  log_info_cont('(1x,A,I6)' ) '--- PARENT_KMAX :', parent_kmax(handling_num)
598  log_info_cont('(1x,A,I6)' ) '--- PARENT_IMAX :', parent_imax(handling_num)
599  log_info_cont('(1x,A,I6)' ) '--- PARENT_JMAX :', parent_jmax(handling_num)
600  log_info_cont('(1x,A,F9.3)') '--- PARENT_DTSEC :', parent_dtsec(handling_num)
601  log_info_cont('(1x,A,I6) ') '--- PARENT_NSTEP :', parent_nstep(handling_num)
602  log_info_cont('(1x,A)' ) 'Informations of Daughter Domain'
603  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_PRC_nprocs :', daughter_prc_nprocs(handling_num)
604  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_PRC_NUM_X :', daughter_prc_num_x(handling_num)
605  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_PRC_NUM_Y :', daughter_prc_num_y(handling_num)
606  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_KMAX :', daughter_kmax(handling_num)
607  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_IMAX :', daughter_imax(handling_num)
608  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_JMAX :', daughter_jmax(handling_num)
609  log_info_cont('(1x,A,F9.3)') '--- DAUGHTER_DTSEC :', daughter_dtsec(handling_num)
610  log_info_cont('(1x,A,I6) ') '--- DAUGHTER_NSTEP :', daughter_nstep(handling_num)
611  log_info_cont('(1x,A,I6) ') 'Limit Num. NCOMM req. :', max_rq
612 
621 
622  call comm_cartesc_nest_setup_nestdown( handling_num )
623 
624  !---------------------------------- end of parent routines
625  endif
626 
627 
628  if( flag_child ) then
629  !-------------------------------------------------
630  if ( .NOT. online_iam_daughter ) then
631  log_error("COMM_CARTESC_NEST_setup",*) '[NEST_setup] Child Flag from launcher is not consistent with namelist!'
632  log_error_cont(*) 'DAUGHTER - domain : ', online_domain_num
633  call prc_abort
634  endif
635 
636  handling_num = 2 !HANDLING_NUM + 1
637  intercomm_id(handling_num) = online_domain_num - 1
638  comm_cartesc_nest_filiation(intercomm_id(handling_num)) = -1
639 
640  intercomm_parent = inter_parent
641  log_info("COMM_CARTESC_NEST_setup",'(1x,A,I2,A)') 'Online Nesting - DAUGHTER [INTERCOMM_ID:', &
642  intercomm_id(handling_num), ' ]'
643  log_info("COMM_CARTESC_NEST_setup",*) 'Online Nesting - INTERCOMM :', intercomm_parent
644 
645  call comm_cartesc_nest_ping( handling_num )
646 
647  call comm_cartesc_nest_parentsize( handling_num )
648 
649  allocate( latlon_catalog(parent_prc_nprocs(handling_num),2,2) )
650  call comm_cartesc_nest_catalogue( handling_num )
651  call mpi_barrier(intercomm_parent, ierr)
652 
654 
664 
665  log_info("COMM_CARTESC_NEST_setup",'(1x,A)' ) 'Informations of Parent Domain'
666  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_nprocs :', parent_prc_nprocs(handling_num)
667  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_NUM_X :', parent_prc_num_x(handling_num)
668  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_NUM_Y :', parent_prc_num_y(handling_num)
669  log_info_cont('(1x,A,I6)' ) '--- PARENT_KMAX :', parent_kmax(handling_num)
670  log_info_cont('(1x,A,I6)' ) '--- PARENT_IMAX :', parent_imax(handling_num)
671  log_info_cont('(1x,A,I6)' ) '--- PARENT_JMAX :', parent_jmax(handling_num)
672  log_info_cont('(1x,A,F9.3)') '--- PARENT_DTSEC :', parent_dtsec(handling_num)
673  log_info_cont('(1x,A,I6)' ) '--- PARENT_NSTEP :', parent_nstep(handling_num)
674  log_info_cont('(1x,A)' ) 'Informations of Daughter Domain [me]'
675  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_PRC_nprocs :', daughter_prc_nprocs(handling_num)
676  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_PRC_NUM_X :', daughter_prc_num_x(handling_num)
677  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_PRC_NUM_Y :', daughter_prc_num_y(handling_num)
678  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_KMAX :', daughter_kmax(handling_num)
679  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_IMAX :', daughter_imax(handling_num)
680  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_JMAX :', daughter_jmax(handling_num)
681  log_info_cont('(1x,A,F9.3)') '--- DAUGHTER_DTSEC :', daughter_dtsec(handling_num)
682  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_NSTEP :', daughter_nstep(handling_num)
683  log_info_cont('(1x,A)' ) 'Informations of Target Tiles'
684  log_info_cont('(1x,A,I6)' ) '--- TILEALL_KA :', tileal_ka(handling_num)
685  log_info_cont('(1x,A,I6)' ) '--- TILEALL_IA :', tileal_ia(handling_num)
686  log_info_cont('(1x,A,I6)' ) '--- TILEALL_JA :', tileal_ja(handling_num)
687  log_info_cont('(1x,A,I6) ') 'Limit Num. NCOMM req. :', max_rq
688 
689  allocate( buffer_2d( parent_ia(handling_num), parent_ja(handling_num) ) )
690  allocate( buffer_3d( parent_ka(handling_num), parent_ia(handling_num), parent_ja(handling_num) ) )
691  allocate( buffer_3df( 0:parent_ka(handling_num), parent_ia(handling_num), parent_ja(handling_num) ) )
692 
693  allocate( recvbuf_3d( parent_ka(handling_num), parent_ia(handling_num), parent_ja(handling_num), max_isu ) )
694  allocate( recvbuf_3df( 0:parent_ka(handling_num), parent_ia(handling_num), parent_ja(handling_num), max_isuf ) )
695 
696  allocate( buffer_ref_lon( tileal_ia(handling_num),tileal_ja(handling_num)) )
697  allocate( buffer_ref_lonuy( tileal_ia(handling_num),tileal_ja(handling_num)) )
698  allocate( buffer_ref_lonxv( tileal_ia(handling_num),tileal_ja(handling_num)) )
699  allocate( buffer_ref_lat( tileal_ia(handling_num),tileal_ja(handling_num)) )
700  allocate( buffer_ref_latuy( tileal_ia(handling_num),tileal_ja(handling_num)) )
701  allocate( buffer_ref_latxv( tileal_ia(handling_num),tileal_ja(handling_num)) )
702 
703  allocate( buffer_ref_cz( tileal_ka(handling_num),tileal_ia(handling_num),tileal_ja(handling_num)) )
704  allocate( buffer_ref_fz(0:tileal_ka(handling_num),tileal_ia(handling_num),tileal_ja(handling_num)) )
705 
706  allocate( buffer_ref_3d( tileal_ka(handling_num),tileal_ia(handling_num),tileal_ja(handling_num)) )
707  allocate( buffer_ref_3df(0:tileal_ka(handling_num),tileal_ia(handling_num),tileal_ja(handling_num)) )
708 
709  allocate( igrd( daughter_ia(handling_num),daughter_ja(handling_num),itp_nh,itp_ng) )
710  allocate( jgrd( daughter_ia(handling_num),daughter_ja(handling_num),itp_nh,itp_ng) )
711  allocate( hfact( daughter_ia(handling_num),daughter_ja(handling_num),itp_nh,itp_ng) )
712  allocate( kgrd(daughter_ka(handling_num),itp_nv,daughter_ia(handling_num),daughter_ja(handling_num),itp_nh,itp_ng) )
713  allocate( vfact(daughter_ka(handling_num),itp_nv,daughter_ia(handling_num),daughter_ja(handling_num),itp_nh,itp_ng) )
714 
715  call comm_cartesc_nest_setup_nestdown( handling_num )
716 
717 
718  ! for scalar points
719  call interp_factor3d( itp_nh, & ! [IN]
720  tileal_ka(handling_num), & ! [IN]
721  khalo+1, & ! [IN]
722  tileal_ka(handling_num)-khalo, & ! [IN]
723  tileal_ia(handling_num), & ! [IN]
724  tileal_ja(handling_num), & ! [IN]
725  buffer_ref_lon(:,:), & ! [IN]
726  buffer_ref_lat(:,:), & ! [IN]
727  buffer_ref_cz(:,:,:), & ! [IN]
728  daughter_ka(handling_num), & ! [IN]
729  datr_ks(handling_num), & ! [IN]
730  datr_ke(handling_num), & ! [IN]
731  daughter_ia(handling_num), & ! [IN]
732  daughter_ja(handling_num), & ! [IN]
733  atmos_grid_cartesc_real_lon(:,:), & ! [IN]
734  atmos_grid_cartesc_real_lat(:,:), & ! [IN]
735  atmos_grid_cartesc_real_cz(:,:,:), & ! [IN]
736  igrd( :,:,:,i_sclr), & ! [OUT]
737  jgrd( :,:,:,i_sclr), & ! [OUT]
738  hfact( :,:,:,i_sclr), & ! [OUT]
739  kgrd(:,:,:,:,:,i_sclr), & ! [OUT]
740  vfact(:,:,:,:,:,i_sclr) ) ! [OUT]
741 
742  ! for z staggered points
743  call interp_factor3d( itp_nh, & ! [IN]
744  tileal_ka(handling_num)+1, & ! [IN]
745  khalo+1, & ! [IN]
746  tileal_ka(handling_num)+1-khalo, & ! [IN]
747  tileal_ia(handling_num), & ! [IN]
748  tileal_ja(handling_num), & ! [IN]
749  buffer_ref_lon(:,:), & ! [IN]
750  buffer_ref_lat(:,:), & ! [IN]
751  buffer_ref_fz(:,:,:), & ! [IN]
752  daughter_ka(handling_num), & ! [IN]
753  datr_ks(handling_num), & ! [IN]
754  datr_ke(handling_num), & ! [IN]
755  daughter_ia(handling_num), & ! [IN]
756  daughter_ja(handling_num), & ! [IN]
757  atmos_grid_cartesc_real_lon(:,:), & ! [IN]
758  atmos_grid_cartesc_real_lat(:,:), & ! [IN]
759  atmos_grid_cartesc_real_fz(1:ka,:,:), & ! [IN]
760  igrd( :,:,:,i_zstg), & ! [OUT]
761  jgrd( :,:,:,i_zstg), & ! [OUT]
762  hfact( :,:,:,i_zstg), & ! [OUT]
763  kgrd(:,:,:,:,:,i_zstg), & ! [OUT]
764  vfact(:,:,:,:,:,i_zstg) ) ! [OUT]
765 
766  ! for x staggered points
767  call interp_factor3d( itp_nh, & ! [IN]
768  tileal_ka(handling_num), & ! [IN]
769  khalo+1, & ! [IN]
770  tileal_ka(handling_num)-khalo, & ! [IN]
771  tileal_ia(handling_num), & ! [IN]
772  tileal_ja(handling_num), & ! [IN]
773  buffer_ref_lonuy(:,:), & ! [IN]
774  buffer_ref_latxv(:,:), & ! [IN]
775  buffer_ref_cz(:,:,:), & ! [IN]
776  daughter_ka(handling_num), & ! [IN]
777  datr_ks(handling_num), & ! [IN]
778  datr_ke(handling_num), & ! [IN]
779  daughter_ia(handling_num), & ! [IN]
780  daughter_ja(handling_num), & ! [IN]
781  atmos_grid_cartesc_real_lonuy(1:ia,1:ja), & ! [IN]
782  atmos_grid_cartesc_real_latuy(1:ia,1:ja), & ! [IN]
783  atmos_grid_cartesc_real_cz(:,:,:), & ! [IN]
784  igrd( :,:,:,i_xstg), & ! [OUT]
785  jgrd( :,:,:,i_xstg), & ! [OUT]
786  hfact( :,:,:,i_xstg), & ! [OUT]
787  kgrd(:,:,:,:,:,i_xstg), & ! [OUT]
788  vfact(:,:,:,:,:,i_xstg) ) ! [OUT]
789 
790  ! for y staggered points
791  call interp_factor3d( itp_nh, & ! [IN]
792  tileal_ka(handling_num), & ! [IN]
793  khalo+1, & ! [IN]
794  tileal_ka(handling_num)-khalo, & ! [IN]
795  tileal_ia(handling_num), & ! [IN]
796  tileal_ja(handling_num), & ! [IN]
797  buffer_ref_lonxv(:,:), & ! [IN]
798  buffer_ref_latxv(:,:), & ! [IN]
799  buffer_ref_cz(:,:,:), & ! [IN]
800  daughter_ka(handling_num), & ! [IN]
801  datr_ks(handling_num), & ! [IN]
802  datr_ke(handling_num), & ! [IN]
803  daughter_ia(handling_num), & ! [IN]
804  daughter_ja(handling_num), & ! [IN]
805  atmos_grid_cartesc_real_lonxv(1:ia,1:ja), & ! [IN]
806  atmos_grid_cartesc_real_latxv(1:ia,1:ja), & ! [IN]
807  atmos_grid_cartesc_real_cz(:,:,:), & ! [IN]
808  igrd( :,:,:,i_ystg), & ! [OUT]
809  jgrd( :,:,:,i_ystg), & ! [OUT]
810  hfact( :,:,:,i_ystg), & ! [OUT]
811  kgrd(:,:,:,:,:,i_ystg), & ! [OUT]
812  vfact(:,:,:,:,:,i_ystg) ) ! [OUT]
813 
814  deallocate( buffer_2d )
815  deallocate( buffer_3d )
816  deallocate( buffer_3df )
817 
818  else
819  online_use_velz = .false.
820  endif
821 
822  !LOG_INFO("COMM_CARTESC_NEST_setup",'(1x,A,I2)') 'Number of Related Domains :', HANDLING_NUM
823  !if ( HANDLING_NUM > 2 ) then
824  ! f( IO_L ) LOG_ERROR("COMM_CARTESC_NEST_setup",*) 'Too much handing domains (up to 2)'
825  ! call PRC_abort
826  !endif
827 
828  endif !--- OFFLINE or NOT
829 
830  endif !--- USE_NESTING
831 
832  return
833  end subroutine comm_cartesc_nest_setup
834 
835  !-----------------------------------------------------------------------------
837  subroutine comm_cartesc_nest_domain_relate( &
838  HANDLE )
839  use scale_prc, only: &
840  prc_myrank, &
841  prc_abort
842  implicit none
843 
844  integer, intent(in) :: HANDLE
845 
846  integer :: x_min, x_max
847  integer :: y_min, y_max
848  logical :: hit(2,2)
849  integer :: p, i, j
850  !---------------------------------------------------------------------------
851 
852  if( .NOT. use_nesting ) return
853 
854  x_min = parent_prc_num_x(handle)
855  x_max = -1
856  y_min = parent_prc_num_y(handle)
857  y_max = -1
858  hit(:,:) = .false.
859  do p = 1, parent_prc_nprocs(handle)
860  if ( ( ( latlon_local(i_min,i_lon) >= latlon_catalog(p,i_min,i_lon) &
861  .AND. latlon_local(i_min,i_lon) <= latlon_catalog(p,i_max,i_lon) ) .OR. &
862  ( latlon_local(i_max,i_lon) >= latlon_catalog(p,i_min,i_lon) &
863  .AND. latlon_local(i_max,i_lon) <= latlon_catalog(p,i_max,i_lon) ) .OR. &
864  ( latlon_catalog(p,i_min,i_lon) >= latlon_local(i_min,i_lon) &
865  .AND. latlon_catalog(p,i_min,i_lon) <= latlon_local(i_max,i_lon) ) .OR. &
866  ( latlon_catalog(p,i_max,i_lon) >= latlon_local(i_min,i_lon) &
867  .AND. latlon_catalog(p,i_max,i_lon) <= latlon_local(i_max,i_lon) ) ) .AND. &
868  ( ( latlon_local(i_min,i_lat) >= latlon_catalog(p,i_min,i_lat) &
869  .AND. latlon_local(i_min,i_lat) <= latlon_catalog(p,i_max,i_lat) ) .OR. &
870  ( latlon_local(i_max,i_lat) >= latlon_catalog(p,i_min,i_lat) &
871  .AND. latlon_local(i_max,i_lat) <= latlon_catalog(p,i_max,i_lat) ) .OR. &
872  ( latlon_catalog(p,i_min,i_lat) >= latlon_local(i_min,i_lat) &
873  .AND. latlon_catalog(p,i_min,i_lat) <= latlon_local(i_max,i_lat) ) .OR. &
874  ( latlon_catalog(p,i_max,i_lat) >= latlon_local(i_min,i_lat) &
875  .AND. latlon_catalog(p,i_max,i_lat) <= latlon_local(i_max,i_lat) ) ) ) then
876  if ( latlon_catalog(p,i_min,i_lon) <= latlon_local(i_min,i_lon) ) hit(i_min,i_lon) = .true.
877  if ( latlon_catalog(p,i_max,i_lon) >= latlon_local(i_max,i_lon) ) hit(i_max,i_lon) = .true.
878  if ( latlon_catalog(p,i_min,i_lat) <= latlon_local(i_min,i_lat) ) hit(i_min,i_lat) = .true.
879  if ( latlon_catalog(p,i_max,i_lat) >= latlon_local(i_max,i_lat) ) hit(i_max,i_lat) = .true.
880  i = mod(p-1, parent_prc_num_x(handle))
881  j = (p-1) / parent_prc_num_x(handle)
882  if ( i < x_min ) x_min = i
883  if ( i > x_max ) x_max = i
884  if ( j < y_min ) y_min = j
885  if ( j > y_max ) y_max = j
886  end if
887  end do
888 
889  if ( .not. ( hit(i_min,i_lon) .and. hit(i_max,i_lon) .and. hit(i_min,i_lat) .and. hit(i_max,i_lat) ) ) then
890  log_error("COMM_CARTESC_NEST_domain_relate",*) 'region of daughter domain is larger than that of parent'
891  log_error_cont(*) ' at rank:', prc_myrank, ' of domain:', online_domain_num
892  log_error_cont(*) 'LON MIN: ',hit(i_min,i_lon), ', LON MAX: ',hit(i_max,i_lon), ', LAT MIN: ',hit(i_min,i_lat), ', LAT MAX: ',hit(i_max,i_lat)
893  log_error_cont('(A,F12.6,1x,F12.6)') 'daughter local (me) MIN-MAX: LON=', &
894  latlon_local(i_min,i_lon), latlon_local(i_max,i_lon)
895  do p = 1, parent_prc_nprocs(handle)
896  log_error_cont('(A,I5,A,F12.6,1x,F12.6)') ' parent (', p,') MIN-MAX: LON=', &
897  latlon_catalog(p,i_min,i_lon) ,latlon_catalog(p,i_max,i_lon)
898  enddo
899  log_error_cont('(A,F12.6,1x,F12.6)') 'daughter local (me): MIN-MAX LAT=', &
900  latlon_local(i_min,i_lat), latlon_local(i_max,i_lat)
901  do p = 1, parent_prc_nprocs(handle)
902  log_error_cont('(A,I5,A,F12.6,1x,F12.6)') ' parent (', p,') MIN-MAX: LAT=', &
903  latlon_catalog(p,i_min,i_lat) ,latlon_catalog(p,i_max,i_lat)
904  enddo
905  call prc_abort
906  end if
907 
908 
909 
910  comm_cartesc_nest_tile_num_x = x_max - x_min + 1
911  comm_cartesc_nest_tile_num_y = y_max - y_min + 1
912 
914 
915  log_info("COMM_CARTESC_NEST_domain_relate",'(1x,A)') 'NEST: target process tile in parent domain'
916  p = 1
919  comm_cartesc_nest_tile_id(p) = x_min + i - 1 + (y_min + j - 1) * parent_prc_num_x(handle)
920  log_info_cont('(1x,A,I4,A,I6)') '(', p, ') target mpi-process:', comm_cartesc_nest_tile_id(p)
921  p = p + 1
922  enddo
923  enddo
924 
925  return
926  end subroutine comm_cartesc_nest_domain_relate
927 
928  !-----------------------------------------------------------------------------
930  ! including definition array size with BND or not in Parent domain
931  subroutine comm_cartesc_nest_domain_shape ( &
932  tilei, &
933  tilej, &
934  cxs, cxe, &
935  cys, cye, &
936  pxs, pxe, &
937  pys, pye, &
938  iloc )
939  implicit none
940 
941  integer, intent(out) :: tilei, tilej
942  integer, intent(out) :: cxs, cxe, cys, cye
943  integer, intent(out) :: pxs, pxe, pys, pye
944  integer, intent(in) :: iloc ! rank number; start from 1
945 
946  integer :: hdl = 1 ! handler number
947  integer :: rank
948  integer :: xloc, yloc
949  integer :: xlocg, ylocg ! location over whole domain
950  !---------------------------------------------------------------------------
951 
952  if( .NOT. use_nesting ) return
953 
954  rank = comm_cartesc_nest_tile_id(iloc)
955  xloc = mod( iloc-1, comm_cartesc_nest_tile_num_x ) + 1
956  yloc = int( real(iloc-1) / real(COMM_CARTESC_NEST_TILE_NUM_X) ) + 1
957  xlocg = mod( rank, offline_parent_prc_num_x ) + 1
958  ylocg = int( real(rank) / real(OFFLINE_PARENT_PRC_NUM_X) ) + 1
959  tilei = parent_imax(hdl)
960  tilej = parent_jmax(hdl)
961 
962  cxs = tilei * (xloc-1) + 1
963  cxe = tilei * xloc
964  cys = tilej * (yloc-1) + 1
965  cye = tilej * yloc
966  pxs = 1
967  pxe = tilei
968  pys = 1
969  pye = tilej
970 
971  if ( xlocg == 1 ) then ! BND_W
972  tilei = tilei + 2
973  pxs = pxs + 2
974  pxe = pxe + 2
975  endif
976  if ( xlocg == offline_parent_prc_num_x ) then ! BND_E
977  tilei = tilei + 2
978  endif
979  if ( ylocg == 1 ) then ! BND_S
980  tilej = tilej + 2
981  pys = pys + 2
982  pye = pye + 2
983  endif
984  if ( ylocg == offline_parent_prc_num_y ) then ! BND_N
985  tilej = tilej + 2
986  endif
987 
988  return
989  end subroutine comm_cartesc_nest_domain_shape
990 
991  !-----------------------------------------------------------------------------
993  subroutine comm_cartesc_nest_parentsize( &
994  HANDLE )
995  use scale_prc, only: &
996  prc_abort, &
997  prc_nprocs, &
998  prc_myrank, &
1000  use scale_prc_cartesc, only: &
1001  prc_num_x, &
1002  prc_num_y
1003  use scale_time, only: &
1004  time_nstep, &
1005  time_dtsec
1006  use scale_comm_cartesc, only: &
1007  comm_bcast
1008  use scale_atmos_hydrometeor, only: &
1009  n_hyd
1010  implicit none
1011 
1012  integer, intent(in) :: HANDLE
1013 
1014  real(RP) :: buffer
1015  integer :: datapack(14)
1016 
1017  integer :: QA_OTHERSIDE
1018  character(len=H_SHORT) :: MP_TYPE_OTHERSIDE
1019 
1020  integer :: ireq1, ireq2, ireq3, ierr1, ierr2, ierr3, ileng
1021  integer :: istatus(mpi_status_size)
1022  integer :: tag
1023  !---------------------------------------------------------------------------
1024 
1025  if( .NOT. use_nesting ) return
1026 
1027  tag = intercomm_id(handle) * 100
1028  ileng = 14
1029 
1030  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
1031 
1032  !##### parent ####
1033 
1034  ! from parent to daughter
1035  datapack( 1) = prc_nprocs
1036  datapack( 2) = prc_num_x
1037  datapack( 3) = prc_num_y
1038  datapack( 4) = kmax
1039  datapack( 5) = imax
1040  datapack( 6) = jmax
1041  datapack( 7) = ks
1042  datapack( 8) = ke
1043  datapack( 9) = is
1044  datapack(10) = ie
1045  datapack(11) = js
1046  datapack(12) = je
1047  datapack(13) = time_nstep
1048  datapack(14) = comm_cartesc_nest_bnd_qa
1049  buffer = time_dtsec
1050 
1051  if ( prc_ismaster ) then
1052  call mpi_isend(datapack, ileng, mpi_integer, prc_myrank, tag, intercomm_daughter, ireq1, ierr1)
1053  call mpi_isend(buffer, 1, mpi_double_precision, prc_myrank, tag+1, intercomm_daughter, ireq2, ierr2)
1054  call mpi_isend(mp_type, h_short, mpi_character, prc_myrank, tag+2, intercomm_daughter, ireq3, ierr3)
1055  call mpi_wait(ireq1, istatus, ierr1)
1056  call mpi_wait(ireq2, istatus, ierr2)
1057  call mpi_wait(ireq3, istatus, ierr3)
1058  endif
1059 
1060  parent_prc_nprocs(handle) = datapack( 1)
1061  parent_prc_num_x(handle) = datapack( 2)
1062  parent_prc_num_y(handle) = datapack( 3)
1063  parent_kmax(handle) = datapack( 4)
1064  parent_imax(handle) = datapack( 5)
1065  parent_jmax(handle) = datapack( 6)
1066  prnt_ks(handle) = datapack( 7)
1067  prnt_ke(handle) = datapack( 8)
1068  prnt_is(handle) = datapack( 9)
1069  prnt_ie(handle) = datapack(10)
1070  prnt_js(handle) = datapack(11)
1071  prnt_je(handle) = datapack(12)
1072  parent_nstep(handle) = datapack(13)
1073  parent_dtsec(handle) = buffer
1074 
1075  ! from daughter to parent
1076  if ( prc_ismaster ) then
1077  call mpi_irecv(datapack, ileng, mpi_integer, prc_myrank, tag+3, intercomm_daughter, ireq1, ierr1)
1078  call mpi_irecv(buffer, 1, mpi_double_precision, prc_myrank, tag+4, intercomm_daughter, ireq2, ierr2)
1079  call mpi_irecv(mp_type_otherside, h_short, mpi_character, prc_myrank, tag+5, intercomm_daughter, ireq3, ierr3)
1080  call mpi_wait(ireq1, istatus, ierr1)
1081  call mpi_wait(ireq2, istatus, ierr2)
1082  call mpi_wait(ireq3, istatus, ierr3)
1083  endif
1084  call comm_bcast(datapack, ileng)
1085  call comm_bcast(buffer)
1086  call comm_bcast(mp_type_otherside)
1087 
1088  daughter_prc_nprocs(handle) = datapack( 1)
1089  daughter_prc_num_x(handle) = datapack( 2)
1090  daughter_prc_num_y(handle) = datapack( 3)
1091  daughter_kmax(handle) = datapack( 4)
1092  daughter_imax(handle) = datapack( 5)
1093  daughter_jmax(handle) = datapack( 6)
1094  datr_ks(handle) = datapack( 7)
1095  datr_ke(handle) = datapack( 8)
1096  datr_is(handle) = datapack( 9)
1097  datr_ie(handle) = datapack(10)
1098  datr_js(handle) = datapack(11)
1099  datr_je(handle) = datapack(12)
1100  daughter_nstep(handle) = datapack(13)
1101  qa_otherside = datapack(14)
1102  daughter_dtsec(handle) = buffer
1103 
1104 
1105  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
1106 
1107  !##### child ####
1108 
1109  ! from parent to daughter
1110  if ( prc_ismaster ) then
1111  call mpi_irecv(datapack, ileng, mpi_integer, prc_myrank, tag, intercomm_parent, ireq1, ierr1)
1112  call mpi_irecv(buffer, 1, mpi_double_precision, prc_myrank, tag+1, intercomm_parent, ireq2, ierr2)
1113  call mpi_irecv(mp_type_otherside, h_short, mpi_character, prc_myrank, tag+2, intercomm_parent, ireq3, ierr3)
1114  call mpi_wait(ireq1, istatus, ierr1)
1115  call mpi_wait(ireq2, istatus, ierr2)
1116  call mpi_wait(ireq3, istatus, ierr3)
1117  endif
1118  call comm_bcast(datapack, ileng)
1119  call comm_bcast(buffer)
1120  call comm_bcast(mp_type_otherside)
1121 
1122  parent_prc_nprocs(handle) = datapack( 1)
1123  parent_prc_num_x(handle) = datapack( 2)
1124  parent_prc_num_y(handle) = datapack( 3)
1125  parent_kmax(handle) = datapack( 4)
1126  parent_imax(handle) = datapack( 5)
1127  parent_jmax(handle) = datapack( 6)
1128  prnt_ks(handle) = datapack( 7)
1129  prnt_ke(handle) = datapack( 8)
1130  prnt_is(handle) = datapack( 9)
1131  prnt_ie(handle) = datapack(10)
1132  prnt_js(handle) = datapack(11)
1133  prnt_je(handle) = datapack(12)
1134  parent_nstep(handle) = datapack(13)
1135  qa_otherside = datapack(14)
1136  parent_dtsec(handle) = buffer
1137 
1138  ! from daughter to parent
1139  datapack( 1) = prc_nprocs
1140  datapack( 2) = prc_num_x
1141  datapack( 3) = prc_num_y
1142  datapack( 4) = kmax
1143  datapack( 5) = imax
1144  datapack( 6) = jmax
1145  datapack( 7) = ks
1146  datapack( 8) = ke
1147  datapack( 9) = is
1148  datapack(10) = ie
1149  datapack(11) = js
1150  datapack(12) = je
1151  datapack(13) = time_nstep
1152  datapack(14) = comm_cartesc_nest_bnd_qa
1153  buffer = time_dtsec
1154 
1155  if ( prc_ismaster ) then
1156  call mpi_isend(datapack, ileng, mpi_integer, prc_myrank, tag+3, intercomm_parent, ireq1, ierr1)
1157  call mpi_isend(buffer, 1, mpi_double_precision, prc_myrank, tag+4, intercomm_parent, ireq2, ierr2)
1158  call mpi_isend(mp_type, h_short, mpi_character, prc_myrank, tag+5, intercomm_parent, ireq3, ierr3)
1159  call mpi_wait(ireq1, istatus, ierr1)
1160  call mpi_wait(ireq2, istatus, ierr2)
1161  call mpi_wait(ireq3, istatus, ierr3)
1162  endif
1163 
1164  daughter_prc_nprocs(handle) = datapack( 1)
1165  daughter_prc_num_x(handle) = datapack( 2)
1166  daughter_prc_num_y(handle) = datapack( 3)
1167  daughter_kmax(handle) = datapack( 4)
1168  daughter_imax(handle) = datapack( 5)
1169  daughter_jmax(handle) = datapack( 6)
1170  datr_ks(handle) = datapack( 7)
1171  datr_ke(handle) = datapack( 8)
1172  datr_is(handle) = datapack( 9)
1173  datr_ie(handle) = datapack(10)
1174  datr_js(handle) = datapack(11)
1175  datr_je(handle) = datapack(12)
1176  daughter_nstep(handle) = datapack(13)
1177  daughter_dtsec(handle) = buffer
1178  else
1179  log_error("COMM_CARTESC_NEST_parentsize",*) '[COMM_CARTESC_NEST_parentsize] internal error'
1180  call prc_abort
1181  endif
1182 
1183  if ( mp_type == mp_type_otherside .and. comm_cartesc_nest_bnd_qa == qa_otherside ) then
1184  online_boundary_diagqhyd = .false.
1185  else
1186  log_info("COMM_CARTESC_NEST_parentsize",*) 'Hydrometeor will be diagnosed'
1187  log_info("COMM_CARTESC_NEST_parentsize",*) 'MP type (remote,local) = ', trim(mp_type_otherside), ", ", trim(mp_type)
1188  log_info("COMM_CARTESC_NEST_parentsize",*) 'Number of QA (remote,local) = ', qa_otherside, comm_cartesc_nest_bnd_qa
1189  comm_cartesc_nest_bnd_qa = n_hyd + 1 ! QV + hydrometeors
1190  online_boundary_diagqhyd = .true.
1191  endif
1192 
1193  return
1194  end subroutine comm_cartesc_nest_parentsize
1195 
1196  !-----------------------------------------------------------------------------
1198  subroutine comm_cartesc_nest_catalogue( &
1199  HANDLE )
1200  use scale_prc, only: &
1201  prc_abort, &
1202  prc_nprocs, &
1203  prc_myrank, &
1204  prc_ismaster
1205  use scale_comm_cartesc, only: &
1206  comm_datatype, &
1207  comm_bcast
1208  use scale_atmos_grid_cartesc_real, only: &
1210  implicit none
1211 
1212  integer, intent(in) :: HANDLE
1213 
1214  integer :: ireq, ierr, ileng
1215  integer :: istatus(mpi_status_size)
1216  integer :: tag
1217  !---------------------------------------------------------------------------
1218 
1219  if( .NOT. use_nesting ) return
1220 
1221  tag = intercomm_id(handle) * 100
1222 
1223  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
1224 
1225  !##### parent ####
1226 
1227  ileng = prc_nprocs * 2 * 2
1228 
1229  if ( prc_ismaster ) then
1230  call mpi_isend(atmos_grid_cartesc_real_domain_catalogue, ileng, comm_datatype, prc_myrank, tag, intercomm_daughter, ireq, ierr)
1231  call mpi_wait(ireq, istatus, ierr)
1232  endif
1233 
1234  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
1235 
1236  !##### child ####
1237 
1238  ileng = parent_prc_nprocs(handle) * 2 * 2
1239 
1240  if ( prc_ismaster ) then
1241  call mpi_irecv(latlon_catalog, ileng, comm_datatype, prc_myrank, tag, intercomm_parent, ireq, ierr)
1242  call mpi_wait(ireq, istatus, ierr)
1243  endif
1244  call comm_bcast( latlon_catalog, parent_prc_nprocs(handle), 2, 2 )
1245 
1246  else
1247  log_error("COMM_CARTESC_NEST_catalogue",*) 'internal error'
1248  call prc_abort
1249  endif
1250 
1251  return
1252  end subroutine comm_cartesc_nest_catalogue
1253 
1254  !-----------------------------------------------------------------------------
1256  subroutine comm_cartesc_nest_ping( &
1257  HANDLE )
1258  use scale_prc, only: &
1259  prc_abort, &
1260  prc_myrank, &
1261  prc_ismaster
1262  use scale_comm_cartesc, only: &
1263  comm_bcast
1264  implicit none
1265 
1266  integer, intent(in) :: HANDLE
1267 
1268  integer :: ping, pong
1269  integer :: ireq1, ireq2, ierr1, ierr2
1270  integer :: istatus(mpi_status_size)
1271  integer :: tag
1272  logical :: ping_error
1273  !---------------------------------------------------------------------------
1274 
1275  if( .NOT. use_nesting ) return
1276 
1277  tag = intercomm_id(handle) * 100
1278  ping_error = .false.
1279 
1280  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
1281 
1282  !##### parent ####
1283 
1284  ping = online_domain_num
1285  pong = 0
1286 
1287  if ( prc_ismaster ) then
1288  call mpi_isend(ping, 1, mpi_integer, prc_myrank, tag+1, intercomm_daughter, ireq1, ierr1)
1289  call mpi_irecv(pong, 1, mpi_integer, prc_myrank, tag+2, intercomm_daughter, ireq2, ierr2)
1290  call mpi_wait(ireq1, istatus, ierr1)
1291  call mpi_wait(ireq2, istatus, ierr2)
1292  endif
1293 
1294  call comm_bcast(pong)
1295 
1296  if ( pong /= intercomm_id(handle)+1 ) ping_error = .true.
1297 
1298  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
1299 
1300  !##### child ####
1301 
1302  ping = online_domain_num
1303  pong = 0
1304 
1305  if ( prc_ismaster ) then
1306  call mpi_isend(ping, 1, mpi_integer, prc_myrank, tag+2, intercomm_parent, ireq1, ierr1)
1307  call mpi_irecv(pong, 1, mpi_integer, prc_myrank, tag+1, intercomm_parent, ireq2, ierr2)
1308  call mpi_wait(ireq1, istatus, ierr1)
1309  call mpi_wait(ireq2, istatus, ierr2)
1310  endif
1311 
1312  call comm_bcast(pong)
1313 
1314  if ( pong /= intercomm_id(handle) ) ping_error = .true.
1315 
1316  else
1317  log_error("COMM_CARTESC_NEST_ping",*) 'internal error'
1318  call prc_abort
1319  endif
1320 
1321  if ( ping_error ) then
1322  log_error("COMM_CARTESC_NEST_ping",*) 'ping destination error'
1323  call prc_abort
1324  endif
1325 
1326  return
1327  end subroutine comm_cartesc_nest_ping
1328 
1329  !-----------------------------------------------------------------------------
1331  subroutine comm_cartesc_nest_setup_nestdown( &
1332  HANDLE )
1333  use scale_prc, only: &
1334  prc_abort, &
1335  prc_myrank, &
1336  prc_ismaster
1337  use scale_comm_cartesc, only: &
1338  comm_world, &
1339  comm_bcast
1340  implicit none
1341 
1342  integer, intent(in) :: HANDLE
1343 
1344  integer, allocatable :: buffer_LIST (:)
1345  integer, allocatable :: buffer_ALLLIST(:)
1346 
1347  integer :: ireq, ierr, ileng
1348  integer :: istatus(mpi_status_size)
1349  integer :: tag, target_rank
1350 
1351  integer :: i, j, k
1352  !---------------------------------------------------------------------------
1353 
1354  if( .NOT. use_nesting ) return
1355 
1356  tag = intercomm_id(handle) * 100
1357 
1358  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
1359 
1360  !##### parent ####
1361 
1362  if ( prc_ismaster ) then
1363  call mpi_irecv(comm_cartesc_nest_tile_allmax_p, 1, mpi_integer, prc_myrank, tag+1, intercomm_daughter, ireq, ierr)
1364  call mpi_wait(ireq, istatus, ierr)
1365  endif
1366  call comm_bcast(comm_cartesc_nest_tile_allmax_p)
1367 
1368  allocate( comm_cartesc_nest_tile_list_p(comm_cartesc_nest_tile_allmax_p,daughter_prc_nprocs(handle)) )
1369  allocate( comm_cartesc_nest_tile_list_yp(comm_cartesc_nest_tile_allmax_p*daughter_prc_nprocs(handle)) )
1370 
1371  ileng = comm_cartesc_nest_tile_allmax_p*daughter_prc_nprocs(handle)
1372  if ( prc_ismaster ) then
1373  call mpi_irecv(comm_cartesc_nest_tile_list_p, ileng, mpi_integer, prc_myrank, tag+2, intercomm_daughter, ireq, ierr)
1374  call mpi_wait(ireq, istatus, ierr)
1375  endif
1376  call comm_bcast(comm_cartesc_nest_tile_list_p, comm_cartesc_nest_tile_allmax_p, daughter_prc_nprocs(handle))
1377 
1378  comm_cartesc_nest_tile_list_yp(:) = -1
1379 
1380  k = 0
1381  do j = 1, daughter_prc_nprocs(handle)
1382  do i = 1, comm_cartesc_nest_tile_allmax_p
1383  if ( comm_cartesc_nest_tile_list_p(i,j) == prc_myrank ) then
1384  k = k + 1
1385  comm_cartesc_nest_tile_list_yp(k) = j - 1 !rank number is started from 1
1386  endif
1387  enddo
1388  enddo
1389  num_yp = k
1390 
1391  log_info("COMM_CARTESC_NEST_setup_nestdown",'(A,I5,A,I5)') "[P] Num YP =",num_yp," Num TILE(MAX) =",comm_cartesc_nest_tile_allmax_p
1392 
1393  if ( prc_ismaster ) then
1394  call mpi_irecv(online_daughter_use_velz, 1, mpi_logical, prc_myrank, tag+3, intercomm_daughter, ireq, ierr)
1395  call mpi_wait(ireq, istatus, ierr)
1396  endif
1397  call comm_bcast(online_daughter_use_velz)
1398 
1399  log_info("COMM_CARTESC_NEST_setup_nestdown",'(1x,A,L2)') 'NEST: ONLINE_DAUGHTER_USE_VELZ =', online_daughter_use_velz
1400 
1401  if ( prc_ismaster ) then
1402  call mpi_irecv(online_daughter_no_rotate, 1, mpi_logical, prc_myrank, tag+4, intercomm_daughter, ireq, ierr)
1403  call mpi_wait(ireq, istatus, ierr)
1404  endif
1405  call comm_bcast(online_daughter_no_rotate)
1406 
1407  if( online_no_rotate .neqv. online_daughter_no_rotate ) then
1408  log_error("COMM_CARTESC_NEST_setup_nestdown",*) 'Flag of NO_ROTATE is not consistent with the child domain'
1409  log_error_cont(*) 'ONLINE_NO_ROTATE = ', online_no_rotate
1410  log_error_cont(*) 'ONLINE_DAUGHTER_NO_ROTATE =', online_daughter_no_rotate
1411  call prc_abort
1412  endif
1413  log_info("COMM_CARTESC_NEST_setup_nestdown",'(1x,A,L2)') 'NEST: ONLINE_DAUGHTER_NO_ROTATE =', online_daughter_no_rotate
1414 
1415  call comm_cartesc_nest_importgrid_nestdown( handle )
1416 
1417  do i = 1, num_yp
1418  target_rank = comm_cartesc_nest_tile_list_yp(i)
1419  call mpi_isend(i, 1, mpi_integer, target_rank, tag+5, intercomm_daughter, ireq, ierr)
1420  call mpi_wait(ireq, istatus, ierr)
1421  enddo
1422 
1423  call mpi_barrier(intercomm_daughter, ierr)
1424 
1425  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
1426 
1427  !##### child ####
1428 
1429  comm_cartesc_nest_tile_all = size( comm_cartesc_nest_tile_id(:) ) ! should be equal to "NEST_TILE_NUM_X*NEST_TILE_NUM_Y"
1430  call mpi_allreduce( comm_cartesc_nest_tile_all, &
1431  comm_cartesc_nest_tile_allmax_d, &
1432  1, &
1433  mpi_integer, &
1434  mpi_max, &
1435  comm_world, &
1436  ierr )
1437  log_info("COMM_CARTESC_NEST_setup_nestdown",'(A,I5,A,I5)') "[D] Num YP =",comm_cartesc_nest_tile_all," Num TILE(MAX) =",comm_cartesc_nest_tile_allmax_d
1438 
1439  if ( prc_ismaster ) then
1440  call mpi_isend(comm_cartesc_nest_tile_allmax_d, 1, mpi_integer, prc_myrank, tag+1, intercomm_parent, ireq, ierr)
1441  call mpi_wait(ireq, istatus, ierr)
1442  endif
1443 
1444  allocate( buffer_list(comm_cartesc_nest_tile_allmax_d) )
1445  allocate( buffer_alllist(comm_cartesc_nest_tile_allmax_d*daughter_prc_nprocs(handle)) )
1446  allocate( comm_cartesc_nest_tile_list_d(comm_cartesc_nest_tile_allmax_d,daughter_prc_nprocs(handle)) )
1447 
1448  do i = 1, comm_cartesc_nest_tile_allmax_d
1449  if ( i <= comm_cartesc_nest_tile_all ) then
1450  buffer_list(i) = comm_cartesc_nest_tile_id(i)
1451  else
1452  buffer_list(i) = -1
1453  endif
1454  enddo
1455 
1456  ileng = comm_cartesc_nest_tile_allmax_d
1457  call mpi_allgather( buffer_list(:), &
1458  ileng, &
1459  mpi_integer, &
1460  buffer_alllist(:), &
1461  ileng, &
1462  mpi_integer, &
1463  comm_world, &
1464  ierr )
1465  k = 1
1466  do j = 1, daughter_prc_nprocs(handle)
1467  do i = 1, comm_cartesc_nest_tile_allmax_d
1468  comm_cartesc_nest_tile_list_d(i,j) = buffer_alllist(k)
1469  k = k + 1
1470  enddo
1471  enddo
1472 
1473  deallocate( buffer_list )
1474  deallocate( buffer_alllist )
1475 
1476  ileng = comm_cartesc_nest_tile_allmax_d*daughter_prc_nprocs(handle)
1477  if ( prc_ismaster ) then
1478  call mpi_isend(comm_cartesc_nest_tile_list_d, ileng, mpi_integer, prc_myrank, tag+2, intercomm_parent, ireq, ierr)
1479  call mpi_wait(ireq, istatus, ierr)
1480  endif
1481 
1482  if ( prc_ismaster ) then
1483  call mpi_isend(online_use_velz, 1, mpi_logical, prc_myrank, tag+3, intercomm_parent, ireq, ierr)
1484  call mpi_wait(ireq, istatus, ierr)
1485  endif
1486 
1487  if ( prc_ismaster ) then
1488  call mpi_isend(online_no_rotate, 1, mpi_logical, prc_myrank, tag+4, intercomm_parent, ireq, ierr)
1489  call mpi_wait(ireq, istatus, ierr)
1490  endif
1491  call comm_bcast(online_daughter_no_rotate)
1492 
1493  call comm_cartesc_nest_importgrid_nestdown( handle )
1494 
1495  do i = 1, comm_cartesc_nest_tile_all
1496  target_rank = comm_cartesc_nest_tile_list_d(i,prc_myrank+1)
1497  call mpi_irecv( call_order(i), 1, mpi_integer, target_rank, tag+5, intercomm_parent, ireq, ierr )
1498  call mpi_wait(ireq, istatus, ierr)
1499  enddo
1500 
1501  call mpi_barrier(intercomm_parent, ierr)
1502  else
1503  log_error("COMM_CARTESC_NEST_setup_nestdown",*) 'internal error'
1504  call prc_abort
1505  endif
1506 
1507  if( num_yp * 16 > max_rq .OR. comm_cartesc_nest_tile_all * 16 > max_rq ) then ! 16 = dyn:5 + qtrc:11
1508  log_error("COMM_CARTESC_NEST_setup_nestdown",*) 'internal error (overflow number of ireq)'
1509  log_error_cont(*) 'NUM_YP x 16 = ', num_yp * 16
1510  log_error_cont(*) 'COMM_CARTESC_NEST_TILE_ALL x 16 = ', comm_cartesc_nest_tile_all * 16
1511  log_error_cont(*) 'max_rq = ', max_rq
1512  call prc_abort
1513  endif
1514 
1515  return
1516  end subroutine comm_cartesc_nest_setup_nestdown
1517 
1518  !-----------------------------------------------------------------------------
1520  subroutine comm_cartesc_nest_importgrid_nestdown( &
1521  HANDLE )
1522  use scale_prc, only: &
1523  prc_myrank, &
1524  prc_abort
1525  use scale_atmos_grid_cartesc_real, only: &
1534  use scale_comm_cartesc, only: &
1536  implicit none
1537 
1538  integer, intent(in) :: HANDLE
1539 
1540  integer :: ierr, ileng
1541  integer :: istatus(mpi_status_size)
1542  integer :: tag, tagbase, target_rank
1543 
1544  integer :: xloc, yloc
1545  integer :: xs, xe
1546  integer :: ys, ye
1547 
1548  real(RP) :: max_ref, max_loc
1549 
1550  integer :: i, k, rq
1551  !---------------------------------------------------------------------------
1552 
1553  if( .NOT. use_nesting ) return
1554 
1555  tagbase = intercomm_id(handle) * 100
1556  rq = 0
1557 
1558  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
1559 
1560  !##### parent [send issue] #####
1561 
1562  do i = 1, num_yp
1563  ! send data to multiple daughter processes
1564  target_rank = comm_cartesc_nest_tile_list_yp(i)
1565 
1566  rq = rq + 1
1567  ileng = parent_ia(handle) * parent_ja(handle)
1568  tag = tagbase + tag_lon
1569  call mpi_isend(atmos_grid_cartesc_real_lon, ileng, comm_datatype, target_rank, tag, intercomm_daughter, ireq_p(rq), ierr)
1570  call mpi_wait(ireq_p(rq), istatus, ierr)
1571 
1572  rq = rq + 1
1573  ileng = parent_ia(handle) * parent_ja(handle)
1574  tag = tagbase + tag_lat
1575  call mpi_isend(atmos_grid_cartesc_real_lat, ileng, comm_datatype, target_rank, tag, intercomm_daughter, ireq_p(rq), ierr)
1576  call mpi_wait(ireq_p(rq), istatus, ierr)
1577 
1578  rq = rq + 1
1579  ileng = parent_ia(handle) * parent_ja(handle)
1580  tag = tagbase + tag_lonuy
1581  call mpi_isend(atmos_grid_cartesc_real_lonuy(1:ia,1:ja), ileng, comm_datatype, target_rank, tag, intercomm_daughter, ireq_p(rq), ierr)
1582  call mpi_wait(ireq_p(rq), istatus, ierr)
1583 
1584  rq = rq + 1
1585  ileng = parent_ia(handle) * parent_ja(handle)
1586  tag = tagbase + tag_latuy
1587  call mpi_isend(atmos_grid_cartesc_real_latuy(1:ia,1:ja), ileng, comm_datatype, target_rank, tag, intercomm_daughter, ireq_p(rq), ierr)
1588  call mpi_wait(ireq_p(rq), istatus, ierr)
1589 
1590  rq = rq + 1
1591  ileng = parent_ia(handle) * parent_ja(handle)
1592  tag = tagbase + tag_lonxv
1593  call mpi_isend(atmos_grid_cartesc_real_lonxv(1:ia,1:ja), ileng, comm_datatype, target_rank, tag, intercomm_daughter, ireq_p(rq), ierr)
1594  call mpi_wait(ireq_p(rq), istatus, ierr)
1595 
1596  rq = rq + 1
1597  ileng = parent_ia(handle) * parent_ja(handle)
1598  tag = tagbase + tag_latxv
1599  call mpi_isend(atmos_grid_cartesc_real_latxv(1:ia,1:ja), ileng, comm_datatype, target_rank, tag, intercomm_daughter, ireq_p(rq), ierr)
1600  call mpi_wait(ireq_p(rq), istatus, ierr)
1601 
1602  rq = rq + 1
1603  ileng = parent_ka(handle) * parent_ia(handle) * parent_ja(handle)
1604  tag = tagbase + tag_cz
1605  call mpi_isend(atmos_grid_cartesc_real_cz, ileng, comm_datatype, target_rank, tag, intercomm_daughter, ireq_p(rq), ierr)
1606  call mpi_wait(ireq_p(rq), istatus, ierr)
1607 
1608  rq = rq + 1
1609  ileng = (parent_ka(handle)+1) * parent_ia(handle) * parent_ja(handle)
1610  tag = tagbase + tag_fz
1611  call mpi_isend(atmos_grid_cartesc_real_fz, ileng, comm_datatype, target_rank, tag, intercomm_daughter, ireq_p(rq), ierr)
1612  call mpi_wait(ireq_p(rq), istatus, ierr)
1613  enddo
1614 
1615  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
1616 
1617  !##### child [recv & wait issue] #####
1618 
1619  do i = 1, comm_cartesc_nest_tile_all
1620  ! receive data from multiple parent tiles
1621  target_rank = comm_cartesc_nest_tile_list_d(i,prc_myrank+1)
1622 
1623  xloc = mod( i-1, comm_cartesc_nest_tile_num_x ) + 1
1624  yloc = int( real(i-1) / real(COMM_CARTESC_NEST_TILE_NUM_X) ) + 1
1625 
1626  xs = parent_imax(handle) * (xloc-1) + 1
1627  xe = parent_imax(handle) * xloc
1628  ys = parent_jmax(handle) * (yloc-1) + 1
1629  ye = parent_jmax(handle) * yloc
1630 
1631  rq = rq + 1
1632  ileng = parent_ia(handle) * parent_ja(handle)
1633  tag = tagbase + tag_lon
1634  call mpi_irecv(buffer_2d, ileng, comm_datatype, target_rank, tag, intercomm_parent, ireq_d(rq), ierr)
1635  call mpi_wait(ireq_d(rq), istatus, ierr)
1636  buffer_ref_lon(xs:xe,ys:ye) = buffer_2d(prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle))
1637 
1638  rq = rq + 1
1639  ileng = parent_ia(handle) * parent_ja(handle)
1640  tag = tagbase + tag_lat
1641  call mpi_irecv(buffer_2d, ileng, comm_datatype, target_rank, tag, intercomm_parent, ireq_d(rq), ierr)
1642  call mpi_wait(ireq_d(rq), istatus, ierr)
1643  buffer_ref_lat(xs:xe,ys:ye) = buffer_2d(prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle))
1644 
1645  rq = rq + 1
1646  ileng = parent_ia(handle) * parent_ja(handle)
1647  tag = tagbase + tag_lonuy
1648  call mpi_irecv(buffer_2d, ileng, comm_datatype, target_rank, tag, intercomm_parent, ireq_d(rq), ierr)
1649  call mpi_wait(ireq_d(rq), istatus, ierr)
1650  buffer_ref_lonuy(xs:xe,ys:ye) = buffer_2d(prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle))
1651 
1652  rq = rq + 1
1653  ileng = parent_ia(handle) * parent_ja(handle)
1654  tag = tagbase + tag_latuy
1655  call mpi_irecv(buffer_2d, ileng, comm_datatype, target_rank, tag, intercomm_parent, ireq_d(rq), ierr)
1656  call mpi_wait(ireq_d(rq), istatus, ierr)
1657  buffer_ref_latuy(xs:xe,ys:ye) = buffer_2d(prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle))
1658 
1659  rq = rq + 1
1660  ileng = parent_ia(handle) * parent_ja(handle)
1661  tag = tagbase + tag_lonxv
1662  call mpi_irecv(buffer_2d, ileng, comm_datatype, target_rank, tag, intercomm_parent, ireq_d(rq), ierr)
1663  call mpi_wait(ireq_d(rq), istatus, ierr)
1664  buffer_ref_lonxv(xs:xe,ys:ye) = buffer_2d(prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle))
1665 
1666  rq = rq + 1
1667  ileng = parent_ia(handle) * parent_ja(handle)
1668  tag = tagbase + tag_latxv
1669  call mpi_irecv(buffer_2d, ileng, comm_datatype, target_rank, tag, intercomm_parent, ireq_d(rq), ierr)
1670  call mpi_wait(ireq_d(rq), istatus, ierr)
1671  buffer_ref_latxv(xs:xe,ys:ye) = buffer_2d(prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle))
1672 
1673  rq = rq + 1
1674  ileng = parent_ka(handle) * parent_ia(handle) * parent_ja(handle)
1675  tag = tagbase + tag_cz
1676  call mpi_irecv(buffer_3d, ileng, comm_datatype, target_rank, tag, intercomm_parent, ireq_d(rq), ierr)
1677  call mpi_wait(ireq_d(rq), istatus, ierr)
1678  do k = 1, parent_ka(handle)
1679  buffer_ref_cz(k,xs:xe,ys:ye) = buffer_3d(k,prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle))
1680  enddo
1681 
1682  rq = rq + 1
1683  ileng = (parent_ka(handle)+1) * parent_ia(handle) * parent_ja(handle)
1684  tag = tagbase + tag_fz
1685  call mpi_irecv(buffer_3df,ileng, comm_datatype, target_rank, tag, intercomm_parent, ireq_d(rq), ierr)
1686  call mpi_wait(ireq_d(rq), istatus, ierr)
1687  do k = 0, parent_ka(handle)
1688  buffer_ref_fz(k,xs:xe,ys:ye) = buffer_3df(k,prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle))
1689  enddo
1690  enddo
1691 
1692  ! check domain compatibility
1693  max_ref = maxval( buffer_ref_fz(:,:,:) )
1694  max_loc = maxval( atmos_grid_cartesc_real_fz(ks-1:ke,:,:) ) ! HALO + 1
1695  if ( max_ref < max_loc ) then
1696  log_error("COMM_CARTESC_NEST_importgrid_nestdown",*) 'REQUESTED DOMAIN IS TOO MUCH BROAD'
1697  log_error_cont(*) '-- VERTICAL direction over the limit'
1698  log_error_cont(*) '-- reference max: ', max_ref
1699  log_error_cont(*) '-- local max: ', max_loc
1700  call prc_abort
1701  endif
1702 
1703  else
1704  log_error("COMM_CARTESC_NEST_importgrid_nestdown",*) 'internal error'
1705  call prc_abort
1706  endif
1707 
1708  return
1709  end subroutine comm_cartesc_nest_importgrid_nestdown
1710 
1711  !-----------------------------------------------------------------------------
1713  subroutine comm_cartesc_nest_nestdown( &
1714  HANDLE, &
1715  BND_QA, &
1716  DENS_send, &
1717  MOMZ_send, &
1718  MOMX_send, &
1719  MOMY_send, &
1720  RHOT_send, &
1721  QTRC_send, &
1722  DENS_recv, &
1723  VELZ_recv, &
1724  VELX_recv, &
1725  VELY_recv, &
1726  POTT_recv, &
1727  QTRC_recv )
1728  use scale_prc, only: &
1729  prc_abort
1730  use scale_comm_cartesc, only: &
1731  comm_vars8, &
1732  comm_wait
1733  use scale_atmos_grid_cartesc_metric, only: &
1735  implicit none
1736 
1737  integer, intent(in) :: HANDLE
1738  integer, intent(in) :: BND_QA
1739  real(RP), intent(in) :: DENS_send(parent_ka (handle),parent_ia (handle),parent_ja (handle))
1740  real(RP), intent(in) :: MOMZ_send(parent_ka (handle),parent_ia (handle),parent_ja (handle))
1741  real(RP), intent(in) :: MOMX_send(parent_ka (handle),parent_ia (handle),parent_ja (handle))
1742  real(RP), intent(in) :: MOMY_send(parent_ka (handle),parent_ia (handle),parent_ja (handle))
1743  real(RP), intent(in) :: RHOT_send(parent_ka (handle),parent_ia (handle),parent_ja (handle))
1744  real(RP), intent(in) :: QTRC_send(parent_ka (handle),parent_ia (handle),parent_ja (handle),bnd_qa)
1745  real(RP), intent(out) :: DENS_recv(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1746  real(RP), intent(out) :: VELZ_recv(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1747  real(RP), intent(out) :: VELX_recv(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1748  real(RP), intent(out) :: VELY_recv(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1749  real(RP), intent(out) :: POTT_recv(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1750  real(RP), intent(out) :: QTRC_recv(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle),bnd_qa)
1751 
1752  real(RP) :: WORK1_send(parent_ka (handle),parent_ia (handle),parent_ja (handle))
1753  real(RP) :: WORK2_send(parent_ka (handle),parent_ia (handle),parent_ja (handle))
1754  real(RP) :: WORK1_recv(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1755  real(RP) :: WORK2_recv(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1756  real(RP) :: U_ll_recv (daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1757  real(RP) :: V_ll_recv (daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1758  real(RP) :: u_on_map, v_on_map
1759 
1760  real(RP) :: dummy(1,1,1)
1761  integer :: tagbase, tagcomm
1762  integer :: isu_tag, isu_tagf
1763 
1764  integer :: ierr
1765  integer :: i, j, k, iq
1766  !---------------------------------------------------------------------------
1767 
1768  if( .NOT. use_nesting ) return
1769 
1770  if ( bnd_qa > i_bndqa ) then
1771  log_error("COMM_CARTESC_NEST_nestdown",*) 'internal error: BND_QA is larger than I_BNDQA'
1772  call prc_abort
1773  endif
1774 
1775  tagcomm = intercomm_id(handle) * order_tag_comm
1776 
1777  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
1778 
1779  !##### parent [send issue] #####
1780 
1781  call prof_rapstart('NEST_total_P', 2)
1782  call prof_rapstart('NEST_pack_P', 2)
1783 
1784  nsend = nsend + 1
1785  log_info("COMM_CARTESC_NEST_nestdown",'(1X,A,I5,A)') "CONeP[P] send( ", nsend, " )"
1786 
1787  ! to keep values at that time by finish of sending process
1788 !OCL XFILL
1789  org_dens(:,:,:) = dens_send(:,:,:)
1790 !OCL XFILL
1791  org_momz(:,:,:) = momz_send(:,:,:)
1792 !OCL XFILL
1793  org_momx(:,:,:) = momx_send(:,:,:)
1794 !OCL XFILL
1795  org_momy(:,:,:) = momy_send(:,:,:)
1796 !OCL XFILL
1797  org_rhot(:,:,:) = rhot_send(:,:,:)
1798  do iq = 1, bnd_qa
1799 !OCL XFILL
1800  org_qtrc(:,:,:,iq) = qtrc_send(:,:,:,iq)
1801  enddo
1802 
1803  !*** request control
1804  !--- do not change the calling order below;
1805  !--- it should be consistent with the order in "COMM_CARTESC_NEST_recvwait_issue"
1806  rq_ctl_p = 0
1807 
1808  if ( .NOT. online_daughter_no_rotate ) then
1809  ! from staggered point to scalar point
1810  do j = 1, parent_ja(handle)
1811  do i = 2, parent_ia(handle)
1812  do k = 1, parent_ka(handle)
1813  work1_send(k,i,j) = ( org_momx(k,i-1,j) + org_momx(k,i,j) ) * 0.5_rp
1814  enddo
1815  enddo
1816  enddo
1817 
1818  do j = 1, parent_ja(handle)
1819  do k = 1, parent_ka(handle)
1820  work1_send(k,1,j) = org_momx(k,1,j)
1821  enddo
1822  enddo
1823 
1824  call comm_vars8( work1_send(:,:,:), 1 )
1825 
1826  do j = 2, parent_ja(handle)
1827  do i = 1, parent_ia(handle)
1828  do k = 1, parent_ka(handle)
1829  work2_send(k,i,j) = ( org_momy(k,i,j-1) + org_momy(k,i,j) ) * 0.5_rp
1830  enddo
1831  enddo
1832  enddo
1833 
1834  do i = 1, parent_ia(handle)
1835  do k = 1, parent_ka(handle)
1836  work2_send(k,i,1) = org_momy(k,i,1)
1837  enddo
1838  enddo
1839 
1840  call comm_vars8( work2_send(:,:,:), 2 )
1841 
1842  call comm_wait ( work1_send(:,:,:), 1, .false. )
1843  call comm_wait ( work2_send(:,:,:), 2, .false. )
1844 
1845  ! rotation from map-projected field to latlon field
1846  do j = 1, parent_ja(handle)
1847  do i = 1, parent_ia(handle)
1848  do k = 1, parent_ka(handle)
1849  u_on_map = work1_send(k,i,j) / org_dens(k,i,j)
1850  v_on_map = work2_send(k,i,j) / org_dens(k,i,j)
1851 
1852  org_u_ll(k,i,j) = u_on_map * rotc(i,j,1) - v_on_map * rotc(i,j,2)
1853  org_v_ll(k,i,j) = u_on_map * rotc(i,j,2) + v_on_map * rotc(i,j,1)
1854  enddo
1855  enddo
1856  enddo
1857  endif
1858 
1859  tagbase = tagcomm + tag_dens*order_tag_var
1860  call comm_cartesc_nest_intercomm_nestdown( org_dens(:,:,:), & ! [IN]
1861  dummy(:,:,:), & ! [OUT]
1862  tagbase, i_sclr, handle, & ! [IN]
1863  isu_tag, isu_tagf, & ! [INOUT]
1864  flag_dens = .true. ) ! [IN]
1865 
1866  tagbase = tagcomm + tag_momz*order_tag_var
1867  if ( online_daughter_use_velz ) then
1868  call comm_cartesc_nest_intercomm_nestdown( org_momz(:,:,:), & ! [IN]
1869  dummy(:,:,:), & ! [OUT]
1870  tagbase, i_zstg, handle, & ! [IN]
1871  isu_tag, isu_tagf ) ! [INOUT]
1872  endif
1873 
1874  tagbase = tagcomm + tag_momx*order_tag_var
1875  if ( online_daughter_no_rotate ) then
1876  call comm_cartesc_nest_intercomm_nestdown( org_momx(:,:,:), & ! [IN]
1877  dummy(:,:,:), & ! [OUT]
1878  tagbase, i_xstg, handle, & ! [IN]
1879  isu_tag, isu_tagf ) ! [INOUT]
1880  else
1881  call comm_cartesc_nest_intercomm_nestdown( org_u_ll(:,:,:), & ! [IN]
1882  dummy(:,:,:), & ! [OUT]
1883  tagbase, i_sclr, handle, & ! [IN]
1884  isu_tag, isu_tagf ) ! [INOUT]
1885  endif
1886 
1887  tagbase = tagcomm + tag_momy*order_tag_var
1888  if ( online_daughter_no_rotate ) then
1889  call comm_cartesc_nest_intercomm_nestdown( org_momy(:,:,:), & ! [IN]
1890  dummy(:,:,:), & ! [OUT]
1891  tagbase, i_ystg, handle, & ! [IN]
1892  isu_tag, isu_tagf ) ! [INOUT]
1893  else
1894  call comm_cartesc_nest_intercomm_nestdown( org_v_ll(:,:,:), & ! [IN]
1895  dummy(:,:,:), & ! [OUT]
1896  tagbase, i_sclr, handle, & ! [IN]
1897  isu_tag, isu_tagf ) ! [INOUT]
1898  endif
1899 
1900  tagbase = tagcomm + tag_rhot*order_tag_var
1901  call comm_cartesc_nest_intercomm_nestdown( org_rhot(:,:,:), & ! [IN]
1902  dummy(:,:,:), & ! [OUT]
1903  tagbase, i_sclr, handle, & ! [IN]
1904  isu_tag, isu_tagf ) ! [INOUT]
1905 
1906  do iq = 1, bnd_qa
1907  tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
1908  call comm_cartesc_nest_intercomm_nestdown( org_qtrc(:,:,:,iq), & ! [IN]
1909  dummy(:,:,:), & ! [OUT]
1910  tagbase, i_sclr, handle, & ! [IN]
1911  isu_tag, isu_tagf ) ! [INOUT]
1912  enddo
1913 
1914  rq_tot_p = rq_ctl_p
1915 
1916  call prof_rapend ('NEST_pack_P', 2)
1917  call prof_rapend ('NEST_total_P', 2)
1918 
1919  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
1920 
1921  !##### child [wait issue] #####
1922 
1923  call prof_rapstart('NEST_total_C', 2)
1924  call prof_rapstart('NEST_wait_C', 2)
1925 
1926  nwait_d = nwait_d + 1
1927  !LOG_INFO("COMM_CARTESC_NEST_nestdown",'(1X,A,I5,A)') "NestIDC [C]: que wait ( ", nwait_d, " )"
1928 
1929  !*** reset issue tag and request control
1930  !--- do not change the calling order below;
1931  !--- it should be consistent with the order in "COMM_CARTESC_NEST_recvwait_issue"
1932  isu_tag = 0
1933  isu_tagf = 0
1934 
1935  call comm_cartesc_nest_waitall( rq_tot_d, ireq_d )
1936 
1937  if ( online_aggressive_comm ) then
1938  ! nothing to do
1939  else
1940  call mpi_barrier(intercomm_parent, ierr)
1941  endif
1942 
1943  call prof_rapend ('NEST_wait_C', 2)
1944  call prof_rapstart('NEST_unpack_C', 2)
1945 
1946  tagbase = tagcomm + tag_dens*order_tag_var
1947  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
1948  work1_recv(:,:,:), & ! [OUT]
1949  tagbase, i_sclr, handle, & ! [IN]
1950  isu_tag, isu_tagf, & ! [INOUT]
1951  flag_dens = .true. ) ! [IN]
1952 !OCL XFILL
1953  do j = 1, daughter_ja(handle)
1954  do i = 1, daughter_ia(handle)
1955  do k = datr_ks(handle), datr_ke(handle)
1956  dens_recv(k,i,j) = work1_recv(k,i,j)
1957  enddo
1958  enddo
1959  enddo
1960 
1961  call comm_vars8( dens_recv, 1 )
1962 
1963  tagbase = tagcomm + tag_momz*order_tag_var
1964  if ( online_use_velz ) then
1965  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
1966  work2_recv(:,:,:), & ! [OUT]
1967  tagbase, i_zstg, handle, & ! [IN]
1968  isu_tag, isu_tagf ) ! [INOUT]
1969 !OCL XFILL
1970  do j = 1, daughter_ja(handle)
1971  do i = 1, daughter_ia(handle)
1972  do k = datr_ks(handle), datr_ke(handle)-1
1973  velz_recv(k,i,j) = work2_recv(k,i,j) / ( work1_recv(k,i,j) + work1_recv(k+1,i,j) ) * 2.0_rp
1974  enddo
1975  enddo
1976  enddo
1977 
1978  do j = 1, daughter_ja(handle)
1979  do i = 1, daughter_ia(handle)
1980  velz_recv(datr_ks(handle)-1,i,j) = 0.0_rp
1981  velz_recv(datr_ke(handle) ,i,j) = 0.0_rp
1982  enddo
1983  enddo
1984  endif
1985 
1986  call comm_wait ( dens_recv, 1, .false. )
1987 
1988  tagbase = tagcomm + tag_momx*order_tag_var
1989  if ( online_no_rotate ) then
1990  ! U_ll_recv receives MOMX
1991  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
1992  work1_recv(:,:,:), & ! [OUT]
1993  tagbase, i_xstg, handle, & ! [IN]
1994  isu_tag, isu_tagf ) ! [INOUT]
1995  else
1996  ! U_ll_recv receives MOMX/DENS
1997  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
1998  u_ll_recv(:,:,:), & ! [OUT]
1999  tagbase, i_sclr, handle, & ! [IN]
2000  isu_tag, isu_tagf ) ! [INOUT]
2001  endif
2002 
2003  tagbase = tagcomm + tag_momy*order_tag_var
2004  if ( online_no_rotate ) then
2005  ! V_ll_recv receives MOMY
2006  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2007  work2_recv(:,:,:), & ! [OUT]
2008  tagbase, i_ystg, handle, & ! [IN]
2009  isu_tag, isu_tagf ) ! [INOUT]
2010  else
2011  ! V_ll_recv receives MOMY/DENS
2012  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2013  v_ll_recv(:,:,:), & ! [OUT]
2014  tagbase, i_sclr, handle, & ! [IN]
2015  isu_tag, isu_tagf ) ! [INOUT]
2016  endif
2017 
2018  if ( online_no_rotate ) then
2019 
2020 !OCL XFILL
2021  do j = 1, daughter_ja(handle)
2022  do i = 1, daughter_ia(handle)-1
2023  do k = datr_ks(handle), datr_ke(handle)
2024  velx_recv(k,i,j) = work1_recv(k,i,j) / ( dens_recv(k,i+1,j) + dens_recv(k,i,j) ) * 2.0_rp
2025  enddo
2026  enddo
2027  enddo
2028 
2029  i = daughter_ia(handle)
2030 !OCL XFILL
2031  do j = 1, daughter_ja(handle)
2032  do k = datr_ks(handle), datr_ke(handle)
2033  velx_recv(k,i,j) = work1_recv(k,i,j) / dens_recv(k,i,j)
2034  enddo
2035  enddo
2036 
2037  call comm_vars8( velx_recv, 2 )
2038 
2039 !OCL XFILL
2040  do j = 1, daughter_ja(handle)-1
2041  do i = 1, daughter_ia(handle)
2042  do k = datr_ks(handle), datr_ke(handle)
2043  vely_recv(k,i,j) = work2_recv(k,i,j) / ( dens_recv(k,i,j+1) + dens_recv(k,i,j) ) * 2.0_rp
2044  enddo
2045  enddo
2046  enddo
2047 
2048  j = daughter_ja(handle)
2049 !OCL XFILL
2050  do i = 1, daughter_ia(handle)
2051  do k = datr_ks(handle), datr_ke(handle)
2052  vely_recv(k,i,j) = work2_recv(k,i,j) / dens_recv(k,i,j)
2053  enddo
2054  enddo
2055 
2056  call comm_vars8( vely_recv, 3 )
2057 
2058  call comm_wait ( velx_recv, 2, .false. )
2059  call comm_wait ( vely_recv, 3, .false. )
2060 
2061  else ! rotate
2062 
2063  ! rotation from latlon field to map-projected field
2064 !OCL XFILL
2065  do j = 1, daughter_ja(handle)
2066  do i = 1, daughter_ia(handle)
2067  do k = datr_ks(handle), datr_ke(handle)
2068  work1_recv(k,i,j) = u_ll_recv(k,i,j) * rotc(i,j,1) + v_ll_recv(k,i,j) * rotc(i,j,2)
2069  work2_recv(k,i,j) = -u_ll_recv(k,i,j) * rotc(i,j,2) + v_ll_recv(k,i,j) * rotc(i,j,1)
2070  enddo
2071  enddo
2072  enddo
2073 
2074  ! from scalar point to staggered point
2075 !OCL XFILL
2076  do j = 1, daughter_ja(handle)
2077  do i = 1, daughter_ia(handle)-1
2078  do k = datr_ks(handle), datr_ke(handle)
2079  velx_recv(k,i,j) = ( work1_recv(k,i+1,j) + work1_recv(k,i,j) ) * 0.5_rp
2080  enddo
2081  enddo
2082  enddo
2083 
2084  i = daughter_ia(handle)
2085 !OCL XFILL
2086  do j = 1, daughter_ja(handle)
2087  do k = datr_ks(handle), datr_ke(handle)
2088  velx_recv(k,i,j) = work1_recv(k,i,j)
2089  enddo
2090  enddo
2091 
2092  call comm_vars8( velx_recv, 2 )
2093 
2094 !OCL XFILL
2095  do j = 1, daughter_ja(handle)-1
2096  do i = 1, daughter_ia(handle)
2097  do k = datr_ks(handle), datr_ke(handle)
2098  vely_recv(k,i,j) = ( work2_recv(k,i,j+1) + work2_recv(k,i,j) ) * 0.5_rp
2099  enddo
2100  enddo
2101  enddo
2102 
2103  j = daughter_ja(handle)
2104 !OCL XFILL
2105  do i = 1, daughter_ia(handle)
2106  do k = datr_ks(handle), datr_ke(handle)
2107  vely_recv(k,i,j) = work2_recv(k,i,j)
2108  enddo
2109  enddo
2110 
2111  call comm_vars8( vely_recv, 3 )
2112 
2113  call comm_wait ( velx_recv, 2, .false. )
2114  call comm_wait ( vely_recv, 3, .false. )
2115 
2116  endif
2117 
2118  tagbase = tagcomm + tag_rhot*order_tag_var
2119  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2120  work1_recv(:,:,:), & ! [OUT]
2121  tagbase, i_sclr, handle, & ! [IN]
2122  isu_tag, isu_tagf ) ! [INOUT]
2123 !OCL XFILL
2124  do j = 1, daughter_ja(handle)
2125  do i = 1, daughter_ia(handle)
2126  do k = datr_ks(handle), datr_ke(handle)
2127  pott_recv(k,i,j) = work1_recv(k,i,j) / dens_recv(k,i,j)
2128  enddo
2129  enddo
2130  enddo
2131 
2132  do iq = 1, bnd_qa
2133  tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
2134  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2135  work1_recv(:,:,:), & ! [OUT]
2136  tagbase, i_sclr, handle, & ! [IN]
2137  isu_tag, isu_tagf ) ! [INOUT]
2138 !OCL XFILL
2139  do j = 1, daughter_ja(handle)
2140  do i = 1, daughter_ia(handle)
2141  do k = datr_ks(handle), datr_ke(handle)
2142  qtrc_recv(k,i,j,iq) = work1_recv(k,i,j)
2143  enddo
2144  enddo
2145  enddo
2146  enddo
2147 
2148  call prof_rapend ('NEST_unpack_C', 2)
2149  call prof_rapend ('NEST_total_C', 2)
2150 
2151  else
2152  log_error("COMM_CARTESC_NEST_nestdown",*) 'internal error'
2153  call prc_abort
2154  endif
2155 
2156  return
2157  end subroutine comm_cartesc_nest_nestdown
2158 
2159  !-----------------------------------------------------------------------------
2161  subroutine comm_cartesc_nest_recvwait_issue( &
2162  HANDLE, &
2163  BND_QA )
2164  use scale_prc, only: &
2165  prc_abort
2166  implicit none
2167 
2168  integer, intent(in) :: HANDLE
2169  integer, intent(in) :: BND_QA
2170 
2171  integer :: isu_tag, isu_tagf
2172  integer :: tagbase, tagcomm
2173  integer :: ierr
2174  integer :: iq
2175  !---------------------------------------------------------------------------
2176 
2177  if( .NOT. use_nesting ) return
2178 
2179  if ( bnd_qa > i_bndqa ) then
2180  log_error("COMM_CARTESC_NEST_recvwait_issue",*) 'internal error: about BND_QA'
2181  call prc_abort
2182  endif
2183 
2184  tagcomm = intercomm_id(handle) * order_tag_comm
2185 
2186  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2187 
2188  !##### parent [wait issue] #####
2189 
2190  call prof_rapstart('NEST_total_P', 2)
2191  call prof_rapstart('NEST_wait_P', 2)
2192 
2193  nwait_p = nwait_p + 1
2194  !LOG_INFO("COMM_CARTESC_NEST_recvwait_issue",'(1X,A,I5,A)') "NestIDC [P]: que wait ( ", nwait_p, " )"
2195 
2196  call comm_cartesc_nest_issuer_of_wait( handle )
2197 
2198  if ( online_aggressive_comm ) then
2199  ! nothing to do
2200  else
2201  call mpi_barrier(intercomm_daughter, ierr)
2202  endif
2203 
2204  call prof_rapend ('NEST_wait_P', 2)
2205  call prof_rapend ('NEST_total_P', 2)
2206 
2207  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2208 
2209  !##### child [receive issue] #####
2210 
2211  call prof_rapstart('NEST_total_C', 2)
2212 
2213  nrecv = nrecv + 1
2214  log_info("COMM_CARTESC_NEST_recvwait_issue",'(1X,A,I5,A)') "NestIDC [C]: que recv ( ", nrecv, " )"
2215 
2216  !*** reset issue tag and request control
2217  !--- do not change the calling order below;
2218  !--- it should be consistent with the order in "COMM_CARTESC_NEST_nestdown"
2219  isu_tag = 0
2220  isu_tagf = 0
2221  rq_ctl_d = 0
2222 
2223  tagbase = tagcomm + tag_dens*order_tag_var
2224  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf )
2225 
2226  tagbase = tagcomm + tag_momz*order_tag_var
2227  if ( online_use_velz ) then
2228  call comm_cartesc_nest_issuer_of_receive( tagbase, i_zstg, handle, isu_tag, isu_tagf )
2229  endif
2230 
2231  tagbase = tagcomm + tag_momx*order_tag_var
2232  if ( online_no_rotate ) then
2233  call comm_cartesc_nest_issuer_of_receive( tagbase, i_xstg, handle, isu_tag, isu_tagf )
2234  else
2235  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf )
2236  endif
2237 
2238  tagbase = tagcomm + tag_momy*order_tag_var
2239  if ( online_no_rotate ) then
2240  call comm_cartesc_nest_issuer_of_receive( tagbase, i_ystg, handle, isu_tag, isu_tagf )
2241  else
2242  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf )
2243  endif
2244 
2245  tagbase = tagcomm + tag_rhot*order_tag_var
2246  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf )
2247 
2248  do iq = 1, bnd_qa
2249  tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
2250  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf )
2251  enddo
2252 
2253  rq_tot_d = rq_ctl_d
2254 
2255  call prof_rapend('NEST_total_C', 2)
2256 
2257  else
2258  log_error("COMM_CARTESC_NEST_recvwait_issue",*) 'internal error'
2259  call prc_abort
2260  endif
2261 
2262  return
2263  end subroutine comm_cartesc_nest_recvwait_issue
2264 
2265  !-----------------------------------------------------------------------------
2267  subroutine comm_cartesc_nest_recv_cancel( &
2268  HANDLE )
2269  use scale_prc, only: &
2270  prc_abort
2271  implicit none
2272 
2273  integer, intent(in) :: HANDLE
2274 
2275  !logical :: flag
2276  !integer :: istatus(MPI_STATUS_SIZE)
2277 
2278  integer :: rq
2279  integer :: ierr
2280  !---------------------------------------------------------------------------
2281 
2282  if( .NOT. use_nesting ) return
2283 
2284  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2285 
2286  !##### parent #####
2287  ! Nothing to do
2288 
2289  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2290 
2291  !##### child #####
2292 
2293  log_info("COMM_CARTESC_NEST_recv_cancel",'(1X,A,I5,A)') "NestIDC [C]: CANCEL recv ( ", nrecv, " )"
2294 
2295  do rq = 1, rq_tot_d
2296  if ( ireq_d(rq) /= mpi_request_null ) then
2297 
2298  call mpi_cancel(ireq_d(rq), ierr)
2299 
2300 ! call MPI_TEST_CANCELLED(istatus, flag, ierr)
2301 ! if ( .NOT. flag ) then
2302 ! LOG_ERROR("COMM_CARTESC_NEST_recv_cancel",*) 'receive actions do not cancelled, req = ', rq
2303 ! endif
2304  endif
2305  enddo
2306 
2307  else
2308  log_error_cont(*) 'internal error'
2309  call prc_abort
2310  endif
2311 
2312  return
2313  end subroutine comm_cartesc_nest_recv_cancel
2314 
2315  !-----------------------------------------------------------------------------
2318  pvar, &
2319  dvar, &
2320  tagbase, &
2321  id_stag, &
2322  HANDLE, &
2323  isu_tag, &
2324  isu_tagf, &
2325  flag_dens )
2326  use scale_prc, only: &
2327  prc_abort
2328  use scale_comm_cartesc, only: &
2330  use scale_interp, only: &
2332  implicit none
2333 
2334  real(RP), intent(in) :: pvar(:,:,:)
2335  real(RP), intent(out) :: dvar(:,:,:)
2336  integer, intent(in) :: tagbase
2337  integer, intent(in) :: id_stag
2338  integer, intent(in) :: HANDLE
2339  integer, intent(inout) :: isu_tag
2340  integer, intent(inout) :: isu_tagf
2341 
2342  logical , intent(in), optional :: flag_dens
2343 
2344  integer :: ileng, tag, target_rank
2345 
2346  integer :: xloc, yloc
2347  integer :: gxs, gxe, gys, gye ! for large domain
2348  integer :: pxs, pxe, pys, pye ! for parent domain
2349  integer :: zs, ze
2350 
2351  integer :: ig, rq, yp
2352  logical :: no_zstag = .true.
2353  logical :: logarithmic = .false.
2354 
2355  integer :: ierr
2356  integer :: i, j
2357  !---------------------------------------------------------------------------
2358 
2359  if( .NOT. use_nesting ) return
2360 
2361  logarithmic = .false.
2362  if ( present(flag_dens) ) then
2363  if( flag_dens ) logarithmic = .true.
2364  endif
2365 
2366  if ( id_stag == i_sclr ) then
2367  no_zstag = .true.
2368  ig = i_sclr
2369  elseif( id_stag == i_zstg ) then
2370  no_zstag = .false.
2371  ig = i_zstg
2372  elseif( id_stag == i_xstg ) then
2373  no_zstag = .true.
2374  ig = i_xstg
2375  elseif( id_stag == i_ystg ) then
2376  no_zstag = .true.
2377  ig = i_ystg
2378  endif
2379 
2380  if ( no_zstag ) then
2381  ileng = (parent_ka(handle) ) * parent_ia(handle) * parent_ja(handle)
2382  else
2383  ileng = (parent_ka(handle)+1) * parent_ia(handle) * parent_ja(handle)
2384  endif
2385 
2386  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2387 
2388  !##### parent #####
2389  rq = rq_ctl_p
2390 
2391  do yp = 1, num_yp
2392  rq = rq + 1
2393 
2394  ! send data to multiple daughter processes
2395  target_rank = comm_cartesc_nest_tile_list_yp(yp)
2396  tag = tagbase + yp
2397 
2398  call mpi_isend( pvar, &
2399  ileng, &
2400  comm_datatype, &
2401  target_rank, &
2402  tag, &
2404  ireq_p(rq), &
2405  ierr )
2406 
2407  dvar(:,:,:) = -1.0_rp ! input as a dummy value
2408  enddo
2409 
2410  rq_ctl_p = rq
2411 
2412  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2413 
2414  !##### child #####
2415  rq = rq_ctl_d
2416 
2417  do yp = 1, comm_cartesc_nest_tile_all
2418  rq = rq + 1
2419 
2420  xloc = mod( yp-1, comm_cartesc_nest_tile_num_x ) + 1
2421  yloc = int( real(yp-1) / real(COMM_CARTESC_NEST_TILE_NUM_X) ) + 1
2422 
2423  gxs = parent_imax(handle) * (xloc-1) + 1
2424  gxe = parent_imax(handle) * xloc
2425  gys = parent_jmax(handle) * (yloc-1) + 1
2426  gye = parent_jmax(handle) * yloc
2427 
2428  pxs = prnt_is(handle)
2429  pxe = prnt_ie(handle)
2430  pys = prnt_js(handle)
2431  pye = prnt_je(handle)
2432 
2433  if ( no_zstag ) then
2434  isu_tag = isu_tag + 1
2435 
2436  zs = 1
2437  ze = parent_ka(handle)
2438 !OCL XFILL
2439  buffer_ref_3d(zs:ze,gxs:gxe,gys:gye) = recvbuf_3d(zs:ze,pxs:pxe,pys:pye,isu_tag)
2440  else
2441  isu_tagf = isu_tagf + 1
2442 
2443  zs = 0
2444  ze = parent_ka(handle)
2445 !OCL XFILL
2446  buffer_ref_3df(zs:ze,gxs:gxe,gys:gye) = recvbuf_3df(zs:ze,pxs:pxe,pys:pye,isu_tagf)
2447  endif
2448 
2449  if ( isu_tag > max_isu .OR. isu_tagf > max_isuf ) then
2450  log_error("COMM_CARTESC_NEST_intercomm_nestdown_3D",*) 'Exceeded maximum issue'
2451  log_error_cont(*) 'isu_tag = ', isu_tag
2452  log_error_cont(*) 'isu_tagf = ', isu_tagf
2453  call prc_abort
2454  endif
2455 
2456  enddo
2457 
2458  rq_ctl_d = rq
2459 
2460  if ( no_zstag ) then
2461  call interp_interp3d( itp_nh, & ! [IN]
2462  tileal_ka(handle), & ! [IN]
2463  tileal_ia(handle), & ! [IN]
2464  tileal_ja(handle), & ! [IN]
2465  daughter_ka(handle), & ! [IN]
2466  datr_ks(handle), & ! [IN]
2467  datr_ke(handle), & ! [IN]
2468  daughter_ia(handle), & ! [IN]
2469  daughter_ja(handle), & ! [IN]
2470  igrd( :,:,:,ig), & ! [IN]
2471  jgrd( :,:,:,ig), & ! [IN]
2472  hfact( :,:,:,ig), & ! [IN]
2473  kgrd(:,:,:,:,:,ig), & ! [IN]
2474  vfact(:,:,:,:,:,ig), & ! [IN]
2475  buffer_ref_3d(:,:,:), & ! [INOUT]
2476  dvar(:,:,:), & ! [OUT]
2477  logwgt = logarithmic ) ! [IN]
2478  else
2479  call interp_interp3d( itp_nh, & ! [IN]
2480  tileal_ka(handle)+1, & ! [IN]
2481  tileal_ia(handle), & ! [IN]
2482  tileal_ja(handle), & ! [IN]
2483  daughter_ka(handle), & ! [IN]
2484  datr_ks(handle), & ! [IN]
2485  datr_ke(handle), & ! [IN]
2486  daughter_ia(handle), & ! [IN]
2487  daughter_ja(handle), & ! [IN]
2488  igrd( :,:,:,ig), & ! [IN]
2489  jgrd( :,:,:,ig), & ! [IN]
2490  hfact( :,:,:,ig), & ! [IN]
2491  kgrd(:,:,:,:,:,ig), & ! [IN]
2492  vfact(:,:,:,:,:,ig), & ! [IN]
2493  buffer_ref_3df(:,:,:), & ! [INOUT]
2494  dvar(:,:,:), & ! [OUT]
2495  logwgt = logarithmic ) ! [IN]
2496  endif
2497 
2498  do j = 1, daughter_ja(handle)
2499  do i = 1, daughter_ia(handle)
2500  dvar( 1:datr_ks(handle)-1,i,j) = 0.0_rp
2501  dvar(datr_ke(handle)+1:daughter_ka(handle) ,i,j) = 0.0_rp
2502  enddo
2503  enddo
2504 
2505  else
2506  log_error("COMM_CARTESC_NEST_intercomm_nestdown_3D",*) 'internal error'
2507  call prc_abort
2508  endif
2509 
2510  return
2512 
2513  !-----------------------------------------------------------------------------
2516  tagbase, &
2517  id_stag, &
2518  HANDLE, &
2519  isu_tag, &
2520  isu_tagf )
2521  use scale_prc, only: &
2522  prc_myrank, &
2523  prc_abort
2524  use scale_comm_cartesc, only: &
2526  implicit none
2527 
2528  integer, intent(in) :: tagbase
2529  integer, intent(in) :: id_stag
2530  integer, intent(in) :: HANDLE
2531  integer, intent(inout) :: isu_tag
2532  integer, intent(inout) :: isu_tagf
2533 
2534  integer :: ierr, ileng
2535  integer :: tag, target_rank
2536 
2537  integer :: ig, rq, yp
2538  logical :: no_zstag = .true.
2539  !---------------------------------------------------------------------------
2540 
2541  if( .NOT. use_nesting ) return
2542 
2543  if ( id_stag == i_sclr ) then
2544  no_zstag = .true.
2545  ig = i_sclr
2546  elseif( id_stag == i_zstg ) then
2547  no_zstag = .false.
2548  ig = i_zstg
2549  elseif( id_stag == i_xstg ) then
2550  no_zstag = .true.
2551  ig = i_xstg
2552  elseif( id_stag == i_ystg ) then
2553  no_zstag = .true.
2554  ig = i_ystg
2555  endif
2556 
2557  if ( no_zstag ) then
2558  ileng = (parent_ka(handle) ) * parent_ia(handle) * parent_ja(handle)
2559  else
2560  ileng = (parent_ka(handle)+1) * parent_ia(handle) * parent_ja(handle)
2561  endif
2562 
2563  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2564 
2565  !##### parent #####
2566  ! nothing to do
2567 
2568  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2569 
2570  !##### child #####
2571  rq = rq_ctl_d
2572 
2573  do yp = 1, comm_cartesc_nest_tile_all
2574  rq = rq + 1
2575 
2576  target_rank = comm_cartesc_nest_tile_list_d(yp,prc_myrank+1)
2577  tag = tagbase + call_order(yp)
2578 
2579  if ( no_zstag ) then
2580  isu_tag = isu_tag + 1
2581 
2582  recvbuf_3d(:,:,:,isu_tag) = 0.0_rp
2583 
2584  call mpi_irecv( recvbuf_3d(:,:,:,isu_tag), &
2585  ileng, &
2586  comm_datatype, &
2587  target_rank, &
2588  tag, &
2589  intercomm_parent, &
2590  ireq_d(rq), &
2591  ierr )
2592  else
2593  isu_tagf = isu_tagf + 1
2594 
2595  recvbuf_3df(:,:,:,isu_tagf) = 0.0_rp
2596 
2597  call mpi_irecv( recvbuf_3df(:,:,:,isu_tagf), &
2598  ileng, &
2599  comm_datatype, &
2600  target_rank, &
2601  tag, &
2602  intercomm_parent, &
2603  ireq_d(rq), &
2604  ierr )
2605  endif
2606 
2607  enddo
2608 
2609  if ( isu_tag > max_isu .OR. isu_tagf > max_isuf ) then
2610  log_error("COMM_CARTESC_NEST_issuer_of_receive_3D",*) 'Exceeded maximum issue'
2611  log_error_cont(*) 'isu_tag = ', isu_tag
2612  log_error_cont(*) 'isu_tagf = ', isu_tagf
2613  call prc_abort
2614  endif
2615 
2616  rq_ctl_d = rq
2617 
2618  else
2619  log_error("COMM_CARTESC_NEST_issuer_of_receive_3D",*) 'internal error'
2620  call prc_abort
2621  endif
2622 
2623  return
2625 
2626  !-----------------------------------------------------------------------------
2629  HANDLE )
2630  use scale_prc, only: &
2631  prc_abort
2632  implicit none
2633 
2634  integer, intent(in) :: HANDLE
2635  !---------------------------------------------------------------------------
2636 
2637  if( .NOT. use_nesting ) return
2638 
2639  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2640 
2641  !##### parent #####
2642  call comm_cartesc_nest_waitall( rq_tot_p, ireq_p(:) )
2643 
2644  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2645 
2646  !##### child #####
2647  ! nothing to do
2648 
2649  else
2650  log_error("COMM_CARTESC_NEST_issuer_of_wait_3D",*) 'internal error'
2651  call prc_abort
2652  endif
2653 
2654  return
2656 
2657  !-----------------------------------------------------------------------------
2659  subroutine comm_cartesc_nest_waitall( &
2660  req_count, &
2661  ireq )
2662  use scale_prc, only: &
2663  prc_abort
2664  implicit none
2665 
2666  integer, intent(in) :: req_count
2667  integer, intent(inout) :: ireq(max_rq)
2668 
2669  integer :: i
2670  integer :: ierr
2671  integer :: istatus(mpi_status_size,req_count)
2672  integer :: req_count2
2673  integer :: ireq2(max_rq)
2674 
2675 ! logical :: flag = .false.
2676 ! integer(8) :: num = 0
2677  !---------------------------------------------------------------------------
2678 
2679  if( .NOT. use_nesting ) return
2680 
2681  req_count2 = 0
2682  do i = 1, req_count
2683  if ( ireq(i) /= mpi_request_null ) then
2684  req_count2 = req_count2 + 1
2685  ireq2(req_count2) = ireq(i)
2686  endif
2687  enddo
2688 
2689  if( req_count2 /= 0 ) call mpi_waitall( req_count2, ireq2(1:req_count2), istatus, ierr )
2690 
2691 ! do while ( .NOT. flag )
2692 ! num = num + 1
2693 ! call MPI_TESTALL( req_count, ireq, flag, istatus, ierr )
2694 !
2695 ! if ( num > ONLINE_WAIT_LIMIT ) then
2696 ! LOG_ERROR("COMM_CARTESC_NEST_waitall",'(1x,A)') 'over the limit of waiting time [NESTCOM]'
2697 ! LOG_ERROR_CONT('(1x,A)') 'over the limit of waiting time [NESTCOM]'
2698 ! call PRC_abort
2699 ! endif
2700 ! enddo
2701 
2702  return
2703  end subroutine comm_cartesc_nest_waitall
2704 
2705  !-----------------------------------------------------------------------------
2707  subroutine comm_cartesc_nest_test( &
2708  HANDLE )
2709  use scale_prc, only: &
2710  prc_abort
2711  implicit none
2712 
2713  integer, intent(in) :: HANDLE
2714 
2715  integer :: istatus(mpi_status_size)
2716  integer :: ierr
2717  logical :: flag
2718  !---------------------------------------------------------------------------
2719 
2720  if( .NOT. use_nesting ) return
2721 
2722  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2723 
2724  !##### parent #####
2725 
2726  call prof_rapstart('NEST_test_P', 2)
2727  if ( rq_ctl_p > 0 ) call mpi_test(ireq_p(1), flag, istatus, ierr)
2728  call prof_rapend('NEST_test_P', 2)
2729 
2730  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2731 
2732  !##### child #####
2733 
2734  call prof_rapstart('NEST_test_C', 2)
2735  if ( rq_ctl_d > 0 ) call mpi_test(ireq_d(1), flag, istatus, ierr)
2736  call prof_rapend('NEST_test_C', 2)
2737 
2738  else
2739  log_error("COMM_CARTESC_NEST_test",*) 'error'
2740  call prc_abort
2741  endif
2742 
2743  return
2744  end subroutine comm_cartesc_nest_test
2745 
2746  !-----------------------------------------------------------------------------
2748  subroutine comm_cartesc_nest_disconnect
2749  use scale_prc, only: &
2751  implicit none
2752 
2753  integer :: ierr
2754  !---------------------------------------------------------------------------
2755 
2756  if( .NOT. use_nesting ) return
2757 
2758  log_info("COMM_CARTESC_NEST_disconnect",'(1x,A)') 'Waiting finish of whole processes'
2759  call mpi_barrier(prc_global_comm_world, ierr)
2760 
2761  if ( online_iam_parent ) then
2762  !LOG_INFO("COMM_CARTESC_NEST_disconnect",'(1x,A)') 'Waiting finish of whole processes as a parent'
2763  !call MPI_BARRIER(INTERCOMM_DAUGHTER, ierr)
2764  call mpi_comm_free(intercomm_daughter, ierr)
2765  log_info("COMM_CARTESC_NEST_disconnect",'(1x,A)') 'Disconnected communication with child'
2766  endif
2767 
2768  if ( online_iam_daughter ) then
2769  !LOG_INFO("COMM_CARTESC_NEST_disconnect",'(1x,A)') 'Waiting finish of whole processes as a child'
2770  !call MPI_BARRIER(INTERCOMM_PARENT, ierr)
2771  call mpi_comm_free(intercomm_parent, ierr)
2772  log_info("COMM_CARTESC_NEST_disconnect",'(1x,A)') 'Disconnected communication with parent'
2773  endif
2774 
2775  return
2776  end subroutine comm_cartesc_nest_disconnect
2777 
2778 end module scale_comm_cartesc_nest
subroutine, public comm_cartesc_nest_nestdown(HANDLE, BND_QA, DENS_send, MOMZ_send, MOMX_send, MOMY_send, RHOT_send, QTRC_send, DENS_recv, VELZ_recv, VELX_recv, VELY_recv, POTT_recv, QTRC_recv)
Boundary data transfer from parent to daughter: nestdown.
integer, dimension(2), public parent_kmax
parent max number in z-direction
module DEBUG
Definition: scale_debug.F90:11
integer, dimension(2), public prnt_je
end index in y-direction in parent
integer, public jmax
of computational cells: y, local
subroutine, public comm_cartesc_nest_disconnect
[finalize: disconnect] Inter-communication
integer, public comm_world
communication world ID
integer, dimension(2), public daughter_ia
daughter max number in x-direction (with halo)
integer, public prc_global_domainid
my domain ID in global communicator
Definition: scale_prc.F90:83
integer, parameter, public khalo
of halo cells: z
integer, dimension(2), public daughter_ja
daughter max number in y-direction (with halo)
subroutine comm_cartesc_nest_issuer_of_receive_3d(tagbase, id_stag, HANDLE, isu_tag, isu_tagf)
[substance of issuer] Inter-communication from parent to daughter: nestdown
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lonuv
longitude at staggered point (uv) [rad,0-2pi]
integer, public comm_cartesc_nest_tile_num_y
parent tile number in y-direction
logical, public online_iam_parent
a flag to say "I am a parent"
integer, dimension(2), public prnt_ke
end index in z-direction in parent
subroutine, public comm_cartesc_nest_domain_relate(HANDLE)
Solve relationship between ParentDomain & Daughter Domain.
integer, public imax
of computational cells: x, local
integer, public time_nstep
total steps [number]
Definition: scale_time.F90:75
integer, dimension(2), public daughter_lkmax
daughter max number in lz-direction
integer, dimension(2), public daughter_ka
daughter max number in z-direction (with halo)
integer, public ia
of whole cells: x, local, with HALO
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_domain_catalogue
domain latlon catalogue [rad]
module Atmosphere Grid CartesianC metirc
integer, public comm_datatype
datatype of variable
integer, dimension(2), public parent_nstep
parent step [number]
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lat
latitude [rad,-pi,pi]
module INTERPOLATION
integer, dimension(2), public datr_is
start index in x-direction in daughter
subroutine, public interp_factor3d(npoints, KA_ref, KS_ref, KE_ref, IA_ref, JA_ref, lon_ref, lat_ref, hgt_ref, KA, KS, KE, IA, JA, lon, lat, hgt, idx_i, idx_j, hfact, idx_k, vfact)
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_fz
geopotential height [m] (wxy)
integer, public ja
of whole cells: y, local, with HALO
module process / cartesC
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
subroutine, public interp_interp3d(npoints, KA_ref, IA_ref, JA_ref, KA, KS, KE, IA, JA, idx_i, idx_j, hfact, idx_k, vfact, val_ref, val, logwgt)
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lon
longitude [rad,0-2pi]
integer, public prc_nprocs
myrank in local communicator
Definition: scale_prc.F90:88
subroutine, public file_open(basename, fid, mode, single, aggregate, rankid, postfix)
Definition: scale_file.F90:464
integer, dimension(2), public daughter_kmax
daughter max number in z-direction
module COMMUNICATION
integer, public comm_cartesc_nest_bnd_qa
number of tracer treated in nesting system
integer, public is
start point of inner domain: x, local
module file
Definition: scale_file.F90:15
integer, public ie
end point of inner domain: x, local
integer, public handling_num
handing number of nesting relation
integer, dimension(2), public parent_ia
parent max number in x-direction (with halo)
module TRACER
module Index
Definition: scale_index.F90:11
integer, dimension(2), public datr_je
end index in y-direction in daughter
subroutine, public comm_cartesc_nest_recv_cancel(HANDLE)
Sub-command for data transfer from parent to daughter: nestdown.
subroutine, public comm_cartesc_nest_setup(QA_MP, MP_TYPE_in, inter_parent, inter_child)
Setup.
module atmosphere / hydrometeor
integer, dimension(10), public comm_cartesc_nest_filiation
index of parent-daughter relation (p>0, d<0)
integer, dimension(2), public datr_js
start index in y-direction in daughter
integer, dimension(2), public parent_jmax
parent max number in y-direction
real(dp), public time_dtsec
time interval of model [sec]
Definition: scale_time.F90:38
integer, dimension(2), public prnt_ie
end index in x-direction in parent
integer, public comm_cartesc_nest_interp_level
horizontal interpolation level
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latuy
latitude at staggered point (uy) [rad,-pi,pi]
module atmosphere / grid / cartesC index
integer, public ke
end point of inner domain: z, local
integer, dimension(2), public daughter_okmax
daughter max number in oz-direction
logical, public online_iam_daughter
a flag to say "I am a daughter"
subroutine, public comm_cartesc_nest_recvwait_issue(HANDLE, BND_QA)
Sub-command for data transfer from parent to daughter: nestdown.
integer function, public io_get_available_fid()
search & get available file ID
Definition: scale_io.F90:313
module PROCESS
Definition: scale_prc.F90:11
integer, public je
end point of inner domain: y, local
integer, dimension(2), public daughter_jmax
daughter max number in y-direction
real(dp), dimension(2), public daughter_dtsec
daughter DT [sec]
integer, dimension(2), public daughter_nstep
daughter steps [number]
integer, dimension(2), public tileal_ja
cells of all tiles in y-direction
module TIME
Definition: scale_time.F90:16
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lonxv
longitude at staggered point (xv) [rad,0-2pi]
integer, public ks
start point of inner domain: z, local
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lonuy
longitude at staggered point (uy) [rad,0-2pi]
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_metric_rotc
rotation coefficient
integer, dimension(2), public parent_ja
parent max number in y-direction (with halo)
integer, dimension(2), public datr_ie
end index in x-direction in daughter
integer, dimension(2), public prnt_is
start index in x-direction in parent
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:89
integer, public kmax
of computational cells: z, local
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
subroutine, public comm_cartesc_nest_test(HANDLE)
[check communication status] Inter-communication
module CONSTANT
Definition: scale_const.F90:11
module Communication CartesianC nesting
integer, public js
start point of inner domain: y, local
subroutine, public interp_setup(weight_order, search_limit)
Setup.
integer, dimension(2), public tileal_ia
cells of all tiles in x-direction
integer, dimension(2), public parent_lkmax
parent max number in lz-direction
subroutine, public comm_cartesc_nest_domain_shape(tilei, tilej, cxs, cxe, cys, cye, pxs, pxe, pys, pye, iloc)
Return shape of ParentDomain at the specified rank (for offline)
integer, parameter, public h_short
Character length (short=16)
Definition: scale_io.F90:43
integer, public prc_num_y
y length of 2D processor topology
integer, dimension(2), public daughter_imax
daughter max number in x-direction
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:157
integer, dimension(:), allocatable, public comm_cartesc_nest_tile_id
parent tile real id
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_cz
geopotential height [m] (zxy)
subroutine comm_cartesc_nest_issuer_of_wait_3d(HANDLE)
[substance of issuer] Inter-communication from parent to daughter: nestdown
module profiler
Definition: scale_prof.F90:11
logical, public prc_ismaster
master process in local communicator?
Definition: scale_prc.F90:90
integer, public comm_cartesc_nest_interp_weight_order
horizontal interpolation weight order
integer, public prc_global_comm_world
global communicator
Definition: scale_prc.F90:78
integer, dimension(2), public prnt_js
start index in y-direction in parent
module Atmosphere GRID CartesC Real(real space)
integer, dimension(2), public datr_ke
end index in z-direction in daughter
module PRECISION
real(dp), dimension(2), public parent_dtsec
parent DT [sec]
integer, dimension(2), public parent_imax
parent max number in x-direction
integer, dimension(2), public parent_ka
parent max number in z-direction (with halo)
integer, public ka
of whole cells: z, local, with HALO
integer, dimension(2), public parent_okmax
parent max number in oz-direction
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latuv
latitude at staggered point (uv) [rad,-pi,pi]
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latxv
latitude at staggered point (xv) [rad,-pi,pi]
subroutine comm_cartesc_nest_intercomm_nestdown_3d(pvar, dvar, tagbase, id_stag, HANDLE, isu_tag, isu_tagf, flag_dens)
Inter-communication from parent to daughter: nestdown.
module STDIO
Definition: scale_io.F90:10
integer, public debug_domain_num
Definition: scale_debug.F90:38
integer, parameter, public n_hyd
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:210
integer, dimension(2), public prnt_ks
start index in z-direction in parent
integer, dimension(2), public datr_ks
start index in z-direction in daughter
integer, public comm_cartesc_nest_tile_num_x
parent tile number in x-direction
subroutine comm_cartesc_nest_waitall(req_count, ireq)
[substance of comm_wait] Inter-communication
integer, dimension(2), public tileal_ka
cells of all tiles in z-direction
integer, public prc_num_x
x length of 2D processor topology