Update version info for release v4.6.1 (#2122)
[WRF.git] / external / io_grib2 / g2lib / gb_info.F
bloba15d632ef1a9d8c25537f1587be60c5b87106677
1       subroutine gb_info(cgrib,lcgrib,listsec0,listsec1,
2      &                    numfields,numlocal,maxlocal,ierr)
3 !$$$  SUBPROGRAM DOCUMENTATION BLOCK
4 !                .      .    .                                       .
5 ! SUBPROGRAM:    gb_info 
6 !   PRGMMR: Gilbert         ORG: W/NP11    DATE: 2000-05-25
8 ! ABSTRACT: This subroutine searches through a GRIB2 message and
9 !   returns the number of gridded fields found in the message and
10 !   the number (and maximum size) of Local Use Sections.
11 !   Also various checks  are performed
12 !   to see if the message is a valid GRIB2 message.
14 ! PROGRAM HISTORY LOG:
15 ! 2000-05-25  Gilbert
17 ! USAGE:    CALL gb_info(cgrib,lcgrib,listsec0,listsec1,
18 !     &                    numfields,numlocal,maxlocal,ierr)
19 !   INPUT ARGUMENT LIST:
20 !     cgrib    - Character array that contains the GRIB2 message
21 !     lcgrib   - Length (in bytes) of GRIB message in array cgrib.
23 !   OUTPUT ARGUMENT LIST:      
24 !     listsec0 - Contains information decoded from GRIB Indicator Section 0.
25 !                Must be dimensioned >= 2.
26 !                listsec0(1)=Discipline-GRIB Master Table Number
27 !                            (see Code Table 0.0)
28 !                listsec0(2)=GRIB Edition Number (currently 2)
29 !                listsec0(3)=Length of GRIB message
30 !     listsec1 - Contains information read from GRIB Identification Section 1.
31 !                Must be dimensioned >= 13.
32 !                listsec1(1)=Id of orginating centre (Common Code Table C-1)
33 !                listsec1(2)=Id of orginating sub-centre (local table)
34 !                listsec1(3)=GRIB Master Tables Version Number (Code Table 1.0)
35 !                listsec1(4)=GRIB Local Tables Version Number 
36 !                listsec1(5)=Significance of Reference Time (Code Table 1.1)
37 !                listsec1(6)=Reference Time - Year (4 digits)
38 !                listsec1(7)=Reference Time - Month
39 !                listsec1(8)=Reference Time - Day
40 !                listsec1(9)=Reference Time - Hour
41 !                listsec1(10)=Reference Time - Minute
42 !                listsec1(11)=Reference Time - Second
43 !                listsec1(12)=Production status of data (Code Table 1.2)
44 !                listsec1(13)=Type of processed data (Code Table 1.3)
45 !     numfields- The number of gridded fieldse found in the GRIB message.
46 !     numlocal - The number of Local Use Sections ( Section 2 ) found in 
47 !                the GRIB message.
48 !     maxlocal-  The size of the largest Local Use Section ( Section 2 ).
49 !                Can be used to ensure that the return array passed
50 !                to subroutine getlocal is dimensioned large enough.
51 !     ierr     - Error return code.
52 !                0 = no error
53 !                1 = Beginning characters "GRIB" not found.
54 !                2 = GRIB message is not Edition 2.
55 !                3 = Could not find Section 1, where expected.
56 !                4 = End string "7777" found, but not where expected.
57 !                5 = End string "7777" not found at end of message.
58 !                6 = Invalid section number found.
60 ! REMARKS: None
62 ! ATTRIBUTES:
63 !   LANGUAGE: Fortran 90
64 !   MACHINE:  IBM SP
66 !$$$
68       character(len=1),intent(in) :: cgrib(lcgrib)
69       integer,intent(in) :: lcgrib
70       integer,intent(out) :: listsec0(3),listsec1(13)
71       integer,intent(out) :: numlocal,numfields,maxlocal,ierr
72       
73       character(len=4),parameter :: grib='GRIB',c7777='7777'
74       character(len=4) :: ctemp
75       integer,parameter :: zero=0,one=1
76       integer,parameter :: mapsec1len=13
77       integer,parameter :: 
78      &        mapsec1(mapsec1len)=(/ 2,2,1,1,1,2,1,1,1,1,1,1,1 /)
79       integer iofst,ibeg,istart
81       ierr=0
82       numlocal=0
83       numfields=0
84       maxlocal=0
86 !  Check for beginning of GRIB message in the first 100 bytes
88       istart=0
89       do j=1,100
90         ctemp=cgrib(j)//cgrib(j+1)//cgrib(j+2)//cgrib(j+3)
91         if (ctemp.eq.grib ) then
92           istart=j
93           exit
94         endif
95       enddo
96       if (istart.eq.0) then
97         print *,'gb_info:  Beginning characters GRIB not found.'
98         ierr=1
99         return
100       endif
102 !  Unpack Section 0 - Indicator Section 
104       iofst=8*(istart+5)
105       call g2lib_gbyte(cgrib,listsec0(1),iofst,8)     ! Discipline
106       iofst=iofst+8
107       call g2lib_gbyte(cgrib,listsec0(2),iofst,8)     ! GRIB edition number
108       iofst=iofst+8
109       iofst=iofst+32
110       call g2lib_gbyte(cgrib,lengrib,iofst,32)        ! Length of GRIB message
111       iofst=iofst+32
112       listsec0(3)=lengrib
113       lensec0=16
114       ipos=istart+lensec0
116 !  Currently handles only GRIB Edition 2.
117 !  
118       if (listsec0(2).ne.2) then
119         print *,'gb_info: can only decode GRIB edition 2.'
120         ierr=2
121         return
122       endif
124 !  Unpack Section 1 - Identification Section
126       call g2lib_gbyte(cgrib,lensec1,iofst,32)        ! Length of Section 1
127       iofst=iofst+32
128       call g2lib_gbyte(cgrib,isecnum,iofst,8)         ! Section number ( 1 )
129       iofst=iofst+8
130       if (isecnum.ne.1) then
131         print *,'gb_info: Could not find section 1.'
132         ierr=3
133         return
134       endif
135       !
136       !   Unpack each input value in array listsec1 into the
137       !   the appropriate number of octets, which are specified in
138       !   corresponding entries in array mapsec1.
139       !
140       do i=1,mapsec1len
141         nbits=mapsec1(i)*8
142         call g2lib_gbyte(cgrib,listsec1(i),iofst,nbits)
143         iofst=iofst+nbits
144       enddo
145       ipos=ipos+lensec1
147 !  Loop through the remaining sections to see if they are valid.
148 !  Also count the number of times Section 2
149 !  and Section 4 appear.
151       do
152         ctemp=cgrib(ipos)//cgrib(ipos+1)//cgrib(ipos+2)//cgrib(ipos+3)
153         if (ctemp.eq.c7777 ) then
154           ipos=ipos+4
155           if (ipos.ne.(istart+lengrib)) then
156             print *,'gb_info: "7777" found, but not where expected.'
157             ierr=4
158             return
159           endif
160           exit
161         endif
162         iofst=(ipos-1)*8
163         call g2lib_gbyte(cgrib,lensec,iofst,32)        ! Get Length of Section
164         iofst=iofst+32
165         call g2lib_gbyte(cgrib,isecnum,iofst,8)         ! Get Section number
166         iofst=iofst+8
167         ipos=ipos+lensec                 ! Update beginning of section pointer
168         if (ipos.gt.(istart+lengrib)) then
169           print *,'gb_info: "7777"  not found at end of GRIB message.'
170           ierr=5
171           return
172         endif
173         if ( isecnum.ge.2.AND.isecnum.le.7 ) then
174            if (isecnum.eq.2) then     ! Local Section 2
175               !   increment counter for total number of local sections found
176               numlocal=numlocal+1
177               lenposs=lensec-5
178               if ( lenposs.gt.maxlocal ) maxlocal=lenposs
179            elseif (isecnum.eq.4) then
180               !   increment counter for total number of fields found
181               numfields=numfields+1
182            endif
183         else
184            print *,'gb_info: Invalid section number found in GRIB',
185      &             ' message: ',isecnum
186            ierr=6
187            return
188         endif
189         
190       enddo
192       return
193       end