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 = 5
 horizontal interpolation level More...
 
integer, public nest_interp_weight_order = 2
 horizontal interpolation weight order More...
 
logical, public use_nesting = .false.
 
logical, public offline = .false.
 
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.
 
logical, public online_boundary_diagqnum = .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
NAMELIST
  • PARAM_NEST
    nametypedefault valuecomment
    LATLON_CATALOGUE_FNAME character(len=H_LONG) 'latlon_domain_catalogue.txt'
    OFFLINE_PARENT_BASENAME character(len=H_LONG) parent file base name
    OFFLINE_PARENT_PRC_NUM_X integer MPI processes in x-direction in parent [for namelist]
    OFFLINE_PARENT_PRC_NUM_Y integer MPI processes in y-direction in parent [for namelist]
    ONLINE_DOMAIN_NUM integer 1
    ONLINE_IAM_PARENT logical .false. a flag to say "I am a parent"
    ONLINE_IAM_DAUGHTER logical .false. a flag to say "I am a daughter"
    ONLINE_USE_VELZ logical .false.
    ONLINE_NO_ROTATE logical .false.
    ONLINE_BOUNDARY_USE_QHYD logical .false.
    ONLINE_BOUNDARY_DIAGQNUM logical .false.
    ONLINE_AGGRESSIVE_COMM logical
    ONLINE_WAIT_LIMIT integer(8) limit times of waiting loop in "NEST_COMM_waitall"
    ONLINE_SPECIFIED_MAXRQ integer 0
    NEST_INTERP_LEVEL integer 5 horizontal interpolation level
    NEST_INTERP_WEIGHT_ORDER integer 2 horizontal interpolation weight order

History Output
No history output

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 278 of file scale_grid_nest.F90.

References 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, gtool_file::filegetshape(), handling_num, scale_atmos_hydrometeor::i_qv, scale_grid_index::ia, 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_fid_nml, scale_stdio::io_get_available_fid(), scale_stdio::io_l, scale_stdio::io_nml, scale_grid_index::ja, scale_grid_index::jhalo, scale_grid_index::khalo, nest_bnd_qa, nest_domain_relate(), nest_filiation, nest_interp_level, nest_interp_weight_order, nest_tile_num_x, nest_tile_num_y, offline, online_boundary_diagqnum, 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_atmos_phy_mp::qa_mp, scale_grid_real::real_cz, 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().

