103 mpp_clock_begin, mpp_clock_end, mpp_clock_sync, &
104 comm_tag_1, comm_tag_2, comm_tag_3, comm_tag_4, &
105 comm_tag_5, comm_tag_6, comm_tag_7, comm_tag_8, &
106 comm_tag_9, comm_tag_10
118 domainug, mpp_get_ug_compute_domains, &
119 mpp_get_ug_domains_index, mpp_get_ug_domain_grid_index, &
121 use constants_mod,
only: pi, radius
126 use stock_constants_mod,
only: istock_top, istock_bottom, istock_side, stock_names, &
128 use gradient_mod,
only: gradient_cubic
130 use fms2_io_mod,
only: fmsnetcdfdomainfile_t,
read_data, get_dimension_size
131 use fms2_io_mod,
only: get_variable_units, dimension_exists
132 use platform_mod,
only: r8_kind, i8_kind, fms_file_len
146 integer,
parameter :: FIRST_ORDER = 1
147 integer,
parameter :: SECOND_ORDER = 2
150 integer,
parameter :: max_fields = 80
157 logical :: debug_stocks = .false.
158 logical :: xgrid_clocks_on = .false.
159 logical :: monotonic_exchange = .false.
167 logical :: do_alltoall = .true.
168 logical :: do_alltoallv = .false.
173 monotonic_exchange,
nsubset, do_alltoall, do_alltoallv, &
179 real(r8_kind),
allocatable,
dimension(:,:) ::
area_atm_model, area_lnd_model, area_ocn_model
181 real(r8_kind),
allocatable,
dimension(:,:) ::
area_atm_sphere, area_lnd_sphere, area_ocn_sphere
216 module procedure put_side2_to_xgrid_ug
224 module procedure get_side2_from_xgrid_ug
225 module procedure get_side1_from_xgrid_ug
278 real(r8_kind) :: area
283 real(r8_kind) :: scale
289 real(r8_kind),
dimension(:,:),
pointer :: dx => null()
290 real(r8_kind),
dimension(:,:),
pointer :: dy => null()
291 real(r8_kind),
dimension(:,:),
pointer :: area => null()
292 real(r8_kind),
dimension(:),
pointer :: edge_w => null()
293 real(r8_kind),
dimension(:),
pointer :: edge_e => null()
294 real(r8_kind),
dimension(:),
pointer :: edge_s => null()
295 real(r8_kind),
dimension(:),
pointer :: edge_n => null()
296 real(r8_kind),
dimension(:,:,:),
pointer :: en1 => null()
297 real(r8_kind),
dimension(:,:,:),
pointer :: en2 => null()
298 real(r8_kind),
dimension(:,:,:),
pointer :: vlon => null()
299 real(r8_kind),
dimension(:,:,:),
pointer :: vlat => null()
305 character(len=3) :: id
307 logical :: on_this_pe
309 integer,
pointer,
dimension(:) :: pelist
313 integer,
pointer,
dimension(:) :: tile =>null()
314 integer,
pointer,
dimension(:) :: is =>null()
315 integer,
pointer,
dimension(:) :: ie =>null()
316 integer,
pointer,
dimension(:) :: js =>null()
317 integer,
pointer,
dimension(:) :: je =>null()
318 integer,
pointer :: is_me =>null()
319 integer,
pointer :: ie_me =>null()
320 integer,
pointer :: js_me =>null()
321 integer,
pointer :: je_me =>null()
330 integer,
pointer :: tile_me
334 real(r8_kind),
pointer,
dimension(:) :: lon =>null()
335 real(r8_kind),
pointer,
dimension(:) :: lat =>null()
336 real(r8_kind),
pointer,
dimension(:,:) :: geolon=>null()
337 real(r8_kind),
pointer,
dimension(:,:) :: geolat=>null()
338 real(r8_kind),
pointer,
dimension(:,:,:) :: frac_area =>null()
339 real(r8_kind),
pointer,
dimension(:,:) :: area =>null()
340 real(r8_kind),
pointer,
dimension(:,:) :: area_inv =>null()
347 integer :: size_repro
356 integer,
pointer :: ls_me =>null()
357 integer,
pointer :: le_me =>null()
358 integer,
pointer,
dimension(:) :: ls =>null(), le =>null()
359 integer,
pointer :: gs_me =>null(), ge_me =>null()
360 integer,
pointer,
dimension(:) :: gs =>null(), ge =>null()
361 integer,
pointer,
dimension(:) :: l_index =>null()
370 real(r8_kind) :: area
381 integer :: i, j, l, k, pos
382 real(r8_kind) :: area
391 integer :: buffer_pos
392 integer,
allocatable :: i(:)
393 integer,
allocatable :: j(:)
394 integer,
allocatable :: g(:)
395 integer,
allocatable :: xLoc(:)
396 integer,
allocatable :: tile(:)
397 real(r8_kind),
allocatable :: di(:)
398 real(r8_kind),
allocatable :: dj(:)
404 integer :: nsend, nrecv
405 integer :: sendsize, recvsize
406 integer,
pointer,
dimension(:) :: unpack_ind=>null()
407 type(overlap_type),
pointer,
dimension(:) :: send=>null()
408 type(overlap_type),
pointer,
dimension(:) :: recv=>null()
418 integer :: me, npes, root_pe
419 logical,
pointer,
dimension(:) :: your1my2 =>null()
422 logical,
pointer,
dimension(:) :: your2my1 =>null()
425 integer,
pointer,
dimension(:) :: your2my1_size=>null()
430 type (
grid_type),
pointer,
dimension(:) :: grids =>null()
435 type(
x1_type),
pointer,
dimension(:) :: x1 =>null()
436 type(
x1_type),
pointer,
dimension(:) :: x1_put =>null()
437 type(
x2_type),
pointer,
dimension(:) :: x2 =>null()
438 type(
x2_type),
pointer,
dimension(:) :: x2_get =>null()
440 integer,
pointer,
dimension(:) :: send_count_repro =>null()
441 integer,
pointer,
dimension(:) :: recv_count_repro =>null()
442 integer :: send_count_repro_tot
443 integer :: recv_count_repro_tot
446 integer,
pointer,
dimension(:) :: ind_get1 =>null()
447 integer,
pointer,
dimension(:) :: ind_put1 =>null()
457 #include<file_version.h>
459 real(r8_kind),
parameter :: eps = 1.0e-10_r8_kind
460 real(r8_kind),
parameter :: large_number = 1.e20_r8_kind
461 logical :: module_is_initialized = .false.
462 integer :: id_put_1_to_xgrid_order_1 = 0
463 integer :: id_put_1_to_xgrid_order_2 = 0
464 integer :: id_get_1_from_xgrid = 0
465 integer :: id_get_1_from_xgrid_repro = 0
466 integer :: id_get_2_from_xgrid = 0
467 integer :: id_put_2_to_xgrid = 0
468 integer :: id_setup_xmap = 0
469 integer :: id_load_xgrid1, id_load_xgrid2, id_load_xgrid3
470 integer :: id_load_xgrid4, id_load_xgrid5
471 integer :: id_load_xgrid, id_set_comm, id_regen, id_conservation_check
475 integer :: nnest=0, tile_nest, tile_parent
476 integer :: is_nest=0, ie_nest=0, js_nest=0, je_nest=0
477 integer :: is_parent=0, ie_parent=0, js_parent=0, je_parent=0
489 module procedure stock_move_ug_3d
511 logical function in_box(i, j, is, ie, js, je)
512 integer,
intent(in) :: i, j, is, ie, js, je
514 in_box = (i>=is) .and. (i<=ie) .and. (j>=js) .and. (j<=je)
523 integer,
intent(out) :: remap_method
526 integer :: iunit, ierr, io, out_unit
528 if (module_is_initialized)
return
529 module_is_initialized = .true.
531 read (input_nml_file, xgrid_nml, iostat=io)
539 if (
mpp_pe() == mpp_root_pe() )
write (iunit,nml=xgrid_nml)
544 'MPP_IO is no longer supported. Please remove use_mpp_io from namelists',&
552 remap_method = first_order
553 if( monotonic_exchange )
call error_mesg(
'xgrid_mod', &
554 'xgrid_nml monotonic_exchange must be .false. when interp_method = first_order', fatal)
555 write(out_unit,*)
"NOTE from xgrid_mod: use first_order conservative exchange"
557 if(monotonic_exchange)
then
558 write(out_unit,*)
"NOTE from xgrid_mod: use monotonic second_order conservative exchange"
560 write(out_unit,*)
"NOTE from xgrid_mod: use second_order conservative exchange"
562 remap_method = second_order
565 ' is not a valid namelist option', fatal)
568 if(xgrid_clocks_on)
then
569 id_put_1_to_xgrid_order_1 =
mpp_clock_id(
"put_1_to_xgrid_order_1", flags=mpp_clock_sync)
570 id_put_1_to_xgrid_order_2 =
mpp_clock_id(
"put_1_to_xgrid_order_2", flags=mpp_clock_sync)
571 id_get_1_from_xgrid =
mpp_clock_id(
"get_1_from_xgrid", flags=mpp_clock_sync)
572 id_get_1_from_xgrid_repro =
mpp_clock_id(
"get_1_from_xgrid_repro", flags=mpp_clock_sync)
573 id_get_2_from_xgrid =
mpp_clock_id(
"get_2_from_xgrid", flags=mpp_clock_sync)
574 id_put_2_to_xgrid =
mpp_clock_id(
"put_2_to_xgrid", flags=mpp_clock_sync)
575 id_setup_xmap =
mpp_clock_id(
"setup_xmap", flags=mpp_clock_sync)
578 id_conservation_check =
mpp_clock_id(
"conservation_check")
593 subroutine load_xgrid (xmap, grid, grid_file, grid1_id, grid_id, tile1, tile2, use_higher_order)
596 character(len=*),
intent(in) :: grid_file
597 character(len=3),
intent(in) :: grid1_id, grid_id
598 integer,
intent(in) :: tile1, tile2
599 logical,
intent(in) :: use_higher_order
601 integer,
pointer,
dimension(:) :: i1=>null(), j1=>null()
602 integer,
pointer,
dimension(:) :: i2=>null(), j2=>null()
603 real(r8_kind),
pointer,
dimension(:) :: di=>null(), dj=>null()
604 real(r8_kind),
pointer,
dimension(:) :: area =>null()
605 integer,
pointer,
dimension(:) :: i1_tmp=>null(), j1_tmp=>null()
606 integer,
pointer,
dimension(:) :: i2_tmp=>null(), j2_tmp=>null()
607 real(r8_kind),
pointer,
dimension(:) :: di_tmp=>null(), dj_tmp=>null()
608 real(r8_kind),
pointer,
dimension(:) :: area_tmp =>null()
609 integer,
pointer,
dimension(:) :: i1_side1=>null(), j1_side1=>null()
610 integer,
pointer,
dimension(:) :: i2_side1=>null(), j2_side1=>null()
611 real(r8_kind),
pointer,
dimension(:) :: di_side1=>null(), dj_side1=>null()
612 real(r8_kind),
pointer,
dimension(:) :: area_side1 =>null()
614 real(r8_kind),
allocatable,
dimension(:,:) :: tmp
615 real(r8_kind),
allocatable,
dimension(:) :: send_buffer, recv_buffer
616 type (
grid_type),
pointer,
save :: grid1 =>null()
617 integer :: l, ll, ll_repro, p, nxgrid, size_prev
619 integer :: size_repro, out_unit
620 logical :: scale_exist = .false.
621 logical :: is_distribute = .false.
622 real(r8_kind),
allocatable,
dimension(:) :: scale
623 real(r8_kind) :: garea
624 integer :: npes, isc, iec, nxgrid_local, pe, nxgrid_local_orig
625 integer :: nxgrid1, nxgrid2, nset1, nset2, ndivs, cur_ind
626 integer :: pos, nsend, nrecv, l1, l2, n, mypos
627 integer :: start(4), nread(4)
629 character(len=128) :: attvalue
630 integer,
dimension(0:xmap%npes-1) :: pelist
631 logical,
dimension(0:xmap%npes-1) :: subset_rootpe
632 integer,
dimension(0:xmap%npes-1) :: nsend1, nsend2, nrecv1, nrecv2
633 integer,
dimension(0:xmap%npes-1) :: send_cnt, recv_cnt
634 integer,
dimension(0:xmap%npes-1) :: send_buffer_pos, recv_buffer_pos
635 integer,
dimension(0:xmap%npes-1) :: ibegin, iend, pebegin, peend
636 integer,
dimension(2*xmap%npes) :: ibuf1, ibuf2
637 integer,
dimension(0:xmap%npes-1) :: pos_x, y2m1_size
638 integer,
allocatable,
dimension(:) :: y2m1_pe
639 integer,
pointer,
save :: iarray(:), jarray(:)
640 integer,
allocatable,
save :: pos_s(:)
641 integer,
pointer,
dimension(:) :: iarray2(:)=>null(), jarray2(:)=>null()
643 integer :: nxgrid1_old
645 type(fmsnetcdffile_t) :: fileobj
647 if(.not.
open_file(fileobj, grid_file,
'read' ))
then
648 call error_mesg(
'xgrid_mod(load_xgrid)',
'Error in opening file '//trim(grid_file), fatal)
651 scale_exist = .false.
652 grid1 => xmap%grids(1)
656 mypos =
mpp_pe()-mpp_root_pe()
658 call mpp_get_current_pelist(pelist)
660 if( npes .NE. pelist(npes-1) - pelist(0) + 1 )
then
661 print*,
"npes =", npes,
", pelist(npes-1)=", pelist(npes-1),
", pelist(0)=", pelist(0)
662 call error_mesg(
'xgrid_mod', .NE.
'npes pelist(npes-1) - pelist(0)', fatal)
665 select case(xmap%version)
668 if (dimension_exists(fileobj,
'i_'//lowercase(grid1_id)//
'X'//lowercase(grid_id)))
then
669 call get_dimension_size(fileobj,
'i_'//lowercase(grid1_id)//
'X'//lowercase(grid_id), nxgrid)
671 if(nxgrid .LE. 0)
return
675 if(nxgrid .LE. 0)
return
679 if(nxgrid > npes)
then
683 if(npes == ndivs)
then
684 p =
mpp_pe()-mpp_root_pe()
687 subset_rootpe(:) = .true.
692 if(pe == pebegin(n))
then
699 subset_rootpe(:) = .false.
702 if(pelist(n) == pebegin(cur_ind))
then
703 subset_rootpe(n) = .true.
705 if(cur_ind == ndivs)
exit
709 is_distribute = .true.
711 is_distribute = .false.
712 isc = 1; iec = nxgrid
717 if(use_higher_order)
then
721 if(scale_exist) nset2 = nset1 + 1
723 call mpp_clock_begin(id_load_xgrid1)
724 if(iec .GE. isc)
then
725 nxgrid_local = iec - isc + 1
726 allocate(i1_tmp(isc:iec), j1_tmp(isc:iec), i2_tmp(isc:iec), j2_tmp(isc:iec), area_tmp(isc:iec) )
727 if(use_higher_order)
allocate(di_tmp(isc:iec), dj_tmp(isc:iec))
731 select case(xmap%version)
733 start(1) = isc; nread(1) = nxgrid_local
734 allocate(tmp(nxgrid_local,1))
735 call read_data(fileobj,
'I_'//grid1_id//
'_'//grid1_id//
'x'//grid_id, tmp, corner=start, edge_lengths=nread)
736 i1_tmp = int(tmp(:,1))
737 call read_data(fileobj,
'J_'//grid1_id//
'_'//grid1_id//
'x'//grid_id, tmp, corner=start, edge_lengths=nread)
738 j1_tmp = int(tmp(:,1))
739 call read_data(fileobj,
'I_'//grid_id//
'_'//grid1_id//
'x'//grid_id, tmp, corner=start, edge_lengths=nread)
740 i2_tmp = int(tmp(:,1))
741 call read_data(fileobj,
'J_'//grid_id//
'_'//grid1_id//
'x'//grid_id, tmp, corner=start, edge_lengths=nread)
742 j2_tmp = int(tmp(:,1))
743 call read_data(fileobj,
'AREA_'//grid1_id//
'x'//grid_id, tmp, corner=start, edge_lengths=nread)
745 if(use_higher_order)
then
746 call read_data(fileobj,
'DI_'//grid1_id//
'x'//grid_id, tmp, corner=start, edge_lengths=nread)
748 call read_data(fileobj,
'DJ_'//grid1_id//
'x'//grid_id, tmp, corner=start, edge_lengths=nread)
753 nread(1) = 2; start(2) = isc; nread(2) = nxgrid_local
754 allocate(tmp(2, isc:iec))
755 call read_data(fileobj,
"tile1_cell", tmp, corner=start, edge_lengths=nread)
756 i1_tmp(isc:iec) = int(tmp(1, isc:iec))
757 j1_tmp(isc:iec) = int(tmp(2, isc:iec))
758 call read_data(fileobj,
"tile2_cell", tmp, corner=start, edge_lengths=nread)
759 i2_tmp(isc:iec) = int(tmp(1, isc:iec))
760 j2_tmp(isc:iec) = int(tmp(2, isc:iec))
761 if(use_higher_order)
then
762 call read_data(fileobj,
"tile1_distance", tmp, corner=start, edge_lengths=nread)
763 di_tmp(isc:iec) = tmp(1, isc:iec)
764 dj_tmp(isc:iec) = tmp(2, isc:iec)
767 start(1) = isc; nread(1) = nxgrid_local
769 allocate(tmp(isc:iec,1) )
771 call read_data(fileobj,
"xgrid_area", tmp(:,1:1), corner=start, edge_lengths=nread)
773 call get_variable_units(fileobj,
"xgrid_area", attvalue)
775 if( trim(attvalue) ==
'm2' )
then
776 garea = 4.0_r8_kind * pi * radius * radius;
777 area_tmp = tmp(:,1)/garea
778 else if( trim(attvalue) ==
'none' )
then
781 call error_mesg(
'xgrid_mod',
'In file '//trim(grid_file)//
', xgrid_area units = '// &
782 trim(attvalue)//
' should be "m2" or "none"', fatal)
787 if(grid1_id ==
'LND' .AND. grid_id ==
'OCN')
then
788 if(variable_exists(fileobj,
"scale"))
then
789 allocate(scale(isc:iec))
790 write(out_unit, *)
"NOTE from load_xgrid(xgrid_mod): field 'scale' exist in the file "// &
791 & trim(grid_file)//
", this field will be read and the exchange grid cell area will be"// &
792 &
" multiplied by scale"
793 call read_data(fileobj,
"scale", tmp, corner=start, edge_lengths=nread)
803 nxgrid_local_orig = nxgrid_local
804 allocate(i1(isc:iec), j1(isc:iec), i2(isc:iec), j2(isc:iec), area(isc:iec) )
805 if(use_higher_order)
allocate(di(isc:iec), dj(isc:iec))
811 if(grid1%tile(p) == tile1)
then
812 if(
in_box_nbr(i1_tmp(l), j1_tmp(l), grid1, p))
then
821 if(grid%tile(p) == tile2)
then
822 if (
in_box_nbr(i2_tmp(l), j2_tmp(l), grid, p))
then
828 area(pos) = area_tmp(l)
829 if(use_higher_order)
then
840 deallocate(i1_tmp, i2_tmp, j1_tmp, j2_tmp, area_tmp)
841 if(use_higher_order)
deallocate( di_tmp, dj_tmp)
843 if(iec .GE. isc)
then
844 nxgrid_local = iec - isc + 1
850 nxgrid_local_orig = 0
855 call mpp_clock_end(id_load_xgrid1)
857 if(is_distribute)
then
861 nsend1(:) = 0; nrecv1(:) = 0
862 nsend2(:) = 0; nrecv2(:) = 0
863 ibuf1(:)= 0; ibuf2(:)= 0
865 call mpp_clock_begin(id_load_xgrid2)
866 if(nxgrid_local>0)
then
867 allocate( send_buffer(nxgrid_local * (nset1+nset2)) )
870 send_buffer_pos(p) = pos
871 if(grid%tile(p) == tile2)
then
874 nsend2(p) = nsend2(p) + 1
875 send_buffer(pos+1) = real(i1(l), r8_kind)
876 send_buffer(pos+2) = real(j1(l), r8_kind)
877 send_buffer(pos+3) = real(i2(l), r8_kind)
878 send_buffer(pos+4) = real(j2(l), r8_kind)
879 send_buffer(pos+5) = area(l)
880 if(use_higher_order)
then
881 send_buffer(pos+6) = di(l)
882 send_buffer(pos+7) = dj(l)
884 if(scale_exist) send_buffer(pos+nset2) = scale(l)
889 if(grid1%tile(p) == tile1)
then
892 nsend1(p) = nsend1(p) + 1
893 send_buffer(pos+1) = real(i1(l), r8_kind)
894 send_buffer(pos+2) = real(j1(l), r8_kind)
895 send_buffer(pos+3) = real(i2(l), r8_kind)
896 send_buffer(pos+4) = real(j2(l), r8_kind)
897 send_buffer(pos+5) = area(l)
898 if(use_higher_order)
then
899 send_buffer(pos+6) = di(l)
900 send_buffer(pos+7) = dj(l)
908 call mpp_clock_end(id_load_xgrid2)
911 call mpp_clock_begin(id_load_xgrid3)
913 if (do_alltoall)
then
915 ibuf1(2*p+1) = nsend1(p)
916 ibuf1(2*p+2) = nsend2(p)
921 p = mod(mypos+npes-n, npes)
922 if(.not. subset_rootpe(p)) cycle
923 call mpp_recv( ibuf2(2*p+1), glen=2, from_pe=pelist(p), block=.false., tag=comm_tag_1)
926 if(nxgrid_local_orig>0)
then
928 p = mod(mypos+n, npes)
929 ibuf1(2*p+1) = nsend1(p)
930 ibuf1(2*p+2) = nsend2(p)
931 call mpp_send( ibuf1(2*p+1), plen=2, to_pe=pelist(p), tag=comm_tag_1)
937 nrecv1(p) = ibuf2(2*p+1)
938 nrecv2(p) = ibuf2(2*p+2)
942 call mpp_clock_end(id_load_xgrid3)
943 call mpp_clock_begin(id_load_xgrid4)
946 recv_buffer_pos(p) = pos
947 pos = pos + nrecv1(p) * nset1 + nrecv2(p) * nset2
951 nxgrid1 = sum(nrecv1)
952 nxgrid2 = sum(nrecv2)
953 if(nxgrid1>0 .OR. nxgrid2>0)
allocate(recv_buffer(nxgrid1*nset1+nxgrid2*nset2))
955 if (do_alltoallv)
then
957 send_cnt(:) = nset1 * nsend1(:) + nset2 * nsend2(:)
958 recv_cnt(:) = nset1 * nrecv1(:) + nset2 * nrecv2(:)
960 call mpp_alltoall(send_buffer, send_cnt, send_buffer_pos, &
961 recv_buffer, recv_cnt, recv_buffer_pos)
964 p = mod(mypos+npes-n, npes)
965 nrecv = nrecv1(p)*nset1+nrecv2(p)*nset2
967 pos = recv_buffer_pos(p)
968 call mpp_recv(recv_buffer(pos+1), glen=nrecv, from_pe=pelist(p), &
969 block=.false., tag=comm_tag_2)
973 p = mod(mypos+n, npes)
974 nsend = nsend1(p)*nset1 + nsend2(p)*nset2
976 pos = send_buffer_pos(p)
977 call mpp_send(send_buffer(pos+1), plen=nsend, to_pe=pelist(p), &
982 call mpp_clock_end(id_load_xgrid4)
984 if( nxgrid_local>0)
then
985 deallocate(i1,j1,i2,j2,area)
988 allocate(i1(nxgrid2), j1(nxgrid2))
989 allocate(i2(nxgrid2), j2(nxgrid2))
990 allocate(area(nxgrid2))
991 allocate(i1_side1(nxgrid1), j1_side1(nxgrid1))
992 allocate(i2_side1(nxgrid1), j2_side1(nxgrid1))
993 allocate(area_side1(nxgrid1))
994 if(use_higher_order)
then
995 if(nxgrid_local>0)
deallocate(di,dj)
996 allocate(di(nxgrid2), dj(nxgrid2))
997 allocate(di_side1(nxgrid1), dj_side1(nxgrid1))
1000 if(nxgrid_local>0)
deallocate(scale)
1001 allocate(scale(nxgrid2))
1008 i1(l2) = int(recv_buffer(pos+1))
1009 j1(l2) = int(recv_buffer(pos+2))
1010 i2(l2) = int(recv_buffer(pos+3))
1011 j2(l2) = int(recv_buffer(pos+4))
1012 area(l2) = recv_buffer(pos+5)
1013 if(use_higher_order)
then
1014 di(l2) = recv_buffer(pos+6)
1015 dj(l2) = recv_buffer(pos+7)
1017 if(scale_exist)scale(l2) = recv_buffer(pos+nset2)
1022 i1_side1(l1) = int(recv_buffer(pos+1))
1023 j1_side1(l1) = int(recv_buffer(pos+2))
1024 i2_side1(l1) = int(recv_buffer(pos+3))
1025 j2_side1(l1) = int(recv_buffer(pos+4))
1026 area_side1(l1) = recv_buffer(pos+5)
1027 if(use_higher_order)
then
1028 di_side1(l1) = recv_buffer(pos+6)
1029 dj_side1(l1) = recv_buffer(pos+7)
1035 if(
allocated(send_buffer))
deallocate(send_buffer)
1036 if(
allocated(recv_buffer))
deallocate(recv_buffer)
1041 i1_side1 => i1; j1_side1 => j1
1042 i2_side1 => i2; j2_side1 => j2
1044 if(use_higher_order)
then
1050 call mpp_clock_begin(id_load_xgrid5)
1053 size_prev = grid%size
1055 if(grid%tile_me == tile2)
then
1057 if (
in_box_me(i2(l), j2(l), grid) )
then
1058 grid%size = grid%size + 1
1060 if( grid1_id .NE.
"ATM" .OR. tile1 .NE. tile_parent .OR. &
1061 .NOT.
in_box(i1(l), j1(l), is_parent, ie_parent, js_parent, je_parent) )
then
1063 lll = grid%l_index((j2(l)-1)*grid%im+i2(l))
1064 grid%area(lll,1) = grid%area(lll,1)+area(l)
1066 grid%area(i2(l),j2(l)) = grid%area(i2(l),j2(l))+area(l)
1070 if(grid1%tile(p) == tile1)
then
1072 xmap%your1my2(p) = .true.
1080 if(grid%size > size_prev)
then
1081 if(size_prev > 0)
then
1082 allocate(x_local(size_prev))
1084 if(
ASSOCIATED(grid%x))
deallocate(grid%x)
1085 allocate( grid%x( grid%size ) )
1086 grid%x(1:size_prev) = x_local
1089 if(
ASSOCIATED(grid%x))
deallocate(grid%x)
1090 allocate( grid%x( grid%size ) )
1091 grid%x%di = 0.0_r8_kind; grid%x%dj = 0.0_r8_kind
1096 if( grid%tile_me == tile2 )
then
1101 grid%x(ll)%i1 = i1(l); grid%x(ll)%i2 = i2(l)
1102 grid%x(ll)%j1 = j1(l); grid%x(ll)%j2 = j2(l)
1104 grid%x(ll)%l2 = grid%l_index((j2(l)-1)*grid%im + i2(l))
1109 grid%x(ll)%tile = tile1
1110 grid%x(ll)%area = area(l)
1111 if(scale_exist)
then
1112 grid%x(ll)%scale = scale(l)
1114 grid%x(ll)%scale = 1.0_r8_kind
1116 if(use_higher_order)
then
1117 grid%x(ll)%di = di(l)
1118 grid%x(ll)%dj = dj(l)
1123 if(grid1%tile(p) == tile1)
then
1125 grid%x(ll)%pe = p + xmap%root_pe
1134 if(grid%id == xmap%grids(
size(xmap%grids(:)))%id)
then
1141 if(grid1%tile_me == tile1)
then
1142 if(
associated(iarray))
then
1143 nxgrid1_old =
size(iarray(:))
1148 allocate(y2m1_pe(nxgrid1))
1149 if(.not. last_grid )
allocate(pos_s(0:xmap%npes-1))
1151 if(nxgrid1_old > 0)
then
1153 y2m1_size(p) = xmap%your2my1_size(p)
1160 if (
in_box_me(i1_side1(l), j1_side1(l), grid1) )
then
1161 if(grid1%is_ug)
then
1162 lll = grid1%l_index((j1_side1(l)-1)*grid1%im+i1_side1(l))
1163 grid1%area(lll,1) = grid1%area(lll,1) + area_side1(l)
1165 grid1%area(i1_side1(l),j1_side1(l)) = grid1%area(i1_side1(l),j1_side1(l))+area_side1(l)
1168 if (grid%tile(p) == tile2)
then
1169 if (
in_box_nbr(i2_side1(l), j2_side1(l), grid, p))
then
1170 xmap%your2my1(p) = .true.
1172 y2m1_size(p) = y2m1_size(p) + 1
1176 size_repro = size_repro + 1
1181 pos_x(p) = pos_x(p-1) + y2m1_size(p-1)
1184 if(.not. last_grid) pos_s(:) = pos_x(:)
1186 if(nxgrid1_old > 0)
then
1187 y2m1_size(:) = xmap%your2my1_size(:)
1190 allocate(iarray(nxgrid1+nxgrid1_old), jarray(nxgrid1+nxgrid1_old))
1193 do n = 1, xmap%your2my1_size(p)
1194 iarray(pos_x(p)+n) = iarray2(pos_s(p)+n)
1195 jarray(pos_x(p)+n) = jarray2(pos_s(p)+n)
1198 deallocate(iarray2, jarray2)
1200 allocate(iarray(nxgrid1), jarray(nxgrid1))
1210 if(y2m1_size(p) > 0)
then
1211 pos = pos_x(p)+y2m1_size(p)
1212 if( i1_side1(l) == iarray(pos) .AND. j1_side1(l) == jarray(pos) )
then
1216 do n = 1, y2m1_size(p)
1218 if(i1_side1(l) == iarray(pos) .AND. j1_side1(l) == jarray(pos))
then
1225 if( (.NOT. found) .OR. monotonic_exchange )
then
1226 y2m1_size(p) = y2m1_size(p)+1
1227 pos = pos_x(p)+y2m1_size(p)
1228 iarray(pos) = i1_side1(l)
1229 jarray(pos) = j1_side1(l)
1232 xmap%your2my1_size(:) = y2m1_size(:)
1235 deallocate(iarray, jarray)
1236 if(
allocated(pos_s))
deallocate(pos_s)
1240 if (grid1%tile_me == tile1 .and. size_repro > 0)
then
1241 ll_repro = grid%size_repro
1242 grid%size_repro = ll_repro + size_repro
1243 if(ll_repro > 0)
then
1244 allocate(x_local(ll_repro))
1245 x_local = grid%x_repro
1246 if(
ASSOCIATED(grid%x_repro))
deallocate(grid%x_repro)
1247 allocate( grid%x_repro(grid%size_repro ) )
1248 grid%x_repro(1:ll_repro) = x_local
1251 if(
ASSOCIATED(grid%x_repro))
deallocate(grid%x_repro)
1252 allocate( grid%x_repro( grid%size_repro ) )
1253 grid%x_repro%di = 0.0_r8_kind; grid%x_repro%dj = 0.0_r8_kind
1256 if (
in_box_me(i1_side1(l),j1_side1(l), grid1) )
then
1257 ll_repro = ll_repro + 1
1258 grid%x_repro(ll_repro)%i1 = i1_side1(l); grid%x_repro(ll_repro)%i2 = i2_side1(l)
1259 grid%x_repro(ll_repro)%j1 = j1_side1(l); grid%x_repro(ll_repro)%j2 = j2_side1(l)
1260 if(grid1%is_ug)
then
1261 grid%x_repro(ll_repro)%l1 = grid1%l_index((j1_side1(l)-1)*grid1%im+i1_side1(l))
1266 grid%x_repro(ll_repro)%tile = tile1
1267 grid%x_repro(ll_repro)%area = area_side1(l)
1268 if(use_higher_order)
then
1269 grid%x_repro(ll_repro)%di = di_side1(l)
1270 grid%x_repro(ll_repro)%dj = dj_side1(l)
1274 if(grid%tile(p) == tile2)
then
1275 if (
in_box_nbr(i2_side1(l), j2_side1(l), grid, p))
then
1276 grid%x_repro(ll_repro)%pe = p + xmap%root_pe
1284 deallocate(i1, j1, i2, j2, area)
1285 if(use_higher_order)
deallocate(di, dj)
1286 if(scale_exist)
deallocate(scale)
1287 if(is_distribute)
then
1288 deallocate(i1_side1, j1_side1, i2_side1, j2_side1, area_side1)
1289 if(use_higher_order)
deallocate(di_side1, dj_side1)
1292 i1=>null(); j1=>null(); i2=>null(); j2=>null()
1293 call mpp_clock_end(id_load_xgrid5)
1297 end subroutine load_xgrid
1306 character(len=3),
intent(in) :: grid_id
1307 character(len=*),
intent(in) :: grid_file
1309 real(r8_kind),
dimension(grid%im) :: lonb
1310 real(r8_kind),
dimension(grid%jm) :: latb
1311 real(r8_kind) :: d2r
1312 integer :: is, ie, js, je
1313 type(fmsnetcdfdomainfile_t) :: fileobj
1315 d2r = pi / 180.0_r8_kind
1317 if(.not.
open_file(fileobj, grid_file,
'read', grid%domain) )
then
1318 call error_mesg(
'xgrid_mod(get_grid_version1)',
'Error in opening file '//trim(grid_file), fatal)
1322 if (
associated(grid%lon))
deallocate(grid%lon)
1323 if (
associated(grid%lat))
deallocate(grid%lat)
1324 allocate(grid%lon(grid%im), grid%lat(grid%jm))
1325 if(grid_id ==
'ATM')
then
1337 else if(grid_id ==
'LND')
then
1340 if(.not.
allocated(area_lnd_model))
then
1341 allocate(area_lnd_model(is:ie, js:je))
1344 if(.not.
allocated(area_lnd_sphere))
then
1345 allocate(area_lnd_sphere(is:ie, js:je))
1348 else if(grid_id ==
'OCN' )
then
1349 if(.not.
allocated(area_ocn_sphere))
then
1350 allocate(area_ocn_sphere(is:ie, js:je))
1355 if(grid_id ==
'LND' .or. grid_id ==
'ATM')
then
1356 grid%lon = lonb * d2r
1357 grid%lat = latb * d2r
1359 grid%is_latlon = .true.
1374 character(len=3),
intent(in) :: grid_id
1375 character(len=*),
intent(in) :: grid_file
1377 real(r8_kind),
allocatable :: tmpx(:,:), tmpy(:,:)
1378 real(r8_kind) :: d2r
1379 integer :: is, ie, js, je, nlon, nlat, i, j
1380 integer :: start(4), nread(4), isc2, iec2, jsc2, jec2
1381 type(fmsnetcdffile_t) :: fileobj
1383 if(.not.
open_file(fileobj, grid_file,
'read') )
then
1384 call error_mesg(
'xgrid_mod(get_grid_version2)',
'Error in opening file '//trim(grid_file), fatal)
1387 d2r = pi / 180.0_r8_kind
1391 call get_dimension_size(fileobj,
"nx", nlon)
1392 call get_dimension_size(fileobj,
"ny", nlat)
1393 if( mod(nlon,2) .NE. 0)
call error_mesg(
'xgrid_mod', &
1394 'flux_exchange_mod: atmos supergrid longitude size can not be divided by 2', fatal)
1395 if( mod(nlat,2) .NE. 0)
call error_mesg(
'xgrid_mod', &
1396 'flux_exchange_mod: atmos supergrid latitude size can not be divided by 2', fatal)
1399 if(nlon .NE. grid%im .OR. nlat .NE. grid%jm)
call error_mesg(
'xgrid_mod', &
1400 'grid size in tile_file does not match the global grid size', fatal)
1402 if( grid_id ==
'LND' .or. grid_id ==
'ATM' .or. grid_id ==
'WAV' )
then
1403 isc2 = 2*grid%is_me-1; iec2 = 2*grid%ie_me+1
1404 jsc2 = 2*grid%js_me-1; jec2 = 2*grid%je_me+1
1405 allocate(tmpx(isc2:iec2, jsc2:jec2) )
1406 allocate(tmpy(isc2:iec2, jsc2:jec2) )
1407 start = 1; nread = 1
1408 start(1) = isc2; nread(1) = iec2 - isc2 + 1
1409 start(2) = jsc2; nread(2) = jec2 - jsc2 + 1
1410 call read_data(fileobj,
'x', tmpx, corner=start, edge_lengths=nread)
1411 call read_data(fileobj,
'y', tmpy, corner=start, edge_lengths=nread)
1412 if(is_lat_lon(tmpx, tmpy) )
then
1413 deallocate(tmpx, tmpy)
1414 start = 1; nread = 1
1415 start(2) = 2; nread(1) = nlon*2+1
1416 allocate(tmpx(nlon*2+1, 1), tmpy(1, nlat*2+1))
1417 call read_data(fileobj,
"x", tmpx, corner=start, edge_lengths=nread)
1418 if (
associated(grid%lon))
deallocate(grid%lon)
1419 if (
associated(grid%lat))
deallocate(grid%lat)
1420 allocate(grid%lon(grid%im), grid%lat(grid%jm))
1422 grid%lon(i) = tmpx(2*i,1) * d2r
1424 start = 1; nread = 1
1425 start(1) = 2; nread(2) = nlat*2+1
1426 call read_data(fileobj,
"y", tmpy, corner=start, edge_lengths=nread)
1428 grid%lat(j) = tmpy(1, 2*j) * d2r
1430 grid%is_latlon = .true.
1432 if (
associated(grid%geolon))
deallocate(grid%geolon)
1433 if (
associated(grid%geolat))
deallocate(grid%geolat)
1434 allocate(grid%geolon(grid%isd_me:grid%ied_me, grid%jsd_me:grid%jed_me))
1435 allocate(grid%geolat(grid%isd_me:grid%ied_me, grid%jsd_me:grid%jed_me))
1436 grid%geolon = 1.0e10_r8_kind
1437 grid%geolat = 1.0e10_r8_kind
1439 do j = grid%js_me,grid%je_me
1440 do i = grid%is_me,grid%ie_me
1441 grid%geolon(i, j) = tmpx(i*2,j*2)*d2r
1442 grid%geolat(i, j) = tmpy(i*2,j*2)*d2r
1447 grid%is_latlon = .false.
1449 deallocate(tmpx, tmpy)
1461 type(fmsnetcdfdomainfile_t),
intent(in) :: fileobj
1462 character(len=*),
intent(in) :: name
1463 real(r8_kind),
intent(out) :: get_area_data(:,:)
1465 if(variable_exists(fileobj, name))
then
1466 call read_data(fileobj, name, get_area_data)
1468 call error_mesg(
'xgrid_mod',
'no field named '//trim(name)//
' in grid file '//trim(fileobj%path)// &
1469 ' Will set data to negative values...', note)
1471 get_area_data = -1.0_r8_kind
1485 type(
domain2d),
intent(in) :: domain
1486 character(len=*),
intent(in) :: grid_file
1487 integer :: is, ie, js, je
1488 type(fmsnetcdffile_t) :: fileobj
1490 if(
allocated(area_ocn_model))
return
1495 allocate(area_ocn_model(is:ie, js:je))
1496 if(ie < is .or. je < js )
return
1498 if(.not.
open_file(fileobj, grid_file,
'read') )
then
1499 call error_mesg(
'xgrid_mod(get_ocean_model_area_elements)',
'Error in opening file '//trim(grid_file), fatal)
1502 if(variable_exists(fileobj,
'AREA_OCN_MODEL') )
then
1503 call read_data(fileobj,
'AREA_OCN_MODEL', area_ocn_model)
1505 deallocate(area_ocn_model)
1516 subroutine setup_xmap(xmap, grid_ids, grid_domains, grid_file, atm_grid, lnd_ug_domain)
1518 character(len=3),
dimension(:),
intent(in ) :: grid_ids
1519 type(
domain2d),
dimension(:),
intent(in ) :: grid_domains
1520 character(len=*),
intent(in ) :: grid_file
1522 type(
domainug),
optional,
intent(in ) :: lnd_ug_domain
1525 integer :: nxgrid_file, i1, i2, i3, tile1, tile2, j
1526 integer :: nxc, nyc, out_unit
1527 type(
grid_type),
pointer :: grid => null()
1528 type(
grid_type),
pointer,
save :: grid1 => null()
1529 real(r8_kind),
dimension(3) :: xxx
1530 real(r8_kind),
dimension(:,:),
allocatable :: check_data
1531 real(r8_kind),
dimension(:,:,:),
allocatable :: check_data_3d
1532 real(r8_kind),
allocatable :: tmp_2d(:,:), tmp_3d(:,:,:)
1533 character(len=FMS_FILE_LEN) :: xgrid_file, xgrid_name
1534 character(len=FMS_FILE_LEN) :: tile_file, mosaic_file
1535 character(len=256) :: mosaic1, mosaic2, contact, xgrid_dimname
1536 character(len=256) :: tile1_name, tile2_name
1537 character(len=256),
allocatable :: tile1_list(:), tile2_list(:)
1538 character(len=FMS_FILE_LEN),
allocatable :: xgrid_filelist(:)
1539 integer :: npes, npes2
1540 integer,
allocatable :: pelist(:)
1542 logical :: use_higher_order = .false.
1543 integer :: lnd_ug_id, l
1544 integer,
allocatable :: grid_index(:)
1545 type(fmsnetcdffile_t) :: gridfileobj, mosaicfileobj, fileobj
1546 type(
grid_type),
allocatable,
target :: grids_tmp(:)
1549 call mpp_clock_begin(id_setup_xmap)
1551 if(
interp_method .ne.
'first_order') use_higher_order = .true.
1556 xmap%root_pe = mpp_root_pe()
1558 if (
associated(xmap%grids))
deallocate(xmap%grids)
1559 allocate( xmap%grids(1:
size(grid_ids(:))) )
1561 if (
associated(xmap%your1my2))
deallocate(xmap%your1my2)
1562 if (
associated(xmap%your2my1))
deallocate(xmap%your2my1)
1563 if (
associated(xmap%your2my1_size))
deallocate(xmap%your2my1_size)
1564 allocate ( xmap%your1my2(0:xmap%npes-1), xmap%your2my1(0:xmap%npes-1) )
1565 allocate ( xmap%your2my1_size(0:xmap%npes-1) )
1567 xmap%your1my2 = .false.; xmap%your2my1 = .false.;
1568 xmap%your2my1_size = 0
1570 if(.not.
open_file(gridfileobj,trim(grid_file),
"read"))
then
1571 call error_mesg(
'xgrid_mod',
'Error when opening file'//trim(grid_file), fatal)
1575 if(variable_exists(gridfileobj,
"AREA_ATMxOCN" ) )
then
1578 else if(variable_exists(gridfileobj,
"ocn_mosaic_file" ) )
then
1581 call error_mesg(
'xgrid_mod',
'both AREA_ATMxOCN and ocn_mosaic_file does not exist in '//trim(grid_file), fatal)
1586 call error_mesg(
'xgrid_mod',
'reading exchange grid information from grid spec file', note)
1588 call error_mesg(
'xgrid_mod',
'reading exchange grid information from mosaic grid file', note)
1593 if(
present(lnd_ug_domain))
then
1594 do g=1,
size(grid_ids(:))
1595 if(grid_ids(g) ==
'LND') lnd_ug_id = g
1599 call mpp_clock_begin(id_load_xgrid)
1603 grids_tmp = xmap%grids
1605 grid1 => xmap%grids(1)
1607 do g=1,
size(grid_ids(:))
1609 grid => grids_tmp(g)
1611 grid%id = grid_ids(g)
1612 grid%domain = grid_domains(g)
1614 if (
associated(grid%is))
deallocate(grid%is)
1615 if (
associated(grid%ie))
deallocate(grid%ie)
1616 if (
associated(grid%js))
deallocate(grid%js)
1617 if (
associated(grid%je))
deallocate(grid%je)
1618 if (
associated(grid%tile))
deallocate(grid%tile)
1619 allocate ( grid%is(0:xmap%npes-1), grid%ie(0:xmap%npes-1) )
1620 allocate ( grid%js(0:xmap%npes-1), grid%je(0:xmap%npes-1) )
1621 allocate ( grid%tile(0:xmap%npes-1) )
1631 select case(xmap%version)
1635 call read_data(gridfileobj, lowercase(grid_ids(g))//
'_mosaic_file', mosaic_file)
1636 if(.not.
open_file(mosaicfileobj,
'INPUT/'//trim(mosaic_file),
"read"))
then
1637 call error_mesg(
'xgrid_mod',
'Error when opening solo mosaic file INPUT/'//trim(mosaic_file), fatal)
1639 call get_dimension_size(mosaicfileobj,
'ntiles', grid%ntile)
1642 if( g == 1 .AND. grid_ids(1) ==
'ATM' )
then
1643 if( .NOT. grid%on_this_pe )
call error_mesg(
'xgrid_mod',
'ATM domain is not defined on some processor' ,fatal)
1646 if( xmap%npes > grid%npes .AND. g == 1 .AND. grid_ids(1) ==
'ATM' )
then
1648 else if(xmap%npes > grid%npes)
then
1654 allocate(grid%pelist(0:npes-1))
1658 call mpp_get_data_domain(grid%domain, grid%isd_me, grid%ied_me, grid%jsd_me, grid%jed_me, &
1659 xsize=grid%nxd_me, ysize=grid%nyd_me)
1662 if( grid%root_pe == xmap%root_pe )
then
1664 ybegin=grid%js(0:npes-1), yend=grid%je(0:npes-1) )
1666 if( xmap%npes > npes .AND. g == 1 .AND. grid_ids(1) ==
'ATM' )
then
1668 ybegin=grid%js(npes:xmap%npes-1), yend=grid%je(npes:xmap%npes-1) )
1672 npes2 = xmap%npes-npes
1674 ybegin=grid%js(0:npes2-1), yend=grid%je(0:npes2-1) )
1675 call mpp_get_compute_domains(grid%domain, xbegin=grid%is(npes2:xmap%npes-1), xend=grid%ie(npes2:xmap%npes-1), &
1676 ybegin=grid%js(npes2:xmap%npes-1), yend=grid%je(npes2:xmap%npes-1) )
1680 if( xmap%npes > grid%npes .AND. g == 1 .AND. grid_ids(1) ==
'ATM' )
then
1684 if( g == 1 .AND. grid_ids(1) ==
'ATM' ) npes = xmap%npes
1686 if(grid%tile(p) > grid%ntile .or. grid%tile(p) < 1)
call error_mesg(
'xgrid_mod', &
1687 'tile id should between 1 and ntile', fatal)
1695 grid%is_me => grid%is(xmap%me-xmap%root_pe); grid%ie_me => grid%ie(xmap%me-xmap%root_pe)
1696 grid%js_me => grid%js(xmap%me-xmap%root_pe); grid%je_me => grid%je(xmap%me-xmap%root_pe)
1697 grid%nxc_me = grid%ie_me - grid%is_me + 1
1698 grid%nyc_me = grid%je_me - grid%js_me + 1
1699 grid%tile_me => grid%tile(xmap%me-xmap%root_pe)
1702 grid%is_ug = .false.
1704 if( g == lnd_ug_id )
then
1706 'does not support unstructured grid for VERSION1 grid' ,fatal)
1708 grid%ug_domain = lnd_ug_domain
1709 if (
associated(grid%ls))
deallocate(grid%ls)
1710 if (
associated(grid%le))
deallocate(grid%le)
1711 if (
associated(grid%gs))
deallocate(grid%gs)
1712 if (
associated(grid%ge))
deallocate(grid%ge)
1713 allocate ( grid%ls(0:xmap%npes-1), grid%le(0:xmap%npes-1) )
1714 allocate ( grid%gs(0:xmap%npes-1), grid%ge(0:xmap%npes-1) )
1719 if(xmap%npes > grid%npes)
then
1722 call mpp_get_ug_compute_domains(grid%ug_domain, begin=grid%ls(0:npes-1),
end=grid%le(0:npes-1) )
1723 call mpp_get_ug_domains_index(grid%ug_domain, grid%gs(0:npes-1), grid%ge(0:npes-1) )
1724 call mpp_get_ug_domain_tile_list(grid%ug_domain, grid%tile(0:npes-1))
1725 grid%ls_me => grid%ls(xmap%me-xmap%root_pe); grid%le_me => grid%le(xmap%me-xmap%root_pe)
1726 grid%gs_me => grid%gs(xmap%me-xmap%root_pe); grid%ge_me => grid%ge(xmap%me-xmap%root_pe)
1727 grid%tile_me => grid%tile(xmap%me-xmap%root_pe)
1728 grid%nxl_me = grid%le_me - grid%ls_me + 1
1729 if (
associated(grid%l_index))
deallocate(grid%l_index)
1730 allocate(grid%l_index(grid%gs_me:grid%ge_me))
1731 allocate(grid_index(grid%ls_me:grid%le_me))
1732 call mpp_get_ug_domain_grid_index(grid%ug_domain, grid_index)
1735 do l = grid%ls_me,grid%le_me
1736 grid%l_index(grid_index(l)) = l
1739 if( grid%on_this_pe )
then
1740 if (
associated(grid%area))
deallocate(grid%area)
1741 if (
associated(grid%area_inv))
deallocate(grid%area_inv)
1742 allocate( grid%area (grid%ls_me:grid%le_me,1) )
1743 allocate( grid%area_inv(grid%ls_me:grid%le_me,1) )
1744 grid%area = 0.0_r8_kind
1748 else if( grid%on_this_pe )
then
1749 if (
associated(grid%area))
deallocate(grid%area)
1750 if (
associated(grid%area_inv))
deallocate(grid%area_inv)
1751 allocate( grid%area (grid%is_me:grid%ie_me, grid%js_me:grid%je_me) )
1752 allocate( grid%area_inv(grid%is_me:grid%ie_me, grid%js_me:grid%je_me) )
1753 grid%area = 0.0_r8_kind
1759 if(.not. grid%is_ug)
then
1760 select case(xmap%version)
1762 if( grid%npes .NE. xmap%npes )
then
1763 call error_mesg(
'xgrid_mod', .NE.
' grid%npes xmap%npes ', fatal)
1767 allocate(pelist(0:xmap%npes-1))
1768 call mpp_get_current_pelist(pelist)
1769 if( grid%on_this_pe )
then
1777 if( g == 1 .AND. grid_ids(1) ==
'ATM' )
then
1779 ie_nest, js_nest, je_nest, is_parent, ie_parent, js_parent, je_parent)
1784 if( use_higher_order .AND. grid%id ==
'ATM')
then
1785 if( nnest > 0 )
call error_mesg(
'xgrid_mod',
'second_order is not supported for nested coupler', fatal)
1786 if( grid%is_latlon )
then
1787 call mpp_modify_domain(grid%domain, grid%domain_with_halo, whalo=1, ehalo=1, shalo=1, nhalo=1)
1788 call mpp_get_data_domain(grid%domain_with_halo, grid%isd_me, grid%ied_me, grid%jsd_me, grid%jed_me, &
1789 xsize=grid%nxd_me, ysize=grid%nyd_me)
1791 if(.NOT.
present(atm_grid))
call error_mesg(
'xgrid_mod', &
1792 'when first grid is "ATM", atm_grid should be present', fatal)
1793 if(grid%is_me-grid%isd_me .NE. 1 .or. grid%ied_me-grid%ie_me .NE. 1 .or. &
1794 grid%js_me-grid%jsd_me .NE. 1 .or. grid%jed_me-grid%je_me .NE. 1 ) &
1795 &
call error_mesg(
'xgrid_mod',
'for non-latlon grid (cubic grid), '//&
1796 &
'the halo size should be 1 in all four direction', fatal)
1797 if(.NOT.(
ASSOCIATED(atm_grid%dx) .AND.
ASSOCIATED(atm_grid%dy) .AND.
ASSOCIATED(atm_grid%edge_w) .AND. &
1798 ASSOCIATED(atm_grid%edge_e) .AND.
ASSOCIATED(atm_grid%edge_s) .AND.
ASSOCIATED(atm_grid%edge_n).AND.&
1799 ASSOCIATED(atm_grid%en1) .AND.
ASSOCIATED(atm_grid%en2) .AND.
ASSOCIATED(atm_grid%vlon) .AND. &
1800 ASSOCIATED(atm_grid%vlat) ) )
call error_mesg(
'xgrid_mod', &
1801 'for non-latlon grid (cubic grid), all the fields in atm_grid data type should be allocated', fatal)
1802 nxc = grid%ie_me - grid%is_me + 1
1803 nyc = grid%je_me - grid%js_me + 1
1804 if(
size(atm_grid%dx,1) .NE. nxc .OR.
size(atm_grid%dx,2) .NE. nyc+1) &
1805 call error_mesg(
'xgrid_mod',
'incorrect dimension size of atm_grid%dx', fatal)
1806 if(
size(atm_grid%dy,1) .NE. nxc+1 .OR.
size(atm_grid%dy,2) .NE. nyc) &
1807 call error_mesg(
'xgrid_mod',
'incorrect dimension sizeof atm_grid%dy', fatal)
1808 if(
size(atm_grid%area,1) .NE. nxc .OR.
size(atm_grid%area,2) .NE. nyc) &
1809 call error_mesg(
'xgrid_mod',
'incorrect dimension size of atm_grid%area', fatal)
1810 if(
size(atm_grid%edge_w(:)) .NE. nyc+1 .OR.
size(atm_grid%edge_e(:)) .NE. nyc+1) &
1811 call error_mesg(
'xgrid_mod',
'incorrect dimension size of atm_grid%edge_w/edge_e', fatal)
1812 if(
size(atm_grid%edge_s(:)) .NE. nxc+1 .OR.
size(atm_grid%edge_n(:)) .NE. nxc+1) &
1813 call error_mesg(
'xgrid_mod',
'incorrect dimension size of atm_grid%edge_s/edge_n', fatal)
1814 if(
size(atm_grid%en1,1) .NE. 3 .OR.
size(atm_grid%en1,2) .NE. nxc .OR.
size(atm_grid%en1,3) .NE. nyc+1) &
1815 call error_mesg(
'xgrid_mod',
'incorrect dimension size of atm_grid%en1', fatal)
1816 if(
size(atm_grid%en2,1) .NE. 3 .OR.
size(atm_grid%en2,2) .NE. nxc+1 .OR.
size(atm_grid%en2,3) .NE. nyc) &
1817 call error_mesg(
'xgrid_mod',
'incorrect dimension size of atm_grid%en2', fatal)
1818 if(
size(atm_grid%vlon,1) .NE. 3 .OR.
size(atm_grid%vlon,2) .NE. nxc .OR.
size(atm_grid%vlon,3) .NE. nyc)&
1819 call error_mesg(
'xgrid_mod',
'incorrect dimension size of atm_grid%vlon', fatal)
1820 if(
size(atm_grid%vlat,1) .NE. 3 .OR.
size(atm_grid%vlat,2) .NE. nxc .OR.
size(atm_grid%vlat,3) .NE. nyc)&
1821 call error_mesg(
'xgrid_mod',
'incorrect dimension size of atm_grid%vlat', fatal)
1822 if (
associated(grid%box%dx))
deallocate(grid%box%dx)
1823 if (
associated(grid%box%dy))
deallocate(grid%box%dy)
1824 if (
associated(grid%box%area))
deallocate(grid%box%area)
1825 if (
associated(grid%box%edge_w))
deallocate(grid%box%edge_w)
1826 if (
associated(grid%box%edge_e))
deallocate(grid%box%edge_e)
1827 if (
associated(grid%box%edge_s))
deallocate(grid%box%edge_s)
1828 if (
associated(grid%box%edge_n))
deallocate(grid%box%edge_n)
1829 if (
associated(grid%box%en1))
deallocate(grid%box%en1)
1830 if (
associated(grid%box%en2))
deallocate(grid%box%en2)
1831 if (
associated(grid%box%vlon))
deallocate(grid%box%vlon)
1832 if (
associated(grid%box%vlat))
deallocate(grid%box%vlat)
1833 allocate(grid%box%dx (grid%is_me:grid%ie_me, grid%js_me:grid%je_me+1 ))
1834 allocate(grid%box%dy (grid%is_me:grid%ie_me+1, grid%js_me:grid%je_me ))
1835 allocate(grid%box%area (grid%is_me:grid%ie_me, grid%js_me:grid%je_me ))
1836 allocate(grid%box%edge_w(grid%js_me:grid%je_me+1))
1837 allocate(grid%box%edge_e(grid%js_me:grid%je_me+1))
1838 allocate(grid%box%edge_s(grid%is_me:grid%ie_me+1))
1839 allocate(grid%box%edge_n(grid%is_me:grid%ie_me+1))
1840 allocate(grid%box%en1 (3, grid%is_me:grid%ie_me, grid%js_me:grid%je_me+1 ))
1841 allocate(grid%box%en2 (3, grid%is_me:grid%ie_me+1, grid%js_me:grid%je_me ))
1842 allocate(grid%box%vlon (3, grid%is_me:grid%ie_me, grid%js_me:grid%je_me ))
1843 allocate(grid%box%vlat (3, grid%is_me:grid%ie_me, grid%js_me:grid%je_me ))
1844 grid%box%dx = atm_grid%dx
1845 grid%box%dy = atm_grid%dy
1846 grid%box%area = atm_grid%area
1847 grid%box%edge_w = atm_grid%edge_w
1848 grid%box%edge_e = atm_grid%edge_e
1849 grid%box%edge_s = atm_grid%edge_s
1850 grid%box%edge_n = atm_grid%edge_n
1851 grid%box%en1 = atm_grid%en1
1852 grid%box%en2 = atm_grid%en2
1853 grid%box%vlon = atm_grid%vlon
1854 grid%box%vlat = atm_grid%vlat
1860 if(grid%on_this_pe)
then
1861 if (
associated(grid%frac_area))
deallocate(grid%frac_area)
1863 allocate( grid%frac_area(grid%ls_me:grid%le_me, 1, grid%km) )
1865 allocate( grid%frac_area(grid%is_me:grid%ie_me, grid%js_me:grid%je_me, grid%km) )
1867 grid%frac_area = 1.0_r8_kind
1871 xmap%grids(g) = grid
1874 select case(xmap%version)
1876 call load_xgrid (xmap, grid, grid_file, grid_ids(1), grid_ids(g), 1, 1, use_higher_order)
1878 select case(grid_ids(1))
1886 call error_mesg(
'xgrid_mod',
'grid_ids(1) should be ATM, LND or WAV', fatal)
1888 select case(grid_ids(g))
1890 xgrid_dimname =
'nfile_'//trim(xgrid_name)//
'Xl'
1891 xgrid_name = trim(xgrid_name)//
'Xl_file'
1893 xgrid_dimname =
'nfile_'//trim(xgrid_name)//
'Xo'
1894 xgrid_name = trim(xgrid_name)//
'Xo_file'
1896 xgrid_dimname =
'nfile_'//trim(xgrid_name)//
'Xw'
1897 xgrid_name = trim(xgrid_name)//
'Xw_file'
1899 call error_mesg(
'xgrid_mod',
'grid_ids(g) should be LND, OCN or WAV', fatal)
1903 call read_data(gridfileobj, lowercase(grid_ids(1))//
'_mosaic_file', mosaic1)
1904 call read_data(gridfileobj, lowercase(grid_ids(g))//
'_mosaic_file', mosaic2)
1906 mosaic1 =
'INPUT/'//trim(mosaic1)
1907 mosaic2 =
'INPUT/'//trim(mosaic2)
1909 allocate(tile1_list(grid1%ntile), tile2_list(grid%ntile) )
1910 if(.not.
open_file(fileobj,mosaic1,
"read"))
then
1911 call error_mesg(
'xgrid_mod(setup_xmap)',
'Error when opening mosaic1 file '//trim(mosaic1), fatal)
1913 call read_data(fileobj,
'gridtiles', tile1_list)
1916 if(.not.
open_file(fileobj,mosaic2,
"read"))
then
1917 call error_mesg(
'xgrid_mod(setup_xmap)',
'Error when opening mosaic2 file '//trim(mosaic2), fatal)
1919 call read_data(fileobj,
'gridtiles', tile2_list)
1922 if(variable_exists(gridfileobj, xgrid_name))
then
1923 call get_dimension_size(gridfileobj, xgrid_dimname, nxgrid_file)
1924 if(nxgrid_file>0)
then
1925 allocate(xgrid_filelist(nxgrid_file))
1926 call read_data(gridfileobj, xgrid_name, xgrid_filelist)
1929 do i = 1, nxgrid_file
1930 xgrid_file =
'INPUT/'//trim(xgrid_filelist(i))
1931 if(.not.
open_file(fileobj,xgrid_file,
"read"))
then
1932 call error_mesg(
'xgrid_mod(setup_xmap)',
'Error when opening xgrid file '// &
1933 & trim(xgrid_file), fatal)
1937 call read_data(fileobj,
"contact", contact)
1938 i1 = index(contact,
":")
1939 i2 = index(contact,
"::")
1940 i3 = index(contact,
":", back=.true. )
1941 if(i1 == 0 .OR. i2 == 0)
call error_mesg(
'xgrid_mod', &
1942 'field contact in file '//trim(xgrid_file)//
' should contains ":" and "::" ', fatal)
1944 'field contact in file '//trim(xgrid_file)//
' should contains two ":"', fatal)
1945 tile1_name = contact(i1+1:i2-1)
1946 tile2_name = contact(i3+1:len_trim(contact))
1947 tile1 = 0; tile2 = 0
1948 do j = 1, grid1%ntile
1949 if( trim(tile1_name) == trim(tile1_list(j)) )
then
1954 do j = 1, grid%ntile
1955 if( tile2_name == tile2_list(j) )
then
1961 if(tile1 == 0)
call error_mesg(
'xgrid_mod', &
1962 trim(tile1_name)//
' is not a tile of mosaic '//trim(mosaic1), fatal)
1963 if(tile2 == 0)
call error_mesg(
'xgrid_mod', &
1964 trim(tile2_name)//
' is not a tile of mosaic '//trim(mosaic2), fatal)
1966 call load_xgrid (xmap, grid, xgrid_file, grid_ids(1), grid_ids(g), tile1, tile2, &
1969 deallocate(xgrid_filelist)
1971 deallocate(tile1_list, tile2_list)
1973 if(grid%on_this_pe)
then
1974 grid%area_inv = 0.0_r8_kind;
1975 where (grid%area>0.0_r8_kind) grid%area_inv = 1.0_r8_kind/grid%area
1980 xmap%grids(g) = grid
1985 call mpp_clock_end(id_load_xgrid)
1987 grid1%area_inv = 0.0_r8_kind;
1988 where (grid1%area>0.0_r8_kind)
1989 grid1%area_inv = 1.0_r8_kind/grid1%area
1992 xmap%your1my2(xmap%me-xmap%root_pe) = .false.
1993 xmap%your2my1(xmap%me-xmap%root_pe) = .false.
1996 if (
associated(xmap%send_count_repro))
deallocate(xmap%send_count_repro)
1997 if (
associated(xmap%recv_count_repro))
deallocate(xmap%recv_count_repro)
1998 allocate( xmap%send_count_repro(0:xmap%npes-1) )
1999 allocate( xmap%recv_count_repro(0:xmap%npes-1) )
2000 xmap%send_count_repro = 0
2001 xmap%recv_count_repro = 0
2002 do g=2,
size(xmap%grids(:))
2004 if(xmap%grids(g)%size >0) &
2005 xmap%send_count_repro(p) = xmap%send_count_repro(p) &
2006 +count(xmap%grids(g)%x (:)%pe==p+xmap%root_pe)
2007 if(xmap%grids(g)%size_repro >0) &
2008 xmap%recv_count_repro(p) = xmap%recv_count_repro(p) &
2009 +count(xmap%grids(g)%x_repro(:)%pe==p+xmap%root_pe)
2012 xmap%send_count_repro_tot = sum(xmap%send_count_repro)
2013 xmap%recv_count_repro_tot = sum(xmap%recv_count_repro)
2015 xmap%send_count_repro_tot = 0
2016 xmap%recv_count_repro_tot = 0
2019 if (
associated(xmap%x1))
deallocate(xmap%x1)
2020 if (
associated(xmap%x2))
deallocate(xmap%x2)
2021 if (
associated(xmap%x1_put))
deallocate(xmap%x1_put)
2022 if (
associated(xmap%x2_get))
deallocate(xmap%x2_get)
2023 allocate( xmap%x1(1:sum(xmap%grids(2:
size(xmap%grids(:)))%size)) )
2024 allocate( xmap%x2(1:sum(xmap%grids(2:
size(xmap%grids(:)))%size)) )
2025 allocate( xmap%x1_put(1:sum(xmap%grids(2:
size(xmap%grids(:)))%size)) )
2026 allocate( xmap%x2_get(1:sum(xmap%grids(2:
size(xmap%grids(:)))%size)) )
2029 if (
associated(xmap%get1))
deallocate(xmap%get1)
2030 if (
associated(xmap%put1))
deallocate(xmap%put1)
2031 allocate(xmap%get1, xmap%put1)
2032 call mpp_clock_begin(id_set_comm)
2034 call set_comm_get1(xmap)
2036 call set_comm_put1(xmap)
2039 if (
associated(xmap%get1_repro))
deallocate(xmap%get1_repro)
2040 allocate(xmap%get1_repro)
2041 call set_comm_get1_repro(xmap)
2044 call mpp_clock_end(id_set_comm)
2046 call mpp_clock_begin(id_regen)
2048 call mpp_clock_end(id_regen)
2050 call mpp_clock_begin(id_conservation_check)
2052 if(lnd_ug_id ==0)
then
2055 allocate(tmp_2d(grid1%is_me:grid1%ie_me, grid1%js_me:grid1%je_me))
2056 tmp_2d = 1.0_r8_kind
2060 write(out_unit,* )
"Checked data is array of constant 1"
2061 write(out_unit,* )grid1%id,
'(',xmap%grids(:)%id,
')=', xxx
2063 if(lnd_ug_id == 0)
then
2064 do g=2,
size(xmap%grids(:))
2065 xxx =
conservation_check(xmap%grids(g)%frac_area*0.0_r8_kind+1.0_r8_kind, xmap%grids(g)%id, xmap )
2066 write( out_unit,* )xmap%grids(g)%id,
'(',xmap%grids(:)%id,
')=', xxx
2069 do g=2,
size(xmap%grids(:))
2070 grid => xmap%grids(g)
2071 allocate(tmp_3d(grid%is_me:grid%ie_me, grid%js_me:grid%je_me,grid%km))
2072 tmp_3d = 1.0_r8_kind
2074 write( out_unit,* )xmap%grids(g)%id,
'(',xmap%grids(:)%id,
')=', xxx
2079 if(grid1%id ==
"ATM")
then
2080 allocate(check_data(
size(grid1%area,1),
size(grid1%area,2)))
2081 call random_number(check_data)
2084 if(lnd_ug_id ==0)
then
2089 write( out_unit,* ) &
2090 "Checked data is array of random number between 0 and 1 using "//trim(
interp_method)
2091 write( out_unit,* )grid1%id,
'(',xmap%grids(:)%id,
')=', xxx
2093 deallocate(check_data)
2094 do g=2,
size(xmap%grids(:))
2095 allocate(check_data_3d(xmap%grids(g)%is_me:xmap%grids(g)%ie_me, &
2096 xmap%grids(g)%js_me:xmap%grids(g)%je_me, grid1%km))
2097 call random_number(check_data_3d)
2098 if(lnd_ug_id ==0)
then
2103 write( out_unit,* )xmap%grids(g)%id,
'(',xmap%grids(:)%id,
')=', xxx
2104 deallocate( check_data_3d)
2107 call mpp_clock_end(id_conservation_check)
2109 call mpp_clock_end(id_setup_xmap)
2118 ie_nest_out, js_nest_out, je_nest_out, is_parent_out, &
2119 ie_parent_out, js_parent_out, je_parent_out) &
2121 type(fmsnetcdffile_t),
intent(in) :: fileobj
2122 integer,
intent(out) :: tile_nest_out, tile_parent_out
2123 integer,
intent(out) :: is_nest_out, ie_nest_out
2124 integer,
intent(out) :: js_nest_out, je_nest_out
2125 integer,
intent(out) :: is_parent_out, ie_parent_out
2126 integer,
intent(out) :: js_parent_out, je_parent_out
2129 integer :: ntiles, ncontacts, n, t1, t2
2130 integer :: nx1_contact, ny1_contact
2131 integer :: nx2_contact, ny2_contact
2132 integer,
allocatable,
dimension(:) :: nx, ny
2133 integer,
allocatable,
dimension(:) :: tile1, tile2
2134 integer,
allocatable,
dimension(:) :: istart1, iend1, jstart1, jend1
2135 integer,
allocatable,
dimension(:) :: istart2, iend2, jstart2, jend2
2137 tile_nest_out = 0; tile_parent_out = 0
2138 is_nest_out = 0; ie_nest_out = 0
2139 js_nest_out = 0; je_nest_out = 0
2140 is_parent_out = 0; ie_parent_out = 0
2141 js_parent_out = 0; je_parent_out = 0
2145 ntiles = get_mosaic_ntiles(fileobj)
2146 if( ntiles == 1 )
return
2148 allocate(nx(ntiles), ny(ntiles))
2149 call get_mosaic_grid_sizes(fileobj, nx, ny)
2151 ncontacts = get_mosaic_ncontacts(fileobj)
2152 if(ncontacts == 0)
return
2153 allocate(tile1(ncontacts), tile2(ncontacts))
2154 allocate(istart1(ncontacts), iend1(ncontacts))
2155 allocate(jstart1(ncontacts), jend1(ncontacts))
2156 allocate(istart2(ncontacts), iend2(ncontacts))
2157 allocate(jstart2(ncontacts), jend2(ncontacts))
2159 call get_mosaic_contact( fileobj, tile1, tile2, istart1, iend1, jstart1, jend1, &
2160 istart2, iend2, jstart2, jend2)
2163 if( tile1(n) == tile2(n) ) cycle
2165 nx1_contact = iend1(n)-istart1(n)+1
2166 ny1_contact = jend1(n)-jstart1(n)+1
2167 nx2_contact = iend2(n)-istart2(n)+1
2168 ny2_contact = jend2(n)-jstart2(n)+1
2172 if( (nx(t1) .NE. nx1_contact .OR. ny(t1) .NE. ny1_contact ) .AND. &
2173 (nx(t2) .NE. nx2_contact .OR. ny(t2) .NE. ny2_contact ) ) cycle
2174 if(nx1_contact == nx2_contact .AND. ny1_contact == ny2_contact)
then
2175 call error_mesg(
'xgrid_mod',
'There is no refinement for the overlapping region', fatal)
2180 call error_mesg(
'xgrid_mod',
'only support one nest region, contact developer' ,fatal)
2182 if(nx2_contact*ny2_contact > nx1_contact*ny1_contact)
then
2183 is_nest_out = istart2(n);
2184 ie_nest_out = iend2(n);
2185 js_nest_out = jstart2(n);
2186 je_nest_out = jend2(n);
2187 tile_nest_out = tile2(n);
2188 is_parent_out = istart1(n);
2189 ie_parent_out = iend1(n);
2190 js_parent_out = jstart1(n);
2191 je_parent_out = jend1(n);
2192 tile_parent_out = tile1(n);
2194 is_nest_out = istart1(n);
2195 ie_nest_out = iend1(n);
2196 js_nest_out = jstart1(n);
2197 je_nest_out = jend1(n);
2198 tile_nest_out = tile1(n);
2199 is_parent_out = istart2(n);
2200 ie_parent_out = iend2(n);
2201 js_parent_out = jstart2(n);
2202 je_parent_out = jend2(n);
2203 tile_parent_out = tile2(n);
2207 deallocate(nx, ny, tile1, tile2)
2208 deallocate(istart1, iend1, jstart1, jend1)
2209 deallocate(istart2, iend2, jstart2, jend2)
2217 subroutine set_comm_get1_repro(xmap)
2219 integer,
dimension(xmap%npes) :: pe_ind, cnt
2220 integer,
dimension(0:xmap%npes-1) :: send_ind, pl
2221 integer :: npes, nsend, nrecv, mypos
2222 integer :: m, p, pos, n, g, l, im, i, j
2223 type(
comm_type),
pointer,
save :: comm => null()
2225 comm => xmap%get1_repro
2229 mypos =
mpp_pe() - mpp_root_pe()
2231 p = mod(mypos+npes-m, npes)
2232 if( xmap%recv_count_repro(p) > 0 )
then
2239 if( nrecv > 0 )
then
2240 if (
associated(comm%recv))
deallocate(comm%recv)
2241 allocate(comm%recv(nrecv))
2245 comm%recv(n)%count = xmap%recv_count_repro(p)
2246 comm%recv(n)%pe = p + xmap%root_pe
2247 comm%recv(n)%buffer_pos = pos
2248 pos = pos + comm%recv(n)%count
2255 mypos =
mpp_pe() - mpp_root_pe()
2257 p = mod(mypos+m, npes)
2258 if( xmap%send_count_repro(p) > 0 )
then
2266 if( nsend > 0 )
then
2267 if (
associated(comm%send))
deallocate(comm%send)
2268 allocate(comm%send(nsend))
2273 comm%send(n)%count = xmap%send_count_repro(p)
2274 comm%send(n)%pe = p + xmap%root_pe
2275 comm%send(n)%buffer_pos = pos
2276 pos = pos + comm%send(n)%count
2277 allocate(comm%send(n)%i(comm%send(n)%count))
2278 allocate(comm%send(n)%j(comm%send(n)%count))
2279 allocate(comm%send(n)%g(comm%send(n)%count))
2280 allocate(comm%send(n)%xLoc(comm%send(n)%count))
2283 do g=2,
size(xmap%grids(:))
2284 im = xmap%grids(g)%im
2285 do l=1,xmap%grids(g)%size
2286 p = xmap%grids(g)%x(l)%pe-xmap%root_pe
2290 i = xmap%grids(g)%x(l)%i2
2291 j = xmap%grids(g)%x(l)%j2
2292 if(xmap%grids(g)%is_ug)
then
2293 comm%send(n)%i(pos) = xmap%grids(g)%l_index((j-1)*im+i)
2294 comm%send(n)%j(pos) = 1
2296 comm%send(n)%i(pos) = xmap%grids(g)%x(l)%i2
2297 comm%send(n)%j(pos) = xmap%grids(g)%x(l)%j2
2299 comm%send(n)%g(pos) = g
2304 if( comm%send(n)%count .NE. cnt(n) )
call error_mesg(
'xgrid_mod', &
2305 .NE.
'comm%send(n)%count cnt(n)', fatal)
2311 do g=2,
size(xmap%grids(:))
2312 do l=1,xmap%grids(g)%size_repro
2313 p = xmap%grids(g)%x_repro(l)%pe-xmap%root_pe
2314 xmap%grids(g)%x_repro(l)%recv_pos = pl(p)
2321 end subroutine set_comm_get1_repro
2324 subroutine set_comm_get1(xmap)
2326 type (
grid_type),
pointer,
save :: grid1 =>null()
2327 integer,
allocatable :: send_size(:)
2328 integer,
allocatable :: recv_size(:)
2329 integer :: max_size, g, npes, l, ll, nset, m
2330 integer :: i1, j1, tile1, p, n, pos, buffer_pos, mypos
2331 integer :: nsend, nrecv, rbuf_size, sbuf_size, msgsize
2333 real(r8_kind),
allocatable :: recv_buf(:), send_buf(:)
2334 real(r8_kind),
allocatable :: diarray(:), djarray(:)
2335 integer,
allocatable :: iarray(:), jarray(:), tarray(:)
2336 integer,
allocatable :: pos_x(:), pelist(:), size_pe(:), pe_side1(:)
2337 integer :: recv_buffer_pos(0:xmap%npes)
2338 integer :: send_buffer_pos(0:xmap%npes)
2339 type(
comm_type),
pointer,
save :: comm => null()
2343 do g=2,
size(xmap%grids(:))
2344 max_size = max_size + xmap%grids(g)%size
2347 grid1 => xmap%grids(1)
2352 allocate(pelist(0:npes-1))
2353 call mpp_get_current_pelist(pelist)
2354 allocate(send_size(0:npes-1))
2355 allocate(recv_size(0:npes-1))
2356 allocate(size_pe(0:npes-1))
2357 allocate(pos_x(0:npes-1))
2362 if(max_size > 0)
then
2363 allocate(pe_side1(max_size))
2364 if (
associated(xmap%ind_get1))
deallocate(xmap%ind_get1)
2365 allocate(xmap%ind_get1(max_size))
2369 do g=2,
size(xmap%grids(:))
2370 do l=1,xmap%grids(g)%size
2371 i1 = xmap%grids(g)%x(l)%i1
2372 j1 = xmap%grids(g)%x(l)%j1
2373 tile1 = xmap%grids(g)%x(l)%tile
2375 if(grid1%tile(p) == tile1)
then
2377 size_pe(p) = size_pe(p) + 1
2382 if( p == npes )
then
2383 call error_mesg(
'xgrid_mod',
'tile is not in grid1%tile(:)', fatal)
2392 pos_x(p) = pos_x(p-1) + size_pe(p-1)
2396 allocate(iarray(max_size))
2397 allocate(jarray(max_size))
2398 allocate(tarray(max_size))
2399 if(monotonic_exchange)
then
2400 allocate(diarray(max_size))
2401 allocate(djarray(max_size))
2406 do g=2,
size(xmap%grids(:))
2407 do l=1,xmap%grids(g)%size
2408 i1 = xmap%grids(g)%x(l)%i1
2409 j1 = xmap%grids(g)%x(l)%j1
2410 tile1 = xmap%grids(g)%x(l)%tile
2415 if(send_size(p) > 0)
then
2416 if( i1 == iarray(pos_x(p)+send_size(p)) .AND. j1 == jarray(pos_x(p)+send_size(p)) &
2417 .AND. tile1 == tarray(pos_x(p)+send_size(p)))
then
2422 do n = 1, send_size(p)
2423 if(i1 == iarray(pos_x(p)+n) .AND. j1 == jarray(pos_x(p)+n) .AND. tile1 == tarray(pos_x(p)+n))
then
2430 if( (.NOT. found) .OR. monotonic_exchange )
then
2431 send_size(p) = send_size(p)+1
2432 pos = pos_x(p)+send_size(p)
2436 if(monotonic_exchange)
then
2437 diarray(pos) = xmap%grids(g)%x(l)%di
2438 djarray(pos) = xmap%grids(g)%x(l)%dj
2442 xmap%ind_get1(ll) = n
2448 pos_x(p) = pos_x(p-1) + send_size(p-1)
2452 do g=2,
size(xmap%grids(:))
2453 do l=1,xmap%grids(g)%size
2456 xmap%ind_get1(ll) = pos_x(p) + xmap%ind_get1(ll)
2461 mypos =
mpp_pe()-mpp_root_pe()
2464 recv_size(:) = xmap%your2my1_size(:)
2465 nsend = count( send_size> 0)
2468 if (
associated(comm%send))
deallocate(comm%send)
2469 allocate(comm%send(nsend))
2470 comm%send(:)%count = 0
2475 send_buffer_pos(p) = pos
2476 pos = pos + send_size(p)
2482 p = mod(mypos+n, npes)
2483 if(send_size(p)>0)
then
2485 allocate(comm%send(pos)%i(send_size(p)))
2486 comm%send(pos)%buffer_pos = send_buffer_pos(p)
2487 comm%send(pos)%count = send_size(p)
2488 comm%send(pos)%pe = pelist(p)
2489 comm%sendsize = comm%sendsize + send_size(p)
2494 if(monotonic_exchange) nset = 5
2495 rbuf_size = sum(recv_size)*nset
2496 sbuf_size = sum(send_size)*nset
2497 if(rbuf_size>0)
allocate(recv_buf(rbuf_size))
2498 if(sbuf_size>0)
allocate(send_buf(sbuf_size))
2502 p = mod(mypos+npes-n, npes)
2503 if(recv_size(p) ==0) cycle
2504 msgsize = recv_size(p)*nset
2505 call mpp_recv(recv_buf(pos+1), glen=msgsize, from_pe=pelist(p), block=.false., tag=comm_tag_4)
2511 pos_x(p) = pos_x(p-1) + size_pe(p-1)
2516 p = mod(mypos+n, npes)
2517 do l = 1, send_size(p)
2518 send_buf(pos+1) = real(iarray(pos_x(p)+l), r8_kind)
2519 send_buf(pos+2) = real(jarray(pos_x(p)+l), r8_kind)
2520 send_buf(pos+3) = real(tarray(pos_x(p)+l), r8_kind)
2521 if(monotonic_exchange)
then
2522 send_buf(pos+4) = diarray(pos_x(p)+l)
2523 send_buf(pos+5) = djarray(pos_x(p)+l)
2531 p = mod(mypos+n, npes)
2532 if(send_size(p) ==0) cycle
2533 msgsize = send_size(p)*nset
2534 call mpp_send(send_buf(pos+1), plen=msgsize, to_pe=pelist(p), tag=comm_tag_4 )
2539 nrecv = count(recv_size>0)
2544 if (
associated(comm%recv))
deallocate(comm%recv)
2545 allocate(comm%recv(nrecv))
2546 comm%recv(:)%count = 0
2550 recv_buffer_pos(p) = buffer_pos
2551 buffer_pos = buffer_pos + recv_size(p)
2556 p = mod(mypos+npes-m, npes)
2557 if(recv_size(p)>0)
then
2559 allocate(comm%recv(pos)%i(recv_size(p)))
2560 allocate(comm%recv(pos)%j(recv_size(p)))
2561 allocate(comm%recv(pos)%tile(recv_size(p)))
2562 comm%recv(pos)%buffer_pos = recv_buffer_pos(p)
2563 comm%recv(pos)%pe = pelist(p)
2564 comm%recv(pos)%count = recv_size(p)
2565 comm%recvsize = comm%recvsize + recv_size(p)
2566 if(monotonic_exchange)
then
2567 allocate(comm%recv(pos)%di(recv_size(p)))
2568 allocate(comm%recv(pos)%dj(recv_size(p)))
2570 if(grid1%is_ug)
then
2571 do n = 1, recv_size(p)
2572 i = int(recv_buf(buffer_pos+1))
2573 j = int(recv_buf(buffer_pos+2))
2574 comm%recv(pos)%i(n) = grid1%l_index((j-1)*grid1%im+i)
2575 comm%recv(pos)%j(n) = 1
2576 comm%recv(pos)%tile(n) = int(recv_buf(buffer_pos+3))
2577 if(monotonic_exchange)
then
2578 comm%recv(pos)%di(n) = recv_buf(buffer_pos+4)
2579 comm%recv(pos)%dj(n) = recv_buf(buffer_pos+5)
2581 buffer_pos = buffer_pos + nset
2584 do n = 1, recv_size(p)
2585 comm%recv(pos)%i(n) = int(recv_buf(buffer_pos+1) )- grid1%is_me + 1
2586 comm%recv(pos)%j(n) = int(recv_buf(buffer_pos+2) )- grid1%js_me + 1
2587 comm%recv(pos)%tile(n) = int(recv_buf(buffer_pos+3))
2588 if(monotonic_exchange)
then
2589 comm%recv(pos)%di(n) = recv_buf(buffer_pos+4)
2590 comm%recv(pos)%dj(n) = recv_buf(buffer_pos+5)
2592 buffer_pos = buffer_pos + nset
2597 if (
associated(comm%unpack_ind))
deallocate(comm%unpack_ind)
2598 allocate(comm%unpack_ind(nrecv))
2601 if(recv_size(p)>0)
then
2604 if(comm%recv(m)%pe == pelist(p))
then
2605 comm%unpack_ind(pos) = m
2614 if(
allocated(send_buf) )
deallocate(send_buf)
2615 if(
allocated(recv_buf) )
deallocate(recv_buf)
2616 if(
allocated(pelist) )
deallocate(pelist)
2617 if(
allocated(pos_x) )
deallocate(pos_x)
2618 if(
allocated(pelist) )
deallocate(pelist)
2619 if(
allocated(iarray) )
deallocate(iarray)
2620 if(
allocated(jarray) )
deallocate(jarray)
2621 if(
allocated(tarray) )
deallocate(tarray)
2622 if(
allocated(size_pe) )
deallocate(size_pe)
2624 end subroutine set_comm_get1
2627 subroutine set_comm_put1(xmap)
2629 type (
grid_type),
pointer,
save :: grid1 =>null()
2630 integer,
allocatable :: send_size(:)
2631 integer,
allocatable :: recv_size(:)
2632 integer :: max_size, g, npes, l, ll, m, mypos
2633 integer :: i1, j1, tile1, p, n, pos, buffer_pos
2634 integer :: nsend, nrecv, msgsize, nset, rbuf_size, sbuf_size
2636 real(r8_kind),
allocatable :: recv_buf(:), send_buf(:)
2637 real(r8_kind),
allocatable :: diarray(:), djarray(:)
2638 integer,
allocatable :: iarray(:), jarray(:), tarray(:)
2639 integer,
allocatable :: pos_x(:), pelist(:), size_pe(:), pe_put1(:)
2640 integer :: recv_buffer_pos(0:xmap%npes)
2641 type(
comm_type),
pointer,
save :: comm => null()
2645 if(nnest == 0 .OR. xmap%grids(1)%id .NE.
'ATM' )
then
2646 comm%nsend = xmap%get1%nrecv
2647 comm%nrecv = xmap%get1%nsend
2648 comm%sendsize = xmap%get1%recvsize
2649 comm%recvsize = xmap%get1%sendsize
2650 comm%send => xmap%get1%recv
2651 comm%recv => xmap%get1%send
2652 xmap%ind_put1 => xmap%ind_get1
2657 do g=2,
size(xmap%grids(:))
2658 max_size = max_size + xmap%grids(g)%size
2660 grid1 => xmap%grids(1)
2664 allocate(pelist(0:npes-1))
2665 call mpp_get_current_pelist(pelist)
2666 allocate(send_size(0:npes-1))
2667 allocate(recv_size(0:npes-1))
2668 allocate(size_pe(0:npes-1))
2669 allocate(pos_x(0:npes-1))
2674 if(max_size > 0)
then
2675 allocate(pe_put1(max_size))
2676 if (
associated(xmap%ind_put1))
deallocate(xmap%ind_put1)
2677 allocate(xmap%ind_put1(max_size))
2681 do g=2,
size(xmap%grids(:))
2682 do l=1,xmap%grids(g)%size
2683 i1 = xmap%grids(g)%x(l)%i1
2684 j1 = xmap%grids(g)%x(l)%j1
2685 tile1 = xmap%grids(g)%x(l)%tile
2687 if(grid1%tile(p) == tile1)
then
2688 if(
in_box(i1, j1, grid1%is(p), grid1%ie(p), grid1%js(p), grid1%je(p)))
then
2689 size_pe(p) = size_pe(p) + 1
2701 pos_x(p) = pos_x(p-1) + size_pe(p-1)
2705 allocate(iarray(max_size))
2706 allocate(jarray(max_size))
2707 allocate(tarray(max_size))
2708 if(monotonic_exchange)
then
2709 allocate(diarray(max_size))
2710 allocate(djarray(max_size))
2715 do g=2,
size(xmap%grids(:))
2716 do l=1,xmap%grids(g)%size
2717 i1 = xmap%grids(g)%x(l)%i1
2718 j1 = xmap%grids(g)%x(l)%j1
2719 tile1 = xmap%grids(g)%x(l)%tile
2724 if(send_size(p) > 0)
then
2725 if( i1 == iarray(pos_x(p)+send_size(p)) .AND. j1 == jarray(pos_x(p)+send_size(p)) &
2726 .AND. tile1 == tarray(pos_x(p)+send_size(p)))
then
2731 do n = 1, send_size(p)
2732 if(i1 == iarray(pos_x(p)+n) .AND. j1 == jarray(pos_x(p)+n) .AND. tile1 == tarray(pos_x(p)+n))
then
2739 if( (.NOT. found) .OR. monotonic_exchange )
then
2740 send_size(p) = send_size(p)+1
2741 pos = pos_x(p)+send_size(p)
2745 if(monotonic_exchange)
then
2746 diarray(pos) = xmap%grids(g)%x(l)%di
2747 djarray(pos) = xmap%grids(g)%x(l)%dj
2751 xmap%ind_put1(ll) = n
2757 pos_x(p) = pos_x(p-1) + send_size(p-1)
2761 do g=2,
size(xmap%grids(:))
2762 do l=1,xmap%grids(g)%size
2763 i1 = xmap%grids(g)%x(l)%i1
2764 j1 = xmap%grids(g)%x(l)%j1
2765 tile1 = xmap%grids(g)%x(l)%tile
2768 xmap%ind_put1(ll) = pos_x(p) + xmap%ind_put1(ll)
2774 mypos =
mpp_pe()-mpp_root_pe()
2776 if (do_alltoall)
then
2777 call mpp_alltoall(send_size, 1, recv_size, 1)
2780 p = mod(mypos+npes-n, npes)
2781 call mpp_recv(recv_size(p), glen=1, from_pe=pelist(p), block=.false., tag=comm_tag_5)
2786 p = mod(mypos+n, npes)
2787 call mpp_send(send_size(p), plen=1, to_pe=pelist(p), tag=comm_tag_5)
2794 nrecv = count( send_size> 0)
2797 if (
associated(comm%recv))
deallocate(comm%recv)
2798 allocate(comm%recv(nrecv))
2799 comm%recv(:)%count = 0
2804 recv_buffer_pos(p) = pos
2805 pos = pos + send_size(p)
2810 p = mod(mypos+npes-n, npes)
2811 if(send_size(p)>0)
then
2813 allocate(comm%recv(pos)%i(send_size(p)))
2814 comm%recv(pos)%buffer_pos = recv_buffer_pos(p)
2815 comm%recv(pos)%count = send_size(p)
2816 comm%recv(pos)%pe = pelist(p)
2817 comm%recvsize = comm%recvsize + send_size(p)
2822 if(monotonic_exchange) nset = 5
2823 rbuf_size = sum(recv_size)*nset
2824 sbuf_size = sum(send_size)*nset
2825 if(rbuf_size>0)
allocate(recv_buf(rbuf_size))
2826 if(sbuf_size>0)
allocate(send_buf(sbuf_size))
2830 p = mod(mypos+npes-n, npes)
2831 if(recv_size(p) ==0) cycle
2832 msgsize = recv_size(p)*nset
2833 call mpp_recv(recv_buf(pos+1), glen=msgsize, from_pe=pelist(p), block=.false., tag=comm_tag_6)
2839 pos_x(p) = pos_x(p-1) + size_pe(p-1)
2844 p = mod(mypos+n, npes)
2845 do l = 1, send_size(p)
2846 send_buf(pos+1) = real(iarray(pos_x(p)+l), r8_kind)
2847 send_buf(pos+2) = real(jarray(pos_x(p)+l), r8_kind)
2848 send_buf(pos+3) = real(tarray(pos_x(p)+l), r8_kind)
2849 if(monotonic_exchange)
then
2850 send_buf(pos+4) = diarray(pos_x(p)+l)
2851 send_buf(pos+5) = djarray(pos_x(p)+l)
2859 p = mod(mypos+n, npes)
2860 if(send_size(p) ==0) cycle
2861 msgsize = send_size(p)*nset
2862 call mpp_send(send_buf(pos+1), plen=msgsize, to_pe=pelist(p), tag=comm_tag_6 )
2867 nsend = count(recv_size>0)
2872 if (
associated(comm%send))
deallocate(comm%send)
2873 allocate(comm%send(nsend))
2874 comm%send(:)%count = 0
2878 p = mod(mypos+npes-m, npes)
2879 if(recv_size(p)>0)
then
2881 allocate(comm%send(pos)%i(recv_size(p)))
2882 allocate(comm%send(pos)%j(recv_size(p)))
2883 allocate(comm%send(pos)%tile(recv_size(p)))
2884 comm%send(pos)%pe = pelist(p)
2885 comm%send(pos)%count = recv_size(p)
2886 comm%sendsize = comm%sendsize + recv_size(p)
2887 if(monotonic_exchange)
then
2888 allocate(comm%send(pos)%di(recv_size(p)))
2889 allocate(comm%send(pos)%dj(recv_size(p)))
2891 do n = 1, recv_size(p)
2892 comm%send(pos)%i(n) = int(recv_buf(buffer_pos+1) )- grid1%is_me + 1
2893 comm%send(pos)%j(n) = int(recv_buf(buffer_pos+2) )- grid1%js_me + 1
2894 comm%send(pos)%tile(n) = int(recv_buf(buffer_pos+3))
2895 if(monotonic_exchange)
then
2896 comm%send(pos)%di(n) = recv_buf(buffer_pos+4)
2897 comm%send(pos)%dj(n) = recv_buf(buffer_pos+5)
2899 buffer_pos = buffer_pos + nset
2906 if(
allocated(send_buf) )
deallocate(send_buf)
2907 if(
allocated(recv_buf) )
deallocate(recv_buf)
2908 if(
allocated(pelist) )
deallocate(pelist)
2909 if(
allocated(pos_x) )
deallocate(pos_x)
2910 if(
allocated(pelist) )
deallocate(pelist)
2911 if(
allocated(iarray) )
deallocate(iarray)
2912 if(
allocated(jarray) )
deallocate(jarray)
2913 if(
allocated(tarray) )
deallocate(tarray)
2914 if(
allocated(size_pe) )
deallocate(size_pe)
2916 end subroutine set_comm_put1
2940 type (xmap_type),
intent(inout) :: xmap
2942 integer :: g, l, k, max_size
2943 integer :: i1, j1, i2, j2, p
2946 logical :: overlap_with_nest
2947 integer :: cnt(xmap%get1%nsend)
2948 integer :: i,j,n,xloc,pos,nsend,m,npes, mypos
2949 integer :: send_ind(0:xmap%npes-1)
2953 do g=2,
size(xmap%grids(:))
2954 max_size = max_size + xmap%grids(g)%size * xmap%grids(g)%km
2957 if (max_size>
size(xmap%x1(:)))
then
2958 if (
associated(xmap%x1))
deallocate(xmap%x1)
2959 if (
associated(xmap%x2))
deallocate(xmap%x2)
2960 allocate( xmap%x1(1:max_size) )
2961 allocate( xmap%x2(1:max_size) )
2965 do g=2,
size(xmap%grids(:))
2966 xmap%grids(g)%first = 1
2967 xmap%grids(g)%last = 0
2972 do g=2,
size(xmap%grids(:))
2973 xmap%grids(g)%first = xmap%size + 1;
2975 do l=1,xmap%grids(g)%size
2976 i1 = xmap%grids(g)%x(l)%i1
2977 j1 = xmap%grids(g)%x(l)%j1
2978 i2 = xmap%grids(g)%x(l)%i2
2979 j2 = xmap%grids(g)%x(l)%j2
2980 tile1 = xmap%grids(g)%x(l)%tile
2982 if(xmap%grids(g)%is_ug)
then
2983 do k=1,xmap%grids(g)%km
2984 lll = xmap%grids(g)%l_index((j2-1)*xmap%grids(g)%im+i2)
2985 if (xmap%grids(g)%frac_area(lll,1,k)/=0.0_r8_kind)
then
2986 xmap%size = xmap%size+1
2987 xmap%x1(xmap%size)%pos = xmap%ind_get1(ll)
2988 xmap%x1(xmap%size)%i = xmap%grids(g)%x(l)%i1
2989 xmap%x1(xmap%size)%j = xmap%grids(g)%x(l)%j1
2990 xmap%x1(xmap%size)%tile = xmap%grids(g)%x(l)%tile
2991 xmap%x1(xmap%size)%area = xmap%grids(g)%x(l)%area &
2992 *xmap%grids(g)%frac_area(lll,1,k)
2993 xmap%x1(xmap%size)%di = xmap%grids(g)%x(l)%di
2994 xmap%x1(xmap%size)%dj = xmap%grids(g)%x(l)%dj
2995 xmap%x2(xmap%size)%i = xmap%grids(g)%x(l)%i2
2996 xmap%x2(xmap%size)%j = xmap%grids(g)%x(l)%j2
2997 xmap%x2(xmap%size)%l = lll
2998 xmap%x2(xmap%size)%k = k
2999 xmap%x2(xmap%size)%area = xmap%grids(g)%x(l)%area * xmap%grids(g)%x(l)%scale
3003 do k=1,xmap%grids(g)%km
3004 if (xmap%grids(g)%frac_area(i2,j2,k)/=0.0_r8_kind)
then
3005 xmap%size = xmap%size+1
3006 xmap%x1(xmap%size)%pos = xmap%ind_get1(ll)
3007 xmap%x1(xmap%size)%i = xmap%grids(g)%x(l)%i1
3008 xmap%x1(xmap%size)%j = xmap%grids(g)%x(l)%j1
3009 xmap%x1(xmap%size)%tile = xmap%grids(g)%x(l)%tile
3010 xmap%x1(xmap%size)%area = xmap%grids(g)%x(l)%area &
3011 *xmap%grids(g)%frac_area(i2,j2,k)
3012 xmap%x1(xmap%size)%di = xmap%grids(g)%x(l)%di
3013 xmap%x1(xmap%size)%dj = xmap%grids(g)%x(l)%dj
3014 xmap%x2(xmap%size)%i = xmap%grids(g)%x(l)%i2
3015 xmap%x2(xmap%size)%j = xmap%grids(g)%x(l)%j2
3016 xmap%x2(xmap%size)%k = k
3017 xmap%x2(xmap%size)%area = xmap%grids(g)%x(l)%area * xmap%grids(g)%x(l)%scale
3022 xmap%grids(g)%last = xmap%size
3026 if (max_size>
size(xmap%x1_put(:)))
then
3027 if (
associated(xmap%x1_put))
deallocate(xmap%x1_put)
3028 allocate( xmap%x1_put(1:max_size) )
3030 if (max_size>
size(xmap%x2_get(:)))
then
3031 if (
associated(xmap%x2_get))
deallocate(xmap%x2_get)
3032 allocate( xmap%x2_get(1:max_size) )
3035 do g=2,
size(xmap%grids(:))
3036 xmap%grids(g)%first_get = 1
3037 xmap%grids(g)%last_get = 0
3043 do g=2,
size(xmap%grids(:))
3044 xmap%grids(g)%first_get = xmap%size_get2 + 1;
3046 do l=1,xmap%grids(g)%size
3047 i1 = xmap%grids(g)%x(l)%i1
3048 j1 = xmap%grids(g)%x(l)%j1
3049 i2 = xmap%grids(g)%x(l)%i2
3050 j2 = xmap%grids(g)%x(l)%j2
3051 tile1 = xmap%grids(g)%x(l)%tile
3053 overlap_with_nest = .false.
3054 if( xmap%grids(1)%id ==
"ATM" .AND. tile1 == tile_parent .AND. &
3055 in_box(i1, j1, is_parent, ie_parent, js_parent, je_parent) ) overlap_with_nest = .true.
3056 if(xmap%grids(g)%is_ug)
then
3057 do k=1,xmap%grids(g)%km
3058 lll = xmap%grids(g)%l_index((j2-1)*xmap%grids(g)%im+i2)
3059 if (xmap%grids(g)%frac_area(lll,1,k)/=0.0_r8_kind)
then
3060 xmap%size_put1 = xmap%size_put1+1
3061 xmap%x1_put(xmap%size_put1)%pos = xmap%ind_put1(ll)
3062 xmap%x1_put(xmap%size_put1)%i = xmap%grids(g)%x(l)%i1
3063 xmap%x1_put(xmap%size_put1)%j = xmap%grids(g)%x(l)%j1
3064 xmap%x1_put(xmap%size_put1)%tile = xmap%grids(g)%x(l)%tile
3065 xmap%x1_put(xmap%size_put1)%area = xmap%grids(g)%x(l)%area &
3066 *xmap%grids(g)%frac_area(lll,1,k)
3067 xmap%x1_put(xmap%size_put1)%di = xmap%grids(g)%x(l)%di
3068 xmap%x1_put(xmap%size_put1)%dj = xmap%grids(g)%x(l)%dj
3069 if( .not. overlap_with_nest)
then
3070 xmap%size_get2 = xmap%size_get2+1
3071 xmap%x2_get(xmap%size_get2)%i = xmap%grids(g)%x(l)%i2
3072 xmap%x2_get(xmap%size_get2)%j = xmap%grids(g)%x(l)%j2
3073 xmap%x2_get(xmap%size_get2)%l = lll
3074 xmap%x2_get(xmap%size_get2)%k = k
3075 xmap%x2_get(xmap%size_get2)%area = xmap%grids(g)%x(l)%area * xmap%grids(g)%x(l)%scale
3076 xmap%x2_get(xmap%size_get2)%pos = xmap%size_put1
3081 do k=1,xmap%grids(g)%km
3082 if (xmap%grids(g)%frac_area(i2,j2,k)/=0.0_r8_kind)
then
3083 xmap%size_put1 = xmap%size_put1+1
3084 xmap%x1_put(xmap%size_put1)%pos = xmap%ind_put1(ll)
3085 xmap%x1_put(xmap%size_put1)%i = xmap%grids(g)%x(l)%i1
3086 xmap%x1_put(xmap%size_put1)%j = xmap%grids(g)%x(l)%j1
3087 xmap%x1_put(xmap%size_put1)%tile = xmap%grids(g)%x(l)%tile
3088 xmap%x1_put(xmap%size_put1)%area = xmap%grids(g)%x(l)%area &
3089 *xmap%grids(g)%frac_area(i2,j2,k)
3090 xmap%x1_put(xmap%size_put1)%di = xmap%grids(g)%x(l)%di
3091 xmap%x1_put(xmap%size_put1)%dj = xmap%grids(g)%x(l)%dj
3092 if( .not. overlap_with_nest)
then
3093 xmap%size_get2 = xmap%size_get2+1
3094 xmap%x2_get(xmap%size_get2)%i = xmap%grids(g)%x(l)%i2
3095 xmap%x2_get(xmap%size_get2)%j = xmap%grids(g)%x(l)%j2
3096 xmap%x2_get(xmap%size_get2)%k = k
3097 xmap%x2_get(xmap%size_get2)%area = xmap%grids(g)%x(l)%area * xmap%grids(g)%x(l)%scale
3098 xmap%x2_get(xmap%size_get2)%pos = xmap%size_put1
3104 xmap%grids(g)%last_get = xmap%size_get2
3109 if (xmap%get1_repro%nsend > 0)
then
3113 mypos =
mpp_pe() - mpp_root_pe()
3116 p = mod(mypos+m, npes)
3117 if( xmap%send_count_repro(p) > 0 )
then
3122 do g=2,
size(xmap%grids(:))
3123 do l=1,xmap%grids(g)%size
3124 p = xmap%grids(g)%x(l)%pe-xmap%root_pe
3128 xmap%get1_repro%send(n)%xLoc(pos) = xloc
3129 if( xmap%grids(g)%is_ug )
then
3130 i = xmap%grids(g)%x(l)%l2
3131 xloc = xloc + count(xmap%grids(g)%frac_area(i,1,:)/=0.0_r8_kind)
3133 i = xmap%grids(g)%x(l)%i2
3134 j = xmap%grids(g)%x(l)%j2
3135 xloc = xloc + count(xmap%grids(g)%frac_area(i,j,:)/=0.0_r8_kind)
3142 end subroutine regen
3173 real(r8_kind),
dimension(:,:,:),
intent(in) :: f
3174 character(len=3),
intent(in) :: grid_id
3175 type (xmap_type),
intent(inout) :: xmap
3178 type(
grid_type),
pointer,
save :: grid =>null()
3180 if (grid_id==xmap%grids(1)%id)
call error_mesg (
'xgrid_mod', &
3181 'set_frac_area called on side 1 grid', fatal)
3182 do g=2,
size(xmap%grids(:))
3183 grid => xmap%grids(g)
3184 if (grid_id==grid%id)
then
3185 if (
size(f,3)/=
size(grid%frac_area,3))
then
3186 if (
associated(grid%frac_area))
deallocate (grid%frac_area)
3187 grid%km =
size(f,3);
3188 allocate( grid%frac_area(grid%is_me:grid%ie_me, grid%js_me:grid%je_me, &
3197 call error_mesg (
'xgrid_mod',
'set_frac_area: could not find grid id', fatal)
3205 real(r8_kind),
dimension(:,:),
intent(in) :: f
3206 character(len=3),
intent(in) :: grid_id
3207 type (xmap_type),
intent(inout) :: xmap
3210 type(
grid_type),
pointer,
save :: grid =>null()
3212 if (grid_id==xmap%grids(1)%id)
call error_mesg (
'xgrid_mod', &
3213 'set_frac_area_ug called on side 1 grid', fatal)
3214 if (grid_id .NE.
'LND' )
call error_mesg (
'xgrid_mod', &
3215 .NE.
'set_frac_area_ug called for grid_id LND', fatal)
3216 do g=2,
size(xmap%grids(:))
3217 grid => xmap%grids(g)
3218 if (grid_id==grid%id)
then
3219 if (
size(f,2)/=
size(grid%frac_area,3))
then
3220 if (
associated(grid%frac_area))
deallocate (grid%frac_area)
3221 grid%km =
size(f,2);
3222 allocate( grid%frac_area(grid%ls_me:grid%le_me, 1, grid%km) )
3224 grid%frac_area(:,1,:) = f(:,:);
3230 call error_mesg (
'xgrid_mod',
'set_frac_area_ug: could not find grid id', fatal)
3248 real(r8_kind),
dimension(:,:),
intent(in) :: d
3249 character(len=3),
intent(in) :: grid_id
3250 real(r8_kind),
dimension(:),
intent(inout) :: x
3251 type (xmap_type),
intent(inout) :: xmap
3252 integer,
intent(in),
optional :: remap_method
3254 logical,
intent(in),
optional :: complete
3256 logical :: is_complete, set_mismatch
3257 integer :: g, method
3258 character(len=2) :: text
3259 integer,
save :: isize=0
3260 integer,
save :: jsize=0
3261 integer,
save :: lsize=0
3262 integer,
save :: xsize=0
3263 integer,
save :: method_saved=0
3264 character(len=3),
save :: grid_id_saved=
""
3265 integer(i8_kind),
dimension(MAX_FIELDS),
save :: d_addrs = -9999_i8_kind
3266 integer(i8_kind),
dimension(MAX_FIELDS),
save :: x_addrs = -9999_i8_kind
3268 if (grid_id==xmap%grids(1)%id)
then
3269 method = first_order
3270 if(
present(remap_method)) method = remap_method
3271 is_complete = .true.
3272 if(
present(complete)) is_complete=complete
3274 if( lsize > max_fields )
then
3275 write( text,
'(i2)' ) max_fields
3276 call error_mesg (
'xgrid_mod',
'MAX_FIELDS='//trim(text)//
' exceeded for group put_side1_to_xgrid', fatal)
3278 d_addrs(lsize) = loc(d)
3279 x_addrs(lsize) = loc(x)
3285 method_saved = method
3286 grid_id_saved = grid_id
3288 set_mismatch = .false.
3289 set_mismatch = set_mismatch .OR. (isize /=
size(d,1))
3290 set_mismatch = set_mismatch .OR. (jsize /=
size(d,2))
3291 set_mismatch = set_mismatch .OR. (xsize /=
size(x(:)))
3292 set_mismatch = set_mismatch .OR. (method_saved /= method)
3293 set_mismatch = set_mismatch .OR. (grid_id_saved /= grid_id)
3294 if(set_mismatch)
then
3295 write( text,
'(i2)' ) lsize
3296 call error_mesg (
'xgrid_mod',
'Incompatible field at count '//text//
' for group put_side1_to_xgrid', fatal )
3300 if(is_complete)
then
3303 if(monotonic_exchange .AND. grid_id ==
'ATM')
then
3304 call put_1_to_xgrid_order_2(d_addrs, x_addrs, xmap, isize, jsize, xsize, lsize)
3305 else if(method == first_order)
then
3306 call put_1_to_xgrid_order_1(d_addrs, x_addrs, xmap, isize, jsize, xsize, lsize)
3308 if(grid_id .NE.
'ATM')
call error_mesg (
'xgrid_mod', &
3309 "second order put_to_xgrid should only be applied to 'ATM' model, "//&
3310 "contact developer", fatal)
3311 call put_1_to_xgrid_order_2(d_addrs, x_addrs, xmap, isize, jsize, xsize, lsize)
3314 d_addrs = -9999_i8_kind
3315 x_addrs = -9999_i8_kind
3326 do g=2,
size(xmap%grids(:))
3327 if (grid_id==xmap%grids(g)%id) &
3328 call error_mesg (
'xgrid_mod', &
3329 'put_to_xgrid expects a 3D side 2 grid', fatal)
3332 call error_mesg (
'xgrid_mod',
'put_to_xgrid: could not find grid id', fatal)
3340 real(r8_kind),
dimension(:,:,:),
intent(in) :: d
3341 character(len=3),
intent(in) :: grid_id
3342 real(r8_kind),
dimension(:),
intent(inout) :: x
3343 type (xmap_type),
intent(inout) :: xmap
3347 if (grid_id==xmap%grids(1)%id) &
3348 call error_mesg (
'xgrid_mod', &
3349 'put_side2_to_xgrid expects a 3D side 2 grid', fatal)
3351 do g=2,
size(xmap%grids(:))
3352 if (grid_id==xmap%grids(g)%id)
then
3353 call put_2_to_xgrid(d, xmap%grids(g), x, xmap)
3358 call error_mesg (
'xgrid_mod',
'put_to_xgrid: could not find grid id', fatal)
3365 real(r8_kind),
dimension(:,:),
intent(out) :: d
3366 character(len=3),
intent(in) :: grid_id
3367 real(r8_kind),
dimension(:),
intent(in) :: x
3368 type (xmap_type),
intent(inout) :: xmap
3369 logical,
intent(in),
optional :: complete
3371 logical :: is_complete, set_mismatch
3373 character(len=2) :: text
3374 integer,
save :: isize=0
3375 integer,
save :: jsize=0
3376 integer,
save :: lsize=0
3377 integer,
save :: xsize=0
3378 character(len=3),
save :: grid_id_saved=
""
3379 integer(i8_kind),
dimension(MAX_FIELDS),
save :: d_addrs = -9999_i8_kind
3380 integer(i8_kind),
dimension(MAX_FIELDS),
save :: x_addrs = -9999_i8_kind
3383 if (grid_id==xmap%grids(1)%id)
then
3384 is_complete = .true.
3385 if(
present(complete)) is_complete=complete
3387 if( lsize > max_fields )
then
3388 write( text,
'(i2)' ) max_fields
3389 call error_mesg (
'xgrid_mod',
'MAX_FIELDS='//trim(text)//
' exceeded for group get_side1_from_xgrid', fatal)
3391 d_addrs(lsize) = loc(d)
3392 x_addrs(lsize) = loc(x)
3398 grid_id_saved = grid_id
3400 set_mismatch = .false.
3401 set_mismatch = set_mismatch .OR. (isize /=
size(d,1))
3402 set_mismatch = set_mismatch .OR. (jsize /=
size(d,2))
3403 set_mismatch = set_mismatch .OR. (xsize /=
size(x(:)))
3404 set_mismatch = set_mismatch .OR. (grid_id_saved /= grid_id)
3405 if(set_mismatch)
then
3406 write( text,
'(i2)' ) lsize
3407 call error_mesg (
'xgrid_mod',
'Incompatible field at count '//text// &
3408 &
' for group get_side1_from_xgrid', fatal )
3412 if(is_complete)
then
3414 call get_1_from_xgrid_repro(d_addrs, x_addrs, xmap, xsize, lsize)
3416 call get_1_from_xgrid(d_addrs, x_addrs, xmap, isize, jsize, xsize, lsize)
3418 d_addrs(1:lsize) = -9999
3419 x_addrs(1:lsize) = -9999
3429 do g=2,
size(xmap%grids(:))
3430 if (grid_id==xmap%grids(g)%id) &
3431 call error_mesg (
'xgrid_mod', &
3432 'get_from_xgrid expects a 3D side 2 grid', fatal)
3435 call error_mesg (
'xgrid_mod',
'get_from_xgrid: could not find grid id', fatal)
3442 real(r8_kind),
dimension(:,:,:),
intent(out) :: d
3443 character(len=3),
intent(in) :: grid_id
3444 real(r8_kind),
dimension(:),
intent(in) :: x
3445 type (xmap_type),
intent(in) :: xmap
3449 if (grid_id==xmap%grids(1)%id) &
3450 call error_mesg (
'xgrid_mod', &
3451 'get_from_xgrid expects a 2D side 1 grid', fatal)
3453 do g=2,
size(xmap%grids(:))
3454 if (grid_id==xmap%grids(g)%id)
then
3455 call get_2_from_xgrid(d, xmap%grids(g), x, xmap)
3460 call error_mesg (
'xgrid_mod',
'get_from_xgrid: could not find grid id', fatal)
3467 subroutine some(xmap, some_arr, grid_id)
3469 character(len=3),
optional,
intent(in) :: grid_id
3470 logical,
dimension(:),
intent(out) :: some_arr
3474 if (.not.
present(grid_id))
then
3476 if(xmap%size > 0)
then
3484 if (grid_id==xmap%grids(1)%id) &
3485 call error_mesg (
'xgrid_mod',
'some expects a side 2 grid id', fatal)
3487 do g=2,
size(xmap%grids(:))
3488 if (grid_id==xmap%grids(g)%id)
then
3490 some_arr(xmap%grids(g)%first:xmap%grids(g)%last) = .true.;
3495 call error_mesg (
'xgrid_mod',
'some could not find grid id', fatal)
3501 subroutine put_2_to_xgrid(d, grid, x, xmap)
3503 real(r8_kind),
dimension(grid%is_me:grid%ie_me, &
grid%js_me:grid%je_me, grid%km),
intent(in) :: d
3504 real(r8_kind),
dimension(:),
intent(inout) :: x
3508 call mpp_clock_begin(id_put_2_to_xgrid)
3510 do l=grid%first,grid%last
3511 x(l) = d(xmap%x2(l)%i,xmap%x2(l)%j,xmap%x2(l)%k)
3514 call mpp_clock_end(id_put_2_to_xgrid)
3515 end subroutine put_2_to_xgrid
3519 subroutine get_2_from_xgrid(d, grid, x, xmap)
3521 real(r8_kind),
dimension(grid%is_me:grid%ie_me, &
grid%js_me:grid%je_me, grid%km),
intent(out) :: d
3522 real(r8_kind),
dimension(:),
intent(in) :: x
3527 call mpp_clock_begin(id_get_2_from_xgrid)
3530 do l=grid%first_get,grid%last_get
3531 d(xmap%x2_get(l)%i,xmap%x2_get(l)%j,xmap%x2_get(l)%k) = &
3532 d(xmap%x2_get(l)%i,xmap%x2_get(l)%j,xmap%x2_get(l)%k) + xmap%x2_get(l)%area*x(xmap%x2_get(l)%pos)
3538 d(:,:,k) = d(:,:,k) * grid%area_inv
3541 call mpp_clock_end(id_get_2_from_xgrid)
3543 end subroutine get_2_from_xgrid
3547 subroutine put_1_to_xgrid_order_1(d_addrs, x_addrs, xmap, isize, jsize, xsize, lsize)
3548 integer(i8_kind),
dimension(:),
intent(in) :: d_addrs
3549 integer(i8_kind),
dimension(:),
intent(in) :: x_addrs
3551 integer,
intent(in) :: isize, jsize, xsize, lsize
3553 integer :: i, j, p, buffer_pos, msgsize
3554 integer :: from_pe, to_pe, pos, n, l, count
3555 integer :: ibegin, istart, iend, start_pos
3556 type (
comm_type),
pointer,
save :: comm =>null()
3557 real(r8_kind) :: recv_buffer(xmap%put1%recvsize*lsize)
3558 real(r8_kind) :: send_buffer(xmap%put1%sendsize*lsize)
3559 real(r8_kind) :: unpack_buffer(xmap%put1%recvsize)
3561 real(r8_kind),
dimension(isize, jsize) :: d
3562 real(r8_kind),
dimension(xsize) :: x
3566 call mpp_clock_begin(id_put_1_to_xgrid_order_1)
3570 do p = 1, comm%nrecv
3571 msgsize = comm%recv(p)%count*lsize
3572 from_pe = comm%recv(p)%pe
3573 buffer_pos = comm%recv(p)%buffer_pos*lsize
3574 call mpp_recv(recv_buffer(buffer_pos+1), glen=msgsize, from_pe = from_pe, block=.false., tag=comm_tag_7)
3579 do p = 1, comm%nsend
3580 msgsize = comm%send(p)%count*lsize
3581 to_pe = comm%send(p)%pe
3585 do n = 1, comm%send(p)%count
3587 i = comm%send(p)%i(n)
3588 j = comm%send(p)%j(n)
3589 send_buffer(pos) = d(i,j)
3592 call mpp_send(send_buffer(buffer_pos+1), plen=msgsize, to_pe = to_pe, tag=comm_tag_7 )
3593 buffer_pos = buffer_pos + msgsize
3599 if( lsize == 1)
then
3601 do l=1,xmap%size_put1
3602 x(l) = recv_buffer(xmap%x1_put(l)%pos)
3610 do p = 1, comm%nrecv
3611 count = comm%recv(p)%count
3612 ibegin = comm%recv(p)%buffer_pos*lsize + 1
3613 istart = ibegin + (l-1)*count
3614 iend = istart + count - 1
3615 pos = comm%recv(p)%buffer_pos
3618 unpack_buffer(pos) = recv_buffer(n)
3621 do i=1,xmap%size_put1
3622 x(i) = unpack_buffer(xmap%x1_put(i)%pos)
3629 call mpp_clock_end(id_put_1_to_xgrid_order_1)
3631 end subroutine put_1_to_xgrid_order_1
3636 subroutine put_1_to_xgrid_order_2(d_addrs, x_addrs, xmap, isize, jsize, xsize, lsize)
3637 integer(i8_kind),
dimension(:),
intent(in) :: d_addrs
3638 integer(i8_kind),
dimension(:),
intent(in) :: x_addrs
3640 integer,
intent(in) :: isize, jsize, xsize, lsize
3643 real(r8_kind),
dimension(0:isize+1, 0:jsize+1, lsize) :: tmp
3644 real(r8_kind),
dimension(isize, jsize, lsize) :: tmpx, tmpy
3645 real(r8_kind),
dimension(isize, jsize, lsize) :: d_bar_max, d_bar_min
3646 real(r8_kind),
dimension(isize, jsize, lsize) :: d_max, d_min
3647 real(r8_kind) :: d_bar
3648 integer :: i, is, ie, j, js, je, ii, jj
3649 integer :: p, l, isd, jsd
3650 type (
grid_type),
pointer,
save :: grid1 =>null()
3651 type (
comm_type),
pointer,
save :: comm =>null()
3652 integer :: buffer_pos, msgsize, from_pe, to_pe, pos, n
3653 integer :: ibegin, count, istart, iend
3654 real(r8_kind) :: recv_buffer(xmap%put1%recvsize*lsize*3)
3655 real(r8_kind) :: send_buffer(xmap%put1%sendsize*lsize*3)
3656 real(r8_kind) :: unpack_buffer(xmap%put1%recvsize*3)
3657 logical :: on_west_edge, on_east_edge, on_south_edge, on_north_edge
3658 real(r8_kind),
dimension(isize, jsize) :: d
3659 real(r8_kind),
dimension(xsize) :: x
3663 call mpp_clock_begin(id_put_1_to_xgrid_order_2)
3664 grid1 => xmap%grids(1)
3666 is = grid1%is_me; ie = grid1%ie_me
3667 js = grid1%js_me; je = grid1%je_me
3673 tmp(:,:,l) = large_number
3675 tmp(1:isize,1:jsize,l) = d(:,:)
3678 if(grid1%is_latlon)
then
3679 call mpp_update_domains(tmp,grid1%domain_with_halo)
3682 tmpy(:,:,l) =
grad_merid_latlon(tmp(:,:,l), grid1%lat, is, ie, js, je, isd, jsd)
3683 tmpx(:,:,l) =
grad_zonal_latlon(tmp(:,:,l), grid1%lon, grid1%lat, is, ie, js, je, isd, jsd)
3686 call mpp_update_domains(tmp,grid1%domain)
3687 on_west_edge = (is==1)
3688 on_east_edge = (ie==grid1%im)
3689 on_south_edge = (js==1)
3690 on_north_edge = (je==grid1%jm)
3694 call gradient_cubic(tmp(:,:,l), grid1%box%dx, grid1%box%dy, grid1%box%area, &
3695 grid1%box%edge_w, grid1%box%edge_e, grid1%box%edge_s, &
3696 grid1%box%edge_n, grid1%box%en1, grid1%box%en2, &
3697 grid1%box%vlon, grid1%box%vlat, tmpx(:,:,l), tmpy(:,:,l), &
3698 on_west_edge, on_east_edge, on_south_edge, on_north_edge)
3705 do p = 1, comm%nrecv
3706 msgsize = comm%recv(p)%count*lsize
3707 buffer_pos = comm%recv(p)%buffer_pos*lsize
3708 if(.NOT. monotonic_exchange)
then
3710 buffer_pos = buffer_pos*3
3712 from_pe = comm%recv(p)%pe
3713 call mpp_recv(recv_buffer(buffer_pos+1), glen=msgsize, from_pe = from_pe, block=.false., tag=comm_tag_8)
3717 if(monotonic_exchange)
then
3722 d_bar_max(i,j,l) = -large_number
3723 d_bar_min(i,j,l) = large_number
3724 d_max(i,j,l) = -large_number
3725 d_min(i,j,l) = large_number
3728 if(tmp(i,j,l) .NE. large_number)
then
3729 if(tmp(i,j,l) > d_bar_max(i,j,l)) d_bar_max(i,j,l) = tmp(i,j,l)
3730 if(tmp(i,j,l) < d_bar_min(i,j,l)) d_bar_min(i,j,l) = tmp(i,j,l)
3741 if(monotonic_exchange)
then
3743 do p = 1, comm%nsend
3744 msgsize = comm%send(p)%count*lsize
3745 to_pe = comm%send(p)%pe
3748 do n = 1, comm%send(p)%count
3750 i = comm%send(p)%i(n)
3751 j = comm%send(p)%j(n)
3752 send_buffer(pos) = d(i,j) + tmpy(i,j,l)*comm%send(p)%dj(n) + tmpx(i,j,l)*comm%send(p)%di(n)
3753 if(send_buffer(pos) > d_max(i,j,l)) d_max(i,j,l) = send_buffer(pos)
3754 if(send_buffer(pos) < d_min(i,j,l)) d_min(i,j,l) = send_buffer(pos)
3759 do p = 1, comm%nsend
3760 msgsize = comm%send(p)%count*lsize
3761 to_pe = comm%send(p)%pe
3765 do n = 1, comm%send(p)%count
3767 i = comm%send(p)%i(n)
3768 j = comm%send(p)%j(n)
3770 if( d_max(i,j,l) > d_bar_max(i,j,l) )
then
3771 send_buffer(pos) = d_bar + ((send_buffer(pos)-d_bar)/(d_max(i,j,l)-d_bar)) * (d_bar_max(i,j,l)-d_bar)
3772 else if( d_min(i,j,l) < d_bar_min(i,j,l) )
then
3773 send_buffer(pos) = d_bar + ((send_buffer(pos)-d_bar)/(d_min(i,j,l)-d_bar)) * (d_bar_min(i,j,l)-d_bar)
3777 call mpp_send(send_buffer(buffer_pos+1), plen=msgsize, to_pe = to_pe, tag=comm_tag_8 )
3778 buffer_pos = buffer_pos + msgsize
3781 do p = 1, comm%nsend
3782 msgsize = comm%send(p)%count*lsize*3
3783 to_pe = comm%send(p)%pe
3787 do n = 1, comm%send(p)%count
3789 i = comm%send(p)%i(n)
3790 j = comm%send(p)%j(n)
3791 send_buffer(pos-2) = d(i,j)
3792 send_buffer(pos-1) = tmpy(i,j,l)
3793 send_buffer(pos ) = tmpx(i,j,l)
3796 call mpp_send(send_buffer(buffer_pos+1), plen=msgsize, to_pe = to_pe, tag=comm_tag_8 )
3797 buffer_pos = buffer_pos + msgsize
3804 if(monotonic_exchange)
then
3805 if( lsize == 1)
then
3807 do l=1,xmap%size_put1
3808 pos = xmap%x1_put(l)%pos
3809 x(l) = recv_buffer(pos)
3815 do p = 1, comm%nsend
3816 count = comm%send(p)%count
3817 ibegin = comm%recv(p)%buffer_pos*lsize + 1
3818 istart = ibegin + (l-1)*count
3819 iend = istart + count - 1
3820 pos = comm%recv(p)%buffer_pos
3823 unpack_buffer(pos) = recv_buffer(n)
3826 do i=1,xmap%size_put1
3827 pos = xmap%x1_put(i)%pos
3828 x(i) = unpack_buffer(pos)
3833 if( lsize == 1)
then
3836 do l=1,xmap%size_put1
3837 pos = xmap%x1_put(l)%pos
3838 x(l) = recv_buffer(3*pos-2) + recv_buffer(3*pos-1)*xmap%x1_put(l)%dj + recv_buffer(3*pos)*xmap%x1_put(l)%di
3847 do p = 1, comm%nrecv
3848 count = comm%recv(p)%count*3
3849 ibegin = comm%recv(p)%buffer_pos*lsize*3 + 1
3850 istart = ibegin + (l-1)*count
3851 iend = istart + count - 1
3852 pos = comm%recv(p)%buffer_pos*3
3855 unpack_buffer(pos) = recv_buffer(n)
3858 do i=1,xmap%size_put1
3859 pos = xmap%x1_put(i)%pos
3860 x(i) = unpack_buffer(3*pos-2) + unpack_buffer(3*pos-1)*xmap%x1_put(i)%dj + unpack_buffer(3*pos) &
3861 & * xmap%x1_put(i)%di
3868 call mpp_clock_end(id_put_1_to_xgrid_order_2)
3870 end subroutine put_1_to_xgrid_order_2
3874 subroutine get_1_from_xgrid(d_addrs, x_addrs, xmap, isize, jsize, xsize, lsize)
3875 integer(i8_kind),
dimension(:),
intent(in) :: d_addrs
3876 integer(i8_kind),
dimension(:),
intent(in) :: x_addrs
3878 integer,
intent(in) :: isize, jsize, xsize, lsize
3880 real(r8_kind),
dimension(xmap%size),
target :: dg(xmap%size, lsize)
3881 integer :: i, j, l, p, n, m
3882 integer :: msgsize, buffer_pos, pos
3883 integer :: istart, iend, count
3884 real(r8_kind) ,
pointer,
save :: dgp =>null()
3885 type(
grid_type) ,
pointer,
save :: grid1 =>null()
3886 type(
comm_type) ,
pointer,
save :: comm =>null()
3889 real(r8_kind) :: recv_buffer(xmap%get1%recvsize*lsize*3)
3890 real(r8_kind) :: send_buffer(xmap%get1%sendsize*lsize*3)
3891 real(r8_kind) :: d(isize,jsize)
3892 real(r8_kind),
dimension(xsize) :: x
3896 call mpp_clock_begin(id_get_1_from_xgrid)
3899 grid1 => xmap%grids(1)
3901 do p = 1, comm%nrecv
3902 recv => comm%recv(p)
3903 msgsize = recv%count*lsize
3904 buffer_pos = recv%buffer_pos*lsize
3905 call mpp_recv(recv_buffer(buffer_pos+1), glen=msgsize, from_pe = recv%pe, block=.false., tag=comm_tag_9)
3913 dgp => dg(xmap%x1(i)%pos,l)
3914 dgp = dgp + xmap%x1(i)%area*x(i)
3922 do p = 1, comm%nsend
3923 send => comm%send(p)
3924 msgsize = send%count*lsize
3926 istart = send%buffer_pos+1
3927 iend = istart + send%count - 1
3931 send_buffer(pos) = dg(n,l)
3934 call mpp_send(send_buffer(buffer_pos+1), plen=msgsize, to_pe = send%pe, tag=comm_tag_9 )
3935 buffer_pos = buffer_pos + msgsize
3948 do p = 1, comm%nrecv
3949 recv => comm%recv(p)
3951 buffer_pos = recv%buffer_pos*lsize
3952 if( recv%pe == xmap%me )
then
3956 pos = buffer_pos + (l-1)*count
3962 d(i,j) = recv_buffer(pos)
3970 do m = 1, comm%nrecv
3971 p = comm%unpack_ind(m)
3972 recv => comm%recv(p)
3973 if( recv%pe == xmap%me )
then
3976 buffer_pos = recv%buffer_pos*lsize
3980 pos = buffer_pos + (l-1)*recv%count
3982 do n = 1, recv%count
3986 d(i,j) = d(i,j) + recv_buffer(pos)
3997 d = d * grid1%area_inv
4000 call mpp_clock_end(id_get_1_from_xgrid)
4002 end subroutine get_1_from_xgrid
4006 subroutine get_1_from_xgrid_repro(d_addrs, x_addrs, xmap, xsize, lsize)
4007 integer(i8_kind),
dimension(:),
intent(in) :: d_addrs
4008 integer(i8_kind),
dimension(:),
intent(in) :: x_addrs
4010 integer,
intent(in) :: xsize, lsize
4012 integer :: g, i, j, k, p, l, n, l2, l3
4013 integer :: msgsize, buffer_pos, pos
4014 type (
grid_type),
pointer,
save :: grid =>null()
4015 type(
comm_type),
pointer,
save :: comm => null()
4018 integer,
dimension(0:xmap%npes-1) :: pl, ml
4019 real(r8_kind) :: recv_buffer(xmap%recv_count_repro_tot*lsize)
4020 real(r8_kind) :: send_buffer(xmap%send_count_repro_tot*lsize)
4021 real(r8_kind) :: d(xmap%grids(1)%is_me:xmap%grids(1)%ie_me, &
4022 xmap%grids(1)%js_me:xmap%grids(1)%je_me)
4023 real(r8_kind),
dimension(xsize) :: x
4027 call mpp_clock_begin(id_get_1_from_xgrid_repro)
4028 comm => xmap%get1_repro
4030 do p = 1, comm%nrecv
4031 recv => comm%recv(p)
4032 msgsize = recv%count*lsize
4033 buffer_pos = recv%buffer_pos*lsize
4034 call mpp_recv(recv_buffer(buffer_pos+1), glen=msgsize, from_pe = recv%pe, block=.false., tag=comm_tag_10)
4035 n = recv%pe -xmap%root_pe
4041 send_buffer(:) = 0.0_r8_kind
4044 do p = 1, comm%nsend
4045 pos = comm%send(p)%buffer_pos*lsize
4046 send => comm%send(p)
4049 do n = 1, send%count
4055 do k =1, xmap%grids(g)%km
4056 if(xmap%grids(g)%frac_area(i,j,k)/=0.0_r8_kind)
then
4058 send_buffer(pos) = send_buffer(pos) + xmap%x1(l2)%area *x(l2)
4066 buffer_pos = comm%send(p)%buffer_pos*lsize
4067 msgsize = comm%send(p)%count*lsize
4068 call mpp_send(send_buffer(buffer_pos+1), plen=msgsize, to_pe=comm%send(p)%pe, tag=comm_tag_10)
4082 do g=2,
size(xmap%grids(:))
4083 grid => xmap%grids(g)
4084 do l3=1,grid%size_repro
4085 i = grid%x_repro(l3)%i1
4086 j = grid%x_repro(l3)%j1
4087 p = grid%x_repro(l3)%pe-xmap%root_pe
4088 pos = pl(p) + (l-1)*ml(p) + grid%x_repro(l3)%recv_pos
4089 d(i,j) = d(i,j) + recv_buffer(pos)
4093 d = d * xmap%grids(1)%area_inv
4098 call mpp_clock_end(id_get_1_from_xgrid_repro)
4100 end subroutine get_1_from_xgrid_repro
4109 real(r8_kind),
dimension(:,:),
intent(in) :: d
4110 character(len=3),
intent(in) :: grid_id
4113 integer,
intent(in),
optional :: remap_method
4116 real(r8_kind),
dimension(xmap%size) :: x_over, x_back
4117 real(r8_kind),
dimension(size(d,1),size(d,2)) :: d1
4118 real(r8_kind),
dimension(:,:,:),
allocatable :: d2
4120 type (
grid_type),
pointer,
save :: grid1 =>null(), grid2 =>null()
4122 grid1 => xmap%grids(1)
4128 call put_to_xgrid (d, grid1%id, x_over, xmap, remap_method)
4129 do g=2,
size(xmap%grids(:))
4130 grid2 => xmap%grids(g)
4131 if(grid2%on_this_pe)
then
4132 allocate (d2(grid2%is_me:grid2%ie_me, grid2%js_me:grid2%je_me, grid2%km) )
4135 if(grid2%on_this_pe)
then
4139 if(
allocated(d2))
deallocate (d2)
4156 real(r8_kind),
dimension(:,:,:),
intent(in) :: d
4157 character(len=3),
intent(in) :: grid_id
4160 integer,
intent(in),
optional :: remap_method
4163 real(r8_kind),
dimension(xmap%size) :: x_over, x_back
4164 real(r8_kind),
dimension(:,: ),
allocatable :: d1
4165 real(r8_kind),
dimension(:,:,:),
allocatable :: d2
4167 type (
grid_type),
pointer,
save :: grid1 =>null(), grid2 =>null()
4169 grid1 => xmap%grids(1)
4171 do g = 2,
size(xmap%grids(:))
4172 grid2 => xmap%grids(g)
4173 if (grid_id==grid2%id)
then
4174 if(grid2%on_this_pe)
then
4179 call put_to_xgrid(0.0_r8_kind * grid2%frac_area, grid2%id, x_over, xmap)
4183 allocate ( d1(
size(grid1%area,1),
size(grid1%area,2)) )
4186 call put_to_xgrid(d1, grid1%id, x_back, xmap,remap_method)
4190 do g = 2,
size(xmap%grids(:))
4191 grid2 => xmap%grids(g)
4192 if(grid2%on_this_pe)
then
4193 allocate ( d2(
size(grid2%frac_area, 1),
size(grid2%frac_area, 2), &
4194 size(grid2%frac_area, 3) ) )
4198 if(
allocated(d2) )
deallocate ( d2 )
4212 real(r8_kind),
dimension(:,:),
intent(in) :: d
4213 character(len=3),
intent(in) :: grid_id
4216 integer,
intent(in),
optional :: remap_method
4218 real(r8_kind),
dimension(xmap%size) :: x_over, x_back
4219 real(r8_kind),
dimension(size(d,1),size(d,2)) :: d1
4220 real(r8_kind),
dimension(:,:,:),
allocatable :: d2
4221 real(r8_kind),
dimension(: ),
allocatable :: d_ug
4222 real(r8_kind),
dimension(:,:),
allocatable :: d2_ug
4224 type (
grid_type),
pointer,
save :: grid1 =>null(), grid2 =>null()
4226 grid1 => xmap%grids(1)
4230 if(grid1%is_ug)
then
4231 allocate(d_ug(grid1%ls_me:grid1%le_me))
4232 call mpp_pass_sg_to_ug(grid1%ug_domain, d, d_ug)
4237 call put_to_xgrid (d, grid1%id, x_over, xmap, remap_method)
4239 do g=2,
size(xmap%grids(:))
4240 grid2 => xmap%grids(g)
4241 if(grid2%is_ug)
then
4242 if(grid2%on_this_pe)
then
4243 allocate (d2_ug(grid2%ls_me:grid2%le_me, grid2%km) )
4247 if(grid2%on_this_pe)
then
4249 sum( grid2%area(:,1) * sum(grid2%frac_area(:,1,:)*d2_ug,dim=2) )
4252 if(
allocated(d2_ug))
deallocate (d2_ug)
4254 if(grid2%on_this_pe)
then
4255 allocate (d2(grid2%is_me:grid2%ie_me, grid2%js_me:grid2%je_me, grid2%km) )
4258 if(grid2%on_this_pe)
then
4260 & + sum( grid2%area * sum(grid2%frac_area*d2,dim=3) )
4263 if(
allocated(d2))
deallocate (d2)
4266 if(grid1%is_ug)
then
4273 if(
allocated(d_ug))
deallocate(d_ug)
4285 real(r8_kind),
dimension(:,:,:),
intent(in) :: d
4286 character(len=3),
intent(in) :: grid_id
4289 integer,
intent(in),
optional :: remap_method
4292 real(r8_kind),
dimension(xmap%size) :: x_over, x_back
4293 real(r8_kind),
dimension(:,: ),
allocatable :: d1, d_ug
4294 real(r8_kind),
dimension(:,:,:),
allocatable :: d2
4296 type (
grid_type),
pointer,
save :: grid1 =>null(), grid2 =>null()
4298 grid1 => xmap%grids(1)
4300 do g = 2,
size(xmap%grids(:))
4301 grid2 => xmap%grids(g)
4302 if (grid_id==grid2%id)
then
4303 if(grid2%on_this_pe)
then
4304 if(grid2%is_ug)
then
4305 allocate(d_ug(grid2%ls_me:grid2%le_me,grid2%km))
4306 call mpp_pass_sg_to_ug(grid2%ug_domain, d, d_ug)
4312 if(grid2%is_ug)
then
4317 if(
allocated(d_ug))
deallocate(d_ug)
4319 if(grid2%is_ug)
then
4320 call put_to_xgrid_ug(0.0_r8_kind * grid2%frac_area(:,1,:), grid2%id, x_over, xmap)
4322 call put_to_xgrid(0.0_r8_kind * grid2%frac_area, grid2%id, x_over, xmap)
4327 allocate ( d1(
size(grid1%area,1),
size(grid1%area,2)) )
4328 if(grid1%is_ug)
then
4334 if(grid1%is_ug)
then
4337 call put_to_xgrid(d1, grid1%id, x_back, xmap,remap_method)
4342 do g = 2,
size(xmap%grids(:))
4343 grid2 => xmap%grids(g)
4344 if(grid2%on_this_pe)
then
4345 allocate ( d2(
size(grid2%frac_area, 1),
size(grid2%frac_area, 2), &
4346 size(grid2%frac_area, 3) ) )
4348 if(grid2%is_ug)
then
4354 if(
allocated(d2) )
deallocate ( d2 )
4365 character(len=3),
intent(in) :: id
4367 real(r8_kind),
dimension(:,:),
intent(out) :: area
4372 do g = 1,
size(xmap%grids(:))
4373 if (id==xmap%grids(g)%id )
then
4374 if(
size(area,1) .NE.
size(xmap%grids(g)%area,1) .OR.
size(area,2) .NE.
size(xmap%grids(g)%area,2) ) &
4375 call error_mesg(
"xgrid_mod",
"size mismatch between area and xmap%grids(g)%area", fatal)
4376 area = xmap%grids(g)%area
4382 if(.not. found)
call error_mesg(
"xgrid_mod", id//
" is not found in xmap%grids id", fatal)
4394 integer,
intent(in) :: isd, jsd
4395 real(r8_kind),
dimension(isd:,jsd:),
intent(in) :: d
4396 real(r8_kind),
dimension(:),
intent(in) :: lon
4397 real(r8_kind),
dimension(:),
intent(in) :: lat
4398 integer,
intent(in) :: is, ie, js, je
4400 real(r8_kind) :: dx, costheta
4401 integer :: i, j, ip1, im1
4407 else if(i==
size(lon(:)) )
then
4410 ip1 = i+1; im1 = i-1
4412 dx = lon(ip1) - lon(im1)
4413 if(abs(dx).lt.eps )
call error_mesg(
'xgrids_mod(grad_zonal_latlon)',
'Improper grid size in lontitude', fatal)
4414 if(dx .gt. pi) dx = dx - 2.0_r8_kind* pi
4415 if(dx .lt. -pi) dx = dx + 2.0_r8_kind* pi
4417 costheta = cos(lat(j))
4418 if(abs(costheta) .lt. eps)
call error_mesg(
'xgrids_mod(grad_zonal_latlon)',
'Improper latitude grid', fatal)
4433 integer,
intent(in) :: isd, jsd
4434 real(r8_kind),
dimension(isd:,jsd:),
intent(in) :: d
4435 real(r8_kind),
dimension(:),
intent(in) :: lat
4436 integer,
intent(in) :: is, ie, js, je
4439 integer :: i, j, jp1, jm1
4445 else if(j ==
size(lat(:)) )
then
4448 jp1 = j+1; jm1 = j-1
4450 dy = lat(jp1) - lat(jm1)
4451 if(abs(dy).lt.eps)
call error_mesg(
'xgrids_mod(grad_merid_latlon)',
'Improper grid size in latitude', fatal)
4462 subroutine get_index_range(xmap, grid_index, is, ie, js, je, km)
4465 integer,
intent(in) :: grid_index
4466 integer,
intent(out) :: is, ie, js, je, km
4468 is = xmap % grids(grid_index) % is_me
4469 ie = xmap % grids(grid_index) % ie_me
4470 js = xmap % grids(grid_index) % js_me
4471 je = xmap % grids(grid_index) % je_me
4472 km = xmap % grids(grid_index) % km
4474 end subroutine get_index_range
4481 subroutine stock_move_3d(from, to, grid_index, stock_data3d, xmap, &
4482 & delta_t, from_side, to_side, radius, verbose, ier)
4492 type(stock_type),
intent(inout),
optional :: from, to
4493 integer,
intent(in) :: grid_index
4494 real(r8_kind),
intent(in) :: stock_data3d(:,:,:)
4496 real(r8_kind),
intent(in) :: delta_t
4497 integer,
intent(in) :: from_side
4498 integer,
intent(in) :: to_side
4499 real(r8_kind),
intent(in) :: radius
4500 character(len=*),
intent(in),
optional :: verbose
4501 integer,
intent(out) :: ier
4503 real(r8_kind) :: from_dq, to_dq
4506 if(grid_index == 1)
then
4512 if(.not.
associated(xmap%grids) )
then
4517 from_dq = delta_t * 4.0_r8_kind * pi * radius**2 * sum( sum(xmap%grids(grid_index)%area * &
4518 & sum(xmap%grids(grid_index)%frac_area * stock_data3d, dim=3), dim=1))
4522 if(
present(to )) to % dq( to_side) = to % dq( to_side) + to_dq
4523 if(
present(from)) from % dq(from_side) = from % dq(from_side) - from_dq
4525 if(
present(verbose).and.debug_stocks)
then
4528 from_dq = from_dq/(4.0_r8_kind*pi*radius**2)
4529 to_dq = to_dq /(4.0_r8_kind*pi*radius**2)
4530 if(
mpp_pe()==mpp_root_pe())
then
4531 write(stocks_file,
'(a,es19.12,a,es19.12,a)') verbose, from_dq,
' [*/m^2]'
4541 subroutine stock_move_2d(from, to, grid_index, stock_data2d, xmap, &
4542 & delta_t, from_side, to_side, radius, verbose, ier)
4551 type(stock_type),
intent(inout),
optional :: from, to
4552 integer,
optional,
intent(in) :: grid_index
4553 real(r8_kind),
intent(in) :: stock_data2d(:,:)
4555 real(r8_kind),
intent(in) :: delta_t
4556 integer,
intent(in) :: from_side
4557 integer,
intent(in) :: to_side
4558 real(r8_kind),
intent(in) :: radius
4559 character(len=*),
intent(in) :: verbose
4560 integer,
intent(out) :: ier
4562 real(r8_kind) :: to_dq, from_dq
4566 if(.not.
associated(xmap%grids) )
then
4571 if( .not.
present(grid_index) .or. grid_index==1 )
then
4574 from_dq = delta_t * 4.0_r8_kind*pi*radius**2 * sum(sum(xmap%grids(1)%area * stock_data2d, dim=1))
4585 if(
present(to )) to % dq( to_side) = to % dq( to_side) + to_dq
4586 if(
present(from)) from % dq(from_side) = from % dq(from_side) - from_dq
4588 if(debug_stocks)
then
4591 from_dq = from_dq/(4.0_r8_kind*pi*radius**2)
4592 to_dq = to_dq /(4.0_r8_kind*pi*radius**2)
4593 if(
mpp_pe()==mpp_root_pe())
then
4594 write(stocks_file,
'(a,es19.12,a,es19.12,a)') verbose, from_dq,
' [*/m^2]'
4598 end subroutine stock_move_2d
4605 subroutine stock_move_ug_3d(from, to, grid_index, stock_ug_data3d, xmap, &
4606 & delta_t, from_side, to_side, radius, verbose, ier)
4616 type(stock_type),
intent(inout),
optional :: from, to
4617 integer,
intent(in) :: grid_index
4618 real(r8_kind),
intent(in) :: stock_ug_data3d(:,:)
4620 real(r8_kind),
intent(in) :: delta_t
4621 integer,
intent(in) :: from_side
4622 integer,
intent(in) :: to_side
4623 real(r8_kind),
intent(in) :: radius
4624 character(len=*),
intent(in),
optional :: verbose
4625 integer,
intent(out) :: ier
4626 real(r8_kind),
dimension(size(stock_ug_data3d,1),size(stock_ug_data3d,2)) :: tmp
4628 real(r8_kind) :: from_dq, to_dq
4631 if(grid_index == 1)
then
4637 if(.not.
associated(xmap%grids) )
then
4642 tmp = xmap%grids(grid_index)%frac_area(:,1,:) * stock_ug_data3d
4643 from_dq = delta_t * 4.0_r8_kind * pi * radius**2 * sum( xmap%grids(grid_index)%area(:,1) * &
4648 if(
present(to )) to % dq( to_side) = to % dq( to_side) + to_dq
4649 if(
present(from)) from % dq(from_side) = from % dq(from_side) - from_dq
4651 if(
present(verbose).and.debug_stocks)
then
4654 from_dq = from_dq/(4.0_r8_kind*pi*radius**2)
4655 to_dq = to_dq /(4.0_r8_kind*pi*radius**2)
4656 if(
mpp_pe()==mpp_root_pe())
then
4657 write(stocks_file,
'(a,es19.12,a,es19.12,a)') verbose, from_dq,
' [*/m^2]'
4661 end subroutine stock_move_ug_3d
4667 subroutine stock_integrate_2d(integrate_data2d, xmap, delta_t, radius, res, ier)
4673 real(r8_kind),
intent(in) :: integrate_data2d(:,:)
4675 real(r8_kind),
intent(in) :: delta_t
4676 real(r8_kind),
intent(in) :: radius
4677 real(r8_kind),
intent(out) :: res
4678 integer,
intent(out) :: ier
4683 if(.not.
associated(xmap%grids) )
then
4688 res = delta_t * 4.0_r8_kind * pi * radius**2 * sum(sum(xmap%grids(1)%area * integrate_data2d, dim=1))
4690 end subroutine stock_integrate_2d
4697 subroutine stock_print(stck, Time, comp_name, index, ref_value, radius, pelist)
4703 type(stock_type),
intent(in) :: stck
4705 character(len=*) :: comp_name
4706 integer,
intent(in) :: index
4707 real(r8_kind),
intent(in) :: ref_value
4708 real(r8_kind),
intent(in) :: radius
4709 integer,
intent(in),
optional :: pelist(:)
4711 integer,
parameter :: initID = -2
4715 real(r8_kind) :: f_value, c_value, planet_area
4716 character(len=80) :: formatString
4717 integer :: iday, isec, hours
4718 integer :: diagID, compInd
4719 integer,
dimension(NELEMS,4),
save :: f_valueDiagID = initid
4720 integer,
dimension(NELEMS,4),
save :: c_valueDiagID = initid
4721 integer,
dimension(NELEMS,4),
save :: fmc_valueDiagID = initid
4723 real(r8_kind) :: diagField
4725 character(len=30) :: field_name, units
4727 f_value = sum(stck % dq)
4728 c_value = ref_value - stck % q_start
4729 if(
present(pelist))
then
4730 call mpp_sum(f_value, pelist=pelist)
4731 call mpp_sum(c_value, pelist=pelist)
4737 if(
mpp_pe() == mpp_root_pe())
then
4739 planet_area = 4.0_r8_kind * pi * radius**2
4740 f_value = f_value / planet_area
4741 c_value = c_value / planet_area
4743 if(comp_name ==
'ATM') compind = 1
4744 if(comp_name ==
'LND') compind = 2
4745 if(comp_name ==
'ICE') compind = 3
4746 if(comp_name ==
'OCN') compind = 4
4749 if(f_valuediagid(index,compind) == initid)
then
4750 field_name = trim(comp_name) // trim(stock_names(index))
4751 field_name = trim(field_name) //
'StocksChange_Flux'
4752 units = trim(stock_units(index))
4757 if(c_valuediagid(index,compind) == initid)
then
4758 field_name = trim(comp_name) // trim(stock_names(index))
4759 field_name = trim(field_name) //
'StocksChange_Comp'
4760 units = trim(stock_units(index))
4765 if(fmc_valuediagid(index,compind) == initid)
then
4766 field_name = trim(comp_name) // trim(stock_names(index))
4767 field_name = trim(field_name) //
'StocksChange_Diff'
4768 units = trim(stock_units(index))
4773 diagid=f_valuediagid(index,compind)
4775 if (diagid > 0) used =
send_data(diagid, diagfield, time = time)
4776 diagid=c_valuediagid(index,compind)
4778 if (diagid > 0) used =
send_data(diagid, diagfield, time)
4779 diagid=fmc_valuediagid(index,compind)
4780 diagfield = f_value-c_value
4781 if (diagid > 0) used =
send_data(diagid, diagfield, time=time)
4785 hours = iday*24 + isec/3600
4786 formatstring =
'(a,a,a,i16,2x,es22.15,2x,es22.15,2x,es22.15)'
4787 write(stocks_file,formatstring) trim(comp_name),stock_names(index),stock_units(index) &
4788 ,hours,f_value,c_value,f_value-c_value
4793 end subroutine stock_print
4798 function is_lat_lon(lon, lat)
4799 real(r8_kind),
dimension(:,:),
intent(in) :: lon, lat
4800 logical :: is_lat_lon
4801 integer :: i, j, nlon, nlat, num
4806 loop_lat:
do j = 1, nlat
4808 if(lat(i,j) .NE. lat(1,j))
then
4809 is_lat_lon = .false.
4816 loop_lon:
do i = 1, nlon
4818 if(lon(i,j) .NE. lon(i,1))
then
4819 is_lat_lon = .false.
4827 if(is_lat_lon) num = 1
4832 is_lat_lon = .false.
4836 end function is_lat_lon
4846 subroutine get_side1_from_xgrid_ug(d, grid_id, x, xmap, complete)
4847 real(r8_kind),
dimension(:),
intent(out) :: d
4848 character(len=3),
intent(in) :: grid_id
4849 real(r8_kind),
dimension(:),
intent(in) :: x
4850 type (xmap_type),
intent(inout) :: xmap
4851 logical,
intent(in),
optional :: complete
4853 logical :: is_complete, set_mismatch
4855 character(len=2) :: text
4856 integer,
save :: isize=0
4857 integer,
save :: lsize=0
4858 integer,
save :: xsize=0
4859 character(len=3),
save :: grid_id_saved=
""
4860 integer(i8_kind),
dimension(MAX_FIELDS),
save :: d_addrs = -9999_i8_kind
4861 integer(i8_kind),
dimension(MAX_FIELDS),
save :: x_addrs = -9999_i8_kind
4864 if (grid_id==xmap%grids(1)%id)
then
4865 is_complete = .true.
4866 if(
present(complete)) is_complete=complete
4868 if( lsize > max_fields )
then
4869 write( text,
'(i2)' ) max_fields
4870 call error_mesg (
'xgrid_mod',
'MAX_FIELDS='//trim(text)//
' exceeded for group get_side1_from_xgrid_ug', fatal)
4872 d_addrs(lsize) = loc(d)
4873 x_addrs(lsize) = loc(x)
4878 grid_id_saved = grid_id
4880 set_mismatch = .false.
4881 set_mismatch = set_mismatch .OR. (isize /=
size(d(:)))
4882 set_mismatch = set_mismatch .OR. (xsize /=
size(x(:)))
4883 set_mismatch = set_mismatch .OR. (grid_id_saved /= grid_id)
4884 if(set_mismatch)
then
4885 write( text,
'(i2)' ) lsize
4886 call error_mesg (
'xgrid_mod',
'Incompatible field at count '//text// &
4887 &
' for group get_side1_from_xgrid_ug', fatal )
4891 if(is_complete)
then
4893 call get_1_from_xgrid_ug_repro(d_addrs, x_addrs, xmap, xsize, lsize)
4895 call get_1_from_xgrid_ug(d_addrs, x_addrs, xmap, isize, xsize, lsize)
4897 d_addrs(1:lsize) = -9999_i8_kind
4898 x_addrs(1:lsize) = -9999_i8_kind
4907 do g=2,
size(xmap%grids(:))
4908 if (grid_id==xmap%grids(g)%id) &
4909 call error_mesg (
'xgrid_mod', &
4910 'get_from_xgrid_ug expects a 3D side 2 grid', fatal)
4913 call error_mesg (
'xgrid_mod',
'get_from_xgrid_ug: could not find grid id', fatal)
4915 end subroutine get_side1_from_xgrid_ug
4929 real(r8_kind),
dimension(:),
intent(in) :: d
4930 character(len=3),
intent(in) :: grid_id
4931 real(r8_kind),
dimension(:),
intent(inout) :: x
4932 type (xmap_type),
intent(inout) :: xmap
4933 logical,
intent(in),
optional :: complete
4935 logical :: is_complete, set_mismatch
4937 character(len=2) :: text
4938 integer,
save :: dsize=0
4939 integer,
save :: lsize=0
4940 integer,
save :: xsize=0
4941 character(len=3),
save :: grid_id_saved=
""
4942 integer(i8_kind),
dimension(MAX_FIELDS),
save :: d_addrs = -9999_i8_kind
4943 integer(i8_kind),
dimension(MAX_FIELDS),
save :: x_addrs = -9999_i8_kind
4945 if (grid_id==xmap%grids(1)%id)
then
4946 is_complete = .true.
4947 if(
present(complete)) is_complete=complete
4949 if( lsize > max_fields )
then
4950 write( text,
'(i2)' ) max_fields
4951 call error_mesg (
'xgrid_mod',
'MAX_FIELDS='//trim(text)//
' exceeded for group put_side1_to_xgrid_ug', fatal)
4953 d_addrs(lsize) = loc(d)
4954 x_addrs(lsize) = loc(x)
4959 grid_id_saved = grid_id
4961 set_mismatch = .false.
4962 set_mismatch = set_mismatch .OR. (dsize /=
size(d(:)))
4963 set_mismatch = set_mismatch .OR. (xsize /=
size(x(:)))
4964 set_mismatch = set_mismatch .OR. (grid_id_saved /= grid_id)
4965 if(set_mismatch)
then
4966 write( text,
'(i2)' ) lsize
4967 call error_mesg (
'xgrid_mod',
'Incompatible field at count '//text// &
4968 &
' for group put_side1_to_xgrid_ug', fatal )
4972 if(is_complete)
then
4973 call put_1_to_xgrid_ug_order_1(d_addrs, x_addrs, xmap, dsize, xsize, lsize)
4974 d_addrs(1:lsize) = -9999_i8_kind
4975 x_addrs(1:lsize) = -9999_i8_kind
4984 do g=2,
size(xmap%grids(:))
4985 if (grid_id==xmap%grids(g)%id) &
4986 call error_mesg (
'xgrid_mod', &
4987 'put_to_xgrid_ug expects a 2D side 2 grid', fatal)
4990 call error_mesg (
'xgrid_mod',
'put_to_xgrid_ug: could not find grid id', fatal)
5003 subroutine put_side2_to_xgrid_ug(d, grid_id, x, xmap)
5004 real(r8_kind),
dimension(:,:),
intent(in) :: d
5005 character(len=3),
intent(in) :: grid_id
5006 real(r8_kind),
dimension(:),
intent(inout) :: x
5007 type (xmap_type),
intent(inout) :: xmap
5011 if (grid_id==xmap%grids(1)%id) &
5012 call error_mesg (
'xgrid_mod', &
5013 'put_to_xgrid_ug expects a 2D side 1 grid', fatal)
5015 do g=2,
size(xmap%grids(:))
5016 if (grid_id==xmap%grids(g)%id)
then
5017 call put_2_to_xgrid_ug(d, xmap%grids(g), x, xmap)
5022 call error_mesg (
'xgrid_mod',
'put_to_xgrid_ug: could not find grid id', fatal)
5024 end subroutine put_side2_to_xgrid_ug
5035 subroutine get_side2_from_xgrid_ug(d, grid_id, x, xmap)
5036 real(r8_kind),
dimension(:,:),
intent(out) :: d
5037 character(len=3),
intent(in) :: grid_id
5038 real(r8_kind),
dimension(:),
intent(in) :: x
5039 type (xmap_type),
intent(in) :: xmap
5043 if (grid_id==xmap%grids(1)%id) &
5044 call error_mesg (
'xgrid_mod', &
5045 'get_from_xgrid_ug expects a 2D side 1 grid', fatal)
5047 do g=2,
size(xmap%grids(:))
5048 if (grid_id==xmap%grids(g)%id)
then
5049 call get_2_from_xgrid_ug(d, xmap%grids(g), x, xmap)
5054 call error_mesg (
'xgrid_mod',
'get_from_xgrid_ug: could not find grid id', fatal)
5056 end subroutine get_side2_from_xgrid_ug
5062 subroutine put_1_to_xgrid_ug_order_1(d_addrs, x_addrs, xmap, dsize, xsize, lsize)
5063 integer(i8_kind),
dimension(:),
intent(in) :: d_addrs
5064 integer(i8_kind),
dimension(:),
intent(in) :: x_addrs
5065 type (xmap_type),
intent(inout) :: xmap
5066 integer,
intent(in) :: dsize, xsize, lsize
5068 integer :: i, p, buffer_pos, msgsize
5069 integer :: from_pe, to_pe, pos, n, l, count
5070 integer :: ibegin, istart, iend, start_pos
5071 type (comm_type),
pointer,
save :: comm =>null()
5072 real(r8_kind) :: recv_buffer(xmap%put1%recvsize*lsize)
5073 real(r8_kind) :: send_buffer(xmap%put1%sendsize*lsize)
5074 real(r8_kind) :: unpack_buffer(xmap%put1%recvsize)
5076 real(r8_kind),
dimension(dsize) :: d
5077 real(r8_kind),
dimension(xsize) :: x
5082 call mpp_clock_begin(id_put_1_to_xgrid_order_1)
5086 do p = 1, comm%nrecv
5087 msgsize = comm%recv(p)%count*lsize
5088 from_pe = comm%recv(p)%pe
5089 buffer_pos = comm%recv(p)%buffer_pos*lsize
5090 call mpp_recv(recv_buffer(buffer_pos+1), glen=msgsize, from_pe = from_pe, block=.false., tag=comm_tag_7)
5095 do p = 1, comm%nsend
5096 msgsize = comm%send(p)%count*lsize
5097 to_pe = comm%send(p)%pe
5101 do n = 1, comm%send(p)%count
5103 lll = comm%send(p)%i(n)
5104 send_buffer(pos) = d(lll)
5107 call mpp_send(send_buffer(buffer_pos+1), plen=msgsize, to_pe = to_pe, tag=comm_tag_7 )
5108 buffer_pos = buffer_pos + msgsize
5114 if( lsize == 1)
then
5116 do l=1,xmap%size_put1
5117 x(l) = recv_buffer(xmap%x1_put(l)%pos)
5125 do p = 1, comm%nrecv
5126 count = comm%recv(p)%count
5127 ibegin = comm%recv(p)%buffer_pos*lsize + 1
5128 istart = ibegin + (l-1)*count
5129 iend = istart + count - 1
5130 pos = comm%recv(p)%buffer_pos
5133 unpack_buffer(pos) = recv_buffer(n)
5136 do i=1,xmap%size_put1
5137 x(i) = unpack_buffer(xmap%x1_put(i)%pos)
5144 call mpp_clock_end(id_put_1_to_xgrid_order_1)
5146 end subroutine put_1_to_xgrid_ug_order_1
5150 subroutine put_2_to_xgrid_ug(d, grid, x, xmap)
5151 type (grid_type),
intent(in) :: grid
5152 real(r8_kind),
dimension(grid%ls_me:grid%le_me, grid%km),
intent(in) :: d
5153 real(r8_kind),
dimension(:),
intent(inout) :: x
5154 type (xmap_type),
intent(in) :: xmap
5157 call mpp_clock_begin(id_put_2_to_xgrid)
5159 do l=grid%first,grid%last
5160 x(l) = d(xmap%x2(l)%l,xmap%x2(l)%k)
5163 call mpp_clock_end(id_put_2_to_xgrid)
5164 end subroutine put_2_to_xgrid_ug
5167 subroutine get_1_from_xgrid_ug(d_addrs, x_addrs, xmap, isize, xsize, lsize)
5168 integer(i8_kind),
dimension(:),
intent(in) :: d_addrs
5169 integer(i8_kind),
dimension(:),
intent(in) :: x_addrs
5170 type (xmap_type),
intent(inout) :: xmap
5171 integer,
intent(in) :: isize, xsize, lsize
5173 real(r8_kind),
dimension(xmap%size),
target :: dg(xmap%size, lsize)
5174 integer :: i, j, l, p, n, m
5175 integer :: msgsize, buffer_pos, pos
5176 integer :: istart, iend, count
5177 real(r8_kind) ,
pointer,
save :: dgp =>null()
5178 type (grid_type) ,
pointer,
save :: grid1 =>null()
5179 type (comm_type) ,
pointer,
save :: comm =>null()
5182 real(r8_kind) :: recv_buffer(xmap%get1%recvsize*lsize*3)
5183 real(r8_kind) :: send_buffer(xmap%get1%sendsize*lsize*3)
5184 real(r8_kind) :: d(isize)
5185 real(r8_kind),
dimension(xsize) :: x
5189 call mpp_clock_begin(id_get_1_from_xgrid)
5192 grid1 => xmap%grids(1)
5194 do p = 1, comm%nrecv
5195 recv => comm%recv(p)
5196 msgsize = recv%count*lsize
5197 buffer_pos = recv%buffer_pos*lsize
5198 call mpp_recv(recv_buffer(buffer_pos+1), glen=msgsize, from_pe = recv%pe, block=.false., tag=comm_tag_9)
5206 dgp => dg(xmap%x1(i)%pos,l)
5207 dgp = dgp + xmap%x1(i)%area*x(i)
5215 do p = 1, comm%nsend
5216 send => comm%send(p)
5217 msgsize = send%count*lsize
5219 istart = send%buffer_pos+1
5220 iend = istart + send%count - 1
5224 send_buffer(pos) = dg(n,l)
5227 call mpp_send(send_buffer(buffer_pos+1), plen=msgsize, to_pe = send%pe, tag=comm_tag_9 )
5228 buffer_pos = buffer_pos + msgsize
5241 do p = 1, comm%nrecv
5242 recv => comm%recv(p)
5244 buffer_pos = recv%buffer_pos*lsize
5245 if( recv%pe == xmap%me )
then
5249 pos = buffer_pos + (l-1)*count
5254 d(i) = recv_buffer(pos)
5262 do m = 1, comm%nrecv
5263 p = comm%unpack_ind(m)
5264 recv => comm%recv(p)
5265 if( recv%pe == xmap%me )
then
5268 buffer_pos = recv%buffer_pos*lsize
5272 pos = buffer_pos + (l-1)*recv%count
5274 do n = 1, recv%count
5277 d(i) = d(i) + recv_buffer(pos)
5288 d = d * grid1%area_inv(:,1)
5291 call mpp_clock_end(id_get_1_from_xgrid)
5293 end subroutine get_1_from_xgrid_ug
5297 subroutine get_1_from_xgrid_ug_repro(d_addrs, x_addrs, xmap, xsize, lsize)
5298 integer(i8_kind),
dimension(:),
intent(in) :: d_addrs
5299 integer(i8_kind),
dimension(:),
intent(in) :: x_addrs
5300 type (xmap_type),
intent(inout) :: xmap
5301 integer,
intent(in) :: xsize, lsize
5303 integer :: g, i, j, k, p, l, n, l2, l3
5304 integer :: msgsize, buffer_pos, pos
5305 type (grid_type),
pointer,
save :: grid =>null()
5306 type(
comm_type),
pointer,
save :: comm => null()
5309 integer,
dimension(0:xmap%npes-1) :: pl, ml
5310 real(r8_kind) :: recv_buffer(xmap%recv_count_repro_tot*lsize)
5311 real(r8_kind) :: send_buffer(xmap%send_count_repro_tot*lsize)
5312 real(r8_kind) :: d(xmap%grids(1)%ls_me:xmap%grids(1)%le_me)
5313 real(r8_kind),
dimension(xsize) :: x
5317 call mpp_clock_begin(id_get_1_from_xgrid_repro)
5318 comm => xmap%get1_repro
5320 do p = 1, comm%nrecv
5321 recv => comm%recv(p)
5322 msgsize = recv%count*lsize
5323 buffer_pos = recv%buffer_pos*lsize
5324 call mpp_recv(recv_buffer(buffer_pos+1), glen=msgsize, from_pe = recv%pe, block=.false., tag=comm_tag_10)
5325 n = recv%pe -xmap%root_pe
5331 send_buffer(:) = 0.0_r8_kind
5334 do p = 1, comm%nsend
5335 pos = comm%send(p)%buffer_pos*lsize
5336 send => comm%send(p)
5339 do n = 1, send%count
5345 do k =1, xmap%grids(g)%km
5346 if(xmap%grids(g)%frac_area(i,j,k)/=0.0_r8_kind)
then
5348 send_buffer(pos) = send_buffer(pos) + xmap%x1(l2)%area *x(l2)
5356 buffer_pos = comm%send(p)%buffer_pos*lsize
5357 msgsize = comm%send(p)%count*lsize
5358 call mpp_send(send_buffer(buffer_pos+1), plen=msgsize, to_pe=comm%send(p)%pe, tag=comm_tag_10)
5372 do g=2,
size(xmap%grids(:))
5373 grid => xmap%grids(g)
5374 do l3=1,grid%size_repro
5375 i = grid%x_repro(l3)%l1
5376 p = grid%x_repro(l3)%pe-xmap%root_pe
5377 pos = pl(p) + (l-1)*ml(p) + grid%x_repro(l3)%recv_pos
5378 d(i) = d(i) + recv_buffer(pos)
5382 d = d * xmap%grids(1)%area_inv(:,1)
5387 call mpp_clock_end(id_get_1_from_xgrid_repro)
5389 end subroutine get_1_from_xgrid_ug_repro
5393 subroutine get_2_from_xgrid_ug(d, grid, x, xmap)
5394 type (grid_type),
intent(in) :: grid
5395 real(r8_kind),
dimension(grid%ls_me:grid%le_me, grid%km),
intent(out) :: d
5396 real(r8_kind),
dimension(:),
intent(in) :: x
5397 type (xmap_type),
intent(in) :: xmap
5401 call mpp_clock_begin(id_get_2_from_xgrid)
5404 do l=grid%first_get,grid%last_get
5405 d(xmap%x2_get(l)%l,xmap%x2_get(l)%k) = &
5406 d(xmap%x2_get(l)%l,xmap%x2_get(l)%k) + xmap%x2_get(l)%area*x(xmap%x2_get(l)%pos)
5412 d(:,k) = d(:,k) * grid%area_inv(:,1)
5415 call mpp_clock_end(id_get_2_from_xgrid)
5417 end subroutine get_2_from_xgrid_ug
5422 integer,
intent(in) :: i, j
5427 g = (j-1)*grid%ni + i
5428 in_box_me = (g>=grid%gs_me) .and. (g<=grid%ge_me)
5430 in_box_me = (i>=grid%is_me) .and. (i<=grid%ie_me) .and. (j>=grid%js_me) .and. (j<=grid%je_me)
5438 integer,
intent(in) :: i, j, p
5443 g = (j-1)*grid%ni + i
5444 in_box_nbr = (g>=grid%gs(p)) .and. (g<=grid%ge(p))
5446 in_box_nbr = (i>=grid%is(p)) .and. (i<=grid%ie(p)) .and. (j>=grid%js(p)) .and. (j<=grid%je(p))
5451 end module xgrid_mod
5454 Register a diagnostic field for a given module.
Send data over to output fields.
Close a netcdf or domain file opened with open_file or open_virtual_file.
Opens a given netcdf or domain file.
Read data from a defined field in a file.
integer function, public check_nml_error(IOSTAT, NML_NAME)
Checks the iostat argument that is returned after reading a namelist and determines if the error code...
subroutine, public write_version_number(version, tag, unit)
Prints to the log file (or a specified unit) the version id string and tag name.
subroutine, public error_mesg(routine, message, level)
Print notes, warnings and error messages; terminates program for warning and error messages....
integer function, public get_mosaic_xgrid_size(fileobj)
return exchange grid size of mosaic xgrid file.
subroutine, public get_mosaic_contact(fileobj, tile1, tile2, istart1, iend1, jstart1, jend1, istart2, iend2, jstart2, jend2)
Get contact information from mosaic_file Example usage: call get_mosaic_contact(mosaic_file,...
integer function, public get_mosaic_ntiles(fileobj)
Get number of tiles in the mosaic_file.
subroutine, public get_mosaic_tile_grid(grid_file, fileobj, domain, tile_count)
Gets the name of a mosaic tile grid file.
integer function, public get_mosaic_ncontacts(fileobj)
Get number of contacts in the mosaic_file.
subroutine, public get_mosaic_grid_sizes(fileobj, nx, ny)
Get grid size of each tile from mosaic_file.
integer function mpp_get_domain_npes(domain)
Set user stack size.
integer function, dimension(size(domain%tile_id(:))) mpp_get_tile_id(domain)
Returns the tile_id on current pe.
integer function mpp_get_current_ntile(domain)
Returns number of tile on current pe.
integer function mpp_get_ntile_count(domain)
Returns number of tiles in mosaic.
subroutine mpp_compute_extent(isg, ieg, ndivs, ibegin, iend, extent)
Computes extents for a grid decomposition with the given indices and divisions.
subroutine mpp_get_tile_list(domain, tiles)
Return the tile_id on current pelist. one-tile-per-pe is assumed.
logical function mpp_domain_is_initialized(domain)
Set user stack size.
subroutine mpp_get_domain_pelist(domain, pelist)
Set user stack size.
integer function mpp_get_domain_root_pe(domain)
Set user stack size.
Broadcasts domain to every pe. Only useful outside the context of it's own pelist.
Deallocate given 1D or 2D domain.
Set up a domain decomposition.
These routines retrieve the axis specifications associated with the compute domains....
Retrieve the entire array of compute domain extents associated with a decomposition.
These routines retrieve the axis specifications associated with the data domains. The domain is a der...
These routines retrieve the axis specifications associated with the global domains....
Global sum of domain-decomposed arrays. mpp_global_sum is used to get the sum of a domain-decomposed...
Modifies the extents (compute, data and global) of a given domain.
Passes data from a structured grid to an unstructured grid Example usage:
Reorganization of distributed global arrays. mpp_redistribute is used to reorganize a distributed ar...
Performs halo updates for a given domain.
One dimensional domain used to manage shared data access between pes.
The domain2D type contains all the necessary information to define the global, compute and data domai...
Domain information for managing data on unstructured grids.
subroutine mpp_sync_self(pelist, check, request, msg_size, msg_type)
This is to check if current PE's outstanding puts are complete but we can't use shmem_fence because w...
integer function stdout()
This function returns the current standard fortran unit numbers for output.
subroutine mpp_set_current_pelist(pelist, no_sync)
Set context pelist.
integer function stdlog()
This function returns the current standard fortran unit numbers for log messages. Log messages,...
integer function mpp_npes()
Returns processor count for current pelist.
integer function mpp_pe()
Returns processor ID.
subroutine mpp_sync(pelist, do_self)
Synchronize PEs in list.
integer function mpp_clock_id(name, flags, grain)
Return an ID for a new or existing clock.
Scatter a vector across all PEs.
Reduction operations. Find the max of scalar a from the PEs in pelist result is also automatically br...
Reduction operations. Find the min of scalar a from the PEs in pelist result is also automatically br...
Recieve data from another PE.
Send data to a receiving PE.
Holds stocks amounts per PE values.
subroutine, public get_time(Time, seconds, days, ticks, err_msg)
Returns days and seconds ( < 86400 ) corresponding to a time. err_msg should be checked for any error...
Type to represent amounts of time. Implemented as seconds and days to allow for larger intervals.
subroutine get_area_elements_fms2_io(fileobj, name, get_area_data)
Read the area elements from NetCDF file.
integer function, public xgrid_count(xmap)
Returns current size of exchange grid variables.
subroutine, public some(xmap, some_arr, grid_id)
Returns logical associating exchange grid cells with given side two grid.
subroutine regen(xmap)
Regenerate/Update the xmap.
integer nsubset
Number of processors to read exchange grid information. Those processors that read the exchange grid ...
subroutine put_side1_to_xgrid(d, grid_id, x, xmap, remap_method, complete)
Scatters data to exchange grid.
integer function get_nest_contact_fms2_io(fileobj, tile_nest_out, tile_parent_out, is_nest_out, ie_nest_out, js_nest_out, je_nest_out, is_parent_out, ie_parent_out, js_parent_out, je_parent_out)
currently we are assuming there is only one nest region
subroutine, public get_ocean_model_area_elements(domain, grid_file)
Read Ocean area element data from netCDF file.
subroutine, public set_frac_area_ug(f, grid_id, xmap)
Changes sub-grid portion areas and/or number.
real(r8_kind) function, dimension(3) conservation_check_ug_side1(d, grid_id, xmap, remap_method)
conservation_check_ug - returns three numbers which are the global sum of a variable (1) on its home ...
subroutine get_side1_from_xgrid(d, grid_id, x, xmap, complete)
subroutine stock_move_3d(from, to, grid_index, stock_data3d, xmap, delta_t, from_side, to_side, radius, verbose, ier)
this version takes rank 3 data, it can be used to compute the flux on anything but the first grid,...
real(r8_kind), dimension(:,:), allocatable, public area_atm_sphere
Area elements based on a the spherical model used by the ICE layer.
logical function in_box_me(i, j, grid)
logical function in_box(i, j, is, ie, js, je)
logical use_mpp_io
use_mpp_io Default = .false. When true, uses mpp_io for IO. When false, uses fms2_io for IO.
integer remapping_method
xgrid nml
subroutine get_grid_version2(grid, grid_id, grid_file)
read the center point of the grid from version 1 grid file. only the grid at the side 1 is needed,...
subroutine get_side2_from_xgrid(d, grid_id, x, xmap)
real(r8_kind) function, dimension(3) conservation_check_ug_side2(d, grid_id, xmap, remap_method)
conservation_check_ug - returns three numbers which are the global sum of a variable (1) on its home ...
real(r8_kind) function, dimension(is:ie, js:je) grad_merid_latlon(d, lat, is, ie, js, je, isd, jsd)
This function is used to calculate the gradient along meridinal direction. Maybe need to setup a limi...
logical make_exchange_reproduce
Set to .true. to make xgrid_mod reproduce answers on different numbers of PEs. This option has a cons...
character(len=64) interp_method
Exchange grid interpolation method. It has two options: "first_order", "second_order".
integer, parameter version2
mosaic grid file
subroutine, public setup_xmap(xmap, grid_ids, grid_domains, grid_file, atm_grid, lnd_ug_domain)
Sets up exchange grid connectivity using grid specification file and processor domain decomposition.
subroutine get_grid_version1(grid, grid_id, grid_file)
read the center point of the grid from version 1 grid file. only the grid at the side 1 is needed,...
subroutine put_side1_to_xgrid_ug(d, grid_id, x, xmap, complete)
Currently only support first order.
subroutine put_side2_to_xgrid(d, grid_id, x, xmap)
Scatters data to exchange grid.
subroutine set_frac_area_sg(f, grid_id, xmap)
Changes sub-grid portion areas and/or number.
logical function in_box_nbr(i, j, grid, p)
integer, parameter version1
grid spec file
real(r8_kind) function, dimension(3) conservation_check_side2(d, grid_id, xmap, remap_method)
conservation_check - returns three numbers which are the global sum of a variable (1) on its home mod...
real(r8_kind) function, dimension(3) conservation_check_side1(d, grid_id, xmap, remap_method)
conservation_check - returns three numbers which are the global sum of a variable (1) on its home mod...
subroutine, public xgrid_init(remap_method)
Initialize the xgrid_mod.
real(r8_kind), dimension(:,:), allocatable, public area_atm_model
Area elements used inside each model.
real(r8_kind) function, dimension(is:ie, js:je) grad_zonal_latlon(d, lon, lat, is, ie, js, je, isd, jsd)
This function is used to calculate the gradient along zonal direction. Maybe need to setup a limit fo...
subroutine, public get_xmap_grid_area(id, xmap, area)
This routine is used to get the grid area of component model with id.
Returns three numbers which are the global sum of a variable.
For an unstructured grid, returns three numbers which are the global sum of a variable (1) on its hom...
Sums data from exchange grid to model grid.
get_from_xgrid for unstructured grids.
Scatters data from model grid onto exchange grid.
put_to_xgrid for unstructured grids.
Sets sub-grid area and numbering in the given exchange grid.
Private type used for exchange grid communication.
Type to hold pointers for grid boxes.
Private type to hold all data needed from given grid for an exchange grid.
Private type for overlap exchange grid data.
Private type for exchange grid data.
Private type for exchange grid data.
Private type for cell indices and data in the exchange grid.
Type for an exchange grid, holds pointers to included grids and any necessary data.