@@ -16,7 +16,7 @@ module MOM_cap_mod
16
16
use MOM_domains, only: MOM_infra_init, MOM_infra_end
17
17
use MOM_file_parser, only: get_param, log_version, param_file_type, close_param_file
18
18
use MOM_get_input, only: get_MOM_input, directories
19
- use MOM_domains, only: pass_var
19
+ use MOM_domains, only: pass_var, pe_here
20
20
use MOM_error_handler, only: MOM_error, FATAL, is_root_pe
21
21
use MOM_grid, only: ocean_grid_type, get_global_grid_size
22
22
use MOM_ocean_model_nuopc, only: ice_ocean_boundary_type
@@ -29,6 +29,7 @@ module MOM_cap_mod
29
29
use MOM_cap_methods, only: med2mod_areacor, state_diagnose
30
30
use MOM_cap_methods, only: ChkErr
31
31
use MOM_ensemble_manager, only: ensemble_manager_init
32
+ use MOM_coms, only: sum_across_PEs
32
33
33
34
#ifdef CESMCOUPLED
34
35
use shr_log_mod, only: shr_log_setLogUnit
@@ -842,6 +843,7 @@ subroutine InitializeRealize(gcomp, importState, exportState, clock, rc)
842
843
type (ocean_grid_type) , pointer :: ocean_grid
843
844
type (ocean_internalstate_wrapper) :: ocean_internalstate
844
845
integer :: npet, ntiles
846
+ integer :: npes ! number of PEs (from FMS).
845
847
integer :: nxg, nyg, cnt
846
848
integer :: isc,iec,jsc,jec
847
849
integer , allocatable :: xb(:),xe(:),yb(:),ye(:),pe(:)
@@ -868,6 +870,8 @@ subroutine InitializeRealize(gcomp, importState, exportState, clock, rc)
868
870
integer :: lsize
869
871
integer :: ig,jg, ni,nj,k
870
872
integer , allocatable :: gindex(:) ! global index space
873
+ integer , allocatable :: gindex_ocn(:) ! global index space for ocean cells (excl. masked cells)
874
+ integer , allocatable :: gindex_elim(:) ! global index space for eliminated cells
871
875
character (len= 128 ) :: fldname
872
876
character (len= 256 ) :: cvalue
873
877
character (len= 256 ) :: frmt ! format specifier for several error msgs
@@ -891,6 +895,11 @@ subroutine InitializeRealize(gcomp, importState, exportState, clock, rc)
891
895
real (ESMF_KIND_R8 ) :: min_areacor_glob(2 )
892
896
real (ESMF_KIND_R8 ) :: max_areacor_glob(2 )
893
897
character (len=* ), parameter :: subname= ' (MOM_cap:InitializeRealize)'
898
+ integer :: niproc, njproc
899
+ integer :: ip, jp, pe_ix
900
+ integer :: num_elim_blocks ! number of blocks to be eliminated
901
+ integer :: num_elim_cells_global, num_elim_cells_local, num_elim_cells_remaining
902
+ integer , allocatable :: cell_mask(:,:)
894
903
real (8 ) :: MPI_Wtime, timeirls
895
904
!- -------------------------------
896
905
@@ -937,19 +946,19 @@ subroutine InitializeRealize(gcomp, importState, exportState, clock, rc)
937
946
rc = ESMF_FAILURE
938
947
call ESMF_LogWrite(subname// ' ntiles must be 1' , ESMF_LOGMSG_ERROR)
939
948
endif
940
- ntiles = mpp_get_domain_npes(ocean_public% domain)
941
- write (tmpstr,' (a,1i6)' ) subname// ' ntiles = ' ,ntiles
949
+ npes = mpp_get_domain_npes(ocean_public% domain)
950
+ write (tmpstr,' (a,1i6)' ) subname// ' npes = ' ,npes
942
951
call ESMF_LogWrite(trim (tmpstr), ESMF_LOGMSG_INFO)
943
952
944
953
!- --------------------------------
945
954
! get start and end indices of each tile and their PET
946
955
!- --------------------------------
947
956
948
- allocate (xb(ntiles ),xe(ntiles ),yb(ntiles ),ye(ntiles ),pe(ntiles ))
957
+ allocate (xb(npes ),xe(npes ),yb(npes ),ye(npes ),pe(npes ))
949
958
call mpp_get_compute_domains(ocean_public% domain, xbegin= xb, xend= xe, ybegin= yb, yend= ye)
950
959
call mpp_get_pelist(ocean_public% domain, pe)
951
960
if (dbug > 1 ) then
952
- do n = 1 ,ntiles
961
+ do n = 1 ,npes
953
962
write (tmpstr,' (a,6i6)' ) subname// ' tiles ' ,n,pe(n),xb(n),xe(n),yb(n),ye(n)
954
963
call ESMF_LogWrite(trim (tmpstr), ESMF_LOGMSG_INFO)
955
964
enddo
@@ -971,17 +980,102 @@ subroutine InitializeRealize(gcomp, importState, exportState, clock, rc)
971
980
call get_global_grid_size(ocean_grid, ni, nj)
972
981
lsize = ( ocean_grid% iec - ocean_grid% isc + 1 ) * ( ocean_grid% jec - ocean_grid% jsc + 1 )
973
982
974
- ! Create the global index space for the computational domain
975
- allocate (gindex(lsize))
976
- k = 0
977
- do j = ocean_grid% jsc, ocean_grid% jec
978
- jg = j + ocean_grid% jdg_offset
979
- do i = ocean_grid% isc, ocean_grid% iec
980
- ig = i + ocean_grid% idg_offset
981
- k = k + 1 ! Increment position within gindex
982
- gindex(k) = ni * (jg - 1 ) + ig
983
+ num_elim_blocks = 0
984
+ num_elim_cells_global = 0
985
+ num_elim_cells_local = 0
986
+ num_elim_cells_remaining = 0
987
+
988
+ ! Compute the number of eliminated blocks (specified in MOM_mask_table)
989
+ if (associated (ocean_grid% Domain% maskmap)) then
990
+ njproc = size (ocean_grid% Domain% maskmap, 1 )
991
+ niproc = size (ocean_grid% Domain% maskmap, 2 )
992
+
993
+ do ip = 1 , niproc
994
+ do jp = 1 , njproc
995
+ if (.not. ocean_grid% Domain% maskmap(jp,ip)) then
996
+ num_elim_blocks = num_elim_blocks+1
997
+ endif
998
+ enddo
983
999
enddo
984
- enddo
1000
+ endif
1001
+
1002
+ ! Apply land block elimination to ESMF gindex
1003
+ ! (Here we assume that each processor gets assigned a single tile. If multi-tile implementation is to be added
1004
+ ! in MOM6 NUOPC cap in the future, below code must be updated accordingly.)
1005
+ if (num_elim_blocks> 0 ) then
1006
+
1007
+ allocate (cell_mask(ni, nj), source= 0 )
1008
+ allocate (gindex_ocn(lsize))
1009
+ k = 0
1010
+ do j = ocean_grid% jsc, ocean_grid% jec
1011
+ jg = j + ocean_grid% jdg_offset
1012
+ do i = ocean_grid% isc, ocean_grid% iec
1013
+ ig = i + ocean_grid% idg_offset
1014
+ k = k + 1 ! Increment position within gindex
1015
+ gindex_ocn(k) = ni * (jg - 1 ) + ig
1016
+ cell_mask(ig, jg) = 1
1017
+ enddo
1018
+ enddo
1019
+ call sum_across_PEs(cell_mask, ni* nj)
1020
+
1021
+ if (maxval (cell_mask) /= 1 ) then
1022
+ call MOM_error(FATAL, " Encountered cells shared by multiple PEs while attempting to determine masked cells." )
1023
+ endif
1024
+
1025
+ num_elim_cells_global = ni * nj - sum (cell_mask)
1026
+ num_elim_cells_local = num_elim_cells_global / npes
1027
+
1028
+ if (pe_here() == pe(npes)) then
1029
+ ! assign all remaining cells to the last PE.
1030
+ num_elim_cells_remaining = num_elim_cells_global - num_elim_cells_local * npes
1031
+ allocate (gindex_elim(num_elim_cells_local+ num_elim_cells_remaining))
1032
+ else
1033
+ allocate (gindex_elim(num_elim_cells_local))
1034
+ endif
1035
+
1036
+ ! Zero-based PE index.
1037
+ pe_ix = pe_here() - pe(1 )
1038
+
1039
+ k = 0
1040
+ do jg = 1 , nj
1041
+ do ig = 1 , ni
1042
+ if (cell_mask(ig, jg) == 0 ) then
1043
+ k = k + 1
1044
+ if (k > pe_ix * num_elim_cells_local .and. &
1045
+ k <= ((pe_ix+1 ) * num_elim_cells_local + num_elim_cells_remaining)) then
1046
+ gindex_elim(k - pe_ix * num_elim_cells_local) = ni * (jg - 1 ) + ig
1047
+ endif
1048
+ endif
1049
+ enddo
1050
+ enddo
1051
+
1052
+ allocate (gindex(lsize + num_elim_cells_local + num_elim_cells_remaining))
1053
+ do k = 1 , lsize
1054
+ gindex(k) = gindex_ocn(k)
1055
+ enddo
1056
+ do k = 1 , num_elim_cells_local + num_elim_cells_remaining
1057
+ gindex(k+ lsize) = gindex_elim(k)
1058
+ enddo
1059
+
1060
+ deallocate (cell_mask)
1061
+ deallocate (gindex_ocn)
1062
+ deallocate (gindex_elim)
1063
+
1064
+ else ! no eliminated land blocks
1065
+
1066
+ ! Create the global index space for the computational domain
1067
+ allocate (gindex(lsize))
1068
+ k = 0
1069
+ do j = ocean_grid% jsc, ocean_grid% jec
1070
+ jg = j + ocean_grid% jdg_offset
1071
+ do i = ocean_grid% isc, ocean_grid% iec
1072
+ ig = i + ocean_grid% idg_offset
1073
+ k = k + 1 ! Increment position within gindex
1074
+ gindex(k) = ni * (jg - 1 ) + ig
1075
+ enddo
1076
+ enddo
1077
+
1078
+ endif
985
1079
986
1080
DistGrid = ESMF_DistGridCreate(arbSeqIndexList= gindex, rc= rc)
987
1081
if (ChkErr(rc,__LINE__,u_FILE_u)) return
@@ -1005,6 +1099,10 @@ subroutine InitializeRealize(gcomp, importState, exportState, clock, rc)
1005
1099
call ESMF_MeshGet(Emesh, spatialDim= spatialDim, numOwnedElements= numOwnedElements, rc= rc)
1006
1100
if (ChkErr(rc,__LINE__,u_FILE_u)) return
1007
1101
1102
+ if (lsize /= numOwnedElements - num_elim_cells_local - num_elim_cells_remaining) then
1103
+ call MOM_error(FATAL, " Discrepancy detected between ESMF mesh and internal MOM6 domain sizes. Check mask table." )
1104
+ endif
1105
+
1008
1106
allocate (ownedElemCoords(spatialDim* numOwnedElements))
1009
1107
allocate (lonMesh(numOwnedElements), lon(numOwnedElements))
1010
1108
allocate (latMesh(numOwnedElements), lat(numOwnedElements))
@@ -1036,7 +1134,7 @@ subroutine InitializeRealize(gcomp, importState, exportState, clock, rc)
1036
1134
end do
1037
1135
1038
1136
eps_omesh = get_eps_omesh(ocean_state)
1039
- do n = 1 ,numOwnedElements
1137
+ do n = 1 ,lsize
1040
1138
diff_lon = abs (mod (lonMesh(n) - lon(n),360.0 ))
1041
1139
if (diff_lon > eps_omesh) then
1042
1140
frmt = " ('ERROR: Difference between ESMF Mesh and MOM6 domain coords is " // &
@@ -1140,11 +1238,11 @@ subroutine InitializeRealize(gcomp, importState, exportState, clock, rc)
1140
1238
1141
1239
! generate delayout and dist_grid
1142
1240
1143
- allocate (deBlockList(2 ,2 ,ntiles ))
1144
- allocate (petMap(ntiles ))
1145
- allocate (deLabelList(ntiles ))
1241
+ allocate (deBlockList(2 ,2 ,npes ))
1242
+ allocate (petMap(npes ))
1243
+ allocate (deLabelList(npes ))
1146
1244
1147
- do n = 1 , ntiles
1245
+ do n = 1 , npes
1148
1246
deLabelList(n) = n
1149
1247
deBlockList(1 ,1 ,n) = xb(n)
1150
1248
deBlockList(1 ,2 ,n) = xe(n)
@@ -1727,7 +1825,7 @@ subroutine ModelAdvance(gcomp, rc)
1727
1825
rpointer_filename = ' rpointer.ocn' // trim (inst_suffix)
1728
1826
1729
1827
write (restartname,' (A,".mom6.r.",I4.4,"-",I2.2,"-",I2.2,"-",I5.5)' ) &
1730
- trim (casename), year, month, day, seconds
1828
+ trim (casename), year, month, day, hour * 3600 + minute * 60 + seconds
1731
1829
call ESMF_LogWrite(" MOM_cap: Writing restart : " // trim (restartname), ESMF_LOGMSG_INFO)
1732
1830
! write restart file(s)
1733
1831
call ocean_model_restart(ocean_state, restartname= restartname, num_rest_files= num_rest_files)
0 commit comments