FV3 Bundle
mosaic.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 mosaic_mod
20 
21 ! <CONTACT EMAIL="Zhi.Liang@noaa.gov">
22 ! Zhi Liang
23 ! </CONTACT>
24 
25 ! <HISTORY SRC="http://www.gfdl.noaa.gov/fms-cgi-bin/cvsweb.cgi/FMS/"/>
26 
27 ! <OVERVIEW>
28 ! <TT>mosaic_mod</TT> implements some utility routines to read mosaic information.
29 ! </OVERVIEW>
30 
31 ! <DESCRIPTION>
32 ! <TT>mosaic_mod</TT> implements some utility routines to read mosaic information.
33 ! The information includes number of tiles and contacts in the mosaic,
34 ! mosaic grid resolution of each tile, mosaic contact information, mosaic exchange
35 ! grid information. Each routine will call a C-version routine to get these information.
36 ! </DESCRIPTION>
37 
38 use fms_mod, only : write_version_number
39 use mpp_mod, only : mpp_error, fatal, mpp_pe, mpp_root_pe
40 use mpp_io_mod, only : mpp_multi
42 use constants_mod, only : pi, radius
43 
44 implicit none
45 private
46 
47 character(len=*), parameter :: &
48  grid_dir = 'INPUT/' ! root directory for all grid files
49 
50 integer, parameter :: &
51  max_name = 256, & ! max length of the variable names
52  max_file = 1024, & ! max length of the file names
53  x_refine = 2, & ! supergrid size/model grid size in x-direction
54  y_refine = 2 ! supergrid size/model grid size in y-direction
55 
56 ! --- public interface
57 
58 
59 public :: get_mosaic_ntiles
60 public :: get_mosaic_ncontacts
61 public :: get_mosaic_grid_sizes
62 public :: get_mosaic_contact
63 public :: get_mosaic_xgrid_size
64 public :: get_mosaic_xgrid
65 public :: calc_mosaic_grid_area
67 public :: is_inside_polygon
68 
69 logical :: module_is_initialized = .true.
70 ! Include variable "version" to be written to log file.
71 #include<file_version.h>
72 
73 contains
74 
75 !#######################################################################
76 
77 ! <SUBROUTINE NAME="mosaic_init">
78 ! <OVERVIEW>
79 ! Initialize the mosaic_mod.
80 ! </OVERVIEW>
81 ! <DESCRIPTION>
82 ! Initialization routine for the mosaic module. It writes the
83 ! version information to the log file.
84 ! </DESCRIPTION>
85 ! <TEMPLATE>
86 ! call mosaic_init ( )
87 ! </TEMPLATE>
88 subroutine mosaic_init()
89 
90  if (module_is_initialized) return
91  module_is_initialized = .true.
92 
93 !--------- write version number and namelist ------------------
94  call write_version_number("MOSAIC_MOD", version)
95 
96 end subroutine mosaic_init
97 ! </SUBROUTINE>
98 
99 !#######################################################################
100 ! <FUNCTION NAME="get_mosaic_xgrid_size">
101 ! <OVERVIEW>
102 ! return exchange grid size of mosaic xgrid file.
103 ! </OVERVIEW>
104 ! <DESCRIPTION>
105 ! return exchange grid size of mosaic xgrid file.
106 ! </DESCRIPTION>
107 ! <TEMPLATE>
108 ! nxgrid = get_mosaic_xgrid_size(xgrid_file)
109 ! </TEMPLATE>
110 ! <IN NAME="xgrid_file" TYPE="character(len=*)">
111 ! The file that contains exchange grid information.
112 ! </IN>
113  function get_mosaic_xgrid_size(xgrid_file)
114  character(len=*), intent(in) :: xgrid_file
115  integer :: get_mosaic_xgrid_size
116 
117  get_mosaic_xgrid_size = dimension_size(xgrid_file, "ncells", no_domain=.true.)
118 
119  return
120 
121  end function get_mosaic_xgrid_size
122 ! </FUNCTION>
123 !#######################################################################
124 ! <SUBROUTINE NAME="get_mosaic_xgrid">
125 ! <OVERVIEW>
126 ! get exchange grid information from mosaic xgrid file.
127 ! </OVERVIEW>
128 ! <DESCRIPTION>
129 ! get exchange grid information from mosaic xgrid file.
130 ! </DESCRIPTION>
131 ! <TEMPLATE>
132 ! call get_mosaic_xgrid(xgrid_file, nxgrid, i1, j1, i2, j2, area)
133 ! </TEMPLATE>
134 ! <IN NAME="xgrid_file" TYPE="character(len=*)">
135 ! The file that contains exchange grid information.
136 ! </IN>
137 ! <INOUT NAME="nxgrid" TYPE="integer">
138 ! number of exchange grid in xgrid_file
139 ! </INOUT>
140 ! <INOUT NAME="i1, j1" TYPE="integer, dimension(:)">
141 ! i and j-index in grid 1 of exchange grid.
142 ! </INOUT>
143 ! <INOUT NAME="i2, j2" TYPE="integer, dimension(:)">
144 ! i and j-index in grid 2 of exchange grid.
145 ! </INOUT>
146 ! <INOUT NAME="area" TYPE="real, dimension(:)">
147 ! area of the exchange grid. The area is scaled to represent unit earth area.
148 ! </INOUT>
149  subroutine get_mosaic_xgrid(xgrid_file, i1, j1, i2, j2, area, ibegin, iend)
150  character(len=*), intent(in) :: xgrid_file
151  integer, intent(inout) :: i1(:), j1(:), i2(:), j2(:)
152  real, intent(inout) :: area(:)
153  integer, optional, intent(in) :: ibegin, iend
154 
155  integer :: start(4), nread(4), istart
156  real, dimension(2, size(i1(:))) :: tile1_cell, tile2_cell
157  integer :: nxgrid, n
158  real :: garea
159  real :: get_global_area;
160 
161  garea = get_global_area();
162 
163  ! When start and nread present, make sure nread(1) is the same as the size of the data
164  if(present(ibegin) .and. present(iend)) then
165  istart = ibegin
166  nxgrid = iend - ibegin + 1
167  if(nxgrid .NE. size(i1(:))) call mpp_error(fatal, .NE."get_mosaic_xgrid: nxgrid size(i1(:))")
168  if(nxgrid .NE. size(j1(:))) call mpp_error(fatal, .NE."get_mosaic_xgrid: nxgrid size(j1(:))")
169  if(nxgrid .NE. size(i2(:))) call mpp_error(fatal, .NE."get_mosaic_xgrid: nxgrid size(i2(:))")
170  if(nxgrid .NE. size(j2(:))) call mpp_error(fatal, .NE."get_mosaic_xgrid: nxgrid size(j2(:))")
171  if(nxgrid .NE. size(area(:))) call mpp_error(fatal, .NE."get_mosaic_xgrid: nxgrid size(area(:))")
172  else
173  istart = 1
174  nxgrid = size(i1(:))
175  endif
176 
177  start = 1; nread = 1
178  start(1) = istart; nread(1) = nxgrid
179  call read_compressed(xgrid_file, 'xgrid_area', area, start=start, nread=nread, threading=mpp_multi)
180  start = 1; nread = 1
181  nread(1) = 2
182  start(2) = istart; nread(2) = nxgrid
183  call read_compressed(xgrid_file, 'tile1_cell', tile1_cell, start=start, nread=nread, threading=mpp_multi)
184  call read_compressed(xgrid_file, 'tile2_cell', tile2_cell, start=start, nread=nread, threading=mpp_multi)
185 
186  do n = 1, nxgrid
187  i1(n) = tile1_cell(1,n)
188  j1(n) = tile1_cell(2,n)
189  i2(n) = tile2_cell(1,n)
190  j2(n) = tile2_cell(2,n)
191  area(n) = area(n)/garea
192  end do
193 
194  return
195 
196  end subroutine get_mosaic_xgrid
197 ! </SUBROUTINE>
198 
199  !###############################################################################
200  ! <SUBROUTINE NAME="get_mosaic_ntiles">
201  ! <OVERVIEW>
202  ! get number of tiles in the mosaic_file.
203  ! </OVERVIEW>
204  ! <DESCRIPTION>
205  ! get number of tiles in the mosaic_file.
206  ! </DESCRIPTION>
207  ! <TEMPLATE>
208  ! ntiles = get_mosaic_ntiles( mosaic_file)
209  ! </TEMPLATE>
210  ! <IN NAME="mosaic_file" TYPE="character(len=*)">
211  ! The file that contains mosaic information.
212  ! </IN>
213  function get_mosaic_ntiles(mosaic_file)
214  character(len=*), intent(in) :: mosaic_file
215  integer :: get_mosaic_ntiles
216 
217  get_mosaic_ntiles = dimension_size(mosaic_file, "ntiles")
218 
219  return
220 
221  end function get_mosaic_ntiles
222 ! </SUBROUTINE>
223 
224  !###############################################################################
225  ! <SUBROUTINE NAME="get_mosaic_ncontacts">
226  ! <OVERVIEW>
227  ! get number of contacts in the mosaic_file.
228  ! </OVERVIEW>
229  ! <DESCRIPTION>
230  ! get number of contacts in the mosaic_file.
231  ! </DESCRIPTION>
232  ! <TEMPLATE>
233  ! ntiles = get_mosaic_ncontacts( mosaic_file)
234  ! </TEMPLATE>
235  ! <IN NAME="mosaic_file" TYPE="character(len=*)">
236  ! The file that contains mosaic information.
237  ! </IN>
238  function get_mosaic_ncontacts( mosaic_file)
239  character(len=*), intent(in) :: mosaic_file
240  integer :: get_mosaic_ncontacts
241 
242  character(len=len_trim(mosaic_file)+1) :: mfile
243  integer :: strlen
244  integer :: read_mosaic_ncontacts
245 
246  if(field_exist(mosaic_file, "contacts") ) then
247  get_mosaic_ncontacts = dimension_size(mosaic_file, "ncontact", no_domain=.true.)
248  else
250  endif
251 
252  return
253 
254  end function get_mosaic_ncontacts
255 ! </SUBROUTINE>
256 
257 
258  !###############################################################################
259  ! <SUBROUTINE NAME="get_mosaic_grid_sizes">
260  ! <OVERVIEW>
261  ! get grid size of each tile from mosaic_file
262  ! </OVERVIEW>
263  ! <DESCRIPTION>
264  ! get grid size of each tile from mosaic_file
265  ! </DESCRIPTION>
266  ! <TEMPLATE>
267  ! call get_mosaic_grid_sizes(mosaic_file, nx, ny)
268  ! </TEMPLATE>
269  ! <IN NAME="mosaic_file" TYPE="character(len=*)">
270  ! The file that contains mosaic information.
271  ! </IN>
272  ! <INOUT NAME="nx" TYPE="integer, dimension(:)">
273  ! List of grid size in x-direction of each tile.
274  ! </INOUT>
275  ! <INOUT NAME="ny" TYPE="integer, dimension(:)">
276  ! List of grid size in y-direction of each tile.
277  ! </INOUT>
278  subroutine get_mosaic_grid_sizes( mosaic_file, nx, ny)
279  character(len=*), intent(in) :: mosaic_file
280  integer, dimension(:), intent(inout) :: nx, ny
281 
282  character(len=MAX_FILE) :: gridfile
283  integer :: ntiles, n
284 
285  ntiles = get_mosaic_ntiles(mosaic_file)
286  if(ntiles .NE. size(nx(:)) .OR. ntiles .NE. size(ny(:)) ) then
287  call mpp_error(fatal, "get_mosaic_grid_sizes: size of nx/ny does not equal to ntiles")
288  endif
289  do n = 1, ntiles
290  call read_data(mosaic_file, 'gridfiles', gridfile, level=n)
291  gridfile = grid_dir//trim(gridfile)
292  nx(n) = dimension_size(gridfile, "nx")
293  ny(n) = dimension_size(gridfile, "ny")
294  if(mod(nx(n),x_refine) .NE. 0) call mpp_error(fatal, "get_mosaic_grid_sizes: nx is not divided by x_refine");
295  if(mod(ny(n),y_refine) .NE. 0) call mpp_error(fatal, "get_mosaic_grid_sizes: ny is not divided by y_refine");
296  nx(n) = nx(n)/x_refine;
297  ny(n) = ny(n)/y_refine;
298  enddo
299 
300  return
301 
302  end subroutine get_mosaic_grid_sizes
303 ! </SUBROUTINE>
304 
305  !###############################################################################
306  ! <SUBROUTINE NAME="get_mosaic_contact">
307  ! <OVERVIEW>
308  ! get contact information from mosaic_file
309  ! </OVERVIEW>
310  ! <DESCRIPTION>
311  ! get contact information from mosaic_file
312  ! </DESCRIPTION>
313  ! <TEMPLATE>
314  ! call get_mosaic_contact(mosaic_file, tile1, tile2, istart1, iend1, jstart1, jend1,
315  ! istart2, iend2, jstart2, jend2)
316  ! </TEMPLATE>
317  ! <IN NAME="mosaic_file" TYPE="character(len=*)">
318  ! The file that contains mosaic information.
319  ! </IN>
320  ! <INOUT NAME="tile1" TYPE="integer, dimension(:)">
321  ! list tile number in tile 1 of each contact.
322  ! </INOUT>
323  ! <INOUT NAME="tile1" TYPE="integer, dimension(:)">
324  ! list tile number in tile 2 of each contact.
325  ! </INOUT>
326  ! <INOUT NAME="istart1" TYPE="integer, dimension(:)">
327  ! list starting i-index in tile 1 of each contact.
328  ! </INOUT>
329  ! <INOUT NAME="iend1" TYPE="integer, dimension(:)">
330  ! list ending i-index in tile 1 of each contact.
331  ! </INOUT>
332  ! <INOUT NAME="jstart1" TYPE="integer, dimension(:)">
333  ! list starting j-index in tile 1 of each contact.
334  ! </INOUT>
335  ! <INOUT NAME="jend1" TYPE="integer, dimension(:)">
336  ! list ending j-index in tile 1 of each contact.
337  ! </INOUT>
338  ! <INOUT NAME="istart2" TYPE="integer, dimension(:)">
339  ! list starting i-index in tile 2 of each contact.
340  ! </INOUT>
341  ! <INOUT NAME="iend2" TYPE="integer, dimension(:)">
342  ! list ending i-index in tile 2 of each contact.
343  ! </INOUT>
344  ! <INOUT NAME="jstart2" TYPE="integer, dimension(:)">
345  ! list starting j-index in tile 2 of each contact.
346  ! </INOUT>
347  ! <INOUT NAME="jend2" TYPE="integer, dimension(:)">
348  ! list ending j-index in tile 2 of each contact.
349  ! </INOUT>
350  subroutine get_mosaic_contact( mosaic_file, tile1, tile2, istart1, iend1, jstart1, jend1, &
351  istart2, iend2, jstart2, jend2)
352  character(len=*), intent(in) :: mosaic_file
353  integer, dimension(:), intent(inout) :: tile1, tile2
354  integer, dimension(:), intent(inout) :: istart1, iend1, jstart1, jend1
355  integer, dimension(:), intent(inout) :: istart2, iend2, jstart2, jend2
356  character(len=MAX_NAME), allocatable :: gridtiles(:)
357  character(len=MAX_NAME) :: contacts
358  character(len=MAX_NAME) :: strlist(8)
359  integer :: ntiles, n, m, ncontacts, nstr, ios
360  integer :: i1_type, j1_type, i2_type, j2_type
361  logical :: found
362 
363  ntiles = get_mosaic_ntiles(mosaic_file)
364  allocate(gridtiles(ntiles))
365  do n = 1, ntiles
366  call read_data(mosaic_file, 'gridtiles', gridtiles(n), level=n)
367  enddo
368 
369  ncontacts = get_mosaic_ncontacts(mosaic_file)
370 
371  do n = 1, ncontacts
372  call read_data(mosaic_file, "contacts", contacts, level=n)
373  nstr = parse_string(contacts, ":", strlist)
374  if(nstr .NE. 4) call mpp_error(fatal, &
375  "mosaic_mod(get_mosaic_contact): number of elements in contact seperated by :/:: should be 4")
376  found = .false.
377  do m = 1, ntiles
378  if(trim(gridtiles(m)) == trim(strlist(2)) ) then !found the tile name
379  found = .true.
380  tile1(n) = m
381  exit
382  endif
383  enddo
384 
385  if(.not.found) call mpp_error(fatal, &
386  "mosaic_mod(get_mosaic_contact):the first tile name specified in contact is not found in tile list")
387 
388  found = .false.
389  do m = 1, ntiles
390  if(trim(gridtiles(m)) == trim(strlist(4)) ) then !found the tile name
391  found = .true.
392  tile2(n) = m
393  exit
394  endif
395  enddo
396 
397  if(.not.found) call mpp_error(fatal, &
398  "mosaic_mod(get_mosaic_contact):the second tile name specified in contact is not found in tile list")
399 
400  call read_data(mosaic_file, "contact_index", contacts, level=n)
401  nstr = parse_string(contacts, ":,", strlist)
402  if(nstr .NE. 8) then
403  if(mpp_pe()==mpp_root_pe()) then
404  print*, "nstr is ", nstr
405  print*, "contacts is ", contacts
406  do m = 1, nstr
407  print*, "strlist is ", trim(strlist(m))
408  enddo
409  endif
410  call mpp_error(fatal, &
411  "mosaic_mod(get_mosaic_contact): number of elements in contact_index seperated by :/, should be 8")
412  endif
413  read(strlist(1), *, iostat=ios) istart1(n)
414  if(ios .NE. 0) call mpp_error(fatal, &
415  "mosaic_mod(get_mosaic_contact): Error in reading istart1")
416  read(strlist(2), *, iostat=ios) iend1(n)
417  if(ios .NE. 0) call mpp_error(fatal, &
418  "mosaic_mod(get_mosaic_contact): Error in reading iend1")
419  read(strlist(3), *, iostat=ios) jstart1(n)
420  if(ios .NE. 0) call mpp_error(fatal, &
421  "mosaic_mod(get_mosaic_contact): Error in reading jstart1")
422  read(strlist(4), *, iostat=ios) jend1(n)
423  if(ios .NE. 0) call mpp_error(fatal, &
424  "mosaic_mod(get_mosaic_contact): Error in reading jend1")
425  read(strlist(5), *, iostat=ios) istart2(n)
426  if(ios .NE. 0) call mpp_error(fatal, &
427  "mosaic_mod(get_mosaic_contact): Error in reading istart2")
428  read(strlist(6), *, iostat=ios) iend2(n)
429  if(ios .NE. 0) call mpp_error(fatal, &
430  "mosaic_mod(get_mosaic_contact): Error in reading iend2")
431  read(strlist(7), *, iostat=ios) jstart2(n)
432  if(ios .NE. 0) call mpp_error(fatal, &
433  "mosaic_mod(get_mosaic_contact): Error in reading jstart2")
434  read(strlist(8), *, iostat=ios) jend2(n)
435  if(ios .NE. 0) call mpp_error(fatal, &
436  "mosaic_mod(get_mosaic_contact): Error in reading jend2")
437 
438  i1_type = transfer_to_model_index(istart1(n), iend1(n), x_refine)
439  j1_type = transfer_to_model_index(jstart1(n), jend1(n), y_refine)
440  i2_type = transfer_to_model_index(istart2(n), iend2(n), x_refine)
441  j2_type = transfer_to_model_index(jstart2(n), jend2(n), y_refine)
442 
443  if( i1_type == 0 .AND. j1_type == 0 ) call mpp_error(fatal, &
444  "mosaic_mod(get_mosaic_contact): istart1==iend1 and jstart1==jend1")
445  if( i2_type == 0 .AND. j2_type == 0 ) call mpp_error(fatal, &
446  "mosaic_mod(get_mosaic_contact): istart2==iend2 and jstart2==jend2")
447  if( i1_type + j1_type .NE. i2_type + j2_type ) call mpp_error(fatal, &
448  "mosaic_mod(get_mosaic_contact): It is not a line or overlap contact")
449 
450  enddo
451 
452  deallocate(gridtiles)
453 
454  end subroutine get_mosaic_contact
455 ! </SUBROUTINE>
456 
457 
458 function transfer_to_model_index(istart, iend, refine_ratio)
459  integer, intent(inout) :: istart, iend
460  integer :: refine_ratio
461  integer :: transfer_to_model_index
462  integer :: istart_in, iend_in
463 
464  istart_in = istart
465  iend_in = iend
466 
467  if( istart_in == iend_in ) then
469  istart = (istart_in + 1)/refine_ratio
470  iend = istart
471  else
473  if( iend_in > istart_in ) then
474  istart = istart_in + 1
475  iend = iend_in
476  else
477  istart = istart_in
478  iend = iend_in + 1
479  endif
480  if( mod(istart, refine_ratio) .NE. 0 .OR. mod(iend,refine_ratio) .NE. 0) call mpp_error(fatal, &
481  "mosaic_mod(transfer_to_model_index): mismatch between refine_ratio and istart/iend")
482  istart = istart/refine_ratio
483  iend = iend/refine_ratio
484 
485  endif
486 
487  return
488 
489 end function transfer_to_model_index
490 
491  !###############################################################################
492  ! <SUBROUTINE NAME="calc_mosaic_grid_area">
493  ! <OVERVIEW>
494  ! calculate grid cell area.
495  ! </OVERVIEW>
496  ! <DESCRIPTION>
497  ! calculate the grid cell area. The purpose of this routine is to make
498  ! sure the consistency between model grid area and exchange grid area.
499  ! </DESCRIPTION>
500  ! <TEMPLATE>
501  ! call calc_mosaic_grid_area(lon, lat, area)
502  ! </TEMPLATE>
503  ! <IN NAME="lon" TYPE="real, dimension(:,:)">
504  ! geographical longitude of grid cell vertices.
505  ! </IN>
506  ! <IN NAME="lat" TYPE="real, dimension(:,:)">
507  ! geographical latitude of grid cell vertices.
508  ! </IN>
509  ! <INOUT NAME="area" TYPE="real, dimension(:,:)">
510  ! grid cell area.
511  ! </INOUT>
512  subroutine calc_mosaic_grid_area(lon, lat, area)
513  real, dimension(:,:), intent(in) :: lon
514  real, dimension(:,:), intent(in) :: lat
515  real, dimension(:,:), intent(inout) :: area
516  integer :: nlon, nlat
517 
518  nlon = size(area,1)
519  nlat = size(area,2)
520  ! make sure size of lon, lat and area are consitency
521  if( size(lon,1) .NE. nlon+1 .OR. size(lat,1) .NE. nlon+1 ) &
522  call mpp_error(fatal, "mosaic_mod: size(lon,1) and size(lat,1) should equal to size(area,1)+1")
523  if( size(lon,2) .NE. nlat+1 .OR. size(lat,2) .NE. nlat+1 ) &
524  call mpp_error(fatal, "mosaic_mod: size(lon,2) and size(lat,2) should equal to size(area,2)+1")
525 
526  call get_grid_area( nlon, nlat, lon, lat, area)
527 
528  end subroutine calc_mosaic_grid_area
529  ! </SUBROUTINE>
530 
531  !###############################################################################
532  ! <SUBROUTINE NAME="calc_mosaic_grid_great_circle_area">
533  ! <OVERVIEW>
534  ! calculate grid cell area using great cirlce algorithm
535  ! </OVERVIEW>
536  ! <DESCRIPTION>
537  ! calculate the grid cell area. The purpose of this routine is to make
538  ! sure the consistency between model grid area and exchange grid area.
539  ! </DESCRIPTION>
540  ! <TEMPLATE>
541  ! call calc_mosaic_grid_great_circle_area(lon, lat, area)
542  ! </TEMPLATE>
543  ! <IN NAME="lon" TYPE="real, dimension(:,:)">
544  ! geographical longitude of grid cell vertices.
545  ! </IN>
546  ! <IN NAME="lat" TYPE="real, dimension(:,:)">
547  ! geographical latitude of grid cell vertices.
548  ! </IN>
549  ! <INOUT NAME="area" TYPE="real, dimension(:,:)">
550  ! grid cell area.
551  ! </INOUT>
552  subroutine calc_mosaic_grid_great_circle_area(lon, lat, area)
553  real, dimension(:,:), intent(in) :: lon
554  real, dimension(:,:), intent(in) :: lat
555  real, dimension(:,:), intent(inout) :: area
556  integer :: nlon, nlat
557 
558 
559  nlon = size(area,1)
560  nlat = size(area,2)
561  ! make sure size of lon, lat and area are consitency
562  if( size(lon,1) .NE. nlon+1 .OR. size(lat,1) .NE. nlon+1 ) &
563  call mpp_error(fatal, "mosaic_mod: size(lon,1) and size(lat,1) should equal to size(area,1)+1")
564  if( size(lon,2) .NE. nlat+1 .OR. size(lat,2) .NE. nlat+1 ) &
565  call mpp_error(fatal, "mosaic_mod: size(lon,2) and size(lat,2) should equal to size(area,2)+1")
566 
567  call get_grid_great_circle_area( nlon, nlat, lon, lat, area)
568 
570  ! </SUBROUTINE>
571 
572  !#####################################################################
573  ! This function check if a point (lon1,lat1) is inside a polygon (lon2(:), lat2(:))
574  ! lon1, lat1, lon2, lat2 are in radians.
575  function is_inside_polygon(lon1, lat1, lon2, lat2 )
576  real, intent(in) :: lon1, lat1
577  real, intent(in) :: lon2(:), lat2(:)
578  logical :: is_inside_polygon
579  real, dimension(size(lon2(:))) :: x2, y2, z2
580  integer :: npts, isinside
581  integer :: inside_a_polygon
582 
583  npts = size(lon2(:))
584 
585  isinside = inside_a_polygon(lon1, lat1, npts, lon2, lat2)
586  if(isinside == 1) then
587  is_inside_polygon = .true.
588  else
589  is_inside_polygon = .false.
590  endif
591 
592  return
593 
594  end function is_inside_polygon
595 
596  function parse_string(string, set, value)
597  character(len=*), intent(in) :: string
598  character(len=*), intent(in) :: set
599  character(len=*), intent(out) :: value(:)
600  integer :: parse_string
601  integer :: nelem, length, first, last
602 
603  nelem = size(value(:))
604  length = len_trim(string)
605 
606  first = 1; last = 0
607  parse_string = 0
608 
609  do while(first .LE. length)
611  if(parse_string>nelem) then
612  call mpp_error(fatal, "mosaic_mod(parse_string) : number of element is greater than size(value(:))")
613  endif
614  last = first - 1 + scan(string(first:length), set)
615  if(last == first-1 ) then ! not found, end of string
616  value(parse_string) = string(first:length)
617  exit
618  else
619  if(last <= first) then
620  call mpp_error(fatal, "mosaic_mod(parse_string) : last <= first")
621  endif
622  value(parse_string) = string(first:(last-1))
623  first = last + 1
624  ! scan to make sure the next is not the character in the set
625  do while (first == last+1)
626  last = first - 1 + scan(string(first:length), set)
627  if(last == first) then
628  first = first+1
629  else
630  exit
631  endif
632  end do
633  endif
634  enddo
635 
636  return
637 
638  end function parse_string
639 
640 
641 
642 end module mosaic_mod
643 
644 
645 #ifdef TEST_MOSAIC
646 program test_mosaic
647 
650 
651 implicit none
652 
653 integer :: ntiles, ncontacts, n
654 integer, allocatable :: tile1(:), tile2(:), nx(:), ny(:)
655 integer, allocatable :: istart1(:), iend1(:), jstart1(:), jend1(:)
656 integer, allocatable :: istart2(:), iend2(:), jstart2(:), jend2(:)
657 character(len=128) :: mosaic_file = "INPUT/mosaic.nc"
658 
659 ntiles = get_mosaic_ntiles(mosaic_file)
660 ncontacts = get_mosaic_ncontacts(mosaic_file)
661 allocate(nx(ntiles), ny(ntiles))
662 allocate(tile1(ncontacts), tile2(ncontacts) )
663 allocate(istart1(ncontacts), iend1(ncontacts), jstart1(ncontacts), jend1(ncontacts) )
664 allocate(istart2(ncontacts), iend2(ncontacts), jstart2(ncontacts), jend2(ncontacts) )
665 
666 call get_mosaic_grid_sizes(mosaic_file, nx, ny )
667 call get_mosaic_contact(mosaic_file, tile1, tile2, istart1, iend1, jstart1, jend1, istart2, iend2, jstart2, jend2)
668 
669 ! print out information
670 
671 print '(a,i3,a,a)', "****** There is ", ntiles, " tiles in ", trim(mosaic_file)
672 do n = 1, ntiles
673  print '(a,i3,a,i3,a,i3)', " tile = ", n, ", nx = ", nx(n), ", ny = ", ny(n)
674 end do
675 
676 print '(a,i3,a,a)', "****** There is ", ncontacts, " contacts in ", trim(mosaic_file)
677 do n = 1, ncontacts
678  print '(a,i3,a,i3,a,i3,a,i4,a,i4,a,i4,a,i4,a,i4,a,i4,a,i4,a,i4)', &
679  "contact=", n, ": tile1=", tile1(n), " tile2=", tile2(n), &
680  " is1=", istart1(n), " ie1=", iend1(n), &
681  " js1=", jstart1(n), " je1=", jend1(n), &
682  " is2=", istart2(n), " ie2=", iend2(n), &
683  " js2=", jstart2(n), " je2=", jend2(n)
684 end do
685 
686 deallocate(tile1, tile2, nx, ny)
687 deallocate(istart1, iend1, jstart1, jend1)
688 deallocate(istart2, iend2, jstart2, jend2)
689 
690 
691 end program test_mosaic
692 #endif
Definition: fms.F90:20
subroutine, public get_mosaic_contact(mosaic_file, tile1, tile2, istart1, iend1, jstart1, jend1, istart2, iend2, jstart2, jend2)
Definition: mosaic.F90:352
real, parameter, public radius
Radius of the Earth [m].
Definition: constants.F90:72
double get_global_area(void)
Definition: read_mosaic.c:404
integer, parameter, public set
integer, parameter, public strlen
logical module_is_initialized
Definition: mosaic.F90:69
subroutine, public calc_mosaic_grid_great_circle_area(lon, lat, area)
Definition: mosaic.F90:553
integer, parameter max_name
Definition: mosaic.F90:50
integer, parameter max_file
Definition: mosaic.F90:50
logical function, public is_inside_polygon(lon1, lat1, lon2, lat2)
Definition: mosaic.F90:576
void get_grid_great_circle_area(const int *nlon, const int *nlat, const double *lon, const double *lat, double *area)
Definition: create_xgrid.c:132
int read_mosaic_ncontacts(const char *mosaic_file)
Definition: read_mosaic.c:630
integer function, public get_mosaic_ntiles(mosaic_file)
Definition: mosaic.F90:214
Definition: mpp.F90:39
void get_grid_area(const int *nlon, const int *nlat, const double *lon, const double *lat, double *area)
Definition: create_xgrid.c:64
subroutine, public get_mosaic_grid_sizes(mosaic_file, nx, ny)
Definition: mosaic.F90:279
subroutine, public get_mosaic_xgrid(xgrid_file, i1, j1, i2, j2, area, ibegin, iend)
Definition: mosaic.F90:150
integer function, public dimension_size(filename, dimname, domain, no_domain)
Definition: fms_io.F90:5044
subroutine, public calc_mosaic_grid_area(lon, lat, area)
Definition: mosaic.F90:513
integer, parameter x_refine
Definition: mosaic.F90:50
integer function transfer_to_model_index(istart, iend, refine_ratio)
Definition: mosaic.F90:459
integer function, public get_mosaic_ncontacts(mosaic_file)
Definition: mosaic.F90:239
int inside_a_polygon(double *lon1, double *lat1, int *npts, double *lon2, double *lat2)
Definition: mosaic_util.c:1357
character(len= *), parameter grid_dir
Definition: mosaic.F90:47
integer function parse_string(string, set, value)
Definition: mosaic.F90:597
integer function, public get_mosaic_xgrid_size(xgrid_file)
Definition: mosaic.F90:114
subroutine mosaic_init()
Definition: mosaic.F90:89
integer, parameter y_refine
Definition: mosaic.F90:50
logical function, public field_exist(file_name, field_name, domain, no_domain)
Definition: fms_io.F90:8298
real(fp), parameter, public pi