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  !
28  public file_tiledata_get_data
29 
30  interface file_tiledata_get_data
31  module procedure file_tiledata_get_data_real
32  module procedure file_tiledata_get_data_int1
33  end interface file_tiledata_get_data
34 
35  !-----------------------------------------------------------------------------
36  !
37  !++ Public parameters & variables
38  !
39  !-----------------------------------------------------------------------------
40  !
41  !++ Private procedure
42  !
43  !-----------------------------------------------------------------------------
44  !
45  !++ Private parameters & variables
46  !
47  !-----------------------------------------------------------------------------
48 contains
49  !-----------------------------------------------------------------------------
51  subroutine file_tiledata_get_info( &
52  TILE_nlim, &
53  TILE_DLAT, TILE_DLON, &
54  DOMAIN_LATS, DOMAIN_LATE, DOMAIN_LONS, DOMAIN_LONE, &
55  catalog_fname, &
56  GLOBAL_IA, &
57  TILE_nmax, &
58  TILE_fname, TILE_hit, &
59  TILE_JS, TILE_JE, TILE_IS, TILE_IE, &
60  nLATH, nLONH, jsh, jeh, ish, ieh, zonal, pole, &
61  single_fname, LATS, LATE, LONS, LONE )
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 
92  call file_tiledata_get_domain_info( tile_dlat, tile_dlon, & ! [IN]
93  domain_lats, domain_late, domain_lons, domain_lone, & ! [IN]
94  domain_js, domain_je, domain_is, domain_ie, & ! [OUT]
95  global_ia ) ! [OUT]
96 
97 
98  if ( catalog_fname /= "" ) then
99  log_newline
100  log_info("FILE_TILEDATA_get_info",*) 'Input catalogue file:', trim(catalog_fname)
101 
102  call file_tiledata_read_catalog_file( tile_nlim, & ! [IN]
103  catalog_fname, & ! [IN]
104  tile_dlat, tile_dlon, & ! [IN]
105  domain_is, global_ia, & ! [IN]
106  tile_nmax, & ! [OUT]
107  tile_fname(:), & ! [OUT]
108  tile_lats(:), tile_late(:), & ! [OUT]
109  tile_lons(:), tile_lone(:) ) ! [OUT]
110  lat_min = minval( tile_lats(1:tile_nmax) )
111  lat_max = maxval( tile_late(1:tile_nmax) )
112 
113  else
114  if ( .not. present(single_fname) ) then
115  log_error("FILE_TILEDATA_get_info",*) "single_fname is required if catalog_fname is empty"
116  call prc_abort
117  end if
118  if ( .not. present(lats) ) then
119  log_error("FILE_TILEDATA_get_info",*) "LATS is required if catalog_fname is empty"
120  call prc_abort
121  end if
122  if ( .not. present(late) ) then
123  log_error("FILE_TILEDATA_get_info",*) "LATE is required if catalog_fname is empty"
124  call prc_abort
125  end if
126  if ( .not. present(lons) ) then
127  log_error("FILE_TILEDATA_get_info",*) "LONS is required if catalog_fname is empty"
128  call prc_abort
129  end if
130  if ( .not. present(lone) ) then
131  log_error("FILE_TILEDATA_get_info",*) "LONE is required if catalog_fname is empty"
132  call prc_abort
133  end if
134 
135  tile_nmax = 1
136  tile_fname(1) = single_fname
137  tile_lats(1) = lats
138  tile_late(1) = late
139  tile_lons(1) = lons
140  tile_lone(1) = lone
141 
142  lat_min = lats
143  lat_max = late
144  end if
145 
146  zonal = ( domain_lone - domain_lons ) / ( 2.0_rp * pi ) > 0.9_rp
147 
148  pole = ( domain_lats < - pi * 0.5_rp + ( domain_late - domain_lats ) * 0.1_rp ) &
149  .or. ( domain_late > pi * 0.5_rp - ( domain_late - domain_lats ) * 0.1_rp )
150 
151  zonal = zonal .or. pole
152 
153  call file_tiledata_get_tile_info( tile_nmax, & ! [IN]
154  domain_js, domain_je, & ! [IN]
155  domain_is, domain_ie, & ! [IN]
156  global_ia, & ! [IN]
157  tile_dlat, tile_dlon, & ! [IN]
158  tile_lats(:), tile_late(:), & ! [IN]
159  tile_lons(:), tile_lone(:), & ! [IN]
160  zonal, & ! [IN]
161  tile_hit(:), & ! [OUT]
162  tile_js(:), tile_je(:), & ! [OUT]
163  tile_is(:), tile_ie(:), & ! [OUT]
164  jsh, jeh, ish, ieh, & ! [OUT]
165  nlath, nlonh ) ! [OUT]
166 
167  return
168  end subroutine file_tiledata_get_info
169 
170  !-----------------------------------------------------------------------------
172  subroutine file_tiledata_get_latlon( &
173  nLAT, nLON, &
174  jsh, ish, &
175  TILE_DLAT, TILE_DLON, &
176  LAT, LON )
177  implicit none
178  integer, intent(in) :: nlat, nlon
179  integer, intent(in) :: jsh, ish
180  real(rp), intent(in) :: tile_dlat, tile_dlon
181  real(rp), intent(out) :: lat(nlat)
182  real(rp), intent(out) :: lon(nlon)
183 
184  integer :: i, j
185 
186  do j = 1, nlat
187  lat(j) = tile_dlat * ( jsh + j - 0.5_rp )
188  end do
189  do i = 1, nlon
190  lon(i) = tile_dlon * ( ish + i - 0.5_rp )
191  end do
192 
193  return
194  end subroutine file_tiledata_get_latlon
195 
196  !-----------------------------------------------------------------------------
198  subroutine file_tiledata_get_data_real( &
199  nLAT, nLON, &
200  dirname, &
201  GLOBAL_IA, &
202  TILE_nmax, &
203  TILE_DLAT, TILE_DLON, &
204  TILE_fname, TILE_hit, &
205  TILE_JS, TILE_JE, &
206  TILE_IS, TILE_IE, &
207  jsh, jeh, ish, ieh, &
208  data_type, &
209  DATA, &
210  min_value, &
211  yrevers, &
212  step )
213  use scale_const, only: &
214  undef => const_undef, &
215  pi => const_pi, &
216  d2r => const_d2r
217  use scale_prc, only: &
218  prc_abort
219  integer, intent(in) :: nlat, nlon
220  character(len=*), intent(in) :: dirname
221  integer, intent(in) :: global_ia
222  integer, intent(in) :: tile_nmax
223  real(rp), intent(in) :: tile_dlat, tile_dlon
224  character(len=*), intent(in) :: tile_fname(:)
225  logical, intent(in) :: tile_hit(:)
226  integer, intent(in) :: tile_js(:)
227  integer, intent(in) :: tile_je(:)
228  integer, intent(in) :: tile_is(:)
229  integer, intent(in) :: tile_ie(:)
230  integer, intent(in) :: jsh
231  integer, intent(in) :: jeh
232  integer, intent(in) :: ish
233  integer, intent(in) :: ieh
234  character(len=*), intent(in) :: data_type
235  real(rp), intent(out) :: data(nlon,nlat)
236  real(rp), intent(in), optional :: min_value
237  logical, intent(in), optional :: yrevers
238  integer, intent(in), optional :: step
239 
240  abstract interface
241  subroutine rd( &
242  jsize, isize, &
243  fname, &
244  TILE_DATA, &
245  yrevers, &
246  step )
247  use scale_precision
248  integer, intent(in) :: jsize
249  integer, intent(in) :: isize
250  character(len=*), intent(in) :: fname
251  real(rp), intent(out) :: tile_data(isize,jsize)
252  logical, intent(in), optional :: yrevers
253  integer, intent(in), optional :: step
254  end subroutine rd
255  end interface
256 
257  procedure(rd), pointer :: read_data
258 
259  real(rp) :: min_value_
260 
261  character(len=H_LONG) :: fname
262  real(rp), allocatable :: tile_data(:,:)
263  integer :: jsize, isize
264  integer :: i, j, ii, jj, t
265 
266  if ( present(min_value) ) then
267  min_value_ = min_value
268  else
269  min_value_ = - abs(undef)
270  end if
271 
272  select case( data_type )
273  case ( "int2", "INT2" )
275  case ( "int4", "INT4" )
277  case ( "real4", "REAL4" )
279  case ( "real8", "REAL8" )
281  case default
282  log_error("FILE_TILEDATA_get_data_real",*) 'data_type is invalid: ', trim(data_type)
283  call prc_abort
284  end select
285 
286  !$omp parallel do
287 !OCL XFILL
288  do j = 1, nlat
289  do i = 1, nlon
290  DATA(i,j) = undef
291  end do
292  end do
293 
294  do t = 1, tile_nmax
295  if ( .not. tile_hit(t) ) cycle
296 
297  fname = trim(dirname) // '/' // trim(tile_fname(t))
298 
299  log_newline
300  log_info("FILE_TILEDATA_get_data_real",*) 'Input data file :', trim(fname)
301  log_info_cont(*) 'Tile (LAT) :', tile_js(t)*tile_dlat/d2r, (tile_je(t)+1)*tile_dlat/d2r
302  log_info_cont(*) ' (LON) :', tile_is(t)*tile_dlon/d2r, (tile_ie(t)+1)*tile_dlon/d2r
303 
304  isize = tile_ie(t) - tile_is(t) + 1
305  jsize = tile_je(t) - tile_js(t) + 1
306 
307  allocate( tile_data(isize,jsize) )
308 
309  call read_data( jsize, isize, & ! [IN]
310  fname, & ! [IN]
311  tile_data(:,:), & ! [OUT]
312  yrevers = yrevers, & ! [IN]
313  step = step ) ! [IN]
314 
315  !$omp parallel do &
316  !$omp private(i,j)
317  do jj = 1, jsize
318  j = tile_js(t) + jj - 1
319  if ( jsh <= j .and. j <= jeh ) then
320  do ii = 1, isize
321  i = tile_is(t) + ii - 1
322  if ( ish <= i .and. i <= ieh ) then
323  if ( tile_data(ii,jj) < min_value_ ) then
324  DATA(i-ish+1,j-jsh+1) = undef
325  else
326  DATA(i-ish+1,j-jsh+1) = tile_data(ii,jj)
327  end if
328  end if
329  i = i - global_ia
330  if ( ish <= i .and. i <= ieh ) then
331  if ( tile_data(ii,jj) < min_value_ ) then
332  DATA(i-ish+1,j-jsh+1) = undef
333  else
334  DATA(i-ish+1,j-jsh+1) = tile_data(ii,jj)
335  end if
336  end if
337  end do
338  end if
339  end do
340 
341  deallocate( tile_data )
342 
343  enddo ! tile loop
344 
345  return
346  end subroutine file_tiledata_get_data_real
347 
348  subroutine file_tiledata_get_data_int1( &
349  nLAT, nLON, &
350  dirname, &
351  GLOBAL_IA, &
352  TILE_nmax, &
353  TILE_DLAT, TILE_DLON, &
354  TILE_fname, TILE_hit, &
355  TILE_JS, TILE_JE, &
356  TILE_IS, TILE_IE, &
357  jsh, jeh, ish, ieh, &
358  data_type, &
359  DATA, &
360  min_value, &
361  yrevers, &
362  step )
363  use scale_const, only: &
364  undef => const_undef, &
365  undef2 => const_undef2, &
366  pi => const_pi, &
367  d2r => const_d2r
368  use scale_prc, only: &
369  prc_abort
370  integer, intent(in) :: nLAT, nLON
371  character(len=*), intent(in) :: dirname
372  integer, intent(in) :: GLOBAL_IA
373  integer, intent(in) :: TILE_nmax
374  real(RP), intent(in) :: TILE_DLAT, TILE_DLON
375  character(len=*), intent(in) :: TILE_fname(:)
376  logical, intent(in) :: TILE_hit(:)
377  integer, intent(in) :: TILE_JS(:)
378  integer, intent(in) :: TILE_JE(:)
379  integer, intent(in) :: TILE_IS(:)
380  integer, intent(in) :: TILE_IE(:)
381  integer, intent(in) :: jsh
382  integer, intent(in) :: jeh
383  integer, intent(in) :: ish
384  integer, intent(in) :: ieh
385  character(len=*), intent(in) :: data_type
386  integer, intent(out) :: DATA(nLON,nLAT)
387 
388  integer, intent(in), optional :: min_value
389  logical, intent(in), optional :: yrevers
390  integer, intent(in), optional :: step
391 
392  abstract interface
393  subroutine rd( &
394  jsize, isize, &
395  fname, &
396  TILE_DATA, &
397  yrevers, &
398  step )
399  use scale_precision
400  integer, intent(in) :: jsize
401  integer, intent(in) :: isize
402  character(len=*), intent(in) :: fname
403  integer, intent(out) :: TILE_DATA(isize,jsize)
404  logical, intent(in), optional :: yrevers
405  integer, intent(in), optional :: step
406  end subroutine rd
407  end interface
408 
409  integer :: min_value_
410 
411  procedure(rd), pointer :: read_data
412 
413  character(len=H_LONG) :: fname
414  integer, allocatable :: TILE_DATA(:,:)
415  integer :: jsize, isize
416  integer :: i, j, ii, jj, t
417 
418  if ( present(min_value) ) then
419  min_value_ = min_value
420  else
421  min_value_ = - abs(undef2)
422  end if
423 
424  select case( data_type )
425  case ( "int1", "INT1" )
427  case ( "int2", "INT2" )
429  case ( "int4", "INT4" )
431  case ( "real4", "REAL4" )
433  case default
434  log_error("FILE_TILEDATA_get_data_int1",*) 'data_type is invalid: ', trim(data_type)
435  call prc_abort
436  end select
437 
438  !$omp parallel do
439 !OCL XFILL
440  do j = 1, nlat
441  do i = 1, nlon
442  DATA(i,j) = undef2
443  end do
444  end do
445 
446  do t = 1, tile_nmax
447  if ( .not. tile_hit(t) ) cycle
448 
449  fname = trim(dirname) // '/' // trim(tile_fname(t))
450 
451  log_newline
452  log_info("FILE_TILEDATA_get_data_int1",*) 'Input data file :', trim(fname)
453  log_info_cont(*) 'Tile (LAT) :', tile_js(t)*tile_dlat/d2r, (tile_je(t)+1)*tile_dlat/d2r
454  log_info_cont(*) ' (LON) :', tile_is(t)*tile_dlon/d2r, (tile_ie(t)+1)*tile_dlon/d2r
455 
456  isize = tile_ie(t) - tile_is(t) + 1
457  jsize = tile_je(t) - tile_js(t) + 1
458 
459  allocate( tile_data(isize,jsize) )
460 
461  call read_data( jsize, isize, & ! [IN]
462  fname, & ! [IN]
463  tile_data(:,:), & ! [OUT]
464  yrevers = yrevers, & ! [IN]
465  step = step ) ! [IN]
466 
467  !$omp parallel do &
468  !$omp private(i,j)
469  do jj = 1, jsize
470  j = tile_js(t) + jj - 1
471  if ( jsh <= j .and. j <= jeh ) then
472  do ii = 1, isize
473  i = tile_is(t) + ii - 1
474  if ( ish <= i .and. i <= ieh ) then
475  if ( tile_data(ii,jj) < min_value_ ) then
476  DATA(i-ish+1,j-jsh+1) = undef2
477  else
478  DATA(i-ish+1,j-jsh+1) = tile_data(ii,jj)
479  end if
480  end if
481  i = i - global_ia
482  if ( ish <= i .and. i <= ieh ) then
483  if ( tile_data(ii,jj) < min_value_ ) then
484  DATA(i-ish+1,j-jsh+1) = undef2
485  else
486  DATA(i-ish+1,j-jsh+1) = tile_data(ii,jj)
487  end if
488  end if
489  end do
490  end if
491  end do
492 
493  deallocate( tile_data )
494 
495  enddo ! tile loop
496 
497  return
498  end subroutine file_tiledata_get_data_int1
499 
500 
501  ! private
502 
503 
504  !-----------------------------------------------------------------------------
506  subroutine file_tiledata_read_catalog_file( &
507  TILE_nlim, &
508  fname, &
509  TILE_DLAT, TILE_DLON, &
510  DOMAIN_IS, GLOBAL_IA, &
511  TILE_nmax, &
512  TILE_fname, &
513  TILE_LATS, TILE_LATE, &
514  TILE_LONS, TILE_LONE )
515  use scale_const, only: &
516  d2r => const_d2r
517  use scale_prc, only: &
518  prc_abort
519  integer, intent(in) :: TILE_nlim
520  character(len=*), intent(in) :: fname
521  real(RP), intent(in) :: TILE_DLAT
522  real(RP), intent(in) :: TILE_DLON
523  integer, intent(in) :: DOMAIN_IS
524  integer, intent(in) :: GLOBAL_IA
525  integer, intent(out) :: TILE_nmax
526  character(len=*), intent(out) :: TILE_fname(:)
527  real(RP), intent(out) :: TILE_LATS(:)
528  real(RP), intent(out) :: TILE_LATE(:)
529  real(RP), intent(out) :: TILE_LONS(:)
530  real(RP), intent(out) :: TILE_LONE(:)
531 
532  character(len=H_LONG) :: name
533  integer :: fid, ierr
534  integer :: index
535  integer :: t
536 
537  fid = io_get_available_fid()
538  call io_get_fname(name, fname)
539  open( fid, &
540  file = name, &
541  form = 'formatted', &
542  status = 'old', &
543  iostat = ierr )
544 
545  if ( ierr /= 0 ) then
546  log_error("FILE_TILEDATA_read_catalog_file",*) 'catalogue file not found! ', trim(name)
547  call prc_abort
548  endif
549 
550  ierr = 0
551  tile_nmax = - 1
552  do while ( ierr == 0 )
553  read(fid,*,iostat=ierr) index, tile_lats(1), tile_late(1), tile_lons(1), tile_lone(1), & ! WEST->EAST
554  tile_fname(1)
555  tile_nmax = tile_nmax + 1
556  end do
557 
558  if ( tile_nmax > tile_nlim ) then
559  log_error('FILE_TILEDATA_read_catalog_file',*) 'TILE_nmax must be >= ', tile_nmax
560  call prc_abort
561  end if
562 
563  rewind(fid)
564  do t = 1, tile_nmax
565  read(fid,*,iostat=ierr) index, tile_lats(t), tile_late(t), & ! South->North
566  tile_lons(t), tile_lone(t), & ! WEST->EAST
567  tile_fname(t)
568 
569  tile_lats(t) = tile_lats(t) * d2r
570  tile_late(t) = tile_late(t) * d2r
571 
572  if ( tile_lons(t) > tile_lone(t) ) tile_lone(t) = tile_lone(t) + 360.0_rp
573  tile_lons(t) = tile_lons(t) * d2r
574  tile_lone(t) = tile_lone(t) * d2r
575  end do
576 
577  close(fid)
578 
579  return
580  end subroutine file_tiledata_read_catalog_file
581 
582  subroutine file_tiledata_get_domain_info( &
583  TILE_DLAT, TILE_DLON, &
584  DOMAIN_LATS, DOMAIN_LATE, DOMAIN_LONS, DOMAIN_LONE, &
585  DOMAIN_JS, DOMAIN_JE, DOMAIN_IS, DOMAIN_IE, &
586  GLOBAL_IA )
587  use scale_const, only: &
588  pi => const_pi
589  real(RP), intent(in) :: TILE_DLAT
590  real(RP), intent(in) :: TILE_DLON
591  real(RP), intent(in) :: DOMAIN_LATS
592  real(RP), intent(in) :: DOMAIN_LATE
593  real(RP), intent(in) :: DOMAIN_LONS
594  real(RP), intent(in) :: DOMAIN_LONE
595  integer, intent(out) :: DOMAIN_JS
596  integer, intent(out) :: DOMAIN_JE
597  integer, intent(out) :: DOMAIN_IS
598  integer, intent(out) :: DOMAIN_IE
599  integer, intent(out) :: GLOBAL_IA
600 
601  domain_js = floor( domain_lats / tile_dlat )
602  domain_je = ceiling( domain_late / tile_dlat )
603  domain_is = floor( domain_lons / tile_dlon )
604  domain_ie = ceiling( domain_lone / tile_dlon )
605 
606  global_ia = nint( 2.0_rp * pi / tile_dlon )
607 
608  return
609  end subroutine file_tiledata_get_domain_info
610 
611  subroutine file_tiledata_get_tile_info( &
612  TILE_nmax, &
613  DOMAIN_JS, DOMAIN_JE, &
614  DOMAIN_IS, DOMAIN_IE, &
615  GLOBAL_IA, &
616  TILE_DLAT, TILE_DLON, &
617  TILE_LATS, TILE_LATE, &
618  TILE_LONS, TILE_LONE, &
619  zonal, &
620  TILE_hit, &
621  TILE_JS, TILE_JE, &
622  TILE_IS, TILE_IE, &
623  jsh, jeh, ish, ieh, &
624  nLAT, nLON )
625  use scale_const, only: &
626  pi => const_pi
627  integer, intent(in) :: TILE_nmax
628  integer, intent(in) :: DOMAIN_JS, DOMAIN_JE, DOMAIN_IS, DOMAIN_IE
629  integer, intent(in) :: GLOBAL_IA
630  real(RP), intent(in) :: TILE_DLAT, TILE_DLON
631  real(RP), intent(in) :: TILE_LATS(:), TILE_LATE(:), TILE_LONS(:), TILE_LONE(:)
632  logical, intent(in) :: zonal
633  logical, intent(out) :: TILE_hit(:)
634  integer, intent(out) :: TILE_JS(:), TILE_JE(:), TILE_IS(:), TILE_IE(:)
635  integer, intent(out) :: jsh, jeh, ish, ieh
636  integer, intent(out) :: nLAT, nLON
637 
638  logical :: hit_lat, hit_lon
639  integer :: nhalo
640  integer :: t
641 
642  nhalo = 2
643 
644  jsh = max( domain_js - nhalo, -floor( 0.5_rp * pi / tile_dlat ) )
645  jeh = min( domain_je + nhalo, floor( 0.5_rp * pi / tile_dlat ) )
646  ish = domain_is - nhalo
647  ieh = domain_ie + nhalo
648 
649  ! data file
650  !$omp parallel do &
651  !$omp private(hit_lat,hit_lon)
652  do t = 1, tile_nmax
653 
654  tile_js(t) = nint( tile_lats(t) / tile_dlat )
655  tile_je(t) = nint( tile_late(t) / tile_dlat ) - 1
656 
657  tile_is(t) = nint( tile_lons(t) / tile_dlon )
658  tile_ie(t) = nint( tile_lone(t) / tile_dlon ) - 1
659 
660  do while ( tile_ie(t) < ish )
661  tile_is(t) = tile_is(t) + global_ia
662  tile_ie(t) = tile_ie(t) + global_ia
663  end do
664  do while ( tile_is(t) - ish >= global_ia )
665  tile_is(t) = tile_is(t) - global_ia
666  tile_ie(t) = tile_ie(t) - global_ia
667  end do
668 
669  if ( ( jsh <= tile_js(t) .AND. tile_js(t) <= jeh ) &
670  .OR. ( jsh <= tile_je(t) .AND. tile_je(t) <= jeh ) &
671  .OR. ( tile_js(t) <= jsh .AND. jsh <= tile_je(t) ) &
672  .OR. ( tile_js(t) <= jeh .AND. jeh <= tile_je(t) ) ) then
673  hit_lat = .true.
674  else
675  hit_lat = .false.
676  endif
677 
678  if ( zonal ) then
679  hit_lon = .true.
680  else if ( ( tile_is(t) <= ieh ) &
681  .OR. ( tile_ie(t) - global_ia >= ish ) ) then
682  hit_lon = .true.
683  else
684  hit_lon = .false.
685  endif
686 
687  tile_hit(t) = ( hit_lat .AND. hit_lon )
688 
689  end do
690 
691  if ( zonal ) then
692  ish = minval(tile_is(1:tile_nmax))
693  ieh = maxval(tile_ie(1:tile_nmax))
694  jsh = max( jsh, minval(tile_js(1:tile_nmax)) )
695  jeh = min( jeh, maxval(tile_je(1:tile_nmax)) )
696  end if
697 
698  nlon = ieh - ish + 1
699  nlat = jeh - jsh + 1
700 
701  return
702  end subroutine file_tiledata_get_tile_info
703 
705  jsize, isize, &
706  fname, &
707  TILE_DATA, &
708  yrevers, &
709  step )
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  integer, intent(in), optional :: step
719 
720  integer(2) :: buf(isize,jsize)
721 
722  character(len=H_LONG) :: name
723  integer :: fid, ierr
724  logical :: yrevers_
725  logical :: lstep
726  integer :: step_
727  integer :: i, j
728 
729  if ( present(yrevers) ) then
730  yrevers_ = yrevers
731  else
732  yrevers_ = .false.
733  end if
734 
735  lstep = present(step)
736  if ( lstep ) then
737  step_ = step
738  inquire( file = fname, number = fid )
739  else
740  step_ = 1
741  fid = -1
742  end if
743 
744  if ( fid < 0 ) then
745  fid = io_get_available_fid()
746  call io_get_fname(name, fname)
747  open( fid, &
748  file = name, &
749  form = 'unformatted', &
750  access = 'direct', &
751  status = 'old', &
752  recl = isize*jsize*2, &
753  iostat = ierr )
754  if ( ierr /= 0 ) then
755  log_error("FILE_TILEDATA_read_data_int2_real",*) 'data file not found!: ', trim(name)
756  call prc_abort
757  endif
758  end if
759 
760 
761  read(fid,rec=step_) buf(:,:)
762  if ( .not. lstep ) close(fid)
763 
764  if ( yrevers_ ) then
765  !$omp parallel do
766 !OCL XFILL
767  do j = 1, jsize
768  do i = 1, isize
769  tile_data(i,j) = buf(i,jsize-j+1)
770  end do
771  end do
772  else
773  !$omp parallel do
774 !OCL XFILL
775  do j = 1, jsize
776  do i = 1, isize
777  tile_data(i,j) = buf(i,j)
778  end do
779  end do
780  end if
781 
782  return
783  end subroutine file_tiledata_read_data_int2_real
784 
786  jsize, isize, &
787  fname, &
788  TILE_DATA, &
789  yrevers, &
790  step )
791  use scale_prc, only: &
792  prc_abort
793  integer, intent(in) :: jsize
794  integer, intent(in) :: isize
795  character(len=*), intent(in) :: fname
796  real(RP), intent(out) :: TILE_DATA(isize,jsize)
797 
798  logical, intent(in), optional :: yrevers
799  integer, intent(in), optional :: step
800 
801  integer(4) :: buf(isize,jsize)
802 
803  character(len=H_LONG) :: name
804  integer :: fid, ierr
805  logical :: yrevers_
806  logical :: lstep
807  integer :: step_
808  integer :: i, j
809 
810  if ( present(yrevers) ) then
811  yrevers_ = yrevers
812  else
813  yrevers_ = .false.
814  end if
815 
816  lstep = present(step)
817  if ( lstep ) then
818  step_ = step
819  inquire( file = fname, number = fid )
820  else
821  step_ = 1
822  fid = -1
823  end if
824 
825  if ( fid < 0 ) then
826  fid = io_get_available_fid()
827  call io_get_fname(name, fname)
828  open( fid, &
829  file = name, &
830  form = 'unformatted', &
831  access = 'direct', &
832  status = 'old', &
833  recl = isize*jsize*4, &
834  iostat = ierr )
835  if ( ierr /= 0 ) then
836  log_error("FILE_TILEDATA_read_data_int4_real",*) 'data file not found!: ', trim(name)
837  call prc_abort
838  end if
839  endif
840 
841  read(fid,rec=step_) buf(:,:)
842  if ( .not. lstep ) close(fid)
843 
844  if ( yrevers_ ) then
845  !$omp parallel do
846 !OCL XFILL
847  do j = 1, jsize
848  do i = 1, isize
849  tile_data(i,j) = buf(i,jsize-j+1)
850  end do
851  end do
852  else
853  !$omp parallel do
854 !OCL XFILL
855  do j = 1, jsize
856  do i = 1, isize
857  tile_data(i,j) = buf(i,j)
858  end do
859  end do
860  end if
861 
862  return
863  end subroutine file_tiledata_read_data_int4_real
864 
866  jsize, isize, &
867  fname, &
868  TILE_DATA, &
869  yrevers, &
870  step )
871  use scale_prc, only: &
872  prc_abort
873  integer, intent(in) :: jsize
874  integer, intent(in) :: isize
875  character(len=*), intent(in) :: fname
876  real(RP), intent(out) :: TILE_DATA(isize,jsize)
877 
878  logical, intent(in), optional :: yrevers
879  integer, intent(in), optional :: step
880 
881  real(4) :: buf(isize,jsize)
882 
883  character(len=H_LONG) :: name
884  integer :: fid, ierr
885  logical :: yrevers_
886  logical :: lstep
887  integer :: step_
888  integer :: i, j
889 
890  if ( present(yrevers) ) then
891  yrevers_ = yrevers
892  else
893  yrevers_ = .false.
894  end if
895 
896  lstep = present(step)
897  if ( lstep ) then
898  step_ = step
899  inquire( file = fname, number = fid )
900  else
901  step_ = 1
902  fid = -1
903  end if
904 
905  if ( fid < 0 ) then
906  fid = io_get_available_fid()
907  call io_get_fname(name, fname)
908  open( fid, &
909  file = name, &
910  form = 'unformatted', &
911  access = 'direct', &
912  status = 'old', &
913  recl = isize*jsize*4, &
914  iostat = ierr )
915  if ( ierr /= 0 ) then
916  log_error("FILE_TILEDATA_read_data_real4_real",*) 'data file not found!: ', trim(name)
917  call prc_abort
918  endif
919  end if
920 
921  read(fid,rec=step_) buf(:,:)
922  if ( .not. lstep ) close(fid)
923 
924  if ( yrevers_ ) then
925  !$omp parallel do
926 !OCL XFILL
927  do j = 1, jsize
928  do i = 1, isize
929  tile_data(i,j) = buf(i,jsize-j+1)
930  end do
931  end do
932  else
933  !$omp parallel do
934 !OCL XFILL
935  do j = 1, jsize
936  do i = 1, isize
937  tile_data(i,j) = buf(i,j)
938  end do
939  end do
940  end if
941 
942  return
944 
946  jsize, isize, &
947  fname, &
948  TILE_DATA, &
949  yrevers, &
950  step )
951  use scale_prc, only: &
952  prc_abort
953  integer, intent(in) :: jsize
954  integer, intent(in) :: isize
955  character(len=*), intent(in) :: fname
956  real(RP), intent(out) :: TILE_DATA(isize,jsize)
957 
958  logical, intent(in), optional :: yrevers
959  integer, intent(in), optional :: step
960 
961  real(8) :: buf(isize,jsize)
962 
963  character(len=H_LONG) :: name
964  integer :: fid, ierr
965  logical :: yrevers_
966  logical :: lstep
967  integer :: step_
968  integer :: i, j
969 
970  if ( present(yrevers) ) then
971  yrevers_ = yrevers
972  else
973  yrevers_ = .false.
974  end if
975 
976  lstep = present(step)
977  if ( lstep ) then
978  step_ = step
979  inquire( file = fname, number = fid )
980  else
981  step_ = 1
982  fid = -1
983  end if
984 
985  if ( fid < 0 ) then
986  fid = io_get_available_fid()
987  call io_get_fname(name, fname)
988  open( fid, &
989  file = name, &
990  form = 'unformatted', &
991  access = 'direct', &
992  status = 'old', &
993  recl = isize*jsize*8, &
994  iostat = ierr )
995  if ( ierr /= 0 ) then
996  log_error("FILE_TILEDATA_read_data_real8_real",*) 'data file not found!: ', trim(name)
997  call prc_abort
998  endif
999  end if
1000 
1001  read(fid,rec=step_) buf(:,:)
1002  if ( .not. lstep ) 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_real8_real
1024 
1025  subroutine file_tiledata_read_data_int1_int( &
1026  jsize, isize, &
1027  fname, &
1028  TILE_DATA, &
1029  yrevers, &
1030  step )
1031  use scale_prc, only: &
1032  prc_abort
1033  integer, intent(in) :: jsize
1034  integer, intent(in) :: isize
1035  character(len=*), intent(in) :: fname
1036  integer, intent(out) :: TILE_DATA(isize,jsize)
1037 
1038  logical, intent(in), optional :: yrevers
1039  integer, intent(in), optional :: step
1040 
1041  integer(1) :: buf(isize,jsize)
1042 
1043  character(len=H_LONG) :: name
1044  integer :: fid, ierr
1045  logical :: yrevers_
1046  logical :: lstep
1047  integer :: step_
1048  integer :: i, j
1049 
1050  if ( present(yrevers) ) then
1051  yrevers_ = yrevers
1052  else
1053  yrevers_ = .false.
1054  end if
1055 
1056  lstep = present(step)
1057  if ( lstep ) then
1058  step_ = step
1059  inquire( file = fname, number = fid )
1060  else
1061  step_ = 1
1062  fid = -1
1063  end if
1064 
1065  if ( fid < 0 ) then
1066  fid = io_get_available_fid()
1067  call io_get_fname(name, fname)
1068  open( fid, &
1069  file = name, &
1070  form = 'unformatted', &
1071  access = 'direct', &
1072  status = 'old', &
1073  recl = isize*jsize*1, &
1074  iostat = ierr )
1075  if ( ierr /= 0 ) then
1076  log_error("FILE_TILEDATA_read_data_int1_int",*) 'data file not found!: ', trim(name)
1077  call prc_abort
1078  endif
1079  end if
1080 
1081  read(fid,rec=step_) buf(:,:)
1082  if ( .not. lstep ) close(fid)
1083 
1084  if ( yrevers_ ) then
1085  !$omp parallel do
1086 !OCL XFILL
1087  do j = 1, jsize
1088  do i = 1, isize
1089  tile_data(i,j) = buf(i,jsize-j+1)
1090  end do
1091  end do
1092  else
1093  !$omp parallel do
1094 !OCL XFILL
1095  do j = 1, jsize
1096  do i = 1, isize
1097  tile_data(i,j) = buf(i,j)
1098  end do
1099  end do
1100  end if
1101 
1102  return
1103  end subroutine file_tiledata_read_data_int1_int
1104 
1105  subroutine file_tiledata_read_data_int2_int( &
1106  jsize, isize, &
1107  fname, &
1108  TILE_DATA, &
1109  yrevers, &
1110  step )
1111  use scale_prc, only: &
1112  prc_abort
1113  integer, intent(in) :: jsize
1114  integer, intent(in) :: isize
1115  character(len=*), intent(in) :: fname
1116  integer, intent(out) :: TILE_DATA(isize,jsize)
1117 
1118  logical, intent(in), optional :: yrevers
1119  integer, intent(in), optional :: step
1120 
1121  integer(2) :: buf(isize,jsize)
1122 
1123  character(len=H_LONG) :: name
1124  integer :: fid, ierr
1125  logical :: yrevers_
1126  logical :: lstep
1127  integer :: step_
1128  integer :: i, j
1129 
1130  if ( present(yrevers) ) then
1131  yrevers_ = yrevers
1132  else
1133  yrevers_ = .false.
1134  end if
1135 
1136  lstep = present(step)
1137  if ( lstep ) then
1138  step_ = step
1139  inquire( file = fname, number = fid )
1140  else
1141  step_ = 1
1142  fid = -1
1143  end if
1144 
1145  if ( fid < 0 ) then
1146  fid = io_get_available_fid()
1147  call io_get_fname(name, fname)
1148  open( fid, &
1149  file = name, &
1150  form = 'unformatted', &
1151  access = 'direct', &
1152  status = 'old', &
1153  recl = isize*jsize*2, &
1154  iostat = ierr )
1155  if ( ierr /= 0 ) then
1156  log_error("FILE_TILEDATA_read_data_int2_int",*) 'data file not found!: ', trim(name)
1157  call prc_abort
1158  endif
1159  end if
1160 
1161  read(fid,rec=step_) buf(:,:)
1162  if ( .not. lstep ) close(fid)
1163 
1164  if ( yrevers_ ) then
1165  !$omp parallel do
1166 !OCL XFILL
1167  do j = 1, jsize
1168  do i = 1, isize
1169  tile_data(i,j) = buf(i,jsize-j+1)
1170  end do
1171  end do
1172  else
1173  !$omp parallel do
1174 !OCL XFILL
1175  do j = 1, jsize
1176  do i = 1, isize
1177  tile_data(i,j) = buf(i,j)
1178  end do
1179  end do
1180  end if
1181 
1182  return
1183  end subroutine file_tiledata_read_data_int2_int
1184 
1185  subroutine file_tiledata_read_data_int4_int( &
1186  jsize, isize, &
1187  fname, &
1188  TILE_DATA, &
1189  yrevers, &
1190  step )
1191  use scale_prc, only: &
1192  prc_abort
1193  integer, intent(in) :: jsize
1194  integer, intent(in) :: isize
1195  character(len=*), intent(in) :: fname
1196  integer, intent(out) :: TILE_DATA(isize,jsize)
1197 
1198  logical, intent(in), optional :: yrevers
1199  integer, intent(in), optional :: step
1200 
1201  integer(4) :: buf(isize,jsize)
1202 
1203  character(len=H_LONG) :: name
1204  integer :: fid, ierr
1205  logical :: yrevers_
1206  logical :: lstep
1207  integer :: step_
1208  integer :: i, j
1209 
1210  if ( present(yrevers) ) then
1211  yrevers_ = yrevers
1212  else
1213  yrevers_ = .false.
1214  end if
1215 
1216  lstep = present(step)
1217  if ( lstep ) then
1218  step_ = step
1219  inquire( file = fname, number = fid )
1220  else
1221  step_ = 1
1222  fid = -1
1223  end if
1224 
1225  if ( fid < 0 ) then
1226  fid = io_get_available_fid()
1227  call io_get_fname(name, fname)
1228  open( fid, &
1229  file = name, &
1230  form = 'unformatted', &
1231  access = 'direct', &
1232  status = 'old', &
1233  recl = isize*jsize*4, &
1234  iostat = ierr )
1235  if ( ierr /= 0 ) then
1236  log_error("FILE_TILEDATA_read_data_int4_int",*) 'data file not found!: ', trim(name)
1237  call prc_abort
1238  endif
1239  end if
1240 
1241  read(fid,rec=step_) buf(:,:)
1242  if ( .not. lstep ) close(fid)
1243 
1244  if ( yrevers_ ) then
1245  !$omp parallel do
1246 !OCL XFILL
1247  do j = 1, jsize
1248  do i = 1, isize
1249  tile_data(i,j) = buf(i,jsize-j+1)
1250  end do
1251  end do
1252  else
1253  !$omp parallel do
1254 !OCL XFILL
1255  do j = 1, jsize
1256  do i = 1, isize
1257  tile_data(i,j) = buf(i,j)
1258  end do
1259  end do
1260  end if
1261 
1262  return
1263  end subroutine file_tiledata_read_data_int4_int
1264 
1265  subroutine file_tiledata_read_data_real4_int( &
1266  jsize, isize, &
1267  fname, &
1268  TILE_DATA, &
1269  yrevers, &
1270  step )
1271  use scale_prc, only: &
1272  prc_abort
1273  integer, intent(in) :: jsize
1274  integer, intent(in) :: isize
1275  character(len=*), intent(in) :: fname
1276  integer, intent(out) :: TILE_DATA(isize,jsize)
1277 
1278  logical, intent(in), optional :: yrevers
1279  integer, intent(in), optional :: step
1280 
1281  real(4) :: buf(isize,jsize)
1282 
1283  character(len=H_LONG) :: name
1284  integer :: fid, ierr
1285  logical :: yrevers_
1286  logical :: lstep
1287  integer :: step_
1288  integer :: i, j
1289 
1290  if ( present(yrevers) ) then
1291  yrevers_ = yrevers
1292  else
1293  yrevers_ = .false.
1294  end if
1295 
1296  lstep = present(step)
1297  if ( lstep ) then
1298  step_ = step
1299  inquire( file = fname, number = fid )
1300  else
1301  step_ = 1
1302  fid = -1
1303  end if
1304 
1305  if ( fid < 0 ) then
1306  fid = io_get_available_fid()
1307  call io_get_fname(name, fname)
1308  open( fid, &
1309  file = name, &
1310  form = 'unformatted', &
1311  access = 'direct', &
1312  status = 'old', &
1313  recl = isize*jsize*4, &
1314  iostat = ierr )
1315  if ( ierr /= 0 ) then
1316  log_error("FILE_TILEDATA_read_data_real4_int",*) 'data file not found!: ', trim(name)
1317  call prc_abort
1318  endif
1319  end if
1320 
1321  read(fid,rec=step_) buf(:,:)
1322  if ( .not. lstep ) close(fid)
1323 
1324  if ( yrevers_ ) then
1325  !$omp parallel do
1326 !OCL XFILL
1327  do j = 1, jsize
1328  do i = 1, isize
1329  tile_data(i,j) = buf(i,jsize-j+1)
1330  end do
1331  end do
1332  else
1333  !$omp parallel do
1334 !OCL XFILL
1335  do j = 1, jsize
1336  do i = 1, isize
1337  tile_data(i,j) = buf(i,j)
1338  end do
1339  end do
1340  end if
1341 
1342  return
1343  end subroutine file_tiledata_read_data_real4_int
1344 
1345 end module scale_file_tiledata
scale_file_tiledata::file_tiledata_read_data_int2_real
subroutine file_tiledata_read_data_int2_real(jsize, isize, fname, TILE_DATA, yrevers, step)
Definition: scale_file_tiledata.F90:710
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_const::const_undef2
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:40
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_file_tiledata::file_tiledata_get_data_int1
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)
Definition: scale_file_tiledata.F90:363
scale_file_tiledata::file_tiledata_get_domain_info
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)
Definition: scale_file_tiledata.F90:587
scale_file_tiledata::file_tiledata_get_info
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
Definition: scale_file_tiledata.F90:62
scale_io::io_get_available_fid
integer function, public io_get_available_fid()
search & get available file ID
Definition: scale_io.F90:373
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_io::io_get_fname
subroutine, public io_get_fname(outstr, instr, rank, ext, len)
generate process specific filename
Definition: scale_io.F90:421
scale_io
module STDIO
Definition: scale_io.F90:10
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_file_tiledata::file_tiledata_read_data_int1_int
subroutine file_tiledata_read_data_int1_int(jsize, isize, fname, TILE_DATA, yrevers, step)
Definition: scale_file_tiledata.F90:1031
scale_file_tiledata::file_tiledata_read_data_int4_real
subroutine file_tiledata_read_data_int4_real(jsize, isize, fname, TILE_DATA, yrevers, step)
Definition: scale_file_tiledata.F90:791
scale_file_tiledata::file_tiledata_read_data_int2_int
subroutine file_tiledata_read_data_int2_int(jsize, isize, fname, TILE_DATA, yrevers, step)
Definition: scale_file_tiledata.F90:1111
scale_file_tiledata::file_tiledata_get_tile_info
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)
Definition: scale_file_tiledata.F90:625
scale_const::const_pi
real(rp), parameter, public const_pi
pi
Definition: scale_const.F90:32
scale_file_tiledata::file_tiledata_get_latlon
subroutine, public file_tiledata_get_latlon(nLAT, nLON, jsh, ish, TILE_DLAT, TILE_DLON, LAT, LON)
get tile data
Definition: scale_file_tiledata.F90:177
scale_file_tiledata::file_tiledata_read_catalog_file
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
Definition: scale_file_tiledata.F90:515
scale_file_tiledata::file_tiledata_read_data_real4_real
subroutine file_tiledata_read_data_real4_real(jsize, isize, fname, TILE_DATA, yrevers, step)
Definition: scale_file_tiledata.F90:871
scale_const::const_d2r
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:33
scale_file_tiledata
module file_tiledata
Definition: scale_file_tiledata.F90:12
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_file_tiledata::file_tiledata_read_data_int4_int
subroutine file_tiledata_read_data_int4_int(jsize, isize, fname, TILE_DATA, yrevers, step)
Definition: scale_file_tiledata.F90:1191
scale_file_tiledata::file_tiledata_read_data_real4_int
subroutine file_tiledata_read_data_real4_int(jsize, isize, fname, TILE_DATA, yrevers, step)
Definition: scale_file_tiledata.F90:1271
scale_file_tiledata::file_tiledata_read_data_real8_real
subroutine file_tiledata_read_data_real8_real(jsize, isize, fname, TILE_DATA, yrevers, step)
Definition: scale_file_tiledata.F90:951