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

module Communication CartesianC nesting More...

Data Types

type  domain_info
 

Functions/Subroutines

subroutine, public comm_cartesc_nest_setup (QA_MP, MP_TYPE_in)
 Setup. More...
 
subroutine, public comm_cartesc_nest_domain_regist_file (dom_id, PARENT_BASENAME, PARENT_PRC_NUM_X, PARENT_PRC_NUM_Y, LATLON_CATALOGUE_FNAME)
 offline setup More...
 
subroutine comm_cartesc_nest_domain_relate (dom_id)
 Solve relationship between ParentDomain & Daughter Domain. More...
 
subroutine, public comm_cartesc_nest_parent_info (dom_id, KMAX, LKMAX, IMAXG, JMAXG, num_tile, tile_id)
 Return infomation of parent domain (for offline) More...
 
subroutine, public comm_cartesc_nest_domain_shape (tilei, tilej, cxs, cxe, cys, cye, pxs, pxe, pys, pye, dom_id, iloc, xstg, ystg)
 Return shape of ParentDomain at the specified rank (for offline) More...
 
subroutine, public comm_cartesc_nest_nestdown_send (DENS_send, MOMZ_send, MOMX_send, MOMY_send, RHOT_send, QTRC_send)
 Boundary data transfer from parent to daughter: nestdown (parent side) More...
 
subroutine, public comm_cartesc_nest_nestdown_recv (DENS_recv, VELZ_recv, VELX_recv, VELY_recv, POTT_recv, QTRC_recv)
 Boundary data transfer from parent to daughter: nestdown (daughter side) More...
 
subroutine, public comm_cartesc_nest_recvwait_issue_send
 Sub-command for data transfer from parent to daughter: nestdown (parent side) More...
 
subroutine, public comm_cartesc_nest_recvwait_issue_recv
 Sub-command for data transfer from parent to daughter: nestdown (daughter side) More...
 
subroutine, public comm_cartesc_nest_recv_cancel_send
 Sub-command for data transfer from parent to daughter: nestdown (parent side) More...
 
subroutine, public comm_cartesc_nest_recv_cancel_recv
 Sub-command for data transfer from parent to daughter: nestdown (daughter side) More...
 
subroutine comm_cartesc_nest_intercomm_nestdown_3d (pvar, dvar, tagbase, id_stag, HANDLE, isu_tag, flag_dens)
 Inter-communication from parent to daughter: nestdown. More...
 
subroutine comm_cartesc_nest_issuer_of_receive_3d (tagbase, id_stag, HANDLE, isu_tag)
 [substance of issuer] Inter-communication from parent to daughter: nestdown More...
 
subroutine comm_cartesc_nest_issuer_of_wait_3d (HANDLE)
 [substance of issuer] Inter-communication from parent to daughter: nestdown More...
 
subroutine comm_cartesc_nest_waitall (req_count, ireq)
 [substance of comm_wait] Inter-communication More...
 
subroutine, public comm_cartesc_nest_test_send
 [check communication status] Inter-communication (parent side) More...
 
subroutine, public comm_cartesc_nest_test_recv
 [check communication status] Inter-communication (daughter side) More...
 
subroutine, public comm_cartesc_nest_finalize
 finalize More...
 

Variables

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_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 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_recv_diagqhyd = .false.
 
logical, public online_send_diagqhyd = .false.
 
integer, public online_recv_qa = 0
 number of tracer received from the parent domain More...
 
integer, public online_send_qa = 0
 number of tracer sent to the daughter domain More...
 
real(dp), public online_parent_dtsec
 parent DT [sec] More...
 
integer, public online_parent_nstep
 parent nsteps More...
 
integer, public online_daughter_nprocs = -1
 

Detailed Description

module Communication CartesianC nesting

Description
Grid module for nesting system
Author
Team SCALE
NAMELIST
  • PARAM_COMM_CARTESC_NEST
    nametypedefault valuecomment
    ONLINE_DOMAIN_NUM integer 1
    ONLINE_IAM_PARENT logical .false. a flag to say "I am a parent"
    ONLINE_IAM_DAUGHTER logical .false. a flag to say "I am a daughter"
    ONLINE_USE_VELZ logical .false.
    ONLINE_NO_ROTATE logical .false.
    ONLINE_BOUNDARY_USE_QHYD logical .false.
    ONLINE_AGGRESSIVE_COMM logical
    ONLINE_WAIT_LIMIT integer(8) limit times of waiting loop in "COMM_CARTESC_NEST_waitall"
    ONLINE_SPECIFIED_MAXRQ integer 0
    COMM_CARTESC_NEST_INTERP_TYPE character(len=H_SHORT) 'LINEAR' ! "LINEAR" or "DIST-WEIGHT"
    COMM_CARTESC_NEST_INTERP_LEVEL integer 5 horizontal interpolation level
    COMM_CARTESC_NEST_INTERP_WEIGHT_ORDER integer 2 horizontal interpolation weight order

History Output
No history output

Function/Subroutine Documentation

◆ comm_cartesc_nest_setup()

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

Setup.

Definition at line 239 of file scale_comm_cartesC_nest.F90.

