diff --git a/src/g2getgb2.F90 b/src/g2getgb2.F90 index 8cbd15fe..314887b3 100644 --- a/src/g2getgb2.F90 +++ b/src/g2getgb2.F90 @@ -1055,6 +1055,7 @@ subroutine getgb2rp2(lugb, idxver, cindex, extract, gribm, leng, iret) call g2_gbytec(cindex, iskp2, mypos, INT4_BITS) ! bytes to skip for section 2 mypos = mypos + INT4_BITS iskp2_8 = iskp2 + mypos = mypos + 32 * INT1_BITS ! skip ahead in the cindex else inc = 12 call g2_gbytec8(cindex, iskip8, mypos, INT8_BITS) ! bytes to skip in file @@ -1062,6 +1063,7 @@ subroutine getgb2rp2(lugb, idxver, cindex, extract, gribm, leng, iret) iskip = int(iskip8, kind(4)) call g2_gbytec8(cindex, iskp2_8, mypos, INT8_BITS) ! bytes to skip for section 2 mypos = mypos + INT8_BITS + mypos = mypos + 36 * INT1_BITS ! skip ahead in the cindex endif if (iskp2_8 .gt. 0) then call bareadl(lugb, iskip8 + iskp2_8, 4_8, lread8, ctemp) @@ -1072,7 +1074,6 @@ subroutine getgb2rp2(lugb, idxver, cindex, extract, gribm, leng, iret) else len2 = 0 endif - mypos = mypos + 32 * INT1_BITS ! skip ahead in the cindex call g2_gbytec(cindex, len1, mypos, INT4_BITS) ! length of section 1 ipos = 44 + len1 mypos = mypos + len1 * INT1_BITS ! skip ahead in the cindex @@ -1173,13 +1174,22 @@ subroutine getgb2rp2(lugb, idxver, cindex, extract, gribm, leng, iret) if (idxver .eq. 1) then call g2_gbytec(cindex, iskip, mypos, INT4_BITS) ! bytes to skip in file mypos = mypos + INT4_BITS + mypos = mypos + 6 * INT4_BITS iskip8 = iskip else call g2_gbytec8(cindex, iskip8, mypos, INT8_BITS) ! bytes to skip in file mypos = mypos + INT8_BITS + mypos = mypos + 2 * INT8_BITS + 4 * INT4_BITS endif - mypos = mypos + 7 * INT4_BITS - call g2_gbytec(cindex, leng, mypos, INT4_BITS) ! length of grib message + + call g2_gbytec8(cindex, leng8, mypos, INT8_BITS) ! length of grib message + leng = int(leng8, kind(4)) +#ifdef LOGGING + write(g2_log_msg, *) ' iskip8 ', iskip8, ' mypos/8 ', mypos/8, & + ' leng ', leng + call g2_log(2) +#endif + if (.not. associated(gribm)) allocate(gribm(leng)) leng8 = leng call bareadl(lugb, iskip8, leng8, lread8, gribm) diff --git a/src/g2index.F90 b/src/g2index.F90 index 55d86ed2..e09f0acd 100644 --- a/src/g2index.F90 +++ b/src/g2index.F90 @@ -523,7 +523,7 @@ subroutine getg2i2(lugi, cbuf, idxver, nlen, nnum, iret) #ifdef LOGGING ! Log results for debugging. - write(g2_log_msg, '(a, i2, a, i1)') 'getg2i2: lugi ', lugi, ' idxver ', idxver + write(g2_log_msg, *) 'getg2i2: lugi ', lugi call g2_log(1) #endif @@ -1300,7 +1300,7 @@ subroutine ix2gb2(lugb, lskip8, idxver, lgrib8, cbuf, numfld, mlen, iret) #ifdef LOGGING ! Log results for debugging. - write(g2_log_msg, '(a, i2, a, i7, a, i1)') 'ix2gb2: lugb ', lugb, ' lskip8 ', lskip8, & + write(g2_log_msg, *) 'ix2gb2: lugb ', lugb, ' lskip8 ', lskip8, & ' idxver ', idxver call g2_log(1) #endif diff --git a/tests/CMakeLists.txt b/tests/CMakeLists.txt index 8de931d3..08a1bdc2 100644 --- a/tests/CMakeLists.txt +++ b/tests/CMakeLists.txt @@ -116,6 +116,7 @@ if(FTP_TEST_FILES) create_test(test_getg2i2r ${kind}) create_test(test_getidx ${kind}) create_test(test_getgb2rp ${kind}) + create_test(test_getgb2rp2 ${kind}) create_test(test_getgb2rp_2 ${kind}) create_test(test_getgb2s ${kind}) create_test(test_getgb2p ${kind}) diff --git a/tests/test_getgb2rp2.F90 b/tests/test_getgb2rp2.F90 new file mode 100644 index 00000000..5a44eaf9 --- /dev/null +++ b/tests/test_getgb2rp2.F90 @@ -0,0 +1,98 @@ +! This is a test program for NCEPLIBS-g2. +! +! This program tests getg2rp2.F90 +! +! Ed Hartnett 5/14/24 +program test_getgb2rp2 + use g2logging + use bacio_module + implicit none + + integer :: lugi + character(len=1), pointer, dimension(:) :: cbuf(:) + integer :: lugb = 3 + integer :: nlen, nnum, iret + logical :: extract + integer (kind = 8) :: leng8 + character(len=1), pointer, dimension(:) :: gribm + integer :: idxver, i + + ! Interfaces are needed due to pointers in the parameter lists. + interface + subroutine getidx2(lugb, lugi, idxver, cindex, nlen, nnum, iret) + integer, intent(in) :: lugb, lugi + integer, intent(inout) :: idxver + character(len = 1), pointer, dimension(:) :: cindex + integer, intent(out) :: nlen, nnum, iret + end subroutine getidx2 + subroutine getgb2rp2(lugb, idxver, cindex, extract, gribm, leng8, iret) + integer, intent(in) :: lugb, idxver + character(len = 1), intent(in) :: cindex(*) + logical, intent(in) :: extract + character(len = 1), pointer, dimension(:) :: gribm + integer(kind = 8), intent(out) :: leng8 + integer, intent(out) :: iret + end subroutine getgb2rp2 + end interface + + print *, 'Testing the getgb2rp() subroutine - expect and ignore error messages during test...' + + do i = 1, 2 + ! Open a real GRIB2 file. + print *, 'Indexing a real GRIB2 file WW3_Regional_US_West_Coast_20220718_0000.grib2...' + call baopenr(lugb, "data/WW3_Regional_US_West_Coast_20220718_0000.grib2", iret) + if (iret .ne. 0) stop 100 + + idxver = i + lugi = 0 + leng8 = 0 + g2_log_level = 1 + + !lugi = lugb ! Force regeneration of index from GRIB2 file. + call getidx2(lugb, lugi, idxver, cbuf, nlen, nnum, iret) + if (iret .ne. 0) stop 101 + if (nnum .ne. 688) stop 102 + if (idxver .eq. 1) then + if (nlen .ne. 137600) then + print *, nlen + stop 103 + endif + else + if (nlen .ne. 145856) then + print *, nlen + stop 104 + endif + endif + print *, 'nlen, nnum: ', nlen, nnum + + ! Extract the whole message. + extract = .false. + nullify(gribm) + call getgb2rp2(lugb, idxver, cbuf, extract, gribm, leng8, iret) + print *, 'leng8 ', leng8 + if (leng8 .ne. 11183) stop 110 + ! Deallocate buffer that got GRIB message. + deallocate(gribm) + + ! Extract just the field (same result). + extract = .true. + nullify(gribm) + call getgb2rp2(lugb, idxver, cbuf, extract, gribm, leng8, iret) + print *, 'leng8 ', leng8 + if (leng8 .ne. 11183) stop 112 + ! Deallocate buffer that got GRIB message. + deallocate(gribm) + + call baclose(lugb, iret) + if (iret .ne. 0) stop 199 + end do + + ! Deallocate the buffer that holds index. + deallocate(cbuf) + + ! Free library memory. + ! call gf_finalize() + + print *, 'SUCCESS!...' + +end program test_getgb2rp2