278 ! flag_parent, &
279 ! flag_child )
280  use gtool_file, only: &
282  use scale_const, only: &
283  d2r => const_d2r
284  use scale_process, only: &
285  prc_mpistop, &
288  use scale_rm_process, only: &
289  prc_has_w, &
290  prc_has_e, &
291  prc_has_s, &
292  prc_has_n
293  use scale_grid_real, only: &
294  real_lonxy, &
295  real_latxy, &
296  my_lon => real_lon, &
297  my_lat => real_lat, &
298  my_lonx => real_lonx, &
299  my_latx => real_latx, &
300  my_lony => real_lony, &
301  my_laty => real_laty, &
302  my_cz => real_cz, &
303  my_fz => real_fz
304  use scale_comm, only: &
305  comm_bcast
306  use scale_interpolation_nest, only: &
307  intrpnest_setup, &
309  use scale_atmos_hydrometeor, only: &
310  i_qv
311  use scale_atmos_phy_mp, only: &
312  qa_mp
313  implicit none
314 
315  integer, intent(in), optional :: inter_parent
316  integer, intent(in), optional :: inter_child
317 
318 
319  character(len=H_LONG) :: LATLON_CATALOGUE_FNAME = 'latlon_domain_catalogue.txt'
320 
321  integer :: ONLINE_SPECIFIED_MAXRQ = 0
322  integer :: i
323  integer :: fid, ierr
324  integer :: parent_id
325  integer, allocatable :: errcodes(:)
326 
327  character(len=2) :: dom_num
328 
329  logical :: flag_parent = .false.
330  logical :: flag_child = .false.
331 
332  integer :: dims(1)
333  logical :: error
334 
335  namelist / param_nest / &
336  latlon_catalogue_fname, &
337  offline_parent_basename, &
338  offline_parent_prc_num_x, &
339  offline_parent_prc_num_y, &
340  online_domain_num, &
341  online_iam_parent, &
342  online_iam_daughter, &
343  online_use_velz, &
344  online_no_rotate, &
345  online_boundary_use_qhyd, &
346  online_boundary_diagqnum, &
347  online_aggressive_comm, &
348  online_wait_limit, &
349  online_specified_maxrq, &
350  nest_interp_level, &
351  nest_interp_weight_order
352 
353  !---------------------------------------------------------------------------
354 
355  if( io_l ) write(io_fid_log,*)
356  if( io_l ) write(io_fid_log,*) '++++++ Module[GRID_NEST] / Categ[ATMOS-RM GRID] / Origin[SCALElib]'
357 
358  if( inter_parent /= mpi_comm_null ) flag_child = .true. ! exist parent, so work as a child
359  if( inter_child /= mpi_comm_null ) flag_parent = .true. ! exist child, so work as a parent
360 
361  offline_parent_basename = ""
362 
363  nwait_p = 0
364  nwait_d = 0
365  nrecv = 0
366  nsend = 0
367 
368  handling_num = 0
369  nest_filiation(:) = 0
370  online_wait_limit = 999999999
371  online_aggressive_comm = .true.
372  interp_search_divnum = 10
373 
374  !--- read namelist
375  rewind(io_fid_conf)
376  read(io_fid_conf,nml=param_nest,iostat=ierr)
377  if( ierr < 0 ) then !--- missing
378  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
379  elseif( ierr > 0 ) then !--- fatal error
380  write(*,*) 'xxx Not appropriate names in namelist PARAM_NEST. Check!'
381  call prc_mpistop
382  endif
383  if( io_nml ) write(io_fid_nml,nml=param_nest)
384 
385  prc_global_domainid = online_domain_num
386 
387  if ( offline_parent_basename /= "" ) then
388 
389  offline = .true.
390  use_nesting = .true.
391 
392  if ( prc_ismaster ) then
393  call filegetshape( dims, offline_parent_basename, "CX", 0 )
394  offline_parent_imax = dims(1)-4
395  call filegetshape( dims, offline_parent_basename, "CY", 0 )
396  offline_parent_jmax = dims(1)-4
397  call filegetshape( dims, offline_parent_basename, "z", 0, error=error )
398  if ( error ) then
399  offline_parent_kmax = 0
400  else
401  offline_parent_kmax = dims(1)
402  end if
403  call filegetshape( dims, offline_parent_basename, "lz", 0, error=error )
404  if ( error ) then
405  offline_parent_lkmax = 0
406  else
407  offline_parent_lkmax = dims(1)
408  end if
409  end if
410  call comm_bcast( offline_parent_imax )
411  call comm_bcast( offline_parent_jmax )
412  call comm_bcast( offline_parent_kmax )
413  call comm_bcast( offline_parent_lkmax )
414  end if
415 
416  if ( online_iam_daughter .or. online_iam_parent ) then
417 
418  if ( offline ) then
419  write(*,*) 'xxx OFFLINE and ONLINE cannot be use at the same time'
420  call prc_mpistop
421  end if
422 
423  use_nesting = .true.
424  call intrpnest_setup( interp_search_divnum, &
425  nest_interp_level, &
426  nest_interp_weight_order, &
427  .false. )
428  else
429  call intrpnest_setup( interp_search_divnum, &
430  nest_interp_level, &
431  nest_interp_weight_order, &
432  .true. )
433  end if
434 
435  itp_nh = int( nest_interp_level )
436  itp_nv = 2
437 
438  debug_domain_num = online_domain_num
439  if( online_specified_maxrq > max_rq ) max_rq = online_specified_maxrq
440  allocate( ireq_p(max_rq) )
441  allocate( ireq_d(max_rq) )
442  allocate( call_order(max_rq) )
443 
444  if( use_nesting ) then
445 
446  if ( offline .OR. online_iam_daughter ) then
447  corner_loc(i_nw,i_lon) = real_lonxy( 0,ja) / d2r
448  corner_loc(i_ne,i_lon) = real_lonxy(ia,ja) / d2r
449  corner_loc(i_sw,i_lon) = real_lonxy( 0, 0) / d2r
450  corner_loc(i_se,i_lon) = real_lonxy(ia, 0) / d2r
451  corner_loc(i_nw,i_lat) = real_latxy( 0,ja) / d2r
452  corner_loc(i_ne,i_lat) = real_latxy(ia,ja) / d2r
453  corner_loc(i_sw,i_lat) = real_latxy( 0, 0) / d2r
454  corner_loc(i_se,i_lat) = real_latxy(ia, 0) / d2r
455 
456  allocate( ncopy(ia,ja,itp_nh,itp_ng) )
457  end if
458 
459  if( offline ) then
460  handling_num = 1
461  parent_prc_num_x(handling_num) = offline_parent_prc_num_x
462  parent_prc_num_y(handling_num) = offline_parent_prc_num_y
463  parent_kmax(handling_num) = offline_parent_kmax
464  parent_imax(handling_num) = offline_parent_imax
465  parent_jmax(handling_num) = offline_parent_jmax
466  parent_lkmax(handling_num) = offline_parent_lkmax
467 
468  parent_prc_nprocs(handling_num) = parent_prc_num_x(handling_num) * parent_prc_num_y(handling_num)
469  allocate( latlon_catalog(parent_prc_nprocs(handling_num),4,2) )
470 
471  !--- read latlon catalogue
472  fid = io_get_available_fid()
473  open( fid, &
474  file = trim(latlon_catalogue_fname), &
475  form = 'formatted', &
476  status = 'old', &
477  iostat = ierr )
478 
479  if ( ierr /= 0 ) then
480  write(*,*) 'xxx [grd_nest/NEST_setup] cannot open latlon-catalogue file!'
481  call prc_mpistop
482  endif
483 
484  do i = 1, parent_prc_nprocs(handling_num)
485  read(fid,'(i8,8f32.24)',iostat=ierr) parent_id, &
486  latlon_catalog(i,i_nw,i_lon), latlon_catalog(i,i_ne,i_lon), & ! LON: NW, NE
487  latlon_catalog(i,i_sw,i_lon), latlon_catalog(i,i_se,i_lon), & ! LON: SW, SE
488  latlon_catalog(i,i_nw,i_lat), latlon_catalog(i,i_ne,i_lat), & ! LAT: NW, NE
489  latlon_catalog(i,i_sw,i_lat), latlon_catalog(i,i_se,i_lat) ! LAT: SW, SE
490  if ( i /= parent_id ) then
491  write(*,*) 'xxx [grd_nest/NEST_setup] internal error: parent mpi id'
492  call prc_mpistop
493  endif
494  if ( ierr /= 0 ) exit
495  enddo
496  close(fid)
497 
498  call nest_domain_relate(handling_num)
499 
500  else ! ONLINE RELATIONSHIP
501 ! if ( present(flag_parent) .AND. present(flag_child) ) then
502 ! if( IO_L ) write(IO_FID_LOG,'(1x,A)') &
503 ! '*** Setup Online Nesting Inter-Domain Communicator (IDC)'
504 ! else
505 ! write(*,*) 'xxx Internal Error:'
506 ! write(*,*) 'xxx The flag_parent and flag_child are needed.'
507 ! write(*,*) ' domain: ', ONLINE_DOMAIN_NUM
508 ! call PRC_MPIstop
509 ! endif
510 
511  if( online_boundary_use_qhyd ) then
512  nest_bnd_qa = qa_mp
513  else if ( i_qv > 0 ) then
514  nest_bnd_qa = 1
515  else
516  nest_bnd_qa = 0
517  endif
518 
519 if( io_l ) write(io_fid_log,*) "flag_parent", flag_parent, "flag_child", flag_child
520 if( io_l ) write(io_fid_log,*) "ONLINE_IAM_PARENT", online_iam_parent, "ONLINE_IAM_DAUGHTER", online_iam_daughter
521 
522  if( flag_parent ) then ! must do first before daughter processes
523  !-------------------------------------------------
524  if ( .NOT. online_iam_parent ) then
525  write(*,*) 'xxx [grd_nest/NEST_setup] Parent Flag from launcher is not consistent with namelist!'
526  write(*,*) 'xxx PARENT - domain : ', online_domain_num
527  call prc_mpistop
528  endif
529 
530  handling_num = 1 !HANDLING_NUM + 1
531  intercomm_id(handling_num) = online_domain_num
532  nest_filiation(intercomm_id(handling_num)) = 1
533 
534  intercomm_daughter = inter_child
535  if( io_l ) write(io_fid_log,'(1x,A,I2,A)') '*** Online Nesting - PARENT [INTERCOMM_ID:', &
536  intercomm_id(handling_num), ' ]'
537  if( io_l ) write(io_fid_log,*) '*** Online Nesting - INTERCOMM :', intercomm_daughter
538 
539  call nest_comm_ping( handling_num )
540 
541  call nest_comm_parentsize( handling_num )
542 
543  call nest_comm_catalogue( handling_num )
544  call mpi_barrier(intercomm_daughter, ierr)
545 
546  parent_ka(handling_num) = parent_kmax(handling_num) + khalo * 2
547  parent_ia(handling_num) = parent_imax(handling_num) + ihalo * 2
548  parent_ja(handling_num) = parent_jmax(handling_num) + jhalo * 2
549  daughter_ka(handling_num) = daughter_kmax(handling_num) + khalo * 2
550  daughter_ia(handling_num) = daughter_imax(handling_num) + ihalo * 2
551  daughter_ja(handling_num) = daughter_jmax(handling_num) + jhalo * 2
552  tileal_ka(handling_num) = 0
553  tileal_ia(handling_num) = 0
554  tileal_ja(handling_num) = 0
555 
556  if ( .NOT. online_no_rotate ) then
557  allocate( u_llp(parent_ka(handling_num), parent_ia(handling_num), parent_ja(handling_num) ) )
558  allocate( v_llp(parent_ka(handling_num), parent_ia(handling_num), parent_ja(handling_num) ) )
559  u_llp(:,:,:) = 0.0_rp
560  v_llp(:,:,:) = 0.0_rp
561  endif
562 
563  if( io_l ) write(io_fid_log,'(1x,A)' ) '*** Informations of Parent Domain [me]'
564  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- PARENT_PRC_nprocs :', parent_prc_nprocs(handling_num)
565  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- PARENT_PRC_NUM_X :', parent_prc_num_x(handling_num)
566  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- PARENT_PRC_NUM_Y :', parent_prc_num_y(handling_num)
567  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- PARENT_KMAX :', parent_kmax(handling_num)
568  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- PARENT_IMAX :', parent_imax(handling_num)
569  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- PARENT_JMAX :', parent_jmax(handling_num)
570  if( io_l ) write(io_fid_log,'(1x,A,F9.3)') '*** --- PARENT_DTSEC :', parent_dtsec(handling_num)
571  if( io_l ) write(io_fid_log,'(1x,A,I6) ') '*** --- PARENT_NSTEP :', parent_nstep(handling_num)
572  if( io_l ) write(io_fid_log,'(1x,A)' ) '*** Informations of Daughter Domain'
573  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- DAUGHTER_PRC_nprocs :', daughter_prc_nprocs(handling_num)
574  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- DAUGHTER_PRC_NUM_X :', daughter_prc_num_x(handling_num)
575  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- DAUGHTER_PRC_NUM_Y :', daughter_prc_num_y(handling_num)
576  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- DAUGHTER_KMAX :', daughter_kmax(handling_num)
577  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- DAUGHTER_IMAX :', daughter_imax(handling_num)
578  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- DAUGHTER_JMAX :', daughter_jmax(handling_num)
579  if( io_l ) write(io_fid_log,'(1x,A,F9.3)') '*** --- DAUGHTER_DTSEC :', daughter_dtsec(handling_num)
580  if( io_l ) write(io_fid_log,'(1x,A,I6) ') '*** --- DAUGHTER_NSTEP :', daughter_nstep(handling_num)
581  if( io_l ) write(io_fid_log,'(1x,A,I6) ') '*** Limit Num. NCOMM req. :', max_rq
582 
583  allocate( org_dens(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
584  allocate( org_momz(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
585  allocate( org_momx(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
586  allocate( org_momy(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
587  allocate( org_rhot(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
588  allocate( org_qtrc(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num),max_bndqa) )
589 
590  call nest_comm_setup_nestdown( handling_num )
591 
592  !---------------------------------- end of parent routines
593  endif
594 
595 
596  if( flag_child ) then
597  !-------------------------------------------------
598  if ( .NOT. online_iam_daughter ) then
599  write(*,*) 'xxx [grd_nest/NEST_setup] Child Flag from launcher is not consistent with namelist!'
600  write(*,*) 'xxx DAUGHTER - domain : ', online_domain_num
601  call prc_mpistop
602  endif
603 
604  handling_num = 2 !HANDLING_NUM + 1
605  intercomm_id(handling_num) = online_domain_num - 1
606  nest_filiation(intercomm_id(handling_num)) = -1
607 
608  intercomm_parent = inter_parent
609  if( io_l ) write(io_fid_log,'(1x,A,I2,A)') '*** Online Nesting - DAUGHTER [INTERCOMM_ID:', &
610  intercomm_id(handling_num), ' ]'
611  if( io_l ) write(io_fid_log,*) '*** Online Nesting - INTERCOMM :', intercomm_parent
612 
613  call nest_comm_ping( handling_num )
614 
615  call nest_comm_parentsize( handling_num )
616 
617  allocate( latlon_catalog(parent_prc_nprocs(handling_num),4,2) )
618  call nest_comm_catalogue( handling_num )
619  call mpi_barrier(intercomm_parent, ierr)
620 
621  call nest_domain_relate( handling_num )
622 
623  parent_ka(handling_num) = parent_kmax(handling_num) + khalo * 2
624  parent_ia(handling_num) = parent_imax(handling_num) + ihalo * 2
625  parent_ja(handling_num) = parent_jmax(handling_num) + jhalo * 2
626  daughter_ka(handling_num) = daughter_kmax(handling_num) + khalo * 2
627  daughter_ia(handling_num) = daughter_imax(handling_num) + ihalo * 2
628  daughter_ja(handling_num) = daughter_jmax(handling_num) + jhalo * 2
629  tileal_ka(handling_num) = parent_ka(handling_num)
630  tileal_ia(handling_num) = parent_imax(handling_num) * nest_tile_num_x
631  tileal_ja(handling_num) = parent_jmax(handling_num) * nest_tile_num_y
632 
633  if( io_l ) write(io_fid_log,'(1x,A)' ) '*** Informations of Parent Domain'
634  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- PARENT_PRC_nprocs :', parent_prc_nprocs(handling_num)
635  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- PARENT_PRC_NUM_X :', parent_prc_num_x(handling_num)
636  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- PARENT_PRC_NUM_Y :', parent_prc_num_y(handling_num)
637  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- PARENT_KMAX :', parent_kmax(handling_num)
638  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- PARENT_IMAX :', parent_imax(handling_num)
639  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- PARENT_JMAX :', parent_jmax(handling_num)
640  if( io_l ) write(io_fid_log,'(1x,A,F9.3)') '*** --- PARENT_DTSEC :', parent_dtsec(handling_num)
641  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- PARENT_NSTEP :', parent_nstep(handling_num)
642  if( io_l ) write(io_fid_log,'(1x,A)' ) '*** Informations of Daughter Domain [me]'
643  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- DAUGHTER_PRC_nprocs :', daughter_prc_nprocs(handling_num)
644  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- DAUGHTER_PRC_NUM_X :', daughter_prc_num_x(handling_num)
645  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- DAUGHTER_PRC_NUM_Y :', daughter_prc_num_y(handling_num)
646  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- DAUGHTER_KMAX :', daughter_kmax(handling_num)
647  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- DAUGHTER_IMAX :', daughter_imax(handling_num)
648  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- DAUGHTER_JMAX :', daughter_jmax(handling_num)
649  if( io_l ) write(io_fid_log,'(1x,A,F9.3)') '*** --- DAUGHTER_DTSEC :', daughter_dtsec(handling_num)
650  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- DAUGHTER_NSTEP :', daughter_nstep(handling_num)
651  if( io_l ) write(io_fid_log,'(1x,A)' ) '*** Informations of Target Tiles'
652  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- TILEALL_KA :', tileal_ka(handling_num)
653  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- TILEALL_IA :', tileal_ia(handling_num)
654  if( io_l ) write(io_fid_log,'(1x,A,I6)' ) '*** --- TILEALL_JA :', tileal_ja(handling_num)
655  if( io_l ) write(io_fid_log,'(1x,A,I6) ') '*** Limit Num. NCOMM req. :', max_rq
656 
657  allocate( buffer_2d( parent_ia(handling_num), parent_ja(handling_num) ) )
658  allocate( buffer_3d( parent_ka(handling_num), parent_ia(handling_num), parent_ja(handling_num) ) )
659  allocate( buffer_3df( 0:parent_ka(handling_num), parent_ia(handling_num), parent_ja(handling_num) ) )
660 
661  allocate( recvbuf_3d( parent_ka(handling_num), parent_ia(handling_num), parent_ja(handling_num), max_isu ) )
662  allocate( recvbuf_3df( 0:parent_ka(handling_num), parent_ia(handling_num), parent_ja(handling_num), max_isuf ) )
663 
664  allocate( buffer_ref_lon( tileal_ia(handling_num),tileal_ja(handling_num)) )
665  allocate( buffer_ref_lonx( tileal_ia(handling_num),tileal_ja(handling_num)) )
666  allocate( buffer_ref_lony( tileal_ia(handling_num),tileal_ja(handling_num)) )
667  allocate( buffer_ref_lat( tileal_ia(handling_num),tileal_ja(handling_num)) )
668  allocate( buffer_ref_latx( tileal_ia(handling_num),tileal_ja(handling_num)) )
669  allocate( buffer_ref_laty( tileal_ia(handling_num),tileal_ja(handling_num)) )
670  allocate( buffer_ref_cz( parent_ka(handling_num),tileal_ia(handling_num),tileal_ja(handling_num)) )
671  allocate( buffer_ref_fz(0:parent_ka(handling_num),tileal_ia(handling_num),tileal_ja(handling_num)) )
672 
673  !allocate( buffer_ref_2D ( TILEAL_IA(HANDLING_NUM),TILEAL_JA(HANDLING_NUM)) )
674  allocate( buffer_ref_3d( parent_ka(handling_num),tileal_ia(handling_num),tileal_ja(handling_num)) )
675  allocate( buffer_ref_3df(0:parent_ka(handling_num),tileal_ia(handling_num),tileal_ja(handling_num)) )
676 
677  allocate( hfact( daughter_ia(handling_num),daughter_ja(handling_num),itp_nh, itp_ng) )
678  allocate( vfact(daughter_ka(handling_num),daughter_ia(handling_num),daughter_ja(handling_num),itp_nh,itp_nv,itp_ng) )
679  allocate( igrd( daughter_ia(handling_num),daughter_ja(handling_num),itp_nh, itp_ng) )
680  allocate( jgrd( daughter_ia(handling_num),daughter_ja(handling_num),itp_nh, itp_ng) )
681  allocate( kgrd(daughter_ka(handling_num),daughter_ia(handling_num),daughter_ja(handling_num),itp_nh,itp_nv,itp_ng) )
682 
683  call nest_comm_setup_nestdown( handling_num )
684 
685 
686  ! for scalar points
687  call intrpnest_interp_fact_llz( hfact(:,:,:,i_sclr), & ! [OUT]
688  vfact(:,:,:,:,:,i_sclr), & ! [OUT]
689  kgrd(:,:,:,:,:,i_sclr), & ! [OUT]
690  igrd(:,:,:,i_sclr), & ! [OUT]
691  jgrd(:,:,:,i_sclr), & ! [OUT]
692  ncopy(:,:,:,i_sclr), & ! [OUT]
693  my_cz(:,:,:), & ! [IN]
694  my_lat(:,:), & ! [IN]
695  my_lon(:,:), & ! [IN]
696  datr_ks(handling_num), & ! [IN]
697  datr_ke(handling_num), & ! [IN]
698  daughter_ia(handling_num), & ! [IN]
699  daughter_ja(handling_num), & ! [IN]
700  buffer_ref_cz(:,:,:), & ! [IN]
701  buffer_ref_lat(:,:), & ! [IN]
702  buffer_ref_lon(:,:), & ! [IN]
703  tileal_ka(handling_num), & ! [IN]
704  tileal_ia(handling_num), & ! [IN]
705  tileal_ja(handling_num) ) ! [IN]
706 
707 
708  ! for z staggered points
709  call intrpnest_interp_fact_llz( hfact(:,:,:,i_zstg), & ! [OUT]
710  vfact(:,:,:,:,:,i_zstg), & ! [OUT]
711  kgrd(:,:,:,:,:,i_zstg), & ! [OUT]
712  igrd(:,:,:,i_zstg), & ! [OUT]
713  jgrd(:,:,:,i_zstg), & ! [OUT]
714  ncopy(:,:,:,i_zstg), & ! [OUT]
715  my_fz(:,:,:), & ! [IN]
716  my_lat(:,:), & ! [IN]
717  my_lon(:,:), & ! [IN]
718  datr_ks(handling_num), & ! [IN]
719  datr_ke(handling_num), & ! [IN]
720  daughter_ia(handling_num), & ! [IN]
721  daughter_ja(handling_num), & ! [IN]
722  buffer_ref_fz(:,:,:), & ! [IN]
723  buffer_ref_lat(:,:), & ! [IN]
724  buffer_ref_lon(:,:), & ! [IN]
725  tileal_ka(handling_num)+1, & ! [IN]
726  tileal_ia(handling_num), & ! [IN]
727  tileal_ja(handling_num) ) ! [IN]
728 
729 
730  ! for x staggered points
731  call intrpnest_interp_fact_llz( hfact(:,:,:,i_xstg), & ! [OUT]
732  vfact(:,:,:,:,:,i_xstg), & ! [OUT]
733  kgrd(:,:,:,:,:,i_xstg), & ! [OUT]
734  igrd(:,:,:,i_xstg), & ! [OUT]
735  jgrd(:,:,:,i_xstg), & ! [OUT]
736  ncopy(:,:,:,i_xstg), & ! [OUT]
737  my_cz(:,:,:), & ! [IN]
738  my_latx(1:ia,:), & ! [IN]
739  my_lonx(1:ia,:), & ! [IN]
740  datr_ks(handling_num), & ! [IN]
741  datr_ke(handling_num), & ! [IN]
742  daughter_ia(handling_num), & ! [IN]
743  daughter_ja(handling_num), & ! [IN]
744  buffer_ref_cz(:,:,:), & ! [IN]
745  buffer_ref_latx(:,:), & ! [IN]
746  buffer_ref_lonx(:,:), & ! [IN]
747  tileal_ka(handling_num), & ! [IN]
748  tileal_ia(handling_num), & ! [IN]
749  tileal_ja(handling_num) ) ! [IN]
750 
751  ! for y staggered points
752  call intrpnest_interp_fact_llz( hfact(:,:,:,i_ystg), & ! [OUT]
753  vfact(:,:,:,:,:,i_ystg), & ! [OUT]
754  kgrd(:,:,:,:,:,i_ystg), & ! [OUT]
755  igrd(:,:,:,i_ystg), & ! [OUT]
756  jgrd(:,:,:,i_ystg), & ! [OUT]
757  ncopy(:,:,:,i_ystg), & ! [OUT]
758  my_cz(:,:,:), & ! [IN]
759  my_laty(:,1:ja), & ! [IN]
760  my_lony(:,1:ja), & ! [IN]
761  datr_ks(handling_num), & ! [IN]
762  datr_ke(handling_num), & ! [IN]
763  daughter_ia(handling_num), & ! [IN]
764  daughter_ja(handling_num), & ! [IN]
765  buffer_ref_cz(:,:,:), & ! [IN]
766  buffer_ref_laty(:,:), & ! [IN]
767  buffer_ref_lony(:,:), & ! [IN]
768  tileal_ka(handling_num), & ! [IN]
769  tileal_ia(handling_num), & ! [IN]
770  tileal_ja(handling_num) ) ! [IN]
771 
772  deallocate( buffer_2d )
773  deallocate( buffer_3d )
774  deallocate( buffer_3df )
775 
776  !---------------------------------- end of daughter routines
777  else
778  online_use_velz = .false.
779  endif
780 
781  !if( IO_L ) write(IO_FID_LOG,'(1x,A,I2)') '*** Number of Related Domains :', HANDLING_NUM
782  !if ( HANDLING_NUM > 2 ) then
783  ! f( IO_L ) write(*,*) 'xxx Too much handing domains (up to 2)'
784  ! call PRC_MPIstop
785  !endif
786 
787  endif !--- OFFLINE or NOT
788 
789  endif !--- USE_NESTING
790 
791  return
module GTOOL_FILE
Definition: gtool_file.f90:17
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
module ATMOSPHERE / Physics Cloud Microphysics
logical, public prc_has_e
real(rp), dimension(:,:,:), allocatable, public real_fz
geopotential height [m] (cell face )
real(rp), dimension(:,:,:), allocatable, public real_cz
geopotential height [m] (cell center)
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)
real(rp), dimension(:,:), allocatable, public real_latx
latitude at staggered point (uy) [rad,-pi,pi]
module GRID (real space)
real(rp), dimension(:,:), allocatable, public real_lonxy
longitude at staggered point (uv) [rad,0-2pi]
module COMMUNICATION
Definition: scale_comm.F90:23
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
module INTERPOLATION (nesting system)
module RM PROCESS
integer, dimension(2), public parent_imax
parent max number in x-direction
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public prc_global_domainid
my domain ID in global communicator
real(rp), dimension(:,:), allocatable, public real_lony
longitude at staggered point (xv) [rad,0-2pi]
logical, public prc_has_w
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 filegetshape(dims, basename, varname, myrank, single, error)
real(rp), dimension(:,:), allocatable, public real_laty
latitude at staggered point (xv) [rad,-pi,pi]
subroutine, public intrpnest_setup(interp_search_divnum, NEST_INTERP_LEVEL, NEST_INTERP_WEIGHT_ORDER, OFFLINE)
Setup.
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 798 of file scale_grid_nest.F90.

References 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(), scale_process::prc_myrank, and use_nesting.

Referenced by nest_setup().

798  use scale_process, only: &
799  prc_myrank, &
801  implicit none
802 
803  integer, intent(in) :: HANDLE
804 
805  logical :: hit = .false.
806  integer, allocatable :: pd_tile_num(:,:)
807 
808  real(RP) :: wid_lon, wid_lat
809  integer :: pd_sw_tile
810  integer :: pd_ne_tile
811  integer :: i, j, k
812  !---------------------------------------------------------------------------
813 
814  if( .not. use_nesting ) then
815  return
816  end if
817 
818  allocate( pd_tile_num(0:parent_prc_nprocs(handle)-1,2) )
819 
820  k = 0 ! MPI process number starts from zero
821  do j = 1, parent_prc_num_y(handle)
822  do i = 1, parent_prc_num_x(handle)
823  pd_tile_num(k,1) = i
824  pd_tile_num(k,2) = j
825  k = k + 1
826  enddo
827  enddo
828 
829  !--- SW search
830  hit = .false.
831  do i = 1, parent_prc_nprocs(handle)
832  wid_lon = abs((latlon_catalog(i,i_sw,i_lon) - latlon_catalog(i,i_se,i_lon)) &
833  / real( PARENT_IMAX(HANDLE)-1, kind=rp )) * 0.8_rp
834  wid_lat = abs((latlon_catalog(i,i_sw,i_lat) - latlon_catalog(i,i_nw,i_lat)) &
835  / real( PARENT_JMAX(HANDLE)-1, kind=rp )) * 0.8_rp
836 
837  if ( corner_loc(i_sw,i_lon) >= min(latlon_catalog(i,i_sw,i_lon),latlon_catalog(i,i_nw,i_lon))-wid_lon .AND. &
838  corner_loc(i_sw,i_lat) >= min(latlon_catalog(i,i_sw,i_lat),latlon_catalog(i,i_se,i_lat))-wid_lat .AND. &
839  corner_loc(i_sw,i_lon) <= max(latlon_catalog(i,i_ne,i_lon),latlon_catalog(i,i_se,i_lon))+wid_lon .AND. &
840  corner_loc(i_sw,i_lat) <= max(latlon_catalog(i,i_ne,i_lat),latlon_catalog(i,i_nw,i_lat))+wid_lat ) then
841 
842  pd_sw_tile = i-1 ! MPI process number starts from zero
843  hit = .true.
844  exit ! exit loop
845  endif
846  enddo
847  if ( .NOT. hit ) then
848  write(*,*) 'xxx [grd_nest/NEST_domain_relate] region of daughter domain is larger than that of parent: SW search'
849  write(*,*) ' at rank:', prc_myrank, ' of domain:', online_domain_num
850  if( io_l ) write(io_fid_log,'(1x,A)') 'xxx region of daughter domain is larger than that of parent: SW search'
851  if( io_l ) write(io_fid_log,*) ' grid width: half width in lat:', wid_lat, ' half width in lon:', wid_lon
852  if( io_l ) write(io_fid_log,'(1x,A,F12.6)') ' daughter local (me): LON=',corner_loc(i_sw,i_lon)
853  do i = 1, parent_prc_nprocs(handle)
854  if( io_l ) write(io_fid_log,'(1x,A,F12.6,1x,F12.6)') ' parent local SW-NE: LON=', &
855  latlon_catalog(i,i_sw,i_lon) ,latlon_catalog(i,i_ne,i_lon)
856  enddo
857  if( io_l ) write(io_fid_log,'(1x,A,F12.6)') ' daughter local (me): LAT=',corner_loc(i_sw,i_lat)
858  do i = 1, parent_prc_nprocs(handle)
859  if( io_l ) write(io_fid_log,'(1x,A,F12.6,1x,F12.6)') ' parent local SW-NE: LAT=', &
860  latlon_catalog(i,i_sw,i_lat) ,latlon_catalog(i,i_ne,i_lat)
861  enddo
862  call prc_mpistop
863  endif
864 
865  !--- NE search
866  hit = .false.
867  do i = parent_prc_nprocs(handle), 1, -1
868  wid_lon = abs((latlon_catalog(i,i_nw,i_lon) - latlon_catalog(i,i_ne,i_lon)) &
869  / real( PARENT_IMAX(HANDLE)-1, kind=rp )) * 0.8_rp
870  wid_lat = abs((latlon_catalog(i,i_se,i_lat) - latlon_catalog(i,i_ne,i_lat)) &
871  / real( PARENT_JMAX(HANDLE)-1, kind=rp )) * 0.8_rp
872 
873  if ( corner_loc(i_ne,i_lon) >= min(latlon_catalog(i,i_sw,i_lon),latlon_catalog(i,i_nw,i_lon))-wid_lon .AND. &
874  corner_loc(i_ne,i_lat) >= min(latlon_catalog(i,i_sw,i_lat),latlon_catalog(i,i_se,i_lat))-wid_lat .AND. &
875  corner_loc(i_ne,i_lon) <= max(latlon_catalog(i,i_ne,i_lon),latlon_catalog(i,i_se,i_lon))+wid_lon .AND. &
876  corner_loc(i_ne,i_lat) <= max(latlon_catalog(i,i_ne,i_lat),latlon_catalog(i,i_nw,i_lat))+wid_lat ) then
877 
878  pd_ne_tile = i-1 ! MPI process number starts from zero
879  hit = .true.
880  exit ! exit loop
881  endif
882  enddo
883  if ( .NOT. hit ) then
884  write(*,*) 'xxx [grd_nest/NEST_domain_relate] region of daughter domain is larger than that of parent: NE search'
885  write(*,*) ' at rank:', prc_myrank, ' of domain:', online_domain_num
886  if( io_l ) write(io_fid_log,'(1x,A)') 'xxx region of daughter domain is larger than that of parent: NE search'
887  if( io_l ) write(io_fid_log,*) ' grid width: half width in lat:', wid_lat, ' half width in lon:', wid_lon
888  if( io_l ) write(io_fid_log,'(1x,A,F12.6)') ' daughter local (me): LON=',corner_loc(i_ne,i_lon)
889  do i = 1, parent_prc_nprocs(handle)
890  if( io_l ) write(io_fid_log,'(1x,A,F12.6,1x,F12.6)') ' parent local SW-NE: LON=', &
891  latlon_catalog(i,i_sw,i_lon) ,latlon_catalog(i,i_ne,i_lon)
892  enddo
893  if( io_l ) write(io_fid_log,'(1x,A,F12.6)') ' daughter local (me): LAT=',corner_loc(i_ne,i_lat)
894  do i = 1, parent_prc_nprocs(handle)
895  if( io_l ) write(io_fid_log,'(1x,A,F12.6,1x,F12.6)') ' parent local SW-NE: LAT=', &
896  latlon_catalog(i,i_sw,i_lat) ,latlon_catalog(i,i_ne,i_lat)
897  enddo
898  call prc_mpistop
899  endif
900 
901  nest_tile_num_x = pd_tile_num(pd_ne_tile,1) - pd_tile_num(pd_sw_tile,1) + 1
902  nest_tile_num_y = pd_tile_num(pd_ne_tile,2) - pd_tile_num(pd_sw_tile,2) + 1
903 
904  allocate( nest_tile_id( nest_tile_num_x*nest_tile_num_y ) )
905 
906  if( io_l ) write(io_fid_log,'(1x,A)') '*** NEST: target process tile in parent domain'
907  k = 1
908  do j = 1, nest_tile_num_y
909  do i = 1, nest_tile_num_x
910  nest_tile_id(k) = pd_sw_tile + (i-1) + parent_prc_num_x(handle)*(j-1)
911  if( io_l ) write(io_fid_log,'(1x,A,I4,A,I6)') ' (', k, ') target mpi-process:', nest_tile_id(k)
912  k = k + 1
913  enddo
914  enddo
915 
916  return
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
integer, public prc_myrank
process num in local communicator
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 930 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::ia, 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::ja, 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_boundary_diagqnum, 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, scale_time::time_nstep, and use_nesting.

Referenced by mod_copytopo::copytopo().

930  implicit none
931 
932  integer, intent(out) :: tilei, tilej
933  integer, intent(out) :: cxs, cxe, cys, cye
934  integer, intent(out) :: pxs, pxe, pys, pye
935  integer, intent(in) :: iloc ! rank number; start from 1
936 
937  integer :: hdl = 1 ! handler number
938  integer :: rank
939  integer :: xloc, yloc
940  integer :: xlocg, ylocg ! location over whole domain
941  !---------------------------------------------------------------------------
942 
943  if( .not. use_nesting ) then
944  return
945  end if
946 
947  rank = nest_tile_id(iloc)
948  xloc = mod( iloc-1, nest_tile_num_x ) + 1
949  yloc = int( real(iloc-1) / real(NEST_TILE_NUM_X) ) + 1
950  xlocg = mod( rank, offline_parent_prc_num_x ) + 1
951  ylocg = int( real(rank) / real(OFFLINE_PARENT_PRC_NUM_X) ) + 1
952  tilei = parent_imax(hdl)
953  tilej = parent_jmax(hdl)
954 
955  cxs = tilei * (xloc-1) + 1
956  cxe = tilei * xloc
957  cys = tilej * (yloc-1) + 1
958  cye = tilej * yloc
959  pxs = 1
960  pxe = tilei
961  pys = 1
962  pye = tilej
963 
964  if ( xlocg == 1 ) then ! BND_W
965  tilei = tilei + 2
966  pxs = pxs + 2
967  pxe = pxe + 2
968  endif
969  if ( xlocg == offline_parent_prc_num_x ) then ! BND_E
970  tilei = tilei + 2
971  endif
972  if ( ylocg == 1 ) then ! BND_S
973  tilej = tilej + 2
974  pys = pys + 2
975  pye = pye + 2
976  endif
977  if ( ylocg == offline_parent_prc_num_y ) then ! BND_N
978  tilej = tilej + 2
979  endif
980 
981  return
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 1695 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(), prnt_ie, prnt_is, prnt_je, prnt_js, prnt_ke, prnt_ks, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), and use_nesting.

Referenced by scale_atmos_boundary::atmos_boundary_update().

1695  use scale_process, only: &
1696  prc_mpistop
1697  use scale_comm, only: &
1698  comm_vars8, &
1699  comm_wait
1700  use scale_gridtrans, only: &
1701  rotc => gtrans_rotc
1702  implicit none
1703 
1704  integer, intent(in) :: HANDLE
1705  integer, intent(in) :: BND_QA
1706 
1707  real(RP), intent(in ) :: ipt_DENS(PARENT_KA(HANDLE),PARENT_IA(HANDLE),PARENT_JA(HANDLE))
1708  real(RP), intent(in ) :: ipt_MOMZ(PARENT_KA(HANDLE),PARENT_IA(HANDLE),PARENT_JA(HANDLE))
1709  real(RP), intent(in ) :: ipt_MOMX(PARENT_KA(HANDLE),PARENT_IA(HANDLE),PARENT_JA(HANDLE))
1710  real(RP), intent(in ) :: ipt_MOMY(PARENT_KA(HANDLE),PARENT_IA(HANDLE),PARENT_JA(HANDLE))
1711  real(RP), intent(in ) :: ipt_RHOT(PARENT_KA(HANDLE),PARENT_IA(HANDLE),PARENT_JA(HANDLE))
1712  real(RP), intent(in ) :: ipt_QTRC(PARENT_KA(HANDLE),PARENT_IA(HANDLE),PARENT_JA(HANDLE),BND_QA)
1713 
1714  real(RP), intent(inout) :: interped_ref_DENS(DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE))
1715  real(RP), intent(inout) :: interped_ref_VELZ(DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE))
1716  real(RP), intent(inout) :: interped_ref_VELX(DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE))
1717  real(RP), intent(inout) :: interped_ref_VELY(DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE))
1718  real(RP), intent(inout) :: interped_ref_POTT(DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE))
1719  real(RP), intent(inout) :: interped_ref_QTRC(DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE),BND_QA)
1720 
1721  real(RP) :: dummy(1,1,1)
1722  real(RP) :: dens (DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE))
1723  real(RP) :: u_lld(DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE))
1724  real(RP) :: v_lld(DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE))
1725  real(RP) :: work1d(DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE))
1726  real(RP) :: work2d(DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE))
1727 
1728  real(RP) :: work1p(PARENT_KA(HANDLE),PARENT_IA(HANDLE),PARENT_JA(HANDLE))
1729  real(RP) :: work2p(PARENT_KA(HANDLE),PARENT_IA(HANDLE),PARENT_JA(HANDLE))
1730 
1731  real(RP) :: u_on_map, v_on_map
1732  integer :: ierr
1733  integer :: tagbase, tagcomm
1734  integer :: isu_tag, isu_tagf
1735  integer :: i, j, k, iq
1736 
1737  integer, parameter :: cosin = 1
1738  integer, parameter :: sine = 2
1739  !---------------------------------------------------------------------------
1740 
1741  if( .not. use_nesting ) then
1742  return
1743  end if
1744 
1745  if ( bnd_qa > i_bndqa ) then
1746  write(*,*) 'xxx [grd_nest/NEST_COMM_nestdown] internal error: BND_QA is larger than I_BNDQA'
1747  call prc_mpistop
1748  endif
1749  if ( bnd_qa > max_bndqa ) then
1750  write(*,*) 'xxx [grd_nest/NEST_COMM_nestdown] internal error: BND_QA is larger than max_bndqa'
1751  call prc_mpistop
1752  endif
1753 
1754  tagcomm = intercomm_id(handle) * order_tag_comm
1755 
1756  if ( nest_filiation( intercomm_id(handle) ) > 0 ) then
1757  !-------------------------------------------------------- parent [send issue]
1758  call prof_rapstart('NEST_total_P', 2)
1759 
1760  nsend = nsend + 1
1761  if( io_l ) write(io_fid_log,'(1X,A,I5,A)') "*** NestIDC [P]: que send ( ", nsend, " )"
1762 
1763  call prof_rapstart('NEST_pack_P', 2)
1764 
1765  ! to keep values at that time by finish of sending process
1766  org_dens(:,:,:) = ipt_dens(:,:,:)
1767  org_momz(:,:,:) = ipt_momz(:,:,:)
1768  org_momx(:,:,:) = ipt_momx(:,:,:)
1769  org_momy(:,:,:) = ipt_momy(:,:,:)
1770  org_rhot(:,:,:) = ipt_rhot(:,:,:)
1771  do iq = 1, bnd_qa
1772  org_qtrc(:,:,:,iq) = ipt_qtrc(:,:,:,iq)
1773  enddo
1774 
1775  !*** request control
1776  !--- do not change the calling order below;
1777  !--- it should be consistent with the order in "NEST_COMM_recvwait_issue"
1778  rq_ctl_p = 0
1779 
1780  if ( .NOT. online_daughter_no_rotate ) then
1781  ! from staggered point to scalar point
1782  do j = 1, parent_ja(handle)
1783  do i = 2, parent_ia(handle)
1784  do k = prnt_ks(handle)-1, prnt_ke(handle)+1
1785  work1p(k,i,j) = ( org_momx(k,i-1,j) + org_momx(k,i,j) ) * 0.5_rp
1786  end do
1787  end do
1788  end do
1789  do j = 1, parent_ja(handle)
1790  do k = prnt_ks(handle)-1, prnt_ke(handle)+1
1791  work1p(k,1,j) = org_momx(k,1,j)
1792  end do
1793  end do
1794  call comm_vars8( work1p(:,:,:), 1 )
1795  do j = 2, parent_ja(handle)
1796  do i = 1, parent_ia(handle)
1797  do k = prnt_ks(handle)-1, prnt_ke(handle)+1
1798  work2p(k,i,j) = ( org_momy(k,i,j-1) + org_momy(k,i,j) ) * 0.5_rp
1799  end do
1800  end do
1801  end do
1802  do i = 1, parent_ia(handle)
1803  do k = prnt_ks(handle)-1, prnt_ke(handle)+1
1804  work2p(k,i,1) = org_momy(k,i,1)
1805  end do
1806  end do
1807  call comm_vars8( work2p(:,:,:), 2 )
1808  call comm_wait ( work1p(:,:,:), 1, .false. )
1809  call comm_wait ( work2p(:,:,:), 2, .false. )
1810 
1811  ! rotation from map-projected field to latlon field
1812  do j = prnt_js(handle), prnt_je(handle)
1813  do i = prnt_is(handle), prnt_ie(handle)
1814  do k = prnt_ks(handle), prnt_ke(handle)
1815  u_on_map = work1p(k,i,j) / org_dens(k,i,j)
1816  v_on_map = work2p(k,i,j) / org_dens(k,i,j)
1817 
1818  u_llp(k,i,j) = u_on_map * rotc(i,j,cosin) - v_on_map * rotc(i,j,sine )
1819  v_llp(k,i,j) = u_on_map * rotc(i,j,sine ) + v_on_map * rotc(i,j,cosin)
1820  enddo
1821  enddo
1822  enddo
1823  end if
1824 
1825  tagbase = tagcomm + tag_dens*order_tag_var
1826  call nest_comm_intercomm_nestdown( org_dens, dummy, tagbase, i_sclr, handle, isu_tag, isu_tagf, .true. )
1827 
1828  tagbase = tagcomm + tag_momz*order_tag_var
1829  if ( online_daughter_use_velz ) then
1830  call nest_comm_intercomm_nestdown( org_momz, dummy, tagbase, i_zstg, handle, isu_tag, isu_tagf )
1831  end if
1832 
1833  tagbase = tagcomm + tag_momx*order_tag_var
1834  if ( online_daughter_no_rotate ) then
1835  call nest_comm_intercomm_nestdown( org_momx, dummy, tagbase, i_xstg, handle, isu_tag, isu_tagf )
1836  else
1837  call nest_comm_intercomm_nestdown( u_llp, dummy, tagbase, i_sclr, handle, isu_tag, isu_tagf )
1838  end if
1839 
1840  tagbase = tagcomm + tag_momy*order_tag_var
1841  if ( online_daughter_no_rotate ) then
1842  call nest_comm_intercomm_nestdown( org_momy, dummy, tagbase, i_ystg, handle, isu_tag, isu_tagf )
1843  else
1844  call nest_comm_intercomm_nestdown( v_llp, dummy, tagbase, i_sclr, handle, isu_tag, isu_tagf )
1845  end if
1846 
1847  tagbase = tagcomm + tag_rhot*order_tag_var
1848  call nest_comm_intercomm_nestdown( org_rhot, dummy, tagbase, i_sclr, handle, isu_tag, isu_tagf )
1849 
1850  do iq = 1, bnd_qa
1851  tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
1852  call nest_comm_intercomm_nestdown( org_qtrc(:,:,:,iq), dummy, tagbase, i_sclr, handle, isu_tag, isu_tagf )
1853  enddo
1854 
1855  rq_tot_p = rq_ctl_p
1856 
1857  call prof_rapend ('NEST_pack_P', 2)
1858 
1859  call prof_rapend ('NEST_total_P', 2)
1860 
1861  elseif ( nest_filiation( intercomm_id(handle) ) < 0 ) then
1862  !-------------------------------------------------------- daughter [wait issue]
1863  call prof_rapstart('NEST_total_C', 2)
1864 
1865  nwait_d = nwait_d + 1
1866  !if( IO_L ) write(IO_FID_LOG,'(1X,A,I5,A)') "*** NestIDC [C]: que wait ( ", nwait_d, " )"
1867 
1868  !*** reset issue tag and request control
1869  !--- do not change the calling order below;
1870  !--- it should be consistent with the order in "NEST_COMM_recvwait_issue"
1871  isu_tag = 0
1872  isu_tagf = 0
1873 
1874  call prof_rapstart('NEST_wait_C', 2)
1875  call nest_comm_waitall( rq_tot_d, ireq_d )
1876 
1877  if ( online_aggressive_comm ) then
1878  ! nothing to do
1879  else
1880  call mpi_barrier(intercomm_parent, ierr)
1881  endif
1882  call prof_rapend ('NEST_wait_C', 2)
1883 
1884  call prof_rapstart('NEST_unpack_C', 2)
1885 
1886  tagbase = tagcomm + tag_dens*order_tag_var
1887  call nest_comm_intercomm_nestdown( dummy, dens, tagbase, i_sclr, handle, isu_tag, isu_tagf, .true. )
1888 !OCL XFILL
1889  do j = 1, daughter_ja(handle)
1890  do i = 1, daughter_ia(handle)
1891  do k = datr_ks(handle), datr_ke(handle)
1892  interped_ref_dens(k,i,j) = dens(k,i,j)
1893  enddo
1894  enddo
1895  enddo
1896 
1897  call comm_vars8( interped_ref_dens, 1 )
1898 
1899  tagbase = tagcomm + tag_momz*order_tag_var
1900  if ( online_use_velz ) then
1901  call nest_comm_intercomm_nestdown( dummy, work1d, tagbase, i_zstg, handle, isu_tag, isu_tagf )
1902 !OCL XFILL
1903  do j = 1, daughter_ja(handle)
1904  do i = 1, daughter_ia(handle)
1905  do k = datr_ks(handle), datr_ke(handle)-1
1906  interped_ref_velz(k,i,j) = work1d(k,i,j) / ( dens(k,i,j) + dens(k+1,i,j) ) * 2.0_rp
1907  enddo
1908  enddo
1909  enddo
1910  end if
1911 
1912  tagbase = tagcomm + tag_momx*order_tag_var
1913  if ( online_no_rotate ) then
1914  ! u_lld receives MOMX
1915  call nest_comm_intercomm_nestdown( dummy, u_lld, tagbase, i_xstg, handle, isu_tag, isu_tagf )
1916  else
1917  ! u_lld receives MOMX/DENS
1918  call nest_comm_intercomm_nestdown( dummy, u_lld, tagbase, i_sclr, handle, isu_tag, isu_tagf )
1919  endif
1920 
1921  tagbase = tagcomm + tag_momy*order_tag_var
1922  if ( online_no_rotate ) then
1923  ! v_lld receives MOMY
1924  call nest_comm_intercomm_nestdown( dummy, v_lld, tagbase, i_ystg, handle, isu_tag, isu_tagf )
1925  else
1926  ! v_lld receives MOMY/DENS
1927  call nest_comm_intercomm_nestdown( dummy, v_lld, tagbase, i_sclr, handle, isu_tag, isu_tagf )
1928  endif
1929 
1930  call comm_wait ( interped_ref_dens, 1, .false. )
1931 
1932  if ( online_no_rotate ) then
1933  do j = 1, daughter_ja(handle)
1934  do i = 1, daughter_ia(handle)-1
1935  do k = datr_ks(handle), datr_ke(handle)
1936  interped_ref_velx(k,i,j) = u_lld(k,i,j) &
1937  / ( interped_ref_dens(k,i+1,j) + interped_ref_dens(k,i,j) ) * 2.0_rp
1938  enddo
1939  enddo
1940  enddo
1941  do j = 1, daughter_ja(handle)
1942  do k = datr_ks(handle), datr_ke(handle)
1943  interped_ref_velx(k,daughter_ia(handle),j) = u_lld(k,daughter_ia(handle),j) &
1944  / interped_ref_dens(k,daughter_ia(handle),j)
1945  enddo
1946  enddo
1947  call comm_vars8( interped_ref_velx, 2 )
1948  do j = 1, daughter_ja(handle)-1
1949  do i = 1, daughter_ia(handle)
1950  do k = datr_ks(handle), datr_ke(handle)
1951  interped_ref_vely(k,i,j) = v_lld(k,i,j) &
1952  / ( interped_ref_dens(k,i,j+1) + interped_ref_dens(k,i,j) ) * 2.0_rp
1953  enddo
1954  enddo
1955  enddo
1956  do i = 1, daughter_ia(handle)
1957  do k = datr_ks(handle), datr_ke(handle)
1958  interped_ref_vely(k,i,daughter_ja(handle)) = v_lld(k,i,daughter_ja(handle)) &
1959  / interped_ref_dens(k,i,daughter_ja(handle))
1960  enddo
1961  enddo
1962  call comm_vars8( interped_ref_vely, 3 )
1963  call comm_wait ( interped_ref_velx, 2, .false. )
1964  call comm_wait ( interped_ref_vely, 3, .false. )
1965  else ! rotate
1966  ! rotation from latlon field to map-projected field
1967  do j = 1, daughter_ja(handle)
1968  do i = 1, daughter_ia(handle)
1969  do k = datr_ks(handle), datr_ke(handle)
1970  work1d(k,i,j) = u_lld(k,i,j) * rotc(i,j,cosin) + v_lld(k,i,j) * rotc(i,j,sine )
1971  work2d(k,i,j) = - u_lld(k,i,j) * rotc(i,j,sine ) + v_lld(k,i,j) * rotc(i,j,cosin)
1972  enddo
1973  enddo
1974  enddo
1975 
1976  ! from scalar point to staggered point
1977  do j = 1, daughter_ja(handle)
1978  do i = 1, daughter_ia(handle)-1
1979  do k = datr_ks(handle), datr_ke(handle)
1980  interped_ref_velx(k,i,j) = ( work1d(k,i+1,j) + work1d(k,i,j) ) * 0.5_rp
1981  end do
1982  end do
1983  end do
1984  do j = 1, daughter_ja(handle)
1985  do k = datr_ks(handle), datr_ke(handle)
1986  interped_ref_velx(k,daughter_ia(handle),j) = work1d(k,daughter_ia(handle),j)
1987  enddo
1988  enddo
1989  call comm_vars8( interped_ref_velx, 2 )
1990  do j = 1, daughter_ja(handle)-1
1991  do i = 1, daughter_ia(handle)
1992  do k = datr_ks(handle), datr_ke(handle)
1993  interped_ref_vely(k,i,j) = ( work2d(k,i,j+1) + work2d(k,i,j) ) * 0.5_rp
1994  end do
1995  end do
1996  end do
1997  do i = 1, daughter_ia(handle)
1998  do k = datr_ks(handle), datr_ke(handle)
1999  interped_ref_vely(k,i,daughter_ja(handle)) = work2d(k,i,daughter_ja(handle))
2000  enddo
2001  enddo
2002  call comm_vars8( interped_ref_vely, 3 )
2003  call comm_wait ( interped_ref_velx, 2, .false. )
2004  call comm_wait ( interped_ref_vely, 3, .false. )
2005  end if
2006 
2007  tagbase = tagcomm + tag_rhot*order_tag_var
2008  call nest_comm_intercomm_nestdown( dummy, work1d, tagbase, i_sclr, handle, isu_tag, isu_tagf )
2009  do j = 1, daughter_ja(handle)
2010  do i = 1, daughter_ia(handle)
2011  do k = datr_ks(handle), datr_ke(handle)
2012  interped_ref_pott(k,i,j) = work1d(k,i,j) / interped_ref_dens(k,i,j)
2013  enddo
2014  enddo
2015  enddo
2016 
2017  do iq = 1, bnd_qa
2018  tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
2019  call nest_comm_intercomm_nestdown( dummy, work1d, tagbase, i_sclr, handle, isu_tag, isu_tagf )
2020  do j = 1, daughter_ja(handle)
2021  do i = 1, daughter_ia(handle)
2022  do k = datr_ks(handle), datr_ke(handle)
2023  interped_ref_qtrc(k,i,j,iq) = work1d(k,i,j)
2024  enddo
2025  enddo
2026  enddo
2027  enddo
2028 
2029  call prof_rapend ('NEST_unpack_C', 2)
2030 
2031  call prof_rapend ('NEST_total_C', 2)
2032  else
2033  write(*,*) 'xxx [grd_nest/NEST_COMM_nestdown] internal error'
2034  call prc_mpistop
2035  endif
2036 
2037  return
subroutine, public prc_mpistop
Abort MPI.
module GRIDTRANS
module COMMUNICATION
Definition: scale_comm.F90:23
module PROCESS
real(rp), dimension(:,:,:), allocatable, public gtrans_rotc
rotation coefficient
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 2045 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(), scale_prof::prof_rapstart(), and use_nesting.

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().