239  use scale_file, only: &
240  file_open, &
241  file_read, &
242  file_get_attribute, &
243  file_get_shape
244  use scale_const, only: &
245  d2r => const_d2r
246  use scale_time, only: &
247  time_nstep, &
248  time_dtsec
249  use scale_prc, only: &
250  prc_abort, &
252  prc_ismaster, &
255  use scale_interp, only: &
256  interp_setup, &
257  interp_factor3d
258  use scale_comm_cartesc, only: &
259  comm_bcast
260  use scale_atmos_grid_cartesc, only: &
265  use scale_atmos_grid_cartesc_real, only: &
276  use scale_atmos_hydrometeor, only: &
278  use scale_mapprojection, only: &
279  mapprojection_lonlat2xy
280  implicit none
281 
282  integer, intent(in) :: QA_MP
283  character(len=*), intent(in) :: MP_TYPE_in
284 
285  character(len=H_SHORT) :: COMM_CARTESC_NEST_INTERP_TYPE = 'LINEAR' ! "LINEAR" or "DIST-WEIGHT"
286  ! LINEAR : bi-linear interpolation
287  ! DIST-WEIGHT: distance-weighted mean of the nearest N-neighbors
288 
289  real(RP), allocatable :: X_ref(:,:)
290  real(RP), allocatable :: Y_ref(:,:)
291 
292  integer :: ONLINE_SPECIFIED_MAXRQ = 0
293  integer :: n, i, j
294  integer :: fid, ierr
295  integer :: parent_id
296 
297  logical :: flag_parent
298  logical :: flag_child
299 
300  integer :: nprocs
301  logical :: parent_periodic_x
302  logical :: parent_periodic_y
303 
304  logical :: error
305 
306  namelist / param_comm_cartesc_nest / &
307  online_domain_num, &
308  online_iam_parent, &
309  online_iam_daughter, &
310  online_use_velz, &
311  online_no_rotate, &
312  online_boundary_use_qhyd, &
313  online_aggressive_comm, &
314  online_wait_limit, &
315  online_specified_maxrq, &
316  comm_cartesc_nest_interp_type, &
317  comm_cartesc_nest_interp_level, &
318  comm_cartesc_nest_interp_weight_order
319 
320  !---------------------------------------------------------------------------
321 
322  log_newline
323  log_info("COMM_CARTESC_NEST_setup",*) 'Setup'
324 
325  flag_child = prc_intercomm_parent /= mpi_comm_null ! exist parent, so work as a child
326  flag_parent = prc_intercomm_child /= mpi_comm_null ! exist child, so work as a parent
327 
328  nwait_p = 0
329  nwait_d = 0
330  nrecv = 0
331  nsend = 0
332 
333  handling_num = 0
334  comm_cartesc_nest_filiation(:) = 0
335  online_wait_limit = 999999999
336  online_aggressive_comm = .true.
337 
338  !--- read namelist
339  rewind(io_fid_conf)
340  read(io_fid_conf,nml=param_comm_cartesc_nest,iostat=ierr)
341  if( ierr < 0 ) then !--- missing
342  log_info("COMM_CARTESC_NEST_setup",*) 'Not found namelist. Default used.'
343  elseif( ierr > 0 ) then !--- fatal error
344  log_error("COMM_CARTESC_NEST_setup",*) 'Not appropriate names in namelist PARAM_COMM_CARTESC_NEST. Check!'
345  call prc_abort
346  endif
347  log_nml(param_comm_cartesc_nest)
348 
349  prc_global_domainid = online_domain_num
350 
351 
352  if ( online_iam_daughter .or. online_iam_parent ) then
353  use_nesting = .true.
354  endif
355 
356  call interp_setup( comm_cartesc_nest_interp_weight_order ) ! [IN]
357 
358  select case ( comm_cartesc_nest_interp_type )
359  case ( 'LINEAR' )
360  itp_nh = 4
361  case ( 'DIST-WEIGHT' )
362  itp_nh = comm_cartesc_nest_interp_level
363  case default
364  log_error("COMM_CARTESC_NEST_setup",*) 'Unsupported type of COMM_CARTESC_NEST_INTERP_TYPE : ', trim(comm_cartesc_nest_interp_type)
365  log_error_cont(*) ' It must be "LINEAR" or "DIST-WEIGHT"'
366  call prc_abort
367  end select
368 
369 
370  latlon_local(i_min,i_lon) = minval(atmos_grid_cartesc_real_lonuv(:,:)) / d2r
371  latlon_local(i_max,i_lon) = maxval(atmos_grid_cartesc_real_lonuv(:,:)) / d2r
372  latlon_local(i_min,i_lat) = minval(atmos_grid_cartesc_real_latuv(:,:)) / d2r
373  latlon_local(i_max,i_lat) = maxval(atmos_grid_cartesc_real_latuv(:,:)) / d2r
374 
375  if ( .not. use_nesting ) return
376 
377 
378  debug_domain_num = online_domain_num
379  if( online_specified_maxrq > max_rq ) max_rq = online_specified_maxrq
380 
381  allocate( ireq_p(max_rq) )
382  allocate( ireq_d(max_rq) )
383  allocate( call_order(max_rq) )
384  ireq_p(:) = mpi_request_null
385  ireq_d(:) = mpi_request_null
386 
387 
388  ! ONLINE_(RECV|SEND)_QA can be modified according to the configuration in the other side
389  ! See COMM_CARTESC_NEST_parentsize
390 
391  if( online_boundary_use_qhyd ) then
392  mp_type = mp_type_in
393  online_recv_qa = qa_mp
394  elseif ( atmos_hydrometeor_dry ) then
395  mp_type = "DRY"
396  online_recv_qa = 0
397  else
398  mp_type = "QV"
399  online_recv_qa = 1
400  endif
401 
402  if ( atmos_hydrometeor_dry ) then
403  mp_type = "DRY"
404  online_send_qa = 0
405  else if ( mp_type_in == "NONE" ) then
406  mp_type = "QV"
407  online_send_qa = 1
408  else
409  mp_type = mp_type_in
410  online_send_qa = qa_mp
411  end if
412 
413  log_info("COMM_CARTESC_NEST_setup",*) "flag_parent", flag_parent, "flag_child", flag_child
414  log_info("COMM_CARTESC_NEST_setup",*) "ONLINE_IAM_PARENT", online_iam_parent, "ONLINE_IAM_DAUGHTER", online_iam_daughter
415 
416  if( flag_parent ) then ! must do first before daughter processes
417  !-------------------------------------------------
418  if ( .NOT. online_iam_parent ) then
419  log_error("COMM_CARTESC_NEST_setup",*) '[NEST_setup] Parent Flag from launcher is not consistent with namelist!'
420  log_error_cont(*) 'PARENT - domain : ', online_domain_num
421  call prc_abort
422  endif
423 
424  handling_num = 1 !HANDLING_NUM + 1
425  intercomm_id(handling_num) = online_domain_num
426  comm_cartesc_nest_filiation(intercomm_id(handling_num)) = 1
427 
428  log_info("COMM_CARTESC_NEST_setup",'(1x,A,I2,A)') 'Online Nesting - PARENT [INTERCOMM_ID:', &
429  intercomm_id(handling_num), ' ]'
430  log_info("COMM_CARTESC_NEST_setup",*) 'Online Nesting - INTERCOMM :', prc_intercomm_child
431 
432  call comm_cartesc_nest_ping( handling_num )
433  call comm_cartesc_nest_parentsize( handling_num )
434  call comm_cartesc_nest_catalogue( handling_num )
435  call mpi_barrier(prc_intercomm_child, ierr)
436 
437  log_info("COMM_CARTESC_NEST_setup",'(1x,A)' ) 'Informations of Daughter Domain'
438  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_nprocs :', online_daughter_nprocs
439  log_info_cont('(1x,A,I6) ') 'Limit Num. NCOMM req. :', max_rq
440 
441  allocate( org_dens(ka,ia,ja) )
442  allocate( org_momz(ka,ia,ja) )
443  allocate( org_momx(ka,ia,ja) )
444  allocate( org_momy(ka,ia,ja) )
445  allocate( org_u_ll(ka,ia,ja) )
446  allocate( org_v_ll(ka,ia,ja) )
447  allocate( org_rhot(ka,ia,ja) )
448  allocate( org_qtrc(ka,ia,ja,max(online_recv_qa,1)) )
449 
450  call comm_cartesc_nest_setup_nestdown( handling_num )
451 
452  !---------------------------------- end of parent routines
453  endif
454 
455 
456  if( flag_child ) then
457  !-------------------------------------------------
458  if ( .NOT. online_iam_daughter ) then
459  log_error("COMM_CARTESC_NEST_setup",*) '[NEST_setup] Child Flag from launcher is not consistent with namelist!'
460  log_error_cont(*) 'DAUGHTER - domain : ', online_domain_num
461  call prc_abort
462  endif
463 
464  handling_num = 2 !HANDLING_NUM + 1
465  intercomm_id(handling_num) = online_domain_num - 1
466  comm_cartesc_nest_filiation(intercomm_id(handling_num)) = -1
467 
468  log_info("COMM_CARTESC_NEST_setup",'(1x,A,I2,A)') 'Online Nesting - DAUGHTER [INTERCOMM_ID:', &
469  intercomm_id(handling_num), ' ]'
470  log_info("COMM_CARTESC_NEST_setup",*) 'Online Nesting - INTERCOMM :', prc_intercomm_parent
471 
472  num_dom = num_dom + 1
473  i_parent = num_dom
474  if ( i_parent > max_dinfo ) then
475  log_error("COMM_CARTESC_NEST_setup",*) 'number of domain exeeds the limit'
476  call prc_abort
477  end if
478 
479  call comm_cartesc_nest_ping( handling_num )
480 
481  call comm_cartesc_nest_parentsize( handling_num )
482 
483  nprocs = dom_info(i_parent)%prc_num_x * dom_info(i_parent)%prc_num_y
484  allocate( dom_info(i_parent)%latlon_catalogue(nprocs,2,2) )
485  call comm_cartesc_nest_catalogue( handling_num )
486  call mpi_barrier(prc_intercomm_parent, ierr)
487 
488  call comm_cartesc_nest_domain_relate( i_parent )
489 
490  tileal_ka = dom_info(i_parent)%KMAX + dom_info(i_parent)%KHALO * 2
491  tileal_ia = dom_info(i_parent)%IMAX * dom_info(i_parent)%tile_num_x
492  tileal_ja = dom_info(i_parent)%JMAX * dom_info(i_parent)%tile_num_y
493 
494  log_info("COMM_CARTESC_NEST_setup",'(1x,A)' ) 'Informations of Parent Domain'
495  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_nprocs :', nprocs
496  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_NUM_X :', dom_info(i_parent)%prc_num_x
497  log_info_cont('(1x,A,I6)' ) '--- PARENT_PRC_NUM_Y :', dom_info(i_parent)%prc_num_y
498  log_info_cont('(1x,A,I6)' ) '--- PARENT_KMAX :', dom_info(i_parent)%KMAX
499  log_info_cont('(1x,A,I6)' ) '--- PARENT_IMAX :', dom_info(i_parent)%IMAX
500  log_info_cont('(1x,A,I6)' ) '--- PARENT_JMAX :', dom_info(i_parent)%JMAX
501  log_info_cont('(1x,A,F9.3)') '--- PARENT_DTSEC :', online_parent_dtsec
502  log_info_cont('(1x,A,I6)' ) '--- PARENT_NSTEP :', online_parent_nstep
503  log_info_cont('(1x,A)' ) 'Informations of Daughter Domain [me]'
504  log_info_cont('(1x,A,F9.3)') '--- DAUGHTER_DTSEC :', time_dtsec
505  log_info_cont('(1x,A,I6)' ) '--- DAUGHTER_NSTEP :', time_nstep
506  log_info_cont('(1x,A)' ) 'Informations of Target Tiles'
507  log_info_cont('(1x,A,I6)' ) '--- TILEALL_KA :', tileal_ka
508  log_info_cont('(1x,A,I6)' ) '--- TILEALL_IA :', tileal_ia
509  log_info_cont('(1x,A,I6)' ) '--- TILEALL_JA :', tileal_ja
510  log_info_cont('(1x,A,I6) ') 'Limit Num. NCOMM req. :', max_rq
511 
512  allocate( buffer_ref_lon(tileal_ia,tileal_ja) )
513  allocate( buffer_ref_lonuy(tileal_ia,tileal_ja) )
514  allocate( buffer_ref_lonxv(tileal_ia,tileal_ja) )
515  allocate( buffer_ref_lat(tileal_ia,tileal_ja) )
516  allocate( buffer_ref_latuy(tileal_ia,tileal_ja) )
517  allocate( buffer_ref_latxv(tileal_ia,tileal_ja) )
518 
519  allocate( buffer_ref_cz(tileal_ka,tileal_ia,tileal_ja) )
520  allocate( buffer_ref_fz(tileal_ka,tileal_ia,tileal_ja) )
521 
522  allocate( buffer_ref_3d(tileal_ka,tileal_ia,tileal_ja) )
523 
524  allocate( igrd(ia,ja,itp_nh,itp_ng) )
525  allocate( jgrd(ia,ja,itp_nh,itp_ng) )
526  allocate( hfact(ia,ja,itp_nh,itp_ng) )
527  allocate( kgrd(ka,2,ia,ja,itp_nh,itp_ng) )
528  allocate( vfact(ka, ia,ja,itp_nh,itp_ng) )
529 
530  call comm_cartesc_nest_setup_nestdown( handling_num )
531 
532 
533  select case ( comm_cartesc_nest_interp_type )
534  case ( 'LINEAR' )
535 
536  allocate( x_ref(tileal_ia,tileal_ja) )
537  allocate( y_ref(tileal_ia,tileal_ja) )
538 
539  ! for scalar points
540  call mapprojection_lonlat2xy( tileal_ia, 1, tileal_ia, &
541  tileal_ja, 1, tileal_ja, &
542  buffer_ref_lon(:,:), & ! [IN]
543  buffer_ref_lat(:,:), & ! [IN]
544  x_ref(:,:), y_ref(:,:) ) ! [OUT]
545  call interp_factor3d( tileal_ka, khalo+1, tileal_ka-khalo, &
546  tileal_ia, tileal_ja, &
547  ka, ks, ke, ia, ja, &
548  x_ref(:,:), y_ref(:,:), & ! [IN]
549  buffer_ref_cz(:,:,:), & ! [IN]
550  atmos_grid_cartesc_cx(:), & ! [IN]
551  atmos_grid_cartesc_cy(:), & ! [IN]
552  atmos_grid_cartesc_real_cz(:,:,:), & ! [IN]
553  igrd( :,:,:,i_sclr), & ! [OUT]
554  jgrd( :,:,:,i_sclr), & ! [OUT]
555  hfact( :,:,:,i_sclr), & ! [OUT]
556  kgrd(:,:,:,:,:,i_sclr), & ! [OUT]
557  vfact(:, :,:,:,i_sclr) ) ! [OUT]
558 
559  ! for z staggered points
560  call interp_factor3d( tileal_ka+1, khalo+1, tileal_ka+1-khalo, &
561  tileal_ia, tileal_ja, &
562  ka, ks, ke, ia, ja, &
563  x_ref(:,:), y_ref(:,:), & ! [IN]
564  buffer_ref_fz(:,:,:), & ! [IN]
565  atmos_grid_cartesc_cx(:), & ! [IN]
566  atmos_grid_cartesc_cy(:), & ! [IN]
567  atmos_grid_cartesc_real_fz(1:ka,:,:), & ! [IN]
568  igrd( :,:,:,i_zstg), & ! [OUT]
569  jgrd( :,:,:,i_zstg), & ! [OUT]
570  hfact( :,:,:,i_zstg), & ! [OUT]
571  kgrd(:,:,:,:,:,i_zstg), & ! [OUT]
572  vfact(:, :,:,:,i_zstg) ) ! [OUT]
573 
574  ! for x staggered points
575  call mapprojection_lonlat2xy( tileal_ia, 1, tileal_ia, &
576  tileal_ja, 1, tileal_ja, &
577  buffer_ref_lonuy(:,:), & ! [IN]
578  buffer_ref_latuy(:,:), & ! [IN]
579  x_ref(:,:), y_ref(:,:) ) ! [OUT]
580  call interp_factor3d( tileal_ka, khalo+1, tileal_ka-khalo, &
581  tileal_ia, tileal_ja, &
582  ka, ks, ke, ia, ja, &
583  x_ref(:,:), y_ref(:,:), & ! [IN]
584  buffer_ref_cz(:,:,:), & ! [IN]
585  atmos_grid_cartesc_fx(1:ia), & ! [IN]
586  atmos_grid_cartesc_cy(:), & ! [IN]
587  atmos_grid_cartesc_real_cz(:,:,:), & ! [IN]
588  igrd( :,:,:,i_xstg), & ! [OUT]
589  jgrd( :,:,:,i_xstg), & ! [OUT]
590  hfact( :,:,:,i_xstg), & ! [OUT]
591  kgrd(:,:,:,:,:,i_xstg), & ! [OUT]
592  vfact(:, :,:,:,i_xstg) ) ! [OUT]
593 
594  ! for y staggered points
595  call mapprojection_lonlat2xy( tileal_ia, 1, tileal_ia, &
596  tileal_ja, 1, tileal_ja, &
597  buffer_ref_lonxv(:,:), & ! [IN]
598  buffer_ref_latxv(:,:), & ! [IN]
599  x_ref(:,:), y_ref(:,:) ) ! [OUT]
600  call interp_factor3d( tileal_ka, khalo+1, tileal_ka-khalo, &
601  tileal_ia, tileal_ja, &
602  ka, ks, ke, ia, ja, &
603  x_ref(:,:), y_ref(:,:), & ! [IN]
604  buffer_ref_cz(:,:,:), & ! [IN]
605  atmos_grid_cartesc_cx(:), & ! [IN]
606  atmos_grid_cartesc_fy(1:ja), & ! [IN]
607  atmos_grid_cartesc_real_cz(:,:,:), & ! [IN]
608  igrd( :,:,:,i_ystg), & ! [OUT]
609  jgrd( :,:,:,i_ystg), & ! [OUT]
610  hfact( :,:,:,i_ystg), & ! [OUT]
611  kgrd(:,:,:,:,:,i_ystg), & ! [OUT]
612  vfact(:, :,:,:,i_ystg) ) ! [OUT]
613 
614  deallocate( x_ref, y_ref )
615 
616  case ( 'DIST-WEIGHT' )
617 
618  ! for scalar points
619  call interp_factor3d( itp_nh, &
620  tileal_ka, khalo+1, tileal_ka-khalo, &
621  tileal_ia, tileal_ja, &
622  ka, ks, ke, ia, js, &
623  buffer_ref_lon(:,:), & ! [IN]
624  buffer_ref_lat(:,:), & ! [IN]
625  buffer_ref_cz(:,:,:), & ! [IN]
626  atmos_grid_cartesc_real_lon(:,:), & ! [IN]
627  atmos_grid_cartesc_real_lat(:,:), & ! [IN]
628  atmos_grid_cartesc_real_cz(:,:,:), & ! [IN]
629  igrd( :,:,:,i_sclr), & ! [OUT]
630  jgrd( :,:,:,i_sclr), & ! [OUT]
631  hfact( :,:,:,i_sclr), & ! [OUT]
632  kgrd(:,:,:,:,:,i_sclr), & ! [OUT]
633  vfact(:, :,:,:,i_sclr) ) ! [OUT]
634 
635  ! for z staggered points
636  call interp_factor3d( itp_nh, &
637  tileal_ka, khalo, tileal_ka-khalo, &
638  tileal_ia, tileal_ja, &
639  ka, ks, ke, ia, ja, &
640  buffer_ref_lon(:,:), & ! [IN]
641  buffer_ref_lat(:,:), & ! [IN]
642  buffer_ref_fz(:,:,:), & ! [IN]
643  atmos_grid_cartesc_real_lon(:,:), & ! [IN]
644  atmos_grid_cartesc_real_lat(:,:), & ! [IN]
645  atmos_grid_cartesc_real_fz(1:ka,:,:), & ! [IN]
646  igrd( :,:,:,i_zstg), & ! [OUT]
647  jgrd( :,:,:,i_zstg), & ! [OUT]
648  hfact( :,:,:,i_zstg), & ! [OUT]
649  kgrd(:,:,:,:,:,i_zstg), & ! [OUT]
650  vfact(:, :,:,:,i_zstg) ) ! [OUT]
651 
652  ! for x staggered points
653  call interp_factor3d( itp_nh, &
654  tileal_ka, khalo+1, tileal_ka-khalo, &
655  tileal_ia, tileal_ja, &
656  ka, ks, ke, ia, ja, &
657  buffer_ref_lonuy(:,:), & ! [IN]
658  buffer_ref_latuy(:,:), & ! [IN]
659  buffer_ref_cz(:,:,:), & ! [IN]
660  atmos_grid_cartesc_real_lonuy(1:ia,1:ja), & ! [IN]
661  atmos_grid_cartesc_real_latuy(1:ia,1:ja), & ! [IN]
662  atmos_grid_cartesc_real_cz(:,:,:), & ! [IN]
663  igrd( :,:,:,i_xstg), & ! [OUT]
664  jgrd( :,:,:,i_xstg), & ! [OUT]
665  hfact( :,:,:,i_xstg), & ! [OUT]
666  kgrd(:,:,:,:,:,i_xstg), & ! [OUT]
667  vfact(:, :,:,:,i_xstg) ) ! [OUT]
668 
669  ! for y staggered points
670  call interp_factor3d( itp_nh, &
671  tileal_ka, khalo+1, tileal_ka-khalo, &
672  tileal_ia, tileal_ja, &
673  ka, ks, ke, ia, ja, &
674  buffer_ref_lonxv(:,:), & ! [IN]
675  buffer_ref_latxv(:,:), & ! [IN]
676  buffer_ref_cz(:,:,:), & ! [IN]
677  atmos_grid_cartesc_real_lonxv(1:ia,1:ja), & ! [IN]
678  atmos_grid_cartesc_real_latxv(1:ia,1:ja), & ! [IN]
679  atmos_grid_cartesc_real_cz(:,:,:), & ! [IN]
680  igrd( :,:,:,i_ystg), & ! [OUT]
681  jgrd( :,:,:,i_ystg), & ! [OUT]
682  hfact( :,:,:,i_ystg), & ! [OUT]
683  kgrd(:,:,:,:,:,i_ystg), & ! [OUT]
684  vfact(:, :,:,:,i_ystg) ) ! [OUT]
685 
686  end select
687 
688  !---------------------------------- end of child routines
689  end if
690 
691  !LOG_INFO("COMM_CARTESC_NEST_setup",'(1x,A,I2)') 'Number of Related Domains :', HANDLING_NUM
692  !if ( HANDLING_NUM > 2 ) then
693  ! f( IO_L ) LOG_ERROR("COMM_CARTESC_NEST_setup",*) 'Too much handing domains (up to 2)'
694  ! call PRC_abort
695  !endif
696 
697  return

References scale_atmos_grid_cartesc::atmos_grid_cartesc_cx, scale_atmos_grid_cartesc::atmos_grid_cartesc_cy, scale_atmos_grid_cartesc::atmos_grid_cartesc_fx, scale_atmos_grid_cartesc::atmos_grid_cartesc_fy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lat, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latuv, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latuy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latxv, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lon, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonuv, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonuy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonxv, scale_atmos_hydrometeor::atmos_hydrometeor_dry, comm_cartesc_nest_domain_relate(), comm_cartesc_nest_filiation, comm_cartesc_nest_interp_level, comm_cartesc_nest_interp_weight_order, scale_const::const_d2r, scale_debug::debug_domain_num, scale_file::file_open(), handling_num, scale_atmos_grid_cartesc_index::ia, scale_interp::interp_setup(), scale_io::io_fid_conf, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::js, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::khalo, scale_atmos_grid_cartesc_index::ks, online_boundary_use_qhyd, online_daughter_nprocs, online_domain_num, online_iam_daughter, online_iam_parent, online_no_rotate, online_parent_dtsec, online_parent_nstep, online_recv_qa, online_send_qa, online_use_velz, scale_prc::prc_abort(), scale_prc::prc_global_domainid, scale_prc::prc_intercomm_child, scale_prc::prc_intercomm_parent, scale_prc::prc_ismaster, scale_time::time_dtsec, scale_time::time_nstep, and use_nesting.

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

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

