NCEPLIBS-g2  3.4.8
ixgb2.F90
Go to the documentation of this file.
1 
4 
55 SUBROUTINE ixgb2(LUGB, LSKIP, LGRIB, CBUF, NUMFLD, MLEN, IRET)
56  USE re_alloc ! NEEDED FOR SUBROUTINE REALLOC
57  implicit none
58 
59  CHARACTER(LEN = 1), POINTER, DIMENSION(:) :: CBUF
60  CHARACTER CVER, CDISC
61  CHARACTER(LEN = 4) :: CTEMP
62  INTEGER LOCLUS, LOCGDS, LENGDS, LOCBMS
63  integer :: indbmp, numsec, next, newsize, mova2i, mbuf, lindex
64  integer :: linmax, ixskp
65  integer :: mxspd, mxskp, mxsgd, mxsdr, mxsbm, mxlus
66  integer :: mxlen, mxds, mxfld, mxbms
67  integer :: init, ixlus, lugb, lskip, lgrib, numfld, mlen, iret
68  integer :: ixsgd, ibread, ibskip, ilndrs, ilnpds, istat, ixds
69  integer :: ixspd, ixfld, ixids, ixlen, ixsbm, ixsdr
70  integer :: lbread, lensec, lensec1
71  parameter(linmax = 5000, init = 50000, next = 10000)
72  parameter(ixskp = 4, ixlus = 8, ixsgd = 12, ixspd = 16, ixsdr = 20, ixsbm = 24, &
73  ixds = 28, ixlen = 36, ixfld = 42, ixids = 44)
74  parameter(mxskp = 4, mxlus = 4, mxsgd = 4, mxspd = 4, mxsdr = 4, mxsbm = 4, &
75  mxds = 4, mxlen = 4, mxfld = 2, mxbms = 6)
76  CHARACTER CBREAD(LINMAX), CINDEX(LINMAX)
77  CHARACTER CIDS(LINMAX), CGDS(LINMAX)
78 
79  loclus = 0
80  iret = 0
81  mlen = 0
82  numfld = 0
83  NULLIFY(cbuf)
84  mbuf = init
85  ALLOCATE(cbuf(mbuf), stat = istat) ! ALLOCATE INITIAL SPACE FOR CBUF
86  IF (istat .NE. 0) THEN
87  iret = 1
88  RETURN
89  ENDIF
90 
91  ! Read sections 0 and 1 for versin number and discipline.
92  ibread = min(lgrib, linmax)
93  CALL baread(lugb, lskip, ibread, lbread, cbread)
94  IF(lbread .NE. ibread) THEN
95  iret = 2
96  RETURN
97  ENDIF
98  IF(cbread(8) .NE. char(2)) THEN ! NOT GRIB EDITION 2
99  iret = 3
100  RETURN
101  ENDIF
102  cver = cbread(8)
103  cdisc = cbread(7)
104  CALL g2_gbytec(cbread, lensec1, 16 * 8, 4 * 8)
105  lensec1 = min(lensec1, ibread)
106  cids(1:lensec1) = cbread(17:16 + lensec1)
107  ibskip = lskip + 16 + lensec1
108 
109  ! Loop through remaining sections creating an index for each field.
110  ibread = max(5, mxbms)
111  DO
112  CALL baread(lugb, ibskip, ibread, lbread, cbread)
113  ctemp = cbread(1)//cbread(2)//cbread(3)//cbread(4)
114  IF (ctemp .EQ. '7777') RETURN ! END OF MESSAGE FOUND
115  IF(lbread .NE. ibread) THEN
116  iret = 2
117  RETURN
118  ENDIF
119  CALL g2_gbytec(cbread, lensec, 0 * 8, 4 * 8)
120  CALL g2_gbytec(cbread, numsec, 4 * 8, 1 * 8)
121 
122  IF (numsec .EQ. 2) THEN ! SAVE LOCAL USE LOCATION
123  loclus = ibskip-lskip
124  ELSEIF (numsec .EQ. 3) THEN ! SAVE GDS INFO
125  lengds = lensec
126  cgds = char(0)
127  CALL baread(lugb, ibskip, lengds, lbread, cgds)
128  IF (lbread .NE. lengds) THEN
129  iret = 2
130  RETURN
131  ENDIF
132  locgds = ibskip-lskip
133  ELSEIF (numsec .EQ. 4) THEN ! FOUND PDS
134  cindex = char(0)
135  CALL g2_sbytec(cindex, lskip, 8 * ixskp, 8 * mxskp) ! BYTES TO SKIP
136  CALL g2_sbytec(cindex, loclus, 8 * ixlus, 8 * mxlus) ! LOCATION OF LOCAL USE
137  CALL g2_sbytec(cindex, locgds, 8 * ixsgd, 8 * mxsgd) ! LOCATION OF GDS
138  CALL g2_sbytec(cindex, ibskip-lskip, 8 * ixspd, 8 * mxspd) ! LOCATION OF PDS
139  CALL g2_sbytec(cindex, lgrib, 8 * ixlen, 8 * mxlen) ! LEN OF GRIB2
140  cindex(41) = cver
141  cindex(42) = cdisc
142  CALL g2_sbytec(cindex, numfld + 1, 8 * ixfld, 8 * mxfld) ! FIELD NUM
143  cindex(ixids + 1:ixids + lensec1) = cids(1:lensec1)
144  lindex = ixids + lensec1
145  cindex(lindex + 1:lindex + lengds) = cgds(1:lengds)
146  lindex = lindex + lengds
147  ilnpds = lensec
148  CALL baread(lugb, ibskip, ilnpds, lbread, cindex(lindex + 1))
149  IF (lbread .NE. ilnpds) THEN
150  iret = 2
151  RETURN
152  ENDIF
153  lindex = lindex + ilnpds
154  ELSEIF (numsec .EQ. 5) THEN ! FOUND DRS
155  CALL g2_sbytec(cindex, ibskip-lskip, 8 * ixsdr, 8 * mxsdr) ! LOCATION OF DRS
156  ilndrs = lensec
157  CALL baread(lugb, ibskip, ilndrs, lbread, cindex(lindex + 1))
158  IF (lbread .NE. ilndrs) THEN
159  iret = 2
160  RETURN
161  ENDIF
162  lindex = lindex + ilndrs
163  ELSEIF (numsec .EQ. 6) THEN ! FOUND BMS
164  indbmp = mova2i(cbread(6))
165  IF (indbmp.LT.254) THEN
166  locbms = ibskip-lskip
167  CALL g2_sbytec(cindex, locbms, 8 * ixsbm, 8 * mxsbm) ! LOC. OF BMS
168  ELSEIF (indbmp.EQ.254) THEN
169  CALL g2_sbytec(cindex, locbms, 8 * ixsbm, 8 * mxsbm) ! LOC. OF BMS
170  ELSEIF (indbmp.EQ.255) THEN
171  CALL g2_sbytec(cindex, ibskip-lskip, 8 * ixsbm, 8 * mxsbm) ! LOC. OF BMS
172  ENDIF
173  cindex(lindex + 1:lindex + mxbms) = cbread(1:mxbms)
174  lindex = lindex + mxbms
175  CALL g2_sbytec(cindex, lindex, 0, 8 * 4) ! NUM BYTES IN INDEX RECORD
176  ELSEIF (numsec .EQ. 7) THEN ! FOUND DATA SECTION
177  CALL g2_sbytec(cindex, ibskip-lskip, 8 * ixds, 8 * mxds) ! LOC. OF DATA SEC.
178  numfld = numfld + 1
179  IF ((lindex + mlen) .GT. mbuf) THEN ! ALLOCATE MORE SPACE IF NECESSARY
180  newsize = max(mbuf + next, mbuf + lindex)
181  CALL realloc(cbuf, mlen, newsize, istat)
182  IF (istat .NE. 0) THEN
183  numfld = numfld-1
184  iret = 4
185  RETURN
186  ENDIF
187  mbuf = newsize
188  ENDIF
189  cbuf(mlen + 1:mlen + lindex) = cindex(1:lindex)
190  mlen = mlen + lindex
191  ELSE ! UNRECOGNIZED SECTION
192  iret = 5
193  RETURN
194  ENDIF
195  ibskip = ibskip + lensec
196  ENDDO
197 END SUBROUTINE ixgb2
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 g2_sbytec(out, in, iskip, nbits)
Put arbitrary size values into a packed bit string, taking the low order bits from each value in the ...
Definition: g2_gbytesc.F90:42
subroutine ixgb2(LUGB, LSKIP, LGRIB, CBUF, NUMFLD, MLEN, IRET)
Generate an index record for each field in a GRIB2 message.
Definition: ixgb2.F90:56
Reallocate memory, preserving contents.
Definition: realloc.f:12