NCEPLIBS-g2  3.4.8
gb_info.F90
Go to the documentation of this file.
1 
5 
53 subroutine gb_info(cgrib, lcgrib, listsec0, listsec1, &
54  numfields, numlocal, maxlocal, ierr)
55  implicit none
56 
57  character(len = 1), intent(in) :: cgrib(lcgrib)
58  integer, intent(in) :: lcgrib
59  integer, intent(out) :: listsec0(3), listsec1(13)
60  integer, intent(out) :: numlocal, numfields, maxlocal, ierr
61 
62  character(len = 4), parameter :: grib = 'GRIB', c7777 = '7777'
63  character(len = 4) :: ctemp
64  integer, parameter :: zero = 0, one = 1
65  integer, parameter :: mapsec1len = 13
66  integer, parameter :: mapsec1(mapsec1len) = (/ 2, 2, 1, 1, 1, 2, 1, 1, 1, 1, 1, 1, 1 /)
67  integer :: iofst, istart
68  integer :: nbits, lensec1, lensec0, lensec, lenposs, lengrib, j
69  integer :: i, ipos, isecnum
70 
71  ierr = 0
72  numlocal = 0
73  numfields = 0
74  maxlocal = 0
75 
76  ! Check for beginning of GRIB message in the first 100 bytes.
77  istart = 0
78  do j = 1, 100
79  ctemp = cgrib(j) // cgrib(j + 1) // cgrib(j + 2) // cgrib(j+3)
80  if (ctemp .eq. grib ) then
81  istart = j
82  exit
83  endif
84  enddo
85  if (istart .eq. 0) then
86  print *, 'gb_info: Beginning characters GRIB not found.'
87  ierr = 1
88  return
89  endif
90 
91  ! Unpack Section 0 - Indicator Section.
92  iofst = 8 * (istart + 5)
93  call g2_gbytec(cgrib, listsec0(1), iofst, 8) ! Discipline
94  iofst = iofst + 8
95  call g2_gbytec(cgrib, listsec0(2), iofst, 8) ! GRIB edition number
96  iofst = iofst+8
97  iofst = iofst + 32
98  call g2_gbytec(cgrib, lengrib, iofst, 32) ! Length of GRIB message
99  iofst = iofst + 32
100  listsec0(3) = lengrib
101  lensec0 = 16
102  ipos = istart + lensec0
103 
104  ! Currently handles only GRIB Edition 2.
105  if (listsec0(2) .ne. 2) then
106  print *, 'gb_info: can only decode GRIB edition 2.'
107  ierr = 2
108  return
109  endif
110 
111  ! Unpack Section 1 - Identification Section.
112  call g2_gbytec(cgrib, lensec1, iofst, 32) ! Length of Section 1
113  iofst = iofst + 32
114  call g2_gbytec(cgrib, isecnum, iofst, 8) ! Section number ( 1 )
115  iofst = iofst + 8
116  if (isecnum .ne. 1) then
117  print *, 'gb_info: Could not find section 1.'
118  ierr = 3
119  return
120  endif
121 
122  ! Unpack each input value in array listsec1 into the the appropriate
123  ! number of octets, which are specified in corresponding entries in
124  ! array mapsec1.
125  do i = 1, mapsec1len
126  nbits = mapsec1(i) * 8
127  call g2_gbytec(cgrib, listsec1(i), iofst, nbits)
128  iofst = iofst + nbits
129  enddo
130  ipos = ipos + lensec1
131 
132  ! Loop through the remaining sections to see if they are valid.
133  ! Also count the number of times Section 2 and Section 4 appear.
134  do
135  ctemp = cgrib(ipos) // cgrib(ipos + 1) // cgrib(ipos + 2) // cgrib(ipos + 3)
136  if (ctemp .eq. c7777 ) then
137  ipos = ipos + 4
138  if (ipos .ne. (istart + lengrib)) then
139  print *, 'gb_info: "7777" found, but not where expected.'
140  ierr = 4
141  return
142  endif
143  exit
144  endif
145  iofst = (ipos - 1) * 8
146  call g2_gbytec(cgrib, lensec, iofst, 32) ! Get Length of Section
147  iofst = iofst + 32
148  call g2_gbytec(cgrib, isecnum, iofst, 8) ! Get Section number
149  iofst = iofst + 8
150  ipos = ipos+lensec ! Update beginning of section pointer
151  if (ipos .gt. (istart + lengrib)) then
152  print *, 'gb_info: "7777" not found at end of GRIB message.'
153  ierr = 5
154  return
155  endif
156  if (isecnum .ge. 2 .AND. isecnum .le. 7) then
157  if (isecnum .eq. 2) then ! Local Section 2
158  ! Increment counter for total number of local sections found.
159  numlocal = numlocal + 1
160  lenposs = lensec - 5
161  if (lenposs .gt. maxlocal) maxlocal = lenposs
162  elseif (isecnum .eq. 4) then
163  ! Increment counter for total number of fields found.
164  numfields = numfields + 1
165  endif
166  else
167  print *, 'gb_info: Invalid section number found in GRIB message: ', isecnum
168  ierr = 6
169  return
170  endif
171  enddo
172 
173 end subroutine gb_info
174 
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 gb_info(cgrib, lcgrib, listsec0, listsec1, numfields, numlocal, maxlocal, ierr)
Find the number of gridded fields and Local Use Sections in a GRIB2 message.
Definition: gb_info.F90:55