2045  use scale_process, only: &
2046  prc_mpistop
2047  implicit none
2048 
2049  integer, intent(in) :: HANDLE
2050  integer, intent(in) :: BND_QA
2051 
2052  integer :: isu_tag, isu_tagf
2053  integer :: tagbase, tagcomm
2054  integer :: ierr
2055  integer :: iq
2056  !---------------------------------------------------------------------------
2057 
2058  if( .not. use_nesting ) then
2059  return
2060  end if
2061 
2062  if ( bnd_qa > i_bndqa ) then
2063  write(*,*) 'xxx [grd_nest/NEST_COMM_recvwait_issue] internal error: about BND_QA'
2064  call prc_mpistop
2065  endif
2066 
2067  tagcomm = intercomm_id(handle) * order_tag_comm
2068 
2069  if ( nest_filiation( intercomm_id(handle) ) > 0 ) then
2070  !-------------------------------------------------------- parent [wait issue]
2071  call prof_rapstart('NEST_total_P', 2)
2072  nwait_p = nwait_p + 1
2073  !if( IO_L ) write(IO_FID_LOG,'(1X,A,I5,A)') "*** NestIDC [P]: que wait ( ", nwait_p, " )"
2074 
2075  call prof_rapstart('NEST_wait_P', 2)
2076  call nest_comm_issuer_of_wait( handle )
2077 
2078  if ( online_aggressive_comm ) then
2079  ! nothing to do
2080  else
2081  call mpi_barrier(intercomm_daughter, ierr)
2082  endif
2083 
2084  call prof_rapend ('NEST_wait_P', 2)
2085 
2086  call prof_rapend('NEST_total_P', 2)
2087  elseif ( nest_filiation( intercomm_id(handle) ) < 0 ) then
2088  !-------------------------------------------------------- daughter [receive issue]
2089  call prof_rapstart('NEST_total_C', 2)
2090  nrecv = nrecv + 1
2091  if( io_l ) write(io_fid_log,'(1X,A,I5,A)') "*** NestIDC [C]: que recv ( ", nrecv, " )"
2092 
2093  !*** reset issue tag and request control
2094  !--- do not change the calling order below;
2095  !--- it should be consistent with the order in "NEST_COMM_nestdown"
2096  isu_tag = 0
2097  isu_tagf = 0
2098  rq_ctl_d = 0
2099 
2100  tagbase = tagcomm + tag_dens*order_tag_var
2101  call nest_comm_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf, .true. )
2102 
2103  tagbase = tagcomm + tag_momz*order_tag_var
2104  if ( online_use_velz ) then
2105  call nest_comm_issuer_of_receive( tagbase, i_zstg, handle, isu_tag, isu_tagf )
2106  end if
2107 
2108  tagbase = tagcomm + tag_momx*order_tag_var
2109  if ( online_no_rotate ) then
2110  call nest_comm_issuer_of_receive( tagbase, i_xstg, handle, isu_tag, isu_tagf )
2111  else
2112  call nest_comm_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf )
2113  endif
2114 
2115  tagbase = tagcomm + tag_momy*order_tag_var
2116  if ( online_no_rotate ) then
2117  call nest_comm_issuer_of_receive( tagbase, i_ystg, handle, isu_tag, isu_tagf )
2118  else
2119  call nest_comm_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf )
2120  endif
2121 
2122  tagbase = tagcomm + tag_rhot*order_tag_var
2123  call nest_comm_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf )
2124 
2125  do iq = 1, bnd_qa
2126  tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
2127  call nest_comm_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf )
2128  enddo
2129 
2130  rq_tot_d = rq_ctl_d
2131 
2132  call prof_rapend('NEST_total_C', 2)
2133  else
2134  write(*,*) 'xxx [grd_nest/NEST_COMM_recvwait_issue] internal error'
2135  call prc_mpistop
2136  endif
2137 
2138  return
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_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 2145 of file scale_grid_nest.F90.

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