◆ comm_cartesc_nest_domain_regist_file()

subroutine, public scale_comm_cartesc_nest::comm_cartesc_nest_domain_regist_file ( integer, intent(out)  dom_id,
character(len=*), intent(in)  PARENT_BASENAME,
integer, intent(in), optional  PARENT_PRC_NUM_X,
integer, intent(in), optional  PARENT_PRC_NUM_Y,
character(len=*), intent(in), optional  LATLON_CATALOGUE_FNAME 
)

offline setup

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

Definition at line 708 of file scale_comm_cartesC_nest.F90.

708  use scale_prc, only: &
709  prc_ismaster, &
710  prc_abort
711  use scale_file, only: &
712  file_open, &
713  file_get_attribute, &
714  file_get_shape, &
715  file_read
716  use scale_comm_cartesc, only: &
717  comm_bcast
718  integer, intent(out) :: dom_id
719 
720  character(len=*), intent(in) :: PARENT_BASENAME
721  integer, intent(in), optional :: PARENT_PRC_NUM_X
722  integer, intent(in), optional :: PARENT_PRC_NUM_Y
723  character(len=*), intent(in), optional :: LATLON_CATALOGUE_FNAME
724 
725  type(domain_info), pointer :: dinfo
726 
727  integer :: nprocs
728  integer :: pnum_x(1), pnum_y(1)
729  integer :: imaxg(1), jmaxg(1)
730  integer :: dims(1), dims2(1)
731  integer :: halos(2)
732  integer :: parent_x, parent_xh
733  integer :: parent_y, parent_yh
734 
735  real(RP), allocatable :: work(:,:), work_uv(:,:), minmax(:,:,:)
736 
737  character(len=H_LONG) :: fname
738  integer :: fid
739  integer :: parent_id
740  logical :: existed, error
741  integer :: ierr
742 
743  integer :: i, j, n
744 
745  do n = 1, num_dom
746  if ( dom_info(n)%basename == parent_basename ) then
747  dom_id = n
748  return
749  end if
750  end do
751 
752  num_dom = num_dom + 1
753  dom_id = num_dom
754  if ( dom_id > max_dinfo ) then
755  log_error("COMM_CARTESC_NEST_domain_regist_file",*) 'number of domains exceed the limit'
756  call prc_abort
757  end if
758  dinfo => dom_info(dom_id)
759  dinfo%basename = parent_basename
760 
761  if ( prc_ismaster ) then
762  call file_open( parent_basename, & ! (in)
763  fid, & ! (out)
764  aggregate = .false., & ! (in)
765  allnodes = .false. ) ! (in)
766 
767  call file_get_attribute( fid, "global", "scale_atmos_grid_cartesC_index_imaxg", &
768  imaxg(:), existed=existed )
769  if ( existed ) then
770  call file_get_attribute( fid, "global", "scale_cartesC_prc_num_x", &
771  pnum_x(:) )
772  dinfo%prc_num_x = pnum_x(1)
773  dinfo%IMAX = imaxg(1) / pnum_x(1)
774 
775  call file_get_attribute( fid, "x", "halo_global", & ! (in)
776  halos(:) ) ! (out)
777  dinfo%IHALO = halos(1)
778  call file_get_attribute( fid, "global", "scale_cartesC_prc_periodic_x", &
779  dinfo%periodic_x )
780 
781  call file_get_attribute( fid, "global", "scale_atmos_grid_cartesC_index_jmaxg", &
782  jmaxg(:) )
783  call file_get_attribute( fid, "global", "scale_cartesC_prc_num_y", &
784  pnum_y(:) )
785  dinfo%prc_num_y = pnum_y(1)
786  dinfo%JMAX = jmaxg(1) / pnum_y(1)
787  call file_get_attribute( fid, "y", "halo_global", & ! (in)
788  halos(:) ) ! (out)
789  dinfo%JHALO = halos(1)
790  call file_get_attribute( fid, "global", "scale_cartesC_prc_periodic_y", &
791  dinfo%periodic_y )
792 
793  else
794  ! for old file (for backward compatibility)
795 
796  if ( present(parent_prc_num_x) .and. present(parent_prc_num_y) ) then
797  dinfo%prc_num_x = parent_prc_num_x
798  dinfo%prc_num_y = parent_prc_num_y
799  else
800  log_error("COMM_CARTESC_NEST_domain_regist_file",*) 'PARENT_PRC_NUM_(X|Y) is needed for files generated by the older version'
801  call prc_abort
802  end if
803 
804  call file_get_shape( fid, "CX", dims(:) )
805  call file_get_shape( fid, "x", dims2(:) )
806  dinfo%IHALO = dims2(1) + ihalo * 2 - dims(1) ! assume IHALO is the same
807  dinfo%IMAX = dims(1) - dinfo%IHALO*2
808 
809  call file_get_shape( fid, "CY", dims(:) )
810  call file_get_shape( fid, "y", dims2(:) )
811  dinfo%JHALO = dims2(1) + jhalo * 2 - dims(1) ! assume JHALO is the same
812  dinfo%JMAX = dims(1) - dinfo%JHALO*2
813 
814  dinfo%periodic_x = .false.
815  dinfo%periodic_y = .false.
816 
817  endif
818 
819  call file_get_attribute( fid, "global", "scale_atmos_grid_cartesC_index_kmax", &
820  dims(:), existed=existed )
821  if ( existed ) then
822  dinfo%KMAX = dims(1)
823  else
824  call file_get_shape( fid, "z", dims(:), error=error )
825  if ( error ) then
826  dinfo%KMAX = 0
827  else
828  dinfo%KMAX = dims(1)
829  endif
830  end if
831  dinfo%KHALO = 0
832 
833  call file_get_attribute( fid, "global", "scale_ocean_grid_cartesC_index_kmax", &
834  dims(:), existed=existed )
835  if ( existed ) then
836  dinfo%OKMAX = dims(1)
837  else
838  call file_get_shape( fid, "oz", dims(:), error=error )
839  if ( error ) then
840  dinfo%OKMAX = 0
841  else
842  dinfo%OKMAX = dims(1)
843  endif
844  end if
845 
846  call file_get_attribute( fid, "global", "scale_land_grid_cartesC_index_kmax", &
847  dims(:), existed=existed )
848  if ( existed ) then
849  dinfo%LKMAX = dims(1)
850  else
851  call file_get_shape( fid, "lz", dims(:), error=error )
852  if ( error ) then
853  dinfo%LKMAX = 0
854  else
855  dinfo%LKMAX = dims(1)
856  endif
857  end if
858 
859  call file_get_attribute( fid, "global", "scale_urban_grid_cartesC_index_kmax", &
860  dims(:), existed=existed )
861  if ( existed ) then
862  dinfo%UKMAX = dims(1)
863  else
864  call file_get_shape( fid, "uz", dims(:), error=error )
865  if ( error ) then
866  dinfo%UKMAX = 0
867  else
868  dinfo%UKMAX = dims(1)
869  endif
870  end if
871 
872  end if ! master node
873 
874 
875  call comm_bcast( dinfo%prc_num_x )
876  call comm_bcast( dinfo%prc_num_y )
877  call comm_bcast( dinfo%KMAX )
878  call comm_bcast( dinfo%OKMAX )
879  call comm_bcast( dinfo%LKMAX )
880  call comm_bcast( dinfo%UKMAX )
881  call comm_bcast( dinfo%IMAX )
882  call comm_bcast( dinfo%JMAX )
883  call comm_bcast( dinfo%IHALO )
884  call comm_bcast( dinfo%JHALO )
885  call comm_bcast( dinfo%periodic_x )
886  call comm_bcast( dinfo%periodic_y )
887 
888  !--- latlon catalogue
889  nprocs = dinfo%prc_num_x * dinfo%prc_num_y
890  allocate( dinfo%latlon_catalogue(nprocs,2,2) )
891 
892  if ( prc_ismaster ) then
893 
894  existed = present(latlon_catalogue_fname)
895  if ( existed ) then
896  existed = latlon_catalogue_fname /= ""
897  end if
898  if ( existed ) then
899  ! read from catalogue file
900 
901  fid = io_get_available_fid()
902  call io_get_fname(fname, latlon_catalogue_fname)
903  open( fid, &
904  file = fname, &
905  form = 'formatted', &
906  status = 'old', &
907  iostat = ierr )
908 
909  if ( ierr /= 0 ) then
910  log_error("COMM_CARTESC_NEST_domain_regist_file",*) 'cannot open latlon-catalogue file!: ', trim(fname)
911  call prc_abort
912  endif
913 
914  do i = 1, nprocs
915  read(fid,'(i8,4f32.24)',iostat=ierr) parent_id, &
916  dinfo%latlon_catalogue(i,i_min,i_lon), dinfo%latlon_catalogue(i,i_max,i_lon), & ! LON: MIN, MAX
917  dinfo%latlon_catalogue(i,i_min,i_lat), dinfo%latlon_catalogue(i,i_max,i_lat) ! LAT: MIN, MAX
918  if ( ierr /= 0 .or. i /= parent_id ) then
919  log_error("COMM_CARTESC_NEST_domain_regist_file",*) 'catalogue file is invalid, ', trim(fname)
920  call prc_abort
921  end if
922  if ( ierr /= 0 ) exit
923  enddo
924  close(fid)
925 
926  else
927  ! read from netcdf file
928 
929  allocate( minmax(nprocs,2,2) )
930 
931  n = 1
932  do j = 1, dinfo%prc_num_y
933  do i = 1, dinfo%prc_num_x
934  call file_open( parent_basename, & ! (in)
935  fid, & ! (out)
936  aggregate = .false., & ! (in)
937  allnodes = .false., & ! (in)
938  rankid = n-1 ) ! (in)
939 
940  call file_get_shape( fid, "xh", dims(:) )
941  parent_xh = dims(1)
942  call file_get_shape( fid, "yh", dims(:) )
943  parent_yh = dims(1)
944  allocate( work_uv( parent_xh, parent_yh ) )
945 
946  if ( dinfo%periodic_x .or. dinfo%periodic_y ) then
947  call file_get_shape( fid, "x", dims(:) )
948  parent_x = dims(1)
949  call file_get_shape( fid, "y", dims(:) )
950  parent_y = dims(1)
951  end if
952 
953  call file_read( fid, "lon_uv", work_uv(:,:) )
954  dinfo%latlon_catalogue(n,i_min,i_lon) = minval( work_uv(:,:) )
955  dinfo%latlon_catalogue(n,i_max,i_lon) = maxval( work_uv(:,:) )
956 
957  if ( i > 1 ) then
958  dinfo%latlon_catalogue(n,i_min,i_lon) = min( dinfo%latlon_catalogue(n,i_min,i_lon), minmax(n-1,i_min,i_lon) )
959  dinfo%latlon_catalogue(n,i_max,i_lon) = max( dinfo%latlon_catalogue(n,i_max,i_lon), minmax(n-1,i_max,i_lon) )
960  else
961  if ( dinfo%periodic_x ) then
962  allocate( work( parent_x, parent_yh ) )
963  call file_read( fid, "lon_xv", work(:,:) )
964  ! This assumes an equally spaced grid
965  work(1,:) = work(1,:) * 2.0_rp - work_uv(1,:)
966  dinfo%latlon_catalogue(n,i_min,i_lon) = min( dinfo%latlon_catalogue(n,i_min,i_lon), minval( work ) )
967  dinfo%latlon_catalogue(n,i_max,i_lon) = max( dinfo%latlon_catalogue(n,i_max,i_lon), maxval( work ) )
968  deallocate( work )
969  end if
970  end if
971  minmax(n,i_min,i_lon) = minval( work_uv(parent_xh,:) )
972  minmax(n,i_max,i_lon) = maxval( work_uv(parent_xh,:) )
973 
974  call file_read( fid, "lat_uv", work_uv(:,:) )
975  dinfo%latlon_catalogue(n,i_min,i_lat) = minval( work_uv(:,:) )
976  dinfo%latlon_catalogue(n,i_max,i_lat) = maxval( work_uv(:,:) )
977 
978  if ( j > 1 ) then
979  dinfo%latlon_catalogue(n,i_min,i_lat) = min( dinfo%latlon_catalogue(n,i_min,i_lat), minmax(n-dinfo%prc_num_x,i_min,i_lat) )
980  dinfo%latlon_catalogue(n,i_max,i_lat) = max( dinfo%latlon_catalogue(n,i_max,i_lat), minmax(n-dinfo%prc_num_x,i_max,i_lat) )
981  else
982  if ( dinfo%periodic_y ) then
983  allocate( work( parent_xh, parent_y ) )
984  call file_read( fid, "lat_uy", work(:,:) )
985  ! This assumes an equally spaced grid
986  work(:,1) = work(:,1) * 2.0_rp - work_uv(:,1)
987  dinfo%latlon_catalogue(n,i_min,i_lat) = min( dinfo%latlon_catalogue(n,i_min,i_lat), minval( work(:,1) ) )
988  dinfo%latlon_catalogue(n,i_max,i_lat) = max( dinfo%latlon_catalogue(n,i_max,i_lat), maxval( work(:,1) ) )
989  deallocate( work )
990  end if
991  end if
992  minmax(n,i_min,i_lat) = minval( work_uv(:,parent_yh) )
993  minmax(n,i_max,i_lat) = maxval( work_uv(:,parent_yh) )
994 
995  deallocate( work_uv )
996 
997  n = n + 1
998  enddo
999  enddo
1000 
1001  deallocate( minmax )
1002 
1003  endif
1004 
1005  end if ! master node
1006 
1007  call comm_bcast( nprocs, 2, 2, dinfo%latlon_catalogue(:,:,:) )
1008 
1009  call comm_cartesc_nest_domain_relate(dom_id)
1010 
1011 

References comm_cartesc_nest_domain_relate(), scale_file::file_open(), scale_atmos_grid_cartesc_index::ihalo, scale_atmos_grid_cartesc_index::imaxg, scale_io::io_get_available_fid(), scale_io::io_get_fname(), scale_atmos_grid_cartesc_index::jhalo, scale_atmos_grid_cartesc_index::jmaxg, scale_prc::prc_abort(), and scale_prc::prc_ismaster.

