SCALE-RM
Data Types | Functions/Subroutines
scale_file_tiledata Module Reference

module file_tiledata More...

Functions/Subroutines

subroutine, public file_tiledata_get_info (TILE_nlim, TILE_DLAT, TILE_DLON, DOMAIN_LATS, DOMAIN_LATE, DOMAIN_LONS, DOMAIN_LONE, catalog_fname, GLOBAL_IA, TILE_nmax, TILE_fname, TILE_hit, TILE_JS, TILE_JE, TILE_IS, TILE_IE, nLATH, nLONH, jsh, jeh, ish, ieh, zonal, pole, single_fname, LATS, LATE, LONS, LONE)
 get tile information More...
 
subroutine, public file_tiledata_get_latlon (nLAT, nLON, jsh, ish, TILE_DLAT, TILE_DLON, LAT, LON)
 get tile data More...
 
subroutine file_tiledata_get_data_int1 (nLAT, nLON, dirname, GLOBAL_IA, TILE_nmax, TILE_DLAT, TILE_DLON, TILE_fname, TILE_hit, TILE_JS, TILE_JE, TILE_IS, TILE_IE, jsh, jeh, ish, ieh, data_type, DATA, min_value, yrevers, step)
 
subroutine file_tiledata_read_catalog_file (TILE_nlim, fname, TILE_DLAT, TILE_DLON, DOMAIN_IS, GLOBAL_IA, TILE_nmax, TILE_fname, TILE_LATS, TILE_LATE, TILE_LONS, TILE_LONE)
 read category file More...
 
subroutine file_tiledata_get_domain_info (TILE_DLAT, TILE_DLON, DOMAIN_LATS, DOMAIN_LATE, DOMAIN_LONS, DOMAIN_LONE, DOMAIN_JS, DOMAIN_JE, DOMAIN_IS, DOMAIN_IE, GLOBAL_IA)
 
subroutine file_tiledata_get_tile_info (TILE_nmax, DOMAIN_JS, DOMAIN_JE, DOMAIN_IS, DOMAIN_IE, GLOBAL_IA, TILE_DLAT, TILE_DLON, TILE_LATS, TILE_LATE, TILE_LONS, TILE_LONE, zonal, TILE_hit, TILE_JS, TILE_JE, TILE_IS, TILE_IE, jsh, jeh, ish, ieh, nLAT, nLON)
 
subroutine file_tiledata_read_data_int2_real (jsize, isize, fname, TILE_DATA, yrevers, step)
 
subroutine file_tiledata_read_data_int4_real (jsize, isize, fname, TILE_DATA, yrevers, step)
 
subroutine file_tiledata_read_data_real4_real (jsize, isize, fname, TILE_DATA, yrevers, step)
 
subroutine file_tiledata_read_data_real8_real (jsize, isize, fname, TILE_DATA, yrevers, step)
 
subroutine file_tiledata_read_data_int1_int (jsize, isize, fname, TILE_DATA, yrevers, step)
 
subroutine file_tiledata_read_data_int2_int (jsize, isize, fname, TILE_DATA, yrevers, step)
 
subroutine file_tiledata_read_data_int4_int (jsize, isize, fname, TILE_DATA, yrevers, step)
 
subroutine file_tiledata_read_data_real4_int (jsize, isize, fname, TILE_DATA, yrevers, step)
 

Detailed Description

module file_tiledata

Description
tile data reader
Author
Team SCALE

Function/Subroutine Documentation

◆ file_tiledata_get_info()

subroutine, public scale_file_tiledata::file_tiledata_get_info ( integer, intent(in)  TILE_nlim,
real(rp), intent(in)  TILE_DLAT,
real(rp), intent(in)  TILE_DLON,
real(rp), intent(in)  DOMAIN_LATS,
real(rp), intent(in)  DOMAIN_LATE,
real(rp), intent(in)  DOMAIN_LONS,
real(rp), intent(in)  DOMAIN_LONE,
character(len=*), intent(in)  catalog_fname,
integer, intent(out)  GLOBAL_IA,
integer, intent(out)  TILE_nmax,
character(len=*), dimension(:), intent(out)  TILE_fname,
logical, dimension(:), intent(out)  TILE_hit,
integer, dimension(:), intent(out)  TILE_JS,
integer, dimension(:), intent(out)  TILE_JE,
integer, dimension(:), intent(out)  TILE_IS,
integer, dimension(:), intent(out)  TILE_IE,
integer, intent(out)  nLATH,
integer, intent(out)  nLONH,
integer, intent(out)  jsh,
integer, intent(out)  jeh,
integer, intent(out)  ish,
integer, intent(out)  ieh,
logical, intent(out)  zonal,
logical, intent(out)  pole,
character(len=*), intent(in), optional  single_fname,
real(rp), intent(in), optional  LATS,
real(rp), intent(in), optional  LATE,
real(rp), intent(in), optional  LONS,
real(rp), intent(in), optional  LONE 
)

get tile information

Definition at line 62 of file scale_file_tiledata.F90.

