NCEPLIBS-g2 4.0.0
Loading...
Searching...
No Matches
g2get.F90
Go to the documentation of this file.
1
4
52subroutine gb_info(cgrib, lcgrib, listsec0, listsec1, &
53 numfields, numlocal, maxlocal, ierr)
54 use g2logging
55 implicit none
56
57 character(len = 1), intent(in) :: cgrib(lcgrib)
58 integer, intent(in) :: lcgrib
59 integer, intent(out) :: listsec0(3), listsec1(13)
60 integer, intent(out) :: numlocal, numfields, maxlocal, ierr
61
62 character(len = 4), parameter :: grib = 'GRIB', c7777 = '7777'
63 character(len = 4) :: ctemp
64 integer, parameter :: zero = 0, one = 1
65 integer, parameter :: mapsec1len = 13
66 integer, parameter :: mapsec1(mapsec1len) = (/ 2, 2, 1, 1, 1, 2, 1, 1, 1, 1, 1, 1, 1 /)
67 integer :: iofst, istart
68 integer :: nbits, lensec1, lensec0, lensec, lenposs, lengrib, j
69 integer (kind = 8) :: lengrib8
70 integer :: i, ipos, isecnum
71
72 interface
73 subroutine g2_gbytec(in, iout, iskip, nbits)
74 character*1, intent(in) :: in(*)
75 integer, intent(inout) :: iout(*)
76 integer, intent(in) :: iskip, nbits
77 end subroutine g2_gbytec
78 subroutine g2_gbytec1(in, siout, iskip, nbits)
79 character*1, intent(in) :: in(*)
80 integer, intent(inout) :: siout
81 integer, intent(in) :: iskip, nbits
82 end subroutine g2_gbytec1
83 subroutine g2_gbytec81(in, siout, iskip, nbits)
84 character*1, intent(in) :: in(*)
85 integer (kind = 8), intent(inout) :: siout
86 integer, intent(in) :: iskip, nbits
87 integer (kind = 8) :: iout(1)
88 end subroutine g2_gbytec81
89 end interface
90
91#ifdef LOGGING
92 write(g2_log_msg, *) 'gb_info: lcgrib ', lcgrib
93 call g2_log(1)
94#endif
95
96 ierr = 0
97 numlocal = 0
98 numfields = 0
99 maxlocal = 0
100
101 ! Check for beginning of GRIB message in the first 100 bytes.
102 istart = 0
103 do j = 1, 100
104 ctemp = cgrib(j) // cgrib(j + 1) // cgrib(j + 2) // cgrib(j+3)
105 if (ctemp .eq. grib ) then
106 istart = j
107 exit
108 endif
109 enddo
110 if (istart .eq. 0) then
111 print *, 'gb_info: Beginning characters GRIB not found.'
112 ierr = 1
113 return
114 endif
115
116 ! Unpack Section 0 - Indicator Section.
117 iofst = 8 * (istart + 5)
118 call g2_gbytec(cgrib, listsec0(1), iofst, 8) ! Discipline
119 iofst = iofst + 8
120 call g2_gbytec(cgrib, listsec0(2), iofst, 8) ! GRIB edition number
121 iofst = iofst + 8
122
123#ifdef LOGGING
124 ! Log results for debugging.
125 write(g2_log_msg, *) 'before getting len: iofst/8 ', iofst/8
126 call g2_log(2)
127#endif
128 ! Bytes 9-16 contain the length of GRIB message. This is an 8-byte
129 ! value, but unfortunately we only have a 4-byte subroutine
130 ! parameter for it.
131 call g2_gbytec81(cgrib, lengrib8, iofst, 64)
132 lengrib = int(lengrib8, kind(4))
133 iofst = iofst + 32
134 iofst = iofst + 32
135 listsec0(3) = lengrib
136 lensec0 = 16
137 ipos = istart + lensec0
138
139 ! The g2 library only handles GRIB2.
140 if (listsec0(2) .ne. 2) then
141 print *, 'gb_info: can only decode GRIB edition 2.'
142 ierr = 2
143 return
144 endif
145
146 ! Unpack Section 1 - Identification Section.
147 call g2_gbytec1(cgrib, lensec1, iofst, 32) ! Length of Section 1
148 iofst = iofst + 32
149 call g2_gbytec1(cgrib, isecnum, iofst, 8) ! Section number ( 1 )
150 iofst = iofst + 8
151 if (isecnum .ne. 1) then
152 print *, 'gb_info: Could not find section 1.'
153 ierr = 3
154 return
155 endif
156
157 ! Unpack each input value in array listsec1 into the the appropriate
158 ! number of octets, which are specified in corresponding entries in
159 ! array mapsec1.
160 do i = 1, mapsec1len
161 nbits = mapsec1(i) * 8
162 call g2_gbytec(cgrib, listsec1(i), iofst, nbits)
163 iofst = iofst + nbits
164 enddo
165 ipos = ipos + lensec1
166
167 ! Loop through the remaining sections to see if they are valid.
168 ! Also count the number of times Section 2 and Section 4 appear.
169 do
170 ctemp = cgrib(ipos) // cgrib(ipos + 1) // cgrib(ipos + 2) // cgrib(ipos + 3)
171 if (ctemp .eq. c7777 ) then
172 ipos = ipos + 4
173 if (ipos .ne. (istart + lengrib)) then
174 print *, 'gb_info: "7777" found, but not where expected.'
175 ierr = 4
176 return
177 endif
178 exit
179 endif
180 iofst = (ipos - 1) * 8
181 call g2_gbytec1(cgrib, lensec, iofst, 32) ! Get Length of Section
182 iofst = iofst + 32
183 call g2_gbytec1(cgrib, isecnum, iofst, 8) ! Get Section number
184 iofst = iofst + 8
185 ipos = ipos+lensec ! Update beginning of section pointer
186 if (ipos .gt. (istart + lengrib)) then
187 print *, 'gb_info: "7777" not found at end of GRIB message.'
188 ierr = 5
189 return
190 endif
191 if (isecnum .ge. 2 .AND. isecnum .le. 7) then
192 if (isecnum .eq. 2) then ! Local Section 2
193 ! Increment counter for total number of local sections found.
194 numlocal = numlocal + 1
195 lenposs = lensec - 5
196 if (lenposs .gt. maxlocal) maxlocal = lenposs
197 elseif (isecnum .eq. 4) then
198 ! Increment counter for total number of fields found.
199 numfields = numfields + 1
200 endif
201 else
202 print *, 'gb_info: Invalid section number found in GRIB message: ', isecnum
203 ierr = 6
204 return
205 endif
206 enddo
207
208end subroutine gb_info
209
282subroutine gribinfo(cgrib, lcgrib, listsec0, listsec1, &
283 numlocal, numfields, maxvals, ierr)
284 implicit none
285
286 character(len = 1), intent(in) :: cgrib(lcgrib)
287 integer, intent(in) :: lcgrib
288 integer, intent(out) :: listsec0(3), listsec1(13), maxvals(7)
289 integer, intent(out) :: numlocal, numfields, ierr
290
291 integer :: i, ipos, isecnum, j, lengrib, lenposs, lensec, lensec0, lensec1
292 integer :: maxcoordlist, maxdeflist, maxdrstmpl, maxgridpts, maxpdstmpl, maxgdstmpl
293 integer :: maxsec2len, nbits, nbyte, ngdpts, numcoord
294 character(len = 4), parameter :: grib = 'GRIB', c7777 = '7777'
295 character(len = 4) :: ctemp
296 integer, parameter :: zero = 0, one = 1
297 integer, parameter :: mapsec1len = 13
298 integer, parameter :: mapsec1(mapsec1len) = (/ 2, 2, 1, 1, 1, 2, 1, 1, &
299 1, 1, 1, 1, 1 /)
300 integer iofst, istart
301 integer (kind = 8) :: lengrib8
302
303 ierr = 0
304 numlocal = 0
305 numfields = 0
306 maxsec2len = 1
307 maxgdstmpl = 1
308 maxdeflist = 1
309 maxpdstmpl = 1
310 maxcoordlist = 1
311 maxdrstmpl = 1
312 maxgridpts = 0
313
314 ! Check for beginning of GRIB message in the first 100 bytes.
315 istart = 0
316 do j = 1, 100
317 ctemp = cgrib(j) // cgrib(j + 1) // cgrib(j + 2) // cgrib(j + 3)
318 if (ctemp .eq. grib) then
319 istart = j
320 exit
321 endif
322 enddo
323 if (istart .eq. 0) then
324 print *, 'gribinfo: Beginning characters GRIB not found.'
325 ierr = 1
326 return
327 endif
328
329 ! Unpack Section 0 - Indicator Section.
330 iofst = 8 * (istart + 5)
331 call g2_gbytec(cgrib, listsec0(1), iofst, 8) ! Discipline
332 iofst = iofst + 8
333 call g2_gbytec(cgrib, listsec0(2), iofst, 8) ! GRIB edition number
334 iofst = iofst + 8
335 ! Bytes 9-16 contain the length of GRIB message. This is an 8-byte
336 ! value, but unfortunately we only have a 4-byte subroutine
337 ! parameter for it.
338 call g2_gbytec81(cgrib, lengrib8, iofst, 64) ! Length of GRIB message
339 lengrib = int(lengrib8, kind(4))
340 iofst = iofst + 32
341 iofst = iofst + 32
342 listsec0(3) = lengrib
343 lensec0 = 16
344 ipos = istart + lensec0
345
346 ! Currently handles only GRIB Edition 2.
347 if (listsec0(2) .ne. 2) then
348 print *, 'gribinfo: can only decode GRIB edition 2.'
349 ierr = 2
350 return
351 endif
352
353 ! Unpack Section 1 - Identification Section.
354 call g2_gbytec1(cgrib, lensec1, iofst, 32) ! Length of Section 1
355 iofst = iofst + 32
356 call g2_gbytec1(cgrib, isecnum, iofst, 8) ! Section number (1)
357 iofst = iofst + 8
358 if (isecnum .ne. 1) then
359 print *, 'gribinfo: Could not find section 1.'
360 ierr = 3
361 return
362 endif
363
364 ! Unpack each input value in array listsec1 into the the appropriate
365 ! number of octets, which are specified in corresponding entries in
366 ! array mapsec1.
367 do i = 1, mapsec1len
368 nbits = mapsec1(i) * 8
369 call g2_gbytec(cgrib, listsec1(i), iofst, nbits)
370 iofst = iofst + nbits
371 enddo
372 ipos = ipos + lensec1
373
374 ! Loop through the remaining sections keeping track of the length of
375 ! each. Also count the number of times Section 2 and Section 4
376 ! appear.
377 do
378 ctemp = cgrib(ipos) // cgrib(ipos + 1) // cgrib(ipos + 2) // &
379 cgrib(ipos + 3)
380 if (ctemp .eq. c7777) then
381 ipos = ipos + 4
382 if (ipos .ne. (istart + lengrib)) then
383 print *, 'gribinfo: "7777" found, but not where expected.'
384 ierr = 4
385 return
386 endif
387 exit
388 endif
389 iofst = (ipos - 1) * 8
390 call g2_gbytec1(cgrib, lensec, iofst, 32) ! Get Length of Section
391 iofst = iofst + 32
392 call g2_gbytec1(cgrib, isecnum, iofst, 8) ! Get Section number
393 iofst = iofst + 8
394 ipos = ipos + lensec ! Update beginning of section pointer
395 if (ipos .gt. (istart + lengrib)) then
396 print *, 'gribinfo: "7777" not found at end of GRIB message.'
397 ierr = 5
398 return
399 endif
400 if (isecnum .eq. 2) then ! Local Section 2
401 ! Increment counter for total number of local sections found
402 ! and determine largest Section 2 in message.
403 numlocal = numlocal + 1
404 lenposs = lensec-5
405 if (lenposs .gt. maxsec2len) maxsec2len = lenposs
406 elseif (isecnum .eq. 3) then
407 iofst = iofst + 8 ! skip source of grid def.
408 call g2_gbytec1(cgrib, ngdpts, iofst, 32) ! Get Num of Grid Points
409 iofst = iofst + 32
410 call g2_gbytec1(cgrib, nbyte, iofst, 8) ! Get Num octets for opt. list
411 iofst = iofst + 8
412 if (ngdpts .gt. maxgridpts) maxgridpts = ngdpts
413 lenposs = lensec - 14
414 if (lenposs .gt. maxgdstmpl) maxgdstmpl = lenposs
415 if (nbyte .ne. 0) then
416 lenposs = lenposs / nbyte
417 if (lenposs .gt. maxdeflist) maxdeflist = lenposs
418 endif
419 elseif (isecnum .eq. 4) then
420 numfields = numfields + 1
421 call g2_gbytec1(cgrib, numcoord, iofst, 16) ! Get Num of Coord Values
422 iofst = iofst + 16
423 if (numcoord .ne. 0) then
424 if (numcoord .gt. maxcoordlist) maxcoordlist = numcoord
425 endif
426 lenposs = lensec - 9
427 if (lenposs.gt.maxpdstmpl) maxpdstmpl = lenposs
428 elseif (isecnum .eq. 5) then
429 lenposs = lensec-11
430 if (lenposs .gt. maxdrstmpl) maxdrstmpl = lenposs
431 endif
432 enddo
433
434 maxvals(1) = maxsec2len
435 maxvals(2) = maxgdstmpl
436 maxvals(3) = maxdeflist
437 maxvals(4) = maxpdstmpl
438 maxvals(5) = maxcoordlist
439 maxvals(6) = maxdrstmpl
440 maxvals(7) = maxgridpts
441end subroutine gribinfo
442
543subroutine getfield(cgrib, lcgrib, ifldnum, igds, igdstmpl, &
544 igdslen, ideflist, idefnum, ipdsnum, ipdstmpl, ipdslen, &
545 coordlist, numcoord, ndpts, idrsnum, idrstmpl, idrslen, &
546 ibmap, bmap, fld, ierr)
547 implicit none
548
549 character(len = 1), intent(in) :: cgrib(lcgrib)
550 integer, intent(in) :: lcgrib, ifldnum
551 integer, intent(out) :: igds(*), igdstmpl(*), ideflist(*)
552 integer, intent(out) :: ipdsnum, ipdstmpl(*)
553 integer, intent(out) :: idrsnum, idrstmpl(*)
554 integer, intent(out) :: ndpts, ibmap, idefnum, numcoord
555 integer, intent(out) :: igdslen, ipdslen, idrslen
556 real, intent(out) :: fld(*), coordlist(*)
557 logical*1, intent(out) :: bmap(*)
558 integer, intent(out) :: ierr
559
560 character(len = 4), parameter :: grib = 'GRIB', c7777 = '7777'
561 character(len = 4) :: ctemp
562 integer :: listsec0(2)
563 integer :: iofst, istart
564 real (kind = 4) :: ieee(1)
565 logical :: have3, have4, have5, have6, have7
566 integer (kind = 8) :: lengrib8
567 integer :: numfld, j, lengrib, lensec0, ipos
568 integer :: lensec, isecnum, jerr, ier, numlocal
569
570 have3 = .false.
571 have4 = .false.
572 have5 = .false.
573 have6 = .false.
574 have7 = .false.
575 ierr = 0
576 numfld = 0
577 numlocal = 0
578
579 ! Check for valid request number.
580 if (ifldnum .le. 0) then
581 print *, 'getfield: Request for field number ' &
582 ,'must be positive.'
583 ierr = 3
584 return
585 endif
586
587 ! Check for beginning of GRIB message in the first 100 bytes.
588 istart = 0
589 do j = 1, 100
590 ctemp = cgrib(j) // cgrib(j + 1) // cgrib(j + 2) // cgrib(j + 3)
591 if (ctemp .eq. grib) then
592 istart = j
593 exit
594 endif
595 enddo
596 if (istart .eq. 0) then
597 print *, 'getfield: Beginning characters GRIB not found.'
598 ierr = 1
599 return
600 endif
601
602 ! Unpack Section 0 - Indicator Section.
603 iofst = 8 * (istart + 5)
604 call g2_gbytec(cgrib, listsec0(1), iofst, 8) ! Discipline
605 iofst = iofst + 8
606 call g2_gbytec(cgrib, listsec0(2), iofst, 8) ! GRIB edition number
607 iofst = iofst + 8
608 ! Bytes 9-16 contain the length of GRIB message. This is an 8-byte
609 ! value, but unfortunately we only have a 4-byte subroutine
610 ! parameter for it.
611 call g2_gbytec81(cgrib, lengrib8, iofst, 64) ! Length of GRIB message
612 lengrib = int(lengrib8, kind(4))
613 iofst = iofst + 32
614 iofst = iofst + 32
615 lensec0 = 16
616 ipos = istart + lensec0
617
618 ! The g2 library only handles GRIB2.
619 if (listsec0(2) .ne. 2) then
620 print *, 'getfield: can only decode GRIB edition 2.'
621 ierr = 2
622 return
623 endif
624
625 ! Loop through the remaining sections keeping track of the length of
626 ! each. Also keep the latest Grid Definition Section info. Unpack
627 ! the requested field number.
628 do
629 ! Check to see if we are at end of GRIB message
630 ctemp = cgrib(ipos) // cgrib(ipos + 1) // cgrib(ipos + 2) // &
631 cgrib(ipos + 3)
632 if (ctemp .eq. c7777) then
633 ipos = ipos + 4
634 ! If end of GRIB message not where expected, issue error
635 if (ipos.ne.(istart + lengrib)) then
636 print *, 'getfield: "7777" found, but not ' &
637 ,'where expected.'
638 ierr = 4
639 return
640 endif
641 exit
642 endif
643 ! Get length of Section and Section number
644 iofst = (ipos - 1) * 8
645 call g2_gbytec1(cgrib, lensec, iofst, 32) ! Get Length of Section
646 iofst = iofst + 32
647 call g2_gbytec1(cgrib, isecnum, iofst, 8) ! Get Section number
648 iofst = iofst + 8
649
650 ! If found Section 3, unpack the GDS info using the appropriate
651 ! template. Save in case this is the latest grid before the
652 ! requested field.
653 if (isecnum .eq. 3) then
654 iofst = iofst - 40 ! reset offset to beginning of section
655 call unpack3(cgrib, lcgrib, iofst, igds, igdstmpl, &
656 igdslen, ideflist, idefnum, jerr)
657 if (jerr .eq. 0) then
658 have3 = .true.
659 else
660 ierr = 10
661 return
662 endif
663 endif
664
665 ! If found Section 4, check to see if this field is the one
666 ! requested.
667 if (isecnum .eq. 4) then
668 numfld = numfld + 1
669 if (numfld .eq. ifldnum) then
670 iofst = iofst - 40 ! reset offset to beginning of section
671 call unpack4(cgrib, lcgrib, iofst, ipdsnum, ipdstmpl, &
672 ipdslen, coordlist, numcoord, jerr)
673 if (jerr .eq. 0) then
674 have4 = .true.
675 else
676 ierr = 11
677 return
678 endif
679 endif
680 endif
681
682 ! If found Section 5, check to see if this field is the one
683 ! requested.
684 if ((isecnum .eq. 5) .and. (numfld .eq. ifldnum)) then
685 iofst = iofst - 40 ! reset offset to beginning of section
686 call unpack5(cgrib, lcgrib, iofst, ndpts, idrsnum, &
687 idrstmpl, idrslen, jerr)
688 if (jerr .eq. 0) then
689 have5 = .true.
690 else
691 ierr = 12
692 return
693 endif
694 endif
695
696 ! If found Section 6, Unpack bitmap. Save in case this is the
697 ! latest bitmap before the requested field.
698 if (isecnum .eq. 6) then
699 iofst = iofst - 40 ! reset offset to beginning of section
700 call unpack6(cgrib, lcgrib, iofst, igds(2), ibmap, bmap, &
701 jerr)
702 if (jerr .eq. 0) then
703 have6 = .true.
704 else
705 ierr = 13
706 return
707 endif
708 endif
709
710 ! If found Section 7, check to see if this field is the one
711 ! requested.
712 if ((isecnum .eq. 7) .and. (numfld .eq. ifldnum)) then
713 if (idrsnum .eq. 0) then
714 call simunpack(cgrib(ipos + 5), lensec - 6, idrstmpl, &
715 ndpts, fld)
716 have7 = .true.
717 elseif (idrsnum .eq. 2 .or. idrsnum .eq. 3) then
718 call comunpack(cgrib(ipos + 5), lensec - 6, lensec, &
719 idrsnum,idrstmpl, ndpts, fld, ier)
720 if (ier .ne. 0) then
721 ierr = 14
722 return
723 endif
724 have7 = .true.
725 elseif (idrsnum .eq. 50) then
726 call simunpack(cgrib(ipos + 5), lensec - 6, idrstmpl, &
727 ndpts - 1, fld(2))
728 ieee = transfer(idrstmpl(5), ieee, 1)
729 call rdieee(ieee, fld(1), 1)
730 have7 = .true.
731 elseif (idrsnum .eq. 40 .or. idrsnum .eq. 40000) then
732 call jpcunpack(cgrib(ipos + 5), lensec - 5, idrstmpl, &
733 ndpts, fld)
734 have7 = .true.
735 elseif (idrsnum .eq. 41 .or. idrsnum .eq. 40010) then
736 call pngunpack(cgrib(ipos + 5), lensec - 5, idrstmpl, &
737 ndpts, fld)
738 have7 = .true.
739#ifdef USE_AEC
740 elseif (idrsnum .eq. 42) then
741 call aecunpack(cgrib(ipos + 5), lensec - 5, idrstmpl, &
742 ndpts, fld)
743 have7 = .true.
744#endif /* USE_AEC */
745 else
746 print *, 'getfield: Data Representation Template ', &
747 idrsnum, ' not yet implemented.'
748 ierr = 9
749 return
750 endif
751 endif
752
753 ! Check to see if we read pass the end of the GRIB message and
754 ! missed the terminator string '7777'.
755 ipos = ipos + lensec ! Update beginning of section pointer
756 if (ipos .gt. (istart + lengrib)) then
757 print *, 'getfield: "7777" not found at end' &
758 ,' of GRIB message.'
759 ierr = 7
760 return
761 endif
762
763 if (have3 .and. have4 .and. have5 .and. have6 .and. have7) &
764 return
765
766 enddo
767
768 ! If exited from above loop, the end of the GRIB message was reached
769 ! before the requested field was found.
770 print *, 'getfield: GRIB message contained ', numlocal, &
771 ' different fields.'
772 print *, 'getfield: The request was for the ', ifldnum, &
773 ' field.'
774 ierr = 6
775
776end subroutine getfield
777
817subroutine unpack3(cgrib, lcgrib, iofst, igds, igdstmpl, &
818 mapgridlen, ideflist, idefnum, ierr)
819 use gridtemplates
820 implicit none
821
822 character(len = 1), intent(in) :: cgrib(lcgrib)
823 integer, intent(in) :: lcgrib
824 integer, intent(inout) :: iofst
825 integer, intent(out) :: igds(*), igdstmpl(*), ideflist(*)
826 integer, intent(out) :: ierr, idefnum
827
828 integer, allocatable :: mapgrid(:)
829 integer :: mapgridlen, ibyttem
830 logical needext
831 integer :: lensec, iret, i, nbits, isign, newmapgridlen
832
833 ierr = 0
834
835 call g2_gbytec1(cgrib, lensec, iofst, 32) ! Get Length of Section
836 iofst = iofst + 32
837 iofst = iofst + 8 ! skip section number
838
839 call g2_gbytec1(cgrib, igds(1), iofst, 8) ! Get source of Grid def.
840 iofst = iofst + 8
841 call g2_gbytec1(cgrib, igds(2), iofst, 32) ! Get number of grid pts.
842 iofst = iofst + 32
843 call g2_gbytec1(cgrib, igds(3), iofst, 8) ! Get num octets for opt. list
844 iofst = iofst + 8
845 call g2_gbytec1(cgrib, igds(4), iofst, 8) ! Get interpret. for opt. list
846 iofst = iofst + 8
847 call g2_gbytec1(cgrib, igds(5), iofst, 16) ! Get Grid Def Template num.
848 iofst = iofst + 16
849 if (igds(1) .eq. 0) then
850 allocate(mapgrid(lensec))
851 ! Get Grid Definition Template
852 call getgridtemplate(igds(5), mapgridlen, mapgrid, needext, &
853 iret)
854 if (iret .ne. 0) then
855 ierr = 5
856 return
857 endif
858 else
859 mapgridlen = 0
860 needext = .false.
861 endif
862
863 ! Unpack each value into array igdstmpl from the the appropriate
864 ! number of octets, which are specified in corresponding entries in
865 ! array mapgrid.
866 ibyttem = 0
867 do i = 1, mapgridlen
868 nbits = iabs(mapgrid(i)) * 8
869 if (mapgrid(i) .ge. 0) then
870 call g2_gbytec1(cgrib, igdstmpl(i), iofst, nbits)
871 else
872 call g2_gbytec1(cgrib, isign, iofst, 1)
873 call g2_gbytec1(cgrib, igdstmpl(i), iofst + 1, nbits-1)
874 if (isign .eq. 1) igdstmpl(i) = -igdstmpl(i)
875 endif
876 iofst = iofst + nbits
877 ibyttem = ibyttem + iabs(mapgrid(i))
878 enddo
879
880 ! Check to see if the Grid Definition Template needs to be
881 ! extended. The number of values in a specific template may vary
882 ! depending on data specified in the "static" part of the template.
883 if (needext) then
884 call extgridtemplate(igds(5), igdstmpl, newmapgridlen, &
885 mapgrid)
886 ! Unpack the rest of the Grid Definition Template
887 do i = mapgridlen + 1, newmapgridlen
888 nbits = iabs(mapgrid(i)) * 8
889 if (mapgrid(i) .ge. 0) then
890 call g2_gbytec1(cgrib, igdstmpl(i), iofst, nbits)
891 else
892 call g2_gbytec1(cgrib, isign, iofst, 1)
893 call g2_gbytec1(cgrib, igdstmpl(i), iofst + 1, nbits - &
894 1)
895 if (isign .eq. 1) igdstmpl(i) = -igdstmpl(i)
896 endif
897 iofst = iofst + nbits
898 ibyttem = ibyttem + iabs(mapgrid(i))
899 enddo
900 mapgridlen = newmapgridlen
901 endif
902
903 ! Unpack optional list of numbers defining number of points in each
904 ! row or column, if included. This is used for non regular grids.
905 if (igds(3) .ne. 0) then
906 nbits = igds(3) * 8
907 idefnum = (lensec - 14 - ibyttem) / igds(3)
908 call g2_gbytesc(cgrib, ideflist, iofst, nbits, 0, idefnum)
909 iofst = iofst + (nbits * idefnum)
910 else
911 idefnum = 0
912 endif
913 if (allocated(mapgrid)) deallocate(mapgrid)
914end subroutine unpack3
915
942subroutine unpack4(cgrib, lcgrib, iofst, ipdsnum, ipdstmpl, &
943 mappdslen, coordlist, numcoord, ierr)
944 use pdstemplates
945 implicit none
946
947 character(len = 1), intent(in) :: cgrib(lcgrib)
948 integer, intent(in) :: lcgrib
949 integer, intent(inout) :: iofst
950 real, intent(out) :: coordlist(*)
951 integer, intent(out) :: ipdsnum, ipdstmpl(*)
952 integer, intent(out) :: ierr, numcoord
953
954 real(4), allocatable :: coordieee(:)
955 integer, allocatable :: mappds(:)
956 integer :: mappdslen
957 logical needext
958 integer :: lensec, iret, i, nbits, isign, newmappdslen
959
960 ierr = 0
961
962 call g2_gbytec1(cgrib, lensec, iofst, 32) ! Get Length of Section
963 iofst = iofst + 32
964 iofst = iofst + 8 ! skip section number
965 allocate(mappds(lensec))
966
967 call g2_gbytec1(cgrib, numcoord, iofst, 16) ! Get num of coordinate values
968 iofst = iofst + 16
969 call g2_gbytec1(cgrib, ipdsnum, iofst, 16) ! Get Prod. Def Template num.
970 iofst = iofst + 16
971
972 ! Get Product Definition Template.
973 call getpdstemplate(ipdsnum, mappdslen, mappds, needext, iret)
974 if (iret.ne.0) then
975 ierr = 5
976 return
977 endif
978
979 ! Unpack each value into array ipdstmpl from the the appropriate
980 ! number of octets, which are specified in corresponding entries in
981 ! array mappds.
982 do i = 1, mappdslen
983 nbits = iabs(mappds(i)) * 8
984 if (mappds(i) .ge. 0) then
985 call g2_gbytec1(cgrib, ipdstmpl(i), iofst, nbits)
986 else
987 call g2_gbytec1(cgrib, isign, iofst, 1)
988 call g2_gbytec1(cgrib, ipdstmpl(i), iofst + 1, nbits-1)
989 if (isign .eq. 1) ipdstmpl(i) = -ipdstmpl(i)
990 endif
991 iofst = iofst + nbits
992 enddo
993
994 ! Check to see if the Product Definition Template needs to be
995 ! extended. The number of values in a specific template may vary
996 ! depending on data specified in the "static" part of the template.
997 if (needext) then
998 call extpdstemplate(ipdsnum, ipdstmpl, newmappdslen, mappds)
999
1000 ! Unpack the rest of the Product Definition Template.
1001 do i = mappdslen + 1, newmappdslen
1002 nbits = iabs(mappds(i)) * 8
1003 if (mappds(i) .ge. 0) then
1004 call g2_gbytec1(cgrib, ipdstmpl(i), iofst, nbits)
1005 else
1006 call g2_gbytec1(cgrib, isign, iofst, 1)
1007 call g2_gbytec1(cgrib, ipdstmpl(i), iofst + 1, nbits-1)
1008 if (isign .eq. 1) ipdstmpl(i) = -ipdstmpl(i)
1009 endif
1010 iofst = iofst + nbits
1011 enddo
1012 mappdslen = newmappdslen
1013 endif
1014
1015 ! Get Optional list of vertical coordinate values after the Product
1016 ! Definition Template, if necessary.
1017 if (numcoord .ne. 0) then
1018 allocate (coordieee(numcoord))
1019 call g2_gbytescr(cgrib, coordieee, iofst, 32, 0, numcoord)
1020 call rdieee(coordieee, coordlist, numcoord)
1021 deallocate (coordieee)
1022 iofst = iofst + (32*numcoord)
1023 endif
1024 if (allocated(mappds)) deallocate(mappds)
1025end subroutine unpack4
1026
1049subroutine unpack5(cgrib, lcgrib, iofst, ndpts, idrsnum, &
1050 idrstmpl, mapdrslen, ierr)
1051 use drstemplates
1052 implicit none
1053
1054 character(len = 1), intent(in) :: cgrib(lcgrib)
1055 integer, intent(in) :: lcgrib
1056 integer, intent(inout) :: iofst
1057 integer, intent(out) :: ndpts, idrsnum, idrstmpl(*)
1058 integer, intent(out) :: ierr
1059
1060 ! integer, allocatable :: mapdrs(:)
1061 integer, allocatable :: mapdrs(:)
1062 integer :: mapdrslen
1063 logical needext
1064 integer :: lensec, i, nbits, isign, newmapdrslen, iret
1065
1066 ierr = 0
1067
1068 call g2_gbytec1(cgrib, lensec, iofst, 32) ! Get Length of Section
1069 iofst = iofst + 32
1070 iofst = iofst + 8 ! skip section number
1071 allocate(mapdrs(lensec))
1072
1073 call g2_gbytec1(cgrib, ndpts, iofst, 32) ! Get num of data points
1074 iofst = iofst + 32
1075 call g2_gbytec1(cgrib, idrsnum, iofst, 16) ! Get Data Rep Template Num.
1076 iofst = iofst + 16
1077
1078 ! Gen Data Representation Template
1079 call getdrstemplate(idrsnum, mapdrslen, mapdrs, needext, iret)
1080 if (iret.ne.0) then
1081 ierr = 7
1082 return
1083 endif
1084
1085 ! Unpack each value into array ipdstmpl from the the appropriate
1086 ! number of octets, which are specified in corresponding entries in
1087 ! array mappds.
1088 do i = 1, mapdrslen
1089 nbits = iabs(mapdrs(i))*8
1090 if (mapdrs(i).ge.0) then
1091 call g2_gbytec1(cgrib, idrstmpl(i), iofst, nbits)
1092 else
1093 call g2_gbytec1(cgrib, isign, iofst, 1)
1094 call g2_gbytec1(cgrib, idrstmpl(i), iofst + 1, nbits-1)
1095 if (isign.eq.1) idrstmpl(i) = -idrstmpl(i)
1096 endif
1097 iofst = iofst + nbits
1098 enddo
1099
1100 ! Check to see if the Data Representation Template needs to be
1101 ! extended. The number of values in a specific template may vary
1102 ! depending on data specified in the "static" part of the template.
1103 if (needext) then
1104 call extdrstemplate(idrsnum, idrstmpl, newmapdrslen, mapdrs)
1105 ! Unpack the rest of the Data Representation Template.
1106 do i = mapdrslen + 1, newmapdrslen
1107 nbits = iabs(mapdrs(i)) * 8
1108 if (mapdrs(i) .ge. 0) then
1109 call g2_gbytec1(cgrib, idrstmpl(i), iofst, nbits)
1110 else
1111 call g2_gbytec1(cgrib, isign, iofst, 1)
1112 call g2_gbytec1(cgrib, idrstmpl(i), iofst + 1, nbits - 1)
1113 if (isign .eq. 1) idrstmpl(i) = -idrstmpl(i)
1114 endif
1115 iofst = iofst + nbits
1116 enddo
1117 mapdrslen = newmapdrslen
1118 endif
1119 if (allocated(mapdrs)) deallocate(mapdrs)
1120end subroutine unpack5
1121
1143subroutine unpack6(cgrib, lcgrib, iofst, ngpts, ibmap, bmap, ierr)
1144 implicit none
1145
1146 character(len = 1), intent(in) :: cgrib(lcgrib)
1147 integer, intent(in) :: lcgrib, ngpts
1148 integer, intent(inout) :: iofst
1149 integer, intent(out) :: ibmap
1150 integer, intent(out) :: ierr
1151 logical*1, intent(out) :: bmap(ngpts)
1152
1153 integer :: intbmap(ngpts)
1154 integer :: j
1155
1156 ierr = 0
1157
1158 iofst = iofst + 32 ! skip Length of Section
1159 iofst = iofst + 8 ! skip section number
1160
1161 call g2_gbytec1(cgrib, ibmap, iofst, 8) ! Get bit-map indicator
1162 iofst = iofst + 8
1163
1164 if (ibmap.eq.0) then ! Unpack bitmap
1165 call g2_gbytesc(cgrib, intbmap, iofst, 1, 0, ngpts)
1166 iofst = iofst + ngpts
1167 do j = 1, ngpts
1168 bmap(j) = .true.
1169 if (intbmap(j).eq.0) bmap(j) = .false.
1170 enddo
1171 elseif (ibmap.eq.254) then ! Use previous bitmap
1172 return
1173 elseif (ibmap.eq.255) then ! No bitmap in message
1174 bmap(1:ngpts) = .true.
1175 else
1176 print *, 'unpack6: Predefined bitmap ', ibmap, &
1177 ' not recognized.'
1178 ierr = 4
1179 endif
1180end subroutine unpack6
1181
1197subroutine getdim(csec3, lcsec3, width, height, iscan)
1198 implicit none
1199
1200 character(len=1),intent(in) :: csec3(*)
1201 integer,intent(in) :: lcsec3
1202 integer,intent(out) :: width,height,iscan
1203
1204 integer,pointer,dimension(:) :: igdstmpl,list_opt
1205 integer :: igds(5)
1206 integer iofst, igdtlen, num_opt, jerr
1207
1208 interface
1209 subroutine gf_unpack3(cgrib, lcgrib, iofst, igds, igdstmpl, &
1210 mapgridlen, ideflist, idefnum, ierr)
1211 character(len = 1), intent(in) :: cgrib(lcgrib)
1212 integer, intent(in) :: lcgrib
1213 integer, intent(inout) :: iofst
1214 integer, pointer, dimension(:) :: igdstmpl, ideflist
1215 integer, intent(out) :: igds(5)
1216 integer, intent(out) :: ierr, idefnum
1217 end subroutine gf_unpack3
1218 end interface
1219
1220 nullify(igdstmpl, list_opt)
1221
1222 iofst = 0 ! set offset to beginning of section
1223 call gf_unpack3(csec3, lcsec3, iofst, igds, igdstmpl, &
1224 igdtlen, list_opt, num_opt, jerr)
1225 if (jerr.eq.0) then
1226 selectcase( igds(5) ) ! Template number
1227 case (0:3) ! Lat/Lon
1228 width = igdstmpl(8)
1229 height = igdstmpl(9)
1230 iscan = igdstmpl(19)
1231 case (10) ! Mercator
1232 width = igdstmpl(8)
1233 height = igdstmpl(9)
1234 iscan = igdstmpl(16)
1235 case (20) ! Polar Stereographic
1236 width = igdstmpl(8)
1237 height = igdstmpl(9)
1238 iscan = igdstmpl(18)
1239 case (30) ! Lambert Conformal
1240 width = igdstmpl(8)
1241 height = igdstmpl(9)
1242 iscan = igdstmpl(18)
1243 case (40:43) ! Gaussian
1244 width = igdstmpl(8)
1245 height = igdstmpl(9)
1246 iscan = igdstmpl(19)
1247 case (90) ! Space View/Orthographic
1248 width = igdstmpl(8)
1249 height = igdstmpl(9)
1250 iscan = igdstmpl(17)
1251 case (110) ! Equatorial Azimuthal
1252 width = igdstmpl(8)
1253 height = igdstmpl(9)
1254 iscan = igdstmpl(16)
1255 case default
1256 width = 0
1257 height = 0
1258 iscan = 0
1259 end select
1260 else
1261 width = 0
1262 height = 0
1263 endif
1264
1265 if (associated(igdstmpl)) deallocate(igdstmpl)
1266 if (associated(list_opt)) deallocate(list_opt)
1267
1268 return
1269end subroutine getdim
1270
1298subroutine getlocal(cgrib, lcgrib, localnum, csec2, lcsec2, ierr)
1299 implicit none
1300
1301 character(len = 1), intent(in) :: cgrib(lcgrib)
1302 integer, intent(in) :: lcgrib, localnum
1303 character(len = 1), intent(out) :: csec2(*)
1304 integer, intent(out) :: lcsec2, ierr
1305
1306 character(len = 4), parameter :: grib = 'GRIB', c7777 = '7777'
1307 character(len = 4) :: ctemp
1308 integer :: listsec0(2)
1309 integer iofst, istart, numlocal
1310 integer :: lengrib, lensec, lensec0, j, ipos, isecnum
1311
1312 ierr = 0
1313 numlocal = 0
1314
1315 ! Check for valid request number.
1316 if (localnum .le. 0) then
1317 print *, 'getlocal: Request for local section must be positive.'
1318 ierr = 3
1319 return
1320 endif
1321
1322 ! Check for beginning of GRIB message in the first 100 bytes
1323 istart = 0
1324 do j = 1, 100
1325 ctemp = cgrib(j)//cgrib(j+1)//cgrib(j+2)//cgrib(j+3)
1326 if (ctemp .eq. grib) then
1327 istart = j
1328 exit
1329 endif
1330 enddo
1331 if (istart .eq. 0) then
1332 print *, 'getlocal: Beginning characters GRIB not found.'
1333 ierr = 1
1334 return
1335 endif
1336
1337 ! Unpack Section 0 - Indicator Section
1338 iofst = 8 * (istart + 5)
1339 call g2_gbytec(cgrib, listsec0(1), iofst, 8) ! Discipline
1340 iofst = iofst + 8
1341 call g2_gbytec(cgrib, listsec0(2), iofst, 8) ! GRIB edition number
1342 iofst = iofst + 8
1343 iofst = iofst + 32
1344 call g2_gbytec1(cgrib, lengrib, iofst, 32) ! Length of GRIB message
1345 iofst = iofst + 32
1346 lensec0 = 16
1347 ipos = istart + lensec0
1348
1349 ! Currently handles only GRIB Edition 2.
1350 if (listsec0(2) .ne. 2) then
1351 print *, 'getlocal: can only decode GRIB edition 2.'
1352 ierr = 2
1353 return
1354 endif
1355
1356 ! Loop through the remaining sections keeping track of the length of
1357 ! each. Also check to see that if the current occurrence of Section
1358 ! 2 is the same as the one requested.
1359 do
1360 ! Check to see if we are at end of GRIB message
1361 ctemp = cgrib(ipos) // cgrib(ipos + 1) // cgrib(ipos + 2) // cgrib(ipos + 3)
1362 if (ctemp .eq. c7777) then
1363 ipos = ipos + 4
1364
1365 ! If end of GRIB message not where expected, issue error
1366 if (ipos .ne. (istart + lengrib)) then
1367 print *, 'getlocal: "7777" found, but not where expected.'
1368 ierr = 4
1369 return
1370 endif
1371 exit
1372 endif
1373
1374 ! Get length of Section and Section number
1375 iofst = (ipos - 1) * 8
1376 call g2_gbytec1(cgrib, lensec, iofst, 32) ! Get Length of Section
1377 iofst = iofst + 32
1378 call g2_gbytec1(cgrib, isecnum, iofst, 8) ! Get Section number
1379 iofst = iofst + 8
1380
1381 ! If found the requested occurrence of Section 2,
1382 ! return the section contents.
1383 if (isecnum .eq. 2) then
1384 numlocal = numlocal + 1
1385 if (numlocal.eq.localnum) then
1386 lcsec2 = lensec - 5
1387 csec2(1:lcsec2) = cgrib(ipos + 5:ipos + lensec - 1)
1388 return
1389 endif
1390 endif
1391
1392 ! Check to see if we read pass the end of the GRIB
1393 ! message and missed the terminator string '7777'.
1394 ipos = ipos + lensec ! Update beginning of section pointer
1395 if (ipos .gt. (istart + lengrib)) then
1396 print *, 'getlocal: "7777" not found at end of GRIB message.'
1397 ierr = 5
1398 return
1399 endif
1400 enddo
1401
1402 ! If exited from above loop, the end of the GRIB message was reached
1403 ! before the requested occurrence of section 2 was found.
1404 print *, 'getlocal: GRIB message contained ', numlocal, ' local sections.'
1405 print *, 'getlocal: The request was for the ', localnum, ' occurrence.'
1406 ierr = 6
1407
1408end subroutine getlocal
1409
1424subroutine getpoly(csec3, lcsec3, jj, kk, mm)
1425 implicit none
1426
1427 character(len = 1), intent(in) :: csec3(*)
1428 integer, intent(in) :: lcsec3
1429 integer, intent(out) :: jj, kk, mm
1430
1431 integer, pointer, dimension(:) :: igdstmpl, list_opt
1432 integer :: igds(5)
1433 integer iofst, igdtlen, num_opt, jerr
1434
1435 interface
1436 subroutine gf_unpack3(cgrib, lcgrib, iofst, igds, igdstmpl, &
1437 mapgridlen, ideflist, idefnum, ierr)
1438 character(len = 1), intent(in) :: cgrib(lcgrib)
1439 integer, intent(in) :: lcgrib
1440 integer, intent(inout) :: iofst
1441 integer, pointer, dimension(:) :: igdstmpl, ideflist
1442 integer, intent(out) :: igds(5)
1443 integer, intent(out) :: ierr, idefnum
1444 end subroutine gf_unpack3
1445 end interface
1446
1447 nullify(igdstmpl, list_opt)
1448
1449 iofst = 0 ! set offset to beginning of section
1450 call gf_unpack3(csec3, lcsec3, iofst, igds, igdstmpl, &
1451 igdtlen, list_opt, num_opt, jerr)
1452 if (jerr.eq.0) then
1453 selectcase( igds(5) ) ! Template number
1454 case (50:53) ! Spherical harmonic coefficients
1455 jj = igdstmpl(1)
1456 kk = igdstmpl(2)
1457 mm = igdstmpl(3)
1458 case default
1459 jj = 0
1460 kk = 0
1461 mm = 0
1462 end select
1463 else
1464 jj = 0
1465 kk = 0
1466 mm = 0
1467 endif
1468
1469 if (associated(igdstmpl)) deallocate(igdstmpl)
1470 if (associated(list_opt)) deallocate(list_opt)
1471
1472end subroutine getpoly
1473
1534subroutine gettemplates(cgrib, lcgrib, ifldnum, igds, igdstmpl, &
1535 igdslen, ideflist, idefnum, ipdsnum, ipdstmpl, ipdslen, coordlist, numcoord, ierr)
1536 implicit none
1537
1538 character(len = 1), intent(in) :: cgrib(lcgrib)
1539 integer, intent(in) :: lcgrib, ifldnum
1540 integer, intent(out) :: igds(*), igdstmpl(*), ideflist(*)
1541 integer, intent(out) :: ipdsnum, ipdstmpl(*)
1542 integer, intent(out) :: idefnum, numcoord
1543 integer, intent(out) :: ierr
1544 real, intent(out) :: coordlist(*)
1545
1546 character(len = 4), parameter :: grib = 'GRIB', c7777 = '7777'
1547 character(len = 4) :: ctemp
1548 integer:: listsec0(2)
1549 integer iofst, istart
1550 logical have3, have4
1551 integer :: igdslen, ipdslen, ipos, isecnum, j, jerr, lengrib, lensec, lensec0, numfld
1552
1553 interface
1554 subroutine g2_gbytec1(in, siout, iskip, nbits)
1555 character*1, intent(in) :: in(*)
1556 integer, intent(inout) :: siout
1557 integer, intent(in) :: iskip, nbits
1558 end subroutine g2_gbytec1
1559 end interface
1560
1561 have3 = .false.
1562 have4 = .false.
1563 ierr = 0
1564 numfld = 0
1565
1566 ! Check for valid request number.
1567 if (ifldnum .le. 0) then
1568 print *, 'gettemplates: Request for field number must be positive.'
1569 ierr = 3
1570 return
1571 endif
1572
1573 ! Check for beginning of GRIB message in the first 100 bytes.
1574 istart = 0
1575 do j = 1, 100
1576 ctemp = cgrib(j) // cgrib(j + 1) // cgrib(j + 2) // cgrib(j + 3)
1577 if (ctemp .eq. grib ) then
1578 istart = j
1579 exit
1580 endif
1581 enddo
1582 if (istart .eq. 0) then
1583 print *, 'gettemplates: Beginning characters GRIB not found.'
1584 ierr = 1
1585 return
1586 endif
1587
1588 ! Unpack Section 0 - Indicator Section.
1589 iofst = 8 * (istart + 5)
1590 call g2_gbytec1(cgrib, listsec0(1), iofst, 8) ! Discipline
1591 iofst = iofst + 8
1592 call g2_gbytec1(cgrib, listsec0(2), iofst, 8) ! GRIB edition number
1593 iofst = iofst + 8
1594 iofst = iofst + 32
1595 call g2_gbytec1(cgrib, lengrib, iofst, 32) ! Length of GRIB message
1596 iofst = iofst + 32
1597 lensec0 = 16
1598 ipos = istart + lensec0
1599
1600 ! Currently handles only GRIB Edition 2.
1601 if (listsec0(2) .ne. 2) then
1602 print *, 'gettemplates: can only decode GRIB edition 2.'
1603 ierr = 2
1604 return
1605 endif
1606
1607 ! Loop through the remaining sections keeping track of the length of
1608 ! each. Also keep the latest Grid Definition Section info. Unpack
1609 ! the requested field number.
1610 do
1611 ! Check to see if we are at end of GRIB message.
1612 ctemp = cgrib(ipos) // cgrib(ipos + 1) // cgrib(ipos + 2) // cgrib(ipos + 3)
1613 if (ctemp .eq. c7777 ) then
1614 ipos = ipos + 4
1615 ! If end of GRIB message not where expected, issue error
1616 if (ipos .ne. (istart + lengrib)) then
1617 print *, 'gettemplates: "7777" found, but not where expected.'
1618 ierr = 4
1619 return
1620 endif
1621 exit
1622 endif
1623 ! Get length of Section and Section number.
1624 iofst = (ipos - 1) * 8
1625 call g2_gbytec1(cgrib, lensec, iofst, 32) ! Get Length of Section
1626 iofst = iofst + 32
1627 call g2_gbytec1(cgrib, isecnum, iofst, 8) ! Get Section number
1628 iofst = iofst + 8
1629 !print *, ' lensec = ', lensec, ' secnum = ', isecnum
1630
1631 ! If found Section 3, unpack the GDS info using the appropriate
1632 ! template. Save in case this is the latest grid before the
1633 ! requested field.
1634 if (isecnum .eq. 3) then
1635 iofst = iofst - 40 ! reset offset to beginning of section
1636 call unpack3(cgrib, lcgrib, iofst, igds, igdstmpl, igdslen, ideflist, idefnum, jerr)
1637 if (jerr .eq. 0) then
1638 have3 = .true.
1639 else
1640 ierr = 10
1641 return
1642 endif
1643 endif
1644
1645 ! If found Section 4, check to see if this field is the one
1646 ! requested.
1647 if (isecnum .eq. 4) then
1648 numfld = numfld + 1
1649 if (numfld .eq. ifldnum) then
1650 iofst = iofst - 40 ! reset offset to beginning of section
1651 call unpack4(cgrib, lcgrib, iofst, ipdsnum, ipdstmpl, ipdslen, coordlist, numcoord, jerr)
1652 if (jerr .eq. 0) then
1653 have4 = .true.
1654 else
1655 ierr = 11
1656 return
1657 endif
1658 endif
1659 endif
1660
1661 ! Check to see if we read pass the end of the GRIB message and
1662 ! missed the terminator string '7777'.
1663 ipos = ipos + lensec ! Update beginning of section pointer
1664 if (ipos .gt. (istart + lengrib)) then
1665 print *, 'gettemplates: "7777" not found at end of GRIB message.'
1666 ierr = 7
1667 return
1668 endif
1669
1670 if (have3 .and. have4) return
1671 enddo
1672
1673 ! If exited from above loop, the end of the GRIB message was reached
1674 ! before the requested field was found.
1675 print *, 'gettemplates: GRIB message contained ', numfld, ' different fields.'
1676 print *, 'gettemplates: The request was for the ', ifldnum, ' field.'
1677 ierr = 6
1678
1679end subroutine gettemplates
subroutine comunpack(cpack, len, lensec, idrsnum, idrstmpl, ndpts, fld, ier)
Unpack a data field that was packed using a complex packing algorithm as defined in the GRIB2 documen...
Definition compack.F90:528
subroutine aecunpack(cpack, len, idrstmpl, ndpts, fld)
Unpack a data field from a AEC code stream as defined in Data Representation Template 5....
Definition g2aec.F90:95
subroutine g2_gbytesc(in, iout, iskip, nbits, nskip, n)
Extract arbitrary size big-endian integer values (up to 32 bits each) from a packed bit string.
Definition g2bytes.F90:120
subroutine rdieee(rieee, a, num)
Copy array of 32-bit IEEE floating point values to local floating point representation.
Definition g2bytes.F90:637
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 g2_gbytec81(in, siout, iskip, nbits)
Extract one arbitrary size big-endian integer value (up to 64 bits) from a packed bit string into a s...
Definition g2bytes.F90:209
subroutine g2_gbytescr(in, rout, iskip, nbits, nskip, n)
Extract big-endian floating-point values (32 bits each) from a packed bit string.
Definition g2bytes.F90:86
subroutine g2_gbytec1(in, siout, iskip, nbits)
Extract one arbitrary size big-endian integer value (up to 32 bits) from a packed bit string into a s...
Definition g2bytes.F90:52
subroutine getfield(cgrib, lcgrib, ifldnum, igds, igdstmpl, igdslen, ideflist, idefnum, ipdsnum, ipdstmpl, ipdslen, coordlist, numcoord, ndpts, idrsnum, idrstmpl, idrslen, ibmap, bmap, fld, ierr)
Return the Grid Definition, Product Definition, Bit-map (if applicable), and the unpacked data for a ...
Definition g2get.F90:547
subroutine gb_info(cgrib, lcgrib, listsec0, listsec1, numfields, numlocal, maxlocal, ierr)
Find the number of gridded fields and Local Use Sections in a GRIB2 message.
Definition g2get.F90:54
subroutine gettemplates(cgrib, lcgrib, ifldnum, igds, igdstmpl, igdslen, ideflist, idefnum, ipdsnum, ipdstmpl, ipdslen, coordlist, numcoord, ierr)
Return the Grid Definition and Product Definition for a data field.
Definition g2get.F90:1536
subroutine unpack6(cgrib, lcgrib, iofst, ngpts, ibmap, bmap, ierr)
This subroutine unpacks Section 6 (Bit-Map Section) starting at octet 6 of that Section.
Definition g2get.F90:1144
subroutine unpack4(cgrib, lcgrib, iofst, ipdsnum, ipdstmpl, mappdslen, coordlist, numcoord, ierr)
This subroutine unpacks Section 4 (Product Definition Section) starting at octet 6 of that Section.
Definition g2get.F90:944
subroutine getlocal(cgrib, lcgrib, localnum, csec2, lcsec2, ierr)
Return the contents of Section 2 (Local Use) from a GRIB2 message.
Definition g2get.F90:1299
subroutine unpack3(cgrib, lcgrib, iofst, igds, igdstmpl, mapgridlen, ideflist, idefnum, ierr)
This subroutine unpacks Section 3 (Grid Definition Section) starting at octet 6 of that Section.
Definition g2get.F90:819
subroutine getdim(csec3, lcsec3, width, height, iscan)
Return the dimensions and scanning mode of a grid definition packed in the Grid Definition Section.
Definition g2get.F90:1198
subroutine getpoly(csec3, lcsec3, jj, kk, mm)
Return the J, K, and M pentagonal resolution parameters specified in a GRIB2 Grid Definition Section ...
Definition g2get.F90:1425
subroutine unpack5(cgrib, lcgrib, iofst, ndpts, idrsnum, idrstmpl, mapdrslen, ierr)
This subroutine unpacks Section 5 (Data Representation Section) starting at octet 6 of that Section.
Definition g2get.F90:1051
subroutine gribinfo(cgrib, lcgrib, listsec0, listsec1, numlocal, numfields, maxvals, ierr)
Find the number of Local Use Sections and gridded fields in a GRIB2 message, and the maximum sizes of...
Definition g2get.F90:284
subroutine jpcunpack(cpack, len, idrstmpl, ndpts, fld)
Unpack a data field from a JPEG2000 code stream as defined in Data Representation Template 5....
Definition g2jpc.F90:96
subroutine pngunpack(cpack, len, idrstmpl, ndpts, fld)
Unpack a data field with PNG, defined in [Data Representation Template 5.40](https://www....
Definition g2png.F90:95
subroutine simunpack(cpack, len, idrstmpl, ndpts, fld)
Unpack a data field that was packed using a simple packing, [Data Representation Template 5....
Definition g2sim.F90:169
subroutine gf_unpack3(cgrib, lcgrib, iofst, igds, igdstmpl, mapgridlen, ideflist, idefnum, ierr)
Unpack Section 3 ([Grid Definition Section] (https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_doc/g...
Definition g2unpack.F90:184
Handles Data Representation Templates used in Section 5.
subroutine getdrstemplate(number, nummap, map, needext, iret)
Return DRS template information for a specified Data Representation Template.
subroutine extdrstemplate(number, list, nummap, map)
Generate the remaining octet map for a given Data Representation Template, if required.
Logging for the g2 library.
Definition g2logging.F90:10
character *120 g2_log_msg
Fill this with the logging message.
Definition g2logging.F90:12
subroutine g2_log(level)
Print a debug message for the g2 library.
Definition g2logging.F90:22
This Fortran module contains info on all the available GRIB2 Grid Definition Templates used in [Secti...
subroutine getgridtemplate(number, nummap, map, needext, iret)
Get the grid template information for a specified Grid Definition Template.
subroutine extgridtemplate(number, list, nummap, map)
Generate the remaining octet map for a given Grid Definition Template, if required.
Information on all GRIB2 Product Definition Templates used in Section 4 - the Product Definition Sect...
subroutine extpdstemplate(number, list, nummap, map)
This subroutine generates the remaining octet map for a given Product Definition Template,...
subroutine getpdstemplate(number, nummap, map, needext, iret)
This subroutine returns PDS template information for a specified Product Definition Template.