NCEPLIBS-g2  3.5.0
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 g2logging
194  use grib_mod
195  implicit none
196 
197  integer, intent(in) :: lugb, lugi, j, jdisc
198  integer, dimension(:) :: jids(*)
199  integer, intent(in) :: jpdtn
200  integer, dimension(:) :: jpdt(*)
201  integer, intent(in) :: jgdtn
202  integer, dimension(:) :: jgdt(*)
203  logical, intent(in) :: unpack
204  integer, intent(in) :: idxver
205  integer, intent(out) :: k
206  type(gribfield), intent(out) :: gfld
207  integer, intent(out) ::iret
208  character(len = 1), pointer, dimension(:) :: cbuf
209  integer :: nnum, nlen, lpos, jk, irgi, irgs
210 
211  ! Declare interfaces (required for cbuf pointer).
212  interface
213  subroutine getidx2(lugb, lugi, idxver, cindex, nlen, nnum, iret)
214  integer, intent(in) :: lugb, lugi, idxver
215  character(len = 1), pointer, dimension(:) :: cindex
216  integer, intent(out) :: nlen, nnum, iret
217  end subroutine getidx2
218  subroutine getgb2s2(cbuf, idxver, nlen, nnum, j, jdisc, jids, jpdtn, jpdt, jgdtn, &
219  jgdt, k, gfld, lpos, iret)
220  import gribfield
221  character(len = 1), intent(in) :: cbuf(nlen)
222  integer, intent(in) :: idxver, nlen, nnum, j, jdisc
223  integer, dimension(:) :: jids(*)
224  integer, intent(in) :: jpdtn
225  integer, dimension(:) :: jpdt(*)
226  integer, intent(in) :: jgdtn
227  integer, dimension(:) :: jgdt(*)
228  integer, intent(out) :: k
229  type(gribfield), intent(out) :: gfld
230  integer, intent(out) :: lpos, iret
231  end subroutine getgb2s2
232  subroutine getgb2l2(lugb, idxver, cindex, gfld, iret)
233  use grib_mod
234  integer, intent(in) :: lugb, idxver
235  character(len = 1), intent(in) :: cindex(*)
236  type(gribfield) :: gfld
237  integer, intent(out) :: iret
238  end subroutine getgb2l2
239  subroutine getgb2r2(lugb, idxver, cindex, gfld, iret)
240  use grib_mod
241  implicit none
242  integer, intent(in) :: lugb, idxver
243  character(len=1), intent(in) :: cindex(*)
244  type(gribfield) :: gfld
245  integer, intent(out) :: iret
246  end subroutine getgb2r2
247  end interface
248 
249 #ifdef LOGGING
250  write(g2_log_msg, '(a, i2, a, i2, a, i5, a, i5, a, l, a, i1)') 'getgb2i2: lugb ', lugb, ' lugi ', lugi, &
251  ' j ', j, ' jdisc ', jdisc, ' unpack ', unpack, ' idxver ', idxver
252  call g2_log(1)
253 #endif
254 
255  ! Fill cbuf with the index records of this file, by recalling them
256  ! from memory, reading them from the index file, or generating them
257  ! from the data file.
258  irgi = 0
259  call getidx2(lugb, lugi, idxver, cbuf, nlen, nnum, irgi)
260  if (irgi .gt. 1) then
261  iret = 96
262  return
263  endif
264 
265  ! Search the index in cbuf for the first message which meets our
266  ! search criteria.
267  call getgb2s2(cbuf, idxver, nlen, nnum, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, jk, &
268  gfld, lpos, irgs)
269  if (irgs .ne. 0) then
270  iret = 99
271  call gf_free(gfld)
272  return
273  endif
274 
275  ! Read local use section for the selected GRIB2 message, if available.
276  call getgb2l2(lugb, idxver, cbuf(lpos), gfld, iret)
277 
278  ! Read and unpack grib record for the selected GRIB2 message.
279  if (unpack) then
280  call getgb2r2(lugb, idxver, cbuf(lpos), gfld, iret)
281  endif
282 
283  ! Return the number of the unpacked field to the caller.
284  k = jk
285 end subroutine getgb2i2
286 
309 subroutine getgb2l(lugb, cindex, gfld, iret)
310  use grib_mod
311  implicit none
312 
313  integer, intent(in) :: lugb
314  character(len = 1), intent(in) :: cindex(*)
315  type(gribfield) :: gfld
316  integer, intent(out) :: iret
317 
318  interface
319  subroutine getgb2l2(lugb, idxver, cindex, gfld, iret)
320  use grib_mod
321  integer, intent(in) :: lugb, idxver
322  character(len = 1), intent(in) :: cindex(*)
323  type(gribfield) :: gfld
324  integer, intent(out) :: iret
325  end subroutine getgb2l2
326  end interface
327 
328  call getgb2l2(lugb, 1, cindex, gfld, iret)
329 
330 end subroutine getgb2l
331 
357 subroutine getgb2l2(lugb, idxver, cindex, gfld, iret)
358  use g2logging
359  use grib_mod
360  implicit none
361 
362  integer, intent(in) :: lugb, idxver
363  character(len = 1), intent(in) :: cindex(*)
364  type(gribfield) :: gfld
365  integer, intent(out) :: iret
366 
367  integer :: lskip, skip2
368  integer (kind = 8) :: lskip8, iskip8, lread8, ilen8, skip28
369  character(len = 1):: csize(4)
370  character(len = 1), allocatable :: ctemp(:)
371  integer :: ilen, iofst, ierr
372  integer :: INT1_BITS, INT2_BITS, INT4_BITS, INT8_BITS
373  parameter(int1_bits = 8, int2_bits = 16, int4_bits = 32, int8_bits = 64)
374  integer :: mypos
375 
376  interface
377  subroutine g2_gbytec1(in, siout, iskip, nbits)
378  character*1, intent(in) :: in(*)
379  integer, intent(inout) :: siout
380  integer, intent(in) :: iskip, nbits
381  end subroutine g2_gbytec1
382  subroutine g2_gbytec81(in, siout, iskip, nbits)
383  character*1, intent(in) :: in(*)
384  integer (kind = 8), intent(inout) :: siout
385  integer, intent(in) :: iskip, nbits
386  integer (kind = 8) :: iout(1)
387  end subroutine g2_gbytec81
388  subroutine gf_unpack2(cgrib, lcgrib, iofst, lencsec2, csec2, ierr)
389  character(len = 1), intent(in) :: cgrib(lcgrib)
390  integer, intent(in) :: lcgrib
391  integer, intent(inout) :: iofst
392  integer, intent(out) :: lencsec2
393  integer, intent(out) :: ierr
394  character(len = 1), pointer, dimension(:) :: csec2
395  end subroutine gf_unpack2
396  end interface
397 
398 
399 #ifdef LOGGING
400  write(g2_log_msg, '(a, i2, a, i1)') 'getgb2l2: lugb ', lugb, ' idxver ', idxver
401  call g2_log(1)
402 #endif
403 
404  ! Get info.
405  nullify(gfld%local)
406  iret = 0
407  mypos = int4_bits ! Skip length of index record.
408 
409  ! These are all 4-byte ints in index format 1, and 8-byte ints in
410  ! index version 2.
411  if (idxver .eq. 1) then
412  ! Read bytes to skip in file before message.
413  call g2_gbytec1(cindex, lskip, mypos, int4_bits)
414  mypos = mypos + int4_bits
415  lskip8 = lskip
416  ! Read bytes to skip in msg before local use.
417  call g2_gbytec1(cindex, skip2, mypos, int4_bits)
418  skip28 = skip2
419  else
420  ! Read bytes to skip in file before message.
421  call g2_gbytec81(cindex, lskip8, mypos, int8_bits)
422  mypos = mypos + int8_bits
423  ! Read bytes to skip in msg before local use.
424  call g2_gbytec81(cindex, skip28, mypos, int8_bits)
425  mypos = mypos + int8_bits
426  endif
427 
428  ! Read and unpack local use section, if present.
429  if (skip28 .ne. 0) then
430  iskip8 = lskip8 + skip28
431 
432  ! Get length of section.
433  call bareadl(lugb, iskip8, 4_8, lread8, csize)
434  call g2_gbytec1(csize, ilen, 0, 32)
435  allocate(ctemp(ilen))
436  ilen8 = ilen
437 
438  ! Read in section.
439  call bareadl(lugb, iskip8, ilen8, lread8, ctemp)
440  if (ilen8 .ne. lread8) then
441  iret = 97
442  deallocate(ctemp)
443  return
444  endif
445  iofst = 0
446  call gf_unpack2(ctemp, ilen, iofst, gfld%locallen, gfld%local, ierr)
447  if (ierr .ne. 0) then
448  iret = 98
449  deallocate(ctemp)
450  return
451  endif
452  deallocate(ctemp)
453  else
454  gfld%locallen = 0
455  endif
456 end subroutine getgb2l2
457 
520 subroutine getgb2p(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
521  extract, k, gribm, leng, iret)
522  use grib_mod
523  implicit none
524 
525  integer, intent(in) :: lugb, lugi, j, jdisc
526  integer, dimension(:) :: jids(*)
527  integer, intent(in) :: jpdtn
528  integer, dimension(:) :: jpdt(*)
529  integer, intent(in) :: jgdtn
530  integer, dimension(:) :: jgdt(*)
531  logical, intent(in) :: extract
532  integer, intent(out) :: k
533  character(len = 1), pointer, dimension(:) :: gribm
534  integer, intent(out) :: leng, iret
535  integer :: idxver
536 
537  integer (kind = 8) :: leng8
538 
539  interface
540  subroutine getgb2p2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
541  extract, idxver, k, gribm, leng8, iret)
542  integer, intent(in) :: lugb, lugi, j, jdisc
543  integer, dimension(:) :: jids(*)
544  integer, intent(in) :: jpdtn
545  integer, dimension(:) :: jpdt(*)
546  integer, intent(in) :: jgdtn
547  integer, dimension(:) :: jgdt(*)
548  logical, intent(in) :: extract
549  integer, intent(inout) :: idxver
550  integer, intent(out) :: k
551  character(len = 1), pointer, dimension(:) :: gribm
552  integer (kind = 8), intent(out) :: leng8
553  integer, intent(out) :: iret
554  end subroutine getgb2p2
555  end interface
556 
557  ! Call the new version of this subroutine, which handles messages > 2 GB.
558  idxver = 1
559  call getgb2p2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
560  extract, idxver, k, gribm, leng8, iret)
561  leng = int(leng8, kind(4))
562 
563 end subroutine getgb2p
564 
654 subroutine getgb2p2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
655  extract, idxver, k, gribm, leng8, iret)
656  use grib_mod
657  implicit none
658 
659  integer, intent(in) :: lugb, lugi, j, jdisc
660  integer, dimension(:) :: jids(*)
661  integer, intent(in) :: jpdtn
662  integer, dimension(:) :: jpdt(*)
663  integer, intent(in) :: jgdtn
664  integer, dimension(:) :: jgdt(*)
665  logical, intent(in) :: extract
666  integer, intent(inout) :: idxver
667  integer, intent(out) :: k
668  character(len = 1), pointer, dimension(:) :: gribm
669  integer (kind = 8), intent(out) :: leng8
670  integer, intent(out) :: iret
671 
672  type(gribfield) :: gfld
673  integer :: irgi, irgs, jk, lpos, mskp, nlen, nmess, nnum
674  character(len = 1), pointer, dimension(:) :: cbuf
675  integer (kind = 8) :: msk1, msk2
676  parameter(msk1 = 32000, msk2 = 4000)
677 
678  interface
679  subroutine getg2i(lugi, cbuf, nlen, nnum, iret)
680  character(len = 1), pointer, dimension(:) :: cbuf
681  integer, intent(in) :: lugi
682  integer, intent(out) :: nlen, nnum, iret
683  end subroutine getg2i
684  subroutine getg2ir(lugb, msk1, msk2, mnum, cbuf, nlen, nnum, &
685  nmess, iret)
686  character(len = 1), pointer, dimension(:) :: cbuf
687  integer, intent(in) :: lugb, msk1, msk2, mnum
688  integer, intent(out) :: nlen, nnum, nmess, iret
689  end subroutine getg2ir
690  subroutine getg2i2(lugi, cbuf, idxver, nlen, nnum, iret)
691  integer, intent(in) :: lugi
692  character(len=1), pointer, dimension(:) :: cbuf
693  integer, intent(out) :: idxver, nlen, nnum, iret
694  end subroutine getg2i2
695  subroutine getg2i2r(lugb, msk1, msk2, mnum, idxver, cbuf, &
696  nlen, nnum, nmess, iret)
697  integer, intent(in) :: lugb
698  integer (kind = 8), intent(in) :: msk1, msk2
699  integer, intent(in) :: mnum, idxver
700  character(len = 1), pointer, dimension(:) :: cbuf
701  integer, intent(out) :: nlen, nnum, nmess, iret
702  end subroutine getg2i2r
703  subroutine getgb2rp2(lugb, idxver, cindex, extract, gribm, leng8, iret)
704  integer, intent(in) :: lugb
705  integer, intent(inout) :: idxver
706  character(len = 1), intent(in) :: cindex(*)
707  logical, intent(in) :: extract
708  character(len = 1), pointer, dimension(:) :: gribm
709  integer (kind = 8), intent(out) :: leng8
710  integer, intent(out) :: iret
711  end subroutine getgb2rp2
712  subroutine getgb2s2(cbuf, idxver, nlen, nnum, j, jdisc, jids, jpdtn, jpdt, jgdtn, &
713  jgdt, k, gfld, lpos, iret)
714  import gribfield
715  character(len = 1), intent(in) :: cbuf(nlen)
716  integer, intent(in) :: idxver, nlen, nnum, j, jdisc
717  integer, dimension(:) :: jids(*)
718  integer, intent(in) :: jpdtn
719  integer, dimension(:) :: jpdt(*)
720  integer, intent(in) :: jgdtn
721  integer, dimension(:) :: jgdt(*)
722  integer, intent(out) :: k
723  type(gribfield), intent(out) :: gfld
724  integer, intent(out) :: lpos, iret
725  end subroutine getgb2s2
726  end interface
727 
728  nullify(gribm)
729 
730  ! Initialize the index information in cbuf.
731  irgi = 0
732  if (lugi .gt. 0) then
733  call getg2i2(lugi, cbuf, idxver, nlen, nnum, irgi)
734  elseif (lugi .le. 0) then
735  mskp = 0
736  call getg2i2r(lugb, msk1, msk2, mskp, idxver, cbuf, nlen, nnum, nmess, irgi)
737  endif
738  if (irgi .gt. 1) then
739  iret = 96
740  return
741  endif
742 
743  ! Find info from index and fill a grib_mod::gribfield variable.
744  call getgb2s2(cbuf, idxver, nlen, nnum, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
745  jk, gfld, lpos, irgs)
746  if (irgs .ne. 0) then
747  iret = 99
748  call gf_free(gfld)
749  return
750  endif
751 
752  ! Extract grib message from file.
753  call getgb2rp2(lugb, idxver, cbuf(lpos:), extract, gribm, leng8, iret)
754 
755  k = jk
756 
757  ! Free cbuf memory allocated in getg2i/getg2ir().
758  if (associated(cbuf)) deallocate(cbuf)
759 
760  call gf_free(gfld)
761 end subroutine getgb2p2
762 
797 subroutine getgb2r(lugb, cindex, gfld, iret)
798  use grib_mod
799  implicit none
800 
801  integer, intent(in) :: lugb
802  character(len=1), intent(in) :: cindex(*)
803  type(gribfield) :: gfld
804  integer, intent(out) :: iret
805 
806  interface
807  subroutine getgb2r2(lugb, idxver, cindex, gfld, iret)
808  use grib_mod
809  implicit none
810  integer, intent(in) :: lugb, idxver
811  character(len=1), intent(in) :: cindex(*)
812  type(gribfield) :: gfld
813  integer, intent(out) :: iret
814  end subroutine getgb2r2
815  end interface
816 
817  call getgb2r2(lugb, 1, cindex, gfld, iret)
818 
819 end subroutine getgb2r
820 
856 subroutine getgb2r2(lugb, idxver, cindex, gfld, iret)
857  use g2logging
858  use grib_mod
859  implicit none
860 
861  integer, intent(in) :: lugb, idxver
862  character(len=1), intent(in) :: cindex(*)
863  type(gribfield) :: gfld
864  integer, intent(out) :: iret
865 
866  integer :: lskip, skip6, skip7
867  integer (kind = 8) :: skip68, skip78
868  character(len=1):: csize(4)
869  character(len=1), allocatable :: ctemp(:)
870  real, pointer, dimension(:) :: newfld
871  integer :: n, j, iskip, iofst, ilen, ierr, idum
872  integer (kind = 8) :: lskip8, lread8, ilen8, iskip8
873  ! Bytes to skip in (version 1 and 2) index record to get to bms.
874  integer :: IXBMS1, IXBMS2
875  parameter(ixbms1 = 24, ixbms2 = 44)
876  ! Bytes to skip in (version 1 and 2) index record to get to data section.
877  integer :: IXDS1, IXDS2
878  parameter(ixds1 = 28, ixds2 = 52)
879  integer :: INT1_BITS, INT2_BITS, INT4_BITS, INT8_BITS
880  parameter(int1_bits = 8, int2_bits = 16, int4_bits = 32, int8_bits = 64)
881 
882  interface
883  subroutine gf_unpack6(cgrib, lcgrib, iofst, ngpts, ibmap, bmap, ierr)
884  character(len=1), intent(in) :: cgrib(lcgrib)
885  integer, intent(in) :: lcgrib, ngpts
886  integer, intent(inout) :: iofst
887  integer, intent(out) :: ibmap
888  integer, intent(out) :: ierr
889  logical*1, pointer, dimension(:) :: bmap
890  end subroutine gf_unpack6
891  subroutine gf_unpack7(cgrib, lcgrib, iofst, igdsnum, igdstmpl, &
892  idrsnum, idrstmpl, ndpts, fld, ierr)
893  character(len=1), intent(in) :: cgrib(lcgrib)
894  integer, intent(in) :: lcgrib, ndpts, idrsnum, igdsnum
895  integer, intent(inout) :: iofst
896  integer, pointer, dimension(:) :: idrstmpl, igdstmpl
897  integer, intent(out) :: ierr
898  real, pointer, dimension(:) :: fld
899  end subroutine gf_unpack7
900  subroutine g2_gbytec(in, iout, iskip, nbits)
901  character*1, intent(in) :: in(*)
902  integer, intent(inout) :: iout(*)
903  integer, intent(in) :: iskip, nbits
904  end subroutine g2_gbytec
905  subroutine g2_gbytec1(in, siout, iskip, nbits)
906  character*1, intent(in) :: in(*)
907  integer, intent(inout) :: siout
908  integer, intent(in) :: iskip, nbits
909  integer (kind = 4) :: iout(1)
910  end subroutine g2_gbytec1
911  subroutine g2_gbytec81(in, siout, iskip, nbits)
912  character*1, intent(in) :: in(*)
913  integer (kind = 8), intent(inout) :: siout
914  integer, intent(in) :: iskip, nbits
915  integer (kind = 8) :: iout(1)
916  end subroutine g2_gbytec81
917  end interface
918 
919 #ifdef LOGGING
920  write(g2_log_msg, *) 'getgb2r2: lugb ', lugb, ' idxver ', idxver
921  call g2_log(1)
922 #endif
923 
924  ! Initialize.
925  nullify(gfld%bmap, gfld%fld)
926  iret = 0
927 
928  ! Get the bytes to skip to reach the local use section. In index
929  ! version 1 this is a 4-byte value, in index version 2 it is an
930  ! 8-byte value. (To reach the local use offset in the index record,
931  ! we skip the first 4 bytes, which is the length of the index
932  ! record.)
933  if (idxver .eq. 1) then
934  call g2_gbytec1(cindex, lskip, int4_bits, int4_bits)
935  lskip8 = lskip
936  else
937  call g2_gbytec81(cindex, lskip8, int4_bits, int8_bits)
938  lskip = int(lskip8, kind(4))
939  endif
940 
941  ! Read the offset to section 6, the BMS section.
942  if (idxver .eq. 1) then
943  call g2_gbytec1(cindex, skip6, ixbms1 * int1_bits, int4_bits)
944  skip68 = skip6
945  else
946  call g2_gbytec81(cindex, skip68, ixbms2 * int1_bits, int8_bits)
947  skip6 = int(skip68, kind(4))
948  endif
949 
950 #ifdef LOGGING
951  write(g2_log_msg, *) ' getgb2r2: skip6', skip6
952  call g2_log(1)
953 #endif
954 
955  ! Read the offset to section 7, the data section.
956  if (idxver .eq. 1) then
957  call g2_gbytec1(cindex, skip7, ixds1 * int1_bits, int4_bits)
958  else
959  call g2_gbytec81(cindex, skip78, ixds2 * int1_bits, int8_bits)
960  skip7 = int(skip78, kind(4))
961  endif
962 
963 #ifdef LOGGING
964  write(g2_log_msg, *) ' getgb2r2: skip7', skip7
965  call g2_log(1)
966 #endif
967 
968  ! Read and unpack bit_map, if present.
969  if (gfld%ibmap .eq. 0 .or. gfld%ibmap .eq. 254) then
970  iskip = lskip + skip6
971  iskip8 = lskip8 + skip6
972 
973  ! Get length of bitmap section.
974  call bareadl(lugb, iskip8, 4_8, lread8, csize)
975  call g2_gbytec1(csize, ilen, 0, int4_bits)
976  allocate(ctemp(ilen))
977  ilen8 = ilen
978 
979  ! Read in bitmap section.
980  call bareadl(lugb, iskip8, ilen8, lread8, ctemp)
981  if (ilen8 .ne. lread8) then
982  iret = 97
983  deallocate(ctemp)
984  return
985  endif
986 
987  ! Unpack bitmap section.
988  iofst = 0
989  call gf_unpack6(ctemp, ilen, iofst, gfld%ngrdpts, idum, gfld%bmap, ierr)
990  if (ierr .ne. 0) then
991  iret = 98
992  deallocate(ctemp)
993  return
994  endif
995  deallocate(ctemp)
996  endif
997 
998  ! Read and unpack data field.
999  iskip = lskip + skip7
1000  iskip8 = lskip8 + skip7
1001 
1002  ! Get length of data section.
1003  call bareadl(lugb, iskip8, 4_8, lread8, csize)
1004  call g2_gbytec1(csize, ilen, 0, int4_bits)
1005  if (ilen .lt. 6) ilen = 6
1006  allocate(ctemp(ilen))
1007  ilen8 = ilen
1008 
1009  ! Read in data section.
1010  call bareadl(lugb, iskip8, ilen8, lread8, ctemp)
1011  if (ilen8 .ne. lread8) then
1012  iret = 97
1013  deallocate(ctemp)
1014  return
1015  endif
1016 
1017  ! Unpack data section.
1018  iofst = 0
1019  call gf_unpack7(ctemp, ilen, iofst, gfld%igdtnum, gfld%igdtmpl, &
1020  gfld%idrtnum, gfld%idrtmpl, gfld%ndpts, gfld%fld, ierr)
1021  if (ierr .ne. 0) then
1022  iret = 98
1023  deallocate(ctemp)
1024  return
1025  endif
1026  deallocate(ctemp)
1027 
1028  ! If bitmap is used with this field, expand data field
1029  ! to grid, if possible.
1030  if (gfld%ibmap .ne. 255 .AND. associated(gfld%bmap)) then
1031  allocate(newfld(gfld%ngrdpts))
1032  n = 1
1033  do j = 1, gfld%ngrdpts
1034  if (gfld%bmap(j)) then
1035  newfld(j) = gfld%fld(n)
1036  n = n+1
1037  else
1038  newfld(j) = 0.0
1039  endif
1040  enddo
1041  deallocate(gfld%fld);
1042  gfld%fld => newfld;
1043  gfld%expanded = .true.
1044  else
1045  gfld%expanded = .true.
1046  endif
1047 end subroutine getgb2r2
1048 
1084 subroutine getgb2rp(lugb, cindex, extract, gribm, leng, iret)
1085  implicit none
1086 
1087  integer, intent(in) :: lugb
1088  character(len = 1), intent(in) :: cindex(*)
1089  logical, intent(in) :: extract
1090  character(len = 1), pointer, dimension(:) :: gribm
1091  integer, intent(out) :: leng, iret
1092 
1093  integer (kind = 8) :: leng8
1094  integer :: idxver
1095 
1096  interface
1097  subroutine getgb2rp2(lugb, idxver, cindex, extract, gribm, leng8, iret)
1098  integer, intent(in) :: lugb
1099  integer, intent(inout) :: idxver
1100  character(len = 1), intent(in) :: cindex(*)
1101  logical, intent(in) :: extract
1102  character(len = 1), pointer, dimension(:) :: gribm
1103  integer (kind = 8), intent(out) :: leng8
1104  integer, intent(out) :: iret
1105  end subroutine getgb2rp2
1106  end interface
1107 
1108  ! Call the legacy version of this function. It will only work with
1109  ! GRIB messages < 2 GB.
1110  idxver = 1
1111  call getgb2rp2(lugb, idxver, cindex, extract, gribm, leng8, iret)
1112  leng = int(leng8, kind(4))
1113 
1114 end subroutine getgb2rp
1115 
1147 subroutine getgb2rp2(lugb, idxver, cindex, extract, gribm, leng8, iret)
1148  use g2logging
1149  implicit none
1150 
1151  integer, intent(in) :: lugb
1152  integer, intent(inout) :: idxver
1153  character(len = 1), intent(in) :: cindex(*)
1154  logical, intent(in) :: extract
1155  character(len = 1), pointer, dimension(:) :: gribm
1156  integer (kind = 8), intent(out) :: leng8
1157  integer, intent(out) :: iret
1158 
1159  integer, parameter :: zero = 0
1160  character(len = 1), allocatable, dimension(:) :: csec2, csec6, csec7
1161  character(len = 4) :: ctemp
1162  integer :: lencur, len0, ibmap = 0, ipos, iskip
1163  integer :: len7 = 0, len8 = 0, len3 = 0, len4 = 0, len5 = 0, len6 = 0, len1 = 0, len2 = 0
1164  integer :: iskp2, iskp6, iskp7
1165  integer (kind = 8) :: iskp2_8
1166  ! Bytes to skip in (version 1 and 2) index record to get to bms.
1167  integer :: IXBMS1, IXBMS2
1168  parameter(ixbms1 = 24, ixbms2 = 44)
1169  ! Bytes to skip in (version 1 and 2) index record to get to data section.
1170  integer :: IXDS1, IXDS2
1171  parameter(ixds1 = 28, ixds2 = 52)
1172  integer :: INT1_BITS, INT2_BITS, INT4_BITS, INT8_BITS
1173  parameter(int1_bits = 8, int2_bits = 16, int4_bits = 32, int8_bits = 64)
1174  integer :: mypos, inc = 0
1175  integer (kind = 8) :: lread8, iskip8, len2_8, len7_8, len6_8, iskp68, iskp78
1176 
1177  interface
1178  subroutine g2_sbytec81(out, sin, iskip, nbits)
1179  character*1, intent(inout) :: out(*)
1180  integer (kind = 8), intent(in) :: sin
1181  integer, intent(in) :: iskip, nbits
1182  end subroutine g2_sbytec81
1183  subroutine g2_gbytec1(in, siout, iskip, nbits)
1184  character*1, intent(in) :: in(*)
1185  integer, intent(inout) :: siout
1186  integer, intent(in) :: iskip, nbits
1187  end subroutine g2_gbytec1
1188  subroutine g2_gbytec81(in, siout, iskip, nbits)
1189  character*1, intent(in) :: in(*)
1190  integer (kind = 8), intent(inout) :: siout
1191  integer, intent(in) :: iskip, nbits
1192  end subroutine g2_gbytec81
1193  end interface
1194 
1195 #ifdef LOGGING
1196  write(g2_log_msg, *) 'getgb2rp2: lugb ', lugb, ' idxver ', idxver, ' extract ', extract
1197  call g2_log(1)
1198 #endif
1199 
1200  iret = 0
1201 
1202  ! Extract grib message from file.
1203  mypos = int4_bits
1204  if (extract) then
1205  len0 = 16
1206  len8 = 4
1207  if (idxver .eq. 1) then
1208  call g2_gbytec1(cindex, iskip, mypos, int4_bits) ! bytes to skip in file
1209  mypos = mypos + int4_bits
1210  iskip8 = iskip
1211  call g2_gbytec1(cindex, iskp2, mypos, int4_bits) ! bytes to skip for section 2
1212  mypos = mypos + int4_bits
1213  iskp2_8 = iskp2
1214  mypos = mypos + 32 * int1_bits ! skip ahead in the cindex
1215  else
1216  inc = 28
1217  call g2_gbytec81(cindex, iskip8, mypos, int8_bits) ! bytes to skip in file
1218  mypos = mypos + int8_bits
1219  iskip = int(iskip8, kind(4))
1220  call g2_gbytec81(cindex, iskp2_8, mypos, int8_bits) ! bytes to skip for section 2
1221  mypos = mypos + int8_bits
1222 
1223  mypos = mypos + 52 * int1_bits ! skip ahead in the cindex
1224  endif
1225 #ifdef LOGGING
1226  write(g2_log_msg, *) 'iskip8', iskip8, 'iskip', iskip, 'mypos/8', mypos/8
1227  call g2_log(2)
1228 #endif
1229 
1230  ! Determine length of local section (section 2).
1231  if (iskp2_8 .gt. 0) then
1232  call bareadl(lugb, iskip8 + iskp2_8, 4_8, lread8, ctemp)
1233  call g2_gbytec1(ctemp, len2, 0, int4_bits) ! length of section 2
1234  allocate(csec2(len2))
1235  len2_8 = len2
1236  call bareadl(lugb, iskip8 + iskp2_8, len2_8, lread8, csec2)
1237  else
1238  len2 = 0
1239  endif
1240 #ifdef LOGGING
1241  write(g2_log_msg, *) 'iskip8 ', iskip8, ' iskp2_8 ', iskp2_8, 'len2', len2, 'mypos/8', mypos/8
1242  call g2_log(2)
1243 #endif
1244 
1245  ! Find the lengths of the sections 1, 3, 4, 5, and 6.
1246  call g2_gbytec1(cindex, len1, mypos, int4_bits) ! length of section 1
1247  mypos = mypos + len1 * int1_bits ! skip ahead in the cindex
1248  call g2_gbytec1(cindex, len3, mypos, int4_bits) ! length of section 3
1249  mypos = mypos + len3 * int1_bits ! skip ahead in the cindex
1250  call g2_gbytec1(cindex, len4, mypos, int4_bits) ! length of section 4
1251  mypos = mypos + len4 * int1_bits ! skip ahead in the cindex
1252  call g2_gbytec1(cindex, len5, mypos, int4_bits) ! length of section 5
1253  mypos = mypos + len5 * int1_bits ! skip ahead in the cindex
1254  call g2_gbytec1(cindex, len6, mypos, int4_bits) ! length of section 6
1255  mypos = mypos + len6 * int1_bits ! skip ahead in the cindex
1256 #ifdef LOGGING
1257  write(g2_log_msg, *) 'len1', len1, 'len3', len3, 'len4', len4, 'len5', len5, 'len6', len6
1258  call g2_log(2)
1259 #endif
1260 
1261  ! Handle the bitmap, if present.
1262  call g2_gbytec1(cindex, ibmap, mypos, int1_bits) ! bitmap indicator
1263  if (ibmap .eq. 254) then
1264  ! Get the bytes to skip for section 6 from the index.
1265  if (idxver .eq. 1) then
1266  call g2_gbytec1(cindex, iskp6, ixbms1 * int1_bits, int4_bits)
1267  else
1268  call g2_gbytec81(cindex, iskp68, ixbms2 * int1_bits, int8_bits)
1269  iskp6 = int(iskp68, kind(4))
1270  endif
1271 #ifdef LOGGING
1272  write(g2_log_msg, *) 'getgb2rp2: iskp6', iskp6
1273  call g2_log(2)
1274 #endif
1275 
1276  ! Read the length of the bitmap section from the data file. (lu, byts to
1277  ! skip, bytes to read, bytes read, buffer for output)
1278  call bareadl(lugb, iskip8 + iskp6, 4_8, lread8, ctemp)
1279  call g2_gbytec1(ctemp, len6, 0, int4_bits) ! length of section 6
1280 #ifdef LOGGING
1281  write(g2_log_msg, *) 'getgb2rp2: len6', len6
1282  call g2_log(2)
1283 #endif
1284  endif
1285 
1286  ! Read the location of section 7 from the index.
1287  if (idxver .eq. 1) then
1288  call g2_gbytec1(cindex, iskp7, ixds1 * int1_bits, int4_bits) ! bytes to skip for section 7
1289  else
1290  call g2_gbytec81(cindex, iskp78, ixds2 * int1_bits, int8_bits) ! bytes to skip for section 7
1291  iskp7 = int(iskp78, kind(4))
1292  endif
1293 #ifdef LOGGING
1294  write(g2_log_msg, *) 'getgb2rp2: iskp7', iskp7, 'IXDS2', ixds2
1295  call g2_log(2)
1296 #endif
1297 
1298  ! Read in the length of section 7 from the data file.
1299  call bareadl(lugb, iskip8 + iskp7, 4_8, lread8, ctemp)
1300  call g2_gbytec1(ctemp, len7, 0, int4_bits) ! length of section 7
1301 
1302  ! Now read in section 7.
1303  allocate(csec7(len7))
1304  len7_8 = len7
1305  call bareadl(lugb, iskip8 + iskp7, len7_8, lread8, csec7)
1306 
1307 #ifdef LOGGING
1308  write(g2_log_msg, *) 'getgb2rp2: len0 ', len0, 'len1', len1, 'len2', len2 , 'len3', len3
1309  call g2_log(2)
1310  write(g2_log_msg, *) 'getgb2rp2: len4', len4, 'len5', len5, 'len5', len5, 'len6', len6, 'len7', len7, 'len8', len8
1311  call g2_log(2)
1312 #endif
1313 
1314  ! Now we know the total length of the grib message.
1315  leng8 = len0 + len1 + len2 + len3 + len4 + len5 + len6 + len7 + len8
1316 
1317 #ifdef LOGGING
1318  write(g2_log_msg, *) 'getgb2rp2: len7 ', len7, 'lread8', lread8, 'calculated leng8', leng8
1319  call g2_log(2)
1320 #endif
1321 
1322  ! Allocate storage for the message.
1323  if (.not. associated(gribm)) allocate(gribm(leng8))
1324 
1325  ! Create Section 0
1326  gribm(1) = 'G'
1327  gribm(2) = 'R'
1328  gribm(3) = 'I'
1329  gribm(4) = 'B'
1330  gribm(5) = char(0)
1331  gribm(6) = char(0)
1332  gribm(7) = cindex(42 + inc) ! discipline
1333  gribm(8) = cindex(41 + inc) ! GRIB version
1334  call g2_sbytec81(gribm, leng8, int8_bits, int8_bits)
1335 
1336 #ifdef LOGGING
1337  write(g2_log_msg, *) 'getgb2rp2: gribm(7) (discipline)', ichar(gribm(7)), &
1338  'gibm(8) (GRIB version)', ichar(gribm(8))
1339  call g2_log(2)
1340 #endif
1341 
1342  ! Copy Section 1 from the index to the message.
1343  gribm(17:16 + len1) = cindex(45 + inc:44 + inc + len1)
1344  lencur = 16 + len1
1345  ipos = 44 + inc + len1
1346 
1347  ! Copy Section 2, if necessary.
1348  if (iskp2_8 .gt. 0) then
1349  gribm(lencur + 1:lencur + len2) = csec2(1:len2)
1350  lencur = lencur + len2
1351  endif
1352 
1353 #ifdef LOGGING
1354  write(g2_log_msg, *) 'getgb2rp2: copied 1, 2'
1355  call g2_log(3)
1356 #endif
1357 
1358  ! Copy Sections 3 through 5 from the index to the message.
1359  gribm(lencur + 1:lencur + len3 + len4 + len5) = cindex(ipos + 1:ipos + len3 + len4 + len5)
1360  lencur = lencur + len3 + len4 + len5
1361  ipos = ipos + len3 + len4 + len5
1362 
1363 #ifdef LOGGING
1364  write(g2_log_msg, *) 'getgb2rp2: copied 3, 4, 5'
1365  call g2_log(3)
1366 #endif
1367 
1368  ! Copy Section 6 from the index to the message.
1369  if (len6 .eq. 6 .and. ibmap .ne. 254) then
1370  gribm(lencur + 1:lencur + len6) = cindex(ipos + 1:ipos + len6)
1371  lencur = lencur + len6
1372  else
1373  if (idxver .eq. 1) then
1374  call g2_gbytec1(cindex, iskp6, ixbms1 * int1_bits, int4_bits) ! bytes to skip for section 6
1375  iskp68 = iskp6
1376  else
1377  call g2_gbytec81(cindex, iskp68, ixbms2 * int1_bits, int8_bits) ! bytes to skip for section 6
1378  endif
1379 #ifdef LOGGING
1380  write(g2_log_msg, *) 'getgb2rp2: iskp68', iskp68
1381  call g2_log(3)
1382 #endif
1383  call bareadl(lugb, iskip8 + iskp68, 4_8, lread8, ctemp)
1384  call g2_gbytec1(ctemp, len6, 0, int4_bits) ! length of section 6
1385 #ifdef LOGGING
1386  write(g2_log_msg, *) 'getgb2rp2: len6', len6
1387  call g2_log(3)
1388 #endif
1389  allocate(csec6(len6))
1390  len6_8 = len6
1391  call bareadl(lugb, iskip8 + iskp68, len6_8, lread8, csec6)
1392  gribm(lencur + 1:lencur + len6) = csec6(1:len6)
1393  lencur = lencur + len6
1394  if (allocated(csec6)) deallocate(csec6)
1395  endif
1396 
1397 #ifdef LOGGING
1398  write(g2_log_msg, *) 'getgb2rp2: copied 6, len7', len7
1399  call g2_log(3)
1400 #endif
1401 
1402  ! Copy Section 7 to the message.
1403  gribm(lencur + 1:lencur + len7) = csec7(1:len7)
1404  lencur = lencur + len7
1405 
1406 #ifdef LOGGING
1407  write(g2_log_msg, *) 'getgb2rp2: copied 7'
1408  call g2_log(3)
1409 #endif
1410 
1411  ! Add Section 8 to the message.
1412  gribm(lencur + 1) = '7'
1413  gribm(lencur + 2) = '7'
1414  gribm(lencur + 3) = '7'
1415  gribm(lencur + 4) = '7'
1416 
1417  ! Clean up.
1418  if (allocated(csec2)) deallocate(csec2)
1419  if (allocated(csec7)) deallocate(csec7)
1420  else ! do not extract field from message : get entire message
1421  if (idxver .eq. 1) then
1422  call g2_gbytec1(cindex, iskip, mypos, int4_bits) ! bytes to skip in file
1423  mypos = mypos + int4_bits
1424  mypos = mypos + 6 * int4_bits
1425  iskip8 = iskip
1426  else
1427  call g2_gbytec81(cindex, iskip8, mypos, int8_bits) ! bytes to skip in file
1428  mypos = mypos + int8_bits
1429  mypos = mypos + 6 * int8_bits
1430  endif
1431 
1432  ! Get the length of the GRIB2 message from the index.
1433  call g2_gbytec8(cindex, leng8, mypos, int8_bits)
1434 #ifdef LOGGING
1435  write(g2_log_msg, *) ' iskip8 ', iskip8, ' mypos/8 ', mypos/8, 'index leng8 ', leng8
1436  call g2_log(2)
1437 #endif
1438 
1439  ! Allocate storage to hold the message.
1440  if (.not. associated(gribm)) allocate(gribm(leng8))
1441 
1442  ! Read the message from the file into buffer gribm.
1443  call bareadl(lugb, iskip8, leng8, lread8, gribm)
1444  if (leng8 .ne. lread8) then
1445  deallocate(gribm)
1446  nullify(gribm)
1447  iret = 97
1448  return
1449  endif
1450 #ifdef LOGGING
1451  write(g2_log_msg, *) ' read message into gribm, lread8', lread8
1452  call g2_log(3)
1453 #endif
1454 
1455  endif
1456 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:178
subroutine g2_sbytec81(out, sin, iskip, nbits)
Put one arbitrary sized (up to 64 bits) scalar into a packed bit string, taking the low order bits fr...
Definition: g2bytes.F90:505
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_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 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:358
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:857
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:1085
subroutine getgb2p2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, extract, idxver, k, gribm, leng8, iret)
Find and extract a GRIB2 message from a file.
Definition: g2getgb2.F90:656
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:310
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:798
subroutine getgb2p(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, extract, k, gribm, leng, iret)
Legacy subroutine to find and extract a GRIB2 message from a file.
Definition: g2getgb2.F90:522
subroutine getgb2rp2(lugb, idxver, cindex, extract, gribm, leng8, iret)
Extract a grib message from a file given the version 1 or 2 index of the requested field.
Definition: g2getgb2.F90:1148
subroutine gf_free(gfld)
Free memory that was used to store array values in derived type grib_mod::gribfield.
Definition: g2gf.F90:585
subroutine getg2i2r(lugb, msk1, msk2, mnum, idxver, cbuf, nlen, nnum, nmess, iret)
Generate a version 1 or 2 index record for each message in a GRIB2 file.
Definition: g2index.F90:649
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:588
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:430
subroutine getg2i2(lugi, cbuf, idxver, nlen, nnum, iret)
Read a version 1 or 2 index file and return its contents.
Definition: g2index.F90:512
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:929
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
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 the declaration of derived type gribfield.
Definition: gribmod.F90:10