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