SCALE-RM
Data Types | Functions/Subroutines | Variables
scale_grid_nest Module Reference

module GRID (nesting system) More...

Functions/Subroutines

subroutine, public nest_setup (inter_parent, inter_child)
 Setup. More...
 
subroutine, public nest_domain_relate (HANDLE)
 Solve relationship between ParentDomain & Daughter Domain. More...
 
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) More...
 
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. More...
 
subroutine, public nest_comm_recvwait_issue (HANDLE, BND_QA)
 Sub-command for data transfer from parent to daughter: nestdown. More...
 
subroutine, public nest_comm_recv_cancel (HANDLE)
 Sub-command for data transfer from parent to daughter: nestdown. More...
 
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. More...
 
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 More...
 
subroutine nest_comm_issuer_of_wait_3d (HANDLE)
 [substance of issuer] Inter-communication from parent to daughter: nestdown More...
 
subroutine nest_comm_waitall (req_count, ireq)
 [substance of comm_wait] Inter-communication More...
 
subroutine, public nest_comm_test (HANDLE)
 [check communication status] Inter-communication More...
 
subroutine, public nest_comm_disconnect ()
 [finalize: disconnect] Inter-communication More...
 

Variables

integer, public intercomm_parent
 
integer, public intercomm_daughter
 
integer, dimension(10), public nest_filiation
 index of parent-daughter relation (p>0, d<0) More...
 
integer, public handling_num
 handing number of nesting relation More...
 
integer, public nest_tile_num_x
 parent tile number in x-direction More...
 
integer, public nest_tile_num_y
 parent tile number in y-direction More...
 
integer, dimension(:), allocatable, public nest_tile_id
 parent tile real id More...
 
integer, dimension(2), public parent_kmax
 parent max number in z-direction More...
 
integer, dimension(2), public parent_imax
 parent max number in x-direction More...
 
integer, dimension(2), public parent_jmax
 parent max number in y-direction More...
 
integer, dimension(2), public parent_ka
 parent max number in z-direction (with halo) More...
 
integer, dimension(2), public parent_ia
 parent max number in x-direction (with halo) More...
 
integer, dimension(2), public parent_ja
 parent max number in y-direction (with halo) More...
 
integer, dimension(2), public parent_lkmax
 parent max number in lz-direction More...
 
real(dp), dimension(2), public parent_dtsec
 parent DT [sec] More...
 
integer, dimension(2), public parent_nstep
 parent step [number] More...
 
integer, dimension(2), public daughter_kmax
 daughter max number in z-direction More...
 
integer, dimension(2), public daughter_imax
 daughter max number in x-direction More...
 
integer, dimension(2), public daughter_jmax
 daughter max number in y-direction More...
 
integer, dimension(2), public daughter_ka
 daughter max number in z-direction (with halo) More...
 
integer, dimension(2), public daughter_ia
 daughter max number in x-direction (with halo) More...
 
integer, dimension(2), public daughter_ja
 daughter max number in y-direction (with halo) More...
 
integer, dimension(2), public daughter_lkmax
 daughter max number in lz-direction More...
 
real(dp), dimension(2), public daughter_dtsec
 daughter DT [sec] More...
 
integer, dimension(2), public daughter_nstep
 daughter steps [number] More...
 
integer, dimension(2), public prnt_ks
 start index in z-direction in parent More...
 
integer, dimension(2), public prnt_ke
 end index in z-direction in parent More...
 
integer, dimension(2), public prnt_is
 start index in x-direction in parent More...
 
integer, dimension(2), public prnt_ie
 end index in x-direction in parent More...
 
integer, dimension(2), public prnt_js
 start index in y-direction in parent More...
 
integer, dimension(2), public prnt_je
 end index in y-direction in parent More...
 
integer, dimension(2), public datr_ks
 start index in z-direction in daughter More...
 
integer, dimension(2), public datr_ke
 end index in z-direction in daughter More...
 
integer, dimension(2), public datr_is
 start index in x-direction in daughter More...
 
integer, dimension(2), public datr_ie
 end index in x-direction in daughter More...
 
integer, dimension(2), public datr_js
 start index in y-direction in daughter More...
 
integer, dimension(2), public datr_je
 end index in y-direction in daughter More...
 
integer, dimension(2), public tileal_ka
 cells of all tiles in z-direction More...
 
integer, dimension(2), public tileal_ia
 cells of all tiles in x-direction More...
 
integer, dimension(2), public tileal_ja
 cells of all tiles in y-direction More...
 
integer, public nest_bnd_qa = 1
 number of tracer treated in nesting system More...
 
integer, public nest_interp_level = 4
 horizontal interpolation level More...
 
logical, public use_nesting = .false.
 
logical, public offline = .true.
 
logical, public online_iam_parent = .false.
 a flag to say "I am a parent" More...
 
logical, public online_iam_daughter = .false.
 a flag to say "I am a daughter" More...
 
integer, public online_domain_num = 1
 
logical, public online_use_velz = .false.
 
logical, public online_no_rotate = .false.
 
