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

module Communication CartesianC nesting More...

Functions/Subroutines

subroutine, public comm_cartesc_nest_setup (QA_MP, MP_TYPE_in, inter_parent, inter_child)
 Setup. More...
 
subroutine, public comm_cartesc_nest_domain_relate (HANDLE)
 Solve relationship between ParentDomain & Daughter Domain. More...
 
subroutine, public comm_cartesc_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 comm_cartesc_nest_nestdown (HANDLE, BND_QA, DENS_send, MOMZ_send, MOMX_send, MOMY_send, RHOT_send, QTRC_send, DENS_recv, VELZ_recv, VELX_recv, VELY_recv, POTT_recv, QTRC_recv)
 Boundary data transfer from parent to daughter: nestdown. More...
 
subroutine, public comm_cartesc_nest_recvwait_issue (HANDLE, BND_QA)
 Sub-command for data transfer from parent to daughter: nestdown. More...
 
subroutine, public comm_cartesc_nest_recv_cancel (HANDLE)
 Sub-command for data transfer from parent to daughter: nestdown. More...
 
subroutine comm_cartesc_nest_intercomm_nestdown_3d (pvar, dvar, tagbase, id_stag, HANDLE, isu_tag, flag_dens)
 Inter-communication from parent to daughter: nestdown. More...
 
subroutine comm_cartesc_nest_issuer_of_receive_3d (tagbase, id_stag, HANDLE, isu_tag)
 [substance of issuer] Inter-communication from parent to daughter: nestdown More...
 
subroutine comm_cartesc_nest_issuer_of_wait_3d (HANDLE)
 [substance of issuer] Inter-communication from parent to daughter: nestdown More...
 
subroutine comm_cartesc_nest_waitall (req_count, ireq)
 [substance of comm_wait] Inter-communication More...
 
subroutine, public comm_cartesc_nest_test (HANDLE)
 [check communication status] Inter-communication More...
 
subroutine, public comm_cartesc_nest_disconnect
 [finalize: disconnect] Inter-communication More...
 

Variables

integer, public intercomm_parent
 
integer, public intercomm_daughter
 
integer, dimension(10), public comm_cartesc_nest_filiation
 index of parent-daughter relation (p>0, d<0) More...
 
integer, public handling_num
 handing number of nesting relation More...
 
integer, public comm_cartesc_nest_tile_num_x
 parent tile number in x-direction More...
 
integer, public comm_cartesc_nest_tile_num_y
 parent tile number in y-direction More...
 
integer, dimension(:), allocatable, public comm_cartesc_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_okmax
 parent max number in oz-direction 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_okmax
 daughter max number in oz-direction 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 comm_cartesc_nest_bnd_qa = 1
 number of tracer treated in nesting system More...
 
integer, public comm_cartesc_nest_interp_level = 5
 horizontal interpolation level More...
 
integer, public comm_cartesc_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_diagqhyd = .false.
 

Detailed Description

module Communication CartesianC nesting

Description
Grid module for nesting system
Author
Team SCALE
NAMELIST
  • PARAM_COMM_CARTESC_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_AGGRESSIVE_COMM logical
    ONLINE_WAIT_LIMIT integer(8) limit times of waiting loop in "COMM_CARTESC_NEST_waitall"
    ONLINE_SPECIFIED_MAXRQ integer 0
    COMM_CARTESC_NEST_INTERP_TYPE character(len=H_SHORT) 'LINEAR' ! "LINEAR" or "DIST-WEIGHT"
    COMM_CARTESC_NEST_INTERP_LEVEL integer 5 horizontal interpolation level
    COMM_CARTESC_NEST_INTERP_WEIGHT_ORDER integer 2 horizontal interpolation weight order

History Output
No history output

Function/Subroutine Documentation

◆ comm_cartesc_nest_setup()

subroutine, public scale_comm_cartesc_nest::comm_cartesc_nest_setup ( integer, intent(in)  QA_MP,
character(len=*), intent(in)  MP_TYPE_in,
integer, intent(in), optional  inter_parent,
integer, intent(in), optional  inter_child 
)

Setup.

Definition at line 261 of file scale_comm_cartesC_nest.F90.

