FMS 2025.01.02-dev
Flexible Modeling System
Loading...
Searching...
No Matches
get_grid_version.inc
1!***********************************************************************
2!* GNU Lesser General Public License
3!*
4!* This file is part of the GFDL Flexible Modeling System (FMS).
5!*
6!* FMS is free software: you can redistribute it and/or modify it under
7!* the terms of the GNU Lesser General Public License as published by
8!* the Free Software Foundation, either version 3 of the License, or (at
9!* your option) any later version.
10!*
11!* FMS is distributed in the hope that it will be useful, but WITHOUT
12!* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13!* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14!* for more details.
15!*
16!* You should have received a copy of the GNU Lesser General Public
17!* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18!***********************************************************************
19
20!> Get global lon and lat of three model (target) grids, with a given file name
21subroutine get_grid_version_1_(grid_file, mod_name, domain, isc, iec, jsc, jec, lon, lat, min_lon, max_lon)
22 integer, parameter :: lkind = fms_get_grid_version_kind_
23
24 character(len=*), intent(in) :: grid_file !< name of grid file
25 character(len=*), intent(in) :: mod_name !< module name
26 type(domain2d), intent(in) :: domain !< 2D domain
27 integer, intent(in) :: isc, iec, jsc, jec
28 real(lkind), dimension(isc:,jsc:), intent(out) :: lon, lat
29 real(lkind), intent(out) :: min_lon, max_lon
30
31 integer :: i, j, siz(4)
32 integer :: nlon, nlat !< size of global lon and lat
33 real(lkind), dimension(:,:,:), allocatable :: lon_vert, lat_vert !< of OCN grid vertices
34 real(lkind), dimension(:), allocatable :: glon, glat !< lon and lat of 1-D grid of atm/lnd
35 logical :: is_new_grid
36 integer :: is, ie, js, je
37 integer :: isd, ied, jsd, jed
38 integer :: isg, ieg, jsg, jeg
39 character(len=3) :: xname, yname
40 integer :: start(2), nread(2)
41 type(FmsNetcdfDomainFile_t) :: fileobj
42 integer :: ndims !< Number of dimensions
43
44 if(.not. open_file(fileobj, grid_file, 'read', domain )) then
45 call mpp_error(fatal, 'data_override_mod(get_grid_version_1): Error in opening file '//trim(grid_file))
46 endif
47
48 call mpp_get_data_domain(domain, isd, ied, jsd, jed)
49 call mpp_get_global_domain(domain, isg, ieg, jsg, jeg)
50
51 select case(mod_name)
52 case('ocn', 'ice')
53 is_new_grid = .false.
54 if(variable_exists(fileobj, 'x_T')) then
55 is_new_grid = .true.
56 else if(variable_exists(fileobj, 'geolon_t')) then
57 is_new_grid = .false.
58 else
59 call mpp_error(fatal,'data_override: both x_T and geolon_t is not in the grid file '//trim(grid_file) )
60 endif
61
62 if(is_new_grid) then
63 ndims = get_variable_num_dimensions(fileobj, 'x_T')
64 call get_variable_size(fileobj, 'x_T', siz(1:ndims))
65 nlon = siz(1); nlat = siz(2)
66 call check_grid_sizes(trim(mod_name)//'_domain ', domain, nlon, nlat)
67 allocate(lon_vert(isc:iec,jsc:jec,4), lat_vert(isc:iec,jsc:jec,4) )
68
69 call read_data(fileobj, 'x_vert_T', lon_vert)
70 call read_data(fileobj, 'y_vert_T', lat_vert)
71
72!2 Global lon and lat of ocean grid cell centers are determined from adjacent vertices
73 lon(:,:) = (lon_vert(:,:,1) + lon_vert(:,:,2) + lon_vert(:,:,3) + lon_vert(:,:,4))*0.25_lkind
74 lat(:,:) = (lat_vert(:,:,1) + lat_vert(:,:,2) + lat_vert(:,:,3) + lat_vert(:,:,4))*0.25_lkind
75 else
76
77 ndims = get_variable_num_dimensions(fileobj, 'geolon_vert_t')
78 call get_variable_size(fileobj, 'geolon_vert_t', siz(1:ndims))
79 nlon = siz(1) - 1; nlat = siz(2) - 1;
80 call check_grid_sizes(trim(mod_name)//'_domain ', domain, nlon, nlat)
81
82 start(1) = isc; nread(1) = iec-isc+2
83 start(2) = jsc; nread(2) = jec-jsc+2
84
85 allocate(lon_vert(isc:iec+1,jsc:jec+1,1))
86 allocate(lat_vert(isc:iec+1,jsc:jec+1,1))
87
88 call read_data(fileobj, 'geolon_vert_t', lon_vert(:,:,1), corner=start, edge_lengths=nread)
89 call read_data(fileobj, 'geolat_vert_t', lat_vert(:,:,1), corner=start, edge_lengths=nread)
90
91 do j = jsc, jec
92 do i = isc, iec
93 lon(i,j) = (lon_vert(i,j,1) + lon_vert(i+1,j,1) + &
94 lon_vert(i+1,j+1,1) + lon_vert(i,j+1,1))*0.25_lkind
95 lat(i,j) = (lat_vert(i,j,1) + lat_vert(i+1,j,1) + &
96 lat_vert(i+1,j+1,1) + lat_vert(i,j+1,1))*0.25_lkind
97 enddo
98 enddo
99 endif
100 deallocate(lon_vert)
101 deallocate(lat_vert)
102 case('atm', 'lnd')
103 if(trim(mod_name) == 'atm') then
104 xname = 'xta'; yname = 'yta'
105 else
106 xname = 'xtl'; yname = 'ytl'
107 endif
108 ndims = get_variable_num_dimensions(fileobj, xname)
109 call get_variable_size(fileobj, xname, siz(1:ndims))
110 nlon = siz(1); allocate(glon(nlon))
111 call read_data(fileobj, xname, glon)
112
113 ndims = get_variable_num_dimensions(fileobj, xname)
114 call get_variable_size(fileobj, yname, siz(1:ndims))
115 nlat = siz(1); allocate(glat(nlat))
116 call read_data(fileobj, yname, glat)
117 call check_grid_sizes(trim(mod_name)//'_domain ', domain, nlon, nlat)
118
119 is = isc - isg + 1; ie = iec - isg + 1
120 js = jsc - jsg + 1; je = jec - jsg + 1
121 do j = js, jec
122 do i = is, ie
123 lon(i,j) = glon(i)
124 lat(i,j) = glat(j)
125 enddo
126 enddo
127 deallocate(glon)
128 deallocate(glat)
129 case default
130 call mpp_error(fatal, "data_override_mod: mod_name should be 'atm', 'ocn', 'ice' or 'lnd' ")
131 end select
132
133 call close_file(fileobj)
134
135 ! convert from degree to radian
136 lon = lon * real(deg_to_rad, lkind)
137 lat = lat* real(deg_to_rad, lkind)
138 min_lon = minval(lon)
139 max_lon = maxval(lon)
140 call mpp_min(min_lon)
141 call mpp_max(max_lon)
142end subroutine get_grid_version_1_
143
144!> Get global lon and lat of three model (target) grids from mosaic.nc.
145!! Currently we assume the refinement ratio is 2 and there is one tile on each pe.
146subroutine get_grid_version_2_(fileobj, mod_name, domain, isc, iec, jsc, jec, lon, lat, min_lon, max_lon, &
147 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 integer, intent(in) :: isc, iec, jsc, jec
154 real(lkind), dimension(isc:,jsc:), intent(out) :: lon, lat
155 real(lkind), intent(out) :: min_lon, max_lon
156 logical, optional, intent(in) :: use_center_grid_points !< Flag indicating whether or not to use the
157 !! centroid values of the supergrid from the
158 !! grid file as opposed to calcuating it by
159 !! taking the average of the four corner points.
160 !! This is only relevant to OCN and ICE grids.
161
162 integer :: i, j, siz(2)
163 integer :: nlon, nlat ! size of global grid
164 integer :: nlon_super, nlat_super ! size of global supergrid.
165 integer :: isd, ied, jsd, jed
166 integer :: isg, ieg, jsg, jeg
167 integer :: isc2, iec2, jsc2, jec2
168 character(len=FMS_PATH_LEN) :: solo_mosaic_file, grid_file
169 real(lkind), allocatable :: tmpx(:,:), tmpy(:,:)
170 logical :: open_solo_mosaic
171 type(FmsNetcdfFile_t) :: mosaicfileobj, tilefileobj
172 integer :: start(2), nread(2)
173 logical :: use_center_grid_points_local
174
175 use_center_grid_points_local = .false.
176 if (present(use_center_grid_points)) use_center_grid_points_local = use_center_grid_points
177
178 if(trim(mod_name) .NE. 'atm' .AND. trim(mod_name) .NE. 'ocn' .AND. &
179 trim(mod_name) .NE. 'ice' .AND. trim(mod_name) .NE. 'lnd' ) call mpp_error(fatal, &
180 "data_override_mod: mod_name should be 'atm', 'ocn', 'ice' or 'lnd' ")
181
182 call mpp_get_data_domain(domain, isd, ied, jsd, jed)
183 call mpp_get_global_domain(domain, isg, ieg, jsg, jeg)
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)
256end subroutine get_grid_version_2_