Referenced by mod_copytopo::copytopo(), mod_realinput_netcdf::parentatmossetupnetcdf(), mod_realinput_netcdf::parentlandsetupnetcdf(), and mod_realinput_netcdf::parentoceansetupnetcdf().

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

◆ comm_cartesc_nest_domain_relate()

subroutine scale_comm_cartesc_nest::comm_cartesc_nest_domain_relate ( integer, intent(in)  dom_id)

Solve relationship between ParentDomain & Daughter Domain.

Parameters
[in]dom_idid number of domain information

Definition at line 1018 of file scale_comm_cartesC_nest.F90.

1018  use scale_prc, only: &
1019  prc_myrank, &
1020  prc_abort
1021  implicit none
1022 
1023  integer, intent(in) :: dom_id
1024 
1025  type(domain_info), pointer :: dinfo
1026 
1027  integer :: nprocs
1028  integer :: x_min, x_max
1029  integer :: y_min, y_max
1030  logical :: hit(2,2)
1031  real(RP) :: dx, dy
1032  integer :: p, i, j
1033  !---------------------------------------------------------------------------
1034 
1035  if ( dom_id < 1 .or. dom_id > num_dom ) then
1036  log_error("COMM_CARTESC_NEST_domain_relate",*) "domain id is invalid: ", dom_id
1037  call prc_abort
1038  end if
1039 
1040  dinfo => dom_info(dom_id)
1041  nprocs = dinfo%prc_num_x * dinfo%prc_num_y
1042 
1043  x_min = dinfo%prc_num_x
1044  x_max = -1
1045  y_min = dinfo%prc_num_y
1046  y_max = -1
1047  hit(:,:) = .false.
1048 
1049  do p = 1, nprocs
1050  dx = ( dinfo%latlon_catalogue(p,i_max,i_lon) - dinfo%latlon_catalogue(p,i_min,i_lon) ) / dinfo%IMAX
1051  dy = ( dinfo%latlon_catalogue(p,i_max,i_lat) - dinfo%latlon_catalogue(p,i_min,i_lat) ) / dinfo%JMAX
1052  if ( ( ( latlon_local(i_min,i_lon) >= dinfo%latlon_catalogue(p,i_min,i_lon) - dx &
1053  .AND. latlon_local(i_min,i_lon) <= dinfo%latlon_catalogue(p,i_max,i_lon) + dx ) .OR. &
1054  ( latlon_local(i_max,i_lon) >= dinfo%latlon_catalogue(p,i_min,i_lon) - dx &
1055  .AND. latlon_local(i_max,i_lon) <= dinfo%latlon_catalogue(p,i_max,i_lon) + dx ) .OR. &
1056  ( dinfo%latlon_catalogue(p,i_min,i_lon) >= latlon_local(i_min,i_lon) - dx &
1057  .AND. dinfo%latlon_catalogue(p,i_min,i_lon) <= latlon_local(i_max,i_lon) + dx ) .OR. &
1058  ( dinfo%latlon_catalogue(p,i_max,i_lon) >= latlon_local(i_min,i_lon) - dx &
1059  .AND. dinfo%latlon_catalogue(p,i_max,i_lon) <= latlon_local(i_max,i_lon) + dx ) ) .AND. &
1060  ( ( latlon_local(i_min,i_lat) >= dinfo%latlon_catalogue(p,i_min,i_lat) - dy &
1061  .AND. latlon_local(i_min,i_lat) <= dinfo%latlon_catalogue(p,i_max,i_lat) + dy ) .OR. &
1062  ( latlon_local(i_max,i_lat) >= dinfo%latlon_catalogue(p,i_min,i_lat) - dy &
1063  .AND. latlon_local(i_max,i_lat) <= dinfo%latlon_catalogue(p,i_max,i_lat) + dy ) .OR. &
1064  ( dinfo%latlon_catalogue(p,i_min,i_lat) >= latlon_local(i_min,i_lat) - dy &
1065  .AND. dinfo%latlon_catalogue(p,i_min,i_lat) <= latlon_local(i_max,i_lat) + dy ) .OR. &
1066  ( dinfo%latlon_catalogue(p,i_max,i_lat) >= latlon_local(i_min,i_lat) - dy &
1067  .AND. dinfo%latlon_catalogue(p,i_max,i_lat) <= latlon_local(i_max,i_lat) + dy ) ) ) then
1068  if ( dinfo%latlon_catalogue(p,i_min,i_lon) <= latlon_local(i_min,i_lon) ) hit(i_min,i_lon) = .true.
1069  if ( dinfo%latlon_catalogue(p,i_max,i_lon) >= latlon_local(i_max,i_lon) ) hit(i_max,i_lon) = .true.
1070  if ( dinfo%latlon_catalogue(p,i_min,i_lat) <= latlon_local(i_min,i_lat) ) hit(i_min,i_lat) = .true.
1071  if ( dinfo%latlon_catalogue(p,i_max,i_lat) >= latlon_local(i_max,i_lat) ) hit(i_max,i_lat) = .true.
1072  i = mod(p-1, dinfo%prc_num_x)
1073  j = (p-1) / dinfo%prc_num_x
1074  if ( i < x_min ) x_min = i
1075  if ( i > x_max ) x_max = i
1076  if ( j < y_min ) y_min = j
1077  if ( j > y_max ) y_max = j
1078  end if
1079  end do
1080 
1081  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
1082  log_error("COMM_CARTESC_NEST_domain_relate",*) 'region of daughter domain is larger than that of parent'
1083  log_error_cont(*) ' at rank:', prc_myrank, ' of domain:', online_domain_num
1084  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)
1085  log_error_cont('(A,F12.6,1x,F12.6)') 'daughter local (me) MIN-MAX: LON=', &
1086  latlon_local(i_min,i_lon), latlon_local(i_max,i_lon)
1087  do p = 1, nprocs
1088  log_error_cont('(A,I5,A,F12.6,1x,F12.6)') ' parent (', p,') MIN-MAX: LON=', &
1089  dinfo%latlon_catalogue(p,i_min,i_lon) ,dinfo%latlon_catalogue(p,i_max,i_lon)
1090  enddo
1091  log_error_cont('(A,F12.6,1x,F12.6)') 'daughter local (me): MIN-MAX LAT=', &
1092  latlon_local(i_min,i_lat), latlon_local(i_max,i_lat)
1093  do p = 1, nprocs
1094  log_error_cont('(A,I5,A,F12.6,1x,F12.6)') ' parent (', p,') MIN-MAX: LAT=', &
1095  dinfo%latlon_catalogue(p,i_min,i_lat) ,dinfo%latlon_catalogue(p,i_max,i_lat)
1096  enddo
1097  call prc_abort
1098  end if
1099 
1100 
1101 
1102  dinfo%tile_num_x = x_max - x_min + 1
1103  dinfo%tile_num_y = y_max - y_min + 1
1104 
1105  allocate( dinfo%tile_id(dinfo%tile_num_x * dinfo%tile_num_y) )
1106 
1107  log_info("COMM_CARTESC_NEST_domain_relate",'(1x,A)') 'NEST: target process tile in parent domain'
1108  p = 1
1109  do j = 1, dinfo%tile_num_y
1110  do i = 1, dinfo%tile_num_x
1111  dinfo%tile_id(p) = x_min + i - 1 + (y_min + j - 1) * dinfo%prc_num_x
1112  log_info_cont('(1x,A,I4,A,I6)') '(', p, ') target mpi-process:', dinfo%tile_id(p)
1113  p = p + 1
1114  enddo
1115  enddo
1116 
1117  return

References online_domain_num, scale_prc::prc_abort(), and scale_prc::prc_myrank.

Referenced by comm_cartesc_nest_domain_regist_file(), and comm_cartesc_nest_setup().

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

◆ comm_cartesc_nest_parent_info()

subroutine, public scale_comm_cartesc_nest::comm_cartesc_nest_parent_info ( integer, intent(in)  dom_id,
integer, intent(out), optional  KMAX,
integer, intent(out), optional  LKMAX,
integer, intent(out), optional  IMAXG,
integer, intent(out), optional  JMAXG,
integer, intent(out), optional  num_tile,
integer, dimension(:), intent(out), optional  tile_id 
)

Return infomation of parent domain (for offline)

Definition at line 1130 of file scale_comm_cartesC_nest.F90.

1130  use scale_prc, only: &
1131  prc_abort
1132  integer, intent(in) :: dom_id
1133 
1134  integer, intent(out), optional :: KMAX
1135  integer, intent(out), optional :: LKMAX
1136  integer, intent(out), optional :: IMAXG
1137  integer, intent(out), optional :: JMAXG
1138  integer, intent(out), optional :: num_tile
1139  integer, intent(out), optional :: tile_id(:)
1140 
1141  integer :: i
1142 
1143  !---------------------------------------------------------------------------
1144 
1145  if ( dom_id < 1 .or. dom_id > num_dom ) then
1146  log_error("COMM_CARTESC_NEST_domina_shape",*) 'domain id is invalid: ', dom_id
1147  call prc_abort
1148  end if
1149 
1150  if ( present(kmax) ) &
1151  kmax = dom_info(dom_id)%KMAX
1152 
1153  if ( present(lkmax) ) &
1154  lkmax = dom_info(dom_id)%LKMAX
1155 
1156  if ( present(imaxg) ) &
1157  imaxg = dom_info(dom_id)%IMAX * dom_info(dom_id)%tile_num_x
1158 
1159  if ( present(jmaxg) ) &
1160  jmaxg = dom_info(dom_id)%JMAX * dom_info(dom_id)%tile_num_y
1161 
1162  if ( present(num_tile) ) &
1163  num_tile = dom_info(dom_id)%tile_num_x * dom_info(dom_id)%tile_num_y
1164 
1165  if ( present(tile_id) ) then
1166  do i = 1, min( size(tile_id), size(dom_info(dom_id)%tile_id) )
1167  tile_id(i) = dom_info(dom_id)%tile_id(i)
1168  end do
1169  end if
1170 
1171  return

References scale_atmos_grid_cartesc_index::imaxg, scale_atmos_grid_cartesc_index::jmaxg, scale_atmos_grid_cartesc_index::kmax, and scale_prc::prc_abort().

Referenced by mod_copytopo::copytopo(), mod_copytopo::copytopo_get_data_scale(), mod_realinput_netcdf::parentatmossetupnetcdf(), mod_realinput_netcdf::parentlandsetupnetcdf(), and mod_realinput_netcdf::parentoceansetupnetcdf().

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)  dom_id,
integer, intent(in)  iloc,
logical, intent(in), optional  xstg,
logical, intent(in), optional  ystg 
)

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

Definition at line 1187 of file scale_comm_cartesC_nest.F90.

1187  use scale_prc, only: &
1188  prc_abort
1189  implicit none
1190 
1191  integer, intent(out) :: tilei, tilej
1192  integer, intent(out) :: cxs, cxe, cys, cye
1193  integer, intent(out) :: pxs, pxe, pys, pye
1194 
1195  integer, intent(in) :: dom_id
1196  integer, intent(in) :: iloc ! rank number; start from 1
1197 
1198  logical, intent(in), optional :: xstg
1199  logical, intent(in), optional :: ystg
1200 
1201  type(domain_info), pointer :: dinfo
1202 
1203  integer :: hdl = 1 ! handler number
1204  integer :: rank
1205  integer :: xloc, yloc
1206  integer :: xlocg, ylocg ! location over whole domain
1207  logical :: xstg_, ystg_
1208  !---------------------------------------------------------------------------
1209 
1210  if ( dom_id < 1 .or. dom_id > num_dom ) then
1211  log_error("COMM_CARTESC_NEST_domina_shape",*) 'domain id is invalid: ', dom_id
1212  call prc_abort
1213  end if
1214 
1215  if ( present(xstg) ) then
1216  xstg_ = xstg
1217  else
1218  xstg_ = .false.
1219  end if
1220  if ( present(ystg) ) then
1221  ystg_ = ystg
1222  else
1223  ystg_ = .false.
1224  end if
1225 
1226  dinfo => dom_info(dom_id)
1227 
1228  rank = dinfo%tile_id(iloc)
1229  xloc = mod( iloc-1, dinfo%tile_num_x ) + 1
1230  yloc = int( real(iloc-1) / real(dinfo%tile_num_x) ) + 1
1231  xlocg = mod( rank, dinfo%prc_num_x ) + 1
1232  ylocg = int( real(rank) / real(dinfo%prc_num_x) ) + 1
1233  tilei = dinfo%IMAX
1234  tilej = dinfo%JMAX
1235 
1236  cxs = tilei * (xloc-1) + 1
1237  cxe = tilei * xloc
1238  cys = tilej * (yloc-1) + 1
1239  cye = tilej * yloc
1240  pxs = 1
1241  pxe = tilei
1242  pys = 1
1243  pye = tilej
1244 
1245  if ( .not. dinfo%periodic_x ) then
1246  if ( xlocg == 1 ) then ! BND_W
1247  tilei = tilei + dinfo%IHALO
1248  pxs = pxs + dinfo%IHALO
1249  pxe = pxe + dinfo%IHALO
1250  endif
1251  if ( xlocg == dinfo%prc_num_x ) then ! BND_E
1252  tilei = tilei + dinfo%IHALO
1253  endif
1254 
1255  if ( xstg_ ) then ! staggarded grid
1256  if ( xlocg == 1 ) then ! BND_W
1257  tilei = tilei + 1
1258  if ( dinfo%IHALO > 0 ) then
1259  pxs = pxs - 1
1260  else
1261  pxe = pxe + 1
1262  end if
1263  else
1264  cxs = cxs + 1
1265  end if
1266  cxe = cxe + 1
1267  end if
1268  end if
1269 
1270  if ( .not. dinfo%periodic_y ) then
1271  if ( ylocg == 1 ) then ! BND_S
1272  tilej = tilej + dinfo%JHALO
1273  pys = pys + dinfo%JHALO
1274  pye = pye + dinfo%JHALO
1275  endif
1276  if ( ylocg == dinfo%prc_num_y ) then ! BND_N
1277  tilej = tilej + dinfo%JHALO
1278  endif
1279 
1280  if ( ystg_ ) then ! staggarded grid
1281  if ( ylocg == 1 ) then ! BND_W
1282  tilej = tilej + 1
1283  if ( dinfo%JHALO > 0 ) then
1284  pys = pys - 1
1285  else
1286  pye = pye + 1
1287  end if
1288  else
1289  cys = cys + 1
1290  end if
1291  cye = cye + 1
1292  end if
1293  end if
1294 
1295  return

References scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_domain_catalogue, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lat, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latuy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latxv, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lon, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonuy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonxv, comm_cartesc_nest_filiation, comm_cartesc_nest_waitall(), scale_comm_cartesc::comm_datatype, scale_comm_cartesc::comm_world, scale_io::h_short, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ihalo, scale_atmos_grid_cartesc_index::imax, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::jhalo, scale_atmos_grid_cartesc_index::jmax, scale_tracer::k, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::khalo, scale_atmos_grid_cartesc_index::kmax, scale_atmos_grid_cartesc_index::ks, scale_atmos_hydrometeor::n_hyd, online_daughter_nprocs, online_domain_num, online_no_rotate, online_parent_dtsec, online_parent_nstep, online_recv_diagqhyd, online_recv_qa, online_send_diagqhyd, online_send_qa, online_use_velz, scale_prc::prc_abort(), scale_prc::prc_intercomm_child, scale_prc::prc_intercomm_parent, scale_prc::prc_ismaster, scale_prc::prc_myrank, scale_prc::prc_nprocs, scale_prc_cartesc::prc_num_x, scale_prc_cartesc::prc_num_y, scale_time::time_dtsec, scale_time::time_nstep, and use_nesting.

Referenced by mod_copytopo::copytopo_get_data_scale(), mod_realinput_netcdf::read2d(), and mod_realinput_netcdf::read3d().

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

◆ comm_cartesc_nest_nestdown_send()

subroutine, public scale_comm_cartesc_nest::comm_cartesc_nest_nestdown_send ( real(rp), dimension(ka,ia,ja), intent(in)  DENS_send,
real(rp), dimension(ka,ia,ja), intent(in)  MOMZ_send,
real(rp), dimension(ka,ia,ja), intent(in)  MOMX_send,
real(rp), dimension(ka,ia,ja), intent(in)  MOMY_send,
real(rp), dimension(ka,ia,ja), intent(in)  RHOT_send,
real(rp), dimension(ka,ia,ja,online_send_qa), intent(in)  QTRC_send 
)

Boundary data transfer from parent to daughter: nestdown (parent side)

Definition at line 2043 of file scale_comm_cartesC_nest.F90.

2043  use scale_prc, only: &
2044  prc_abort, &
2047  use scale_comm_cartesc, only: &
2048  comm_vars8, &
2049  comm_wait
2050  use scale_atmos_grid_cartesc_metric, only: &
2052  implicit none
2053 
2054  real(RP), intent(in) :: DENS_send(KA,IA,JA)
2055  real(RP), intent(in) :: MOMZ_send(KA,IA,JA)
2056  real(RP), intent(in) :: MOMX_send(KA,IA,JA)
2057  real(RP), intent(in) :: MOMY_send(KA,IA,JA)
2058  real(RP), intent(in) :: RHOT_send(KA,IA,JA)
2059  real(RP), intent(in) :: QTRC_send(KA,IA,JA,ONLINE_SEND_QA)
2060 
2061  integer, parameter :: HANDLE = 1
2062 
2063  real(RP) :: WORK1_send(KA,IA,JA)
2064  real(RP) :: WORK2_send(KA,IA,JA)
2065  real(RP) :: u_on_map, v_on_map
2066 
2067  real(RP) :: dummy(1,1,1)
2068  integer :: tagbase, tagcomm
2069  integer :: isu_tag
2070 
2071  integer :: ierr
2072  integer :: i, j, k, iq
2073  !---------------------------------------------------------------------------
2074 
2075  if( .NOT. use_nesting ) return
2076 
2077  tagcomm = intercomm_id(handle) * order_tag_comm
2078 
2079  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2080 
2081  !##### parent [send issue] #####
2082 
2083  call prof_rapstart('NEST_total_P', 2)
2084  call prof_rapstart('NEST_pack_P', 2)
2085 
2086  nsend = nsend + 1
2087  log_info("COMM_CARTESC_NEST_nestdown",'(1X,A,I5,A)') "CONeP[P] send( ", nsend, " )"
2088 
2089  ! to keep values at that time by finish of sending process
2090 !OCL XFILL
2091  org_dens(:,:,:) = dens_send(:,:,:)
2092 !OCL XFILL
2093  org_momz(:,:,:) = momz_send(:,:,:)
2094 !OCL XFILL
2095  org_momx(:,:,:) = momx_send(:,:,:)
2096 !OCL XFILL
2097  org_momy(:,:,:) = momy_send(:,:,:)
2098 !OCL XFILL
2099  org_rhot(:,:,:) = rhot_send(:,:,:)
2100  do iq = 1, online_send_qa
2101 !OCL XFILL
2102  org_qtrc(:,:,:,iq) = qtrc_send(:,:,:,iq)
2103  enddo
2104 
2105  !*** request control
2106  !--- do not change the calling order below;
2107  !--- it should be consistent with the order in "COMM_CARTESC_NEST_recvwait_issue"
2108  rq_ctl_p = 0
2109 
2110  if ( .NOT. online_daughter_no_rotate ) then
2111  ! from staggered point to scalar point
2112  do j = 1, ja
2113  do i = 2, ia
2114  do k = 1, ka
2115  work1_send(k,i,j) = ( org_momx(k,i-1,j) + org_momx(k,i,j) ) * 0.5_rp
2116  enddo
2117  enddo
2118  enddo
2119 
2120  do j = 1, ja
2121  do k = 1, ka
2122  work1_send(k,1,j) = org_momx(k,1,j)
2123  enddo
2124  enddo
2125 
2126  call comm_vars8( work1_send(:,:,:), 1 )
2127 
2128  do j = 2, ja
2129  do i = 1, ia
2130  do k = 1, ka
2131  work2_send(k,i,j) = ( org_momy(k,i,j-1) + org_momy(k,i,j) ) * 0.5_rp
2132  enddo
2133  enddo
2134  enddo
2135 
2136  do i = 1, ia
2137  do k = 1, ka
2138  work2_send(k,i,1) = org_momy(k,i,1)
2139  enddo
2140  enddo
2141 
2142  call comm_vars8( work2_send(:,:,:), 2 )
2143 
2144  call comm_wait ( work1_send(:,:,:), 1, .false. )
2145  call comm_wait ( work2_send(:,:,:), 2, .false. )
2146 
2147  ! rotation from map-projected field to latlon field
2148  do j = 1, ja
2149  do i = 1, ia
2150  do k = 1, ka
2151  u_on_map = work1_send(k,i,j) / org_dens(k,i,j)
2152  v_on_map = work2_send(k,i,j) / org_dens(k,i,j)
2153 
2154  org_u_ll(k,i,j) = u_on_map * rotc(i,j,1) - v_on_map * rotc(i,j,2)
2155  org_v_ll(k,i,j) = u_on_map * rotc(i,j,2) + v_on_map * rotc(i,j,1)
2156  enddo
2157  enddo
2158  enddo
2159  endif
2160 
2161  tagbase = tagcomm + tag_dens*order_tag_var
2162  call comm_cartesc_nest_intercomm_nestdown( org_dens(:,:,:), & ! [IN]
2163  dummy(:,:,:), & ! [OUT]
2164  tagbase, i_sclr, handle, & ! [IN]
2165  isu_tag, & ! [INOUT]
2166  flag_dens = .true. ) ! [IN]
2167 
2168  tagbase = tagcomm + tag_momz*order_tag_var
2169  if ( online_daughter_use_velz ) then
2170  call comm_cartesc_nest_intercomm_nestdown( org_momz(:,:,:), & ! [IN]
2171  dummy(:,:,:), & ! [OUT]
2172  tagbase, i_zstg, handle, & ! [IN]
2173  isu_tag ) ! [INOUT]
2174  endif
2175 
2176  tagbase = tagcomm + tag_momx*order_tag_var
2177  if ( online_daughter_no_rotate ) then
2178  call comm_cartesc_nest_intercomm_nestdown( org_momx(:,:,:), & ! [IN]
2179  dummy(:,:,:), & ! [OUT]
2180  tagbase, i_xstg, handle, & ! [IN]
2181  isu_tag ) ! [INOUT]
2182  else
2183  call comm_cartesc_nest_intercomm_nestdown( org_u_ll(:,:,:), & ! [IN]
2184  dummy(:,:,:), & ! [OUT]
2185  tagbase, i_sclr, handle, & ! [IN]
2186  isu_tag ) ! [INOUT]
2187  endif
2188 
2189  tagbase = tagcomm + tag_momy*order_tag_var
2190  if ( online_daughter_no_rotate ) then
2191  call comm_cartesc_nest_intercomm_nestdown( org_momy(:,:,:), & ! [IN]
2192  dummy(:,:,:), & ! [OUT]
2193  tagbase, i_ystg, handle, & ! [IN]
2194  isu_tag ) ! [INOUT]
2195  else
2196  call comm_cartesc_nest_intercomm_nestdown( org_v_ll(:,:,:), & ! [IN]
2197  dummy(:,:,:), & ! [OUT]
2198  tagbase, i_sclr, handle, & ! [IN]
2199  isu_tag ) ! [INOUT]
2200  endif
2201 
2202  tagbase = tagcomm + tag_rhot*order_tag_var
2203  call comm_cartesc_nest_intercomm_nestdown( org_rhot(:,:,:), & ! [IN]
2204  dummy(:,:,:), & ! [OUT]
2205  tagbase, i_sclr, handle, & ! [IN]
2206  isu_tag ) ! [INOUT]
2207 
2208  do iq = 1, online_send_qa
2209  tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
2210  call comm_cartesc_nest_intercomm_nestdown( org_qtrc(:,:,:,iq), & ! [IN]
2211  dummy(:,:,:), & ! [OUT]
2212  tagbase, i_sclr, handle, & ! [IN]
2213  isu_tag ) ! [INOUT]
2214  enddo
2215 
2216  rq_tot_p = rq_ctl_p
2217 
2218  call prof_rapend ('NEST_pack_P', 2)
2219  call prof_rapend ('NEST_total_P', 2)
2220 
2221  else
2222  log_error("COMM_CARTESC_NEST_nestdown_send",*) 'internal error'
2223  call prc_abort
2224  endif
2225 
2226  return

References scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_rotc, comm_cartesc_nest_filiation, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ja, scale_tracer::k, scale_atmos_grid_cartesc_index::ka, online_send_qa, scale_prc::prc_abort(), scale_prc::prc_intercomm_child, scale_prc::prc_intercomm_parent, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), and use_nesting.

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_send().

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

◆ comm_cartesc_nest_nestdown_recv()

subroutine, public scale_comm_cartesc_nest::comm_cartesc_nest_nestdown_recv ( real(rp), dimension(ka,ia,ja), intent(out)  DENS_recv,
real(rp), dimension(ka,ia,ja), intent(out)  VELZ_recv,
real(rp), dimension(ka,ia,ja), intent(out)  VELX_recv,
real(rp), dimension(ka,ia,ja), intent(out)  VELY_recv,
real(rp), dimension(ka,ia,ja), intent(out)  POTT_recv,
real(rp), dimension(ka,ia,ja,online_recv_qa), intent(out)  QTRC_recv 
)

Boundary data transfer from parent to daughter: nestdown (daughter side)

Definition at line 2237 of file scale_comm_cartesC_nest.F90.

