FV3 Bundle
grid.F90
Go to the documentation of this file.
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 module grid_mod
20 
21 use mpp_mod, only : mpp_root_pe
22 use constants_mod, only : pi, radius
23 use fms_mod, only : uppercase, lowercase, field_exist, field_size, read_data, &
24  error_mesg, string, fatal, note
28 
29 ! the following two use statement are only needed for define_cube_mosaic
30 use mpp_domains_mod, only : domain2d, mpp_define_mosaic, mpp_get_compute_domain, &
33 
34 implicit none;private
35 
36 ! ==== public interfaces =====================================================
37 ! grid dimension inquiry subroutines
38 public :: get_grid_ntiles ! returns number of tiles
39 public :: get_grid_size ! returns horizontal sizes of the grid
40 ! grid geometry inquiry subroutines
41 public :: get_grid_cell_centers
42 public :: get_grid_cell_vertices
43 ! grid area inquiry subroutines
44 public :: get_grid_cell_area
45 public :: get_grid_comp_area
46 ! decompose cubed sphere domains -- probably does not belong here, but it should
47 ! be in some place available for component models
48 public :: define_cube_mosaic
49 ! ==== end of public interfaces ==============================================
50 
51 interface get_grid_size
52  module procedure get_grid_size_for_all_tiles
53  module procedure get_grid_size_for_one_tile
54 end interface
55 
57  module procedure get_grid_cell_vertices_1d
58  module procedure get_grid_cell_vertices_2d
59  module procedure get_grid_cell_vertices_ug
60 end interface
61 
63  module procedure get_grid_cell_centers_1d
64  module procedure get_grid_cell_centers_2d
65  module procedure get_grid_cell_centers_ug
66 end interface
67 
69  module procedure get_grid_cell_area_sg
70  module procedure get_grid_cell_area_ug
71 end interface get_grid_cell_area
72 
74  module procedure get_grid_comp_area_sg
75  module procedure get_grid_comp_area_ug
76 end interface get_grid_comp_area
77 
78 ! ==== module constants ======================================================
79 character(len=*), parameter :: &
80  module_name = 'grid_mod'
81 
82 ! Include variable "version" to be written to log file.
83 #include<file_version.h>
84 
85 character(len=*), parameter :: &
86  grid_dir = 'INPUT/', & ! root directory for all grid files
87  grid_file = 'INPUT/grid_spec.nc' ! name of the grid spec file
88 
89 integer, parameter :: &
90  max_name = 256, & ! max length of the variable names
91  max_file = 1024, & ! max length of the file names
92  version_0 = 0, &
93  version_1 = 1, &
94  version_2 = 2
95 
96 integer, parameter :: bufsize = 1048576 ! This is used to control memory usage in get_grid_comp_area
97  ! We may change this to a namelist variable is needed.
98 
99 ! ==== module variables ======================================================
100 integer :: grid_version = -1
101 logical :: great_circle_algorithm = .false.
102 logical :: first_call = .true.
103 
104 
105 contains
106 
107 function get_grid_version()
108  integer :: get_grid_version
109 
110  if(first_call) then
112  first_call = .false.
113  endif
114 
115  if(grid_version<0) then
116  if(field_exist(grid_file, 'geolon_t')) then
118  else if(field_exist(grid_file, 'x_T')) then
120  else if(field_exist(grid_file, 'ocn_mosaic_file') ) then
122  else
123  call error_mesg(module_name//'/get_grid_version',&
124  'Can''t determine the version of the grid spec: none of "x_T", "geolon_t", or "ocn_mosaic_file" exist in file "'//trim(grid_file)//'"', &
125  fatal )
126  endif
127  endif
129 end function get_grid_version
130 
131 
132 ! ============================================================================
133 ! returns number of tiles for a given component
134 ! ============================================================================
135 subroutine get_grid_ntiles(component,ntiles)
136  character(len=*) :: component
137  integer, intent(out) :: ntiles
138 
139  ! local vars
140  character(len=MAX_FILE) :: component_mosaic
141 
142  select case (get_grid_version())
143  case(version_0,version_1)
144  ntiles = 1
145  case(version_2)
146  call read_data(grid_file,trim(lowercase(component))//'_mosaic_file',component_mosaic)
147  ntiles = get_mosaic_ntiles(grid_dir//trim(component_mosaic))
148  end select
149 end subroutine get_grid_ntiles
150 
151 
152 ! ============================================================================
153 ! returns size of the grid for each of the tiles
154 ! ============================================================================
155 subroutine get_grid_size_for_all_tiles(component,nx,ny)
156  character(len=*) :: component
157  integer, intent(inout) :: nx(:),ny(:)
158 
159  ! local vars
160  integer :: siz(4) ! for the size of external fields
161  character(len=MAX_NAME) :: varname1, varname2
162  character(len=MAX_FILE) :: component_mosaic
163 
164  varname1 = 'AREA_'//trim(uppercase(component))
165  varname2 = trim(lowercase(component))//'_mosaic_file'
166 
167  select case (get_grid_version())
168  case(version_0,version_1)
169  call field_size(grid_file, varname1, siz)
170  nx(1) = siz(1); ny(1)=siz(2)
171  case(version_2) ! mosaic file
172  call read_data(grid_file,varname2, component_mosaic)
173  call get_mosaic_grid_sizes(grid_dir//trim(component_mosaic),nx,ny)
174  end select
175 end subroutine get_grid_size_for_all_tiles
176 
177 
178 ! ============================================================================
179 ! returns size of the grid for one of the tiles
180 ! ============================================================================
181 subroutine get_grid_size_for_one_tile(component,tile,nx,ny)
182  character(len=*) :: component
183  integer, intent(in) :: tile
184  integer, intent(inout) :: nx,ny
185 
186  ! local vars
187  integer, allocatable :: nnx(:), nny(:)
188  integer :: ntiles
189 
190  call get_grid_ntiles(component, ntiles)
191  if(tile>0.and.tile<=ntiles) then
192  allocate(nnx(ntiles),nny(ntiles))
193  call get_grid_size_for_all_tiles(component,nnx,nny)
194  nx = nnx(tile); ny = nny(tile)
195  deallocate(nnx,nny)
196  else
197  call error_mesg('get_grid_size',&
198  'requested tile index '//trim(string(tile))//' is out of bounds (1:'//trim(string(ntiles))//')',&
199  fatal)
200  endif
201 end subroutine get_grid_size_for_one_tile
202 
203 ! ============================================================================
204 ! return grid cell area for the specified model component and tile
205 ! ============================================================================
206 subroutine get_grid_cell_area_sg(component, tile, cellarea, domain)
207  character(len=*), intent(in) :: component
208  integer , intent(in) :: tile
209  real , intent(inout) :: cellarea(:,:)
210  type(domain2d) , intent(in), optional :: domain
211 
212  ! local vars
213  integer :: nlon, nlat
214  real, allocatable :: glonb(:,:), glatb(:,:)
215 
216  select case(get_grid_version())
217  case(version_0,version_1)
218  select case(trim(component))
219  case('LND')
220  call read_data(grid_file, 'AREA_LND_CELL', cellarea, &
221  no_domain=.not.present(domain), domain=domain)
222  case('ATM','OCN')
223  call read_data(grid_file, 'AREA_'//trim(uppercase(component)),cellarea,&
224  no_domain=.not.present(domain),domain=domain)
225  case default
226  call error_mesg(module_name//'/get_grid_cell_area',&
227  'Illegal component name "'//trim(component)//'": must be one of ATM, LND, or OCN',&
228  fatal)
229  end select
230  ! convert area to m2
231  cellarea = cellarea*4.*pi*radius**2
232  case(version_2)
233  if (present(domain)) then
234  call mpp_get_compute_domain(domain,xsize=nlon,ysize=nlat)
235  else
236  call get_grid_size(component,tile,nlon,nlat)
237  endif
238  allocate(glonb(nlon+1,nlat+1),glatb(nlon+1,nlat+1))
239  call get_grid_cell_vertices(component, tile, glonb, glatb, domain)
240  if (great_circle_algorithm) then
241  call calc_mosaic_grid_great_circle_area(glonb*pi/180.0, glatb*pi/180.0, cellarea)
242  else
243  call calc_mosaic_grid_area(glonb*pi/180.0, glatb*pi/180.0, cellarea)
244  end if
245  deallocate(glonb,glatb)
246  end select
247 
248 end subroutine get_grid_cell_area_sg
249 
250 ! ============================================================================
251 ! get the area of the component per grid cell
252 ! ============================================================================
253 subroutine get_grid_comp_area_sg(component,tile,area,domain)
254  character(len=*) :: component
255  integer, intent(in) :: tile
256  real, intent(inout) :: area(:,:)
257  type(domain2d), intent(in), optional :: domain
258  ! local vars
259  integer :: n_xgrid_files ! number of exchange grid files in the mosaic
260  integer :: siz(4), nxgrid
261  integer :: i,j,m,n
262  integer, allocatable :: i1(:), j1(:), i2(:), j2(:)
263  real, allocatable :: xgrid_area(:)
264  real, allocatable :: rmask(:,:)
265  character(len=MAX_NAME) :: &
266  xgrid_name, & ! name of the variable holding xgrid names
267  tile_name, & ! name of the tile
268  xgrid_file, & ! name of the current xgrid file
269  mosaic_name,& ! name of the mosaic
270  mosaic_file,&
271  tilefile
272  character(len=4096) :: attvalue
273  character(len=MAX_NAME), allocatable :: nest_tile_name(:)
274  character(len=MAX_NAME) :: varname1, varname2
275  integer :: is,ie,js,je ! boundaries of our domain
276  integer :: i0, j0 ! offsets for x and y, respectively
277  integer :: num_nest_tile, ntiles
278  logical :: is_nest
279  integer :: found_xgrid_files ! how many xgrid files we actually found in the grid spec
280  integer :: ibegin, iend, bsize, l
281 
282  select case (get_grid_version())
283  case(version_0,version_1)
284  select case(component)
285  case('ATM')
286  call read_data(grid_file,'AREA_ATM',area, no_domain=.not.present(domain),domain=domain)
287  case('OCN')
288  allocate(rmask(size(area,1),size(area,2)))
289  call read_data(grid_file,'AREA_OCN',area, no_domain=.not.present(domain),domain=domain)
290  call read_data(grid_file,'wet', rmask,no_domain=.not.present(domain),domain=domain)
291  area = area*rmask
292  deallocate(rmask)
293  case('LND')
294  call read_data(grid_file,'AREA_LND',area,no_domain=.not.present(domain),domain=domain)
295  case default
296  call error_mesg(module_name//'/get_grid_comp_area',&
297  'Illegal component name "'//trim(component)//'": must be one of ATM, LND, or OCN',&
298  fatal)
299  end select
300  case(version_2) ! mosaic gridspec
301  select case (component)
302  case ('ATM')
303  ! just read the grid cell area and return
304  call get_grid_cell_area(component,tile,area)
305  return
306  case ('LND')
307  xgrid_name = 'aXl_file'
308  call read_data(grid_file, 'lnd_mosaic', mosaic_name)
309  tile_name = trim(mosaic_name)//'_tile'//char(tile+ichar('0'))
310  case ('OCN')
311  xgrid_name = 'aXo_file'
312  call read_data(grid_file, 'ocn_mosaic', mosaic_name)
313  tile_name = trim(mosaic_name)//'_tile'//char(tile+ichar('0'))
314  case default
315  call error_mesg(module_name//'/get_grid_comp_area',&
316  'Illegal component name "'//trim(component)//'": must be one of ATM, LND, or OCN',&
317  fatal)
318  end select
319  ! get the boundaries of the requested domain
320  if(present(domain)) then
321  call mpp_get_compute_domain(domain,is,ie,js,je)
322  i0 = 1-is ; j0=1-js
323  else
324  call get_grid_size(component,tile,ie,je)
325  is = 1 ; i0 = 0
326  js = 1 ; j0 = 0
327  endif
328  if (size(area,1)/=ie-is+1.or.size(area,2)/=je-js+1) &
329  call error_mesg(module_name//'/get_grid_comp_area',&
330  'size of the output argument "area" is not consistent with the domain',fatal)
331 
332  ! find the nest tile
333  call read_data(grid_file, 'atm_mosaic', mosaic_name)
334  call read_data(grid_file,'atm_mosaic_file',mosaic_file)
335  mosaic_file = grid_dir//trim(mosaic_file)
336  ntiles = get_mosaic_ntiles(trim(mosaic_file))
337  allocate(nest_tile_name(ntiles))
338  num_nest_tile = 0
339  do n = 1, ntiles
340  call read_data(mosaic_file, 'gridfiles', tilefile, level=n)
341  tilefile = grid_dir//trim(tilefile)
342  if( get_global_att_value(tilefile, "nest_grid", attvalue) ) then
343  if(trim(attvalue) == "TRUE") then
344  num_nest_tile = num_nest_tile + 1
345  nest_tile_name(num_nest_tile) = trim(mosaic_name)//'_tile'//char(n+ichar('0'))
346  else if(trim(attvalue) .NE. "FALSE") then
347  call error_mesg(module_name//'/get_grid_comp_area', 'value of global attribute nest_grid in file'// &
348  trim(tilefile)//' should be TRUE of FALSE', fatal)
349  endif
350  end if
351  end do
352  area(:,:) = 0.
353  if(field_exist(grid_file,xgrid_name)) then
354  ! get the number of the exchange-grid files
355  call field_size(grid_file,xgrid_name,siz)
356  n_xgrid_files = siz(2)
357  found_xgrid_files = 0
358  ! loop through all exchange grid files
359  do n = 1, n_xgrid_files
360  ! get the name of the current exchange grid file
361  call read_data(grid_file,xgrid_name,xgrid_file,level=n)
362  ! skip the rest of the loop if the name of the current tile isn't found
363  ! in the file name, but check this only if there is more than 1 tile
364  if(n_xgrid_files>1) then
365  if(index(xgrid_file,trim(tile_name))==0) cycle
366  endif
367  found_xgrid_files = found_xgrid_files + 1
368  !---make sure the atmosphere grid is not a nested grid
369  is_nest = .false.
370  do m = 1, num_nest_tile
371  if(index(xgrid_file, trim(nest_tile_name(m))) .NE. 0) then
372  is_nest = .true.
373  exit
374  end if
375  end do
376  if(is_nest) cycle
377 
378  ! finally read the exchange grid
379  nxgrid = get_mosaic_xgrid_size(grid_dir//xgrid_file)
380  if(nxgrid < bufsize) then
381  allocate(i1(nxgrid), j1(nxgrid), i2(nxgrid), j2(nxgrid), xgrid_area(nxgrid))
382  else
383  allocate(i1(bufsize), j1(bufsize), i2(bufsize), j2(bufsize), xgrid_area(bufsize))
384  endif
385  ibegin = 1
386  do l = 1,nxgrid,bufsize
387  bsize = min(bufsize, nxgrid-l+1)
388  iend = ibegin + bsize - 1
389  call get_mosaic_xgrid(grid_dir//xgrid_file, i1(1:bsize), j1(1:bsize), i2(1:bsize), j2(1:bsize), &
390  xgrid_area(1:bsize), ibegin, iend)
391  ! and sum the exchange grid areas
392  do m = 1, bsize
393  i = i2(m); j = j2(m)
394  if (i<is.or.i>ie) cycle
395  if (j<js.or.j>je) cycle
396  area(i+i0,j+j0) = area(i+i0,j+j0) + xgrid_area(m)
397  end do
398  ibegin = iend + 1
399  enddo
400  deallocate(i1, j1, i2, j2, xgrid_area)
401  enddo
402  if (found_xgrid_files == 0) &
403  call error_mesg('get_grid_comp_area', 'no xgrid files were found for component '&
404  //trim(component)//' (mosaic name is '//trim(mosaic_name)//')', fatal)
405 
406  endif
407  deallocate(nest_tile_name)
408  end select ! version
409  ! convert area to m2
410  area = area*4.*pi*radius**2
411 end subroutine get_grid_comp_area_sg
412 
413 !======================================================================
414 subroutine get_grid_cell_area_ug(component, tile, cellarea, SG_domain, UG_domain)
415  character(len=*), intent(in) :: component
416  integer , intent(in) :: tile
417  real , intent(inout) :: cellarea(:)
418  type(domain2d) , intent(in) :: SG_domain
419  type(domainUG) , intent(in) :: UG_domain
420  integer :: is, ie, js, je
421  real, allocatable :: SG_area(:,:)
422 
423  call mpp_get_compute_domain(sg_domain, is, ie, js, je)
424  allocate(sg_area(is:ie, js:je))
425  call get_grid_cell_area_sg(component, tile, sg_area, sg_domain)
426  call mpp_pass_sg_to_ug(ug_domain, sg_area, cellarea)
427  deallocate(sg_area)
428 
429 end subroutine get_grid_cell_area_ug
430 
431 subroutine get_grid_comp_area_ug(component, tile, area, SG_domain, UG_domain)
432  character(len=*), intent(in) :: component
433  integer , intent(in) :: tile
434  real , intent(inout) :: area(:)
435  type(domain2d) , intent(in) :: SG_domain
436  type(domainUG) , intent(in) :: UG_domain
437  integer :: is, ie, js, je
438  real, allocatable :: SG_area(:,:)
439 
440  call mpp_get_compute_domain(sg_domain, is, ie, js, je)
441  allocate(sg_area(is:ie, js:je))
442  call get_grid_comp_area_sg(component, tile, sg_area, sg_domain)
443  call mpp_pass_sg_to_ug(ug_domain, sg_area, area)
444  deallocate(sg_area)
445 
446 end subroutine get_grid_comp_area_ug
447 
448 
449 ! ============================================================================
450 ! returns arrays of global grid cell boundaries for given model component and
451 ! mosaic tile number.
452 ! NOTE that in case of non-lat-lon grid the returned coordinates may have be not so
453 ! meaningful, by the very nature of such grids. But presumably these 1D coordinate
454 ! arrays are good enough for diag axis and such.
455 ! ============================================================================
456 subroutine get_grid_cell_vertices_1d(component, tile, glonb, glatb)
457  character(len=*), intent(in) :: component
458  integer, intent(in) :: tile
459  real, intent(inout) :: glonb(:),glatb(:)
460 
461  integer :: nlon, nlat
462  integer :: start(4), nread(4)
463  real, allocatable :: tmp(:,:), x_vert_t(:,:,:), y_vert_t(:,:,:)
464  character(len=MAX_FILE) :: filename1, filename2
465 
466  call get_grid_size_for_one_tile(component, tile, nlon, nlat)
467  if (size(glonb(:))/=nlon+1) &
468  call error_mesg ( module_name//'/get_grid_cell_vertices_1D',&
469  'Size of argument "glonb" is not consistent with the grid size',fatal)
470  if (size(glatb(:))/=nlat+1) &
471  call error_mesg ( module_name//'/get_grid_cell_vertices_1D',&
472  'Size of argument "glatb" is not consistent with the grid size',fatal)
473  if(trim(component) .NE. 'ATM' .AND. component .NE. 'LND' .AND. component .NE. 'OCN') then
474  call error_mesg(module_name//'/get_grid_cell_vertices_1D',&
475  'Illegal component name "'//trim(component)//'": must be one of ATM, LND, or OCN',&
476  fatal)
477  endif
478 
479  select case(get_grid_version())
480  case(version_0)
481  select case(trim(component))
482  case('ATM','LND')
483  call read_data(grid_file, 'xb'//lowercase(component(1:1)), glonb, no_domain=.true.)
484  call read_data(grid_file, 'yb'//lowercase(component(1:1)), glatb, no_domain=.true.)
485  case('OCN')
486  call read_data(grid_file, "gridlon_vert_t", glonb, no_domain=.true.)
487  call read_data(grid_file, "gridlat_vert_t", glatb, no_domain=.true.)
488  end select
489  case(version_1)
490  select case(trim(component))
491  case('ATM','LND')
492  call read_data(grid_file, 'xb'//lowercase(component(1:1)), glonb, no_domain=.true.)
493  call read_data(grid_file, 'yb'//lowercase(component(1:1)), glatb, no_domain=.true.)
494  case('OCN')
495  allocate (x_vert_t(nlon,1,2), y_vert_t(1,nlat,2) )
496  start = 1; nread = 1
497  nread(1) = nlon; nread(2) = 1; start(3) = 1
498  call read_data(grid_file, "x_vert_T", x_vert_t(:,:,1), start, nread, no_domain=.true.)
499  nread(1) = nlon; nread(2) = 1; start(3) = 2
500  call read_data(grid_file, "x_vert_T", x_vert_t(:,:,2), start, nread, no_domain=.true.)
501 
502  nread(1) = 1; nread(2) = nlat; start(3) = 1
503  call read_data(grid_file, "y_vert_T", y_vert_t(:,:,1), start, nread, no_domain=.true.)
504  nread(1) = 1; nread(2) = nlat; start(3) = 4
505  call read_data(grid_file, "y_vert_T", y_vert_t(:,:,2), start, nread, no_domain=.true.)
506  glonb(1:nlon) = x_vert_t(1:nlon,1,1)
507  glonb(nlon+1) = x_vert_t(nlon,1,2)
508  glatb(1:nlat) = y_vert_t(1,1:nlat,1)
509  glatb(nlat+1) = y_vert_t(1,nlat,2)
510  deallocate(x_vert_t, y_vert_t)
511  end select
512  case(version_2)
513  ! get the name of the mosaic file for the component
514  call read_data(grid_file, trim(lowercase(component))//'_mosaic_file', filename1)
515  filename1=grid_dir//trim(filename1)
516  ! get the name of the grid file for the component and tile
517  call read_data(filename1, 'gridfiles', filename2, level=tile)
518  filename2 = grid_dir//trim(filename2)
519 
520  start = 1; nread = 1
521  nread(1) = 2*nlon+1
522  allocate( tmp(2*nlon+1,1) )
523  call read_data(filename2, "x", tmp, start, nread, no_domain=.true.)
524  glonb(1:nlon+1) = tmp(1:2*nlon+1:2,1)
525  deallocate(tmp)
526  allocate(tmp(1,2*nlat+1))
527 
528  start = 1; nread = 1
529  nread(2) = 2*nlat+1
530  call read_data(filename2, "y", tmp, start, nread, no_domain=.true.)
531  glatb(1:nlat+1) = tmp(1,1:2*nlat+1:2)
532  deallocate(tmp)
533  end select
534 
535 end subroutine get_grid_cell_vertices_1d
536 
537 ! ============================================================================
538 ! returns cell vertices for the specified model component and mosaic tile number
539 ! ============================================================================
540 subroutine get_grid_cell_vertices_2d(component, tile, lonb, latb, domain)
541  character(len=*), intent(in) :: component
542  integer, intent(in) :: tile
543  real, intent(inout) :: lonb(:,:),latb(:,:)
544  type(domain2d), optional, intent(in) :: domain
545 
546  ! local vars
547  character(len=MAX_FILE) :: filename1, filename2
548  integer :: nlon, nlat
549  integer :: i,j
550  real, allocatable :: buffer(:), tmp(:,:), x_vert_t(:,:,:), y_vert_t(:,:,:)
551  integer :: is,ie,js,je ! boundaries of our domain
552  integer :: i0,j0 ! offsets for coordinates
553  integer :: isg, jsg
554  integer :: start(4), nread(4)
555 
556  call get_grid_size_for_one_tile(component, tile, nlon, nlat)
557  if (present(domain)) then
558  call mpp_get_compute_domain(domain,is,ie,js,je)
559  else
560  is = 1 ; ie = nlon
561  js = 1 ; je = nlat
562  !--- domain normally should be present
563  call error_mesg ( module_name//'/get_grid_cell_vertices',&
564  'domain is not present, global data will be read', note)
565  endif
566  i0 = -is+1; j0 = -js+1
567 
568  ! verify that lonb and latb sizes are consistent with the size of domain
569  if (size(lonb,1)/=ie-is+2.or.size(lonb,2)/=je-js+2) &
570  call error_mesg ( module_name//'/get_grid_cell_vertices',&
571  'Size of argument "lonb" is not consistent with the domain size',fatal)
572  if (size(latb,1)/=ie-is+2.or.size(latb,2)/=je-js+2) &
573  call error_mesg ( module_name//'/get_grid_cell_vertices',&
574  'Size of argument "latb" is not consistent with the domain size',fatal)
575  if(trim(component) .NE. 'ATM' .AND. component .NE. 'LND' .AND. component .NE. 'OCN') then
576  call error_mesg(module_name//'/get_grid_cell_vertices',&
577  'Illegal component name "'//trim(component)//'": must be one of ATM, LND, or OCN',&
578  fatal)
579  endif
580 
581  select case(get_grid_version())
582  case(version_0)
583  select case(component)
584  case('ATM','LND')
585  allocate(buffer(max(nlon,nlat)+1))
586  ! read coordinates of grid cell vertices
587  call read_data(grid_file, 'xb'//lowercase(component(1:1)), buffer(1:nlon+1), no_domain=.true.)
588  do j = js, je+1
589  do i = is, ie+1
590  lonb(i+i0,j+j0) = buffer(i)
591  enddo
592  enddo
593  call read_data(grid_file, 'yb'//lowercase(component(1:1)), buffer(1:nlat+1), no_domain=.true.)
594  do j = js, je+1
595  do i = is, ie+1
596  latb(i+i0,j+j0) = buffer(j)
597  enddo
598  enddo
599  deallocate(buffer)
600  case('OCN')
601  if (present(domain)) then
602  start = 1; nread = 1
603  start(1) = is; start(2) = js
604  nread(1) = ie-is+2; nread(2) = je-js+2
605  call read_data(grid_file, 'geolon_vert_t', lonb, start, nread, no_domain=.true. )
606  call read_data(grid_file, 'geolat_vert_t', latb, start, nread, no_domain=.true. )
607  else
608  call read_data(grid_file, 'geolon_vert_t', lonb, no_domain=.true. )
609  call read_data(grid_file, 'geolat_vert_t', latb, no_domain=.true. )
610  endif
611  end select
612  case(version_1)
613  select case(component)
614  case('ATM','LND')
615  allocate(buffer(max(nlon,nlat)+1))
616  ! read coordinates of grid cell vertices
617  call read_data(grid_file, 'xb'//lowercase(component(1:1)), buffer(1:nlon+1), no_domain=.true.)
618  do j = js, je+1
619  do i = is, ie+1
620  lonb(i+i0,j+j0) = buffer(i)
621  enddo
622  enddo
623  call read_data(grid_file, 'yb'//lowercase(component(1:1)), buffer(1:nlat+1), no_domain=.true.)
624  do j = js, je+1
625  do i = is, ie+1
626  latb(i+i0,j+j0) = buffer(j)
627  enddo
628  enddo
629  deallocate(buffer)
630  case('OCN')
631  nlon=ie-is+1; nlat=je-js+1
632  allocate (x_vert_t(nlon,nlat,4), y_vert_t(nlon,nlat,4) )
633  call read_data(grid_file, 'x_vert_T', x_vert_t, no_domain=.not.present(domain), domain=domain )
634  call read_data(grid_file, 'y_vert_T', y_vert_t, no_domain=.not.present(domain), domain=domain )
635  lonb(1:nlon,1:nlat) = x_vert_t(1:nlon,1:nlat,1)
636  lonb(nlon+1,1:nlat) = x_vert_t(nlon,1:nlat,2)
637  lonb(1:nlon,nlat+1) = x_vert_t(1:nlon,nlat,4)
638  lonb(nlon+1,nlat+1) = x_vert_t(nlon,nlat,3)
639  latb(1:nlon,1:nlat) = y_vert_t(1:nlon,1:nlat,1)
640  latb(nlon+1,1:nlat) = y_vert_t(nlon,1:nlat,2)
641  latb(1:nlon,nlat+1) = y_vert_t(1:nlon,nlat,4)
642  latb(nlon+1,nlat+1) = y_vert_t(nlon,nlat,3)
643  deallocate(x_vert_t, y_vert_t)
644  end select
645  case(version_2)
646  ! get the name of the mosaic file for the component
647  call read_data(grid_file, trim(lowercase(component))//'_mosaic_file', filename1)
648  filename1=grid_dir//trim(filename1)
649  ! get the name of the grid file for the component and tile
650  call read_data(filename1, 'gridfiles', filename2, level=tile)
651  filename2 = grid_dir//trim(filename2)
652  if(PRESENT(domain)) then
653  call mpp_get_global_domain(domain, xbegin=isg, ybegin=jsg)
654  start = 1; nread = 1
655  start(1) = 2*(is-isg+1) - 1; nread(1) = 2*(ie-is)+3
656  start(2) = 2*(js-jsg+1) - 1; nread(2) = 2*(je-js)+3
657  allocate(tmp(nread(1), nread(2)) )
658  call read_data(filename2, 'x', tmp, start, nread, no_domain=.true.)
659  do j = 1, je-js+2
660  do i = 1, ie-is+2
661  lonb(i,j) = tmp(2*i-1,2*j-1)
662  enddo
663  enddo
664  call read_data(filename2, 'y', tmp, start, nread, no_domain=.true.)
665  do j = 1, je-js+2
666  do i = 1, ie-is+2
667  latb(i,j) = tmp(2*i-1,2*j-1)
668  enddo
669  enddo
670  else
671  allocate(tmp(2*nlon+1,2*nlat+1))
672  call read_data(filename2, 'x', tmp, no_domain=.true.)
673  do j = js, je+1
674  do i = is, ie+1
675  lonb(i+i0,j+j0) = tmp(2*i-1,2*j-1)
676  end do
677  end do
678  call read_data(filename2, 'y', tmp, no_domain=.true.)
679  do j = js, je+1
680  do i = is, ie+1
681  latb(i+i0,j+j0) = tmp(2*i-1,2*j-1)
682  end do
683  end do
684  endif
685  deallocate(tmp)
686  end select
687 
688 end subroutine get_grid_cell_vertices_2d
689 
690 
691 subroutine get_grid_cell_vertices_ug(component, tile, lonb, latb, SG_domain, UG_domain)
692  character(len=*), intent(in) :: component
693  integer, intent(in) :: tile
694  real, intent(inout) :: lonb(:,:),latb(:,:) ! The second dimension is 4
695  type(domain2d) , intent(in) :: SG_domain
696  type(domainUG) , intent(in) :: UG_domain
697  integer :: is, ie, js, je, i, j
698  real, allocatable :: SG_lonb(:,:), SG_latb(:,:), tmp(:,:,:)
699 
700  call mpp_get_compute_domain(sg_domain, is, ie, js, je)
701  allocate(sg_lonb(is:ie+1, js:je+1))
702  allocate(sg_latb(is:ie+1, js:je+1))
703  allocate(tmp(is:ie,js:je,4))
704  call get_grid_cell_vertices_2d(component, tile, sg_lonb, sg_latb, sg_domain)
705  do j = js, je
706  do i = is, ie
707  tmp(i,j,1) = sg_lonb(i,j)
708  tmp(i,j,2) = sg_lonb(i+1,j)
709  tmp(i,j,3) = sg_lonb(i+1,j+1)
710  tmp(i,j,4) = sg_lonb(i,j+1)
711  enddo
712  enddo
713  call mpp_pass_sg_to_ug(ug_domain, tmp, lonb)
714  do j = js, je
715  do i = is, ie
716  tmp(i,j,1) = sg_latb(i,j)
717  tmp(i,j,2) = sg_latb(i+1,j)
718  tmp(i,j,3) = sg_latb(i+1,j+1)
719  tmp(i,j,4) = sg_latb(i,j+1)
720  enddo
721  enddo
722  call mpp_pass_sg_to_ug(ug_domain, tmp, latb)
723 
724 
725  deallocate(sg_lonb, sg_latb, tmp)
726 
727 end subroutine get_grid_cell_vertices_ug
728 
729 ! ============================================================================
730 ! returns global coordinate arrays fro given model component and mosaic tile number
731 ! NOTE that in case of non-lat-lon grid those coordinates may have be not so
732 ! meaningful, by the very nature of such grids. But presumably these 1D coordinate
733 ! arrays are good enough for diag axis and such.
734 ! ============================================================================
735 subroutine get_grid_cell_centers_1d(component, tile, glon, glat)
736  character(len=*), intent(in) :: component
737  integer, intent(in) :: tile
738  real, intent(inout) :: glon(:),glat(:)
739  integer :: nlon, nlat
740  integer :: start(4), nread(4)
741  real, allocatable :: tmp(:,:)
742  character(len=MAX_FILE) :: filename1, filename2
743 
744  call get_grid_size_for_one_tile(component, tile, nlon, nlat)
745  if (size(glon(:))/=nlon) &
746  call error_mesg ( module_name//'/get_grid_cell_centers_1D',&
747  'Size of argument "glon" is not consistent with the grid size',fatal)
748  if (size(glat(:))/=nlat) &
749  call error_mesg ( module_name//'/get_grid_cell_centers_1D',&
750  'Size of argument "glat" is not consistent with the grid size',fatal)
751  if(trim(component) .NE. 'ATM' .AND. component .NE. 'LND' .AND. component .NE. 'OCN') then
752  call error_mesg(module_name//'/get_grid_cell_centers_1D',&
753  'Illegal component name "'//trim(component)//'": must be one of ATM, LND, or OCN',&
754  fatal)
755  endif
756 
757  select case(get_grid_version())
758  case(version_0)
759  select case(trim(component))
760  case('ATM','LND')
761  call read_data(grid_file, 'xt'//lowercase(component(1:1)), glon, no_domain=.true.)
762  call read_data(grid_file, 'yt'//lowercase(component(1:1)), glat, no_domain=.true.)
763  case('OCN')
764  call read_data(grid_file, "gridlon_t", glon, no_domain=.true.)
765  call read_data(grid_file, "gridlat_t", glat, no_domain=.true.)
766  end select
767  case(version_1)
768  select case(trim(component))
769  case('ATM','LND')
770  call read_data(grid_file, 'xt'//lowercase(component(1:1)), glon, no_domain=.true.)
771  call read_data(grid_file, 'yt'//lowercase(component(1:1)), glat, no_domain=.true.)
772  case('OCN')
773  call read_data(grid_file, "grid_x_T", glon, no_domain=.true.)
774  call read_data(grid_file, "grid_y_T", glat, no_domain=.true.)
775  end select
776  case(version_2)
777  ! get the name of the mosaic file for the component
778  call read_data(grid_file, trim(lowercase(component))//'_mosaic_file', filename1)
779  filename1=grid_dir//trim(filename1)
780  ! get the name of the grid file for the component and tile
781  call read_data(filename1, 'gridfiles', filename2, level=tile)
782  filename2 = grid_dir//trim(filename2)
783 
784  start = 1; nread = 1
785  nread(1) = 2*nlon+1; start(2) = 2
786  allocate( tmp(2*nlon+1,1) )
787  call read_data(filename2, "x", tmp, start, nread, no_domain=.true.)
788  glon(1:nlon) = tmp(2:2*nlon:2,1)
789  deallocate(tmp)
790  allocate(tmp(1, 2*nlat+1))
791 
792  start = 1; nread = 1
793  nread(2) = 2*nlat+1; start(1) = 2
794  call read_data(filename2, "y", tmp, start, nread, no_domain=.true.)
795  glat(1:nlat) = tmp(1,2:2*nlat:2)
796  deallocate(tmp)
797  end select
798 
799 
800 end subroutine get_grid_cell_centers_1d
801 
802 ! ============================================================================
803 ! returns grid cell centers for specified model component and mosaic tile number
804 ! ============================================================================
805 subroutine get_grid_cell_centers_2d(component, tile, lon, lat, domain)
806  character(len=*), intent(in) :: component
807  integer, intent(in) :: tile
808  real, intent(inout) :: lon(:,:),lat(:,:)
809  type(domain2d), intent(in), optional :: domain
810  ! local vars
811  character(len=MAX_NAME) :: varname
812  character(len=MAX_FILE) :: filename1, filename2
813  integer :: nlon, nlat
814  integer :: i,j
815  real, allocatable :: buffer(:),tmp(:,:)
816  integer :: is,ie,js,je ! boundaries of our domain
817  integer :: i0,j0 ! offsets for coordinates
818  integer :: isg, jsg
819  integer :: start(4), nread(4)
820 
821  call get_grid_size_for_one_tile(component, tile, nlon, nlat)
822  if (present(domain)) then
823  call mpp_get_compute_domain(domain,is,ie,js,je)
824  else
825  is = 1 ; ie = nlon
826  js = 1 ; je = nlat
827  !--- domain normally should be present
828  call error_mesg ( module_name//'/get_grid_cell_centers',&
829  'domain is not present, global data will be read', note)
830  endif
831  i0 = -is+1; j0 = -js+1
832 
833  ! verify that lon and lat sizes are consistent with the size of domain
834  if (size(lon,1)/=ie-is+1.or.size(lon,2)/=je-js+1) &
835  call error_mesg ( module_name//'/get_grid_cell_centers',&
836  'Size of array "lon" is not consistent with the domain size',&
837  fatal )
838  if (size(lat,1)/=ie-is+1.or.size(lat,2)/=je-js+1) &
839  call error_mesg ( module_name//'/get_grid_cell_centers',&
840  'Size of array "lat" is not consistent with the domain size',&
841  fatal )
842  if(trim(component) .NE. 'ATM' .AND. component .NE. 'LND' .AND. component .NE. 'OCN') then
843  call error_mesg(module_name//'/get_grid_cell_vertices',&
844  'Illegal component name "'//trim(component)//'": must be one of ATM, LND, or OCN',&
845  fatal)
846  endif
847 
848  select case(get_grid_version())
849  case(version_0)
850  select case (trim(component))
851  case('ATM','LND')
852  allocate(buffer(max(nlon,nlat)))
853  ! read coordinates of grid cell vertices
854  call read_data(grid_file, 'xt'//lowercase(component(1:1)), buffer(1:nlon), no_domain=.true.)
855  do j = js,je
856  do i = is,ie
857  lon(i+i0,j+j0) = buffer(i)
858  enddo
859  enddo
860  call read_data(grid_file, 'yt'//lowercase(component(1:1)), buffer(1:nlat), no_domain=.true.)
861  do j = js,je
862  do i = is,ie
863  lat(i+i0,j+j0) = buffer(j)
864  enddo
865  enddo
866  deallocate(buffer)
867  case('OCN')
868  call read_data(grid_file, 'geolon_t', lon, no_domain=.not.present(domain), domain=domain )
869  call read_data(grid_file, 'geolat_t', lat, no_domain=.not.present(domain), domain=domain )
870  end select
871  case(version_1)
872  select case(trim(component))
873  case('ATM','LND')
874  allocate(buffer(max(nlon,nlat)))
875  ! read coordinates of grid cell vertices
876  call read_data(grid_file, 'xt'//lowercase(component(1:1)), buffer(1:nlon), no_domain=.true.)
877  do j = js,je
878  do i = is,ie
879  lon(i+i0,j+j0) = buffer(i)
880  enddo
881  enddo
882  call read_data(grid_file, 'yt'//lowercase(component(1:1)), buffer(1:nlat), no_domain=.true.)
883  do j = js,je
884  do i = is,ie
885  lat(i+i0,j+j0) = buffer(j)
886  enddo
887  enddo
888  deallocate(buffer)
889  case('OCN')
890  call read_data(grid_file, 'x_T', lon, no_domain=.not.present(domain), domain=domain )
891  call read_data(grid_file, 'y_T', lat, no_domain=.not.present(domain), domain=domain )
892  end select
893  case(version_2) ! mosaic grid file
894  ! get the name of the mosaic file for the component
895  call read_data(grid_file, trim(lowercase(component))//'_mosaic_file', filename1)
896  filename1=grid_dir//trim(filename1)
897  ! get the name of the grid file for the component and tile
898  call read_data(filename1, 'gridfiles', filename2, level=tile)
899  filename2 = grid_dir//trim(filename2)
900  if(PRESENT(domain)) then
901  call mpp_get_global_domain(domain, xbegin=isg, ybegin=jsg)
902  start = 1; nread = 1
903  start(1) = 2*(is-isg+1) - 1; nread(1) = 2*(ie-is)+3
904  start(2) = 2*(js-jsg+1) - 1; nread(2) = 2*(je-js)+3
905  allocate(tmp(nread(1), nread(2)))
906  call read_data(filename2, 'x', tmp, start, nread, no_domain=.true.)
907  do j = 1, je-js+1
908  do i = 1, ie-is+1
909  lon(i,j) = tmp(2*i,2*j)
910  enddo
911  enddo
912  call read_data(filename2, 'y', tmp, start, nread, no_domain=.true.)
913  do j = 1, je-js+1
914  do i = 1, ie-is+1
915  lat(i,j) = tmp(2*i,2*j)
916  enddo
917  enddo
918  else
919  allocate(tmp(2*nlon+1,2*nlat+1))
920  call read_data(filename2, 'x', tmp, no_domain=.true.)
921  do j = js,je
922  do i = is,ie
923  lon(i+i0,j+j0) = tmp(2*i,2*j)
924  end do
925  end do
926  call read_data(filename2, 'y', tmp, no_domain=.true.)
927  do j = js,je
928  do i = is,ie
929  lat(i+i0,j+j0) = tmp(2*i,2*j)
930  end do
931  end do
932  deallocate(tmp)
933  endif
934  end select
935 
936 end subroutine get_grid_cell_centers_2d
937 
938 subroutine get_grid_cell_centers_ug(component, tile, lon, lat, SG_domain, UG_domain)
939  character(len=*), intent(in) :: component
940  integer, intent(in) :: tile
941  real, intent(inout) :: lon(:),lat(:)
942  type(domain2d) , intent(in) :: SG_domain
943  type(domainUG) , intent(in) :: UG_domain
944  integer :: is, ie, js, je
945  real, allocatable :: SG_lon(:,:), SG_lat(:,:)
946 
947  call mpp_get_compute_domain(sg_domain, is, ie, js, je)
948  allocate(sg_lon(is:ie, js:je))
949  allocate(sg_lat(is:ie, js:je))
950  call get_grid_cell_centers_2d(component, tile, sg_lon, sg_lat, sg_domain)
951  call mpp_pass_sg_to_ug(ug_domain, sg_lon, lon)
952  call mpp_pass_sg_to_ug(ug_domain, sg_lat, lat)
953  deallocate(sg_lon, sg_lat)
954 
955 end subroutine get_grid_cell_centers_ug
956 
957 ! ============================================================================
958 ! given a model component, a layout, and (optionally) a halo size, returns a
959 ! domain for current processor
960 ! ============================================================================
961 ! this subroutine probably does not belong in the grid_mod
962 subroutine define_cube_mosaic ( component, domain, layout, halo, maskmap )
963  character(len=*) , intent(in) :: component
964  type(domain2d) , intent(inout) :: domain
965  integer , intent(in) :: layout(2)
966  integer, optional, intent(in) :: halo
967  logical, optional, intent(in) :: maskmap(:,:,:)
968 
969  ! ---- local constants
970 
971  ! ---- local vars
972  character(len=MAX_NAME) :: varname
973  character(len=MAX_FILE) :: mosaic_file
974  integer :: ntiles ! number of tiles
975  integer :: ncontacts ! number of contacts between mosaic tiles
976  integer :: n
977  integer :: ng, pe_pos, npes ! halo size
978  integer, allocatable :: nlon(:), nlat(:), global_indices(:,:)
979  integer, allocatable :: pe_start(:), pe_end(:), layout_2d(:,:)
980  integer, allocatable :: tile1(:),tile2(:)
981  integer, allocatable :: is1(:),ie1(:),js1(:),je1(:)
982  integer, allocatable :: is2(:),ie2(:),js2(:),je2(:)
983 
984  call get_grid_ntiles(component,ntiles)
985  allocate(nlon(ntiles), nlat(ntiles))
986  allocate(global_indices(4,ntiles))
987  allocate(pe_start(ntiles),pe_end(ntiles))
988  allocate(layout_2d(2,ntiles))
989  call get_grid_size(component,nlon,nlat)
990 
991  pe_pos = mpp_root_pe()
992  do n = 1, ntiles
993  global_indices(:,n) = (/ 1, nlon(n), 1, nlat(n) /)
994  layout_2d(:,n) = layout
995  if(present(maskmap)) then
996  npes = count(maskmap(:,:,n))
997  else
998  npes = layout(1)*layout(2)
999  endif
1000  pe_start(n) = pe_pos
1001  pe_end(n) = pe_pos + npes - 1
1002  pe_pos = pe_end(n) + 1
1003  enddo
1004 
1005  varname=trim(lowercase(component))//'_mosaic_file'
1006  call read_data(grid_file,varname,mosaic_file)
1007  mosaic_file = grid_dir//mosaic_file
1008 
1009  ! get the contact information from mosaic file
1010  ncontacts = get_mosaic_ncontacts(mosaic_file)
1011  allocate(tile1(ncontacts),tile2(ncontacts))
1012  allocate(is1(ncontacts),ie1(ncontacts),js1(ncontacts),je1(ncontacts))
1013  allocate(is2(ncontacts),ie2(ncontacts),js2(ncontacts),je2(ncontacts))
1014  call get_mosaic_contact(mosaic_file, tile1, tile2, &
1015  is1, ie1, js1, je1, is2, ie2, js2, je2)
1016 
1017  ng = 0
1018  if(present(halo)) ng = halo
1019  ! create the domain2d variable
1020  call mpp_define_mosaic ( global_indices, layout_2d, domain, &
1021  ntiles, ncontacts, tile1, tile2, &
1022  is1, ie1, js1, je1, &
1023  is2, ie2, js2, je2, &
1024  pe_start=pe_start, pe_end=pe_end, symmetry=.true., &
1025  shalo = ng, nhalo = ng, whalo = ng, ehalo = ng, &
1026  maskmap = maskmap, &
1027  name = trim(component)//'Cubic-Sphere Grid' )
1028 
1029  deallocate(nlon,nlat,global_indices,pe_start,pe_end,layout_2d)
1030  deallocate(tile1,tile2)
1031  deallocate(is1,ie1,js1,je1)
1032  deallocate(is2,ie2,js2,je2)
1033 
1034 end subroutine define_cube_mosaic
1035 
1036 end module grid_mod
Definition: fms.F90:20
subroutine get_grid_cell_vertices_ug(component, tile, lonb, latb, SG_domain, UG_domain)
Definition: grid.F90:692
subroutine, public get_mosaic_contact(mosaic_file, tile1, tile2, istart1, iend1, jstart1, jend1, istart2, iend2, jstart2, jend2)
Definition: mosaic.F90:352
logical function, public get_great_circle_algorithm()
Definition: fms_io.F90:8566
real, parameter, public radius
Radius of the Earth [m].
Definition: constants.F90:72
integer, parameter version_0
Definition: grid.F90:89
subroutine, public define_cube_mosaic(component, domain, layout, halo, maskmap)
Definition: grid.F90:963
character(len= *), parameter grid_dir
Definition: grid.F90:85
subroutine get_grid_size_for_one_tile(component, tile, nx, ny)
Definition: grid.F90:182
integer, parameter bufsize
Definition: grid.F90:96
logical first_call
Definition: grid.F90:102
subroutine, public calc_mosaic_grid_great_circle_area(lon, lat, area)
Definition: mosaic.F90:553
character(len= *), parameter grid_file
Definition: grid.F90:85
subroutine get_grid_cell_vertices_2d(component, tile, lonb, latb, domain)
Definition: grid.F90:541
subroutine get_grid_cell_area_ug(component, tile, cellarea, SG_domain, UG_domain)
Definition: grid.F90:415
subroutine get_grid_cell_area_sg(component, tile, cellarea, domain)
Definition: grid.F90:207
integer function, public get_mosaic_ntiles(mosaic_file)
Definition: mosaic.F90:214
integer, parameter max_file
Definition: grid.F90:89
Definition: mpp.F90:39
integer grid_version
Definition: grid.F90:100
subroutine get_grid_cell_centers_1d(component, tile, glon, glat)
Definition: grid.F90:736
subroutine, public get_mosaic_grid_sizes(mosaic_file, nx, ny)
Definition: mosaic.F90:279
integer, parameter version_2
Definition: grid.F90:89
subroutine, public get_mosaic_xgrid(xgrid_file, i1, j1, i2, j2, area, ibegin, iend)
Definition: mosaic.F90:150
subroutine, public calc_mosaic_grid_area(lon, lat, area)
Definition: mosaic.F90:513
subroutine get_grid_comp_area_sg(component, tile, area, domain)
Definition: grid.F90:254
integer, parameter max_name
Definition: grid.F90:89
int field_exist(const char *file, const char *name)
Definition: read_mosaic.c:69
logical great_circle_algorithm
Definition: grid.F90:101
integer, parameter version_1
Definition: grid.F90:89
subroutine, public get_grid_ntiles(component, ntiles)
Definition: grid.F90:136
integer function, public get_mosaic_ncontacts(mosaic_file)
Definition: mosaic.F90:239
character(len= *), parameter module_name
Definition: grid.F90:79
subroutine get_grid_comp_area_ug(component, tile, area, SG_domain, UG_domain)
Definition: grid.F90:432
subroutine get_grid_cell_centers_2d(component, tile, lon, lat, domain)
Definition: grid.F90:806
#define max(a, b)
Definition: mosaic_util.h:33
integer function get_grid_version()
Definition: grid.F90:108
integer function, public get_mosaic_xgrid_size(xgrid_file)
Definition: mosaic.F90:114
subroutine get_grid_size_for_all_tiles(component, nx, ny)
Definition: grid.F90:156
subroutine get_grid_cell_vertices_1d(component, tile, glonb, glatb)
Definition: grid.F90:457
#define min(a, b)
Definition: mosaic_util.h:32
subroutine get_grid_cell_centers_ug(component, tile, lon, lat, SG_domain, UG_domain)
Definition: grid.F90:939
subroutine, public error_mesg(routine, message, level)
Definition: fms.F90:529
real(fp), parameter, public pi