Referenced by scale_atmos_boundary::atmos_boundary_finalize().

2145  use scale_process, only: &
2146  prc_mpistop
2147  implicit none
2148 
2149  integer, intent(in) :: HANDLE
2150 
2151  integer :: ierr
2152  !logical :: flag
2153  !integer :: istatus(MPI_STATUS_SIZE)
2154 
2155  integer :: rq
2156  !---------------------------------------------------------------------------
2157 
2158  if( .not. use_nesting ) then
2159  return
2160  end if
2161 
2162  if ( nest_filiation( intercomm_id(handle) ) > 0 ) then
2163  !-------------------------------------------------------- parent
2164  !--- Nothing to do
2165 
2166  elseif ( nest_filiation( intercomm_id(handle) ) < 0 ) then
2167  !-------------------------------------------------------- daughter [receive issue]
2168  if( io_l ) write(io_fid_log,'(1X,A,I5,A)') "*** NestIDC [C]: CANCEL recv ( ", nrecv, " )"
2169  do rq = 1, rq_tot_d
2170  if ( ireq_d(rq) /= mpi_request_null ) then
2171  call mpi_cancel(ireq_d(rq), ierr)
2172  !call MPI_TEST_CANCELLED(istatus, flag, ierr)
2173  !if ( .NOT. flag ) then
2174  ! write(*,*) 'xxx receive actions do not cancelled, req = ', rq
2175  !endif
2176  endif
2177  enddo
2178 
2179  else
2180  write(*,*) 'xxx [grd_nest/NEST_COMM_recv_cancel] internal error'
2181  call prc_mpistop
2182  endif
2183 
2184  return
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_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 2198 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, nest_filiation, nest_tile_num_x, parent_ia, parent_imax, parent_ja, parent_jmax, parent_ka, scale_process::prc_mpistop(), prnt_ie, prnt_is, prnt_je, prnt_js, and use_nesting.

