Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

adding getgb2i2() to support new index format for files > 2 GB #646

Merged
merged 8 commits into from
Mar 20, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
129 changes: 111 additions & 18 deletions src/g2getgb2.F90
Original file line number Diff line number Diff line change
Expand Up @@ -3,11 +3,93 @@
!> file.
!> @author Ed Hartnett @date Mar 6, 2024

!> Find and unpack a GRIB2 message in a file. It reads a GRIB index
!> file (or optionally the GRIB file itself) to get the index buffer
!> (i.e. table of contents) for the GRIB file.
!> This is a legacy version of getgb2i2(). It finds and unpacks a
!> GRIB2 message in a file, using an version 1 index record which is
!> either found in memory, read from an index file, or generated.
!>
!> Find in the index buffer a reference to the GRIB field requested.
!> Users of this routine must:
!> 1. include a use grib_mod
!> 2. call gf_free() on the @ref grib_mod::gribfield parameter
!> 3. free library memory with gf_finalize()
!>
!> For more details, see getgb2i2().
!>
!> @param[in] lugb integer unit of the unblocked grib data file.
!> File must be opened with [baopen() or baopenr()]
!> (https://noaa-emc.github.io/NCEPLIBS-bacio/) before calling
!> this routine.
!> @param[in] lugi integer unit of the unblocked grib index file.
!> If nonzero, file must be opened with [baopen() or baopenr()]
!> (https://noaa-emc.github.io/NCEPLIBS-bacio/) before
!> calling this routine. lugi may be:
!> - > 0 read index from index file lugi, if index doesn"t already exist.
!> - = 0 to get index buffer from the grib file, if index
!> doesn"t already exist.
!> - < 0 force reread of index from index file abs(lugi).
!> - = lugb force regeneration of index from GRIB2 file lugb.
!> @param[in] j integer number of fields to skip (0 to search from
!> beginning).
!> @param[in] jdisc GRIB2 discipline number of requested field:
!> --1 accept any discipline
!> - 0 meteorological products
!> - 1 hydrological products
!> - 2 land surface products
!> - 3 space products
!> - 10 oceanographic products
!> @param[in] JIDS integer array of values in the identification section
!> (=-9999 for wildcard).
!> @param[in] jpdtn integer product definition template number (n)
!> (if = -1, don't bother matching pdt - accept any)
!> @param[in] JPDT integer array of values defining the product definition
!> template 4.n of the field for which to search (=-9999 for wildcard)
!> @param[in] jgdtn integer grid definition template number (m)
!> (if = -1, don't bother matching gdt - accept any )
!> @param[in] JGDT integer array of values defining the grid definition
!> template 3.m of the field for which to search (=-9999 for wildcard)
!> @param[in] unpack logical value indicating whether to unpack bitmap/data
!> - .TRUE. unpack bitmap and data values
!> - .FALSE. do not unpack bitmap and data values
!> @param[out] k integer field number unpacked
!> @param[out] gfld derived type @ref grib_mod::gribfield.
!> @param[out] iret integer return code
!> - 0 all ok
!> - 96 error reading index
!> - 97 error reading grib file
!> - 99 request not found
!> - other gf_getfld grib2 unpacker return code
!>
!> @author Mark Iredell, Ed Hartnett @date 1994-04-01
subroutine getgb2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
unpack, k, gfld, iret)
use grib_mod
implicit none

integer, intent(in) :: lugb, lugi, j, jdisc
integer, dimension(:) :: jids(*)
integer, intent(in) :: jpdtn
integer, dimension(:) :: jpdt(*)
integer, intent(in) :: jgdtn
integer, dimension(:) :: jgdt(*)
logical, intent(in) :: unpack
integer, intent(out) :: k
type(gribfield), intent(out) :: gfld
integer, intent(out) ::iret
integer :: idxver = 1

call getgb2i2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
unpack, idxver, k, gfld, iret)
end subroutine getgb2

!> Find and unpack a GRIB2 message in a file, using an version 1 or 2
!> index record which is either found in memory, read from an index
!> file, or generated.
!>
!> This subroutine is similar to getgb2(), but allows the caller to
!> specify the index version of the generated index (if one is
!> generated). Use index version 2 for all new code, as it handles all
!> GRIB2 files, including files > 2 GB. Use index version 1 for
!> backward compatibility, but it does not work past the first 2 GB of
!> the file.
!>
!> The GRIB field request specifies the number of fields to skip and
!> the unpacked identification section, grid definition template and
Expand All @@ -28,7 +110,7 @@
!> Derived type @ref grib_mod::gribfield contains pointers to many
!> arrays of data. Users must free this memory by calling gf_free().
!>
!> This subroutine calls getidx(), which allocates memory and stores
!> This subroutine calls getidx2(), which allocates memory and stores
!> the resulting pointers in an array that is a Fortran "save"
!> variable. The result is that the memory will not be freed by the
!> library and cannot be reached by the caller. To free this memory
Expand Down Expand Up @@ -95,6 +177,7 @@
!> @param[in] unpack logical value indicating whether to unpack bitmap/data
!> - .TRUE. unpack bitmap and data values
!> - .FALSE. do not unpack bitmap and data values
!> @param[in] idxver Index version of the cindex buffer.
!> @param[out] k integer field number unpacked
!> @param[out] gfld derived type @ref grib_mod::gribfield.
!> @param[out] iret integer return code
Expand All @@ -104,19 +187,25 @@
!> - 99 request not found
!> - other gf_getfld grib2 unpacker return code
!>
!> @author Mark Iredell @date 1994-04-01
subroutine getgb2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
unpack, k, gfld, iret)
!> @author Ed Hartnett @date 2024-03-19
subroutine getgb2i2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
unpack, idxver, k, gfld, iret)
use grib_mod
implicit none

