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 = .false.
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  allocate( pd_tile_num(0:parent_prc_nprocs(handle)-1,2) )
770 
771  k = 0 ! MPI process number starts from zero
772  do j = 1, parent_prc_num_y(handle)
773  do i = 1, parent_prc_num_x(handle)
774  pd_tile_num(k,1) = i
775  pd_tile_num(k,2) = j
776  k = k + 1
777  enddo
778  enddo
779 
780  !--- SW search
781  hit = .false.
782  do i = 1, parent_prc_nprocs(handle)
783  wid_lon = abs((latlon_catalog(i,i_sw,i_lon) - latlon_catalog(i,i_se,i_lon)) &
784  / real( parent_imax(handle)-1, kind=rp )) * 0.8_rp
785  wid_lat = abs((latlon_catalog(i,i_sw,i_lat) - latlon_catalog(i,i_nw,i_lat)) &
786  / real( parent_jmax(handle)-1, kind=rp )) * 0.8_rp
787 
788  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. &
789  corner_loc(i_sw,i_lat) >= min(latlon_catalog(i,i_sw,i_lat),latlon_catalog(i,i_se,i_lat))-wid_lat .AND. &
790  corner_loc(i_sw,i_lon) <= max(latlon_catalog(i,i_ne,i_lon),latlon_catalog(i,i_se,i_lon))+wid_lon .AND. &
791  corner_loc(i_sw,i_lat) <= max(latlon_catalog(i,i_ne,i_lat),latlon_catalog(i,i_nw,i_lat))+wid_lat ) then
792 
793  pd_sw_tile = i-1 ! MPI process number starts from zero
794  hit = .true.
795  exit ! exit loop
796  endif
797  enddo
798  if ( .NOT. hit ) then
799  write(*,*) 'xxx region of daughter domain is larger than that of parent: SW search'
800  write(*,*) ' at rank:', prc_myrank, ' of domain:', online_domain_num
801  if( io_l ) write(io_fid_log,'(1x,A)') 'xxx region of daughter domain is larger than that of parent: SW search'
802  if( io_l ) write(io_fid_log,*) ' grid width: half width in lat:', wid_lat, ' half width in lon:', wid_lon
803  if( io_l ) write(io_fid_log,'(1x,A,F12.6)') ' daughter local (me): LON=',corner_loc(i_sw,i_lon)
804  do i = 1, parent_prc_nprocs(handle)
805  if( io_l ) write(io_fid_log,'(1x,A,F12.6,1x,F12.6)') ' parent local SW-NE: LON=', &
806  latlon_catalog(i,i_sw,i_lon) ,latlon_catalog(i,i_ne,i_lon)
807  enddo
808  if( io_l ) write(io_fid_log,'(1x,A,F12.6)') ' daughter local (me): LAT=',corner_loc(i_sw,i_lat)
809  do i = 1, parent_prc_nprocs(handle)
810  if( io_l ) write(io_fid_log,'(1x,A,F12.6,1x,F12.6)') ' parent local SW-NE: LAT=', &
811  latlon_catalog(i,i_sw,i_lat) ,latlon_catalog(i,i_ne,i_lat)
812  enddo
813  call prc_mpistop
814  endif
815 
816  !--- NE search
817  hit = .false.
818  do i = parent_prc_nprocs(handle), 1, -1
819  wid_lon = abs((latlon_catalog(i,i_nw,i_lon) - latlon_catalog(i,i_ne,i_lon)) &
820  / real( parent_imax(handle)-1, kind=rp )) * 0.8_rp
821  wid_lat = abs((latlon_catalog(i,i_se,i_lat) - latlon_catalog(i,i_ne,i_lat)) &
822  / real( parent_jmax(handle)-1, kind=rp )) * 0.8_rp
823 
824  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. &
825  corner_loc(i_ne,i_lat) >= min(latlon_catalog(i,i_sw,i_lat),latlon_catalog(i,i_se,i_lat))-wid_lat .AND. &
826  corner_loc(i_ne,i_lon) <= max(latlon_catalog(i,i_ne,i_lon),latlon_catalog(i,i_se,i_lon))+wid_lon .AND. &
827  corner_loc(i_ne,i_lat) <= max(latlon_catalog(i,i_ne,i_lat),latlon_catalog(i,i_nw,i_lat))+wid_lat ) then
828 
829  pd_ne_tile = i-1 ! MPI process number starts from zero
830  hit = .true.
831  exit ! exit loop
832  endif
833  enddo
834  if ( .NOT. hit ) then
835  write(*,*) 'xxx region of daughter domain is larger than that of parent: NE search'
836  write(*,*) ' at rank:', prc_myrank, ' of domain:', online_domain_num
837  if( io_l ) write(io_fid_log,'(1x,A)') 'xxx region of daughter domain is larger than that of parent: NE search'
838  if( io_l ) write(io_fid_log,*) ' grid width: half width in lat:', wid_lat, ' half width in lon:', wid_lon
839  if( io_l ) write(io_fid_log,'(1x,A,F12.6)') ' daughter local (me): LON=',corner_loc(i_ne,i_lon)
840  do i = 1, parent_prc_nprocs(handle)
841  if( io_l ) write(io_fid_log,'(1x,A,F12.6,1x,F12.6)') ' parent local SW-NE: LON=', &
842  latlon_catalog(i,i_sw,i_lon) ,latlon_catalog(i,i_ne,i_lon)
843  enddo
844  if( io_l ) write(io_fid_log,'(1x,A,F12.6)') ' daughter local (me): LAT=',corner_loc(i_ne,i_lat)
845  do i = 1, parent_prc_nprocs(handle)
846  if( io_l ) write(io_fid_log,'(1x,A,F12.6,1x,F12.6)') ' parent local SW-NE: LAT=', &
847  latlon_catalog(i,i_sw,i_lat) ,latlon_catalog(i,i_ne,i_lat)
848  enddo
849  call prc_mpistop
850  endif
851 
852  nest_tile_num_x = pd_tile_num(pd_ne_tile,1) - pd_tile_num(pd_sw_tile,1) + 1
853  nest_tile_num_y = pd_tile_num(pd_ne_tile,2) - pd_tile_num(pd_sw_tile,2) + 1
854 
856 
857  if( io_l ) write(io_fid_log,'(1x,A)') '*** NEST: target process tile in parent domain'
858  k = 1
859  do j = 1, nest_tile_num_y
860  do i = 1, nest_tile_num_x
861  nest_tile_id(k) = pd_sw_tile + (i-1) + parent_prc_num_x(handle)*(j-1)
862  if( io_l ) write(io_fid_log,'(1x,A,I4,A,I6)') ' (', k, ') target mpi-process:', nest_tile_id(k)
863  k = k + 1
864  enddo
865  enddo
866 
867  return
868  end subroutine nest_domain_relate
869 
870  !-----------------------------------------------------------------------------
872  ! including definition array size with BND or not in Parent domain
873  subroutine nest_domain_shape ( &
874  tilei, & ! [out]
875  tilej, & ! [out]
876  cxs, cxe, & ! [out]
877  cys, cye, & ! [out]
878  pxs, pxe, & ! [out]
879  pys, pye, & ! [out]
880  iloc ) ! [in ]
881  implicit none
882 
883  integer, intent(out) :: tilei, tilej
884  integer, intent(out) :: cxs, cxe, cys, cye
885  integer, intent(out) :: pxs, pxe, pys, pye
886  integer, intent(in) :: iloc ! rank number; start from 1
887 
888  integer :: hdl = 1 ! handler number
889  integer :: rank
890  integer :: xloc, yloc
891  integer :: xlocg, ylocg ! location over whole domain
892  !---------------------------------------------------------------------------
893 
894  rank = nest_tile_id(iloc)
895  xloc = mod( iloc-1, nest_tile_num_x ) + 1
896  yloc = int( real(iloc-1) / real(NEST_TILE_NUM_X) ) + 1
897  xlocg = mod( rank, offline_parent_prc_num_x ) + 1
898  ylocg = int( real(rank) / real(OFFLINE_PARENT_PRC_NUM_X) ) + 1
899  tilei = parent_imax(hdl)
900  tilej = parent_jmax(hdl)
901 
902  cxs = tilei * (xloc-1) + 1
903  cxe = tilei * xloc
904  cys = tilej * (yloc-1) + 1
905  cye = tilej * yloc
906  pxs = 1
907  pxe = tilei
908  pys = 1
909  pye = tilej
910 
911  if ( xlocg == 1 ) then ! BND_W
912  tilei = tilei + 2
913  pxs = pxs + 2
914  pxe = pxe + 2
915  endif
916  if ( xlocg == offline_parent_prc_num_x ) then ! BND_E
917  tilei = tilei + 2
918  endif
919  if ( ylocg == 1 ) then ! BND_S
920  tilej = tilej + 2
921  pys = pys + 2
922  pye = pye + 2
923  endif
924  if ( ylocg == offline_parent_prc_num_y ) then ! BND_N
925  tilej = tilej + 2
926  endif
927 
928  return
929  end subroutine nest_domain_shape
930 
931  !-----------------------------------------------------------------------------
933  subroutine nest_comm_parentsize( &
934  HANDLE )
935  use scale_process, only: &
936  prc_mpistop, &
937  prc_nprocs, &
938  prc_myrank, &
940  use scale_rm_process, only: &
941  prc_num_x, &
942  prc_num_y
943  use scale_time, only: &
944  time_nstep, &
945  time_dtsec
946  use scale_comm, only: &
947  comm_bcast
948  implicit none
949 
950  integer, intent(in) :: HANDLE
951 
952  real(RP) :: buffer
953  integer :: datapack(14)
954  integer :: QA_OTHERSIDE
955  integer :: ireq1, ireq2, ierr1, ierr2, ileng
956  integer :: istatus(mpi_status_size)
957  integer :: tag
958  !---------------------------------------------------------------------------
959 
960  tag = intercomm_id(handle) * 100
961  ileng = 14
962 
963  if ( nest_filiation( intercomm_id(handle) ) > 0 ) then !--- parent
964  ! from parent to daughter
965  datapack( 1) = prc_nprocs
966  datapack( 2) = prc_num_x
967  datapack( 3) = prc_num_y
968  datapack( 4) = kmax
969  datapack( 5) = imax
970  datapack( 6) = jmax
971  datapack( 7) = ks
972  datapack( 8) = ke
973  datapack( 9) = is
974  datapack(10) = ie
975  datapack(11) = js
976  datapack(12) = je
977  datapack(13) = time_nstep
978  datapack(14) = nest_bnd_qa
979  buffer = time_dtsec
980 
981  if ( prc_ismaster ) then
982  call mpi_isend(datapack, ileng, mpi_integer, prc_myrank, tag, intercomm_daughter, ireq1, ierr1)
983  call mpi_isend(buffer, 1, mpi_double_precision, prc_myrank, tag+1, intercomm_daughter, ireq2, ierr2)
984  call mpi_wait(ireq1, istatus, ierr1)
985  call mpi_wait(ireq2, istatus, ierr2)
986  endif
987 
988  parent_prc_nprocs(handle) = datapack( 1)
989  parent_prc_num_x(handle) = datapack( 2)
990  parent_prc_num_y(handle) = datapack( 3)
991  parent_kmax(handle) = datapack( 4)
992  parent_imax(handle) = datapack( 5)
993  parent_jmax(handle) = datapack( 6)
994  prnt_ks(handle) = datapack( 7)
995  prnt_ke(handle) = datapack( 8)
996  prnt_is(handle) = datapack( 9)
997  prnt_ie(handle) = datapack(10)
998  prnt_js(handle) = datapack(11)
999  prnt_je(handle) = datapack(12)
1000  parent_nstep(handle) = datapack(13)
1001  parent_dtsec(handle) = buffer
1002 
1003  ! from daughter to parent
1004  if ( prc_ismaster ) then
1005  call mpi_irecv(datapack, ileng, mpi_integer, prc_myrank, tag+2, intercomm_daughter, ireq1, ierr1)
1006  call mpi_irecv(buffer, 1, mpi_double_precision, prc_myrank, tag+3, intercomm_daughter, ireq2, ierr2)
1007  call mpi_wait(ireq1, istatus, ierr1)
1008  call mpi_wait(ireq2, istatus, ierr2)
1009  endif
1010  call comm_bcast(datapack, ileng)
1011  call comm_bcast(buffer)
1012 
1013  daughter_prc_nprocs(handle) = datapack( 1)
1014  daughter_prc_num_x(handle) = datapack( 2)
1015  daughter_prc_num_y(handle) = datapack( 3)
1016  daughter_kmax(handle) = datapack( 4)
1017  daughter_imax(handle) = datapack( 5)
1018  daughter_jmax(handle) = datapack( 6)
1019  datr_ks(handle) = datapack( 7)
1020  datr_ke(handle) = datapack( 8)
1021  datr_is(handle) = datapack( 9)
1022  datr_ie(handle) = datapack(10)
1023  datr_js(handle) = datapack(11)
1024  datr_je(handle) = datapack(12)
1025  daughter_nstep(handle) = datapack(13)
1026  qa_otherside = datapack(14)
1027  daughter_dtsec(handle) = buffer
1028 
1029 
1030  elseif ( nest_filiation( intercomm_id(handle) ) < 0 ) then !--- daughter
1031  ! from parent to daughter
1032  if ( prc_ismaster ) then
1033  call mpi_irecv(datapack, ileng, mpi_integer, prc_myrank, tag, intercomm_parent, ireq1, ierr1)
1034  call mpi_irecv(buffer, 1, mpi_double_precision, prc_myrank, tag+1, intercomm_parent, ireq2, ierr2)
1035  call mpi_wait(ireq1, istatus, ierr1)
1036  call mpi_wait(ireq2, istatus, ierr2)
1037  endif
1038  call comm_bcast(datapack, ileng)
1039  call comm_bcast(buffer)
1040 
1041  parent_prc_nprocs(handle) = datapack( 1)
1042  parent_prc_num_x(handle) = datapack( 2)
1043  parent_prc_num_y(handle) = datapack( 3)
1044  parent_kmax(handle) = datapack( 4)
1045  parent_imax(handle) = datapack( 5)
1046  parent_jmax(handle) = datapack( 6)
1047  prnt_ks(handle) = datapack( 7)
1048  prnt_ke(handle) = datapack( 8)
1049  prnt_is(handle) = datapack( 9)
1050  prnt_ie(handle) = datapack(10)
1051  prnt_js(handle) = datapack(11)
1052  prnt_je(handle) = datapack(12)
1053  parent_nstep(handle) = datapack(13)
1054  qa_otherside = datapack(14)
1055  parent_dtsec(handle) = buffer
1056 
1057  ! from daughter to parent
1058  datapack( 1) = prc_nprocs
1059  datapack( 2) = prc_num_x
1060  datapack( 3) = prc_num_y
1061  datapack( 4) = kmax
1062  datapack( 5) = imax
1063  datapack( 6) = jmax
1064  datapack( 7) = ks
1065  datapack( 8) = ke
1066  datapack( 9) = is
1067  datapack(10) = ie
1068  datapack(11) = js
1069  datapack(12) = je
1070  datapack(13) = time_nstep
1071  datapack(14) = nest_bnd_qa
1072  buffer = time_dtsec
1073 
1074  if ( prc_ismaster ) then
1075  call mpi_isend(datapack, ileng, mpi_integer, prc_myrank, tag+2, intercomm_parent, ireq1, ierr1)
1076  call mpi_isend(buffer, 1, mpi_double_precision, prc_myrank, tag+3, intercomm_parent, ireq2, ierr2)
1077  call mpi_wait(ireq1, istatus, ierr1)
1078  call mpi_wait(ireq2, istatus, ierr2)
1079  endif
1080 
1081  daughter_prc_nprocs(handle) = datapack( 1)
1082  daughter_prc_num_x(handle) = datapack( 2)
1083  daughter_prc_num_y(handle) = datapack( 3)
1084  daughter_kmax(handle) = datapack( 4)
1085  daughter_imax(handle) = datapack( 5)
1086  daughter_jmax(handle) = datapack( 6)
1087  datr_ks(handle) = datapack( 7)
1088  datr_ke(handle) = datapack( 8)
1089  datr_is(handle) = datapack( 9)
1090  datr_ie(handle) = datapack(10)
1091  datr_js(handle) = datapack(11)
1092  datr_je(handle) = datapack(12)
1093  daughter_nstep(handle) = datapack(13)
1094  daughter_dtsec(handle) = buffer
1095  else
1096  write(*,*) 'xxx internal error [parentsize: nest/grid]'
1097  call prc_mpistop
1098  endif
1099 
1100  if( qa_otherside /= nest_bnd_qa ) then
1101  write(*,*) 'xxx ERROR: NUMBER of QA are not matched! [parentsize: nest/grid]'
1102  write(*,*) 'xxx check a flag of ONLINE_BOUNDARY_USE_QHYD.', qa_otherside, nest_bnd_qa
1103  call prc_mpistop
1104  endif
1105 
1106  return
1107  end subroutine nest_comm_parentsize
1108 
1109  !-----------------------------------------------------------------------------
1111  subroutine nest_comm_catalogue( &
1112  HANDLE )
1113  use scale_process, only: &
1114  prc_mpistop, &
1115  prc_nprocs, &
1116  prc_myrank, &
1117  prc_ismaster
1118  use scale_grid_real, only: &
1120  use scale_comm, only: &
1121  comm_datatype, &
1122  comm_bcast
1123  implicit none
1124 
1125  integer, intent(in) :: HANDLE
1126 
1127  integer :: ireq, ierr, ileng
1128  integer :: istatus(mpi_status_size)
1129  integer :: tag
1130  !---------------------------------------------------------------------------
1131 
1132  tag = intercomm_id(handle) * 100
1133 
1134  if ( nest_filiation( intercomm_id(handle) ) > 0 ) then !--- parent
1135  ileng = prc_nprocs * 4 * 2
1136  if ( prc_ismaster ) then
1137  call mpi_isend(real_domain_catalogue, ileng, comm_datatype, prc_myrank, tag, intercomm_daughter, ireq, ierr)
1138  call mpi_wait(ireq, istatus, ierr)
1139  endif
1140 
1141  elseif ( nest_filiation( intercomm_id(handle) ) < 0 ) then !--- daughter
1142  ileng = parent_prc_nprocs(handle) * 4 * 2
1143  if ( prc_ismaster ) then
1144  call mpi_irecv(latlon_catalog, ileng, comm_datatype, prc_myrank, tag, intercomm_parent, ireq, ierr)
1145  call mpi_wait(ireq, istatus, ierr)
1146  endif
1147  call comm_bcast( latlon_catalog, parent_prc_nprocs(handle), 4, 2 )
1148  else
1149  if( io_l ) write(*,*) 'xxx internal error [nest/grid]'
1150  call prc_mpistop
1151  endif
1152 
1153  return
1154  end subroutine nest_comm_catalogue
1155 
1156  !-----------------------------------------------------------------------------
1158  subroutine nest_comm_ping( &
1159  HANDLE )
1160  use scale_process, only: &
1161  prc_mpistop, &
1162  prc_myrank, &
1163  prc_ismaster
1164  use scale_comm, only: &
1165  comm_bcast
1166  implicit none
1167 
1168  integer, intent(in) :: HANDLE
1169 
1170  integer :: ping, pong
1171  integer :: ireq1, ireq2, ierr1, ierr2
1172  integer :: istatus(mpi_status_size)
1173  integer :: tag
1174  logical :: ping_error
1175  !---------------------------------------------------------------------------
1176 
1177  tag = intercomm_id(handle) * 100
1178  ping_error = .false.
1179 
1180  if ( nest_filiation( intercomm_id(handle) ) > 0 ) then !--- parent
1181  ping = online_domain_num
1182  pong = 0
1183 
1184  if ( prc_ismaster ) then
1185  call mpi_isend(ping, 1, mpi_integer, prc_myrank, tag+1, intercomm_daughter, ireq1, ierr1)
1186  call mpi_irecv(pong, 1, mpi_integer, prc_myrank, tag+2, intercomm_daughter, ireq2, ierr2)
1187  call mpi_wait(ireq1, istatus, ierr1)
1188  call mpi_wait(ireq2, istatus, ierr2)
1189  endif
1190 
1191  call comm_bcast(pong)
1192 
1193  if ( pong /= intercomm_id(handle)+1 ) ping_error = .true.
1194 
1195  elseif ( nest_filiation( intercomm_id(handle) ) < 0 ) then !--- daughter
1196  ping = online_domain_num
1197  pong = 0
1198 
1199  if ( prc_ismaster ) then
1200  call mpi_isend(ping, 1, mpi_integer, prc_myrank, tag+2, intercomm_parent, ireq1, ierr1)
1201  call mpi_irecv(pong, 1, mpi_integer, prc_myrank, tag+1, intercomm_parent, ireq2, ierr2)
1202  call mpi_wait(ireq1, istatus, ierr1)
1203  call mpi_wait(ireq2, istatus, ierr2)
1204  endif
1205 
1206  call comm_bcast(pong)
1207 
1208  if ( pong /= intercomm_id(handle) ) ping_error = .true.
1209 
1210  else
1211  if( io_l ) write(*,*) 'xxx internal error [nest/grid]'
1212  call prc_mpistop
1213  endif
1214 
1215  if ( ping_error ) then
1216  if( io_l ) write(*,*) 'xxx ping destination error [nest/grid]'
1217  call prc_mpistop
1218  endif
1219 
1220  return
1221  end subroutine nest_comm_ping
1222 
1223  !-----------------------------------------------------------------------------
1225  subroutine nest_comm_setup_nestdown( &
1226  HANDLE )
1227  use scale_process, only: &
1228  prc_mpistop, &
1229  prc_myrank, &
1230  prc_ismaster
1231  use scale_grid_real, only: &
1233  use scale_comm, only: &
1234  comm_world, &
1235  comm_bcast
1236  implicit none
1237 
1238  integer, intent(in) :: HANDLE
1239 
1240  integer, allocatable :: buffer_LIST(:)
1241  integer, allocatable :: buffer_ALLLIST(:)
1242 
1243  integer :: ireq, ierr, ileng
1244  integer :: istatus(mpi_status_size)
1245  integer :: tag, target_rank
1246 
1247  integer :: i, j, k
1248  !---------------------------------------------------------------------------
1249 
1250  tag = intercomm_id(handle) * 100
1251 
1252  if ( nest_filiation( intercomm_id(handle) ) > 0 ) then
1253  !--------------------------------------------------- parent
1254 
1255  if ( prc_ismaster ) then
1256  call mpi_irecv(nest_tile_allmax_p, 1, mpi_integer, prc_myrank, tag+1, intercomm_daughter, ireq, ierr)
1257  call mpi_wait(ireq, istatus, ierr)
1258  endif
1259  call comm_bcast(nest_tile_allmax_p)
1260 
1261  allocate( nest_tile_list_p(nest_tile_allmax_p,daughter_prc_nprocs(handle)) )
1262  allocate( nest_tile_list_yp(nest_tile_allmax_p*daughter_prc_nprocs(handle)) )
1263 
1264  ileng = nest_tile_allmax_p*daughter_prc_nprocs(handle)
1265  if ( prc_ismaster ) then
1266  call mpi_irecv(nest_tile_list_p, ileng, mpi_integer, prc_myrank, tag+2, intercomm_daughter, ireq, ierr)
1267  call mpi_wait(ireq, istatus, ierr)
1268  endif
1269  call comm_bcast(nest_tile_list_p, nest_tile_allmax_p, daughter_prc_nprocs(handle))
1270 
1271  nest_tile_list_yp(:) = -1
1272 
1273  k = 0
1274  do j = 1, daughter_prc_nprocs(handle)
1275  do i = 1, nest_tile_allmax_p
1276  if ( nest_tile_list_p(i,j) == prc_myrank ) then
1277  k = k + 1
1278  nest_tile_list_yp(k) = j - 1 !rank number is started from 1
1279  endif
1280  enddo
1281  enddo
1282  num_yp = k
1283 
1284  if( io_l ) write(io_fid_log,'(A,I5,A,I5)') "[P] Num YP =",num_yp," Num TILE(MAX) =",nest_tile_allmax_p
1285 
1286  if ( prc_ismaster ) then
1287  call mpi_irecv(online_daughter_use_velz, 1, mpi_logical, prc_myrank, tag+3, intercomm_daughter, ireq, ierr)
1288  call mpi_wait(ireq, istatus, ierr)
1289  endif
1290  call comm_bcast(online_daughter_use_velz)
1291 
1292  if( io_l ) write(io_fid_log,'(1x,A,L2)') '*** NEST: ONLINE_DAUGHTER_USE_VELZ =', online_daughter_use_velz
1293 
1294  if ( prc_ismaster ) then
1295  call mpi_irecv(online_daughter_no_rotate, 1, mpi_logical, prc_myrank, tag+4, intercomm_daughter, ireq, ierr)
1296  call mpi_wait(ireq, istatus, ierr)
1297  endif
1298  call comm_bcast(online_daughter_no_rotate)
1299 
1300  if( online_no_rotate .neqv. online_daughter_no_rotate ) then
1301  write(*,*) 'xxx Flag of NO_ROTATE is not consistent with the child domain'
1302  if( io_l ) write(io_fid_log,*) 'xxx ONLINE_NO_ROTATE = ', online_no_rotate
1303  if( io_l ) write(io_fid_log,*) 'xxx ONLINE_DAUGHTER_NO_ROTATE =', online_daughter_no_rotate
1304  call prc_mpistop
1305  endif
1306  if( io_l ) write(io_fid_log,'(1x,A,L2)') '*** NEST: ONLINE_DAUGHTER_NO_ROTATE =', online_daughter_no_rotate
1307 
1308  call nest_comm_importgrid_nestdown( handle )
1309 
1310  do i = 1, num_yp
1311  target_rank = nest_tile_list_yp(i)
1312  call mpi_isend(i, 1, mpi_integer, target_rank, tag+5, intercomm_daughter, ireq, ierr)
1313  call mpi_wait(ireq, istatus, ierr)
1314  enddo
1315 
1316  call mpi_barrier(intercomm_daughter, ierr)
1317 
1318  elseif ( nest_filiation( intercomm_id(handle) ) < 0 ) then
1319  !--------------------------------------------------- daughter
1320 
1321  nest_tile_all = size( nest_tile_id(:) ) ! should be equal to "NEST_TILE_NUM_X*NEST_TILE_NUM_Y"
1322  call mpi_allreduce( nest_tile_all, &
1323  nest_tile_allmax_d, &
1324  1, &
1325  mpi_integer, &
1326  mpi_max, &
1327  comm_world, &
1328  ierr )
1329  if( io_l ) write(io_fid_log,'(A,I5,A,I5)') "[D] Num YP =",nest_tile_all," Num TILE(MAX) =",nest_tile_allmax_d
1330 
1331  if ( prc_ismaster ) then
1332  call mpi_isend(nest_tile_allmax_d, 1, mpi_integer, prc_myrank, tag+1, intercomm_parent, ireq, ierr)
1333  call mpi_wait(ireq, istatus, ierr)
1334  endif
1335 
1336  allocate( buffer_list(nest_tile_allmax_d) )
1337  allocate( buffer_alllist(nest_tile_allmax_d*daughter_prc_nprocs(handle)) )
1338  allocate( nest_tile_list_d(nest_tile_allmax_d,daughter_prc_nprocs(handle)) )
1339 
1340  do i = 1, nest_tile_allmax_d
1341  if ( i <= nest_tile_all ) then
1342  buffer_list(i) = nest_tile_id(i)
1343  else
1344  buffer_list(i) = -1
1345  endif
1346  enddo
1347 
1348  ileng = nest_tile_allmax_d
1349  call mpi_allgather( buffer_list(:), &
1350  ileng, &
1351  mpi_integer, &
1352  buffer_alllist(:), &
1353  ileng, &
1354  mpi_integer, &
1355  comm_world, &
1356  ierr )
1357  k = 1
1358  do j = 1, daughter_prc_nprocs(handle)
1359  do i = 1, nest_tile_allmax_d
1360  nest_tile_list_d(i,j) = buffer_alllist(k)
1361  k = k + 1
1362  enddo
1363  enddo
1364 
1365  deallocate( buffer_list )
1366  deallocate( buffer_alllist )
1367 
1368  ileng = nest_tile_allmax_d*daughter_prc_nprocs(handle)
1369  if ( prc_ismaster ) then
1370  call mpi_isend(nest_tile_list_d, ileng, mpi_integer, prc_myrank, tag+2, intercomm_parent, ireq, ierr)
1371  call mpi_wait(ireq, istatus, ierr)
1372  endif
1373 
1374  if ( prc_ismaster ) then
1375  call mpi_isend(online_use_velz, 1, mpi_logical, prc_myrank, tag+3, intercomm_parent, ireq, ierr)
1376  call mpi_wait(ireq, istatus, ierr)
1377  endif
1378 
1379  if ( prc_ismaster ) then
1380  call mpi_isend(online_no_rotate, 1, mpi_logical, prc_myrank, tag+4, intercomm_parent, ireq, ierr)
1381  call mpi_wait(ireq, istatus, ierr)
1382  endif
1383  call comm_bcast(online_daughter_no_rotate)
1384 
1385  call nest_comm_importgrid_nestdown( handle )
1386 
1387  do i = 1, nest_tile_all
1388  target_rank = nest_tile_list_d(i,prc_myrank+1)
1389  call mpi_irecv( call_order(i), 1, mpi_integer, target_rank, tag+5, intercomm_parent, ireq, ierr )
1390  call mpi_wait(ireq, istatus, ierr)
1391  enddo
1392 
1393  call mpi_barrier(intercomm_parent, ierr)
1394  else
1395  !---------------------------------------------------
1396  if( io_l ) write(*,*) 'xxx internal error [nest/grid]'
1397  call prc_mpistop
1398  endif
1399 
1400  if( num_yp * 16 > max_rq .OR. nest_tile_all * 16 > max_rq ) then ! 16 = dyn:5 + qtrc:11
1401  write(*,*) 'xxx internal error (overflow number of ireq) [nest/grid]'
1402  write(*,*) ' NUM_YP x 16 = ', num_yp * 16
1403  write(*,*) ' NEST_TILE_ALL x 16 = ', nest_tile_all * 16
1404  write(*,*) ' max_rq = ', max_rq
1405  call prc_mpistop
1406  endif
1407 
1408  return
1409  end subroutine nest_comm_setup_nestdown
1410 
1411  !-----------------------------------------------------------------------------
1413  subroutine nest_comm_importgrid_nestdown( &
1414  HANDLE )
1415  use scale_process, only: &
1416  prc_myrank, &
1417  prc_nprocs, &
1418  prc_mpistop
1419  use scale_grid_real, only: &
1420  real_lon, &
1421  real_lat, &
1422  real_lonx, &
1423  real_lony, &
1424  real_latx, &
1425  real_laty, &
1426  real_cz, &
1427  real_fz
1428  use scale_comm, only: &
1430  implicit none
1431 
1432  integer, intent(in) :: HANDLE
1433 
1434  integer :: ierr, ileng
1435  integer :: istatus(mpi_status_size)
1436  integer :: tag, tagbase, target_rank
1437 
1438  integer :: xloc, yloc
1439  integer :: xs, xe
1440  integer :: ys, ye
1441 
1442  real(RP) :: max_ref, max_loc
1443 
1444  integer :: i, k, rq
1445  !---------------------------------------------------------------------------
1446 
1447  tagbase = intercomm_id(handle) * 100
1448  rq = 0
1449 
1450  if ( nest_filiation( intercomm_id(handle) ) > 0 ) then
1451  !--------------------------------------------------- parent
1452  do i = 1, num_yp
1453  ! send data to multiple daughter processes
1454  target_rank = nest_tile_list_yp(i)
1455 
1456  rq = rq + 1
1457  ileng = parent_ia(handle) * parent_ja(handle)
1458  tag = tagbase + tag_lon
1459  call mpi_isend(real_lon, ileng, comm_datatype, target_rank, tag, intercomm_daughter, ireq_p(rq), ierr)
1460  call mpi_wait(ireq_p(rq), istatus, ierr)
1461 
1462  rq = rq + 1
1463  ileng = parent_ia(handle) * parent_ja(handle)
1464  tag = tagbase + tag_lat
1465  call mpi_isend(real_lat, ileng, comm_datatype, target_rank, tag, intercomm_daughter, ireq_p(rq), ierr)
1466  call mpi_wait(ireq_p(rq), istatus, ierr)
1467 
1468  rq = rq + 1
1469  ileng = parent_ia(handle) * parent_ja(handle)
1470  tag = tagbase + tag_lonx
1471  call mpi_isend(real_lonx, ileng, comm_datatype, target_rank, tag, intercomm_daughter, ireq_p(rq), ierr)
1472  call mpi_wait(ireq_p(rq), istatus, ierr)
1473 
1474  rq = rq + 1
1475  ileng = parent_ia(handle) * parent_ja(handle)
1476  tag = tagbase + tag_latx
1477  call mpi_isend(real_latx, ileng, comm_datatype, target_rank, tag, intercomm_daughter, ireq_p(rq), ierr)
1478  call mpi_wait(ireq_p(rq), istatus, ierr)
1479 
1480  rq = rq + 1
1481  ileng = parent_ia(handle) * parent_ja(handle)
1482  tag = tagbase + tag_lony
1483  call mpi_isend(real_lony, ileng, comm_datatype, target_rank, tag, intercomm_daughter, ireq_p(rq), ierr)
1484  call mpi_wait(ireq_p(rq), istatus, ierr)
1485 
1486  rq = rq + 1
1487  ileng = parent_ia(handle) * parent_ja(handle)
1488  tag = tagbase + tag_laty
1489  call mpi_isend(real_laty, ileng, comm_datatype, target_rank, tag, intercomm_daughter, ireq_p(rq), ierr)
1490  call mpi_wait(ireq_p(rq), istatus, ierr)
1491 
1492  rq = rq + 1
1493  ileng = parent_ka(handle) * parent_ia(handle) * parent_ja(handle)
1494  tag = tagbase + tag_cz
1495  call mpi_isend(real_cz, ileng, comm_datatype, target_rank, tag, intercomm_daughter, ireq_p(rq), ierr)
1496  call mpi_wait(ireq_p(rq), istatus, ierr)
1497 
1498  rq = rq + 1
1499  ileng = (parent_ka(handle)+1) * parent_ia(handle) * parent_ja(handle)
1500  tag = tagbase + tag_fz
1501  call mpi_isend(real_fz, ileng, comm_datatype, target_rank, tag, intercomm_daughter, ireq_p(rq), ierr)
1502  call mpi_wait(ireq_p(rq), istatus, ierr)
1503  enddo
1504 
1505 
1506  elseif ( nest_filiation( intercomm_id(handle) ) < 0 ) then
1507  !--------------------------------------------------- daughter
1508  do i = 1, nest_tile_all
1509  ! receive data from multiple parent tiles
1510  target_rank = nest_tile_list_d(i,prc_myrank+1)
1511 
1512  xloc = mod( i-1, nest_tile_num_x ) + 1
1513  yloc = int( real(i-1) / real(NEST_TILE_NUM_X) ) + 1
1514 
1515  xs = parent_imax(handle) * (xloc-1) + 1
1516  xe = parent_imax(handle) * xloc
1517  ys = parent_jmax(handle) * (yloc-1) + 1
1518  ye = parent_jmax(handle) * yloc
1519 
1520  rq = rq + 1
1521  ileng = parent_ia(handle) * parent_ja(handle)
1522  tag = tagbase + tag_lon
1523  call mpi_irecv(buffer_2d, ileng, comm_datatype, target_rank, tag, intercomm_parent, ireq_d(rq), ierr)
1524  call mpi_wait(ireq_d(rq), istatus, ierr)
1525  buffer_ref_lon(xs:xe,ys:ye) = buffer_2d(prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle))
1526 
1527  rq = rq + 1
1528  ileng = parent_ia(handle) * parent_ja(handle)
1529  tag = tagbase + tag_lat
1530  call mpi_irecv(buffer_2d, ileng, comm_datatype, target_rank, tag, intercomm_parent, ireq_d(rq), ierr)
1531  call mpi_wait(ireq_d(rq), istatus, ierr)
1532  buffer_ref_lat(xs:xe,ys:ye) = buffer_2d(prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle))
1533 
1534  rq = rq + 1
1535  ileng = parent_ia(handle) * parent_ja(handle)
1536  tag = tagbase + tag_lonx
1537  call mpi_irecv(buffer_2d, ileng, comm_datatype, target_rank, tag, intercomm_parent, ireq_d(rq), ierr)
1538  call mpi_wait(ireq_d(rq), istatus, ierr)
1539  buffer_ref_lonx(xs:xe,ys:ye) = buffer_2d(prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle))
1540 
1541  rq = rq + 1
1542  ileng = parent_ia(handle) * parent_ja(handle)
1543  tag = tagbase + tag_latx
1544  call mpi_irecv(buffer_2d, ileng, comm_datatype, target_rank, tag, intercomm_parent, ireq_d(rq), ierr)
1545  call mpi_wait(ireq_d(rq), istatus, ierr)
1546  buffer_ref_latx(xs:xe,ys:ye) = buffer_2d(prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle))
1547 
1548  rq = rq + 1
1549  ileng = parent_ia(handle) * parent_ja(handle)
1550  tag = tagbase + tag_lony
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_lony(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_laty
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_laty(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_ka(handle) * parent_ia(handle) * parent_ja(handle)
1564  tag = tagbase + tag_cz
1565  call mpi_irecv(buffer_3d, ileng, comm_datatype, target_rank, tag, intercomm_parent, ireq_d(rq), ierr)
1566  call mpi_wait(ireq_d(rq), istatus, ierr)
1567  do k = 1, parent_ka(handle)
1568  buffer_ref_cz(k,xs:xe,ys:ye) = buffer_3d(k,prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle))
1569  enddo
1570 
1571  rq = rq + 1
1572  ileng = (parent_ka(handle)+1) * parent_ia(handle) * parent_ja(handle)
1573  tag = tagbase + tag_fz
1574  call mpi_irecv(buffer_3df,ileng, comm_datatype, target_rank, tag, intercomm_parent, ireq_d(rq), ierr)
1575  call mpi_wait(ireq_d(rq), istatus, ierr)
1576  do k = 0, parent_ka(handle)
1577  buffer_ref_fz(k,xs:xe,ys:ye) = buffer_3df(k,prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle))
1578  enddo
1579  enddo
1580 
1581  ! check domain compatibility
1582  max_ref = maxval( buffer_ref_fz(:,:,:) )
1583  max_loc = maxval( real_fz(ks-1:ke,:,:) ) ! HALO + 1
1584  if ( max_ref < max_loc ) then
1585  write(*,*) 'xxx ERROR: REQUESTED DOMAIN IS TOO MUCH BROAD'
1586  write(*,*) 'xxx -- VERTICAL direction over the limit'
1587  write(*,*) 'xxx -- reference max: ', max_ref
1588  write(*,*) 'xxx -- local max: ', max_loc
1589  call prc_mpistop
1590  endif
1591 
1592  else
1593  !---------------------------------------------------
1594  write(*,*) 'xxx internal error [nest/grid]'
1595  call prc_mpistop
1596  endif
1597 
1598  return
1599  end subroutine nest_comm_importgrid_nestdown
1600 
1601  !-----------------------------------------------------------------------------
1603  subroutine nest_comm_nestdown( &
1604  HANDLE, & ! [in ]
1605  BND_QA, & ! [in ]
1606  ipt_DENS, & ! [in ]
1607  ipt_MOMZ, & ! [in ]
1608  ipt_MOMX, & ! [in ]
1609  ipt_MOMY, & ! [in ]
1610  ipt_RHOT, & ! [in ]
1611  ipt_QTRC, & ! [in ]
1612  interped_ref_DENS, & ! [inout]
1613  interped_ref_VELZ, & ! [inout]
1614  interped_ref_VELX, & ! [inout]
1615  interped_ref_VELY, & ! [inout]
1616  interped_ref_POTT, & ! [inout]
1617  interped_ref_QTRC ) ! [inout]
1618  use scale_process, only: &
1619  prc_myrank, &
1620  prc_nprocs, &
1621  prc_mpistop
1622  use scale_grid_real, only: &
1624  use scale_comm, only: &
1625  comm_vars8, &
1626  comm_wait
1627  use scale_gridtrans, only: &
1628  rotc => gtrans_rotc
1629  implicit none
1630 
1631  integer, intent(in) :: HANDLE
1632  integer, intent(in) :: BND_QA
1633 
1634  real(RP), intent(in ) :: ipt_DENS(parent_ka(handle),parent_ia(handle),parent_ja(handle))
1635  real(RP), intent(in ) :: ipt_MOMZ(parent_ka(handle),parent_ia(handle),parent_ja(handle))
1636  real(RP), intent(in ) :: ipt_MOMX(parent_ka(handle),parent_ia(handle),parent_ja(handle))
1637  real(RP), intent(in ) :: ipt_MOMY(parent_ka(handle),parent_ia(handle),parent_ja(handle))
1638  real(RP), intent(in ) :: ipt_RHOT(parent_ka(handle),parent_ia(handle),parent_ja(handle))
1639  real(RP), intent(in ) :: ipt_QTRC(parent_ka(handle),parent_ia(handle),parent_ja(handle),bnd_qa)
1640 
1641  real(RP), intent(inout) :: interped_ref_DENS(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1642  real(RP), intent(inout) :: interped_ref_VELZ(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1643  real(RP), intent(inout) :: interped_ref_VELX(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1644  real(RP), intent(inout) :: interped_ref_VELY(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1645  real(RP), intent(inout) :: interped_ref_POTT(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1646  real(RP), intent(inout) :: interped_ref_QTRC(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle),bnd_qa)
1647 
1648  real(RP) :: dummy(1,1,1)
1649  real(RP) :: dens (daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1650  real(RP) :: u_lld(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1651  real(RP) :: v_lld(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1652  real(RP) :: work1d(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1653  real(RP) :: work2d(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1654 
1655  real(RP) :: work1p(parent_ka(handle),parent_ia(handle),parent_ja(handle))
1656  real(RP) :: work2p(parent_ka(handle),parent_ia(handle),parent_ja(handle))
1657 
1658  real(RP) :: u_on_map, v_on_map
1659  integer :: ierr
1660  integer :: tagbase, tagcomm
1661  integer :: isu_tag, isu_tagf
1662  integer :: i, j, k, iq
1663 
1664  integer, parameter :: cosin = 1
1665  integer, parameter :: sine = 2
1666  !---------------------------------------------------------------------------
1667 
1668  if ( bnd_qa > i_bndqa ) then
1669  if( io_l ) write(*,*) 'xxx internal error: BND_QA is larger than I_BNDQA [nest/grid]'
1670  call prc_mpistop
1671  endif
1672  if ( bnd_qa > max_bndqa ) then
1673  if( io_l ) write(*,*) 'xxx internal error: BND_QA is larger than max_bndqa [nest/grid]'
1674  call prc_mpistop
1675  endif
1676 
1677  tagcomm = intercomm_id(handle) * order_tag_comm
1678 
1679  if ( nest_filiation( intercomm_id(handle) ) > 0 ) then
1680  !-------------------------------------------------------- parent [send issue]
1681  call prof_rapstart('NEST_total_P', 2)
1682 
1683  nsend = nsend + 1
1684  if( io_l ) write(io_fid_log,'(1X,A,I5,A)') "*** NestIDC [P]: que send ( ", nsend, " )"
1685 
1686  call prof_rapstart('NEST_pack_P', 2)
1687 
1688  ! to keep values at that time by finish of sending process
1689  org_dens(:,:,:) = ipt_dens(:,:,:)
1690  org_momz(:,:,:) = ipt_momz(:,:,:)
1691  org_momx(:,:,:) = ipt_momx(:,:,:)
1692  org_momy(:,:,:) = ipt_momy(:,:,:)
1693  org_rhot(:,:,:) = ipt_rhot(:,:,:)
1694  do iq = 1, bnd_qa
1695  org_qtrc(:,:,:,iq) = ipt_qtrc(:,:,:,iq)
1696  enddo
1697 
1698  !*** request control
1699  !--- do not change the calling order below;
1700  !--- it should be consistent with the order in "NEST_COMM_recvwait_issue"
1701  rq_ctl_p = 0
1702 
1703  if ( .NOT. online_daughter_no_rotate ) then
1704  ! from staggered point to scalar point
1705  do j = 1, parent_ja(handle)
1706  do i = 2, parent_ia(handle)
1707  do k = prnt_ks(handle)-1, prnt_ke(handle)+1
1708  work1p(k,i,j) = ( org_momx(k,i-1,j) + org_momx(k,i,j) ) * 0.5_rp
1709  end do
1710  end do
1711  end do
1712  do j = 1, parent_ja(handle)
1713  do k = prnt_ks(handle)-1, prnt_ke(handle)+1
1714  work1p(k,1,j) = org_momx(k,1,j)
1715  end do
1716  end do
1717  call comm_vars8( work1p(:,:,:), 1 )
1718  do j = 2, parent_ja(handle)
1719  do i = 1, parent_ia(handle)
1720  do k = prnt_ks(handle)-1, prnt_ke(handle)+1
1721  work2p(k,i,j) = ( org_momy(k,i,j-1) + org_momy(k,i,j) ) * 0.5_rp
1722  end do
1723  end do
1724  end do
1725  do i = 1, parent_ia(handle)
1726  do k = prnt_ks(handle)-1, prnt_ke(handle)+1
1727  work2p(k,i,1) = org_momy(k,i,1)
1728  end do
1729  end do
1730  call comm_vars8( work2p(:,:,:), 2 )
1731  call comm_wait ( work1p(:,:,:), 1, .false. )
1732  call comm_wait ( work2p(:,:,:), 2, .false. )
1733 
1734  ! rotation from map-projected field to latlon field
1735  do j = prnt_js(handle), prnt_je(handle)
1736  do i = prnt_is(handle), prnt_ie(handle)
1737  do k = prnt_ks(handle), prnt_ke(handle)
1738  u_on_map = work1p(k,i,j) / org_dens(k,i,j)
1739  v_on_map = work2p(k,i,j) / org_dens(k,i,j)
1740 
1741  u_llp(k,i,j) = u_on_map * rotc(i,j,cosin) - v_on_map * rotc(i,j,sine )
1742  v_llp(k,i,j) = u_on_map * rotc(i,j,sine ) + v_on_map * rotc(i,j,cosin)
1743  enddo
1744  enddo
1745  enddo
1746  end if
1747 
1748  tagbase = tagcomm + tag_dens*order_tag_var
1749  call nest_comm_intercomm_nestdown( org_dens, dummy, tagbase, i_sclr, handle, isu_tag, isu_tagf, .true. )
1750 
1751  tagbase = tagcomm + tag_momz*order_tag_var
1752  if ( online_daughter_use_velz ) then
1753  call nest_comm_intercomm_nestdown( org_momz, dummy, tagbase, i_zstg, handle, isu_tag, isu_tagf )
1754  end if
1755 
1756  tagbase = tagcomm + tag_momx*order_tag_var
1757  if ( online_daughter_no_rotate ) then
1758  call nest_comm_intercomm_nestdown( org_momx, dummy, tagbase, i_xstg, handle, isu_tag, isu_tagf )
1759  else
1760  call nest_comm_intercomm_nestdown( u_llp, dummy, tagbase, i_sclr, handle, isu_tag, isu_tagf )
1761  end if
1762 
1763  tagbase = tagcomm + tag_momy*order_tag_var
1764  if ( online_daughter_no_rotate ) then
1765  call nest_comm_intercomm_nestdown( org_momy, dummy, tagbase, i_ystg, handle, isu_tag, isu_tagf )
1766  else
1767  call nest_comm_intercomm_nestdown( v_llp, dummy, tagbase, i_sclr, handle, isu_tag, isu_tagf )
1768  end if
1769 
1770  tagbase = tagcomm + tag_rhot*order_tag_var
1771  call nest_comm_intercomm_nestdown( org_rhot, dummy, tagbase, i_sclr, handle, isu_tag, isu_tagf )
1772 
1773  do iq = 1, bnd_qa
1774  tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
1775  call nest_comm_intercomm_nestdown( org_qtrc(:,:,:,iq), dummy, tagbase, i_sclr, handle, isu_tag, isu_tagf )
1776  enddo
1777 
1778  rq_tot_p = rq_ctl_p
1779 
1780  call prof_rapend ('NEST_pack_P', 2)
1781 
1782  call prof_rapend ('NEST_total_P', 2)
1783 
1784  elseif ( nest_filiation( intercomm_id(handle) ) < 0 ) then
1785  !-------------------------------------------------------- daughter [wait issue]
1786  call prof_rapstart('NEST_total_C', 2)
1787 
1788  nwait_d = nwait_d + 1
1789  !if( IO_L ) write(IO_FID_LOG,'(1X,A,I5,A)') "*** NestIDC [C]: que wait ( ", nwait_d, " )"
1790 
1791  !*** reset issue tag and request control
1792  !--- do not change the calling order below;
1793  !--- it should be consistent with the order in "NEST_COMM_recvwait_issue"
1794  isu_tag = 0
1795  isu_tagf = 0
1796 
1797  call prof_rapstart('NEST_wait_C', 2)
1798  call nest_comm_waitall( rq_tot_d, ireq_d )
1799 
1800  if ( online_aggressive_comm ) then
1801  ! nothing to do
1802  else
1803  call mpi_barrier(intercomm_parent, ierr)
1804  endif
1805  call prof_rapend ('NEST_wait_C', 2)
1806 
1807  call prof_rapstart('NEST_unpack_C', 2)
1808 
1809  tagbase = tagcomm + tag_dens*order_tag_var
1810  call nest_comm_intercomm_nestdown( dummy, dens, tagbase, i_sclr, handle, isu_tag, isu_tagf, .true. )
1811 !OCL XFILL
1812  do j = 1, daughter_ja(handle)
1813  do i = 1, daughter_ia(handle)
1814  do k = datr_ks(handle), datr_ke(handle)
1815  interped_ref_dens(k,i,j) = dens(k,i,j)
1816  enddo
1817  enddo
1818  enddo
1819 
1820  call comm_vars8( interped_ref_dens, 1 )
1821 
1822  tagbase = tagcomm + tag_momz*order_tag_var
1823  if ( online_use_velz ) then
1824  call nest_comm_intercomm_nestdown( dummy, work1d, tagbase, i_zstg, handle, isu_tag, isu_tagf )
1825 !OCL XFILL
1826  do j = 1, daughter_ja(handle)
1827  do i = 1, daughter_ia(handle)
1828  do k = datr_ks(handle), datr_ke(handle)-1
1829  interped_ref_velz(k,i,j) = work1d(k,i,j) / ( dens(k,i,j) + dens(k+1,i,j) ) * 2.0_rp
1830  enddo
1831  enddo
1832  enddo
1833  end if
1834 
1835  tagbase = tagcomm + tag_momx*order_tag_var
1836  if ( online_no_rotate ) then
1837  ! u_lld receives MOMX
1838  call nest_comm_intercomm_nestdown( dummy, u_lld, tagbase, i_xstg, handle, isu_tag, isu_tagf )
1839  else
1840  ! u_lld receives MOMX/DENS
1841  call nest_comm_intercomm_nestdown( dummy, u_lld, tagbase, i_sclr, handle, isu_tag, isu_tagf )
1842  endif
1843 
1844  tagbase = tagcomm + tag_momy*order_tag_var
1845  if ( online_no_rotate ) then
1846  ! v_lld receives MOMY
1847  call nest_comm_intercomm_nestdown( dummy, v_lld, tagbase, i_ystg, handle, isu_tag, isu_tagf )
1848  else
1849  ! v_lld receives MOMY/DENS
1850  call nest_comm_intercomm_nestdown( dummy, v_lld, tagbase, i_sclr, handle, isu_tag, isu_tagf )
1851  endif
1852 
1853  call comm_wait ( interped_ref_dens, 1, .false. )
1854 
1855  if ( online_no_rotate ) then
1856  do j = 1, daughter_ja(handle)
1857  do i = 1, daughter_ia(handle)-1
1858  do k = datr_ks(handle), datr_ke(handle)
1859  interped_ref_velx(k,i,j) = u_lld(k,i,j) &
1860  / ( interped_ref_dens(k,i+1,j) + interped_ref_dens(k,i,j) ) * 2.0_rp
1861  enddo
1862  enddo
1863  enddo
1864  do j = 1, daughter_ja(handle)
1865  do k = datr_ks(handle), datr_ke(handle)
1866  interped_ref_velx(k,daughter_ia(handle),j) = u_lld(k,daughter_ia(handle),j) &
1867  / interped_ref_dens(k,daughter_ia(handle),j)
1868  enddo
1869  enddo
1870  call comm_vars8( interped_ref_velx, 2 )
1871  do j = 1, daughter_ja(handle)-1
1872  do i = 1, daughter_ia(handle)
1873  do k = datr_ks(handle), datr_ke(handle)
1874  interped_ref_vely(k,i,j) = v_lld(k,i,j) &
1875  / ( interped_ref_dens(k,i,j+1) + interped_ref_dens(k,i,j) ) * 2.0_rp
1876  enddo
1877  enddo
1878  enddo
1879  do i = 1, daughter_ia(handle)
1880  do k = datr_ks(handle), datr_ke(handle)
1881  interped_ref_vely(k,i,daughter_ja(handle)) = v_lld(k,i,daughter_ja(handle)) &
1882  / interped_ref_dens(k,i,daughter_ja(handle))
1883  enddo
1884  enddo
1885  call comm_vars8( interped_ref_vely, 3 )
1886  call comm_wait ( interped_ref_velx, 2, .false. )
1887  call comm_wait ( interped_ref_vely, 3, .false. )
1888  else ! rotate
1889  ! rotation from latlon field to map-projected field
1890  do j = 1, daughter_ja(handle)
1891  do i = 1, daughter_ia(handle)
1892  do k = datr_ks(handle), datr_ke(handle)
1893  work1d(k,i,j) = u_lld(k,i,j) * rotc(i,j,cosin) + v_lld(k,i,j) * rotc(i,j,sine )
1894  work2d(k,i,j) = - u_lld(k,i,j) * rotc(i,j,sine ) + v_lld(k,i,j) * rotc(i,j,cosin)
1895  enddo
1896  enddo
1897  enddo
1898 
1899  ! from scalar point to staggered point
1900  do j = 1, daughter_ja(handle)
1901  do i = 1, daughter_ia(handle)-1
1902  do k = datr_ks(handle), datr_ke(handle)
1903  interped_ref_velx(k,i,j) = ( work1d(k,i+1,j) + work1d(k,i,j) ) * 0.5_rp
1904  end do
1905  end do
1906  end do
1907  do j = 1, daughter_ja(handle)
1908  do k = datr_ks(handle), datr_ke(handle)
1909  interped_ref_velx(k,daughter_ia(handle),j) = work1d(k,daughter_ia(handle),j)
1910  enddo
1911  enddo
1912  call comm_vars8( interped_ref_velx, 2 )
1913  do j = 1, daughter_ja(handle)-1
1914  do i = 1, daughter_ia(handle)
1915  do k = datr_ks(handle), datr_ke(handle)
1916  interped_ref_vely(k,i,j) = ( work2d(k,i,j+1) + work2d(k,i,j) ) * 0.5_rp
1917  end do
1918  end do
1919  end do
1920  do i = 1, daughter_ia(handle)
1921  do k = datr_ks(handle), datr_ke(handle)
1922  interped_ref_vely(k,i,daughter_ja(handle)) = work2d(k,i,daughter_ja(handle))
1923  enddo
1924  enddo
1925  call comm_vars8( interped_ref_vely, 3 )
1926  call comm_wait ( interped_ref_velx, 2, .false. )
1927  call comm_wait ( interped_ref_vely, 3, .false. )
1928  end if
1929 
1930  tagbase = tagcomm + tag_rhot*order_tag_var
1931  call nest_comm_intercomm_nestdown( dummy, work1d, tagbase, i_sclr, handle, isu_tag, isu_tagf )
1932  do j = 1, daughter_ja(handle)
1933  do i = 1, daughter_ia(handle)
1934  do k = datr_ks(handle), datr_ke(handle)
1935  interped_ref_pott(k,i,j) = work1d(k,i,j) / interped_ref_dens(k,i,j)
1936  enddo
1937  enddo
1938  enddo
1939 
1940  do iq = 1, bnd_qa
1941  tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
1942  call nest_comm_intercomm_nestdown( dummy, work1d, tagbase, i_sclr, handle, isu_tag, isu_tagf )
1943  do j = 1, daughter_ja(handle)
1944  do i = 1, daughter_ia(handle)
1945  do k = datr_ks(handle), datr_ke(handle)
1946  interped_ref_qtrc(k,i,j,iq) = work1d(k,i,j)
1947  enddo
1948  enddo
1949  enddo
1950  enddo
1951 
1952  call prof_rapend ('NEST_unpack_C', 2)
1953 
1954  call prof_rapend ('NEST_total_C', 2)
1955  else
1956  write(*,*) 'xxx internal error [nestdown: nest/grid]'
1957  call prc_mpistop
1958  endif
1959 
1960  return
1961  end subroutine nest_comm_nestdown
1962 
1963  !-----------------------------------------------------------------------------
1965  subroutine nest_comm_recvwait_issue( &
1966  HANDLE, & ! [in]
1967  BND_QA ) ! [in]
1968  use scale_process, only: &
1969  prc_mpistop
1970  implicit none
1971 
1972  integer, intent(in) :: HANDLE
1973  integer, intent(in) :: BND_QA
1974 
1975  integer :: isu_tag, isu_tagf
1976  integer :: tagbase, tagcomm
1977  integer :: ierr
1978  integer :: iq
1979  !---------------------------------------------------------------------------
1980 
1981  if ( bnd_qa > i_bndqa ) then
1982  write(*,*) 'xxx internal error: about BND_QA [nest/grid]'
1983  call prc_mpistop
1984  endif
1985 
1986  tagcomm = intercomm_id(handle) * order_tag_comm
1987 
1988  if ( nest_filiation( intercomm_id(handle) ) > 0 ) then
1989  !-------------------------------------------------------- parent [wait issue]
1990  call prof_rapstart('NEST_total_P', 2)
1991  nwait_p = nwait_p + 1
1992  !if( IO_L ) write(IO_FID_LOG,'(1X,A,I5,A)') "*** NestIDC [P]: que wait ( ", nwait_p, " )"
1993 
1994  call prof_rapstart('NEST_wait_P', 2)
1995  call nest_comm_issuer_of_wait( handle )
1996 
1997  if ( online_aggressive_comm ) then
1998  ! nothing to do
1999  else
2000  call mpi_barrier(intercomm_daughter, ierr)
2001  endif
2002 
2003  call prof_rapend ('NEST_wait_P', 2)
2004 
2005  call prof_rapend('NEST_total_P', 2)
2006  elseif ( nest_filiation( intercomm_id(handle) ) < 0 ) then
2007  !-------------------------------------------------------- daughter [receive issue]
2008  call prof_rapstart('NEST_total_C', 2)
2009  nrecv = nrecv + 1
2010  if( io_l ) write(io_fid_log,'(1X,A,I5,A)') "*** NestIDC [C]: que recv ( ", nrecv, " )"
2011 
2012  !*** reset issue tag and request control
2013  !--- do not change the calling order below;
2014  !--- it should be consistent with the order in "NEST_COMM_nestdown"
2015  isu_tag = 0
2016  isu_tagf = 0
2017  rq_ctl_d = 0
2018 
2019  tagbase = tagcomm + tag_dens*order_tag_var
2020  call nest_comm_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf, .true. )
2021 
2022  tagbase = tagcomm + tag_momz*order_tag_var
2023  if ( online_use_velz ) then
2024  call nest_comm_issuer_of_receive( tagbase, i_zstg, handle, isu_tag, isu_tagf )
2025  end if
2026 
2027  tagbase = tagcomm + tag_momx*order_tag_var
2028  if ( online_no_rotate ) then
2029  call nest_comm_issuer_of_receive( tagbase, i_xstg, handle, isu_tag, isu_tagf )
2030  else
2031  call nest_comm_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf )
2032  endif
2033 
2034  tagbase = tagcomm + tag_momy*order_tag_var
2035  if ( online_no_rotate ) then
2036  call nest_comm_issuer_of_receive( tagbase, i_ystg, handle, isu_tag, isu_tagf )
2037  else
2038  call nest_comm_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf )
2039  endif
2040 
2041  tagbase = tagcomm + tag_rhot*order_tag_var
2042  call nest_comm_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf )
2043 
2044  do iq = 1, bnd_qa
2045  tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
2046  call nest_comm_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf )
2047  enddo
2048 
2049  rq_tot_d = rq_ctl_d
2050 
2051  call prof_rapend('NEST_total_C', 2)
2052  else
2053  write(*,*) 'xxx internal error [issue: nest/grid]'
2054  call prc_mpistop
2055  endif
2056 
2057  return
2058  end subroutine nest_comm_recvwait_issue
2059 
2060  !-----------------------------------------------------------------------------
2062  subroutine nest_comm_recv_cancel( &
2063  HANDLE ) ! [in]
2064  use scale_process, only: &
2065  prc_myrank, &
2066  prc_nprocs, &
2067  prc_mpistop
2068  implicit none
2069 
2070  integer, intent(in) :: HANDLE
2071 
2072  integer :: ierr
2073  !integer :: istatus(MPI_STATUS_SIZE)
2074 
2075  integer :: rq
2076  logical :: flag
2077  !---------------------------------------------------------------------------
2078 
2079  if ( nest_filiation( intercomm_id(handle) ) > 0 ) then
2080  !-------------------------------------------------------- parent
2081  !--- Nothing to do
2082 
2083  elseif ( nest_filiation( intercomm_id(handle) ) < 0 ) then
2084  !-------------------------------------------------------- daughter [receive issue]
2085  if( io_l ) write(io_fid_log,'(1X,A,I5,A)') "*** NestIDC [C]: CANCEL recv ( ", nrecv, " )"
2086  do rq = 1, rq_tot_d
2087  if ( ireq_d(rq) /= mpi_request_null ) then
2088  call mpi_cancel(ireq_d(rq), ierr)
2089  !call MPI_TEST_CANCELLED(istatus, flag, ierr)
2090  !if ( .NOT. flag ) then
2091  ! write(IO_FID_LOG,*) 'XXX ERROR: receive actions do not cancelled, req = ', rq
2092  !endif
2093  endif
2094  enddo
2095 
2096  else
2097  write(*,*) 'xxx internal error [cancel: nest/grid]'
2098  call prc_mpistop
2099  endif
2100 
2101  return
2102  end subroutine nest_comm_recv_cancel
2103 
2104  !-----------------------------------------------------------------------------
2106  subroutine nest_comm_intercomm_nestdown_3d( &
2107  pvar, & ! [in ]
2108  dvar, & ! [out]
2109  tagbase, & ! [in ]
2110  id_stag, & ! [in ]
2111  HANDLE, & ! [in ]
2112  isu_tag, & ! [inout]
2113  isu_tagf, & ! [inout]
2114  flag_dens ) ! [in ]: optional
2115  use scale_process, only: &
2116  prc_myrank, &
2117  prc_nprocs, &
2118  prc_mpistop
2119  use scale_comm, only: &
2121  use scale_interpolation_nest, only: &
2123  implicit none
2124 
2125  real(RP), intent(in) :: pvar(:,:,:)
2126  real(RP), intent(out) :: dvar(:,:,:)
2127  integer , intent(in) :: tagbase
2128  integer , intent(in) :: id_stag
2129  integer , intent(in) :: HANDLE
2130  integer , intent(inout) :: isu_tag
2131  integer , intent(inout) :: isu_tagf
2132 
2133  logical , intent(in), optional :: flag_dens
2134 
2135  integer :: ierr, ileng
2136  integer :: tag, target_rank
2137 
2138  integer :: xloc, yloc
2139  integer :: xs, xe
2140  integer :: ys, ye
2141 
2142  integer :: i, j, k
2143  integer :: ig, rq, yp
2144  logical :: no_zstag = .true.
2145  logical :: logarithmic = .false.
2146  !---------------------------------------------------------------------------
2147  logarithmic = .false.
2148  if ( present(flag_dens) ) then
2149  if ( flag_dens ) then
2150  logarithmic = .true.
2151  endif
2152  endif
2153 
2154  if ( id_stag == i_sclr ) then
2155  no_zstag = .true.
2156  ig = i_sclr
2157  elseif ( id_stag == i_zstg ) then
2158  no_zstag = .false.
2159  ig = i_zstg
2160  elseif ( id_stag == i_xstg ) then
2161  no_zstag = .true.
2162  ig = i_xstg
2163  elseif ( id_stag == i_ystg ) then
2164  no_zstag = .true.
2165  ig = i_ystg
2166  endif
2167 
2168  if ( no_zstag ) then
2169  ileng = parent_ka(handle) * parent_ia(handle) * parent_ja(handle)
2170  else
2171  ileng = (parent_ka(handle)+1) * parent_ia(handle) * parent_ja(handle)
2172  endif
2173 
2174  if ( nest_filiation( intercomm_id(handle) ) > 0 ) then
2175  !--------------------------------------------------- parent
2176 
2177  rq = rq_ctl_p
2178 
2179  do yp = 1, num_yp
2180  rq = rq + 1
2181 
2182  ! send data to multiple daughter processes
2183  target_rank = nest_tile_list_yp(yp)
2184  tag = tagbase + yp
2185  call mpi_isend(pvar, ileng, comm_datatype, target_rank, tag, intercomm_daughter, ireq_p(rq), ierr)
2186 
2187  dvar(:,:,:) = -1.0_rp ! input as a dummy value
2188  enddo
2189 
2190  rq_ctl_p = rq
2191 
2192  elseif ( nest_filiation( intercomm_id(handle) ) < 0 ) then
2193  !--------------------------------------------------- daughter
2194 
2195  rq = rq_ctl_d
2196 
2197  do yp = 1, nest_tile_all ! YP Loop
2198  rq = rq + 1
2199 
2200  xloc = mod( yp-1, nest_tile_num_x ) + 1
2201  yloc = int( real(yp-1) / real(NEST_TILE_NUM_X) ) + 1
2202 
2203  xs = parent_imax(handle) * (xloc-1) + 1
2204  xe = parent_imax(handle) * xloc
2205  ys = parent_jmax(handle) * (yloc-1) + 1
2206  ye = parent_jmax(handle) * yloc
2207 
2208  if ( no_zstag ) then
2209  isu_tag = isu_tag + 1
2210 
2211  if ( .NOT. logarithmic ) then
2212 !OCL XFILL
2213  ! linear interpolation
2214  do k = 1, parent_ka(handle)
2215  buffer_ref_3d(k,xs:xe,ys:ye) &
2216  = recvbuf_3d(k,prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle),isu_tag)
2217  enddo
2218  else
2219 !OCL XFILL
2220  ! logarithmic weighted interpolation
2221  do k = 1, parent_ka(handle)
2222  buffer_ref_3d(k,xs:xe,ys:ye) &
2223  = log( recvbuf_3d(k,prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle),isu_tag) )
2224  enddo
2225  endif
2226  else
2227  isu_tagf = isu_tagf + 1
2228 !OCL XFILL
2229  do k = 0, parent_ka(handle)
2230  buffer_ref_3df(k,xs:xe,ys:ye) &
2231  = recvbuf_3df(k,prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle),isu_tagf)
2232  enddo
2233  endif
2234 
2235  if ( isu_tag > max_isu .OR. isu_tagf > max_isuf ) then
2236  write(*,*) 'xxx Exceeded maximum issue [intercomm: nest/grid]'
2237  write(*,*) 'xxx isu_tag = ', isu_tag
2238  write(*,*) 'xxx isu_tagf = ', isu_tagf
2239  call prc_mpistop
2240  endif
2241 
2242  enddo ! YP Loop
2243  rq_ctl_d = rq
2244 
2245 !OCL XFILL
2246  dvar(:,:,:) = 0.0_rp
2247 
2248  if ( no_zstag ) then
2249  call intrpnest_interp_3d( dvar, &
2250  buffer_ref_3d, &
2251  hfact(:,:,:,ig), &
2252  vfact(:,:,:,:,:,ig), &
2253  kgrd(:,:,:,:,:,ig), &
2254  igrd(:,:,:,ig), &
2255  jgrd(:,:,:,ig), &
2256  daughter_ia(handle), &
2257  daughter_ja(handle), &
2258  datr_ks(handle), &
2259  datr_ke(handle), &
2260  logarithmic )
2261  else
2262  ! linear interpolation (z-staggered)
2263  call intrpnest_interp_3d( dvar, &
2264  buffer_ref_3df, &
2265  hfact(:,:,:,ig), &
2266  vfact(:,:,:,:,:,ig), &
2267  kgrd(:,:,:,:,:,ig), &
2268  igrd(:,:,:,ig), &
2269  jgrd(:,:,:,ig), &
2270  daughter_ia(handle), &
2271  daughter_ja(handle), &
2272  datr_ks(handle), &
2273  datr_ke(handle), &
2274  logarithmic )
2275  endif
2276  else
2277  !---------------------------------------------------
2278  write(*,*) 'xxx internal error [nest/grid]'
2279  call prc_mpistop
2280  endif
2281 
2282  return
2283  end subroutine nest_comm_intercomm_nestdown_3d
2284 
2285  !-----------------------------------------------------------------------------
2287  subroutine nest_comm_issuer_of_receive_3d( &
2288  tagbase, & ! [in ]
2289  id_stag, & ! [in ]
2290  HANDLE, & ! [in ]
2291  isu_tag, & ! [inout]
2292  isu_tagf, & ! [inout]
2293  flag_dens ) ! [in ]: optional
2294  use scale_process, only: &
2295  prc_myrank, &
2296  prc_nprocs, &
2297  prc_mpistop
2298  use scale_comm, only: &
2300  implicit none
2301 
2302  integer , intent(in) :: tagbase
2303  integer , intent(in) :: id_stag
2304  integer , intent(in) :: HANDLE
2305  integer , intent(inout) :: isu_tag
2306  integer , intent(inout) :: isu_tagf
2307  logical , intent(in), optional :: flag_dens
2308 
2309  integer :: ierr, ileng
2310  integer :: tag, target_rank
2311 
2312  integer :: ig, rq, yp
2313  logical :: no_zstag = .true.
2314  logical :: logarithmic = .false.
2315  !---------------------------------------------------------------------------
2316 
2317  logarithmic = .false.
2318  if ( present(flag_dens) ) then
2319  if ( flag_dens ) then
2320  logarithmic = .true.
2321  endif
2322  endif
2323 
2324  if ( id_stag == i_sclr ) then
2325  no_zstag = .true.
2326  ig = i_sclr
2327  elseif ( id_stag == i_zstg ) then
2328  no_zstag = .false.
2329  ig = i_zstg
2330  elseif ( id_stag == i_xstg ) then
2331  no_zstag = .true.
2332  ig = i_xstg
2333  elseif ( id_stag == i_ystg ) then
2334  no_zstag = .true.
2335  ig = i_ystg
2336  endif
2337 
2338  if ( no_zstag ) then
2339  ileng = parent_ka(handle) * parent_ia(handle) * parent_ja(handle)
2340  else
2341  ileng = (parent_ka(handle)+1) * parent_ia(handle) * parent_ja(handle)
2342  endif
2343 
2344  if ( nest_filiation( intercomm_id(handle) ) > 0 ) then
2345  !--------------------------------------------------- parent
2346  ! nothing to do
2347  ! rq = rq_ctl_p
2348  ! rq_ctl_p = rq
2349  elseif ( nest_filiation( intercomm_id(handle) ) < 0 ) then
2350  !--------------------------------------------------- daughter
2351 
2352  rq = rq_ctl_d
2353 
2354  do yp = 1, nest_tile_all ! YP Loop
2355  rq = rq + 1
2356 
2357  target_rank = nest_tile_list_d(yp,prc_myrank+1)
2358 
2359  tag = tagbase + call_order(yp)
2360  if ( no_zstag ) then
2361  isu_tag = isu_tag + 1
2362  recvbuf_3d(:,:,:,isu_tag) = 0.0_rp
2363  call mpi_irecv( recvbuf_3d(:,:,:,isu_tag), &
2364  ileng, &
2365  comm_datatype, &
2366  target_rank, &
2367  tag, &
2368  intercomm_parent, &
2369  ireq_d(rq), &
2370  ierr )
2371  else
2372  isu_tagf = isu_tagf + 1
2373  recvbuf_3df(:,:,:,isu_tagf) = 0.0_rp
2374  call mpi_irecv( recvbuf_3df(:,:,:,isu_tagf), &
2375  ileng, &
2376  comm_datatype, &
2377  target_rank, &
2378  tag, &
2379  intercomm_parent, &
2380  ireq_d(rq), &
2381  ierr )
2382  endif
2383 
2384  enddo ! YP Loop
2385 
2386  if ( isu_tag > max_isu .OR. isu_tagf > max_isuf ) then
2387  write(*,*) 'xxx Exceeded maximum issue [receive: nest/grid]'
2388  write(*,*) 'xxx isu_tag = ', isu_tag
2389  write(*,*) 'xxx isu_tagf = ', isu_tagf
2390  call prc_mpistop
2391  endif
2392 
2393  rq_ctl_d = rq
2394 
2395  else
2396  !---------------------------------------------------
2397  write(*,*) 'xxx internal error [receive: nest/grid]'
2398  call prc_mpistop
2399  endif
2400 
2401  return
2402  end subroutine nest_comm_issuer_of_receive_3d
2403 
2404  !-----------------------------------------------------------------------------
2406  subroutine nest_comm_issuer_of_wait_3d( &
2407  HANDLE ) ! [in]
2408  use scale_process, only: &
2409  prc_mpistop
2410  implicit none
2411 
2412  integer, intent(in) :: HANDLE
2413  !---------------------------------------------------------------------------
2414 
2415  if ( nest_filiation( intercomm_id(handle) ) > 0 ) then
2416  !--------------------------------------------------- parent
2417 
2418  call nest_comm_waitall( rq_tot_p, ireq_p )
2419 
2420  elseif ( nest_filiation( intercomm_id(handle) ) < 0 ) then
2421  !--------------------------------------------------- daughter
2422  ! nothing to do
2423  else
2424  !---------------------------------------------------
2425  write(*,*) 'xxx internal error [wait: nest/grid]'
2426  call prc_mpistop
2427  endif
2428 
2429  return
2430  end subroutine nest_comm_issuer_of_wait_3d
2431 
2432  !-----------------------------------------------------------------------------
2434  subroutine nest_comm_waitall( &
2435  req_count, & ! [in]
2436  ireq ) ! [in]
2437  use scale_process, only: &
2438  prc_mpistop
2439  implicit none
2440 
2441  integer, intent(in) :: req_count
2442  integer, intent(inout) :: ireq(max_rq)
2443 
2444  integer :: i
2445  integer :: ierr
2446  integer :: istatus(mpi_status_size,req_count)
2447  integer :: req_count2
2448  integer :: ireq2(max_rq)
2449 
2450 ! logical :: flag
2451 ! integer(8) :: num
2452  !---------------------------------------------------------------------------
2453 ! num = 0
2454 ! flag = .false.
2455 
2456  req_count2 = 0
2457  do i=1, req_count
2458  if (ireq(i) /= mpi_request_null) then
2459  req_count2 = req_count2 + 1
2460  ireq2(req_count2) = ireq(i)
2461  endif
2462  enddo
2463  if ( req_count2 /= 0 ) call mpi_waitall( req_count2, ireq2, istatus, ierr )
2464 
2465 ! do while ( .NOT. flag )
2466 ! num = num + 1
2467 ! call MPI_TESTALL( req_count, ireq, flag, istatus, ierr )
2468 
2469 ! if ( num > ONLINE_WAIT_LIMIT ) then
2470 ! if( IO_L ) write(IO_FID_LOG,'(1x,A)') '*** ERROR: over the limit of waiting time [NESTCOM]'
2471 ! write(*,'(1x,A)') '*** ERROR: over the limit of waiting time [NESTCOM]'
2472 ! call PRC_MPIstop
2473 ! endif
2474 ! enddo
2475 
2476  return
2477  end subroutine nest_comm_waitall
2478 
2479  !-----------------------------------------------------------------------------
2481  subroutine nest_comm_test( &
2482  HANDLE ) ! [in]
2483  use scale_process, only: &
2484  prc_mpistop
2485  implicit none
2486 
2487  integer, intent(in) :: HANDLE
2488 
2489  integer :: istatus(mpi_status_size)
2490  integer :: ierr
2491  logical :: flag
2492  !---------------------------------------------------------------------------
2493 
2494  if ( nest_filiation( intercomm_id(handle) ) > 0 ) then
2495  !--------------------------------------------------- parent
2496  call prof_rapstart('NEST_test_P', 2)
2497  if ( rq_ctl_p > 0 ) call mpi_test(ireq_p(1), flag, istatus, ierr)
2498  call prof_rapend('NEST_test_P', 2)
2499 
2500  elseif ( nest_filiation( intercomm_id(handle) ) < 0 ) then
2501  !--------------------------------------------------- daughter
2502  call prof_rapstart('NEST_test_C', 2)
2503  if ( rq_ctl_d > 0 ) call mpi_test(ireq_d(1), flag, istatus, ierr)
2504  call prof_rapend('NEST_test_C', 2)
2505  else
2506  !---------------------------------------------------
2507  write(*,*) 'xxx internal error [test: nest/grid]'
2508  call prc_mpistop
2509  endif
2510 
2511  return
2512  end subroutine nest_comm_test
2513 
2514  !-----------------------------------------------------------------------------
2516  subroutine nest_comm_disconnect ( )
2517  use scale_process, only: &
2519  implicit none
2520 
2521  integer :: ierr
2522  !---------------------------------------------------------------------------
2523 
2524  if( io_l ) write(io_fid_log,'(1x,A)') '*** Waiting finish of whole processes'
2525  call mpi_barrier(prc_global_comm_world, ierr)
2526 
2527  if ( online_iam_parent ) then
2528  !if( IO_L ) write(IO_FID_LOG,'(1x,A)') '*** Waiting finish of whole processes as a parent'
2529  !call MPI_BARRIER(INTERCOMM_DAUGHTER, ierr)
2530  call mpi_comm_free(intercomm_daughter, ierr)
2531  if( io_l ) write(io_fid_log,'(1x,A)') '*** Disconnected communication with child'
2532  endif
2533 
2534  if ( online_iam_daughter ) then
2535  !if( IO_L ) write(IO_FID_LOG,'(1x,A)') '*** Waiting finish of whole processes as a child'
2536  !call MPI_BARRIER(INTERCOMM_PARENT, ierr)
2537  call mpi_comm_free(intercomm_parent, ierr)
2538  if( io_l ) write(io_fid_log,'(1x,A)') '*** Disconnected communication with parent'
2539  endif
2540 
2541  return
2542  end subroutine nest_comm_disconnect
2543 
2544 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)