FMS  2025.04
Flexible Modeling System
get_grid_version.inc
1 !***********************************************************************
2 !* Apache License 2.0
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* Licensed under the Apache License, Version 2.0 (the "License");
7 !* you may not use this file except in compliance with the License.
8 !* You may obtain a copy of the License at
9 !*
10 !* http://www.apache.org/licenses/LICENSE-2.0
11 !*
12 !* FMS is distributed in the hope that it will be useful, but WITHOUT
13 !* WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied;
14 !* without even the implied warranty of MERCHANTABILITY or FITNESS FOR A
15 !* PARTICULAR PURPOSE. See the License for the specific language
16 !* governing permissions and limitations under the License.
17 !***********************************************************************
18 
19 !> Get global lon and lat of three model (target) grids, with a given file name
20 subroutine get_grid_version_1_(grid_file, mod_name, domain, lon, lat, min_lon, max_lon)
21  integer, parameter :: lkind = fms_get_grid_version_kind_
22 
23  character(len=*), intent(in) :: grid_file !< name of grid file
24  character(len=*), intent(in) :: mod_name !< module name
25  type(domain2d), intent(in) :: domain !< 2D domain
26  real(lkind), intent(out), dimension(:,:), allocatable :: lon, lat
27  real(lkind), intent(out) :: min_lon, max_lon
28 
29  integer :: i, j, siz(4)
30  integer :: nlon, nlat !< size of global lon and lat
31  real(lkind), dimension(:,:,:), allocatable :: lon_vert, lat_vert !< of OCN grid vertices
32  real(lkind), dimension(:), allocatable :: glon, glat !< lon and lat of 1-D grid of atm/lnd
33  logical :: is_new_grid
34  integer :: is, ie, js, je
35  integer :: isg, ieg, jsg, jeg
36  integer :: isc, iec, jsc, jec
37  character(len=3) :: xname, yname
38  integer :: start(2), nread(2)
39  type(FmsNetcdfDomainFile_t) :: fileobj
40  integer :: ndims !< Number of dimensions
41 
42  if(.not. open_file(fileobj, grid_file, 'read', domain )) then
43  call mpp_error(fatal, 'data_override_mod(get_grid_version_1): Error in opening file '//trim(grid_file))
44  endif
45 
46  call mpp_get_global_domain(domain, isg, ieg, jsg, jeg)
47  call mpp_get_compute_domain(domain, isc, iec, jsc, jec)
48 
49  allocate(lon(isc:iec, jsc:jec))
50  allocate(lat(isc:iec, jsc:jec))
51 
52  select case(mod_name)
53  case('ocn', 'ice')
54  is_new_grid = .false.
55  if(variable_exists(fileobj, 'x_T')) then
56  is_new_grid = .true.
57  else if(variable_exists(fileobj, 'geolon_t')) then
58  is_new_grid = .false.
59  else
60  call mpp_error(fatal,'data_override: both x_T and geolon_t is not in the grid file '//trim(grid_file) )
61  endif
62 
63  if(is_new_grid) then
64  ndims = get_variable_num_dimensions(fileobj, 'x_T')
65  call get_variable_size(fileobj, 'x_T', siz(1:ndims))
66  nlon = siz(1); nlat = siz(2)
67  call check_grid_sizes(trim(mod_name)//'_domain ', domain, nlon, nlat)
68  allocate(lon_vert(isc:iec,jsc:jec,4), lat_vert(isc:iec,jsc:jec,4) )
69 
70  call read_data(fileobj, 'x_vert_T', lon_vert)
71  call read_data(fileobj, 'y_vert_T', lat_vert)
72 
73 !2 Global lon and lat of ocean grid cell centers are determined from adjacent vertices
74  lon(:,:) = (lon_vert(:,:,1) + lon_vert(:,:,2) + lon_vert(:,:,3) + lon_vert(:,:,4))*0.25_lkind
75  lat(:,:) = (lat_vert(:,:,1) + lat_vert(:,:,2) + lat_vert(:,:,3) + lat_vert(:,:,4))*0.25_lkind
76  else
77 
78  ndims = get_variable_num_dimensions(fileobj, 'geolon_vert_t')
79  call get_variable_size(fileobj, 'geolon_vert_t', siz(1:ndims))
80  nlon = siz(1) - 1; nlat = siz(2) - 1;
81  call check_grid_sizes(trim(mod_name)//'_domain ', domain, nlon, nlat)
82 
83  start(1) = isc; nread(1) = iec-isc+2
84  start(2) = jsc; nread(2) = jec-jsc+2
85 
86  allocate(lon_vert(isc:iec+1,jsc:jec+1,1))
87  allocate(lat_vert(isc:iec+1,jsc:jec+1,1))
88 
89  call read_data(fileobj, 'geolon_vert_t', lon_vert(:,:,1), corner=start, edge_lengths=nread)
90  call read_data(fileobj, 'geolat_vert_t', lat_vert(:,:,1), corner=start, edge_lengths=nread)
91 
92  do j = jsc, jec
93  do i = isc, iec
94  lon(i,j) = (lon_vert(i,j,1) + lon_vert(i+1,j,1) + &
95  lon_vert(i+1,j+1,1) + lon_vert(i,j+1,1))*0.25_lkind
96  lat(i,j) = (lat_vert(i,j,1) + lat_vert(i+1,j,1) + &
97  lat_vert(i+1,j+1,1) + lat_vert(i,j+1,1))*0.25_lkind
98  enddo
99  enddo
100  endif
101  deallocate(lon_vert)
102  deallocate(lat_vert)
103  case('atm', 'lnd')
104  if(trim(mod_name) == 'atm') then
105  xname = 'xta'; yname = 'yta'
106  else
107  xname = 'xtl'; yname = 'ytl'
108  endif
109  ndims = get_variable_num_dimensions(fileobj, xname)
110  call get_variable_size(fileobj, xname, siz(1:ndims))
111  nlon = siz(1); allocate(glon(nlon))
112  call read_data(fileobj, xname, glon)
113 
114  ndims = get_variable_num_dimensions(fileobj, xname)
115  call get_variable_size(fileobj, yname, siz(1:ndims))
116  nlat = siz(1); allocate(glat(nlat))
117  call read_data(fileobj, yname, glat)
118  call check_grid_sizes(trim(mod_name)//'_domain ', domain, nlon, nlat)
119 
120  is = isc - isg + 1; ie = iec - isg + 1
121  js = jsc - jsg + 1; je = jec - jsg + 1
122  do j = js, jec
123  do i = is, ie
124  lon(i,j) = glon(i)
125  lat(i,j) = glat(j)
126  enddo
127  enddo
128  deallocate(glon)
129  deallocate(glat)
130  case default
131  call mpp_error(fatal, "data_override_mod: mod_name should be 'atm', 'ocn', 'ice' or 'lnd' ")
132  end select
133 
134  call close_file(fileobj)
135 
136  ! convert from degree to radian
137  lon = lon * real(deg_to_rad, lkind)
138  lat = lat* real(deg_to_rad, lkind)
139  min_lon = minval(lon)
140  max_lon = maxval(lon)
141  call mpp_min(min_lon)
142  call mpp_max(max_lon)
143 end subroutine get_grid_version_1_
144 
145 !> Get global lon and lat of three model (target) grids from mosaic.nc.
146 !! Currently we assume the refinement ratio is 2 and there is one tile on each pe.
147 subroutine get_grid_version_2_(fileobj, mod_name, domain, lon, lat, min_lon, max_lon, use_center_grid_points)
148  integer, parameter :: lkind = fms_get_grid_version_kind_
149 
150  type(FmsNetcdfFile_t), intent(in) :: fileobj !< file object for grid file
151  character(len=*), intent(in) :: mod_name !< module name
152  type(domain2d), intent(in) :: domain !< 2D domain
153  real(lkind), intent(out), dimension(:,:), allocatable :: lon, lat
154  real(lkind), intent(out) :: min_lon, max_lon
155  logical, optional, intent(in) :: use_center_grid_points !< Flag indicating whether or not to use the
156  !! centroid values of the supergrid from the
157  !! grid file as opposed to calcuating it by
158  !! taking the average of the four corner points.
159  !! This is only relevant to OCN and ICE grids.
160 
161  integer :: i, j, siz(2)
162  integer :: nlon, nlat ! size of global grid
163  integer :: nlon_super, nlat_super ! size of global supergrid.
164  integer :: isc, iec, jsc, jec
165  integer :: isc2, iec2, jsc2, jec2
166  character(len=FMS_PATH_LEN) :: solo_mosaic_file, grid_file
167  real(lkind), allocatable :: tmpx(:,:), tmpy(:,:)
168  logical :: open_solo_mosaic
169  type(FmsNetcdfFile_t) :: mosaicfileobj, tilefileobj
170  integer :: start(2), nread(2)
171  logical :: use_center_grid_points_local
172 
173  use_center_grid_points_local = .false.
174  if (present(use_center_grid_points)) use_center_grid_points_local = use_center_grid_points
175 
176  if(trim(mod_name) .NE. 'atm' .AND. trim(mod_name) .NE. 'ocn' .AND. &
177  trim(mod_name) .NE. 'ice' .AND. trim(mod_name) .NE. 'lnd' ) call mpp_error(fatal, &
178  "data_override_mod: mod_name should be 'atm', 'ocn', 'ice' or 'lnd' ")
179 
180  call mpp_get_compute_domain(domain, isc, iec, jsc, jec)
181 
182  allocate(lon(isc:iec, jsc:jec))
183  allocate(lat(isc:iec, jsc:jec))
184 
185  ! get the grid file to read
186 
187  if(variable_exists(fileobj, trim(mod_name)//'_mosaic_file' )) then
188  call read_data(fileobj, trim(mod_name)//'_mosaic_file', solo_mosaic_file)
189 
190  solo_mosaic_file = 'INPUT/'//trim(solo_mosaic_file)
191  if(.not. open_file(mosaicfileobj, solo_mosaic_file, 'read')) then
192  call mpp_error(fatal, 'data_override_mod(get_grid_version_2): Error in opening solo mosaic file '// &
193  & trim(solo_mosaic_file))
194  endif
195  open_solo_mosaic=.true.
196  else
197  mosaicfileobj = fileobj
198  open_solo_mosaic = .false.
199  end if
200 
201  call get_mosaic_tile_grid(grid_file, mosaicfileobj, domain)
202 
203  if(.not. open_file(tilefileobj, grid_file, 'read')) then
204  call mpp_error(fatal, 'data_override_mod(get_grid_version_2): Error in opening tile file '//trim(grid_file))
205  endif
206 
207  call get_variable_size(tilefileobj, 'area', siz)
208  nlon_super = siz(1); nlat_super = siz(2)
209  if( mod(nlon_super,2) .NE. 0) call mpp_error(fatal, &
210  'data_override_mod: '//trim(mod_name)//' supergrid longitude size can not be divided by 2')
211  if( mod(nlat_super,2) .NE. 0) call mpp_error(fatal, &
212  'data_override_mod: '//trim(mod_name)//' supergrid latitude size can not be divided by 2')
213  nlon = nlon_super/2;
214  nlat = nlat_super/2;
215  call check_grid_sizes(trim(mod_name)//'_domain ', domain, nlon, nlat)
216  isc2 = 2*isc-1; iec2 = 2*iec+1
217  jsc2 = 2*jsc-1; jec2 = 2*jec+1
218 
219  start(1) = isc2; nread(1) = iec2-isc2+1
220  start(2) = jsc2; nread(2) = jec2-jsc2+1
221 
222  allocate(tmpx(isc2:iec2, jsc2:jec2), tmpy(isc2:iec2, jsc2:jec2) )
223 
224  call read_data( tilefileobj, 'x', tmpx, corner=start,edge_lengths=nread)
225  call read_data( tilefileobj, 'y', tmpy, corner=start,edge_lengths=nread)
226 
227  ! copy data onto model grid
228  if(trim(mod_name) == 'atm' .OR. trim(mod_name) == 'lnd' .OR. use_center_grid_points_local) then
229  do j = jsc, jec
230  do i = isc, iec
231  lon(i,j) = tmpx(i*2,j*2)
232  lat(i,j) = tmpy(i*2,j*2)
233  end do
234  end do
235  else
236  do j = jsc, jec
237  do i = isc, iec
238  lon(i,j) = (tmpx(i*2-1,j*2-1)+tmpx(i*2+1,j*2-1)+tmpx(i*2+1,j*2+1)+tmpx(i*2-1,j*2+1))*0.25_lkind
239  lat(i,j) = (tmpy(i*2-1,j*2-1)+tmpy(i*2+1,j*2-1)+tmpy(i*2+1,j*2+1)+tmpy(i*2-1,j*2+1))*0.25_lkind
240  end do
241  end do
242  endif
243 
244  ! convert to radian
245  lon = lon * real(deg_to_rad, lkind)
246  lat = lat * real(deg_to_rad, lkind)
247 
248  deallocate(tmpx, tmpy)
249  min_lon = minval(lon)
250  max_lon = maxval(lon)
251  call mpp_min(min_lon)
252  call mpp_max(max_lon)
253 
254  call close_file(tilefileobj)
255  if(open_solo_mosaic) call close_file(mosaicfileobj)
256 end subroutine get_grid_version_2_