integer, intent(in) :: lugb, lugi, j, jdisc, jpdtn, jgdtn
integer, dimension(:) :: jids(*), jpdt(*), jgdt(*)
integer, intent(in) :: lugb, lugi, j, jdisc
integer, dimension(:) :: jids(*)
integer, intent(in) :: jpdtn
integer, dimension(:) :: jpdt(*)
integer, intent(in) :: jgdtn
integer, dimension(:) :: jgdt(*)
logical, intent(in) :: unpack
integer, intent(out) :: k, iret
integer, intent(in) :: idxver
integer, intent(out) :: k
type(gribfield), intent(out) :: gfld
integer, intent(out) ::iret
character(len = 1), pointer, dimension(:) :: cbuf
integer :: nnum, nlen, lpos, jk, irgi, irgs, idxver
integer :: nnum, nlen, lpos, jk, irgi, irgs

! Declare interfaces (required for cbuf pointer).
interface
Expand Down Expand Up @@ -156,16 +245,18 @@ subroutine getgb2r2(lugb, idxver, cindex, gfld, iret)
end subroutine getgb2r2
end interface

! Determine whether index buffer needs to be initialized.
! Fill cbuf with the index records of this file, by recalling them
! from memory, reading them from the index file, or generating them
! from the data file.
irgi = 0
idxver = 2
call getidx2(lugb, lugi, idxver, cbuf, nlen, nnum, irgi)
if (irgi .gt. 1) then
iret = 96
return
endif

! search index buffer.
! Search the index in cbuf for the first message which meets our
! search criteria.
call getgb2s2(cbuf, idxver, nlen, nnum, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, jk, &
gfld, lpos, irgs)
if (irgs .ne. 0) then
Expand All @@ -174,15 +265,17 @@ end subroutine getgb2r2
return
endif

! Read local use section, if available.
! Read local use section for the selected GRIB2 message, if available.
call getgb2l2(lugb, idxver, cbuf(lpos), gfld, iret)

! Read and unpack grib record.
! Read and unpack grib record for the selected GRIB2 message.
if (unpack) then
call getgb2r2(lugb, idxver, cbuf(lpos), gfld, iret)
endif

! Return the number of the unpacked field to the caller.
k = jk
end subroutine getgb2
end subroutine getgb2i2

!> Read and unpack a local use section from a GRIB2 index record
!> (index format 1) and GRIB2 file.
Expand Down
5 changes: 3 additions & 2 deletions src/g2index.F90
Original file line number Diff line number Diff line change
Expand Up @@ -360,8 +360,9 @@ end subroutine getg2i2r
return
endif