logical, public online_boundary_use_qhyd = .false.
 

Detailed Description

module GRID (nesting system)

Description
Grid module for nesting system
Author
Team SCALE
History
  • 2014-07-28 (R.Yoshida) [new]
  • 2014-09-05 (R.Yoshida) [add] online communication system

Function/Subroutine Documentation

◆ nest_setup()

subroutine, public scale_grid_nest::nest_setup ( integer, intent(in), optional  inter_parent,
integer, intent(in), optional  inter_child 
)

Setup.

Parameters
[in]inter_childmetadata files for lat-lon domain for all processes

Definition at line 275 of file scale_grid_nest.F90.

References scale_comm::comm_world, scale_const::const_d2r, datr_ke, datr_ks, daughter_dtsec, daughter_ia, daughter_imax, daughter_ja, daughter_jmax, daughter_ka, daughter_kmax, daughter_nstep, scale_debug::debug_domain_num, handling_num, scale_tracer::i_qv, scale_grid_index::ia, scale_grid_index::ie, scale_grid_index::ihalo, intercomm_daughter, intercomm_parent, scale_interpolation_nest::intrpnest_interp_fact_llz(), scale_interpolation_nest::intrpnest_setup(), scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_get_available_fid(), scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::is, scale_grid_index::ja, scale_grid_index::je, scale_grid_index::jhalo, scale_grid_index::js, scale_grid_index::khalo, nest_bnd_qa, nest_domain_relate(), nest_filiation, nest_interp_level, nest_tile_num_x, nest_tile_num_y, offline, online_boundary_use_qhyd, online_domain_num, online_iam_daughter, online_iam_parent, online_no_rotate, online_use_velz, parent_dtsec, parent_ia, parent_imax, parent_ja, parent_jmax, parent_ka, parent_kmax, parent_lkmax, parent_nstep, scale_process::prc_global_domainid, scale_rm_process::prc_has_e, scale_rm_process::prc_has_n, scale_rm_process::prc_has_s, scale_rm_process::prc_has_w, scale_process::prc_ismaster, scale_process::prc_mpistop(), scale_process::prc_nprocs, scale_tracer::qa, scale_grid_real::real_cz, scale_grid_real::real_domain_catalogue, scale_grid_real::real_fz, scale_grid_real::real_lat, scale_grid_real::real_latx, scale_grid_real::real_latxy, scale_grid_real::real_laty, scale_grid_real::real_lon, scale_grid_real::real_lonx, scale_grid_real::real_lonxy, scale_grid_real::real_lony, tileal_ia, tileal_ja, tileal_ka, and use_nesting.