2198  use scale_process, only: &
2199  prc_mpistop
2200  use scale_comm, only: &
2201  comm_datatype
2202  use scale_interpolation_nest, only: &
2204  implicit none
2205 
2206  real(RP), intent(in) :: pvar(:,:,:)
2207  real(RP), intent(out) :: dvar(:,:,:)
2208  integer , intent(in) :: tagbase
2209  integer , intent(in) :: id_stag
2210  integer , intent(in) :: HANDLE
2211  integer , intent(inout) :: isu_tag
2212  integer , intent(inout) :: isu_tagf
2213 
2214  logical , intent(in), optional :: flag_dens
2215 
2216  integer :: ierr, ileng
2217  integer :: tag, target_rank
2218 
2219  integer :: xloc, yloc
2220  integer :: xs, xe
2221  integer :: ys, ye
2222 
2223  integer :: k
2224  integer :: ig, rq, yp
2225  logical :: no_zstag = .true.
2226  logical :: logarithmic = .false.
2227  !---------------------------------------------------------------------------
2228 
2229  if( .not. use_nesting ) then
2230  return
2231  end if
2232 
2233  logarithmic = .false.
2234  if ( present(flag_dens) ) then
2235  if ( flag_dens ) then
2236  logarithmic = .true.
2237  endif
2238  endif
2239 
2240  if ( id_stag == i_sclr ) then
2241  no_zstag = .true.
2242  ig = i_sclr
2243  elseif ( id_stag == i_zstg ) then
2244  no_zstag = .false.
2245  ig = i_zstg
2246  elseif ( id_stag == i_xstg ) then
2247  no_zstag = .true.
2248  ig = i_xstg
2249  elseif ( id_stag == i_ystg ) then
2250  no_zstag = .true.
2251  ig = i_ystg
2252  endif
2253 
2254  if ( no_zstag ) then
2255  ileng = parent_ka(handle) * parent_ia(handle) * parent_ja(handle)
2256  else
2257  ileng = (parent_ka(handle)+1) * parent_ia(handle) * parent_ja(handle)
2258  endif
2259 
2260  if ( nest_filiation( intercomm_id(handle) ) > 0 ) then
2261  !--------------------------------------------------- parent
2262 
2263  rq = rq_ctl_p
2264 
2265  do yp = 1, num_yp
2266  rq = rq + 1
2267 
2268  ! send data to multiple daughter processes
2269  target_rank = nest_tile_list_yp(yp)
2270  tag = tagbase + yp
2271  call mpi_isend(pvar, ileng, comm_datatype, target_rank, tag, intercomm_daughter, ireq_p(rq), ierr)
2272 
2273  dvar(:,:,:) = -1.0_rp ! input as a dummy value
2274  enddo
2275 
2276  rq_ctl_p = rq
2277 
2278  elseif ( nest_filiation( intercomm_id(handle) ) < 0 ) then
2279  !--------------------------------------------------- daughter
2280 
2281  rq = rq_ctl_d
2282 
2283  do yp = 1, nest_tile_all ! YP Loop
2284  rq = rq + 1
2285 
2286  xloc = mod( yp-1, nest_tile_num_x ) + 1
2287  yloc = int( real(yp-1) / real(NEST_TILE_NUM_X) ) + 1
2288 
2289  xs = parent_imax(handle) * (xloc-1) + 1
2290  xe = parent_imax(handle) * xloc
2291  ys = parent_jmax(handle) * (yloc-1) + 1
2292  ye = parent_jmax(handle) * yloc
2293 
2294  if ( no_zstag ) then
2295  isu_tag = isu_tag + 1
2296 
2297  if ( .NOT. logarithmic ) then
2298 !OCL XFILL
2299  ! linear interpolation
2300  do k = 1, parent_ka(handle)
2301  buffer_ref_3d(k,xs:xe,ys:ye) &
2302  = recvbuf_3d(k,prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle),isu_tag)
2303  enddo
2304  else
2305 !OCL XFILL
2306  ! logarithmic weighted interpolation
2307  do k = 1, parent_ka(handle)
2308  buffer_ref_3d(k,xs:xe,ys:ye) &
2309  = log( recvbuf_3d(k,prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle),isu_tag) )
2310  enddo
2311  endif
2312  else
2313  isu_tagf = isu_tagf + 1
2314 !OCL XFILL
2315  do k = 0, parent_ka(handle)
2316  buffer_ref_3df(k,xs:xe,ys:ye) &
2317  = recvbuf_3df(k,prnt_is(handle):prnt_ie(handle),prnt_js(handle):prnt_je(handle),isu_tagf)
2318  enddo
2319  endif
2320 
2321  if ( isu_tag > max_isu .OR. isu_tagf > max_isuf ) then
2322  write(*,*) 'xxx [grd_nest/NEST_COMM_intercomm_nestdown_3D] Exceeded maximum issue'
2323  write(*,*) 'xxx isu_tag = ', isu_tag
2324  write(*,*) 'xxx isu_tagf = ', isu_tagf
2325  call prc_mpistop
2326  endif
2327 
2328  enddo ! YP Loop
2329  rq_ctl_d = rq
2330 
2331 !OCL XFILL
2332  dvar(:,:,:) = 0.0_rp
2333 
2334  if ( no_zstag ) then
2335  call intrpnest_interp_3d( dvar, &
2336  buffer_ref_3d, &
2337  hfact(:,:,:,ig), &
2338  vfact(:,:,:,:,:,ig), &
2339  kgrd(:,:,:,:,:,ig), &
2340  igrd(:,:,:,ig), &
2341  jgrd(:,:,:,ig), &
2342  daughter_ia(handle), &
2343  daughter_ja(handle), &
2344  datr_ks(handle), &
2345  datr_ke(handle), &
2346  logarithmic )
2347  else
2348  ! linear interpolation (z-staggered)
2349  call intrpnest_interp_3d( dvar, &
2350  buffer_ref_3df, &
2351  hfact(:,:,:,ig), &
2352  vfact(:,:,:,:,:,ig), &
2353  kgrd(:,:,:,:,:,ig), &
2354  igrd(:,:,:,ig), &
2355  jgrd(:,:,:,ig), &
2356  daughter_ia(handle), &
2357  daughter_ja(handle), &
2358  datr_ks(handle), &
2359  datr_ke(handle), &
2360  logarithmic )
2361  endif
2362  else
2363  !---------------------------------------------------
2364  write(*,*) 'xxx [grd_nest/NEST_COMM_intercomm_nestdown_3D] internal error'
2365  call prc_mpistop
2366  endif
2367 
2368  return
subroutine, public prc_mpistop
Abort MPI.
module COMMUNICATION
Definition: scale_comm.F90:23
module PROCESS
procedure(intrpnest_intfc_interp_3d), pointer, public intrpnest_interp_3d
module INTERPOLATION (nesting system)
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 2380 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 use_nesting.