62  use scale_const, only: &
63  pi => const_pi
64  use scale_prc, only: &
65  prc_abort
66  integer, intent(in) :: TILE_nlim
67  real(RP), intent(in) :: TILE_DLAT, TILE_DLON
68  real(RP), intent(in) :: DOMAIN_LATS, DOMAIN_LATE, DOMAIN_LONS, DOMAIN_LONE
69  character(len=*), intent(in) :: catalog_fname
70  integer, intent(out) :: GLOBAL_IA
71  integer, intent(out) :: TILE_nmax
72  logical, intent(out) :: TILE_hit(:)
73  character(len=*), intent(out) :: TILE_fname(:)
74  integer, intent(out) :: TILE_JS(:), TILE_JE(:), TILE_IS(:), TILE_IE(:)
75  integer, intent(out) :: nLATH, nLONH
76  integer, intent(out) :: jsh, jeh, ish, ieh
77  logical, intent(out) :: zonal, pole
78 
79  character(len=*), intent(in), optional :: single_fname
80  real(RP), intent(in), optional :: LATS
81  real(RP), intent(in), optional :: LATE
82  real(RP), intent(in), optional :: LONS
83  real(RP), intent(in), optional :: LONE
84 
85  real(RP) :: TILE_LATS(TILE_nlim), TILE_LATE(TILE_nlim)
86  real(RP) :: TILE_LONS(TILE_nlim), TILE_LONE(TILE_nlim)
87  real(RP) :: LAT_MIN, LAT_MAX
88 
89  integer :: DOMAIN_JS, DOMAIN_JE, DOMAIN_IS, DOMAIN_IE
90 
91  call file_tiledata_get_domain_info( tile_dlat, tile_dlon, & ! [IN]
92  domain_lats, domain_late, domain_lons, domain_lone, & ! [IN]
93  domain_js, domain_je, domain_is, domain_ie, & ! [OUT]
94  global_ia ) ! [OUT]
95 
96 
97  if ( catalog_fname /= "" ) then
98  log_info("FILE_TILEDATA_get_info",*) 'Input catalogue file:', trim(catalog_fname)
99 
100  call file_tiledata_read_catalog_file( tile_nlim, & ! [IN]
101  catalog_fname, & ! [IN]
102  tile_dlat, tile_dlon, & ! [IN]
103  domain_is, global_ia, & ! [IN]
104  tile_nmax, & ! [OUT]
105  tile_fname(:), & ! [OUT]
106  tile_lats(:), tile_late(:), & ! [OUT]
107  tile_lons(:), tile_lone(:) ) ! [OUT]
108  lat_min = minval( tile_lats(1:tile_nmax) )
109  lat_max = maxval( tile_late(1:tile_nmax) )
110 
111  else
112  if ( .not. present(single_fname) ) then
113  log_error("FILE_TILEDATA_get_info",*) "single_fname is required if catalog_fname is empty"
114  call prc_abort
115  end if
116  if ( .not. present(lats) ) then
117  log_error("FILE_TILEDATA_get_info",*) "LATS is required if catalog_fname is empty"
118  call prc_abort
119  end if
120  if ( .not. present(late) ) then
121  log_error("FILE_TILEDATA_get_info",*) "LATE is required if catalog_fname is empty"
122  call prc_abort
123  end if
124  if ( .not. present(lons) ) then
125  log_error("FILE_TILEDATA_get_info",*) "LONS is required if catalog_fname is empty"
126  call prc_abort
127  end if
128  if ( .not. present(lone) ) then
129  log_error("FILE_TILEDATA_get_info",*) "LONE is required if catalog_fname is empty"
130  call prc_abort
131  end if
132 
133  tile_nmax = 1
134  tile_fname(1) = single_fname
135  tile_lats(1) = lats
136  tile_late(1) = late
137  tile_lons(1) = lons
138  tile_lone(1) = lone
139 
140  lat_min = lats
141  lat_max = late
142  end if
143 
144  zonal = ( domain_lone - domain_lons ) / ( 2.0_rp * pi ) > 0.9_rp
145 
146  pole = ( domain_lats < - pi * 0.5_rp + ( domain_late - domain_lats ) * 0.1_rp ) &
147  .or. ( domain_late > pi * 0.5_rp - ( domain_late - domain_lats ) * 0.1_rp )
148 
149  zonal = zonal .or. pole
150 
151  call file_tiledata_get_tile_info( tile_nmax, & ! [IN]
152  domain_js, domain_je, & ! [IN]
153  domain_is, domain_ie, & ! [IN]
154  global_ia, & ! [IN]
155  tile_dlat, tile_dlon, & ! [IN]
156  tile_lats(:), tile_late(:), & ! [IN]
157  tile_lons(:), tile_lone(:), & ! [IN]
158  zonal, & ! [IN]
159  tile_hit(:), & ! [OUT]
160  tile_js(:), tile_je(:), & ! [OUT]
161  tile_is(:), tile_ie(:), & ! [OUT]
162  jsh, jeh, ish, ieh, & ! [OUT]
163  nlath, nlonh ) ! [OUT]
164 
165  return

References scale_const::const_pi, file_tiledata_get_domain_info(), file_tiledata_get_tile_info(), file_tiledata_read_catalog_file(), and scale_prc::prc_abort().

Referenced by mod_cnv2d::cnv2d_tile_init(), and mod_cnvlanduse::cnvlanduse().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ file_tiledata_get_latlon()

subroutine, public scale_file_tiledata::file_tiledata_get_latlon ( integer, intent(in)  nLAT,
integer, intent(in)  nLON,
integer, intent(in)  jsh,
integer, intent(in)  ish,
real(rp), intent(in)  TILE_DLAT,
real(rp), intent(in)  TILE_DLON,
real(rp), dimension(nlat), intent(out)  LAT,
real(rp), dimension(nlon), intent(out)  LON 
)

get tile data

Definition at line 175 of file scale_file_tiledata.F90.

175  implicit none
176  integer, intent(in) :: nLAT, nLON
177  integer, intent(in) :: jsh, ish
178  real(RP), intent(in) :: TILE_DLAT, TILE_DLON
179  real(RP), intent(out) :: LAT(nLAT)
180  real(RP), intent(out) :: LON(nLON)
181 
182  integer :: i, j
183 
184  do j = 1, nlat
185  lat(j) = tile_dlat * ( jsh + j - 0.5_rp )
186  end do
187  do i = 1, nlon
188  lon(i) = tile_dlon * ( ish + i - 0.5_rp )
189  end do
190 
191  return

References scale_const::const_d2r, scale_const::const_pi, scale_const::const_undef, file_tiledata_read_data_int2_real(), file_tiledata_read_data_int4_real(), file_tiledata_read_data_real4_real(), file_tiledata_read_data_real8_real(), and scale_prc::prc_abort().

Referenced by mod_cnv2d::cnv2d_tile_init(), and mod_cnvlanduse::cnvlanduse().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ file_tiledata_get_data_int1()