261  use scale_file, only: &
262  file_open, &
263  file_get_attribute, &
264  file_get_shape
265  use scale_const, only: &
266  d2r => const_d2r
267  use scale_prc, only: &
268  prc_abort, &
271  use scale_interp, only: &
272  interp_setup, &
273  interp_factor3d
274  use scale_comm_cartesc, only: &
275  comm_bcast
276  use scale_atmos_grid_cartesc, only: &
281  use scale_atmos_grid_cartesc_real, only: &
292  use scale_atmos_hydrometeor, only: &
294  use scale_mapprojection, only: &
295  mapprojection_lonlat2xy
296  implicit none
297 
298  integer, intent(in) :: QA_MP
299  character(len=*), intent(in) :: MP_TYPE_in
300  integer, intent(in), optional :: inter_parent
301  integer, intent(in), optional :: inter_child
302 
303  character(len=H_SHORT) :: COMM_CARTESC_NEST_INTERP_TYPE = 'LINEAR' ! "LINEAR" or "DIST-WEIGHT"
304  ! LINEAR : bi-linear interpolation
305  ! DIST-WEIGHT: distance-weighted mean of the nearest N-neighbors
306 
307 
308  character(len=H_LONG) :: LATLON_CATALOGUE_FNAME = 'latlon_domain_catalogue.txt'
309 
310  real(RP), allocatable :: X_ref(:,:)
311  real(RP), allocatable :: Y_ref(:,:)
312 
313  integer :: ONLINE_SPECIFIED_MAXRQ = 0
314  integer :: i
315  integer :: fid, ierr
316  integer :: parent_id
317 
318  logical :: flag_parent = .false.
319  logical :: flag_child = .false.
320 
321  integer :: imaxg(1), jmaxg(1)
322  integer :: pnum_x(1), pnum_y(1)
323  integer :: dims(1)
324  logical :: error, existed
325 
326  namelist / param_comm_cartesc_nest / &
327  latlon_catalogue_fname, &
328  offline_parent_basename, &
329  offline_parent_prc_num_x, &
330  offline_parent_prc_num_y, &
331  online_domain_num, &
332  online_iam_parent, &
333  online_iam_daughter, &
334  online_use_velz, &
335  online_no_rotate, &
336  online_boundary_use_qhyd, &
337  online_aggressive_comm, &
338  online_wait_limit, &
339  online_specified_maxrq, &
340  comm_cartesc_nest_interp_type, &
341  comm_cartesc_nest_interp_level, &
342  comm_cartesc_nest_interp_weight_order
343 
344  !---------------------------------------------------------------------------
345 
346  log_newline
347  log_info("COMM_CARTESC_NEST_setup",*) 'Setup'
348 
349  if ( present(inter_parent) ) then
350  if( inter_parent /= mpi_comm_null ) flag_child = .true. ! exist parent, so work as a child
351  endif
352  if ( present(inter_child) ) then
353  if( inter_child /= mpi_comm_null ) flag_parent = .true. ! exist child, so work as a parent
354  endif
355 
356  offline_parent_basename = ""
357 
358  nwait_p = 0
359  nwait_d = 0
360  nrecv = 0
361  nsend = 0
362 
363  handling_num = 0
364  comm_cartesc_nest_filiation(:) = 0
365  online_wait_limit = 999999999
366  online_aggressive_comm = .true.
367 
368  !--- read namelist
369  rewind(io_fid_conf)
370  read(io_fid_conf,nml=param_comm_cartesc_nest,iostat=ierr)
371  if( ierr < 0 ) then !--- missing
372  log_info("COMM_CARTESC_NEST_setup",*) 'Not found namelist. Default used.'
373  elseif( ierr > 0 ) then !--- fatal error
374  log_error("COMM_CARTESC_NEST_setup",*) 'Not appropriate names in namelist PARAM_COMM_CARTESC_NEST. Check!'
375  call prc_abort
376  endif
377  log_nml(param_comm_cartesc_nest)
378 
379  prc_global_domainid = online_domain_num
380 
381  if ( offline_parent_basename /= "" ) then
382 
383  offline = .true.
384  use_nesting = .true.
385 
386  if ( prc_ismaster ) then
387  call file_open( offline_parent_basename, & ! (in)
388  fid, & ! (out)
389  aggregate = .false. ) ! (in)
390 
391  call file_get_attribute( fid, "global", "scale_atmos_grid_cartesC_index_imaxg", &
392  imaxg(:), existed=existed )
393  if ( existed ) then
394  call file_get_attribute( fid, "global", "scale_cartesC_prc_num_x", &
395  pnum_x(:), existed=existed )
396  end if
397  if ( existed ) then
398  offline_parent_imax = imaxg(1) / pnum_x(1)
399  else
400  ! for old file
401  call file_get_shape( fid, "CX", dims(:) )
402  offline_parent_imax = dims(1)-ihalo*2
403  end if
404 
405  call file_get_attribute( fid, "global", "scale_atmos_grid_cartesC_index_jmaxg", &
406  jmaxg(:), existed=existed )
407  if ( existed ) then
408  call file_get_attribute( fid, "global", "scale_cartesC_prc_num_y", &
409  pnum_y(:), existed=existed )
410  end if
411  if ( existed ) then
412  offline_parent_jmax = jmaxg(1) / pnum_y(1)
413  else
414  ! for old file
415  call file_get_shape( fid, "CY", dims(:) )
416  offline_parent_jmax = dims(1)-jhalo*2
417  end if
418 
419  call file_get_attribute( fid, "global", "scale_atmos_grid_cartesC_index_kmax", &
420  dims(:), existed=existed )
421  if ( existed ) then
422  offline_parent_kmax = dims(1)
423  else
424  call file_get_shape( fid, "z", dims(:), error=error )
425  if ( error ) then
426  offline_parent_kmax = 0
427  else
428  offline_parent_kmax = dims(1)
429  endif
430  end if
431 
432  call file_get_attribute( fid, "global", "scale_ocean_grid_cartesC_index_kmax", &
433  dims(:), existed=existed )
434  if ( existed ) then
435  offline_parent_okmax = dims(1)
436  else
437  ! for old file
438  call file_get_shape( fid, "oz", dims(:), error=error )
439  if ( error ) then
440  offline_parent_okmax = 0
441  else
442  offline_parent_okmax = dims(1)
443  endif
444  end if
445 
446  call file_get_attribute( fid, "global", "scale_land_grid_cartesC_index_kmax", &
447  dims(:), existed=existed )
448  if ( existed ) then
449  offline_parent_lkmax = dims(1)
450  else
451  ! for old file
452  call file_get_shape( fid, "lz", dims(:), error=error )
453  if ( error ) then
454  offline_parent_lkmax = 0
455  else
456  offline_parent_lkmax = dims(1)
457  endif
458  end if
459 
460  endif
461  call comm_bcast( offline_parent_imax )
462  call comm_bcast( offline_parent_jmax )
463  call comm_bcast( offline_parent_kmax )
464  call comm_bcast( offline_parent_okmax )
465  call comm_bcast( offline_parent_lkmax )
466  endif
467 
468  if ( online_iam_daughter .or. online_iam_parent ) then
469 
470  if ( offline ) then
471  log_error("COMM_CARTESC_NEST_setup",*) 'OFFLINE and ONLINE cannot be use at the same time'
472  call prc_abort
473  endif
474 
475  use_nesting = .true.
476  endif
477 
478  call interp_setup( comm_cartesc_nest_interp_weight_order ) ! [IN]
479 
480  select case ( comm_cartesc_nest_interp_type )
481  case ( 'LINEAR' )
482  itp_nh = 4
483  case ( 'DIST-WEIGHT' )
484  itp_nh = comm_cartesc_nest_interp_level
485  case default
486  log_error("COMM_CARTESC_NEST_setup",*) 'Unsupported type of COMM_CARTESC_NEST_INTERP_TYPE : ', trim(comm_cartesc_nest_interp_type)
487  log_error_cont(*) ' It must be "LINEAR" or "DIST-WEIGHT"'
488  call prc_abort
489  end select
490 
491  debug_domain_num = online_domain_num
492  if( online_specified_maxrq > max_rq ) max_rq = online_specified_maxrq
493 
494  allocate( ireq_p(max_rq) )
495  allocate( ireq_d(max_rq) )
496  allocate( call_order(max_rq) )
497  ireq_p(:) = mpi_request_null
498  ireq_d(:) = mpi_request_null
499 
500  if ( use_nesting ) then
501 
502  if ( offline .OR. online_iam_daughter ) then
503  latlon_local(i_min,i_lon) = minval(atmos_grid_cartesc_real_lonuv(:,:)) / d2r
504  latlon_local(i_max,i_lon) = maxval(atmos_grid_cartesc_real_lonuv(:,:)) / d2r
505  latlon_local(i_min,i_lat) = minval(atmos_grid_cartesc_real_latuv(:,:)) / d2r
506  latlon_local(i_max,i_lat) = maxval(atmos_grid_cartesc_real_latuv(:,:)) / d2r
507  endif
508 
509  if ( offline ) then
510 
511  handling_num = 1
512  parent_prc_num_x(handling_num) = offline_parent_prc_num_x
513  parent_prc_num_y(handling_num) = offline_parent_prc_num_y
514  parent_kmax(handling_num) = offline_parent_kmax
515  parent_imax(handling_num) = offline_parent_imax
516  parent_jmax(handling_num) = offline_parent_jmax
517  parent_okmax(handling_num) = offline_parent_okmax
518  parent_lkmax(handling_num) = offline_parent_lkmax
519 
520  parent_prc_nprocs(handling_num) = parent_prc_num_x(handling_num) * parent_prc_num_y(handling_num)
521  allocate( latlon_catalog(parent_prc_nprocs(handling_num),2,2) )
522 
523  !--- read latlon catalogue
524  fid = io_get_available_fid()
525  open( fid, &
526  file = trim(latlon_catalogue_fname), &
527  form = 'formatted', &
528  status = 'old', &
529  iostat = ierr )
530 
531  if ( ierr /= 0 ) then
532  log_error("COMM_CARTESC_NEST_setup",*) 'cannot open latlon-catalogue file!'
533  call prc_abort
534  endif
535 
536  do i = 1, parent_prc_nprocs(handling_num)
537  read(fid,'(i8,4f32.24)',iostat=ierr) parent_id, &
538  latlon_catalog(i,i_min,i_lon), latlon_catalog(i,i_max,i_lon), & ! LON: MIN, MAX
539  latlon_catalog(i,i_min,i_lat), latlon_catalog(i,i_max,i_lat) ! LAT: MIN, MAX
540  if ( i /= parent_id ) then
541  log_error("COMM_CARTESC_NEST_setup",*) 'internal error: parent mpi id'
542  call prc_abort
543  endif
544  if ( ierr /= 0 ) exit
545  enddo
546  close(fid)
547 
548  call comm_cartesc_nest_domain_relate(handling_num)
549 
550  else ! ONLINE RELATIONSHIP
551 ! if ( present(flag_parent) .AND. present(flag_child) ) then
552 ! LOG_INFO("COMM_CARTESC_NEST_setup",'(1x,A)') &
553 ! '*** Setup Online Nesting Inter-Domain Communicator (IDC)'
554 ! else
555 ! LOG_ERROR("COMM_CARTESC_NEST_setup",*) 'Internal Error:'
556 ! LOG_ERROR_CONT(*) 'The flag_parent and flag_child are needed.'
557 ! LOG_WARN("COMM_CARTESC_NEST_setup",*) ' domain: ', ONLINE_DOMAIN_NUM
558 ! call PRC_abort
559 ! endif
560 
561  if( online_boundary_use_qhyd ) then
562  mp_type = mp_type_in
563  comm_cartesc_nest_bnd_qa = qa_mp
564  elseif ( atmos_hydrometeor_dry ) then
565  mp_type = "DRY"
566  comm_cartesc_nest_bnd_qa = 0
567  else
568  mp_type = "QV"
569  comm_cartesc_nest_bnd_qa = 1
570  endif
571 
572  log_info("COMM_CARTESC_NEST_setup",*) "flag_parent", flag_parent, "flag_child", flag_child
573  log_info("COMM_CARTESC_NEST_setup",*) "ONLINE_IAM_PARENT", online_iam_parent, "ONLINE_IAM_DAUGHTER", online_iam_daughter
574 
575  if( flag_parent ) then ! must do first before daughter processes
576  !-------------------------------------------------
577  if ( .NOT. online_iam_parent ) then
578  log_error("COMM_CARTESC_NEST_setup",*) '[NEST_setup] Parent Flag from launcher is not consistent with namelist!'
579  log_error_cont(*) 'PARENT - domain : ', online_domain_num
580  call prc_abort
581  endif
582 
583  handling_num = 1 !HANDLING_NUM + 1
584  intercomm_id(handling_num) = online_domain_num
585  comm_cartesc_nest_filiation(intercomm_id(handling_num)) = 1
586 
587  intercomm_daughter = inter_child
588  log_info("COMM_CARTESC_NEST_setup",'(1x,A,I2,A)') 'Online Nesting - PARENT [INTERCOMM_ID:', &
589  intercomm_id(handling_num), ' ]'
590  log_info("COMM_CARTESC_NEST_setup",*) 'Online Nesting - INTERCOMM :', intercomm_daughter
591 
592  call comm_cartesc_nest_ping( handling_num )
593 
594  call comm_cartesc_nest_parentsize( handling_num )
595 
596  call comm_cartesc_nest_catalogue( handling_num )
597  call mpi_barrier(intercomm_daughter, ierr)
598 
599  parent_ka(handling_num) = parent_kmax(handling_num) + khalo * 2
600  parent_ia(handling_num) = parent_imax(handling_num) + ihalo * 2
601  parent_ja(handling_num) = parent_jmax(handling_num) + jhalo * 2
602  daughter_ka(handling_num) = daughter_kmax(handling_num) + khalo * 2
603  daughter_ia(handling_num) = daughter_imax(handling_num) + ihalo * 2
604  daughter_ja(handling_num) = daughter_jmax(handling_num) + jhalo * 2
605  tileal_ka(handling_num) = 0
606  tileal_ia(handling_num) = 0
607  tileal_ja(handling_num) = 0
608 
609  log_info("COMM_CARTESC_NEST_setup",'(1x,A)' ) 'Informations of Parent Domain [me]'
610  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_nprocs :', parent_prc_nprocs(handling_num)
611  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_NUM_X :', parent_prc_num_x(handling_num)
612  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_NUM_Y :', parent_prc_num_y(handling_num)
613  log_info_cont('(1x,A,I6)' ) '--- PARENT_KMAX :', parent_kmax(handling_num)
614  log_info_cont('(1x,A,I6)' ) '--- PARENT_IMAX :', parent_imax(handling_num)
615  log_info_cont('(1x,A,I6)' ) '--- PARENT_JMAX :', parent_jmax(handling_num)
616  log_info_cont('(1x,A,F9.3)') '--- PARENT_DTSEC :', parent_dtsec(handling_num)
617  log_info_cont('(1x,A,I6) ') '--- PARENT_NSTEP :', parent_nstep(handling_num)
618  log_info_cont('(1x,A)' ) 'Informations of Daughter Domain'
619  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_PRC_nprocs :', daughter_prc_nprocs(handling_num)
620  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_PRC_NUM_X :', daughter_prc_num_x(handling_num)
621  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_PRC_NUM_Y :', daughter_prc_num_y(handling_num)
622  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_KMAX :', daughter_kmax(handling_num)
623  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_IMAX :', daughter_imax(handling_num)
624  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_JMAX :', daughter_jmax(handling_num)
625  log_info_cont('(1x,A,F9.3)') '--- DAUGHTER_DTSEC :', daughter_dtsec(handling_num)
626  log_info_cont('(1x,A,I6) ') '--- DAUGHTER_NSTEP :', daughter_nstep(handling_num)
627  log_info_cont('(1x,A,I6) ') 'Limit Num. NCOMM req. :', max_rq
628 
629  allocate( org_dens(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
630  allocate( org_momz(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
631  allocate( org_momx(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
632  allocate( org_momy(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
633  allocate( org_u_ll(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
634  allocate( org_v_ll(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
635  allocate( org_rhot(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
636  allocate( org_qtrc(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num),max(comm_cartesc_nest_bnd_qa,1)) )
637 
638  call comm_cartesc_nest_setup_nestdown( handling_num )
639 
640  !---------------------------------- end of parent routines
641  endif
642 
643 
644  if( flag_child ) then
645  !-------------------------------------------------
646  if ( .NOT. online_iam_daughter ) then
647  log_error("COMM_CARTESC_NEST_setup",*) '[NEST_setup] Child Flag from launcher is not consistent with namelist!'
648  log_error_cont(*) 'DAUGHTER - domain : ', online_domain_num
649  call prc_abort
650  endif
651 
652  handling_num = 2 !HANDLING_NUM + 1
653  intercomm_id(handling_num) = online_domain_num - 1
654  comm_cartesc_nest_filiation(intercomm_id(handling_num)) = -1
655 
656  intercomm_parent = inter_parent
657  log_info("COMM_CARTESC_NEST_setup",'(1x,A,I2,A)') 'Online Nesting - DAUGHTER [INTERCOMM_ID:', &
658  intercomm_id(handling_num), ' ]'
659  log_info("COMM_CARTESC_NEST_setup",*) 'Online Nesting - INTERCOMM :', intercomm_parent
660 
661  call comm_cartesc_nest_ping( handling_num )
662 
663  call comm_cartesc_nest_parentsize( handling_num )
664 
665  allocate( latlon_catalog(parent_prc_nprocs(handling_num),2,2) )
666  call comm_cartesc_nest_catalogue( handling_num )
667  call mpi_barrier(intercomm_parent, ierr)
668 
669  call comm_cartesc_nest_domain_relate( handling_num )
670 
671  parent_ka(handling_num) = parent_kmax(handling_num) + khalo * 2
672  parent_ia(handling_num) = parent_imax(handling_num) + ihalo * 2
673  parent_ja(handling_num) = parent_jmax(handling_num) + jhalo * 2
674  daughter_ka(handling_num) = daughter_kmax(handling_num) + khalo * 2
675  daughter_ia(handling_num) = daughter_imax(handling_num) + ihalo * 2
676  daughter_ja(handling_num) = daughter_jmax(handling_num) + jhalo * 2
677  tileal_ka(handling_num) = parent_ka(handling_num)
678  tileal_ia(handling_num) = parent_imax(handling_num) * comm_cartesc_nest_tile_num_x
679  tileal_ja(handling_num) = parent_jmax(handling_num) * comm_cartesc_nest_tile_num_y
680 
681  log_info("COMM_CARTESC_NEST_setup",'(1x,A)' ) 'Informations of Parent Domain'
682  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_nprocs :', parent_prc_nprocs(handling_num)
683  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_NUM_X :', parent_prc_num_x(handling_num)
684  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_NUM_Y :', parent_prc_num_y(handling_num)
685  log_info_cont('(1x,A,I6)' ) '--- PARENT_KMAX :', parent_kmax(handling_num)
686  log_info_cont('(1x,A,I6)' ) '--- PARENT_IMAX :', parent_imax(handling_num)
687  log_info_cont('(1x,A,I6)' ) '--- PARENT_JMAX :', parent_jmax(handling_num)
688  log_info_cont('(1x,A,F9.3)') '--- PARENT_DTSEC :', parent_dtsec(handling_num)
689  log_info_cont('(1x,A,I6)' ) '--- PARENT_NSTEP :', parent_nstep(handling_num)
690  log_info_cont('(1x,A)' ) 'Informations of Daughter Domain [me]'
691  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_PRC_nprocs :', daughter_prc_nprocs(handling_num)
692  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_PRC_NUM_X :', daughter_prc_num_x(handling_num)
693  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_PRC_NUM_Y :', daughter_prc_num_y(handling_num)
694  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_KMAX :', daughter_kmax(handling_num)
695  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_IMAX :', daughter_imax(handling_num)
696  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_JMAX :', daughter_jmax(handling_num)
697  log_info_cont('(1x,A,F9.3)') '--- DAUGHTER_DTSEC :', daughter_dtsec(handling_num)
698  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_NSTEP :', daughter_nstep(handling_num)
699  log_info_cont('(1x,A)' ) 'Informations of Target Tiles'
700  log_info_cont('(1x,A,I6)' ) '--- TILEALL_KA :', tileal_ka(handling_num)
701  log_info_cont('(1x,A,I6)' ) '--- TILEALL_IA :', tileal_ia(handling_num)
702  log_info_cont('(1x,A,I6)' ) '--- TILEALL_JA :', tileal_ja(handling_num)
703  log_info_cont('(1x,A,I6) ') 'Limit Num. NCOMM req. :', max_rq
704 
705  allocate( recvbuf_3d( parent_ka(handling_num), parent_ia(handling_num), parent_ja(handling_num), max_isu ) )
706 
707  allocate( buffer_ref_lon( tileal_ia(handling_num),tileal_ja(handling_num)) )
708  allocate( buffer_ref_lonuy( tileal_ia(handling_num),tileal_ja(handling_num)) )
709  allocate( buffer_ref_lonxv( tileal_ia(handling_num),tileal_ja(handling_num)) )
710  allocate( buffer_ref_lat( tileal_ia(handling_num),tileal_ja(handling_num)) )
711  allocate( buffer_ref_latuy( tileal_ia(handling_num),tileal_ja(handling_num)) )
712  allocate( buffer_ref_latxv( tileal_ia(handling_num),tileal_ja(handling_num)) )
713 
714  allocate( buffer_ref_cz(tileal_ka(handling_num),tileal_ia(handling_num),tileal_ja(handling_num)) )
715  allocate( buffer_ref_fz(tileal_ka(handling_num),tileal_ia(handling_num),tileal_ja(handling_num)) )
716 
717  allocate( buffer_ref_3d(tileal_ka(handling_num),tileal_ia(handling_num),tileal_ja(handling_num)) )
718 
719  allocate( igrd( daughter_ia(handling_num),daughter_ja(handling_num),itp_nh,itp_ng) )
720  allocate( jgrd( daughter_ia(handling_num),daughter_ja(handling_num),itp_nh,itp_ng) )
721  allocate( hfact( daughter_ia(handling_num),daughter_ja(handling_num),itp_nh,itp_ng) )
722  allocate( kgrd(daughter_ka(handling_num),2,daughter_ia(handling_num),daughter_ja(handling_num),itp_nh,itp_ng) )
723  allocate( vfact(daughter_ka(handling_num), daughter_ia(handling_num),daughter_ja(handling_num),itp_nh,itp_ng) )
724 
725  call comm_cartesc_nest_setup_nestdown( handling_num )
726 
727 
728  select case ( comm_cartesc_nest_interp_type )
729  case ( 'LINEAR' )
730 
731  allocate( x_ref(tileal_ia(handling_num),tileal_ja(handling_num)) )
732  allocate( y_ref(tileal_ia(handling_num),tileal_ja(handling_num)) )
733 
734  ! for scalar points
735  call mapprojection_lonlat2xy( tileal_ia(handling_num), 1, tileal_ia(handling_num), &
736  tileal_ja(handling_num), 1, tileal_ja(handling_num), &
737  buffer_ref_lon(:,:), & ! [IN]
738  buffer_ref_lat(:,:), & ! [IN]
739  x_ref(:,:), y_ref(:,:) ) ! [OUT]
740  call interp_factor3d( tileal_ka(handling_num), & ! [IN]
741  khalo+1, & ! [IN]
742  tileal_ka(handling_num)-khalo, & ! [IN]
743  tileal_ia(handling_num), & ! [IN]
744  tileal_ja(handling_num), & ! [IN]
745  daughter_ka(handling_num), & ! [IN]
746  datr_ks(handling_num), & ! [IN]
747  datr_ke(handling_num), & ! [IN]
748  daughter_ia(handling_num), & ! [IN]
749  daughter_ja(handling_num), & ! [IN]
750  x_ref(:,:), y_ref(:,:), & ! [IN]
751  buffer_ref_cz(:,:,:), & ! [IN]
752  atmos_grid_cartesc_cx(:), & ! [IN]
753  atmos_grid_cartesc_cy(:), & ! [IN]
754  atmos_grid_cartesc_real_cz(:,:,:), & ! [IN]
755  igrd( :,:,:,i_sclr), & ! [OUT]
756  jgrd( :,:,:,i_sclr), & ! [OUT]
757  hfact( :,:,:,i_sclr), & ! [OUT]
758  kgrd(:,:,:,:,:,i_sclr), & ! [OUT]
759  vfact(:, :,:,:,i_sclr) ) ! [OUT]
760 
761  ! for z staggered points
762  call interp_factor3d( tileal_ka(handling_num)+1, & ! [IN]
763  khalo+1, & ! [IN]
764  tileal_ka(handling_num)+1-khalo, & ! [IN]
765  tileal_ia(handling_num), & ! [IN]
766  tileal_ja(handling_num), & ! [IN]
767  daughter_ka(handling_num), & ! [IN]
768  datr_ks(handling_num), & ! [IN]
769  datr_ke(handling_num), & ! [IN]
770  daughter_ia(handling_num), & ! [IN]
771  daughter_ja(handling_num), & ! [IN]
772  x_ref(:,:), y_ref(:,:), & ! [IN]
773  buffer_ref_fz(:,:,:), & ! [IN]
774  atmos_grid_cartesc_cx(:), & ! [IN]
775  atmos_grid_cartesc_cy(:), & ! [IN]
776  atmos_grid_cartesc_real_fz(1:ka,:,:), & ! [IN]
777  igrd( :,:,:,i_zstg), & ! [OUT]
778  jgrd( :,:,:,i_zstg), & ! [OUT]
779  hfact( :,:,:,i_zstg), & ! [OUT]
780  kgrd(:,:,:,:,:,i_zstg), & ! [OUT]
781  vfact(:, :,:,:,i_zstg) ) ! [OUT]
782 
783  ! for x staggered points
784  call mapprojection_lonlat2xy( tileal_ia(handling_num), 1, tileal_ia(handling_num), &
785  tileal_ja(handling_num), 1, tileal_ja(handling_num), &
786  buffer_ref_lonuy(:,:), & ! [IN]
787  buffer_ref_latuy(:,:), & ! [IN]
788  x_ref(:,:), y_ref(:,:) ) ! [OUT]
789  call interp_factor3d( tileal_ka(handling_num), & ! [IN]
790  khalo+1, & ! [IN]
791  tileal_ka(handling_num)-khalo, & ! [IN]
792  tileal_ia(handling_num), & ! [IN]
793  tileal_ja(handling_num), & ! [IN]
794  daughter_ka(handling_num), & ! [IN]
795  datr_ks(handling_num), & ! [IN]
796  datr_ke(handling_num), & ! [IN]
797  daughter_ia(handling_num), & ! [IN]
798  daughter_ja(handling_num), & ! [IN]
799  x_ref(:,:), y_ref(:,:), & ! [IN]
800  buffer_ref_cz(:,:,:), & ! [IN]
801  atmos_grid_cartesc_fx(1:ia), & ! [IN]
802  atmos_grid_cartesc_cy(:), & ! [IN]
803  atmos_grid_cartesc_real_cz(:,:,:), & ! [IN]
804  igrd( :,:,:,i_xstg), & ! [OUT]
805  jgrd( :,:,:,i_xstg), & ! [OUT]
806  hfact( :,:,:,i_xstg), & ! [OUT]
807  kgrd(:,:,:,:,:,i_xstg), & ! [OUT]
808  vfact(:, :,:,:,i_xstg) ) ! [OUT]
809 
810  ! for y staggered points
811  call mapprojection_lonlat2xy( tileal_ia(handling_num), 1, tileal_ia(handling_num), &
812  tileal_ja(handling_num), 1, tileal_ja(handling_num), &
813  buffer_ref_lonxv(:,:), & ! [IN]
814  buffer_ref_latxv(:,:), & ! [IN]
815  x_ref(:,:), y_ref(:,:) ) ! [OUT]
816  call interp_factor3d( tileal_ka(handling_num), & ! [IN]
817  khalo+1, & ! [IN]
818  tileal_ka(handling_num)-khalo, & ! [IN]
819  tileal_ia(handling_num), & ! [IN]
820  tileal_ja(handling_num), & ! [IN]
821  daughter_ka(handling_num), & ! [IN]
822  datr_ks(handling_num), & ! [IN]
823  datr_ke(handling_num), & ! [IN]
824  daughter_ia(handling_num), & ! [IN]
825  daughter_ja(handling_num), & ! [IN]
826  x_ref(:,:), y_ref(:,:), & ! [IN]
827  buffer_ref_cz(:,:,:), & ! [IN]
828  atmos_grid_cartesc_cx(:), & ! [IN]
829  atmos_grid_cartesc_fy(1:ja), & ! [IN]
830  atmos_grid_cartesc_real_cz(:,:,:), & ! [IN]
831  igrd( :,:,:,i_ystg), & ! [OUT]
832  jgrd( :,:,:,i_ystg), & ! [OUT]
833  hfact( :,:,:,i_ystg), & ! [OUT]
834  kgrd(:,:,:,:,:,i_ystg), & ! [OUT]
835  vfact(:, :,:,:,i_ystg) ) ! [OUT]
836 
837  deallocate( x_ref, y_ref )
838 
839  case ( 'DIST-WEIGHT' )
840 
841  ! for scalar points
842  call interp_factor3d( itp_nh, & ! [IN]
843  tileal_ka(handling_num), & ! [IN]
844  khalo+1, & ! [IN]
845  tileal_ka(handling_num)-khalo, & ! [IN]
846  tileal_ia(handling_num), & ! [IN]
847  tileal_ja(handling_num), & ! [IN]
848  daughter_ka(handling_num), & ! [IN]
849  datr_ks(handling_num), & ! [IN]
850  datr_ke(handling_num), & ! [IN]
851  daughter_ia(handling_num), & ! [IN]
852  daughter_ja(handling_num), & ! [IN]
853  buffer_ref_lon(:,:), & ! [IN]
854  buffer_ref_lat(:,:), & ! [IN]
855  buffer_ref_cz(:,:,:), & ! [IN]
856  atmos_grid_cartesc_real_lon(:,:), & ! [IN]
857  atmos_grid_cartesc_real_lat(:,:), & ! [IN]
858  atmos_grid_cartesc_real_cz(:,:,:), & ! [IN]
859  igrd( :,:,:,i_sclr), & ! [OUT]
860  jgrd( :,:,:,i_sclr), & ! [OUT]
861  hfact( :,:,:,i_sclr), & ! [OUT]
862  kgrd(:,:,:,:,:,i_sclr), & ! [OUT]
863  vfact(:, :,:,:,i_sclr) ) ! [OUT]
864 
865  ! for z staggered points
866  call interp_factor3d( itp_nh, & ! [IN]
867  tileal_ka(handling_num), & ! [IN]
868  khalo, & ! [IN]
869  tileal_ka(handling_num)-khalo, & ! [IN]
870  tileal_ia(handling_num), & ! [IN]
871  tileal_ja(handling_num), & ! [IN]
872  daughter_ka(handling_num), & ! [IN]
873  datr_ks(handling_num), & ! [IN]
874  datr_ke(handling_num), & ! [IN]
875  daughter_ia(handling_num), & ! [IN]
876  daughter_ja(handling_num), & ! [IN]
877  buffer_ref_lon(:,:), & ! [IN]
878  buffer_ref_lat(:,:), & ! [IN]
879  buffer_ref_fz(:,:,:), & ! [IN]
880  atmos_grid_cartesc_real_lon(:,:), & ! [IN]
881  atmos_grid_cartesc_real_lat(:,:), & ! [IN]
882  atmos_grid_cartesc_real_fz(1:ka,:,:), & ! [IN]
883  igrd( :,:,:,i_zstg), & ! [OUT]
884  jgrd( :,:,:,i_zstg), & ! [OUT]
885  hfact( :,:,:,i_zstg), & ! [OUT]
886  kgrd(:,:,:,:,:,i_zstg), & ! [OUT]
887  vfact(:, :,:,:,i_zstg) ) ! [OUT]
888 
889  ! for x staggered points
890  call interp_factor3d( itp_nh, & ! [IN]
891  tileal_ka(handling_num), & ! [IN]
892  khalo+1, & ! [IN]
893  tileal_ka(handling_num)-khalo, & ! [IN]
894  tileal_ia(handling_num), & ! [IN]
895  tileal_ja(handling_num), & ! [IN]
896  daughter_ka(handling_num), & ! [IN]
897  datr_ks(handling_num), & ! [IN]
898  datr_ke(handling_num), & ! [IN]
899  daughter_ia(handling_num), & ! [IN]
900  daughter_ja(handling_num), & ! [IN]
901  buffer_ref_lonuy(:,:), & ! [IN]
902  buffer_ref_latuy(:,:), & ! [IN]
903  buffer_ref_cz(:,:,:), & ! [IN]
904  atmos_grid_cartesc_real_lonuy(1:ia,1:ja), & ! [IN]
905  atmos_grid_cartesc_real_latuy(1:ia,1:ja), & ! [IN]
906  atmos_grid_cartesc_real_cz(:,:,:), & ! [IN]
907  igrd( :,:,:,i_xstg), & ! [OUT]
908  jgrd( :,:,:,i_xstg), & ! [OUT]
909  hfact( :,:,:,i_xstg), & ! [OUT]
910  kgrd(:,:,:,:,:,i_xstg), & ! [OUT]
911  vfact(:, :,:,:,i_xstg) ) ! [OUT]
912 
913  ! for y staggered points
914  call interp_factor3d( itp_nh, & ! [IN]
915  tileal_ka(handling_num), & ! [IN]
916  khalo+1, & ! [IN]
917  tileal_ka(handling_num)-khalo, & ! [IN]
918  tileal_ia(handling_num), & ! [IN]
919  tileal_ja(handling_num), & ! [IN]
920  daughter_ka(handling_num), & ! [IN]
921  datr_ks(handling_num), & ! [IN]
922  datr_ke(handling_num), & ! [IN]
923  daughter_ia(handling_num), & ! [IN]
924  daughter_ja(handling_num), & ! [IN]
925  buffer_ref_lonxv(:,:), & ! [IN]
926  buffer_ref_latxv(:,:), & ! [IN]
927  buffer_ref_cz(:,:,:), & ! [IN]
928  atmos_grid_cartesc_real_lonxv(1:ia,1:ja), & ! [IN]
929  atmos_grid_cartesc_real_latxv(1:ia,1:ja), & ! [IN]
930  atmos_grid_cartesc_real_cz(:,:,:), & ! [IN]
931  igrd( :,:,:,i_ystg), & ! [OUT]
932  jgrd( :,:,:,i_ystg), & ! [OUT]
933  hfact( :,:,:,i_ystg), & ! [OUT]
934  kgrd(:,:,:,:,:,i_ystg), & ! [OUT]
935  vfact(:, :,:,:,i_ystg) ) ! [OUT]
936 
937  end select
938 
939 
940  else
941  online_use_velz = .false.
942  endif
943 
944  !LOG_INFO("COMM_CARTESC_NEST_setup",'(1x,A,I2)') 'Number of Related Domains :', HANDLING_NUM
945  !if ( HANDLING_NUM > 2 ) then
946  ! f( IO_L ) LOG_ERROR("COMM_CARTESC_NEST_setup",*) 'Too much handing domains (up to 2)'
947  ! call PRC_abort
948  !endif
949 
950  endif !--- OFFLINE or NOT
951 
952  endif !--- USE_NESTING
953 
954  return

References scale_atmos_grid_cartesc::atmos_grid_cartesc_cx, scale_atmos_grid_cartesc::atmos_grid_cartesc_cy, scale_atmos_grid_cartesc::atmos_grid_cartesc_fx, scale_atmos_grid_cartesc::atmos_grid_cartesc_fy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lat, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latuv, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latuy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latxv, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lon, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonuv, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonuy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonxv, scale_atmos_hydrometeor::atmos_hydrometeor_dry, comm_cartesc_nest_bnd_qa, comm_cartesc_nest_domain_relate(), comm_cartesc_nest_filiation, comm_cartesc_nest_interp_level, comm_cartesc_nest_interp_weight_order, comm_cartesc_nest_tile_num_x, comm_cartesc_nest_tile_num_y, 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, scale_file::file_open(), handling_num, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ihalo, scale_atmos_grid_cartesc_index::imaxg, intercomm_daughter, intercomm_parent, scale_interp::interp_setup(), scale_io::io_fid_conf, scale_io::io_get_available_fid(), scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::jhalo, scale_atmos_grid_cartesc_index::jmaxg, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::khalo, offline, online_boundary_use_qhyd, online_domain_num, online_iam_daughter, online_iam_parent, online_no_rotate, online_use_velz, parent_dtsec, parent_ia, parent_imax, parent_ja, parent_jmax, parent_ka, parent_kmax, parent_lkmax, parent_nstep, parent_okmax, scale_prc::prc_abort(), scale_prc::prc_global_domainid, scale_prc::prc_ismaster, tileal_ia, tileal_ja, tileal_ka, and use_nesting.

Referenced by mod_rm_driver::rm_driver(), and mod_rm_prep::rm_prep().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ comm_cartesc_nest_domain_relate()

subroutine, public scale_comm_cartesc_nest::comm_cartesc_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 961 of file scale_comm_cartesC_nest.F90.

961  use scale_prc, only: &
962  prc_myrank, &
963  prc_abort
964  implicit none
965 
966  integer, intent(in) :: HANDLE
967 
968  integer :: x_min, x_max
969  integer :: y_min, y_max
970  logical :: hit(2,2)
971  integer :: p, i, j
972  !---------------------------------------------------------------------------
973 
974  if( .NOT. use_nesting ) return
975 
976  x_min = parent_prc_num_x(handle)
977  x_max = -1
978  y_min = parent_prc_num_y(handle)
979  y_max = -1
980  hit(:,:) = .false.
981  do p = 1, parent_prc_nprocs(handle)
982  if ( ( ( latlon_local(i_min,i_lon) >= latlon_catalog(p,i_min,i_lon) &
983  .AND. latlon_local(i_min,i_lon) <= latlon_catalog(p,i_max,i_lon) ) .OR. &
984  ( latlon_local(i_max,i_lon) >= latlon_catalog(p,i_min,i_lon) &
985  .AND. latlon_local(i_max,i_lon) <= latlon_catalog(p,i_max,i_lon) ) .OR. &
986  ( latlon_catalog(p,i_min,i_lon) >= latlon_local(i_min,i_lon) &
987  .AND. latlon_catalog(p,i_min,i_lon) <= latlon_local(i_max,i_lon) ) .OR. &
988  ( latlon_catalog(p,i_max,i_lon) >= latlon_local(i_min,i_lon) &
989  .AND. latlon_catalog(p,i_max,i_lon) <= latlon_local(i_max,i_lon) ) ) .AND. &
990  ( ( latlon_local(i_min,i_lat) >= latlon_catalog(p,i_min,i_lat) &
991  .AND. latlon_local(i_min,i_lat) <= latlon_catalog(p,i_max,i_lat) ) .OR. &
992  ( latlon_local(i_max,i_lat) >= latlon_catalog(p,i_min,i_lat) &
993  .AND. latlon_local(i_max,i_lat) <= latlon_catalog(p,i_max,i_lat) ) .OR. &
994  ( latlon_catalog(p,i_min,i_lat) >= latlon_local(i_min,i_lat) &
995  .AND. latlon_catalog(p,i_min,i_lat) <= latlon_local(i_max,i_lat) ) .OR. &
996  ( latlon_catalog(p,i_max,i_lat) >= latlon_local(i_min,i_lat) &
997  .AND. latlon_catalog(p,i_max,i_lat) <= latlon_local(i_max,i_lat) ) ) ) then
998  if ( latlon_catalog(p,i_min,i_lon) <= latlon_local(i_min,i_lon) ) hit(i_min,i_lon) = .true.
999  if ( latlon_catalog(p,i_max,i_lon) >= latlon_local(i_max,i_lon) ) hit(i_max,i_lon) = .true.
1000  if ( latlon_catalog(p,i_min,i_lat) <= latlon_local(i_min,i_lat) ) hit(i_min,i_lat) = .true.
1001  if ( latlon_catalog(p,i_max,i_lat) >= latlon_local(i_max,i_lat) ) hit(i_max,i_lat) = .true.
1002  i = mod(p-1, parent_prc_num_x(handle))
1003  j = (p-1) / parent_prc_num_x(handle)
1004  if ( i < x_min ) x_min = i
1005  if ( i > x_max ) x_max = i
1006  if ( j < y_min ) y_min = j
1007  if ( j > y_max ) y_max = j
1008  end if
1009  end do
1010 
1011  if ( .not. ( hit(i_min,i_lon) .and. hit(i_max,i_lon) .and. hit(i_min,i_lat) .and. hit(i_max,i_lat) ) ) then
1012  log_error("COMM_CARTESC_NEST_domain_relate",*) 'region of daughter domain is larger than that of parent'
1013  log_error_cont(*) ' at rank:', prc_myrank, ' of domain:', online_domain_num
1014  log_error_cont(*) 'LON MIN: ',hit(i_min,i_lon), ', LON MAX: ',hit(i_max,i_lon), ', LAT MIN: ',hit(i_min,i_lat), ', LAT MAX: ',hit(i_max,i_lat)
1015  log_error_cont('(A,F12.6,1x,F12.6)') 'daughter local (me) MIN-MAX: LON=', &
1016  latlon_local(i_min,i_lon), latlon_local(i_max,i_lon)
1017  do p = 1, parent_prc_nprocs(handle)
1018  log_error_cont('(A,I5,A,F12.6,1x,F12.6)') ' parent (', p,') MIN-MAX: LON=', &
1019  latlon_catalog(p,i_min,i_lon) ,latlon_catalog(p,i_max,i_lon)
1020  enddo
1021  log_error_cont('(A,F12.6,1x,F12.6)') 'daughter local (me): MIN-MAX LAT=', &
1022  latlon_local(i_min,i_lat), latlon_local(i_max,i_lat)
1023  do p = 1, parent_prc_nprocs(handle)
1024  log_error_cont('(A,I5,A,F12.6,1x,F12.6)') ' parent (', p,') MIN-MAX: LAT=', &
1025  latlon_catalog(p,i_min,i_lat) ,latlon_catalog(p,i_max,i_lat)
1026  enddo
1027  call prc_abort
1028  end if
1029 
1030 
1031 
1032  comm_cartesc_nest_tile_num_x = x_max - x_min + 1
1033  comm_cartesc_nest_tile_num_y = y_max - y_min + 1
1034 
1035  allocate( comm_cartesc_nest_tile_id( comm_cartesc_nest_tile_num_x * comm_cartesc_nest_tile_num_y ) )
1036 
1037  log_info("COMM_CARTESC_NEST_domain_relate",'(1x,A)') 'NEST: target process tile in parent domain'
1038  p = 1
1039  do j = 1, comm_cartesc_nest_tile_num_y
1040  do i = 1, comm_cartesc_nest_tile_num_x
1041  comm_cartesc_nest_tile_id(p) = x_min + i - 1 + (y_min + j - 1) * parent_prc_num_x(handle)
1042  log_info_cont('(1x,A,I4,A,I6)') '(', p, ') target mpi-process:', comm_cartesc_nest_tile_id(p)
1043  p = p + 1
1044  enddo
1045  enddo
1046 
1047  return

References comm_cartesc_nest_tile_id, comm_cartesc_nest_tile_num_x, comm_cartesc_nest_tile_num_y, online_domain_num, scale_prc::prc_abort(), scale_prc::prc_myrank, and use_nesting.

Referenced by comm_cartesc_nest_setup().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ comm_cartesc_nest_domain_shape()

subroutine, public scale_comm_cartesc_nest::comm_cartesc_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 1061 of file scale_comm_cartesC_nest.F90.

1061  implicit none
1062 
1063  integer, intent(out) :: tilei, tilej
1064  integer, intent(out) :: cxs, cxe, cys, cye
1065  integer, intent(out) :: pxs, pxe, pys, pye
1066  integer, intent(in) :: iloc ! rank number; start from 1
1067 
1068  integer :: hdl = 1 ! handler number
1069  integer :: rank
1070  integer :: xloc, yloc
1071  integer :: xlocg, ylocg ! location over whole domain
1072  !---------------------------------------------------------------------------
1073 
1074  if( .NOT. use_nesting ) return
1075 
1076  rank = comm_cartesc_nest_tile_id(iloc)
1077  xloc = mod( iloc-1, comm_cartesc_nest_tile_num_x ) + 1
1078  yloc = int( real(iloc-1) / real(comm_cartesc_nest_tile_num_x) ) + 1
1079  xlocg = mod( rank, offline_parent_prc_num_x ) + 1
1080  ylocg = int( real(rank) / real(offline_parent_prc_num_x) ) + 1
1081  tilei = parent_imax(hdl)
1082  tilej = parent_jmax(hdl)
1083 
1084  cxs = tilei * (xloc-1) + 1
1085  cxe = tilei * xloc
1086  cys = tilej * (yloc-1) + 1
1087  cye = tilej * yloc
1088  pxs = 1
1089  pxe = tilei
1090  pys = 1
1091  pye = tilej
1092 
1093  if ( xlocg == 1 ) then ! BND_W
1094  tilei = tilei + 2
1095  pxs = pxs + 2
1096  pxe = pxe + 2
1097  endif
1098  if ( xlocg == offline_parent_prc_num_x ) then ! BND_E
1099  tilei = tilei + 2
1100  endif
1101  if ( ylocg == 1 ) then ! BND_S
1102  tilej = tilej + 2
1103  pys = pys + 2
1104  pye = pye + 2
1105  endif
1106  if ( ylocg == offline_parent_prc_num_y ) then ! BND_N
1107  tilej = tilej + 2
1108  endif
1109 
1110  return

References scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_domain_catalogue, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lat, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latuy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latxv, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lon, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonuy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonxv, comm_cartesc_nest_bnd_qa, comm_cartesc_nest_filiation, comm_cartesc_nest_tile_id, comm_cartesc_nest_tile_num_x, comm_cartesc_nest_waitall(), scale_comm_cartesc::comm_datatype, scale_comm_cartesc::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_io::h_short, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::imax, intercomm_daughter, intercomm_parent, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::jmax, scale_atmos_grid_cartesc_index::js, scale_tracer::k, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::kmax, scale_atmos_grid_cartesc_index::ks, scale_atmos_hydrometeor::n_hyd, online_boundary_diagqhyd, 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_prc::prc_abort(), scale_prc::prc_ismaster, scale_prc::prc_myrank, scale_prc::prc_nprocs, scale_prc_cartesc::prc_num_x, scale_prc_cartesc::prc_num_y, prnt_ie, prnt_is, prnt_je, prnt_js, prnt_ke, prnt_ks, scale_time::time_dtsec, scale_time::time_nstep, and use_nesting.

Referenced by mod_copytopo::copytopo_get_data_scale().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ comm_cartesc_nest_nestdown()

subroutine, public scale_comm_cartesc_nest::comm_cartesc_nest_nestdown ( integer, intent(in)  HANDLE,
integer, intent(in)  BND_QA,
real(rp), dimension(parent_ka (handle),parent_ia (handle),parent_ja (handle)), intent(in)  DENS_send,
real(rp), dimension(parent_ka (handle),parent_ia (handle),parent_ja (handle)), intent(in)  MOMZ_send,
real(rp), dimension(parent_ka (handle),parent_ia (handle),parent_ja (handle)), intent(in)  MOMX_send,
real(rp), dimension(parent_ka (handle),parent_ia (handle),parent_ja (handle)), intent(in)  MOMY_send,
real(rp), dimension(parent_ka (handle),parent_ia (handle),parent_ja (handle)), intent(in)  RHOT_send,
real(rp), dimension(parent_ka (handle),parent_ia (handle),parent_ja (handle),bnd_qa), intent(in)  QTRC_send,
real(rp), dimension(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle)), intent(out)  DENS_recv,
real(rp), dimension(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle)), intent(out)  VELZ_recv,
real(rp), dimension(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle)), intent(out)  VELX_recv,
real(rp), dimension(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle)), intent(out)  VELY_recv,
real(rp), dimension(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle)), intent(out)  POTT_recv,
real(rp), dimension(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle),bnd_qa), intent(out)  QTRC_recv 
)

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 1868 of file scale_comm_cartesC_nest.F90.

1868  use scale_prc, only: &
1869  prc_abort
1870  use scale_comm_cartesc, only: &
1871  comm_vars8, &
1872  comm_wait
1873  use scale_atmos_grid_cartesc_metric, only: &
1875  implicit none
1876 
1877  integer, intent(in) :: HANDLE
1878  integer, intent(in) :: BND_QA
1879  real(RP), intent(in) :: DENS_send(PARENT_KA (HANDLE),PARENT_IA (HANDLE),PARENT_JA (HANDLE))
1880  real(RP), intent(in) :: MOMZ_send(PARENT_KA (HANDLE),PARENT_IA (HANDLE),PARENT_JA (HANDLE))
1881  real(RP), intent(in) :: MOMX_send(PARENT_KA (HANDLE),PARENT_IA (HANDLE),PARENT_JA (HANDLE))
1882  real(RP), intent(in) :: MOMY_send(PARENT_KA (HANDLE),PARENT_IA (HANDLE),PARENT_JA (HANDLE))
1883  real(RP), intent(in) :: RHOT_send(PARENT_KA (HANDLE),PARENT_IA (HANDLE),PARENT_JA (HANDLE))
1884  real(RP), intent(in) :: QTRC_send(PARENT_KA (HANDLE),PARENT_IA (HANDLE),PARENT_JA (HANDLE),BND_QA)
1885  real(RP), intent(out) :: DENS_recv(DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE))
1886  real(RP), intent(out) :: VELZ_recv(DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE))
1887  real(RP), intent(out) :: VELX_recv(DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE))
1888  real(RP), intent(out) :: VELY_recv(DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE))
1889  real(RP), intent(out) :: POTT_recv(DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE))
1890  real(RP), intent(out) :: QTRC_recv(DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE),BND_QA)
1891 
1892  real(RP) :: WORK1_send(PARENT_KA (HANDLE),PARENT_IA (HANDLE),PARENT_JA (HANDLE))
1893  real(RP) :: WORK2_send(PARENT_KA (HANDLE),PARENT_IA (HANDLE),PARENT_JA (HANDLE))
1894  real(RP) :: WORK1_recv(DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE))
1895  real(RP) :: WORK2_recv(DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE))
1896  real(RP) :: U_ll_recv (DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE))
1897  real(RP) :: V_ll_recv (DAUGHTER_KA(HANDLE),DAUGHTER_IA(HANDLE),DAUGHTER_JA(HANDLE))
1898  real(RP) :: u_on_map, v_on_map
1899 
1900  real(RP) :: dummy(1,1,1)
1901  integer :: tagbase, tagcomm
1902  integer :: isu_tag
1903 
1904  integer :: ierr
1905  integer :: i, j, k, iq
1906  !---------------------------------------------------------------------------
1907 
1908  if( .NOT. use_nesting ) return
1909 
1910  if ( bnd_qa > i_bndqa ) then
1911  log_error("COMM_CARTESC_NEST_nestdown",*) 'internal error: BND_QA is larger than I_BNDQA'
1912  call prc_abort
1913  endif
1914 
1915  tagcomm = intercomm_id(handle) * order_tag_comm
1916 
1917  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
1918 
1919  !##### parent [send issue] #####
1920 
1921  call prof_rapstart('NEST_total_P', 2)
1922  call prof_rapstart('NEST_pack_P', 2)
1923 
1924  nsend = nsend + 1
1925  log_info("COMM_CARTESC_NEST_nestdown",'(1X,A,I5,A)') "CONeP[P] send( ", nsend, " )"
1926 
1927  ! to keep values at that time by finish of sending process
1928 !OCL XFILL
1929  org_dens(:,:,:) = dens_send(:,:,:)
1930 !OCL XFILL
1931  org_momz(:,:,:) = momz_send(:,:,:)
1932 !OCL XFILL
1933  org_momx(:,:,:) = momx_send(:,:,:)
1934 !OCL XFILL
1935  org_momy(:,:,:) = momy_send(:,:,:)
1936 !OCL XFILL
1937  org_rhot(:,:,:) = rhot_send(:,:,:)
1938  do iq = 1, bnd_qa
1939 !OCL XFILL
1940  org_qtrc(:,:,:,iq) = qtrc_send(:,:,:,iq)
1941  enddo
1942 
1943  !*** request control
1944  !--- do not change the calling order below;
1945  !--- it should be consistent with the order in "COMM_CARTESC_NEST_recvwait_issue"
1946  rq_ctl_p = 0
1947 
1948  if ( .NOT. online_daughter_no_rotate ) then
1949  ! from staggered point to scalar point
1950  do j = 1, parent_ja(handle)
1951  do i = 2, parent_ia(handle)
1952  do k = 1, parent_ka(handle)
1953  work1_send(k,i,j) = ( org_momx(k,i-1,j) + org_momx(k,i,j) ) * 0.5_rp
1954  enddo
1955  enddo
1956  enddo
1957 
1958  do j = 1, parent_ja(handle)
1959  do k = 1, parent_ka(handle)
1960  work1_send(k,1,j) = org_momx(k,1,j)
1961  enddo
1962  enddo
1963 
1964  call comm_vars8( work1_send(:,:,:), 1 )
1965 
1966  do j = 2, parent_ja(handle)
1967  do i = 1, parent_ia(handle)
1968  do k = 1, parent_ka(handle)
1969  work2_send(k,i,j) = ( org_momy(k,i,j-1) + org_momy(k,i,j) ) * 0.5_rp
1970  enddo
1971  enddo
1972  enddo
1973 
1974  do i = 1, parent_ia(handle)
1975  do k = 1, parent_ka(handle)
1976  work2_send(k,i,1) = org_momy(k,i,1)
1977  enddo
1978  enddo
1979 
1980  call comm_vars8( work2_send(:,:,:), 2 )
1981 
1982  call comm_wait ( work1_send(:,:,:), 1, .false. )
1983  call comm_wait ( work2_send(:,:,:), 2, .false. )
1984 
1985  ! rotation from map-projected field to latlon field
1986  do j = 1, parent_ja(handle)
1987  do i = 1, parent_ia(handle)
1988  do k = 1, parent_ka(handle)
1989  u_on_map = work1_send(k,i,j) / org_dens(k,i,j)
1990  v_on_map = work2_send(k,i,j) / org_dens(k,i,j)
1991 
1992  org_u_ll(k,i,j) = u_on_map * rotc(i,j,1) - v_on_map * rotc(i,j,2)
1993  org_v_ll(k,i,j) = u_on_map * rotc(i,j,2) + v_on_map * rotc(i,j,1)
1994  enddo
1995  enddo
1996  enddo
1997  endif
1998 
1999  tagbase = tagcomm + tag_dens*order_tag_var
2000  call comm_cartesc_nest_intercomm_nestdown( org_dens(:,:,:), & ! [IN]
2001  dummy(:,:,:), & ! [OUT]
2002  tagbase, i_sclr, handle, & ! [IN]
2003  isu_tag, & ! [INOUT]
2004  flag_dens = .true. ) ! [IN]
2005 
2006  tagbase = tagcomm + tag_momz*order_tag_var
2007  if ( online_daughter_use_velz ) then
2008  call comm_cartesc_nest_intercomm_nestdown( org_momz(:,:,:), & ! [IN]
2009  dummy(:,:,:), & ! [OUT]
2010  tagbase, i_zstg, handle, & ! [IN]
2011  isu_tag ) ! [INOUT]
2012  endif
2013 
2014  tagbase = tagcomm + tag_momx*order_tag_var
2015  if ( online_daughter_no_rotate ) then
2016  call comm_cartesc_nest_intercomm_nestdown( org_momx(:,:,:), & ! [IN]
2017  dummy(:,:,:), & ! [OUT]
2018  tagbase, i_xstg, handle, & ! [IN]
2019  isu_tag ) ! [INOUT]
2020  else
2021  call comm_cartesc_nest_intercomm_nestdown( org_u_ll(:,:,:), & ! [IN]
2022  dummy(:,:,:), & ! [OUT]
2023  tagbase, i_sclr, handle, & ! [IN]
2024  isu_tag ) ! [INOUT]
2025  endif
2026 
2027  tagbase = tagcomm + tag_momy*order_tag_var
2028  if ( online_daughter_no_rotate ) then
2029  call comm_cartesc_nest_intercomm_nestdown( org_momy(:,:,:), & ! [IN]
2030  dummy(:,:,:), & ! [OUT]
2031  tagbase, i_ystg, handle, & ! [IN]
2032  isu_tag ) ! [INOUT]
2033  else
2034  call comm_cartesc_nest_intercomm_nestdown( org_v_ll(:,:,:), & ! [IN]
2035  dummy(:,:,:), & ! [OUT]
2036  tagbase, i_sclr, handle, & ! [IN]
2037  isu_tag ) ! [INOUT]
2038  endif
2039 
2040  tagbase = tagcomm + tag_rhot*order_tag_var
2041  call comm_cartesc_nest_intercomm_nestdown( org_rhot(:,:,:), & ! [IN]
2042  dummy(:,:,:), & ! [OUT]
2043  tagbase, i_sclr, handle, & ! [IN]
2044  isu_tag ) ! [INOUT]
2045 
2046  do iq = 1, bnd_qa
2047  tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
2048  call comm_cartesc_nest_intercomm_nestdown( org_qtrc(:,:,:,iq), & ! [IN]
2049  dummy(:,:,:), & ! [OUT]
2050  tagbase, i_sclr, handle, & ! [IN]
2051  isu_tag ) ! [INOUT]
2052  enddo
2053 
2054  rq_tot_p = rq_ctl_p
2055 
2056  call prof_rapend ('NEST_pack_P', 2)
2057  call prof_rapend ('NEST_total_P', 2)
2058 
2059  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2060 
2061  !##### child [wait issue] #####
2062 
2063  call prof_rapstart('NEST_total_C', 2)
2064  call prof_rapstart('NEST_wait_C', 2)
2065 
2066  nwait_d = nwait_d + 1
2067  !LOG_INFO("COMM_CARTESC_NEST_nestdown",'(1X,A,I5,A)') "NestIDC [C]: que wait ( ", nwait_d, " )"
2068 
2069  !*** reset issue tag and request control
2070  !--- do not change the calling order below;
2071  !--- it should be consistent with the order in "COMM_CARTESC_NEST_recvwait_issue"
2072  isu_tag = 0
2073 
2074  call comm_cartesc_nest_waitall( rq_tot_d, ireq_d )
2075 
2076  if ( online_aggressive_comm ) then
2077  ! nothing to do
2078  else
2079  call mpi_barrier(intercomm_parent, ierr)
2080  endif
2081 
2082  call prof_rapend ('NEST_wait_C', 2)
2083  call prof_rapstart('NEST_unpack_C', 2)
2084 
2085  tagbase = tagcomm + tag_dens*order_tag_var
2086  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2087  work1_recv(:,:,:), & ! [OUT]
2088  tagbase, i_sclr, handle, & ! [IN]
2089  isu_tag, & ! [INOUT]
2090  flag_dens = .true. ) ! [IN]
2091 !OCL XFILL
2092  do j = 1, daughter_ja(handle)
2093  do i = 1, daughter_ia(handle)
2094  do k = datr_ks(handle), datr_ke(handle)
2095  dens_recv(k,i,j) = work1_recv(k,i,j)
2096  enddo
2097  enddo
2098  enddo
2099 
2100  call comm_vars8( dens_recv, 1 )
2101 
2102  tagbase = tagcomm + tag_momz*order_tag_var
2103  if ( online_use_velz ) then
2104  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2105  work2_recv(:,:,:), & ! [OUT]
2106  tagbase, i_zstg, handle, & ! [IN]
2107  isu_tag ) ! [INOUT]
2108 !OCL XFILL
2109  do j = 1, daughter_ja(handle)
2110  do i = 1, daughter_ia(handle)
2111  do k = datr_ks(handle), datr_ke(handle)-1
2112  velz_recv(k,i,j) = work2_recv(k,i,j) / ( work1_recv(k,i,j) + work1_recv(k+1,i,j) ) * 2.0_rp
2113  enddo
2114  enddo
2115  enddo
2116 
2117  do j = 1, daughter_ja(handle)
2118  do i = 1, daughter_ia(handle)
2119  velz_recv(datr_ks(handle)-1,i,j) = 0.0_rp
2120  velz_recv(datr_ke(handle) ,i,j) = 0.0_rp
2121  enddo
2122  enddo
2123  endif
2124 
2125  call comm_wait ( dens_recv, 1, .false. )
2126 
2127  tagbase = tagcomm + tag_momx*order_tag_var
2128  if ( online_no_rotate ) then
2129  ! U_ll_recv receives MOMX
2130  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2131  work1_recv(:,:,:), & ! [OUT]
2132  tagbase, i_xstg, handle, & ! [IN]
2133  isu_tag ) ! [INOUT]
2134  else
2135  ! U_ll_recv receives MOMX/DENS
2136  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2137  u_ll_recv(:,:,:), & ! [OUT]
2138  tagbase, i_sclr, handle, & ! [IN]
2139  isu_tag ) ! [INOUT]
2140  endif
2141 
2142  tagbase = tagcomm + tag_momy*order_tag_var
2143  if ( online_no_rotate ) then
2144  ! V_ll_recv receives MOMY
2145  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2146  work2_recv(:,:,:), & ! [OUT]
2147  tagbase, i_ystg, handle, & ! [IN]
2148  isu_tag ) ! [INOUT]
2149  else
2150  ! V_ll_recv receives MOMY/DENS
2151  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2152  v_ll_recv(:,:,:), & ! [OUT]
2153  tagbase, i_sclr, handle, & ! [IN]
2154  isu_tag ) ! [INOUT]
2155  endif
2156 
2157  if ( online_no_rotate ) then
2158 
2159 !OCL XFILL
2160  do j = 1, daughter_ja(handle)
2161  do i = 1, daughter_ia(handle)-1
2162  do k = datr_ks(handle), datr_ke(handle)
2163  velx_recv(k,i,j) = work1_recv(k,i,j) / ( dens_recv(k,i+1,j) + dens_recv(k,i,j) ) * 2.0_rp
2164  enddo
2165  enddo
2166  enddo
2167 
2168  i = daughter_ia(handle)
2169 !OCL XFILL
2170  do j = 1, daughter_ja(handle)
2171  do k = datr_ks(handle), datr_ke(handle)
2172  velx_recv(k,i,j) = work1_recv(k,i,j) / dens_recv(k,i,j)
2173  enddo
2174  enddo
2175 
2176  call comm_vars8( velx_recv, 2 )
2177 
2178 !OCL XFILL
2179  do j = 1, daughter_ja(handle)-1
2180  do i = 1, daughter_ia(handle)
2181  do k = datr_ks(handle), datr_ke(handle)
2182  vely_recv(k,i,j) = work2_recv(k,i,j) / ( dens_recv(k,i,j+1) + dens_recv(k,i,j) ) * 2.0_rp
2183  enddo
2184  enddo
2185  enddo
2186 
2187  j = daughter_ja(handle)
2188 !OCL XFILL
2189  do i = 1, daughter_ia(handle)
2190  do k = datr_ks(handle), datr_ke(handle)
2191  vely_recv(k,i,j) = work2_recv(k,i,j) / dens_recv(k,i,j)
2192  enddo
2193  enddo
2194 
2195  call comm_vars8( vely_recv, 3 )
2196 
2197  call comm_wait ( velx_recv, 2, .false. )
2198  call comm_wait ( vely_recv, 3, .false. )
2199 
2200  else ! rotate
2201 
2202  ! rotation from latlon field to map-projected field
2203 !OCL XFILL
2204  do j = 1, daughter_ja(handle)
2205  do i = 1, daughter_ia(handle)
2206  do k = datr_ks(handle), datr_ke(handle)
2207  work1_recv(k,i,j) = u_ll_recv(k,i,j) * rotc(i,j,1) + v_ll_recv(k,i,j) * rotc(i,j,2)
2208  work2_recv(k,i,j) = -u_ll_recv(k,i,j) * rotc(i,j,2) + v_ll_recv(k,i,j) * rotc(i,j,1)
2209  enddo
2210  enddo
2211  enddo
2212 
2213  ! from scalar point to staggered point
2214 !OCL XFILL
2215  do j = 1, daughter_ja(handle)
2216  do i = 1, daughter_ia(handle)-1
2217  do k = datr_ks(handle), datr_ke(handle)
2218  velx_recv(k,i,j) = ( work1_recv(k,i+1,j) + work1_recv(k,i,j) ) * 0.5_rp
2219  enddo
2220  enddo
2221  enddo
2222 
2223  i = daughter_ia(handle)
2224 !OCL XFILL
2225  do j = 1, daughter_ja(handle)
2226  do k = datr_ks(handle), datr_ke(handle)
2227  velx_recv(k,i,j) = work1_recv(k,i,j)
2228  enddo
2229  enddo
2230 
2231  call comm_vars8( velx_recv, 2 )
2232 
2233 !OCL XFILL
2234  do j = 1, daughter_ja(handle)-1
2235  do i = 1, daughter_ia(handle)
2236  do k = datr_ks(handle), datr_ke(handle)
2237  vely_recv(k,i,j) = ( work2_recv(k,i,j+1) + work2_recv(k,i,j) ) * 0.5_rp
2238  enddo
2239  enddo
2240  enddo
2241 
2242  j = daughter_ja(handle)
2243 !OCL XFILL
2244  do i = 1, daughter_ia(handle)
2245  do k = datr_ks(handle), datr_ke(handle)
2246  vely_recv(k,i,j) = work2_recv(k,i,j)
2247  enddo
2248  enddo
2249 
2250  call comm_vars8( vely_recv, 3 )
2251 
2252  call comm_wait ( velx_recv, 2, .false. )
2253  call comm_wait ( vely_recv, 3, .false. )
2254 
2255  endif
2256 
2257  tagbase = tagcomm + tag_rhot*order_tag_var
2258  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2259  work1_recv(:,:,:), & ! [OUT]
2260  tagbase, i_sclr, handle, & ! [IN]
2261  isu_tag ) ! [INOUT]
2262 !OCL XFILL
2263  do j = 1, daughter_ja(handle)
2264  do i = 1, daughter_ia(handle)
2265  do k = datr_ks(handle), datr_ke(handle)
2266  pott_recv(k,i,j) = work1_recv(k,i,j) / dens_recv(k,i,j)
2267  enddo
2268  enddo
2269  enddo
2270 
2271  do iq = 1, bnd_qa
2272  tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
2273  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2274  work1_recv(:,:,:), & ! [OUT]
2275  tagbase, i_sclr, handle, & ! [IN]
2276  isu_tag ) ! [INOUT]
2277 !OCL XFILL
2278  do j = 1, daughter_ja(handle)
2279  do i = 1, daughter_ia(handle)
2280  do k = datr_ks(handle), datr_ke(handle)
2281  qtrc_recv(k,i,j,iq) = work1_recv(k,i,j)
2282  enddo
2283  enddo
2284  enddo
2285  enddo
2286 
2287  call prof_rapend ('NEST_unpack_C', 2)
2288  call prof_rapend ('NEST_total_C', 2)
2289 
2290  else
2291  log_error("COMM_CARTESC_NEST_nestdown",*) 'internal error'
2292  call prc_abort
2293  endif
2294 
2295  return

References scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_rotc, comm_cartesc_nest_filiation, comm_cartesc_nest_waitall(), datr_ke, datr_ks, daughter_ia, daughter_ja, intercomm_parent, scale_tracer::k, online_no_rotate, online_use_velz, parent_ia, parent_ja, parent_ka, scale_prc::prc_abort(), scale_prof::prof_rapend(), scale_prof::prof_rapstart(), and use_nesting.

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_send().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ comm_cartesc_nest_recvwait_issue()

subroutine, public scale_comm_cartesc_nest::comm_cartesc_nest_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 2303 of file scale_comm_cartesC_nest.F90.

2303  use scale_prc, only: &
2304  prc_abort
2305  implicit none
2306 
2307  integer, intent(in) :: HANDLE
2308  integer, intent(in) :: BND_QA
2309 
2310  integer :: isu_tag
2311  integer :: tagbase, tagcomm
2312  integer :: ierr
2313  integer :: iq
2314  !---------------------------------------------------------------------------
2315 
2316  if( .NOT. use_nesting ) return
2317 
2318  if ( bnd_qa > i_bndqa ) then
2319  log_error("COMM_CARTESC_NEST_recvwait_issue",*) 'internal error: about BND_QA'
2320  call prc_abort
2321  endif
2322 
2323  tagcomm = intercomm_id(handle) * order_tag_comm
2324 
2325  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2326 
2327  !##### parent [wait issue] #####
2328 
2329  call prof_rapstart('NEST_total_P', 2)
2330  call prof_rapstart('NEST_wait_P', 2)
2331 
2332  nwait_p = nwait_p + 1
2333  !LOG_INFO("COMM_CARTESC_NEST_recvwait_issue",'(1X,A,I5,A)') "NestIDC [P]: que wait ( ", nwait_p, " )"
2334 
2335  call comm_cartesc_nest_issuer_of_wait( handle )
2336 
2337  if ( online_aggressive_comm ) then
2338  ! nothing to do
2339  else
2340  call mpi_barrier(intercomm_daughter, ierr)
2341  endif
2342 
2343  call prof_rapend ('NEST_wait_P', 2)
2344  call prof_rapend ('NEST_total_P', 2)
2345 
2346  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2347 
2348  !##### child [receive issue] #####
2349 
2350  call prof_rapstart('NEST_total_C', 2)
2351 
2352  nrecv = nrecv + 1
2353  log_info("COMM_CARTESC_NEST_recvwait_issue",'(1X,A,I5,A)') "NestIDC [C]: que recv ( ", nrecv, " )"
2354 
2355  !*** reset issue tag and request control
2356  !--- do not change the calling order below;
2357  !--- it should be consistent with the order in "COMM_CARTESC_NEST_nestdown"
2358  isu_tag = 0
2359  rq_ctl_d = 0
2360 
2361  tagbase = tagcomm + tag_dens*order_tag_var
2362  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag )
2363 
2364  tagbase = tagcomm + tag_momz*order_tag_var
2365  if ( online_use_velz ) then
2366  call comm_cartesc_nest_issuer_of_receive( tagbase, i_zstg, handle, isu_tag )
2367  endif
2368 
2369  tagbase = tagcomm + tag_momx*order_tag_var
2370  if ( online_no_rotate ) then
2371  call comm_cartesc_nest_issuer_of_receive( tagbase, i_xstg, handle, isu_tag )
2372  else
2373  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag )
2374  endif
2375 
2376  tagbase = tagcomm + tag_momy*order_tag_var
2377  if ( online_no_rotate ) then
2378  call comm_cartesc_nest_issuer_of_receive( tagbase, i_ystg, handle, isu_tag )
2379  else
2380  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag )
2381  endif
2382 
2383  tagbase = tagcomm + tag_rhot*order_tag_var
2384  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag )
2385 
2386  do iq = 1, bnd_qa
2387  tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
2388  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag )
2389  enddo
2390 
2391  rq_tot_d = rq_ctl_d
2392 
2393  call prof_rapend('NEST_total_C', 2)
2394 
2395  else
2396  log_error("COMM_CARTESC_NEST_recvwait_issue",*) 'internal error'
2397  call prc_abort
2398  endif
2399 
2400  return