Referenced by mod_rm_driver::scalerm(), and mod_rm_prep::scalerm_prep().

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, &
338  online_domain_num, &
339  online_iam_parent, &
340  online_iam_daughter, &
341  online_use_velz, &
342  online_no_rotate, &
343  online_boundary_use_qhyd, &
344  online_aggressive_comm, &
345  online_wait_limit, &
346  online_specified_maxrq, &
347  nest_interp_level
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 
379  prc_global_domainid = online_domain_num
380 
381  call intrpnest_setup ( interp_search_divnum, nest_interp_level, offline )
382  itp_nh = int( nest_interp_level )
383  itp_nv = 2
384 
385  debug_domain_num = online_domain_num
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 
453  call nest_domain_relate(handling_num)
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 
499  parent_ka(handling_num) = parent_kmax(handling_num) + khalo * 2
500  parent_ia(handling_num) = parent_imax(handling_num) + ihalo * 2
501  parent_ja(handling_num) = parent_jmax(handling_num) + jhalo * 2
502  daughter_ka(handling_num) = daughter_kmax(handling_num) + khalo * 2
503  daughter_ia(handling_num) = daughter_imax(handling_num) + ihalo * 2
504  daughter_ja(handling_num) = daughter_jmax(handling_num) + jhalo * 2
505  tileal_ka(handling_num) = 0
506  tileal_ia(handling_num) = 0
507  tileal_ja(handling_num) = 0
508 
509  if ( .NOT. online_no_rotate ) then
510  allocate( u_llp(parent_ka(handling_num), parent_ia(handling_num), parent_ja(handling_num) ) )
511  allocate( v_llp(parent_ka(handling_num), parent_ia(handling_num), parent_ja(handling_num) ) )
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 
536  allocate( org_dens(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
537  allocate( org_momz(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
538  allocate( org_momx(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
539  allocate( org_momy(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
540  allocate( org_rhot(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
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 
574  call nest_domain_relate( handling_num )
575 
576  parent_ka(handling_num) = parent_kmax(handling_num) + khalo * 2
577  parent_ia(handling_num) = parent_imax(handling_num) + ihalo * 2
578  parent_ja(handling_num) = parent_jmax(handling_num) + jhalo * 2
579  daughter_ka(handling_num) = daughter_kmax(handling_num) + khalo * 2
580  daughter_ia(handling_num) = daughter_imax(handling_num) + ihalo * 2
581  daughter_ja(handling_num) = daughter_jmax(handling_num) + jhalo * 2
582  tileal_ka(handling_num) = parent_ka(handling_num)
583  tileal_ia(handling_num) = parent_imax(handling_num) * nest_tile_num_x
584  tileal_ja(handling_num) = parent_jmax(handling_num) * nest_tile_num_y
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
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
logical, public prc_ismaster
master process in local communicator?
logical, public prc_has_n
subroutine, public prc_mpistop
Abort MPI.
integer, dimension(2), public parent_jmax
parent max number in y-direction
logical, public prc_has_e
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)
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:35
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, parameter, public i_sw
integer, public ia
of x whole cells (local, with HALO)
real(rp), dimension(:,:), allocatable, public real_latx
latitude at staggered point (uy) [rad,-pi,pi]
module GRID (real space)
integer, public i_qv
integer, public comm_world
communication world ID
Definition: scale_comm.F90:118
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
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
integer, parameter, public khalo
of halo cells: z
real(rp), dimension(:,:,:), allocatable, public real_domain_catalogue
domain latlon catalogue [rad]
module INTERPOLATION (nesting system)
module RM PROCESS
integer, public ie
end point of inner domain: x, local
integer, dimension(2), public parent_imax
parent max number in x-direction
real(rp), dimension(:,:), allocatable, public real_lon
longitude [rad,0-2pi]
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]
logical, public prc_has_w
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]
subroutine, public intrpnest_setup(interp_search_divnum, NEST_INTERP_LEVEL, OFFLINE)
Setup.
integer, public ihalo
of halo cells: x
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)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ nest_domain_relate()

subroutine, public scale_grid_nest::nest_domain_relate ( integer, intent(in)  HANDLE)

Solve relationship between ParentDomain & Daughter Domain.

Parameters
[in]handleid number of nesting relation in this process target

Definition at line 751 of file scale_grid_nest.F90.

References scale_const::const_eps, scale_stdio::io_fid_log, scale_stdio::io_l, nest_tile_id, nest_tile_num_x, nest_tile_num_y, online_domain_num, scale_process::prc_mpistop(), and scale_process::prc_myrank.

Referenced by nest_setup().

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 
855  allocate( nest_tile_id( nest_tile_num_x*nest_tile_num_y ) )
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
subroutine, public prc_mpistop
Abort MPI.
integer, dimension(2), public parent_jmax
parent max number in y-direction
integer, parameter, public i_sw
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
integer, public prc_myrank
process num in local communicator
integer, dimension(2), public parent_imax
parent max number in x-direction
real(rp), public const_eps
small number
Definition: scale_const.F90:36
integer, parameter, public rp
Here is the call graph for this function:
Here is the caller graph for this function:

◆ nest_domain_shape()

subroutine, public scale_grid_nest::nest_domain_shape ( integer, intent(out)  tilei,
integer, intent(out)  tilej,
integer, intent(out)  cxs,
integer, intent(out)  cxe,
integer, intent(out)  cys,
integer, intent(out)  cye,
integer, intent(out)  pxs,
integer, intent(out)  pxe,
integer, intent(out)  pys,
integer, intent(out)  pye,
integer, intent(in)  iloc 
)

Return shape of ParentDomain at the specified rank (for offline)

Definition at line 881 of file scale_grid_nest.F90.

References scale_comm::comm_datatype, scale_comm::comm_world, datr_ie, datr_is, datr_je, datr_js, datr_ke, datr_ks, daughter_dtsec, daughter_imax, daughter_jmax, daughter_kmax, daughter_nstep, scale_grid_index::ie, scale_grid_index::imax, intercomm_daughter, intercomm_parent, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::is, scale_grid_index::je, scale_grid_index::jmax, scale_grid_index::js, scale_grid_index::ke, scale_grid_index::kmax, scale_grid_index::ks, nest_bnd_qa, nest_filiation, nest_tile_id, nest_tile_num_x, online_domain_num, online_no_rotate, online_use_velz, parent_dtsec, parent_ia, parent_imax, parent_ja, parent_jmax, parent_ka, parent_kmax, parent_nstep, scale_process::prc_ismaster, scale_process::prc_mpistop(), scale_process::prc_myrank, scale_process::prc_nprocs, scale_rm_process::prc_num_x, scale_rm_process::prc_num_y, prnt_ie, prnt_is, prnt_je, prnt_js, prnt_ke, prnt_ks, scale_grid_real::real_cz, scale_grid_real::real_domain_catalogue, scale_grid_real::real_fz, scale_grid_real::real_lat, scale_grid_real::real_latx, scale_grid_real::real_laty, scale_grid_real::real_lon, scale_grid_real::real_lonx, scale_grid_real::real_lony, scale_time::time_dtsec, and scale_time::time_nstep.

Referenced by mod_copytopo::copytopo().

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
integer, dimension(2), public parent_jmax
parent max number in y-direction
integer, dimension(2), public parent_imax
parent max number in x-direction
Here is the call graph for this function:
Here is the caller graph for this function:

◆ nest_comm_nestdown()

subroutine, public scale_grid_nest::nest_comm_nestdown ( integer, intent(in)  HANDLE,
integer, intent(in)  BND_QA,
real(rp), dimension(parent_ka(handle),parent_ia(handle),parent_ja(handle)), intent(in)  ipt_DENS,
real(rp), dimension(parent_ka(handle),parent_ia(handle),parent_ja(handle)), intent(in)  ipt_MOMZ,
real(rp), dimension(parent_ka(handle),parent_ia(handle),parent_ja(handle)), intent(in)  ipt_MOMX,
real(rp), dimension(parent_ka(handle),parent_ia(handle),parent_ja(handle)), intent(in)  ipt_MOMY,
real(rp), dimension(parent_ka(handle),parent_ia(handle),parent_ja(handle)), intent(in)  ipt_RHOT,
real(rp), dimension(parent_ka(handle),parent_ia(handle),parent_ja(handle),bnd_qa), intent(in)  ipt_QTRC,
real(rp), dimension(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle)), intent(inout)  interped_ref_DENS,
real(rp), dimension(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle)), intent(inout)  interped_ref_VELZ,
real(rp), dimension(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle)), intent(inout)  interped_ref_VELX,
real(rp), dimension(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle)), intent(inout)  interped_ref_VELY,
real(rp), dimension(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle)), intent(inout)  interped_ref_POTT,
real(rp), dimension(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle),bnd_qa), intent(inout)  interped_ref_QTRC 
)

