NCEPLIBS-g2  3.4.7
getlocal.F90
Go to the documentation of this file.
1 
5 
34 subroutine getlocal(cgrib, lcgrib, localnum, csec2, lcsec2, ierr)
35  implicit none
36 
37  character(len = 1), intent(in) :: cgrib(lcgrib)
38  integer, intent(in) :: lcgrib, localnum
39  character(len = 1), intent(out) :: csec2(*)
40  integer, intent(out) :: lcsec2, ierr
41 
42  character(len = 4), parameter :: grib = 'GRIB', c7777 = '7777'
43  character(len = 4) :: ctemp
44  integer :: listsec0(2)
45  integer iofst, istart, numlocal
46  integer :: lengrib, lensec, lensec0, j, ipos, isecnum
47 
48  ierr = 0
49  numlocal = 0
50 
51  ! Check for valid request number.
52  if (localnum .le. 0) then
53  print *, 'getlocal: Request for local section must be positive.'
54  ierr = 3
55  return
56  endif
57 
58  ! Check for beginning of GRIB message in the first 100 bytes
59  istart = 0
60  do j = 1, 100
61  ctemp = cgrib(j)//cgrib(j+1)//cgrib(j+2)//cgrib(j+3)
62  if (ctemp .eq. grib) then
63  istart = j
64  exit
65  endif
66  enddo
67  if (istart .eq. 0) then
68  print *, 'getlocal: Beginning characters GRIB not found.'
69  ierr = 1
70  return
71  endif
72 
73  ! Unpack Section 0 - Indicator Section
74  iofst = 8 * (istart + 5)
75  call g2_gbytec(cgrib, listsec0(1), iofst, 8) ! Discipline
76  iofst = iofst + 8
77  call g2_gbytec(cgrib, listsec0(2), iofst, 8) ! GRIB edition number
78  iofst = iofst + 8
79  iofst = iofst + 32
80  call g2_gbytec(cgrib, lengrib, iofst, 32) ! Length of GRIB message
81  iofst = iofst + 32
82  lensec0 = 16
83  ipos = istart + lensec0
84 
85  ! Currently handles only GRIB Edition 2.
86  if (listsec0(2) .ne. 2) then
87  print *, 'getlocal: can only decode GRIB edition 2.'
88  ierr = 2
89  return
90  endif
91 
92  ! Loop through the remaining sections keeping track of the length of
93  ! each. Also check to see that if the current occurrence of Section
94  ! 2 is the same as the one requested.
95  do
96  ! Check to see if we are at end of GRIB message
97  ctemp = cgrib(ipos) // cgrib(ipos + 1) // cgrib(ipos + 2) // cgrib(ipos + 3)
98  if (ctemp .eq. c7777) then
99  ipos = ipos + 4
100 
101  ! If end of GRIB message not where expected, issue error
102  if (ipos .ne. (istart + lengrib)) then
103  print *, 'getlocal: "7777" found, but not where expected.'
104  ierr = 4
105  return
106  endif
107  exit
108  endif
109 
110  ! Get length of Section and Section number
111  iofst = (ipos - 1) * 8
112  call g2_gbytec(cgrib, lensec, iofst, 32) ! Get Length of Section
113  iofst = iofst + 32
114  call g2_gbytec(cgrib, isecnum, iofst, 8) ! Get Section number
115  iofst = iofst + 8
116 
117  ! If found the requested occurrence of Section 2,
118  ! return the section contents.
119  if (isecnum .eq. 2) then
120  numlocal = numlocal + 1
121  if (numlocal.eq.localnum) then
122  lcsec2 = lensec - 5
123  csec2(1:lcsec2) = cgrib(ipos + 5:ipos + lensec - 1)
124  return
125  endif
126  endif
127 
128  ! Check to see if we read pass the end of the GRIB
129  ! message and missed the terminator string '7777'.
130  ipos = ipos + lensec ! Update beginning of section pointer
131  if (ipos .gt. (istart + lengrib)) then
132  print *, 'getlocal: "7777" not found at end of GRIB message.'
133  ierr = 5
134  return
135  endif
136  enddo
137 
138  ! If exited from above loop, the end of the GRIB message was reached
139  ! before the requested occurrence of section 2 was found.
140  print *, 'getlocal: GRIB message contained ', numlocal, ' local sections.'
141  print *, 'getlocal: The request was for the ', localnum, ' occurrence.'
142  ierr = 6
143 
144 end subroutine getlocal
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 getlocal(cgrib, lcgrib, localnum, csec2, lcsec2, ierr)
This subroutine returns the contents of Section 2 from a GRIB2 message.
Definition: getlocal.F90:35