NCEPLIBS-g2  3.4.9
g2getgb2.F90
Go to the documentation of this file.
1 
5 
62 subroutine getgb2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
63  unpack, k, gfld, iret)
64  use grib_mod
65  implicit none
66 
67  integer, intent(in) :: lugb, lugi, j, jdisc
68  integer, dimension(:) :: jids(*)
69  integer, intent(in) :: jpdtn
70  integer, dimension(:) :: jpdt(*)
71  integer, intent(in) :: jgdtn
72  integer, dimension(:) :: jgdt(*)
73  logical, intent(in) :: unpack
74  integer, intent(out) :: k
75  type(gribfield), intent(out) :: gfld
76  integer, intent(out) ::iret
77  integer :: idxver = 1
78 
79  call getgb2i2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
80  unpack, idxver, k, gfld, iret)
81 end subroutine getgb2
82 
191 subroutine getgb2i2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
192  unpack, idxver, k, gfld, iret)
193  use grib_mod
194  implicit none
195 
196  integer, intent(in) :: lugb, lugi, j, jdisc
197  integer, dimension(:) :: jids(*)
198  integer, intent(in) :: jpdtn
199  integer, dimension(:) :: jpdt(*)
200  integer, intent(in) :: jgdtn
201  integer, dimension(:) :: jgdt(*)
202  logical, intent(in) :: unpack
203  integer, intent(in) :: idxver
204  integer, intent(out) :: k
205  type(gribfield), intent(out) :: gfld
206  integer, intent(out) ::iret
207  character(len = 1), pointer, dimension(:) :: cbuf
208  integer :: nnum, nlen, lpos, jk, irgi, irgs
209 
210  ! Declare interfaces (required for cbuf pointer).
211  interface
212  subroutine getidx2(lugb, lugi, idxver, cindex, nlen, nnum, iret)
213  integer, intent(in) :: lugb, lugi, idxver
214  character(len = 1), pointer, dimension(:) :: cindex
215  integer, intent(out) :: nlen, nnum, iret
216  end subroutine getidx2
217  subroutine getgb2s2(cbuf, idxver, nlen, nnum, j, jdisc, jids, jpdtn, jpdt, jgdtn, &
218  jgdt, k, gfld, lpos, iret)
219  import gribfield
220  character(len = 1), intent(in) :: cbuf(nlen)
221  integer, intent(in) :: idxver, nlen, nnum, j, jdisc
222  integer, dimension(:) :: jids(*)
223  integer, intent(in) :: jpdtn
224  integer, dimension(:) :: jpdt(*)
225  integer, intent(in) :: jgdtn
226  integer, dimension(:) :: jgdt(*)
227  integer, intent(out) :: k
228  type(gribfield), intent(out) :: gfld
229  integer, intent(out) :: lpos, iret
230  end subroutine getgb2s2
231  subroutine getgb2l2(lugb, idxver, cindex, gfld, iret)
232  use grib_mod
233  integer, intent(in) :: lugb, idxver
234  character(len = 1), intent(in) :: cindex(*)
235  type(gribfield) :: gfld
236  integer, intent(out) :: iret
237  end subroutine getgb2l2
238  subroutine getgb2r2(lugb, idxver, cindex, gfld, iret)
239  use grib_mod
240  implicit none
241  integer, intent(in) :: lugb, idxver
242  character(len=1), intent(in) :: cindex(*)
243  type(gribfield) :: gfld
244  integer, intent(out) :: iret
245  end subroutine getgb2r2
246  end interface
247 
248  ! Fill cbuf with the index records of this file, by recalling them
249  ! from memory, reading them from the index file, or generating them
250  ! from the data file.
251  irgi = 0
252  call getidx2(lugb, lugi, idxver, cbuf, nlen, nnum, irgi)
253  if (irgi .gt. 1) then
254  iret = 96
255  return
256  endif
257 
258  ! Search the index in cbuf for the first message which meets our
259  ! search criteria.
260  call getgb2s2(cbuf, idxver, nlen, nnum, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, jk, &
261  gfld, lpos, irgs)
262  if (irgs .ne. 0) then
263  iret = 99
264  call gf_free(gfld)
265  return
266  endif
267 
268  ! Read local use section for the selected GRIB2 message, if available.
269  call getgb2l2(lugb, idxver, cbuf(lpos), gfld, iret)
270 
271  ! Read and unpack grib record for the selected GRIB2 message.
272  if (unpack) then
273  call getgb2r2(lugb, idxver, cbuf(lpos), gfld, iret)
274  endif
275 
276  ! Return the number of the unpacked field to the caller.
277  k = jk
278 end subroutine getgb2i2
279 
302 subroutine getgb2l(lugb, cindex, gfld, iret)
303  use grib_mod
304  implicit none
305 
306  integer, intent(in) :: lugb
307  character(len = 1), intent(in) :: cindex(*)
308  type(gribfield) :: gfld
309  integer, intent(out) :: iret
310 
311  interface
312  subroutine getgb2l2(lugb, idxver, cindex, gfld, iret)
313  use grib_mod
314  integer, intent(in) :: lugb, idxver
315  character(len = 1), intent(in) :: cindex(*)
316  type(gribfield) :: gfld
317  integer, intent(out) :: iret
318  end subroutine getgb2l2
319  end interface
320 
321  call getgb2l2(lugb, 1, cindex, gfld, iret)
322 
323 end subroutine getgb2l
324 
350 subroutine getgb2l2(lugb, idxver, cindex, gfld, iret)
351  use grib_mod
352  implicit none
353 
354  integer, intent(in) :: lugb, idxver
355  character(len = 1), intent(in) :: cindex(*)
356  type(gribfield) :: gfld
357  integer, intent(out) :: iret
358 
359  integer :: lskip, skip2
360  integer (kind = 8) :: lskip8, iskip8, lread8, ilen8
361  character(len = 1):: csize(4)
362  character(len = 1), allocatable :: ctemp(:)
363  integer :: ilen, iofst, ierr
364  integer :: INT1_BITS, INT2_BITS, INT4_BITS, INT8_BITS
365  parameter(int1_bits = 8, int2_bits = 16, int4_bits = 32, int8_bits = 64)
366  integer :: mypos
367 
368  interface
369  subroutine gf_unpack2(cgrib, lcgrib, iofst, lencsec2, csec2, ierr)
370  character(len = 1), intent(in) :: cgrib(lcgrib)
371  integer, intent(in) :: lcgrib
372  integer, intent(inout) :: iofst
373  integer, intent(out) :: lencsec2
374  integer, intent(out) :: ierr
375  character(len = 1), pointer, dimension(:) :: csec2
376  end subroutine gf_unpack2
377  end interface
378 
379  ! Get info.
380  nullify(gfld%local)
381  iret = 0
382  mypos = int4_bits
383  if (idxver .eq. 1) then
384  call g2_gbytec(cindex, lskip, mypos, int4_bits)
385  mypos = mypos + int4_bits
386  lskip8 = lskip
387  else
388  call g2_gbytec8(cindex, lskip8, mypos, int8_bits)
389  mypos = mypos + int8_bits
390  endif
391  call g2_gbytec(cindex, skip2, mypos, int4_bits)
392 
393  ! Read and unpack local use section, if present.
394  if (skip2 .ne. 0) then
395  iskip8 = lskip8 + skip2
396 
397  ! Get length of section.
398  call bareadl(lugb, iskip8, 4_8, lread8, csize)
399  call g2_gbytec(csize, ilen, 0, 32)
400  allocate(ctemp(ilen))
401  ilen8 = ilen
402 
403  ! Read in section.
404  call bareadl(lugb, iskip8, ilen8, lread8, ctemp)
405  if (ilen8 .ne. lread8) then
406  iret = 97
407  deallocate(ctemp)
408  return
409  endif
410  iofst = 0
411  call gf_unpack2(ctemp, ilen, iofst, gfld%locallen, gfld%local, ierr)
412  if (ierr .ne. 0) then
413  iret = 98
414  deallocate(ctemp)
415  return
416  endif
417  deallocate(ctemp)
418  else
419  gfld%locallen = 0
420  endif
421 end subroutine getgb2l2
422 
510 subroutine getgb2p(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
511  extract, k, gribm, leng, iret)
512  use grib_mod
513  implicit none
514 
515  integer, intent(in) :: lugb, lugi, j, jdisc, jpdtn, jgdtn
516  integer, dimension(:) :: jids(*), jpdt(*), jgdt(*)
517  logical, intent(in) :: extract
518  integer, intent(out) :: k, iret, leng
519  character(len = 1), pointer, dimension(:) :: gribm
520 
521  type(gribfield) :: gfld
522  integer :: msk1, irgi, irgs, jk, lpos, msk2, mskp, nlen, nmess, nnum
523 
524  character(len = 1), pointer, dimension(:) :: cbuf
525  parameter(msk1 = 32000, msk2 = 4000)
526 
527  ! Declare interfaces (required for cbuf pointer).
528  interface
529  subroutine getg2i(lugi, cbuf, nlen, nnum, iret)
530  character(len = 1), pointer, dimension(:) :: cbuf
531  integer, intent(in) :: lugi
532  integer, intent(out) :: nlen, nnum, iret
533  end subroutine getg2i
534  subroutine getg2ir(lugb, msk1, msk2, mnum, cbuf, nlen, nnum, &
535  nmess, iret)
536  character(len = 1), pointer, dimension(:) :: cbuf
537  integer, intent(in) :: lugb, msk1, msk2, mnum
538  integer, intent(out) :: nlen, nnum, nmess, iret
539  end subroutine getg2ir
540  subroutine getgb2rp(lugb, cindex, extract, gribm, leng, iret)
541  integer, intent(in) :: lugb
542  character(len = 1), intent(in) :: cindex(*)
543  logical, intent(in) :: extract
544  integer, intent(out) :: leng, iret
545  character(len = 1), pointer, dimension(:) :: gribm
546  end subroutine getgb2rp
547  end interface
548 
549  ! Initialize the index information in cbuf.
550  irgi = 0
551  if (lugi .gt. 0) then
552  call getg2i(lugi, cbuf, nlen, nnum, irgi)
553  elseif (lugi .le. 0) then
554  mskp = 0
555  call getg2ir(lugb, msk1, msk2, mskp, cbuf, nlen, nnum, nmess, irgi)
556  endif
557  if (irgi .gt. 1) then
558  iret = 96
559  return
560  endif
561 
562  ! Find info from index and fill a grib_mod::gribfield variable.
563  call getgb2s(cbuf, nlen, nnum, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
564  jk, gfld, lpos, irgs)
565  if (irgs .ne. 0) then
566  iret = 99
567  call gf_free(gfld)
568  return
569  endif
570 
571  ! Extract grib message from file.
572  nullify(gribm)
573  call getgb2rp(lugb, cbuf(lpos:), extract, gribm, leng, iret)
574 
575  k = jk
576 
577  ! Free cbuf memory allocated in getg2i/getg2ir().
578  if (associated(cbuf)) deallocate(cbuf)
579 
580  call gf_free(gfld)
581 end subroutine getgb2p
582 
617 subroutine getgb2r(lugb, cindex, gfld, iret)
618  use grib_mod
619  implicit none
620 
621  integer, intent(in) :: lugb
622  character(len=1), intent(in) :: cindex(*)
623  type(gribfield) :: gfld
624  integer, intent(out) :: iret
625 
626  interface
627  subroutine getgb2r2(lugb, idxver, cindex, gfld, iret)
628  use grib_mod
629  implicit none
630  integer, intent(in) :: lugb, idxver
631  character(len=1), intent(in) :: cindex(*)
632  type(gribfield) :: gfld
633  integer, intent(out) :: iret
634  end subroutine getgb2r2
635  end interface
636 
637  call getgb2r2(lugb, 1, cindex, gfld, iret)
638 
639 end subroutine getgb2r
640 
676 subroutine getgb2r2(lugb, idxver, cindex, gfld, iret)
677  use grib_mod
678  implicit none
679 
680  integer, intent(in) :: lugb, idxver
681  character(len=1), intent(in) :: cindex(*)
682  type(gribfield) :: gfld
683  integer, intent(out) :: iret
684 
685  integer :: lskip, skip6, skip7
686  character(len=1):: csize(4)
687  character(len=1), allocatable :: ctemp(:)
688  real, pointer, dimension(:) :: newfld
689  integer :: n, j, iskip, iofst, ilen, ierr, idum
690  integer :: inc
691  integer (kind = 8) :: lskip8, lread8, ilen8, iskip8
692  integer :: INT1_BITS, INT2_BITS, INT4_BITS, INT8_BITS
693  parameter(int1_bits = 8, int2_bits = 16, int4_bits = 32, int8_bits = 64)
694 
695  interface
696  subroutine gf_unpack6(cgrib, lcgrib, iofst, ngpts, ibmap, bmap, ierr)
697  character(len=1), intent(in) :: cgrib(lcgrib)
698  integer, intent(in) :: lcgrib, ngpts
699  integer, intent(inout) :: iofst
700  integer, intent(out) :: ibmap
701  integer, intent(out) :: ierr
702  logical*1, pointer, dimension(:) :: bmap
703  end subroutine gf_unpack6
704  subroutine gf_unpack7(cgrib, lcgrib, iofst, igdsnum, igdstmpl, &
705  idrsnum, idrstmpl, ndpts, fld, ierr)
706  character(len=1), intent(in) :: cgrib(lcgrib)
707  integer, intent(in) :: lcgrib, ndpts, idrsnum, igdsnum
708  integer, intent(inout) :: iofst
709  integer, pointer, dimension(:) :: idrstmpl, igdstmpl
710  integer, intent(out) :: ierr
711  real, pointer, dimension(:) :: fld
712  end subroutine gf_unpack7
713  end interface
714 
715  ! Get info.
716  nullify(gfld%bmap, gfld%fld)
717  iret = 0
718  inc = 0
719  if (idxver .eq. 1) then
720  call g2_gbytec(cindex, lskip, int4_bits, int4_bits)
721  lskip8 = lskip
722  else
723  inc = 4
724  call g2_gbytec8(cindex, lskip8, int4_bits, int8_bits)
725  lskip = int(lskip8, kind(4))
726  endif
727  call g2_gbytec(cindex, skip6, (24 + inc) * int1_bits, int4_bits)
728  call g2_gbytec(cindex, skip7, (28 + inc) * int1_bits, int4_bits)
729 
730  ! Read and unpack bit_map, if present.
731  if (gfld%ibmap .eq. 0 .or. gfld%ibmap .eq. 254) then
732  iskip = lskip + skip6
733  iskip8 = lskip8 + skip6
734 
735  ! get length of section.
736  call bareadl(lugb, iskip8, 4_8, lread8, csize)
737  call g2_gbytec(csize, ilen, 0, 32)
738  allocate(ctemp(ilen))
739  ilen8 = ilen
740 
741  ! read in section.
742  call bareadl(lugb, iskip8, ilen8, lread8, ctemp)
743  if (ilen8 .ne. lread8) then
744  iret = 97
745  deallocate(ctemp)
746  return
747  endif
748  iofst = 0
749  call gf_unpack6(ctemp, ilen, iofst, gfld%ngrdpts, idum, gfld%bmap, ierr)
750  if (ierr .ne. 0) then
751  iret = 98
752  deallocate(ctemp)
753  return
754  endif
755  deallocate(ctemp)
756  endif
757 
758  ! Read and unpack data field.
759  iskip = lskip + skip7
760  iskip8 = lskip8 + skip7
761 
762  ! Get length of section.
763  call bareadl(lugb, iskip8, 4_8, lread8, csize)
764  call g2_gbytec(csize, ilen, 0, 32)
765  if (ilen .lt. 6) ilen = 6
766  allocate(ctemp(ilen))
767  ilen8 = ilen
768 
769  ! Read in section.
770  call bareadl(lugb, iskip8, ilen8, lread8, ctemp)
771  if (ilen8 .ne. lread8) then
772  iret = 97
773  deallocate(ctemp)
774  return
775  endif
776  iofst = 0
777  call gf_unpack7(ctemp, ilen, iofst, gfld%igdtnum, gfld%igdtmpl, &
778  gfld%idrtnum, gfld%idrtmpl, gfld%ndpts, gfld%fld, ierr)
779  if (ierr .ne. 0) then
780  iret = 98
781  deallocate(ctemp)
782  return
783  endif
784  deallocate(ctemp)
785 
786  ! If bitmap is used with this field, expand data field
787  ! to grid, if possible.
788  if (gfld%ibmap .ne. 255 .AND. associated(gfld%bmap)) then
789  allocate(newfld(gfld%ngrdpts))
790  n = 1
791  do j = 1, gfld%ngrdpts
792  if (gfld%bmap(j)) then
793  newfld(j) = gfld%fld(n)
794  n = n+1
795  else
796  newfld(j) = 0.0
797  endif
798  enddo
799  deallocate(gfld%fld);
800  gfld%fld => newfld;
801  gfld%expanded = .true.
802  else
803  gfld%expanded = .true.
804  endif
805 end subroutine getgb2r2
806 
842 subroutine getgb2rp(lugb, cindex, extract, gribm, leng, iret)
843  implicit none
844 
845  integer, intent(in) :: lugb
846  character(len = 1), intent(in) :: cindex(*)
847  logical, intent(in) :: extract
848  character(len = 1), pointer, dimension(:) :: gribm
849  integer, intent(out) :: leng, iret
850 
851  interface
852  subroutine getgb2rp2(lugb, idxver, cindex, extract, gribm, leng, iret)
853  integer, intent(in) :: lugb, idxver
854  character(len = 1), intent(in) :: cindex(*)
855  logical, intent(in) :: extract
856  character(len = 1), pointer, dimension(:) :: gribm
857  integer, intent(out) :: leng, iret
858  end subroutine getgb2rp2
859  end interface
860 
861  call getgb2rp2(lugb, 1, cindex, extract, gribm, leng, iret)
862 
863 end subroutine getgb2rp
864 
896 subroutine getgb2rp2(lugb, idxver, cindex, extract, gribm, leng, iret)
897  implicit none
898 
899  integer, intent(in) :: lugb, idxver
900  character(len = 1), intent(in) :: cindex(*)
901  logical, intent(in) :: extract
902  integer, intent(out) :: leng, iret
903  character(len = 1), pointer, dimension(:) :: gribm
904 
905  integer, parameter :: zero = 0
906  character(len = 1), allocatable, dimension(:) :: csec2, csec6, csec7
907  character(len = 4) :: ctemp
908  integer :: lencur, len0, ibmap = 0, ipos, iskip
909  integer :: len7 = 0, len8 = 0, len3 = 0, len4 = 0, len5 = 0, len6 = 0, len1 = 0, len2 = 0
910  integer :: iskp2, iskp6, iskp7
911  integer :: INT1_BITS, INT2_BITS, INT4_BITS, INT8_BITS
912  parameter(int1_bits = 8, int2_bits = 16, int4_bits = 32, int8_bits = 64)
913  integer :: mypos, inc = 0
914  integer (kind = 8) :: lread8, iskip8, leng8, len2_8, len7_8, len6_8
915 
916  iret = 0
917 
918  ! Extract grib message from file.
919  mypos = int4_bits
920  if (extract) then
921  len0 = 16
922  len8 = 4
923  if (idxver .eq. 1) then
924  call g2_gbytec(cindex, iskip, mypos, int4_bits) ! bytes to skip in file
925  mypos = mypos + int4_bits
926  iskip8 = iskip
927  else
928  inc = 4
929  call g2_gbytec8(cindex, iskip8, mypos, int8_bits) ! bytes to skip in file
930  mypos = mypos + int8_bits
931  iskip = int(iskip8, kind(4))
932  endif
933  call g2_gbytec(cindex, iskp2, mypos, int4_bits) ! bytes to skip for section 2
934  mypos = mypos + int4_bits
935  if (iskp2 .gt. 0) then
936  call bareadl(lugb, iskip8 + iskp2, 4_8, lread8, ctemp)
937  call g2_gbytec(ctemp, len2, 0, int4_bits) ! length of section 2
938  allocate(csec2(len2))
939  len2_8 = len2
940  call bareadl(lugb, iskip8 + iskp2, len2_8, lread8, csec2)
941  else
942  len2 = 0
943  endif
944  mypos = mypos + 32 * int1_bits ! skip ahead in the cindex
945  call g2_gbytec(cindex, len1, mypos, int4_bits) ! length of section 1
946  ipos = 44 + len1
947  mypos = mypos + len1 * int1_bits ! skip ahead in the cindex
948  call g2_gbytec(cindex, len3, mypos, int4_bits) ! length of section 3
949  ipos = ipos + len3
950  mypos = mypos + len3 * int1_bits ! skip ahead in the cindex
951  call g2_gbytec(cindex, len4, mypos, int4_bits) ! length of section 4
952  ipos = ipos + len4
953  mypos = mypos + len4 * int1_bits ! skip ahead in the cindex
954  call g2_gbytec(cindex, len5, mypos, int4_bits) ! length of section 5
955  ipos = ipos + len5
956  mypos = mypos + len5 * int1_bits ! skip ahead in the cindex
957  call g2_gbytec(cindex, len6, mypos, int4_bits) ! length of section 6
958  ipos = ipos + 5
959  mypos = mypos + len6 * int1_bits ! skip ahead in the cindex
960  call g2_gbytec(cindex, ibmap, mypos, int1_bits) ! bitmap indicator
961  if (ibmap .eq. 254) then
962  call g2_gbytec(cindex, iskp6, (24 + inc) * int1_bits, int4_bits) ! bytes to skip for section 6
963  !call baread(lugb, iskip + iskp6, 4, lread, ctemp)
964  call bareadl(lugb, iskip8 + iskp6, 4_8, lread8, ctemp)
965  call g2_gbytec(ctemp, len6, 0, int4_bits) ! length of section 6
966  endif
967 
968  ! read in section 7 from file
969  call g2_gbytec(cindex, iskp7, (28 + inc) * int1_bits, int4_bits) ! bytes to skip for section 7
970  !call baread(lugb, iskip + iskp7, 4, lread, ctemp)
971  call bareadl(lugb, iskip8 + iskp7, 4_8, lread8, ctemp)
972  call g2_gbytec(ctemp, len7, 0, int4_bits) ! length of section 7
973  allocate(csec7(len7))
974  !call baread(lugb, iskip + iskp7, len7, lread, csec7)
975  len7_8 = len7
976  call bareadl(lugb, iskip8 + iskp7, len7_8, lread8, csec7)
977 
978  leng = len0 + len1 + len2 + len3 + len4 + len5 + len6 + len7 + len8
979  if (.not. associated(gribm)) allocate(gribm(leng))
980 
981  ! Create Section 0
982  gribm(1) = 'G'
983  gribm(2) = 'R'
984  gribm(3) = 'I'
985  gribm(4) = 'B'
986  gribm(5) = char(0)
987  gribm(6) = char(0)
988  gribm(7) = cindex(42 + inc)
989  gribm(8) = cindex(41 + inc)
990  gribm(9) = char(0)
991  gribm(10) = char(0)
992  gribm(11) = char(0)
993  gribm(12) = char(0)
994  call g2_sbytec(gribm, leng, 12*8, int4_bits)
995 
996  ! Copy Section 1
997  gribm(17:16 + len1) = cindex(45 + inc:44 + inc + len1)
998  lencur = 16 + inc + len1
999  ipos = 44 + inc + len1
1000 
1001  ! Copy Section 2, if necessary
1002  if (iskp2 .gt. 0) then
1003  gribm(lencur + 1:lencur + len2) = csec2(1:len2)
1004  lencur = lencur + len2
1005  endif
1006 
1007  ! Copy Sections 3 through 5
1008  gribm(lencur + 1:lencur + len3 + len4 + len5) = cindex(ipos + 1:ipos + len3 + len4 + len5)
1009  lencur = lencur + len3 + len4 + len5
1010  ipos = ipos + len3 + len4 + len5
1011 
1012  ! Copy Section 6
1013  if (len6 .eq. 6 .and. ibmap .ne. 254) then
1014  gribm(lencur + 1:lencur + len6) = cindex(ipos + 1:ipos + len6)
1015  lencur = lencur + len6
1016  else
1017  call g2_gbytec(cindex, iskp6, (24 + inc) * 8, int4_bits) ! bytes to skip for section 6
1018  call bareadl(lugb, iskip8 + iskp6, 4_8, lread8, ctemp)
1019  call g2_gbytec(ctemp, len6, 0, int4_bits) ! length of section 6
1020  allocate(csec6(len6))
1021  len6_8 = len6
1022  call bareadl(lugb, iskip8 + iskp6, len6_8, lread8, csec6)
1023  gribm(lencur + 1:lencur + len6) = csec6(1:len6)
1024  lencur = lencur + len6
1025  if (allocated(csec6)) deallocate(csec6)
1026  endif
1027 
1028  ! Copy Section 7
1029  gribm(lencur + 1:lencur + len7) = csec7(1:len7)
1030  lencur = lencur + len7
1031 
1032  ! Section 8
1033  gribm(lencur + 1) = '7'
1034  gribm(lencur + 2) = '7'
1035  gribm(lencur + 3) = '7'
1036  gribm(lencur + 4) = '7'
1037 
1038  ! clean up
1039  if (allocated(csec2)) deallocate(csec2)
1040  if (allocated(csec7)) deallocate(csec7)
1041  else ! do not extract field from message : get entire message
1042  if (idxver .eq. 1) then
1043  call g2_gbytec(cindex, iskip, mypos, int4_bits) ! bytes to skip in file
1044  mypos = mypos + int4_bits
1045  iskip8 = iskip
1046  else
1047  call g2_gbytec8(cindex, iskip8, mypos, int8_bits) ! bytes to skip in file
1048  mypos = mypos + int8_bits
1049  endif
1050  mypos = mypos + 7 * int4_bits
1051  call g2_gbytec(cindex, leng, mypos, int4_bits) ! length of grib message
1052  if (.not. associated(gribm)) allocate(gribm(leng))
1053  leng8 = leng
1054  call bareadl(lugb, iskip8, leng8, lread8, gribm)
1055  if (leng8 .ne. lread8) then
1056  deallocate(gribm)
1057  nullify(gribm)
1058  iret = 97
1059  return
1060  endif
1061  endif
1062 end subroutine getgb2rp2
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_gbytec8(in, iout, iskip, nbits)
Extract one arbitrary sized (up to 64-bits) values from a packed bit string, right justifying each va...
Definition: g2bytes.F90:152
subroutine g2_sbytec(out, in, iskip, nbits)
Put one arbitrary sized (up to 32 bits) value from an integer array, into a packed bit string,...
Definition: g2bytes.F90:238
subroutine getgb2i2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, unpack, idxver, k, gfld, iret)
Find and unpack a GRIB2 message in a file, using an version 1 or 2 index record which is either found...
Definition: g2getgb2.F90:193
subroutine getgb2l2(lugb, idxver, cindex, gfld, iret)
Read and unpack a local use section from a GRIB2 index record (index format 1 or 2) and GRIB2 file.
Definition: g2getgb2.F90:351
subroutine getgb2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, unpack, k, gfld, iret)
This is a legacy version of getgb2i2().
Definition: g2getgb2.F90:64
subroutine getgb2r2(lugb, idxver, cindex, gfld, iret)
Read and unpack sections 6 and 7 from a GRIB2 message using a version 1 or version 2 index record.
Definition: g2getgb2.F90:677
subroutine getgb2rp(lugb, cindex, extract, gribm, leng, iret)
Extract a grib message from a file given the index (index format 1) of the requested field.
Definition: g2getgb2.F90:843
subroutine getgb2l(lugb, cindex, gfld, iret)
Read and unpack a local use section from a GRIB2 index record (index format 1) and GRIB2 file.
Definition: g2getgb2.F90:303
subroutine getgb2rp2(lugb, idxver, cindex, extract, gribm, leng, iret)
Extract a grib message from a file given the version 1 or 2 index of the requested field.
Definition: g2getgb2.F90:897
subroutine getgb2r(lugb, cindex, gfld, iret)
Read and unpack sections 6 and 7 from a GRIB2 message using a version 1 index record.
Definition: g2getgb2.F90:618
subroutine getgb2p(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, extract, k, gribm, leng, iret)
Find and extract a GRIB2 message from a file.
Definition: g2getgb2.F90:512
subroutine gf_free(gfld)
Free memory that was used to store array values in derived type grib_mod::gribfield.
Definition: g2gf.F90:585
subroutine getgb2s(cbuf, nlen, nnum, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, k, gfld, lpos, iret)
Find information about a GRIB field from the index and fill a grib_mod::gribfield.
Definition: g2index.F90:786
subroutine getg2ir(lugb, msk1, msk2, mnum, cbuf, nlen, nnum, nmess, iret)
Generate a version 1 index record for each message in a GRIB2 file.
Definition: g2index.F90:566
subroutine getidx2(lugb, lugi, idxver, cindex, nlen, nnum, iret)
Find, read or generate a version 1 or 2 GRIB2 index for a GRIB2 file (which may be > 2 GB).
Definition: g2index.F90:266
subroutine getg2i(lugi, cbuf, nlen, nnum, iret)
Read a version 1 index file and return its contents.
Definition: g2index.F90:422
subroutine getgb2s2(cbuf, idxver, nlen, nnum, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, k, gfld, lpos, iret)
Find information about a GRIB field from the index and fill a grib_mod::gribfield.
Definition: g2index.F90:899
subroutine gf_unpack2(cgrib, lcgrib, iofst, lencsec2, csec2, ierr)
Unpack Section 2 (Local Use Section) of a GRIB2 message.
Definition: g2unpack.F90:100
subroutine gf_unpack7(cgrib, lcgrib, iofst, igdsnum, igdstmpl, idrsnum, idrstmpl, ndpts, fld, ierr)
Unpack Section 7 (Data Section) of a GRIB2 message.
Definition: g2unpack.F90:650
subroutine gf_unpack6(cgrib, lcgrib, iofst, ngpts, ibmap, bmap, ierr)
Unpack Section 6 (Bit-Map Section) of a GRIB2 message, starting at octet 6 of that Section.
Definition: g2unpack.F90:576
This Fortran module contains the declaration of derived type gribfield.
Definition: gribmod.F90:10