SCALE-RM
scale_file_grads.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 #include "scalelib.h"
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_io
18  use scale_prof
19  !-----------------------------------------------------------------------------
20  implicit none
21  private
22  !-----------------------------------------------------------------------------
23  !
24  !++ Public procedure
25  !
26  public :: file_grads_finalize
27  public :: file_grads_open
28  public :: file_grads_varid
29  public :: file_grads_varcheck
30  public :: file_grads_isoned
31  public :: file_grads_get_shape
32  public :: file_grads_read
33  public :: file_grads_close
34 
35  interface file_grads_get_shape
36  module procedure file_grads_get_shape_name
37  module procedure file_grads_get_shape_id
38  end interface file_grads_get_shape
39 
40  interface file_grads_read
41  module procedure file_grads_read_1d_name
42  module procedure file_grads_read_1d_id
43  module procedure file_grads_read_2d_name
44  module procedure file_grads_read_2d_id
45  module procedure file_grads_read_3d_name
46  module procedure file_grads_read_3d_id
47  end interface file_grads_read
48  !-----------------------------------------------------------------------------
49  !
50  !++ Public parameters & variables
51  !
52  !-----------------------------------------------------------------------------
53  !
54  !++ Private procedure
55  !
56  !-----------------------------------------------------------------------------
57  !
58  !++ Private parameters & variables
59  !
60  integer, parameter :: nmls_max = 10
61  integer, parameter :: vars_max = 100
62  integer, parameter :: lvars_max = 1000
63  type t_var
64  character(len=H_SHORT) :: name
65  character(len=H_LONG) :: fname
66  character(len=H_SHORT) :: dtype
67  real(RP) :: swpoint
68  real(RP) :: dd
69  integer :: lnum
70  real(RP), allocatable :: lvars(:)
71  integer :: startrec
72  integer :: totalrec
73  real(SP) :: missval
74  integer :: nx
75  integer :: ny
76  integer :: nz
77  logical :: yrev
78  integer :: endian ! 0: little, 1: big
79  character(len=H_SHORT) :: bintype
80  end type t_var
81  type t_nml
82  character(len=H_LONG) :: fname
83  integer :: nx, ny, nz
84  type(t_var), allocatable :: vars(:)
85  integer :: nvars
86  end type t_nml
87  type t_file
88  character(len=H_LONG) :: fname
89  character(len=H_SHORT) :: postfix
90  integer :: fid
91  end type t_file
92  type(t_nml), save :: nmls(nmls_max)
93  integer :: nnmls = 0
94  type(t_file) :: files(vars_max)
95  integer :: nfiles = 0
96 
97  !-----------------------------------------------------------------------------
98 contains
99  !-----------------------------------------------------------------------------
101  subroutine file_grads_open( &
102  file_name, &
103  file_id )
104  use scale_prc, only: &
105  prc_abort
106  use scale_const, only: &
107  undef => const_undef
108  implicit none
109  character(len=*), intent(in) :: file_name
110  integer, intent(out) :: file_id
111 
112  character(len=H_SHORT) :: name ! up to 16 characters
113  character(len=H_SHORT) :: dtype ! 'linear','levels','map'
114  character(len=H_LONG) :: fname ! head of file name
115  real(rp) :: swpoint ! start point (south-west point) for "linear"
116  real(rp) :: dd ! dlon,dlat for "linear"
117  integer :: lnum ! number of data
118  real(rp) :: lvars(lvars_max) ! values for "levels"
119  integer :: startrec ! record position
120  integer :: totalrec ! total record number per one time
121  real(sp) :: missval ! missing value
122  integer :: nx ! optional
123  integer :: ny ! optional
124  integer :: nz ! optional
125  character(len=H_SHORT) :: fendian ! option for "map"
126  logical :: yrev ! option for "map", if yrev=.true., order of data is NW to SE.
127  character(len=H_SHORT) :: bintype ! binary type: 'int?' or 'real?' ?=2,4, or 8
128 
129  namelist /grads_dims/ &
130  nx, &
131  ny, &
132  nz
133 
134  namelist /grads_item/ &
135  name, & ! necessary
136  dtype, & ! necessary
137  fname, & ! necessary except for linear data
138  swpoint, & ! for linear data
139  dd, & ! for linear data
140  lnum, & ! for levels data
141  lvars, & ! for levels data
142  startrec, & ! for map data
143  totalrec, & ! for map data
144  missval, & ! option
145  nx, & ! option
146  ny, & ! option
147  nz, & ! option
148  yrev, & ! option
149  bintype ! option
150 ! fendian ! option
151 
152  character(len=H_LONG) :: dirname
153 
154  integer :: fid
155  integer :: nvars
156  integer :: ierr
157  integer :: n
158  !---------------------------------------------------------------------------
159 
160  call io_get_fname(fname, file_name)
161 
162  log_newline
163  log_info("FILE_GrADS_open",*) 'open namelist file :', trim(fname)
164 
165  ! check exist
166  do n = 1, nnmls
167  if ( nmls(n)%fname == fname ) then
168  ! alread read
169  file_id = n
170  return
171  end if
172  end do
173 
174 
175  fid = io_get_available_fid()
176  !--- open namelist
177  open( fid, &
178  file = fname, &
179  form = 'formatted', &
180  status = 'old', &
181  action = 'read', &
182  iostat = ierr )
183  if ( ierr /= 0 ) then
184  log_error("FILE_GrADS_open",*) 'Input file is not found! ', trim(fname)
185  call prc_abort
186  endif
187 
188  call check_oldnamelist( fid )
189 
190  !--- read namelist dims
191  read(fid,nml=grads_dims,iostat=ierr)
192  if( ierr /= 0 ) then !--- missing or fatal error
193  log_error("FILE_GrADS_open",*) 'Not appropriate names in GrADS_DIMS in ', trim(fname),'. Check!'
194  call prc_abort
195  endif
196  log_nml(grads_dims)
197 
198 
199  nnmls = nnmls + 1
200  if ( nnmls > nmls_max ) then
201  log_error("FILE_GrADS_open",*) 'Number of GrADS file to be open is exceeded the maximum', nmls_max
202  call prc_abort
203  end if
204  file_id = nnmls
205 
206  nmls(file_id)%fname = file_name
207  nmls(file_id)%nx = nx
208  nmls(file_id)%ny = ny
209  nmls(file_id)%nz = nz
210 
211 
212 
213  !--- count the number of variables
214  rewind(fid)
215  nvars = 0
216  do n = 1, vars_max
217  read(fid, nml=grads_item, iostat=ierr)
218  if( ierr > 0 )then
219  log_error("FILE_GrADS_open",*) 'Not appropriate names in GrADS_ITEM in ', trim(file_name),'. Check!'
220  call prc_abort
221  else if( ierr < 0 )then
222  exit
223  endif
224  nvars = nvars + 1
225  enddo
226 
227  if ( nvars > vars_max ) then
228  log_error("FILE_GRADS_open",*) 'The number of grads vars exceeds the limit! ', &
229  nvars, ' > ', vars_max
230  call prc_abort
231  endif
232 
233  nmls(file_id)%nvars = nvars
234  allocate( nmls(file_id)%vars(nvars) )
235 
236  n = index( file_name, '/', back=.true. )
237  if ( n > 0 ) then
238  dirname = file_name(1:n)
239  else
240  dirname = ""
241  end if
242 
243  !--- read information of the variables
244  rewind(fid)
245  do n = 1, nvars
246 
247  ! set default
248  name = ''
249  dtype = ''
250  fname = ''
251  swpoint = undef
252  dd = undef
253  lnum = -1
254  lvars(:) = undef
255  startrec = -1
256  totalrec = -1
257  nx = nmls(file_id)%nx
258  ny = nmls(file_id)%ny
259  nz = nmls(file_id)%nz
260  yrev = .false.
261  fendian = 'big'
262  missval = undef
263  bintype = 'real4'
264 
265  ! read namelist
266  read(fid, nml=grads_item, iostat=ierr)
267  if( ierr /= 0 ) exit
268 
269  nmls(file_id)%vars(n)%name = upcase(name)
270  if ( fname(1:1) == "/" ) then
271  nmls(file_id)%vars(n)%fname = fname
272  else
273  nmls(file_id)%vars(n)%fname = trim(dirname) // fname
274  end if
275  nmls(file_id)%vars(n)%dtype = dtype
276  nmls(file_id)%vars(n)%swpoint = swpoint
277  nmls(file_id)%vars(n)%dd = dd
278  nmls(file_id)%vars(n)%lnum = lnum
279  if ( lnum > 0 ) then
280  if ( lnum > lvars_max ) then
281  log_error("FILE_GrADS_open",*) 'lnum exceeds the limit', lvars_max
282  call prc_abort
283  end if
284  allocate( nmls(file_id)%vars(n)%lvars(lnum) )
285  nmls(file_id)%vars(n)%lvars(:) = lvars(1:lnum)
286  end if
287  nmls(file_id)%vars(n)%startrec = startrec
288  nmls(file_id)%vars(n)%totalrec = totalrec
289  nmls(file_id)%vars(n)%missval = missval
290  nmls(file_id)%vars(n)%nx = nx
291  nmls(file_id)%vars(n)%ny = ny
292  nmls(file_id)%vars(n)%nz = nz
293  nmls(file_id)%vars(n)%yrev = yrev
294  if ( fendian == "big" ) then
295  nmls(file_id)%vars(n)%endian = 1
296  else
297  nmls(file_id)%vars(n)%endian = 0
298  end if
299  nmls(file_id)%vars(n)%bintype = bintype
300 
301  end do
302 
303  close( fid )
304 
305  end subroutine file_grads_open
306 
307  !-----------------------------------------------------------------------------
308  subroutine file_grads_varid( &
309  file_id, &
310  var_name, &
311  var_id )
312  use scale_prc, only: &
313  prc_abort
314  implicit none
315  integer, intent(in) :: file_id
316  character(len=*), intent(in) :: var_name
317  integer, intent(out) :: var_id
318 
319  character(len=len(var_name)) :: vname
320  integer :: n
321 
322  if ( file_id < 0 ) then
323  log_error("FILE_GrADS_varid",*) 'file_id is invalid: ', file_id
324  call prc_abort
325  end if
326 
327  vname = upcase(var_name)
328 
329  var_id = -1
330  do n = 1, nmls(file_id)%nvars
331  if ( nmls(file_id)%vars(n)%name == vname ) then
332  var_id = n
333  return
334  end if
335  end do
336 
337  return
338  end subroutine file_grads_varid
339 
340  subroutine file_grads_varcheck( &
341  file_id, &
342  var_name, &
343  exist )
344  use scale_prc, only: &
345  prc_abort
346  implicit none
347  integer, intent(in) :: file_id
348  character(len=*), intent(in) :: var_name
349  logical, intent(out) :: exist
350 
351  integer :: var_id
352 
353  if ( file_id < 0 ) then
354  log_error("FILE_GrADS_varcheck",*) 'file_id is invalid: ', file_id
355  call prc_abort
356  end if
357 
358  exist = .true.
359 
360  call file_grads_varid( file_id, var_name, & ! (in)
361  var_id ) ! (out)
362  if ( var_id < 0 ) then
363  exist = .false.
364  end if
365 
366  return
367  end subroutine file_grads_varcheck
368 
369  function file_grads_isoned( &
370  file_id, &
371  var_id ) &
372  result(ret)
373  use scale_prc, only: &
374  prc_abort
375  implicit none
376  integer, intent(in) :: file_id
377  integer, intent(in) :: var_id
378  logical :: ret
379 
380  if ( file_id < 0 ) then
381  log_error("FILE_GrADS_isOneD",*) 'file_id is invalid: ', file_id
382  call prc_abort
383  end if
384  if ( var_id < 0 ) then
385  log_error("FILE_GrADS_isOneD",*) 'var_id is invalid: ', var_id
386  call prc_abort
387  end if
388 
389  select case( nmls(file_id)%vars(var_id)%dtype )
390  case ('linear', 'levels')
391  ret = .true.
392  case default
393  ret = .false.
394  end select
395 
396  return
397  end function file_grads_isoned
398 
399  !-----------------------------------------------------------------------------
400  subroutine file_grads_get_shape_name( &
401  file_id, &
402  var_name, &
403  shape )
404  use scale_prc, only: &
405  prc_abort
406  implicit none
407  integer, intent(in) :: file_id
408  character(len=*), intent(in) :: var_name
409  integer, intent(out) :: shape(:)
410 
411  integer :: var_id
412 
413  call file_grads_varid( file_id, var_name, & ! (in)
414  var_id ) ! (out)
415 
416  if ( var_id < 0 ) then
417  log_error("FILE_GrADS_get_shape",*) 'variable "', trim(var_name), ' is not founed in file "', trim(nmls(file_id)%fname), '"'
418  call prc_abort
419  end if
420 
421  call file_grads_get_shape_id( file_id, var_id, & ! (in)
422  shape(:) ) ! (out)
423 
424  return
425  end subroutine file_grads_get_shape_name
426  !-----------------------------------------------------------------------------
427  subroutine file_grads_get_shape_id( &
428  file_id, &
429  var_id, &
430  shape )
431  implicit none
432  integer, intent(in) :: file_id
433  integer, intent(in) :: var_id
434  integer, intent(out) :: shape(:)
435 
436  intrinsic :: size
437 
438  if ( file_grads_isoned( file_id, var_id ) ) then
439  if ( nmls(file_id)%vars(var_id)%dtype == "levels" ) then
440  shape(1) = nmls(file_id)%vars(var_id)%lnum
441  else
442  shape(1) = -1
443  end if
444  else if ( size(shape) == 2 ) then
445  shape(1) = nmls(file_id)%vars(var_id)%nx
446  shape(2) = nmls(file_id)%vars(var_id)%ny
447  else
448  shape(1) = nmls(file_id)%vars(var_id)%nz
449  shape(2) = nmls(file_id)%vars(var_id)%nx
450  shape(3) = nmls(file_id)%vars(var_id)%ny
451  end if
452 
453  return
454  end subroutine file_grads_get_shape_id
455 
456  !-----------------------------------------------------------------------------
457  subroutine file_grads_read_1d_name( &
458  file_id, &
459  var_name, &
460  var, &
461  step, &
462  start, &
463  count, &
464  postfix )
465  use scale_prc, only: &
466  prc_abort
467  implicit none
468  integer, intent(in) :: file_id
469  character(len=*), intent(in) :: var_name
470  real(RP), intent(out) :: var(:)
471  integer, intent(in), optional :: step
472  integer, intent(in), optional :: start(1)
473  integer, intent(in), optional :: count(1)
474  character(len=*), intent(in), optional :: postfix
475 
476  integer :: var_id
477  !---------------------------------------------------------------------------
478 
479  call file_grads_varid( file_id, var_name, & ! (in)
480  var_id ) ! (out)
481 
482  if ( var_id < 0 ) then
483  log_error("FILE_GrADS_read_1D_name",*) 'variable "', trim(var_name), ' is not founed in file "', trim(nmls(file_id)%fname), '"'
484  call prc_abort
485  end if
486 
487 
488  call file_grads_read_1d_id( file_id, var_id, & ! (in)
489  var(:), & ! (out)
490  step = step, & ! (in)
491  start = start, & ! (in)
492  count = count, & ! (in)
493  postfix = postfix ) ! (in)
494 
495  return
496  end subroutine file_grads_read_1d_name
497  !-----------------------------------------------------------------------------
498  subroutine file_grads_read_1d_id( &
499  file_id, &
500  var_id, &
501  var, &
502  step, &
503  start, &
504  count, &
505  postfix )
506  implicit none
507  integer, intent(in) :: file_id
508  integer, intent(in) :: var_id
509  real(RP), intent(out) :: var(:)
510  integer, intent(in), optional :: step
511  integer, intent(in), optional :: start(1)
512  integer, intent(in), optional :: count(1)
513  character(len=*), intent(in), optional :: postfix
514 
515  logical :: exist
516  integer :: vid
517 
518  intrinsic :: size
519  !---------------------------------------------------------------------------
520 
521  if ( file_id < 0 ) then
522  log_error("FILE_GrADS_read_1D_vid",*) 'file_id is invalid: ', file_id
523  end if
524  if ( var_id < 0 ) then
525  log_error("FILE_GrADS_read_1D_vid",*) 'var_id is invalid: ', var_id
526  end if
527 
528  call file_grads_read_data( nmls(file_id)%vars(var_id), & ! (in)
529  1, size(var), shape(var), & ! (in)
530  var(:), & ! (out)
531  step = step, & ! (int)
532  start = start, & ! (int)
533  count = count, & ! (int)
534  postfix = postfix ) ! (in)
535 
536  return
537  end subroutine file_grads_read_1d_id
538 
539  !-----------------------------------------------------------------------------
540  subroutine file_grads_read_2d_name( &
541  file_id, &
542  var_name, &
543  var, &
544  step, &
545  start, &
546  count, &
547  postfix )
548  use scale_prc, only: &
549  prc_abort
550  implicit none
551  integer, intent(in) :: file_id
552  character(len=*), intent(in) :: var_name
553  real(RP), intent(out) :: var(:,:)
554  integer, intent(in), optional :: step
555  integer, intent(in), optional :: start(2)
556  integer, intent(in), optional :: count(2)
557  character(len=*), intent(in), optional :: postfix
558 
559  integer :: var_id
560  !---------------------------------------------------------------------------
561 
562  call file_grads_varid( file_id, var_name, & ! (in)
563  var_id ) ! (out)
564 
565  if ( var_id < 0 ) then
566  log_error("FILE_GrADS_read_2D_name",*) 'variable "', trim(var_name), ' is not founed in file "', trim(nmls(file_id)%fname), '"'
567  call prc_abort
568  end if
569 
570  call file_grads_read_2d_id( file_id, var_id, & ! (in)
571  var(:,:), & ! (out)
572  step = step, & ! (in)
573  start = start, & ! (in)
574  count = count, & ! (in)
575  postfix = postfix ) ! (in)
576 
577  return
578  end subroutine file_grads_read_2d_name
579  !-----------------------------------------------------------------------------
580  subroutine file_grads_read_2d_id( &
581  file_id, &
582  var_id, &
583  var, &
584  step, &
585  start, &
586  count, &
587  postfix )
588  implicit none
589  integer, intent(in) :: file_id
590  integer, intent(in) :: var_id
591  real(RP), intent(out) :: var(:,:)
592  integer, intent(in), optional :: step
593  integer, intent(in), optional :: start(2)
594  integer, intent(in), optional :: count(2)
595  character(len=*), intent(in), optional :: postfix
596 
597  integer :: vid
598  intrinsic :: size
599  !---------------------------------------------------------------------------
600 
601  if ( file_id < 0 ) then
602  log_error("FILE_GrADS_read_2D_vid",*) 'file_id is invalid: ', file_id
603  end if
604  if ( var_id < 0 ) then
605  log_error("FILE_GrADS_read_2D_vid",*) 'var_id is invalid: ', var_id
606  end if
607 
608  call file_grads_read_data( nmls(file_id)%vars(var_id), & ! (in)
609  2, size(var), shape(var), & ! (in)
610  var(:,:), & ! (out)
611  step = step, & ! (in)
612  start = start, & ! (in)
613  count = count, & ! (in)
614  postfix = postfix ) ! (in)
615 
616  return
617  end subroutine file_grads_read_2d_id
618 
619  !-----------------------------------------------------------------------------
620  subroutine file_grads_read_3d_name( &
621  file_id, &
622  var_name, &
623  var, &
624  step, &
625  start, &
626  count, &
627  postfix )
628  use scale_prc, only: &
629  prc_abort
630  implicit none
631  integer, intent(in) :: file_id
632  character(len=*), intent(in) :: var_name
633  real(RP), intent(out) :: var(:,:,:)
634  integer, intent(in), optional :: step
635  integer, intent(in), optional :: start(3)
636  integer, intent(in), optional :: count(3)
637  character(len=*), intent(in), optional :: postfix
638 
639  integer :: var_id
640  !---------------------------------------------------------------------------
641 
642  call file_grads_varid( file_id, var_name, & ! (in)
643  var_id ) ! (out)
644 
645  if ( var_id < 0 ) then
646  log_error("FILE_GrADS_read_3D_name",*) 'variable "', trim(var_name), ' is not founed in file "', trim(nmls(file_id)%fname), '"'
647  call prc_abort
648  end if
649 
650  call file_grads_read_3d_id( file_id, var_id, & ! (in)
651  var(:,:,:), & ! (out)
652  step = step, & ! (in)
653  start = start, & ! (in)
654  count = count, & ! (in)
655  postfix = postfix ) ! (in)
656 
657  return
658  end subroutine file_grads_read_3d_name
659  !-----------------------------------------------------------------------------
660  subroutine file_grads_read_3d_id( &
661  file_id, &
662  var_id, &
663  var, &
664  step, &
665  start, &
666  count, &
667  postfix )
668  implicit none
669  integer, intent(in) :: file_id
670  integer, intent(in) :: var_id
671  real(RP), intent(out) :: var(:,:,:)
672  integer, intent(in), optional :: step
673  integer, intent(in), optional :: start(3)
674  integer, intent(in), optional :: count(3)
675  character(len=*), intent(in), optional :: postfix
676 
677  integer :: vid
678  intrinsic :: size
679  !---------------------------------------------------------------------------
680 
681  if ( file_id < 0 ) then
682  log_error("FILE_GrADS_read_3D_vid",*) 'file_id is invalid: ', file_id
683  end if
684  if ( var_id < 0 ) then
685  log_error("FILE_GrADS_read_3D_vid",*) 'var_id is invalid: ', var_id
686  end if
687 
688  call file_grads_read_data( nmls(file_id)%vars(var_id), & ! (in)
689  3, size(var), shape(var), & ! (in)
690  var(:,:,:), & ! (out)
691  step = step, & ! (in)
692  start = start, & ! (in)
693  count = count, & ! (in)
694  postfix = postfix ) ! (in)
695 
696  return
697  end subroutine file_grads_read_3d_id
698 
699  subroutine file_grads_finalize
700  integer :: n
701 
702  do n = 1, nnmls
703  call file_grads_close(n)
704  end do
705  nnmls = 0
706 
707  return
708  end subroutine file_grads_finalize
709 
710  !-----------------------------------------------------------------------------
711  subroutine file_grads_close( &
712  file_id )
713  implicit none
714  integer, intent(in) :: file_id
715 
716  integer :: n, m
717 
718  if ( file_id < 0 ) return
719 
720  do n = 1, nmls(file_id)%nvars
721  do m = 1, nfiles
722  if ( files(m)%fname == nmls(file_id)%vars(n)%fname ) then
723  if ( files(m)%fid > 0 ) then
724  close( files(m)%fid )
725  files(m)%fid = -1
726  files(m)%postfix = ""
727  end if
728  exit
729  end if
730  end do
731  if ( nmls(file_id)%vars(n)%lnum > 0 ) deallocate( nmls(file_id)%vars(n)%lvars )
732  nmls(file_id)%vars(n)%lnum = -1
733  end do
734  deallocate( nmls(file_id)%vars )
735  nmls(file_id)%fname = ""
736  nmls(file_id)%nvars = 0
737 
738  return
739  end subroutine file_grads_close
740 
741  !-----------------------------------------------------------------------------
742  ! private
743  !-----------------------------------------------------------------------------
744 
745  subroutine file_grads_read_data( &
746  var_info, &
747  ndims, n, &
748  shape, &
749  var, &
750  step, &
751  start, &
752  count, &
753  postfix )
754  use scale_prc, only: &
755  prc_abort
756  use scale_const, only: &
757  undef => const_undef, &
758  eps => const_eps
759  implicit none
760  type(t_var), intent(in) :: var_info
761  integer, intent(in) :: ndims
762  integer, intent(in) :: n
763  integer, intent(in) :: shape(ndims)
764  real(rp), intent(out) :: var(n)
765  integer, intent(in), optional :: step
766  integer, intent(in), optional :: start(ndims)
767  integer, intent(in), optional :: count(ndims)
768  character(len=*), intent(in), optional :: postfix
769 
770  integer :: fid
771  character(len=H_LONG) :: gfile
772  real(sp) :: buf(var_info%nx,var_info%ny)
773 
774  integer :: step_
775  integer :: start_(3)
776  integer :: count_(3)
777  character(len=H_SHORT) :: postfix_
778 
779  integer(8) :: pos
780 
781  integer :: nxy, nz, ka
782  integer :: isize, ierr
783  integer :: i, j, k, ii, jj, kk
784 
785  abstract interface
786  subroutine rd( fid, pos, nx, cz, k, sx, cx, cy, yrev, var, ierr )
787  use scale_precision
788  integer, intent(in) :: fid
789  integer(8), intent(in) :: pos
790  integer, intent(in) :: nx, cz, k, sx, cx, cy
791  logical, intent(in) :: yrev
792  real(rp), intent(out) :: var(:)
793  integer, intent(out) :: ierr
794  end subroutine rd
795  end interface
796  procedure(rd), pointer :: read_data
797 
798  select case( var_info%dtype )
799  case("linear")
800 
801  if ( ndims > 1 ) then
802  log_error("FILE_GrADS_read_data",*) '"linear" is invalid for dtype of 2D or 2D var!'
803  end if
804 
805  if( var_info%swpoint == undef .or.var_info%dd == undef ) then
806  log_error("FILE_GrADS_read_data",*) '"swpoint" and "dd" are required for linear data! ', var_info%swpoint
807  call prc_abort
808  endif
809 
810  if ( present(start) ) then
811  start_(1) = start(1)
812  else
813  start_(1) = 1
814  end if
815  if ( present(count) ) then
816  count_(1) = count(1)
817  else
818  count_(1) = n - start_(1) + 1
819  end if
820 
821  do i = 1, count_(1)
822  ii = i + start_(1) - 2 ! from 0
823  var(i) = var_info%swpoint + ii * var_info%dd
824  end do
825 
826  case("levels")
827 
828  if ( ndims > 1 ) then
829  log_error("FILE_GrADS_read_data",*) '"levels" is invalid for dtype of 2D or 2D var!'
830  end if
831 
832  if ( var_info%lnum < 0 )then
833  log_error("FILE_GrADS_read_data",*) '"lnum" is required for levels data! '
834  call prc_abort
835  endif
836  if ( var_info%lnum .ne. n ) then
837  log_error("FILE_GrADS_read_data",*) '"lnum" and size of var are not the same', var_info%lnum, n
838  call prc_abort
839  end if
840 
841  if ( present(start) ) then
842  start_(1) = start(1)
843  else
844  start_(1) = 1
845  end if
846  if ( present(count) ) then
847  count_(1) = count(1)
848  else
849  count_(1) = var_info%lnum - start_(1) + 1
850  end if
851 
852  do k = 1, count_(1)
853  kk = k + start_(1) - 1
854  if ( kk > var_info%lnum ) exit
855  var(k) = var_info%lvars(kk)
856  if( var(k) == undef ) then
857  log_error("FILE_GrADS_read_data",*) '"lvars" must be specified for levels data! '
858  call prc_abort
859  endif
860  end do
861 
862  case("map")
863 
864  if ( present(postfix) ) then
865  postfix_ = postfix
866  else
867  postfix_ = ""
868  end if
869  if ( present(step) ) then
870  step_ = step
871  else
872  step_ = 1
873  end if
874 
875  if ( ndims == 1 ) then
876  log_error("FILE_GrADS_read_data",*) '"map" is invalid for dtype of 1D var!'
877  end if
878 
879  if( var_info%startrec < 0 .or. var_info%totalrec < 0 )then
880  log_error("FILE_GrADS_read_data",*) '"startrec" and "totalrec" are required for map data! ', var_info%startrec, var_info%totalrec
881  call prc_abort
882  endif
883  if( var_info%fname == "" )then
884  log_error("FILE_GrADS_read_data",*) '"fname" is required for map data!'
885  call prc_abort
886  endif
887 
888  ! get file_id
889  fid = -1
890  do i = 1, nfiles
891  if ( files(i)%fname == var_info%fname ) then
892  fid = i
893  exit
894  end if
895  end do
896  if ( fid < 0 ) then
897  nfiles = nfiles + 1
898  if ( nfiles > vars_max ) then
899  log_error("FILE_GrADS_read_data",*) 'The number of files exceeds the limit', vars_max
900  call prc_abort
901  end if
902  fid = nfiles
903  files(fid)%fname = var_info%fname
904  files(fid)%postfix = ""
905  files(fid)%fid = -1
906  end if
907 
908  call io_get_fname(gfile, trim(var_info%fname)//trim(postfix_)//'.grd')
909 
910  if ( files(fid)%postfix == postfix_ .and. files(fid)%fid > 0 ) then
911  fid = files(fid)%fid
912  else
913  if ( files(fid)%fid > 0 ) close( files(fid)%fid )
914  files(fid)%fid = io_get_available_fid()
915  files(fid)%postfix = postfix_
916  fid = files(fid)%fid
917  open( fid, &
918  file = gfile, &
919  form = 'unformatted', &
920  access = 'stream', &
921  status = 'old', &
922  iostat = ierr )
923  if ( ierr /= 0 ) then
924  log_error("FILE_GrADS_read_data",*) 'Failed to open the grads data file! ', trim(gfile)
925  call prc_abort
926  end if
927  end if
928 
929  if ( ndims == 2 ) then
930  nz = 1
931  ka = 1
932  start_(1) = 1
933  if ( present(start) ) then
934  start_(2:3) = start(:)
935  else
936  start_(2:3) = 1
937  end if
938  count_(1) = 1
939  if ( present(count) ) then
940  count_(2:3) = count(:)
941  else
942  count_(2:3) = (/ var_info%nx, var_info%ny /)
943  end if
944  if ( shape(1) .ne. count_(2) .or. shape(2) .ne. count_(3) ) then
945  log_error("FILE_GrADS_read_data",*) 'nx or ny is different', shape(:), count_(2:3)
946  call prc_abort
947  end if
948  else
949  nz = var_info%nz
950  ka = shape(1)
951  if ( present(start) ) then
952  start_(:) = start(:)
953  else
954  start_(:) = 1
955  end if
956  if ( present(count) ) then
957  count_(:) = count(:)
958  else
959  count_(:) = (/ nz, var_info%nx, var_info%ny /)
960  end if
961  if ( ka < count_(1) ) then
962  log_error("FILE_GrADS_read_data",*) 'size is too small for the 1st dimension'
963  call prc_abort
964  end if
965  if ( shape(2) .ne. count_(2) .or. shape(3) .ne. count_(3) ) then
966  log_error("FILE_GrADS_read_data",*) 'nx or ny is different', shape(2:), count_(2:)
967  call prc_abort
968  end if
969  count_(1) = min( count_(1), nz )
970  end if
971  nxy = count_(2) * count_(3)
972 
973 
974  if ( n .ne. nxy * ka ) then
975  log_error("FILE_GrADS_read_data",*) 'size of var is not consitent with namelist info! ', n, ka, count_(2:3)
976  call prc_abort
977  end if
978 
979  select case ( var_info%bintype )
980  case ( 'int1' )
981  read_data => read_data_int1
982  isize = 1
983  case ( 'int2' )
984  read_data => read_data_int2
985  isize = 2
986  case ( 'int4' )
987  read_data => read_data_int4
988  isize = 4
989  case ( 'real4' )
990  read_data => read_data_real4
991  isize = 4
992  case ( 'int8' )
993  read_data => read_data_int8
994  isize = 8
995  case ( 'real8' )
996  read_data => read_data_real8
997  isize = 8
998  case default
999  log_error("FILE_GrADS_read_data",*) 'bintype is invalid for ', trim(var_info%name)
1000  end select
1001 
1002  do k = 1, count_(1)
1003  kk = k + start_(1) - 2 ! from 0
1004  pos = ( var_info%totalrec * ( step_ - 1 ) + var_info%startrec - 1 + kk ) * var_info%ny
1005  if ( var_info%yrev ) then
1006  pos = pos + var_info%ny - ( start_(3) + count_(3) - 1 )
1007  else
1008  pos = pos + start_(3) - 1
1009  end if
1010  pos = pos * var_info%nx * isize + 1
1011  call read_data( fid, pos, ka, var_info%nx, k, start_(2), count_(2), count_(3), var_info%yrev, var(:), ierr )
1012  if ( ierr /= 0 ) then
1013  log_error("FILE_GrADS_read_data",*) 'Failed to read data! ', trim(var_info%name), ', k=',k,', step=',step_, ' in ', trim(gfile)
1014  log_error_cont(*) 'pos=', pos
1015  call prc_abort
1016  end if
1017  end do
1018  if ( var_info%missval .ne. undef ) then
1019  !$omp parallel do &
1020  !$omp private(ii)
1021  do i = 1, nxy
1022  do k = 1, count_(1)
1023  ii = (i-1) * ka + k
1024  if ( abs( var(ii) - var_info%missval ) < eps ) var(ii) = undef
1025  end do
1026  end do
1027  end if
1028  !$omp parallel do &
1029  !$omp private(ii)
1030  do i = 1, nxy
1031  do k = count_(1)+1, ka
1032  ii = (i-1) * ka + k
1033  var(ii) = undef
1034  end do
1035  end do
1036 
1037  end select
1038 
1039  return
1040  end subroutine file_grads_read_data
1041 
1042  subroutine read_data_int1( fid, pos, ka, nx, k, sx, cx, cy, yrev, var, ierr )
1043  integer, intent(in) :: fid
1044  integer(8), intent(in) :: pos
1045  integer, intent(in) :: ka, nx, k, sx, cx, cy
1046  logical, intent(in) :: yrev
1047 
1048  real(RP), intent(out) :: var(:)
1049  integer, intent(out) :: ierr
1050 
1051  integer(1) :: buf(nx,cy)
1052  integer :: i, j
1053 
1054  read(fid, pos=pos, iostat=ierr) buf(:,:)
1055  if ( ierr /= 0 ) return
1056 
1057  if ( yrev ) then
1058  !$omp parallel do
1059  do j = 1, cy
1060  do i = 1, cx
1061  var(k+(i-1)*ka+(j-1)*cx*ka) = buf(sx+i-1,cy-j+1)
1062  end do
1063  end do
1064  else
1065  !$omp parallel do
1066  do j = 1, cy
1067  do i = 1, cx
1068  var(k+(i-1)*ka+(j-1)*cx*ka) = buf(sx+i-1,j)
1069  end do
1070  end do
1071  end if
1072 
1073  return
1074  end subroutine read_data_int1
1075 
1076  subroutine read_data_int2( fid, pos, ka, nx, k, sx, cx, cy, yrev, var, ierr )
1077  integer, intent(in) :: fid
1078  integer(8), intent(in) :: pos
1079  integer, intent(in) :: ka, nx, k, sx, cx, cy
1080  logical, intent(in) :: yrev
1081 
1082  real(RP), intent(out) :: var(:)
1083  integer, intent(out) :: ierr
1084 
1085  integer(2) :: buf(nx,cy)
1086  integer :: i, j
1087 
1088  read(fid, pos=pos, iostat=ierr) buf(:,:)
1089  if ( ierr /= 0 ) return
1090 
1091  if ( yrev ) then
1092  !$omp parallel do
1093  do j = 1, cy
1094  do i = 1, cx
1095  var(k+(i-1)*ka+(j-1)*cx*ka) = buf(sx+i-1,cy-j+1)
1096  end do
1097  end do
1098  else
1099  !$omp parallel do
1100  do j = 1, cy
1101  do i = 1, cx
1102  var(k+(i-1)*ka+(j-1)*cx*ka) = buf(sx+i-1,j)
1103  end do
1104  end do
1105  end if
1106 
1107  return
1108  end subroutine read_data_int2
1109 
1110  subroutine read_data_int4( fid, pos, ka, nx, k, sx, cx, cy, yrev, var, ierr )
1111  integer, intent(in) :: fid
1112  integer(8), intent(in) :: pos
1113  integer, intent(in) :: ka, nx, k, sx, cx, cy
1114  logical, intent(in) :: yrev
1115 
1116  real(RP), intent(out) :: var(:)
1117  integer, intent(out) :: ierr
1118 
1119  integer(4) :: buf(nx,cy)
1120  integer :: i, j
1121 
1122  read(fid, pos=pos, iostat=ierr) buf(:,:)
1123  if ( ierr /= 0 ) return
1124 
1125  if ( yrev ) then
1126  !$omp parallel do
1127  do j = 1, cy
1128  do i = 1, cx
1129  var(k+(i-1)*ka+(j-1)*cx*ka) = buf(sx+i-1,cy-j+1)
1130  end do
1131  end do
1132  else
1133  !$omp parallel do
1134  do j = 1, cy
1135  do i = 1, cx
1136  var(k+(i-1)*ka+(j-1)*cx*ka) = buf(sx+i-1,j)
1137  end do
1138  end do
1139  end if
1140 
1141  return
1142  end subroutine read_data_int4
1143 
1144  subroutine read_data_real4( fid, pos, ka, nx, k, sx, cx, cy, yrev, var, ierr )
1145  integer, intent(in) :: fid
1146  integer(8), intent(in) :: pos
1147  integer, intent(in) :: ka, nx, k, sx, cx, cy
1148  logical, intent(in) :: yrev
1149 
1150  real(RP), intent(out) :: var(:)
1151  integer, intent(out) :: ierr
1152 
1153  real(4) :: buf(nx,cy)
1154  integer :: i, j
1155 
1156  read(fid, pos=pos, iostat=ierr) buf(:,:)
1157  if ( ierr /= 0 ) return
1158 
1159  if ( yrev ) then
1160  !$omp parallel do
1161  do j = 1, cy
1162  do i = 1, cx
1163  var(k+(i-1)*ka+(j-1)*cx*ka) = buf(sx+i-1,cy-j+1)
1164  end do
1165  end do
1166  else
1167  !$omp parallel do
1168  do j = 1, cy
1169  do i = 1, cx
1170  var(k+(i-1)*ka+(j-1)*cx*ka) = buf(sx+i-1,j)
1171  end do
1172  end do
1173  end if
1174 
1175  return
1176  end subroutine read_data_real4
1177 
1178  subroutine read_data_int8( fid, pos, ka, nx, k, sx, cx, cy, yrev, var, ierr )
1179  integer, intent(in) :: fid
1180  integer(8), intent(in) :: pos
1181  integer, intent(in) :: ka, nx, k, sx, cx, cy
1182  logical, intent(in) :: yrev
1183 
1184  real(RP), intent(out) :: var(:)
1185  integer, intent(out) :: ierr
1186 
1187  integer(8) :: buf(nx,cy)
1188  integer :: i, j
1189 
1190  read(fid, pos=pos, iostat=ierr) buf(:,:)
1191  if ( ierr /= 0 ) return
1192 
1193  if ( yrev ) then
1194  !$omp parallel do
1195  do j = 1, cy
1196  do i = 1, cx
1197  var(k+(i-1)*ka+(j-1)*cx*ka) = buf(sx+i-1,cy-j+1)
1198  end do
1199  end do
1200  else
1201  !$omp parallel do
1202  do j = 1, cy
1203  do i = 1, cx
1204  var(k+(i-1)*ka+(j-1)*cx*ka) = buf(sx+i-1,j)
1205  end do
1206  end do
1207  end if
1208 
1209  return
1210  end subroutine read_data_int8
1211 
1212  subroutine read_data_real8( fid, pos, ka, nx, k, sx, cx, cy, yrev, var, ierr )
1213  integer, intent(in) :: fid
1214  integer(8), intent(in) :: pos
1215  integer, intent(in) :: ka, nx, k, sx, cx, cy
1216  logical, intent(in) :: yrev
1217 
1218  real(RP), intent(out) :: var(:)
1219  integer, intent(out) :: ierr
1220 
1221  real(8) :: buf(nx,cy)
1222  integer :: i, j
1223 
1224  read(fid, pos=pos, iostat=ierr) buf(:,:)
1225  if ( ierr /= 0 ) return
1226 
1227  if ( yrev ) then
1228  !$omp parallel do
1229  do j = 1, cy
1230  do i = 1, cx
1231  var(k+(i-1)*ka+(j-1)*cx*ka) = buf(sx+i-1,cy-j+1)
1232  end do
1233  end do
1234  else
1235  !$omp parallel do
1236  do j = 1, cy
1237  do i = 1, cx
1238  var(k+(i-1)*ka+(j-1)*cx*ka) = buf(sx+i-1,j)
1239  end do
1240  end do
1241  end if
1242 
1243  return
1244  end subroutine read_data_real8
1245 
1246  subroutine check_oldnamelist( fid )
1247  use scale_prc, only: &
1248  prc_abort
1249  implicit none
1250  integer, intent(in) :: fid
1251 
1252  integer :: ierr
1253  logical :: dummy
1254 
1255  namelist /nml_grads_grid/ dummy
1256  namelist /grdvar/ dummy
1257 
1258  read(fid, nml=nml_grads_grid, iostat=ierr)
1259  if( ierr > 0 )then
1260  log_error("check_oldnamelist",*) 'The old namelist "nml_grads_grid" is found.'
1261  log_error_cont(*) 'Use "GrADS_DIMS" instead.'
1262  call prc_abort
1263  endif
1264  rewind(fid)
1265 
1266  read(fid, nml=grdvar, iostat=ierr)
1267  if( ierr > 0 )then
1268  log_error("check_oldnamelist",*) 'The old namelist "grdvar" is found.'
1269  log_error_cont(*) 'Use "GrADS_ITEM" instead.'
1270  call prc_abort
1271  endif
1272  rewind(fid)
1273 
1274  return
1275  end subroutine check_oldnamelist
1276 
1277  ! convert string to upcase
1278  function upcase( str )
1279  character(len=*), intent(in) :: str
1280  character(len=len(str)) :: upcase
1281  integer :: i
1282  do i = 1, len_trim(str)
1283  if ( str(i:i) >= 'a' .and. str(i:i) <= 'z' ) then
1284  upcase(i:i) = char( ichar(str(i:i)) - 32 )
1285  else
1286  upcase(i:i) = str(i:i)
1287  end if
1288  end do
1289  do i = len_trim(str)+1, len(str)
1290  upcase(i:i) = ' '
1291  end do
1292  end function upcase
1293 
1294 end module scale_file_grads
scale_precision::sp
integer, parameter, public sp
Definition: scale_precision.F90:31
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_file_grads::file_grads_isoned
logical function, public file_grads_isoned(file_id, var_id)
Definition: scale_file_grads.F90:373
scale_file_grads::file_grads_open
subroutine, public file_grads_open(file_name, file_id)
Open.
Definition: scale_file_grads.F90:104
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_file_grads::file_grads_read_1d_id
subroutine file_grads_read_1d_id(file_id, var_id, var, step, start, count, postfix)
Definition: scale_file_grads.F90:506
scale_file_grads::upcase
character(len=len(str)) function upcase(str)
Definition: scale_file_grads.F90:1279
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:35
scale_file_grads::file_grads_varcheck
subroutine, public file_grads_varcheck(file_id, var_name, exist)
Definition: scale_file_grads.F90:344
scale_file_grads::file_grads_varid
subroutine, public file_grads_varid(file_id, var_name, var_id)
Definition: scale_file_grads.F90:312
scale_io::io_get_available_fid
integer function, public io_get_available_fid()
search & get available file ID
Definition: scale_io.F90:373
scale_file_grads::file_grads_finalize
subroutine, public file_grads_finalize
Definition: scale_file_grads.F90:700
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_grads::file_grads_read_3d_id
subroutine file_grads_read_3d_id(file_id, var_id, var, step, start, count, postfix)
Definition: scale_file_grads.F90:668
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_file_grads
module file_grads
Definition: scale_file_grads.F90:11
scale_file_grads::file_grads_get_shape_id
subroutine file_grads_get_shape_id(file_id, var_id, shape)
Definition: scale_file_grads.F90:431
scale_file_grads::file_grads_close
subroutine, public file_grads_close(file_id)
Definition: scale_file_grads.F90:713
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_file_grads::file_grads_get_shape_name
subroutine file_grads_get_shape_name(file_id, var_name, shape)
Definition: scale_file_grads.F90:404
scale_file_grads::read_data_int1
subroutine read_data_int1(fid, pos, ka, nx, k, sx, cx, cy, yrev, var, ierr)
Definition: scale_file_grads.F90:1043
scale_file_grads::file_grads_read_2d_id
subroutine file_grads_read_2d_id(file_id, var_id, var, step, start, count, postfix)
Definition: scale_file_grads.F90:588