38 use fms_mod,
only : write_version_number
47 character(len=*),
parameter :: &
50 integer,
parameter :: &
71 #include<file_version.h> 94 call write_version_number(
"MOSAIC_MOD", version)
114 character(len=*),
intent(in) :: xgrid_file
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
155 integer :: start(4), nread(4), istart
156 real,
dimension(2, size(i1(:))) :: tile1_cell, tile2_cell
164 if(
present(ibegin) .and.
present(iend))
then 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(:))")
178 start(1) = istart; nread(1) = nxgrid
179 call read_compressed(xgrid_file,
'xgrid_area', area, start=start, nread=nread, threading=mpp_multi)
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)
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
214 character(len=*),
intent(in) :: mosaic_file
239 character(len=*),
intent(in) :: mosaic_file
242 character(len=len_trim(mosaic_file)+1) :: mfile
279 character(len=*),
intent(in) :: mosaic_file
280 integer,
dimension(:),
intent(inout) :: nx, ny
282 character(len=MAX_FILE) :: gridfile
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")
290 call read_data(mosaic_file,
'gridfiles', gridfile, level=n)
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");
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
364 allocate(gridtiles(ntiles))
366 call read_data(mosaic_file,
'gridtiles', gridtiles(n), level=n)
372 call read_data(mosaic_file,
"contacts", contacts, level=n)
375 "mosaic_mod(get_mosaic_contact): number of elements in contact seperated by :/:: should be 4")
378 if(trim(gridtiles(m)) == trim(strlist(2)) )
then 386 "mosaic_mod(get_mosaic_contact):the first tile name specified in contact is not found in tile list")
390 if(trim(gridtiles(m)) == trim(strlist(4)) )
then 398 "mosaic_mod(get_mosaic_contact):the second tile name specified in contact is not found in tile list")
400 call read_data(mosaic_file,
"contact_index", contacts, level=n)
403 if(mpp_pe()==mpp_root_pe())
then 404 print*,
"nstr is ", nstr
405 print*,
"contacts is ", contacts
407 print*,
"strlist is ", trim(strlist(m))
411 "mosaic_mod(get_mosaic_contact): number of elements in contact_index seperated by :/, should be 8")
413 read(strlist(1), *, iostat=ios) istart1(n)
415 "mosaic_mod(get_mosaic_contact): Error in reading istart1")
416 read(strlist(2), *, iostat=ios) iend1(n)
418 "mosaic_mod(get_mosaic_contact): Error in reading iend1")
419 read(strlist(3), *, iostat=ios) jstart1(n)
421 "mosaic_mod(get_mosaic_contact): Error in reading jstart1")
422 read(strlist(4), *, iostat=ios) jend1(n)
424 "mosaic_mod(get_mosaic_contact): Error in reading jend1")
425 read(strlist(5), *, iostat=ios) istart2(n)
427 "mosaic_mod(get_mosaic_contact): Error in reading istart2")
428 read(strlist(6), *, iostat=ios) iend2(n)
430 "mosaic_mod(get_mosaic_contact): Error in reading iend2")
431 read(strlist(7), *, iostat=ios) jstart2(n)
433 "mosaic_mod(get_mosaic_contact): Error in reading jstart2")
434 read(strlist(8), *, iostat=ios) jend2(n)
436 "mosaic_mod(get_mosaic_contact): Error in reading jend2")
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")
452 deallocate(gridtiles)
459 integer,
intent(inout) :: istart, iend
460 integer :: refine_ratio
462 integer :: istart_in, iend_in
467 if( istart_in == iend_in )
then 469 istart = (istart_in + 1)/refine_ratio
473 if( iend_in > istart_in )
then 474 istart = istart_in + 1
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
513 real,
dimension(:,:),
intent(in) :: lon
514 real,
dimension(:,:),
intent(in) :: lat
515 real,
dimension(:,:),
intent(inout) :: area
516 integer :: nlon, nlat
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")
553 real,
dimension(:,:),
intent(in) :: lon
554 real,
dimension(:,:),
intent(in) :: lat
555 real,
dimension(:,:),
intent(inout) :: area
556 integer :: nlon, nlat
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")
576 real,
intent(in) :: lon1, lat1
577 real,
intent(in) :: lon2(:), lat2(:)
579 real,
dimension(size(lon2(:))) :: x2, y2, z2
580 integer :: npts, isinside
586 if(isinside == 1)
then 598 character(len=*),
intent(in) ::
set 599 character(len=*),
intent(out) ::
value(:)
601 integer :: nelem, length, first, last
603 nelem =
size(value(:))
609 do while(first .LE. length)
612 call mpp_error(fatal,
"mosaic_mod(parse_string) : number of element is greater than size(value(:))")
614 last = first - 1 + scan(
string(first:length),
set)
615 if(last == first-1 )
then 619 if(last <= first)
then 620 call mpp_error(fatal,
"mosaic_mod(parse_string) : last <= first")
625 do while (first == last+1)
626 last = first - 1 + scan(
string(first:length),
set)
627 if(last == first)
then 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" 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) )
667 call get_mosaic_contact(mosaic_file, tile1, tile2, istart1, iend1, jstart1, jend1, istart2, iend2, jstart2, jend2)
671 print
'(a,i3,a,a)',
"****** There is ", ntiles,
" tiles in ", trim(mosaic_file)
673 print
'(a,i3,a,i3,a,i3)',
" tile = ", n,
", nx = ", nx(n),
", ny = ", ny(n)
676 print
'(a,i3,a,a)',
"****** There is ", ncontacts,
" contacts in ", trim(mosaic_file)
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)
686 deallocate(tile1, tile2, nx, ny)
687 deallocate(istart1, iend1, jstart1, jend1)
688 deallocate(istart2, iend2, jstart2, jend2)
691 end program test_mosaic
subroutine, public get_mosaic_contact(mosaic_file, tile1, tile2, istart1, iend1, jstart1, jend1, istart2, iend2, jstart2, jend2)
real, parameter, public radius
Radius of the Earth [m].
double get_global_area(void)
integer, parameter, public set
integer, parameter, public strlen
logical module_is_initialized
subroutine, public calc_mosaic_grid_great_circle_area(lon, lat, area)
integer, parameter max_name
integer, parameter max_file
logical function, public is_inside_polygon(lon1, lat1, lon2, lat2)
void get_grid_great_circle_area(const int *nlon, const int *nlat, const double *lon, const double *lat, double *area)
int read_mosaic_ncontacts(const char *mosaic_file)
integer function, public get_mosaic_ntiles(mosaic_file)
void get_grid_area(const int *nlon, const int *nlat, const double *lon, const double *lat, double *area)
subroutine, public get_mosaic_grid_sizes(mosaic_file, nx, ny)
subroutine, public get_mosaic_xgrid(xgrid_file, i1, j1, i2, j2, area, ibegin, iend)
integer function, public dimension_size(filename, dimname, domain, no_domain)
subroutine, public calc_mosaic_grid_area(lon, lat, area)
integer, parameter x_refine
integer function transfer_to_model_index(istart, iend, refine_ratio)
integer function, public get_mosaic_ncontacts(mosaic_file)
int inside_a_polygon(double *lon1, double *lat1, int *npts, double *lon2, double *lat2)
character(len= *), parameter grid_dir
integer function parse_string(string, set, value)
integer function, public get_mosaic_xgrid_size(xgrid_file)
integer, parameter y_refine
logical function, public field_exist(file_name, field_name, domain, no_domain)
real(fp), parameter, public pi