subroutine scale_file_tiledata::file_tiledata_get_data_int1 ( integer, intent(in)  nLAT,
integer, intent(in)  nLON,
character(len=*), intent(in)  dirname,
integer, intent(in)  GLOBAL_IA,
integer, intent(in)  TILE_nmax,
real(rp), intent(in)  TILE_DLAT,
real(rp), intent(in)  TILE_DLON,
character(len=*), dimension(:), intent(in)  TILE_fname,
logical, dimension(:), intent(in)  TILE_hit,
integer, dimension(:), intent(in)  TILE_JS,
integer, dimension(:), intent(in)  TILE_JE,
integer, dimension(:), intent(in)  TILE_IS,
integer, dimension(:), intent(in)  TILE_IE,
integer, intent(in)  jsh,
integer, intent(in)  jeh,
integer, intent(in)  ish,
integer, intent(in)  ieh,
character(len=*), intent(in)  data_type,
integer, dimension(nlon,nlat), intent(out)  DATA,
integer, intent(in), optional  min_value,
logical, intent(in), optional  yrevers,
integer, intent(in), optional  step 
)

Definition at line 361 of file scale_file_tiledata.F90.

361  use scale_const, only: &
362  undef => const_undef, &
363  undef2 => const_undef2, &
364  pi => const_pi, &
365  d2r => const_d2r
366  use scale_prc, only: &
367  prc_abort
368  integer, intent(in) :: nLAT, nLON
369  character(len=*), intent(in) :: dirname
370  integer, intent(in) :: GLOBAL_IA
371  integer, intent(in) :: TILE_nmax
372  real(RP), intent(in) :: TILE_DLAT, TILE_DLON
373  character(len=*), intent(in) :: TILE_fname(:)
374  logical, intent(in) :: TILE_hit(:)
375  integer, intent(in) :: TILE_JS(:)
376  integer, intent(in) :: TILE_JE(:)
377  integer, intent(in) :: TILE_IS(:)
378  integer, intent(in) :: TILE_IE(:)
379  integer, intent(in) :: jsh
380  integer, intent(in) :: jeh
381  integer, intent(in) :: ish
382  integer, intent(in) :: ieh
383  character(len=*), intent(in) :: data_type
384  integer, intent(out) :: DATA(nLON,nLAT)
385 
386  integer, intent(in), optional :: min_value
387  logical, intent(in), optional :: yrevers
388  integer, intent(in), optional :: step
389 
390  abstract interface
391  subroutine rd( &
392  jsize, isize, &
393  fname, &
394  TILE_DATA, &
395  yrevers, &
396  step )
397  use scale_precision
398  integer, intent(in) :: jsize
399  integer, intent(in) :: isize
400  character(len=*), intent(in) :: fname
401  integer, intent(out) :: TILE_DATA(isize,jsize)
402  logical, intent(in), optional :: yrevers
403  integer, intent(in), optional :: step
404  end subroutine rd
405  end interface
406 
407  integer :: min_value_
408 
409  procedure(rd), pointer :: read_data
410 
411  character(len=H_LONG) :: fname
412  integer, allocatable :: TILE_DATA(:,:)
413  integer :: jsize, isize
414  integer :: i, j, ii, jj, t
415 
416  if ( present(min_value) ) then
417  min_value_ = min_value
418  else
419  min_value_ = - abs(undef2)
420  end if
421 
422  select case( data_type )
423  case ( "int1", "INT1" )
424  read_data => file_tiledata_read_data_int1_int
425  case ( "int2", "INT2" )
426  read_data => file_tiledata_read_data_int2_int
427  case ( "int4", "INT4" )
428  read_data => file_tiledata_read_data_int4_int
429  case ( "real4", "REAL4" )
430  read_data => file_tiledata_read_data_real4_int
431  case default
432  log_error("FILE_TILEDATA_get_data_int1",*) 'data_type is invalid: ', trim(data_type)
433  call prc_abort
434  end select
435 
436  !$omp parallel do
437 !OCL XFILL
438  do j = 1, nlat
439  do i = 1, nlon
440  DATA(i,j) = undef2
441  end do
442  end do
443 
444  do t = 1, tile_nmax
445  if ( .not. tile_hit(t) ) cycle
446 
447  fname = trim(dirname) // '/' // trim(tile_fname(t))
448 
449  log_newline
450  log_info("FILE_TILEDATA_get_data_int1",*) 'Input data file :', trim(fname)
451  log_info_cont(*) 'Tile (LAT) :', tile_js(t)*tile_dlat/d2r, (tile_je(t)+1)*tile_dlat/d2r
452  log_info_cont(*) ' (LON) :', tile_is(t)*tile_dlon/d2r, (tile_ie(t)+1)*tile_dlon/d2r
453 
454  isize = tile_ie(t) - tile_is(t) + 1
455  jsize = tile_je(t) - tile_js(t) + 1
456 
457  allocate( tile_data(isize,jsize) )
458 
459  call read_data( jsize, isize, & ! [IN]
460  fname, & ! [IN]
461  tile_data(:,:), & ! [OUT]
462  yrevers = yrevers, & ! [IN]
463  step = step ) ! [IN]
464 
465  !$omp parallel do &
466  !$omp private(i,j)
467  do jj = 1, jsize
468  j = tile_js(t) + jj - 1
469  if ( jsh <= j .and. j <= jeh ) then
470  do ii = 1, isize
471  i = tile_is(t) + ii - 1
472  if ( ish <= i .and. i <= ieh ) then
473  if ( tile_data(ii,jj) < min_value_ ) then
474  DATA(i-ish+1,j-jsh+1) = undef2
475  else
476  DATA(i-ish+1,j-jsh+1) = tile_data(ii,jj)
477  end if
478  end if
479  i = i - global_ia
480  if ( ish <= i .and. i <= ieh ) then
481  if ( tile_data(ii,jj) < min_value_ ) then
482  DATA(i-ish+1,j-jsh+1) = undef2
483  else
484  DATA(i-ish+1,j-jsh+1) = tile_data(ii,jj)
485  end if
486  end if
487  end do
488  end if
489  end do
490 
491  deallocate( tile_data )
492 
493  enddo ! tile loop
494 
495  return

References scale_const::const_d2r, scale_const::const_pi, scale_const::const_undef, scale_const::const_undef2, file_tiledata_read_data_int1_int(), file_tiledata_read_data_int2_int(), file_tiledata_read_data_int4_int(), file_tiledata_read_data_real4_int(), and scale_prc::prc_abort().

Here is the call graph for this function:

◆ file_tiledata_read_catalog_file()