2380  use scale_process, only: &
2381  prc_myrank, &
2382  prc_mpistop
2383  use scale_comm, only: &
2384  comm_datatype
2385  implicit none
2386 
2387  integer , intent(in) :: tagbase
2388  integer , intent(in) :: id_stag
2389  integer , intent(in) :: HANDLE
2390  integer , intent(inout) :: isu_tag
2391  integer , intent(inout) :: isu_tagf
2392  logical , intent(in), optional :: flag_dens
2393 
2394  integer :: ierr, ileng
2395  integer :: tag, target_rank
2396 
2397  integer :: ig, rq, yp
2398  logical :: no_zstag = .true.
2399  logical :: logarithmic = .false.
2400  !---------------------------------------------------------------------------
2401 
2402  if( .not. use_nesting ) then
2403  return
2404  end if
2405 
2406  logarithmic = .false.
2407  if ( present(flag_dens) ) then
2408  if ( flag_dens ) then
2409  logarithmic = .true.
2410  endif
2411  endif
2412 
2413  if ( id_stag == i_sclr ) then
2414  no_zstag = .true.
2415  ig = i_sclr
2416  elseif ( id_stag == i_zstg ) then
2417  no_zstag = .false.
2418  ig = i_zstg
2419  elseif ( id_stag == i_xstg ) then
2420  no_zstag = .true.
2421  ig = i_xstg
2422  elseif ( id_stag == i_ystg ) then
2423  no_zstag = .true.
2424  ig = i_ystg
2425  endif
2426 
2427  if ( no_zstag ) then
2428  ileng = parent_ka(handle) * parent_ia(handle) * parent_ja(handle)
2429  else
2430  ileng = (parent_ka(handle)+1) * parent_ia(handle) * parent_ja(handle)
2431  endif
2432 
2433  if ( nest_filiation( intercomm_id(handle) ) > 0 ) then
2434  !--------------------------------------------------- parent
2435  ! nothing to do
2436  ! rq = rq_ctl_p
2437  ! rq_ctl_p = rq
2438  elseif ( nest_filiation( intercomm_id(handle) ) < 0 ) then
2439  !--------------------------------------------------- daughter
2440 
2441  rq = rq_ctl_d
2442 
2443  do yp = 1, nest_tile_all ! YP Loop
2444  rq = rq + 1
2445 
2446  target_rank = nest_tile_list_d(yp,prc_myrank+1)
2447 
2448  tag = tagbase + call_order(yp)
2449  if ( no_zstag ) then
2450  isu_tag = isu_tag + 1
2451  recvbuf_3d(:,:,:,isu_tag) = 0.0_rp
2452  call mpi_irecv( recvbuf_3d(:,:,:,isu_tag), &
2453  ileng, &
2454  comm_datatype, &
2455  target_rank, &
2456  tag, &
2457  intercomm_parent, &
2458  ireq_d(rq), &
2459  ierr )
2460  else
2461  isu_tagf = isu_tagf + 1
2462  recvbuf_3df(:,:,:,isu_tagf) = 0.0_rp
2463  call mpi_irecv( recvbuf_3df(:,:,:,isu_tagf), &
2464  ileng, &
2465  comm_datatype, &
2466  target_rank, &
2467  tag, &
2468  intercomm_parent, &
2469  ireq_d(rq), &
2470  ierr )
2471  endif
2472 
2473  enddo ! YP Loop
2474 
2475  if ( isu_tag > max_isu .OR. isu_tagf > max_isuf ) then
2476  write(*,*) 'xxx [grd_nest/NEST_COMM_issuer_of_receive_3D] Exceeded maximum issue'
2477  write(*,*) 'xxx isu_tag = ', isu_tag
2478  write(*,*) 'xxx isu_tagf = ', isu_tagf
2479  call prc_mpistop
2480  endif
2481 
2482  rq_ctl_d = rq
2483 
2484  else
2485  !---------------------------------------------------
2486  write(*,*) 'xxx [grd_nest/NEST_COMM_issuer_of_receive_3D] internal error'
2487  call prc_mpistop
2488  endif
2489 
2490  return
subroutine, public prc_mpistop
Abort MPI.
module COMMUNICATION
Definition: scale_comm.F90:23
module PROCESS
integer, public prc_myrank
process num 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 2497 of file scale_grid_nest.F90.

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

