FMS  2025.04
Flexible Modeling System
mosaic2.F90
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 !> @defgroup mosaic2_mod mosaic2_mod
19 !> @ingroup mosaic2
20 !> @brief Implements some utility routines to read mosaic information.
21 !!
22 !> @author Zhi Liang
23 !!
24 !! Implements some utility routines to read mosaic information.
25 !! The information includes number of tiles and contacts in the mosaic,
26 !! mosaic grid resolution of each tile, mosaic contact information, mosaic exchange
27 !! grid information. Each routine will call a C-version routine to get these information.
28 
29 !> @addtogroup mosaic2_mod
30 !> @{
31 module mosaic2_mod
32 
33 !use fms_mod, only : write_version_number
34 
35 use mpp_mod, only : mpp_error, fatal, mpp_pe, mpp_root_pe
36 use mpp_domains_mod, only : domain2d, mpp_get_current_ntile, mpp_get_tile_id
37 use constants_mod, only : pi, radius
38 use fms2_io_mod, only : fmsnetcdffile_t, open_file, close_file, get_dimension_size
39 use fms2_io_mod, only : read_data, variable_exists
40 use platform_mod, only : r4_kind, r8_kind, fms_path_len
41 
42 implicit none
43 private
44 
45 character(len=*), parameter :: &
46  grid_dir = 'INPUT/' !> root directory for all grid files
47 
48 integer, parameter :: &
49  MAX_NAME = 256, & !> max length of the variable names
50  x_refine = 2, & !> supergrid size/model grid size in x-direction
51  y_refine = 2 !> supergrid size/model grid size in y-direction
52 
53 ! --- public interface
54 
55 
56 public :: get_mosaic_ntiles
57 public :: get_mosaic_ncontacts
58 public :: get_mosaic_grid_sizes
59 public :: get_mosaic_contact
60 public :: get_mosaic_xgrid_size
61 public :: get_mosaic_xgrid
62 public :: get_mosaic_tile_grid
63 public :: calc_mosaic_grid_area
65 public :: is_inside_polygon
66 
67 
69  module procedure get_mosaic_xgrid_r4
70  module procedure get_mosaic_xgrid_r8
71 end interface get_mosaic_xgrid
72 
74  module procedure calc_mosaic_grid_area_r4
75  module procedure calc_mosaic_grid_area_r8
76 end interface calc_mosaic_grid_area
77 
79  module procedure calc_mosaic_grid_great_circle_area_r4
80  module procedure calc_mosaic_grid_great_circle_area_r8
82 
84  module procedure is_inside_polygon_r4
85  module procedure is_inside_polygon_r8
86 end interface is_inside_polygon
87 
88 
89 logical :: module_is_initialized = .true.
90 ! Include variable "version" to be written to log file.
91 #include<file_version.h>
92 
93 contains
94 
95 !#######################################################################
96 
97 !> @brief Initialize the mosaic_mod.
98 !> Initialization routine for the mosaic module. It writes the
99 !! version information to the log file.
100 !! <br>Example usage:
101 !! call mosaic_init ( )
102 subroutine mosaic_init()
103 
104  if (module_is_initialized) return
105  module_is_initialized = .true.
106 
107 !--------- write version number and namelist ------------------
108 ! call write_version_number("MOSAIC_MOD", version)
109 
110 end subroutine mosaic_init
111 
112 !#######################################################################
113 !> @brief return exchange grid size of mosaic xgrid file.
114 !!
115 !> <br>Example usage:
116 !! nxgrid = get_mosaic_xgrid_size(xgrid_file)
117  function get_mosaic_xgrid_size(fileobj)
118  type(fmsnetcdffile_t), intent(in) :: fileobj !> File that contains exchange grid information.
119  integer :: get_mosaic_xgrid_size
120 
121  call get_dimension_size(fileobj, "ncells", get_mosaic_xgrid_size)
122 
123  return
124 
125  end function get_mosaic_xgrid_size
126  !###############################################################################
127  !> Get number of tiles in the mosaic_file.
128  !> @param fileobj mosaic file object
129  !> @return Number of tiles in given file
130  !> <br>Example usage:
131  !! ntiles = get_mosaic_ntiles( mosaic_file)
132  !!
133  function get_mosaic_ntiles(fileobj)
134  type(fmsnetcdffile_t), intent(in) :: fileobj
135  integer :: get_mosaic_ntiles
136 
137  call get_dimension_size(fileobj, "ntiles", get_mosaic_ntiles)
138 
139  return
140 
141  end function get_mosaic_ntiles
142 
143  !###############################################################################
144  !> Get number of contacts in the mosaic_file.
145  !> @param fileobj mosaic file object
146  !> @return number of contacts in a given file
147  !> <br>Example usage:
148  !! ntiles = get_mosaic_ncontacts( mosaic_file)
149  function get_mosaic_ncontacts(fileobj)
150  type(fmsnetcdffile_t), intent(in) :: fileobj
151  integer :: get_mosaic_ncontacts
152 
153  if(variable_exists(fileobj, "contacts") ) then
154  call get_dimension_size(fileobj, "ncontact", get_mosaic_ncontacts)
155  else
157  endif
158 
159  return
160 
161  end function get_mosaic_ncontacts
162 
163 
164  !###############################################################################
165  !> Get grid size of each tile from mosaic_file
166  !> @param fileobj mosaic file object
167  !> @param[inout] nx List of grid size in x-direction of each tile
168  !> @param[inout] ny List of grid size in y-direction of each tile
169  !> <br>Example usage:
170  !! call get_mosaic_grid_sizes(mosaic_file, nx, ny)
171  subroutine get_mosaic_grid_sizes( fileobj, nx, ny)
172  type(fmsnetcdffile_t), intent(in) :: fileobj
173  integer, dimension(:), intent(inout) :: nx, ny
174 
175  character(len=FMS_PATH_LEN) :: gridfile
176  integer :: ntiles, n
177  type(fmsnetcdffile_t) :: gridobj
178 
179  ntiles = get_mosaic_ntiles(fileobj)
180  if(ntiles .NE. size(nx(:)) .OR. ntiles .NE. size(ny(:)) ) then
181  call mpp_error(fatal, "get_mosaic_grid_sizes: size of nx/ny does not equal to ntiles")
182  endif
183 
184  do n = 1, ntiles
185  call read_data(fileobj, 'gridfiles', gridfile, corner=n)
186  gridfile = grid_dir//trim(gridfile)
187 
188  if(.not. open_file(gridobj, gridfile, 'read')) then
189  call mpp_error(fatal, 'mosaic_mod(get_mosaic_grid_sizes):Error in opening file '//trim(gridfile))
190  endif
191  call get_dimension_size(gridobj, "nx", nx(n))
192  call get_dimension_size(gridobj, "ny", ny(n))
193  call close_file(gridobj)
194  if(mod(nx(n),x_refine) .NE. 0) call mpp_error(fatal, "get_mosaic_grid_sizes: nx is not divided by x_refine");
195  if(mod(ny(n),y_refine) .NE. 0) call mpp_error(fatal, "get_mosaic_grid_sizes: ny is not divided by y_refine");
196  nx(n) = nx(n)/x_refine;
197  ny(n) = ny(n)/y_refine;
198  enddo
199 
200  return
201 
202  end subroutine get_mosaic_grid_sizes
203 
204  !###############################################################################
205  !> Get contact information from mosaic_file
206  !> <br>Example usage:
207  !! call get_mosaic_contact(mosaic_file, tile1, tile2, istart1, iend1, jstart1, jend1,
208  !! istart2, iend2, jstart2, jend2)
209  !> @param fileobj mosaic file object
210  !> @param[inout] tile1 list of tile numbers in tile 1 of each contact
211  !> @param[inout] tile2 list of tile numbers in tile 2 of each contact
212  !> @param[inout] istart1 list of starting i-index in tile 1 of each contact
213  !> @param[inout] iend1 list of ending i-index in tile 1 of each contact
214  !> @param[inout] jstart1 list of starting j-index in tile 1 of each contact
215  !> @param[inout] jend1 list of ending j-index in tile 1 of each contact
216  !> @param[inout] istart2 list of starting i-index in tile 2 of each contact
217  !> @param[inout] iend2 list of ending i-index in tile 2 of each contact
218  !> @param[inout] jstart2 list of starting j-index in tile 2 of each contact
219  !> @param[inout] jend2 list of ending j-index in tile 2 of each contact
220  subroutine get_mosaic_contact( fileobj, tile1, tile2, istart1, iend1, jstart1, jend1, &
221  istart2, iend2, jstart2, jend2)
222  type(fmsnetcdffile_t), intent(in) :: fileobj
223  integer, dimension(:), intent(inout) :: tile1, tile2
224  integer, dimension(:), intent(inout) :: istart1, iend1, jstart1, jend1
225  integer, dimension(:), intent(inout) :: istart2, iend2, jstart2, jend2
226  character(len=MAX_NAME), allocatable :: gridtiles(:)
227  character(len=MAX_NAME), allocatable :: contacts(:)
228  character(len=MAX_NAME), allocatable :: contacts_index(:)
229  character(len=MAX_NAME) :: strlist(8)
230  integer :: ntiles, n, m, ncontacts, nstr, ios
231  integer :: i1_type, j1_type, i2_type, j2_type
232  logical :: found
233 
234  ntiles = get_mosaic_ntiles(fileobj)
235  allocate(gridtiles(ntiles))
236  if(mpp_pe()==mpp_root_pe()) then
237  do n = 1, ntiles
238  do m = 1,max_name
239  gridtiles(n)(m:m) = " "
240  enddo
241  enddo
242  endif
243  call read_data(fileobj, 'gridtiles', gridtiles)
244 
245  ncontacts = get_mosaic_ncontacts(fileobj)
246  if(ncontacts>0) then
247  allocate(contacts(ncontacts), contacts_index(ncontacts))
248  if(mpp_pe()==mpp_root_pe()) then
249  do n = 1, ncontacts
250  do m = 1,max_name
251  contacts(n)(m:m) = " "
252  contacts_index(n)(m:m) = " "
253  enddo
254  enddo
255  endif
256  call read_data(fileobj, "contacts", contacts)
257  call read_data(fileobj, "contact_index", contacts_index)
258  endif
259  do n = 1, ncontacts
260  nstr = parse_string(contacts(n), ":", strlist)
261  if(nstr .NE. 4) call mpp_error(fatal, &
262  "mosaic_mod(get_mosaic_contact): number of elements in contact separated by :/:: should be 4")
263  found = .false.
264  do m = 1, ntiles
265  if(trim(gridtiles(m)) == trim(strlist(2)) ) then !found the tile name
266  found = .true.
267  tile1(n) = m
268  exit
269  endif
270  enddo
271 
272  if(.not.found) call mpp_error(fatal, &
273  "mosaic_mod(get_mosaic_contact):the first tile name specified in contact is not found in tile list")
274 
275  found = .false.
276  do m = 1, ntiles
277  if(trim(gridtiles(m)) == trim(strlist(4)) ) then !found the tile name
278  found = .true.
279  tile2(n) = m
280  exit
281  endif
282  enddo
283 
284  if(.not.found) call mpp_error(fatal, &
285  "mosaic_mod(get_mosaic_contact):the second tile name specified in contact is not found in tile list")
286 
287  nstr = parse_string(contacts_index(n), ":,", strlist)
288  if(nstr .NE. 8) then
289  if(mpp_pe()==mpp_root_pe()) then
290  print*, "nstr is ", nstr
291  print*, "contacts is ", contacts_index(n)
292  do m = 1, nstr
293  print*, "strlist is ", trim(strlist(m))
294  enddo
295  endif
296  call mpp_error(fatal, &
297  "mosaic_mod(get_mosaic_contact): number of elements in contact_index separated by :/, should be 8")
298  endif
299  read(strlist(1), *, iostat=ios) istart1(n)
300  if(ios .NE. 0) call mpp_error(fatal, &
301  "mosaic_mod(get_mosaic_contact): Error in reading istart1")
302  read(strlist(2), *, iostat=ios) iend1(n)
303  if(ios .NE. 0) call mpp_error(fatal, &
304  "mosaic_mod(get_mosaic_contact): Error in reading iend1")
305  read(strlist(3), *, iostat=ios) jstart1(n)
306  if(ios .NE. 0) call mpp_error(fatal, &
307  "mosaic_mod(get_mosaic_contact): Error in reading jstart1")
308  read(strlist(4), *, iostat=ios) jend1(n)
309  if(ios .NE. 0) call mpp_error(fatal, &
310  "mosaic_mod(get_mosaic_contact): Error in reading jend1")
311  read(strlist(5), *, iostat=ios) istart2(n)
312  if(ios .NE. 0) call mpp_error(fatal, &
313  "mosaic_mod(get_mosaic_contact): Error in reading istart2")
314  read(strlist(6), *, iostat=ios) iend2(n)
315  if(ios .NE. 0) call mpp_error(fatal, &
316  "mosaic_mod(get_mosaic_contact): Error in reading iend2")
317  read(strlist(7), *, iostat=ios) jstart2(n)
318  if(ios .NE. 0) call mpp_error(fatal, &
319  "mosaic_mod(get_mosaic_contact): Error in reading jstart2")
320  read(strlist(8), *, iostat=ios) jend2(n)
321  if(ios .NE. 0) call mpp_error(fatal, &
322  "mosaic_mod(get_mosaic_contact): Error in reading jend2")
323 
324  i1_type = transfer_to_model_index(istart1(n), iend1(n), x_refine)
325  j1_type = transfer_to_model_index(jstart1(n), jend1(n), y_refine)
326  i2_type = transfer_to_model_index(istart2(n), iend2(n), x_refine)
327  j2_type = transfer_to_model_index(jstart2(n), jend2(n), y_refine)
328 
329  if( i1_type == 0 .AND. j1_type == 0 ) call mpp_error(fatal, &
330  "mosaic_mod(get_mosaic_contact): istart1==iend1 and jstart1==jend1")
331  if( i2_type == 0 .AND. j2_type == 0 ) call mpp_error(fatal, &
332  "mosaic_mod(get_mosaic_contact): istart2==iend2 and jstart2==jend2")
333  if( i1_type + j1_type .NE. i2_type + j2_type ) call mpp_error(fatal, &
334  "mosaic_mod(get_mosaic_contact): It is not a line or overlap contact")
335 
336  enddo
337 
338  deallocate(gridtiles)
339  if(ncontacts>0) deallocate(contacts, contacts_index)
340 
341  end subroutine get_mosaic_contact
342 
343 
344 function transfer_to_model_index(istart, iend, refine_ratio)
345  integer, intent(inout) :: istart, iend
346  integer :: refine_ratio
347  integer :: transfer_to_model_index
348  integer :: istart_in, iend_in
349 
350  istart_in = istart
351  iend_in = iend
352 
353  if( istart_in == iend_in ) then
354  transfer_to_model_index = 0
355  istart = (istart_in + 1)/refine_ratio
356  iend = istart
357  else
358  transfer_to_model_index = 1
359  if( iend_in > istart_in ) then
360  istart = istart_in + 1
361  iend = iend_in
362  else
363  istart = istart_in
364  iend = iend_in + 1
365  endif
366  if( mod(istart, refine_ratio) .NE. 0 .OR. mod(iend,refine_ratio) .NE. 0) call mpp_error(fatal, &
367  "mosaic_mod(transfer_to_model_index): mismatch between refine_ratio and istart/iend")
368  istart = istart/refine_ratio
369  iend = iend/refine_ratio
370 
371  endif
372 
373  return
374 
375 end function transfer_to_model_index
376 !#####################################################################
377 function parse_string(string, set, sval)
378  character(len=*), intent(in) :: string
379  character(len=*), intent(in) :: set
380  character(len=*), intent(out) :: sval(:)
381  integer :: parse_string
382  integer :: nelem, length, first, last
383 
384  nelem = size(sval(:))
385  length = len_trim(string)
386 
387  first = 1; last = 0
388  parse_string = 0
389 
390  do while(first .LE. length)
391  parse_string = parse_string + 1
392  if(parse_string>nelem) then
393  call mpp_error(fatal, "mosaic_mod(parse_string) : number of element is greater than size(sval(:))")
394  endif
395  last = first - 1 + scan(string(first:length), set)
396  if(last == first-1 ) then ! not found, end of string
397  sval(parse_string) = string(first:length)
398  exit
399  else
400  if(last <= first) then
401  call mpp_error(fatal, "mosaic_mod(parse_string) : last <= first")
402  endif
403  sval(parse_string) = string(first:(last-1))
404  first = last + 1
405  ! scan to make sure the next is not the character in the set
406  do while (first == last+1)
407  last = first - 1 + scan(string(first:length), set)
408  if(last == first) then
409  first = first+1
410  else
411  exit
412  endif
413  end do
414  endif
415  enddo
416 
417  return
418 
419 end function parse_string
420 
421  !#############################################################################
422  !> Gets the name of a mosaic tile grid file
423  !> @param[out] grid_file name of grid file
424  !> @param fileobj mosaic file object
425  !> @param domain current domain
426  !> @param tile_count optional count of tiles
427 subroutine get_mosaic_tile_grid(grid_file, fileobj, domain, tile_count)
428  character(len=*), intent(out) :: grid_file
429  type(fmsnetcdffile_t), intent(in) :: fileobj
430  type(domain2d), intent(in) :: domain
431  integer, intent(in), optional :: tile_count
432  integer :: tile, ntileme
433  integer, dimension(:), allocatable :: tile_id
434  character(len=256), allocatable :: filelist(:)
435  integer :: ntiles
436 
437  ntiles = get_mosaic_ntiles(fileobj)
438  allocate(filelist(ntiles))
439  tile = 1
440  if(present(tile_count)) tile = tile_count
441  ntileme = mpp_get_current_ntile(domain)
442  allocate(tile_id(ntileme))
443  tile_id = mpp_get_tile_id(domain)
444  call read_data(fileobj, "gridfiles", filelist)
445  grid_file = 'INPUT/'//trim(filelist(tile_id(tile)))
446  deallocate(tile_id, filelist)
447 
448 end subroutine get_mosaic_tile_grid
449 
450 #include "mosaic2_r4.fh"
451 #include "mosaic2_r8.fh"
452 
453 end module mosaic2_mod
454 !> @}
455 ! close documentation grouping
Close a netcdf or domain file opened with open_file or open_virtual_file.
Definition: fms2_io.F90:165
Opens a given netcdf or domain file.
Definition: fms2_io.F90:121
Read data from a defined field in a file.
Definition: fms2_io.F90:291
integer function, public get_mosaic_xgrid_size(fileobj)
return exchange grid size of mosaic xgrid file.
Definition: mosaic2.F90:118
subroutine, public get_mosaic_contact(fileobj, tile1, tile2, istart1, iend1, jstart1, jend1, istart2, iend2, jstart2, jend2)
Get contact information from mosaic_file Example usage: call get_mosaic_contact(mosaic_file,...
Definition: mosaic2.F90:222
integer function, public get_mosaic_ntiles(fileobj)
Get number of tiles in the mosaic_file.
Definition: mosaic2.F90:134
subroutine, public get_mosaic_tile_grid(grid_file, fileobj, domain, tile_count)
Gets the name of a mosaic tile grid file.
Definition: mosaic2.F90:428
subroutine mosaic_init()
Initialize the mosaic_mod. Initialization routine for the mosaic module. It writes the version inform...
Definition: mosaic2.F90:103
integer function, public get_mosaic_ncontacts(fileobj)
Get number of contacts in the mosaic_file.
Definition: mosaic2.F90:150
subroutine, public get_mosaic_grid_sizes(fileobj, nx, ny)
Get grid size of each tile from mosaic_file.
Definition: mosaic2.F90:172
integer function, dimension(size(domain%tile_id(:))) mpp_get_tile_id(domain)
Returns the tile_id on current pe.
integer function mpp_get_current_ntile(domain)
Returns number of tile on current pe.
The domain2D type contains all the necessary information to define the global, compute and data domai...
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:406
Error handler.
Definition: mpp.F90:381