References comm_cartesc_nest_filiation, intercomm_daughter, online_no_rotate, online_use_velz, scale_prc::prc_abort(), scale_prof::prof_rapend(), scale_prof::prof_rapstart(), and use_nesting.

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_finalize(), mod_atmos_bnd_driver::atmos_boundary_driver_send(), and mod_atmos_bnd_driver::atmos_boundary_set_file().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ comm_cartesc_nest_recv_cancel()

subroutine, public scale_comm_cartesc_nest::comm_cartesc_nest_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 2407 of file scale_comm_cartesC_nest.F90.

2407  use scale_prc, only: &
2408  prc_abort
2409  implicit none
2410 
2411  integer, intent(in) :: HANDLE
2412 
2413  !logical :: flag
2414  !integer :: istatus(MPI_STATUS_SIZE)
2415 
2416  integer :: rq
2417  integer :: ierr
2418  !---------------------------------------------------------------------------
2419 
2420  if( .NOT. use_nesting ) return
2421 
2422  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2423 
2424  !##### parent #####
2425  ! Nothing to do
2426 
2427  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2428 
2429  !##### child #####
2430 
2431  log_info("COMM_CARTESC_NEST_recv_cancel",'(1X,A,I5,A)') "NestIDC [C]: CANCEL recv ( ", nrecv, " )"
2432 
2433  do rq = 1, rq_tot_d
2434  if ( ireq_d(rq) /= mpi_request_null ) then
2435 
2436  call mpi_cancel(ireq_d(rq), ierr)
2437 
2438 ! call MPI_TEST_CANCELLED(istatus, flag, ierr)
2439 ! if ( .NOT. flag ) then
2440 ! LOG_ERROR("COMM_CARTESC_NEST_recv_cancel",*) 'receive actions do not cancelled, req = ', rq
2441 ! endif
2442  endif
2443  enddo
2444 
2445  else
2446  log_error_cont(*) 'internal error'
2447  call prc_abort
2448  endif
2449 
2450  return