2497  use scale_process, only: &
2498  prc_mpistop
2499  implicit none
2500 
2501  integer, intent(in) :: HANDLE
2502  !---------------------------------------------------------------------------
2503 
2504  if( .not. use_nesting ) then
2505  return
2506  end if
2507 
2508  if ( nest_filiation( intercomm_id(handle) ) > 0 ) then
2509  !--------------------------------------------------- parent
2510 
2511  call nest_comm_waitall( rq_tot_p, ireq_p )
2512 
2513  elseif ( nest_filiation( intercomm_id(handle) ) < 0 ) then
2514  !--------------------------------------------------- daughter
2515  ! nothing to do
2516  else
2517  !---------------------------------------------------
2518  write(*,*) 'xxx [grd_nest/NEST_COMM_issuer_of_wait_3D] internal error'
2519  call prc_mpistop
2520  endif
2521 
2522  return
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 2530 of file scale_grid_nest.F90.

References scale_process::prc_mpistop(), and use_nesting.

Referenced by nest_comm_issuer_of_wait_3d(), and nest_comm_nestdown().

2530  use scale_process, only: &
2531  prc_mpistop
2532  implicit none
2533 
2534  integer, intent(in) :: req_count
2535  integer, intent(inout) :: ireq(max_rq)
2536 
2537  integer :: i
2538  integer :: ierr
2539  integer :: istatus(MPI_STATUS_SIZE,req_count)
2540  integer :: req_count2
2541  integer :: ireq2(max_rq)
2542 
2543 ! logical :: flag
2544 ! integer(8) :: num
2545  !---------------------------------------------------------------------------
2546 ! num = 0
2547 ! flag = .false.
2548 
2549  if( .not. use_nesting ) then
2550  return
2551  end if
2552 
2553  req_count2 = 0
2554  do i=1, req_count
2555  if (ireq(i) /= mpi_request_null) then
2556  req_count2 = req_count2 + 1
2557  ireq2(req_count2) = ireq(i)
2558  endif
2559  enddo
2560  if ( req_count2 /= 0 ) call mpi_waitall( req_count2, ireq2, istatus, ierr )
2561 
2562 ! do while ( .NOT. flag )
2563 ! num = num + 1
2564 ! call MPI_TESTALL( req_count, ireq, flag, istatus, ierr )
2565 
2566 ! if ( num > ONLINE_WAIT_LIMIT ) then
2567 ! if( IO_L ) write(IO_FID_LOG,'(1x,A)') '*** ERROR: over the limit of waiting time [NESTCOM]'
2568 ! write(*,'(1x,A)') '*** ERROR: over the limit of waiting time [NESTCOM]'
2569 ! call PRC_MPIstop
2570 ! endif
2571 ! enddo
2572 
2573  return
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 2580 of file scale_grid_nest.F90.

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

