62 subroutine getgb2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
63 unpack, k, gfld, iret)
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
76 integer,
intent(out) ::iret
79 call getgb2i2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
80 unpack, idxver, k, gfld, iret)
191 subroutine getgb2i2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
192 unpack, idxver, k, gfld, iret)
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
207 integer,
intent(out) ::iret
208 character(len = 1),
pointer,
dimension(:) :: cbuf
209 integer :: nnum, nlen, lpos, jk, irgi, irgs
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
218 subroutine getgb2s2(cbuf, idxver, nlen, nnum, j, jdisc, jids, jpdtn, jpdt, jgdtn, &
219 jgdt, k, gfld, lpos, iret)
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
230 integer,
intent(out) :: lpos, iret
232 subroutine getgb2l2(lugb, idxver, cindex, gfld, iret)
234 integer,
intent(in) :: lugb, idxver
235 character(len = 1),
intent(in) :: cindex(*)
237 integer,
intent(out) :: iret
239 subroutine getgb2r2(lugb, idxver, cindex, gfld, iret)
242 integer,
intent(in) :: lugb, idxver
243 character(len=1),
intent(in) :: cindex(*)
245 integer,
intent(out) :: iret
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
259 call getidx2(lugb, lugi, idxver, cbuf, nlen, nnum, irgi)
260 if (irgi .gt. 1)
then
267 call getgb2s2(cbuf, idxver, nlen, nnum, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, jk, &
269 if (irgs .ne. 0)
then
276 call getgb2l2(lugb, idxver, cbuf(lpos), gfld, iret)
280 call getgb2r2(lugb, idxver, cbuf(lpos), gfld, iret)
313 integer,
intent(in) :: lugb
314 character(len = 1),
intent(in) :: cindex(*)
316 integer,
intent(out) :: iret
319 subroutine getgb2l2(lugb, idxver, cindex, gfld, iret)
321 integer,
intent(in) :: lugb, idxver
322 character(len = 1),
intent(in) :: cindex(*)
324 integer,
intent(out) :: iret
328 call getgb2l2(lugb, 1, cindex, gfld, iret)
357 subroutine getgb2l2(lugb, idxver, cindex, gfld, iret)
362 integer,
intent(in) :: lugb, idxver
363 character(len = 1),
intent(in) :: cindex(*)
365 integer,
intent(out) :: iret
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)
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
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)
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
400 write(
g2_log_msg,
'(a, i2, a, i1)')
'getgb2l2: lugb ', lugb,
' idxver ', idxver
411 if (idxver .eq. 1)
then
413 call g2_gbytec1(cindex, lskip, mypos, int4_bits)
414 mypos = mypos + int4_bits
417 call g2_gbytec1(cindex, skip2, mypos, int4_bits)
422 mypos = mypos + int8_bits
425 mypos = mypos + int8_bits
429 if (skip28 .ne. 0)
then
430 iskip8 = lskip8 + skip28
433 call bareadl(lugb, iskip8, 4_8, lread8, csize)
435 allocate(ctemp(ilen))
439 call bareadl(lugb, iskip8, ilen8, lread8, ctemp)
440 if (ilen8 .ne. lread8)
then
446 call gf_unpack2(ctemp, ilen, iofst, gfld%locallen, gfld%local, ierr)
447 if (ierr .ne. 0)
then
520 subroutine getgb2p(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
521 extract, k, gribm, leng, iret)
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
537 integer (kind = 8) :: leng8
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
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))
654 subroutine getgb2p2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
655 extract, idxver, k, gribm, leng8, iret)
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
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)
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
684 subroutine getg2ir(lugb, msk1, msk2, mnum, cbuf, nlen, nnum, &
686 character(len = 1),
pointer,
dimension(:) :: cbuf
687 integer,
intent(in) :: lugb, msk1, msk2, mnum
688 integer,
intent(out) :: nlen, nnum, nmess, iret
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
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
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
712 subroutine getgb2s2(cbuf, idxver, nlen, nnum, j, jdisc, jids, jpdtn, jpdt, jgdtn, &
713 jgdt, k, gfld, lpos, iret)
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
724 integer,
intent(out) :: lpos, iret
732 if (lugi .gt. 0)
then
733 call getg2i2(lugi, cbuf, idxver, nlen, nnum, irgi)
734 elseif (lugi .le. 0)
then
736 call getg2i2r(lugb, msk1, msk2, mskp, idxver, cbuf, nlen, nnum, nmess, irgi)
738 if (irgi .gt. 1)
then
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
753 call getgb2rp2(lugb, idxver, cbuf(lpos:), extract, gribm, leng8, iret)
758 if (
associated(cbuf))
deallocate(cbuf)
801 integer,
intent(in) :: lugb
802 character(len=1),
intent(in) :: cindex(*)
804 integer,
intent(out) :: iret
807 subroutine getgb2r2(lugb, idxver, cindex, gfld, iret)
810 integer,
intent(in) :: lugb, idxver
811 character(len=1),
intent(in) :: cindex(*)
813 integer,
intent(out) :: iret
817 call getgb2r2(lugb, 1, cindex, gfld, iret)
856 subroutine getgb2r2(lugb, idxver, cindex, gfld, iret)
861 integer,
intent(in) :: lugb, idxver
862 character(len=1),
intent(in) :: cindex(*)
864 integer,
intent(out) :: iret
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
874 integer :: IXBMS1, IXBMS2
875 parameter(ixbms1 = 24, ixbms2 = 44)
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)
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
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
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
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)
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)
920 write(
g2_log_msg, *)
'getgb2r2: lugb ', lugb,
' idxver ', idxver
925 nullify(gfld%bmap, gfld%fld)
933 if (idxver .eq. 1)
then
934 call g2_gbytec1(cindex, lskip, int4_bits, int4_bits)
937 call g2_gbytec81(cindex, lskip8, int4_bits, int8_bits)
938 lskip = int(lskip8, kind(4))
942 if (idxver .eq. 1)
then
943 call g2_gbytec1(cindex, skip6, ixbms1 * int1_bits, int4_bits)
946 call g2_gbytec81(cindex, skip68, ixbms2 * int1_bits, int8_bits)
947 skip6 = int(skip68, kind(4))
951 write(
g2_log_msg, *)
' getgb2r2: skip6', skip6
956 if (idxver .eq. 1)
then
957 call g2_gbytec1(cindex, skip7, ixds1 * int1_bits, int4_bits)
959 call g2_gbytec81(cindex, skip78, ixds2 * int1_bits, int8_bits)
960 skip7 = int(skip78, kind(4))
964 write(
g2_log_msg, *)
' getgb2r2: skip7', skip7
969 if (gfld%ibmap .eq. 0 .or. gfld%ibmap .eq. 254)
then
970 iskip = lskip + skip6
971 iskip8 = lskip8 + skip6
974 call bareadl(lugb, iskip8, 4_8, lread8, csize)
976 allocate(ctemp(ilen))
980 call bareadl(lugb, iskip8, ilen8, lread8, ctemp)
981 if (ilen8 .ne. lread8)
then
989 call gf_unpack6(ctemp, ilen, iofst, gfld%ngrdpts, idum, gfld%bmap, ierr)
990 if (ierr .ne. 0)
then
999 iskip = lskip + skip7
1000 iskip8 = lskip8 + skip7
1003 call bareadl(lugb, iskip8, 4_8, lread8, csize)
1005 if (ilen .lt. 6) ilen = 6
1006 allocate(ctemp(ilen))
1010 call bareadl(lugb, iskip8, ilen8, lread8, ctemp)
1011 if (ilen8 .ne. lread8)
then
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
1030 if (gfld%ibmap .ne. 255 .AND.
associated(gfld%bmap))
then
1031 allocate(newfld(gfld%ngrdpts))
1033 do j = 1, gfld%ngrdpts
1034 if (gfld%bmap(j))
then
1035 newfld(j) = gfld%fld(n)
1041 deallocate(gfld%fld);
1043 gfld%expanded = .true.
1045 gfld%expanded = .true.
1084 subroutine getgb2rp(lugb, cindex, extract, gribm, leng, iret)
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
1093 integer (kind = 8) :: leng8
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
1111 call getgb2rp2(lugb, idxver, cindex, extract, gribm, leng8, iret)
1112 leng = int(leng8, kind(4))
1147 subroutine getgb2rp2(lugb, idxver, cindex, extract, gribm, leng8, iret)
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
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
1167 integer :: IXBMS1, IXBMS2
1168 parameter(ixbms1 = 24, ixbms2 = 44)
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
1179 character*1,
intent(inout) :: out(*)
1180 integer (kind = 8),
intent(in) :: sin
1181 integer,
intent(in) :: iskip, nbits
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
1189 character*1,
intent(in) :: in(*)
1190 integer (kind = 8),
intent(inout) :: siout
1191 integer,
intent(in) :: iskip, nbits
1196 write(
g2_log_msg, *)
'getgb2rp2: lugb ', lugb,
' idxver ', idxver,
' extract ', extract
1207 if (idxver .eq. 1)
then
1208 call g2_gbytec1(cindex, iskip, mypos, int4_bits)
1209 mypos = mypos + int4_bits
1211 call g2_gbytec1(cindex, iskp2, mypos, int4_bits)
1212 mypos = mypos + int4_bits
1214 mypos = mypos + 32 * int1_bits
1217 call g2_gbytec81(cindex, iskip8, mypos, int8_bits)
1218 mypos = mypos + int8_bits
1219 iskip = int(iskip8, kind(4))
1220 call g2_gbytec81(cindex, iskp2_8, mypos, int8_bits)
1221 mypos = mypos + int8_bits
1223 mypos = mypos + 52 * int1_bits
1226 write(
g2_log_msg, *)
'iskip8', iskip8,
'iskip', iskip,
'mypos/8', mypos/8
1231 if (iskp2_8 .gt. 0)
then
1232 call bareadl(lugb, iskip8 + iskp2_8, 4_8, lread8, ctemp)
1234 allocate(csec2(len2))
1236 call bareadl(lugb, iskip8 + iskp2_8, len2_8, lread8, csec2)
1241 write(
g2_log_msg, *)
'iskip8 ', iskip8,
' iskp2_8 ', iskp2_8,
'len2', len2,
'mypos/8', mypos/8
1246 call g2_gbytec1(cindex, len1, mypos, int4_bits)
1247 mypos = mypos + len1 * int1_bits
1248 call g2_gbytec1(cindex, len3, mypos, int4_bits)
1249 mypos = mypos + len3 * int1_bits
1250 call g2_gbytec1(cindex, len4, mypos, int4_bits)
1251 mypos = mypos + len4 * int1_bits
1252 call g2_gbytec1(cindex, len5, mypos, int4_bits)
1253 mypos = mypos + len5 * int1_bits
1254 call g2_gbytec1(cindex, len6, mypos, int4_bits)
1255 mypos = mypos + len6 * int1_bits
1257 write(
g2_log_msg, *)
'len1', len1,
'len3', len3,
'len4', len4,
'len5', len5,
'len6', len6
1262 call g2_gbytec1(cindex, ibmap, mypos, int1_bits)
1263 if (ibmap .eq. 254)
then
1265 if (idxver .eq. 1)
then
1266 call g2_gbytec1(cindex, iskp6, ixbms1 * int1_bits, int4_bits)
1268 call g2_gbytec81(cindex, iskp68, ixbms2 * int1_bits, int8_bits)
1269 iskp6 = int(iskp68, kind(4))
1272 write(
g2_log_msg, *)
'getgb2rp2: iskp6', iskp6
1278 call bareadl(lugb, iskip8 + iskp6, 4_8, lread8, ctemp)
1287 if (idxver .eq. 1)
then
1288 call g2_gbytec1(cindex, iskp7, ixds1 * int1_bits, int4_bits)
1290 call g2_gbytec81(cindex, iskp78, ixds2 * int1_bits, int8_bits)
1291 iskp7 = int(iskp78, kind(4))
1294 write(
g2_log_msg, *)
'getgb2rp2: iskp7', iskp7,
'IXDS2', ixds2
1299 call bareadl(lugb, iskip8 + iskp7, 4_8, lread8, ctemp)
1303 allocate(csec7(len7))
1305 call bareadl(lugb, iskip8 + iskp7, len7_8, lread8, csec7)
1308 write(
g2_log_msg, *)
'getgb2rp2: len0 ', len0,
'len1', len1,
'len2', len2 ,
'len3', len3
1310 write(
g2_log_msg, *)
'getgb2rp2: len4', len4,
'len5', len5,
'len5', len5,
'len6', len6,
'len7', len7,
'len8', len8
1315 leng8 = len0 + len1 + len2 + len3 + len4 + len5 + len6 + len7 + len8
1318 write(
g2_log_msg, *)
'getgb2rp2: len7 ', len7,
'lread8', lread8,
'calculated leng8', leng8
1323 if (.not.
associated(gribm))
allocate(gribm(leng8))
1332 gribm(7) = cindex(42 + inc)
1333 gribm(8) = cindex(41 + inc)
1334 call g2_sbytec81(gribm, leng8, int8_bits, int8_bits)
1337 write(
g2_log_msg, *)
'getgb2rp2: gribm(7) (discipline)', ichar(gribm(7)), &
1338 'gibm(8) (GRIB version)', ichar(gribm(8))
1343 gribm(17:16 + len1) = cindex(45 + inc:44 + inc + len1)
1345 ipos = 44 + inc + len1
1348 if (iskp2_8 .gt. 0)
then
1349 gribm(lencur + 1:lencur + len2) = csec2(1:len2)
1350 lencur = lencur + len2
1354 write(
g2_log_msg, *)
'getgb2rp2: copied 1, 2'
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
1364 write(
g2_log_msg, *)
'getgb2rp2: copied 3, 4, 5'
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
1373 if (idxver .eq. 1)
then
1374 call g2_gbytec1(cindex, iskp6, ixbms1 * int1_bits, int4_bits)
1377 call g2_gbytec81(cindex, iskp68, ixbms2 * int1_bits, int8_bits)
1380 write(
g2_log_msg, *)
'getgb2rp2: iskp68', iskp68
1383 call bareadl(lugb, iskip8 + iskp68, 4_8, lread8, ctemp)
1389 allocate(csec6(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)
1398 write(
g2_log_msg, *)
'getgb2rp2: copied 6, len7', len7
1403 gribm(lencur + 1:lencur + len7) = csec7(1:len7)
1404 lencur = lencur + len7
1412 gribm(lencur + 1) =
'7'
1413 gribm(lencur + 2) =
'7'
1414 gribm(lencur + 3) =
'7'
1415 gribm(lencur + 4) =
'7'
1418 if (
allocated(csec2))
deallocate(csec2)
1419 if (
allocated(csec7))
deallocate(csec7)
1421 if (idxver .eq. 1)
then
1422 call g2_gbytec1(cindex, iskip, mypos, int4_bits)
1423 mypos = mypos + int4_bits
1424 mypos = mypos + 6 * int4_bits
1427 call g2_gbytec81(cindex, iskip8, mypos, int8_bits)
1428 mypos = mypos + int8_bits
1429 mypos = mypos + 6 * int8_bits
1433 call g2_gbytec8(cindex, leng8, mypos, int8_bits)
1435 write(
g2_log_msg, *)
' iskip8 ', iskip8,
' mypos/8 ', mypos/8,
'index leng8 ', leng8
1440 if (.not.
associated(gribm))
allocate(gribm(leng8))
1443 call bareadl(lugb, iskip8, leng8, lread8, gribm)
1444 if (leng8 .ne. lread8)
then
1451 write(
g2_log_msg, *)
' read message into gribm, lread8', lread8
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...
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...
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...
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...
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...
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...
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.
subroutine getgb2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, unpack, k, gfld, iret)
This is a legacy version of getgb2i2().
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.
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.
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.
subroutine getgb2l(lugb, cindex, gfld, iret)
Read and unpack a local use section from a GRIB2 index record (index format 1) and GRIB2 file.
subroutine getgb2r(lugb, cindex, gfld, iret)
Read and unpack sections 6 and 7 from a GRIB2 message using a version 1 index record.
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.
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.
subroutine gf_free(gfld)
Free memory that was used to store array values in derived type grib_mod::gribfield.
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.
subroutine getg2ir(lugb, msk1, msk2, mnum, cbuf, nlen, nnum, nmess, iret)
Generate a version 1 index record for each message in a GRIB2 file.
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).
subroutine getg2i(lugi, cbuf, nlen, nnum, iret)
Read a version 1 index file and return its contents.
subroutine getg2i2(lugi, cbuf, idxver, nlen, nnum, iret)
Read a version 1 or 2 index file and return its contents.
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.
subroutine gf_unpack2(cgrib, lcgrib, iofst, lencsec2, csec2, ierr)
Unpack Section 2 (Local Use Section) of a GRIB2 message.
subroutine gf_unpack7(cgrib, lcgrib, iofst, igdsnum, igdstmpl, idrsnum, idrstmpl, ndpts, fld, ierr)
Unpack Section 7 (Data Section) of a GRIB2 message.
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.
Logging for the g2 library.
character *120 g2_log_msg
Fill this with the logging message.
subroutine g2_log(level)
Print a debug message for the g2 library.
This Fortran module contains the declaration of derived type gribfield.