29 public :: file_grads_get_shape
30 public :: file_grads_read
33 interface file_grads_get_shape
36 end interface file_grads_get_shape
38 interface file_grads_read
39 module procedure file_grads_read_1d_name
41 module procedure file_grads_read_2d_name
43 module procedure file_grads_read_3d_name
45 end interface file_grads_read
58 integer,
parameter :: nmls_max = 10
59 integer,
parameter :: vars_max = 100
60 integer,
parameter :: lvars_max = 1000
62 character(len=H_SHORT) :: name
63 character(len=H_LONG) :: fname
64 character(len=H_SHORT) :: dtype
68 real(RP),
allocatable :: lvars(:)
77 character(len=H_SHORT) :: bintype
80 character(len=H_LONG) :: fname
82 type(t_var),
allocatable :: vars(:)
86 character(len=H_LONG) :: fname
87 character(len=H_SHORT) :: postfix
90 type(t_nml),
save :: nmls(nmls_max)
92 type(t_file) :: files(vars_max)
107 character(len=*),
intent(in) :: file_name
108 integer,
intent(out) :: file_id
110 character(len=H_SHORT) :: name
111 character(len=H_SHORT) :: dtype
112 character(len=H_LONG) :: fname
116 real(
rp) :: lvars(lvars_max)
123 character(len=H_SHORT) :: fendian
125 character(len=H_SHORT) :: bintype
127 namelist /grads_dims/ &
132 namelist /grads_item/ &
150 character(len=H_LONG) :: dirname
159 log_info(
"FILE_GrADS_open",*)
'open namelist file :', trim(file_name)
163 if ( nmls(n)%fname == file_name )
then
174 file = trim(file_name), &
175 form =
'formatted', &
179 if ( ierr /= 0 )
then
180 log_error(
"FILE_GrADS_open",*)
'Input file is not found! ', trim(file_name)
184 call check_oldnamelist( fid )
187 read(fid,nml=grads_dims,iostat=ierr)
189 log_error(
"FILE_GrADS_open",*)
'Not appropriate names in GrADS_DIMS in ', trim(file_name),
'. Check!'
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
209 read(fid, nml=grads_item, iostat=ierr)
211 log_error(
"FILE_GrADS_open",*)
'Not appropriate names in GrADS_ITEM in ', trim(file_name),
'. Check!'
213 else if( ierr < 0 )
then
219 if ( nvars > vars_max )
then
220 log_error(
"FILE_GRADS_open",*)
'The number of grads vars exceeds the limit! ', &
221 nvars,
' > ', vars_max
225 nmls(file_id)%nvars = nvars
226 allocate( nmls(file_id)%vars(nvars) )
228 n = index( file_name,
'/', back=.true. )
230 dirname = file_name(1:n)
249 nx = nmls(file_id)%nx
250 ny = nmls(file_id)%ny
251 nz = nmls(file_id)%nz
258 read(fid, nml=grads_item, iostat=ierr)
261 nmls(file_id)%vars(n)%name =
upcase(name)
262 if ( fname(1:1) ==
"/" )
then
263 nmls(file_id)%vars(n)%fname = fname
265 nmls(file_id)%vars(n)%fname = trim(dirname) // fname
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
272 allocate( nmls(file_id)%vars(n)%lvars(lnum) )
273 nmls(file_id)%vars(n)%lvars(:) = lvars(1:lnum)
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
285 nmls(file_id)%vars(n)%endian = 0
287 nmls(file_id)%vars(n)%bintype = bintype
303 integer,
intent(in) :: file_id
304 character(len=*),
intent(in) :: var_name
305 integer,
intent(out) :: var_id
307 character(len=len(var_name)) :: vname
310 if ( file_id < 0 )
then
311 log_error(
"FILE_GrADS_varid",*)
'file_id is invalid: ', file_id
318 do n = 1, nmls(file_id)%nvars
319 if ( nmls(file_id)%vars(n)%name == vname )
then
335 integer,
intent(in) :: file_id
336 integer,
intent(in) :: var_id
339 if ( file_id < 0 )
then
340 log_error(
"FILE_GrADS_isOneD",*)
'file_id is invalid: ', file_id
343 if ( var_id < 0 )
then
344 log_error(
"FILE_GrADS_isOneD",*)
'var_id is invalid: ', var_id
348 select case( nmls(file_id)%vars(var_id)%dtype )
349 case (
'linear',
'levels')
366 integer,
intent(in) :: file_id
367 character(len=*),
intent(in) :: var_name
368 integer,
intent(out) :: shape(:)
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),
'"'
391 integer,
intent(in) :: file_id
392 integer,
intent(in) :: var_id
393 integer,
intent(out) :: shape(:)
398 if ( nmls(file_id)%vars(var_id)%dtype ==
"levels" )
then
399 shape(1) = nmls(file_id)%vars(var_id)%lnum
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
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
416 subroutine file_grads_read_1d_name( &
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
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),
'"'
449 end subroutine file_grads_read_1d_name
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
470 if ( file_id < 0 )
then
471 log_error(
"FILE_GrADS_read_1D_vid",*)
'file_id is invalid: ', file_id
473 if ( var_id < 0 )
then
474 log_error(
"FILE_GrADS_read_1D_vid",*)
'var_id is invalid: ', var_id
477 call file_grads_read_data( nmls(file_id)%vars(var_id), &
487 subroutine file_grads_read_2d_name( &
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
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),
'"'
519 end subroutine file_grads_read_2d_name
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
538 if ( file_id < 0 )
then
539 log_error(
"FILE_GrADS_read_2D_vid",*)
'file_id is invalid: ', file_id
541 if ( var_id < 0 )
then
542 log_error(
"FILE_GrADS_read_2D_vid",*)
'var_id is invalid: ', var_id
545 call file_grads_read_data( nmls(file_id)%vars(var_id), &
555 subroutine file_grads_read_3d_name( &
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
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),
'"'
587 end subroutine file_grads_read_3d_name
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
606 if ( file_id < 0 )
then
607 log_error(
"FILE_GrADS_read_3D_vid",*)
'file_id is invalid: ', file_id
609 if ( var_id < 0 )
then
610 log_error(
"FILE_GrADS_read_3D_vid",*)
'var_id is invalid: ', var_id
613 call file_grads_read_data( nmls(file_id)%vars(var_id), &
616 step = step, postfix = postfix )
625 integer,
intent(in) :: file_id
629 if ( file_id < 0 )
return
631 do n = 1, nmls(file_id)%nvars
633 if ( files(m)%fname == nmls(file_id)%vars(n)%fname )
then
634 if ( files(m)%fid > 0 )
then
635 close( files(m)%fid )
637 files(m)%postfix =
""
642 if ( nmls(file_id)%vars(n)%lnum > 0 )
deallocate( nmls(file_id)%vars(n)%lvars )
643 nmls(file_id)%vars(n)%lnum = -1
645 deallocate( nmls(file_id)%vars )
646 nmls(file_id)%fname =
""
647 nmls(file_id)%nvars = 0
656 subroutine file_grads_read_data( &
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
676 character(len=H_LONG) :: gfile
677 real(
sp) :: buf(var_info%nx,var_info%ny)
680 character(len=H_SHORT) :: postfix_
683 integer :: irecl, isize, ierr
687 subroutine rd( fid, irecl, nx, ny, nz, k, yrev, var, ierr )
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
696 procedure(rd),
pointer :: read_data
698 select case( var_info%dtype )
701 if ( ndims > 1 )
then
702 log_error(
"FILE_GrADS_read_data",*)
'"linear" is invalid for dtype of 2D or 2D var!'
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
711 var(i) = var_info%swpoint + (i-1) * var_info%dd
716 if ( ndims > 1 )
then
717 log_error(
"FILE_GrADS_read_data",*)
'"levels" is invalid for dtype of 2D or 2D var!'
720 if ( var_info%lnum < 0 )
then
721 log_error(
"FILE_GrADS_read_data",*)
'"lnum" is required for levels data! '
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
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! '
739 if (
present(postfix) )
then
744 if (
present(step) )
then
750 if ( ndims == 1 )
then
751 log_error(
"FILE_GrADS_read_data",*)
'"map" is invalid for dtype of 1D var!'
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
758 if( var_info%fname ==
"" )
then
759 log_error(
"FILE_GrADS_read_data",*)
'"fname" is required for map data!'
766 if ( files(i)%fname == var_info%fname )
then
774 files(fid)%fname = var_info%fname
775 files(fid)%postfix =
""
779 gfile = trim(var_info%fname)//trim(postfix_)//
'.grd'
781 if ( files(fid)%postfix == postfix_ .and. files(fid)%fid > 0 )
then
784 if ( files(fid)%fid > 0 )
close( files(fid)%fid )
786 files(fid)%postfix = postfix_
788 select case ( var_info%bintype )
793 case (
'real4',
'int4' )
795 case (
'real8',
'int8' )
798 log_error(
"FILE_GrADS_read_data",*)
'bintype is invalid for ', trim(var_info%name)
800 irecl = var_info%nx * var_info%ny * isize
803 form =
'unformatted', &
808 if ( ierr /= 0 )
then
809 log_error(
"FILE_GrADS_read_data",*)
'Failed to open the grads data file! ', trim(gfile)
814 if ( ndims == 2 )
then
819 nxy = var_info%nx * var_info%ny
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
826 select case ( var_info%bintype )
830 read_data => read_data_int2
832 read_data => read_data_int4
834 read_data => read_data_real4
836 read_data => read_data_int8
838 read_data => read_data_real8
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
850 if ( var_info%missval .ne. undef )
then
853 if ( abs( var(i) - var_info%missval ) < eps ) var(i) = undef
860 end subroutine file_grads_read_data
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
867 real(RP),
intent(out) :: var(:)
868 integer,
intent(out) :: ierr
870 integer(1) :: buf(nx,ny)
873 read(fid, rec=irecl, iostat=ierr) buf(:,:)
874 if ( ierr /= 0 )
return
880 var(k+(i-1)*nz+(j-1)*nx*nz) = buf(i,ny-j+1)
887 var(k+(i-1)*nz+(j-1)*nx*nz) = buf(i,j)
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
900 real(RP),
intent(out) :: var(:)
901 integer,
intent(out) :: ierr
903 integer(2) :: buf(nx,ny)
906 read(fid, rec=irecl, iostat=ierr) buf(:,:)
907 if ( ierr /= 0 )
return
913 var(k+(i-1)*nz+(j-1)*nx*nz) = buf(i,ny-j+1)
920 var(k+(i-1)*nz+(j-1)*nx*nz) = buf(i,j)
926 end subroutine read_data_int2
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
933 real(RP),
intent(out) :: var(:)
934 integer,
intent(out) :: ierr
936 integer(4) :: buf(nx,ny)
939 read(fid, rec=irecl, iostat=ierr) buf(:,:)
940 if ( ierr /= 0 )
return
946 var(k+(i-1)*nz+(j-1)*nx*nz) = buf(i,ny-j+1)
953 var(k+(i-1)*nz+(j-1)*nx*nz) = buf(i,j)
959 end subroutine read_data_int4
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
966 real(RP),
intent(out) :: var(:)
967 integer,
intent(out) :: ierr
969 real(4) :: buf(nx,ny)
972 read(fid, rec=irecl, iostat=ierr) buf(:,:)
973 if ( ierr /= 0 )
return
979 var(k+(i-1)*nz+(j-1)*nx*nz) = buf(i,ny-j+1)
986 var(k+(i-1)*nz+(j-1)*nx*nz) = buf(i,j)
992 end subroutine read_data_real4
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
999 real(RP),
intent(out) :: var(:)
1000 integer,
intent(out) :: ierr
1002 integer(8) :: buf(nx,ny)
1005 read(fid, rec=irecl, iostat=ierr) buf(:,:)
1006 if ( ierr /= 0 )
return
1012 var(k+(i-1)*nz+(j-1)*nx*nz) = buf(i,ny-j+1)
1019 var(k+(i-1)*nz+(j-1)*nx*nz) = buf(i,j)
1025 end subroutine read_data_int8
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
1032 real(RP),
intent(out) :: var(:)
1033 integer,
intent(out) :: ierr
1035 real(8) :: buf(nx,ny)
1038 read(fid, rec=irecl, iostat=ierr) buf(:,:)
1039 if ( ierr /= 0 )
return
1045 var(k+(i-1)*nz+(j-1)*nx*nz) = buf(i,ny-j+1)
1052 var(k+(i-1)*nz+(j-1)*nx*nz) = buf(i,j)
1058 end subroutine read_data_real8
1060 subroutine check_oldnamelist( fid )
1064 integer,
intent(in) :: fid
1069 namelist /nml_grads_grid/ dummy
1070 namelist /grdvar/ dummy
1072 read(fid, nml=nml_grads_grid, iostat=ierr)
1074 log_error(
"check_oldnamelist",*)
'The old namelist "nml_grads_grid" is found.'
1075 log_error_cont(*)
'Use "GrADS_DIMS" instead.'
1080 read(fid, nml=grdvar, iostat=ierr)
1082 log_error(
"check_oldnamelist",*)
'The old namelist "grdvar" is found.'
1083 log_error_cont(*)
'Use "GrADS_ITEM" instead.'
1089 end subroutine check_oldnamelist
1093 character(len=*),
intent(in) :: str
1094 character(len=len(str)) ::
upcase
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 )
1103 do i = len_trim(str)+1, len(str)