Created a tag for the 2012 HWRF baseline tests.
[WPS-merge.git] / hwrf-baseline-20120103-1354 / ungrib / src / ngl / g2 / gdt2gds.f
blob705952245ed15e366e343c4d7cab78944465f9a7
1 subroutine gdt2gds(igds,igdstmpl,idefnum,ideflist,kgds,
2 & igrid,iret)
3 C$$$ SUBPROGRAM DOCUMENTATION BLOCK
4 C . . . .
5 C SUBPROGRAM: gdt2gds
6 C PRGMMR: Gilbert ORG: W/NP11 DATE: 2003-06-17
8 C ABSTRACT: This routine converts grid information from a GRIB2
9 C Grid Description Section as well as its
10 C Grid Definition Template to GRIB1 GDS info. In addition,
11 C a check is made to determine if the grid is an NCEP
12 C predefined grid.
14 C PROGRAM HISTORY LOG:
15 C 2003-06-17 Gilbert
16 C 2004-04-27 Gilbert - Added support for gaussian grids.
17 C 2007-04-16 Vuong - Added Curvilinear Orthogonal grids.
18 C 2007-05-29 Vuong - Added Rotate Lat/Lon E-grid (203)
20 C USAGE: CALL gdt2gds(igds,igdstmpl,idefnum,ideflist,kgds,igrid,iret)
21 C INPUT ARGUMENT LIST:
22 C igds() - Contains information read from the appropriate GRIB Grid
23 C Definition Section 3 for the field being returned.
24 C Must be dimensioned >= 5.
25 C igds(1)=Source of grid definition (see Code Table 3.0)
26 C igds(2)=Number of grid points in the defined grid.
27 C igds(3)=Number of octets needed for each
28 C additional grid points definition.
29 C Used to define number of
30 C points in each row ( or column ) for
31 C non-regular grids.
32 C = 0, if using regular grid.
33 C igds(4)=Interpretation of list for optional points
34 C definition. (Code Table 3.11)
35 C igds(5)=Grid Definition Template Number (Code Table 3.1)
36 C igdstmpl() - Grid Definition Template values for GDT 3.igds(5)
37 C idefnum - The number of entries in array ideflist.
38 C i.e. number of rows ( or columns )
39 C for which optional grid points are defined.
40 C ideflist() - Optional integer array containing
41 C the number of grid points contained in each row (or column).
43 C OUTPUT ARGUMENT LIST: (INCLUDING WORK ARRAYS)
44 C kgds() - GRIB1 GDS as described in w3fi63 format.
45 C igrid - NCEP predefined GRIB1 grid number
46 C set to 255, if not NCEP grid
47 C iret - Error return value:
48 C 0 = Successful
49 C 1 = Unrecognized GRIB2 GDT number 3.igds(5)
51 C REMARKS: LIST CAVEATS, OTHER HELPFUL HINTS OR INFORMATION
53 C ATTRIBUTES:
54 C LANGUAGE: INDICATE EXTENSIONS, COMPILER OPTIONS
55 C MACHINE: IBM SP
57 C$$$
59 integer,intent(in) :: idefnum
60 integer,intent(in) :: igds(*),igdstmpl(*),ideflist(*)
61 integer,intent(out) :: kgds(*),igrid,iret
63 integer :: kgds72(200),kgds71(200),idum(200),jdum(200)
65 iret=0
66 if (igds(5).eq.0) then ! Lat/Lon grid
67 kgds(1)=0
68 kgds(2)=igdstmpl(8) ! Ni
69 kgds(3)=igdstmpl(9) ! Nj
70 kgds(4)=igdstmpl(12)/1000 ! Lat of 1st grid point
71 kgds(5)=igdstmpl(13)/1000 ! Long of 1st grid point
72 kgds(6)=0 ! resolution and component flags
73 if (igdstmpl(1)==2 ) kgds(6)=64
74 if ( btest(igdstmpl(14),4).OR.btest(igdstmpl(14),5) )
75 & kgds(6)=kgds(6)+128
76 if ( btest(igdstmpl(14),3) ) kgds(6)=kgds(6)+8
77 kgds(7)=igdstmpl(15)/1000 ! Lat of last grid point
78 kgds(8)=igdstmpl(16)/1000 ! Long of last grid point
79 kgds(9)=igdstmpl(17)/1000 ! Di
80 kgds(10)=igdstmpl(18)/1000 ! Dj
81 kgds(11)=igdstmpl(19) ! Scanning mode
82 kgds(12)=0
83 kgds(13)=0
84 kgds(14)=0
85 kgds(15)=0
86 kgds(16)=0
87 kgds(17)=0
88 kgds(18)=0
89 kgds(19)=0
90 kgds(20)=255
91 kgds(21)=0
92 kgds(22)=0
94 ! Process irreg grid stuff, if necessary
96 if ( idefnum.ne.0 ) then
97 if ( igdstmpl(8).eq.-1 ) then
98 kgds(2)=65535
99 kgds(9)=65535
100 endif
101 if ( igdstmpl(9).eq.-1 ) then
102 kgds(3)=65535
103 kgds(10)=65535
104 endif
105 kgds(19)=0
106 kgds(20)=33
107 if ( kgds(1).eq.1.OR.kgds(1).eq.3 ) kgds(20)=43
108 kgds(21)=igds(2) ! num of grid points
109 do j=1,idefnum
110 kgds(21+j)=ideflist(j)
111 enddo
112 endif
113 elseif (igds(5).eq.10) then ! Mercator grid
114 kgds(1)=1 ! Grid Definition Template number
115 kgds(2)=igdstmpl(8) ! Ni
116 kgds(3)=igdstmpl(9) ! Nj
117 kgds(4)=igdstmpl(10)/1000 ! Lat of 1st grid point
118 kgds(5)=igdstmpl(11)/1000 ! Long of 1st grid point
119 kgds(6)=0 ! resolution and component flags
120 if (igdstmpl(1)==2 ) kgds(6)=64
121 if ( btest(igdstmpl(12),4).OR.btest(igdstmpl(12),5) )
122 & kgds(6)=kgds(6)+128
123 if ( btest(igdstmpl(12),3) ) kgds(6)=kgds(6)+8
124 kgds(7)=igdstmpl(14)/1000 ! Lat of last grid point
125 kgds(8)=igdstmpl(15)/1000 ! Long of last grid point
126 kgds(9)=igdstmpl(13)/1000 ! Lat intersects earth
127 kgds(10)=0
128 kgds(11)=igdstmpl(16) ! Scanning mode
129 kgds(12)=igdstmpl(18)/1000 ! Di
130 kgds(13)=igdstmpl(19)/1000 ! Dj
131 kgds(14)=0
132 kgds(15)=0
133 kgds(16)=0
134 kgds(17)=0
135 kgds(18)=0
136 kgds(19)=0
137 kgds(20)=255
138 kgds(21)=0
139 kgds(22)=0
140 elseif (igds(5).eq.30) then ! Lambert Conformal Grid
141 kgds(1)=3
142 kgds(2)=igdstmpl(8) ! Nx
143 kgds(3)=igdstmpl(9) ! Ny
144 kgds(4)=igdstmpl(10)/1000 ! Lat of 1st grid point
145 kgds(5)=igdstmpl(11)/1000 ! Long of 1st grid point
146 kgds(6)=0 ! resolution and component flags
147 if (igdstmpl(1)==2 ) kgds(6)=64
148 if ( btest(igdstmpl(12),4).OR.btest(igdstmpl(12),5) )
149 & kgds(6)=kgds(6)+128
150 if ( btest(igdstmpl(12),3) ) kgds(6)=kgds(6)+8
151 kgds(7)=igdstmpl(14)/1000 ! Lon of orientation
152 kgds(8)=igdstmpl(15)/1000 ! Dx
153 kgds(9)=igdstmpl(16)/1000 ! Dy
154 kgds(10)=igdstmpl(17) ! Projection Center Flag
155 kgds(11)=igdstmpl(18) ! Scanning mode
156 kgds(12)=igdstmpl(19)/1000 ! Lat in 1
157 kgds(13)=igdstmpl(20)/1000 ! Lat in 2
158 kgds(14)=igdstmpl(21)/1000 ! Lat of S. Pole of projection
159 kgds(15)=igdstmpl(22)/1000 ! Lon of S. Pole of projection
160 kgds(16)=0
161 kgds(17)=0
162 kgds(18)=0
163 kgds(19)=0
164 kgds(20)=255
165 kgds(21)=0
166 kgds(22)=0
167 elseif (igds(5).eq.40) then ! Gaussian Lat/Lon grid
168 kgds(1)=4
169 kgds(2)=igdstmpl(8) ! Ni
170 kgds(3)=igdstmpl(9) ! Nj
171 kgds(4)=igdstmpl(12)/1000 ! Lat of 1st grid point
172 kgds(5)=igdstmpl(13)/1000 ! Long of 1st grid point
173 kgds(6)=0 ! resolution and component flags
174 if (igdstmpl(1)==2 ) kgds(6)=64
175 if ( btest(igdstmpl(14),4).OR.btest(igdstmpl(14),5) )
176 & kgds(6)=kgds(6)+128
177 if ( btest(igdstmpl(14),3) ) kgds(6)=kgds(6)+8
178 kgds(7)=igdstmpl(15)/1000 ! Lat of last grid point
179 kgds(8)=igdstmpl(16)/1000 ! Long of last grid point
180 kgds(9)=igdstmpl(17)/1000 ! Di
181 kgds(10)=igdstmpl(18) ! N - Number of parallels
182 kgds(11)=igdstmpl(19) ! Scanning mode
183 kgds(12)=0
184 kgds(13)=0
185 kgds(14)=0
186 kgds(15)=0
187 kgds(16)=0
188 kgds(17)=0
189 kgds(18)=0
190 kgds(19)=0
191 kgds(20)=255
192 kgds(21)=0
193 kgds(22)=0
194 elseif (igds(5).eq.20) then ! Polar Stereographic Grid
195 kgds(1)=5
196 kgds(2)=igdstmpl(8) ! Nx
197 kgds(3)=igdstmpl(9) ! Ny
198 kgds(4)=igdstmpl(10)/1000 ! Lat of 1st grid point
199 kgds(5)=igdstmpl(11)/1000 ! Long of 1st grid point
200 kgds(6)=0 ! resolution and component flags
201 if (igdstmpl(1)==2 ) kgds(6)=64
202 if ( btest(igdstmpl(12),4).OR.btest(igdstmpl(12),5) )
203 & kgds(6)=kgds(6)+128
204 if ( btest(igdstmpl(12),3) ) kgds(6)=kgds(6)+8
205 kgds(7)=igdstmpl(14)/1000 ! Lon of orientation
206 kgds(8)=igdstmpl(15)/1000 ! Dx
207 kgds(9)=igdstmpl(16)/1000 ! Dy
208 kgds(10)=igdstmpl(17) ! Projection Center Flag
209 kgds(11)=igdstmpl(18) ! Scanning mode
210 kgds(12)=0
211 kgds(13)=0
212 kgds(14)=0
213 kgds(15)=0
214 kgds(16)=0
215 kgds(17)=0
216 kgds(18)=0
217 kgds(19)=0
218 kgds(20)=255
219 kgds(21)=0
220 kgds(22)=0
221 elseif (igds(5).eq.204) then ! Curvilinear Orthogonal
222 kgds(1)=204
223 kgds(2)=igdstmpl(8) ! Ni
224 kgds(3)=igdstmpl(9) ! Nj
225 kgds(4)=0
226 kgds(5)=0
227 kgds(6)=0 ! resolution and component flags
228 if (igdstmpl(1)==2 ) kgds(6)=64
229 if ( btest(igdstmpl(14),4).OR.btest(igdstmpl(14),5) )
230 & kgds(6)=kgds(6)+128
231 if ( btest(igdstmpl(14),3) ) kgds(6)=kgds(6)+8
232 kgds(7)=0
233 kgds(8)=0
234 kgds(9)=0
235 kgds(10)=0
236 kgds(11)=igdstmpl(19) ! Scanning mode
237 kgds(12)=0
238 kgds(13)=0
239 kgds(14)=0
240 kgds(15)=0
241 kgds(16)=0
242 kgds(17)=0
243 kgds(18)=0
244 kgds(19)=0
245 kgds(20)=255
246 kgds(21)=0
247 kgds(22)=0
249 ! Process irreg grid stuff, if necessary
251 if ( idefnum.ne.0 ) then
252 if ( igdstmpl(8).eq.-1 ) then
253 kgds(2)=65535
254 kgds(9)=65535
255 endif
256 if ( igdstmpl(9).eq.-1 ) then
257 kgds(3)=65535
258 kgds(10)=65535
259 endif
260 kgds(19)=0
261 kgds(20)=33
262 if ( kgds(1).eq.1.OR.kgds(1).eq.3 ) kgds(20)=43
263 kgds(21)=igds(2) ! num of grid points
264 do j=1,idefnum
265 kgds(21+j)=ideflist(j)
266 enddo
267 endif
268 elseif (igds(5).eq.32768) then ! Rotate Lat/Lon grid
269 kgds(1)=0 ! Arakawa Staggerred E/B grid
270 kgds(2)=igdstmpl(8) ! Ni
271 kgds(3)=igdstmpl(9) ! Nj
272 kgds(4)=igdstmpl(12)/1000 ! Lat of 1st grid point
273 kgds(5)=igdstmpl(13)/1000 ! Long of 1st grid point
274 kgds(6)=0 ! resolution and component flags
275 if (igdstmpl(1)==2 ) kgds(6)=64
276 if ( btest(igdstmpl(14),4).OR.btest(igdstmpl(14),5) )
277 & kgds(6)=kgds(6)+128
278 if ( btest(igdstmpl(14),3) ) kgds(6)=kgds(6)+8
279 kgds(7)=igdstmpl(15)/1000 ! Lat of last grid point
280 kgds(8)=igdstmpl(16)/1000 ! Long of last grid point
281 kgds(9)=igdstmpl(17)/1000 ! Di
282 kgds(10)=igdstmpl(18)/1000 ! Dj
283 kgds(11)=igdstmpl(19) ! Scanning mode
284 kgds(12)=0
285 kgds(13)=0
286 kgds(14)=0
287 kgds(15)=0
288 kgds(16)=0
289 kgds(17)=0
290 kgds(18)=0
291 kgds(19)=0
292 kgds(20)=255
293 kgds(21)=0
294 kgds(22)=0
296 ! Process irreg grid stuff, if necessary
298 if ( idefnum.ne.0 ) then
299 if ( igdstmpl(8).eq.-1 ) then
300 kgds(2)=65535
301 kgds(9)=65535
302 endif
303 if ( igdstmpl(9).eq.-1 ) then
304 kgds(3)=65535
305 kgds(10)=65535
306 endif
307 kgds(19)=0
308 kgds(20)=33
309 if ( kgds(1).eq.1.OR.kgds(1).eq.3 ) kgds(20)=43
310 kgds(21)=igds(2) ! num of grid points
311 do j=1,idefnum
312 kgds(21+j)=ideflist(j)
313 enddo
314 endif
315 elseif (igds(5).eq.32769) then ! Rotate Lat/Lon grid
316 kgds(1)=0 ! Arakawa Staggerred for Non-E Stagger grid
317 kgds(2)=igdstmpl(8) ! Ni
318 kgds(3)=igdstmpl(9) ! Nj
319 kgds(4)=igdstmpl(12)/1000 ! Lat of 1st grid point
320 kgds(5)=igdstmpl(13)/1000 ! Long of 1st grid point
321 kgds(6)=0 ! resolution and component flags
322 if (igdstmpl(1)==2 ) kgds(6)=64
323 if ( btest(igdstmpl(14),4).OR.btest(igdstmpl(14),5) )
324 & kgds(6)=kgds(6)+128
325 if ( btest(igdstmpl(14),3) ) kgds(6)=kgds(6)+8
326 kgds(7)=igdstmpl(15)/1000 ! Lat of last grid point
327 kgds(8)=igdstmpl(16)/1000 ! Long of last grid point
328 kgds(9)=igdstmpl(17)/1000 ! Di
329 kgds(10)=igdstmpl(18)/1000 ! Dj
330 kgds(11)=igdstmpl(19) ! Scanning mode
331 kgds(12)=igdstmpl(20)/1000
332 kgds(13)=igdstmpl(21)/1000
333 kgds(14)=0
334 kgds(15)=0
335 kgds(16)=0
336 kgds(17)=0
337 kgds(18)=0
338 kgds(19)=0
339 kgds(20)=255
340 kgds(21)=0
341 kgds(22)=0
342 else
343 Print *,'gdt2gds: Unrecognized GRIB2 GDT = 3.',igds(5)
344 iret=1
345 kgds(1:22)=0
346 return
347 endif
349 ! Can we determine NCEP grid number ?
351 igrid=255
352 do j=254,1,-1
353 !do j=225,225
354 kgds71=0
355 kgds72=0
356 call w3fi71(j,kgds71,ierr)
357 if ( ierr.ne.0 ) cycle
358 ! convert W to E for longitudes
359 if ( kgds71(3).eq.0 ) then ! lat/lon
360 if ( kgds71(7).lt.0 ) kgds71(7)=360000+kgds71(7)
361 if ( kgds71(10).lt.0 ) kgds71(10)=360000+kgds71(10)
362 elseif ( kgds71(3).eq.1 ) then ! mercator
363 if ( kgds71(7).lt.0 ) kgds71(7)=360000+kgds71(7)
364 if ( kgds71(10).lt.0 ) kgds71(10)=360000+kgds71(10)
365 elseif ( kgds71(3).eq.3 ) then ! lambert conformal
366 if ( kgds71(7).lt.0 ) kgds71(7)=360000+kgds71(7)
367 if ( kgds71(9).lt.0 ) kgds71(9)=360000+kgds71(9)
368 if ( kgds71(18).lt.0 ) kgds71(18)=360000+kgds71(18)
369 elseif ( kgds71(3).eq.4 ) then ! Guassian lat/lon
370 if ( kgds71(7).lt.0 ) kgds71(7)=360000+kgds71(7)
371 if ( kgds71(10).lt.0 ) kgds71(10)=360000+kgds71(10)
372 elseif ( kgds71(3).eq.5 ) then ! polar stereographic
373 if ( kgds71(7).lt.0 ) kgds71(7)=360000+kgds71(7)
374 if ( kgds71(9).lt.0 ) kgds71(9)=360000+kgds71(9)
375 endif
376 call r63w72(idum,kgds,jdum,kgds72)
377 if ( kgds72(3).eq.3 ) kgds72(14)=0 ! lambert conformal fix
378 if ( kgds72(3).eq.1 ) kgds72(15:18)=0 ! mercator fix
379 if ( kgds72(3).eq.5 ) kgds72(14:18)=0 ! polar str fix
380 c print *,' kgds71(',j,')= ', kgds71(1:30)
381 c print *,' kgds72 = ', kgds72(1:30)
382 if ( all(kgds71.eq.kgds72) ) then
383 igrid=j
384 return
385 endif
386 enddo
388 return