Boundary data transfer from parent to daughter: nestdown.

Parameters
[in]handleid number of nesting relation in this process target
[in]bnd_qanum of tracer

Definition at line 1618 of file scale_grid_nest.F90.

References datr_ke, datr_ks, daughter_ia, daughter_ja, scale_gridtrans::gtrans_rotc, intercomm_parent, scale_stdio::io_fid_log, scale_stdio::io_l, nest_comm_waitall(), nest_filiation, online_no_rotate, online_use_velz, parent_ia, parent_ja, scale_process::prc_mpistop(), scale_process::prc_myrank, scale_process::prc_nprocs, prnt_ie, prnt_is, prnt_je, prnt_js, prnt_ke, prnt_ks, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), and scale_grid_real::real_domain_catalogue.

Referenced by scale_atmos_boundary::atmos_boundary_update().

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
subroutine, public prc_mpistop
Abort MPI.
module GRIDTRANS
module GRID (real space)
module COMMUNICATION
Definition: scale_comm.F90:23
module PROCESS
real(rp), dimension(:,:,:), allocatable, public gtrans_rotc
rotation coefficient
integer, public prc_myrank
process num in local communicator
real(rp), dimension(:,:,:), allocatable, public real_domain_catalogue
domain latlon catalogue [rad]
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:132
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:178
integer, public prc_nprocs
myrank in local communicator
Here is the call graph for this function:
Here is the caller graph for this function:

◆ nest_comm_recvwait_issue()

subroutine, public scale_grid_nest::nest_comm_recvwait_issue ( integer, intent(in)  HANDLE,
integer, intent(in)  BND_QA 
)

Sub-command for data transfer from parent to daughter: nestdown.

Parameters
[in]handleid number of nesting relation in this process target
[in]bnd_qanum of tracer in online-nesting

Definition at line 1968 of file scale_grid_nest.F90.

References intercomm_daughter, scale_stdio::io_fid_log, scale_stdio::io_l, nest_filiation, online_no_rotate, online_use_velz, scale_process::prc_mpistop(), scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

Referenced by scale_atmos_boundary::atmos_boundary_finalize(), scale_atmos_boundary::atmos_boundary_resume_file(), scale_atmos_boundary::atmos_boundary_resume_online(), and scale_atmos_boundary::atmos_boundary_update().

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
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:132
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:178
Here is the call graph for this function:
Here is the caller graph for this function:

◆ nest_comm_recv_cancel()

subroutine, public scale_grid_nest::nest_comm_recv_cancel ( integer, intent(in)  HANDLE)

Sub-command for data transfer from parent to daughter: nestdown.

Parameters
[in]handleid number of nesting relation in this process target

Definition at line 2064 of file scale_grid_nest.F90.

References scale_stdio::io_fid_log, scale_stdio::io_l, nest_filiation, scale_process::prc_mpistop(), scale_process::prc_myrank, and scale_process::prc_nprocs.

Referenced by scale_atmos_boundary::atmos_boundary_finalize().

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
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
integer, public prc_myrank
process num in local communicator
integer, public prc_nprocs
myrank in local communicator
Here is the call graph for this function:
Here is the caller graph for this function:

◆ nest_comm_intercomm_nestdown_3d()

subroutine scale_grid_nest::nest_comm_intercomm_nestdown_3d ( real(rp), dimension(:,:,:), intent(in)  pvar,
real(rp), dimension(:,:,:), intent(out)  dvar,
integer, intent(in)  tagbase,
integer, intent(in)  id_stag,
integer, intent(in)  HANDLE,
integer, intent(inout)  isu_tag,
integer, intent(inout)  isu_tagf,
logical, intent(in), optional  flag_dens 
)

