NCEPLIBS-g2 4.0.0
Loading...
Searching...
No Matches
prgrib2.F90
Go to the documentation of this file.
1
4
23subroutine prlevel(ipdtn, ipdtmpl, labbrev)
24 implicit none
25
26 integer, intent(in) :: ipdtn
27 integer, intent(in) :: ipdtmpl(*)
28 character(len = 40), intent(out) :: labbrev
29 integer :: ipos
30 character(len = 10) :: tmpval1, tmpval2
31
32 labbrev(1:40) = " "
33
34 ! Use the template number to determine which element of the product
35 ! template array holds the "Type of first fixed surface" value
36 ! (which will be stored in ipos).
37 selectcase (ipdtn)
38 case(0:15)
39 ipos = 10
40 case(40:43)
41 ipos = 11
42 case(44:47)
43 ipos = 16
44 case(48)
45 ipos = 21
46 case(50:51)
47 ipos = 10
48 case(52)
49 ipos = 13
50 case(91)
51 ipos = 10
52 case default
53 ipos = 10
54 end select
55
56 ! Construct level description based on "Type of first fixed surface"
57 ! value and some other values from the template.
58 if (ipdtmpl(ipos) .eq. 100 .and. & ! Pressure Level
59 ipdtmpl(ipos + 3) .eq. 255) then
60 call frmt(tmpval1, ipdtmpl(ipos + 2), ipdtmpl(ipos + 1) + 2)
61 labbrev = trim(tmpval1)//" mb"
62 elseif (ipdtmpl(ipos) .eq. 100 .and. & ! Pressure Layer
63 ipdtmpl(ipos + 3) .eq. 100) then
64 call frmt(tmpval1, ipdtmpl(ipos + 2), ipdtmpl(ipos + 1) + 2)
65 call frmt(tmpval2, ipdtmpl(ipos + 5), ipdtmpl(ipos + 4) + 2)
66 labbrev = trim(tmpval1)//" - "//trim(tmpval2)//" mb"
67 elseif (ipdtmpl(ipos) .eq. 101) then ! Mean Sea Level
68 labbrev(1:30) = " Mean Sea Level "
69 elseif (ipdtmpl(ipos) .eq. 102 .and. & ! Altitude above MSL
70 ipdtmpl(ipos + 3) .eq. 255) then
71 call frmt(tmpval1, ipdtmpl(ipos + 2), ipdtmpl(ipos + 1))
72 labbrev = trim(tmpval1)//" m above MSL"
73 elseif (ipdtmpl(ipos) .eq. 103 .and. & ! Height above Ground
74 ipdtmpl(ipos + 3) .eq. 255) then
75 call frmt(tmpval1, ipdtmpl(ipos + 2), ipdtmpl(ipos + 1))
76 labbrev = trim(tmpval1)//" m above ground"
77 elseif (ipdtmpl(ipos) .eq. 103 .and. & ! Height above Ground
78 ipdtmpl(ipos + 3) .eq. 103) then
79 call frmt(tmpval1, ipdtmpl(ipos + 2), ipdtmpl(ipos + 1))
80 call frmt(tmpval2, ipdtmpl(ipos + 5), ipdtmpl(ipos + 4))
81 labbrev = trim(tmpval1)//" - "//trim(tmpval2)//" m AGL"
82 elseif (ipdtmpl(ipos) .eq. 104 .and. & ! Sigma Level
83 ipdtmpl(ipos + 3) .eq. 255) then
84 call frmt(tmpval1, ipdtmpl(ipos + 2), ipdtmpl(ipos + 1))
85 labbrev = trim(tmpval1)//" sigma"
86 elseif (ipdtmpl(ipos) .eq. 104 .and. & ! Sigma Layer
87 ipdtmpl(ipos + 3) .eq. 104) then
88 call frmt(tmpval1, ipdtmpl(ipos + 2), ipdtmpl(ipos + 1))
89 call frmt(tmpval2, ipdtmpl(ipos + 5), ipdtmpl(ipos + 4))
90 labbrev = trim(tmpval1)//" - "//trim(tmpval2)//" sigma"
91 elseif (ipdtmpl(ipos) .eq. 105 .and. & ! Hybrid Level
92 ipdtmpl(ipos + 3) .eq. 255) then
93 call frmt(tmpval1, ipdtmpl(ipos + 2), ipdtmpl(ipos + 1))
94 labbrev = trim(tmpval1)//" hybrid lvl"
95 elseif (ipdtmpl(ipos).eq.105 .and. & ! Hybrid Level
96 ipdtmpl(ipos + 3).eq.105) then
97 call frmt(tmpval1, ipdtmpl(ipos + 2), ipdtmpl(ipos + 1))
98 call frmt(tmpval2, ipdtmpl(ipos + 5), ipdtmpl(ipos + 4))
99 labbrev = trim(tmpval1)//" - "//trim(tmpval2)//" hybrid lvl"
100 elseif (ipdtmpl(ipos) .eq. 106 .and. & ! Depth Below Land Sfc
101 ipdtmpl(ipos + 3) .eq. 255) then
102 call frmt(tmpval1, ipdtmpl(ipos + 2), ipdtmpl(ipos + 1))
103 labbrev = trim(tmpval1)//" m below land"
104 elseif (ipdtmpl(ipos).eq.106 .and. & ! Depth Below Land Sfc Layer
105 ipdtmpl(ipos + 3).eq.106) then
106 call frmt(tmpval1, ipdtmpl(ipos + 2), ipdtmpl(ipos + 1))
107 call frmt(tmpval2, ipdtmpl(ipos + 5), ipdtmpl(ipos + 4))
108 labbrev = trim(tmpval1)//" - "//trim(tmpval2)//" m DBLY"
109 elseif (ipdtmpl(ipos) .eq. 107) then ! Isentrophic (theta) level (THEL)
110 labbrev(1:30) = " Isentropic level"
111 elseif (ipdtmpl(ipos).eq.108 .and. & ! Press Diff from Ground Layer
112 ipdtmpl(ipos + 3).eq.108) then
113 call frmt(tmpval1, ipdtmpl(ipos + 2), ipdtmpl(ipos + 1) + 2)
114 call frmt(tmpval2, ipdtmpl(ipos + 5), ipdtmpl(ipos + 4) + 2)
115 labbrev = trim(tmpval1)//" - "//trim(tmpval2)//" mb SPDY"
116 elseif (ipdtmpl(ipos) .eq. 110) then ! Layer between two hybrid levels (HYBY)
117 labbrev(1:30) = " Layer bet 2-hyb lvl"
118 elseif (ipdtmpl(ipos).eq.109 .and. & ! Potential Vorticity Sfc
119 ipdtmpl(ipos + 3).eq.255) then
120 call frmt(tmpval1, ipdtmpl(ipos + 2), ipdtmpl(ipos + 1)-6)
121 labbrev = trim(tmpval1)//" pv surface"
122 elseif (ipdtmpl(ipos) .eq. 111) then ! Eta Level (EtaL)
123 labbrev(1:30) = " Eta level"
124 elseif (ipdtmpl(ipos) .eq. 114) then ! Layer between two isentropic levels (THEY)
125 labbrev(1:30) = " Layer bet. 2-isent."
126 elseif (ipdtmpl(ipos) .eq. 117) then ! Mixed layer depth
127 labbrev(1:30) = " Mixed layer depth"
128 elseif (ipdtmpl(ipos) .eq. 120) then ! Layer between two Eta levels (EtaY)
129 labbrev(1:30) = " Layer bet. 2-Eta lvl"
130 elseif (ipdtmpl(ipos) .eq. 121) then ! Layer between two isobaric surface (IBYH)
131 labbrev(1:30) = " Layer bet. 2-isob."
132 elseif (ipdtmpl(ipos) .eq. 125) then ! Specified height level above ground (HGLH)
133 labbrev(1:30) = " Specified height lvl"
134 elseif (ipdtmpl(ipos) .eq. 126) then ! Isobaric Level (ISBP)
135 labbrev(1:30) = " Isobaric level"
136 elseif (ipdtmpl(ipos) .eq. 160) then ! Depth below sea level
137 labbrev(1:30) = " Depth below sea lvl"
138 elseif (ipdtmpl(ipos) .eq. 170) then ! Ionospheric D-region
139 labbrev(1:30) = " Ionospheric D-region lvl"
140 elseif (ipdtmpl(ipos) .eq. 1) then ! Ground/Water Surface
141 labbrev(1:30) = " Surface "
142 elseif (ipdtmpl(ipos) .eq. 2) then ! Cloud base level (CBL)
143 labbrev(1:30) = " Cloud base lvl"
144 elseif (ipdtmpl(ipos) .eq. 3) then ! Cloud top level (CTL)
145 labbrev(1:30) = " Cloud top lvl"
146 elseif (ipdtmpl(ipos) .eq. 4) then ! Freezing Level
147 labbrev(1:30) = " 0 Deg Isotherm"
148 elseif (ipdtmpl(ipos) .eq. 5) then ! Level of adiabatic condensation lifted
149 labbrev(1:30) = " Level of adiabatic" ! from the surface
150 elseif (ipdtmpl(ipos) .eq. 6) then ! Max Wind Level
151 labbrev(1:30) = " Max wind lvl"
152 elseif (ipdtmpl(ipos) .eq. 7) then ! Tropopause
153 labbrev(1:30) = " Tropopause"
154 elseif (ipdtmpl(ipos) .eq. 8) then ! Nominal top of Atmosphere
155 labbrev(1:30) = " Nom. top"
156 elseif (ipdtmpl(ipos) .eq. 9) then ! Sea bottom
157 labbrev(1:30) = " Sea Bottom"
158 elseif (ipdtmpl(ipos) .eq. 10) then ! Entire Atmosphere (EATM)
159 labbrev(1:30) = " Entire Atmosphere"
160 elseif (ipdtmpl(ipos) .eq. 11) then ! Cumulonimbus Base
161 labbrev(1:30) = " Cumulonimbus Base"
162 elseif (ipdtmpl(ipos) .eq. 12) then ! Cumulonimbus Top
163 labbrev(1:30) = " Cumulonimbus Top"
164 elseif (ipdtmpl(ipos) .eq. 20) then ! Isothermal level
165 labbrev(1:30) = " Isothermal level"
166 elseif (ipdtmpl(ipos) .eq. 200) then ! Entire Atmosphere (EATM)
167 labbrev(1:30) = " Entire Atmosphere"
168 elseif (ipdtmpl(ipos) .eq. 201) then ! Entire ocean (EOCN)
169 labbrev(1:30) = " Entire ocean"
170 elseif (ipdtmpl(ipos) .eq. 204) then ! Highest tropospheric freezing level (HTFL)
171 labbrev(1:30) = " Highest Frz. lvl"
172 elseif (ipdtmpl(ipos) .eq. 206) then ! Grid scale cloud bottom level (GCBL)
173 labbrev(1:30) = " Grid scale cloud bl"
174 elseif (ipdtmpl(ipos) .eq. 207) then ! Grid scale cloud top level (GCTL)
175 labbrev(1:30) = " Grid scale cloud tl"
176 elseif (ipdtmpl(ipos) .eq. 209) then ! Boundary layer cloud bottom level (BCBL)
177 labbrev(1:30) = " Boundary layer cbl"
178 elseif (ipdtmpl(ipos) .eq. 210) then ! Boundary layer cloud top level (BCTL)
179 labbrev(1:30) = " Boundary layer ctl"
180 elseif (ipdtmpl(ipos) .eq. 211) then ! Boundary layer cloud layer (BCY)
181 labbrev(1:30) = " Boundary layer cl"
182 elseif (ipdtmpl(ipos) .eq. 212) then ! Low cloud bottom level (LCBL)
183 labbrev(1:30) = " Low cloud bot. lvl"
184 elseif (ipdtmpl(ipos) .eq. 213) then ! Low cloud top level (LCTL)
185 labbrev(1:30) = " Low cloud top lvl"
186 elseif (ipdtmpl(ipos) .eq. 214) then ! Low cloud layer (LCY)
187 labbrev(1:30) = " Low cloud layer"
188 elseif (ipdtmpl(ipos) .eq. 215) then ! Cloud ceiling (CEIL)
189 labbrev(1:30) = " Cloud ceiling"
190 elseif (ipdtmpl(ipos) .eq. 220) then ! Planetary Boundary Layer (PBLRI)
191 labbrev(1:30) = " Planetary boundary"
192 elseif (ipdtmpl(ipos) .eq. 221) then ! Layer Between Two Hybrid Levels (HYBY)
193 labbrev(1:30) = " Layer 2 Hybrid lvl "
194 elseif (ipdtmpl(ipos) .eq. 222) then ! Middle cloud bottom level (MCBL)
195 labbrev(1:30) = " Mid. cloud bot. lvl"
196 elseif (ipdtmpl(ipos) .eq. 223) then ! Middle cloud top level (MCTL)
197 labbrev(1:30) = " Mid. cloud top lvl"
198 elseif (ipdtmpl(ipos) .eq. 224) then ! Middle cloud layer (MCY)
199 labbrev(1:30) = " Middle cloud layer"
200 elseif (ipdtmpl(ipos) .eq. 232) then ! High cloud bottom level (HCBL)
201 labbrev(1:30) = " High cloud bot. lvl"
202 elseif (ipdtmpl(ipos) .eq. 233) then ! High cloud top level (HCTL)
203 labbrev(1:30) = " High cloud top lvl"
204 elseif (ipdtmpl(ipos) .eq. 234) then ! High cloud layer (HCY)
205 labbrev(1:30) = " High cloud layer"
206 elseif (ipdtmpl(ipos) .eq. 235) then ! Ocean isotherm level (OITL)
207 labbrev(1:30) = " Ocean Isotherm lvl"
208 elseif (ipdtmpl(ipos) .eq. 236) then ! Layer between two depth below ocean sfc (OLYR)
209 labbrev(1:30) = " Layer 2-depth below"
210 elseif (ipdtmpl(ipos) .eq. 237) then ! Bottom of Ocean mixed layer (OBML)
211 labbrev(1:30) = " Bot. Ocean mix. lyr"
212 elseif (ipdtmpl(ipos) .eq. 238) then ! Bottom of Ocean iisothermal layer (OBIL)
213 labbrev(1:30) = " Bot. Ocean iso. lyr"
214 elseif (ipdtmpl(ipos) .eq. 239) then ! Layer ocean surface and 26C ocean
215 labbrev(1:30) = " layer ocean sfc 26C" ! isothermal level (S26CY)
216 elseif (ipdtmpl(ipos) .eq. 240) then ! Ocean Mixed Layer
217 labbrev(1:30) = " Ocean Mixed Layer"
218 elseif (ipdtmpl(ipos) .eq. 241) then ! Ordered Sequence of Data
219 labbrev(1:30) = " Order Seq. Of Data"
220 elseif (ipdtmpl(ipos) .eq. 242) then ! Convective cloud bottom level (CCBL)
221 labbrev(1:30) = " Con. cloud bot. lvl"
222 elseif (ipdtmpl(ipos) .eq. 243) then ! Convective cloud top level (CCTL)
223 labbrev(1:30) = " Con. cloud top lvl"
224 elseif (ipdtmpl(ipos) .eq. 244) then ! Convective cloud layer (CCY)
225 labbrev(1:30) = " Conv. cloud layer"
226 elseif (ipdtmpl(ipos) .eq. 245) then ! Lowest level of the wet bulb zero (LLTW)
227 labbrev(1:30) = " Lowest lvl wet bulb"
228 elseif (ipdtmpl(ipos) .eq. 246) then ! Maximum equiv. potential temp. level (MTHE)
229 labbrev(1:30) = " Max. equi. potential"
230 elseif (ipdtmpl(ipos) .eq. 247) then ! Equilibrium level (EHLT)
231 labbrev(1:30) = " Equilibrium level"
232 elseif (ipdtmpl(ipos) .eq. 248) then ! Shallow convective cloud bottom level (SCBL)
233 labbrev(1:30) = " Shallow con. cld bl"
234 elseif (ipdtmpl(ipos) .eq. 249) then ! Shallow convective cloud top level (SCTL)
235 labbrev(1:30) = " Shallow con. cld tl"
236 elseif (ipdtmpl(ipos) .eq. 251) then ! Deep convective cloud bottom level (DCBL)
237 labbrev(1:30) = " Deep conv. cld bl"
238 elseif (ipdtmpl(ipos) .eq. 252) then ! Deep convective cloud top level (DCTL)
239 labbrev(1:30) = " Deep conv. cld tl"
240 elseif (ipdtmpl(ipos) .eq. 253) then ! Lowest bottom level of supercooled
241 labbrev(1:30) = " Lowest bot. lvl sup" ! liquid water layer (LBLSW)
242 elseif (ipdtmpl(ipos) .eq. 254) then ! Highest top level of supercooled
243 labbrev(1:30) = " highest top lvl sup" ! liquid water layer (HBLSW)
244 else
245 write(labbrev, fmt = '(1x,I4," (Unknown Lvl)")') ipdtmpl(ipos)
246 endif
247
248end subroutine prlevel
249
257subroutine frmt(cval, ival, iscal)
258 implicit none
259
260 character(len = 10), intent(out) :: cval
261 integer, intent(in) :: ival, iscal
262
263 real :: rval
264 integer :: newscal
265 character(len = 7) :: cformat
266
267 if (iscal .eq. 0) then
268 write(cval, '(I0)') ival
269 else
270 newscal = -1 * iscal
271 rval = real(ival) * (10.0**newscal)
272 if (rval .eq. real(nint(rval))) then
273 write(cval, '(1X,I0)') nint(rval)
274 else
275 write(cformat, fmt = '("(f0.",I1,")")') iabs(iscal)
276 write(cval, fmt = cformat) rval
277 endif
278 endif
279
280 return
281end subroutine frmt
282
296subroutine prvtime(ipdtn, ipdtmpl, listsec1, tabbrev)
297 implicit none
298
299 integer, intent(in) :: ipdtn
300 integer, intent(in) :: ipdtmpl(*), listsec1(*)
301 character(len = 110), intent(out) :: tabbrev
302
303 character(len = 16) :: reftime, endtime
304 character(len = 12) :: tmpval2
305 character(len = 12) :: tmpval
306 character(len = 10) :: tunit
307 integer, dimension(200) :: ipos, ipos2
308 integer :: is, itemp, itemp2, iunit, iuni2t2, iunit2, iutpos, iutpos2, j
309
310 data ipos /7*0, 16, 23, 17, 19, 18, 32, 31, 27*0, 17, 20, 0, 0, 22, &
311 25, 43*0, 23, 109*0/
312
313 data ipos2 /7*0, 26, 33, 27, 29, 28, 42, 41, 27*0, 22, 30, 0, 0, 32, &
314 35, 43*0, 33, 109*0/
315
316 tabbrev(1:100) = " "
317
318 ! Determine unit of time range.
319 if ((ipdtn .ge. 0 .and. ipdtn .le. 15) .or. ipdtn .eq. 32 &
320 .or. ipdtn .eq. 50 .or. ipdtn .eq. 51 &
321 .or. ipdtn .eq. 91) then
322 iutpos = 8
323 elseif (ipdtn .ge. 40 .and. ipdtn .le. 43) then
324 iutpos = 9
325 elseif (ipdtn .ge. 44 .and. ipdtn .le. 47) then
326 iutpos = 14
327 elseif (ipdtn .eq. 48) then
328 iutpos = 19
329 elseif (ipdtn .eq. 52) then
330 iutpos = 11
331 else
332 iutpos = 8
333 endif
334
335 ! Determine first unit of time range.
336 selectcase(ipdtmpl(iutpos))
337 case (0)
338 tunit = "minute"
339 iunit = 1
340 case (1)
341 tunit = "hour"
342 iunit = 1
343 case (2)
344 tunit = "day"
345 iunit = 1
346 case (3)
347 tunit = "month"
348 iunit = 1
349 case (4)
350 tunit = "year"
351 iunit = 1
352 case (10)
353 tunit = "hour"
354 iunit = 3
355 case (11)
356 tunit = "hour"
357 iunit = 6
358 case default
359 tunit = "hour"
360 iunit = 1
361 end select
362
363 ! Determine second unit of time range.
364 if (ipdtn .eq. 0) then
365 iunit2 = 1
366 iutpos2 = 0
367 else
368 iutpos2 = ipos2(ipdtn)
369 if (iutpos2 .gt. 0) then
370 selectcase(ipdtmpl(iutpos2))
371 case (0)
372 iunit2 = 1
373 case (1)
374 iunit2 = 1
375 case (2)
376 iunit2 = 1
377 case (3)
378 iuni2t2 = 1
379 case (4)
380 iunit2 = 1
381 case (10)
382 iunit2 = 3
383 case (11)
384 iunit2 = 6
385 case default
386 iunit2 = 1
387 end select
388 endif
389 endif
390
391 write(reftime, fmt = '(i4,3i2.2,":",i2.2,":",i2.2)') (listsec1(j), j = 6, 11)
392 itemp = abs(ipdtmpl(iutpos + 1)) * iunit
393 write(tmpval, '(I0)') itemp
394 write(tabbrev, fmt = '("valid at ", i4)') ipdtmpl(iutpos + 1)
395
396 ! Determine Reference Time: Year, Month, Day, Hour, Minute, Second.
397 if ((ipdtn .ge. 0 .and. ipdtn .le. 7) .or. ipdtn .eq. 15 &
398 .or. ipdtn .eq. 20 .or. (ipdtn .ge. 30 .and. ipdtn .le. 32) &
399 .or. ipdtn .eq. 40 .or. ipdtn .eq. 41 .or. ipdtn .eq. 44 &
400 .or. ipdtn .eq. 45 .or. ipdtn .eq. 48 .or. &
401 (ipdtn .ge. 50 .and. ipdtn .le. 52)) then ! Point in time
402 tabbrev = "valid " // trim(tmpval) // " " // trim(tunit) // " after " // reftime
403 else
404 is = ipos(ipdtn) ! Continuous time interval
405 write(endtime, fmt = '(i4,3i2.2,":",i2.2,":",i2.2)') (ipdtmpl(j), j = is, is + 5)
406 itemp2 = abs(ipdtmpl(iutpos2 + 1)) * iunit2
407 itemp2 = itemp + itemp2
408 write(tmpval2, '(I0)') itemp2
409 if (ipdtn .eq. 8 .and. ipdtmpl(9) .lt. 0) then
410 tabbrev = "(" // trim(tmpval) // " -" &
411 // trim(tmpval2) // ") valid " // trim(tmpval) // &
412 " " // trim(tunit) // " before " &
413 // reftime // " to " //endtime
414 elseif ((ipdtn .ge. 8 .and. ipdtn .le. 14) .or. &
415 (ipdtn .ge. 42 .and. ipdtn .le. 47) .or. &
416 ipdtn .eq. 91) then ! Continuous time interval
417 tabbrev = "(" // trim(tmpval) // " -" &
418 // trim(tmpval2) // " hr) valid " // trim(tmpval) // &
419 " " // trim(tunit) // " after " &
420 // reftime // " to " // endtime
421 endif
422 endif
423
424 return
425end subroutine prvtime
subroutine frmt(cval, ival, iscal)
Format the level description.
Definition prgrib2.F90:258
subroutine prlevel(ipdtn, ipdtmpl, labbrev)
Print level information to a character array, given the GRIB2 Product Definition Template information...
Definition prgrib2.F90:24
subroutine prvtime(ipdtn, ipdtmpl, listsec1, tabbrev)
Convert date and time from GRIB2 info to string output.
Definition prgrib2.F90:297