Referenced by scale_atmos_boundary::atmos_boundary_update().

2580  use scale_process, only: &
2581  prc_mpistop
2582  implicit none
2583 
2584  integer, intent(in) :: HANDLE
2585 
2586  integer :: istatus(MPI_STATUS_SIZE)
2587  integer :: ierr
2588  logical :: flag
2589  !---------------------------------------------------------------------------
2590 
2591  if( .not. use_nesting ) then
2592  return
2593  end if
2594 
2595  if ( nest_filiation( intercomm_id(handle) ) > 0 ) then
2596  !--------------------------------------------------- parent
2597  call prof_rapstart('NEST_test_P', 2)
2598  if ( rq_ctl_p > 0 ) call mpi_test(ireq_p(1), flag, istatus, ierr)
2599  call prof_rapend('NEST_test_P', 2)
2600 
2601  elseif ( nest_filiation( intercomm_id(handle) ) < 0 ) then
2602  !--------------------------------------------------- daughter
2603  call prof_rapstart('NEST_test_C', 2)
2604  if ( rq_ctl_d > 0 ) call mpi_test(ireq_d(1), flag, istatus, ierr)
2605  call prof_rapend('NEST_test_C', 2)
2606  else
2607  !---------------------------------------------------
2608  write(*,*) 'xxx [grd_nest/NEST_COMM_test] error'
2609  call prc_mpistop
2610  endif
2611 
2612  return
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_disconnect()

subroutine, public scale_grid_nest::nest_comm_disconnect ( )

[finalize: disconnect] Inter-communication

Definition at line 2618 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, scale_process::prc_global_comm_world, and use_nesting.

Referenced by mod_atmos_driver::atmos_driver_finalize().

2618  use scale_process, only: &
2620  implicit none
2621 
2622  integer :: ierr
2623  !---------------------------------------------------------------------------
2624 
2625  if( .not. use_nesting ) then
2626  return
2627  end if
2628 
2629  if( io_l ) write(io_fid_log,'(1x,A)') '*** Waiting finish of whole processes'
2630  call mpi_barrier(prc_global_comm_world, ierr)
2631 
2632  if ( online_iam_parent ) then
2633  !if( IO_L ) write(IO_FID_LOG,'(1x,A)') '*** Waiting finish of whole processes as a parent'
2634  !call MPI_BARRIER(INTERCOMM_DAUGHTER, ierr)
2635  call mpi_comm_free(intercomm_daughter, ierr)
2636  if( io_l ) write(io_fid_log,'(1x,A)') '*** Disconnected communication with child'
2637  endif
2638 
2639  if ( online_iam_daughter ) then
2640  !if( IO_L ) write(IO_FID_LOG,'(1x,A)') '*** Waiting finish of whole processes as a child'
2641  !call MPI_BARRIER(INTERCOMM_PARENT, ierr)
2642  call mpi_comm_free(intercomm_parent, ierr)
2643  if( io_l ) write(io_fid_log,'(1x,A)') '*** Disconnected communication with parent'
2644  endif
2645 
2646  return
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 = 5

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 = 5

◆ nest_interp_weight_order

integer, public scale_grid_nest::nest_interp_weight_order = 2

horizontal interpolation weight order

Definition at line 104 of file scale_grid_nest.F90.

Referenced by nest_setup().

104  integer, public :: NEST_INTERP_WEIGHT_ORDER = 2

◆ use_nesting

logical, public scale_grid_nest::use_nesting = .false.

◆ offline

logical, public scale_grid_nest::offline = .false.

Definition at line 107 of file scale_grid_nest.F90.

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

107  logical, public :: OFFLINE = .false.

◆ online_iam_parent

logical, public scale_grid_nest::online_iam_parent = .false.

a flag to say "I am a parent"

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_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 109 of file scale_grid_nest.F90.

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

109  logical, public :: ONLINE_IAM_DAUGHTER = .false.

◆ online_domain_num

integer, public scale_grid_nest::online_domain_num = 1

Definition at line 110 of file scale_grid_nest.F90.

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

110  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 112 of file scale_grid_nest.F90.

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

112  logical, public :: ONLINE_NO_ROTATE = .false.

◆ online_boundary_use_qhyd

logical, public scale_grid_nest::online_boundary_use_qhyd = .false.

Definition at line 113 of file scale_grid_nest.F90.

Referenced by nest_setup().

113  logical, public :: ONLINE_BOUNDARY_USE_QHYD = .false.

◆ online_boundary_diagqnum

logical, public scale_grid_nest::online_boundary_diagqnum = .false.

Definition at line 114 of file scale_grid_nest.F90.

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

114  logical, public :: ONLINE_BOUNDARY_DIAGQNUM = .false.