Inter-communication from parent to daughter: nestdown.

Parameters
[in]pvarvariable from parent domain (PARENT_KA,PARENT_IA,PARENT_JA / 1,1,1)
[out]dvarvariable to daughter domain (1,1,1 / MY_KA,MY_IA,MY_JA)
[in]tagbasecommunication tag of the variable
[in]id_stagid of staggered grid option
[in]handleid number of nesting relation in this process target
[in,out]isu_tagtag for receive buffer
[in,out]isu_tagftag for receive buffer
[in]flag_densflag of logarithmic interpolation for density

Definition at line 2115 of file scale_grid_nest.F90.

References scale_comm::comm_datatype, datr_ke, datr_ks, daughter_ia, daughter_ja, intercomm_daughter, scale_interpolation_nest::intrpnest_interp_3d, dc_log::log(), nest_filiation, nest_tile_num_x, parent_ia, parent_imax, parent_ja, parent_jmax, parent_ka, scale_process::prc_mpistop(), scale_process::prc_myrank, scale_process::prc_nprocs, prnt_ie, prnt_is, prnt_je, and prnt_js.

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
integer, public comm_datatype
datatype of variable
Definition: scale_comm.F90:117
subroutine, public prc_mpistop
Abort MPI.
integer, dimension(2), public parent_jmax
parent max number in y-direction
module COMMUNICATION
Definition: scale_comm.F90:23
module PROCESS
subroutine, public log(type, message)
Definition: dc_log.f90:133
procedure(intrpnest_intfc_interp_3d), pointer, public intrpnest_interp_3d
integer, public prc_myrank
process num in local communicator
module INTERPOLATION (nesting system)
integer, dimension(2), public parent_imax
parent max number in x-direction
integer, public prc_nprocs
myrank in local communicator
Here is the call graph for this function:

◆ nest_comm_issuer_of_receive_3d()

subroutine scale_grid_nest::nest_comm_issuer_of_receive_3d ( integer, intent(in)  tagbase,
integer, intent(in)  id_stag,
integer, intent(in)  HANDLE,
integer, intent(inout)  isu_tag,
integer, intent(inout)  isu_tagf,
logical, intent(in), optional  flag_dens 
)

[substance of issuer] Inter-communication from parent to daughter: nestdown

Parameters
[in]tagbasecommunication tag of the variable
[in]id_stagid of staggered grid option
[in]handleid number of nesting relation in this process target
[in,out]isu_tagtag for receive buffer
[in,out]isu_tagftag for receive buffer
[in]flag_densflag of logarithmic interpolation for density

Definition at line 2294 of file scale_grid_nest.F90.

References scale_comm::comm_datatype, intercomm_parent, nest_filiation, parent_ia, parent_ja, parent_ka, scale_process::prc_mpistop(), scale_process::prc_myrank, and scale_process::prc_nprocs.

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
integer, public comm_datatype
datatype of variable
Definition: scale_comm.F90:117
subroutine, public prc_mpistop
Abort MPI.
module COMMUNICATION
Definition: scale_comm.F90:23
module PROCESS
integer, public prc_myrank
process num in local communicator
integer, public prc_nprocs
myrank in local communicator
Here is the call graph for this function:

◆ nest_comm_issuer_of_wait_3d()

subroutine scale_grid_nest::nest_comm_issuer_of_wait_3d ( integer, intent(in)  HANDLE)

[substance of issuer] Inter-communication from parent to daughter: nestdown

Parameters
[in]handleid number of nesting relation in this process target

Definition at line 2408 of file scale_grid_nest.F90.

References nest_comm_waitall(), nest_filiation, and scale_process::prc_mpistop().

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
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
Here is the call graph for this function:

◆ nest_comm_waitall()

subroutine scale_grid_nest::nest_comm_waitall ( integer, intent(in)  req_count,
integer, dimension(max_rq), intent(inout)  ireq 
)

[substance of comm_wait] Inter-communication

Definition at line 2437 of file scale_grid_nest.F90.

References scale_process::prc_mpistop().

Referenced by nest_comm_issuer_of_wait_3d(), and nest_comm_nestdown().

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
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
Here is the call graph for this function:
Here is the caller graph for this function:

◆ nest_comm_test()

subroutine, public scale_grid_nest::nest_comm_test ( integer, intent(in)  HANDLE)

[check communication status] Inter-communication

Parameters
[in]handleid number of nesting relation in this process target

Definition at line 2483 of file scale_grid_nest.F90.

References nest_filiation, scale_process::prc_mpistop(), scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

Referenced by scale_atmos_boundary::atmos_boundary_update().

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
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:132
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:178
Here is the call graph for this function:
Here is the caller graph for this function:

◆ nest_comm_disconnect()

subroutine, public scale_grid_nest::nest_comm_disconnect ( )

[finalize: disconnect] Inter-communication

Definition at line 2517 of file scale_grid_nest.F90.

References intercomm_daughter, intercomm_parent, scale_stdio::io_fid_log, scale_stdio::io_l, online_iam_daughter, online_iam_parent, and scale_process::prc_global_comm_world.

