NCEPLIBS-g2 4.0.0
Loading...
Searching...
No Matches
g2index.F90
Go to the documentation of this file.
1
4
20subroutine g2_create_index(lugb, lugi, idxver, filename, iret)
21 implicit none
22
23 integer, intent(in) :: lugb, lugi, idxver
24 character*(*) :: filename
25 integer, intent(out) :: iret
26
27 integer (kind = 8) :: msk1, msk2
28 parameter(msk1 = 32000_8, msk2 = 4000_8)
29 character(len=1), pointer, dimension(:) :: cbuf
30 integer :: numtot, nnum, nlen, mnum, kw
31 integer :: irgi, iw, nmess
32
33 interface
34 subroutine getg2i2r(lugb, msk1, msk2, mnum, idxver, cbuf, &
35 nlen, nnum, nmess, iret)
36 integer, intent(in) :: lugb
37 integer (kind = 8), intent(in) :: msk1, msk2
38 integer, intent(in) :: mnum, idxver
39 character(len = 1), pointer, dimension(:) :: cbuf
40 integer, intent(out) :: nlen, nnum, nmess, iret
41 end subroutine getg2i2r
42 end interface
43
44 ! Assume success.
45 iret = 0
46 numtot = 0
47 nlen = 0
48
49 ! Generate index records for all messages in file, or until memory
50 ! runs out.
51 mnum = 0
52 call getg2i2r(lugb, msk1, msk2, mnum, idxver, cbuf, &
53 nlen, nnum, nmess, irgi)
54 if (irgi .gt. 1 .or. nnum .eq. 0 .or. nlen .eq. 0) then
55 iret = 92
56 if (associated(cbuf)) then
57 deallocate(cbuf)
58 nullify(cbuf)
59 endif
60 return
61 endif
62 numtot = numtot + nnum
63 mnum = mnum + nmess
64
65 ! Write headers.
66 call g2_write_index_headers(lugi, nlen, numtot, idxver, filename)
67 iw = 162
68
69 ! Write the index data we have so far.
70 call bawrite(lugi, iw, nlen, kw, cbuf)
71 iw = iw + nlen
72
73 ! Extend index file if index buffer length too large to hold in memory.
74 if (irgi .eq. 1) then
75 do while (irgi .eq. 1 .and. nnum .gt. 0)
76 if (associated(cbuf)) then
77 deallocate(cbuf)
78 nullify(cbuf)
79 endif
80 call getg2i2r(11, msk1, msk2, mnum, idxver, cbuf, &
81 nlen, nnum, nmess, irgi)
82 if (irgi .le. 1 .and. nnum .gt. 0) then
83 numtot = numtot + nnum
84 mnum = mnum + nmess
85 call bawrite(lugi, iw, nlen, kw, cbuf)
86 iw = iw + nlen
87 endif
88 enddo
89 ! Go back and overwrite headers with new info.
90 call g2_write_index_headers(lugi, iw, numtot, idxver, filename)
91 endif
92 deallocate(cbuf)
93
94end subroutine g2_create_index
95
106subroutine g2_write_index_headers(lugi, nlen, nnum, idxver, filename)
107 implicit none
108
109 integer, intent(in) :: lugi, nlen, nnum, idxver
110 character, intent(in) :: filename*(*)
111
112 character cd8*8, ct10*10, hostname*15
113#ifdef __GFORTRAN__
114 integer istat
115#else
116 character hostnam*15
117 integer hostnm
118#endif
119 character chead(2)*81
120 integer :: kw
121
122 ! fill first 81-byte header
123 call date_and_time(cd8, ct10)
124 chead(1) = '!GFHDR!'
125 chead(1)(9:10) = ' 1'
126 chead(1)(12:14) = ' 1'
127 write(chead(1)(16:20),'(i5)') 162
128 chead(1)(22:31) = cd8(1:4) // '-' // cd8(5:6) // '-' // cd8(7:8)
129 chead(1)(33:40) = ct10(1:2) // ':' // ct10(3:4) // ':' // ct10(5:6)
130 chead(1)(42:47) = 'GB2IX1'
131 chead(1)(49:54) = ' '
132#ifdef __GFORTRAN__
133 istat = hostnm(hostname)
134 if (istat .eq. 0) then
135 chead(1)(56:70) = '0000'
136 else
137 chead(1)(56:70) = '0001'
138 endif
139#else
140 chead(1)(56:70) = hostnam(hostname)
141#endif
142 chead(1)(72:80) = 'grb2index'
143 chead(1)(81:81) = char(10)
144
145 ! fill second 81-byte header
146 if (idxver .eq. 1) then
147 chead(2) = 'IX1FORM:'
148 else
149 chead(2) = 'IX2FORM:'
150 endif
151 write(chead(2)(9:38),'(3i10)') 162, nlen, nnum
152 chead(2)(41:80) = filename
153 chead(2)(81:81) = char(10)
154
155 ! write headers at beginning of index file
156 call bawrite(lugi, 0, 162, kw, chead)
157
158 return
159end subroutine g2_write_index_headers
160
201subroutine getidx(lugb, lugi, cindex, nlen, nnum, iret)
202 implicit none
203
204 integer, intent(in) :: lugb, lugi
205 character(len = 1), pointer, dimension(:) :: cindex
206 integer, intent(out) :: nlen, nnum, iret
207 integer :: idxver
208
209 interface
210 subroutine getidx2(lugb, lugi, idxver, cindex, nlen, nnum, iret)
211 integer, intent(in) :: lugb, lugi
212 integer, intent(inout) :: idxver
213 character(len = 1), pointer, dimension(:) :: cindex
214 integer, intent(out) :: nlen, nnum, iret
215 end subroutine getidx2
216 end interface
217
218 ! When getidx() is called, always use index version 1. Call
219 ! getidx2() for a chance to set the index version.
220 idxver = 1
221 call getidx2(lugb, lugi, idxver, cindex, nlen, nnum, iret)
222
223 ! If index version 2 is being used, return error.
224 if (iret .eq. 0 .and. idxver .eq. 2) iret = 97
225
226end subroutine getidx
227
269subroutine getidx2(lugb, lugi, idxver, cindex, nlen, nnum, iret)
270 use g2logging
271 implicit none
272
273 integer, intent(in) :: lugb, lugi
274 integer, intent(inout) :: idxver
275 character(len = 1), pointer, dimension(:) :: cindex
276 integer, intent(out) :: nlen, nnum, iret
277
278 integer, parameter :: maxidx = 10000
279 integer (kind = 8), parameter :: msk1 = 32000_8, msk2 = 4000_8
280 integer :: lux
281 integer :: irgi, mskp, nmess, i
282
283 type gindex
284 integer :: nlen
285 integer :: nnum
286 integer :: idxver
287 character(len = 1), pointer, dimension(:) :: cbuf
288 end type gindex
289
290 type(gindex), save :: idxlist(10000)
291
292 data lux/0/
293
294 ! declare interfaces (required for cbuf pointer)
295 interface
296 subroutine getg2i2(lugi, cbuf, idxver, nlen, nnum, iret)
297 integer, intent(in) :: lugi
298 character(len=1), pointer, dimension(:) :: cbuf
299 integer, intent(out) :: idxver, nlen, nnum, iret
300 end subroutine getg2i2
301 subroutine getg2i2r(lugb, msk1, msk2, mnum, idxver, cbuf, &
302 nlen, nnum, nmess, iret)
303 integer, intent(in) :: lugb
304 integer (kind = 8), intent(in) :: msk1, msk2
305 integer, intent(in) :: mnum, idxver
306 character(len = 1), pointer, dimension(:) :: cbuf
307 integer, intent(out) :: nlen, nnum, nmess, iret
308 end subroutine getg2i2r
309 end interface
310
311#ifdef LOGGING
312 ! Log results for debugging.
313 write(g2_log_msg, '(a, i2, a, i2, a, i1)') 'getidx2: lugb ', lugb, ' lugi ', lugi, &
314 ' idxver ', idxver
315 call g2_log(1)
316#endif
317
318 ! Free all associated memory and exit.
319 if (lugb .eq. 0) then
320 !print *, 'getidx: Freeing all memory'
321 do i = 1, 10000
322 if (associated(idxlist(i)%cbuf)) then
323 !print *, 'deallocating ', loc(idxlist(i)%cbuf)
324 deallocate(idxlist(i)%cbuf)
325 nullify(idxlist(i)%cbuf)
326 endif
327 end do
328 iret = 0
329 return
330 endif
331
332 ! Determine whether index buffer needs to be initialized.
333 lux = 0
334 iret = 0
335 if (lugb .le. 0 .or. lugb .gt. 9999) then
336 print *, ' file unit number out of range'
337 print *, ' use unit numbers in range: 0 - 9999 '
338 iret = 90
339 return
340 endif
341
342 ! Force regeneration of index from GRIB2 file.
343 if (lugi .eq. lugb) then
344 if (associated(idxlist(lugb)%cbuf)) &
345 deallocate(idxlist(lugb)%cbuf)
346 !print *, 'Force regeneration'
347 nullify(idxlist(lugb)%cbuf)
348 idxlist(lugb)%nlen = 0
349 idxlist(lugb)%nnum = 0
350 lux = 0
351 endif
352
353 ! Force re-read of index from indexfile.
354 if (lugi .lt. 0) then
355 ! associated with unit abs(lugi)
356 if (associated(idxlist(lugb)%cbuf)) &
357 deallocate(idxlist(lugb)%cbuf)
358 !print *, 'Force re-read'
359 nullify(idxlist(lugb)%cbuf)
360 idxlist(lugb)%nlen = 0
361 idxlist(lugb)%nnum = 0
362 lux = abs(lugi)
363 endif
364
365 ! Check if index already exists in memory.
366 if (associated(idxlist(lugb)%cbuf)) then
367 !print *, 'Index exists in memory!'
368 cindex => idxlist(lugb)%cbuf
369 nlen = idxlist(lugb)%nlen
370 nnum = idxlist(lugb)%nnum
371 idxver = idxlist(lugb)%idxver
372 return
373 endif
374
375 ! Either read index from index file, or generate it from the GRIB2
376 ! file. This is an index for all messages in the file, each message
377 ! gets an index record, all stuffed into idxlist(lugbb)%cbuf.
378 irgi = 0
379 if (lux .gt. 0) then
380 call getg2i2(lux, idxlist(lugb)%cbuf, idxver, nlen, nnum, irgi)
381 elseif (lux .le. 0) then
382 mskp = 0
383 call getg2i2r(lugb, msk1, msk2, mskp, idxver, idxlist(lugb)%cbuf, &
384 nlen, nnum, nmess, irgi)
385 endif
386
387 ! Handle errors.
388 if (irgi .ne. 0) then
389 nlen = 0
390 nnum = 0
391 print *, ' error reading index file '
392 iret = 96
393 return
394 endif
395
396 ! Fill these values.
397 cindex => idxlist(lugb)%cbuf
398 idxlist(lugb)%nlen = nlen
399 idxlist(lugb)%nnum = nnum
400 idxlist(lugb)%idxver = idxver
401
402end subroutine getidx2
403
433subroutine getg2i(lugi, cbuf, nlen, nnum, iret)
434 implicit none
435
436 integer, intent(in) :: lugi
437 character(len=1), pointer, dimension(:) :: cbuf
438 integer, intent(out) :: nlen, nnum, iret
439 integer :: idxver
440
441 interface
442 subroutine getg2i2(lugi, cbuf, idxver, nlen, nnum, iret)
443 integer, intent(in) :: lugi
444 character(len=1), pointer, dimension(:) :: cbuf
445 integer, intent(out) :: idxver, nlen, nnum, iret
446 end subroutine getg2i2
447 end interface
448
449 ! Call the version of this subroutine that handles version 1 and
450 ! version 2. If getg2i() is called on a version 2 index, return
451 ! error.
452 call getg2i2(lugi, cbuf, idxver, nlen, nnum, iret)
453 if (idxver .eq. 2) iret = 5
454
455end subroutine getg2i
456
515subroutine getg2i2(lugi, cbuf, idxver, nlen, nnum, iret)
516 use g2logging
517 implicit none
518
519 integer, intent(in) :: lugi
520 character(len=1), pointer, dimension(:) :: cbuf
521 integer, intent(out) :: idxver, nlen, nnum, iret
522
523 character chead*162
524 integer :: ios, istat, lbuf, lhead, nskp
525
526#ifdef LOGGING
527 write(g2_log_msg, *) 'getg2i2: lugi ', lugi
528 call g2_log(1)
529#endif
530
531 nullify(cbuf)
532 nlen = 0
533 nnum = 0
534 iret = 4
535 call baread(lugi, 0, 162, lhead, chead)
536#ifdef LOGGING
537 write(g2_log_msg, *) 'lhead', lhead
538 call g2_log(2)
539#endif
540 if (lhead .eq. 162 .and. chead(42:47) .eq. 'GB2IX1') then
541 read(chead(82:162), '(2x, i1, 5x, 3i10, 2x, a40)', iostat = ios) idxver, nskp, nlen, nnum
542#ifdef LOGGING
543 write(g2_log_msg, *) 'ios', ios, 'idxver', idxver, 'nskp', nskp, 'nlen', nlen, 'nnum', nnum
544 call g2_log(2)
545#endif
546 if (ios .eq. 0) then
547 allocate(cbuf(nlen), stat = istat) ! Allocate space for cbuf.
548 if (istat .ne. 0) then
549 iret = 2
550 return
551 endif
552 iret = 0
553 call baread(lugi, nskp, nlen, lbuf, cbuf)
554 if (lbuf .ne. nlen) iret = 3
555 endif
556 endif
557end subroutine getg2i2
558
591subroutine getg2ir(lugb, msk1, msk2, mnum, cbuf, nlen, nnum, nmess, iret)
592 use re_alloc ! needed for subroutine realloc
593 implicit none
594
595 integer, intent(in) :: lugb
596 integer, intent(in) :: msk1, msk2
597 integer, intent(in) :: mnum
598 character(len = 1), pointer, dimension(:) :: cbuf
599 integer, intent(out) :: nlen, nnum, nmess, iret
600
601 integer (kind = 8) :: msk1_8, msk2_8
602
603 interface
604 subroutine getg2i2r(lugb, msk1, msk2, mnum, idxver, cbuf, &
605 nlen, nnum, nmess, iret)
606 integer, intent(in) :: lugb
607 integer (kind = 8), intent(in) :: msk1, msk2
608 integer, intent(in) :: mnum, idxver
609 character(len = 1), pointer, dimension(:) :: cbuf
610 integer, intent(out) :: nlen, nnum, nmess, iret
611 end subroutine getg2i2r
612 end interface
613
614 msk1_8 = msk1
615 msk2_8 = msk2
616 call getg2i2r(lugb, msk1_8, msk2_8, mnum, 1, cbuf, nlen, nnum, nmess, iret)
617end subroutine getg2ir
618
652subroutine getg2i2r(lugb, msk1, msk2, mnum, idxver, cbuf, nlen, nnum, nmess, iret)
653 use g2logging
654 use re_alloc ! needed for subroutine realloc
655 implicit none
656
657 integer, intent(in) :: lugb
658 integer (kind = 8), intent(in) :: msk1, msk2
659 integer, intent(in) :: mnum, idxver
660 character(len = 1), pointer, dimension(:) :: cbuf
661 integer, intent(out) :: nlen, nnum, nmess, iret
662
663 character(len = 1), pointer, dimension(:) :: cbuftmp
664 integer :: nbytes, newsize, next, numfld, m, mbuf
665 integer (kind = 8) :: iseek, lskip, lgrib
666 integer :: istat, init, iret1, lgrib4
667 parameter(init = 50000, next = 10000)
668
669 interface ! required for cbuf pointer
670 subroutine ix2gb2(lugb, lskip8, idxver, lgrib8, cbuf, numfld, mlen, iret)
671 integer, intent(in) :: lugb
672 integer (kind = 8), intent(in) :: lskip8
673 integer, intent(in) :: idxver
674 integer (kind = 8), intent(in) :: lgrib8
675 character(len = 1), pointer, intent(inout), dimension(:) :: cbuf
676 integer, intent(out) :: numfld, mlen, iret
677 end subroutine ix2gb2
678 end interface
679
680#ifdef LOGGING
681 ! Log results for debugging.
682 write(g2_log_msg, *) 'getg2i2r: lugb ', lugb, ' msk1 ', msk1, ' msk2 ', msk2, 'idxver', idxver
683 call g2_log(1)
684#endif
685
686 ! Initialize.
687 iret = 0
688 nullify(cbuf)
689 mbuf = init
690 allocate(cbuf(mbuf), stat = istat) ! allocate initial space for cbuf.
691 if (istat .ne. 0) then
692 iret = 2
693 return
694 endif
695
696 ! Search for first grib message.
697 iseek = 0_8
698 call skgb8(lugb, iseek, msk1, lskip, lgrib)
699 do m = 1, mnum
700 if (lgrib .gt. 0) then
701 iseek = lskip + lgrib
702 call skgb8(lugb, iseek, msk2, lskip, lgrib)
703 endif
704 enddo
705
706 ! Get index records for every grib message found.
707 nlen = 0
708 nnum = 0
709 nmess = mnum
710 do while (iret .eq. 0 .and. lgrib .gt. 0)
711 lgrib4 = int(lgrib, kind(4))
712 call ix2gb2(lugb, lskip, idxver, lgrib, cbuftmp, numfld, nbytes, iret1)
713 if (iret1 .ne. 0) print *, ' SAGT ', numfld, nbytes, iret1
714 if (nbytes + nlen .gt. mbuf) then ! Allocate more space, if necessary.
715 newsize = max(mbuf + next, mbuf + nbytes)
716 call realloc(cbuf, nlen, newsize, istat)
717 if (istat .ne. 0) then
718 iret = 1
719 return
720 endif
721 mbuf = newsize
722 endif
723
724 ! If index records were returned in cbuftmp from ixgb2,
725 ! copy cbuftmp into cbuf, then deallocate cbuftmp when done.
726 if (associated(cbuftmp)) then
727 cbuf(nlen + 1 : nlen + nbytes) = cbuftmp(1 : nbytes)
728 deallocate(cbuftmp, stat = istat)
729 if (istat .ne. 0) then
730 print *, ' deallocating cbuftmp ... ', istat
731 iret = 3
732 return
733 endif
734 nullify(cbuftmp)
735 nnum = nnum + numfld
736 nlen = nlen + nbytes
737 nmess = nmess + 1
738 endif
739
740 ! Look for next grib message.
741 iseek = lskip + lgrib
742 call skgb8(lugb, iseek, msk2, lskip, lgrib)
743 enddo
744end subroutine getg2i2r
745
818subroutine getgb2s(cbuf, nlen, nnum, j, jdisc, jids, jpdtn, jpdt, jgdtn, &
819 jgdt, k, gfld, lpos, iret)
820 use grib_mod
821 implicit none
822
823 character(len = 1), intent(in) :: cbuf(nlen)
824 integer, intent(in) :: nlen, nnum, j, jdisc
825 integer, dimension(:) :: jids(*)
826 integer, intent(in) :: jpdtn
827 integer, dimension(:) :: jpdt(*)
828 integer, intent(in) :: jgdtn
829 integer, dimension(:) :: jgdt(*)
830 integer, intent(out) :: k
831 type(gribfield), intent(out) :: gfld
832 integer, intent(out) :: lpos, iret
833
834 interface
835 subroutine getgb2s2(cbuf, idxver, nlen, nnum, j, jdisc, jids, jpdtn, jpdt, jgdtn, &
836 jgdt, k, gfld, lpos, iret)
837 import gribfield
838 character(len = 1), intent(in) :: cbuf(nlen)
839 integer, intent(in) :: idxver, nlen, nnum, j, jdisc
840 integer, dimension(:) :: jids(*)
841 integer, intent(in) :: jpdtn
842 integer, dimension(:) :: jpdt(*)
843 integer, intent(in) :: jgdtn
844 integer, dimension(:) :: jgdt(*)
845 integer, intent(out) :: k
846 type(gribfield), intent(out) :: gfld
847 integer, intent(out) :: lpos, iret
848 end subroutine getgb2s2
849 end interface
850
851 ! When getgb2s() is called, always use index version 1. Call
852 ! getgb2s2() to handle version 1 or 2.
853 call getgb2s2(cbuf, 1, nlen, nnum, j, jdisc, jids, jpdtn, jpdt, jgdtn, &
854 jgdt, k, gfld, lpos, iret)
855
856end subroutine getgb2s
857
931subroutine getgb2s2(cbuf, idxver, nlen, nnum, j, jdisc, jids, jpdtn, jpdt, jgdtn, &
932 jgdt, k, gfld, lpos, iret)
933 use grib_mod
934 use g2logging
935 implicit none
936
937 character(len = 1), intent(in) :: cbuf(nlen)
938 integer, intent(in) :: idxver, nlen, nnum, j, jdisc
939 integer, dimension(:) :: jids(*)
940 integer, intent(in) :: jpdtn
941 integer, dimension(:) :: jpdt(*)
942 integer, intent(in) :: jgdtn
943 integer, dimension(:) :: jgdt(*)
944 integer, intent(out) :: k
945 type(gribfield), intent(out) :: gfld
946 integer, intent(out) :: lpos, iret
947
948 integer :: kgds(5)
949 logical :: match1, match3, match4
950 integer :: i, icnd, inlen, iof, ipos, jpos, lsec1, lsec3, lsec4, lsec5, numgdt, numpdt, inc
951
952 interface
953 subroutine g2_gbytec1(in, siout, iskip, nbits)
954 character*1, intent(in) :: in(*)
955 integer, intent(inout) :: siout
956 integer, intent(in) :: iskip, nbits
957 end subroutine g2_gbytec1
958 subroutine gf_unpack1(cgrib, lcgrib, iofst, ids, idslen, ierr)
959 character(len = 1), intent(in) :: cgrib(lcgrib)
960 integer, intent(in) :: lcgrib
961 integer, intent(inout) :: iofst
962 integer, pointer, dimension(:) :: ids
963 integer, intent(out) :: ierr, idslen
964 end subroutine gf_unpack1
965 subroutine gf_unpack3(cgrib, lcgrib, iofst, igds, igdstmpl, &
966 mapgridlen, ideflist, idefnum, ierr)
967 character(len = 1), intent(in) :: cgrib(lcgrib)
968 integer, intent(in) :: lcgrib
969 integer, intent(inout) :: iofst
970 integer, pointer, dimension(:) :: igdstmpl, ideflist
971 integer, intent(out) :: igds(5)
972 integer, intent(out) :: ierr, idefnum
973 end subroutine gf_unpack3
974 subroutine gf_unpack4(cgrib, lcgrib, iofst, ipdsnum, ipdstmpl, &
975 mappdslen, coordlist, numcoord, ierr)
976 character(len = 1), intent(in) :: cgrib(lcgrib)
977 integer, intent(in) :: lcgrib
978 integer, intent(inout) :: iofst
979 real, pointer, dimension(:) :: coordlist
980 integer, pointer, dimension(:) :: ipdstmpl
981 integer, intent(out) :: ipdsnum
982 integer, intent(out) :: ierr, numcoord
983 end subroutine gf_unpack4
984 subroutine gf_unpack5(cgrib, lcgrib, iofst, ndpts, idrsnum, &
985 idrstmpl, mapdrslen, ierr)
986 character(len = 1), intent(in) :: cgrib(lcgrib)
987 integer, intent(in) :: lcgrib
988 integer, intent(inout) :: iofst
989 integer, intent(out) :: ndpts, idrsnum
990 integer, pointer, dimension(:) :: idrstmpl
991 integer, intent(out) :: ierr
992 end subroutine gf_unpack5
993 end interface
994
995#ifdef LOGGING
996 ! Log results for debugging.
997 write(g2_log_msg, *) 'getgb2s2: idxver ', idxver, ' nlen ', nlen, &
998 ' nnum ', nnum, ' j ', j, ' jdisc ', jdisc
999 call g2_log(1)
1000#endif
1001
1002 ! Initialize.
1003 k = 0
1004 lpos = 0
1005 iret = 1
1006 ipos = 0
1007 nullify(gfld%idsect, gfld%local)
1008 nullify(gfld%list_opt, gfld%igdtmpl, gfld%ipdtmpl)
1009 nullify(gfld%coord_list, gfld%idrtmpl, gfld%bmap, gfld%fld)
1010 if (idxver .eq. 1) then
1011 inc = 0
1012 else
1013 ! Add the extra 24 bytes in the version 2 index record, starting
1014 ! at byte 9.
1015 inc = 28
1016 endif
1017
1018 ! Search for request.
1019 do while (iret .ne. 0 .and. k .lt. nnum)
1020 k = k + 1
1021 ! Get length of current index record.
1022 call g2_gbytec1(cbuf, inlen, ipos * 8, 4 * 8)
1023 if (k .le. j) then ! skip this index
1024 ipos = ipos + inlen
1025 cycle
1026 endif
1027
1028 ! Check if grib2 discipline is a match.
1029 call g2_gbytec1(cbuf, gfld%discipline, (ipos + inc + 41) * 8, 1 * 8)
1030 if (jdisc .ne. -1 .and. jdisc .ne. gfld%discipline) then
1031 ipos = ipos + inlen
1032 cycle
1033 endif
1034
1035 ! Check if identification section is a match.
1036 match1 = .false.
1037 ! Get length of ids.
1038 call g2_gbytec1(cbuf, lsec1, (ipos + inc + 44) * 8, 4 * 8)
1039 iof = 0
1040 call gf_unpack1(cbuf(ipos + inc + 45), lsec1, iof, gfld%idsect, gfld%idsectlen, icnd)
1041 if (icnd .eq. 0) then
1042 match1 = .true.
1043 do i = 1, gfld%idsectlen
1044 if (jids(i) .ne. -9999 .and. jids(i) .ne. gfld%idsect(i)) then
1045 match1 = .false.
1046 exit
1047 endif
1048 enddo
1049 endif
1050 if (.not. match1) then
1051 deallocate(gfld%idsect)
1052 ipos = ipos + inlen
1053 cycle
1054 endif
1055
1056 ! Check if grid definition template is a match.
1057 jpos = ipos + 44 + inc + lsec1
1058 match3 = .false.
1059 call g2_gbytec1(cbuf, lsec3, jpos * 8, 4 * 8) ! get length of gds
1060 if (jgdtn .eq. -1) then
1061 match3 = .true.
1062 else
1063 call g2_gbytec1(cbuf, numgdt, (jpos + 12) * 8, 2 * 8) ! get gdt template no.
1064 if (jgdtn .eq. numgdt) then
1065 iof = 0
1066 call gf_unpack3(cbuf(jpos + 1), lsec3, iof, kgds, gfld%igdtmpl, &
1067 gfld%igdtlen, gfld%list_opt, gfld%num_opt, icnd)
1068 if (icnd .eq. 0) then
1069 match3 = .true.
1070 do i = 1, gfld%igdtlen
1071 if (jgdt(i) .ne. -9999 .and. jgdt(i).ne.gfld%igdtmpl(i)) then
1072 match3 = .false.
1073 exit
1074 endif
1075 enddo
1076 endif
1077 endif
1078 endif
1079 if (.not. match3) then
1080 if (associated(gfld%idsect)) deallocate(gfld%idsect)
1081 if (associated(gfld%igdtmpl)) deallocate(gfld%igdtmpl)
1082 if (associated(gfld%list_opt)) deallocate(gfld%list_opt)
1083 ipos = ipos + inlen
1084 cycle
1085 else
1086 gfld%griddef = kgds(1)
1087 gfld%ngrdpts = kgds(2)
1088 gfld%numoct_opt = kgds(3)
1089 gfld%interp_opt = kgds(4)
1090 gfld%igdtnum = kgds(5)
1091 endif
1092
1093 ! Check if product definition template is a match.
1094 jpos = jpos + lsec3
1095 match4 = .false.
1096
1097 ! Get length of pds.
1098 call g2_gbytec1(cbuf, lsec4, jpos * 8, 4 * 8)
1099 if (jpdtn .eq. -1) then
1100 match4 = .true.
1101 else
1102 ! Get pdt template no.
1103 call g2_gbytec1(cbuf, numpdt, (jpos + 7) * 8, 2 * 8)
1104 if (jpdtn .eq. numpdt) then
1105 iof = 0
1106 call gf_unpack4(cbuf(jpos + 1), lsec4, iof, gfld%ipdtnum, &
1107 gfld%ipdtmpl, gfld%ipdtlen, gfld%coord_list, gfld%num_coord, icnd)
1108 if (icnd .eq. 0) then
1109 match4 = .true.
1110 do i = 1, gfld%ipdtlen
1111 if (jpdt(i) .ne. -9999 .and. jpdt(i) .ne. gfld%ipdtmpl(i)) then
1112 match4 = .false.
1113 exit
1114 endif
1115 enddo
1116 endif
1117 endif
1118 endif
1119 if (.not. match4) then
1120 if (associated(gfld%idsect)) deallocate(gfld%idsect)
1121 if (associated(gfld%ipdtmpl)) deallocate(gfld%ipdtmpl)
1122 if (associated(gfld%coord_list)) deallocate(gfld%coord_list)
1123 endif
1124
1125 ! If request is found set values for derived type gfld and return.
1126 if (match1 .and. match3 .and. match4) then
1127 lpos = ipos + 1
1128 call g2_gbytec1(cbuf, gfld%version, (ipos + inc + 40) * 8, 1 * 8)
1129 call g2_gbytec1(cbuf, gfld%ifldnum, (ipos + inc + 42) * 8, 2 * 8)
1130 gfld%unpacked = .false.
1131 jpos = ipos + 44 + inc + lsec1
1132 if (jgdtn .eq. -1) then ! unpack gds, if not done before
1133 iof = 0
1134 call gf_unpack3(cbuf(jpos + 1), lsec3, iof, kgds, gfld%igdtmpl, &
1135 gfld%igdtlen, gfld%list_opt, gfld%num_opt, icnd)
1136 gfld%griddef = kgds(1)
1137 gfld%ngrdpts = kgds(2)
1138 gfld%numoct_opt = kgds(3)
1139 gfld%interp_opt = kgds(4)
1140 gfld%igdtnum = kgds(5)
1141 endif
1142 jpos = jpos + lsec3
1143 if (jpdtn .eq. -1) then ! unpack pds, if not done before
1144 iof = 0
1145 call gf_unpack4(cbuf(jpos + 1), lsec4, iof, gfld%ipdtnum, &
1146 gfld%ipdtmpl, gfld%ipdtlen, gfld%coord_list, gfld%num_coord, icnd)
1147 endif
1148 jpos = jpos + lsec4
1149 call g2_gbytec1(cbuf, lsec5, jpos * 8, 4 * 8) ! get length of drs
1150 iof = 0
1151 call gf_unpack5(cbuf(jpos + 1), lsec5, iof, gfld%ndpts, &
1152 gfld%idrtnum, gfld%idrtmpl, gfld%idrtlen, icnd)
1153 jpos = jpos + lsec5
1154 call g2_gbytec1(cbuf, gfld%ibmap, (jpos + 5) * 8, 1 * 8) ! get ibmap
1155 iret = 0
1156 else ! pdt did not match
1157 ipos = ipos + inlen
1158 endif
1159 enddo
1160end subroutine getgb2s2
1161
1195subroutine ixgb2(lugb, lskip, lgrib, cbuf, numfld, mlen, iret)
1196 use re_alloc ! needed for subroutine realloc
1197 implicit none
1198
1199 integer :: lugb, lskip, lgrib
1200 character(len = 1), pointer, dimension(:) :: cbuf
1201 integer :: numfld, mlen, iret
1202 integer (kind = 8) :: lskip8, lgrib8
1203
1204 interface
1205 subroutine ix2gb2(lugb, lskip8, idxver, lgrib8, cbuf, numfld, mlen, iret)
1206 integer, intent(in) :: lugb
1207 integer (kind = 8), intent(in) :: lskip8
1208 integer, intent(in) :: idxver
1209 integer (kind = 8), intent(in) :: lgrib8
1210 character(len = 1), pointer, intent(inout), dimension(:) :: cbuf
1211 integer, intent(out) :: numfld, mlen, iret
1212 end subroutine ix2gb2
1213 end interface
1214
1215 ! Always use index version 1 from this subroutine.
1216 lskip8 = lskip
1217 lgrib8 = lgrib
1218 call ix2gb2(lugb, lskip8, 1, lgrib8, cbuf, numfld, mlen, iret)
1219 lgrib = int(lgrib8, 4)
1220end subroutine ixgb2
1221
1259subroutine ix2gb2(lugb, lskip8, idxver, lgrib8, cbuf, numfld, mlen, iret)
1260 use re_alloc ! needed for subroutine realloc
1261 use g2logging
1262 implicit none
1263
1264 ! Subroutine parameters.
1265 integer, intent(in) :: lugb
1266 integer (kind = 8), intent(in) :: lskip8
1267 integer, intent(in) :: idxver
1268 integer (kind = 8), intent(in) :: lgrib8
1269 character(len = 1), pointer, intent(inout), dimension(:) :: cbuf
1270 integer, intent(out) :: numfld, mlen, iret
1271
1272 character cver, cdisc
1273 character(len = 4) :: ctemp
1274 integer (kind = 8) :: loclus8, locgds8, locbms8
1275 integer locbms
1276 integer :: indbmp, numsec, newsize, g2_mova2i, mbuf, lindex
1277 integer :: lskip
1278 integer :: ilndrs, ilnpds, istat
1279 integer (kind = 8) :: ibread8, lbread8, ibskip8, lengds8
1280 integer (kind = 8) :: ilnpds8, ilndrs8
1281 integer :: lensec, lensec1
1282 integer :: mypos, inc
1283
1284 ! Parameters.
1285 ! Size of the internal char buffers used in this subroutine.
1286 integer :: LINMAX
1287 parameter(linmax = 5000)
1288 ! Initial size of buffer to hold index record in cbuf array.
1289 integer :: INIT
1290 parameter(init = 50000)
1291 ! If more space is needed in cbuf, it will get this much more.
1292 integer :: NEXT
1293 parameter(next = 10000)
1294 ! Number of bytes of the BMS section put into index record.
1295 integer :: MXBMS
1296 parameter(mxbms = 6)
1297 integer :: IXDS1, IXDS2
1298 parameter(ixds1 = 28, ixds2 = 52)
1299 ! Bytes to skip in (version 1) index record to get to section 0.
1300 integer :: IXIDS
1301 parameter(ixids = 44)
1302 integer :: IXSDR
1303 parameter(ixsdr = 20)
1304 ! Bytes to skip in (version 1 and 2) index record to get to bms.
1305 integer :: IXBMS1, IXBMS2
1306 parameter(ixbms1 = 24, ixbms2 = 44)
1307 ! Sizes of integers in bits.
1308 integer :: INT1_BITS, INT2_BITS, INT4_BITS, INT8_BITS
1309 parameter(int1_bits = 8, int2_bits = 16, int4_bits = 32, int8_bits = 64)
1310 ! Location of bytes to drs field in index version 1 and 2.
1311 integer :: IXDRS1, IXDRS2
1312 parameter(ixdrs1 = 20, ixdrs2 = 36)
1313
1314 ! Buffers.
1315 character cbread(LINMAX), cindex(LINMAX)
1316 character cids(LINMAX), cgds(LINMAX)
1317
1318 interface
1319 subroutine g2_gbytec1(in, siout, iskip, nbits)
1320 character*1, intent(in) :: in(*)
1321 integer, intent(inout) :: siout
1322 integer, intent(in) :: iskip, nbits
1323 end subroutine g2_gbytec1
1324 subroutine g2_sbytec81(out, sin, iskip, nbits)
1325 character*1, intent(inout) :: out(*)
1326 integer (kind = 8), intent(in) :: sin
1327 integer, intent(in) :: iskip, nbits
1328 end subroutine g2_sbytec81
1329 subroutine g2_sbytec1(out, in, iskip, nbits)
1330 character*1, intent(inout) :: out(*)
1331 integer, intent(in) :: in
1332 integer, intent(in) :: iskip, nbits
1333 end subroutine g2_sbytec1
1334 end interface
1335
1336#ifdef LOGGING
1337 ! Log results for debugging.
1338 write(g2_log_msg, *) 'ix2gb2: lugb ', lugb, ' lskip8 ', lskip8, ' idxver ', idxver
1339 call g2_log(1)
1340#endif
1341
1342 ! Are we using index version 1 (legacy), or version 2 (introduced to
1343 ! handle files > 2 GB).
1344 if (idxver .eq. 1) then
1345 inc = 0
1346 else
1347 ! Add the extra bytes in the version 2 index record, starting at
1348 ! byte 9. This is because some values early in the index record
1349 ! changed from 4-byte ints to 8-byte ints. This is the total
1350 ! extra bytes that were added to the beginning of the index
1351 ! record in version 2.
1352 inc = 28
1353 endif
1354
1355 ! Initialize values and allocate buffer (at the user-provided cbuf
1356 ! pointer) where the index data will be written. When subroutine is
1357 ! complete, cbuf will hold either version 1 or 2 index data.
1358 loclus8 = 0
1359 iret = 0
1360 mlen = 0
1361 numfld = 0
1362 nullify(cbuf)
1363 mbuf = init
1364 allocate(cbuf(mbuf), stat = istat) ! allocate initial space for cbuf
1365 if (istat .ne. 0) then
1366 iret = 1
1367 return
1368 endif
1369
1370 ! Read up to 5000 bytes from the file into buffer cbread. lskip8
1371 ! must be set to the beginning of a GRIB2 message in the file.
1372 ibread8 = min(lgrib8, linmax)
1373 call bareadl(lugb, lskip8, ibread8, lbread8, cbread)
1374 if (lbread8 .ne. ibread8) then
1375 iret = 2
1376 return
1377 endif
1378
1379 ! Check GRIB version from section 0, must be 2.
1380 if (cbread(8) .ne. char(2)) then
1381 iret = 3
1382 return
1383 endif
1384
1385 ! Remember the GRIB version and discipline.
1386 cver = cbread(8)
1387 cdisc = cbread(7)
1388
1389 ! Read the length of section 1 from the file data buffer.
1390 call g2_gbytec1(cbread, lensec1, 16 * int1_bits, int4_bits)
1391 lensec1 = min(lensec1, int(ibread8, kind(lensec1)))
1392
1393 ! Copy section 1 values into cids array.
1394 cids(1:lensec1) = cbread(17:16 + lensec1)
1395 ! do i = 1, lensec1
1396 ! print *, i, ichar(cids(i))
1397 ! end do
1398
1399 ! Skip past section 1 in the data buffer.
1400 ibskip8 = lskip8 + 16_8 + int(lensec1, kind(8))
1401
1402 ! Loop through remaining sections creating an index for each
1403 ! field. This overwrites the cbread data buffer.
1404 ibread8 = max(5, mxbms)
1405 do
1406 ! Read 6 bytes from file into cbread buffer. (Why 6? Should this
1407 ! be 5?)
1408 call bareadl(lugb, ibskip8, ibread8, lbread8, cbread)
1409
1410 ! Check if the first 4 bytes are '7777', indicating end of
1411 ! message. If we found end of message, we are done.
1412 ctemp = cbread(1)//cbread(2)//cbread(3)//cbread(4)
1413 if (ctemp .eq. '7777') return ! end of message found
1414
1415 ! If we did not find end of message, check that we read 6 bytes.
1416 if (lbread8 .ne. ibread8) then
1417 iret = 2
1418 return
1419 endif
1420
1421 ! Read the 4-byte section length, and then the 1-byte section
1422 ! number.
1423 call g2_gbytec1(cbread, lensec, 0, int4_bits)
1424 call g2_gbytec1(cbread, numsec, int4_bits, int1_bits)
1425
1426 ! Based on the section number, generate index data for each
1427 ! section.
1428 if (numsec .eq. 2) then
1429 ! Save the location of the local use section in the message.
1430 loclus8 = ibskip8 - lskip8
1431 elseif (numsec .eq. 3) then
1432 ! For the GDS section, read the whole section into the cgds
1433 ! buffer.
1434 lengds8 = lensec
1435 cgds = char(0)
1436 call bareadl(lugb, ibskip8, lengds8, lbread8, cgds)
1437 if (lbread8 .ne. lengds8) then
1438 iret = 2
1439 return
1440 endif
1441 ! Remember the GDS location in the message.
1442 locgds8 = ibskip8 - lskip8
1443 elseif (numsec .eq. 4) then
1444 ! Having found the PDS, we write the beginning of the index
1445 ! record into the cindex buffer.
1446 cindex = char(0)
1447
1448 ! The first 4 bytes are the length of the index record. We
1449 ! don't know that yet, so skip it for now,
1450 mypos = int4_bits
1451
1452 ! Write the beginning of the index record. It contains bytes
1453 ! to skip in the file to reach the message, and offsets to
1454 ! each section within the message. Index version 1 uses 4-byte
1455 ! ints for these values, index version 2 uses 8-byte ints.
1456 if (idxver .eq. 1) then
1457 inc = 0
1458 lskip = int(lskip8, kind(4))
1459 call g2_sbytec1(cindex, lskip, mypos, int4_bits) ! bytes to skip
1460 mypos = mypos + int4_bits
1461 call g2_sbytec1(cindex, int(loclus8, kind(4)), mypos, int4_bits) ! location of local use
1462 mypos = mypos + int4_bits
1463 call g2_sbytec1(cindex, int(locgds8, kind(4)), mypos, int4_bits) ! location of gds
1464 mypos = mypos + int4_bits
1465 call g2_sbytec1(cindex, int(ibskip8 - lskip8, kind(4)), mypos, int4_bits) ! location of pds
1466#ifdef LOGGING
1467 write(g2_log_msg, *) ' writing pds location to index: mypos/8 ', mypos/8, &
1468 ' loc ', int(ibskip8 - lskip8, kind(4))
1469 call g2_log(4)
1470#endif
1471 mypos = mypos + int4_bits * 4 ! skip ahead in cbuf
1472 else
1473 inc = 28
1474 call g2_sbytec81(cindex, lskip8, mypos, int8_bits) ! bytes to skip
1475 mypos = mypos + int8_bits
1476 call g2_sbytec81(cindex, loclus8, mypos, int8_bits) ! location of local use
1477 mypos = mypos + int8_bits
1478 call g2_sbytec81(cindex, locgds8, mypos, int8_bits) ! location of gds
1479 mypos = mypos + int8_bits
1480 call g2_sbytec81(cindex, ibskip8 - lskip8, mypos, int8_bits) ! location of pds
1481 mypos = mypos + int8_bits * 4 ! skip ahead in cbuf
1482 endif
1483
1484#ifdef LOGGING
1485 write(g2_log_msg, *) ' writing total len to index: mypos/8 ', mypos/8, lgrib8
1486 call g2_log(4)
1487#endif
1488 ! These ints are the same size in index version 1 and 2. The
1489 ! mypos variable contains the proper offset, which is
1490 call g2_sbytec81(cindex, lgrib8, mypos, int8_bits) ! length of grib2
1491 mypos = mypos + int8_bits
1492 cindex((mypos / 8) + 1) = cver
1493 mypos = mypos + int1_bits
1494 cindex((mypos / 8) + 1) = cdisc
1495 mypos = mypos + int1_bits
1496 call g2_sbytec1(cindex, numfld + 1, mypos, int2_bits) ! field num
1497 mypos = mypos + int2_bits
1498
1499 ! Copy the section 1 values into the cindex buffer.
1500 cindex(ixids + 1 + inc:ixids + lensec1 + inc) = cids(1:lensec1)
1501 lindex = ixids + lensec1 + inc
1502
1503 ! Copy the GDS values into the cindex buffer.
1504 cindex(lindex + 1:lindex + lengds8) = cgds(1:lengds8)
1505 lindex = lindex + int(lengds8, kind(lindex))
1506
1507 ! Now read the PDS values from the file directly into cindex.
1508 ilnpds = lensec
1509 ilnpds8 = ilnpds
1510 call bareadl(lugb, ibskip8, ilnpds8, lbread8, cindex(lindex + 1))
1511 if (lbread8 .ne. ilnpds8) then
1512 iret = 2
1513 return
1514 endif
1515 lindex = lindex + ilnpds
1516 mypos = mypos + ilnpds
1517#ifdef LOGGING
1518 write(g2_log_msg, *) ' after writing pds location to index: mypos/8 ', mypos/8
1519 call g2_log(4)
1520#endif
1521 elseif (numsec .eq. 5) then
1522 ! Write the byte offset to the DRS section into the cindex buffer.
1523 !mypos = (IXSDR + inc) * INT1_BITS
1524#ifdef LOGGING
1525 write(g2_log_msg, *) ' before writing drs to index: ibskip8 - lskip8 ', ibskip8 - lskip8, ixdrs2
1526 call g2_log(4)
1527#endif
1528 ! Write the bytes to skip to the drs section into the index record.
1529 if (idxver .eq. 1) then
1530 call g2_sbytec1(cindex, int(ibskip8 - lskip8, kind(4)), ixdrs1 * int1_bits, int4_bits)
1531 else
1532 call g2_sbytec81(cindex, ibskip8 - lskip8, ixdrs2 * int1_bits, int8_bits) ! location of drs
1533 endif
1534
1535 ! Read the DRS section directly into the cindex buffer.
1536 ilndrs = lensec
1537 ilndrs8 = ilndrs
1538 call bareadl(lugb, ibskip8, ilndrs8, lbread8, cindex(lindex + 1))
1539 if (lbread8 .ne. ilndrs8) then
1540 iret = 2
1541 return
1542 endif
1543 lindex = lindex + ilndrs
1544 elseif (numsec .eq. 6) then
1545 ! Write the location of the BMS section in the message into
1546 ! the cindex buffer.
1547 indbmp = g2_mova2i(cbread(6))
1548#ifdef LOGGING
1549 write(g2_log_msg, *) ' section 6: indbmp', indbmp
1550 call g2_log(4)
1551#endif
1552
1553 if (indbmp .lt. 254) then
1554 if (idxver .eq. 1) then
1555 locbms = int(ibskip8 - lskip8, kind(4))
1556 call g2_sbytec1(cindex, locbms, ixbms1 * int1_bits, int4_bits) ! loc. of bms
1557 else
1558 locbms8 = ibskip8 - lskip8
1559 call g2_sbytec81(cindex, locbms8, ixbms2 * int1_bits, int8_bits) ! loc. of bms
1560 endif
1561#ifdef LOGGING
1562 write(g2_log_msg, *) ' section 6: locbms', locbms, 'locbms8', locbms8
1563 call g2_log(4)
1564#endif
1565 elseif (indbmp .eq. 254) then
1566 if (idxver .eq. 1) then
1567 call g2_sbytec1(cindex, locbms, ixbms1 * int1_bits, int4_bits) ! loc. of bms
1568 else
1569 call g2_sbytec81(cindex, locbms8, ixbms2 * int1_bits, int8_bits) ! loc. of bms
1570 endif
1571 elseif (indbmp .eq. 255) then
1572 if (idxver .eq. 1) then
1573 call g2_sbytec1(cindex, int(ibskip8 - lskip8, kind(4)), ixbms1 * int1_bits, int4_bits) ! loc. of bms
1574 else
1575 call g2_sbytec81(cindex, ibskip8 - lskip8, ixbms2 * int1_bits, int8_bits) ! loc. of bms
1576 endif
1577 endif
1578
1579 ! Copy 6 bytes of the BMS from data buffer to the cindex buffer.
1580 cindex(lindex + 1:lindex + mxbms) = cbread(1:mxbms)
1581 lindex = lindex + mxbms
1582
1583 ! The size of the index record is now known, so write it to
1584 ! the cindex buffer.
1585 call g2_sbytec1(cindex, lindex, 0, int4_bits) ! num bytes in index record
1586 elseif (numsec .eq. 7) then ! found data section
1587
1588#ifdef LOGGING
1589 write(g2_log_msg, *) ' writing offset to the data in cindex: ibskip8 - lskip8 ', ibskip8 - lskip8, &
1590 'IXDS2', ixds2
1591 call g2_log(4)
1592#endif
1593 ! Write the offset to the data section in the cindex buffer.
1594 if (idxver .eq. 1) then
1595 call g2_sbytec1(cindex, int(ibskip8 - lskip8, kind(4)), ixds1 * int1_bits, int4_bits)
1596 else
1597 call g2_sbytec81(cindex, ibskip8 - lskip8, ixds2 * int1_bits, int8_bits)
1598 endif
1599
1600 ! Increment the field count.
1601 numfld = numfld + 1
1602
1603 ! Allocate more space in cbuf if necessary. The index record
1604 ! will be copied there.
1605 if (lindex + mlen .gt. mbuf) then
1606 newsize = max(mbuf + next, mbuf + lindex)
1607 call realloc(cbuf, mlen, newsize, istat)
1608 if (istat .ne. 0) then
1609 numfld = numfld - 1
1610 iret = 4
1611 return
1612 endif
1613 mbuf = newsize
1614 endif
1615
1616 ! Copy the index record into cbuf.
1617 cbuf(mlen + 1:mlen + lindex) = cindex(1:lindex)
1618 mlen = mlen + lindex
1619 else
1620 ! Unrecognized section.
1621 iret = 5
1622 return
1623 endif
1624
1625 ! Skip past this section in the data buffer.
1626 ibskip8 = ibskip8 + lensec
1627 enddo ! next section
1628end subroutine ix2gb2
1629
1636subroutine gf_finalize(iret)
1637 implicit none
1638
1639 integer, intent(out) :: iret
1640 character(len = 1), pointer, dimension(:) :: cindex
1641 integer :: nlen, nnum
1642
1643 ! Declare interfaces (required for cbuf pointer).
1644 interface
1645 subroutine getidx(lugb, lugi, cbuf, nlen, nnum, irgi)
1646 character(len = 1), pointer, dimension(:) :: cbuf
1647 integer, intent(in) :: lugb, lugi
1648 integer, intent(out) :: nlen, nnum, irgi
1649 end subroutine getidx
1650 end interface
1651
1652 ! Call getidx with 0 for the first parameter, ensuring that the
1653 ! internal memory is freed.
1654 call getidx(0, 0, cindex, nlen, nnum, iret)
1655
1656end subroutine gf_finalize
1657
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_sbytec1(out, in, iskip, nbits)
Put one arbitrary sized (up to 32 bits) values from an integer scalar into a packed bit string,...
Definition g2bytes.F90:337
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 g2_create_index(lugb, lugi, idxver, filename, iret)
Create a version 1 or 2 index file for a GRIB2 file.
Definition g2index.F90:21
subroutine ix2gb2(lugb, lskip8, idxver, lgrib8, cbuf, numfld, mlen, iret)
Generate a version 1 or 2 index record for each field in a GRIB2 message.
Definition g2index.F90:1260
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:820
subroutine getidx(lugb, lugi, cindex, nlen, nnum, iret)
Find, read or generate a version 1 GRIB2 index for a GRIB2 file (which must be < 2 GB).
Definition g2index.F90:202
subroutine g2_write_index_headers(lugi, nlen, nnum, idxver, filename)
Write index headers.
Definition g2index.F90:107
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:653
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:592
subroutine ixgb2(lugb, lskip, lgrib, cbuf, numfld, mlen, iret)
Generate a version 1 index record for each field in a GRIB2 message.
Definition g2index.F90:1196
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:270
subroutine gf_finalize(iret)
Free all memory associated with the library.
Definition g2index.F90:1637
subroutine getg2i(lugi, cbuf, nlen, nnum, iret)
Read a version 1 index file and return its contents.
Definition g2index.F90:434
subroutine getg2i2(lugi, cbuf, idxver, nlen, nnum, iret)
Read a version 1 or 2 index file and return its contents.
Definition g2index.F90:516
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:933
subroutine gf_unpack1(cgrib, lcgrib, iofst, ids, idslen, ierr)
Unpack Section 1 ([Identification Section] (https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_doc/gr...
Definition g2unpack.F90:42
subroutine gf_unpack5(cgrib, lcgrib, iofst, ndpts, idrsnum, idrstmpl, mapdrslen, ierr)
Unpack Section 5 ([Data Representation Section] (https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_d...
Definition g2unpack.F90:467
subroutine gf_unpack4(cgrib, lcgrib, iofst, ipdsnum, ipdstmpl, mappdslen, coordlist, numcoord, ierr)
Unpack Section 4 ([Product Definition Section] (https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_do...
Definition g2unpack.F90:334
subroutine gf_unpack3(cgrib, lcgrib, iofst, igds, igdstmpl, mapgridlen, ideflist, idefnum, ierr)
Unpack Section 3 ([Grid Definition Section] (https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_doc/g...
Definition g2unpack.F90:184
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
Reallocate memory, preserving contents.
Definition realloc.F90:12
subroutine skgb8(lugb, iseek8, mseek8, lskip8, lgrib8)
Search a file for the next GRIB1 or GRIB2 message.
Definition skgb.F90:53