NCEPLIBS-g2  3.4.8
getgb2.F90
Go to the documentation of this file.
1 
4 
112 SUBROUTINE getgb2(LUGB, LUGI, J, JDISC, JIDS, JPDTN, JPDT, JGDTN, JGDT, &
113  UNPACK, K, GFLD, IRET)
114  USE grib_mod
115  implicit none
116 
117  INTEGER, INTENT(IN) :: LUGB, LUGI, J, JDISC, JPDTN, JGDTN
118  INTEGER, DIMENSION(:) :: JIDS(*), JPDT(*), JGDT(*)
119  LOGICAL, INTENT(IN) :: UNPACK
120  INTEGER, INTENT(OUT) :: K, IRET
121  TYPE(gribfield), INTENT(OUT) :: GFLD
122  CHARACTER(LEN = 1), POINTER, DIMENSION(:) :: CBUF
123  integer :: nnum, nlen, lpos, jk, irgi, irgs
124 
125  ! Declare interfaces (required for cbuf pointer).
126  INTERFACE
127  SUBROUTINE getidx(LUGB, LUGI, CBUF, NLEN, NNUM, IRGI)
128  CHARACTER(LEN = 1), POINTER, DIMENSION(:) :: CBUF
129  INTEGER, INTENT(IN) :: LUGB, LUGI
130  INTEGER, INTENT(OUT) :: NLEN, NNUM, IRGI
131  END SUBROUTINE getidx
132  END INTERFACE
133 
134  ! Determine whether index buffer needs to be initialized.
135  irgi = 0
136  CALL getidx(lugb, lugi, cbuf, nlen, nnum, irgi)
137  IF (irgi .GT. 1) THEN
138  iret = 96
139  RETURN
140  ENDIF
141 
142  ! Search index buffer.
143  CALL getgb2s(cbuf, nlen, nnum, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, jk, &
144  gfld, lpos, irgs)
145  IF (irgs .NE. 0) THEN
146  iret = 99
147  CALL gf_free(gfld)
148  RETURN
149  ENDIF
150 
151  ! Read local use section, if available.
152  CALL getgb2l(lugb, cbuf(lpos), gfld, iret)
153 
154  ! Read and unpack grib record.
155  IF (unpack) THEN
156  CALL getgb2r(lugb, cbuf(lpos), gfld, iret)
157  ENDIF
158  k = jk
159 END SUBROUTINE getgb2
subroutine getgb2(LUGB, LUGI, J, JDISC, JIDS, JPDTN, JPDT, JGDTN, JGDT, UNPACK, K, GFLD, IRET)
Find and unpack a GRIB2 message in a file.
Definition: getgb2.F90:114
subroutine getgb2l(LUGB, CINDEX, GFLD, IRET)
This subroutine reads and unpacks a local use section from a GRIB2 message.
Definition: getgb2l.F90:37
subroutine getgb2r(LUGB, CINDEX, GFLD, IRET)
Read and unpack sections 6 and 7 from a GRIB2 message.
Definition: getgb2r.F90:40
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: getgb2s.F90:80
subroutine getidx(lugb, lugi, cindex, nlen, nnum, iret)
Find, read or generate a GRIB2 index for the GRIB2 file associated with unit lugb.
Definition: getidx.F90:44
subroutine gf_free(gfld)
Free memory that was used to store array values in derived type grib_mod::gribfield.
Definition: gf_free.F90:12
This Fortran module contains the declaration of derived type gribfield.
Definition: gribmod.F90:10