SCALE-RM
scale_file_tiledata.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 #include "scalelib.h"
13  !-----------------------------------------------------------------------------
14  !
15  !++ used modules
16  !
17  use scale_precision
18  use scale_io
19  !-----------------------------------------------------------------------------
20  implicit none
21  private
22  !-----------------------------------------------------------------------------
23  !
24  !++ Public procedure
25  !
27  public file_tiledata_get_data
28 
29  interface file_tiledata_get_data
30  module procedure file_tiledata_get_data_real
31  module procedure file_tiledata_get_data_int1
32  end interface file_tiledata_get_data
33 
34  !-----------------------------------------------------------------------------
35  !
36  !++ Public parameters & variables
37  !
38  !-----------------------------------------------------------------------------
39  !
40  !++ Private procedure
41  !
42  !-----------------------------------------------------------------------------
43  !
44  !++ Private parameters & variables
45  !
46  !-----------------------------------------------------------------------------
47 contains
48  !-----------------------------------------------------------------------------
50  subroutine file_tiledata_get_info( &
51  TILE_nlim, &
52  TILE_DLAT, TILE_DLON, &
53  DOMAIN_LATS, DOMAIN_LATE, DOMAIN_LONS, DOMAIN_LONE, &
54  catalog_fname, &
55  GLOBAL_IA, &
56  TILE_nmax, &
57  TILE_fname, TILE_hit, &
58  TILE_JS, TILE_JE, TILE_IS, TILE_IE, &
59  nLATH, nLONH, jsh, jeh, ish, ieh, &
60  single_fname, LATS, LATE, LONS, LONE )
61  use scale_prc, only: &
62  prc_abort
63  integer, intent(in) :: TILE_nlim
64  real(RP), intent(in) :: TILE_DLAT, TILE_DLON
65  real(RP), intent(in) :: DOMAIN_LATS, DOMAIN_LATE, DOMAIN_LONS, DOMAIN_LONE
66  character(len=*), intent(in) :: catalog_fname
67  integer, intent(out) :: GLOBAL_IA
68  integer, intent(out) :: TILE_nmax
69  logical, intent(out) :: TILE_hit(:)
70  character(len=*), intent(out) :: TILE_fname(:)
71  integer, intent(out) :: TILE_JS(:), TILE_JE(:), TILE_IS(:), TILE_IE(:)
72  integer, intent(out) :: nLATH, nLONH
73  integer, intent(out) :: jsh, jeh, ish, ieh
74 
75  character(len=*), intent(in), optional :: single_fname
76  real(RP), intent(in), optional :: LATS
77  real(RP), intent(in), optional :: LATE
78  real(RP), intent(in), optional :: LONS
79  real(RP), intent(in), optional :: LONE
80 
81  real(RP) :: TILE_LATS(tile_nlim), TILE_LATE(tile_nlim)
82  real(RP) :: TILE_LONS(tile_nlim), TILE_LONE(tile_nlim)
83 
84  integer :: DOMAIN_JS, DOMAIN_JE, DOMAIN_IS, DOMAIN_IE
85 
86  call file_tiledata_get_domain_info( tile_dlat, tile_dlon, & ! [IN]
87  domain_lats, domain_late, domain_lons, domain_lone, & ! [IN]
88  domain_js, domain_je, domain_is, domain_ie, & ! [OUT]
89  global_ia ) ! [OUT]
90 
91 
92  if ( catalog_fname /= "" ) then
93  log_info("FILE_TILEDATA_get_info",*) 'Input catalogue file:', trim(catalog_fname)
94 
95  call file_tiledata_read_catalog_file( tile_nlim, & ! [IN]
96  catalog_fname, & ! [IN]
97  tile_dlat, tile_dlon, & ! [IN]
98  domain_is, global_ia, & ! [IN]
99  tile_nmax, & ! [OUT]
100  tile_fname(:), & ! [OUT]
101  tile_lats(:), tile_late(:), & ! [OUT]
102  tile_lons(:), tile_lone(:) ) ! [OUT]
103 
104  else
105  if ( .not. present(single_fname) ) then
106  log_error("FILE_TILEDATA_get_info",*) "single_fname is required if catalog_fname is empty"
107  call prc_abort
108  end if
109  if ( .not. present(lats) ) then
110  log_error("FILE_TILEDATA_get_info",*) "LATS is required if catalog_fname is empty"
111  call prc_abort
112  end if
113  if ( .not. present(late) ) then
114  log_error("FILE_TILEDATA_get_info",*) "LATE is required if catalog_fname is empty"
115  call prc_abort
116  end if
117  if ( .not. present(lons) ) then
118  log_error("FILE_TILEDATA_get_info",*) "LONS is required if catalog_fname is empty"
119  call prc_abort
120  end if
121  if ( .not. present(lone) ) then
122  log_error("FILE_TILEDATA_get_info",*) "LONE is required if catalog_fname is empty"
123  call prc_abort
124  end if
125 
126  tile_nmax = 1
127  tile_fname(1) = single_fname
128  tile_lats(1) = lats
129  tile_late(1) = late
130  tile_lons(1) = lons
131  tile_lone(1) = lone
132  end if
133 
134  call file_tiledata_get_tile_info( tile_nmax, & ! [IN]
135  domain_js, domain_je, & ! [IN]
136  domain_is, domain_ie, & ! [IN]
137  global_ia, & ! [IN]
138  tile_dlat, tile_dlon, & ! [IN]
139  tile_lats(:), tile_late(:), & ! [IN]
140  tile_lons(:), tile_lone(:), & ! [IN]
141  tile_hit(:), & ! [OUT]
142  tile_js(:), tile_je(:), & ! [OUT]
143  tile_is(:), tile_ie(:), & ! [OUT]
144  jsh, jeh, ish, ieh, & ! [OUT]
145  nlath, nlonh ) ! [OUT]
146 
147  return
148  end subroutine file_tiledata_get_info
149 
150  !-----------------------------------------------------------------------------
152  subroutine file_tiledata_get_data_real( &
153  nLATH, nLONH, &
154  dirname, &
155  GLOBAL_IA, &
156  TILE_nmax, &
157  TILE_DLAT, TILE_DLON, &
158  TILE_fname, TILE_hit, &
159  TILE_JS, TILE_JE, &
160  TILE_IS, TILE_IE, &
161  jsh, jeh, ish, ieh, &
162  data_type, &
163  DATA, LATH, LONH, &
164  min_value, &
165  yrevers )
166  use scale_const, only: &
167  undef => const_undef, &
168  pi => const_pi, &
169  d2r => const_d2r
170  use scale_prc, only: &
171  prc_abort
172  integer, intent(in) :: nLATH
173  integer, intent(in) :: nLONH
174  character(len=*), intent(in) :: dirname
175  integer, intent(in) :: GLOBAL_IA
176  integer, intent(in) :: TILE_nmax
177  real(RP), intent(in) :: TILE_DLAT
178  real(RP), intent(in) :: TILE_DLON
179  character(len=*), intent(in) :: TILE_fname(:)
180  logical, intent(in) :: TILE_hit(:)
181  integer, intent(in) :: TILE_JS(:)
182  integer, intent(in) :: TILE_JE(:)
183  integer, intent(in) :: TILE_IS(:)
184  integer, intent(in) :: TILE_IE(:)
185  integer, intent(in) :: jsh
186  integer, intent(in) :: jeh
187  integer, intent(in) :: ish
188  integer, intent(in) :: ieh
189  character(len=*), intent(in) :: data_type
190  real(RP), intent(out) :: DATA(nlonh,nlath)
191  real(RP), intent(out) :: LATH(nlonh,nlath)
192  real(RP), intent(out) :: LONH(nlonh,nlath)
193  real(RP), intent(in), optional :: min_value
194  logical, intent(in), optional :: yrevers
195 
196  abstract interface
197  subroutine rd( &
198  jsize, isize, &
199  fname, &
200  TILE_DATA, &
201  yrevers )
202  use scale_precision
203  integer, intent(in) :: jsize
204  integer, intent(in) :: isize
205  character(len=*), intent(in) :: fname
206  real(RP), intent(out) :: TILE_DATA(isize,jsize)
207  logical, intent(in), optional :: yrevers
208  end subroutine rd
209  end interface
210 
211  procedure(rd), pointer :: read_data
212 
213  real(RP) :: min_value_
214 
215  character(len=H_LONG) :: fname
216  real(RP), allocatable :: TILE_DATA(:,:)
217  integer :: jsize, isize
218  integer :: i, j, ii, jj, t
219 
220  if ( present(min_value) ) then
221  min_value_ = min_value
222  else
223  min_value_ = - abs(undef)
224  end if
225 
226  select case( data_type )
227  case ( "int2", "INT2" )
229  case ( "int4", "INT4" )
231  case ( "real4", "REAL4" )
233  case ( "real8", "REAL8" )
235  case default
236  log_error("FILE_TILEDATA_get_data_real",*) 'data_type is invalid: ', trim(data_type)
237  call prc_abort
238  end select
239 
240  !$omp parallel do
241 !OCL XFILL
242  do j = 1, nlath
243  do i = 1, nlonh
244  DATA(i,j) = undef
245  lath(i,j) = undef
246  lonh(i,j) = undef
247  end do
248  end do
249 
250  do t = 1, tile_nmax
251  if ( .not. tile_hit(t) ) cycle
252 
253  fname = trim(dirname) // '/' // trim(tile_fname(t))
254 
255  log_newline
256  log_info("FILE_TILEDATA_get_data_real",*) 'Input data file :', trim(fname)
257  log_info_cont(*) 'Tile (LAT) :', tile_js(t)*tile_dlat/d2r, (tile_je(t)+1)*tile_dlat/d2r
258  log_info_cont(*) ' (LON) :', tile_is(t)*tile_dlon/d2r, (tile_ie(t)+1)*tile_dlon/d2r
259 
260  isize = tile_ie(t) - tile_is(t) + 1
261  jsize = tile_je(t) - tile_js(t) + 1
262 
263  allocate( tile_data(isize,jsize) )
264 
265  call read_data( jsize, isize, & ! [IN]
266  fname, & ! [IN]
267  tile_data(:,:), & ! [OUT]
268  yrevers = yrevers ) ! [IN]
269 
270  !$omp parallel do &
271  !$omp private(i,j)
272  do jj = 1, jsize
273  j = tile_js(t) + jj - 1
274  if ( jsh <= j .and. j <= jeh ) then
275  do ii = 1, isize
276  i = tile_is(t) + ii - 1
277  if ( ish <= i .and. i <= ieh ) then
278  if ( tile_data(ii,jj) < min_value_ ) then
279  DATA(i-ish+1,j-jsh+1) = undef
280  else
281  DATA(i-ish+1,j-jsh+1) = tile_data(ii,jj)
282  end if
283  lath(i-ish+1,j-jsh+1) = tile_dlat * ( tile_js(t) + jj - 1 + 0.5_rp )
284  lonh(i-ish+1,j-jsh+1) = tile_dlon * ( tile_is(t) + ii - 1 + 0.5_rp )
285  end if
286  i = i - global_ia
287  if ( ish <= i .and. i <= ieh ) then
288  if ( tile_data(ii,jj) < min_value_ ) then
289  DATA(i-ish+1,j-jsh+1) = undef
290  else
291  DATA(i-ish+1,j-jsh+1) = tile_data(ii,jj)
292  end if
293  lath(i-ish+1,j-jsh+1) = tile_dlat * ( tile_js(t) + jj - 1 + 0.5_rp )
294  lonh(i-ish+1,j-jsh+1) = tile_dlon * ( tile_is(t) + ii - 1 + 0.5_rp ) - 2.0 * pi
295  end if
296  end do
297  end if
298  end do
299 
300  deallocate( tile_data )
301 
302  enddo ! tile loop
303 
304  return
305  end subroutine file_tiledata_get_data_real
306 
307  subroutine file_tiledata_get_data_int1( &
308  nLATH, nLONH, &
309  dirname, &
310  GLOBAL_IA, &
311  TILE_nmax, &
312  TILE_DLAT, TILE_DLON, &
313  TILE_fname, TILE_hit, &
314  TILE_JS, TILE_JE, &
315  TILE_IS, TILE_IE, &
316  jsh, jeh, ish, ieh, &
317  data_type, &
318  DATA, LATH, LONH, &
319  min_value, &
320  yrevers )
321  use scale_const, only: &
322  undef2 => const_undef2, &
323  pi => const_pi, &
324  d2r => const_d2r
325  use scale_prc, only: &
326  prc_abort
327  integer, intent(in) :: nLATH
328  integer, intent(in) :: nLONH
329  character(len=*), intent(in) :: dirname
330  integer, intent(in) :: GLOBAL_IA
331  integer, intent(in) :: TILE_nmax
332  real(RP), intent(in) :: TILE_DLAT
333  real(RP), intent(in) :: TILE_DLON
334  character(len=*), intent(in) :: TILE_fname(:)
335  logical, intent(in) :: TILE_hit(:)
336  integer, intent(in) :: TILE_JS(:)
337  integer, intent(in) :: TILE_JE(:)
338  integer, intent(in) :: TILE_IS(:)
339  integer, intent(in) :: TILE_IE(:)
340  integer, intent(in) :: jsh
341  integer, intent(in) :: jeh
342  integer, intent(in) :: ish
343  integer, intent(in) :: ieh
344  character(len=*), intent(in) :: data_type
345  integer, intent(out) :: DATA(nlonh,nlath)
346  real(RP), intent(out) :: LATH(nlonh,nlath)
347  real(RP), intent(out) :: LONH(nlonh,nlath)
348 
349  integer, intent(in), optional :: min_value
350  logical, intent(in), optional :: yrevers
351 
352  abstract interface
353  subroutine rd( &
354  jsize, isize, &
355  fname, &
356  TILE_DATA, &
357  yrevers )
358  use scale_precision
359  integer, intent(in) :: jsize
360  integer, intent(in) :: isize
361  character(len=*), intent(in) :: fname
362  integer, intent(out) :: TILE_DATA(isize,jsize)
363  logical, intent(in), optional :: yrevers
364  end subroutine rd
365  end interface
366 
367  integer :: min_value_
368 
369  procedure(rd), pointer :: read_data
370 
371  character(len=H_LONG) :: fname
372  integer, allocatable :: TILE_DATA(:,:)
373  integer :: jsize, isize
374  integer :: i, j, ii, jj, t
375 
376  if ( present(min_value) ) then
377  min_value_ = min_value
378  else
379  min_value_ = - abs(undef2)
380  end if
381 
382  select case( data_type )
383  case ( "int1", "INT1" )
385  case ( "int2", "INT2" )
387  case ( "int4", "INT4" )
389  case ( "real4", "REAL4" )
391  case default
392  log_error("FILE_TILEDATA_get_data_int1",*) 'data_type is invalid: ', trim(data_type)
393  call prc_abort
394  end select
395 
396  !$omp parallel do
397 !OCL XFILL
398  do j = 1, nlath
399  do i = 1, nlonh
400  DATA(i,j) = - 1
401  end do
402  end do
403 
404  do t = 1, tile_nmax
405  if ( .not. tile_hit(t) ) cycle
406 
407  fname = trim(dirname) // '/' // trim(tile_fname(t))
408 
409  log_newline
410  log_info("FILE_TILEDATA_get_data_int1",*) 'Input data file :', trim(fname)
411  log_info_cont(*) 'Tile (LAT) :', tile_js(t)*tile_dlat/d2r, (tile_je(t)+1)*tile_dlat/d2r
412  log_info_cont(*) ' (LON) :', tile_is(t)*tile_dlon/d2r, (tile_ie(t)+1)*tile_dlon/d2r
413 
414  isize = tile_ie(t) - tile_is(t) + 1
415  jsize = tile_je(t) - tile_js(t) + 1
416 
417  allocate( tile_data(isize,jsize) )
418 
419  call read_data( jsize, isize, & ! [IN]
420  fname, & ! [IN]
421  tile_data(:,:), & ! [OUT]
422  yrevers = yrevers ) ! [IN]
423 
424  !$omp parallel do &
425  !$omp private(i,j)
426  do jj = 1, jsize
427  do ii = 1, isize
428  i = tile_is(t) + ii - 1
429  j = tile_js(t) + jj - 1
430  if ( jsh <= j .and. j <= jeh ) then
431  if ( ish <= i .and. i <= ieh ) then
432  if ( tile_data(ii,jj) < min_value_ ) then
433  DATA(i-ish+1,j-jsh+1) = undef2
434  else
435  DATA(i-ish+1,j-jsh+1) = tile_data(ii,jj)
436  end if
437  lath(i-ish+1,j-jsh+1) = tile_dlat * ( tile_js(t) + jj - 1 + 0.5_rp )
438  lonh(i-ish+1,j-jsh+1) = tile_dlon * ( tile_is(t) + ii - 1 + 0.5_rp )
439  end if
440  i = i - global_ia
441  if ( ish <= i .and. i <= ieh ) then
442  if ( tile_data(ii,jj) < min_value_ ) then
443  DATA(i-ish+1,j-jsh+1) = undef2
444  else
445  DATA(i-ish+1,j-jsh+1) = tile_data(ii,jj)
446  end if
447  lath(i-ish+1,j-jsh+1) = tile_dlat * ( tile_js(t) + jj - 1 + 0.5_rp )
448  lonh(i-ish+1,j-jsh+1) = tile_dlon * ( tile_is(t) + ii - 1 + 0.5_rp ) - 2.0 * pi
449  end if
450  end if
451  end do
452  end do
453 
454  deallocate( tile_data )
455 
456  enddo ! tile loop
457 
458  return
459  end subroutine file_tiledata_get_data_int1
460 
461 
462  ! private
463 
464 
465  !-----------------------------------------------------------------------------
467  subroutine file_tiledata_read_catalog_file( &
468  TILE_nlim, &
469  fname, &
470  TILE_DLAT, TILE_DLON, &
471  DOMAIN_IS, GLOBAL_IA, &
472  TILE_nmax, &
473  TILE_fname, &
474  TILE_LATS, TILE_LATE, &
475  TILE_LONS, TILE_LONE )
476  use scale_const, only: &
477  d2r => const_d2r
478  use scale_prc, only: &
479  prc_abort
480  integer, intent(in) :: TILE_nlim
481  character(len=*), intent(in) :: fname
482  real(RP), intent(in) :: TILE_DLAT
483  real(RP), intent(in) :: TILE_DLON
484  integer, intent(in) :: DOMAIN_IS
485  integer, intent(in) :: GLOBAL_IA
486  integer, intent(out) :: TILE_nmax
487  character(len=*), intent(out) :: TILE_fname(:)
488  real(RP), intent(out) :: TILE_LATS(:)
489  real(RP), intent(out) :: TILE_LATE(:)
490  real(RP), intent(out) :: TILE_LONS(:)
491  real(RP), intent(out) :: TILE_LONE(:)
492 
493  integer :: fid, ierr
494  integer :: index
495  integer :: t
496 
497  fid = io_get_available_fid()
498  open( fid, &
499  file = fname, &
500  form = 'formatted', &
501  status = 'old', &
502  iostat = ierr )
503 
504  if ( ierr /= 0 ) then
505  log_error("FILE_TILEDATA_read_catalog_file",*) 'catalogue file not found! ', trim(fname)
506  call prc_abort
507  endif
508 
509  do t = 1, tile_nlim
510  read(fid,*,iostat=ierr) index, tile_lats(t), tile_late(t), & ! South->North
511  tile_lons(t), tile_lone(t), & ! WEST->EAST
512  tile_fname(t)
513 
514  if ( ierr /= 0 ) exit
515 
516  tile_lats(t) = tile_lats(t) * d2r
517  tile_late(t) = tile_late(t) * d2r
518 
519  if ( tile_lons(t) > tile_lone(t) ) tile_lone(t) = tile_lone(t) + 360.0_rp
520  tile_lons(t) = tile_lons(t) * d2r
521  tile_lone(t) = tile_lone(t) * d2r
522  end do
523 
524  tile_nmax = t - 1
525 
526  close(fid)
527 
528  return
529  end subroutine file_tiledata_read_catalog_file
530 
531  subroutine file_tiledata_get_domain_info( &
532  TILE_DLAT, TILE_DLON, &
533  DOMAIN_LATS, DOMAIN_LATE, DOMAIN_LONS, DOMAIN_LONE, &
534  DOMAIN_JS, DOMAIN_JE, DOMAIN_IS, DOMAIN_IE, &
535  GLOBAL_IA )
536  use scale_const, only: &
537  pi => const_pi
538  real(RP), intent(in) :: TILE_DLAT
539  real(RP), intent(in) :: TILE_DLON
540  real(RP), intent(in) :: DOMAIN_LATS
541  real(RP), intent(in) :: DOMAIN_LATE
542  real(RP), intent(in) :: DOMAIN_LONS
543  real(RP), intent(in) :: DOMAIN_LONE
544  integer, intent(out) :: DOMAIN_JS
545  integer, intent(out) :: DOMAIN_JE
546  integer, intent(out) :: DOMAIN_IS
547  integer, intent(out) :: DOMAIN_IE
548  integer, intent(out) :: GLOBAL_IA
549 
550  domain_js = floor( domain_lats / tile_dlat )
551  domain_je = ceiling( domain_late / tile_dlat )
552  domain_is = floor( domain_lons / tile_dlon )
553  domain_ie = ceiling( domain_lone / tile_dlon )
554 
555  global_ia = nint( 2.0_rp * pi / tile_dlon )
556 
557  return
558  end subroutine file_tiledata_get_domain_info
559 
560  subroutine file_tiledata_get_tile_info( &
561  TILE_nmax, &
562  DOMAIN_JS, DOMAIN_JE, &
563  DOMAIN_IS, DOMAIN_IE, &
564  GLOBAL_IA, &
565  TILE_DLAT, TILE_DLON, &
566  TILE_LATS, TILE_LATE, &
567  TILE_LONS, TILE_LONE, &
568  TILE_hit, &
569  TILE_JS, TILE_JE, &
570  TILE_IS, TILE_IE, &
571  jsh, jeh, ish, ieh, &
572  nLATH, nLONH )
573  use scale_const, only: &
574  pi => const_pi
575  integer, intent(in) :: TILE_nmax
576  integer, intent(in) :: DOMAIN_JS, DOMAIN_JE, DOMAIN_IS, DOMAIN_IE
577  integer, intent(in) :: GLOBAL_IA
578  real(RP), intent(in) :: TILE_DLAT, TILE_DLON
579  real(RP), intent(in) :: TILE_LATS(:), TILE_LATE(:), TILE_LONS(:), TILE_LONE(:)
580  logical, intent(out) :: TILE_hit(:)
581  integer, intent(out) :: TILE_JS(:), TILE_JE(:), TILE_IS(:), TILE_IE(:)
582  integer, intent(out) :: jsh, jeh, ish, ieh
583  integer, intent(out) :: nLATH, nLONH
584 
585  logical :: hit_lat, hit_lon
586  integer :: nhalo
587  integer :: t
588 
589  nhalo = 2
590 
591  jsh = max( domain_js - nhalo, -floor( 0.5_rp * pi / tile_dlat ) )
592  jeh = min( domain_je + nhalo, floor( 0.5_rp * pi / tile_dlat ) )
593  ish = domain_is - nhalo
594  ieh = domain_ie + nhalo
595 
596  nlonh = ieh - ish + 1
597  nlath = jeh - jsh + 1
598 
599  ! data file
600  !$omp parallel do &
601  !$omp private(hit_lat,hit_lon)
602  do t = 1, tile_nmax
603 
604  tile_js(t) = nint( tile_lats(t) / tile_dlat )
605  tile_je(t) = nint( tile_late(t) / tile_dlat ) - 1
606 
607  tile_is(t) = nint( tile_lons(t) / tile_dlon )
608  tile_ie(t) = nint( tile_lone(t) / tile_dlon ) - 1
609 
610  do while ( tile_ie(t) < domain_is )
611  tile_is(t) = tile_is(t) + global_ia
612  tile_ie(t) = tile_ie(t) + global_ia
613  end do
614  do while ( tile_is(t) - domain_is >= global_ia )
615  tile_is(t) = tile_is(t) - global_ia
616  tile_ie(t) = tile_ie(t) - global_ia
617  end do
618 
619  if ( ( jsh <= tile_js(t) .AND. tile_js(t) <= jeh ) &
620  .OR. ( jsh <= tile_je(t) .AND. tile_je(t) <= jeh ) &
621  .OR. ( tile_js(t) <= jsh .AND. jsh <= tile_je(t) ) &
622  .OR. ( tile_js(t) <= jeh .AND. jeh <= tile_je(t) ) ) then
623  hit_lat = .true.
624  else
625  hit_lat = .false.
626  endif
627 
628  if ( ( tile_is(t) <= ieh ) &
629  .OR. ( ish <= tile_ie(t) - global_ia ) ) then
630  hit_lon = .true.
631  else
632  hit_lon = .false.
633  endif
634 
635  tile_hit(t) = ( hit_lat .AND. hit_lon )
636  end do
637 
638  return
639  end subroutine file_tiledata_get_tile_info
640 
642  jsize, isize, &
643  fname, &
644  TILE_DATA, &
645  yrevers )
646  use scale_prc, only: &
647  prc_abort
648  integer, intent(in) :: jsize
649  integer, intent(in) :: isize
650  character(len=*), intent(in) :: fname
651  real(RP), intent(out) :: TILE_DATA(isize,jsize)
652 
653  logical, intent(in), optional :: yrevers
654 
655  integer(2) :: buf(isize,jsize)
656 
657  integer :: fid, ierr
658  logical :: yrevers_
659  integer :: i, j
660 
661  if ( present(yrevers) ) then
662  yrevers_ = yrevers
663  else
664  yrevers_ = .false.
665  end if
666 
667  fid = io_get_available_fid()
668  open( fid, &
669  file = fname, &
670  form = 'unformatted', &
671  access = 'direct', &
672  status = 'old', &
673  recl = isize*jsize*2, &
674  iostat = ierr )
675 
676  if ( ierr /= 0 ) then
677  log_error("FILE_TILEDATA_read_data_int2_real",*) 'data file not found!: ', trim(fname)
678  call prc_abort
679  endif
680 
681  read(fid,rec=1) buf(:,:)
682  close(fid)
683 
684  if ( yrevers_ ) then
685  !$omp parallel do
686 !OCL XFILL
687  do j = 1, jsize
688  do i = 1, isize
689  tile_data(i,j) = buf(i,jsize-j+1)
690  end do
691  end do
692  else
693  !$omp parallel do
694 !OCL XFILL
695  do j = 1, jsize
696  do i = 1, isize
697  tile_data(i,j) = buf(i,j)
698  end do
699  end do
700  end if
701 
702  return
703  end subroutine file_tiledata_read_data_int2_real
704 
706  jsize, isize, &
707  fname, &
708  TILE_DATA, &
709  yrevers )
710  use scale_prc, only: &
711  prc_abort
712  integer, intent(in) :: jsize
713  integer, intent(in) :: isize
714  character(len=*), intent(in) :: fname
715  real(RP), intent(out) :: TILE_DATA(isize,jsize)
716 
717  logical, intent(in), optional :: yrevers
718 
719  integer(4) :: buf(isize,jsize)
720 
721  integer :: fid, ierr
722  logical :: yrevers_
723  integer :: i, j
724 
725  if ( present(yrevers) ) then
726  yrevers_ = yrevers
727  else
728  yrevers_ = .false.
729  end if
730 
731  fid = io_get_available_fid()
732  open( fid, &
733  file = fname, &
734  form = 'unformatted', &
735  access = 'direct', &
736  status = 'old', &
737  recl = isize*jsize*4, &
738  iostat = ierr )
739 
740  if ( ierr /= 0 ) then
741  log_error("FILE_TILEDATA_read_data_int4_real",*) 'data file not found!: ', trim(fname)
742  call prc_abort
743  endif
744 
745  read(fid,rec=1) buf(:,:)
746  close(fid)
747 
748  if ( yrevers_ ) then
749  !$omp parallel do
750 !OCL XFILL
751  do j = 1, jsize
752  do i = 1, isize
753  tile_data(i,j) = buf(i,jsize-j+1)
754  end do
755  end do
756  else
757  !$omp parallel do
758 !OCL XFILL
759  do j = 1, jsize
760  do i = 1, isize
761  tile_data(i,j) = buf(i,j)
762  end do
763  end do
764  end if
765 
766  return
767  end subroutine file_tiledata_read_data_int4_real
768 
770  jsize, isize, &
771  fname, &
772  TILE_DATA, &
773  yrevers )
774  use scale_prc, only: &
775  prc_abort
776  integer, intent(in) :: jsize
777  integer, intent(in) :: isize
778  character(len=*), intent(in) :: fname
779  real(RP), intent(out) :: TILE_DATA(isize,jsize)
780 
781  logical, intent(in), optional :: yrevers
782 
783  real(4) :: buf(isize,jsize)
784 
785  integer :: fid, ierr
786  logical :: yrevers_
787  integer :: i, j
788 
789  if ( present(yrevers) ) then
790  yrevers_ = yrevers
791  else
792  yrevers_ = .false.
793  end if
794 
795  fid = io_get_available_fid()
796  open( fid, &
797  file = fname, &
798  form = 'unformatted', &
799  access = 'direct', &
800  status = 'old', &
801  recl = isize*jsize*4, &
802  iostat = ierr )
803 
804  if ( ierr /= 0 ) then
805  log_error("FILE_TILEDATA_read_data_real4_real",*) 'data file not found!: ', trim(fname)
806  call prc_abort
807  endif
808 
809  read(fid,rec=1) buf(:,:)
810  close(fid)
811 
812  if ( yrevers_ ) then
813  !$omp parallel do
814 !OCL XFILL
815  do j = 1, jsize
816  do i = 1, isize
817  tile_data(i,j) = buf(i,jsize-j+1)
818  end do
819  end do
820  else
821  !$omp parallel do
822 !OCL XFILL
823  do j = 1, jsize
824  do i = 1, isize
825  tile_data(i,j) = buf(i,j)
826  end do
827  end do
828  end if
829 
830  return
832 
834  jsize, isize, &
835  fname, &
836  TILE_DATA, &
837  yrevers )
838  use scale_prc, only: &
839  prc_abort
840  integer, intent(in) :: jsize
841  integer, intent(in) :: isize
842  character(len=*), intent(in) :: fname
843  real(RP), intent(out) :: TILE_DATA(isize,jsize)
844 
845  logical, intent(in), optional :: yrevers
846 
847  real(8) :: buf(isize,jsize)
848 
849  integer :: fid, ierr
850  logical :: yrevers_
851  integer :: i, j
852 
853  if ( present(yrevers) ) then
854  yrevers_ = yrevers
855  else
856  yrevers_ = .false.
857  end if
858 
859  fid = io_get_available_fid()
860  open( fid, &
861  file = fname, &
862  form = 'unformatted', &
863  access = 'direct', &
864  status = 'old', &
865  recl = isize*jsize*8, &
866  iostat = ierr )
867 
868  if ( ierr /= 0 ) then
869  log_error("FILE_TILEDATA_read_data_real8_real",*) 'data file not found!: ', trim(fname)
870  call prc_abort
871  endif
872 
873  read(fid,rec=1) buf(:,:)
874  close(fid)
875 
876  if ( yrevers_ ) then
877  !$omp parallel do
878 !OCL XFILL
879  do j = 1, jsize
880  do i = 1, isize
881  tile_data(i,j) = buf(i,jsize-j+1)
882  end do
883  end do
884  else
885  !$omp parallel do
886 !OCL XFILL
887  do j = 1, jsize
888  do i = 1, isize
889  tile_data(i,j) = buf(i,j)
890  end do
891  end do
892  end if
893 
894  return
896 
898  jsize, isize, &
899  fname, &
900  TILE_DATA, &
901  yrevers )
902  use scale_prc, only: &
903  prc_abort
904  integer, intent(in) :: jsize
905  integer, intent(in) :: isize
906  character(len=*), intent(in) :: fname
907  integer, intent(out) :: TILE_DATA(isize,jsize)
908 
909  logical, intent(in), optional :: yrevers
910 
911  integer(1) :: buf(isize,jsize)
912 
913  integer :: fid, ierr
914  logical :: yrevers_
915  integer :: i, j
916 
917  if ( present(yrevers) ) then
918  yrevers_ = yrevers
919  else
920  yrevers_ = .false.
921  end if
922 
923  fid = io_get_available_fid()
924  open( fid, &
925  file = fname, &
926  form = 'unformatted', &
927  access = 'direct', &
928  status = 'old', &
929  recl = isize*jsize*1, &
930  iostat = ierr )
931 
932  if ( ierr /= 0 ) then
933  log_error("FILE_TILEDATA_read_data_int1_int",*) 'data file not found!: ', trim(fname)
934  call prc_abort
935  endif
936 
937  read(fid,rec=1) buf(:,:)
938  close(fid)
939 
940  if ( yrevers_ ) then
941  !$omp parallel do
942 !OCL XFILL
943  do j = 1, jsize
944  do i = 1, isize
945  tile_data(i,j) = buf(i,jsize-j+1)
946  end do
947  end do
948  else
949  !$omp parallel do
950 !OCL XFILL
951  do j = 1, jsize
952  do i = 1, isize
953  tile_data(i,j) = buf(i,j)
954  end do
955  end do
956  end if
957 
958  return
959  end subroutine file_tiledata_read_data_int1_int
960 
962  jsize, isize, &
963  fname, &
964  TILE_DATA, &
965  yrevers )
966  use scale_prc, only: &
967  prc_abort
968  integer, intent(in) :: jsize
969  integer, intent(in) :: isize
970  character(len=*), intent(in) :: fname
971  integer, intent(out) :: TILE_DATA(isize,jsize)
972 
973  logical, intent(in), optional :: yrevers
974 
975  integer(2) :: buf(isize,jsize)
976 
977  integer :: fid, ierr
978  logical :: yrevers_
979  integer :: i, j
980 
981  if ( present(yrevers) ) then
982  yrevers_ = yrevers
983  else
984  yrevers_ = .false.
985  end if
986 
987  fid = io_get_available_fid()
988  open( fid, &
989  file = fname, &
990  form = 'unformatted', &
991  access = 'direct', &
992  status = 'old', &
993  recl = isize*jsize*2, &
994  iostat = ierr )
995 
996  if ( ierr /= 0 ) then
997  log_error("FILE_TILEDATA_read_data_int2_int",*) 'data file not found!: ', trim(fname)
998  call prc_abort
999  endif
1000 
1001  read(fid,rec=1) buf(:,:)
1002  close(fid)
1003 
1004  if ( yrevers_ ) then
1005  !$omp parallel do
1006 !OCL XFILL
1007  do j = 1, jsize
1008  do i = 1, isize
1009  tile_data(i,j) = buf(i,jsize-j+1)
1010  end do
1011  end do
1012  else
1013  !$omp parallel do
1014 !OCL XFILL
1015  do j = 1, jsize
1016  do i = 1, isize
1017  tile_data(i,j) = buf(i,j)
1018  end do
1019  end do
1020  end if
1021 
1022  return
1023  end subroutine file_tiledata_read_data_int2_int
1024 
1025  subroutine file_tiledata_read_data_int4_int( &
1026  jsize, isize, &
1027  fname, &
1028  TILE_DATA, &
1029  yrevers )
1030  use scale_prc, only: &
1031  prc_abort
1032  integer, intent(in) :: jsize
1033  integer, intent(in) :: isize
1034  character(len=*), intent(in) :: fname
1035  integer, intent(out) :: TILE_DATA(isize,jsize)
1036 
1037  logical, intent(in), optional :: yrevers
1038 
1039  integer(4) :: buf(isize,jsize)
1040 
1041  integer :: fid, ierr
1042  logical :: yrevers_
1043  integer :: i, j
1044 
1045  if ( present(yrevers) ) then
1046  yrevers_ = yrevers
1047  else
1048  yrevers_ = .false.
1049  end if
1050 
1051  fid = io_get_available_fid()
1052  open( fid, &
1053  file = fname, &
1054  form = 'unformatted', &
1055  access = 'direct', &
1056  status = 'old', &
1057  recl = isize*jsize*4, &
1058  iostat = ierr )
1059 
1060  if ( ierr /= 0 ) then
1061  log_error("FILE_TILEDATA_read_data_int4_int",*) 'data file not found!: ', trim(fname)
1062  call prc_abort
1063  endif
1064 
1065  read(fid,rec=1) buf(:,:)
1066  close(fid)
1067 
1068  if ( yrevers_ ) then
1069  !$omp parallel do
1070 !OCL XFILL
1071  do j = 1, jsize
1072  do i = 1, isize
1073  tile_data(i,j) = buf(i,jsize-j+1)
1074  end do
1075  end do
1076  else
1077  !$omp parallel do
1078 !OCL XFILL
1079  do j = 1, jsize
1080  do i = 1, isize
1081  tile_data(i,j) = buf(i,j)
1082  end do
1083  end do
1084  end if
1085 
1086  return
1087  end subroutine file_tiledata_read_data_int4_int
1088 
1089  subroutine file_tiledata_read_data_real4_int( &
1090  jsize, isize, &
1091  fname, &
1092  TILE_DATA, &
1093  yrevers )
1094  use scale_prc, only: &
1095  prc_abort
1096  integer, intent(in) :: jsize
1097  integer, intent(in) :: isize
1098  character(len=*), intent(in) :: fname
1099  integer, intent(out) :: TILE_DATA(isize,jsize)
1100 
1101  logical, intent(in), optional :: yrevers
1102 
1103  real(4) :: buf(isize,jsize)
1104 
1105  integer :: fid, ierr
1106  logical :: yrevers_
1107  integer :: i, j
1108 
1109  if ( present(yrevers) ) then
1110  yrevers_ = yrevers
1111  else
1112  yrevers_ = .false.
1113  end if
1114 
1115  fid = io_get_available_fid()
1116  open( fid, &
1117  file = fname, &
1118  form = 'unformatted', &
1119  access = 'direct', &
1120  status = 'old', &
1121  recl = isize*jsize*4, &
1122  iostat = ierr )
1123 
1124  if ( ierr /= 0 ) then
1125  log_error("FILE_TILEDATA_read_data_real4_int",*) 'data file not found!: ', trim(fname)
1126  call prc_abort
1127  endif
1128 
1129  read(fid,rec=1) buf(:,:)
1130  close(fid)
1131 
1132  if ( yrevers_ ) then
1133  !$omp parallel do
1134 !OCL XFILL
1135  do j = 1, jsize
1136  do i = 1, isize
1137  tile_data(i,j) = buf(i,jsize-j+1)
1138  end do
1139  end do
1140  else
1141  !$omp parallel do
1142 !OCL XFILL
1143  do j = 1, jsize
1144  do i = 1, isize
1145  tile_data(i,j) = buf(i,j)
1146  end do
1147  end do
1148  end if
1149 
1150  return
1151  end subroutine file_tiledata_read_data_real4_int
1152 
1153 end module scale_file_tiledata
subroutine file_tiledata_get_data_int1(nLATH, nLONH, 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, LATH, LONH, min_value, yrevers)
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, TILE_hit, TILE_JS, TILE_JE, TILE_IS, TILE_IE, jsh, jeh, ish, ieh, nLATH, nLONH)
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
subroutine file_tiledata_read_data_real4_real(jsize, isize, fname, TILE_DATA, yrevers)
real(rp), public const_undef
Definition: scale_const.F90:41
subroutine file_tiledata_get_data_real(nLATH, nLONH, 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, LATH, LONH, min_value, yrevers)
get tile data
subroutine file_tiledata_read_data_real4_int(jsize, isize, fname, TILE_DATA, yrevers)
subroutine file_tiledata_read_data_int4_real(jsize, isize, fname, TILE_DATA, yrevers)
module file_tiledata
integer function, public io_get_available_fid()
search & get available file ID
Definition: scale_io.F90:313
module PROCESS
Definition: scale_prc.F90:11
subroutine file_tiledata_read_data_int2_int(jsize, isize, fname, TILE_DATA, yrevers)
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:38
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module CONSTANT
Definition: scale_const.F90:11
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_read_data_int4_int(jsize, isize, fname, TILE_DATA, yrevers)
subroutine file_tiledata_read_data_int1_int(jsize, isize, fname, TILE_DATA, yrevers)
module PRECISION
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, single_fname, LATS, LATE, LONS, LONE)
get tile information
real(rp), public const_pi
pi
Definition: scale_const.F90:31
module STDIO
Definition: scale_io.F90:10
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
subroutine file_tiledata_read_data_int2_real(jsize, isize, fname, TILE_DATA, yrevers)
subroutine file_tiledata_read_data_real8_real(jsize, isize, fname, TILE_DATA, yrevers)