Referenced by mod_atmos_driver::atmos_driver_finalize().

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
integer, public prc_global_comm_world
global communicator
module PROCESS
Here is the caller graph for this function:

Variable Documentation

◆ intercomm_parent

integer, public scale_grid_nest::intercomm_parent

Definition at line 55 of file scale_grid_nest.F90.

Referenced by nest_comm_disconnect(), nest_comm_issuer_of_receive_3d(), nest_comm_nestdown(), nest_domain_shape(), and nest_setup().

55  integer, public :: intercomm_parent ! inter-communicator to parent

◆ intercomm_daughter

integer, public scale_grid_nest::intercomm_daughter

Definition at line 56 of file scale_grid_nest.F90.

Referenced by nest_comm_disconnect(), nest_comm_intercomm_nestdown_3d(), nest_comm_recvwait_issue(), nest_domain_shape(), and nest_setup().

56  integer, public :: intercomm_daughter ! inter-communicator to daughter

◆ nest_filiation

integer, dimension(10), public scale_grid_nest::nest_filiation

index of parent-daughter relation (p>0, d<0)

Definition at line 58 of file scale_grid_nest.F90.

Referenced by nest_comm_intercomm_nestdown_3d(), nest_comm_issuer_of_receive_3d(), nest_comm_issuer_of_wait_3d(), nest_comm_nestdown(), nest_comm_recv_cancel(), nest_comm_recvwait_issue(), nest_comm_test(), nest_domain_shape(), and nest_setup().

58  integer, public :: nest_filiation(10)

◆ handling_num

integer, public scale_grid_nest::handling_num

handing number of nesting relation

Definition at line 59 of file scale_grid_nest.F90.

Referenced by nest_setup().

59  integer, public :: handling_num

◆ nest_tile_num_x

integer, public scale_grid_nest::nest_tile_num_x

◆ nest_tile_num_y

integer, public scale_grid_nest::nest_tile_num_y

◆ nest_tile_id

integer, dimension(:), allocatable, public scale_grid_nest::nest_tile_id

◆ parent_kmax

integer, dimension(2), public scale_grid_nest::parent_kmax

parent max number in z-direction

Definition at line 64 of file scale_grid_nest.F90.

Referenced by nest_domain_shape(), nest_setup(), and mod_realinput_scale::parentatomsetupscale().

64  integer, public :: parent_kmax(2)

◆ parent_imax

integer, dimension(2), public scale_grid_nest::parent_imax

◆ parent_jmax

integer, dimension(2), public scale_grid_nest::parent_jmax

◆ parent_ka

integer, dimension(2), public scale_grid_nest::parent_ka

parent max number in z-direction (with halo)

Definition at line 67 of file scale_grid_nest.F90.

Referenced by scale_atmos_boundary::atmos_boundary_update(), nest_comm_intercomm_nestdown_3d(), nest_comm_issuer_of_receive_3d(), nest_domain_shape(), and nest_setup().

67  integer, public :: parent_ka(2)

◆ parent_ia

integer, dimension(2), public scale_grid_nest::parent_ia

parent max number in x-direction (with halo)

Definition at line 68 of file scale_grid_nest.F90.

Referenced by scale_atmos_boundary::atmos_boundary_update(), nest_comm_intercomm_nestdown_3d(), nest_comm_issuer_of_receive_3d(), nest_comm_nestdown(), nest_domain_shape(), and nest_setup().

68  integer, public :: parent_ia(2)

◆ parent_ja

integer, dimension(2), public scale_grid_nest::parent_ja

parent max number in y-direction (with halo)

Definition at line 69 of file scale_grid_nest.F90.

Referenced by scale_atmos_boundary::atmos_boundary_update(), nest_comm_intercomm_nestdown_3d(), nest_comm_issuer_of_receive_3d(), nest_comm_nestdown(), nest_domain_shape(), and nest_setup().

69  integer, public :: parent_ja(2)

◆ parent_lkmax

integer, dimension(2), public scale_grid_nest::parent_lkmax

parent max number in lz-direction

Definition at line 70 of file scale_grid_nest.F90.

Referenced by nest_setup(), and mod_realinput_scale::parentlandsetupscale().

70  integer, public :: parent_lkmax(2)

◆ parent_dtsec

real(dp), dimension(2), public scale_grid_nest::parent_dtsec

parent DT [sec]

Definition at line 71 of file scale_grid_nest.F90.

Referenced by scale_atmos_boundary::atmos_boundary_resume_file(), nest_domain_shape(), and nest_setup().

71  real(DP), public :: parent_dtsec(2)

◆ parent_nstep

integer, dimension(2), public scale_grid_nest::parent_nstep

parent step [number]

Definition at line 72 of file scale_grid_nest.F90.

Referenced by scale_atmos_boundary::atmos_boundary_resume_online(), nest_domain_shape(), and nest_setup().

72  integer, public :: parent_nstep(2)

◆ daughter_kmax

integer, dimension(2), public scale_grid_nest::daughter_kmax

