updated top-level README and version_decl for V4.4.2 (#1795)
[WRF.git] / external / io_grib2 / g2lib / getlocal.F
blobc101a5ea5af28f9a8394f521fd18234763e16895
1       subroutine getlocal(cgrib,lcgrib,localnum,csec2,lcsec2,ierr)
2 !$$$  SUBPROGRAM DOCUMENTATION BLOCK
3 !                .      .    .                                       .
4 ! SUBPROGRAM:    getlocal 
5 !   PRGMMR: Gilbert         ORG: W/NP11    DATE: 2000-05-25
7 ! ABSTRACT: This subroutine returns the contents of Section 2 ( Local 
8 !   Use Section ) from a GRIB2 message.  Since there can be multiple
9 !   occurrences of Section 2 within a GRIB message, the calling routine
10 !   indicates which occurrence is being requested with the localnum argument.
12 ! PROGRAM HISTORY LOG:
13 ! 2000-05-25  Gilbert
15 ! USAGE:    CALL getlocal(cgrib,lcgrib,localnum,csec2,lcsec2,ierr)
16 !   INPUT ARGUMENT LIST:
17 !     cgrib    - Character array that contains the GRIB2 message
18 !     lcgrib   - Length (in bytes) of GRIB message in array cgrib.
19 !     localnum - The nth occurrence of Section 2 requested.
21 !   OUTPUT ARGUMENT LIST:      
22 !     csec2    - Character array containing information read from 
23 !                Section 2.
24 !                The dimension of this array can be obtained in advance
25 !                from argument maxlocal, which is returned from subroutine 
26 !                gb_info.
27 !     lcsec2   - Number of bytes of character array csec2 read from
28 !                Section 2.
29 !     ierr     - Error return code.
30 !                0 = no error
31 !                1 = Beginning characters "GRIB" not found.
32 !                2 = GRIB message is not Edition 2.
33 !                3 = The section 2 request number was not positive.
34 !                4 = End string "7777" found, but not where expected.
35 !                5 = End string "7777" not found at end of message.
36 !                6 = GRIB message did not contain the requested number of
37 !                    Local Use Sections.
39 ! REMARKS: Note that subroutine gb_info can be used to first determine
40 !          how many Local Use sections exist in a given GRIB message.
42 ! ATTRIBUTES:
43 !   LANGUAGE: Fortran 90
44 !   MACHINE:  IBM SP
46 !$$$
48       character(len=1),intent(in) :: cgrib(lcgrib)
49       integer,intent(in) :: lcgrib,localnum
50       character(len=1),intent(out) :: csec2(*)
51       integer,intent(out) :: lcsec2,ierr
52       
53       character(len=4),parameter :: grib='GRIB',c7777='7777'
54       character(len=4) :: ctemp
55       integer :: listsec0(2)
56       integer iofst,ibeg,istart,numlocal
58       ierr=0
59       numlocal=0
61 !  Check for valid request number
62 !  
63       if (localnum.le.0) then
64         print *,'getlocal: Request for local section must be positive.'
65         ierr=3
66         return
67       endif
69 !  Check for beginning of GRIB message in the first 100 bytes
71       istart=0
72       do j=1,100
73         ctemp=cgrib(j)//cgrib(j+1)//cgrib(j+2)//cgrib(j+3)
74         if (ctemp.eq.grib ) then
75           istart=j
76           exit
77         endif
78       enddo
79       if (istart.eq.0) then
80         print *,'getlocal:  Beginning characters GRIB not found.'
81         ierr=1
82         return
83       endif
85 !  Unpack Section 0 - Indicator Section 
87       iofst=8*(istart+5)
88       call g2lib_gbyte(cgrib,listsec0(1),iofst,8)     ! Discipline
89       iofst=iofst+8
90       call g2lib_gbyte(cgrib,listsec0(2),iofst,8)     ! GRIB edition number
91       iofst=iofst+8
92       iofst=iofst+32
93       call g2lib_gbyte(cgrib,lengrib,iofst,32)        ! Length of GRIB message
94       iofst=iofst+32
95       lensec0=16
96       ipos=istart+lensec0
98 !  Currently handles only GRIB Edition 2.
99 !  
100       if (listsec0(2).ne.2) then
101         print *,'getlocal: can only decode GRIB edition 2.'
102         ierr=2
103         return
104       endif
106 !  Loop through the remaining sections keeping track of the 
107 !  length of each.  Also check to see that if the current occurrence
108 !  of Section 2 is the same as the one requested.
110       do
111         !    Check to see if we are at end of GRIB message
112         ctemp=cgrib(ipos)//cgrib(ipos+1)//cgrib(ipos+2)//cgrib(ipos+3)
113         if (ctemp.eq.c7777 ) then
114           ipos=ipos+4
115           !    If end of GRIB message not where expected, issue error
116           if (ipos.ne.(istart+lengrib)) then
117             print *,'getlocal: "7777" found, but not where expected.'
118             ierr=4
119             return
120           endif
121           exit
122         endif
123         !     Get length of Section and Section number
124         iofst=(ipos-1)*8
125         call g2lib_gbyte(cgrib,lensec,iofst,32)        ! Get Length of Section
126         iofst=iofst+32
127         call g2lib_gbyte(cgrib,isecnum,iofst,8)         ! Get Section number
128         iofst=iofst+8
129         !   If found the requested occurrence of Section 2,
130         !   return the section contents.
131         if (isecnum.eq.2) then
132           numlocal=numlocal+1
133           if (numlocal.eq.localnum) then
134             lcsec2=lensec-5
135             csec2(1:lcsec2)=cgrib(ipos+5:ipos+lensec-1)
136             return
137           endif
138         endif
139         !   Check to see if we read pass the end of the GRIB
140         !   message and missed the terminator string '7777'.
141         ipos=ipos+lensec                 ! Update beginning of section pointer
142         if (ipos.gt.(istart+lengrib)) then
143           print *,'getlocal: "7777"  not found at end of GRIB message.'
144           ierr=5
145           return
146         endif
147         
148       enddo
151 !  If exited from above loop, the end of the GRIB message was reached
152 !  before the requested occurrence of section 2 was found.
154       print *,'getlocal: GRIB message contained ',numlocal,
155      &        ' local sections.'
156       print *,'getlocal: The request was for the ',localnum,
157      &        ' occurrence.'
158       ierr=6
160       return
161       end