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