subroutine scale_file_tiledata::file_tiledata_read_catalog_file ( integer, intent(in)  TILE_nlim,
character(len=*), intent(in)  fname,
real(rp), intent(in)  TILE_DLAT,
real(rp), intent(in)  TILE_DLON,
integer, intent(in)  DOMAIN_IS,
integer, intent(in)  GLOBAL_IA,
integer, intent(out)  TILE_nmax,
character(len=*), dimension(:), intent(out)  TILE_fname,
real(rp), dimension(:), intent(out)  TILE_LATS,
real(rp), dimension(:), intent(out)  TILE_LATE,
real(rp), dimension(:), intent(out)  TILE_LONS,
real(rp), dimension(:), intent(out)  TILE_LONE 
)

read category file

Definition at line 513 of file scale_file_tiledata.F90.

513  use scale_const, only: &
514  d2r => const_d2r
515  use scale_prc, only: &
516  prc_abort
517  integer, intent(in) :: TILE_nlim
518  character(len=*), intent(in) :: fname
519  real(RP), intent(in) :: TILE_DLAT
520  real(RP), intent(in) :: TILE_DLON
521  integer, intent(in) :: DOMAIN_IS
522  integer, intent(in) :: GLOBAL_IA
523  integer, intent(out) :: TILE_nmax
524  character(len=*), intent(out) :: TILE_fname(:)
525  real(RP), intent(out) :: TILE_LATS(:)
526  real(RP), intent(out) :: TILE_LATE(:)
527  real(RP), intent(out) :: TILE_LONS(:)
528  real(RP), intent(out) :: TILE_LONE(:)
529 
530  integer :: fid, ierr
531  integer :: index
532  integer :: t
533 
534  fid = io_get_available_fid()
535  open( fid, &
536  file = fname, &
537  form = 'formatted', &
538  status = 'old', &
539  iostat = ierr )
540 
541  if ( ierr /= 0 ) then
542  log_error("FILE_TILEDATA_read_catalog_file",*) 'catalogue file not found! ', trim(fname)
543  call prc_abort
544  endif
545 
546  ierr = 0
547  tile_nmax = - 1
548  do while ( ierr == 0 )
549  read(fid,*,iostat=ierr) index, tile_lats(1), tile_late(1), tile_lons(1), tile_lone(1), & ! WEST->EAST
550  tile_fname(1)
551  tile_nmax = tile_nmax + 1
552  end do
553 
554  if ( tile_nmax > tile_nlim ) then
555  log_error('FILE_TILEDATA_read_catalog_file',*) 'TILE_nmax must be >= ', tile_nmax
556  call prc_abort
557  end if
558 
559  rewind(fid)
560  do t = 1, tile_nmax
561  read(fid,*,iostat=ierr) index, tile_lats(t), tile_late(t), & ! South->North
562  tile_lons(t), tile_lone(t), & ! WEST->EAST
563  tile_fname(t)
564 
565  tile_lats(t) = tile_lats(t) * d2r
566  tile_late(t) = tile_late(t) * d2r
567 
568  if ( tile_lons(t) > tile_lone(t) ) tile_lone(t) = tile_lone(t) + 360.0_rp
569  tile_lons(t) = tile_lons(t) * d2r
570  tile_lone(t) = tile_lone(t) * d2r
571  end do
572 
573  close(fid)
574 
575  return

References scale_const::const_d2r, scale_io::io_get_available_fid(), and scale_prc::prc_abort().

Referenced by file_tiledata_get_info().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ file_tiledata_get_domain_info()

subroutine scale_file_tiledata::file_tiledata_get_domain_info ( real(rp), intent(in)  TILE_DLAT,
real(rp), intent(in)  TILE_DLON,
real(rp), intent(in)  DOMAIN_LATS,
real(rp), intent(in)  DOMAIN_LATE,
real(rp), intent(in)  DOMAIN_LONS,
real(rp), intent(in)  DOMAIN_LONE,
integer, intent(out)  DOMAIN_JS,
integer, intent(out)  DOMAIN_JE,
integer, intent(out)  DOMAIN_IS,
integer, intent(out)  DOMAIN_IE,
integer, intent(out)  GLOBAL_IA 
)

Definition at line 583 of file scale_file_tiledata.F90.

583  use scale_const, only: &
584  pi => const_pi
585  real(RP), intent(in) :: TILE_DLAT
586  real(RP), intent(in) :: TILE_DLON
587  real(RP), intent(in) :: DOMAIN_LATS
588  real(RP), intent(in) :: DOMAIN_LATE
589  real(RP), intent(in) :: DOMAIN_LONS
590  real(RP), intent(in) :: DOMAIN_LONE
591  integer, intent(out) :: DOMAIN_JS
592  integer, intent(out) :: DOMAIN_JE
593  integer, intent(out) :: DOMAIN_IS
594  integer, intent(out) :: DOMAIN_IE
595  integer, intent(out) :: GLOBAL_IA
596 
597  domain_js = floor( domain_lats / tile_dlat )
598  domain_je = ceiling( domain_late / tile_dlat )
599  domain_is = floor( domain_lons / tile_dlon )
600  domain_ie = ceiling( domain_lone / tile_dlon )
601 
602  global_ia = nint( 2.0_rp * pi / tile_dlon )
603 
604  return

References scale_const::const_pi.

Referenced by file_tiledata_get_info().

Here is the caller graph for this function:

◆ file_tiledata_get_tile_info()

subroutine scale_file_tiledata::file_tiledata_get_tile_info ( integer, intent(in)  TILE_nmax,
integer, intent(in)  DOMAIN_JS,
integer, intent(in)  DOMAIN_JE,
integer, intent(in)  DOMAIN_IS,
integer, intent(in)  DOMAIN_IE,
integer, intent(in)  GLOBAL_IA,
real(rp), intent(in)  TILE_DLAT,
real(rp), intent(in)  TILE_DLON,
real(rp), dimension(:), intent(in)  TILE_LATS,
real(rp), dimension(:), intent(in)  TILE_LATE,
real(rp), dimension(:), intent(in)  TILE_LONS,
real(rp), dimension(:), intent(in)  TILE_LONE,
logical, intent(in)  zonal,
logical, dimension(:), intent(out)  TILE_hit,
integer, dimension(:), intent(out)  TILE_JS,
integer, dimension(:), intent(out)  TILE_JE,
integer, dimension(:), intent(out)  TILE_IS,
integer, dimension(:), intent(out)  TILE_IE,
integer, intent(out)  jsh,
integer, intent(out)  jeh,
integer, intent(out)  ish,
integer, intent(out)  ieh,
integer, intent(out)  nLAT,
integer, intent(out)  nLON 
)

