FMS  2025.04
Flexible Modeling System
grid2.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 !> @file
20 
21 !> @brief return grid cell area for the specified model component and tile
22 subroutine get_grid_cell_area_sg_(component, tile, cellarea, domain)
23  character(len=*), intent(in) :: component !< Component model (atm, lnd, ocn)
24  integer , intent(in) :: tile !< Tile number
25  real(kind=fms_mos_kind_) , intent(inout) :: cellarea(:,:) !< Cell area
26  type(domain2d) , intent(in), optional :: domain !< Domain
27 
28  ! local vars
29  integer :: nlon, nlat
30  real(kind=r8_kind), allocatable :: glonb(:,:), glatb(:,:)
31  real(kind=r8_kind), allocatable :: cellarea8(:,:)
32 
33  call init_checks("get_grid_cell_area")
34  allocate(cellarea8(size(cellarea,1),size(cellarea,2)))
35 
36  select case(grid_version)
37  case(version_geolon_t,version_x_t)
38  select case(trim(component))
39  case('LND')
40  call read_data(gridfileobj, 'AREA_LND_CELL', cellarea8)
41  case('ATM','OCN')
42  call read_data(gridfileobj, 'AREA_'//trim(uppercase(component)),cellarea8)
43  case default
44  call mpp_error(fatal, module_name//'/get_grid_cell_area'//&
45  'Illegal component name "'//trim(component)//'": must be one of ATM, LND, or OCN')
46  end select
47  ! convert area to m2
48  cellarea = real( cellarea8*4.0_r8_kind*pi*radius**2, fms_mos_kind_)
49  case(version_ocn_mosaic_file, version_gridfiles)
50  if (present(domain)) then
51  call mpp_get_compute_domain(domain,xsize=nlon,ysize=nlat)
52  else
53  call get_grid_size(component,tile,nlon,nlat)
54  endif
55  allocate(glonb(nlon+1,nlat+1),glatb(nlon+1,nlat+1))
56  call get_grid_cell_vertices(component, tile, glonb, glatb, domain)
57  if (great_circle_algorithm) then
58  call calc_mosaic_grid_great_circle_area(glonb*pi/180.0_r8_kind, glatb*pi/180_r8_kind, cellarea8)
59  cellarea=real(cellarea8,fms_mos_kind_)
60  else
61  call calc_mosaic_grid_area(glonb*pi/180.0_r8_kind, glatb*pi/180_r8_kind, cellarea8)
62  cellarea=real(cellarea8,fms_mos_kind_)
63  end if
64  deallocate(glonb,glatb)
65  end select
66 
67  deallocate(cellarea8)
68 
69 end subroutine get_grid_cell_area_sg_
70 
71 !> @brief get the area of the component per grid cell
72 subroutine get_grid_comp_area_sg_(component,tile,area,domain)
73  character(len=*) :: component !< Component model (atm, lnd, ocn)
74  integer, intent(in) :: tile !< Tile number
75  real(kind=fms_mos_kind_), intent(inout) :: area(:,:) !< Area of grid cell
76  type(domain2d), intent(in), optional :: domain !< Domain
77  ! local vars
78  integer :: n_xgrid_files ! number of exchange grid files in the mosaic
79  integer :: siz(2), nxgrid
80  integer :: i,j,m,n
81  integer, allocatable :: i1(:), j1(:), i2(:), j2(:)
82  real(kind=r8_kind), allocatable :: xgrid_area(:)
83  real(kind=r8_kind), allocatable :: rmask(:,:)
84  character(len=MAX_NAME) :: &
85  xgrid_name, & ! name of the variable holding xgrid names
86  tile_name, & ! name of the tile
87  mosaic_name ! name of the mosaic
88  character(len=FMS_PATH_LEN) :: &
89  tilefile, & ! name of current tile file
90  xgrid_file ! name of the current xgrid file
91  character(len=4096) :: attvalue
92  character(len=MAX_NAME), allocatable :: nest_tile_name(:)
93  integer :: is,ie,js,je ! boundaries of our domain
94  integer :: i0, j0 ! offsets for x and y, respectively
95  integer :: num_nest_tile, ntiles
96  logical :: is_nest
97  integer :: found_xgrid_files ! how many xgrid files we actually found in the grid spec
98  integer :: ibegin, iend, bsize, l
99  type(FmsNetcdfFile_t) :: tilefileobj, xgrid_fileobj
100 
101  real(r8_kind),allocatable :: area8(:,:)
102 
103  call init_checks("get_grid_comp_area")
104  allocate(area8(size(area,1),size(area,2)))
105 
106  select case (grid_version )
107  case(version_geolon_t,version_x_t)
108  select case(component)
109  case('ATM')
110  call read_data(gridfileobj,'AREA_ATM',area8)
111  case('OCN')
112  allocate(rmask(size(area8,1),size(area8,2)))
113  call read_data(gridfileobj,'AREA_OCN',area8)
114  call read_data(gridfileobj,'wet', rmask)
115  area = real(area8*rmask, fms_mos_kind_)
116  deallocate(rmask)
117  case('LND')
118  call read_data(gridfileobj,'AREA_LND',area8)
119  case default
120  call mpp_error(fatal, module_name//'/get_grid_comp_area'//&
121  'Illegal component name "'//trim(component)//'": must be one of ATM, LND, or OCN')
122  end select
123  case(version_ocn_mosaic_file, version_gridfiles) ! mosaic gridspec
124  select case (component)
125  case ('ATM')
126  ! just read the grid cell area and return
127  call get_grid_cell_area(component,tile,area8)
128  area = real(area8, fms_mos_kind_)
129  return
130  case ('LND')
131  xgrid_name = 'aXl_file'
132  call read_data(gridfileobj, 'lnd_mosaic', mosaic_name)
133  tile_name = trim(mosaic_name)//'_tile'//char(tile+ichar('0'))
134  case ('OCN')
135  xgrid_name = 'aXo_file'
136  call read_data(gridfileobj, 'ocn_mosaic', mosaic_name)
137  tile_name = trim(mosaic_name)//'_tile'//char(tile+ichar('0'))
138  case default
139  call mpp_error(fatal, module_name//'/get_grid_comp_area'//&
140  'Illegal component name "'//trim(component)//'": must be one of ATM, LND, or OCN')
141  end select
142  ! get the boundaries of the requested domain
143  if(present(domain)) then
144  call mpp_get_compute_domain(domain,is,ie,js,je)
145  i0 = 1-is ; j0=1-js
146  else
147  call get_grid_size(component,tile,ie,je)
148  is = 1 ; i0 = 0
149  js = 1 ; j0 = 0
150  endif
151  if (size(area8,1)/=ie-is+1.or.size(area8,2)/=je-js+1) &
152  call mpp_error(fatal, module_name//'/get_grid_comp_area '//&
153  'size of the output argument "area" is not consistent with the domain')
154 
155  ! find the nest tile
156  call read_data(gridfileobj, 'atm_mosaic', mosaic_name)
157  call get_grid_ntiles('atm', ntiles)
158  allocate(nest_tile_name(ntiles))
159  num_nest_tile = 0
160  do n = 1, ntiles
161  tilefile = read_file_name(mosaic_fileobj(1), 'gridfiles', n)
162  call open_grid_file(tilefileobj, grid_dir//tilefile)
163  if (global_att_exists(tilefileobj, "nest_grid")) then
164  call get_global_attribute(tilefileobj, "nest_grid", attvalue)
165  if(trim(attvalue) == "TRUE") then
166  num_nest_tile = num_nest_tile + 1
167  nest_tile_name(num_nest_tile) = trim(mosaic_name)//'_tile'//char(n+ichar('0'))
168  else if(trim(attvalue) .NE. "FALSE") then
169  call mpp_error(fatal,module_name//'/get_grid_comp_area value of global attribute nest_grid in file'//&
170  trim(tilefile)//' should be TRUE or FALSE')
171  endif
172  end if
173  call close_file(tilefileobj)
174  end do
175  area8(:,:) = 0.0_r8_kind
176  if(variable_exists(gridfileobj,xgrid_name)) then
177  ! get the number of the exchange-grid files
178  call get_variable_size(gridfileobj,xgrid_name,siz)
179  n_xgrid_files = siz(2)
180  found_xgrid_files = 0
181  ! loop through all exchange grid files
182  do n = 1, n_xgrid_files
183  ! get the name of the current exchange grid file
184  xgrid_file = read_file_name(gridfileobj,xgrid_name,n)
185  call open_grid_file(xgrid_fileobj, grid_dir//xgrid_file)
186  ! skip the rest of the loop if the name of the current tile isn't found
187  ! in the file name, but check this only if there is more than 1 tile
188  if(n_xgrid_files>1) then
189  if(index(xgrid_file,trim(tile_name))==0) cycle
190  endif
191  found_xgrid_files = found_xgrid_files + 1
192  !---make sure the atmosphere grid is not a nested grid
193  is_nest = .false.
194  do m = 1, num_nest_tile
195  if(index(xgrid_file, trim(nest_tile_name(m))) .NE. 0) then
196  is_nest = .true.
197  exit
198  end if
199  end do
200  if(is_nest) cycle
201 
202  ! finally read the exchange grid
203  nxgrid = get_mosaic_xgrid_size(xgrid_fileobj)
204  if(nxgrid < bufsize) then
205  allocate(i1(nxgrid), j1(nxgrid), i2(nxgrid), j2(nxgrid), xgrid_area(nxgrid))
206  else
207  allocate(i1(bufsize), j1(bufsize), i2(bufsize), j2(bufsize), xgrid_area(bufsize))
208  endif
209  ibegin = 1
210  do l = 1,nxgrid,bufsize
211  bsize = min(bufsize, nxgrid-l+1)
212  iend = ibegin + bsize - 1
213  call get_mosaic_xgrid(xgrid_fileobj, i1(1:bsize), j1(1:bsize), i2(1:bsize), j2(1:bsize), &
214  xgrid_area(1:bsize), ibegin, iend)
215  ! and sum the exchange grid areas
216  do m = 1, bsize
217  i = i2(m); j = j2(m)
218  if (i<is.or.i>ie) cycle
219  if (j<js.or.j>je) cycle
220  area8(i+i0,j+j0) = area8(i+i0,j+j0) + xgrid_area(m)
221  end do
222  ibegin = iend + 1
223  enddo
224  deallocate(i1, j1, i2, j2, xgrid_area)
225  call close_file(xgrid_fileobj)
226  enddo
227  if (found_xgrid_files == 0) &
228  call mpp_error(fatal, 'get_grid_comp_area no xgrid files were found for component '&
229  //trim(component)//' (mosaic name is '//trim(mosaic_name)//')')
230 
231  endif
232  deallocate(nest_tile_name)
233  end select ! version
234  ! convert area to m2
235  area = real(area8*4.0_r8_kind*pi*radius**2, fms_mos_kind_)
236 
237  deallocate(area8)
238 
239 end subroutine get_grid_comp_area_sg_
240 
241 !> @brief return grid cell area for the specified model component and tile on an
242 !! unstructured domain
243 subroutine get_grid_cell_area_ug_(component, tile, cellarea, SG_domain, UG_domain)
244  character(len=*), intent(in) :: component !< Component model (atm, lnd, ocn)
245  integer , intent(in) :: tile !< Tile number
246  real(kind=fms_mos_kind_), intent(inout) :: cellarea(:) !< Cell area
247  type(domain2d) , intent(in) :: SG_domain !< Structured Domain
248  type(domainUG) , intent(in) :: UG_domain !< Unstructured Domain
249  integer :: is, ie, js, je
250  real(kind=fms_mos_kind_), allocatable :: sg_area(:,:)
251 
252  call init_checks("get_grid_cell_area")
253  call mpp_get_compute_domain(sg_domain, is, ie, js, je)
254  allocate(sg_area(is:ie, js:je))
255  call get_grid_cell_area(component, tile, sg_area, sg_domain)
256  call mpp_pass_sg_to_ug(ug_domain, sg_area, cellarea)
257  deallocate(sg_area)
258 end subroutine get_grid_cell_area_ug_
259 
260 !> @brief get the area of the component per grid cell for an unstructured domain
261 subroutine get_grid_comp_area_ug_(component, tile, area, SG_domain, UG_domain)
262  character(len=*), intent(in) :: component !< Component model (atm, lnd, ocn)
263  integer , intent(in) :: tile !< Tile number
264  real(kind=fms_mos_kind_), intent(inout) :: area(:) !< Area of the component
265  type(domain2d) , intent(in) :: SG_domain !< Structured domain
266  type(domainUG) , intent(in) :: UG_domain !< Unstructured domain
267  integer :: is, ie, js, je
268  real(kind=fms_mos_kind_), allocatable :: sg_area(:,:)
269 
270  call init_checks("get_grid_comp_area")
271  call mpp_get_compute_domain(sg_domain, is, ie, js, je)
272  allocate(sg_area(is:ie, js:je))
273  call get_grid_comp_area(component, tile, sg_area, sg_domain)
274  call mpp_pass_sg_to_ug(ug_domain, sg_area, area)
275  deallocate(sg_area)
276 
277 end subroutine get_grid_comp_area_ug_
278 
279 !> @brief returns arrays of global grid cell boundaries for given model component and
280 !! mosaic tile number.
281 subroutine get_grid_cell_vertices_1d_(component, tile, glonb, glatb)
282  character(len=*), intent(in) :: component !< Component model (atm, lnd, ocn)
283  integer, intent(in) :: tile !< Tile number
284  real(kind=fms_mos_kind_),intent(inout) :: glonb(:),glatb(:) !< Grid cell vertices
285 
286  integer :: nlon, nlat
287  integer :: start(4), nread(4)
288  real(kind=fms_mos_kind_), allocatable :: tmp(:,:), x_vert_t(:,:,:), y_vert_t(:,:,:)
289  character(len=FMS_PATH_LEN) :: tilefile
290  type(FmsNetcdfFile_t) :: tilefileobj
291 
292  call init_checks("get_grid_cell_vertices")
293  call get_grid_size_for_one_tile(component, tile, nlon, nlat)
294  if (size(glonb(:))/=nlon+1) &
295  call mpp_error (fatal, module_name//'/get_grid_cell_vertices_1D '//&
296  'Size of argument "glonb" is not consistent with the grid size')
297  if (size(glatb(:))/=nlat+1) &
298  call mpp_error (fatal, module_name//'/get_grid_cell_vertices_1D '//&
299  'Size of argument "glatb" is not consistent with the grid size')
300  if(trim(component) .NE. 'ATM' .AND. component .NE. 'LND' .AND. component .NE. 'OCN') then
301  call mpp_error(fatal, module_name//'/get_grid_cell_vertices_1D '//&
302  'Illegal component name "'//trim(component)//'": must be one of ATM, LND, or OCN')
303  endif
304 
305  select case(grid_version)
306  case(version_geolon_t)
307  select case(trim(component))
308  case('ATM','LND')
309  call read_data(gridfileobj, 'xb'//lowercase(component(1:1)), glonb)
310  call read_data(gridfileobj, 'yb'//lowercase(component(1:1)), glatb)
311  case('OCN')
312  call read_data(gridfileobj, "gridlon_vert_t", glonb)
313  call read_data(gridfileobj, "gridlat_vert_t", glatb)
314  end select
315  case(version_x_t)
316  select case(trim(component))
317  case('ATM','LND')
318  call read_data(gridfileobj, 'xb'//lowercase(component(1:1)), glonb)
319  call read_data(gridfileobj, 'yb'//lowercase(component(1:1)), glatb)
320  case('OCN')
321  allocate (x_vert_t(nlon,1,2), y_vert_t(1,nlat,2) )
322  start = 1; nread = 1
323  nread(1) = nlon; nread(2) = 1; start(3) = 1
324  call read_data(gridfileobj, "x_vert_T", x_vert_t(:,:,1), corner=start, edge_lengths=nread)
325  nread(1) = nlon; nread(2) = 1; start(3) = 2
326  call read_data(gridfileobj, "x_vert_T", x_vert_t(:,:,2), corner=start, edge_lengths=nread)
327 
328  nread(1) = 1; nread(2) = nlat; start(3) = 1
329  call read_data(gridfileobj, "y_vert_T", y_vert_t(:,:,1), corner=start, edge_lengths=nread)
330  nread(1) = 1; nread(2) = nlat; start(3) = 4
331  call read_data(gridfileobj, "y_vert_T", y_vert_t(:,:,2), corner=start, edge_lengths=nread)
332  glonb(1:nlon) = x_vert_t(1:nlon,1,1)
333  glonb(nlon+1) = x_vert_t(nlon,1,2)
334  glatb(1:nlat) = y_vert_t(1,1:nlat,1)
335  glatb(nlat+1) = y_vert_t(1,nlat,2)
336  deallocate(x_vert_t, y_vert_t)
337  end select
338  case(version_ocn_mosaic_file, version_gridfiles)
339  ! get the name of the grid file for the component and tile
340  tilefile = read_file_name(mosaic_fileobj(get_component_number(trim(component))), 'gridfiles',tile)
341  call open_grid_file(tilefileobj, grid_dir//tilefile)
342 
343  start = 1; nread = 1
344  nread(1) = 2*nlon+1
345  allocate( tmp(2*nlon+1,1) )
346  call read_data(tilefileobj, "x", tmp, corner=start, edge_lengths=nread)
347  glonb(1:nlon+1) = tmp(1:2*nlon+1:2,1)
348  deallocate(tmp)
349  allocate(tmp(1,2*nlat+1))
350 
351  start = 1; nread = 1
352  nread(2) = 2*nlat+1
353  call read_data(tilefileobj, "y", tmp, corner=start, edge_lengths=nread)
354  glatb(1:nlat+1) = tmp(1,1:2*nlat+1:2)
355  deallocate(tmp)
356  call close_file(tilefileobj)
357  end select
358 end subroutine get_grid_cell_vertices_1d_
359 
360 !> @brief returns cell vertices for the specified model component and mosaic tile number
361 subroutine get_grid_cell_vertices_2d_(component, tile, lonb, latb, domain)
362  character(len=*), intent(in) :: component !< Component model (atm, lnd, ocn)
363  integer, intent(in) :: tile !< Tile number
364  real(kind=fms_mos_kind_), intent(inout) :: lonb(:,:),latb(:,:) !< Cell vertices
365  type(domain2d), optional, intent(in) :: domain !< Domain
366 
367  ! local vars
368  integer :: nlon, nlat
369  integer :: i,j
370  real(kind=fms_mos_kind_), allocatable :: buffer(:), tmp(:,:), x_vert_t(:,:,:), y_vert_t(:,:,:)
371  integer :: is,ie,js,je ! boundaries of our domain
372  integer :: i0,j0 ! offsets for coordinates
373  integer :: isg, jsg
374  integer :: start(4), nread(4)
375  character(len=FMS_PATH_LEN) :: tilefile
376  type(FmsNetcdfFile_t) :: tilefileobj
377 
378  call init_checks("get_grid_cell_vertices")
379  call get_grid_size_for_one_tile(component, tile, nlon, nlat)
380 
381  if (present(domain)) then
382  call mpp_get_compute_domain(domain,is,ie,js,je)
383  else
384  is = 1 ; ie = nlon
385  js = 1 ; je = nlat
386  !--- domain normally should be present
387  call mpp_error (note, module_name//'/get_grid_cell_vertices '//&
388  'domain is not present, global data will be read')
389  endif
390  i0 = -is+1; j0 = -js+1
391 
392  ! verify that lonb and latb sizes are consistent with the size of domain
393  if (size(lonb,1)/=ie-is+2.or.size(lonb,2)/=je-js+2) &
394  call mpp_error (fatal, module_name//'/get_grid_cell_vertices '//&
395  'Size of argument "lonb" is not consistent with the domain size')
396  if (size(latb,1)/=ie-is+2.or.size(latb,2)/=je-js+2) &
397  call mpp_error (fatal, module_name//'/get_grid_cell_vertices '//&
398  'Size of argument "latb" is not consistent with the domain size')
399  if(trim(component) .NE. 'ATM' .AND. component .NE. 'LND' .AND. component .NE. 'OCN') then
400  call mpp_error(fatal, module_name//'/get_grid_cell_vertices '//&
401  'Illegal component name "'//trim(component)//'": must be one of ATM, LND, or OCN')
402  endif
403 
404  !! use lonb, latb as r4
405  select case(grid_version)
406  case(version_geolon_t)
407  select case(component)
408  case('ATM','LND')
409  allocate(buffer(max(nlon,nlat)+1))
410  ! read coordinates of grid cell vertices
411  call read_data(gridfileobj, 'xb'//lowercase(component(1:1)), buffer(1:nlon+1))
412  do j = js, je+1
413  do i = is, ie+1
414  lonb(i+i0,j+j0) = buffer(i)
415  enddo
416  enddo
417  call read_data(gridfileobj, 'yb'//lowercase(component(1:1)), buffer(1:nlat+1))
418  do j = js, je+1
419  do i = is, ie+1
420  latb(i+i0,j+j0) = buffer(j)
421  enddo
422  enddo
423  deallocate(buffer)
424  case('OCN')
425  if (present(domain)) then
426  start = 1; nread = 1
427  start(1) = is; start(2) = js
428  nread(1) = ie-is+2; nread(2) = je-js+2
429  call read_data(gridfileobj, "geolon_vert_t", lonb, corner=start, edge_lengths=nread)
430  call read_data(gridfileobj, "geolat_vert_t", latb, corner=start, edge_lengths=nread)
431  else
432  call read_data(gridfileobj, "geolon_vert_t", lonb)
433  call read_data(gridfileobj, "geolat_vert_t", latb)
434  endif
435  end select
436  case(version_x_t)
437  select case(component)
438  case('ATM','LND')
439  allocate(buffer(max(nlon,nlat)+1))
440  ! read coordinates of grid cell vertices
441  call read_data(gridfileobj, 'xb'//lowercase(component(1:1)), buffer(1:nlon+1))
442  do j = js, je+1
443  do i = is, ie+1
444  lonb(i+i0,j+j0) = buffer(i)
445  enddo
446  enddo
447  call read_data(gridfileobj, 'yb'//lowercase(component(1:1)), buffer(1:nlat+1))
448  do j = js, je+1
449  do i = is, ie+1
450  latb(i+i0,j+j0) = buffer(j)
451  enddo
452  enddo
453  deallocate(buffer)
454  case('OCN')
455  nlon=ie-is+1; nlat=je-js+1
456  allocate (x_vert_t(nlon,nlat,4), y_vert_t(nlon,nlat,4) )
457  call read_data(gridfileobj, 'x_vert_T', x_vert_t)
458  call read_data(gridfileobj, 'y_vert_T', y_vert_t)
459  lonb(1:nlon,1:nlat) = x_vert_t(1:nlon,1:nlat,1)
460  lonb(nlon+1,1:nlat) = x_vert_t(nlon,1:nlat,2)
461  lonb(1:nlon,nlat+1) = x_vert_t(1:nlon,nlat,4)
462  lonb(nlon+1,nlat+1) = x_vert_t(nlon,nlat,3)
463  latb(1:nlon,1:nlat) = y_vert_t(1:nlon,1:nlat,1)
464  latb(nlon+1,1:nlat) = y_vert_t(nlon,1:nlat,2)
465  latb(1:nlon,nlat+1) = y_vert_t(1:nlon,nlat,4)
466  latb(nlon+1,nlat+1) = y_vert_t(nlon,nlat,3)
467  deallocate(x_vert_t, y_vert_t)
468  end select
469  case(version_ocn_mosaic_file, version_gridfiles)
470  ! get the name of the grid file for the component and tile
471  tilefile = read_file_name(mosaic_fileobj(get_component_number(trim(component))), 'gridfiles',tile)
472  call open_grid_file(tilefileobj, grid_dir//tilefile)
473  if(PRESENT(domain)) then
474  call mpp_get_global_domain(domain, xbegin=isg, ybegin=jsg)
475  start = 1; nread = 1
476  start(1) = 2*(is-isg+1) - 1; nread(1) = 2*(ie-is)+3
477  start(2) = 2*(js-jsg+1) - 1; nread(2) = 2*(je-js)+3
478  allocate(tmp(nread(1), nread(2)) )
479  call read_data(tilefileobj, "x", tmp, corner=start, edge_lengths=nread)
480  do j = 1, je-js+2
481  do i = 1, ie-is+2
482  lonb(i,j) = tmp(2*i-1,2*j-1)
483  enddo
484  enddo
485  call read_data(tilefileobj, "y", tmp, corner=start, edge_lengths=nread)
486  do j = 1, je-js+2
487  do i = 1, ie-is+2
488  latb(i,j) = tmp(2*i-1,2*j-1)
489  enddo
490  enddo
491  else
492  allocate(tmp(2*nlon+1,2*nlat+1))
493  call read_data(tilefileobj, "x", tmp)
494  do j = js, je+1
495  do i = is, ie+1
496  lonb(i+i0,j+j0) = tmp(2*i-1,2*j-1)
497  end do
498  end do
499  call read_data(tilefileobj, "y", tmp)
500  do j = js, je+1
501  do i = is, ie+1
502  latb(i+i0,j+j0) = tmp(2*i-1,2*j-1)
503  end do
504  end do
505  endif
506  deallocate(tmp)
507  call close_file(tilefileobj)
508  end select ! end grid_version
509  end subroutine get_grid_cell_vertices_2d_
510 
511 !> @brief returns cell vertices for the specified model component and mosaic tile number for
512 !! an unstructured domain
513 subroutine get_grid_cell_vertices_ug_(component, tile, lonb, latb, SG_domain, UG_domain)
514  character(len=*), intent(in) :: component !< Component model (atm, lnd, ocn)
515  integer, intent(in) :: tile !< Tile number
516  real(kind=fms_mos_kind_), intent(inout) :: lonb(:,:),latb(:,:) ! The second dimension is 4
517  type(domain2d) , intent(in) :: SG_domain !< Structured domain
518  type(domainUG) , intent(in) :: UG_domain !< Unstructured domain
519  integer :: is, ie, js, je, i, j
520  real(kind=fms_mos_kind_), allocatable :: sg_lonb(:,:), sg_latb(:,:), tmp(:,:,:)
521 
522  call init_checks("get_grid_cell_vertices")
523  call mpp_get_compute_domain(sg_domain, is, ie, js, je)
524  allocate(sg_lonb(is:ie+1, js:je+1))
525  allocate(sg_latb(is:ie+1, js:je+1))
526  allocate(tmp(is:ie,js:je,4))
527  call get_grid_cell_vertices(component, tile, sg_lonb, sg_latb, sg_domain)
528  do j = js, je
529  do i = is, ie
530  tmp(i,j,1) = sg_lonb(i,j)
531  tmp(i,j,2) = sg_lonb(i+1,j)
532  tmp(i,j,3) = sg_lonb(i+1,j+1)
533  tmp(i,j,4) = sg_lonb(i,j+1)
534  enddo
535  enddo
536  call mpp_pass_sg_to_ug(ug_domain, tmp, lonb)
537  do j = js, je
538  do i = is, ie
539  tmp(i,j,1) = sg_latb(i,j)
540  tmp(i,j,2) = sg_latb(i+1,j)
541  tmp(i,j,3) = sg_latb(i+1,j+1)
542  tmp(i,j,4) = sg_latb(i,j+1)
543  enddo
544  enddo
545  call mpp_pass_sg_to_ug(ug_domain, tmp, latb)
546 
547 
548  deallocate(sg_lonb, sg_latb, tmp)
549 end subroutine get_grid_cell_vertices_ug_
550 
551 !> @brief returns grid cell centers given model component and mosaic tile number
552 subroutine get_grid_cell_centers_1d_(component, tile, glon, glat)
553  character(len=*), intent(in) :: component !< Component model (atm, lnd, ocn)
554  integer, intent(in) :: tile !< Tile number
555  real(kind=fms_mos_kind_), intent(inout) :: glon(:),glat(:) !< Grid cell centers
556 
557  integer :: nlon, nlat
558  integer :: start(4), nread(4)
559  real(kind=fms_mos_kind_), allocatable :: tmp(:,:)
560  character(len=FMS_PATH_LEN) :: tilefile
561  type(FmsNetcdfFile_t) :: tilefileobj
562 
563  call init_checks("get_grid_cell_centers")
564  call get_grid_size_for_one_tile(component, tile, nlon, nlat)
565  if (size(glon(:))/=nlon) &
566  call mpp_error (fatal, module_name//'/get_grid_cell_centers_1D '//&
567  'Size of argument "glon" is not consistent with the grid size')
568  if (size(glat(:))/=nlat) &
569  call mpp_error (fatal, module_name//'/get_grid_cell_centers_1D '//&
570  'Size of argument "glat" is not consistent with the grid size')
571  if(trim(component) .NE. 'ATM' .AND. component .NE. 'LND' .AND. component .NE. 'OCN') then
572  call mpp_error(fatal, module_name//'/get_grid_cell_centers_1D '//&
573  'Illegal component name "'//trim(component)//'": must be one of ATM, LND, or OCN')
574  endif
575 
576  select case(grid_version)
577  case(version_geolon_t)
578  select case(trim(component))
579  case('ATM','LND')
580  call read_data(gridfileobj, 'xt'//lowercase(component(1:1)), glon)
581  call read_data(gridfileobj, 'yt'//lowercase(component(1:1)), glat)
582  case('OCN')
583  call read_data(gridfileobj, "gridlon_t", glon)
584  call read_data(gridfileobj, "gridlat_t", glat)
585  end select
586  case(version_x_t)
587  select case(trim(component))
588  case('ATM','LND')
589  call read_data(gridfileobj, 'xt'//lowercase(component(1:1)), glon)
590  call read_data(gridfileobj, 'yt'//lowercase(component(1:1)), glat)
591  case('OCN')
592  call read_data(gridfileobj, "grid_x_T", glon)
593  call read_data(gridfileobj, "grid_y_T", glat)
594  end select
595  case(version_ocn_mosaic_file, version_gridfiles)
596  ! get the name of the grid file for the component and tile
597  tilefile = read_file_name(mosaic_fileobj(get_component_number(trim(component))), 'gridfiles',tile)
598  call open_grid_file(tilefileobj, grid_dir//tilefile)
599 
600  start = 1; nread = 1
601  nread(1) = 2*nlon+1; start(2) = 2
602  allocate( tmp(2*nlon+1,1) )
603  call read_data(tilefileobj, "x", tmp, corner=start, edge_lengths=nread)
604  glon(1:nlon) = tmp(2:2*nlon:2,1)
605  deallocate(tmp)
606  allocate(tmp(1, 2*nlat+1))
607 
608  start = 1; nread = 1
609  nread(2) = 2*nlat+1; start(1) = 2
610  call read_data(tilefileobj, "y", tmp, corner=start, edge_lengths=nread)
611  glat(1:nlat) = tmp(1,2:2*nlat:2)
612  deallocate(tmp)
613  call close_file(tilefileobj)
614  end select
615 end subroutine get_grid_cell_centers_1d_
616 
617 !> @brief returns grid cell centers given model component and mosaic tile number
618 subroutine get_grid_cell_centers_2d_(component, tile, lon, lat, domain)
619  character(len=*), intent(in) :: component !< Component model (atm, lnd, ocn)
620  integer, intent(in) :: tile !< Tile number
621  real(kind=fms_mos_kind_), intent(inout) :: lon(:,:),lat(:,:) !< Grid cell centers
622  type(domain2d), intent(in), optional :: domain !< Domain
623  ! local vars
624  integer :: nlon, nlat
625  integer :: i,j
626  real(kind=fms_mos_kind_), allocatable :: buffer(:),tmp(:,:)
627  integer :: is,ie,js,je ! boundaries of our domain
628  integer :: i0,j0 ! offsets for coordinates
629  integer :: isg, jsg
630  integer :: start(4), nread(4)
631  character(len=FMS_PATH_LEN) :: tilefile
632  type(FmsNetcdfFile_t) :: tilefileobj
633 
634  call init_checks("get_grid_cell_centers")
635  call get_grid_size_for_one_tile(component, tile, nlon, nlat)
636  if (present(domain)) then
637  call mpp_get_compute_domain(domain,is,ie,js,je)
638  else
639  is = 1 ; ie = nlon
640  js = 1 ; je = nlat
641  !--- domain normally should be present
642  call mpp_error (note, module_name//'/get_grid_cell_centers '//&
643  'domain is not present, global data will be read')
644  endif
645  i0 = -is+1; j0 = -js+1
646 
647  ! verify that lon and lat sizes are consistent with the size of domain
648  if (size(lon,1)/=ie-is+1.or.size(lon,2)/=je-js+1) &
649  call mpp_error (fatal, module_name//'/get_grid_cell_centers '//&
650  'Size of array "lon" is not consistent with the domain size')
651  if (size(lat,1)/=ie-is+1.or.size(lat,2)/=je-js+1) &
652  call mpp_error (fatal, module_name//'/get_grid_cell_centers '//&
653  'Size of array "lat" is not consistent with the domain size')
654  if(trim(component) .NE. 'ATM' .AND. component .NE. 'LND' .AND. component .NE. 'OCN') then
655  call mpp_error(fatal, module_name//'/get_grid_cell_vertices '//&
656  'Illegal component name "'//trim(component)//'": must be one of ATM, LND, or OCN')
657  endif
658 
659  select case(grid_version)
660  case(version_geolon_t)
661  select case (trim(component))
662  case('ATM','LND')
663  allocate(buffer(max(nlon,nlat)))
664  ! read coordinates of grid cell vertices
665  call read_data(gridfileobj, 'xt'//lowercase(component(1:1)), buffer(1:nlon))
666  do j = js,je
667  do i = is,ie
668  lon(i+i0,j+j0) = buffer(i)
669  enddo
670  enddo
671  call read_data(gridfileobj, 'yt'//lowercase(component(1:1)), buffer(1:nlat))
672  do j = js,je
673  do i = is,ie
674  lat(i+i0,j+j0) = buffer(j)
675  enddo
676  enddo
677  deallocate(buffer)
678  case('OCN')
679  call read_data(gridfileobj, 'geolon_t', lon)
680  call read_data(gridfileobj, 'geolat_t', lat)
681  end select
682  case(version_x_t)
683  select case(trim(component))
684  case('ATM','LND')
685  allocate(buffer(max(nlon,nlat)))
686  ! read coordinates of grid cell vertices
687  call read_data(gridfileobj, 'xt'//lowercase(component(1:1)), buffer(1:nlon))
688  do j = js,je
689  do i = is,ie
690  lon(i+i0,j+j0) = buffer(i)
691  enddo
692  enddo
693  call read_data(gridfileobj, 'yt'//lowercase(component(1:1)), buffer(1:nlat))
694  do j = js,je
695  do i = is,ie
696  lat(i+i0,j+j0) = buffer(j)
697  enddo
698  enddo
699  deallocate(buffer)
700  case('OCN')
701  call read_data(gridfileobj, 'x_T', lon)
702  call read_data(gridfileobj, 'y_T', lat)
703  end select
704  case(version_ocn_mosaic_file, version_gridfiles) ! mosaic grid file
705  ! get the name of the grid file for the component and tile
706  tilefile = read_file_name(mosaic_fileobj(get_component_number(trim(component))), 'gridfiles',tile)
707  call open_grid_file(tilefileobj, grid_dir//tilefile)
708 
709  if(PRESENT(domain)) then
710  call mpp_get_global_domain(domain, xbegin=isg, ybegin=jsg)
711  start = 1; nread = 1
712  start(1) = 2*(is-isg+1) - 1; nread(1) = 2*(ie-is)+3
713  start(2) = 2*(js-jsg+1) - 1; nread(2) = 2*(je-js)+3
714  allocate(tmp(nread(1), nread(2)))
715  call read_data(tilefileobj, "x", tmp, corner=start, edge_lengths=nread)
716  do j = 1, je-js+1
717  do i = 1, ie-is+1
718  lon(i,j) = tmp(2*i,2*j)
719  enddo
720  enddo
721  call read_data(tilefileobj, "y", tmp, corner=start, edge_lengths=nread)
722  do j = 1, je-js+1
723  do i = 1, ie-is+1
724  lat(i,j) = tmp(2*i,2*j)
725  enddo
726  enddo
727  else
728  allocate(tmp(2*nlon+1,2*nlat+1))
729  call read_data(tilefileobj, 'x', tmp)
730  do j = js,je
731  do i = is,ie
732  lon(i+i0,j+j0) = tmp(2*i,2*j)
733  end do
734  end do
735  call read_data(tilefileobj, 'y', tmp)
736  do j = js,je
737  do i = is,ie
738  lat(i+i0,j+j0) = tmp(2*i,2*j)
739  end do
740  end do
741  deallocate(tmp)
742  endif
743  call close_file(tilefileobj)
744  end select
745 end subroutine get_grid_cell_centers_2d_
746 
747 !> @brief returns grid cell centers given model component and mosaic tile number
748 !! for unstructured domain
749 subroutine get_grid_cell_centers_ug_(component, tile, lon, lat, SG_domain, UG_domain)
750  character(len=*), intent(in) :: component !< Component model (atm, lnd, ocn)
751  integer, intent(in) :: tile !< Tile number
752  real(kind=fms_mos_kind_), intent(inout) :: lon(:),lat(:) !< Grid cell centers
753  type(domain2d) , intent(in) :: SG_domain !< Structured domain
754  type(domainUG) , intent(in) :: UG_domain !< Unstructured domain
755  integer :: is, ie, js, je
756  real(kind=fms_mos_kind_), allocatable :: sg_lon(:,:), sg_lat(:,:)
757 
758  call init_checks("get_grid_cell_centers")
759  call mpp_get_compute_domain(sg_domain, is, ie, js, je)
760  allocate(sg_lon(is:ie, js:je))
761  allocate(sg_lat(is:ie, js:je))
762  call get_grid_cell_centers(component, tile, sg_lon, sg_lat, sg_domain)
763  call mpp_pass_sg_to_ug(ug_domain, sg_lon, lon)
764  call mpp_pass_sg_to_ug(ug_domain, sg_lat, lat)
765  deallocate(sg_lon, sg_lat)
766 end subroutine get_grid_cell_centers_ug_
767 
768 !> @}
769 ! close documentation grouping