65 real(
rp),
allocatable :: latlon_catalogue(:,:,:)
66 integer,
allocatable :: tile_id(:)
69 character(len=FILE_HLONG) :: basename
101 private :: comm_cartesc_nest_ping
102 private :: comm_cartesc_nest_parentsize
103 private :: comm_cartesc_nest_catalogue
104 private :: comm_cartesc_nest_setup_nestdown
105 private :: comm_cartesc_nest_importgrid_nestdown
106 private :: comm_cartesc_nest_intercomm_nestdown
107 private :: comm_cartesc_nest_issuer_of_receive
108 private :: comm_cartesc_nest_issuer_of_wait
110 interface comm_cartesc_nest_intercomm_nestdown
112 end interface comm_cartesc_nest_intercomm_nestdown
114 interface comm_cartesc_nest_issuer_of_receive
116 end interface comm_cartesc_nest_issuer_of_receive
118 interface comm_cartesc_nest_issuer_of_wait
120 end interface comm_cartesc_nest_issuer_of_wait
126 integer,
private,
parameter :: max_dinfo = 3
127 integer,
private :: num_dom = 0
128 type(
domain_info),
private,
target :: dom_info(max_dinfo)
129 integer,
private :: i_parent = -1
131 real(
rp),
private :: latlon_local (4,2)
133 integer,
private :: comm_cartesc_nest_tile_all
134 integer,
private :: comm_cartesc_nest_tile_allmax_p
135 integer,
private :: comm_cartesc_nest_tile_allmax_d
136 integer,
private,
allocatable :: comm_cartesc_nest_tile_list_p(:,:)
137 integer,
private,
allocatable :: comm_cartesc_nest_tile_list_d(:,:)
138 integer,
private,
allocatable :: comm_cartesc_nest_tile_list_yp(:)
139 integer,
private :: num_yp
141 integer(8),
private :: online_wait_limit
142 logical,
private :: online_daughter_use_velz
143 logical,
private :: online_daughter_no_rotate
144 logical,
private :: online_aggressive_comm
146 integer,
private :: tileal_ka
147 integer,
private :: tileal_ia
148 integer,
private :: tileal_ja
150 integer,
parameter :: i_lon = 1
151 integer,
parameter :: i_lat = 2
153 integer,
parameter :: i_min = 1
154 integer,
parameter :: i_max = 2
156 integer,
parameter :: i_sclr = 1
157 integer,
parameter :: i_zstg = 2
158 integer,
parameter :: i_xstg = 3
159 integer,
parameter :: i_ystg = 4
161 integer,
parameter :: itp_ng = 4
162 integer,
private :: itp_nh = 4
164 integer,
parameter :: tag_lon = 1
165 integer,
parameter :: tag_lat = 2
166 integer,
parameter :: tag_lonuy = 3
167 integer,
parameter :: tag_latuy = 4
168 integer,
parameter :: tag_lonxv = 5
169 integer,
parameter :: tag_latxv = 6
170 integer,
parameter :: tag_cz = 7
171 integer,
parameter :: tag_fz = 8
173 integer,
parameter :: tag_dens = 1
174 integer,
parameter :: tag_momz = 2
175 integer,
parameter :: tag_momx = 3
176 integer,
parameter :: tag_momy = 4
177 integer,
parameter :: tag_rhot = 5
178 integer,
parameter :: tag_qx = 6
180 integer,
parameter :: order_tag_comm = 100000
181 integer,
parameter :: order_tag_var = 1000
185 integer,
private,
parameter :: interp_search_divnum = 10
187 integer,
private :: intercomm_id(2)
189 integer,
private :: max_isu
190 integer,
private :: max_rq = 1000
191 integer,
private :: rq_ctl_p
192 integer,
private :: rq_ctl_d
193 integer,
private :: rq_tot_p
194 integer,
private :: rq_tot_d
195 integer,
private,
allocatable :: ireq_p(:)
196 integer,
private,
allocatable :: ireq_d(:)
197 integer,
private,
allocatable :: call_order(:)
199 real(
rp),
private,
allocatable :: recvbuf_3d(:,:,:,:)
201 real(
rp),
private,
allocatable :: buffer_ref_lon (:,:)
202 real(
rp),
private,
allocatable :: buffer_ref_lonuy(:,:)
203 real(
rp),
private,
allocatable :: buffer_ref_lonxv(:,:)
204 real(
rp),
private,
allocatable :: buffer_ref_lat (:,:)
205 real(
rp),
private,
allocatable :: buffer_ref_latuy(:,:)
206 real(
rp),
private,
allocatable :: buffer_ref_latxv(:,:)
207 real(
rp),
private,
allocatable :: buffer_ref_cz (:,:,:)
208 real(
rp),
private,
allocatable :: buffer_ref_fz (:,:,:)
211 real(
rp),
private,
allocatable :: buffer_ref_3d (:,:,:)
213 real(
rp),
private,
allocatable :: org_dens(:,:,:)
214 real(
rp),
private,
allocatable :: org_momz(:,:,:)
215 real(
rp),
private,
allocatable :: org_momx(:,:,:)
216 real(
rp),
private,
allocatable :: org_momy(:,:,:)
217 real(
rp),
private,
allocatable :: org_u_ll(:,:,:)
218 real(
rp),
private,
allocatable :: org_v_ll(:,:,:)
219 real(
rp),
private,
allocatable :: org_rhot(:,:,:)
220 real(
rp),
private,
allocatable :: org_qtrc(:,:,:,:)
222 integer,
private,
allocatable :: igrd (:,:,:,:)
223 integer,
private,
allocatable :: jgrd (:,:,:,:)
224 real(
rp),
private,
allocatable :: hfact(:,:,:,:)
225 integer,
private,
allocatable :: kgrd (:,:,:,:,:,:)
226 real(
rp),
private,
allocatable :: vfact(:, :,:,:,:)
228 integer(8),
private :: nwait_p, nwait_d, nrecv, nsend
230 character(len=H_SHORT) :: mp_type
242 file_get_attribute, &
279 mapprojection_lonlat2xy
282 integer,
intent(in) :: qa_mp
283 character(len=*),
intent(in) :: mp_type_in
285 character(len=H_SHORT) :: comm_cartesc_nest_interp_type =
'LINEAR'
289 real(
rp),
allocatable :: x_ref(:,:)
290 real(
rp),
allocatable :: y_ref(:,:)
292 integer :: online_specified_maxrq = 0
297 logical :: flag_parent
298 logical :: flag_child
301 logical :: parent_periodic_x
302 logical :: parent_periodic_y
306 namelist / param_comm_cartesc_nest / &
313 online_aggressive_comm, &
315 online_specified_maxrq, &
316 comm_cartesc_nest_interp_type, &
323 log_info(
"COMM_CARTESC_NEST_setup",*)
'Setup'
335 online_wait_limit = 999999999
336 online_aggressive_comm = .true.
340 read(
io_fid_conf,nml=param_comm_cartesc_nest,iostat=ierr)
342 log_info(
"COMM_CARTESC_NEST_setup",*)
'Not found namelist. Default used.'
343 elseif( ierr > 0 )
then
344 log_error(
"COMM_CARTESC_NEST_setup",*)
'Not appropriate names in namelist PARAM_COMM_CARTESC_NEST. Check!'
347 log_nml(param_comm_cartesc_nest)
358 select case ( comm_cartesc_nest_interp_type )
361 case (
'DIST-WEIGHT' )
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"'
379 if( online_specified_maxrq > max_rq ) max_rq = online_specified_maxrq
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
405 else if ( mp_type_in ==
"NONE" )
then
413 log_info(
"COMM_CARTESC_NEST_setup",*)
"flag_parent", flag_parent,
"flag_child", flag_child
416 if( flag_parent )
then
419 log_error(
"COMM_CARTESC_NEST_setup",*)
'[NEST_setup] Parent Flag from launcher is not consistent with namelist!'
428 log_info(
"COMM_CARTESC_NEST_setup",
'(1x,A,I2,A)')
'Online Nesting - PARENT [INTERCOMM_ID:', &
437 log_info(
"COMM_CARTESC_NEST_setup",
'(1x,A)' )
'Informations of Daughter Domain'
439 log_info_cont(
'(1x,A,I6) ')
'Limit Num. NCOMM req. :', max_rq
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) )
456 if( flag_child )
then
459 log_error(
"COMM_CARTESC_NEST_setup",*)
'[NEST_setup] Child Flag from launcher is not consistent with namelist!'
468 log_info(
"COMM_CARTESC_NEST_setup",
'(1x,A,I2,A)')
'Online Nesting - DAUGHTER [INTERCOMM_ID:', &
472 num_dom = num_dom + 1
474 if ( i_parent > max_dinfo )
then
475 log_error(
"COMM_CARTESC_NEST_setup",*)
'number of domain exeeds the limit'
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) )
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
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
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
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) )
519 allocate( buffer_ref_cz(tileal_ka,tileal_ia,tileal_ja) )
520 allocate( buffer_ref_fz(tileal_ka,tileal_ia,tileal_ja) )
522 allocate( buffer_ref_3d(tileal_ka,tileal_ia,tileal_ja) )
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) )
533 select case ( comm_cartesc_nest_interp_type )
536 allocate( x_ref(tileal_ia,tileal_ja) )
537 allocate( y_ref(tileal_ia,tileal_ja) )
540 call mapprojection_lonlat2xy( tileal_ia, 1, tileal_ia, &
541 tileal_ja, 1, tileal_ja, &
542 buffer_ref_lon(:,:), &
543 buffer_ref_lat(:,:), &
544 x_ref(:,:), y_ref(:,:) )
545 call interp_factor3d( tileal_ka,
khalo+1, tileal_ka-
khalo, &
546 tileal_ia, tileal_ja, &
548 x_ref(:,:), y_ref(:,:), &
549 buffer_ref_cz(:,:,:), &
553 igrd( :,:,:,i_sclr), &
554 jgrd( :,:,:,i_sclr), &
555 hfact( :,:,:,i_sclr), &
556 kgrd(:,:,:,:,:,i_sclr), &
557 vfact(:, :,:,:,i_sclr) )
560 call interp_factor3d( tileal_ka+1,
khalo+1, tileal_ka+1-
khalo, &
561 tileal_ia, tileal_ja, &
563 x_ref(:,:), y_ref(:,:), &
564 buffer_ref_fz(:,:,:), &
568 igrd( :,:,:,i_zstg), &
569 jgrd( :,:,:,i_zstg), &
570 hfact( :,:,:,i_zstg), &
571 kgrd(:,:,:,:,:,i_zstg), &
572 vfact(:, :,:,:,i_zstg) )
575 call mapprojection_lonlat2xy( tileal_ia, 1, tileal_ia, &
576 tileal_ja, 1, tileal_ja, &
577 buffer_ref_lonuy(:,:), &
578 buffer_ref_latuy(:,:), &
579 x_ref(:,:), y_ref(:,:) )
580 call interp_factor3d( tileal_ka,
khalo+1, tileal_ka-
khalo, &
581 tileal_ia, tileal_ja, &
583 x_ref(:,:), y_ref(:,:), &
584 buffer_ref_cz(:,:,:), &
588 igrd( :,:,:,i_xstg), &
589 jgrd( :,:,:,i_xstg), &
590 hfact( :,:,:,i_xstg), &
591 kgrd(:,:,:,:,:,i_xstg), &
592 vfact(:, :,:,:,i_xstg) )
595 call mapprojection_lonlat2xy( tileal_ia, 1, tileal_ia, &
596 tileal_ja, 1, tileal_ja, &
597 buffer_ref_lonxv(:,:), &
598 buffer_ref_latxv(:,:), &
599 x_ref(:,:), y_ref(:,:) )
600 call interp_factor3d( tileal_ka,
khalo+1, tileal_ka-
khalo, &
601 tileal_ia, tileal_ja, &
603 x_ref(:,:), y_ref(:,:), &
604 buffer_ref_cz(:,:,:), &
608 igrd( :,:,:,i_ystg), &
609 jgrd( :,:,:,i_ystg), &
610 hfact( :,:,:,i_ystg), &
611 kgrd(:,:,:,:,:,i_ystg), &
612 vfact(:, :,:,:,i_ystg) )
614 deallocate( x_ref, y_ref )
616 case (
'DIST-WEIGHT' )
619 call interp_factor3d( itp_nh, &
621 tileal_ia, tileal_ja, &
623 buffer_ref_lon(:,:), &
624 buffer_ref_lat(:,:), &
625 buffer_ref_cz(:,:,:), &
629 igrd( :,:,:,i_sclr), &
630 jgrd( :,:,:,i_sclr), &
631 hfact( :,:,:,i_sclr), &
632 kgrd(:,:,:,:,:,i_sclr), &
633 vfact(:, :,:,:,i_sclr) )
636 call interp_factor3d( itp_nh, &
638 tileal_ia, tileal_ja, &
640 buffer_ref_lon(:,:), &
641 buffer_ref_lat(:,:), &
642 buffer_ref_fz(:,:,:), &
646 igrd( :,:,:,i_zstg), &
647 jgrd( :,:,:,i_zstg), &
648 hfact( :,:,:,i_zstg), &
649 kgrd(:,:,:,:,:,i_zstg), &
650 vfact(:, :,:,:,i_zstg) )
653 call interp_factor3d( itp_nh, &
655 tileal_ia, tileal_ja, &
657 buffer_ref_lonuy(:,:), &
658 buffer_ref_latuy(:,:), &
659 buffer_ref_cz(:,:,:), &
663 igrd( :,:,:,i_xstg), &
664 jgrd( :,:,:,i_xstg), &
665 hfact( :,:,:,i_xstg), &
666 kgrd(:,:,:,:,:,i_xstg), &
667 vfact(:, :,:,:,i_xstg) )
670 call interp_factor3d( itp_nh, &
672 tileal_ia, tileal_ja, &
674 buffer_ref_lonxv(:,:), &
675 buffer_ref_latxv(:,:), &
676 buffer_ref_cz(:,:,:), &
680 igrd( :,:,:,i_ystg), &
681 jgrd( :,:,:,i_ystg), &
682 hfact( :,:,:,i_ystg), &
683 kgrd(:,:,:,:,:,i_ystg), &
684 vfact(:, :,:,:,i_ystg) )
707 LATLON_CATALOGUE_FNAME )
713 file_get_attribute, &
718 integer,
intent(out) :: dom_id
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
728 integer :: pnum_x(1), pnum_y(1)
730 integer :: dims(1), dims2(1)
732 integer :: parent_x, parent_xh
733 integer :: parent_y, parent_yh
735 real(
rp),
allocatable :: work(:,:), work_uv(:,:), minmax(:,:,:)
737 character(len=H_LONG) :: fname
740 logical :: existed, error
746 if ( dom_info(n)%basename == parent_basename )
then
752 num_dom = num_dom + 1
754 if ( dom_id > max_dinfo )
then
755 log_error(
"COMM_CARTESC_NEST_domain_regist_file",*)
'number of domains exceed the limit'
758 dinfo => dom_info(dom_id)
759 dinfo%basename = parent_basename
764 aggregate = .false., &
767 call file_get_attribute( fid,
"global",
"scale_atmos_grid_cartesC_index_imaxg", &
768 imaxg(:), existed=existed )
770 call file_get_attribute( fid,
"global",
"scale_cartesC_prc_num_x", &
772 dinfo%prc_num_x = pnum_x(1)
773 dinfo%IMAX =
imaxg(1) / pnum_x(1)
775 call file_get_attribute( fid,
"x",
"halo_global", &
777 dinfo%IHALO = halos(1)
778 call file_get_attribute( fid,
"global",
"scale_cartesC_prc_periodic_x", &
781 call file_get_attribute( fid,
"global",
"scale_atmos_grid_cartesC_index_jmaxg", &
783 call file_get_attribute( fid,
"global",
"scale_cartesC_prc_num_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", &
789 dinfo%JHALO = halos(1)
790 call file_get_attribute( fid,
"global",
"scale_cartesC_prc_periodic_y", &
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
800 log_error(
"COMM_CARTESC_NEST_domain_regist_file",*)
'PARENT_PRC_NUM_(X|Y) is needed for files generated by the older version'
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)
807 dinfo%IMAX = dims(1) - dinfo%IHALO*2
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)
812 dinfo%JMAX = dims(1) - dinfo%JHALO*2
814 dinfo%periodic_x = .false.
815 dinfo%periodic_y = .false.
819 call file_get_attribute( fid,
"global",
"scale_atmos_grid_cartesC_index_kmax", &
820 dims(:), existed=existed )
824 call file_get_shape( fid,
"z", dims(:), error=error )
833 call file_get_attribute( fid,
"global",
"scale_ocean_grid_cartesC_index_kmax", &
834 dims(:), existed=existed )
836 dinfo%OKMAX = dims(1)
838 call file_get_shape( fid,
"oz", dims(:), error=error )
842 dinfo%OKMAX = dims(1)
846 call file_get_attribute( fid,
"global",
"scale_land_grid_cartesC_index_kmax", &
847 dims(:), existed=existed )
849 dinfo%LKMAX = dims(1)
851 call file_get_shape( fid,
"lz", dims(:), error=error )
855 dinfo%LKMAX = dims(1)
859 call file_get_attribute( fid,
"global",
"scale_urban_grid_cartesC_index_kmax", &
860 dims(:), existed=existed )
862 dinfo%UKMAX = dims(1)
864 call file_get_shape( fid,
"uz", dims(:), error=error )
868 dinfo%UKMAX = dims(1)
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 )
889 nprocs = dinfo%prc_num_x * dinfo%prc_num_y
890 allocate( dinfo%latlon_catalogue(nprocs,2,2) )
894 existed =
present(latlon_catalogue_fname)
896 existed = latlon_catalogue_fname /=
""
905 form =
'formatted', &
909 if ( ierr /= 0 )
then
910 log_error(
"COMM_CARTESC_NEST_domain_regist_file",*)
'cannot open latlon-catalogue file!: ', trim(fname)
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), &
917 dinfo%latlon_catalogue(i,i_min,i_lat), dinfo%latlon_catalogue(i,i_max,i_lat)
918 if ( ierr /= 0 .or. i /= parent_id )
then
919 log_error(
"COMM_CARTESC_NEST_domain_regist_file",*)
'catalogue file is invalid, ', trim(fname)
922 if ( ierr /= 0 )
exit
929 allocate( minmax(nprocs,2,2) )
932 do j = 1, dinfo%prc_num_y
933 do i = 1, dinfo%prc_num_x
936 aggregate = .false., &
937 allnodes = .false., &
940 call file_get_shape( fid,
"xh", dims(:) )
942 call file_get_shape( fid,
"yh", dims(:) )
944 allocate( work_uv( parent_xh, parent_yh ) )
946 if ( dinfo%periodic_x .or. dinfo%periodic_y )
then
947 call file_get_shape( fid,
"x", dims(:) )
949 call file_get_shape( fid,
"y", dims(:) )
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(:,:) )
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) )
961 if ( dinfo%periodic_x )
then
962 allocate( work( parent_x, parent_yh ) )
963 call file_read( fid,
"lon_xv", work(:,:) )
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 ) )
971 minmax(n,i_min,i_lon) = minval( work_uv(parent_xh,:) )
972 minmax(n,i_max,i_lon) = maxval( work_uv(parent_xh,:) )
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(:,:) )
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) )
982 if ( dinfo%periodic_y )
then
983 allocate( work( parent_xh, parent_y ) )
984 call file_read( fid,
"lat_uy", work(:,:) )
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) ) )
992 minmax(n,i_min,i_lat) = minval( work_uv(:,parent_yh) )
993 minmax(n,i_max,i_lat) = maxval( work_uv(:,parent_yh) )
995 deallocate( work_uv )
1001 deallocate( minmax )
1007 call comm_bcast( nprocs, 2, 2, dinfo%latlon_catalogue(:,:,:) )
1023 integer,
intent(in) :: dom_id
1028 integer :: x_min, x_max
1029 integer :: y_min, y_max
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
1040 dinfo => dom_info(dom_id)
1041 nprocs = dinfo%prc_num_x * dinfo%prc_num_y
1043 x_min = dinfo%prc_num_x
1045 y_min = dinfo%prc_num_y
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
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'
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)
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)
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)
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)
1102 dinfo%tile_num_x = x_max - x_min + 1
1103 dinfo%tile_num_y = y_max - y_min + 1
1105 allocate( dinfo%tile_id(dinfo%tile_num_x * dinfo%tile_num_y) )
1107 log_info(
"COMM_CARTESC_NEST_domain_relate",
'(1x,A)')
'NEST: target process tile in parent domain'
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)
1132 integer,
intent(in) :: dom_id
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(:)
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
1150 if (
present(
kmax) ) &
1151 kmax = dom_info(dom_id)%KMAX
1153 if (
present(lkmax) ) &
1154 lkmax = dom_info(dom_id)%LKMAX
1156 if (
present(
imaxg) ) &
1157 imaxg = dom_info(dom_id)%IMAX * dom_info(dom_id)%tile_num_x
1159 if (
present(
jmaxg) ) &
1160 jmaxg = dom_info(dom_id)%JMAX * dom_info(dom_id)%tile_num_y
1162 if (
present(num_tile) ) &
1163 num_tile = dom_info(dom_id)%tile_num_x * dom_info(dom_id)%tile_num_y
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)
1191 integer,
intent(out) :: tilei, tilej
1192 integer,
intent(out) :: cxs, cxe, cys, cye
1193 integer,
intent(out) :: pxs, pxe, pys, pye
1195 integer,
intent(in) :: dom_id
1196 integer,
intent(in) :: iloc
1198 logical,
intent(in),
optional :: xstg
1199 logical,
intent(in),
optional :: ystg
1205 integer :: xloc, yloc
1206 integer :: xlocg, ylocg
1207 logical :: xstg_, ystg_
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
1215 if (
present(xstg) )
then
1220 if (
present(ystg) )
then
1226 dinfo => dom_info(dom_id)
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
1236 cxs = tilei * (xloc-1) + 1
1238 cys = tilej * (yloc-1) + 1
1245 if ( .not. dinfo%periodic_x )
then
1246 if ( xlocg == 1 )
then
1247 tilei = tilei + dinfo%IHALO
1248 pxs = pxs + dinfo%IHALO
1249 pxe = pxe + dinfo%IHALO
1251 if ( xlocg == dinfo%prc_num_x )
then
1252 tilei = tilei + dinfo%IHALO
1256 if ( xlocg == 1 )
then
1258 if ( dinfo%IHALO > 0 )
then
1270 if ( .not. dinfo%periodic_y )
then
1271 if ( ylocg == 1 )
then
1272 tilej = tilej + dinfo%JHALO
1273 pys = pys + dinfo%JHALO
1274 pye = pye + dinfo%JHALO
1276 if ( ylocg == dinfo%prc_num_y )
then
1277 tilej = tilej + dinfo%JHALO
1281 if ( ylocg == 1 )
then
1283 if ( dinfo%JHALO > 0 )
then
1300 subroutine comm_cartesc_nest_parentsize( &
1321 integer,
intent(in) :: handle
1324 integer,
parameter :: ileng = 10
1325 integer :: datapack(ileng)
1327 integer :: qa_otherside
1328 character(len=H_SHORT) :: mp_type_otherside
1330 integer :: ireq1, ireq2, ireq3, ireq4, ireq5, ireq6
1331 integer :: ierr1, ierr2, ierr3, ierr4, ierr5, ierr6
1332 integer :: istatus(mpi_status_size)
1338 tag = intercomm_id(handle) * 100
1347 datapack( 4) =
khalo
1349 datapack( 6) =
ihalo
1351 datapack( 8) =
jhalo
1366 call mpi_wait(ireq1, istatus, ierr1)
1367 call mpi_wait(ireq2, istatus, ierr2)
1368 call mpi_wait(ireq3, istatus, ierr3)
1370 call mpi_wait(ireq4, istatus, ierr4)
1371 call mpi_wait(ireq5, istatus, ierr5)
1372 call mpi_wait(ireq6, istatus, ierr6)
1376 call comm_bcast(qa_otherside)
1377 call comm_bcast(mp_type_otherside)
1379 if ( mp_type == mp_type_otherside .and.
online_send_qa == qa_otherside )
then
1381 else if ( mp_type_otherside ==
"DRY" .or. mp_type ==
"DRY" )
then
1384 else if ( mp_type_otherside ==
"QV" .or. mp_type ==
"QV" )
then
1388 log_info(
"COMM_CARTESC_NEST_parentsize",*)
'Hydrometeor will be diagnosed on children side'
1389 log_info(
"COMM_CARTESC_NEST_parentsize",*)
'MP type (remote,local) = ', trim(mp_type_otherside),
", ", trim(mp_type)
1390 log_info(
"COMM_CARTESC_NEST_parentsize",*)
'Number of QA (remote,local) = ', qa_otherside,
online_send_qa
1411 call mpi_wait(ireq4, istatus, ierr4)
1412 call mpi_wait(ireq5, istatus, ierr5)
1413 call mpi_wait(ireq6, istatus, ierr6)
1415 call mpi_wait(ireq1, istatus, ierr1)
1416 call mpi_wait(ireq2, istatus, ierr2)
1417 call mpi_wait(ireq3, istatus, ierr3)
1420 call comm_bcast(ileng, datapack)
1421 call comm_bcast(buffer)
1422 call comm_bcast(mp_type_otherside)
1424 dom_info(i_parent)%prc_num_x = datapack( 1)
1425 dom_info(i_parent)%prc_num_y = datapack( 2)
1426 dom_info(i_parent)%KMAX = datapack( 3)
1427 dom_info(i_parent)%KHALO = datapack( 4)
1428 dom_info(i_parent)%IMAX = datapack( 5)
1429 dom_info(i_parent)%IHALO = datapack( 6)
1430 dom_info(i_parent)%JMAX = datapack( 7)
1431 dom_info(i_parent)%JHALO = datapack( 8)
1433 qa_otherside = datapack(10)
1436 if ( mp_type == mp_type_otherside .and.
online_recv_qa == qa_otherside )
then
1438 else if ( mp_type ==
"DRY" )
then
1441 else if ( mp_type ==
"QV" )
then
1445 log_info(
"COMM_CARTESC_NEST_parentsize",*)
'Hydrometeor will be diagnosed on this side'
1446 log_info(
"COMM_CARTESC_NEST_parentsize",*)
'MP type (remote,local) = ', trim(mp_type_otherside),
", ", trim(mp_type)
1447 log_info(
"COMM_CARTESC_NEST_parentsize",*)
'Number of QA (remote,local) = ', qa_otherside,
online_recv_qa
1453 log_error(
"COMM_CARTESC_NEST_parentsize",*)
'[COMM_CARTESC_NEST_parentsize] internal error'
1458 end subroutine comm_cartesc_nest_parentsize
1462 subroutine comm_cartesc_nest_catalogue( &
1478 integer,
intent(in) :: handle
1481 integer :: ireq, ierr, ileng
1482 integer :: istatus(mpi_status_size)
1488 tag = intercomm_id(handle) * 100
1498 call mpi_wait(ireq, istatus, ierr)
1505 nprocs = dom_info(i_parent)%prc_num_x * dom_info(i_parent)%prc_num_y
1506 ileng = nprocs * 2 * 2
1510 call mpi_wait(ireq, istatus, ierr)
1512 call comm_bcast( nprocs, 2, 2, dom_info(i_parent)%latlon_catalogue )
1515 log_error(
"COMM_CARTESC_NEST_catalogue",*)
'internal error'
1520 end subroutine comm_cartesc_nest_catalogue
1524 subroutine comm_cartesc_nest_ping( &
1536 integer,
intent(in) :: handle
1538 integer :: ping, pong
1539 integer :: ireq1, ireq2, ierr1, ierr2
1540 integer :: istatus(mpi_status_size)
1542 logical :: ping_error
1547 tag = intercomm_id(handle) * 100
1548 ping_error = .false.
1560 call mpi_wait(ireq1, istatus, ierr1)
1561 call mpi_wait(ireq2, istatus, ierr2)
1564 call comm_bcast(pong)
1566 if ( pong /= intercomm_id(handle)+1 ) ping_error = .true.
1578 call mpi_wait(ireq1, istatus, ierr1)
1579 call mpi_wait(ireq2, istatus, ierr2)
1582 call comm_bcast(pong)
1584 if ( pong /= intercomm_id(handle) ) ping_error = .true.
1587 log_error(
"COMM_CARTESC_NEST_ping",*)
'internal error'
1591 if ( ping_error )
then
1592 log_error(
"COMM_CARTESC_NEST_ping",*)
'ping destination error'
1597 end subroutine comm_cartesc_nest_ping
1601 subroutine comm_cartesc_nest_setup_nestdown( &
1616 integer,
intent(in) :: handle
1618 integer,
allocatable :: buffer_list (:)
1619 integer,
allocatable :: buffer_alllist(:)
1621 integer :: parent_ka
1622 integer :: parent_ia
1623 integer :: parent_ja
1625 integer :: ireq, ierr, ileng
1626 integer :: istatus(mpi_status_size)
1627 integer :: tag, target_rank
1634 tag = intercomm_id(handle) * 100
1640 if ( prc_ismaster )
then
1641 call mpi_irecv(comm_cartesc_nest_tile_allmax_p, 1, mpi_integer, prc_myrank, tag+1, prc_intercomm_child, ireq, ierr)
1642 call mpi_wait(ireq, istatus, ierr)
1644 call comm_bcast(comm_cartesc_nest_tile_allmax_p)
1650 if ( prc_ismaster )
then
1651 call mpi_irecv(comm_cartesc_nest_tile_list_p, ileng, mpi_integer, prc_myrank, tag+2, prc_intercomm_child, ireq, ierr)
1652 call mpi_wait(ireq, istatus, ierr)
1656 comm_cartesc_nest_tile_list_yp(:) = -1
1660 do i = 1, comm_cartesc_nest_tile_allmax_p
1661 if ( comm_cartesc_nest_tile_list_p(i,j) == prc_myrank )
then
1663 comm_cartesc_nest_tile_list_yp(
k) = j - 1
1669 log_info(
"COMM_CARTESC_NEST_setup_nestdown",
'(A,I5,A,I5)')
"[P] Num YP =",num_yp,
" Num TILE(MAX) =",comm_cartesc_nest_tile_allmax_p
1671 if ( prc_ismaster )
then
1672 call mpi_irecv(online_daughter_use_velz, 1, mpi_logical, prc_myrank, tag+3, prc_intercomm_child, ireq, ierr)
1673 call mpi_wait(ireq, istatus, ierr)
1675 call comm_bcast(online_daughter_use_velz)
1677 log_info(
"COMM_CARTESC_NEST_setup_nestdown",
'(1x,A,L2)')
'NEST: ONLINE_DAUGHTER_USE_VELZ =', online_daughter_use_velz
1679 if ( prc_ismaster )
then
1680 call mpi_irecv(online_daughter_no_rotate, 1, mpi_logical, prc_myrank, tag+4, prc_intercomm_child, ireq, ierr)
1681 call mpi_wait(ireq, istatus, ierr)
1683 call comm_bcast(online_daughter_no_rotate)
1686 log_error(
"COMM_CARTESC_NEST_setup_nestdown",*)
'Flag of NO_ROTATE is not consistent with the child domain'
1688 log_error_cont(*)
'ONLINE_DAUGHTER_NO_ROTATE =', online_daughter_no_rotate
1691 log_info(
"COMM_CARTESC_NEST_setup_nestdown",
'(1x,A,L2)')
'NEST: ONLINE_DAUGHTER_NO_ROTATE =', online_daughter_no_rotate
1693 call comm_cartesc_nest_importgrid_nestdown( handle )
1696 target_rank = comm_cartesc_nest_tile_list_yp(i)
1697 call mpi_isend(i, 1, mpi_integer, target_rank, tag+5, prc_intercomm_child, ireq, ierr)
1698 call mpi_wait(ireq, istatus, ierr)
1701 call mpi_barrier(prc_intercomm_child, ierr)
1707 comm_cartesc_nest_tile_all =
size( dom_info(i_parent)%TILE_ID(:) )
1708 call mpi_allreduce( comm_cartesc_nest_tile_all, &
1709 comm_cartesc_nest_tile_allmax_d, &
1715 log_info(
"COMM_CARTESC_NEST_setup_nestdown",
'(A,I5,A,I5)')
"[D] Num YP =",comm_cartesc_nest_tile_all,
" Num TILE(MAX) =",comm_cartesc_nest_tile_allmax_d
1717 if ( prc_ismaster )
then
1718 call mpi_isend(comm_cartesc_nest_tile_allmax_d, 1, mpi_integer, prc_myrank, tag+1, prc_intercomm_parent, ireq, ierr)
1719 call mpi_wait(ireq, istatus, ierr)
1722 parent_ka = dom_info(i_parent)%KMAX + dom_info(i_parent)%KHALO * 2
1723 parent_ia = dom_info(i_parent)%IMAX + dom_info(i_parent)%IHALO * 2
1724 parent_ja = dom_info(i_parent)%JMAX + dom_info(i_parent)%JHALO * 2
1728 max_isu = comm_cartesc_nest_tile_all * max_isu
1729 allocate( recvbuf_3d( parent_ka, parent_ia, parent_ja, max_isu ) )
1731 allocate( buffer_list(comm_cartesc_nest_tile_allmax_d) )
1732 allocate( buffer_alllist(comm_cartesc_nest_tile_allmax_d*
prc_nprocs) )
1733 allocate( comm_cartesc_nest_tile_list_d(comm_cartesc_nest_tile_allmax_d,
prc_nprocs) )
1735 do i = 1, comm_cartesc_nest_tile_allmax_d
1736 if ( i <= comm_cartesc_nest_tile_all )
then
1737 buffer_list(i) = dom_info(i_parent)%TILE_ID(i)
1743 ileng = comm_cartesc_nest_tile_allmax_d
1744 call mpi_allgather( buffer_list(:), &
1747 buffer_alllist(:), &
1754 do i = 1, comm_cartesc_nest_tile_allmax_d
1755 comm_cartesc_nest_tile_list_d(i,j) = buffer_alllist(
k)
1760 deallocate( buffer_list )
1761 deallocate( buffer_alllist )
1763 ileng = comm_cartesc_nest_tile_allmax_d *
prc_nprocs
1764 if ( prc_ismaster )
then
1765 call mpi_isend(comm_cartesc_nest_tile_list_d, ileng, mpi_integer, prc_myrank, tag+2, prc_intercomm_parent, ireq, ierr)
1766 call mpi_wait(ireq, istatus, ierr)
1769 if ( prc_ismaster )
then
1770 call mpi_isend(
online_use_velz, 1, mpi_logical, prc_myrank, tag+3, prc_intercomm_parent, ireq, ierr)
1771 call mpi_wait(ireq, istatus, ierr)
1774 if ( prc_ismaster )
then
1775 call mpi_isend(
online_no_rotate, 1, mpi_logical, prc_myrank, tag+4, prc_intercomm_parent, ireq, ierr)
1776 call mpi_wait(ireq, istatus, ierr)
1778 call comm_bcast(online_daughter_no_rotate)
1780 call comm_cartesc_nest_importgrid_nestdown( handle )
1782 do i = 1, comm_cartesc_nest_tile_all
1783 target_rank = comm_cartesc_nest_tile_list_d(i,prc_myrank+1)
1784 call mpi_irecv( call_order(i), 1, mpi_integer, target_rank, tag+5, prc_intercomm_parent, ireq, ierr )
1785 call mpi_wait(ireq, istatus, ierr)
1788 call mpi_barrier(prc_intercomm_parent, ierr)
1790 log_error(
"COMM_CARTESC_NEST_setup_nestdown",*)
'internal error'
1794 if( num_yp * 16 > max_rq .OR. comm_cartesc_nest_tile_all * 16 > max_rq )
then
1795 log_error(
"COMM_CARTESC_NEST_setup_nestdown",*)
'internal error (overflow number of ireq)'
1796 log_error_cont(*)
'NUM_YP x 16 = ', num_yp * 16
1797 log_error_cont(*)
'COMM_CARTESC_NEST_TILE_ALL x 16 = ', comm_cartesc_nest_tile_all * 16
1798 log_error_cont(*)
'max_rq = ', max_rq
1803 end subroutine comm_cartesc_nest_setup_nestdown
1807 subroutine comm_cartesc_nest_importgrid_nestdown( &
1827 integer,
intent(in) :: handle
1829 integer :: parent_ka
1830 integer :: parent_ia, parent_is, parent_ie, parent_imax
1831 integer :: parent_ja, parent_js, parent_je, parent_jmax
1833 integer :: ierr, ileng
1834 integer :: istatus(mpi_status_size)
1835 integer :: tag, tagbase, target_rank
1836 integer :: rq_str, rq_end, rq_tot
1838 integer :: xloc, yloc
1842 real(
rp) :: max_ref, max_loc
1844 real(
rp),
allocatable :: sendbuf_2d(:,:,:)
1845 real(
rp),
allocatable :: sendbuf_3d(:,:,:,:)
1846 real(
rp),
allocatable :: recvbuf_2d(:,:,:)
1853 tagbase = intercomm_id(handle) * 100
1860 allocate( sendbuf_2d(
ia,
ja, 4 ) )
1861 allocate( sendbuf_3d(
ka,
ia,
ja, 1 ) )
1865 target_rank = comm_cartesc_nest_tile_list_yp(i)
1871 tag = tagbase + tag_lon
1876 tag = tagbase + tag_lat
1882 tag = tagbase + tag_lonuy
1888 tag = tagbase + tag_latuy
1894 tag = tagbase + tag_lonxv
1900 tag = tagbase + tag_latxv
1905 tag = tagbase + tag_cz
1911 tag = tagbase + tag_fz
1912 call mpi_isend(sendbuf_3d(:,:,:,1), ileng,
comm_datatype, target_rank, tag,
prc_intercomm_child, ireq_p(rq), ierr)
1915 rq_tot = rq_end - rq_str + 1
1920 deallocate( sendbuf_2d )
1921 deallocate( sendbuf_3d )
1927 parent_ka = dom_info(i_parent)%KMAX + dom_info(i_parent)%KHALO * 2
1929 parent_imax = dom_info(i_parent)%IMAX
1930 parent_ia = parent_imax + dom_info(i_parent)%IHALO * 2
1931 parent_is = dom_info(i_parent)%IHALO + 1
1932 parent_ie = parent_imax + dom_info(i_parent)%IHALO
1934 parent_jmax = dom_info(i_parent)%JMAX
1935 parent_ja = parent_jmax + dom_info(i_parent)%JHALO * 2
1936 parent_js = dom_info(i_parent)%JHALO + 1
1937 parent_je = parent_jmax + dom_info(i_parent)%JHALO
1939 allocate( recvbuf_2d( parent_ia, parent_ja, 6 ) )
1941 do i = 1, comm_cartesc_nest_tile_all
1943 target_rank = comm_cartesc_nest_tile_list_d(i,
prc_myrank+1)
1945 xloc = mod( i-1, dom_info(i_parent)%tile_num_x ) + 1
1946 yloc = int( real(i-1) / real(dom_info(i_parent)%tile_num_x) ) + 1
1948 xs = parent_imax * (xloc-1) + 1
1949 xe = parent_imax * xloc
1950 ys = parent_jmax * (yloc-1) + 1
1951 ye = parent_jmax * yloc
1956 ileng = parent_ia * parent_ja
1957 tag = tagbase + tag_lon
1961 ileng = parent_ia * parent_ja
1962 tag = tagbase + tag_lat
1966 ileng = parent_ia * parent_ja
1967 tag = tagbase + tag_lonuy
1971 ileng = parent_ia * parent_ja
1972 tag = tagbase + tag_latuy
1976 ileng = parent_ia * parent_ja
1977 tag = tagbase + tag_lonxv
1981 ileng = parent_ia * parent_ja
1982 tag = tagbase + tag_latxv
1986 ileng = parent_ka * parent_ia * parent_ja
1987 tag = tagbase + tag_cz
1988 call mpi_irecv(recvbuf_3d(:,:,:,tag_cz), ileng,
comm_datatype, target_rank, tag,
prc_intercomm_parent, ireq_d(rq), ierr)
1991 ileng = parent_ka * parent_ia * parent_ja
1992 tag = tagbase + tag_fz
1993 call mpi_irecv(recvbuf_3d(:,:,:,tag_fz), ileng,
comm_datatype, target_rank, tag,
prc_intercomm_parent, ireq_d(rq), ierr)
1996 rq_tot = rq_end - rq_str + 1
2000 buffer_ref_lon(xs:xe,ys:ye) = recvbuf_2d(parent_is:parent_ie,parent_js:parent_je,tag_lon )
2001 buffer_ref_lat(xs:xe,ys:ye) = recvbuf_2d(parent_is:parent_ie,parent_js:parent_je,tag_lat )
2002 buffer_ref_lonuy(xs:xe,ys:ye) = recvbuf_2d(parent_is:parent_ie,parent_js:parent_je,tag_lonuy)
2003 buffer_ref_latuy(xs:xe,ys:ye) = recvbuf_2d(parent_is:parent_ie,parent_js:parent_je,tag_latuy)
2004 buffer_ref_lonxv(xs:xe,ys:ye) = recvbuf_2d(parent_is:parent_ie,parent_js:parent_je,tag_lonxv)
2005 buffer_ref_latxv(xs:xe,ys:ye) = recvbuf_2d(parent_is:parent_ie,parent_js:parent_je,tag_latxv)
2008 buffer_ref_cz(
k,xs:xe,ys:ye) = recvbuf_3d(
k,parent_is:parent_ie,parent_js:parent_je,tag_cz)
2009 buffer_ref_fz(
k,xs:xe,ys:ye) = recvbuf_3d(
k,parent_is:parent_ie,parent_js:parent_je,tag_fz)
2014 max_ref = maxval( buffer_ref_fz(:,:,:) )
2016 if ( max_ref < max_loc )
then
2017 log_error(
"COMM_CARTESC_NEST_importgrid_nestdown",*)
'REQUESTED DOMAIN IS TOO MUCH BROAD'
2018 log_error_cont(*)
'-- VERTICAL direction over the limit'
2019 log_error_cont(*)
'-- reference max: ', max_ref
2020 log_error_cont(*)
'-- local max: ', max_loc
2024 deallocate( recvbuf_2d )
2027 log_error(
"COMM_CARTESC_NEST_importgrid_nestdown",*)
'internal error'
2032 end subroutine comm_cartesc_nest_importgrid_nestdown
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)
2061 integer,
parameter :: handle = 1
2065 real(
rp) :: u_on_map, v_on_map
2067 real(
rp) :: dummy(1,1,1)
2068 integer :: tagbase, tagcomm
2072 integer :: i, j,
k, iq
2077 tagcomm = intercomm_id(handle) * order_tag_comm
2087 log_info(
"COMM_CARTESC_NEST_nestdown",
'(1X,A,I5,A)')
"CONeP[P] send( ", nsend,
" )"
2091 org_dens(:,:,:) = dens_send(:,:,:)
2093 org_momz(:,:,:) = momz_send(:,:,:)
2095 org_momx(:,:,:) = momx_send(:,:,:)
2097 org_momy(:,:,:) = momy_send(:,:,:)
2099 org_rhot(:,:,:) = rhot_send(:,:,:)
2102 org_qtrc(:,:,:,iq) = qtrc_send(:,:,:,iq)
2110 if ( .NOT. online_daughter_no_rotate )
then
2115 work1_send(
k,i,j) = ( org_momx(
k,i-1,j) + org_momx(
k,i,j) ) * 0.5_rp
2122 work1_send(
k,1,j) = org_momx(
k,1,j)
2126 call comm_vars8( work1_send(:,:,:), 1 )
2131 work2_send(
k,i,j) = ( org_momy(
k,i,j-1) + org_momy(
k,i,j) ) * 0.5_rp
2138 work2_send(
k,i,1) = org_momy(
k,i,1)
2142 call comm_vars8( work2_send(:,:,:), 2 )
2144 call comm_wait ( work1_send(:,:,:), 1, .false. )
2145 call comm_wait ( work2_send(:,:,:), 2, .false. )
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)
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)
2161 tagbase = tagcomm + tag_dens*order_tag_var
2162 call comm_cartesc_nest_intercomm_nestdown( org_dens(:,:,:), &
2164 tagbase, i_sclr, handle, &
2166 flag_dens = .true. )
2168 tagbase = tagcomm + tag_momz*order_tag_var
2169 if ( online_daughter_use_velz )
then
2170 call comm_cartesc_nest_intercomm_nestdown( org_momz(:,:,:), &
2172 tagbase, i_zstg, handle, &
2176 tagbase = tagcomm + tag_momx*order_tag_var
2177 if ( online_daughter_no_rotate )
then
2178 call comm_cartesc_nest_intercomm_nestdown( org_momx(:,:,:), &
2180 tagbase, i_xstg, handle, &
2183 call comm_cartesc_nest_intercomm_nestdown( org_u_ll(:,:,:), &
2185 tagbase, i_sclr, handle, &
2189 tagbase = tagcomm + tag_momy*order_tag_var
2190 if ( online_daughter_no_rotate )
then
2191 call comm_cartesc_nest_intercomm_nestdown( org_momy(:,:,:), &
2193 tagbase, i_ystg, handle, &
2196 call comm_cartesc_nest_intercomm_nestdown( org_v_ll(:,:,:), &
2198 tagbase, i_sclr, handle, &
2202 tagbase = tagcomm + tag_rhot*order_tag_var
2203 call comm_cartesc_nest_intercomm_nestdown( org_rhot(:,:,:), &
2205 tagbase, i_sclr, handle, &
2209 tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
2210 call comm_cartesc_nest_intercomm_nestdown( org_qtrc(:,:,:,iq), &
2212 tagbase, i_sclr, handle, &
2222 log_error(
"COMM_CARTESC_NEST_nestdown_send",*)
'internal error'
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)
2255 integer,
parameter :: handle = 2
2261 real(
rp) :: u_on_map, v_on_map
2263 real(
rp) :: dummy(1,1,1)
2264 integer :: tagbase, tagcomm
2268 integer :: i, j,
k, iq
2273 tagcomm = intercomm_id(handle) * order_tag_comm
2282 nwait_d = nwait_d + 1
2292 if ( online_aggressive_comm )
then
2301 tagbase = tagcomm + tag_dens*order_tag_var
2302 call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), &
2303 work1_recv(:,:,:), &
2304 tagbase, i_sclr, handle, &
2306 flag_dens = .true. )
2311 dens_recv(
k,i,j) = work1_recv(
k,i,j)
2316 call comm_vars8( dens_recv, 1 )
2318 tagbase = tagcomm + tag_momz*order_tag_var
2320 call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), &
2321 work2_recv(:,:,:), &
2322 tagbase, i_zstg, handle, &
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
2335 velz_recv(
ks-1,i,j) = 0.0_rp
2336 velz_recv(
ke ,i,j) = 0.0_rp
2341 call comm_wait ( dens_recv, 1, .false. )
2343 tagbase = tagcomm + tag_momx*order_tag_var
2346 call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), &
2347 work1_recv(:,:,:), &
2348 tagbase, i_xstg, handle, &
2352 call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), &
2354 tagbase, i_sclr, handle, &
2358 tagbase = tagcomm + tag_momy*order_tag_var
2361 call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), &
2362 work2_recv(:,:,:), &
2363 tagbase, i_ystg, handle, &
2367 call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), &
2369 tagbase, i_sclr, handle, &
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
2387 velx_recv(
k,
ia,j) = work1_recv(
k,
ia,j) / dens_recv(
k,
ia,j)
2391 call comm_vars8( velx_recv, 2 )
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
2405 vely_recv(
k,i,
ja) = work2_recv(
k,i,
ja) / dens_recv(
k,i,
ja)
2409 call comm_vars8( vely_recv, 3 )
2411 call comm_wait ( velx_recv, 2, .false. )
2412 call comm_wait ( vely_recv, 3, .false. )
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)
2432 velx_recv(
k,i,j) = ( work1_recv(
k,i+1,j) + work1_recv(
k,i,j) ) * 0.5_rp
2440 velx_recv(
k,
ia,j) = work1_recv(
k,
ia,j)
2444 call comm_vars8( velx_recv, 2 )
2450 vely_recv(
k,i,j) = ( work2_recv(
k,i,j+1) + work2_recv(
k,i,j) ) * 0.5_rp
2458 vely_recv(
k,i,
ja) = work2_recv(
k,i,
ja)
2462 call comm_vars8( vely_recv, 3 )
2464 call comm_wait ( velx_recv, 2, .false. )
2465 call comm_wait ( vely_recv, 3, .false. )
2469 tagbase = tagcomm + tag_rhot*order_tag_var
2470 call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), &
2471 work1_recv(:,:,:), &
2472 tagbase, i_sclr, handle, &
2478 pott_recv(
k,i,j) = work1_recv(
k,i,j) / dens_recv(
k,i,j)
2484 tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
2485 call comm_cartesc_nest_intercomm_nestdown( dummy(:,:,:), &
2486 work1_recv(:,:,:), &
2487 tagbase, i_sclr, handle, &
2493 qtrc_recv(
k,i,j,iq) = work1_recv(
k,i,j)
2503 log_error(
"COMM_CARTESC_NEST_nestdown_recv",*)
'internal error'
2518 integer,
parameter :: handle = 1
2521 integer :: tagbase, tagcomm
2528 tagcomm = intercomm_id(handle) * order_tag_comm
2537 nwait_p = nwait_p + 1
2540 call comm_cartesc_nest_issuer_of_wait( handle )
2542 if ( online_aggressive_comm )
then
2552 log_error(
"COMM_CARTESC_NEST_recvwait_issue_send",*)
'internal error'
2567 integer,
parameter :: handle = 2
2570 integer :: tagbase, tagcomm
2577 tagcomm = intercomm_id(handle) * order_tag_comm
2586 log_info(
"COMM_CARTESC_NEST_recvwait_issue_recv",
'(1X,A,I5,A)')
"NestIDC [C]: que recv ( ", nrecv,
" )"
2594 tagbase = tagcomm + tag_dens*order_tag_var
2595 call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag )
2597 tagbase = tagcomm + tag_momz*order_tag_var
2599 call comm_cartesc_nest_issuer_of_receive( tagbase, i_zstg, handle, isu_tag )
2602 tagbase = tagcomm + tag_momx*order_tag_var
2604 call comm_cartesc_nest_issuer_of_receive( tagbase, i_xstg, handle, isu_tag )
2606 call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag )
2609 tagbase = tagcomm + tag_momy*order_tag_var
2611 call comm_cartesc_nest_issuer_of_receive( tagbase, i_ystg, handle, isu_tag )
2613 call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag )
2616 tagbase = tagcomm + tag_rhot*order_tag_var
2617 call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag )
2620 tagbase = tagcomm + (tag_qx*10+iq)*order_tag_var
2621 call comm_cartesc_nest_issuer_of_receive( tagbase, i_sclr, handle, isu_tag )
2629 log_error(
"COMM_CARTESC_NEST_recvwait_issue_recv",*)
'internal error'
2643 integer,
parameter :: handle = 1
2660 log_error(
"COMM_CARTESC_NEST_recv_cancel_send",*)
'internal error'
2674 integer,
parameter :: handle = 2
2689 log_info(
"COMM_CARTESC_NEST_recv_cancel_recv",
'(1X,A,I5,A)')
"NestIDC [C]: CANCEL recv ( ", nrecv,
" )"
2692 if ( ireq_d(rq) /= mpi_request_null )
then
2694 call mpi_cancel(ireq_d(rq), ierr)
2704 log_error(
"COMM_CARTESC_NEST_recv_cancel_recv",*)
'internal error'
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
2741 logical ,
intent(in),
optional :: flag_dens
2743 integer :: tile_num_x
2745 integer :: ileng, tag, target_rank
2747 integer :: parent_KA
2749 integer :: xloc, yloc
2750 integer :: gxs, gxe, gys, gye
2751 integer :: pxs, pxe, pys, pye
2754 integer :: ig, rq, yp
2757 logical :: logarithmic
2765 logarithmic = .false.
2767 if (
present(flag_dens) )
then
2768 if( flag_dens )
then
2769 logarithmic = .true.
2774 if ( id_stag == i_sclr )
then
2777 elseif( id_stag == i_zstg )
then
2780 elseif( id_stag == i_xstg )
then
2783 elseif( id_stag == i_ystg )
then
2799 target_rank = comm_cartesc_nest_tile_list_yp(yp)
2802 call mpi_isend( pvar, &
2811 dvar(:,:,:) = -1.0_rp
2820 parent_ka = dom_info(i_parent)%KMAX + dom_info(i_parent)%KHALO * 2
2822 * ( dom_info(i_parent)%IMAX + dom_info(i_parent)%IHALO * 2 ) &
2823 * ( dom_info(i_parent)%JMAX + dom_info(i_parent)%JHALO * 2 )
2825 tile_num_x = dom_info(i_parent)%tile_num_x
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
2837 do yp = 1, comm_cartesc_nest_tile_all
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
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
2848 isu_tag = isu_tag + 1
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
2857 buffer_ref_3d(zs:ze,gxs:gxe,gys:gye) = recvbuf_3d(zs:ze,pxs:pxe,pys:pye,isu_tag)
2863 if ( no_zstag )
then
2866 tileal_ia, tileal_ja, &
2871 kgrd(:,:,:,:,:,ig), &
2872 vfact(:, :,:,:,ig), &
2873 buffer_ref_cz(:,:,:), &
2875 buffer_ref_3d(:,:,:), &
2878 logwgt = logarithmic )
2883 tileal_ia, tileal_ja, &
2888 kgrd(:,:,:,:,:,ig), &
2889 vfact(:, :,:,:,ig), &
2890 buffer_ref_fz(:,:,:), &
2892 buffer_ref_3d(:,:,:), &
2895 logwgt = logarithmic )
2900 dvar( 1:
ks-1,i,j) = 0.0_rp
2901 dvar(
ke+1:
ka ,i,j) = 0.0_rp
2906 log_error(
"COMM_CARTESC_NEST_intercomm_nestdown_3D",*)
'internal error'
2928 integer,
intent(in) :: tagbase
2929 integer,
intent(in) :: id_stag
2930 integer,
intent(in) :: HANDLE
2931 integer,
intent(inout) :: isu_tag
2933 integer :: ierr, ileng
2934 integer :: tag, target_rank
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 )
2956 do yp = 1, comm_cartesc_nest_tile_all
2959 target_rank = comm_cartesc_nest_tile_list_d(yp,
prc_myrank+1)
2960 tag = tagbase + call_order(yp)
2962 isu_tag = isu_tag + 1
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
2970 recvbuf_3d(:,:,:,isu_tag) = 0.0_rp
2972 call mpi_irecv( recvbuf_3d(:,:,:,isu_tag), &
2986 log_error(
"COMM_CARTESC_NEST_issuer_of_receive_3D",*)
'internal error'
3001 integer,
intent(in) :: HANDLE
3017 log_error(
"COMM_CARTESC_NEST_issuer_of_wait_3D",*)
'internal error'
3033 integer,
intent(in) :: req_count
3034 integer,
intent(inout) :: ireq(max_rq)
3038 integer :: istatus(MPI_STATUS_SIZE,req_count)
3039 integer :: req_count2
3040 integer :: ireq2(max_rq)
3050 if ( ireq(i) /= mpi_request_null )
then
3051 req_count2 = req_count2 + 1
3052 ireq2(req_count2) = ireq(i)
3056 if( req_count2 /= 0 )
call mpi_waitall( req_count2, ireq2(1:req_count2), istatus, ierr )
3079 integer,
parameter :: handle = 1
3081 integer :: istatus(mpi_status_size)
3093 if ( rq_ctl_p > 0 )
call mpi_test(ireq_p(1), flag, istatus, ierr)
3097 log_error(
"COMM_CARTESC_NEST_test_send",*)
'error'
3111 integer,
parameter :: handle = 2
3113 integer :: istatus(mpi_status_size)
3125 if ( rq_ctl_d > 0 )
call mpi_test(ireq_d(1), flag, istatus, ierr)
3129 log_error(
"COMM_CARTESC_NEST_test_recv",*)
'error'
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 )
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 )
3153 if (
allocated( ireq_p ) )
deallocate( ireq_p )
3154 if (
allocated( ireq_d ) )
deallocate( ireq_d )
3156 if (
allocated( call_order ) )
deallocate( call_order )
3157 if (
allocated( recvbuf_3d ) )
deallocate( recvbuf_3d )
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 )
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 )
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 )