NCEPLIBS-g2  3.4.7
gribinfo.F90
Go to the documentation of this file.
1 
5 
78 subroutine gribinfo(cgrib, lcgrib, listsec0, listsec1, &
79  numlocal, numfields, maxvals, ierr)
80  implicit none
81 
82  character(len = 1), intent(in) :: cgrib(lcgrib)
83  integer, intent(in) :: lcgrib
84  integer, intent(out) :: listsec0(3), listsec1(13), maxvals(7)
85  integer, intent(out) :: numlocal, numfields, ierr
86 
87  integer :: i, ipos, isecnum, j, lengrib, lenposs, lensec, lensec0, lensec1
88  integer :: maxcoordlist, maxdeflist, maxdrstmpl, maxgridpts, maxpdstmpl, maxgdstmpl
89  integer :: maxsec2len, nbits, nbyte, ngdpts, numcoord
90  character(len = 4), parameter :: grib = 'GRIB', c7777 = '7777'
91  character(len = 4) :: ctemp
92  integer, parameter :: zero = 0, one = 1
93  integer, parameter :: mapsec1len = 13
94  integer, parameter :: mapsec1(mapsec1len) = (/ 2, 2, 1, 1, 1, 2, 1, 1, &
95  1, 1, 1, 1, 1 /)
96  integer iofst, istart
97 
98  ierr = 0
99  numlocal = 0
100  numfields = 0
101  maxsec2len = 1
102  maxgdstmpl = 1
103  maxdeflist = 1
104  maxpdstmpl = 1
105  maxcoordlist = 1
106  maxdrstmpl = 1
107  maxgridpts = 0
108 
109  ! Check for beginning of GRIB message in the first 100 bytes
110  istart = 0
111  do j = 1, 100
112  ctemp = cgrib(j) // cgrib(j + 1) // cgrib(j + 2) // cgrib(j + 3)
113  if (ctemp .eq. grib) then
114  istart = j
115  exit
116  endif
117  enddo
118  if (istart .eq. 0) then
119  print *, 'gribinfo: Beginning characters GRIB not found.'
120  ierr = 1
121  return
122  endif
123 
124  ! Unpack Section 0 - Indicator Section.
125  iofst = 8 * (istart + 5)
126  call g2_gbytec(cgrib, listsec0(1), iofst, 8) ! Discipline
127  iofst = iofst + 8
128  call g2_gbytec(cgrib, listsec0(2), iofst, 8) ! GRIB edition number
129  iofst = iofst + 8
130  iofst = iofst + 32
131  call g2_gbytec(cgrib, lengrib, iofst, 32) ! Length of GRIB message
132  iofst = iofst + 32
133  listsec0(3) = lengrib
134  lensec0 = 16
135  ipos = istart + lensec0
136 
137  ! Currently handles only GRIB Edition 2.
138  if (listsec0(2) .ne. 2) then
139  print *, 'gribinfo: can only decode GRIB edition 2.'
140  ierr = 2
141  return
142  endif
143 
144  ! Unpack Section 1 - Identification Section.
145  call g2_gbytec(cgrib, lensec1, iofst, 32) ! Length of Section 1
146  iofst = iofst + 32
147  call g2_gbytec(cgrib, isecnum, iofst, 8) ! Section number (1)
148  iofst = iofst + 8
149  if (isecnum .ne. 1) then
150  print *, 'gribinfo: Could not find section 1.'
151  ierr = 3
152  return
153  endif
154 
155  ! Unpack each input value in array listsec1 into the the appropriate
156  ! number of octets, which are specified in corresponding entries in
157  ! array mapsec1.
158  do i = 1, mapsec1len
159  nbits = mapsec1(i) * 8
160  call g2_gbytec(cgrib, listsec1(i), iofst, nbits)
161  iofst = iofst + nbits
162  enddo
163  ipos = ipos + lensec1
164 
165  ! Loop through the remaining sections keeping track of the length of
166  ! each. Also count the number of times Section 2 and Section 4
167  ! appear.
168  do
169  ctemp = cgrib(ipos) // cgrib(ipos + 1) // cgrib(ipos + 2) // &
170  cgrib(ipos + 3)
171  if (ctemp .eq. c7777) then
172  ipos = ipos + 4
173  if (ipos .ne. (istart + lengrib)) then
174  print *, 'gribinfo: "7777" found, but not where expected.'
175  ierr = 4
176  return
177  endif
178  exit
179  endif
180  iofst = (ipos - 1) * 8
181  call g2_gbytec(cgrib, lensec, iofst, 32) ! Get Length of Section
182  iofst = iofst + 32
183  call g2_gbytec(cgrib, isecnum, iofst, 8) ! Get Section number
184  iofst = iofst + 8
185  ipos = ipos + lensec ! Update beginning of section pointer
186  if (ipos .gt. (istart + lengrib)) then
187  print *, 'gribinfo: "7777" not found at end of GRIB message.'
188  ierr = 5
189  return
190  endif
191  if (isecnum .eq. 2) then ! Local Section 2
192  ! Increment counter for total number of local sections found
193  ! and determine largest Section 2 in message.
194  numlocal = numlocal + 1
195  lenposs = lensec-5
196  if (lenposs .gt. maxsec2len) maxsec2len = lenposs
197  elseif (isecnum .eq. 3) then
198  iofst = iofst + 8 ! skip source of grid def.
199  call g2_gbytec(cgrib, ngdpts, iofst, 32) ! Get Num of Grid Points
200  iofst = iofst + 32
201  call g2_gbytec(cgrib, nbyte, iofst, 8) ! Get Num octets for opt. list
202  iofst = iofst + 8
203  if (ngdpts .gt. maxgridpts) maxgridpts = ngdpts
204  lenposs = lensec - 14
205  if (lenposs .gt. maxgdstmpl) maxgdstmpl = lenposs
206  if (nbyte .ne. 0) then
207  lenposs = lenposs / nbyte
208  if (lenposs .gt. maxdeflist) maxdeflist = lenposs
209  endif
210  elseif (isecnum .eq. 4) then
211  numfields = numfields + 1
212  call g2_gbytec(cgrib, numcoord, iofst, 16) ! Get Num of Coord Values
213  iofst = iofst + 16
214  if (numcoord .ne. 0) then
215  if (numcoord .gt. maxcoordlist) maxcoordlist = numcoord
216  endif
217  lenposs = lensec - 9
218  if (lenposs.gt.maxpdstmpl) maxpdstmpl = lenposs
219  elseif (isecnum .eq. 5) then
220  lenposs = lensec-11
221  if (lenposs .gt. maxdrstmpl) maxdrstmpl = lenposs
222  endif
223  enddo
224 
225  maxvals(1) = maxsec2len
226  maxvals(2) = maxgdstmpl
227  maxvals(3) = maxdeflist
228  maxvals(4) = maxpdstmpl
229  maxvals(5) = maxcoordlist
230  maxvals(6) = maxdrstmpl
231  maxvals(7) = maxgridpts
232 end subroutine gribinfo
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 gribinfo(cgrib, lcgrib, listsec0, listsec1, numlocal, numfields, maxvals, ierr)
Find the number of Local Use Sections and gridded fields in a GRIB2 message, and the maximum sizes of...
Definition: gribinfo.F90:80