daughter max number in z-direction

Definition at line 74 of file scale_grid_nest.F90.

Referenced by nest_domain_shape(), and nest_setup().

74  integer, public :: daughter_kmax(2)

◆ daughter_imax

integer, dimension(2), public scale_grid_nest::daughter_imax

daughter max number in x-direction

Definition at line 75 of file scale_grid_nest.F90.

Referenced by nest_domain_shape(), and nest_setup().

75  integer, public :: daughter_imax(2)

◆ daughter_jmax

integer, dimension(2), public scale_grid_nest::daughter_jmax

daughter max number in y-direction

Definition at line 76 of file scale_grid_nest.F90.

Referenced by nest_domain_shape(), and nest_setup().

76  integer, public :: daughter_jmax(2)

◆ daughter_ka

integer, dimension(2), public scale_grid_nest::daughter_ka

daughter max number in z-direction (with halo)

Definition at line 77 of file scale_grid_nest.F90.

Referenced by scale_atmos_boundary::atmos_boundary_update(), and nest_setup().

77  integer, public :: daughter_ka(2)

◆ daughter_ia

integer, dimension(2), public scale_grid_nest::daughter_ia

daughter max number in x-direction (with halo)

Definition at line 78 of file scale_grid_nest.F90.

Referenced by scale_atmos_boundary::atmos_boundary_update(), nest_comm_intercomm_nestdown_3d(), nest_comm_nestdown(), and nest_setup().

78  integer, public :: daughter_ia(2)

◆ daughter_ja

integer, dimension(2), public scale_grid_nest::daughter_ja

daughter max number in y-direction (with halo)

Definition at line 79 of file scale_grid_nest.F90.

Referenced by scale_atmos_boundary::atmos_boundary_update(), nest_comm_intercomm_nestdown_3d(), nest_comm_nestdown(), and nest_setup().

79  integer, public :: daughter_ja(2)

◆ daughter_lkmax

integer, dimension(2), public scale_grid_nest::daughter_lkmax

daughter max number in lz-direction

Definition at line 80 of file scale_grid_nest.F90.

80  integer, public :: daughter_lkmax(2)

◆ daughter_dtsec

real(dp), dimension(2), public scale_grid_nest::daughter_dtsec

daughter DT [sec]

Definition at line 81 of file scale_grid_nest.F90.

Referenced by nest_domain_shape(), and nest_setup().

81  real(DP), public :: daughter_dtsec(2)

◆ daughter_nstep

integer, dimension(2), public scale_grid_nest::daughter_nstep

daughter steps [number]

Definition at line 82 of file scale_grid_nest.F90.

Referenced by nest_domain_shape(), and nest_setup().

82  integer, public :: daughter_nstep(2)

◆ prnt_ks

integer, dimension(2), public scale_grid_nest::prnt_ks

start index in z-direction in parent

Definition at line 84 of file scale_grid_nest.F90.

Referenced by nest_comm_nestdown(), and nest_domain_shape().

84  integer, public :: prnt_ks(2)

◆ prnt_ke

integer, dimension(2), public scale_grid_nest::prnt_ke

end index in z-direction in parent

Definition at line 85 of file scale_grid_nest.F90.

Referenced by nest_comm_nestdown(), and nest_domain_shape().

85  integer, public :: prnt_ke(2)

◆ prnt_is

integer, dimension(2), public scale_grid_nest::prnt_is

start index in x-direction in parent

Definition at line 86 of file scale_grid_nest.F90.

Referenced by nest_comm_intercomm_nestdown_3d(), nest_comm_nestdown(), and nest_domain_shape().

86  integer, public :: prnt_is(2)

◆ prnt_ie

integer, dimension(2), public scale_grid_nest::prnt_ie

end index in x-direction in parent

Definition at line 87 of file scale_grid_nest.F90.

Referenced by nest_comm_intercomm_nestdown_3d(), nest_comm_nestdown(), and nest_domain_shape().

87  integer, public :: prnt_ie(2)

◆ prnt_js

integer, dimension(2), public scale_grid_nest::prnt_js

start index in y-direction in parent

Definition at line 88 of file scale_grid_nest.F90.

Referenced by nest_comm_intercomm_nestdown_3d(), nest_comm_nestdown(), and nest_domain_shape().

88  integer, public :: prnt_js(2)

◆ prnt_je

integer, dimension(2), public scale_grid_nest::prnt_je

end index in y-direction in parent

Definition at line 89 of file scale_grid_nest.F90.

Referenced by nest_comm_intercomm_nestdown_3d(), nest_comm_nestdown(), and nest_domain_shape().

89  integer, public :: prnt_je(2)

◆ datr_ks

integer, dimension(2), public scale_grid_nest::datr_ks

start index in z-direction in daughter

Definition at line 91 of file scale_grid_nest.F90.

Referenced by nest_comm_intercomm_nestdown_3d(), nest_comm_nestdown(), nest_domain_shape(), and nest_setup().

91  integer, public :: datr_ks(2)

◆ datr_ke

