NCEPLIBS-g2  3.4.8
getg2ir.F90
Go to the documentation of this file.
1 
4 
37 SUBROUTINE getg2ir(LUGB, MSK1, MSK2, MNUM, CBUF, NLEN, NNUM, NMESS, IRET)
38  USE re_alloc ! NEEDED FOR SUBROUTINE REALLOC
39  implicit none
40 
41  CHARACTER(LEN = 1), POINTER, DIMENSION(:) :: CBUF
42  INTEGER, INTENT(IN) :: LUGB, MSK1, MSK2, MNUM
43  INTEGER, INTENT(OUT) :: NLEN, NNUM, NMESS, IRET
44  CHARACTER(LEN = 1), POINTER, DIMENSION(:) :: CBUFTMP
45  integer :: nbytes, newsize, next, numfld, m, mbuf, lskip, lgrib
46  integer :: istat, iseek, init, iret1
47  parameter(init = 50000, next = 10000)
48 
49  INTERFACE ! REQUIRED FOR CBUF POINTER
50  SUBROUTINE ixgb2(LUGB, LSKIP, LGRIB, CBUF, NUMFLD, MLEN, IRET)
51  INTEGER, INTENT(IN) :: LUGB, LSKIP, LGRIB
52  CHARACTER(LEN = 1), POINTER, DIMENSION(:) :: CBUF
53  INTEGER, INTENT(OUT) :: NUMFLD, MLEN, IRET
54  END SUBROUTINE ixgb2
55  END INTERFACE
56 
57  ! Initialize.
58  iret = 0
59  NULLIFY(cbuf)
60  mbuf = init
61  ALLOCATE(cbuf(mbuf), stat = istat) ! Allocate initial space for cbuf.
62  IF (istat .NE. 0) THEN
63  iret = 2
64  RETURN
65  ENDIF
66 
67  ! Search for first grib message.
68  iseek = 0
69  CALL skgb(lugb, iseek, msk1, lskip, lgrib)
70  DO m = 1, mnum
71  IF(lgrib.GT.0) THEN
72  iseek = lskip + lgrib
73  CALL skgb(lugb, iseek, msk2, lskip, lgrib)
74  ENDIF
75  ENDDO
76 
77  ! Get index records for every grib message found.
78  nlen = 0
79  nnum = 0
80  nmess = mnum
81  DO WHILE(iret .EQ. 0 .AND. lgrib .GT. 0)
82  CALL ixgb2(lugb, lskip, lgrib, cbuftmp, numfld, nbytes, iret1)
83  IF (iret1 .NE. 0) print *, ' SAGT ', numfld, nbytes, iret1
84  IF((nbytes + nlen) .GT. mbuf) THEN ! Allocate more space, if necessary.
85  newsize = max(mbuf + next, mbuf + nbytes)
86  CALL realloc(cbuf, nlen, newsize, istat)
87  IF (istat .NE. 0) THEN
88  iret = 1
89  RETURN
90  ENDIF
91  mbuf = newsize
92  ENDIF
93 
94  ! If index records were returned in cbuftmp from ixgb2,
95  ! copy cbuftmp into cbuf, then deallocate cbuftmp when done.
96  IF (ASSOCIATED(cbuftmp)) THEN
97  cbuf(nlen + 1 : nlen + nbytes) = cbuftmp(1 : nbytes)
98  DEALLOCATE(cbuftmp, stat = istat)
99  IF (istat.NE.0) THEN
100  print *, ' deallocating cbuftmp ... ', istat
101  iret = 3
102  RETURN
103  ENDIF
104  NULLIFY(cbuftmp)
105  nnum = nnum + numfld
106  nlen = nlen + nbytes
107  nmess = nmess + 1
108  ENDIF
109 
110  ! Look for next grib message.
111  iseek = lskip + lgrib
112  CALL skgb(lugb, iseek, msk2, lskip, lgrib)
113  ENDDO
114 END SUBROUTINE getg2ir
subroutine getg2ir(LUGB, MSK1, MSK2, MNUM, CBUF, NLEN, NNUM, NMESS, IRET)
Generate an index record for a message in a GRIB2 file.
Definition: getg2ir.F90:38
subroutine ixgb2(LUGB, LSKIP, LGRIB, CBUF, NUMFLD, MLEN, IRET)
Generate an index record for each field in a GRIB2 message.
Definition: ixgb2.F90:56
Reallocate memory, preserving contents.
Definition: realloc.f:12
subroutine skgb(lugb, iseek, mseek, lskip, lgrib)
Search a file for the next GRIB1 or GRIB2 message.
Definition: skgb.F90:21