Definition at line 621 of file scale_file_tiledata.F90.

621  use scale_const, only: &
622  pi => const_pi
623  integer, intent(in) :: TILE_nmax
624  integer, intent(in) :: DOMAIN_JS, DOMAIN_JE, DOMAIN_IS, DOMAIN_IE
625  integer, intent(in) :: GLOBAL_IA
626  real(RP), intent(in) :: TILE_DLAT, TILE_DLON
627  real(RP), intent(in) :: TILE_LATS(:), TILE_LATE(:), TILE_LONS(:), TILE_LONE(:)
628  logical, intent(in) :: zonal
629  logical, intent(out) :: TILE_hit(:)
630  integer, intent(out) :: TILE_JS(:), TILE_JE(:), TILE_IS(:), TILE_IE(:)
631  integer, intent(out) :: jsh, jeh, ish, ieh
632  integer, intent(out) :: nLAT, nLON
633 
634  logical :: hit_lat, hit_lon
635  integer :: nhalo
636  integer :: t
637 
638  nhalo = 2
639 
640  jsh = max( domain_js - nhalo, -floor( 0.5_rp * pi / tile_dlat ) )
641  jeh = min( domain_je + nhalo, floor( 0.5_rp * pi / tile_dlat ) )
642  ish = domain_is - nhalo
643  ieh = domain_ie + nhalo
644 
645  ! data file
646  !$omp parallel do &
647  !$omp private(hit_lat,hit_lon)
648  do t = 1, tile_nmax
649 
650  tile_js(t) = nint( tile_lats(t) / tile_dlat )
651  tile_je(t) = nint( tile_late(t) / tile_dlat ) - 1
652 
653  tile_is(t) = nint( tile_lons(t) / tile_dlon )
654  tile_ie(t) = nint( tile_lone(t) / tile_dlon ) - 1
655 
656  do while ( tile_ie(t) < ish )
657  tile_is(t) = tile_is(t) + global_ia
658  tile_ie(t) = tile_ie(t) + global_ia
659  end do
660  do while ( tile_is(t) - ish >= global_ia )
661  tile_is(t) = tile_is(t) - global_ia
662  tile_ie(t) = tile_ie(t) - global_ia
663  end do
664 
665  if ( ( jsh <= tile_js(t) .AND. tile_js(t) <= jeh ) &
666  .OR. ( jsh <= tile_je(t) .AND. tile_je(t) <= jeh ) &
667  .OR. ( tile_js(t) <= jsh .AND. jsh <= tile_je(t) ) &
668  .OR. ( tile_js(t) <= jeh .AND. jeh <= tile_je(t) ) ) then
669  hit_lat = .true.
670  else
671  hit_lat = .false.
672  endif
673 
674  if ( zonal ) then
675  hit_lon = .true.
676  else if ( ( tile_is(t) <= ieh ) &
677  .OR. ( tile_ie(t) - global_ia >= ish ) ) then
678  hit_lon = .true.
679  else
680  hit_lon = .false.
681  endif
682 
683  tile_hit(t) = ( hit_lat .AND. hit_lon )
684 
685  end do
686 
687  if ( zonal ) then
688  ish = minval(tile_is(1:tile_nmax))
689  ieh = maxval(tile_ie(1:tile_nmax))
690  jsh = max( jsh, minval(tile_js(1:tile_nmax)) )
691  jeh = min( jeh, maxval(tile_je(1:tile_nmax)) )
692  end if
693 
694  nlon = ieh - ish + 1
695  nlat = jeh - jsh + 1
696 
697  return

References scale_const::const_pi.

Referenced by file_tiledata_get_info().

Here is the caller graph for this function:

◆ file_tiledata_read_data_int2_real()

subroutine scale_file_tiledata::file_tiledata_read_data_int2_real ( integer, intent(in)  jsize,
integer, intent(in)  isize,
character(len=*), intent(in)  fname,
real(rp), dimension(isize,jsize), intent(out)  TILE_DATA,
logical, intent(in), optional  yrevers,
integer, intent(in), optional  step 
)

Definition at line 706 of file scale_file_tiledata.F90.

706  use scale_prc, only: &
707  prc_abort
708  integer, intent(in) :: jsize
709  integer, intent(in) :: isize
710  character(len=*), intent(in) :: fname
711  real(RP), intent(out) :: TILE_DATA(isize,jsize)
712 
713  logical, intent(in), optional :: yrevers
714  integer, intent(in), optional :: step
715 
716  integer(2) :: buf(isize,jsize)
717 
718  integer :: fid, ierr
719  logical :: yrevers_
720  logical :: lstep
721  integer :: step_
722  integer :: i, j
723 
724  if ( present(yrevers) ) then
725  yrevers_ = yrevers
726  else
727  yrevers_ = .false.
728  end if
729 
730  lstep = present(step)
731  if ( lstep ) then
732  step_ = step
733  inquire( file = fname, number = fid )
734  else
735  step_ = 1
736  fid = -1
737  end if
738 
739  if ( fid < 0 ) then
740  fid = io_get_available_fid()
741  open( fid, &
742  file = fname, &
743  form = 'unformatted', &
744  access = 'direct', &
745  status = 'old', &
746  recl = isize*jsize*2, &
747  iostat = ierr )
748  if ( ierr /= 0 ) then
749  log_error("FILE_TILEDATA_read_data_int2_real",*) 'data file not found!: ', trim(fname)
750  call prc_abort
751  endif
752  end if
753 
754 
755  read(fid,rec=step_) buf(:,:)
756  if ( .not. lstep ) close(fid)
757 
758  if ( yrevers_ ) then
759  !$omp parallel do
760 !OCL XFILL
761  do j = 1, jsize
762  do i = 1, isize
763  tile_data(i,j) = buf(i,jsize-j+1)
764  end do
765  end do
766  else
767  !$omp parallel do
768 !OCL XFILL
769  do j = 1, jsize
770  do i = 1, isize
771  tile_data(i,j) = buf(i,j)
772  end do
773  end do
774  end if
775 
776  return