integer, dimension(2), public scale_grid_nest::datr_ke

end index in z-direction in daughter

Definition at line 92 of file scale_grid_nest.F90.

Referenced by nest_comm_intercomm_nestdown_3d(), nest_comm_nestdown(), nest_domain_shape(), and nest_setup().

92  integer, public :: datr_ke(2)

◆ datr_is

integer, dimension(2), public scale_grid_nest::datr_is

start index in x-direction in daughter

Definition at line 93 of file scale_grid_nest.F90.

Referenced by nest_domain_shape().

93  integer, public :: datr_is(2)

◆ datr_ie

integer, dimension(2), public scale_grid_nest::datr_ie

end index in x-direction in daughter

Definition at line 94 of file scale_grid_nest.F90.

Referenced by nest_domain_shape().

94  integer, public :: datr_ie(2)

◆ datr_js

integer, dimension(2), public scale_grid_nest::datr_js

start index in y-direction in daughter

Definition at line 95 of file scale_grid_nest.F90.

Referenced by nest_domain_shape().

95  integer, public :: datr_js(2)

◆ datr_je

integer, dimension(2), public scale_grid_nest::datr_je

end index in y-direction in daughter

Definition at line 96 of file scale_grid_nest.F90.

Referenced by nest_domain_shape().

96  integer, public :: datr_je(2)

◆ tileal_ka

integer, dimension(2), public scale_grid_nest::tileal_ka

cells of all tiles in z-direction

Definition at line 98 of file scale_grid_nest.F90.

Referenced by nest_setup().

98  integer, public :: tileal_ka(2)

◆ tileal_ia

integer, dimension(2), public scale_grid_nest::tileal_ia

cells of all tiles in x-direction

Definition at line 99 of file scale_grid_nest.F90.

Referenced by nest_setup().

99  integer, public :: tileal_ia(2)

◆ tileal_ja

integer, dimension(2), public scale_grid_nest::tileal_ja

cells of all tiles in y-direction

Definition at line 100 of file scale_grid_nest.F90.

Referenced by nest_setup().

100  integer, public :: tileal_ja(2)

◆ nest_bnd_qa

integer, public scale_grid_nest::nest_bnd_qa = 1

number of tracer treated in nesting system

Definition at line 102 of file scale_grid_nest.F90.

Referenced by scale_atmos_boundary::atmos_boundary_finalize(), scale_atmos_boundary::atmos_boundary_resume_file(), scale_atmos_boundary::atmos_boundary_update(), nest_domain_shape(), and nest_setup().

102  integer, public :: nest_bnd_qa = 1

◆ nest_interp_level

integer, public scale_grid_nest::nest_interp_level = 4

horizontal interpolation level

Definition at line 103 of file scale_grid_nest.F90.

Referenced by mod_copytopo::copytopo(), nest_setup(), and mod_realinput::parentatomsetup().

103  integer, public :: nest_interp_level = 4

◆ use_nesting

logical, public scale_grid_nest::use_nesting = .false.

Definition at line 105 of file scale_grid_nest.F90.

Referenced by scale_atmos_boundary::atmos_boundary_setup(), and nest_setup().

105  logical, public :: use_nesting = .false.

◆ offline

logical, public scale_grid_nest::offline = .true.

Definition at line 106 of file scale_grid_nest.F90.

Referenced by scale_atmos_boundary::atmos_boundary_setup(), and nest_setup().

106  logical, public :: offline = .true.

◆ online_iam_parent

logical, public scale_grid_nest::online_iam_parent = .false.

a flag to say "I am a parent"

Definition at line 107 of file scale_grid_nest.F90.

Referenced by scale_atmos_boundary::atmos_boundary_setup(), nest_comm_disconnect(), and nest_setup().

107  logical, public :: online_iam_parent = .false.

◆ online_iam_daughter

logical, public scale_grid_nest::online_iam_daughter = .false.

a flag to say "I am a daughter"

Definition at line 108 of file scale_grid_nest.F90.

Referenced by scale_atmos_boundary::atmos_boundary_setup(), nest_comm_disconnect(), and nest_setup().

108  logical, public :: online_iam_daughter = .false.

◆ online_domain_num

integer, public scale_grid_nest::online_domain_num = 1

Definition at line 109 of file scale_grid_nest.F90.

Referenced by nest_domain_relate(), nest_domain_shape(), and nest_setup().

109  integer, public :: online_domain_num = 1

◆ online_use_velz

logical, public scale_grid_nest::online_use_velz = .false.

◆ online_no_rotate

logical, public scale_grid_nest::online_no_rotate = .false.

Definition at line 111 of file scale_grid_nest.F90.

Referenced by nest_comm_nestdown(), nest_comm_recvwait_issue(), nest_domain_shape(), and nest_setup().

111  logical, public :: online_no_rotate = .false.

◆ online_boundary_use_qhyd

logical, public scale_grid_nest::online_boundary_use_qhyd = .false.

Definition at line 112 of file scale_grid_nest.F90.

Referenced by nest_setup().

112  logical, public :: online_boundary_use_qhyd = .false.