ungrib build
[WPS.git] / ungrib / src / ngl / g2 / gf_unpack3.f
blob3ed3268e24381a14123a9a6cf2617b86a2f6ebca
1 subroutine gf_unpack3(cgrib,lcgrib,iofst,igds,igdstmpl,
2 & mapgridlen,ideflist,idefnum,ierr)
3 !$$$ SUBPROGRAM DOCUMENTATION BLOCK
4 ! . . . .
5 ! SUBPROGRAM: gf_unpack3
6 ! PRGMMR: Gilbert ORG: W/NP11 DATE: 2000-05-26
8 ! ABSTRACT: This subroutine unpacks Section 3 (Grid 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_unpack3(cgrib,lcgrib,lensec,iofst,igds,igdstmpl,
18 ! & mapgridlen,ideflist,idefnum,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 3.
24 ! OUTPUT ARGUMENT LIST:
25 ! iofst - Bit offset at the end of Section 3, returned.
26 ! igds - Contains information read from the appropriate GRIB Grid
27 ! Definition Section 3 for the field being returned.
28 ! Must be dimensioned >= 5.
29 ! igds(1)=Source of grid definition (see Code Table 3.0)
30 ! igds(2)=Number of grid points in the defined grid.
31 ! igds(3)=Number of octets needed for each
32 ! additional grid points definition.
33 ! Used to define number of
34 ! points in each row ( or column ) for
35 ! non-regular grids.
36 ! = 0, if using regular grid.
37 ! igds(4)=Interpretation of list for optional points
38 ! definition. (Code Table 3.11)
39 ! igds(5)=Grid Definition Template Number (Code Table 3.1)
40 ! igdstmpl - Pointer to integer array containing the data values for
41 ! the specified Grid Definition
42 ! Template ( NN=igds(5) ). Each element of this integer
43 ! array contains an entry (in the order specified) of Grid
44 ! Defintion Template 3.NN
45 ! mapgridlen- Number of elements in igdstmpl(). i.e. number of entries
46 ! in Grid Defintion Template 3.NN ( NN=igds(5) ).
47 ! ideflist - (Used if igds(3) .ne. 0) Pointer to integer array containing
48 ! the number of grid points contained in each row ( or column ).
49 ! (part of Section 3)
50 ! idefnum - (Used if igds(3) .ne. 0) The number of entries
51 ! in array ideflist. i.e. number of rows ( or columns )
52 ! for which optional grid points are defined.
53 ! ierr - Error return code.
54 ! 0 = no error
55 ! 5 = "GRIB" message contains an undefined Grid Definition
56 ! Template.
57 ! 6 = memory allocation error
59 ! REMARKS: Uses Fortran 90 module gridtemplates and module re_alloc.
61 ! ATTRIBUTES:
62 ! LANGUAGE: Fortran 90
63 ! MACHINE: IBM SP
65 !$$$
67 use gridtemplates
68 use re_alloc ! needed for subroutine realloc
70 character(len=1),intent(in) :: cgrib(lcgrib)
71 integer,intent(in) :: lcgrib
72 integer,intent(inout) :: iofst
73 integer,pointer,dimension(:) :: igdstmpl,ideflist
74 integer,intent(out) :: igds(5)
75 integer,intent(out) :: ierr,idefnum
77 integer,allocatable :: mapgrid(:)
78 integer :: mapgridlen,ibyttem
79 logical needext
81 ierr=0
82 nullify(igdstmpl,ideflist)
84 call gbyte(cgrib,lensec,iofst,32) ! Get Length of Section
85 iofst=iofst+32
86 iofst=iofst+8 ! skip section number
88 call gbyte(cgrib,igds(1),iofst,8) ! Get source of Grid def.
89 iofst=iofst+8
90 call gbyte(cgrib,igds(2),iofst,32) ! Get number of grid pts.
91 iofst=iofst+32
92 call gbyte(cgrib,igds(3),iofst,8) ! Get num octets for opt. list
93 iofst=iofst+8
94 call gbyte(cgrib,igds(4),iofst,8) ! Get interpret. for opt. list
95 iofst=iofst+8
96 call gbyte(cgrib,igds(5),iofst,16) ! Get Grid Def Template num.
97 iofst=iofst+16
98 ! if (igds(1).eq.0) then
99 if (igds(1).eq.0.OR.igds(1).eq.255) then ! FOR ECMWF TEST ONLY
100 allocate(mapgrid(lensec))
101 ! Get Grid Definition Template
102 call getgridtemplate(igds(5),mapgridlen,mapgrid,needext,
103 & iret)
104 if (iret.ne.0) then
105 ierr=5
106 if( allocated(mapgrid) ) deallocate(mapgrid)
107 return
108 endif
109 else
110 ! igdstmpl=-1
111 mapgridlen=0
112 needext=.false.
113 endif
115 ! Unpack each value into array igdstmpl from the
116 ! the appropriate number of octets, which are specified in
117 ! corresponding entries in array mapgrid.
119 istat=0
120 if (mapgridlen.gt.0) allocate(igdstmpl(mapgridlen),stat=istat)
121 if (istat.ne.0) then
122 ierr=6
123 nullify(igdstmpl)
124 if( allocated(mapgrid) ) deallocate(mapgrid)
125 return
126 endif
127 ibyttem=0
128 do i=1,mapgridlen
129 nbits=iabs(mapgrid(i))*8
130 if ( mapgrid(i).ge.0 ) then
131 call gbyte(cgrib,igdstmpl(i),iofst,nbits)
132 else
133 call gbyte(cgrib,isign,iofst,1)
134 call gbyte(cgrib,igdstmpl(i),iofst+1,nbits-1)
135 if (isign.eq.1) igdstmpl(i)=-igdstmpl(i)
136 endif
137 iofst=iofst+nbits
138 ibyttem=ibyttem+iabs(mapgrid(i))
139 enddo
141 ! Check to see if the Grid Definition Template needs to be
142 ! extended.
143 ! The number of values in a specific template may vary
144 ! depending on data specified in the "static" part of the
145 ! template.
147 if ( needext ) then
148 call extgridtemplate(igds(5),igdstmpl,newmapgridlen,mapgrid)
149 ! Unpack the rest of the Grid Definition Template
150 call realloc(igdstmpl,mapgridlen,newmapgridlen,istat)
151 do i=mapgridlen+1,newmapgridlen
152 nbits=iabs(mapgrid(i))*8
153 if ( mapgrid(i).ge.0 ) then
154 call gbyte(cgrib,igdstmpl(i),iofst,nbits)
155 else
156 call gbyte(cgrib,isign,iofst,1)
157 call gbyte(cgrib,igdstmpl(i),iofst+1,nbits-1)
158 if (isign.eq.1) igdstmpl(i)=-igdstmpl(i)
159 endif
160 iofst=iofst+nbits
161 ibyttem=ibyttem+iabs(mapgrid(i))
162 enddo
163 mapgridlen=newmapgridlen
164 endif
165 if( allocated(mapgrid) ) deallocate(mapgrid)
167 ! Unpack optional list of numbers defining number of points
168 ! in each row or column, if included. This is used for non regular
169 ! grids.
171 if ( igds(3).ne.0 ) then
172 nbits=igds(3)*8
173 idefnum=(lensec-14-ibyttem)/igds(3)
174 istat=0
175 if (idefnum.gt.0) allocate(ideflist(idefnum),stat=istat)
176 if (istat.ne.0) then
177 ierr=6
178 nullify(ideflist)
179 return
180 endif
181 call gbytes(cgrib,ideflist,iofst,nbits,0,idefnum)
182 iofst=iofst+(nbits*idefnum)
183 else
184 idefnum=0
185 nullify(ideflist)
186 endif
188 return ! End of Section 3 processing