References scale_io::io_get_available_fid(), and scale_prc::prc_abort().

Referenced by file_tiledata_get_latlon().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ file_tiledata_read_data_int4_real()

subroutine scale_file_tiledata::file_tiledata_read_data_int4_real ( integer, intent(in)  jsize,
integer, intent(in)  isize,
character(len=*), intent(in)  fname,
real(rp), dimension(isize,jsize), intent(out)  TILE_DATA,
logical, intent(in), optional  yrevers,
integer, intent(in), optional  step 
)

Definition at line 785 of file scale_file_tiledata.F90.

785  use scale_prc, only: &
786  prc_abort
787  integer, intent(in) :: jsize
788  integer, intent(in) :: isize
789  character(len=*), intent(in) :: fname
790  real(RP), intent(out) :: TILE_DATA(isize,jsize)
791 
792  logical, intent(in), optional :: yrevers
793  integer, intent(in), optional :: step
794 
795  integer(4) :: buf(isize,jsize)
796 
797  integer :: fid, ierr
798  logical :: yrevers_
799  logical :: lstep
800  integer :: step_
801  integer :: i, j
802 
803  if ( present(yrevers) ) then
804  yrevers_ = yrevers
805  else
806  yrevers_ = .false.
807  end if
808 
809  lstep = present(step)
810  if ( lstep ) then
811  step_ = step
812  inquire( file = fname, number = fid )
813  else
814  step_ = 1
815  fid = -1
816  end if
817 
818  if ( fid < 0 ) then
819  fid = io_get_available_fid()
820  open( fid, &
821  file = fname, &
822  form = 'unformatted', &
823  access = 'direct', &
824  status = 'old', &
825  recl = isize*jsize*4, &
826  iostat = ierr )
827  if ( ierr /= 0 ) then
828  log_error("FILE_TILEDATA_read_data_int4_real",*) 'data file not found!: ', trim(fname)
829  call prc_abort
830  end if
831  endif
832 
833  read(fid,rec=step_) buf(:,:)
834  if ( .not. lstep ) close(fid)
835 
836  if ( yrevers_ ) then
837  !$omp parallel do
838 !OCL XFILL
839  do j = 1, jsize
840  do i = 1, isize
841  tile_data(i,j) = buf(i,jsize-j+1)
842  end do
843  end do
844  else
845  !$omp parallel do
846 !OCL XFILL
847  do j = 1, jsize
848  do i = 1, isize
849  tile_data(i,j) = buf(i,j)
850  end do
851  end do
852  end if
853 
854  return

References scale_io::io_get_available_fid(), and scale_prc::prc_abort().

Referenced by file_tiledata_get_latlon().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ file_tiledata_read_data_real4_real()

subroutine scale_file_tiledata::file_tiledata_read_data_real4_real ( integer, intent(in)  jsize,
integer, intent(in)  isize,
character(len=*), intent(in)  fname,
real(rp), dimension(isize,jsize), intent(out)  TILE_DATA,
logical, intent(in), optional  yrevers,
integer, intent(in), optional  step 
)

Definition at line 863 of file scale_file_tiledata.F90.

863  use scale_prc, only: &
864  prc_abort
865  integer, intent(in) :: jsize
866  integer, intent(in) :: isize
867  character(len=*), intent(in) :: fname
868  real(RP), intent(out) :: TILE_DATA(isize,jsize)
869 
870  logical, intent(in), optional :: yrevers
871  integer, intent(in), optional :: step
872 
873  real(4) :: buf(isize,jsize)
874 
875  integer :: fid, ierr
876  logical :: yrevers_
877  logical :: lstep
878  integer :: step_
879  integer :: i, j
880 
881  if ( present(yrevers) ) then
882  yrevers_ = yrevers
883  else
884  yrevers_ = .false.
885  end if
886 
887  lstep = present(step)
888  if ( lstep ) then
889  step_ = step
890  inquire( file = fname, number = fid )
891  else
892  step_ = 1
893  fid = -1
894  end if
895 
896  if ( fid < 0 ) then
897  fid = io_get_available_fid()
898  open( fid, &
899  file = fname, &
900  form = 'unformatted', &
901  access = 'direct', &
902  status = 'old', &
903  recl = isize*jsize*4, &
904  iostat = ierr )
905  if ( ierr /= 0 ) then
906  log_error("FILE_TILEDATA_read_data_real4_real",*) 'data file not found!: ', trim(fname)
907  call prc_abort
908  endif
909  end if
910 
911  read(fid,rec=step_) buf(:,:)
912  if ( .not. lstep ) close(fid)
913 
914  if ( yrevers_ ) then
915  !$omp parallel do
916 !OCL XFILL
917  do j = 1, jsize
918  do i = 1, isize
919  tile_data(i,j) = buf(i,jsize-j+1)
920  end do
921  end do
922  else
923  !$omp parallel do
924 !OCL XFILL
925  do j = 1, jsize
926  do i = 1, isize
927  tile_data(i,j) = buf(i,j)
928  end do
929  end do
930  end if
931 
932  return

References scale_io::io_get_available_fid(), and scale_prc::prc_abort().

Referenced by file_tiledata_get_latlon().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ file_tiledata_read_data_real8_real()

subroutine scale_file_tiledata::file_tiledata_read_data_real8_real ( integer, intent(in)  jsize,
integer, intent(in)  isize,
character(len=*), intent(in)  fname,
real(rp), dimension(isize,jsize), intent(out)  TILE_DATA,
logical, intent(in), optional  yrevers,
integer, intent(in), optional  step 
)

Definition at line 941 of file scale_file_tiledata.F90.

