NCEPLIBS-g2 4.0.0
Loading...
Searching...
No Matches
cnv12.F90
Go to the documentation of this file.
1
4
47subroutine cnv12(ifl1, ifl2, ipack, usemiss, imiss, uvvect, table_ver)
48
49 use params
50 use params_ecmwf
51 integer, intent(in) :: ifl1, ifl2, ipack
52 logical, intent(in) :: usemiss, uvvect
53
54 parameter(maxpts = 40000000, msk1 = 32000)
55 CHARACTER(len = 1), allocatable, dimension(:) :: cgrib, cgribin
56 integer KPDS(200), KGDS(200), KPTR(200)
57 integer LPDS(200), LGDS(200), KENS(200), LENS(200)
58 integer KPROB(2), KCLUST(16), KMEMBR(80)
59 real XPROB(2)
60 real, allocatable, dimension(:) :: FLD
61 real, allocatable, dimension(:) :: FLDV
62 real, allocatable, dimension(:) :: coordlist
63 integer :: listsec0(2) = (/0, 2/), imiss
64 integer :: listsec1(13) = (/7, 0, 2, 1, 1, 0, 0, 0, 0, 0, 0, 0, 0/)
65 integer :: ideflist(MAXPTS), idefnum
66 integer :: igds(5) = (/0, 0, 0, 0, 0/), igdstmpl(200), ipdstmpl(200)
67 integer :: ipdstmplv(200)
68 integer :: idrstmpl(200), idrstmplv(200)
69 integer :: currlen = 0, table_ver
70 integer, parameter :: mingrib = 500
71 logical :: ensemble, ecmwf
72 Logical*1, allocatable, dimension(:) :: bmp, bmpv
73 !
74 icnd = 0
75 ifli1 = 0
76 allocate(fld(maxpts))
77 allocate(coordlist(maxpts))
78 allocate(bmp(maxpts))
79 listsec1(3) = table_ver
80 !
81 iseek = 0
82 do
83 call skgb(ifl1, iseek, msk1, lskip, lgrib)
84 if (lgrib .eq. 0) exit ! end loop at EOF or problem
85 if (lgrib.gt.currlen) then
86 if (allocated(cgribin)) deallocate(cgribin)
87 allocate(cgribin(lgrib), stat = is)
88 currlen = lgrib
89 lcgrib = lgrib*2
90 if (lcgrib .lt. mingrib) lcgrib = mingrib
91 if (allocated(cgrib)) deallocate(cgrib)
92 allocate(cgrib(lcgrib), stat = is)
93 endif
94 call baread(ifl1, lskip, lgrib, lengrib, cgribin)
95 if (lgrib .eq. lengrib) then
96 call w3fi63(cgribin, kpds, kgds, bmp, fld, kptr, iret)
97 numpts = kptr(10)
98 if (iret .ne. 0) then
99 print *, ' cnvgrib: Error unpacking GRIB field.', iret
100 iseek = lskip+lgrib
101 cycle
102 endif
103 else
104 print *, ' cnvgrib: IO Error on input GRIB file.'
105 cycle
106 endif
107 iseek = lskip+lgrib
108 !print *, 'kpds:', kpds(1:28)
109 !print *, 'kpds:', kpds(1:45)
110 if ((kpds(5) .eq. 34).AND.uvvect) cycle ! V-comp already processed with U
111 listsec1(1) = kpds(1)
112 listsec1(2) = kpds(23)
113 listsec1(5) = 1
114 if (kpds(16) .eq. 1) listsec1(5) = 0
115 listsec1(6) = ((kpds(21)-1)*100)+kpds(8)
116 listsec1(7) = kpds(9)
117 listsec1(8) = kpds(10)
118 listsec1(9) = kpds(11)
119 listsec1(10) = kpds(12)
120 listsec1(13) = 1
121 if (kpds(16) .eq. 1) listsec1(13) = 0
122 ensemble = .false.
123 if ((kpds(23) .eq. 2) .or. &
124 (kptr(3).gt.28 .and. kpds(19) .eq. 2 .and. &
125 (kpds(5) .eq. 191 .or. kpds(5) .eq. 192))) then ! ensemble forecast
126 ensemble = .true.
127 endif
128 if (ensemble) then ! ensemble forecast
129 call g2_gbytec(cgribin(9), ilast, 0, 24)
130 call pdseup(kens, kprob, xprob, kclust, kmembr, ilast, cgribin(9))
131 if (kens(2) .eq. 1) listsec1(13) = 3
132 if (kens(2) .eq. 2 .OR. kens(2) .eq. 3) listsec1(13) = 4
133 if (kens(2) .eq. 5) listsec1(13) = 5
134 endif
135 ecmwf = .false.
136 if (kpds(1) .eq. 98) ecmwf = .true.
137 if (ecmwf) then ! treat ecmwf data conversion seperately
138 call param_ecmwf_g1_to_g2(kpds(5), kpds(19), listsec0(1), idum, &
139 jdum) ! set discipline
140 else
141 if (ensemble.and.(kpds(5) .eq. 191 .or. kpds(5) .eq. 192).and. &
142 kpds(19) .eq. 2) then
143 !kprob(1) = 61
144 call param_g1_to_g2(kprob(1), kpds(19), listsec0(1), idum, &
145 jdum) ! set discipline
146 else
147 call param_g1_to_g2(kpds(5), kpds(19), listsec0(1), idum, &
148 jdum) ! set discipline
149 endif
150 endif
151 call gribcreate(cgrib, lcgrib, listsec0, listsec1, ierr)
152 if (ierr .ne. 0) then
153 write(6, *) ' ERROR creating new GRIB2 field = ', ierr
154 cycle
155 endif
156 !
157 ! convert grid info
158 call gds2gdt(kgds, igds, igdstmpl, idefnum, ideflist, ierr)
159 if (ierr .ne. 0) then
160 cycle
161 endif
162 if (listsec1(1) .eq. 7) igdstmpl(1) = 6 ! FOR NWS/NCEP
163 if ((listsec1(1) .eq. 7 .and. igds(5) .eq. 20 & ! For Snow Cover Analysis
164 .and. kpds(2) .eq. 25) .and. & ! Polar Stereographic Grid
165 (kpds(5) .eq. 91 .or. kpds(5) .eq. 238)) then
166 igdstmpl(1) = 2
167 end if
168 call addgrid(cgrib, lcgrib, igds, igdstmpl, 200, ideflist, &
169 idefnum, ierr)
170 if (ierr .ne. 0) then
171 write(6, *) ' ERROR adding GRIB2 grid = ', ierr
172 cycle
173 endif
174
175 ! set PDS Template
176 if (ensemble) then ! ensemble forecast
177 call pds2pdtens(kpds, kens, kprob, xprob, kclust, kmembr, &
178 ipdsnum, ipdstmpl, numcoord, coordlist, ierr)
179 else
180 call pds2pdt(kpds, ipdsnum, ipdstmpl, numcoord, coordlist, ierr)
181 endif
182 if (ierr .ne. 0) then
183 cycle
184 endif
185
186 ! set bitmap flag
187 idrstmpl = 0
188 if (btest(kpds(4), 6)) then
189 ibmap = 0
190 !fld = pack(fld, mask = bmp(1:numpts))
191 !itemp = count(bmp(1:numpts))
192 !numpts = itemp
193 !
194 ! convert bitmap to "missing" values, if requested.
195 !
196 if ((usemiss) .AND. (ipack .eq. 2 .OR. ipack .eq. 31 .OR. &
197 ipack .eq. 32)) then
198 ibmap = 255
199 rmiss = minval(fld(1:numpts))
200 if (rmiss .lt. -9999.0) then
201 rmiss = rmiss*10.0
202 else
203 rmiss = -9999.0
204 endif
205 do i = 1, numpts
206 if (.NOT. bmp(i)) then
207 fld(i) = rmiss
208 bmp(i) = .true.
209 endif
210 enddo
211 idrstmpl(7) = imiss ! Missing value management
212 call mkieee(rmiss, idrstmpl(8), 1)
213 endif
214 else
215 ibmap = 255
216 idrstmpl(7) = 0 ! No missing values
217 endif
218
219 ! Set DRT info (packing info)
220 if (ipack .eq. 0) then
221 idrsnum = 0
222 elseif (ipack .eq. 2) then
223 idrsnum = 2
224 idrstmpl(6) = 1 ! general group split
225 elseif (ipack .eq. 31 .OR. ipack .eq. 32) then
226 idrsnum = ipack/10
227 idrstmpl(6) = 1 ! general group split
228 idrstmpl(17) = mod(ipack, 10) ! order of s.d.
229 elseif (ipack .eq. 40 .OR. ipack .eq. 41 .OR. &
230 ipack .eq. 40000 .OR. ipack .eq. 40010) then
231 idrsnum = ipack
232 idrstmpl(6) = 0
233 idrstmpl(7) = 255
234 !idrstmpl(6) = 1
235 !idrstmpl(7) = 15
236 else
237 idrsnum = 3
238 idrstmpl(17) = 1 ! order of s.d.
239 idrstmpl(6) = 1 ! general group split
240 if (kpds(5) .eq. 61) idrsnum = 2
241 endif
242 idrstmpl(2) = kptr(19) ! binary scale
243 idrstmpl(3) = kpds(22) ! decimal scale
244 !idrstmpl(2) = -4 ! binary scale
245 !idrstmpl(3) = 0 ! decimal scale
246 call addfield(cgrib, lcgrib, ipdsnum, ipdstmpl, 200, &
247 coordlist, numcoord, idrsnum, idrstmpl, 200, &
248 fld, numpts, ibmap, bmp, ierr)
249 ! print *, 'done with addfield'
250 if (ierr .ne. 0) then
251 write(6, *) ' ERROR adding GRIB2 field = ', ierr
252 cycle
253 endif
254
255 if ((kpds(5) .eq. 33) .AND. uvvect) then
256 if (.not.allocated(fldv)) allocate(fldv(maxpts))
257 if (.not.allocated(bmpv)) allocate(bmpv(maxpts))
258 lgds = kgds
259 lens = kens
260 lpds = kpds
261 lpds(22) = -1
262 lpds(5) = 34
263 jsrch = 0
264 CALL getgbe(ifl1, ifli1, maxpts, jsrch, lpds, lgds, lens, numptso, &
265 jsrch, kpds, kgds, kens, bmpv, fldv, icnd)
266 if (icnd .ne. 0) then
267 write(6, *) ' ERROR READING/UNPACKING GRIB1 V = ', icnd
268 exit
269 endif
270 ipdstmplv = ipdstmpl
271 if (ecmwf) then ! treat ecmwf data conversion seperately
272 ! print *, ' param_ecmwf call 2'
273 call param_ecmwf_g1_to_g2(kpds(5), kpds(19), idum, &
274 ipdstmplv(1), ipdstmplv(2))
275 ! print *, ' done with call 2'
276 else
277 call param_g1_to_g2(kpds(5), kpds(19), idum, ipdstmplv(1), &
278 ipdstmplv(2))
279 endif
280 ! set bitmap flag
281 idrstmplv = 0
282 if (btest(kpds(4), 6)) then
283 !fldv = pack(fldv, mask = bmpv(1:numpts))
284 if (any(bmp(1:igds(2)) .NEQV. bmpv(1:igds(2)))) then
285 !print *, 'SAGT: BITMAP different'
286 ibmap = 0
287 ! convert bitmap to "missing" values, if requested.
288 if ((usemiss) .AND. (ipack .eq. 2 .OR. ipack .eq. 31 .OR. &
289 ipack .eq. 32)) then
290 ibmap = 255
291 rmiss = minval(fldv(1:numpts))
292 if (rmiss .lt. -9999.0) then
293 rmiss = rmiss*10.0
294 else
295 rmiss = -9999.0
296 endif
297 do i = 1, numpts
298 if (.NOT. bmpv(i)) then
299 fldv(i) = rmiss
300 bmpv(i) = .true.
301 endif
302 enddo
303 idrstmplv(7) = imiss ! Missing values management
304 call mkieee(rmiss, idrstmplv(8), 1)
305 endif
306 else
307 !print *, 'SAGT: BITMAP SAME'
308 ibmap = 254
309 endif
310 else
311 ibmap = 255
312 idrstmplv(7) = 0 ! No missing values
313 endif
314 ! Set DRT info (packing info)
315 if (ipack .eq. 0) then
316 idrsnum = 0
317 elseif (ipack .eq. 2) then
318 idrsnum = 2
319 idrstmplv(6) = 1 ! general group split
320 elseif (ipack .eq. 31 .OR. ipack .eq. 32) then
321 idrsnum = ipack/10
322 idrstmplv(6) = 1 ! general group split
323 idrstmplv(17) = mod(ipack, 10) ! order of s.d.
324 elseif (ipack .eq. 40 .OR. ipack .eq. 41 .OR. &
325 ipack .eq. 40000 .OR. ipack .eq. 40010) then
326 idrsnum = ipack
327 idrstmplv(6) = 0
328 idrstmplv(7) = 255
329 !idrstmplv(6) = 1
330 !idrstmplv(7) = 15
331 else
332 idrsnum = 3
333 idrstmplv(17) = 1 ! order of s.d.
334 idrstmplv(6) = 1 ! general group split
335 if (kpds(5) .eq. 61) idrsnum = 2
336 endif
337 idrstmplv(2) = kptr(19) ! binary scale
338 idrstmplv(3) = kpds(22) ! decimal scale
339 !idrstmplv(2) = -4 ! binary scale
340 !idrstmplv(3) = 0 ! decimal scale
341 call addfield(cgrib, lcgrib, ipdsnum, ipdstmplv, 200, &
342 coordlist, numcoord, idrsnum, idrstmplv, 200, &
343 fldv, numpts, ibmap, bmpv, ierr)
344 if (ierr .ne. 0) then
345 write(6, *) ' ERROR adding second GRIB2 field = ', ierr
346 cycle
347 endif
348 endif
349 ! End GRIB2 field
350 call gribend(cgrib, lcgrib, lengrib, ierr)
351 if (ierr .ne. 0) then
352 write(6, *) ' ERROR ending new GRIB2 message = ', ierr
353 cycle
354 endif
355 ! print *, ' writing ', lengrib, ' bytes...'
356 call wryte(ifl2, lengrib, cgrib)
357
358 enddo
359
360 if (allocated(cgribin)) deallocate(cgribin)
361 if (allocated(cgrib)) deallocate(cgrib)
362 if (allocated(fld)) deallocate(fld)
363 if (allocated(fldv)) deallocate(fldv)
364 if (allocated(coordlist)) deallocate(coordlist)
365 if (allocated(bmp)) deallocate(bmp)
366 if (allocated(bmpv)) deallocate(bmpv)
367
368 return
369end subroutine cnv12
370
406subroutine gds2gdt(kgds,igds,igdstmpl,idefnum,ideflist,iret)
407
408 integer,intent(in) :: kgds(*)
409 integer,intent(out) :: igds(*),igdstmpl(*),ideflist(*)
410 integer,intent(out) :: idefnum,iret
411
412 iret=0
413 if (kgds(1).eq.0) then ! Lat/Lon grid
414 idefnum=0
415 igds(1)=0 ! grid def specfied in template
416 igds(2)=kgds(2)*kgds(3) ! num of grid points
417 igds(3)=0 ! octets for further grid definition
418 igds(4)=0 ! interpretation of optional list
419 igds(5)=0 ! Grid Definition Template number
420 if (btest(kgds(6),6)) then ! shape of Earth
421 igdstmpl(1)=2
422 else
423 igdstmpl(1)=0
424 endif
425 igdstmpl(2)=0
426 igdstmpl(3)=0
427 igdstmpl(4)=0
428 igdstmpl(5)=0
429 igdstmpl(6)=0
430 igdstmpl(7)=0
431 igdstmpl(8)=kgds(2) !Ni
432 igdstmpl(9)=kgds(3) !Nj
433 igdstmpl(10)=0
434 igdstmpl(11)=0
435 igdstmpl(12)=kgds(4)*1000 ! Lat of 1st grid point
436 if (kgds(5).lt.0) then ! Lon of 1st grid point
437 igdstmpl(13)=(360000+kgds(5))*1000 ! convert W to E
438 else
439 igdstmpl(13)=kgds(5)*1000
440 endif
441 igdstmpl(14)=0 ! Resolution and Component flags
442 if (btest(kgds(6),7)) igdstmpl(14)=48
443 if (btest(kgds(6),3)) igdstmpl(14)=igdstmpl(14)+8
444 igdstmpl(15)=kgds(7)*1000 ! Lat of last grid point
445 if (kgds(8).lt.0) then ! Lon of last grid point
446 igdstmpl(16)=(360000+kgds(8))*1000 ! convert W to E
447 else
448 igdstmpl(16)=kgds(8)*1000
449 endif
450 igdstmpl(17)=kgds(9)*1000 ! Di
451 igdstmpl(18)=kgds(10)*1000 ! Dj
452 igdstmpl(19)=kgds(11) ! Scanning mode
453 if (kgds(20).ne.255) then ! irregular grid (eg WAFS)
454 igds(2)=kgds(21) ! num of grid points
455 !idefnum=kgds(19)
456 if (kgds(2).eq.65535) idefnum=kgds(3)
457 if (kgds(3).eq.65535) idefnum=kgds(2)
458 imax=0
459 do j=1,idefnum
460 ideflist(j)=kgds(21+j)
461 if (ideflist(j).gt.imax) imax=ideflist(j)
462 enddo
463 igds(3)=1 ! octets for further grid definition
464 if (imax.gt.255) igds(3)=2
465 if (imax.gt.65535) igds(3)=3
466 if (imax.gt.16777215) igds(3)=4
467 igds(4)=1 ! interpretation of optional list
468 igdstmpl(8)=-1
469 igdstmpl(17)=-1
470 endif
471 elseif (kgds(1).eq.1) then ! Mercator grid
472 idefnum=0
473 igds(1)=0 ! grid def specfied in template
474 igds(2)=kgds(2)*kgds(3) ! num of grid points
475 igds(3)=0 ! octets for further grid definition
476 igds(4)=0 ! interpretation of optional list
477 igds(5)=10 ! Grid Definition Template number
478 if (btest(kgds(6),6)) then ! shape of Earth
479 igdstmpl(1)=2
480 else
481 igdstmpl(1)=0
482 endif
483 igdstmpl(2)=0
484 igdstmpl(3)=0
485 igdstmpl(4)=0
486 igdstmpl(5)=0
487 igdstmpl(6)=0
488 igdstmpl(7)=0
489 igdstmpl(8)=kgds(2) ! Ni
490 igdstmpl(9)=kgds(3) ! Nj
491 igdstmpl(10)=kgds(4)*1000 ! Lat of 1st grid point
492 if (kgds(5).lt.0) then ! Lon of 1st grid point
493 igdstmpl(11)=(360000+kgds(5))*1000 ! convert W to E
494 else
495 igdstmpl(11)=kgds(5)*1000
496 endif
497 igdstmpl(12)=0 ! Resolution and Component flags
498 if (btest(kgds(6),7)) igdstmpl(12)=48
499 if (btest(kgds(6),3)) igdstmpl(12)=igdstmpl(12)+8
500 igdstmpl(13)=kgds(9)*1000 ! Lat intersects earth
501 igdstmpl(14)=kgds(7)*1000 ! Lat of last grid point
502 if (kgds(8).lt.0) then ! Lon of last grid point
503 igdstmpl(15)=(360000+kgds(8))*1000 ! convert W to E
504 else
505 igdstmpl(15)=kgds(8)*1000
506 endif
507 igdstmpl(16)=kgds(11) ! Scanning mode
508 igdstmpl(17)=0 ! Orientation of grid
509 igdstmpl(18)=kgds(12)*1000 ! Di
510 igdstmpl(19)=kgds(13)*1000 ! Dj
511 elseif (kgds(1).eq.3) then ! Lambert Conformal Grid
512 idefnum=0
513 igds(1)=0 ! grid def specfied in template
514 igds(2)=kgds(2)*kgds(3) ! num of grid points
515 igds(3)=0 ! octets for further grid definition
516 igds(4)=0 ! interpretation of optional list
517 igds(5)=30 ! Grid Definition Template number
518 if (btest(kgds(6),6)) then ! shape of Earth
519 igdstmpl(1)=2
520 else
521 igdstmpl(1)=0
522 endif
523 igdstmpl(2)=0
524 igdstmpl(3)=0
525 igdstmpl(4)=0
526 igdstmpl(5)=0
527 igdstmpl(6)=0
528 igdstmpl(7)=0
529 igdstmpl(8)=kgds(2) ! Nx
530 igdstmpl(9)=kgds(3) ! Ny
531 igdstmpl(10)=kgds(4)*1000 ! Lat of 1st grid point
532 if (kgds(5).lt.0) then ! Lon of 1st grid point
533 igdstmpl(11)=(360000+kgds(5))*1000 ! convert W to E
534 else
535 igdstmpl(11)=kgds(5)*1000
536 endif
537 igdstmpl(12)=0 ! Resolution and Component flags
538 if (btest(kgds(6),7)) igdstmpl(12)=48
539 if (btest(kgds(6),3)) igdstmpl(12)=igdstmpl(12)+8
540 igdstmpl(13)=kgds(12)*1000 ! Lat where Dx and Dy specified
541 if (kgds(7).lt.0) then ! Lon of orientation
542 igdstmpl(14)=(360000+kgds(7))*1000 ! convert W to E
543 else
544 igdstmpl(14)=kgds(7)*1000
545 endif
546 igdstmpl(15)=kgds(8)*1000 ! Dx
547 igdstmpl(16)=kgds(9)*1000 ! Dy
548 igdstmpl(17)=kgds(10) ! Projection Center Flag
549 igdstmpl(18)=kgds(11) ! Scanning mode
550 igdstmpl(19)=kgds(12)*1000 ! Latin 1
551 igdstmpl(20)=kgds(13)*1000 ! Latin 2
552 igdstmpl(21)=kgds(14)*1000 ! Lat of S. Pole of projection
553 if (kgds(15).lt.0) then ! Lon of S. Pole of projection
554 igdstmpl(22)=(360000+kgds(15))*1000 ! convert W to E
555 else
556 igdstmpl(22)=kgds(15)*1000
557 endif
558 elseif (kgds(1).eq.4) then ! Gaussian Lat/Lon grid
559 idefnum=0
560 igds(1)=0 ! grid def specfied in template
561 igds(2)=kgds(2)*kgds(3) ! num of grid points
562 igds(3)=0 ! octets for further grid definition
563 igds(4)=0 ! interpretation of optional list
564 igds(5)=40 ! Grid Definition Template number
565 if (btest(kgds(6),6)) then ! shape of Earth
566 igdstmpl(1)=2
567 else
568 igdstmpl(1)=0
569 endif
570 igdstmpl(2)=0
571 igdstmpl(3)=0
572 igdstmpl(4)=0
573 igdstmpl(5)=0
574 igdstmpl(6)=0
575 igdstmpl(7)=0
576 igdstmpl(8)=kgds(2) !Ni
577 igdstmpl(9)=kgds(3) !Nj
578 igdstmpl(10)=0
579 igdstmpl(11)=0
580 igdstmpl(12)=kgds(4)*1000 ! Lat of 1st grid point
581 if (kgds(5).lt.0) then ! Lon of 1st grid point
582 igdstmpl(13)=(360000+kgds(5))*1000 ! convert W to E
583 else
584 igdstmpl(13)=kgds(5)*1000
585 endif
586 igdstmpl(14)=0 ! Resolution and Component flags
587 if (btest(kgds(6),7)) igdstmpl(14)=48
588 if (btest(kgds(6),3)) igdstmpl(14)=igdstmpl(14)+8
589 igdstmpl(15)=kgds(7)*1000 ! Lat of last grid point
590 if (kgds(8).lt.0) then ! Lon of last grid point
591 igdstmpl(16)=(360000+kgds(8))*1000 ! convert W to E
592 else
593 igdstmpl(16)=kgds(8)*1000
594 endif
595 igdstmpl(17)=kgds(9)*1000 ! Di
596 igdstmpl(18)=kgds(10) ! D - Number of parallels
597 igdstmpl(19)=kgds(11) ! Scanning mode
598 elseif (kgds(1).eq.5) then ! Polar Stereographic Grid
599 idefnum=0
600 igds(1)=0 ! grid def specfied in template
601 igds(2)=kgds(2)*kgds(3) ! num of grid points
602 igds(3)=0 ! octets for further grid definition
603 igds(4)=0 ! interpretation of optional list
604 igds(5)=20 ! Grid Definition Template number
605 if (btest(kgds(6),6)) then ! shape of Earth
606 igdstmpl(1)=2
607 else
608 igdstmpl(1)=0
609 endif
610 igdstmpl(2)=0
611 igdstmpl(3)=0
612 igdstmpl(4)=0
613 igdstmpl(5)=0
614 igdstmpl(6)=0
615 igdstmpl(7)=0
616 igdstmpl(8)=kgds(2) ! Nx
617 igdstmpl(9)=kgds(3) ! Ny
618 igdstmpl(10)=kgds(4)*1000 ! Lat of 1st grid point
619 if (kgds(5).lt.0) then ! Lon of 1st grid point
620 igdstmpl(11)=(360000+kgds(5))*1000 ! convert W to E
621 else
622 igdstmpl(11)=kgds(5)*1000
623 endif
624 igdstmpl(12)=0 ! Resolution and Component flags
625 if (btest(kgds(6),7)) igdstmpl(12)=48
626 if (btest(kgds(6),3)) igdstmpl(12)=igdstmpl(12)+8
627 igdstmpl(13)=60000000 ! Lat where Dx and Dy specified
628 if (btest(kgds(10),7)) igdstmpl(13)=-60000000
629 if (kgds(7).lt.0) then ! Lon of orientation
630 igdstmpl(14)=(360000+kgds(7))*1000 ! convert W to E
631 else
632 igdstmpl(14)=kgds(7)*1000
633 endif
634 igdstmpl(15)=kgds(8)*1000 ! Dx
635 igdstmpl(16)=kgds(9)*1000 ! Dy
636 igdstmpl(17)=kgds(10) ! Projection Center Flag
637 igdstmpl(18)=kgds(11) ! Scanning mode
638 elseif (kgds(1).eq.204) then ! Curivilinear Orthogonal Grid (Used by RTOFS)
639 idefnum=0
640 igds(1)=0 ! grid def specfied in template
641 igds(2)=kgds(2)*kgds(3) ! num of grid points
642 igds(3)=0 ! octets for further grid definition
643 igds(4)=0 ! interpretation of optional list
644 igds(5)=204 ! Grid Definition Template number
645 if (btest(kgds(6),6)) then ! shape of Earth
646 igdstmpl(1)=2
647 else
648 igdstmpl(1)=0
649 endif
650 igdstmpl(2)=0
651 igdstmpl(3)=0
652 igdstmpl(4)=0
653 igdstmpl(5)=0
654 igdstmpl(6)=0
655 igdstmpl(7)=0
656 igdstmpl(8)=kgds(2) !Ni - No of points along x-grid direction
657 igdstmpl(9)=kgds(3) !Nj - No of points along y-grid direction
658 igdstmpl(10)=0
659 igdstmpl(11)=0
660 igdstmpl(12)=0
661 igdstmpl(13)=0
662 igdstmpl(14)=0 ! Resolution and Component flags
663 if (btest(kgds(6),7)) igdstmpl(14)=48
664 if (btest(kgds(6),3)) igdstmpl(14)=igdstmpl(14)+8
665 igdstmpl(15)=0
666 igdstmpl(16)=0
667 igdstmpl(17)=0
668 igdstmpl(18)=0
669 igdstmpl(19)=kgds(11) ! Scanning mode
670 elseif (kgds(1).eq.203) then ! Rot Lat/Lon grid (Arakawa)
671 idefnum=0
672 igds(1)=0 ! grid def specfied in template
673 igds(2)=kgds(2)*kgds(3) ! num of grid points
674 igds(3)=0 ! octets for further grid definition
675 igds(4)=0 ! interpretation of optional list
676 igds(5)=32768 ! Grid Definition Template number
677 if (btest(kgds(6),6)) then ! shape of Earth
678 igdstmpl(1)=2
679 else
680 igdstmpl(1)=0
681 endif
682 igdstmpl(2)=0
683 igdstmpl(3)=0
684 igdstmpl(4)=0
685 igdstmpl(5)=0
686 igdstmpl(6)=0
687 igdstmpl(7)=0
688 igdstmpl(8)=kgds(2) !Ni
689 igdstmpl(9)=kgds(3) !Nj
690 igdstmpl(10)=0
691 igdstmpl(11)=0
692 igdstmpl(12)=kgds(4)*1000 ! Lat of 1st grid point
693 if (kgds(5).lt.0) then ! Lon of 1st grid point
694 igdstmpl(13)=(360000+kgds(5))*1000 ! convert W to E
695 else
696 igdstmpl(13)=kgds(5)*1000
697 endif
698 igdstmpl(14)=0 ! Resolution and Component flags
699 if (btest(kgds(6),7)) igdstmpl(14)=48
700 if (btest(kgds(6),3)) igdstmpl(14)=igdstmpl(14)+8
701 igdstmpl(15)=kgds(7)*1000 ! Lat of last grid point
702 if (kgds(8).lt.0) then ! Lon of last grid point
703 igdstmpl(16)=(360000+kgds(8))*1000 ! convert W to E
704 else
705 igdstmpl(16)=kgds(8)*1000
706 endif
707 igdstmpl(17)=kgds(9)*1000 ! Di
708 igdstmpl(18)=kgds(10)*1000 ! Dj
709 igdstmpl(19)=kgds(11) ! Scanning mode
710 elseif (kgds(1).eq.205) then ! Rot Lat/Lon for Non-E Stagger grid (Arakawa)
711 idefnum=0
712 igds(1)=0 ! grid def specfied in template
713 igds(2)=kgds(2)*kgds(3) ! num of grid points
714 igds(3)=0 ! octets for further grid definition
715 igds(4)=0 ! interpretation of optional list
716 igds(5)=32769 ! Grid Definition Template number
717 if (btest(kgds(6),6)) then ! shape of Earth
718 igdstmpl(1)=2
719 else
720 igdstmpl(1)=0
721 endif
722 igdstmpl(2)=0
723 igdstmpl(3)=0
724 igdstmpl(4)=0
725 igdstmpl(5)=0
726 igdstmpl(6)=0
727 igdstmpl(7)=0
728 igdstmpl(8)=kgds(2) !Ni
729 igdstmpl(9)=kgds(3) !Nj
730 igdstmpl(10)=0
731 igdstmpl(11)=0
732 igdstmpl(12)=kgds(4)*1000 ! Lat of 1st grid point
733 if (kgds(5).lt.0) then ! Lon of 1st grid point
734 igdstmpl(13)=(360000+kgds(5))*1000 ! convert W to E
735 else
736 igdstmpl(13)=kgds(5)*1000
737 endif
738 igdstmpl(14)=0 ! Resolution and Component flags
739 if (btest(kgds(6),7)) igdstmpl(14)=48
740 if (btest(kgds(6),3)) igdstmpl(14)=igdstmpl(14)+8
741 igdstmpl(15)=kgds(7)*1000 ! Lat of last grid point
742 if (kgds(8).lt.0) then ! Lon of last grid point
743 igdstmpl(16)=(360000+kgds(8))*1000 ! convert W to E
744 else
745 igdstmpl(16)=kgds(8)*1000
746 endif
747 igdstmpl(17)=kgds(9)*1000 ! Di
748 igdstmpl(18)=kgds(10)*1000 ! Dj
749 igdstmpl(19)=kgds(11) ! Scanning mode
750 igdstmpl(20)=kgds(12)*1000
751 igdstmpl(21)=kgds(13)*1000
752 else
753 print *,'gds2gdt: Unrecognized GRIB1 Grid type = ',kgds(1)
754 iret=1
755 endif
756 return
757end subroutine gds2gdt
758
788subroutine pds2pdt(kpds,ipdsnum,ipdstmpl,numcoord,coordlist, &
789 iret)
790
791 use params
792 use params_ecmwf
793
794 integer,intent(in) :: kpds(*)
795 integer,intent(out) :: ipdstmpl(*)
796 real,intent(out) :: coordlist(*)
797 integer,intent(out) :: ipdsnum,numcoord,iret
798
799 integer :: idat(8),jdat(8)
800 real :: rinc(5)
801 logical :: ecmwf
802
803 iret=0
804 numcoord=0
805 ecmwf=.false.
806
807 if (kpds(1).eq.98) ecmwf=.true.
808 !
809 ! Special check for WAFS products for parameters (Max Icing, TP and CAT)
810 ! to PDT 4.15
811 !
812 if ((kpds(2).eq.96 .AND. kpds(3).eq.45 .AND. &
813 kpds(16).eq.10) .AND. &
814 (kpds(19).eq.140 .AND. (kpds(5).ge.168 .AND. &
815 kpds(5).le.173))) then
816 ipdsnum=15
817 ! get GRIB2 parameter category and number from GRIB1
818 ! parameter number
819 if (ecmwf) then ! treat ecmwf data conversion seperately
820 call param_ecmwf_g1_to_g2(kpds(5),kpds(19),idum, &
821 ipdstmpl(1),ipdstmpl(2))
822 else
823 call param_g1_to_g2(kpds(5),kpds(19),idum,ipdstmpl(1), &
824 ipdstmpl(2))
825 endif
826 ipdstmpl(3)=2
827 ipdstmpl(4)=0
828 ipdstmpl(5)=kpds(2)
829 ipdstmpl(6)=0
830 ipdstmpl(7)=0
831 ipdstmpl(8)=kpds(13)
832 if (kpds(13).eq.254) ipdstmpl(8)=13
833 ipdstmpl(9)=kpds(14)
834 call cnvlevel(kpds(6),kpds(7),ipdstmpl)
835 if (kpds(5).eq.168.or.kpds(5).eq.170.or. & ! statistical process WAFS-ICAO
836 kpds(5).eq.172) ipdstmpl(16)=0 ! for Mean Icing, CT, CAT
837 if (kpds(5).eq.169.or.kpds(5).eq.171.or. & ! statistical process WAFS-ICAO
838 kpds(5).eq.173) ipdstmpl(16)=2 ! for MAX Icing, CT, CAT
839 ipdstmpl(17)=3 ! Neighbor interpolation
840 ! Output values is set to nearest input values
841 ipdstmpl(18)=1 ! Number of data point (grid 45)
842 elseif (kpds(16).eq.0.or.kpds(16).eq.1.or.kpds(16).eq.10) then
843 ipdsnum=0
844 ! get GRIB2 parameter category and number from GRIB1
845 ! parameter number
846 if (ecmwf) then ! treat ecmwf data conversion seperately
847 call param_ecmwf_g1_to_g2(kpds(5),kpds(19),idum, &
848 ipdstmpl(1),ipdstmpl(2))
849 else
850 call param_g1_to_g2(kpds(5),kpds(19),idum,ipdstmpl(1), &
851 ipdstmpl(2))
852 endif
853 if (kpds(16).eq.1) then
854 ipdstmpl(3)=0
855 else
856 ipdstmpl(3)=2
857 endif
858 ipdstmpl(4)=0
859 ipdstmpl(5)=kpds(2)
860 ipdstmpl(6)=0
861 ipdstmpl(7)=0
862 ipdstmpl(8)=kpds(13)
863 if (kpds(13).eq.254) ipdstmpl(8)=13
864 !if (kpds(16).eq.10) then
865 ! ipdstmpl(9)=(kpds(14)*256)+kpds(15)
866 !else
867 ipdstmpl(9)=kpds(14)
868 !endif
869 call cnvlevel(kpds(6),kpds(7),ipdstmpl)
870 if (kpds(2).eq.96 .AND. kpds(3).eq.45 .AND. &
871 kpds(16).eq.10) then
872 if (kpds(5).eq.174) ipdstmpl(10) = 10
873 if (kpds(5).eq.179) ipdstmpl(10) = 11
874 if (kpds(5).eq.180) ipdstmpl(10) = 12
875 end if
876 elseif (kpds(16).ge.2.AND.kpds(16).le.5) then
877 ipdsnum=8
878 ! get GRIB2 parameter category and number from GRIB1
879 ! parameter number
880 if (ecmwf) then ! treat ecmwf data conversion seperately
881 call param_ecmwf_g1_to_g2(kpds(5),kpds(19),idum, &
882 ipdstmpl(1),ipdstmpl(2))
883 else
884 call param_g1_to_g2(kpds(5),kpds(19),idum,ipdstmpl(1), &
885 ipdstmpl(2))
886 endif
887 ipdstmpl(3)=2
888 ipdstmpl(4)=0
889 ipdstmpl(5)=kpds(2)
890 ipdstmpl(6)=0
891 ipdstmpl(7)=0
892 ipdstmpl(8)=kpds(13)
893 if (kpds(13).eq.254) ipdstmpl(8)=13
894 ipdstmpl(9)=kpds(14)
895 call cnvlevel(kpds(6),kpds(7),ipdstmpl)
896 ! calculate ending time using initial ref-time, idat,
897 ! and increment rinc.
898 idat=0
899 idat(1)=((kpds(21)-1)*100)+kpds(8)
900 idat(2)=kpds(9)
901 idat(3)=kpds(10)
902 idat(4)=-500 ! EST
903 idat(5)=kpds(11)
904 idat(6)=kpds(12)
905 rinc=0.0
906 if ( ipdstmpl(8).eq.0 ) then
907 rinc(3)=kpds(15)
908 elseif ( ipdstmpl(8).eq.1 ) then
909 rinc(2)=kpds(15)
910 elseif ( ipdstmpl(8).eq.2 ) then
911 rinc(1)=kpds(15)
912 elseif ( ipdstmpl(8).eq.10 ) then
913 rinc(2)=kpds(15) * 3
914 elseif ( ipdstmpl(8).eq.11 ) then
915 rinc(2)=kpds(15) * 6
916 elseif ( ipdstmpl(8).eq.12 ) then
917 rinc(2)=kpds(15) * 12
918 elseif ( ipdstmpl(8).eq.13 ) then
919 rinc(4)=kpds(15)
920 endif
921 call w3movdat(rinc,idat,jdat) ! calculate end date/time
922 ipdstmpl(16)=jdat(1) ! year of end time
923 ipdstmpl(17)=jdat(2) ! month of end time
924 ipdstmpl(18)=jdat(3) ! day of end time
925 ipdstmpl(19)=jdat(5) ! hour of end time
926 ipdstmpl(20)=jdat(6) ! minute of end time
927 ipdstmpl(21)=jdat(7) ! second of end time
928 ipdstmpl(22)=1 ! # of time ranges
929 ipdstmpl(23)=kpds(20) ! # of values missing
930 if (kpds(16).eq.2) then ! statistical process
931 ipdstmpl(24)=255
932 elseif (kpds(16).eq.3) then
933 ipdstmpl(24)=0
934 elseif (kpds(16).eq.4) then
935 ipdstmpl(24)=1
936 elseif (kpds(16).eq.5) then
937 ipdstmpl(24)=4
938 endif
939 ipdstmpl(25)=2
940 if (kpds(19).eq.129 .AND. &
941 (kpds(5).eq.235 .or. kpds(5).eq.236 .or. &
942 kpds(5).eq.237 .or. kpds(5).eq.238 .or. &
943 kpds(5).eq.239 .or. kpds(5).eq.253 .or. &
944 kpds(5).eq.254 )) then
945 ipdstmpl(24)=2
946 endif
947 ipdstmpl(26)=kpds(13)
948 if (kpds(13).eq.254) ipdstmpl(26)=13
949 ipdstmpl(27)=kpds(15)-kpds(14)
950 ipdstmpl(28)=255
951 ipdstmpl(29)=0
952 elseif (kpds(16).eq.7) then
953 ipdsnum=8
954 ! get GRIB2 parameter category and number from GRIB1
955 ! parameter number
956 if (ecmwf) then ! treat ecmwf data conversion seperately
957 call param_ecmwf_g1_to_g2(kpds(5),kpds(19),idum, &
958 ipdstmpl(1),ipdstmpl(2))
959 else
960 call param_g1_to_g2(kpds(5),kpds(19),idum,ipdstmpl(1), &
961 ipdstmpl(2))
962 endif
963 ipdstmpl(3)=2
964 ipdstmpl(4)=0
965 ipdstmpl(5)=kpds(2)
966 ipdstmpl(6)=0
967 ipdstmpl(7)=0
968 ipdstmpl(8)=kpds(13)
969 if (kpds(13).eq.254) ipdstmpl(8)=13
970 ipdstmpl(9)= - kpds(14)
971 call cnvlevel(kpds(6),kpds(7),ipdstmpl)
972 ! calculate ending time using initial ref-time, idat,
973 ! and increment rinc.
974 idat=0
975 idat(1)=((kpds(21)-1)*100)+kpds(8)
976 idat(2)=kpds(9)
977 idat(3)=kpds(10)
978 idat(4)=-500 ! EST
979 idat(5)=kpds(11)
980 idat(6)=kpds(12)
981 rinc=0.0
982 if ( ipdstmpl(8).eq.0 ) then
983 rinc(3)=kpds(15)
984 elseif ( ipdstmpl(8).eq.1 ) then
985 rinc(2)=kpds(15)
986 elseif ( ipdstmpl(8).eq.2 ) then
987 rinc(1)=kpds(15)
988 elseif ( ipdstmpl(8).eq.10 ) then
989 rinc(2)=kpds(15) * 3
990 elseif ( ipdstmpl(8).eq.11 ) then
991 rinc(2)=kpds(15) * 6
992 elseif ( ipdstmpl(8).eq.12 ) then
993 rinc(2)=kpds(15) * 12
994 elseif ( ipdstmpl(8).eq.13 ) then
995 rinc(4)=kpds(15)
996 endif
997 call w3movdat(rinc,idat,jdat) ! calculate end date/time
998 ipdstmpl(16)=jdat(1) ! year of end time
999 ipdstmpl(17)=jdat(2) ! month of end time
1000 ipdstmpl(18)=jdat(3) ! day of end time
1001 ipdstmpl(19)=jdat(5) ! hour of end time
1002 ipdstmpl(20)=jdat(6) ! minute of end time
1003 ipdstmpl(21)=jdat(7) ! second of end time
1004 ipdstmpl(22)=1 ! # of time ranges
1005 ipdstmpl(23)=kpds(20) ! # of values missing
1006 ipdstmpl(24)=0
1007 ipdstmpl(25)=2
1008 ipdstmpl(26)=kpds(13)
1009 if (kpds(13).eq.254) ipdstmpl(26)=13
1010 ipdstmpl(27)=kpds(15) + kpds(14)
1011 ipdstmpl(28)=255
1012 ipdstmpl(29)=0
1013 elseif (kpds(16).eq.51) then
1014 ipdsnum=8
1015 ! get GRIB2 parameter category and number from GRIB1
1016 ! parameter number
1017 if (ecmwf) then ! treat ecmwf data conversion seperately
1018 call param_ecmwf_g1_to_g2(kpds(5),kpds(19),idum, &
1019 ipdstmpl(1),ipdstmpl(2))
1020 else
1021 call param_g1_to_g2(kpds(5),kpds(19),idum,ipdstmpl(1), &
1022 ipdstmpl(2))
1023 endif
1024 ipdstmpl(3)=2
1025 ipdstmpl(4)=0
1026 ipdstmpl(5)=kpds(2)
1027 ipdstmpl(6)=0
1028 ipdstmpl(7)=0
1029 ipdstmpl(8)=kpds(13)
1030 if (kpds(13).eq.254) ipdstmpl(8)=13
1031 ipdstmpl(9)=kpds(14)
1032 call cnvlevel(kpds(6),kpds(7),ipdstmpl)
1033 ! calculate ending time using initial ref-time, idat,
1034 ! and increment rinc.
1035 idat=0
1036 idat(1)=((kpds(21)-1)*100)+kpds(8)
1037 idat(2)=kpds(9)
1038 idat(3)=kpds(10)
1039 idat(4)=-500 ! EST
1040 idat(5)=kpds(11)
1041 idat(6)=kpds(12)
1042 rinc=0.0
1043 if ( ipdstmpl(8).eq.0 ) then
1044 rinc(3)=kpds(15)
1045 elseif ( ipdstmpl(8).eq.1 ) then
1046 rinc(2)=kpds(15)
1047 elseif ( ipdstmpl(8).eq.2 ) then
1048 rinc(1)=kpds(15)
1049 elseif ( ipdstmpl(8).eq.10 ) then
1050 rinc(2)=kpds(15) * 3
1051 elseif ( ipdstmpl(8).eq.11 ) then
1052 rinc(2)=kpds(15) * 6
1053 elseif ( ipdstmpl(8).eq.12 ) then
1054 rinc(2)=kpds(15) * 12
1055 elseif ( ipdstmpl(8).eq.13 ) then
1056 rinc(4)=kpds(15)
1057 endif
1058 call w3movdat(rinc,idat,jdat) ! calculate end date/time
1059 ipdstmpl(16)=jdat(1) ! year of end time
1060 ipdstmpl(17)=jdat(2) ! month of end time
1061 ipdstmpl(18)=jdat(3) ! day of end time
1062 ipdstmpl(19)=jdat(5) ! hour of end time
1063 ipdstmpl(20)=jdat(6) ! minute of end time
1064 ipdstmpl(21)=jdat(7) ! second of end time
1065 ipdstmpl(22)=1 ! # of time ranges
1066 ipdstmpl(23)=kpds(20) ! # of values missing
1067 ipdstmpl(24)=51 ! Climatological Mean
1068 ipdstmpl(25)=2
1069 ipdstmpl(26)=kpds(13)
1070 if (kpds(13).eq.254) ipdstmpl(26)=13
1071 ipdstmpl(27)=kpds(15)-kpds(14)
1072 ipdstmpl(28)=255
1073 ipdstmpl(29)=0
1074 else
1075 print *,' Unrecognized Time Range Indicator = ',kpds(16)
1076 print *,'pds2pdt: Couldn:t construct PDS Template '
1077 iret=1
1078 endif
1079
1080 return
1081end subroutine pds2pdt
1082
1098subroutine cnvlevel(ltype,lval,ipdstmpl)
1099
1100 integer,intent(in) :: ltype,lval
1101 integer,intent(inout) :: ipdstmpl(*)
1102
1103 ipdstmpl(10)=ltype
1104 ipdstmpl(11)=0
1105 ipdstmpl(12)=0
1106 ipdstmpl(13)=255
1107 ipdstmpl(14)=0
1108 ipdstmpl(15)=0
1109
1110 if (ltype.eq.100) then
1111 ipdstmpl(12)=lval*100
1112 elseif (ltype.eq.101) then
1113 ipdstmpl(10)=100
1114 ipdstmpl(12)=(lval/256)*1000
1115 ipdstmpl(13)=100
1116 ipdstmpl(15)=mod(lval,256)*1000
1117 elseif (ltype.eq.102) then
1118 ipdstmpl(10)=101
1119 elseif (ltype.eq.103) then
1120 ipdstmpl(10)=102
1121 ipdstmpl(12)=lval
1122 elseif (ltype.eq.104) then
1123 ipdstmpl(10)=102
1124 ipdstmpl(12)=lval/256
1125 ipdstmpl(13)=102
1126 ipdstmpl(15)=mod(lval,256)
1127 elseif (ltype.eq.105) then
1128 ipdstmpl(10)=103
1129 ipdstmpl(12)=lval
1130 elseif (ltype.eq.106) then
1131 ipdstmpl(10)=103
1132 ipdstmpl(12)=(lval/256)*100
1133 ipdstmpl(13)=103
1134 ipdstmpl(15)=mod(lval,256)*100
1135 elseif (ltype.eq.107) then
1136 ipdstmpl(10)=104
1137 ipdstmpl(11)=4
1138 ipdstmpl(12)=lval
1139 elseif (ltype.eq.108) then
1140 ipdstmpl(10)=104
1141 ipdstmpl(11)=2
1142 ipdstmpl(12)=lval/256
1143 ipdstmpl(13)=104
1144 ipdstmpl(14)=2
1145 ipdstmpl(15)=mod(lval,256)
1146 elseif (ltype.eq.109) then
1147 ipdstmpl(10)=105
1148 ipdstmpl(12)=lval
1149 elseif (ltype.eq.110) then
1150 ipdstmpl(10)=105
1151 ipdstmpl(12)=lval/256
1152 ipdstmpl(13)=105
1153 ipdstmpl(15)=mod(lval,256)
1154 elseif (ltype.eq.111) then
1155 ipdstmpl(10)=106
1156 ipdstmpl(11)=2
1157 ipdstmpl(12)=lval
1158 elseif (ltype.eq.112) then
1159 ipdstmpl(10)=106
1160 ipdstmpl(11)=2
1161 ipdstmpl(12)=lval/256
1162 ipdstmpl(13)=106
1163 ipdstmpl(14)=2
1164 ipdstmpl(15)=mod(lval,256)
1165 elseif (ltype.eq.113) then
1166 ipdstmpl(10)=107
1167 ipdstmpl(12)=lval
1168 elseif (ltype.eq.114) then
1169 ipdstmpl(10)=107
1170 ipdstmpl(12)=475+(lval/256)
1171 ipdstmpl(13)=107
1172 ipdstmpl(15)=475+mod(lval,256)
1173 elseif (ltype.eq.115) then
1174 ipdstmpl(10)=108
1175 ipdstmpl(12)=lval*100
1176 elseif (ltype.eq.116) then
1177 ipdstmpl(10)=108
1178 ipdstmpl(12)=(lval/256)*100
1179 ipdstmpl(13)=108
1180 ipdstmpl(15)=mod(lval,256)*100
1181 elseif (ltype.eq.117) then
1182 ipdstmpl(10)=109
1183 ipdstmpl(11)=9
1184 ipdstmpl(12)=lval
1185 if ( btest(lval,15) ) then
1186 ipdstmpl(12)=-1*mod(lval,32768)
1187 endif
1188 elseif (ltype.eq.119) then
1189 ipdstmpl(10)=111
1190 ipdstmpl(11)=4
1191 ipdstmpl(12)=lval
1192 elseif (ltype.eq.120) then
1193 ipdstmpl(10)=111
1194 ipdstmpl(11)=2
1195 ipdstmpl(12)=lval/256
1196 ipdstmpl(13)=111
1197 ipdstmpl(14)=2
1198 ipdstmpl(15)=mod(lval,256)
1199 elseif (ltype.eq.121) then
1200 ipdstmpl(10)=100
1201 ipdstmpl(12)=(1100+(lval/256))*100
1202 ipdstmpl(13)=100
1203 ipdstmpl(15)=(1100+mod(lval,256))*100
1204 elseif (ltype.eq.125) then
1205 ipdstmpl(10)=103
1206 ipdstmpl(11)=2
1207 ipdstmpl(12)=lval
1208 elseif (ltype.eq.126) then
1209 ipdstmpl(10)=100
1210 ipdstmpl(12)=lval
1211 elseif (ltype.eq.128) then
1212 ipdstmpl(10)=104
1213 ipdstmpl(11)=3
1214 ipdstmpl(12)=1100+(lval/256)
1215 ipdstmpl(13)=104
1216 ipdstmpl(14)=3
1217 ipdstmpl(15)=1100+mod(lval,256)
1218 elseif (ltype.eq.141) then
1219 ipdstmpl(10)=100
1220 ipdstmpl(12)=(lval/256)*100
1221 ipdstmpl(13)=100
1222 ipdstmpl(15)=(1100+mod(lval,256))*100
1223 elseif (ltype.eq.160) then
1224 ipdstmpl(10)=160
1225 ipdstmpl(12)=lval
1226 elseif (ltype.gt.99.AND.ltype.lt.200) then
1227 print *,'cnvlevel: GRIB1 Level ',ltype,' not recognized.'
1228 ipdstmpl(10)=255
1229 elseif (ltype.eq.235) then
1230 ipdstmpl(10)=ltype
1231 ipdstmpl(12)=lval
1232 elseif (ltype.eq.236) then
1233 ipdstmpl(10)=ltype
1234 ipdstmpl(12)=lval/256
1235 ipdstmpl(13)=ltype
1236 ipdstmpl(15)=mod(lval,256)
1237 elseif (ltype.ge.237.AND.ltype.le.239) then
1238 ipdstmpl(10)=ltype
1239 ipdstmpl(12)=lval
1240 endif
1241
1242 return
1243end subroutine cnvlevel
1244
1275subroutine pds2pdtens(kpds,kens,kprob,xprob,kclust,kmember, &
1276 ipdsnum,ipdstmpl,numcoord,coordlist, &
1277 iret)
1278
1279 use params
1280
1281 integer,intent(in) :: kpds(*),kens(*),kprob(*),kclust(*)
1282 integer,intent(in) :: kmember(*)
1283 real,intent(in) :: xprob(*)
1284 integer,intent(out) :: ipdstmpl(*)
1285 real,intent(out) :: coordlist(*)
1286 integer,intent(out) :: ipdsnum,numcoord,iret
1287
1288 integer :: idat(8),jdat(8)
1289 real :: rinc(5)
1290
1291 iret=0
1292 numcoord=0
1293 if (kens(2).eq.1.or.kens(2).eq.2.or.kens(2).eq.3) then
1294 ! individual ensemble fcst...
1295 if (kpds(16).eq.0.or.kpds(16).eq.1.or.kpds(16).eq.10) then
1296 ! At specific point in time...
1297 ipdsnum=1
1298 ! get GRIB2 parameter category and number from GRIB1
1299 ! parameter number
1300 call param_g1_to_g2(kpds(5),kpds(19),idum,ipdstmpl(1), &
1301 ipdstmpl(2))
1302 ipdstmpl(3)=4
1303 ipdstmpl(4)=0
1304 ipdstmpl(5)=kpds(2)
1305 ipdstmpl(6)=0
1306 ipdstmpl(7)=0
1307 ipdstmpl(8)=kpds(13)
1308 if (kpds(13).eq.254) ipdstmpl(8)=13
1309 !if (kpds(16).eq.10) then
1310 ! ipdstmpl(9)=(kpds(14)*256)+kpds(15)
1311 !else
1312 ipdstmpl(9)=kpds(14)
1313 !endif
1314 call cnvlevel(kpds(6),kpds(7),ipdstmpl)
1315 if (kens(2).eq.1) then
1316 ! if (kens(3).eq.1) ipdstmpl(16)=0
1317 ! if (kens(3).eq.2) ipdstmpl(16)=1
1318 ipdstmpl(16)=kens(3)-1
1319 ipdstmpl(17)=0
1320 elseif (kens(2).eq.2) then
1321 ipdstmpl(16)=2
1322 ipdstmpl(17)=kens(3)
1323 elseif (kens(2).eq.3) then
1324 ipdstmpl(16)=3
1325 ipdstmpl(17)=kens(3)
1326 endif
1327 ipdstmpl(18)=10
1328 elseif (kpds(16).ge.2.AND.kpds(16).le.5) then
1329 ! over time range...
1330 ipdsnum=11
1331 ! get GRIB2 parameter category and number from GRIB1
1332 ! parameter number
1333 call param_g1_to_g2(kpds(5),kpds(19),idum,ipdstmpl(1), &
1334 ipdstmpl(2))
1335 ipdstmpl(3)=4
1336 ipdstmpl(4)=0
1337 ipdstmpl(5)=kpds(2)
1338 ipdstmpl(6)=0
1339 ipdstmpl(7)=0
1340 ipdstmpl(8)=kpds(13)
1341 if (kpds(13).eq.254) ipdstmpl(8)=13
1342 ipdstmpl(9)=kpds(14)
1343 call cnvlevel(kpds(6),kpds(7),ipdstmpl)
1344 !ipdstmpl(9)=kpds(15)
1345 if (kens(2).eq.1) then
1346 ! if (kens(3).eq.1) ipdstmpl(16)=0
1347 ! if (kens(3).eq.2) ipdstmpl(16)=1
1348 ipdstmpl(16)=kens(3)-1
1349 ipdstmpl(17)=0
1350 elseif (kens(2).eq.2) then
1351 ipdstmpl(16)=2
1352 ipdstmpl(17)=kens(3)
1353 elseif (kens(2).eq.3) then
1354 ipdstmpl(16)=3
1355 ipdstmpl(17)=kens(3)
1356 endif
1357 ipdstmpl(18)=10
1358 ! calculate ending time using initial ref-time, idat,
1359 ! and increment rinc.
1360 idat=0
1361 idat(1)=((kpds(21)-1)*100)+kpds(8)
1362 idat(2)=kpds(9)
1363 idat(3)=kpds(10)
1364 idat(4)=-500 ! EST
1365 idat(5)=kpds(11)
1366 idat(6)=kpds(12)
1367 rinc=0
1368 if ( ipdstmpl(8).eq.0 ) then
1369 rinc(3)=kpds(15)
1370 elseif ( ipdstmpl(8).eq.1 ) then
1371 rinc(2)=kpds(15)
1372 elseif ( ipdstmpl(8).eq.2 ) then
1373 rinc(1)=kpds(15)
1374 elseif ( ipdstmpl(8).eq.10 ) then
1375 rinc(2)=kpds(15) * 3
1376 elseif ( ipdstmpl(8).eq.11 ) then
1377 rinc(2)=kpds(15) * 6
1378 elseif ( ipdstmpl(8).eq.12 ) then
1379 rinc(2)=kpds(15) * 12
1380 elseif ( ipdstmpl(8).eq.13 ) then
1381 rinc(4)=kpds(15)
1382 endif
1383 call w3movdat(rinc,idat,jdat) ! calculate end date/time
1384 ipdstmpl(19)=jdat(1) ! year of end time
1385 ipdstmpl(20)=jdat(2) ! month of end time
1386 ipdstmpl(21)=jdat(3) ! day of end time
1387 ipdstmpl(22)=jdat(5) ! hour of end time
1388 ipdstmpl(23)=jdat(6) ! minute of end time
1389 ipdstmpl(24)=jdat(7) ! second of end time
1390 ipdstmpl(25)=1
1391 ipdstmpl(26)=0
1392 if (kpds(16).eq.2) then
1393 ipdstmpl(27)=255
1394 if (kpds(5).eq.15) ipdstmpl(27)=2
1395 if (kpds(5).eq.16) ipdstmpl(27)=3
1396 elseif (kpds(16).eq.3) then
1397 ipdstmpl(27)=0
1398 elseif (kpds(16).eq.4) then
1399 ipdstmpl(27)=1
1400 elseif (kpds(16).eq.5) then
1401 ipdstmpl(27)=4
1402 endif
1403 ipdstmpl(28)=2
1404 ipdstmpl(29)=kpds(13)
1405 if (kpds(13).eq.254) ipdstmpl(29)=13
1406 ipdstmpl(30)=kpds(15)-kpds(14)
1407 ipdstmpl(31)=255
1408 ipdstmpl(32)=0
1409 elseif (kpds(16).eq.51) then
1410 ! over time range...
1411 ipdsnum=11
1412 ! get GRIB2 parameter category and number from GRIB1
1413 ! parameter number
1414 call param_g1_to_g2(kpds(5),kpds(19),idum,ipdstmpl(1), &
1415 ipdstmpl(2))
1416 ipdstmpl(3)=4
1417 ipdstmpl(4)=0
1418 ipdstmpl(5)=kpds(2)
1419 ipdstmpl(6)=0
1420 ipdstmpl(7)=0
1421 ipdstmpl(8)=kpds(13)
1422 if (kpds(13).eq.254) ipdstmpl(8)=13
1423 ipdstmpl(9)=kpds(14)
1424 call cnvlevel(kpds(6),kpds(7),ipdstmpl)
1425 !ipdstmpl(9)=kpds(15)
1426 if (kens(2).eq.1) then
1427 ! if (kens(3).eq.1) ipdstmpl(16)=0
1428 ! if (kens(3).eq.2) ipdstmpl(16)=1
1429 ipdstmpl(16)=kens(3)-1
1430 ipdstmpl(17)=0
1431 elseif (kens(2).eq.2) then
1432 ipdstmpl(16)=2
1433 ipdstmpl(17)=kens(3)
1434 elseif (kens(2).eq.3) then
1435 ipdstmpl(16)=3
1436 ipdstmpl(17)=kens(3)
1437 endif
1438 ipdstmpl(18)=10
1439 ! calculate ending time using initial ref-time, idat,
1440 ! and increment rinc.
1441 idat=0
1442 idat(1)=((kpds(21)-1)*100)+kpds(8)
1443 idat(2)=kpds(9)
1444 idat(3)=kpds(10)
1445 idat(4)=-500 ! EST
1446 idat(5)=kpds(11)
1447 idat(6)=kpds(12)
1448 rinc=0
1449 if ( ipdstmpl(8).eq.0 ) then
1450 rinc(3)=kpds(15)
1451 elseif ( ipdstmpl(8).eq.1 ) then
1452 rinc(2)=kpds(15)
1453 elseif ( ipdstmpl(8).eq.2 ) then
1454 rinc(1)=kpds(15)
1455 elseif ( ipdstmpl(8).eq.10 ) then
1456 rinc(2)=kpds(15) * 3
1457 elseif ( ipdstmpl(8).eq.11 ) then
1458 rinc(2)=kpds(15) * 6
1459 elseif ( ipdstmpl(8).eq.12 ) then
1460 rinc(2)=kpds(15) * 12
1461 elseif ( ipdstmpl(8).eq.13 ) then
1462 rinc(4)=kpds(15)
1463 endif
1464 call w3movdat(rinc,idat,jdat) ! calculate end date/time
1465 ipdstmpl(19)=jdat(1) ! year of end time
1466 ipdstmpl(20)=jdat(2) ! month of end time
1467 ipdstmpl(21)=jdat(3) ! day of end time
1468 ipdstmpl(22)=jdat(5) ! hour of end time
1469 ipdstmpl(23)=jdat(6) ! minute of end time
1470 ipdstmpl(24)=jdat(7) ! second of end time
1471 ipdstmpl(25)=1
1472 ipdstmpl(26)=0
1473 ipdstmpl(27)=51
1474 ipdstmpl(28)=2
1475 ipdstmpl(29)=kpds(13)
1476 if (kpds(13).eq.254) ipdstmpl(29)=13
1477 ipdstmpl(30)=kpds(15)-kpds(14)
1478 ipdstmpl(31)=255
1479 ipdstmpl(32)=0
1480 else
1481 print *,' Unrecognized Time Range Ind for ensembles = ', &
1482 kpds(16),kens(2)
1483 print *,'pds2pdtens: Couldn:t construct PDS Template '
1484 iret=1
1485 endif
1486
1487 elseif (kens(2).eq.5) then ! WHOLE or CLUSTERENSEMBLE type
1488 if (kpds(5).eq.191.OR.kpds(5).eq.192) then ! probs
1489 if (kpds(16).eq.0.or.kpds(16).eq.1.or.kpds(16).eq.10) then
1490 ! At specific point in time...
1491 ipdsnum=5
1492 ! get GRIB2 parameter category and number from GRIB1
1493 ! parameter number
1494 call param_g1_to_g2(kprob(1),kpds(19),idum,ipdstmpl(1), &
1495 ipdstmpl(2))
1496 ipdstmpl(3)=5
1497 ipdstmpl(4)=0
1498 ipdstmpl(5)=kpds(2)
1499 ipdstmpl(6)=0
1500 ipdstmpl(7)=0
1501 ipdstmpl(8)=kpds(13)
1502 if (kpds(13).eq.254) ipdstmpl(8)=13
1503 !if (kpds(16).eq.10) then
1504 ! ipdstmpl(9)=(kpds(14)*256)+kpds(15)
1505 !else
1506 ipdstmpl(9)=kpds(14)
1507 !endif
1508 call cnvlevel(kpds(6),kpds(7),ipdstmpl)
1509 ipdstmpl(16)=0 !?
1510 ipdstmpl(17)=kclust(1) !?
1511 ipdstmpl(18)=kprob(2)-1
1512 if (ipdstmpl(18).eq.0.OR.ipdstmpl(18).eq.2) then
1513 ipdstmpl(19)=3
1514 ipdstmpl(20)=nint(xprob(1)*1000.0)
1515 else
1516 ipdstmpl(19)=0
1517 ipdstmpl(20)=0
1518 endif
1519 if (ipdstmpl(18).eq.1.OR.ipdstmpl(18).eq.2) then
1520 ipdstmpl(21)=3
1521 ipdstmpl(22)=nint(xprob(2)*1000.0)
1522 else
1523 ipdstmpl(21)=0
1524 ipdstmpl(22)=0
1525 endif
1526 elseif (kpds(16).ge.2.AND.kpds(16).le.5) then
1527 ! over time range...
1528 ipdsnum=9
1529 ! get GRIB2 parameter category and number from GRIB1
1530 ! parameter number
1531 call param_g1_to_g2(kprob(1),kpds(19),idum,ipdstmpl(1), &
1532 ipdstmpl(2))
1533 ipdstmpl(3)=5
1534 ipdstmpl(4)=0
1535 ipdstmpl(5)=kpds(2)
1536 ipdstmpl(6)=0
1537 ipdstmpl(7)=0
1538 ipdstmpl(8)=kpds(13)
1539 if (kpds(13).eq.254) ipdstmpl(8)=13
1540 ipdstmpl(9)=kpds(14)
1541 call cnvlevel(kpds(6),kpds(7),ipdstmpl)
1542 !ipdstmpl(9)=kpds(15)
1543 ipdstmpl(16)=0 !?
1544 ipdstmpl(17)=kclust(1) !?
1545 ipdstmpl(18)=kprob(2)-1
1546 if (ipdstmpl(18).eq.0.OR.ipdstmpl(18).eq.2) then
1547 ipdstmpl(19)=3
1548 ipdstmpl(20)=nint(xprob(1)*1000.0)
1549 else
1550 ipdstmpl(19)=0
1551 ipdstmpl(20)=0
1552 endif
1553 if (ipdstmpl(18).eq.1.OR.ipdstmpl(18).eq.2) then
1554 ipdstmpl(21)=3
1555 ipdstmpl(22)=nint(xprob(2)*1000.0)
1556 else
1557 ipdstmpl(21)=0
1558 ipdstmpl(22)=0
1559 endif
1560 ! calculate ending time using initial ref-time, idat,
1561 ! and increment rinc.
1562 idat=0
1563 idat(1)=((kpds(21)-1)*100)+kpds(8)
1564 idat(2)=kpds(9)
1565 idat(3)=kpds(10)
1566 idat(4)=-500 ! EST
1567 idat(5)=kpds(11)
1568 idat(6)=kpds(12)
1569 rinc=0
1570 if ( ipdstmpl(8).eq.0 ) then
1571 rinc(3)=kpds(15)
1572 elseif ( ipdstmpl(8).eq.1 ) then
1573 rinc(2)=kpds(15)
1574 elseif ( ipdstmpl(8).eq.2 ) then
1575 rinc(1)=kpds(15)
1576 elseif ( ipdstmpl(8).eq.10 ) then
1577 rinc(2)=kpds(15) * 3
1578 elseif ( ipdstmpl(8).eq.11 ) then
1579 rinc(2)=kpds(15) * 6
1580 elseif ( ipdstmpl(8).eq.12 ) then
1581 rinc(2)=kpds(15) * 12
1582 elseif ( ipdstmpl(8).eq.13 ) then
1583 rinc(4)=kpds(15)
1584 endif
1585 call w3movdat(rinc,idat,jdat) ! calculate end date/time
1586 ipdstmpl(23)=jdat(1) ! year of end time
1587 ipdstmpl(24)=jdat(2) ! month of end time
1588 ipdstmpl(25)=jdat(3) ! day of end time
1589 ipdstmpl(26)=jdat(5) ! hour of end time
1590 ipdstmpl(27)=jdat(6) ! minute of end time
1591 ipdstmpl(28)=jdat(7) ! second of end time
1592 ipdstmpl(29)=1
1593 ipdstmpl(30)=0
1594 if (kpds(16).eq.2) then
1595 ipdstmpl(31)=255
1596 if (kpds(5).eq.15) ipdstmpl(31)=2
1597 if (kpds(5).eq.16) ipdstmpl(31)=3
1598 elseif (kpds(16).eq.3) then
1599 ipdstmpl(31)=0
1600 elseif (kpds(16).eq.4) then
1601 ipdstmpl(31)=1
1602 elseif (kpds(16).eq.5) then
1603 ipdstmpl(31)=4
1604 endif
1605 ipdstmpl(32)=2
1606 ipdstmpl(33)=kpds(13)
1607 if (kpds(13).eq.254) ipdstmpl(33)=13
1608 ipdstmpl(34)=kpds(15)-kpds(14)
1609 ipdstmpl(35)=255
1610 ipdstmpl(36)=0
1611 elseif (kpds(16).eq.51) then
1612 ! over time range...
1613 ipdsnum=9
1614 ! get GRIB2 parameter category and number from GRIB1
1615 ! parameter number
1616 call param_g1_to_g2(kprob(1),kpds(19),idum,ipdstmpl(1), &
1617 ipdstmpl(2))
1618 ipdstmpl(3)=5
1619 ipdstmpl(4)=0
1620 ipdstmpl(5)=kpds(2)
1621 ipdstmpl(6)=0
1622 ipdstmpl(7)=0
1623 ipdstmpl(8)=kpds(13)
1624 if (kpds(13).eq.254) ipdstmpl(8)=13
1625 ipdstmpl(9)=kpds(14)
1626 call cnvlevel(kpds(6),kpds(7),ipdstmpl)
1627 !ipdstmpl(9)=kpds(15)
1628 ipdstmpl(16)=0 !?
1629 ipdstmpl(17)=kclust(1) !?
1630 ipdstmpl(18)=kprob(2)-1
1631 if (ipdstmpl(18).eq.0.OR.ipdstmpl(18).eq.2) then
1632 ipdstmpl(19)=3
1633 ipdstmpl(20)=nint(xprob(1)*1000.0)
1634 else
1635 ipdstmpl(19)=0
1636 ipdstmpl(20)=0
1637 endif
1638 if (ipdstmpl(18).eq.1.OR.ipdstmpl(18).eq.2) then
1639 ipdstmpl(21)=3
1640 ipdstmpl(22)=nint(xprob(2)*1000.0)
1641 else
1642 ipdstmpl(21)=0
1643 ipdstmpl(22)=0
1644 endif
1645 ! calculate ending time using initial ref-time, idat,
1646 ! and increment rinc.
1647 idat=0
1648 idat(1)=((kpds(21)-1)*100)+kpds(8)
1649 idat(2)=kpds(9)
1650 idat(3)=kpds(10)
1651 idat(4)=-500 ! EST
1652 idat(5)=kpds(11)
1653 idat(6)=kpds(12)
1654 rinc=0
1655 if ( ipdstmpl(8).eq.0 ) then
1656 rinc(3)=kpds(15)
1657 elseif ( ipdstmpl(8).eq.1 ) then
1658 rinc(2)=kpds(15)
1659 elseif ( ipdstmpl(8).eq.2 ) then
1660 rinc(1)=kpds(15)
1661 elseif ( ipdstmpl(8).eq.10 ) then
1662 rinc(2)=kpds(15) * 3
1663 elseif ( ipdstmpl(8).eq.11 ) then
1664 rinc(2)=kpds(15) * 6
1665 elseif ( ipdstmpl(8).eq.12 ) then
1666 rinc(2)=kpds(15) * 12
1667 elseif ( ipdstmpl(8).eq.13 ) then
1668 rinc(4)=kpds(15)
1669 endif
1670 call w3movdat(rinc,idat,jdat) ! calculate end date/time
1671 ipdstmpl(23)=jdat(1) ! year of end time
1672 ipdstmpl(24)=jdat(2) ! month of end time
1673 ipdstmpl(25)=jdat(3) ! day of end time
1674 ipdstmpl(26)=jdat(5) ! hour of end time
1675 ipdstmpl(27)=jdat(6) ! minute of end time
1676 ipdstmpl(28)=jdat(7) ! second of end time
1677 ipdstmpl(29)=1
1678 ipdstmpl(30)=0
1679 ipdstmpl(31)=51
1680 ipdstmpl(32)=2
1681 ipdstmpl(33)=kpds(13)
1682 if (kpds(13).eq.254) ipdstmpl(33)=13
1683 ipdstmpl(34)=kpds(15)-kpds(14)
1684 ipdstmpl(35)=255
1685 ipdstmpl(36)=0
1686 else
1687 print *,' Unrecognized Time Range Ind for Probs = ', &
1688 kpds(16),kens(2)
1689 print *,'pds2pdtens: Couldn:t construct PDS Template '
1690 iret=10
1691 endif
1692 else ! Non-probablility Whole Ens Fcst
1693 if (kpds(16).eq.0.or.kpds(16).eq.1.or.kpds(16).eq.10) then
1694 ! At specific point in time...
1695 ipdsnum=2
1696 ! get GRIB2 parameter category and number from GRIB1
1697 ! parameter number
1698 call param_g1_to_g2(kpds(5),kpds(19),idum,ipdstmpl(1), &
1699 ipdstmpl(2))
1700 ipdstmpl(3)=4
1701 ipdstmpl(4)=0
1702 ipdstmpl(5)=kpds(2)
1703 ipdstmpl(6)=0
1704 ipdstmpl(7)=0
1705 ipdstmpl(8)=kpds(13)
1706 if (kpds(13).eq.254) ipdstmpl(8)=13
1707 !if (kpds(16).eq.10) then
1708 ! ipdstmpl(9)=(kpds(14)*256)+kpds(15)
1709 !else
1710 ipdstmpl(9)=kpds(14)
1711 !endif
1712 call cnvlevel(kpds(6),kpds(7),ipdstmpl)
1713 if (kens(4).eq.1) then
1714 ipdstmpl(16)=0
1715 elseif (kens(4).eq.2) then
1716 ipdstmpl(16)=1
1717 elseif (kens(4).eq.11) then
1718 ipdstmpl(16)=2
1719 elseif (kens(4).eq.12) then
1720 ipdstmpl(16)=3
1721 endif
1722 ipdstmpl(17)=kclust(1)
1723 elseif (kpds(16).ge.2.AND.kpds(16).le.5) then
1724 ! over time range...
1725 ipdsnum=12
1726 ! get GRIB2 parameter category and number from GRIB1
1727 ! parameter number
1728 call param_g1_to_g2(kpds(5),kpds(19),idum,ipdstmpl(1), &
1729 ipdstmpl(2))
1730 ipdstmpl(3)=4
1731 ipdstmpl(4)=0
1732 ipdstmpl(5)=kpds(2)
1733 ipdstmpl(6)=0
1734 ipdstmpl(7)=0
1735 ipdstmpl(8)=kpds(13)
1736 if (kpds(13).eq.254) ipdstmpl(8)=13
1737 ipdstmpl(9)=kpds(14)
1738 call cnvlevel(kpds(6),kpds(7),ipdstmpl)
1739 !ipdstmpl(9)=kpds(15)
1740 if (kens(4).eq.1) then
1741 ipdstmpl(16)=0
1742 elseif (kens(4).eq.2) then
1743 ipdstmpl(16)=1
1744 elseif (kens(4).eq.11) then
1745 ipdstmpl(16)=2
1746 elseif (kens(4).eq.12) then
1747 ipdstmpl(16)=3
1748 endif
1749 ipdstmpl(17)=kclust(1)
1750 ! calculate ending time using initial ref-time, idat,
1751 ! and increment rinc.
1752 idat=0
1753 idat(1)=((kpds(21)-1)*100)+kpds(8)
1754 idat(2)=kpds(9)
1755 idat(3)=kpds(10)
1756 idat(4)=-500 ! EST
1757 idat(5)=kpds(11)
1758 idat(6)=kpds(12)
1759 rinc=0
1760 if ( ipdstmpl(8).eq.0 ) then
1761 rinc(3)=kpds(15)
1762 elseif ( ipdstmpl(8).eq.1 ) then
1763 rinc(2)=kpds(15)
1764 elseif ( ipdstmpl(8).eq.2 ) then
1765 rinc(1)=kpds(15)
1766 elseif ( ipdstmpl(8).eq.10 ) then
1767 rinc(2)=kpds(15) * 3
1768 elseif ( ipdstmpl(8).eq.11 ) then
1769 rinc(2)=kpds(15) * 6
1770 elseif ( ipdstmpl(8).eq.12 ) then
1771 rinc(2)=kpds(15) * 12
1772 elseif ( ipdstmpl(8).eq.13 ) then
1773 rinc(4)=kpds(15)
1774 endif
1775 call w3movdat(rinc,idat,jdat) ! calculate end date/time
1776 ipdstmpl(18)=jdat(1) ! year of end time
1777 ipdstmpl(19)=jdat(2) ! month of end time
1778 ipdstmpl(20)=jdat(3) ! day of end time
1779 ipdstmpl(21)=jdat(5) ! hour of end time
1780 ipdstmpl(22)=jdat(6) ! minute of end time
1781 ipdstmpl(23)=jdat(7) ! second of end time
1782 ipdstmpl(24)=1
1783 ipdstmpl(25)=0
1784 if (kpds(16).eq.2) then
1785 ipdstmpl(26)=255
1786 if (kpds(5).eq.15) ipdstmpl(26)=2
1787 if (kpds(5).eq.16) ipdstmpl(26)=3
1788 elseif (kpds(16).eq.3) then
1789 ipdstmpl(26)=0
1790 elseif (kpds(16).eq.4) then
1791 ipdstmpl(26)=1
1792 elseif (kpds(16).eq.5) then
1793 ipdstmpl(26)=4
1794 endif
1795 ipdstmpl(27)=2
1796 ipdstmpl(28)=kpds(13)
1797 if (kpds(13).eq.254) ipdstmpl(28)=13
1798 ipdstmpl(29)=kpds(15)-kpds(14)
1799 ipdstmpl(30)=255
1800 ipdstmpl(31)=0
1801 elseif (kpds(16).eq.51) then
1802 ! over time range...
1803 ipdsnum=12
1804 ! get GRIB2 parameter category and number from GRIB1
1805 ! parameter number
1806 call param_g1_to_g2(kpds(5),kpds(19),idum,ipdstmpl(1), &
1807 ipdstmpl(2))
1808 ipdstmpl(3)=4
1809 ipdstmpl(4)=0
1810 ipdstmpl(5)=kpds(2)
1811 ipdstmpl(6)=0
1812 ipdstmpl(7)=0
1813 ipdstmpl(8)=kpds(13)
1814 if (kpds(13).eq.254) ipdstmpl(8)=13
1815 ipdstmpl(9)=kpds(14)
1816 call cnvlevel(kpds(6),kpds(7),ipdstmpl)
1817 !ipdstmpl(9)=kpds(15)
1818 if (kens(4).eq.1) then
1819 ipdstmpl(16)=0
1820 elseif (kens(4).eq.2) then
1821 ipdstmpl(16)=1
1822 elseif (kens(4).eq.11) then
1823 ipdstmpl(16)=2
1824 elseif (kens(4).eq.12) then
1825 ipdstmpl(16)=3
1826 endif
1827 ipdstmpl(17)=kclust(1)
1828 ! calculate ending time using initial ref-time, idat,
1829 ! and increment rinc.
1830 idat=0
1831 idat(1)=((kpds(21)-1)*100)+kpds(8)
1832 idat(2)=kpds(9)
1833 idat(3)=kpds(10)
1834 idat(4)=-500 ! EST
1835 idat(5)=kpds(11)
1836 idat(6)=kpds(12)
1837 rinc=0
1838 if ( ipdstmpl(8).eq.0 ) then
1839 rinc(3)=kpds(15)
1840 elseif ( ipdstmpl(8).eq.1 ) then
1841 rinc(2)=kpds(15)
1842 elseif ( ipdstmpl(8).eq.2 ) then
1843 rinc(1)=kpds(15)
1844 elseif ( ipdstmpl(8).eq.10 ) then
1845 rinc(2)=kpds(15) * 3
1846 elseif ( ipdstmpl(8).eq.11 ) then
1847 rinc(2)=kpds(15) * 6
1848 elseif ( ipdstmpl(8).eq.12 ) then
1849 rinc(2)=kpds(15) * 12
1850 elseif ( ipdstmpl(8).eq.13 ) then
1851 rinc(4)=kpds(15)
1852 endif
1853 call w3movdat(rinc,idat,jdat) ! calculate end date/time
1854 ipdstmpl(18)=jdat(1) ! year of end time
1855 ipdstmpl(19)=jdat(2) ! month of end time
1856 ipdstmpl(20)=jdat(3) ! day of end time
1857 ipdstmpl(21)=jdat(5) ! hour of end time
1858 ipdstmpl(22)=jdat(6) ! minute of end time
1859 ipdstmpl(23)=jdat(7) ! second of end time
1860 ipdstmpl(24)=1
1861 ipdstmpl(25)=0
1862 ipdstmpl(26)=51
1863 ipdstmpl(27)=2
1864 ipdstmpl(28)=kpds(13)
1865 if (kpds(13).eq.254) ipdstmpl(28)=13
1866 ipdstmpl(29)=kpds(15)-kpds(14)
1867 ipdstmpl(30)=255
1868 ipdstmpl(31)=0
1869 endif
1870 endif
1871 else
1872 print *,' Unrecognized Ensemble type = ',kens(2)
1873 print *,'pds2pdtens: Couldn:t construct PDS Template '
1874 iret=2
1875 endif
1876
1877 return
1878end subroutine pds2pdtens
subroutine pds2pdtens(kpds, kens, kprob, xprob, kclust, kmember, ipdsnum, ipdstmpl, numcoord, coordlist, iret)
This routine converts a GRIB1 PDS (Section 1) that includes NCEP ensemble PDS extensions to a GRIB2 P...
Definition cnv12.F90:1278
subroutine cnvlevel(ltype, lval, ipdstmpl)
This routine converts a GRIB1 Level type and Level value to GRIB2 values and fills in the appropriate...
Definition cnv12.F90:1099
subroutine gds2gdt(kgds, igds, igdstmpl, idefnum, ideflist, iret)
This routine converts a GRIB1 GDS (in format specfied in w3fi63.f) to necessary info for a GRIB2 Grid...
Definition cnv12.F90:407
subroutine pds2pdt(kpds, ipdsnum, ipdstmpl, numcoord, coordlist, iret)
This routine converts a GRIB1 PDS (Section 1) info to a GRIB2 PDS (Section 4) info with appropriate P...
Definition cnv12.F90:790
subroutine cnv12(ifl1, ifl2, ipack, usemiss, imiss, uvvect, table_ver)
This subroutine converts every GRIB1 field in a file to a GRIB2 field.
Definition cnv12.F90:48
subroutine g2_gbytec(in, iout, iskip, nbits)
Extract one arbitrary size big-endian value (up to 32 bits) from a packed bit string into one element...
Definition g2bytes.F90:21
subroutine mkieee(a, rieee, num)
Copy an array of real to an array of 32-bit IEEE floating points.
Definition g2bytes.F90:685
subroutine addgrid(cgrib, lcgrib, igds, igdstmpl, igdstmplen, ideflist, idefnum, ierr)
Add a Grid Definition Section (Section 3) to a GRIB2 message.
Definition g2create.F90:682
subroutine gribcreate(cgrib, lcgrib, listsec0, listsec1, ierr)
Initialize a new GRIB2 message and pack GRIB2 sections 0 (Indicator) and 1 (Identification).
Definition g2create.F90:54
subroutine addfield(cgrib, lcgrib, ipdsnum, ipdstmpl, ipdstmplen, coordlist, numcoord, idrsnum, idrstmpl, idrstmplen, fld, ngrdpts, ibmap, bmap, ierr)
Pack up Sections 4 through 7 for a field and add them to a GRIB2 message.
Definition g2create.F90:199
subroutine gribend(cgrib, lcgrib, lengrib, ierr)
Finalize a GRIB2 message after all grids and fields have been added.
This Fortran Module contains info on all the available ECMWF GRIB Parameters.
subroutine param_ecmwf_g1_to_g2(g1val, g1ver, g2disc, g2cat, g2num)
This subroutine returns the corresponding GRIB2 Discipline Category and Number for a given GRIB1 para...
This Fortran Module contains info on all the available GRIB Parameters, and their GRIB1 and GRIB2 cod...
Definition params.F90:12
subroutine param_g1_to_g2(g1val, g1ver, g2disc, g2cat, g2num)
This subroutine returns the corresponding GRIB2 Discipline Category and Number for a given GRIB1 para...
Definition params.F90:1084
subroutine skgb(lugb, iseek, mseek, lskip, lgrib)
Search a file for the next GRIB1 or GRIB2 message.
Definition skgb.F90:21