NCEPLIBS-g2  3.4.5
getlocal.f
Go to the documentation of this file.
1 
6 
34 
35  subroutine getlocal(cgrib,lcgrib,localnum,csec2,lcsec2,ierr)
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,ibeg,istart,numlocal
46 
47  ierr=0
48  numlocal=0
49 !
50 ! Check for valid request number
51 !
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 !
60  istart=0
61  do j=1,100
62  ctemp=cgrib(j)//cgrib(j+1)//cgrib(j+2)//cgrib(j+3)
63  if (ctemp.eq.grib ) then
64  istart=j
65  exit
66  endif
67  enddo
68  if (istart.eq.0) then
69  print *,'getlocal: Beginning characters GRIB not found.'
70  ierr=1
71  return
72  endif
73 !
74 ! Unpack Section 0 - Indicator Section
75 !
76  iofst=8*(istart+5)
77  call g2_gbytec(cgrib,listsec0(1),iofst,8) ! Discipline
78  iofst=iofst+8
79  call g2_gbytec(cgrib,listsec0(2),iofst,8) ! GRIB edition number
80  iofst=iofst+8
81  iofst=iofst+32
82  call g2_gbytec(cgrib,lengrib,iofst,32) ! Length of GRIB message
83  iofst=iofst+32
84  lensec0=16
85  ipos=istart+lensec0
86 !
87 ! Currently handles only GRIB Edition 2.
88 !
89  if (listsec0(2).ne.2) then
90  print *,'getlocal: can only decode GRIB edition 2.'
91  ierr=2
92  return
93  endif
94 !
95 ! Loop through the remaining sections keeping track of the
96 ! length of each. Also check to see that if the current occurrence
97 ! of Section 2 is the same as the one requested.
98 !
99  do
100  ! Check to see if we are at end of GRIB message
101  ctemp=cgrib(ipos)//cgrib(ipos+1)//cgrib(ipos+2)//cgrib(ipos+3)
102  if (ctemp.eq.c7777 ) then
103  ipos=ipos+4
104  ! If end of GRIB message not where expected, issue error
105  if (ipos.ne.(istart+lengrib)) then
106  print *,'getlocal: "7777" found, but not where expected.'
107  ierr=4
108  return
109  endif
110  exit
111  endif
112  ! Get length of Section and Section number
113  iofst=(ipos-1)*8
114  call g2_gbytec(cgrib,lensec,iofst,32) ! Get Length of Section
115  iofst=iofst+32
116  call g2_gbytec(cgrib,isecnum,iofst,8) ! Get Section number
117  iofst=iofst+8
118  ! If found the requested occurrence of Section 2,
119  ! return the section contents.
120  if (isecnum.eq.2) then
121  numlocal=numlocal+1
122  if (numlocal.eq.localnum) then
123  lcsec2=lensec-5
124  csec2(1:lcsec2)=cgrib(ipos+5:ipos+lensec-1)
125  return
126  endif
127  endif
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 
137  enddo
138 
139 !
140 ! If exited from above loop, the end of the GRIB message was reached
141 ! before the requested occurrence of section 2 was found.
142 !
143  print *,'getlocal: GRIB message contained ',numlocal,
144  & ' local sections.'
145  print *,'getlocal: The request was for the ',localnum,
146  & ' occurrence.'
147  ierr=6
148 
149  return
150  end
151 
152 
153 
154 
155 
156 
157 
getlocal
subroutine getlocal(cgrib, lcgrib, localnum, csec2, lcsec2, ierr)
This subroutine returns the contents of Section 2 from a GRIB2 message.
Definition: getlocal.f:36
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