NCEPLIBS-g2  3.4.7
getgb2r.F90
Go to the documentation of this file.
1 
4 
39 SUBROUTINE getgb2r(LUGB, CINDEX, GFLD, IRET)
40  use grib_mod
41  implicit none
42 
43  INTEGER, INTENT(IN) :: LUGB
44  CHARACTER(LEN=1), INTENT(IN) :: CINDEX(*)
45  TYPE(gribfield) :: GFLD
46  INTEGER, INTENT(OUT) :: IRET
47 
48  INTEGER :: LSKIP, SKIP6, SKIP7
49  CHARACTER(LEN=1):: CSIZE(4)
50  CHARACTER(LEN=1), ALLOCATABLE :: CTEMP(:)
51  real, pointer, dimension(:) :: newfld
52  integer :: n, lread, j, iskip, iofst, ilen, ierr, idum
53 
54  interface
55  subroutine gf_unpack6(cgrib, lcgrib, iofst, ngpts, ibmap, bmap, ierr)
56  character(len=1), intent(in) :: cgrib(lcgrib)
57  integer, intent(in) :: lcgrib, ngpts
58  integer, intent(inout) :: iofst
59  integer, intent(out) :: ibmap
60  integer, intent(out) :: ierr
61  logical*1, pointer, dimension(:) :: bmap
62  end subroutine gf_unpack6
63  subroutine gf_unpack7(cgrib, lcgrib, iofst, igdsnum, igdstmpl, &
64  idrsnum, idrstmpl, ndpts, fld, ierr)
65  character(len=1), intent(in) :: cgrib(lcgrib)
66  integer, intent(in) :: lcgrib, ndpts, idrsnum, igdsnum
67  integer, intent(inout) :: iofst
68  integer, pointer, dimension(:) :: idrstmpl, igdstmpl
69  integer, intent(out) :: ierr
70  real, pointer, dimension(:) :: fld
71  end subroutine gf_unpack7
72  end interface
73 
74  ! Get info.
75  NULLIFY(gfld%bmap, gfld%fld)
76  iret = 0
77  CALL g2_gbytec(cindex, lskip, 4*8, 4*8)
78  CALL g2_gbytec(cindex, skip6, 24*8, 4*8)
79  CALL g2_gbytec(cindex, skip7, 28*8, 4*8)
80 
81  ! Read and unpack bit_map, if present.
82  IF (gfld%ibmap .eq. 0 .OR. gfld%ibmap .eq. 254) THEN
83  iskip = lskip + skip6
84 
85  ! Get length of section.
86  CALL baread(lugb, iskip, 4, lread, csize)
87  CALL g2_gbytec(csize, ilen, 0, 32)
88  ALLOCATE(ctemp(ilen))
89 
90  ! Read in section.
91  CALL baread(lugb, iskip, ilen, lread, ctemp)
92  IF (ilen .NE. lread) THEN
93  iret = 97
94  DEALLOCATE(ctemp)
95  RETURN
96  ENDIF
97  iofst = 0
98  CALL gf_unpack6(ctemp, ilen, iofst, gfld%ngrdpts, idum, gfld%bmap, ierr)
99  IF (ierr .NE. 0) THEN
100  iret = 98
101  DEALLOCATE(ctemp)
102  RETURN
103  ENDIF
104  DEALLOCATE(ctemp)
105  ENDIF
106 
107  ! Read and unpack data field.
108  iskip = lskip + skip7
109 
110  ! Get length of section.
111  CALL baread(lugb, iskip, 4, lread, csize)
112  CALL g2_gbytec(csize, ilen, 0, 32)
113  if (ilen .lt. 6) ilen = 6
114  ALLOCATE(ctemp(ilen))
115 
116  ! Read in section.
117  CALL baread(lugb, iskip, ilen, lread, ctemp)
118  IF (ilen .NE. lread) THEN
119  iret = 97
120  DEALLOCATE(ctemp)
121  RETURN
122  ENDIF
123  iofst = 0
124  CALL gf_unpack7(ctemp, ilen, iofst, gfld%igdtnum, gfld%igdtmpl, &
125  gfld%idrtnum, gfld%idrtmpl, gfld%ndpts, gfld%fld, ierr)
126  IF (ierr .NE. 0) THEN
127  iret = 98
128  DEALLOCATE(ctemp)
129  RETURN
130  ENDIF
131  DEALLOCATE(ctemp)
132 
133  ! If bitmap is used with this field, expand data field
134  ! to grid, if possible.
135  if (gfld%ibmap .ne. 255 .AND. associated(gfld%bmap)) then
136  allocate(newfld(gfld%ngrdpts))
137  n = 1
138  do j = 1, gfld%ngrdpts
139  if (gfld%bmap(j)) then
140  newfld(j) = gfld%fld(n)
141  n = n+1
142  else
143  newfld(j) = 0.0
144  endif
145  enddo
146  deallocate(gfld%fld);
147  gfld%fld => newfld;
148  gfld%expanded = .true.
149  else
150  gfld%expanded = .true.
151  endif
152 END SUBROUTINE getgb2r
subroutine g2_gbytec(in, iout, iskip, nbits)
Extract arbitrary size values from a packed bit string, right justifying each value in the unpacked a...
Definition: g2_gbytesc.F90:20
subroutine getgb2r(LUGB, CINDEX, GFLD, IRET)
Read and unpack sections 6 and 7 from a GRIB2 message.
Definition: getgb2r.F90:40
subroutine gf_unpack6(cgrib, lcgrib, iofst, ngpts, ibmap, bmap, ierr)
Unpack Section 6 (Bit-Map Section) of a GRIB2 message, starting at octet 6 of that Section.
Definition: gf_unpack6.F90:30
subroutine gf_unpack7(cgrib, lcgrib, iofst, igdsnum, igdstmpl, idrsnum, idrstmpl, ndpts, fld, ierr)
Unpack Section 7 (Data Section) of a GRIB2 message.
Definition: gf_unpack7.F90:42
This Fortran module contains the declaration of derived type gribfield.
Definition: gribmod.F90:10