NCEPLIBS-g2  3.4.9
g2get.F90
Go to the documentation of this file.
1 
4 
52 subroutine gb_info(cgrib, lcgrib, listsec0, listsec1, &
53  numfields, numlocal, maxlocal, ierr)
54  implicit none
55 
56  character(len = 1), intent(in) :: cgrib(lcgrib)
57  integer, intent(in) :: lcgrib
58  integer, intent(out) :: listsec0(3), listsec1(13)
59  integer, intent(out) :: numlocal, numfields, maxlocal, ierr
60 
61  character(len = 4), parameter :: grib = 'GRIB', c7777 = '7777'
62  character(len = 4) :: ctemp
63  integer, parameter :: zero = 0, one = 1
64  integer, parameter :: mapsec1len = 13
65  integer, parameter :: mapsec1(mapsec1len) = (/ 2, 2, 1, 1, 1, 2, 1, 1, 1, 1, 1, 1, 1 /)
66  integer :: iofst, istart
67  integer :: nbits, lensec1, lensec0, lensec, lenposs, lengrib, j
68  integer :: i, ipos, isecnum
69 
70  ierr = 0
71  numlocal = 0
72  numfields = 0
73  maxlocal = 0
74 
75  ! Check for beginning of GRIB message in the first 100 bytes.
76  istart = 0
77  do j = 1, 100
78  ctemp = cgrib(j) // cgrib(j + 1) // cgrib(j + 2) // cgrib(j+3)
79  if (ctemp .eq. grib ) then
80  istart = j
81  exit
82  endif
83  enddo
84  if (istart .eq. 0) then
85  print *, 'gb_info: Beginning characters GRIB not found.'
86  ierr = 1
87  return
88  endif
89 
90  ! Unpack Section 0 - Indicator Section.
91  iofst = 8 * (istart + 5)
92  call g2_gbytec(cgrib, listsec0(1), iofst, 8) ! Discipline
93  iofst = iofst + 8
94  call g2_gbytec(cgrib, listsec0(2), iofst, 8) ! GRIB edition number
95  iofst = iofst+8
96  iofst = iofst + 32
97  call g2_gbytec(cgrib, lengrib, iofst, 32) ! Length of GRIB message
98  iofst = iofst + 32
99  listsec0(3) = lengrib
100  lensec0 = 16
101  ipos = istart + lensec0
102 
103  ! Currently handles only GRIB Edition 2.
104  if (listsec0(2) .ne. 2) then
105  print *, 'gb_info: can only decode GRIB edition 2.'
106  ierr = 2
107  return
108  endif
109 
110  ! Unpack Section 1 - Identification Section.
111  call g2_gbytec(cgrib, lensec1, iofst, 32) ! Length of Section 1
112  iofst = iofst + 32
113  call g2_gbytec(cgrib, isecnum, iofst, 8) ! Section number ( 1 )
114  iofst = iofst + 8
115  if (isecnum .ne. 1) then
116  print *, 'gb_info: Could not find section 1.'
117  ierr = 3
118  return
119  endif
120 
121  ! Unpack each input value in array listsec1 into the the appropriate
122  ! number of octets, which are specified in corresponding entries in
123  ! array mapsec1.
124  do i = 1, mapsec1len
125  nbits = mapsec1(i) * 8
126  call g2_gbytec(cgrib, listsec1(i), iofst, nbits)
127  iofst = iofst + nbits
128  enddo
129  ipos = ipos + lensec1
130 
131  ! Loop through the remaining sections to see if they are valid.
132  ! Also count the number of times Section 2 and Section 4 appear.
133  do
134  ctemp = cgrib(ipos) // cgrib(ipos + 1) // cgrib(ipos + 2) // cgrib(ipos + 3)
135  if (ctemp .eq. c7777 ) then
136  ipos = ipos + 4
137  if (ipos .ne. (istart + lengrib)) then
138  print *, 'gb_info: "7777" found, but not where expected.'
139  ierr = 4
140  return
141  endif
142  exit
143  endif
144  iofst = (ipos - 1) * 8
145  call g2_gbytec(cgrib, lensec, iofst, 32) ! Get Length of Section
146  iofst = iofst + 32
147  call g2_gbytec(cgrib, isecnum, iofst, 8) ! Get Section number
148  iofst = iofst + 8
149  ipos = ipos+lensec ! Update beginning of section pointer
150  if (ipos .gt. (istart + lengrib)) then
151  print *, 'gb_info: "7777" not found at end of GRIB message.'
152  ierr = 5
153  return
154  endif
155  if (isecnum .ge. 2 .AND. isecnum .le. 7) then
156  if (isecnum .eq. 2) then ! Local Section 2
157  ! Increment counter for total number of local sections found.
158  numlocal = numlocal + 1
159  lenposs = lensec - 5
160  if (lenposs .gt. maxlocal) maxlocal = lenposs
161  elseif (isecnum .eq. 4) then
162  ! Increment counter for total number of fields found.
163  numfields = numfields + 1
164  endif
165  else
166  print *, 'gb_info: Invalid section number found in GRIB message: ', isecnum
167  ierr = 6
168  return
169  endif
170  enddo
171 
172 end subroutine gb_info
173 
246 subroutine gribinfo(cgrib, lcgrib, listsec0, listsec1, &
247  numlocal, numfields, maxvals, ierr)
248  implicit none
249 
250  character(len = 1), intent(in) :: cgrib(lcgrib)
251  integer, intent(in) :: lcgrib
252  integer, intent(out) :: listsec0(3), listsec1(13), maxvals(7)
253  integer, intent(out) :: numlocal, numfields, ierr
254 
255  integer :: i, ipos, isecnum, j, lengrib, lenposs, lensec, lensec0, lensec1
256  integer :: maxcoordlist, maxdeflist, maxdrstmpl, maxgridpts, maxpdstmpl, maxgdstmpl
257  integer :: maxsec2len, nbits, nbyte, ngdpts, numcoord
258  character(len = 4), parameter :: grib = 'GRIB', c7777 = '7777'
259  character(len = 4) :: ctemp
260  integer, parameter :: zero = 0, one = 1
261  integer, parameter :: mapsec1len = 13
262  integer, parameter :: mapsec1(mapsec1len) = (/ 2, 2, 1, 1, 1, 2, 1, 1, &
263  1, 1, 1, 1, 1 /)
264  integer iofst, istart
265 
266  ierr = 0
267  numlocal = 0
268  numfields = 0
269  maxsec2len = 1
270  maxgdstmpl = 1
271  maxdeflist = 1
272  maxpdstmpl = 1
273  maxcoordlist = 1
274  maxdrstmpl = 1
275  maxgridpts = 0
276 
277  ! Check for beginning of GRIB message in the first 100 bytes
278  istart = 0
279  do j = 1, 100
280  ctemp = cgrib(j) // cgrib(j + 1) // cgrib(j + 2) // cgrib(j + 3)
281  if (ctemp .eq. grib) then
282  istart = j
283  exit
284  endif
285  enddo
286  if (istart .eq. 0) then
287  print *, 'gribinfo: Beginning characters GRIB not found.'
288  ierr = 1
289  return
290  endif
291 
292  ! Unpack Section 0 - Indicator Section.
293  iofst = 8 * (istart + 5)
294  call g2_gbytec(cgrib, listsec0(1), iofst, 8) ! Discipline
295  iofst = iofst + 8
296  call g2_gbytec(cgrib, listsec0(2), iofst, 8) ! GRIB edition number
297  iofst = iofst + 8
298  iofst = iofst + 32
299  call g2_gbytec(cgrib, lengrib, iofst, 32) ! Length of GRIB message
300  iofst = iofst + 32
301  listsec0(3) = lengrib
302  lensec0 = 16
303  ipos = istart + lensec0
304 
305  ! Currently handles only GRIB Edition 2.
306  if (listsec0(2) .ne. 2) then
307  print *, 'gribinfo: can only decode GRIB edition 2.'
308  ierr = 2
309  return
310  endif
311 
312  ! Unpack Section 1 - Identification Section.
313  call g2_gbytec(cgrib, lensec1, iofst, 32) ! Length of Section 1
314  iofst = iofst + 32
315  call g2_gbytec(cgrib, isecnum, iofst, 8) ! Section number (1)
316  iofst = iofst + 8
317  if (isecnum .ne. 1) then
318  print *, 'gribinfo: Could not find section 1.'
319  ierr = 3
320  return
321  endif
322 
323  ! Unpack each input value in array listsec1 into the the appropriate
324  ! number of octets, which are specified in corresponding entries in
325  ! array mapsec1.
326  do i = 1, mapsec1len
327  nbits = mapsec1(i) * 8
328  call g2_gbytec(cgrib, listsec1(i), iofst, nbits)
329  iofst = iofst + nbits
330  enddo
331  ipos = ipos + lensec1
332 
333  ! Loop through the remaining sections keeping track of the length of
334  ! each. Also count the number of times Section 2 and Section 4
335  ! appear.
336  do
337  ctemp = cgrib(ipos) // cgrib(ipos + 1) // cgrib(ipos + 2) // &
338  cgrib(ipos + 3)
339  if (ctemp .eq. c7777) then
340  ipos = ipos + 4
341  if (ipos .ne. (istart + lengrib)) then
342  print *, 'gribinfo: "7777" found, but not where expected.'
343  ierr = 4
344  return
345  endif
346  exit
347  endif
348  iofst = (ipos - 1) * 8
349  call g2_gbytec(cgrib, lensec, iofst, 32) ! Get Length of Section
350  iofst = iofst + 32
351  call g2_gbytec(cgrib, isecnum, iofst, 8) ! Get Section number
352  iofst = iofst + 8
353  ipos = ipos + lensec ! Update beginning of section pointer
354  if (ipos .gt. (istart + lengrib)) then
355  print *, 'gribinfo: "7777" not found at end of GRIB message.'
356  ierr = 5
357  return
358  endif
359  if (isecnum .eq. 2) then ! Local Section 2
360  ! Increment counter for total number of local sections found
361  ! and determine largest Section 2 in message.
362  numlocal = numlocal + 1
363  lenposs = lensec-5
364  if (lenposs .gt. maxsec2len) maxsec2len = lenposs
365  elseif (isecnum .eq. 3) then
366  iofst = iofst + 8 ! skip source of grid def.
367  call g2_gbytec(cgrib, ngdpts, iofst, 32) ! Get Num of Grid Points
368  iofst = iofst + 32
369  call g2_gbytec(cgrib, nbyte, iofst, 8) ! Get Num octets for opt. list
370  iofst = iofst + 8
371  if (ngdpts .gt. maxgridpts) maxgridpts = ngdpts
372  lenposs = lensec - 14
373  if (lenposs .gt. maxgdstmpl) maxgdstmpl = lenposs
374  if (nbyte .ne. 0) then
375  lenposs = lenposs / nbyte
376  if (lenposs .gt. maxdeflist) maxdeflist = lenposs
377  endif
378  elseif (isecnum .eq. 4) then
379  numfields = numfields + 1
380  call g2_gbytec(cgrib, numcoord, iofst, 16) ! Get Num of Coord Values
381  iofst = iofst + 16
382  if (numcoord .ne. 0) then
383  if (numcoord .gt. maxcoordlist) maxcoordlist = numcoord
384  endif
385  lenposs = lensec - 9
386  if (lenposs.gt.maxpdstmpl) maxpdstmpl = lenposs
387  elseif (isecnum .eq. 5) then
388  lenposs = lensec-11
389  if (lenposs .gt. maxdrstmpl) maxdrstmpl = lenposs
390  endif
391  enddo
392 
393  maxvals(1) = maxsec2len
394  maxvals(2) = maxgdstmpl
395  maxvals(3) = maxdeflist
396  maxvals(4) = maxpdstmpl
397  maxvals(5) = maxcoordlist
398  maxvals(6) = maxdrstmpl
399  maxvals(7) = maxgridpts
400 end subroutine gribinfo
401 
502 subroutine getfield(cgrib, lcgrib, ifldnum, igds, igdstmpl, &
503  igdslen, ideflist, idefnum, ipdsnum, ipdstmpl, ipdslen, &
504  coordlist, numcoord, ndpts, idrsnum, idrstmpl, idrslen, &
505  ibmap, bmap, fld, ierr)
506 
507  implicit none
508 
509  character(len = 1), intent(in) :: cgrib(lcgrib)
510  integer, intent(in) :: lcgrib, ifldnum
511  integer, intent(out) :: igds(*), igdstmpl(*), ideflist(*)
512  integer, intent(out) :: ipdsnum, ipdstmpl(*)
513  integer, intent(out) :: idrsnum, idrstmpl(*)
514  integer, intent(out) :: ndpts, ibmap, idefnum, numcoord
515  integer, intent(out) :: ierr
516  logical*1, intent(out) :: bmap(*)
517  real, intent(out) :: fld(*), coordlist(*)
518 
519  character(len = 4), parameter :: grib = 'GRIB', c7777 = '7777'
520  character(len = 4) :: ctemp
521  integer :: listsec0(2)
522  integer :: iofst, istart
523  real (kind = 4) :: ieee(1)
524  logical :: have3, have4, have5, have6, have7
525 
526  !implicit none additions
527  integer, intent(out) :: igdslen, ipdslen, idrslen
528  integer :: numfld, j, lengrib, lensec0, ipos
529  integer :: lensec, isecnum, jerr, ier, numlocal
530 
531  have3 = .false.
532  have4 = .false.
533  have5 = .false.
534  have6 = .false.
535  have7 = .false.
536  ierr = 0
537  numfld = 0
538  numlocal = 0
539 
540  ! Check for valid request number
541  if (ifldnum .le. 0) then
542  print *, 'getfield: Request for field number ' &
543  ,'must be positive.'
544  ierr = 3
545  return
546  endif
547 
548  ! Check for beginning of GRIB message in the first 100 bytes
549  istart = 0
550  do j = 1, 100
551  ctemp = cgrib(j) // cgrib(j + 1) // cgrib(j + 2) // cgrib(j + 3)
552  if (ctemp .eq. grib) then
553  istart = j
554  exit
555  endif
556  enddo
557  if (istart .eq. 0) then
558  print *, 'getfield: Beginning characters GRIB not found.'
559  ierr = 1
560  return
561  endif
562 
563  ! Unpack Section 0 - Indicator Section
564  iofst = 8 * (istart + 5)
565  call g2_gbytec(cgrib, listsec0(1), iofst, 8) ! Discipline
566  iofst = iofst + 8
567  call g2_gbytec(cgrib, listsec0(2), iofst, 8) ! GRIB edition number
568  iofst = iofst + 8
569  iofst = iofst + 32
570  call g2_gbytec(cgrib, lengrib, iofst, 32) ! Length of GRIB message
571  iofst = iofst + 32
572  lensec0 = 16
573  ipos = istart + lensec0
574 
575  ! Currently handles only GRIB Edition 2.
576  if (listsec0(2) .ne. 2) then
577  print *, 'getfield: can only decode GRIB edition 2.'
578  ierr = 2
579  return
580  endif
581 
582  ! Loop through the remaining sections keeping track of the length of
583  ! each. Also keep the latest Grid Definition Section info. Unpack
584  ! the requested field number.
585  do
586  ! Check to see if we are at end of GRIB message
587  ctemp = cgrib(ipos) // cgrib(ipos + 1) // cgrib(ipos + 2) // &
588  cgrib(ipos + 3)
589  if (ctemp .eq. c7777) then
590  ipos = ipos + 4
591  ! If end of GRIB message not where expected, issue error
592  if (ipos.ne.(istart + lengrib)) then
593  print *, 'getfield: "7777" found, but not ' &
594  ,'where expected.'
595  ierr = 4
596  return
597  endif
598  exit
599  endif
600  ! Get length of Section and Section number
601  iofst = (ipos - 1) * 8
602  call g2_gbytec(cgrib, lensec, iofst, 32) ! Get Length of Section
603  iofst = iofst + 32
604  call g2_gbytec(cgrib, isecnum, iofst, 8) ! Get Section number
605  iofst = iofst + 8
606 
607  ! If found Section 3, unpack the GDS info using the appropriate
608  ! template. Save in case this is the latest grid before the
609  ! requested field.
610  if (isecnum .eq. 3) then
611  iofst = iofst - 40 ! reset offset to beginning of section
612  call unpack3(cgrib, lcgrib, iofst, igds, igdstmpl, &
613  igdslen, ideflist, idefnum, jerr)
614  if (jerr .eq. 0) then
615  have3 = .true.
616  else
617  ierr = 10
618  return
619  endif
620  endif
621 
622  ! If found Section 4, check to see if this field is the one
623  ! requested.
624  if (isecnum .eq. 4) then
625  numfld = numfld + 1
626  if (numfld .eq. ifldnum) then
627  iofst = iofst - 40 ! reset offset to beginning of section
628  call unpack4(cgrib, lcgrib, iofst, ipdsnum, ipdstmpl, &
629  ipdslen, coordlist, numcoord, jerr)
630  if (jerr .eq. 0) then
631  have4 = .true.
632  else
633  ierr = 11
634  return
635  endif
636  endif
637  endif
638 
639  ! If found Section 5, check to see if this field is the one
640  ! requested.
641  if ((isecnum .eq. 5) .and. (numfld .eq. ifldnum)) then
642  iofst = iofst - 40 ! reset offset to beginning of section
643  call unpack5(cgrib, lcgrib, iofst, ndpts, idrsnum, &
644  idrstmpl, idrslen, jerr)
645  if (jerr .eq. 0) then
646  have5 = .true.
647  else
648  ierr = 12
649  return
650  endif
651  endif
652 
653  ! If found Section 6, Unpack bitmap. Save in case this is the
654  ! latest bitmap before the requested field.
655  if (isecnum .eq. 6) then
656  iofst = iofst - 40 ! reset offset to beginning of section
657  call unpack6(cgrib, lcgrib, iofst, igds(2), ibmap, bmap, &
658  jerr)
659  if (jerr .eq. 0) then
660  have6 = .true.
661  else
662  ierr = 13
663  return
664  endif
665  endif
666 
667  ! If found Section 7, check to see if this field is the one
668  ! requested.
669  if ((isecnum .eq. 7) .and. (numfld .eq. ifldnum)) then
670  if (idrsnum .eq. 0) then
671  call simunpack(cgrib(ipos + 5), lensec - 6, idrstmpl, &
672  ndpts, fld)
673  have7 = .true.
674  elseif (idrsnum .eq. 2 .or. idrsnum .eq. 3) then
675  call comunpack(cgrib(ipos + 5), lensec - 6, lensec, &
676  idrsnum,idrstmpl, ndpts, fld, ier)
677  if (ier .ne. 0) then
678  ierr = 14
679  return
680  endif
681  have7 = .true.
682  elseif (idrsnum .eq. 50) then
683  call simunpack(cgrib(ipos + 5), lensec - 6, idrstmpl, &
684  ndpts - 1, fld(2))
685  ieee = transfer(idrstmpl(5), ieee, 1)
686  call rdieee(ieee, fld(1), 1)
687  have7 = .true.
688  elseif (idrsnum .eq. 40 .or. idrsnum .eq. 40000) then
689  call jpcunpack(cgrib(ipos + 5), lensec - 5, idrstmpl, &
690  ndpts, fld)
691  have7 = .true.
692  elseif (idrsnum .eq. 41 .or. idrsnum .eq. 40010) then
693  call pngunpack(cgrib(ipos + 5), lensec - 5, idrstmpl, &
694  ndpts, fld)
695  have7 = .true.
696  else
697  print *, 'getfield: Data Representation Template ', &
698  idrsnum, ' not yet implemented.'
699  ierr = 9
700  return
701  endif
702  endif
703 
704  ! Check to see if we read pass the end of the GRIB message and
705  ! missed the terminator string '7777'.
706  ipos = ipos + lensec ! Update beginning of section pointer
707  if (ipos .gt. (istart + lengrib)) then
708  print *, 'getfield: "7777" not found at end' &
709  ,' of GRIB message.'
710  ierr = 7
711  return
712  endif
713 
714  if (have3 .and. have4 .and. have5 .and. have6 .and. have7) &
715  return
716 
717  enddo
718 
719  ! If exited from above loop, the end of the GRIB message was reached
720  ! before the requested field was found.
721  print *, 'getfield: GRIB message contained ', numlocal, &
722  ' different fields.'
723  print *, 'getfield: The request was for the ', ifldnum, &
724  ' field.'
725  ierr = 6
726 
727 end subroutine getfield
728 
739 
770 subroutine unpack3(cgrib, lcgrib, iofst, igds, igdstmpl, &
771  mapgridlen, ideflist, idefnum, ierr)
772 
773  use gridtemplates
774  implicit none
775 
776  character(len = 1), intent(in) :: cgrib(lcgrib)
777  integer, intent(in) :: lcgrib
778  integer, intent(inout) :: iofst
779  integer, intent(out) :: igds(*), igdstmpl(*), ideflist(*)
780  integer, intent(out) :: ierr, idefnum
781 
782  integer, allocatable :: mapgrid(:)
783  integer :: mapgridlen, ibyttem
784  logical needext
785 
786  !implicit none additions
787  integer :: lensec, iret, i, nbits, isign, newmapgridlen
788 
789  ierr = 0
790 
791  call g2_gbytec(cgrib, lensec, iofst, 32) ! Get Length of Section
792  iofst = iofst + 32
793  iofst = iofst + 8 ! skip section number
794 
795  call g2_gbytec(cgrib, igds(1), iofst, 8) ! Get source of Grid def.
796  iofst = iofst + 8
797  call g2_gbytec(cgrib, igds(2), iofst, 32) ! Get number of grid pts.
798  iofst = iofst + 32
799  call g2_gbytec(cgrib, igds(3), iofst, 8) ! Get num octets for opt. list
800  iofst = iofst + 8
801  call g2_gbytec(cgrib, igds(4), iofst, 8) ! Get interpret. for opt. list
802  iofst = iofst + 8
803  call g2_gbytec(cgrib, igds(5), iofst, 16) ! Get Grid Def Template num.
804  iofst = iofst + 16
805  if (igds(1) .eq. 0) then
806  ! if (igds(1).eq.0.OR.igds(1).eq.255) then ! FOR ECMWF TEST ONLY
807  allocate(mapgrid(lensec))
808  ! Get Grid Definition Template
809  call getgridtemplate(igds(5), mapgridlen, mapgrid, needext, &
810  iret)
811  if (iret .ne. 0) then
812  ierr = 5
813  return
814  endif
815  else
816  ! igdstmpl = -1
817  mapgridlen = 0
818  needext = .false.
819  endif
820 
821  ! Unpack each value into array igdstmpl from the the appropriate
822  ! number of octets, which are specified in corresponding entries in
823  ! array mapgrid.
824  ibyttem = 0
825  do i = 1, mapgridlen
826  nbits = iabs(mapgrid(i)) * 8
827  if (mapgrid(i) .ge. 0) then
828  call g2_gbytec(cgrib, igdstmpl(i), iofst, nbits)
829  else
830  call g2_gbytec(cgrib, isign, iofst, 1)
831  call g2_gbytec(cgrib, igdstmpl(i), iofst + 1, nbits-1)
832  if (isign .eq. 1) igdstmpl(i) = -igdstmpl(i)
833  endif
834  iofst = iofst + nbits
835  ibyttem = ibyttem + iabs(mapgrid(i))
836  enddo
837 
838  ! Check to see if the Grid Definition Template needs to be
839  ! extended. The number of values in a specific template may vary
840  ! depending on data specified in the "static" part of the template.
841  if (needext) then
842  call extgridtemplate(igds(5), igdstmpl, newmapgridlen, &
843  mapgrid)
844  ! Unpack the rest of the Grid Definition Template
845  do i = mapgridlen + 1, newmapgridlen
846  nbits = iabs(mapgrid(i)) * 8
847  if (mapgrid(i) .ge. 0) then
848  call g2_gbytec(cgrib, igdstmpl(i), iofst, nbits)
849  else
850  call g2_gbytec(cgrib, isign, iofst, 1)
851  call g2_gbytec(cgrib, igdstmpl(i), iofst + 1, nbits - &
852  1)
853  if (isign .eq. 1) igdstmpl(i) = -igdstmpl(i)
854  endif
855  iofst = iofst + nbits
856  ibyttem = ibyttem + iabs(mapgrid(i))
857  enddo
858  mapgridlen = newmapgridlen
859  endif
860 
861  ! Unpack optional list of numbers defining number of points in each
862  ! row or column, if included. This is used for non regular grids.
863  if (igds(3) .ne. 0) then
864  nbits = igds(3) * 8
865  idefnum = (lensec - 14 - ibyttem) / igds(3)
866  call g2_gbytesc(cgrib, ideflist, iofst, nbits, 0, idefnum)
867  iofst = iofst + (nbits * idefnum)
868  else
869  idefnum = 0
870  endif
871  if (allocated(mapgrid)) deallocate(mapgrid)
872 end subroutine unpack3
873 
900 subroutine unpack4(cgrib, lcgrib, iofst, ipdsnum, ipdstmpl, &
901  mappdslen, coordlist, numcoord, ierr)
902 
903  use pdstemplates
904  implicit none
905 
906  character(len = 1), intent(in) :: cgrib(lcgrib)
907  integer, intent(in) :: lcgrib
908  integer, intent(inout) :: iofst
909  real, intent(out) :: coordlist(*)
910  integer, intent(out) :: ipdsnum, ipdstmpl(*)
911  integer, intent(out) :: ierr, numcoord
912 
913  real(4), allocatable :: coordieee(:)
914  integer, allocatable :: mappds(:)
915  integer :: mappdslen
916  logical needext
917 
918  !implicit none additions
919  integer :: lensec, iret, i, nbits, isign, newmappdslen
920 
921  ierr = 0
922 
923  call g2_gbytec(cgrib, lensec, iofst, 32) ! Get Length of Section
924  iofst = iofst + 32
925  iofst = iofst + 8 ! skip section number
926  allocate(mappds(lensec))
927 
928  call g2_gbytec(cgrib, numcoord, iofst, 16) ! Get num of coordinate values
929  iofst = iofst + 16
930  call g2_gbytec(cgrib, ipdsnum, iofst, 16) ! Get Prod. Def Template num.
931  iofst = iofst + 16
932  ! Get Product Definition Template.
933  call getpdstemplate(ipdsnum, mappdslen, mappds, needext, iret)
934  if (iret.ne.0) then
935  ierr = 5
936  return
937  endif
938 
939  ! Unpack each value into array ipdstmpl from the the appropriate
940  ! number of octets, which are specified in corresponding entries in
941  ! array mappds.
942  do i = 1, mappdslen
943  nbits = iabs(mappds(i))*8
944  if (mappds(i).ge.0) then
945  call g2_gbytec(cgrib, ipdstmpl(i), iofst, nbits)
946  else
947  call g2_gbytec(cgrib, isign, iofst, 1)
948  call g2_gbytec(cgrib, ipdstmpl(i), iofst + 1, nbits-1)
949  if (isign.eq.1) ipdstmpl(i) = -ipdstmpl(i)
950  endif
951  iofst = iofst + nbits
952  enddo
953 
954  ! Check to see if the Product Definition Template needs to be
955  ! extended. The number of values in a specific template may vary
956  ! depending on data specified in the "static" part of the template.
957  if (needext) then
958  call extpdstemplate(ipdsnum, ipdstmpl, newmappdslen, mappds)
959 
960  ! Unpack the rest of the Product Definition Template
961  do i = mappdslen + 1, newmappdslen
962  nbits = iabs(mappds(i))*8
963  if (mappds(i).ge.0) then
964  call g2_gbytec(cgrib, ipdstmpl(i), iofst, nbits)
965  else
966  call g2_gbytec(cgrib, isign, iofst, 1)
967  call g2_gbytec(cgrib, ipdstmpl(i), iofst + 1, nbits-1)
968  if (isign.eq.1) ipdstmpl(i) = -ipdstmpl(i)
969  endif
970  iofst = iofst + nbits
971  enddo
972  mappdslen = newmappdslen
973  endif
974 
975  ! Get Optional list of vertical coordinate values after the Product
976  ! Definition Template, if necessary.
977  if (numcoord .ne. 0) then
978  allocate (coordieee(numcoord))
979  call g2_gbytescr(cgrib, coordieee, iofst, 32, 0, numcoord)
980  call rdieee(coordieee, coordlist, numcoord)
981  deallocate (coordieee)
982  iofst = iofst + (32*numcoord)
983  endif
984  if (allocated(mappds)) deallocate(mappds)
985 end subroutine unpack4
986 
1009 subroutine unpack5(cgrib, lcgrib, iofst, ndpts, idrsnum, &
1010  idrstmpl, mapdrslen, ierr)
1011 
1012  use drstemplates
1013  implicit none
1014 
1015  character(len = 1), intent(in) :: cgrib(lcgrib)
1016  integer, intent(in) :: lcgrib
1017  integer, intent(inout) :: iofst
1018  integer, intent(out) :: ndpts, idrsnum, idrstmpl(*)
1019  integer, intent(out) :: ierr
1020 
1021  ! integer, allocatable :: mapdrs(:)
1022  integer, allocatable :: mapdrs(:)
1023  integer :: mapdrslen
1024  logical needext
1025 
1026  !implicit none additions
1027  integer :: lensec, i, nbits, isign, newmapdrslen, iret
1028 
1029  ierr = 0
1030 
1031  call g2_gbytec(cgrib, lensec, iofst, 32) ! Get Length of Section
1032  iofst = iofst + 32
1033  iofst = iofst + 8 ! skip section number
1034  allocate(mapdrs(lensec))
1035 
1036  call g2_gbytec(cgrib, ndpts, iofst, 32) ! Get num of data points
1037  iofst = iofst + 32
1038  call g2_gbytec(cgrib, idrsnum, iofst, 16) ! Get Data Rep Template Num.
1039  iofst = iofst + 16
1040  ! Gen Data Representation Template
1041  call getdrstemplate(idrsnum, mapdrslen, mapdrs, needext, iret)
1042  if (iret.ne.0) then
1043  ierr = 7
1044  return
1045  endif
1046 
1047  ! Unpack each value into array ipdstmpl from the the appropriate
1048  ! number of octets, which are specified in corresponding entries in
1049  ! array mappds.
1050  do i = 1, mapdrslen
1051  nbits = iabs(mapdrs(i))*8
1052  if (mapdrs(i).ge.0) then
1053  call g2_gbytec(cgrib, idrstmpl(i), iofst, nbits)
1054  else
1055  call g2_gbytec(cgrib, isign, iofst, 1)
1056  call g2_gbytec(cgrib, idrstmpl(i), iofst + 1, nbits-1)
1057  if (isign.eq.1) idrstmpl(i) = -idrstmpl(i)
1058  endif
1059  iofst = iofst + nbits
1060  enddo
1061 
1062  ! Check to see if the Data Representation Template needs to be
1063  ! extended. The number of values in a specific template may vary
1064  ! depending on data specified in the "static" part of the template.
1065  if (needext) then
1066  call extdrstemplate(idrsnum, idrstmpl, newmapdrslen, mapdrs)
1067  ! Unpack the rest of the Data Representation Template
1068  do i = mapdrslen + 1, newmapdrslen
1069  nbits = iabs(mapdrs(i))*8
1070  if (mapdrs(i).ge.0) then
1071  call g2_gbytec(cgrib, idrstmpl(i), iofst, nbits)
1072  else
1073  call g2_gbytec(cgrib, isign, iofst, 1)
1074  call g2_gbytec(cgrib, idrstmpl(i), iofst + 1, nbits - 1)
1075  if (isign.eq.1) idrstmpl(i) = -idrstmpl(i)
1076  endif
1077  iofst = iofst + nbits
1078  enddo
1079  mapdrslen = newmapdrslen
1080  endif
1081  if (allocated(mapdrs)) deallocate(mapdrs)
1082 end subroutine unpack5
1083 
1105 subroutine unpack6(cgrib, lcgrib, iofst, ngpts, ibmap, bmap, ierr)
1106  implicit none
1107 
1108  character(len = 1), intent(in) :: cgrib(lcgrib)
1109  integer, intent(in) :: lcgrib, ngpts
1110  integer, intent(inout) :: iofst
1111  integer, intent(out) :: ibmap
1112  integer, intent(out) :: ierr
1113  logical*1, intent(out) :: bmap(ngpts)
1114 
1115  integer :: intbmap(ngpts)
1116 
1117  !implicit none additions
1118  integer :: j
1119 
1120  ierr = 0
1121 
1122  iofst = iofst + 32 ! skip Length of Section
1123  iofst = iofst + 8 ! skip section number
1124 
1125  call g2_gbytec(cgrib, ibmap, iofst, 8) ! Get bit-map indicator
1126  iofst = iofst + 8
1127 
1128  if (ibmap.eq.0) then ! Unpack bitmap
1129  call g2_gbytesc(cgrib, intbmap, iofst, 1, 0, ngpts)
1130  iofst = iofst + ngpts
1131  do j = 1, ngpts
1132  bmap(j) = .true.
1133  if (intbmap(j).eq.0) bmap(j) = .false.
1134  enddo
1135  elseif (ibmap.eq.254) then ! Use previous bitmap
1136  return
1137  elseif (ibmap.eq.255) then ! No bitmap in message
1138  bmap(1:ngpts) = .true.
1139  else
1140  print *, 'unpack6: Predefined bitmap ', ibmap, &
1141  ' not recognized.'
1142  ierr = 4
1143  endif
1144 end subroutine unpack6
1145 
1161 subroutine getdim(csec3,lcsec3,width,height,iscan)
1162  implicit none
1163 
1164  character(len=1),intent(in) :: csec3(*)
1165  integer,intent(in) :: lcsec3
1166  integer,intent(out) :: width,height,iscan
1167 
1168  integer,pointer,dimension(:) :: igdstmpl,list_opt
1169  integer :: igds(5)
1170  integer iofst,igdtlen,num_opt,jerr
1171 
1172  interface
1173  subroutine gf_unpack3(cgrib,lcgrib,iofst,igds,igdstmpl, &
1174  mapgridlen,ideflist,idefnum,ierr)
1175  character(len=1),intent(in) :: cgrib(lcgrib)
1176  integer,intent(in) :: lcgrib
1177  integer,intent(inout) :: iofst
1178  integer,pointer,dimension(:) :: igdstmpl,ideflist
1179  integer,intent(out) :: igds(5)
1180  integer,intent(out) :: ierr,idefnum
1181  end subroutine gf_unpack3
1182  end interface
1183 
1184  nullify(igdstmpl,list_opt)
1185 
1186  iofst=0 ! set offset to beginning of section
1187  call gf_unpack3(csec3,lcsec3,iofst,igds,igdstmpl, &
1188  igdtlen,list_opt,num_opt,jerr)
1189  if (jerr.eq.0) then
1190  selectcase( igds(5) ) ! Template number
1191  case (0:3) ! Lat/Lon
1192  width=igdstmpl(8)
1193  height=igdstmpl(9)
1194  iscan=igdstmpl(19)
1195  case (10) ! Mercator
1196  width=igdstmpl(8)
1197  height=igdstmpl(9)
1198  iscan=igdstmpl(16)
1199  case (20) ! Polar Stereographic
1200  width=igdstmpl(8)
1201  height=igdstmpl(9)
1202  iscan=igdstmpl(18)
1203  case (30) ! Lambert Conformal
1204  width=igdstmpl(8)
1205  height=igdstmpl(9)
1206  iscan=igdstmpl(18)
1207  case (40:43) ! Gaussian
1208  width=igdstmpl(8)
1209  height=igdstmpl(9)
1210  iscan=igdstmpl(19)
1211  case (90) ! Space View/Orthographic
1212  width=igdstmpl(8)
1213  height=igdstmpl(9)
1214  iscan=igdstmpl(17)
1215  case (110) ! Equatorial Azimuthal
1216  width=igdstmpl(8)
1217  height=igdstmpl(9)
1218  iscan=igdstmpl(16)
1219  case default
1220  width=0
1221  height=0
1222  iscan=0
1223  end select
1224  else
1225  width=0
1226  height=0
1227  endif
1228 
1229  if (associated(igdstmpl)) deallocate(igdstmpl)
1230  if (associated(list_opt)) deallocate(list_opt)
1231 
1232  return
1233 end subroutine getdim
1234 
1262 subroutine getlocal(cgrib, lcgrib, localnum, csec2, lcsec2, ierr)
1263  implicit none
1264 
1265  character(len = 1), intent(in) :: cgrib(lcgrib)
1266  integer, intent(in) :: lcgrib, localnum
1267  character(len = 1), intent(out) :: csec2(*)
1268  integer, intent(out) :: lcsec2, ierr
1269 
1270  character(len = 4), parameter :: grib = 'GRIB', c7777 = '7777'
1271  character(len = 4) :: ctemp
1272  integer :: listsec0(2)
1273  integer iofst, istart, numlocal
1274  integer :: lengrib, lensec, lensec0, j, ipos, isecnum
1275 
1276  ierr = 0
1277  numlocal = 0
1278 
1279  ! Check for valid request number.
1280  if (localnum .le. 0) then
1281  print *, 'getlocal: Request for local section must be positive.'
1282  ierr = 3
1283  return
1284  endif
1285 
1286  ! Check for beginning of GRIB message in the first 100 bytes
1287  istart = 0
1288  do j = 1, 100
1289  ctemp = cgrib(j)//cgrib(j+1)//cgrib(j+2)//cgrib(j+3)
1290  if (ctemp .eq. grib) then
1291  istart = j
1292  exit
1293  endif
1294  enddo
1295  if (istart .eq. 0) then
1296  print *, 'getlocal: Beginning characters GRIB not found.'
1297  ierr = 1
1298  return
1299  endif
1300 
1301  ! Unpack Section 0 - Indicator Section
1302  iofst = 8 * (istart + 5)
1303  call g2_gbytec(cgrib, listsec0(1), iofst, 8) ! Discipline
1304  iofst = iofst + 8
1305  call g2_gbytec(cgrib, listsec0(2), iofst, 8) ! GRIB edition number
1306  iofst = iofst + 8
1307  iofst = iofst + 32
1308  call g2_gbytec(cgrib, lengrib, iofst, 32) ! Length of GRIB message
1309  iofst = iofst + 32
1310  lensec0 = 16
1311  ipos = istart + lensec0
1312 
1313  ! Currently handles only GRIB Edition 2.
1314  if (listsec0(2) .ne. 2) then
1315  print *, 'getlocal: can only decode GRIB edition 2.'
1316  ierr = 2
1317  return
1318  endif
1319 
1320  ! Loop through the remaining sections keeping track of the length of
1321  ! each. Also check to see that if the current occurrence of Section
1322  ! 2 is the same as the one requested.
1323  do
1324  ! Check to see if we are at end of GRIB message
1325  ctemp = cgrib(ipos) // cgrib(ipos + 1) // cgrib(ipos + 2) // cgrib(ipos + 3)
1326  if (ctemp .eq. c7777) then
1327  ipos = ipos + 4
1328 
1329  ! If end of GRIB message not where expected, issue error
1330  if (ipos .ne. (istart + lengrib)) then
1331  print *, 'getlocal: "7777" found, but not where expected.'
1332  ierr = 4
1333  return
1334  endif
1335  exit
1336  endif
1337 
1338  ! Get length of Section and Section number
1339  iofst = (ipos - 1) * 8
1340  call g2_gbytec(cgrib, lensec, iofst, 32) ! Get Length of Section
1341  iofst = iofst + 32
1342  call g2_gbytec(cgrib, isecnum, iofst, 8) ! Get Section number
1343  iofst = iofst + 8
1344 
1345  ! If found the requested occurrence of Section 2,
1346  ! return the section contents.
1347  if (isecnum .eq. 2) then
1348  numlocal = numlocal + 1
1349  if (numlocal.eq.localnum) then
1350  lcsec2 = lensec - 5
1351  csec2(1:lcsec2) = cgrib(ipos + 5:ipos + lensec - 1)
1352  return
1353  endif
1354  endif
1355 
1356  ! Check to see if we read pass the end of the GRIB
1357  ! message and missed the terminator string '7777'.
1358  ipos = ipos + lensec ! Update beginning of section pointer
1359  if (ipos .gt. (istart + lengrib)) then
1360  print *, 'getlocal: "7777" not found at end of GRIB message.'
1361  ierr = 5
1362  return
1363  endif
1364  enddo
1365 
1366  ! If exited from above loop, the end of the GRIB message was reached
1367  ! before the requested occurrence of section 2 was found.
1368  print *, 'getlocal: GRIB message contained ', numlocal, ' local sections.'
1369  print *, 'getlocal: The request was for the ', localnum, ' occurrence.'
1370  ierr = 6
1371 
1372 end subroutine getlocal
1373 
1388 subroutine getpoly(csec3,lcsec3,jj,kk,mm)
1389  implicit none
1390 
1391  character(len=1),intent(in) :: csec3(*)
1392  integer,intent(in) :: lcsec3
1393  integer,intent(out) :: jj,kk,mm
1394 
1395  integer,pointer,dimension(:) :: igdstmpl,list_opt
1396  integer :: igds(5)
1397  integer iofst,igdtlen,num_opt,jerr
1398 
1399  interface
1400  subroutine gf_unpack3(cgrib,lcgrib,iofst,igds,igdstmpl, &
1401  mapgridlen,ideflist,idefnum,ierr)
1402  character(len=1),intent(in) :: cgrib(lcgrib)
1403  integer,intent(in) :: lcgrib
1404  integer,intent(inout) :: iofst
1405  integer,pointer,dimension(:) :: igdstmpl,ideflist
1406  integer,intent(out) :: igds(5)
1407  integer,intent(out) :: ierr,idefnum
1408  end subroutine gf_unpack3
1409  end interface
1410 
1411  nullify(igdstmpl,list_opt)
1412 
1413  iofst=0 ! set offset to beginning of section
1414  call gf_unpack3(csec3,lcsec3,iofst,igds,igdstmpl, &
1415  igdtlen,list_opt,num_opt,jerr)
1416  if (jerr.eq.0) then
1417  selectcase( igds(5) ) ! Template number
1418  case (50:53) ! Spherical harmonic coefficients
1419  jj=igdstmpl(1)
1420  kk=igdstmpl(2)
1421  mm=igdstmpl(3)
1422  case default
1423  jj=0
1424  kk=0
1425  mm=0
1426  end select
1427  else
1428  jj=0
1429  kk=0
1430  mm=0
1431  endif
1432 
1433  if (associated(igdstmpl)) deallocate(igdstmpl)
1434  if (associated(list_opt)) deallocate(list_opt)
1435 
1436 end subroutine getpoly
1437 
1498 subroutine gettemplates(cgrib, lcgrib, ifldnum, igds, igdstmpl, &
1499  igdslen, ideflist, idefnum, ipdsnum, ipdstmpl, ipdslen, coordlist, numcoord, ierr)
1500  implicit none
1501 
1502  character(len = 1), intent(in) :: cgrib(lcgrib)
1503  integer, intent(in) :: lcgrib, ifldnum
1504  integer, intent(out) :: igds(*), igdstmpl(*), ideflist(*)
1505  integer, intent(out) :: ipdsnum, ipdstmpl(*)
1506  integer, intent(out) :: idefnum, numcoord
1507  integer, intent(out) :: ierr
1508  real, intent(out) :: coordlist(*)
1509 
1510  character(len = 4), parameter :: grib = 'GRIB', c7777 = '7777'
1511  character(len = 4) :: ctemp
1512  integer:: listsec0(2)
1513  integer iofst, istart
1514  logical have3, have4
1515  integer :: igdslen, ipdslen, ipos, isecnum, j, jerr, lengrib, lensec, lensec0, numfld
1516 
1517  have3 = .false.
1518  have4 = .false.
1519  ierr = 0
1520  numfld = 0
1521 
1522  ! Check for valid request number.
1523  if (ifldnum .le. 0) then
1524  print *, 'gettemplates: Request for field number must be positive.'
1525  ierr = 3
1526  return
1527  endif
1528 
1529  ! Check for beginning of GRIB message in the first 100 bytes.
1530  istart = 0
1531  do j = 1, 100
1532  ctemp = cgrib(j) // cgrib(j + 1) // cgrib(j + 2) // cgrib(j + 3)
1533  if (ctemp .eq. grib ) then
1534  istart = j
1535  exit
1536  endif
1537  enddo
1538  if (istart .eq. 0) then
1539  print *, 'gettemplates: Beginning characters GRIB not found.'
1540  ierr = 1
1541  return
1542  endif
1543 
1544  ! Unpack Section 0 - Indicator Section.
1545  iofst = 8 * (istart + 5)
1546  call g2_gbytec(cgrib, listsec0(1), iofst, 8) ! Discipline
1547  iofst = iofst + 8
1548  call g2_gbytec(cgrib, listsec0(2), iofst, 8) ! GRIB edition number
1549  iofst = iofst + 8
1550  iofst = iofst + 32
1551  call g2_gbytec(cgrib, lengrib, iofst, 32) ! Length of GRIB message
1552  iofst = iofst + 32
1553  lensec0 = 16
1554  ipos = istart + lensec0
1555 
1556  ! Currently handles only GRIB Edition 2.
1557  if (listsec0(2) .ne. 2) then
1558  print *, 'gettemplates: can only decode GRIB edition 2.'
1559  ierr = 2
1560  return
1561  endif
1562 
1563  ! Loop through the remaining sections keeping track of the length of
1564  ! each. Also keep the latest Grid Definition Section info. Unpack
1565  ! the requested field number.
1566  do
1567  ! Check to see if we are at end of GRIB message.
1568  ctemp = cgrib(ipos) // cgrib(ipos + 1) // cgrib(ipos + 2) // cgrib(ipos + 3)
1569  if (ctemp .eq. c7777 ) then
1570  ipos = ipos+4
1571  ! If end of GRIB message not where expected, issue error
1572  if (ipos .ne. (istart + lengrib)) then
1573  print *, 'gettemplates: "7777" found, but not where expected.'
1574  ierr = 4
1575  return
1576  endif
1577  exit
1578  endif
1579  ! Get length of Section and Section number.
1580  iofst = (ipos - 1) * 8
1581  call g2_gbytec(cgrib, lensec, iofst, 32) ! Get Length of Section
1582  iofst = iofst + 32
1583  call g2_gbytec(cgrib, isecnum, iofst, 8) ! Get Section number
1584  iofst = iofst + 8
1585  !print *, ' lensec = ', lensec, ' secnum = ', isecnum
1586 
1587  ! If found Section 3, unpack the GDS info using the appropriate
1588  ! template. Save in case this is the latest grid before the
1589  ! requested field.
1590  if (isecnum .eq. 3) then
1591  iofst = iofst - 40 ! reset offset to beginning of section
1592  call unpack3(cgrib, lcgrib, iofst, igds, igdstmpl, igdslen, ideflist, idefnum, jerr)
1593  if (jerr .eq. 0) then
1594  have3 = .true.
1595  else
1596  ierr = 10
1597  return
1598  endif
1599  endif
1600 
1601  ! If found Section 4, check to see if this field is the one
1602  ! requested.
1603  if (isecnum .eq. 4) then
1604  numfld = numfld + 1
1605  if (numfld .eq. ifldnum) then
1606  iofst = iofst - 40 ! reset offset to beginning of section
1607  call unpack4(cgrib, lcgrib, iofst, ipdsnum, ipdstmpl, ipdslen, coordlist, numcoord, jerr)
1608  if (jerr .eq. 0) then
1609  have4 = .true.
1610  else
1611  ierr = 11
1612  return
1613  endif
1614  endif
1615  endif
1616 
1617  ! Check to see if we read pass the end of the GRIB message and
1618  ! missed the terminator string '7777'.
1619  ipos = ipos + lensec ! Update beginning of section pointer
1620  if (ipos .gt. (istart + lengrib)) then
1621  print *, 'gettemplates: "7777" not found at end of GRIB message.'
1622  ierr = 7
1623  return
1624  endif
1625 
1626  if (have3.and.have4) return
1627  enddo
1628 
1629  ! If exited from above loop, the end of the GRIB message was reached
1630  ! before the requested field was found.
1631  print *, 'gettemplates: GRIB message contained ', numfld, ' different fields.'
1632  print *, 'gettemplates: The request was for the ', ifldnum, ' field.'
1633  ierr = 6
1634 
1635 end 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 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:94
subroutine rdieee(rieee, a, num)
Copy array of 32-bit IEEE floating point values to local floating point representation.
Definition: g2bytes.F90:492
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_gbytescr(in, rout, iskip, nbits, nskip, n)
Extract big-endian floating-point values (32 bits each) from a packed bit string.
Definition: g2bytes.F90:68
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:506
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:1500
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:1106
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:902
subroutine getlocal(cgrib, lcgrib, localnum, csec2, lcsec2, ierr)
Return the contents of Section 2 (Local Use) from a GRIB2 message.
Definition: g2get.F90:1263
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:772
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:1162
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:1389
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:1011
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:248
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:177
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:168
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) of a GRIB2 message, starting at octet 6 of that Section.
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.
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.