NCEPLIBS-g2  3.4.5
ixgb2.f
Go to the documentation of this file.
1 C> @file
2 C> @brief This subroutine generates an index record for each field in
3 C> a grib2 message.
4 C> @author Mark Iredell @date 1995-10-31
5 C>
6 
7 C> This subroutine generates an index record for each field in
8 C> a grib2 message. The index records are written to index buffer
9 C> pointed to by cbuf.
10 C> - byte 001 - 004 length of index record
11 C> - byte 005 - 008 bytes to skip in data file before grib message
12 C> - byte 009 - 012 bytes to skip in message before lus (local use)
13 C> set = 0, if no local use section in grib2 message.
14 C> - byte 013 - 016 bytes to skip in message before gds
15 C> - byte 017 - 020 bytes to skip in message before pds
16 C> - byte 021 - 024 bytes to skip in message before drs
17 C> - byte 025 - 028 bytes to skip in message before bms
18 C> - byte 029 - 032 bytes to skip in message before data section
19 C> - byte 033 - 040 bytes total in the message
20 C> - byte 041 - 041 grib version number (currently 2)
21 C> - byte 042 - 042 message discipline
22 C> - byte 043 - 044 field number within grib2 message
23 C> - byte 045 - ii identification section (ids)
24 C> - byte ii+1 - jj grid definition section (gds)
25 C> - byte jj+1 - kk product definition section (pds)
26 C> - byte kk+1 - ll the data representation section (drs)
27 C> - byte ll+1 - ll+6 first 6 bytes of the bit map section (bms)
28 C>
29 C> Program history log:
30 C> - 1995-10-31 Mark Iredell
31 C> - 1996-10-31 Mark Iredell augmented optional definitions to byte 320.
32 C> - 2001-12-10 Stephen Gilbert modified from ixgb to create grib2 indexes.
33 C> - 2002-01-31 Stephen Gilbert added identification section to index record.
34 C>
35 C> @param[in] LUGB integer unit of the unblocked grib file.
36 C> @param[in] LSKIP integer number of bytes to skip before grib message.
37 C> @param[in] LGRIB integer number of bytes in grib message.
38 C> @param[out] CBUF character*1 pointer to a buffer that contains
39 C> index records users should free memory that cbuf points to
40 C> using deallocate(cbuf) when cbuf is no longer needed.
41 C> @param[out] NUMFLD integer number of index records created.
42 C> @param[out] MLEN integer total length of all index records.
43 C> @param[out] IRET integer return code
44 C> - 0 all ok
45 C> - 1 not enough memory available to hold full index buffer
46 C> - 2 i/o error in read
47 C> - 3 grib message is not edition 2
48 C> - 4 not enough memory to allocate extent to index buffer
49 C> - 5 unidentified grib section encountered
50 C>
51 C> @author Mark Iredell @date 1995-10-31
52 C>
53 
54  SUBROUTINE ixgb2(LUGB,LSKIP,LGRIB,CBUF,NUMFLD,MLEN,IRET)
55 
56  USE re_alloc ! NEEDED FOR SUBROUTINE REALLOC
57  CHARACTER(LEN=1),POINTER,DIMENSION(:) :: CBUF
58  parameter(linmax=5000,init=50000,next=10000)
59  parameter(ixskp=4,ixlus=8,ixsgd=12,ixspd=16,ixsdr=20,ixsbm=24,
60  & ixds=28,ixlen=36,ixfld=42,ixids=44)
61  parameter(mxskp=4,mxlus=4,mxsgd=4,mxspd=4,mxsdr=4,mxsbm=4,
62  & mxds=4,mxlen=4,mxfld=2,mxbms=6)
63  CHARACTER CBREAD(LINMAX),CINDEX(LINMAX)
64  CHARACTER CVER,CDISC
65  CHARACTER CIDS(LINMAX),CGDS(LINMAX),CBMS(6)
66  CHARACTER(LEN=4) :: CTEMP
67  INTEGER LOCLUS,LOCGDS,LENGDS,LOCBMS
68 C - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69  loclus=0
70  iret=0
71  mlen=0
72  numfld=0
73  IF (ASSOCIATED(cbuf)) NULLIFY(cbuf)
74  mbuf=init
75  ALLOCATE(cbuf(mbuf),stat=istat) ! ALLOCATE INITIAL SPACE FOR CBUF
76  IF (istat.NE.0) THEN
77  iret=1
78  RETURN
79  ENDIF
80 C - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81 C READ SECTIONS 0 AND 1 FOR VERSIN NUMBER AND DISCIPLINE
82  ibread=min(lgrib,linmax)
83  CALL baread(lugb,lskip,ibread,lbread,cbread)
84  IF(lbread.NE.ibread) THEN
85  iret=2
86  RETURN
87  ENDIF
88  IF(cbread(8).NE.char(2)) THEN ! NOT GRIB EDITION 2
89  iret=3
90  RETURN
91  ENDIF
92  cver=cbread(8)
93  cdisc=cbread(7)
94  CALL g2_gbytec(cbread,lensec1,16*8,4*8)
95  lensec1=min(lensec1,ibread)
96  cids(1:lensec1)=cbread(17:16+lensec1)
97  ibskip=lskip+16+lensec1
98 C - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
99 C LOOP THROUGH REMAINING SECTIONS CREATING AN INDEX FOR EACH FIELD
100  ibread=max(5,mxbms)
101  DO
102  CALL baread(lugb,ibskip,ibread,lbread,cbread)
103  ctemp=cbread(1)//cbread(2)//cbread(3)//cbread(4)
104  IF (ctemp.EQ.'7777') RETURN ! END OF MESSAGE FOUND
105  IF(lbread.NE.ibread) THEN
106  iret=2
107  RETURN
108  ENDIF
109  CALL g2_gbytec(cbread,lensec,0*8,4*8)
110  CALL g2_gbytec(cbread,numsec,4*8,1*8)
111 
112  IF (numsec.EQ.2) THEN ! SAVE LOCAL USE LOCATION
113  loclus=ibskip-lskip
114  ELSEIF (numsec.EQ.3) THEN ! SAVE GDS INFO
115  lengds=lensec
116  cgds=char(0)
117  CALL baread(lugb,ibskip,lengds,lbread,cgds)
118  IF(lbread.NE.lengds) THEN
119  iret=2
120  RETURN
121  ENDIF
122  locgds=ibskip-lskip
123  ELSEIF (numsec.EQ.4) THEN ! FOUND PDS
124  cindex=char(0)
125  CALL g2_sbytec(cindex,lskip,8*ixskp,8*mxskp) ! BYTES TO SKIP
126  CALL g2_sbytec(cindex,loclus,8*ixlus,8*mxlus) ! LOCATION OF LOCAL USE
127  CALL g2_sbytec(cindex,locgds,8*ixsgd,8*mxsgd) ! LOCATION OF GDS
128  CALL g2_sbytec(cindex,ibskip-lskip,8*ixspd,8*mxspd) ! LOCATION OF PDS
129  CALL g2_sbytec(cindex,lgrib,8*ixlen,8*mxlen) ! LEN OF GRIB2
130  cindex(41)=cver
131  cindex(42)=cdisc
132  CALL g2_sbytec(cindex,numfld+1,8*ixfld,8*mxfld) ! FIELD NUM
133  cindex(ixids+1:ixids+lensec1)=cids(1:lensec1)
134  lindex=ixids+lensec1
135  cindex(lindex+1:lindex+lengds)=cgds(1:lengds)
136  lindex=lindex+lengds
137  ilnpds=lensec
138  CALL baread(lugb,ibskip,ilnpds,lbread,cindex(lindex+1))
139  IF(lbread.NE.ilnpds) THEN
140  iret=2
141  RETURN
142  ENDIF
143  ! CINDEX(LINDEX+1:LINDEX+ILNPDS)=CBREAD(1:ILNPDS)
144  lindex=lindex+ilnpds
145  ELSEIF (numsec.EQ.5) THEN ! FOUND DRS
146  CALL g2_sbytec(cindex,ibskip-lskip,8*ixsdr,8*mxsdr) ! LOCATION OF DRS
147  ilndrs=lensec
148  CALL baread(lugb,ibskip,ilndrs,lbread,cindex(lindex+1))
149  IF(lbread.NE.ilndrs) THEN
150  iret=2
151  RETURN
152  ENDIF
153  ! CINDEX(LINDEX+1:LINDEX+ILNDRS)=CBREAD(1:ILNDRS)
154  lindex=lindex+ilndrs
155  ELSEIF (numsec.EQ.6) THEN ! FOUND BMS
156  indbmp=mova2i(cbread(6))
157  IF ( indbmp.LT.254 ) THEN
158  locbms=ibskip-lskip
159  CALL g2_sbytec(cindex,locbms,8*ixsbm,8*mxsbm) ! LOC. OF BMS
160  ELSEIF ( indbmp.EQ.254 ) THEN
161  CALL g2_sbytec(cindex,locbms,8*ixsbm,8*mxsbm) ! LOC. OF BMS
162  ELSEIF ( indbmp.EQ.255 ) THEN
163  CALL g2_sbytec(cindex,ibskip-lskip,8*ixsbm,8*mxsbm) ! LOC. OF BMS
164  ENDIF
165  cindex(lindex+1:lindex+mxbms)=cbread(1:mxbms)
166  lindex=lindex+mxbms
167  CALL g2_sbytec(cindex,lindex,0,8*4) ! NUM BYTES IN INDEX RECORD
168  ELSEIF (numsec.EQ.7) THEN ! FOUND DATA SECTION
169  CALL g2_sbytec(cindex,ibskip-lskip,8*ixds,8*mxds) ! LOC. OF DATA SEC.
170  numfld=numfld+1
171  IF ((lindex+mlen).GT.mbuf) THEN ! ALLOCATE MORE SPACE IF
172  ! NECESSARY
173  newsize=max(mbuf+next,mbuf+lindex)
174  CALL realloc(cbuf,mlen,newsize,istat)
175  IF ( istat .NE. 0 ) THEN
176  numfld=numfld-1
177  iret=4
178  RETURN
179  ENDIF
180  mbuf=newsize
181  ENDIF
182  cbuf(mlen+1:mlen+lindex)=cindex(1:lindex)
183  mlen=mlen+lindex
184  ELSE ! UNRECOGNIZED SECTION
185  iret=5
186  RETURN
187  ENDIF
188  ibskip=ibskip+lensec
189  ENDDO
190 
191 C - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
192  RETURN
193  END
re_alloc
This module contains three subroutines to reorganize the integer, real and character data in memory i...
Definition: realloc.f:14
re_alloc::realloc
Definition: realloc.f:16
ixgb2
subroutine ixgb2(LUGB, LSKIP, LGRIB, CBUF, NUMFLD, MLEN, IRET)
This subroutine generates an index record for each field in a grib2 message.
Definition: ixgb2.f:55
g2_sbytec
subroutine g2_sbytec(OUT, IN, ISKIP, NBYTE)
This subrountine is to put arbitrary size values into a packed bit string, taking the low order bits ...
Definition: g2_gbytesc.f:39
g2_gbytec
subroutine g2_gbytec(IN, IOUT, ISKIP, NBYTE)
This subrountine is to extract arbitrary size values from a packed bit string, right justifying each ...
Definition: g2_gbytesc.f:20