2237  use scale_prc, only: &
2238  prc_abort, &
2241  use scale_comm_cartesc, only: &
2242  comm_vars8, &
2243  comm_wait
2244  use scale_atmos_grid_cartesc_metric, only: &
2246  implicit none
2247 
2248  real(RP), intent(out) :: DENS_recv(KA,IA,JA)
2249  real(RP), intent(out) :: VELZ_recv(KA,IA,JA)
2250  real(RP), intent(out) :: VELX_recv(KA,IA,JA)
2251  real(RP), intent(out) :: VELY_recv(KA,IA,JA)
2252  real(RP), intent(out) :: POTT_recv(KA,IA,JA)
2253  real(RP), intent(out) :: QTRC_recv(KA,IA,JA,ONLINE_RECV_QA)
2254 
2255  integer, parameter :: HANDLE = 2
2256 
2257  real(RP) :: WORK1_recv(KA,IA,JA)
2258  real(RP) :: WORK2_recv(KA,IA,JA)
2259  real(RP) :: U_ll_recv (KA,IA,JA)
2260  real(RP) :: V_ll_recv (KA,IA,JA)
2261  real(RP) :: u_on_map, v_on_map
2262 
2263  real(RP) :: dummy(1,1,1)
2264  integer :: tagbase, tagcomm
2265  integer :: isu_tag
2266 
2267  integer :: ierr
2268  integer :: i, j, k, iq
2269  !---------------------------------------------------------------------------
2270 
2271  if( .NOT. use_nesting ) return
2272 
2273  tagcomm = intercomm_id(handle) * order_tag_comm
2274 
2275  if( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2276 
2277  !##### child [wait issue] #####
2278 
2279  call prof_rapstart('NEST_total_C', 2)
2280  call prof_rapstart('NEST_wait_C', 2)
2281 
2282  nwait_d = nwait_d + 1
2283  !LOG_INFO("COMM_CARTESC_NEST_nestdown",'(1X,A,I5,A)') "NestIDC [C]: que wait ( ", nwait_d, " )"
2284 
2285  !*** reset issue tag and request control
2286  !--- do not change the calling order below;
2287  !--- it should be consistent with the order in "COMM_CARTESC_NEST_recvwait_issue"
2288  isu_tag = 0
2289 
2290  call comm_cartesc_nest_waitall( rq_tot_d, ireq_d )
2291 
2292  if ( online_aggressive_comm ) then
2293  ! nothing to do
2294  else
2295  call mpi_barrier(prc_intercomm_parent, ierr)
2296  endif
2297 
2298  call prof_rapend ('NEST_wait_C', 2)
2299  call prof_rapstart('NEST_unpack_C', 2)
2300 
2301  tagbase = tagcomm + tag_dens*order_tag_var
2302  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2303  work1_recv(:,:,:), & ! [OUT]
2304  tagbase, i_sclr, handle, & ! [IN]
2305  isu_tag, & ! [INOUT]
2306  flag_dens = .true. ) ! [IN]
2307 !OCL XFILL
2308  do j = 1, ja
2309  do i = 1, ia
2310  do k = ks, ke
2311  dens_recv(k,i,j) = work1_recv(k,i,j)
2312  enddo
2313  enddo
2314  enddo
2315 
2316  call comm_vars8( dens_recv, 1 )
2317 
2318  tagbase = tagcomm + tag_momz*order_tag_var
2319  if ( online_use_velz ) then
2320  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2321  work2_recv(:,:,:), & ! [OUT]
2322  tagbase, i_zstg, handle, & ! [IN]
2323  isu_tag ) ! [INOUT]
2324 !OCL XFILL
2325  do j = 1, ja
2326  do i = 1, ia
2327  do k = ks, ke-1
2328  velz_recv(k,i,j) = work2_recv(k,i,j) / ( work1_recv(k,i,j) + work1_recv(k+1,i,j) ) * 2.0_rp
2329  enddo
2330  enddo
2331  enddo
2332 
2333  do j = 1, ja
2334  do i = 1, ia
2335  velz_recv(ks-1,i,j) = 0.0_rp
2336  velz_recv(ke ,i,j) = 0.0_rp
2337  enddo
2338  enddo
2339  endif
2340 
2341  call comm_wait ( dens_recv, 1, .false. )
2342 
2343  tagbase = tagcomm + tag_momx*order_tag_var
2344  if ( online_no_rotate ) then
2345  ! U_ll_recv receives MOMX
2346  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2347  work1_recv(:,:,:), & ! [OUT]
2348  tagbase, i_xstg, handle, & ! [IN]
2349  isu_tag ) ! [INOUT]
2350  else
2351  ! U_ll_recv receives MOMX/DENS
2352  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2353  u_ll_recv(:,:,:), & ! [OUT]
2354  tagbase, i_sclr, handle, & ! [IN]
2355  isu_tag ) ! [INOUT]
2356  endif
2357 
2358  tagbase = tagcomm + tag_momy*order_tag_var
2359  if ( online_no_rotate ) then
2360  ! V_ll_recv receives MOMY
2361  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2362  work2_recv(:,:,:), & ! [OUT]
2363  tagbase, i_ystg, handle, & ! [IN]
2364  isu_tag ) ! [INOUT]
2365  else
2366  ! V_ll_recv receives MOMY/DENS
2367  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2368  v_ll_recv(:,:,:), & ! [OUT]
2369  tagbase, i_sclr, handle, & ! [IN]
2370  isu_tag ) ! [INOUT]
2371  endif
2372 
2373  if ( online_no_rotate ) then
2374 
2375 !OCL XFILL
2376  do j = 1, ja
2377  do i = 1, ia-1
2378  do k = ks, ke
2379  velx_recv(k,i,j) = work1_recv(k,i,j) / ( dens_recv(k,i+1,j) + dens_recv(k,i,j) ) * 2.0_rp
2380  enddo
2381  enddo
2382  enddo
2383 
2384 !OCL XFILL
2385  do j = 1, ja
2386  do k = ks, ke
2387  velx_recv(k,ia,j) = work1_recv(k,ia,j) / dens_recv(k,ia,j)
2388  enddo
2389  enddo
2390 
2391  call comm_vars8( velx_recv, 2 )
2392 
2393 !OCL XFILL
2394  do j = 1, ja-1
2395  do i = 1, ia
2396  do k = ks, ke
2397  vely_recv(k,i,j) = work2_recv(k,i,j) / ( dens_recv(k,i,j+1) + dens_recv(k,i,j) ) * 2.0_rp
2398  enddo
2399  enddo
2400  enddo
2401 
2402 !OCL XFILL
2403  do i = 1, ia
2404  do k = ks, ke
2405  vely_recv(k,i,ja) = work2_recv(k,i,ja) / dens_recv(k,i,ja)
2406  enddo
2407  enddo
2408 
2409  call comm_vars8( vely_recv, 3 )
2410 
2411  call comm_wait ( velx_recv, 2, .false. )
2412  call comm_wait ( vely_recv, 3, .false. )
2413 
2414  else ! rotate
2415 
2416  ! rotation from latlon field to map-projected field
2417 !OCL XFILL
2418  do j = 1, ja
2419  do i = 1, ia
2420  do k = ks, ke
2421  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)
2422  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)
2423  enddo
2424  enddo
2425  enddo
2426 
2427  ! from scalar point to staggered point
2428 !OCL XFILL
2429  do j = 1, ja
2430  do i = 1, ia-1
2431  do k = ks, ke
2432  velx_recv(k,i,j) = ( work1_recv(k,i+1,j) + work1_recv(k,i,j) ) * 0.5_rp
2433  enddo
2434  enddo
2435  enddo
2436 
2437 !OCL XFILL
2438  do j = 1, ja
2439  do k = ks, ke
2440  velx_recv(k,ia,j) = work1_recv(k,ia,j)
2441  enddo
2442  enddo
2443 
2444  call comm_vars8( velx_recv, 2 )
2445 
2446 !OCL XFILL
2447  do j = 1, ja-1
2448  do i = 1, ia
2449  do k = ks, ke
2450  vely_recv(k,i,j) = ( work2_recv(k,i,j+1) + work2_recv(k,i,j) ) * 0.5_rp
2451  enddo
2452  enddo
2453  enddo
2454 
2455 !OCL XFILL
2456  do i = 1, ia
2457  do k = ks, ke
2458  vely_recv(k,i,ja) = work2_recv(k,i,ja)
2459  enddo
2460  enddo
2461 
2462  call comm_vars8( vely_recv, 3 )
2463 
2464  call comm_wait ( velx_recv, 2, .false. )
2465  call comm_wait ( vely_recv, 3, .false. )
2466 
2467  endif
2468 
2469  tagbase = tagcomm + tag_rhot*order_tag_var
2470  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2471  work1_recv(:,:,:), & ! [OUT]
2472  tagbase, i_sclr, handle, & ! [IN]
2473  isu_tag ) ! [INOUT]
2474 !OCL XFILL
2475  do j = 1, ja
2476  do i = 1, ia
2477  do k = ks, ke
2478  pott_recv(k,i,j) = work1_recv(k,i,j) / dens_recv(k,i,j)
2479  enddo
2480  enddo
2481  enddo
2482 
2483  do iq = 1, online_recv_qa
2484  tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
2485  call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), & ! [IN]
2486  work1_recv(:,:,:), & ! [OUT]
2487  tagbase, i_sclr, handle, & ! [IN]
2488  isu_tag ) ! [INOUT]
2489 !OCL XFILL
2490  do j = 1, ja
2491  do i = 1, ia
2492  do k = ks, ke
2493  qtrc_recv(k,i,j,iq) = work1_recv(k,i,j)
2494  enddo
2495  enddo
2496  enddo
2497  enddo
2498 
2499  call prof_rapend ('NEST_unpack_C', 2)
2500  call prof_rapend ('NEST_total_C', 2)
2501 
2502  else
2503  log_error("COMM_CARTESC_NEST_nestdown_recv",*) 'internal error'
2504  call prc_abort
2505  endif
2506 
2507  return

References scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_rotc, comm_cartesc_nest_filiation, comm_cartesc_nest_waitall(), scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ja, scale_tracer::k, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, online_no_rotate, online_recv_qa, online_use_velz, scale_prc::prc_abort(), scale_prc::prc_intercomm_child, scale_prc::prc_intercomm_parent, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), and use_nesting.

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_send().

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

◆ comm_cartesc_nest_recvwait_issue_send()

subroutine, public scale_comm_cartesc_nest::comm_cartesc_nest_recvwait_issue_send

Sub-command for data transfer from parent to daughter: nestdown (parent side)

Definition at line 2513 of file scale_comm_cartesC_nest.F90.

2513  use scale_prc, only: &
2514  prc_abort, &
2516  implicit none
2517 
2518  integer, parameter :: HANDLE = 1
2519 
2520  integer :: isu_tag
2521  integer :: tagbase, tagcomm
2522  integer :: ierr
2523  integer :: iq
2524  !---------------------------------------------------------------------------
2525 
2526  if( .NOT. use_nesting ) return
2527 
2528  tagcomm = intercomm_id(handle) * order_tag_comm
2529 
2530  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2531 
2532  !##### parent [wait issue] #####
2533 
2534  call prof_rapstart('NEST_total_P', 2)
2535  call prof_rapstart('NEST_wait_P', 2)
2536 
2537  nwait_p = nwait_p + 1
2538  !LOG_INFO("COMM_CARTESC_NEST_recvwait_issue",'(1X,A,I5,A)') "NestIDC [P]: que wait ( ", nwait_p, " )"
2539 
2540  call comm_cartesc_nest_issuer_of_wait( handle )
2541 
2542  if ( online_aggressive_comm ) then
2543  ! nothing to do
2544  else
2545  call mpi_barrier(prc_intercomm_child, ierr)
2546  endif
2547 
2548  call prof_rapend ('NEST_wait_P', 2)
2549  call prof_rapend ('NEST_total_P', 2)
2550 
2551  else
2552  log_error("COMM_CARTESC_NEST_recvwait_issue_send",*) 'internal error'
2553  call prc_abort
2554  endif
2555 
2556  return

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

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

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

◆ comm_cartesc_nest_recvwait_issue_recv()

subroutine, public scale_comm_cartesc_nest::comm_cartesc_nest_recvwait_issue_recv

Sub-command for data transfer from parent to daughter: nestdown (daughter side)

Definition at line 2562 of file scale_comm_cartesC_nest.F90.

2562  use scale_prc, only: &
2563  prc_abort, &
2565  implicit none
2566 
2567  integer, parameter :: HANDLE = 2
2568 
2569  integer :: isu_tag
2570  integer :: tagbase, tagcomm
2571  integer :: ierr
2572  integer :: iq
2573  !---------------------------------------------------------------------------
2574 
2575  if( .NOT. use_nesting ) return
2576 
2577  tagcomm = intercomm_id(handle) * order_tag_comm
2578 
2579  if( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2580 
2581  !##### child [receive issue] #####
2582 
2583  call prof_rapstart('NEST_total_C', 2)
2584 
2585  nrecv = nrecv + 1
2586  log_info("COMM_CARTESC_NEST_recvwait_issue_recv",'(1X,A,I5,A)') "NestIDC [C]: que recv ( ", nrecv, " )"
2587 
2588  !*** reset issue tag and request control
2589  !--- do not change the calling order below;
2590  !--- it should be consistent with the order in "COMM_CARTESC_NEST_nestdown"
2591  isu_tag = 0
2592  rq_ctl_d = 0
2593 
2594  tagbase = tagcomm + tag_dens*order_tag_var
2595  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag )
2596 
2597  tagbase = tagcomm + tag_momz*order_tag_var
2598  if ( online_use_velz ) then
2599  call comm_cartesc_nest_issuer_of_receive( tagbase, i_zstg, handle, isu_tag )
2600  endif
2601 
2602  tagbase = tagcomm + tag_momx*order_tag_var
2603  if ( online_no_rotate ) then
2604  call comm_cartesc_nest_issuer_of_receive( tagbase, i_xstg, handle, isu_tag )
2605  else
2606  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag )
2607  endif
2608 
2609  tagbase = tagcomm + tag_momy*order_tag_var
2610  if ( online_no_rotate ) then
2611  call comm_cartesc_nest_issuer_of_receive( tagbase, i_ystg, handle, isu_tag )
2612  else
2613  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag )
2614  endif
2615 
2616  tagbase = tagcomm + tag_rhot*order_tag_var
2617  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag )
2618 
2619  do iq = 1, online_recv_qa
2620  tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
2621  call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag )
2622  enddo
2623 
2624  rq_tot_d = rq_ctl_d
2625 
2626  call prof_rapend('NEST_total_C', 2)
2627 
2628  else
2629  log_error("COMM_CARTESC_NEST_recvwait_issue_recv",*) 'internal error'
2630  call prc_abort
2631  endif
2632 
2633  return

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

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

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

◆ comm_cartesc_nest_recv_cancel_send()

subroutine, public scale_comm_cartesc_nest::comm_cartesc_nest_recv_cancel_send

Sub-command for data transfer from parent to daughter: nestdown (parent side)

Definition at line 2639 of file scale_comm_cartesC_nest.F90.

2639  use scale_prc, only: &
2640  prc_abort
2641  implicit none
2642 
2643  integer, parameter :: HANDLE = 1
2644 
2645  !logical :: flag
2646  !integer :: istatus(MPI_STATUS_SIZE)
2647 
2648  integer :: rq
2649  integer :: ierr
2650  !---------------------------------------------------------------------------
2651 
2652  if( .NOT. use_nesting ) return
2653 
2654  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2655 
2656  !##### parent #####
2657  ! Nothing to do
2658 
2659  else
2660  log_error("COMM_CARTESC_NEST_recv_cancel_send",*) 'internal error'
2661  call prc_abort
2662  endif
2663 
2664  return

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

Here is the call graph for this function:

◆ comm_cartesc_nest_recv_cancel_recv()

subroutine, public scale_comm_cartesc_nest::comm_cartesc_nest_recv_cancel_recv

Sub-command for data transfer from parent to daughter: nestdown (daughter side)

Definition at line 2670 of file scale_comm_cartesC_nest.F90.

