updated top-level README and version_decl for V4.4.2 (#1795)
[WRF.git] / external / io_grib2 / g2lib / gf_unpack4.F
blob6b1a0ff8403d63fac9bc76c0b2517c1f1356ff63
1       subroutine gf_unpack4(cgrib,lcgrib,iofst,ipdsnum,ipdstmpl,
2      &                      mappdslen,coordlist,numcoord,ierr)
3 !$$$  SUBPROGRAM DOCUMENTATION BLOCK
4 !                .      .    .                                       .
5 ! SUBPROGRAM:    gf_unpack4 
6 !   PRGMMR: Gilbert         ORG: W/NP11    DATE: 2000-05-26
8 ! ABSTRACT: This subroutine unpacks Section 4 (Product Definition Section)
9 !   starting at octet 6 of that Section.  
11 ! PROGRAM HISTORY LOG:
12 ! 2000-05-26  Gilbert
13 ! 2002-01-24  Gilbert  - Changed to dynamically allocate arrays
14 !                        and to pass pointers to those arrays through
15 !                        the argument list.
17 ! USAGE:    CALL gf_unpack4(cgrib,lcgrib,iofst,ipdsnum,ipdstmpl,mappdslen,
18 !    &                   coordlist,numcoord,ierr)
19 !   INPUT ARGUMENT LIST:
20 !     cgrib    - Character array that contains the GRIB2 message
21 !     lcgrib   - Length (in bytes) of GRIB message array cgrib.
22 !     iofst    - Bit offset of the beginning of Section 4.
24 !   OUTPUT ARGUMENT LIST:      
25 !     iofst    - Bit offset of the end of Section 4, returned.
26 !     ipdsnum  - Product Definition Template Number ( see Code Table 4.0)
27 !     ipdstmpl - Pointer to integer array containing the data values for 
28 !                the specified Product Definition
29 !                Template ( N=ipdsnum ).  Each element of this integer
30 !                array contains an entry (in the order specified) of Product
31 !                Defintion Template 4.N
32 !     mappdslen- Number of elements in ipdstmpl().  i.e. number of entries
33 !                in Product Defintion Template 4.N  ( N=ipdsnum ).
34 !     coordlist- Pointer to real array containing floating point values 
35 !                intended to document
36 !                the vertical discretisation associated to model data
37 !                on hybrid coordinate vertical levels.  (part of Section 4)
38 !     numcoord - number of values in array coordlist.
39 !     ierr     - Error return code.
40 !                0 = no error
41 !                5 = "GRIB" message contains an undefined Product Definition
42 !                    Template.
43 !                6 = memory allocation error
45 ! REMARKS: Uses Fortran 90 module pdstemplates and module re_alloc.
47 ! ATTRIBUTES:
48 !   LANGUAGE: Fortran 90
49 !   MACHINE:  IBM SP
51 !$$$
53       use pdstemplates
54       use re_alloc        !  needed for subroutine realloc
56       character(len=1),intent(in) :: cgrib(lcgrib)
57       integer,intent(in) :: lcgrib
58       integer,intent(inout) :: iofst
59       real,pointer,dimension(:) :: coordlist
60       integer,pointer,dimension(:) :: ipdstmpl
61       integer,intent(out) :: ipdsnum
62       integer,intent(out) :: ierr,numcoord
64       real(4),allocatable :: coordieee(:)
65       integer,allocatable :: mappds(:)
66       integer :: mappdslen
67       logical needext
69       ierr=0
70       nullify(ipdstmpl,coordlist)
72       call g2lib_gbyte(cgrib,lensec,iofst,32)        ! Get Length of Section
73       iofst=iofst+32
74       iofst=iofst+8     ! skip section number
75       allocate(mappds(lensec))
77       call g2lib_gbyte(cgrib,numcoord,iofst,16)    ! Get num of coordinate values
78       iofst=iofst+16
79       call g2lib_gbyte(cgrib,ipdsnum,iofst,16)    ! Get Prod. Def Template num.
80       iofst=iofst+16
81       !   Get Product Definition Template
82       call getpdstemplate(ipdsnum,mappdslen,mappds,needext,iret)
83       if (iret.ne.0) then
84         ierr=5
85         if( allocated(mappds) ) deallocate(mappds)
86         return
87       endif
88       !
89       !   Unpack each value into array ipdstmpl from the
90       !   the appropriate number of octets, which are specified in
91       !   corresponding entries in array mappds.
92       !
93       istat=0
94       if (mappdslen.gt.0) allocate(ipdstmpl(mappdslen),stat=istat)
95       if (istat.ne.0) then
96          ierr=6
97          nullify(ipdstmpl)
98          if( allocated(mappds) ) deallocate(mappds)
99          return
100       endif
101       do i=1,mappdslen
102         nbits=iabs(mappds(i))*8
103         if ( mappds(i).ge.0 ) then
104           call g2lib_gbyte(cgrib,ipdstmpl(i),iofst,nbits)
105         else
106           call g2lib_gbyte(cgrib,isign,iofst,1)
107           call g2lib_gbyte(cgrib,ipdstmpl(i),iofst+1,nbits-1)
108           if (isign.eq.1) ipdstmpl(i)=-ipdstmpl(i)
109         endif
110         iofst=iofst+nbits
111       enddo
112       !
113       !   Check to see if the Product Definition Template needs to be
114       !   extended.
115       !   The number of values in a specific template may vary
116       !   depending on data specified in the "static" part of the
117       !   template.
118       !
119       if ( needext ) then
120         call extpdstemplate(ipdsnum,ipdstmpl,newmappdslen,mappds)
121         call realloc(ipdstmpl,mappdslen,newmappdslen,istat)
122         !   Unpack the rest of the Product Definition Template
123         do i=mappdslen+1,newmappdslen
124           nbits=iabs(mappds(i))*8
125           if ( mappds(i).ge.0 ) then
126             call g2lib_gbyte(cgrib,ipdstmpl(i),iofst,nbits)
127           else
128             call g2lib_gbyte(cgrib,isign,iofst,1)
129             call g2lib_gbyte(cgrib,ipdstmpl(i),iofst+1,nbits-1)
130             if (isign.eq.1) ipdstmpl(i)=-ipdstmpl(i)
131           endif
132           iofst=iofst+nbits
133         enddo
134         mappdslen=newmappdslen
135       endif
136       if( allocated(mappds) ) deallocate(mappds)
137       !
138       !   Get Optional list of vertical coordinate values
139       !   after the Product Definition Template, if necessary.
140       !
141       nullify(coordlist)
142       if ( numcoord .ne. 0 ) then
143          allocate (coordieee(numcoord),stat=istat1)
144          allocate(coordlist(numcoord),stat=istat)
145          if ((istat1+istat).ne.0) then
146             ierr=6
147             nullify(coordlist)
148             if( allocated(coordieee) ) deallocate(coordieee)
149             return
150          endif
151         call g2lib_gbytes(cgrib,coordieee,iofst,32,0,numcoord)
152         call rdieee(coordieee,coordlist,numcoord)
153         deallocate (coordieee)
154         iofst=iofst+(32*numcoord)
155       endif
156       
157       return    ! End of Section 4 processing
158       end