941  use scale_prc, only: &
942  prc_abort
943  integer, intent(in) :: jsize
944  integer, intent(in) :: isize
945  character(len=*), intent(in) :: fname
946  real(RP), intent(out) :: TILE_DATA(isize,jsize)
947 
948  logical, intent(in), optional :: yrevers
949  integer, intent(in), optional :: step
950 
951  real(8) :: buf(isize,jsize)
952 
953  integer :: fid, ierr
954  logical :: yrevers_
955  logical :: lstep
956  integer :: step_
957  integer :: i, j
958 
959  if ( present(yrevers) ) then
960  yrevers_ = yrevers
961  else
962  yrevers_ = .false.
963  end if
964 
965  lstep = present(step)
966  if ( lstep ) then
967  step_ = step
968  inquire( file = fname, number = fid )
969  else
970  step_ = 1
971  fid = -1
972  end if
973 
974  if ( fid < 0 ) then
975  fid = io_get_available_fid()
976  open( fid, &
977  file = fname, &
978  form = 'unformatted', &
979  access = 'direct', &
980  status = 'old', &
981  recl = isize*jsize*8, &
982  iostat = ierr )
983  if ( ierr /= 0 ) then
984  log_error("FILE_TILEDATA_read_data_real8_real",*) 'data file not found!: ', trim(fname)
985  call prc_abort
986  endif
987  end if
988 
989  read(fid,rec=step_) buf(:,:)
990  if ( .not. lstep ) close(fid)
991 
992  if ( yrevers_ ) then
993  !$omp parallel do
994 !OCL XFILL
995  do j = 1, jsize
996  do i = 1, isize
997  tile_data(i,j) = buf(i,jsize-j+1)
998  end do
999  end do
1000  else
1001  !$omp parallel do
1002 !OCL XFILL
1003  do j = 1, jsize
1004  do i = 1, isize
1005  tile_data(i,j) = buf(i,j)
1006  end do
1007  end do
1008  end if
1009 
1010  return

References scale_io::io_get_available_fid(), and scale_prc::prc_abort().

Referenced by file_tiledata_get_latlon().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ file_tiledata_read_data_int1_int()

subroutine scale_file_tiledata::file_tiledata_read_data_int1_int ( integer, intent(in)  jsize,
integer, intent(in)  isize,
character(len=*), intent(in)  fname,
integer, dimension(isize,jsize), intent(out)  TILE_DATA,
logical, intent(in), optional  yrevers,
integer, intent(in), optional  step 
)

Definition at line 1019 of file scale_file_tiledata.F90.

1019  use scale_prc, only: &
1020  prc_abort
1021  integer, intent(in) :: jsize
1022  integer, intent(in) :: isize
1023  character(len=*), intent(in) :: fname
1024  integer, intent(out) :: TILE_DATA(isize,jsize)
1025 
1026  logical, intent(in), optional :: yrevers
1027  integer, intent(in), optional :: step
1028 
1029  integer(1) :: buf(isize,jsize)
1030 
1031  integer :: fid, ierr
1032  logical :: yrevers_
1033  logical :: lstep
1034  integer :: step_
1035  integer :: i, j
1036 
1037  if ( present(yrevers) ) then
1038  yrevers_ = yrevers
1039  else
1040  yrevers_ = .false.
1041  end if
1042 
1043  lstep = present(step)
1044  if ( lstep ) then
1045  step_ = step
1046  inquire( file = fname, number = fid )
1047  else
1048  step_ = 1
1049  fid = -1
1050  end if
1051 
1052  if ( fid < 0 ) then
1053  fid = io_get_available_fid()
1054  open( fid, &
1055  file = fname, &
1056  form = 'unformatted', &
1057  access = 'direct', &
1058  status = 'old', &
1059  recl = isize*jsize*1, &
1060  iostat = ierr )
1061  if ( ierr /= 0 ) then
1062  log_error("FILE_TILEDATA_read_data_int1_int",*) 'data file not found!: ', trim(fname)
1063  call prc_abort
1064  endif
1065  end if
1066 
1067  read(fid,rec=step_) buf(:,:)
1068  if ( .not. lstep ) close(fid)
1069 
1070  if ( yrevers_ ) then
1071  !$omp parallel do
1072 !OCL XFILL
1073  do j = 1, jsize
1074  do i = 1, isize
1075  tile_data(i,j) = buf(i,jsize-j+1)
1076  end do
1077  end do
1078  else
1079  !$omp parallel do
1080 !OCL XFILL
1081  do j = 1, jsize
1082  do i = 1, isize
1083  tile_data(i,j) = buf(i,j)
1084  end do
1085  end do
1086  end if
1087 
1088  return

References scale_io::io_get_available_fid(), and scale_prc::prc_abort().

Referenced by file_tiledata_get_data_int1().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ file_tiledata_read_data_int2_int()

subroutine scale_file_tiledata::file_tiledata_read_data_int2_int ( integer, intent(in)  jsize,
integer, intent(in)  isize,
character(len=*), intent(in)  fname,
integer, dimension(isize,jsize), intent(out)  TILE_DATA,
logical, intent(in), optional  yrevers,
integer, intent(in), optional  step 
)

Definition at line 1097 of file scale_file_tiledata.F90.

1097  use scale_prc, only: &
1098  prc_abort
1099  integer, intent(in) :: jsize
1100  integer, intent(in) :: isize
1101  character(len=*), intent(in) :: fname
1102  integer, intent(out) :: TILE_DATA(isize,jsize)
1103 
1104  logical, intent(in), optional :: yrevers
1105  integer, intent(in), optional :: step
1106 
1107  integer(2) :: buf(isize,jsize)
1108 
1109  integer :: fid, ierr
1110  logical :: yrevers_
1111  logical :: lstep
1112  integer :: step_
1113  integer :: i, j
1114 
1115  if ( present(yrevers) ) then
1116  yrevers_ = yrevers
1117  else
1118  yrevers_ = .false.
1119  end if
1120 
1121  lstep = present(step)
1122  if ( lstep ) then
1123  step_ = step
1124  inquire( file = fname, number = fid )
1125  else
1126  step_ = 1
1127  fid = -1
1128  end if
1129 
1130  if ( fid < 0 ) then
1131  fid = io_get_available_fid()
1132  open( fid, &
1133  file = fname, &
1134  form = 'unformatted', &
1135  access = 'direct', &
1136  status = 'old', &
1137  recl = isize*jsize*2, &
1138  iostat = ierr )
1139  if ( ierr /= 0 ) then
1140  log_error("FILE_TILEDATA_read_data_int2_int",*) 'data file not found!: ', trim(fname)
1141  call prc_abort
1142  endif
1143  end if
1144 
1145  read(fid,rec=step_) buf(:,:)
1146  if ( .not. lstep ) close(fid)
1147 
1148  if ( yrevers_ ) then
1149  !$omp parallel do
1150 !OCL XFILL
1151  do j = 1, jsize
1152  do i = 1, isize
1153  tile_data(i,j) = buf(i,jsize-j+1)
1154  end do
1155  end do
1156  else
1157  !$omp parallel do
1158 !OCL XFILL
1159  do j = 1, jsize
1160  do i = 1, isize
1161  tile_data(i,j) = buf(i,j)
1162  end do
1163  end do
1164  end if
1165 
1166  return