2670  use scale_prc, only: &
2671  prc_abort
2672  implicit none
2673 
2674  integer, parameter :: HANDLE = 2
2675 
2676  !logical :: flag
2677  !integer :: istatus(MPI_STATUS_SIZE)
2678 
2679  integer :: rq
2680  integer :: ierr
2681  !---------------------------------------------------------------------------
2682 
2683  if( .NOT. use_nesting ) return
2684 
2685  if( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2686 
2687  !##### child #####
2688 
2689  log_info("COMM_CARTESC_NEST_recv_cancel_recv",'(1X,A,I5,A)') "NestIDC [C]: CANCEL recv ( ", nrecv, " )"
2690 
2691  do rq = 1, rq_tot_d
2692  if ( ireq_d(rq) /= mpi_request_null ) then
2693 
2694  call mpi_cancel(ireq_d(rq), ierr)
2695 
2696 ! call MPI_TEST_CANCELLED(istatus, flag, ierr)
2697 ! if ( .NOT. flag ) then
2698 ! LOG_ERROR("COMM_CARTESC_NEST_recv_cancel_recv",*) 'receive actions do not cancelled, req = ', rq
2699 ! endif
2700  endif
2701  enddo
2702 
2703  else
2704  log_error("COMM_CARTESC_NEST_recv_cancel_recv",*) 'internal error'
2705  call prc_abort
2706  endif
2707 
2708  return

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

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_finalize().

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

◆ comm_cartesc_nest_intercomm_nestdown_3d()

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

Inter-communication from parent to daughter: nestdown.

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

Definition at line 2721 of file scale_comm_cartesC_nest.F90.

2721  use scale_prc, only: &
2722  prc_abort, &
2725  use scale_comm_cartesc, only: &
2727  use scale_interp, only: &
2730  real_cz => atmos_grid_cartesc_real_cz, &
2731  real_fz => atmos_grid_cartesc_real_fz
2732  implicit none
2733 
2734  real(RP), intent(in) :: pvar(:,:,:)
2735  real(RP), intent(out) :: dvar(:,:,:)
2736  integer, intent(in) :: tagbase
2737  integer, intent(in) :: id_stag
2738  integer, intent(in) :: HANDLE
2739  integer, intent(inout) :: isu_tag
2740 
2741  logical , intent(in), optional :: flag_dens
2742 
2743  integer :: tile_num_x
2744 
2745  integer :: ileng, tag, target_rank
2746 
2747  integer :: parent_KA
2748 
2749  integer :: xloc, yloc
2750  integer :: gxs, gxe, gys, gye ! for large domain
2751  integer :: pxs, pxe, pys, pye ! for parent domain
2752  integer :: zs, ze
2753 
2754  integer :: ig, rq, yp
2755  logical :: no_zstag
2756  logical :: spline
2757  logical :: logarithmic
2758 
2759  integer :: ierr
2760  integer :: i, j
2761  !---------------------------------------------------------------------------
2762 
2763  if( .NOT. use_nesting ) return
2764 
2765  logarithmic = .false.
2766  spline = .false.
2767  if ( present(flag_dens) ) then
2768  if( flag_dens ) then
2769  logarithmic = .true.
2770  spline = .true.
2771  end if
2772  endif
2773 
2774  if ( id_stag == i_sclr ) then
2775  no_zstag = .true.
2776  ig = i_sclr
2777  elseif( id_stag == i_zstg ) then
2778  no_zstag = .false.
2779  ig = i_zstg
2780  elseif( id_stag == i_xstg ) then
2781  no_zstag = .true.
2782  ig = i_xstg
2783  elseif( id_stag == i_ystg ) then
2784  no_zstag = .true.
2785  ig = i_ystg
2786  endif
2787 
2788  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2789  !##### parent #####
2790 
2791  ileng = ka * ia * ja
2792 
2793  rq = rq_ctl_p
2794 
2795  do yp = 1, num_yp
2796  rq = rq + 1
2797 
2798  ! send data to multiple daughter processes
2799  target_rank = comm_cartesc_nest_tile_list_yp(yp)
2800  tag = tagbase + yp
2801 
2802  call mpi_isend( pvar, &
2803  ileng, &
2804  comm_datatype, &
2805  target_rank, &
2806  tag, &
2808  ireq_p(rq), &
2809  ierr )
2810 
2811  dvar(:,:,:) = -1.0_rp ! input as a dummy value
2812  enddo
2813 
2814  rq_ctl_p = rq
2815 
2816  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2817 
2818  !##### child #####
2819 
2820  parent_ka = dom_info(i_parent)%KMAX + dom_info(i_parent)%KHALO * 2
2821  ileng = parent_ka &
2822  * ( dom_info(i_parent)%IMAX + dom_info(i_parent)%IHALO * 2 ) &
2823  * ( dom_info(i_parent)%JMAX + dom_info(i_parent)%JHALO * 2 )
2824 
2825  tile_num_x = dom_info(i_parent)%tile_num_x
2826 
2827  zs = 1
2828  ze = parent_ka
2829 
2830  pxs = dom_info(i_parent)%IHALO + 1
2831  pxe = dom_info(i_parent)%IMAX + dom_info(i_parent)%IHALO
2832  pys = dom_info(i_parent)%JHALO + 1
2833  pye = dom_info(i_parent)%JMAX + dom_info(i_parent)%JHALO
2834 
2835  rq = rq_ctl_d
2836 
2837  do yp = 1, comm_cartesc_nest_tile_all
2838  rq = rq + 1
2839 
2840  xloc = mod( yp-1, dom_info(i_parent)%TILE_NUM_X ) + 1
2841  yloc = int( real(yp-1) / real(dom_info(i_parent)%TILE_NUM_X) ) + 1
2842 
2843  gxs = dom_info(i_parent)%IMAX * (xloc-1) + 1
2844  gxe = dom_info(i_parent)%IMAX * xloc
2845  gys = dom_info(i_parent)%JMAX * (yloc-1) + 1
2846  gye = dom_info(i_parent)%JMAX * yloc
2847 
2848  isu_tag = isu_tag + 1
2849 
2850  if ( isu_tag > max_isu ) then
2851  log_error("COMM_CARTESC_NEST_intercomm_nestdown_3D",*) 'Exceeded maximum issue'
2852  log_error_cont(*) 'isu_tag = ', isu_tag
2853  call prc_abort
2854  endif
2855 
2856 !OCL XFILL
2857  buffer_ref_3d(zs:ze,gxs:gxe,gys:gye) = recvbuf_3d(zs:ze,pxs:pxe,pys:pye,isu_tag)
2858 
2859  enddo
2860 
2861  rq_ctl_d = rq
2862 
2863  if ( no_zstag ) then
2864  call interp_interp3d( itp_nh, &
2865  tileal_ka, khalo+1, tileal_ka-khalo, &
2866  tileal_ia, tileal_ja, &
2867  ka, ks, ke, ia, ja, &
2868  igrd( :,:,:,ig), & ! [IN]
2869  jgrd( :,:,:,ig), & ! [IN]
2870  hfact( :,:,:,ig), & ! [IN]
2871  kgrd(:,:,:,:,:,ig), & ! [IN]
2872  vfact(:, :,:,:,ig), & ! [IN]
2873  buffer_ref_cz(:,:,:), & ! [IN]
2874  real_cz(:,:,:), & ! [IN]
2875  buffer_ref_3d(:,:,:), & ! [IN]
2876  dvar(:,:,:), & ! [OUT]
2877  spline = spline, & ! [IN, optional]
2878  logwgt = logarithmic ) ! [IN, optional]
2879 
2880  else
2881  call interp_interp3d( itp_nh, &
2882  tileal_ka, khalo, tileal_ka-khalo, &
2883  tileal_ia, tileal_ja, &
2884  ka, ks, ke, ia, ja, &
2885  igrd( :,:,:,ig), & ! [IN]
2886  jgrd( :,:,:,ig), & ! [IN]
2887  hfact( :,:,:,ig), & ! [IN]
2888  kgrd(:,:,:,:,:,ig), & ! [IN]
2889  vfact(:, :,:,:,ig), & ! [IN]
2890  buffer_ref_fz(:,:,:), & ! [IN]
2891  real_fz(1:,:,:), & ! [IN]
2892  buffer_ref_3d(:,:,:), & ! [IN]
2893  dvar(:,:,:), & ! [OUT]
2894  spline = spline, & ! [IN, optional]
2895  logwgt = logarithmic ) ! [IN, optional]
2896  endif
2897 
2898  do j = 1, ja
2899  do i = 1, ia
2900  dvar( 1:ks-1,i,j) = 0.0_rp
2901  dvar(ke+1:ka ,i,j) = 0.0_rp
2902  enddo
2903  enddo
2904 
2905  else
2906  log_error("COMM_CARTESC_NEST_intercomm_nestdown_3D",*) 'internal error'
2907  call prc_abort
2908  endif
2909 
2910  return

References scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz, comm_cartesc_nest_filiation, scale_comm_cartesc::comm_datatype, scale_atmos_grid_cartesc_index::ia, scale_interp::interp_interp3d(), scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::khalo, scale_atmos_grid_cartesc_index::ks, scale_prc::prc_abort(), scale_prc::prc_intercomm_child, scale_prc::prc_intercomm_parent, and use_nesting.

Here is the call graph for this function:

◆ comm_cartesc_nest_issuer_of_receive_3d()

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

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

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

Definition at line 2920 of file scale_comm_cartesC_nest.F90.

2920  use scale_prc, only: &
2921  prc_myrank, &
2922  prc_abort, &
2924  use scale_comm_cartesc, only: &
2926  implicit none
2927 
2928  integer, intent(in) :: tagbase
2929  integer, intent(in) :: id_stag
2930  integer, intent(in) :: HANDLE
2931  integer, intent(inout) :: isu_tag
2932 
2933  integer :: ierr, ileng
2934  integer :: tag, target_rank
2935 
2936  integer :: rq, yp
2937  !---------------------------------------------------------------------------
2938 
2939  if( .NOT. use_nesting ) return
2940 
2941  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
2942 
2943  !##### parent #####
2944  ! nothing to do
2945 
2946  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
2947 
2948  !##### child #####
2949 
2950  ileng = ( dom_info(i_parent)%KMAX + dom_info(i_parent)%KHALO * 2 ) &
2951  * ( dom_info(i_parent)%IMAX + dom_info(i_parent)%IHALO * 2 ) &
2952  * ( dom_info(i_parent)%JMAX + dom_info(i_parent)%JHALO * 2 )
2953 
2954  rq = rq_ctl_d
2955 
2956  do yp = 1, comm_cartesc_nest_tile_all
2957  rq = rq + 1
2958 
2959  target_rank = comm_cartesc_nest_tile_list_d(yp,prc_myrank+1)
2960  tag = tagbase + call_order(yp)
2961 
2962  isu_tag = isu_tag + 1
2963 
2964  if ( isu_tag > max_isu ) then
2965  log_error("COMM_CARTESC_NEST_issuer_of_receive_3D",*) 'Exceeded maximum issue'
2966  log_error_cont(*) 'isu_tag = ', isu_tag
2967  call prc_abort
2968  endif
2969 
2970  recvbuf_3d(:,:,:,isu_tag) = 0.0_rp
2971 
2972  call mpi_irecv( recvbuf_3d(:,:,:,isu_tag), &
2973  ileng, &
2974  comm_datatype, &
2975  target_rank, &
2976  tag, &
2978  ireq_d(rq), &
2979  ierr )
2980 
2981  enddo
2982 
2983  rq_ctl_d = rq
2984 
2985  else
2986  log_error("COMM_CARTESC_NEST_issuer_of_receive_3D",*) 'internal error'
2987  call prc_abort
2988  endif
2989 
2990  return

References comm_cartesc_nest_filiation, scale_comm_cartesc::comm_datatype, scale_prc::prc_abort(), scale_prc::prc_intercomm_parent, scale_prc::prc_myrank, and use_nesting.

Here is the call graph for this function:

◆ comm_cartesc_nest_issuer_of_wait_3d()

subroutine scale_comm_cartesc_nest::comm_cartesc_nest_issuer_of_wait_3d ( integer, intent(in)  HANDLE)

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

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

Definition at line 2997 of file scale_comm_cartesC_nest.F90.

2997  use scale_prc, only: &
2998  prc_abort
2999  implicit none
3000 
3001  integer, intent(in) :: HANDLE
3002  !---------------------------------------------------------------------------
3003 
3004  if( .NOT. use_nesting ) return
3005 
3006  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
3007 
3008  !##### parent #####
3009  call comm_cartesc_nest_waitall( rq_tot_p, ireq_p(:) )
3010 
3011  elseif( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
3012 
3013  !##### child #####
3014  ! nothing to do
3015 
3016  else
3017  log_error("COMM_CARTESC_NEST_issuer_of_wait_3D",*) 'internal error'
3018  call prc_abort
3019  endif
3020 
3021  return

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

Here is the call graph for this function:

◆ comm_cartesc_nest_waitall()

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

[substance of comm_wait] Inter-communication

Definition at line 3029 of file scale_comm_cartesC_nest.F90.

3029  use scale_prc, only: &
3030  prc_abort
3031  implicit none
3032 
3033  integer, intent(in) :: req_count
3034  integer, intent(inout) :: ireq(max_rq)
3035 
3036  integer :: i
3037  integer :: ierr
3038  integer :: istatus(MPI_STATUS_SIZE,req_count)
3039  integer :: req_count2
3040  integer :: ireq2(max_rq)
3041 
3042 ! logical :: flag = .false.
3043 ! integer(8) :: num = 0
3044  !---------------------------------------------------------------------------
3045 
3046  if( .NOT. use_nesting ) return
3047 
3048  req_count2 = 0
3049  do i = 1, req_count
3050  if ( ireq(i) /= mpi_request_null ) then
3051  req_count2 = req_count2 + 1
3052  ireq2(req_count2) = ireq(i)
3053  endif
3054  enddo
3055 
3056  if( req_count2 /= 0 ) call mpi_waitall( req_count2, ireq2(1:req_count2), istatus, ierr )
3057 
3058 ! do while ( .NOT. flag )
3059 ! num = num + 1
3060 ! call MPI_TESTALL( req_count, ireq, flag, istatus, ierr )
3061 !
3062 ! if ( num > ONLINE_WAIT_LIMIT ) then
3063 ! LOG_ERROR("COMM_CARTESC_NEST_waitall",'(1x,A)') 'over the limit of waiting time [NESTCOM]'
3064 ! LOG_ERROR_CONT('(1x,A)') 'over the limit of waiting time [NESTCOM]'
3065 ! call PRC_abort
3066 ! endif
3067 ! enddo
3068 
3069  return

References scale_prc::prc_abort(), and use_nesting.

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

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

◆ comm_cartesc_nest_test_send()

subroutine, public scale_comm_cartesc_nest::comm_cartesc_nest_test_send

[check communication status] Inter-communication (parent side)

Definition at line 3075 of file scale_comm_cartesC_nest.F90.

3075  use scale_prc, only: &
3076  prc_abort
3077  implicit none
3078 
3079  integer, parameter :: HANDLE = 1
3080 
3081  integer :: istatus(MPI_STATUS_SIZE)
3082  integer :: ierr
3083  logical :: flag
3084  !---------------------------------------------------------------------------
3085 
3086  if( .NOT. use_nesting ) return
3087 
3088  if ( comm_cartesc_nest_filiation( intercomm_id(handle) ) > 0 ) then
3089 
3090  !##### parent #####
3091 
3092  call prof_rapstart('NEST_test_P', 2)
3093  if ( rq_ctl_p > 0 ) call mpi_test(ireq_p(1), flag, istatus, ierr)
3094  call prof_rapend('NEST_test_P', 2)
3095 
3096  else
3097  log_error("COMM_CARTESC_NEST_test_send",*) 'error'
3098  call prc_abort
3099  endif
3100 
3101  return

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

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_send().

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

◆ comm_cartesc_nest_test_recv()

subroutine, public scale_comm_cartesc_nest::comm_cartesc_nest_test_recv

[check communication status] Inter-communication (daughter side)

Definition at line 3107 of file scale_comm_cartesC_nest.F90.

3107  use scale_prc, only: &
3108  prc_abort
3109  implicit none
3110 
3111  integer, parameter :: HANDLE = 2
3112 
3113  integer :: istatus(MPI_STATUS_SIZE)
3114  integer :: ierr
3115  logical :: flag
3116  !---------------------------------------------------------------------------
3117 
3118  if( .NOT. use_nesting ) return
3119 
3120  if( comm_cartesc_nest_filiation( intercomm_id(handle) ) < 0 ) then
3121 
3122  !##### child #####
3123 
3124  call prof_rapstart('NEST_test_C', 2)
3125  if ( rq_ctl_d > 0 ) call mpi_test(ireq_d(1), flag, istatus, ierr)
3126  call prof_rapend('NEST_test_C', 2)
3127 
3128  else
3129  log_error("COMM_CARTESC_NEST_test_recv",*) 'error'
3130  call prc_abort
3131  endif
3132 
3133  return

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

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_update().

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

◆ comm_cartesc_nest_finalize()

subroutine, public scale_comm_cartesc_nest::comm_cartesc_nest_finalize

finalize

Definition at line 3139 of file scale_comm_cartesC_nest.F90.

3139  implicit none
3140 
3141  integer :: i
3142 
3143  do i = 1, num_dom
3144  if ( allocated( dom_info(i)%latlon_catalogue ) ) deallocate( dom_info(i)%latlon_catalogue )
3145  if ( allocated( dom_info(i)%tile_id ) ) deallocate( dom_info(i)%tile_id )
3146  num_dom = 0
3147  end do
3148 
3149  if ( allocated( comm_cartesc_nest_tile_list_p ) ) deallocate( comm_cartesc_nest_tile_list_p )
3150  if ( allocated( comm_cartesc_nest_tile_list_d ) ) deallocate( comm_cartesc_nest_tile_list_d )
3151  if ( allocated( comm_cartesc_nest_tile_list_yp ) ) deallocate( comm_cartesc_nest_tile_list_yp )
3152 
3153  if ( allocated( ireq_p ) ) deallocate( ireq_p )
3154  if ( allocated( ireq_d ) ) deallocate( ireq_d )
3155 
3156  if ( allocated( call_order ) ) deallocate( call_order )
3157  if ( allocated( recvbuf_3d ) ) deallocate( recvbuf_3d )
3158 
3159  if ( allocated( buffer_ref_lon ) ) deallocate( buffer_ref_lon )
3160  if ( allocated( buffer_ref_lonuy ) ) deallocate( buffer_ref_lonuy )
3161  if ( allocated( buffer_ref_lonxv ) ) deallocate( buffer_ref_lonxv )
3162  if ( allocated( buffer_ref_lat ) ) deallocate( buffer_ref_lat )
3163  if ( allocated( buffer_ref_latuy ) ) deallocate( buffer_ref_latuy )
3164  if ( allocated( buffer_ref_latxv ) ) deallocate( buffer_ref_latxv )
3165  if ( allocated( buffer_ref_cz ) ) deallocate( buffer_ref_cz )
3166  if ( allocated( buffer_ref_fz ) ) deallocate( buffer_ref_fz )
3167  if ( allocated( buffer_ref_3d ) ) deallocate( buffer_ref_3d )
3168 
3169  if ( allocated( org_dens ) ) deallocate( org_dens )
3170  if ( allocated( org_momz ) ) deallocate( org_momz )
3171  if ( allocated( org_momx ) ) deallocate( org_momx )
3172  if ( allocated( org_momy ) ) deallocate( org_momy )
3173  if ( allocated( org_u_ll ) ) deallocate( org_u_ll )
3174  if ( allocated( org_v_ll ) ) deallocate( org_v_ll )
3175  if ( allocated( org_rhot ) ) deallocate( org_rhot )
3176  if ( allocated( org_qtrc ) ) deallocate( org_qtrc )
3177 
3178  if ( allocated( igrd ) ) deallocate( igrd )
3179  if ( allocated( jgrd ) ) deallocate( jgrd )
3180  if ( allocated( hfact ) ) deallocate( hfact )
3181  if ( allocated( kgrd ) ) deallocate( kgrd )
3182  if ( allocated( vfact ) ) deallocate( vfact )
3183 
3184  return

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

Here is the caller graph for this function:

Variable Documentation

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

74  integer, public :: HANDLING_NUM

Referenced by comm_cartesc_nest_setup().

◆ comm_cartesc_nest_interp_level

integer, public scale_comm_cartesc_nest::comm_cartesc_nest_interp_level = 5

horizontal interpolation level

Definition at line 76 of file scale_comm_cartesC_nest.F90.

76  integer, public :: COMM_CARTESC_NEST_INTERP_LEVEL = 5

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

◆ comm_cartesc_nest_interp_weight_order

integer, public scale_comm_cartesc_nest::comm_cartesc_nest_interp_weight_order = 2

horizontal interpolation weight order

Definition at line 77 of file scale_comm_cartesC_nest.F90.

77  integer, public :: COMM_CARTESC_NEST_INTERP_WEIGHT_ORDER = 2

Referenced by comm_cartesc_nest_setup().

◆ use_nesting

logical, public scale_comm_cartesc_nest::use_nesting = .false.

◆ online_iam_parent

logical, public scale_comm_cartesc_nest::online_iam_parent = .false.

a flag to say "I am a parent"

Definition at line 80 of file scale_comm_cartesC_nest.F90.

80  logical, public :: ONLINE_IAM_PARENT = .false.

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

◆ online_iam_daughter

logical, public scale_comm_cartesc_nest::online_iam_daughter = .false.

a flag to say "I am a daughter"

Definition at line 81 of file scale_comm_cartesC_nest.F90.

81  logical, public :: ONLINE_IAM_DAUGHTER = .false.

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

◆ online_domain_num

integer, public scale_comm_cartesc_nest::online_domain_num = 1

Definition at line 82 of file scale_comm_cartesC_nest.F90.

82  integer, public :: ONLINE_DOMAIN_NUM = 1

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

◆ online_use_velz

logical, public scale_comm_cartesc_nest::online_use_velz = .false.

◆ online_no_rotate

logical, public scale_comm_cartesc_nest::online_no_rotate = .false.

Definition at line 84 of file scale_comm_cartesC_nest.F90.

84  logical, public :: ONLINE_NO_ROTATE = .false.

Referenced by comm_cartesc_nest_domain_shape(), comm_cartesc_nest_nestdown_recv(), comm_cartesc_nest_recvwait_issue_recv(), and comm_cartesc_nest_setup().

◆ online_boundary_use_qhyd

logical, public scale_comm_cartesc_nest::online_boundary_use_qhyd = .false.

Definition at line 85 of file scale_comm_cartesC_nest.F90.

85  logical, public :: ONLINE_BOUNDARY_USE_QHYD = .false.

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

◆ online_recv_diagqhyd

logical, public scale_comm_cartesc_nest::online_recv_diagqhyd = .false.

Definition at line 87 of file scale_comm_cartesC_nest.F90.

87  logical, public :: ONLINE_RECV_DIAGQHYD = .false.

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

◆ online_send_diagqhyd

logical, public scale_comm_cartesc_nest::online_send_diagqhyd = .false.

Definition at line 88 of file scale_comm_cartesC_nest.F90.

88  logical, public :: ONLINE_SEND_DIAGQHYD = .false.

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

◆ online_recv_qa

integer, public scale_comm_cartesc_nest::online_recv_qa = 0

◆ online_send_qa

integer, public scale_comm_cartesc_nest::online_send_qa = 0

number of tracer sent to the daughter domain

Definition at line 90 of file scale_comm_cartesC_nest.F90.

90  integer, public :: ONLINE_SEND_QA = 0

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

◆ online_parent_dtsec

real(dp), public scale_comm_cartesc_nest::online_parent_dtsec

◆ online_parent_nstep

integer, public scale_comm_cartesc_nest::online_parent_nstep

◆ online_daughter_nprocs

integer, public scale_comm_cartesc_nest::online_daughter_nprocs = -1

Definition at line 95 of file scale_comm_cartesC_nest.F90.

95  integer, public :: ONLINE_DAUGHTER_nprocs = -1

Referenced by comm_cartesc_nest_domain_shape(), and comm_cartesc_nest_setup().

scale_comm_cartesc::comm_datatype
integer, public comm_datatype
datatype of variable
Definition: scale_comm_cartesC.F90:105
scale_atmos_grid_cartesc_index::ke
integer, public ke
end point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:52
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
mod_atmos_phy_mp_vars::qa_mp
integer, public qa_mp
Definition: mod_atmos_phy_mp_vars.F90:78
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_cz
geopotential height [m] (zxy)
Definition: scale_atmos_grid_cartesC_real.F90:39
scale_atmos_grid_cartesc_index::ihalo
integer, public ihalo
Definition: scale_atmos_grid_cartesC_index.F90:44
scale_file::file_open
subroutine, public file_open(basename, fid, mode, single, allnodes, aggregate, rankid, postfix)
Definition: scale_file.F90:536
scale_atmos_grid_cartesc_index::ka
integer, public ka
Definition: scale_atmos_grid_cartesC_index.F90:47
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonxv
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lonxv
longitude at staggered point (xv) [rad,0-2pi]
Definition: scale_atmos_grid_cartesC_real.F90:51
scale_interp
module INTERPOLATION
Definition: scale_interp.F90:12
scale_atmos_grid_cartesc_metric
module Atmosphere Grid CartesianC metirc
Definition: scale_atmos_grid_cartesC_metric.F90:12
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_atmos_grid_cartesc::atmos_grid_cartesc_fx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fx
face coordinate [m]: x, local
Definition: scale_atmos_grid_cartesC.F90:58
scale_prc::prc_myrank
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:91
scale_atmos_grid_cartesc_index::khalo
integer, parameter, public khalo
Definition: scale_atmos_grid_cartesC_index.F90:43
scale_atmos_hydrometeor::atmos_hydrometeor_dry
logical, public atmos_hydrometeor_dry
Definition: scale_atmos_hydrometeor.F90:114
scale_atmos_grid_cartesc_index::imaxg
integer, public imaxg
Definition: scale_atmos_grid_cartesC_index.F90:72
scale_atmos_grid_cartesc_real
module Atmosphere GRID CartesC Real(real space)
Definition: scale_atmos_grid_cartesC_real.F90:11
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lon
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lon
longitude [rad,0-2pi]
Definition: scale_atmos_grid_cartesC_real.F90:49
scale_file
module file
Definition: scale_file.F90:15
scale_atmos_grid_cartesc_index::jmaxg
integer, public jmaxg
Definition: scale_atmos_grid_cartesC_index.F90:73
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_mapprojection
module Map projection
Definition: scale_mapprojection.F90:12
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latuv
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latuv
latitude at staggered point (uv) [rad,-pi,pi]
Definition: scale_atmos_grid_cartesC_real.F90:56
scale_atmos_grid_cartesc_index::kmax
integer, public kmax
Definition: scale_atmos_grid_cartesC_index.F90:36
scale_atmos_grid_cartesc::atmos_grid_cartesc_fy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fy
face coordinate [m]: y, local
Definition: scale_atmos_grid_cartesC.F90:59
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latuy
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latuy
latitude at staggered point (uy) [rad,-pi,pi]
Definition: scale_atmos_grid_cartesC_real.F90:54
scale_atmos_grid_cartesc_index::jhalo
integer, public jhalo
Definition: scale_atmos_grid_cartesC_index.F90:45
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
scale_time
module TIME
Definition: scale_time.F90:11
scale_interp::interp_interp3d
subroutine, public interp_interp3d(npoints, KA_ref, KS_ref, KE_ref, IA_ref, JA_ref, KA, KS, KE, IA, JA, idx_i, idx_j, hfact, idx_k, vfact, hgt_ref, hgt, val_ref, val, spline, logwgt, threshold_undef, wsum, val2)
Definition: scale_interp.F90:1470
scale_prc::prc_global_domainid
integer, public prc_global_domainid
my domain ID in global communicator
Definition: scale_prc.F90:85
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_latxv
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latxv
latitude at staggered point (xv) [rad,-pi,pi]
Definition: scale_atmos_grid_cartesC_real.F90:55
scale_atmos_grid_cartesc_index::ks
integer, public ks
start point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:51
scale_prc::prc_intercomm_child
integer, public prc_intercomm_child
communicator between this rank and child domain
Definition: scale_prc.F90:97
scale_prc::prc_intercomm_parent
integer, public prc_intercomm_parent
communicator between this rank and parent domain
Definition: scale_prc.F90:96
scale_time::time_dtsec
real(dp), public time_dtsec
time interval of model [sec]
Definition: scale_time.F90:33
scale_comm_cartesc
module COMMUNICATION
Definition: scale_comm_cartesC.F90:11
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonuv
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lonuv
longitude at staggered point (uv) [rad,0-2pi]
Definition: scale_atmos_grid_cartesC_real.F90:52
scale_atmos_grid_cartesc_index::js
integer, public js
start point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:55
scale_atmos_grid_cartesc::atmos_grid_cartesc_cy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cy
center coordinate [m]: y, local
Definition: scale_atmos_grid_cartesC.F90:57
scale_interp::interp_setup
subroutine, public interp_setup(weight_order, search_limit)
Setup.
Definition: scale_interp.F90:90
scale_const::const_d2r
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:33
scale_time::time_nstep
integer, public time_nstep
total steps [number]
Definition: scale_time.F90:74
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lat
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lat
latitude [rad,-pi,pi]
Definition: scale_atmos_grid_cartesC_real.F90:53
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_fz
geopotential height [m] (wxy)
Definition: scale_atmos_grid_cartesC_real.F90:43
scale_atmos_grid_cartesc
module atmosphere / grid / cartesC
Definition: scale_atmos_grid_cartesC.F90:12
scale_atmos_grid_cartesc::atmos_grid_cartesc_cx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cx
center coordinate [m]: x, local
Definition: scale_atmos_grid_cartesC.F90:56
scale_prc::prc_ismaster
logical, public prc_ismaster
master process in local communicator?
Definition: scale_prc.F90:92
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lonuy
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lonuy
longitude at staggered point (uy) [rad,0-2pi]
Definition: scale_atmos_grid_cartesC_real.F90:50
scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_rotc
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_metric_rotc
rotation coefficient
Definition: scale_atmos_grid_cartesC_metric.F90:36