References comm_cartesc_nest_filiation, scale_prc::prc_abort(), and use_nesting.

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_finalize().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ comm_cartesc_nest_intercomm_nestdown_3d()

subroutine scale_comm_cartesc_nest::comm_cartesc_nest_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,
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]flag_densflag of logarithmic interpolation for density

Definition at line 2463 of file scale_comm_cartesC_nest.F90.

2463  use scale_prc, only: &
2464  prc_abort
2465  use scale_comm_cartesc, only: &
2466  comm_datatype
2467  use scale_interp, only: &
2470  real_cz => atmos_grid_cartesc_real_cz, &
2471  real_fz => atmos_grid_cartesc_real_fz
2472  implicit none
2473 
2474  real(RP), intent(in) :: pvar(:,:,:)
2475  real(RP), intent(out) :: dvar(:,:,:)
2476  integer, intent(in) :: tagbase
2477  integer, intent(in) :: id_stag
2478  integer, intent(in) :: HANDLE
2479  integer, intent(inout) :: isu_tag
2480 
2481  logical , intent(in), optional :: flag_dens
2482 
2483  integer :: ileng, tag, target_rank
2484 
2485  integer :: xloc, yloc
2486  integer :: gxs, gxe, gys, gye ! for large domain
2487  integer :: pxs, pxe, pys, pye ! for parent domain
2488  integer :: zs, ze
2489 
2490  integer :: ig, rq, yp
2491  logical :: no_zstag
2492  logical :: spline
2493  logical :: logarithmic
2494 
2495  integer :: ierr
2496  integer :: i, j
2497  !---------------------------------------------------------------------------
2498 
2499  if( .NOT. use_nesting ) return
2500 
2501  logarithmic = .false.
2502  spline = .false.
2503  if ( present(flag_dens) ) then
2504  if( flag_dens ) then
2505  logarithmic = .true.
2506  spline = .true.
2507  end if
2508  endif
2509 
2510  if ( id_stag == i_sclr ) then
2511  no_zstag = .true.
2512  ig = i_sclr
2513  elseif( id_stag == i_zstg ) then
2514  no_zstag = .false.
2515  ig = i_zstg
2516  elseif( id_stag == i_xstg ) then
2517  no_zstag = .true.
2518  ig = i_xstg
2519  elseif( id_stag == i_ystg ) then
2520  no_zstag = .true.
2521  ig = i_ystg
2522  endif
2523 
2524  ileng = parent_ka(handle) * parent_ia(handle) * parent_ja(handle)
2525 
2526  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2527 
2528  !##### parent #####
2529  rq = rq_ctl_p
2530 
2531  do yp = 1, num_yp
2532  rq = rq + 1
2533 
2534  ! send data to multiple daughter processes
2535  target_rank = comm_cartesc_nest_tile_list_yp(yp)
2536  tag = tagbase + yp
2537 
2538  call mpi_isend( pvar, &
2539  ileng, &
2540  comm_datatype, &
2541  target_rank, &
2542  tag, &
2543  intercomm_daughter, &
2544  ireq_p(rq), &
2545  ierr )
2546 
2547  dvar(:,:,:) = -1.0_rp ! input as a dummy value
2548  enddo
2549 
2550  rq_ctl_p = rq
2551 
2552  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2553 
2554  !##### child #####
2555  rq = rq_ctl_d
2556 
2557  do yp = 1, comm_cartesc_nest_tile_all
2558  rq = rq + 1
2559 
2560  xloc = mod( yp-1, comm_cartesc_nest_tile_num_x ) + 1
2561  yloc = int( real(yp-1) / real(comm_cartesc_nest_tile_num_x) ) + 1
2562 
2563  gxs = parent_imax(handle) * (xloc-1) + 1
2564  gxe = parent_imax(handle) * xloc
2565  gys = parent_jmax(handle) * (yloc-1) + 1
2566  gye = parent_jmax(handle) * yloc
2567 
2568  pxs = prnt_is(handle)
2569  pxe = prnt_ie(handle)
2570  pys = prnt_js(handle)
2571  pye = prnt_je(handle)
2572 
2573  isu_tag = isu_tag + 1
2574 
2575  zs = 1
2576  ze = parent_ka(handle)
2577 !OCL XFILL
2578  buffer_ref_3d(zs:ze,gxs:gxe,gys:gye) = recvbuf_3d(zs:ze,pxs:pxe,pys:pye,isu_tag)
2579 
2580  if ( isu_tag > max_isu ) then
2581  log_error("COMM_CARTESC_NEST_intercomm_nestdown_3D",*) 'Exceeded maximum issue'
2582  log_error_cont(*) 'isu_tag = ', isu_tag
2583  call prc_abort
2584  endif
2585 
2586  enddo
2587 
2588  rq_ctl_d = rq
2589 
2590  if ( no_zstag ) then
2591  call interp_interp3d( itp_nh, & ! [IN]
2592  tileal_ka(handle), & ! [IN]
2593  khalo+1, & ! [IN]
2594  tileal_ka(handle)-khalo, & ! [IN]
2595  tileal_ia(handle), & ! [IN]
2596  tileal_ja(handle), & ! [IN]
2597  daughter_ka(handle), & ! [IN]
2598  datr_ks(handle), & ! [IN]
2599  datr_ke(handle), & ! [IN]
2600  daughter_ia(handle), & ! [IN]
2601  daughter_ja(handle), & ! [IN]
2602  igrd( :,:,:,ig), & ! [IN]
2603  jgrd( :,:,:,ig), & ! [IN]
2604  hfact( :,:,:,ig), & ! [IN]
2605  kgrd(:,:,:,:,:,ig), & ! [IN]
2606  vfact(:, :,:,:,ig), & ! [IN]
2607  buffer_ref_cz(:,:,:), & ! [IN]
2608  real_cz(:,:,:), & ! [IN]
2609  buffer_ref_3d(:,:,:), & ! [IN]
2610  dvar(:,:,:), & ! [OUT]
2611  spline = spline, & ! [IN, optional]
2612  logwgt = logarithmic ) ! [IN, optional]
2613 
2614  else
2615  call interp_interp3d( itp_nh, & ! [IN]
2616  tileal_ka(handle) , & ! [IN]
2617  khalo, & ! [IN]
2618  tileal_ka(handle)-khalo, & ! [IN]
2619  tileal_ia(handle), & ! [IN]
2620  tileal_ja(handle), & ! [IN]
2621  daughter_ka(handle), & ! [IN]
2622  datr_ks(handle), & ! [IN]
2623  datr_ke(handle), & ! [IN]
2624  daughter_ia(handle), & ! [IN]
2625  daughter_ja(handle), & ! [IN]
2626  igrd( :,:,:,ig), & ! [IN]
2627  jgrd( :,:,:,ig), & ! [IN]
2628  hfact( :,:,:,ig), & ! [IN]
2629  kgrd(:,:,:,:,:,ig), & ! [IN]
2630  vfact(:, :,:,:,ig), & ! [IN]
2631  buffer_ref_fz(:,:,:), & ! [IN]
2632  real_fz(1:,:,:), & ! [IN]
2633  buffer_ref_3d(:,:,:), & ! [IN]
2634  dvar(:,:,:), & ! [OUT]
2635  spline = spline, & ! [IN, optional]
2636  logwgt = logarithmic ) ! [IN, optional]
2637  endif
2638 
2639  do j = 1, daughter_ja(handle)
2640  do i = 1, daughter_ia(handle)
2641  dvar( 1:datr_ks(handle)-1,i,j) = 0.0_rp
2642  dvar(datr_ke(handle)+1:daughter_ka(handle) ,i,j) = 0.0_rp
2643  enddo
2644  enddo
2645 
2646  else
2647  log_error("COMM_CARTESC_NEST_intercomm_nestdown_3D",*) 'internal error'
2648  call prc_abort
2649  endif
2650 
2651  return

References scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz, comm_cartesc_nest_filiation, comm_cartesc_nest_tile_num_x, scale_comm_cartesc::comm_datatype, datr_ke, datr_ks, daughter_ia, daughter_ja, daughter_ka, intercomm_daughter, scale_interp::interp_interp3d(), scale_atmos_grid_cartesc_index::khalo, parent_ia, parent_imax, parent_ja, parent_jmax, parent_ka, scale_prc::prc_abort(), prnt_ie, prnt_is, prnt_je, prnt_js, tileal_ia, tileal_ja, tileal_ka, and use_nesting.

Here is the call graph for this function:

◆ comm_cartesc_nest_issuer_of_receive_3d()

subroutine scale_comm_cartesc_nest::comm_cartesc_nest_issuer_of_receive_3d ( integer, intent(in)  tagbase,
integer, intent(in)  id_stag,
integer, intent(in)  HANDLE,
integer, intent(inout)  isu_tag 
)

[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

Definition at line 2661 of file scale_comm_cartesC_nest.F90.

2661  use scale_prc, only: &
2662  prc_myrank, &
2663  prc_abort
2664  use scale_comm_cartesc, only: &
2665  comm_datatype
2666  implicit none
2667 
2668  integer, intent(in) :: tagbase
2669  integer, intent(in) :: id_stag
2670  integer, intent(in) :: HANDLE
2671  integer, intent(inout) :: isu_tag
2672 
2673  integer :: ierr, ileng
2674  integer :: tag, target_rank
2675 
2676  integer :: rq, yp
2677  !---------------------------------------------------------------------------
2678 
2679  if( .NOT. use_nesting ) return
2680 
2681  ileng = parent_ka(handle) * parent_ia(handle) * parent_ja(handle)
2682 
2683  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2684 
2685  !##### parent #####
2686  ! nothing to do
2687 
2688  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2689 
2690  !##### child #####
2691  rq = rq_ctl_d
2692 
2693  do yp = 1, comm_cartesc_nest_tile_all
2694  rq = rq + 1
2695 
2696  target_rank = comm_cartesc_nest_tile_list_d(yp,prc_myrank+1)
2697  tag = tagbase + call_order(yp)
2698 
2699  isu_tag = isu_tag + 1
2700 
2701  recvbuf_3d(:,:,:,isu_tag) = 0.0_rp
2702 
2703  call mpi_irecv( recvbuf_3d(:,:,:,isu_tag), &
2704  ileng, &
2705  comm_datatype, &
2706  target_rank, &
2707  tag, &
2708  intercomm_parent, &
2709  ireq_d(rq), &
2710  ierr )
2711 
2712  enddo
2713 
2714  if ( isu_tag > max_isu ) then
2715  log_error("COMM_CARTESC_NEST_issuer_of_receive_3D",*) 'Exceeded maximum issue'
2716  log_error_cont(*) 'isu_tag = ', isu_tag
2717  call prc_abort
2718  endif
2719 
2720  rq_ctl_d = rq
2721 
2722  else
2723  log_error("COMM_CARTESC_NEST_issuer_of_receive_3D",*) 'internal error'
2724  call prc_abort
2725  endif
2726 
2727  return

References comm_cartesc_nest_filiation, scale_comm_cartesc::comm_datatype, intercomm_parent, parent_ia, parent_ja, parent_ka, scale_prc::prc_abort(), scale_prc::prc_myrank, and use_nesting.

Here is the call graph for this function:

◆ comm_cartesc_nest_issuer_of_wait_3d()

subroutine scale_comm_cartesc_nest::comm_cartesc_nest_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 2734 of file scale_comm_cartesC_nest.F90.

2734  use scale_prc, only: &
2735  prc_abort
2736  implicit none
2737 
2738  integer, intent(in) :: HANDLE
2739  !---------------------------------------------------------------------------
2740 
2741  if( .NOT. use_nesting ) return
2742 
2743  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2744 
2745  !##### parent #####
2746  call comm_cartesc_nest_waitall( rq_tot_p, ireq_p(:) )
2747 
2748  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2749 
2750  !##### child #####
2751  ! nothing to do
2752 
2753  else
2754  log_error("COMM_CARTESC_NEST_issuer_of_wait_3D",*) 'internal error'
2755  call prc_abort
2756  endif
2757 
2758  return

References comm_cartesc_nest_filiation, comm_cartesc_nest_waitall(), scale_prc::prc_abort(), and use_nesting.

Here is the call graph for this function:

◆ comm_cartesc_nest_waitall()

subroutine scale_comm_cartesc_nest::comm_cartesc_nest_waitall ( integer, intent(in)  req_count,
integer, dimension(max_rq), intent(inout)  ireq 
)

[substance of comm_wait] Inter-communication

Definition at line 2766 of file scale_comm_cartesC_nest.F90.

2766  use scale_prc, only: &
2767  prc_abort
2768  implicit none
2769 
2770  integer, intent(in) :: req_count
2771  integer, intent(inout) :: ireq(max_rq)
2772 
2773  integer :: i
2774  integer :: ierr
2775  integer :: istatus(MPI_STATUS_SIZE,req_count)
2776  integer :: req_count2
2777  integer :: ireq2(max_rq)
2778 
2779 ! logical :: flag = .false.
2780 ! integer(8) :: num = 0
2781  !---------------------------------------------------------------------------
2782 
2783  if( .NOT. use_nesting ) return
2784 
2785  req_count2 = 0
2786  do i = 1, req_count
2787  if ( ireq(i) /= mpi_request_null ) then
2788  req_count2 = req_count2 + 1
2789  ireq2(req_count2) = ireq(i)
2790  endif
2791  enddo
2792 
2793  if( req_count2 /= 0 ) call mpi_waitall( req_count2, ireq2(1:req_count2), istatus, ierr )
2794 
2795 ! do while ( .NOT. flag )
2796 ! num = num + 1
2797 ! call MPI_TESTALL( req_count, ireq, flag, istatus, ierr )
2798 !
2799 ! if ( num > ONLINE_WAIT_LIMIT ) then
2800 ! LOG_ERROR("COMM_CARTESC_NEST_waitall",'(1x,A)') 'over the limit of waiting time [NESTCOM]'
2801 ! LOG_ERROR_CONT('(1x,A)') 'over the limit of waiting time [NESTCOM]'
2802 ! call PRC_abort
2803 ! endif
2804 ! enddo
2805 
2806  return

References scale_prc::prc_abort(), and use_nesting.

Referenced by comm_cartesc_nest_domain_shape(), comm_cartesc_nest_issuer_of_wait_3d(), and comm_cartesc_nest_nestdown().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ comm_cartesc_nest_test()

subroutine, public scale_comm_cartesc_nest::comm_cartesc_nest_test ( integer, intent(in)  HANDLE)

[check communication status] Inter-communication

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

Definition at line 2813 of file scale_comm_cartesC_nest.F90.

2813  use scale_prc, only: &
2814  prc_abort
2815  implicit none
2816 
2817  integer, intent(in) :: HANDLE
2818 
2819  integer :: istatus(MPI_STATUS_SIZE)
2820  integer :: ierr
2821  logical :: flag
2822  !---------------------------------------------------------------------------
2823 
2824  if( .NOT. use_nesting ) return
2825 
2826  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2827 
2828  !##### parent #####
2829 
2830  call prof_rapstart('NEST_test_P', 2)
2831  if ( rq_ctl_p > 0 ) call mpi_test(ireq_p(1), flag, istatus, ierr)
2832  call prof_rapend('NEST_test_P', 2)
2833 
2834  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2835 
2836  !##### child #####
2837 
2838  call prof_rapstart('NEST_test_C', 2)
2839  if ( rq_ctl_d > 0 ) call mpi_test(ireq_d(1), flag, istatus, ierr)
2840  call prof_rapend('NEST_test_C', 2)
2841 
2842  else
2843  log_error("COMM_CARTESC_NEST_test",*) 'error'
2844  call prc_abort
2845  endif
2846 
2847  return

References comm_cartesc_nest_filiation, scale_prc::prc_abort(), scale_prof::prof_rapend(), scale_prof::prof_rapstart(), and use_nesting.

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_send(), and mod_atmos_bnd_driver::atmos_boundary_driver_update().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ comm_cartesc_nest_disconnect()

subroutine, public scale_comm_cartesc_nest::comm_cartesc_nest_disconnect

[finalize: disconnect] Inter-communication

Definition at line 2853 of file scale_comm_cartesC_nest.F90.

2853  use scale_prc, only: &
2855  implicit none
2856 
2857  integer :: ierr
2858  !---------------------------------------------------------------------------
2859 
2860  if( .NOT. use_nesting ) return
2861 
2862  log_info("COMM_CARTESC_NEST_disconnect",'(1x,A)') 'Waiting finish of whole processes'
2863  call mpi_barrier(prc_global_comm_world, ierr)
2864 
2865  if ( online_iam_parent ) then
2866  !LOG_INFO("COMM_CARTESC_NEST_disconnect",'(1x,A)') 'Waiting finish of whole processes as a parent'
2867  !call MPI_BARRIER(INTERCOMM_DAUGHTER, ierr)
2868  call mpi_comm_free(intercomm_daughter, ierr)
2869  log_info("COMM_CARTESC_NEST_disconnect",'(1x,A)') 'Disconnected communication with child'
2870  endif
2871 
2872  if ( online_iam_daughter ) then
2873  !LOG_INFO("COMM_CARTESC_NEST_disconnect",'(1x,A)') 'Waiting finish of whole processes as a child'
2874  !call MPI_BARRIER(INTERCOMM_PARENT, ierr)
2875  call mpi_comm_free(intercomm_parent, ierr)
2876  log_info("COMM_CARTESC_NEST_disconnect",'(1x,A)') 'Disconnected communication with parent'
2877  endif
2878 
2879  return

References intercomm_daughter, intercomm_parent, online_iam_daughter, online_iam_parent, scale_prc::prc_global_comm_world, and use_nesting.

Referenced by mod_atmos_driver::atmos_driver_finalize().

Here is the caller graph for this function:

Variable Documentation

◆ intercomm_parent

integer, public scale_comm_cartesc_nest::intercomm_parent

Definition at line 45 of file scale_comm_cartesC_nest.F90.

45  integer, public :: INTERCOMM_PARENT ! inter-communicator to parent

Referenced by comm_cartesc_nest_disconnect(), comm_cartesc_nest_domain_shape(), comm_cartesc_nest_issuer_of_receive_3d(), comm_cartesc_nest_nestdown(), and comm_cartesc_nest_setup().

◆ intercomm_daughter

integer, public scale_comm_cartesc_nest::intercomm_daughter

Definition at line 46 of file scale_comm_cartesC_nest.F90.

46  integer, public :: INTERCOMM_DAUGHTER ! inter-communicator to daughter

Referenced by comm_cartesc_nest_disconnect(), comm_cartesc_nest_domain_shape(), comm_cartesc_nest_intercomm_nestdown_3d(), comm_cartesc_nest_recvwait_issue(), and comm_cartesc_nest_setup().

◆ comm_cartesc_nest_filiation

integer, dimension(10), public scale_comm_cartesc_nest::comm_cartesc_nest_filiation

◆ handling_num

integer, public scale_comm_cartesc_nest::handling_num

handing number of nesting relation

Definition at line 49 of file scale_comm_cartesC_nest.F90.

49  integer, public :: HANDLING_NUM

Referenced by comm_cartesc_nest_setup().

◆ comm_cartesc_nest_tile_num_x

integer, public scale_comm_cartesc_nest::comm_cartesc_nest_tile_num_x

◆ comm_cartesc_nest_tile_num_y

integer, public scale_comm_cartesc_nest::comm_cartesc_nest_tile_num_y

◆ comm_cartesc_nest_tile_id

integer, dimension(:), allocatable, public scale_comm_cartesc_nest::comm_cartesc_nest_tile_id

◆ parent_kmax

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

parent max number in z-direction

Definition at line 54 of file scale_comm_cartesC_nest.F90.

54  integer, public :: PARENT_KMAX(2)

Referenced by comm_cartesc_nest_domain_shape(), comm_cartesc_nest_setup(), and mod_realinput_scale::parentatmossetupscale().

◆ parent_imax

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

◆ parent_jmax

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

◆ parent_ka

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

◆ parent_ia

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

◆ parent_ja

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

◆ parent_okmax

integer, dimension(2), public scale_comm_cartesc_nest::parent_okmax

parent max number in oz-direction

Definition at line 60 of file scale_comm_cartesC_nest.F90.

60  integer, public :: PARENT_OKMAX(2)

Referenced by comm_cartesc_nest_setup().

◆ parent_lkmax

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

parent max number in lz-direction

Definition at line 61 of file scale_comm_cartesC_nest.F90.

61  integer, public :: PARENT_LKMAX(2)

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

◆ parent_dtsec

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

parent DT [sec]

Definition at line 62 of file scale_comm_cartesC_nest.F90.

62  real(DP), public :: PARENT_DTSEC(2)

Referenced by mod_atmos_bnd_driver::atmos_boundary_set_file(), comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_setup().

◆ parent_nstep

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

parent step [number]

Definition at line 63 of file scale_comm_cartesC_nest.F90.

63  integer, public :: PARENT_NSTEP(2)

Referenced by mod_atmos_bnd_driver::atmos_boundary_set_online(), comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_setup().

◆ daughter_kmax

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

daughter max number in z-direction

Definition at line 65 of file scale_comm_cartesC_nest.F90.

65  integer, public :: DAUGHTER_KMAX(2)

Referenced by comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_setup().

◆ daughter_imax

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

daughter max number in x-direction

Definition at line 66 of file scale_comm_cartesC_nest.F90.

66  integer, public :: DAUGHTER_IMAX(2)

Referenced by comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_setup().

◆ daughter_jmax

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

daughter max number in y-direction

Definition at line 67 of file scale_comm_cartesC_nest.F90.

67  integer, public :: DAUGHTER_JMAX(2)

Referenced by comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_setup().

◆ daughter_ka

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

daughter max number in z-direction (with halo)

Definition at line 68 of file scale_comm_cartesC_nest.F90.

68  integer, public :: DAUGHTER_KA(2)

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_send(), comm_cartesc_nest_intercomm_nestdown_3d(), and comm_cartesc_nest_setup().

◆ daughter_ia

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

daughter max number in x-direction (with halo)

Definition at line 69 of file scale_comm_cartesC_nest.F90.

69  integer, public :: DAUGHTER_IA(2)

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_send(), comm_cartesc_nest_intercomm_nestdown_3d(), comm_cartesc_nest_nestdown(), and comm_cartesc_nest_setup().

◆ daughter_ja

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

daughter max number in y-direction (with halo)

Definition at line 70 of file scale_comm_cartesC_nest.F90.

70  integer, public :: DAUGHTER_JA(2)

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_send(), comm_cartesc_nest_intercomm_nestdown_3d(), comm_cartesc_nest_nestdown(), and comm_cartesc_nest_setup().

◆ daughter_okmax

integer, dimension(2), public scale_comm_cartesc_nest::daughter_okmax

daughter max number in oz-direction

Definition at line 71 of file scale_comm_cartesC_nest.F90.

71  integer, public :: DAUGHTER_OKMAX(2)

◆ daughter_lkmax

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

daughter max number in lz-direction

Definition at line 72 of file scale_comm_cartesC_nest.F90.

72  integer, public :: DAUGHTER_LKMAX(2)

◆ daughter_dtsec

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

daughter DT [sec]

Definition at line 73 of file scale_comm_cartesC_nest.F90.

73  real(DP), public :: DAUGHTER_DTSEC(2)

Referenced by comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_setup().

◆ daughter_nstep

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

daughter steps [number]

Definition at line 74 of file scale_comm_cartesC_nest.F90.

74  integer, public :: DAUGHTER_NSTEP(2)

Referenced by comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_setup().

◆ prnt_ks

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

start index in z-direction in parent

Definition at line 76 of file scale_comm_cartesC_nest.F90.

76  integer, public :: PRNT_KS(2)

Referenced by comm_cartesc_nest_domain_shape().

◆ prnt_ke

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

end index in z-direction in parent

Definition at line 77 of file scale_comm_cartesC_nest.F90.

77  integer, public :: PRNT_KE(2)

Referenced by comm_cartesc_nest_domain_shape().

◆ prnt_is

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

start index in x-direction in parent

Definition at line 78 of file scale_comm_cartesC_nest.F90.

78  integer, public :: PRNT_IS(2)

Referenced by comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_intercomm_nestdown_3d().

◆ prnt_ie

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

end index in x-direction in parent

Definition at line 79 of file scale_comm_cartesC_nest.F90.

79  integer, public :: PRNT_IE(2)

Referenced by comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_intercomm_nestdown_3d().

◆ prnt_js

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

start index in y-direction in parent

Definition at line 80 of file scale_comm_cartesC_nest.F90.

80  integer, public :: PRNT_JS(2)

Referenced by comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_intercomm_nestdown_3d().

◆ prnt_je

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

end index in y-direction in parent

Definition at line 81 of file scale_comm_cartesC_nest.F90.

81  integer, public :: PRNT_JE(2)

Referenced by comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_intercomm_nestdown_3d().

◆ datr_ks

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

start index in z-direction in daughter

Definition at line 83 of file scale_comm_cartesC_nest.F90.

83  integer, public :: DATR_KS(2)

Referenced by comm_cartesc_nest_domain_shape(), comm_cartesc_nest_intercomm_nestdown_3d(), comm_cartesc_nest_nestdown(), and comm_cartesc_nest_setup().

◆ datr_ke

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

end index in z-direction in daughter

Definition at line 84 of file scale_comm_cartesC_nest.F90.

84  integer, public :: DATR_KE(2)

Referenced by comm_cartesc_nest_domain_shape(), comm_cartesc_nest_intercomm_nestdown_3d(), comm_cartesc_nest_nestdown(), and comm_cartesc_nest_setup().

◆ datr_is

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

start index in x-direction in daughter

Definition at line 85 of file scale_comm_cartesC_nest.F90.

85  integer, public :: DATR_IS(2)

Referenced by comm_cartesc_nest_domain_shape().

◆ datr_ie

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

end index in x-direction in daughter

Definition at line 86 of file scale_comm_cartesC_nest.F90.

86  integer, public :: DATR_IE(2)

Referenced by comm_cartesc_nest_domain_shape().

◆ datr_js

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

start index in y-direction in daughter

Definition at line 87 of file scale_comm_cartesC_nest.F90.

87  integer, public :: DATR_JS(2)

Referenced by comm_cartesc_nest_domain_shape().

◆ datr_je

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

end index in y-direction in daughter

Definition at line 88 of file scale_comm_cartesC_nest.F90.

88  integer, public :: DATR_JE(2)

Referenced by comm_cartesc_nest_domain_shape().

◆ tileal_ka

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

cells of all tiles in z-direction

Definition at line 90 of file scale_comm_cartesC_nest.F90.

90  integer, public :: TILEAL_KA(2)

Referenced by comm_cartesc_nest_intercomm_nestdown_3d(), and comm_cartesc_nest_setup().

◆ tileal_ia

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

cells of all tiles in x-direction

Definition at line 91 of file scale_comm_cartesC_nest.F90.

91  integer, public :: TILEAL_IA(2)

Referenced by comm_cartesc_nest_intercomm_nestdown_3d(), and comm_cartesc_nest_setup().

◆ tileal_ja

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

cells of all tiles in y-direction

Definition at line 92 of file scale_comm_cartesC_nest.F90.

92  integer, public :: TILEAL_JA(2)

Referenced by comm_cartesc_nest_intercomm_nestdown_3d(), and comm_cartesc_nest_setup().

◆ comm_cartesc_nest_bnd_qa

integer, public scale_comm_cartesc_nest::comm_cartesc_nest_bnd_qa = 1

◆ comm_cartesc_nest_interp_level

integer, public scale_comm_cartesc_nest::comm_cartesc_nest_interp_level = 5

horizontal interpolation level

Definition at line 95 of file scale_comm_cartesC_nest.F90.

95  integer, public :: COMM_CARTESC_NEST_INTERP_LEVEL = 5

Referenced by comm_cartesc_nest_setup(), mod_copytopo::copytopo(), mod_realinput::realinput_atmos(), and mod_realinput::realinput_surface().

◆ comm_cartesc_nest_interp_weight_order

integer, public scale_comm_cartesc_nest::comm_cartesc_nest_interp_weight_order = 2

horizontal interpolation weight order

Definition at line 96 of file scale_comm_cartesC_nest.F90.

96  integer, public :: COMM_CARTESC_NEST_INTERP_WEIGHT_ORDER = 2

Referenced by comm_cartesc_nest_setup().

◆ use_nesting

logical, public scale_comm_cartesc_nest::use_nesting = .false.

◆ offline

logical, public scale_comm_cartesc_nest::offline = .false.

Definition at line 99 of file scale_comm_cartesC_nest.F90.

99  logical, public :: OFFLINE = .false.

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_setup(), and comm_cartesc_nest_setup().

◆ online_iam_parent

logical, public scale_comm_cartesc_nest::online_iam_parent = .false.

a flag to say "I am a parent"

Definition at line 100 of file scale_comm_cartesC_nest.F90.

100  logical, public :: ONLINE_IAM_PARENT = .false.

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_setup(), comm_cartesc_nest_disconnect(), and comm_cartesc_nest_setup().

◆ online_iam_daughter

logical, public scale_comm_cartesc_nest::online_iam_daughter = .false.

a flag to say "I am a daughter"

Definition at line 101 of file scale_comm_cartesC_nest.F90.

101  logical, public :: ONLINE_IAM_DAUGHTER = .false.

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_setup(), comm_cartesc_nest_disconnect(), and comm_cartesc_nest_setup().

◆ online_domain_num

integer, public scale_comm_cartesc_nest::online_domain_num = 1

Definition at line 102 of file scale_comm_cartesC_nest.F90.

102  integer, public :: ONLINE_DOMAIN_NUM = 1

Referenced by comm_cartesc_nest_domain_relate(), comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_setup().

◆ online_use_velz

logical, public scale_comm_cartesc_nest::online_use_velz = .false.

◆ online_no_rotate

logical, public scale_comm_cartesc_nest::online_no_rotate = .false.

Definition at line 104 of file scale_comm_cartesC_nest.F90.

104  logical, public :: ONLINE_NO_ROTATE = .false.

Referenced by comm_cartesc_nest_domain_shape(), comm_cartesc_nest_nestdown(), comm_cartesc_nest_recvwait_issue(), and comm_cartesc_nest_setup().

◆ online_boundary_use_qhyd

logical, public scale_comm_cartesc_nest::online_boundary_use_qhyd = .false.

Definition at line 105 of file scale_comm_cartesC_nest.F90.

105  logical, public :: ONLINE_BOUNDARY_USE_QHYD = .false.

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_setup(), and comm_cartesc_nest_setup().

◆ online_boundary_diagqhyd

logical, public scale_comm_cartesc_nest::online_boundary_diagqhyd = .false.

Definition at line 106 of file scale_comm_cartesC_nest.F90.

106  logical, public :: ONLINE_BOUNDARY_DIAGQHYD = .false.

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_send(), mod_atmos_bnd_driver::atmos_boundary_driver_setup(), and comm_cartesc_nest_domain_shape().

scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_cz
geopotential height [m] (zxy)
Definition: scale_atmos_grid_cartesC_real.F90:38
scale_atmos_grid_cartesc_index::ihalo
integer, public ihalo
Definition: scale_atmos_grid_cartesC_index.F90:44
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonxv
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lonxv
longitude at staggered point (xv) [rad,0-2pi]
Definition: scale_atmos_grid_cartesC_real.F90:50
scale_interp
module INTERPOLATION
Definition: scale_interp.F90:12
scale_atmos_grid_cartesc_metric
module Atmosphere Grid CartesianC metirc
Definition: scale_atmos_grid_cartesC_metric.F90:12
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_atmos_grid_cartesc::atmos_grid_cartesc_fx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fx
face coordinate [m]: x, local
Definition: scale_atmos_grid_cartesC.F90:57
scale_prc::prc_myrank
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:90
scale_atmos_grid_cartesc_index::khalo
integer, parameter, public khalo
Definition: scale_atmos_grid_cartesC_index.F90:43
scale_atmos_hydrometeor::atmos_hydrometeor_dry
logical, public atmos_hydrometeor_dry
Definition: scale_atmos_hydrometeor.F90:97
scale_atmos_grid_cartesc_real
module Atmosphere GRID CartesC Real(real space)
Definition: scale_atmos_grid_cartesC_real.F90:11
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lon
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lon
longitude [rad,0-2pi]
Definition: scale_atmos_grid_cartesC_real.F90:48
scale_file
module file
Definition: scale_file.F90:15
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_mapprojection
module Map projection
Definition: scale_mapprojection.F90:12
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latuv
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latuv
latitude at staggered point (uv) [rad,-pi,pi]
Definition: scale_atmos_grid_cartesC_real.F90:55
scale_atmos_grid_cartesc::atmos_grid_cartesc_fy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fy
face coordinate [m]: y, local
Definition: scale_atmos_grid_cartesC.F90:58
scale_file::file_open
subroutine, public file_open(basename, fid, mode, single, aggregate, rankid, postfix)
Definition: scale_file.F90:487
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latuy
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latuy
latitude at staggered point (uy) [rad,-pi,pi]
Definition: scale_atmos_grid_cartesC_real.F90:53
scale_atmos_grid_cartesc_index::jhalo
integer, public jhalo
Definition: scale_atmos_grid_cartesC_index.F90:45
scale_interp::interp_interp3d
subroutine, public interp_interp3d(npoints, KA_ref, KS_ref, KE_ref, IA_ref, JA_ref, KA, KS, KE, IA, JA, idx_i, idx_j, hfact, idx_k, vfact, hgt_ref, hgt, val_ref, val, spline, logwgt, threshold_undef, wsum, val2)
Definition: scale_interp.F90:1322
scale_prc::prc_global_domainid
integer, public prc_global_domainid
my domain ID in global communicator
Definition: scale_prc.F90:84
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latxv
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latxv
latitude at staggered point (xv) [rad,-pi,pi]
Definition: scale_atmos_grid_cartesC_real.F90:54
scale_comm_cartesc
module COMMUNICATION
Definition: scale_comm_cartesC.F90:11
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonuv
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lonuv
longitude at staggered point (uv) [rad,0-2pi]
Definition: scale_atmos_grid_cartesC_real.F90:51
scale_atmos_grid_cartesc::atmos_grid_cartesc_cy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cy
center coordinate [m]: y, local
Definition: scale_atmos_grid_cartesC.F90:56
scale_interp::interp_setup
subroutine, public interp_setup(weight_order, search_limit)
Setup.
Definition: scale_interp.F90:88
scale_const::const_d2r
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lat
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lat
latitude [rad,-pi,pi]
Definition: scale_atmos_grid_cartesC_real.F90:52
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_fz
geopotential height [m] (wxy)
Definition: scale_atmos_grid_cartesC_real.F90:42
scale_prc::prc_global_comm_world
integer, public prc_global_comm_world
global communicator
Definition: scale_prc.F90:79
scale_atmos_grid_cartesc
module atmosphere / grid / cartesC
Definition: scale_atmos_grid_cartesC.F90:12
scale_atmos_grid_cartesc::atmos_grid_cartesc_cx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cx
center coordinate [m]: x, local
Definition: scale_atmos_grid_cartesC.F90:55
scale_prc::prc_ismaster
logical, public prc_ismaster
master process in local communicator?
Definition: scale_prc.F90:91
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonuy
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lonuy
longitude at staggered point (uy) [rad,0-2pi]
Definition: scale_atmos_grid_cartesC_real.F90:49
scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_rotc
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_metric_rotc
rotation coefficient
Definition: scale_atmos_grid_cartesC_metric.F90:35