! Either read index record from index file, or generate it from the
! GRIB2 file.
! Either read index from index file, or generate it from the GRIB2
! file. This is an index for all messages in the file, each message
! gets an index record, all stuffed into idxlist(lugbb)%cbuf.
irgi = 0
if (lux .gt. 0) then
call getg2i2(lux, idxlist(lugb)%cbuf, idxver, nlen, nnum, irgi)
Expand Down
77 changes: 58 additions & 19 deletions tests/test_files.F90
Original file line number Diff line number Diff line change
Expand Up @@ -12,31 +12,70 @@ program test_files
integer, dimension(200) :: jids, jpdt, jgdt
logical :: unpack = .true.
type(gribfield) :: gfld
integer :: iret
integer :: expected_idsect(13) = (/ 7, 0, 2, 1, 1, 2021, 11, 30, 0, 0, 0, 0, 1 /)
integer :: expected_igdtmpl(19) = (/ 6, 0, 0, 0, 0, 0, 0, 241, 151, 0, 0, 50000000, &
210000000, 48, 25000000, 250000000, 166667, 166667, 0 /)
integer :: expected_ipdtmpl(15) = (/ 2, 1, 2, 0, 11, 0, 0, 1, 0, 1, 0, 1, 255, 0, 0 /)
integer :: expected_idrtmpl(7) = (/ 1092616192, 0, 2, 11, 0, 0, 255 /)
integer :: i, idxver, iret, j, k

print *, 'Testing reading GRIB2 file gdaswave.t00z.wcoast.0p16.f000.grib2...'

! Open the file.
call baopenr(lugb, 'gdaswave.t00z.wcoast.0p16.f000.grib2', iret)
if (iret .ne. 0) stop 2
do j = 1, 2
idxver = j
! Open the file.
call baopenr(lugb, 'gdaswave.t00z.wcoast.0p16.f000.grib2', iret)
if (iret .ne. 0) stop 2

! Learn about the file.
jids = -9999
jpdt = -9999
jgdt = -9999
call getgb2(lugb, 0, jskp, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
unpack, jskp, gfld, iret)
if (iret .ne. 0) stop 3
! Learn about the file.
jids = -9999
jpdt = -9999
jgdt = -9999
call getgb2i2(lugb, 0, jskp, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
unpack, idxver, k, gfld, iret)
if (iret .ne. 0) stop 3
! Make sure we found message 1 in the file.
if (k .ne. 1) stop 4
if (gfld%version .ne. 2 .or. gfld%discipline .ne. 0 .or. gfld%idsectlen .ne. 13 .or. &
gfld%locallen .ne. 0 .or. gfld%ifldnum .ne. 1 .or. gfld%griddef .ne. 0 .or. &
gfld%ngrdpts .ne. 36391 .or. gfld%numoct_opt .ne. 0 .or. gfld%interp_opt .ne. 0 .or. &
gfld%num_opt .ne. 0 .or. gfld%igdtnum .ne. 0 .or. gfld%ipdtnum .ne. 0 .or. &
gfld%ipdtlen .ne. 15 .or. gfld%ndpts .ne. 11041 .or. gfld%idrtnum .ne. 40 .or. &
gfld%idrtlen .ne. 7 .or. gfld%unpacked .neqv. .false. .or. gfld%expanded .neqv. .true. .or. &
gfld%ibmap .ne. 0) stop 10
do i = 1, 13
if (gfld%idsect(i) .ne. expected_idsect(i)) stop 20
end do
do i = 1, 19
if (gfld%igdtmpl(i) .ne. expected_igdtmpl(i)) stop 30
end do
do i = 1, 15
if (gfld%ipdtmpl(i) .ne. expected_ipdtmpl(i)) then
print *, 'got gfld%ipdtmpl', gfld%ipdtmpl
print *, 'expected ', expected_ipdtmpl
stop 40
endif
end do
do i = 1, 7
if (gfld%idrtmpl(i) .ne. expected_idrtmpl(i)) stop 50
end do
! do i = 1, 100
! print *, gfld%fld(1)
! end do

! Close the file.
call baclose(lugb, iret)
if (iret .ne. 0) stop 200
! Close the file.
call baclose(lugb, iret)
if (iret .ne. 0) stop 200

! Free memory.
call gf_free(gfld)

end do

! Free internal library memory.
call gf_finalize(iret)
if (iret .ne. 0) stop 5

! Free memory.
call gf_free(gfld)
call gf_finalize(iret)
if (iret .ne. 0) stop 5

print *, 'OK!'
print *, 'SUCCESS!'
end program test_files
4 changes: 2 additions & 2 deletions tests/test_getgb2.F90
Original file line number Diff line number Diff line change
Expand Up @@ -47,8 +47,8 @@ program test_getgb2
if (iret .ne. 0) stop 3

! Read a field from the test file.
call getgb2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
unpack, k, gfld, iret)
call getgb2i2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
unpack, 2, k, gfld, iret)
if (iret .ne. 0) stop 4

! Check results.
Expand Down
Loading