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, isu_tagf, flag_dens)
 Inter-communication from parent to daughter: nestdown. More...
 
subroutine comm_cartesc_nest_issuer_of_receive_3d (tagbase, id_stag, HANDLE, isu_tag, isu_tagf)
 [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_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.

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

Definition at line 268 of file scale_comm_cartesC_nest.F90.

References 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, intercomm_daughter, intercomm_parent, scale_interp::interp_factor3d(), 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::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().

268  use scale_file, only: &
269  file_open, &
270  file_get_attribute, &
271  file_get_shape
272  use scale_const, only: &
273  d2r => const_d2r
274  use scale_prc, only: &
275  prc_abort, &
278  use scale_interp, only: &
279  interp_setup, &
281  use scale_comm_cartesc, only: &
282  comm_bcast
283  use scale_atmos_grid_cartesc_real, only: &
294  use scale_atmos_hydrometeor, only: &
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 
304  character(len=H_LONG) :: latlon_catalogue_fname = 'latlon_domain_catalogue.txt'
305 
306  integer :: online_specified_maxrq = 0
307  integer :: i
308  integer :: fid, ierr
309  integer :: parent_id
310 
311  logical :: flag_parent = .false.
312  logical :: flag_child = .false.
313 
314  integer :: imaxg(1), jmaxg(1)
315  integer :: pnum_x(1), pnum_y(1)
316  integer :: dims(1)
317  logical :: error, existed
318 
319  namelist / param_comm_cartesc_nest / &
320  latlon_catalogue_fname, &
321  offline_parent_basename, &
322  offline_parent_prc_num_x, &
323  offline_parent_prc_num_y, &
324  online_domain_num, &
325  online_iam_parent, &
326  online_iam_daughter, &
327  online_use_velz, &
328  online_no_rotate, &
329  online_boundary_use_qhyd, &
330  online_aggressive_comm, &
331  online_wait_limit, &
332  online_specified_maxrq, &
333  comm_cartesc_nest_interp_level, &
334  comm_cartesc_nest_interp_weight_order
335 
336  !---------------------------------------------------------------------------
337 
338  log_newline
339  log_info("COMM_CARTESC_NEST_setup",*) 'Setup'
340 
341  if ( present(inter_parent) ) then
342  if( inter_parent /= mpi_comm_null ) flag_child = .true. ! exist parent, so work as a child
343  endif
344  if ( present(inter_child) ) then
345  if( inter_child /= mpi_comm_null ) flag_parent = .true. ! exist child, so work as a parent
346  endif
347 
348  offline_parent_basename = ""
349 
350  nwait_p = 0
351  nwait_d = 0
352  nrecv = 0
353  nsend = 0
354 
355  handling_num = 0
356  comm_cartesc_nest_filiation(:) = 0
357  online_wait_limit = 999999999
358  online_aggressive_comm = .true.
359 
360  !--- read namelist
361  rewind(io_fid_conf)
362  read(io_fid_conf,nml=param_comm_cartesc_nest,iostat=ierr)
363  if( ierr < 0 ) then !--- missing
364  log_info("COMM_CARTESC_NEST_setup",*) 'Not found namelist. Default used.'
365  elseif( ierr > 0 ) then !--- fatal error
366  log_error("COMM_CARTESC_NEST_setup",*) 'Not appropriate names in namelist PARAM_COMM_CARTESC_NEST. Check!'
367  call prc_abort
368  endif
369  log_nml(param_comm_cartesc_nest)
370 
371  prc_global_domainid = online_domain_num
372 
373  if ( offline_parent_basename /= "" ) then
374 
375  offline = .true.
376  use_nesting = .true.
377 
378  if ( prc_ismaster ) then
379  call file_open( offline_parent_basename, & ! (in)
380  fid, & ! (out)
381  aggregate = .false. ) ! (in)
382 
383  call file_get_attribute( fid, "global", "scale_atmos_grid_cartesC_index_imaxg", &
384  imaxg(:), existed=existed )
385  if ( existed ) then
386  call file_get_attribute( fid, "global", "scale_cartesC_prc_num_x", &
387  pnum_x(:), existed=existed )
388  end if
389  if ( existed ) then
390  offline_parent_imax = imaxg(1) / pnum_x(1)
391  else
392  ! for old file
393  call file_get_shape( fid, "CX", dims(:) )
394  offline_parent_imax = dims(1)-ihalo*2
395  end if
396 
397  call file_get_attribute( fid, "global", "scale_atmos_grid_cartesC_index_jmaxg", &
398  jmaxg(:), existed=existed )
399  if ( existed ) then
400  call file_get_attribute( fid, "global", "scale_cartesC_prc_num_y", &
401  pnum_y(:), existed=existed )
402  end if
403  if ( existed ) then
404  offline_parent_jmax = jmaxg(1) / pnum_y(1)
405  else
406  ! for old file
407  call file_get_shape( fid, "CY", dims(:) )
408  offline_parent_jmax = dims(1)-jhalo*2
409  end if
410 
411  call file_get_attribute( fid, "global", "scale_atmos_grid_cartesC_index_kmax", &
412  dims(:), existed=existed )
413  if ( existed ) then
414  offline_parent_kmax = dims(1)
415  else
416  call file_get_shape( fid, "z", dims(:), error=error )
417  if ( error ) then
418  offline_parent_kmax = 0
419  else
420  offline_parent_kmax = dims(1)
421  endif
422  end if
423 
424  call file_get_attribute( fid, "global", "scale_ocean_grid_cartesC_index_kmax", &
425  dims(:), existed=existed )
426  if ( existed ) then
427  offline_parent_okmax = dims(1)
428  else
429  ! for old file
430  call file_get_shape( fid, "oz", dims(:), error=error )
431  if ( error ) then
432  offline_parent_okmax = 0
433  else
434  offline_parent_okmax = dims(1)
435  endif
436  end if
437 
438  call file_get_attribute( fid, "global", "scale_land_grid_cartesC_index_kmax", &
439  dims(:), existed=existed )
440  if ( existed ) then
441  offline_parent_lkmax = dims(1)
442  else
443  ! for old file
444  call file_get_shape( fid, "lz", dims(:), error=error )
445  if ( error ) then
446  offline_parent_lkmax = 0
447  else
448  offline_parent_lkmax = dims(1)
449  endif
450  end if
451 
452  endif
453  call comm_bcast( offline_parent_imax )
454  call comm_bcast( offline_parent_jmax )
455  call comm_bcast( offline_parent_kmax )
456  call comm_bcast( offline_parent_okmax )
457  call comm_bcast( offline_parent_lkmax )
458  endif
459 
460  if ( online_iam_daughter .or. online_iam_parent ) then
461 
462  if ( offline ) then
463  log_error("COMM_CARTESC_NEST_setup",*) 'OFFLINE and ONLINE cannot be use at the same time'
464  call prc_abort
465  endif
466 
467  use_nesting = .true.
468  endif
469 
470  call interp_setup( comm_cartesc_nest_interp_weight_order ) ! [IN]
471 
472  itp_nh = comm_cartesc_nest_interp_level
473  itp_nv = 2
474 
475  debug_domain_num = online_domain_num
476  if( online_specified_maxrq > max_rq ) max_rq = online_specified_maxrq
477 
478  allocate( ireq_p(max_rq) )
479  allocate( ireq_d(max_rq) )
480  allocate( call_order(max_rq) )
481  ireq_p(:) = mpi_request_null
482  ireq_d(:) = mpi_request_null
483 
484  if ( use_nesting ) then
485 
486  if ( offline .OR. online_iam_daughter ) then
487  latlon_local(i_min,i_lon) = minval(atmos_grid_cartesc_real_lonuv(:,:)) / d2r
488  latlon_local(i_max,i_lon) = maxval(atmos_grid_cartesc_real_lonuv(:,:)) / d2r
489  latlon_local(i_min,i_lat) = minval(atmos_grid_cartesc_real_latuv(:,:)) / d2r
490  latlon_local(i_max,i_lat) = maxval(atmos_grid_cartesc_real_latuv(:,:)) / d2r
491  endif
492 
493  if ( offline ) then
494 
495  handling_num = 1
496  parent_prc_num_x(handling_num) = offline_parent_prc_num_x
497  parent_prc_num_y(handling_num) = offline_parent_prc_num_y
498  parent_kmax(handling_num) = offline_parent_kmax
499  parent_imax(handling_num) = offline_parent_imax
500  parent_jmax(handling_num) = offline_parent_jmax
501  parent_okmax(handling_num) = offline_parent_okmax
502  parent_lkmax(handling_num) = offline_parent_lkmax
503 
504  parent_prc_nprocs(handling_num) = parent_prc_num_x(handling_num) * parent_prc_num_y(handling_num)
505  allocate( latlon_catalog(parent_prc_nprocs(handling_num),2,2) )
506 
507  !--- read latlon catalogue
508  fid = io_get_available_fid()
509  open( fid, &
510  file = trim(latlon_catalogue_fname), &
511  form = 'formatted', &
512  status = 'old', &
513  iostat = ierr )
514 
515  if ( ierr /= 0 ) then
516  log_error("COMM_CARTESC_NEST_setup",*) 'cannot open latlon-catalogue file!'
517  call prc_abort
518  endif
519 
520  do i = 1, parent_prc_nprocs(handling_num)
521  read(fid,'(i8,4f32.24)',iostat=ierr) parent_id, &
522  latlon_catalog(i,i_min,i_lon), latlon_catalog(i,i_max,i_lon), & ! LON: MIN, MAX
523  latlon_catalog(i,i_min,i_lat), latlon_catalog(i,i_max,i_lat) ! LAT: MIN, MAX
524  if ( i /= parent_id ) then
525  log_error("COMM_CARTESC_NEST_setup",*) 'internal error: parent mpi id'
526  call prc_abort
527  endif
528  if ( ierr /= 0 ) exit
529  enddo
530  close(fid)
531 
532  call comm_cartesc_nest_domain_relate(handling_num)
533 
534  else ! ONLINE RELATIONSHIP
535 ! if ( present(flag_parent) .AND. present(flag_child) ) then
536 ! LOG_INFO("COMM_CARTESC_NEST_setup",'(1x,A)') &
537 ! '*** Setup Online Nesting Inter-Domain Communicator (IDC)'
538 ! else
539 ! LOG_ERROR("COMM_CARTESC_NEST_setup",*) 'Internal Error:'
540 ! LOG_ERROR_CONT(*) 'The flag_parent and flag_child are needed.'
541 ! LOG_WARN("COMM_CARTESC_NEST_setup",*) ' domain: ', ONLINE_DOMAIN_NUM
542 ! call PRC_abort
543 ! endif
544 
545  if( online_boundary_use_qhyd ) then
546  mp_type = mp_type_in
547  comm_cartesc_nest_bnd_qa = qa_mp
548  elseif ( atmos_hydrometeor_dry ) then
549  mp_type = "DRY"
550  comm_cartesc_nest_bnd_qa = 0
551  else
552  mp_type = "QV"
553  comm_cartesc_nest_bnd_qa = 1
554  endif
555 
556  log_info("COMM_CARTESC_NEST_setup",*) "flag_parent", flag_parent, "flag_child", flag_child
557  log_info("COMM_CARTESC_NEST_setup",*) "ONLINE_IAM_PARENT", online_iam_parent, "ONLINE_IAM_DAUGHTER", online_iam_daughter
558 
559  if( flag_parent ) then ! must do first before daughter processes
560  !-------------------------------------------------
561  if ( .NOT. online_iam_parent ) then
562  log_error("COMM_CARTESC_NEST_setup",*) '[NEST_setup] Parent Flag from launcher is not consistent with namelist!'
563  log_error_cont(*) 'PARENT - domain : ', online_domain_num
564  call prc_abort
565  endif
566 
567  handling_num = 1 !HANDLING_NUM + 1
568  intercomm_id(handling_num) = online_domain_num
569  comm_cartesc_nest_filiation(intercomm_id(handling_num)) = 1
570 
571  intercomm_daughter = inter_child
572  log_info("COMM_CARTESC_NEST_setup",'(1x,A,I2,A)') 'Online Nesting - PARENT [INTERCOMM_ID:', &
573  intercomm_id(handling_num), ' ]'
574  log_info("COMM_CARTESC_NEST_setup",*) 'Online Nesting - INTERCOMM :', intercomm_daughter
575 
576  call comm_cartesc_nest_ping( handling_num )
577 
578  call comm_cartesc_nest_parentsize( handling_num )
579 
580  call comm_cartesc_nest_catalogue( handling_num )
581  call mpi_barrier(intercomm_daughter, ierr)
582 
583  parent_ka(handling_num) = parent_kmax(handling_num) + khalo * 2
584  parent_ia(handling_num) = parent_imax(handling_num) + ihalo * 2
585  parent_ja(handling_num) = parent_jmax(handling_num) + jhalo * 2
586  daughter_ka(handling_num) = daughter_kmax(handling_num) + khalo * 2
587  daughter_ia(handling_num) = daughter_imax(handling_num) + ihalo * 2
588  daughter_ja(handling_num) = daughter_jmax(handling_num) + jhalo * 2
589  tileal_ka(handling_num) = 0
590  tileal_ia(handling_num) = 0
591  tileal_ja(handling_num) = 0
592 
593  log_info("COMM_CARTESC_NEST_setup",'(1x,A)' ) 'Informations of Parent Domain [me]'
594  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_nprocs :', parent_prc_nprocs(handling_num)
595  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_NUM_X :', parent_prc_num_x(handling_num)
596  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_NUM_Y :', parent_prc_num_y(handling_num)
597  log_info_cont('(1x,A,I6)' ) '--- PARENT_KMAX :', parent_kmax(handling_num)
598  log_info_cont('(1x,A,I6)' ) '--- PARENT_IMAX :', parent_imax(handling_num)
599  log_info_cont('(1x,A,I6)' ) '--- PARENT_JMAX :', parent_jmax(handling_num)
600  log_info_cont('(1x,A,F9.3)') '--- PARENT_DTSEC :', parent_dtsec(handling_num)
601  log_info_cont('(1x,A,I6) ') '--- PARENT_NSTEP :', parent_nstep(handling_num)
602  log_info_cont('(1x,A)' ) 'Informations of Daughter Domain'
603  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_PRC_nprocs :', daughter_prc_nprocs(handling_num)
604  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_PRC_NUM_X :', daughter_prc_num_x(handling_num)
605  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_PRC_NUM_Y :', daughter_prc_num_y(handling_num)
606  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_KMAX :', daughter_kmax(handling_num)
607  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_IMAX :', daughter_imax(handling_num)
608  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_JMAX :', daughter_jmax(handling_num)
609  log_info_cont('(1x,A,F9.3)') '--- DAUGHTER_DTSEC :', daughter_dtsec(handling_num)
610  log_info_cont('(1x,A,I6) ') '--- DAUGHTER_NSTEP :', daughter_nstep(handling_num)
611  log_info_cont('(1x,A,I6) ') 'Limit Num. NCOMM req. :', max_rq
612 
613  allocate( org_dens(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
614  allocate( org_momz(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
615  allocate( org_momx(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
616  allocate( org_momy(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
617  allocate( org_u_ll(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
618  allocate( org_v_ll(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
619  allocate( org_rhot(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num)) )
620  allocate( org_qtrc(parent_ka(handling_num),parent_ia(handling_num),parent_ja(handling_num),max(comm_cartesc_nest_bnd_qa,1)) )
621 
622  call comm_cartesc_nest_setup_nestdown( handling_num )
623 
624  !---------------------------------- end of parent routines
625  endif
626 
627 
628  if( flag_child ) then
629  !-------------------------------------------------
630  if ( .NOT. online_iam_daughter ) then
631  log_error("COMM_CARTESC_NEST_setup",*) '[NEST_setup] Child Flag from launcher is not consistent with namelist!'
632  log_error_cont(*) 'DAUGHTER - domain : ', online_domain_num
633  call prc_abort
634  endif
635 
636  handling_num = 2 !HANDLING_NUM + 1
637  intercomm_id(handling_num) = online_domain_num - 1
638  comm_cartesc_nest_filiation(intercomm_id(handling_num)) = -1
639 
640  intercomm_parent = inter_parent
641  log_info("COMM_CARTESC_NEST_setup",'(1x,A,I2,A)') 'Online Nesting - DAUGHTER [INTERCOMM_ID:', &
642  intercomm_id(handling_num), ' ]'
643  log_info("COMM_CARTESC_NEST_setup",*) 'Online Nesting - INTERCOMM :', intercomm_parent
644 
645  call comm_cartesc_nest_ping( handling_num )
646 
647  call comm_cartesc_nest_parentsize( handling_num )
648 
649  allocate( latlon_catalog(parent_prc_nprocs(handling_num),2,2) )
650  call comm_cartesc_nest_catalogue( handling_num )
651  call mpi_barrier(intercomm_parent, ierr)
652 
653  call comm_cartesc_nest_domain_relate( handling_num )
654 
655  parent_ka(handling_num) = parent_kmax(handling_num) + khalo * 2
656  parent_ia(handling_num) = parent_imax(handling_num) + ihalo * 2
657  parent_ja(handling_num) = parent_jmax(handling_num) + jhalo * 2
658  daughter_ka(handling_num) = daughter_kmax(handling_num) + khalo * 2
659  daughter_ia(handling_num) = daughter_imax(handling_num) + ihalo * 2
660  daughter_ja(handling_num) = daughter_jmax(handling_num) + jhalo * 2
661  tileal_ka(handling_num) = parent_ka(handling_num)
662  tileal_ia(handling_num) = parent_imax(handling_num) * comm_cartesc_nest_tile_num_x
663  tileal_ja(handling_num) = parent_jmax(handling_num) * comm_cartesc_nest_tile_num_y
664 
665  log_info("COMM_CARTESC_NEST_setup",'(1x,A)' ) 'Informations of Parent Domain'
666  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_nprocs :', parent_prc_nprocs(handling_num)
667  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_NUM_X :', parent_prc_num_x(handling_num)
668  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_NUM_Y :', parent_prc_num_y(handling_num)
669  log_info_cont('(1x,A,I6)' ) '--- PARENT_KMAX :', parent_kmax(handling_num)
670  log_info_cont('(1x,A,I6)' ) '--- PARENT_IMAX :', parent_imax(handling_num)
671  log_info_cont('(1x,A,I6)' ) '--- PARENT_JMAX :', parent_jmax(handling_num)
672  log_info_cont('(1x,A,F9.3)') '--- PARENT_DTSEC :', parent_dtsec(handling_num)
673  log_info_cont('(1x,A,I6)' ) '--- PARENT_NSTEP :', parent_nstep(handling_num)
674  log_info_cont('(1x,A)' ) 'Informations of Daughter Domain [me]'
675  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_PRC_nprocs :', daughter_prc_nprocs(handling_num)
676  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_PRC_NUM_X :', daughter_prc_num_x(handling_num)
677  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_PRC_NUM_Y :', daughter_prc_num_y(handling_num)
678  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_KMAX :', daughter_kmax(handling_num)
679  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_IMAX :', daughter_imax(handling_num)
680  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_JMAX :', daughter_jmax(handling_num)
681  log_info_cont('(1x,A,F9.3)') '--- DAUGHTER_DTSEC :', daughter_dtsec(handling_num)
682  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_NSTEP :', daughter_nstep(handling_num)
683  log_info_cont('(1x,A)' ) 'Informations of Target Tiles'
684  log_info_cont('(1x,A,I6)' ) '--- TILEALL_KA :', tileal_ka(handling_num)
685  log_info_cont('(1x,A,I6)' ) '--- TILEALL_IA :', tileal_ia(handling_num)
686  log_info_cont('(1x,A,I6)' ) '--- TILEALL_JA :', tileal_ja(handling_num)
687  log_info_cont('(1x,A,I6) ') 'Limit Num. NCOMM req. :', max_rq
688 
689  allocate( buffer_2d( parent_ia(handling_num), parent_ja(handling_num) ) )
690  allocate( buffer_3d( parent_ka(handling_num), parent_ia(handling_num), parent_ja(handling_num) ) )
691  allocate( buffer_3df( 0:parent_ka(handling_num), parent_ia(handling_num), parent_ja(handling_num) ) )
692 
693  allocate( recvbuf_3d( parent_ka(handling_num), parent_ia(handling_num), parent_ja(handling_num), max_isu ) )
694  allocate( recvbuf_3df( 0:parent_ka(handling_num), parent_ia(handling_num), parent_ja(handling_num), max_isuf ) )
695 
696  allocate( buffer_ref_lon( tileal_ia(handling_num),tileal_ja(handling_num)) )
697  allocate( buffer_ref_lonuy( tileal_ia(handling_num),tileal_ja(handling_num)) )
698  allocate( buffer_ref_lonxv( tileal_ia(handling_num),tileal_ja(handling_num)) )
699  allocate( buffer_ref_lat( tileal_ia(handling_num),tileal_ja(handling_num)) )
700  allocate( buffer_ref_latuy( tileal_ia(handling_num),tileal_ja(handling_num)) )
701  allocate( buffer_ref_latxv( tileal_ia(handling_num),tileal_ja(handling_num)) )
702 
703  allocate( buffer_ref_cz( tileal_ka(handling_num),tileal_ia(handling_num),tileal_ja(handling_num)) )
704  allocate( buffer_ref_fz(0:tileal_ka(handling_num),tileal_ia(handling_num),tileal_ja(handling_num)) )
705 
706  allocate( buffer_ref_3d( tileal_ka(handling_num),tileal_ia(handling_num),tileal_ja(handling_num)) )
707  allocate( buffer_ref_3df(0:tileal_ka(handling_num),tileal_ia(handling_num),tileal_ja(handling_num)) )
708 
709  allocate( igrd( daughter_ia(handling_num),daughter_ja(handling_num),itp_nh,itp_ng) )
710  allocate( jgrd( daughter_ia(handling_num),daughter_ja(handling_num),itp_nh,itp_ng) )
711  allocate( hfact( daughter_ia(handling_num),daughter_ja(handling_num),itp_nh,itp_ng) )
712  allocate( kgrd(daughter_ka(handling_num),itp_nv,daughter_ia(handling_num),daughter_ja(handling_num),itp_nh,itp_ng) )
713  allocate( vfact(daughter_ka(handling_num),itp_nv,daughter_ia(handling_num),daughter_ja(handling_num),itp_nh,itp_ng) )
714 
715  call comm_cartesc_nest_setup_nestdown( handling_num )
716 
717 
718  ! for scalar points
719  call interp_factor3d( itp_nh, & ! [IN]
720  tileal_ka(handling_num), & ! [IN]
721  khalo+1, & ! [IN]
722  tileal_ka(handling_num)-khalo, & ! [IN]
723  tileal_ia(handling_num), & ! [IN]
724  tileal_ja(handling_num), & ! [IN]
725  buffer_ref_lon(:,:), & ! [IN]
726  buffer_ref_lat(:,:), & ! [IN]
727  buffer_ref_cz(:,:,:), & ! [IN]
728  daughter_ka(handling_num), & ! [IN]
729  datr_ks(handling_num), & ! [IN]
730  datr_ke(handling_num), & ! [IN]
731  daughter_ia(handling_num), & ! [IN]
732  daughter_ja(handling_num), & ! [IN]
733  atmos_grid_cartesc_real_lon(:,:), & ! [IN]
734  atmos_grid_cartesc_real_lat(:,:), & ! [IN]
735  atmos_grid_cartesc_real_cz(:,:,:), & ! [IN]
736  igrd( :,:,:,i_sclr), & ! [OUT]
737  jgrd( :,:,:,i_sclr), & ! [OUT]
738  hfact( :,:,:,i_sclr), & ! [OUT]
739  kgrd(:,:,:,:,:,i_sclr), & ! [OUT]
740  vfact(:,:,:,:,:,i_sclr) ) ! [OUT]
741 
742  ! for z staggered points
743  call interp_factor3d( itp_nh, & ! [IN]
744  tileal_ka(handling_num)+1, & ! [IN]
745  khalo+1, & ! [IN]
746  tileal_ka(handling_num)+1-khalo, & ! [IN]
747  tileal_ia(handling_num), & ! [IN]
748  tileal_ja(handling_num), & ! [IN]
749  buffer_ref_lon(:,:), & ! [IN]
750  buffer_ref_lat(:,:), & ! [IN]
751  buffer_ref_fz(:,:,:), & ! [IN]
752  daughter_ka(handling_num), & ! [IN]
753  datr_ks(handling_num), & ! [IN]
754  datr_ke(handling_num), & ! [IN]
755  daughter_ia(handling_num), & ! [IN]
756  daughter_ja(handling_num), & ! [IN]
757  atmos_grid_cartesc_real_lon(:,:), & ! [IN]
758  atmos_grid_cartesc_real_lat(:,:), & ! [IN]
759  atmos_grid_cartesc_real_fz(1:ka,:,:), & ! [IN]
760  igrd( :,:,:,i_zstg), & ! [OUT]
761  jgrd( :,:,:,i_zstg), & ! [OUT]
762  hfact( :,:,:,i_zstg), & ! [OUT]
763  kgrd(:,:,:,:,:,i_zstg), & ! [OUT]
764  vfact(:,:,:,:,:,i_zstg) ) ! [OUT]
765 
766  ! for x staggered points
767  call interp_factor3d( itp_nh, & ! [IN]
768  tileal_ka(handling_num), & ! [IN]
769  khalo+1, & ! [IN]
770  tileal_ka(handling_num)-khalo, & ! [IN]
771  tileal_ia(handling_num), & ! [IN]
772  tileal_ja(handling_num), & ! [IN]
773  buffer_ref_lonuy(:,:), & ! [IN]
774  buffer_ref_latxv(:,:), & ! [IN]
775  buffer_ref_cz(:,:,:), & ! [IN]
776  daughter_ka(handling_num), & ! [IN]
777  datr_ks(handling_num), & ! [IN]
778  datr_ke(handling_num), & ! [IN]
779  daughter_ia(handling_num), & ! [IN]
780  daughter_ja(handling_num), & ! [IN]
781  atmos_grid_cartesc_real_lonuy(1:ia,1:ja), & ! [IN]
782  atmos_grid_cartesc_real_latuy(1:ia,1:ja), & ! [IN]
783  atmos_grid_cartesc_real_cz(:,:,:), & ! [IN]
784  igrd( :,:,:,i_xstg), & ! [OUT]
785  jgrd( :,:,:,i_xstg), & ! [OUT]
786  hfact( :,:,:,i_xstg), & ! [OUT]
787  kgrd(:,:,:,:,:,i_xstg), & ! [OUT]
788  vfact(:,:,:,:,:,i_xstg) ) ! [OUT]
789 
790  ! for y staggered points
791  call interp_factor3d( itp_nh, & ! [IN]
792  tileal_ka(handling_num), & ! [IN]
793  khalo+1, & ! [IN]
794  tileal_ka(handling_num)-khalo, & ! [IN]
795  tileal_ia(handling_num), & ! [IN]
796  tileal_ja(handling_num), & ! [IN]
797  buffer_ref_lonxv(:,:), & ! [IN]
798  buffer_ref_latxv(:,:), & ! [IN]
799  buffer_ref_cz(:,:,:), & ! [IN]
800  daughter_ka(handling_num), & ! [IN]
801  datr_ks(handling_num), & ! [IN]
802  datr_ke(handling_num), & ! [IN]
803  daughter_ia(handling_num), & ! [IN]
804  daughter_ja(handling_num), & ! [IN]
805  atmos_grid_cartesc_real_lonxv(1:ia,1:ja), & ! [IN]
806  atmos_grid_cartesc_real_latxv(1:ia,1:ja), & ! [IN]
807  atmos_grid_cartesc_real_cz(:,:,:), & ! [IN]
808  igrd( :,:,:,i_ystg), & ! [OUT]
809  jgrd( :,:,:,i_ystg), & ! [OUT]
810  hfact( :,:,:,i_ystg), & ! [OUT]
811  kgrd(:,:,:,:,:,i_ystg), & ! [OUT]
812  vfact(:,:,:,:,:,i_ystg) ) ! [OUT]
813 
814  deallocate( buffer_2d )
815  deallocate( buffer_3d )
816  deallocate( buffer_3df )
817 
818  else
819  online_use_velz = .false.
820  endif
821 
822  !LOG_INFO("COMM_CARTESC_NEST_setup",'(1x,A,I2)') 'Number of Related Domains :', HANDLING_NUM
823  !if ( HANDLING_NUM > 2 ) then
824  ! f( IO_L ) LOG_ERROR("COMM_CARTESC_NEST_setup",*) 'Too much handing domains (up to 2)'
825  ! call PRC_abort
826  !endif
827 
828  endif !--- OFFLINE or NOT
829 
830  endif !--- USE_NESTING
831 
832  return
integer, public prc_global_domainid
my domain ID in global communicator
Definition: scale_prc.F90:83
integer, parameter, public khalo
of halo cells: z
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lonuv
longitude at staggered point (uv) [rad,0-2pi]
integer, public ia
of whole cells: x, local, with HALO
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lat
latitude [rad,-pi,pi]
module INTERPOLATION
subroutine, public interp_factor3d(npoints, KA_ref, KS_ref, KE_ref, IA_ref, JA_ref, lon_ref, lat_ref, hgt_ref, KA, KS, KE, IA, JA, lon, lat, hgt, idx_i, idx_j, hfact, idx_k, vfact)
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_fz
geopotential height [m] (wxy)
integer, public ja
of whole cells: y, local, with HALO
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lon
longitude [rad,0-2pi]
subroutine, public file_open(basename, fid, mode, single, aggregate, rankid, postfix)
Definition: scale_file.F90:464
integer, public jmaxg
of computational cells: y, global
module COMMUNICATION
module file
Definition: scale_file.F90:15
module atmosphere / hydrometeor
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latuy
latitude at staggered point (uy) [rad,-pi,pi]
integer, parameter, public i_min
[index] minute
integer function, public io_get_available_fid()
search & get available file ID
Definition: scale_io.F90:313
module PROCESS
Definition: scale_prc.F90:11
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lonxv
longitude at staggered point (xv) [rad,0-2pi]
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lonuy
longitude at staggered point (uy) [rad,0-2pi]
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module CONSTANT
Definition: scale_const.F90:11
subroutine, public interp_setup(weight_order, search_limit)
Setup.
integer, public imaxg
of computational cells: x, global
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_cz
geopotential height [m] (zxy)
logical, public prc_ismaster
master process in local communicator?
Definition: scale_prc.F90:90
module Atmosphere GRID CartesC Real(real space)
integer, public ka
of whole cells: z, local, with HALO
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latuv
latitude at staggered point (uv) [rad,-pi,pi]
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latxv
latitude at staggered point (xv) [rad,-pi,pi]
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 839 of file scale_comm_cartesC_nest.F90.

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

839  use scale_prc, only: &
840  prc_myrank, &
841  prc_abort
842  implicit none
843 
844  integer, intent(in) :: handle
845 
846  integer :: x_min, x_max
847  integer :: y_min, y_max
848  logical :: hit(2,2)
849  integer :: p, i, j
850  !---------------------------------------------------------------------------
851 
852  if( .NOT. use_nesting ) return
853 
854  x_min = parent_prc_num_x(handle)
855  x_max = -1
856  y_min = parent_prc_num_y(handle)
857  y_max = -1
858  hit(:,:) = .false.
859  do p = 1, parent_prc_nprocs(handle)
860  if ( ( ( latlon_local(i_min,i_lon) >= latlon_catalog(p,i_min,i_lon) &
861  .AND. latlon_local(i_min,i_lon) <= latlon_catalog(p,i_max,i_lon) ) .OR. &
862  ( latlon_local(i_max,i_lon) >= latlon_catalog(p,i_min,i_lon) &
863  .AND. latlon_local(i_max,i_lon) <= latlon_catalog(p,i_max,i_lon) ) .OR. &
864  ( latlon_catalog(p,i_min,i_lon) >= latlon_local(i_min,i_lon) &
865  .AND. latlon_catalog(p,i_min,i_lon) <= latlon_local(i_max,i_lon) ) .OR. &
866  ( latlon_catalog(p,i_max,i_lon) >= latlon_local(i_min,i_lon) &
867  .AND. latlon_catalog(p,i_max,i_lon) <= latlon_local(i_max,i_lon) ) ) .AND. &
868  ( ( latlon_local(i_min,i_lat) >= latlon_catalog(p,i_min,i_lat) &
869  .AND. latlon_local(i_min,i_lat) <= latlon_catalog(p,i_max,i_lat) ) .OR. &
870  ( latlon_local(i_max,i_lat) >= latlon_catalog(p,i_min,i_lat) &
871  .AND. latlon_local(i_max,i_lat) <= latlon_catalog(p,i_max,i_lat) ) .OR. &
872  ( latlon_catalog(p,i_min,i_lat) >= latlon_local(i_min,i_lat) &
873  .AND. latlon_catalog(p,i_min,i_lat) <= latlon_local(i_max,i_lat) ) .OR. &
874  ( latlon_catalog(p,i_max,i_lat) >= latlon_local(i_min,i_lat) &
875  .AND. latlon_catalog(p,i_max,i_lat) <= latlon_local(i_max,i_lat) ) ) ) then
876  if ( latlon_catalog(p,i_min,i_lon) <= latlon_local(i_min,i_lon) ) hit(i_min,i_lon) = .true.
877  if ( latlon_catalog(p,i_max,i_lon) >= latlon_local(i_max,i_lon) ) hit(i_max,i_lon) = .true.
878  if ( latlon_catalog(p,i_min,i_lat) <= latlon_local(i_min,i_lat) ) hit(i_min,i_lat) = .true.
879  if ( latlon_catalog(p,i_max,i_lat) >= latlon_local(i_max,i_lat) ) hit(i_max,i_lat) = .true.
880  i = mod(p-1, parent_prc_num_x(handle))
881  j = (p-1) / parent_prc_num_x(handle)
882  if ( i < x_min ) x_min = i
883  if ( i > x_max ) x_max = i
884  if ( j < y_min ) y_min = j
885  if ( j > y_max ) y_max = j
886  end if
887  end do
888 
889  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
890  log_error("COMM_CARTESC_NEST_domain_relate",*) 'region of daughter domain is larger than that of parent'
891  log_error_cont(*) ' at rank:', prc_myrank, ' of domain:', online_domain_num
892  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)
893  log_error_cont('(A,F12.6,1x,F12.6)') 'daughter local (me) MIN-MAX: LON=', &
894  latlon_local(i_min,i_lon), latlon_local(i_max,i_lon)
895  do p = 1, parent_prc_nprocs(handle)
896  log_error_cont('(A,I5,A,F12.6,1x,F12.6)') ' parent (', p,') MIN-MAX: LON=', &
897  latlon_catalog(p,i_min,i_lon) ,latlon_catalog(p,i_max,i_lon)
898  enddo
899  log_error_cont('(A,F12.6,1x,F12.6)') 'daughter local (me): MIN-MAX LAT=', &
900  latlon_local(i_min,i_lat), latlon_local(i_max,i_lat)
901  do p = 1, parent_prc_nprocs(handle)
902  log_error_cont('(A,I5,A,F12.6,1x,F12.6)') ' parent (', p,') MIN-MAX: LAT=', &
903  latlon_catalog(p,i_min,i_lat) ,latlon_catalog(p,i_max,i_lat)
904  enddo
905  call prc_abort
906  end if
907 
908 
909 
910  comm_cartesc_nest_tile_num_x = x_max - x_min + 1
911  comm_cartesc_nest_tile_num_y = y_max - y_min + 1
912 
913  allocate( comm_cartesc_nest_tile_id( comm_cartesc_nest_tile_num_x * comm_cartesc_nest_tile_num_y ) )
914 
915  log_info("COMM_CARTESC_NEST_domain_relate",'(1x,A)') 'NEST: target process tile in parent domain'
916  p = 1
917  do j = 1, comm_cartesc_nest_tile_num_y
918  do i = 1, comm_cartesc_nest_tile_num_x
919  comm_cartesc_nest_tile_id(p) = x_min + i - 1 + (y_min + j - 1) * parent_prc_num_x(handle)
920  log_info_cont('(1x,A,I4,A,I6)') '(', p, ') target mpi-process:', comm_cartesc_nest_tile_id(p)
921  p = p + 1
922  enddo
923  enddo
924 
925  return
integer, parameter, public i_min
[index] minute
module PROCESS
Definition: scale_prc.F90:11
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:89
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
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 939 of file scale_comm_cartesC_nest.F90.

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

939  implicit none
940 
941  integer, intent(out) :: tilei, tilej
942  integer, intent(out) :: cxs, cxe, cys, cye
943  integer, intent(out) :: pxs, pxe, pys, pye
944  integer, intent(in) :: iloc ! rank number; start from 1
945 
946  integer :: hdl = 1 ! handler number
947  integer :: rank
948  integer :: xloc, yloc
949  integer :: xlocg, ylocg ! location over whole domain
950  !---------------------------------------------------------------------------
951 
952  if( .NOT. use_nesting ) return
953 
954  rank = comm_cartesc_nest_tile_id(iloc)
955  xloc = mod( iloc-1, comm_cartesc_nest_tile_num_x ) + 1
956  yloc = int( real(iloc-1) / real(COMM_CARTESC_NEST_TILE_NUM_X) ) + 1
957  xlocg = mod( rank, offline_parent_prc_num_x ) + 1
958  ylocg = int( real(rank) / real(OFFLINE_PARENT_PRC_NUM_X) ) + 1
959  tilei = parent_imax(hdl)
960  tilej = parent_jmax(hdl)
961 
962  cxs = tilei * (xloc-1) + 1
963  cxe = tilei * xloc
964  cys = tilej * (yloc-1) + 1
965  cye = tilej * yloc
966  pxs = 1
967  pxe = tilei
968  pys = 1
969  pye = tilej
970 
971  if ( xlocg == 1 ) then ! BND_W
972  tilei = tilei + 2
973  pxs = pxs + 2
974  pxe = pxe + 2
975  endif
976  if ( xlocg == offline_parent_prc_num_x ) then ! BND_E
977  tilei = tilei + 2
978  endif
979  if ( ylocg == 1 ) then ! BND_S
980  tilej = tilej + 2
981  pys = pys + 2
982  pye = pye + 2
983  endif
984  if ( ylocg == offline_parent_prc_num_y ) then ! BND_N
985  tilej = tilej + 2
986  endif
987 
988  return
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 1728 of file scale_comm_cartesC_nest.F90.

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

1728  use scale_prc, only: &
1729  prc_abort
1730  use scale_comm_cartesc, only: &
1731  comm_vars8, &
1732  comm_wait
1733  use scale_atmos_grid_cartesc_metric, only: &
1735  implicit none
1736 
1737  integer, intent(in) :: handle
1738  integer, intent(in) :: bnd_qa
1739  real(RP), intent(in) :: dens_send(parent_ka (handle),parent_ia (handle),parent_ja (handle))
1740  real(RP), intent(in) :: momz_send(parent_ka (handle),parent_ia (handle),parent_ja (handle))
1741  real(RP), intent(in) :: momx_send(parent_ka (handle),parent_ia (handle),parent_ja (handle))
1742  real(RP), intent(in) :: momy_send(parent_ka (handle),parent_ia (handle),parent_ja (handle))
1743  real(RP), intent(in) :: rhot_send(parent_ka (handle),parent_ia (handle),parent_ja (handle))
1744  real(RP), intent(in) :: qtrc_send(parent_ka (handle),parent_ia (handle),parent_ja (handle),bnd_qa)
1745  real(RP), intent(out) :: dens_recv(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1746  real(RP), intent(out) :: velz_recv(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1747  real(RP), intent(out) :: velx_recv(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1748  real(RP), intent(out) :: vely_recv(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1749  real(RP), intent(out) :: pott_recv(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1750  real(RP), intent(out) :: qtrc_recv(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle),bnd_qa)
1751 
1752  real(RP) :: work1_send(parent_ka (handle),parent_ia (handle),parent_ja (handle))
1753  real(RP) :: work2_send(parent_ka (handle),parent_ia (handle),parent_ja (handle))
1754  real(RP) :: work1_recv(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1755  real(RP) :: work2_recv(daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1756  real(RP) :: u_ll_recv (daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1757  real(RP) :: v_ll_recv (daughter_ka(handle),daughter_ia(handle),daughter_ja(handle))
1758  real(RP) :: u_on_map, v_on_map
1759 
1760  real(RP) :: dummy(1,1,1)
1761  integer :: tagbase, tagcomm
1762  integer :: isu_tag, isu_tagf
1763 
1764  integer :: ierr
1765  integer :: i, j, k, iq
1766  !---------------------------------------------------------------------------
1767 
1768  if( .NOT. use_nesting ) return
1769 
1770  if ( bnd_qa > i_bndqa ) then
1771  log_error("COMM_CARTESC_NEST_nestdown",*) 'internal error: BND_QA is larger than I_BNDQA'
1772  call prc_abort
1773  endif
1774 
1775  tagcomm = intercomm_id(handle) * order_tag_comm
1776 
1777  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
1778 
1779  !##### parent [send issue] #####
1780 
1781  call prof_rapstart('NEST_total_P', 2)
1782  call prof_rapstart('NEST_pack_P', 2)
1783 
1784  nsend = nsend + 1
1785  log_info("COMM_CARTESC_NEST_nestdown",'(1X,A,I5,A)') "CONeP[P] send( ", nsend, " )"
1786 
1787  ! to keep values at that time by finish of sending process
1788 !OCL XFILL
1789  org_dens(:,:,:) = dens_send(:,:,:)
1790 !OCL XFILL
1791  org_momz(:,:,:) = momz_send(:,:,:)
1792 !OCL XFILL
1793  org_momx(:,:,:) = momx_send(:,:,:)
1794 !OCL XFILL
1795  org_momy(:,:,:) = momy_send(:,:,:)
1796 !OCL XFILL
1797  org_rhot(:,:,:) = rhot_send(:,:,:)
1798  do iq = 1, bnd_qa
1799 !OCL XFILL
1800  org_qtrc(:,:,:,iq) = qtrc_send(:,:,:,iq)
1801  enddo
1802 
1803  !*** request control
1804  !--- do not change the calling order below;
1805  !--- it should be consistent with the order in "COMM_CARTESC_NEST_recvwait_issue"
1806  rq_ctl_p = 0
1807 
1808  if ( .NOT. online_daughter_no_rotate ) then
1809  ! from staggered point to scalar point
1810  do j = 1, parent_ja(handle)
1811  do i = 2, parent_ia(handle)
1812  do k = 1, parent_ka(handle)
1813  work1_send(k,i,j) = ( org_momx(k,i-1,j) + org_momx(k,i,j) ) * 0.5_rp
1814  enddo
1815  enddo
1816  enddo
1817 
1818  do j = 1, parent_ja(handle)
1819  do k = 1, parent_ka(handle)
1820  work1_send(k,1,j) = org_momx(k,1,j)
1821  enddo
1822  enddo
1823 
1824  call comm_vars8( work1_send(:,:,:), 1 )
1825 
1826  do j = 2, parent_ja(handle)
1827  do i = 1, parent_ia(handle)
1828  do k = 1, parent_ka(handle)
1829  work2_send(k,i,j) = ( org_momy(k,i,j-1) + org_momy(k,i,j) ) * 0.5_rp
1830  enddo
1831  enddo
1832  enddo
1833 
1834  do i = 1, parent_ia(handle)
1835  do k = 1, parent_ka(handle)
1836  work2_send(k,i,1) = org_momy(k,i,1)
1837  enddo
1838  enddo
1839 
1840  call comm_vars8( work2_send(:,:,:), 2 )
1841 
1842  call comm_wait ( work1_send(:,:,:), 1, .false. )
1843  call comm_wait ( work2_send(:,:,:), 2, .false. )
1844 
1845  ! rotation from map-projected field to latlon field
1846  do j = 1, parent_ja(handle)
1847  do i = 1, parent_ia(handle)
1848  do k = 1, parent_ka(handle)
1849  u_on_map = work1_send(k,i,j) / org_dens(k,i,j)
1850  v_on_map = work2_send(k,i,j) / org_dens(k,i,j)
1851 
1852  org_u_ll(k,i,j) = u_on_map * rotc(i,j,1) - v_on_map * rotc(i,j,2)
1853  org_v_ll(k,i,j) = u_on_map * rotc(i,j,2) + v_on_map * rotc(i,j,1)
1854  enddo
1855  enddo
1856  enddo
1857  endif
1858 
1859  tagbase = tagcomm + tag_dens*order_tag_var
1860  call comm_cartesc_nest_intercomm_nestdown( org_dens(:,:,:), & ! [IN]
1861  dummy(:,:,:), & ! [OUT]
1862  tagbase, i_sclr, handle, & ! [IN]
1863  isu_tag, isu_tagf, & ! [INOUT]
1864  flag_dens = .true. ) ! [IN]
1865 
1866  tagbase = tagcomm + tag_momz*order_tag_var
1867  if ( online_daughter_use_velz ) then
1868  call comm_cartesc_nest_intercomm_nestdown( org_momz(:,:,:), & ! [IN]
1869  dummy(:,:,:), & ! [OUT]
1870  tagbase, i_zstg, handle, & ! [IN]
1871  isu_tag, isu_tagf ) ! [INOUT]
1872  endif
1873 
1874  tagbase = tagcomm + tag_momx*order_tag_var
1875  if ( online_daughter_no_rotate ) then
1876  call comm_cartesc_nest_intercomm_nestdown( org_momx(:,:,:), & ! [IN]
1877  dummy(:,:,:), & ! [OUT]
1878  tagbase, i_xstg, handle, & ! [IN]
1879  isu_tag, isu_tagf ) ! [INOUT]
1880  else
1881  call comm_cartesc_nest_intercomm_nestdown( org_u_ll(:,:,:), & ! [IN]
1882  dummy(:,:,:), & ! [OUT]
1883  tagbase, i_sclr, handle, & ! [IN]
1884  isu_tag, isu_tagf ) ! [INOUT]
1885  endif
1886 
1887  tagbase = tagcomm + tag_momy*order_tag_var
1888  if ( online_daughter_no_rotate ) then
1889  call comm_cartesc_nest_intercomm_nestdown( org_momy(:,:,:), & ! [IN]
1890  dummy(:,:,:), & ! [OUT]
1891  tagbase, i_ystg, handle, & ! [IN]
1892  isu_tag, isu_tagf ) ! [INOUT]
1893  else
1894  call comm_cartesc_nest_intercomm_nestdown( org_v_ll(:,:,:), & ! [IN]
1895  dummy(:,:,:), & ! [OUT]
1896  tagbase, i_sclr, handle, & ! [IN]
1897  isu_tag, isu_tagf ) ! [INOUT]
1898  endif
1899 
1900  tagbase = tagcomm + tag_rhot*order_tag_var
1901  call comm_cartesc_nest_intercomm_nestdown( org_rhot(:,:,:), & ! [IN]
1902  dummy(:,:,:), & ! [OUT]
1903  tagbase, i_sclr, handle, & ! [IN]
1904  isu_tag, isu_tagf ) ! [INOUT]
1905 
1906  do iq = 1, bnd_qa
1907  tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
1908  call comm_cartesc_nest_intercomm_nestdown( org_qtrc(:,:,:,iq), & ! [IN]
1909  dummy(:,:,:), & ! [OUT]
1910  tagbase, i_sclr, handle, & ! [IN]
1911  isu_tag, isu_tagf ) ! [INOUT]
1912  enddo
1913 
1914  rq_tot_p = rq_ctl_p
1915 
1916  call prof_rapend ('NEST_pack_P', 2)
1917  call prof_rapend ('NEST_total_P', 2)
1918 
1919  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
1920 
1921  !##### child [wait issue] #####
1922 
1923  call prof_rapstart('NEST_total_C', 2)
1924  call prof_rapstart('NEST_wait_C', 2)
1925 
1926  nwait_d = nwait_d + 1
1927  !LOG_INFO("COMM_CARTESC_NEST_nestdown",'(1X,A,I5,A)') "NestIDC [C]: que wait ( ", nwait_d, " )"
1928 
1929  !*** reset issue tag and request control
1930  !--- do not change the calling order below;
1931  !--- it should be consistent with the order in "COMM_CARTESC_NEST_recvwait_issue"
1932  isu_tag = 0
1933  isu_tagf = 0
1934 
1935  call comm_cartesc_nest_waitall( rq_tot_d, ireq_d )
1936 
1937  if ( online_aggressive_comm ) then
1938  ! nothing to do
1939  else
1940  call mpi_barrier(intercomm_parent, ierr)
1941  endif
1942 
1943  call prof_rapend ('NEST_wait_C', 2)
1944  call prof_rapstart('NEST_unpack_C', 2)
1945 
1946  tagbase = tagcomm + tag_dens*order_tag_var
1947  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
1948  work1_recv(:,:,:), & ! [OUT]
1949  tagbase, i_sclr, handle, & ! [IN]
1950  isu_tag, isu_tagf, & ! [INOUT]
1951  flag_dens = .true. ) ! [IN]
1952 !OCL XFILL
1953  do j = 1, daughter_ja(handle)
1954  do i = 1, daughter_ia(handle)
1955  do k = datr_ks(handle), datr_ke(handle)
1956  dens_recv(k,i,j) = work1_recv(k,i,j)
1957  enddo
1958  enddo
1959  enddo
1960 
1961  call comm_vars8( dens_recv, 1 )
1962 
1963  tagbase = tagcomm + tag_momz*order_tag_var
1964  if ( online_use_velz ) then
1965  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
1966  work2_recv(:,:,:), & ! [OUT]
1967  tagbase, i_zstg, handle, & ! [IN]
1968  isu_tag, isu_tagf ) ! [INOUT]
1969 !OCL XFILL
1970  do j = 1, daughter_ja(handle)
1971  do i = 1, daughter_ia(handle)
1972  do k = datr_ks(handle), datr_ke(handle)-1
1973  velz_recv(k,i,j) = work2_recv(k,i,j) / ( work1_recv(k,i,j) + work1_recv(k+1,i,j) ) * 2.0_rp
1974  enddo
1975  enddo
1976  enddo
1977 
1978  do j = 1, daughter_ja(handle)
1979  do i = 1, daughter_ia(handle)
1980  velz_recv(datr_ks(handle)-1,i,j) = 0.0_rp
1981  velz_recv(datr_ke(handle) ,i,j) = 0.0_rp
1982  enddo
1983  enddo
1984  endif
1985 
1986  call comm_wait ( dens_recv, 1, .false. )
1987 
1988  tagbase = tagcomm + tag_momx*order_tag_var
1989  if ( online_no_rotate ) then
1990  ! U_ll_recv receives MOMX
1991  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
1992  work1_recv(:,:,:), & ! [OUT]
1993  tagbase, i_xstg, handle, & ! [IN]
1994  isu_tag, isu_tagf ) ! [INOUT]
1995  else
1996  ! U_ll_recv receives MOMX/DENS
1997  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
1998  u_ll_recv(:,:,:), & ! [OUT]
1999  tagbase, i_sclr, handle, & ! [IN]
2000  isu_tag, isu_tagf ) ! [INOUT]
2001  endif
2002 
2003  tagbase = tagcomm + tag_momy*order_tag_var
2004  if ( online_no_rotate ) then
2005  ! V_ll_recv receives MOMY
2006  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2007  work2_recv(:,:,:), & ! [OUT]
2008  tagbase, i_ystg, handle, & ! [IN]
2009  isu_tag, isu_tagf ) ! [INOUT]
2010  else
2011  ! V_ll_recv receives MOMY/DENS
2012  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2013  v_ll_recv(:,:,:), & ! [OUT]
2014  tagbase, i_sclr, handle, & ! [IN]
2015  isu_tag, isu_tagf ) ! [INOUT]
2016  endif
2017 
2018  if ( online_no_rotate ) then
2019 
2020 !OCL XFILL
2021  do j = 1, daughter_ja(handle)
2022  do i = 1, daughter_ia(handle)-1
2023  do k = datr_ks(handle), datr_ke(handle)
2024  velx_recv(k,i,j) = work1_recv(k,i,j) / ( dens_recv(k,i+1,j) + dens_recv(k,i,j) ) * 2.0_rp
2025  enddo
2026  enddo
2027  enddo
2028 
2029  i = daughter_ia(handle)
2030 !OCL XFILL
2031  do j = 1, daughter_ja(handle)
2032  do k = datr_ks(handle), datr_ke(handle)
2033  velx_recv(k,i,j) = work1_recv(k,i,j) / dens_recv(k,i,j)
2034  enddo
2035  enddo
2036 
2037  call comm_vars8( velx_recv, 2 )
2038 
2039 !OCL XFILL
2040  do j = 1, daughter_ja(handle)-1
2041  do i = 1, daughter_ia(handle)
2042  do k = datr_ks(handle), datr_ke(handle)
2043  vely_recv(k,i,j) = work2_recv(k,i,j) / ( dens_recv(k,i,j+1) + dens_recv(k,i,j) ) * 2.0_rp
2044  enddo
2045  enddo
2046  enddo
2047 
2048  j = daughter_ja(handle)
2049 !OCL XFILL
2050  do i = 1, daughter_ia(handle)
2051  do k = datr_ks(handle), datr_ke(handle)
2052  vely_recv(k,i,j) = work2_recv(k,i,j) / dens_recv(k,i,j)
2053  enddo
2054  enddo
2055 
2056  call comm_vars8( vely_recv, 3 )
2057 
2058  call comm_wait ( velx_recv, 2, .false. )
2059  call comm_wait ( vely_recv, 3, .false. )
2060 
2061  else ! rotate
2062 
2063  ! rotation from latlon field to map-projected field
2064 !OCL XFILL
2065  do j = 1, daughter_ja(handle)
2066  do i = 1, daughter_ia(handle)
2067  do k = datr_ks(handle), datr_ke(handle)
2068  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)
2069  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)
2070  enddo
2071  enddo
2072  enddo
2073 
2074  ! from scalar point to staggered point
2075 !OCL XFILL
2076  do j = 1, daughter_ja(handle)
2077  do i = 1, daughter_ia(handle)-1
2078  do k = datr_ks(handle), datr_ke(handle)
2079  velx_recv(k,i,j) = ( work1_recv(k,i+1,j) + work1_recv(k,i,j) ) * 0.5_rp
2080  enddo
2081  enddo
2082  enddo
2083 
2084  i = daughter_ia(handle)
2085 !OCL XFILL
2086  do j = 1, daughter_ja(handle)
2087  do k = datr_ks(handle), datr_ke(handle)
2088  velx_recv(k,i,j) = work1_recv(k,i,j)
2089  enddo
2090  enddo
2091 
2092  call comm_vars8( velx_recv, 2 )
2093 
2094 !OCL XFILL
2095  do j = 1, daughter_ja(handle)-1
2096  do i = 1, daughter_ia(handle)
2097  do k = datr_ks(handle), datr_ke(handle)
2098  vely_recv(k,i,j) = ( work2_recv(k,i,j+1) + work2_recv(k,i,j) ) * 0.5_rp
2099  enddo
2100  enddo
2101  enddo
2102 
2103  j = daughter_ja(handle)
2104 !OCL XFILL
2105  do i = 1, daughter_ia(handle)
2106  do k = datr_ks(handle), datr_ke(handle)
2107  vely_recv(k,i,j) = work2_recv(k,i,j)
2108  enddo
2109  enddo
2110 
2111  call comm_vars8( vely_recv, 3 )
2112 
2113  call comm_wait ( velx_recv, 2, .false. )
2114  call comm_wait ( vely_recv, 3, .false. )
2115 
2116  endif
2117 
2118  tagbase = tagcomm + tag_rhot*order_tag_var
2119  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2120  work1_recv(:,:,:), & ! [OUT]
2121  tagbase, i_sclr, handle, & ! [IN]
2122  isu_tag, isu_tagf ) ! [INOUT]
2123 !OCL XFILL
2124  do j = 1, daughter_ja(handle)
2125  do i = 1, daughter_ia(handle)
2126  do k = datr_ks(handle), datr_ke(handle)
2127  pott_recv(k,i,j) = work1_recv(k,i,j) / dens_recv(k,i,j)
2128  enddo
2129  enddo
2130  enddo
2131 
2132  do iq = 1, bnd_qa
2133  tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
2134  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2135  work1_recv(:,:,:), & ! [OUT]
2136  tagbase, i_sclr, handle, & ! [IN]
2137  isu_tag, isu_tagf ) ! [INOUT]
2138 !OCL XFILL
2139  do j = 1, daughter_ja(handle)
2140  do i = 1, daughter_ia(handle)
2141  do k = datr_ks(handle), datr_ke(handle)
2142  qtrc_recv(k,i,j,iq) = work1_recv(k,i,j)
2143  enddo
2144  enddo
2145  enddo
2146  enddo
2147 
2148  call prof_rapend ('NEST_unpack_C', 2)
2149  call prof_rapend ('NEST_total_C', 2)
2150 
2151  else
2152  log_error("COMM_CARTESC_NEST_nestdown",*) 'internal error'
2153  call prc_abort
2154  endif
2155 
2156  return
module Atmosphere Grid CartesianC metirc
module COMMUNICATION
module PROCESS
Definition: scale_prc.F90:11
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_metric_rotc
rotation coefficient
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
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 2164 of file scale_comm_cartesC_nest.F90.

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_update(), and mod_atmos_bnd_driver::atmos_boundary_set_file().

2164  use scale_prc, only: &
2165  prc_abort
2166  implicit none
2167 
2168  integer, intent(in) :: handle
2169  integer, intent(in) :: bnd_qa
2170 
2171  integer :: isu_tag, isu_tagf
2172  integer :: tagbase, tagcomm
2173  integer :: ierr
2174  integer :: iq
2175  !---------------------------------------------------------------------------
2176 
2177  if( .NOT. use_nesting ) return
2178 
2179  if ( bnd_qa > i_bndqa ) then
2180  log_error("COMM_CARTESC_NEST_recvwait_issue",*) 'internal error: about BND_QA'
2181  call prc_abort
2182  endif
2183 
2184  tagcomm = intercomm_id(handle) * order_tag_comm
2185 
2186  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2187 
2188  !##### parent [wait issue] #####
2189 
2190  call prof_rapstart('NEST_total_P', 2)
2191  call prof_rapstart('NEST_wait_P', 2)
2192 
2193  nwait_p = nwait_p + 1
2194  !LOG_INFO("COMM_CARTESC_NEST_recvwait_issue",'(1X,A,I5,A)') "NestIDC [P]: que wait ( ", nwait_p, " )"
2195 
2196  call comm_cartesc_nest_issuer_of_wait( handle )
2197 
2198  if ( online_aggressive_comm ) then
2199  ! nothing to do
2200  else
2201  call mpi_barrier(intercomm_daughter, ierr)
2202  endif
2203 
2204  call prof_rapend ('NEST_wait_P', 2)
2205  call prof_rapend ('NEST_total_P', 2)
2206 
2207  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2208 
2209  !##### child [receive issue] #####
2210 
2211  call prof_rapstart('NEST_total_C', 2)
2212 
2213  nrecv = nrecv + 1
2214  log_info("COMM_CARTESC_NEST_recvwait_issue",'(1X,A,I5,A)') "NestIDC [C]: que recv ( ", nrecv, " )"
2215 
2216  !*** reset issue tag and request control
2217  !--- do not change the calling order below;
2218  !--- it should be consistent with the order in "COMM_CARTESC_NEST_nestdown"
2219  isu_tag = 0
2220  isu_tagf = 0
2221  rq_ctl_d = 0
2222 
2223  tagbase = tagcomm + tag_dens*order_tag_var
2224  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf )
2225 
2226  tagbase = tagcomm + tag_momz*order_tag_var
2227  if ( online_use_velz ) then
2228  call comm_cartesc_nest_issuer_of_receive( tagbase, i_zstg, handle, isu_tag, isu_tagf )
2229  endif
2230 
2231  tagbase = tagcomm + tag_momx*order_tag_var
2232  if ( online_no_rotate ) then
2233  call comm_cartesc_nest_issuer_of_receive( tagbase, i_xstg, handle, isu_tag, isu_tagf )
2234  else
2235  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf )
2236  endif
2237 
2238  tagbase = tagcomm + tag_momy*order_tag_var
2239  if ( online_no_rotate ) then
2240  call comm_cartesc_nest_issuer_of_receive( tagbase, i_ystg, handle, isu_tag, isu_tagf )
2241  else
2242  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf )
2243  endif
2244 
2245  tagbase = tagcomm + tag_rhot*order_tag_var
2246  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf )
2247 
2248  do iq = 1, bnd_qa
2249  tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
2250  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag, isu_tagf )
2251  enddo
2252 
2253  rq_tot_d = rq_ctl_d
2254 
2255  call prof_rapend('NEST_total_C', 2)
2256 
2257  else
2258  log_error("COMM_CARTESC_NEST_recvwait_issue",*) 'internal error'
2259  call prc_abort
2260  endif
2261 
2262  return
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
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 2269 of file scale_comm_cartesC_nest.F90.

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

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_finalize().

2269  use scale_prc, only: &
2270  prc_abort
2271  implicit none
2272 
2273  integer, intent(in) :: handle
2274 
2275  !logical :: flag
2276  !integer :: istatus(MPI_STATUS_SIZE)
2277 
2278  integer :: rq
2279  integer :: ierr
2280  !---------------------------------------------------------------------------
2281 
2282  if( .NOT. use_nesting ) return
2283 
2284  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2285 
2286  !##### parent #####
2287  ! Nothing to do
2288 
2289  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2290 
2291  !##### child #####
2292 
2293  log_info("COMM_CARTESC_NEST_recv_cancel",'(1X,A,I5,A)') "NestIDC [C]: CANCEL recv ( ", nrecv, " )"
2294 
2295  do rq = 1, rq_tot_d
2296  if ( ireq_d(rq) /= mpi_request_null ) then
2297 
2298  call mpi_cancel(ireq_d(rq), ierr)
2299 
2300 ! call MPI_TEST_CANCELLED(istatus, flag, ierr)
2301 ! if ( .NOT. flag ) then
2302 ! LOG_ERROR("COMM_CARTESC_NEST_recv_cancel",*) 'receive actions do not cancelled, req = ', rq
2303 ! endif
2304  endif
2305  enddo
2306 
2307  else
2308  log_error_cont(*) 'internal error'
2309  call prc_abort
2310  endif
2311 
2312  return
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
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,
integer, intent(inout)  isu_tagf,
logical, intent(in), optional  flag_dens 
)

Inter-communication from parent to daughter: nestdown.

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

Definition at line 2326 of file scale_comm_cartesC_nest.F90.

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

2326  use scale_prc, only: &
2327  prc_abort
2328  use scale_comm_cartesc, only: &
2330  use scale_interp, only: &
2332  implicit none
2333 
2334  real(RP), intent(in) :: pvar(:,:,:)
2335  real(RP), intent(out) :: dvar(:,:,:)
2336  integer, intent(in) :: tagbase
2337  integer, intent(in) :: id_stag
2338  integer, intent(in) :: handle
2339  integer, intent(inout) :: isu_tag
2340  integer, intent(inout) :: isu_tagf
2341 
2342  logical , intent(in), optional :: flag_dens
2343 
2344  integer :: ileng, tag, target_rank
2345 
2346  integer :: xloc, yloc
2347  integer :: gxs, gxe, gys, gye ! for large domain
2348  integer :: pxs, pxe, pys, pye ! for parent domain
2349  integer :: zs, ze
2350 
2351  integer :: ig, rq, yp
2352  logical :: no_zstag = .true.
2353  logical :: logarithmic = .false.
2354 
2355  integer :: ierr
2356  integer :: i, j
2357  !---------------------------------------------------------------------------
2358 
2359  if( .NOT. use_nesting ) return
2360 
2361  logarithmic = .false.
2362  if ( present(flag_dens) ) then
2363  if( flag_dens ) logarithmic = .true.
2364  endif
2365 
2366  if ( id_stag == i_sclr ) then
2367  no_zstag = .true.
2368  ig = i_sclr
2369  elseif( id_stag == i_zstg ) then
2370  no_zstag = .false.
2371  ig = i_zstg
2372  elseif( id_stag == i_xstg ) then
2373  no_zstag = .true.
2374  ig = i_xstg
2375  elseif( id_stag == i_ystg ) then
2376  no_zstag = .true.
2377  ig = i_ystg
2378  endif
2379 
2380  if ( no_zstag ) then
2381  ileng = (parent_ka(handle) ) * parent_ia(handle) * parent_ja(handle)
2382  else
2383  ileng = (parent_ka(handle)+1) * parent_ia(handle) * parent_ja(handle)
2384  endif
2385 
2386  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2387 
2388  !##### parent #####
2389  rq = rq_ctl_p
2390 
2391  do yp = 1, num_yp
2392  rq = rq + 1
2393 
2394  ! send data to multiple daughter processes
2395  target_rank = comm_cartesc_nest_tile_list_yp(yp)
2396  tag = tagbase + yp
2397 
2398  call mpi_isend( pvar, &
2399  ileng, &
2400  comm_datatype, &
2401  target_rank, &
2402  tag, &
2403  intercomm_daughter, &
2404  ireq_p(rq), &
2405  ierr )
2406 
2407  dvar(:,:,:) = -1.0_rp ! input as a dummy value
2408  enddo
2409 
2410  rq_ctl_p = rq
2411 
2412  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2413 
2414  !##### child #####
2415  rq = rq_ctl_d
2416 
2417  do yp = 1, comm_cartesc_nest_tile_all
2418  rq = rq + 1
2419 
2420  xloc = mod( yp-1, comm_cartesc_nest_tile_num_x ) + 1
2421  yloc = int( real(yp-1) / real(COMM_CARTESC_NEST_TILE_NUM_X) ) + 1
2422 
2423  gxs = parent_imax(handle) * (xloc-1) + 1
2424  gxe = parent_imax(handle) * xloc
2425  gys = parent_jmax(handle) * (yloc-1) + 1
2426  gye = parent_jmax(handle) * yloc
2427 
2428  pxs = prnt_is(handle)
2429  pxe = prnt_ie(handle)
2430  pys = prnt_js(handle)
2431  pye = prnt_je(handle)
2432 
2433  if ( no_zstag ) then
2434  isu_tag = isu_tag + 1
2435 
2436  zs = 1
2437  ze = parent_ka(handle)
2438 !OCL XFILL
2439  buffer_ref_3d(zs:ze,gxs:gxe,gys:gye) = recvbuf_3d(zs:ze,pxs:pxe,pys:pye,isu_tag)
2440  else
2441  isu_tagf = isu_tagf + 1
2442 
2443  zs = 0
2444  ze = parent_ka(handle)
2445 !OCL XFILL
2446  buffer_ref_3df(zs:ze,gxs:gxe,gys:gye) = recvbuf_3df(zs:ze,pxs:pxe,pys:pye,isu_tagf)
2447  endif
2448 
2449  if ( isu_tag > max_isu .OR. isu_tagf > max_isuf ) then
2450  log_error("COMM_CARTESC_NEST_intercomm_nestdown_3D",*) 'Exceeded maximum issue'
2451  log_error_cont(*) 'isu_tag = ', isu_tag
2452  log_error_cont(*) 'isu_tagf = ', isu_tagf
2453  call prc_abort
2454  endif
2455 
2456  enddo
2457 
2458  rq_ctl_d = rq
2459 
2460  if ( no_zstag ) then
2461  call interp_interp3d( itp_nh, & ! [IN]
2462  tileal_ka(handle), & ! [IN]
2463  tileal_ia(handle), & ! [IN]
2464  tileal_ja(handle), & ! [IN]
2465  daughter_ka(handle), & ! [IN]
2466  datr_ks(handle), & ! [IN]
2467  datr_ke(handle), & ! [IN]
2468  daughter_ia(handle), & ! [IN]
2469  daughter_ja(handle), & ! [IN]
2470  igrd( :,:,:,ig), & ! [IN]
2471  jgrd( :,:,:,ig), & ! [IN]
2472  hfact( :,:,:,ig), & ! [IN]
2473  kgrd(:,:,:,:,:,ig), & ! [IN]
2474  vfact(:,:,:,:,:,ig), & ! [IN]
2475  buffer_ref_3d(:,:,:), & ! [INOUT]
2476  dvar(:,:,:), & ! [OUT]
2477  logwgt = logarithmic ) ! [IN]
2478  else
2479  call interp_interp3d( itp_nh, & ! [IN]
2480  tileal_ka(handle)+1, & ! [IN]
2481  tileal_ia(handle), & ! [IN]
2482  tileal_ja(handle), & ! [IN]
2483  daughter_ka(handle), & ! [IN]
2484  datr_ks(handle), & ! [IN]
2485  datr_ke(handle), & ! [IN]
2486  daughter_ia(handle), & ! [IN]
2487  daughter_ja(handle), & ! [IN]
2488  igrd( :,:,:,ig), & ! [IN]
2489  jgrd( :,:,:,ig), & ! [IN]
2490  hfact( :,:,:,ig), & ! [IN]
2491  kgrd(:,:,:,:,:,ig), & ! [IN]
2492  vfact(:,:,:,:,:,ig), & ! [IN]
2493  buffer_ref_3df(:,:,:), & ! [INOUT]
2494  dvar(:,:,:), & ! [OUT]
2495  logwgt = logarithmic ) ! [IN]
2496  endif
2497 
2498  do j = 1, daughter_ja(handle)
2499  do i = 1, daughter_ia(handle)
2500  dvar( 1:datr_ks(handle)-1,i,j) = 0.0_rp
2501  dvar(datr_ke(handle)+1:daughter_ka(handle) ,i,j) = 0.0_rp
2502  enddo
2503  enddo
2504 
2505  else
2506  log_error("COMM_CARTESC_NEST_intercomm_nestdown_3D",*) 'internal error'
2507  call prc_abort
2508  endif
2509 
2510  return
integer, public comm_datatype
datatype of variable
module INTERPOLATION
subroutine, public interp_interp3d(npoints, KA_ref, IA_ref, JA_ref, KA, KS, KE, IA, JA, idx_i, idx_j, hfact, idx_k, vfact, val_ref, val, logwgt)
module COMMUNICATION
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
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,
integer, intent(inout)  isu_tagf 
)

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

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

Definition at line 2521 of file scale_comm_cartesC_nest.F90.

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.

2521  use scale_prc, only: &
2522  prc_myrank, &
2523  prc_abort
2524  use scale_comm_cartesc, only: &
2526  implicit none
2527 
2528  integer, intent(in) :: tagbase
2529  integer, intent(in) :: id_stag
2530  integer, intent(in) :: handle
2531  integer, intent(inout) :: isu_tag
2532  integer, intent(inout) :: isu_tagf
2533 
2534  integer :: ierr, ileng
2535  integer :: tag, target_rank
2536 
2537  integer :: ig, rq, yp
2538  logical :: no_zstag = .true.
2539  !---------------------------------------------------------------------------
2540 
2541  if( .NOT. use_nesting ) return
2542 
2543  if ( id_stag == i_sclr ) then
2544  no_zstag = .true.
2545  ig = i_sclr
2546  elseif( id_stag == i_zstg ) then
2547  no_zstag = .false.
2548  ig = i_zstg
2549  elseif( id_stag == i_xstg ) then
2550  no_zstag = .true.
2551  ig = i_xstg
2552  elseif( id_stag == i_ystg ) then
2553  no_zstag = .true.
2554  ig = i_ystg
2555  endif
2556 
2557  if ( no_zstag ) then
2558  ileng = (parent_ka(handle) ) * parent_ia(handle) * parent_ja(handle)
2559  else
2560  ileng = (parent_ka(handle)+1) * parent_ia(handle) * parent_ja(handle)
2561  endif
2562 
2563  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2564 
2565  !##### parent #####
2566  ! nothing to do
2567 
2568  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2569 
2570  !##### child #####
2571  rq = rq_ctl_d
2572 
2573  do yp = 1, comm_cartesc_nest_tile_all
2574  rq = rq + 1
2575 
2576  target_rank = comm_cartesc_nest_tile_list_d(yp,prc_myrank+1)
2577  tag = tagbase + call_order(yp)
2578 
2579  if ( no_zstag ) then
2580  isu_tag = isu_tag + 1
2581 
2582  recvbuf_3d(:,:,:,isu_tag) = 0.0_rp
2583 
2584  call mpi_irecv( recvbuf_3d(:,:,:,isu_tag), &
2585  ileng, &
2586  comm_datatype, &
2587  target_rank, &
2588  tag, &
2589  intercomm_parent, &
2590  ireq_d(rq), &
2591  ierr )
2592  else
2593  isu_tagf = isu_tagf + 1
2594 
2595  recvbuf_3df(:,:,:,isu_tagf) = 0.0_rp
2596 
2597  call mpi_irecv( recvbuf_3df(:,:,:,isu_tagf), &
2598  ileng, &
2599  comm_datatype, &
2600  target_rank, &
2601  tag, &
2602  intercomm_parent, &
2603  ireq_d(rq), &
2604  ierr )
2605  endif
2606 
2607  enddo
2608 
2609  if ( isu_tag > max_isu .OR. isu_tagf > max_isuf ) then
2610  log_error("COMM_CARTESC_NEST_issuer_of_receive_3D",*) 'Exceeded maximum issue'
2611  log_error_cont(*) 'isu_tag = ', isu_tag
2612  log_error_cont(*) 'isu_tagf = ', isu_tagf
2613  call prc_abort
2614  endif
2615 
2616  rq_ctl_d = rq
2617 
2618  else
2619  log_error("COMM_CARTESC_NEST_issuer_of_receive_3D",*) 'internal error'
2620  call prc_abort
2621  endif
2622 
2623  return
integer, public comm_datatype
datatype of variable
module COMMUNICATION
module PROCESS
Definition: scale_prc.F90:11
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:89
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
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 2630 of file scale_comm_cartesC_nest.F90.

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

2630  use scale_prc, only: &
2631  prc_abort
2632  implicit none
2633 
2634  integer, intent(in) :: handle
2635  !---------------------------------------------------------------------------
2636 
2637  if( .NOT. use_nesting ) return
2638 
2639  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2640 
2641  !##### parent #####
2642  call comm_cartesc_nest_waitall( rq_tot_p, ireq_p(:) )
2643 
2644  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2645 
2646  !##### child #####
2647  ! nothing to do
2648 
2649  else
2650  log_error("COMM_CARTESC_NEST_issuer_of_wait_3D",*) 'internal error'
2651  call prc_abort
2652  endif
2653 
2654  return
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
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 2662 of file scale_comm_cartesC_nest.F90.

References scale_prc::prc_abort(), and use_nesting.

Referenced by comm_cartesc_nest_issuer_of_wait_3d(), and comm_cartesc_nest_nestdown().

2662  use scale_prc, only: &
2663  prc_abort
2664  implicit none
2665 
2666  integer, intent(in) :: req_count
2667  integer, intent(inout) :: ireq(max_rq)
2668 
2669  integer :: i
2670  integer :: ierr
2671  integer :: istatus(mpi_status_size,req_count)
2672  integer :: req_count2
2673  integer :: ireq2(max_rq)
2674 
2675 ! logical :: flag = .false.
2676 ! integer(8) :: num = 0
2677  !---------------------------------------------------------------------------
2678 
2679  if( .NOT. use_nesting ) return
2680 
2681  req_count2 = 0
2682  do i = 1, req_count
2683  if ( ireq(i) /= mpi_request_null ) then
2684  req_count2 = req_count2 + 1
2685  ireq2(req_count2) = ireq(i)
2686  endif
2687  enddo
2688 
2689  if( req_count2 /= 0 ) call mpi_waitall( req_count2, ireq2(1:req_count2), istatus, ierr )
2690 
2691 ! do while ( .NOT. flag )
2692 ! num = num + 1
2693 ! call MPI_TESTALL( req_count, ireq, flag, istatus, ierr )
2694 !
2695 ! if ( num > ONLINE_WAIT_LIMIT ) then
2696 ! LOG_ERROR("COMM_CARTESC_NEST_waitall",'(1x,A)') 'over the limit of waiting time [NESTCOM]'
2697 ! LOG_ERROR_CONT('(1x,A)') 'over the limit of waiting time [NESTCOM]'
2698 ! call PRC_abort
2699 ! endif
2700 ! enddo
2701 
2702  return
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
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 2709 of file scale_comm_cartesC_nest.F90.

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

2709  use scale_prc, only: &
2710  prc_abort
2711  implicit none
2712 
2713  integer, intent(in) :: handle
2714 
2715  integer :: istatus(mpi_status_size)
2716  integer :: ierr
2717  logical :: flag
2718  !---------------------------------------------------------------------------
2719 
2720  if( .NOT. use_nesting ) return
2721 
2722  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2723 
2724  !##### parent #####
2725 
2726  call prof_rapstart('NEST_test_P', 2)
2727  if ( rq_ctl_p > 0 ) call mpi_test(ireq_p(1), flag, istatus, ierr)
2728  call prof_rapend('NEST_test_P', 2)
2729 
2730  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2731 
2732  !##### child #####
2733 
2734  call prof_rapstart('NEST_test_C', 2)
2735  if ( rq_ctl_d > 0 ) call mpi_test(ireq_d(1), flag, istatus, ierr)
2736  call prof_rapend('NEST_test_C', 2)
2737 
2738  else
2739  log_error("COMM_CARTESC_NEST_test",*) 'error'
2740  call prc_abort
2741  endif
2742 
2743  return
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
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 2749 of file scale_comm_cartesC_nest.F90.

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

2749  use scale_prc, only: &
2751  implicit none
2752 
2753  integer :: ierr
2754  !---------------------------------------------------------------------------
2755 
2756  if( .NOT. use_nesting ) return
2757 
2758  log_info("COMM_CARTESC_NEST_disconnect",'(1x,A)') 'Waiting finish of whole processes'
2759  call mpi_barrier(prc_global_comm_world, ierr)
2760 
2761  if ( online_iam_parent ) then
2762  !LOG_INFO("COMM_CARTESC_NEST_disconnect",'(1x,A)') 'Waiting finish of whole processes as a parent'
2763  !call MPI_BARRIER(INTERCOMM_DAUGHTER, ierr)
2764  call mpi_comm_free(intercomm_daughter, ierr)
2765  log_info("COMM_CARTESC_NEST_disconnect",'(1x,A)') 'Disconnected communication with child'
2766  endif
2767 
2768  if ( online_iam_daughter ) then
2769  !LOG_INFO("COMM_CARTESC_NEST_disconnect",'(1x,A)') 'Waiting finish of whole processes as a child'
2770  !call MPI_BARRIER(INTERCOMM_PARENT, ierr)
2771  call mpi_comm_free(intercomm_parent, ierr)
2772  log_info("COMM_CARTESC_NEST_disconnect",'(1x,A)') 'Disconnected communication with parent'
2773  endif
2774 
2775  return
module PROCESS
Definition: scale_prc.F90:11
integer, public prc_global_comm_world
global communicator
Definition: scale_prc.F90:78
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.

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

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

◆ intercomm_daughter

integer, public scale_comm_cartesc_nest::intercomm_daughter

Definition at line 46 of file scale_comm_cartesC_nest.F90.

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

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

◆ 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.

Referenced by comm_cartesc_nest_setup().

49  integer, public :: handling_num

◆ 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.

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

54  integer, public :: parent_kmax(2)

◆ 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.

Referenced by comm_cartesc_nest_setup().

60  integer, public :: parent_okmax(2)

◆ 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.

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

61  integer, public :: parent_lkmax(2)

◆ 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.

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

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

◆ 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.

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

63  integer, public :: parent_nstep(2)

◆ 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.

Referenced by comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_setup().

65  integer, public :: daughter_kmax(2)

◆ 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.

Referenced by comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_setup().

66  integer, public :: daughter_imax(2)

◆ 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.

Referenced by comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_setup().

67  integer, public :: daughter_jmax(2)

◆ 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.

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

68  integer, public :: daughter_ka(2)

◆ 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.

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

69  integer, public :: daughter_ia(2)

◆ 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.

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

70  integer, public :: daughter_ja(2)

◆ 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.

Referenced by comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_setup().

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

◆ 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.

Referenced by comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_setup().

74  integer, public :: daughter_nstep(2)

◆ 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.

Referenced by comm_cartesc_nest_domain_shape().

76  integer, public :: prnt_ks(2)

◆ 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.

Referenced by comm_cartesc_nest_domain_shape().

77  integer, public :: prnt_ke(2)

◆ 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.

Referenced by comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_intercomm_nestdown_3d().

78  integer, public :: prnt_is(2)

◆ 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.

Referenced by comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_intercomm_nestdown_3d().

79  integer, public :: prnt_ie(2)

◆ 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.

Referenced by comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_intercomm_nestdown_3d().

80  integer, public :: prnt_js(2)

◆ 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.

Referenced by comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_intercomm_nestdown_3d().

81  integer, public :: prnt_je(2)

◆ 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.

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

83  integer, public :: datr_ks(2)

◆ 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.

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

84  integer, public :: datr_ke(2)

◆ 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.

Referenced by comm_cartesc_nest_domain_shape().

85  integer, public :: datr_is(2)

◆ 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.

Referenced by comm_cartesc_nest_domain_shape().

86  integer, public :: datr_ie(2)

◆ 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.

Referenced by comm_cartesc_nest_domain_shape().

87  integer, public :: datr_js(2)

◆ 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.

Referenced by comm_cartesc_nest_domain_shape().

88  integer, public :: datr_je(2)

◆ 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.

Referenced by comm_cartesc_nest_intercomm_nestdown_3d(), and comm_cartesc_nest_setup().

90  integer, public :: tileal_ka(2)

◆ 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.

Referenced by comm_cartesc_nest_intercomm_nestdown_3d(), and comm_cartesc_nest_setup().

91  integer, public :: tileal_ia(2)

◆ 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.

Referenced by comm_cartesc_nest_intercomm_nestdown_3d(), and comm_cartesc_nest_setup().

92  integer, public :: tileal_ja(2)

◆ 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.

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

95  integer, public :: comm_cartesc_nest_interp_level = 5

◆ 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.

Referenced by comm_cartesc_nest_setup().

96  integer, public :: comm_cartesc_nest_interp_weight_order = 2

◆ 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.

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

99  logical, public :: offline = .false.

◆ 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.

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

100  logical, public :: online_iam_parent = .false.

◆ 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.

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

101  logical, public :: online_iam_daughter = .false.

◆ online_domain_num

integer, public scale_comm_cartesc_nest::online_domain_num = 1

Definition at line 102 of file scale_comm_cartesC_nest.F90.

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

102  integer, public :: online_domain_num = 1

◆ 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.

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

104  logical, public :: online_no_rotate = .false.

◆ 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.

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

105  logical, public :: online_boundary_use_qhyd = .false.

◆ online_boundary_diagqhyd

logical, public scale_comm_cartesc_nest::online_boundary_diagqhyd = .false.

Definition at line 106 of file scale_comm_cartesC_nest.F90.

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

106  logical, public :: online_boundary_diagqhyd = .false.