References scale_io::io_get_available_fid(), and scale_prc::prc_abort().

Referenced by file_tiledata_get_data_int1().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ file_tiledata_read_data_int4_int()

subroutine scale_file_tiledata::file_tiledata_read_data_int4_int ( integer, intent(in)  jsize,
integer, intent(in)  isize,
character(len=*), intent(in)  fname,
integer, dimension(isize,jsize), intent(out)  TILE_DATA,
logical, intent(in), optional  yrevers,
integer, intent(in), optional  step 
)

Definition at line 1175 of file scale_file_tiledata.F90.

1175  use scale_prc, only: &
1176  prc_abort
1177  integer, intent(in) :: jsize
1178  integer, intent(in) :: isize
1179  character(len=*), intent(in) :: fname
1180  integer, intent(out) :: TILE_DATA(isize,jsize)
1181 
1182  logical, intent(in), optional :: yrevers
1183  integer, intent(in), optional :: step
1184 
1185  integer(4) :: buf(isize,jsize)
1186 
1187  integer :: fid, ierr
1188  logical :: yrevers_
1189  logical :: lstep
1190  integer :: step_
1191  integer :: i, j
1192 
1193  if ( present(yrevers) ) then
1194  yrevers_ = yrevers
1195  else
1196  yrevers_ = .false.
1197  end if
1198 
1199  lstep = present(step)
1200  if ( lstep ) then
1201  step_ = step
1202  inquire( file = fname, number = fid )
1203  else
1204  step_ = 1
1205  fid = -1
1206  end if
1207 
1208  if ( fid < 0 ) then
1209  fid = io_get_available_fid()
1210  open( fid, &
1211  file = fname, &
1212  form = 'unformatted', &
1213  access = 'direct', &
1214  status = 'old', &
1215  recl = isize*jsize*4, &
1216  iostat = ierr )
1217  if ( ierr /= 0 ) then
1218  log_error("FILE_TILEDATA_read_data_int4_int",*) 'data file not found!: ', trim(fname)
1219  call prc_abort
1220  endif
1221  end if
1222 
1223  read(fid,rec=step_) buf(:,:)
1224  if ( .not. lstep ) close(fid)
1225 
1226  if ( yrevers_ ) then
1227  !$omp parallel do
1228 !OCL XFILL
1229  do j = 1, jsize
1230  do i = 1, isize
1231  tile_data(i,j) = buf(i,jsize-j+1)
1232  end do
1233  end do
1234  else
1235  !$omp parallel do
1236 !OCL XFILL
1237  do j = 1, jsize
1238  do i = 1, isize
1239  tile_data(i,j) = buf(i,j)
1240  end do
1241  end do
1242  end if
1243 
1244  return

References scale_io::io_get_available_fid(), and scale_prc::prc_abort().

Referenced by file_tiledata_get_data_int1().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ file_tiledata_read_data_real4_int()

subroutine scale_file_tiledata::file_tiledata_read_data_real4_int ( integer, intent(in)  jsize,
integer, intent(in)  isize,
character(len=*), intent(in)  fname,
integer, dimension(isize,jsize), intent(out)  TILE_DATA,
logical, intent(in), optional  yrevers,
integer, intent(in), optional  step 
)

Definition at line 1253 of file scale_file_tiledata.F90.

1253  use scale_prc, only: &
1254  prc_abort
1255  integer, intent(in) :: jsize
1256  integer, intent(in) :: isize
1257  character(len=*), intent(in) :: fname
1258  integer, intent(out) :: TILE_DATA(isize,jsize)
1259 
1260  logical, intent(in), optional :: yrevers
1261  integer, intent(in), optional :: step
1262 
1263  real(4) :: buf(isize,jsize)
1264 
1265  integer :: fid, ierr
1266  logical :: yrevers_
1267  logical :: lstep
1268  integer :: step_
1269  integer :: i, j
1270 
1271  if ( present(yrevers) ) then
1272  yrevers_ = yrevers
1273  else
1274  yrevers_ = .false.
1275  end if
1276 
1277  lstep = present(step)
1278  if ( lstep ) then
1279  step_ = step
1280  inquire( file = fname, number = fid )
1281  else
1282  step_ = 1
1283  fid = -1
1284  end if
1285 
1286  if ( fid < 0 ) then
1287  fid = io_get_available_fid()
1288  open( fid, &
1289  file = fname, &
1290  form = 'unformatted', &
1291  access = 'direct', &
1292  status = 'old', &
1293  recl = isize*jsize*4, &
1294  iostat = ierr )
1295  if ( ierr /= 0 ) then
1296  log_error("FILE_TILEDATA_read_data_real4_int",*) 'data file not found!: ', trim(fname)
1297  call prc_abort
1298  endif
1299  end if
1300 
1301  read(fid,rec=step_) buf(:,:)
1302  if ( .not. lstep ) close(fid)
1303 
1304  if ( yrevers_ ) then
1305  !$omp parallel do
1306 !OCL XFILL
1307  do j = 1, jsize
1308  do i = 1, isize
1309  tile_data(i,j) = buf(i,jsize-j+1)
1310  end do
1311  end do
1312  else
1313  !$omp parallel do
1314 !OCL XFILL
1315  do j = 1, jsize
1316  do i = 1, isize
1317  tile_data(i,j) = buf(i,j)
1318  end do
1319  end do
1320  end if
1321 
1322  return

References scale_io::io_get_available_fid(), and scale_prc::prc_abort().

Referenced by file_tiledata_get_data_int1().

Here is the call graph for this function:
Here is the caller graph for this function:
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_const
module CONSTANT
Definition: scale_const.F90:11