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