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