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 tests #173

Merged
merged 22 commits into from
Sep 27, 2023
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
10 changes: 5 additions & 5 deletions src/w3fi73.f
Original file line number Diff line number Diff line change
Expand Up @@ -12,11 +12,11 @@
C> @param[in] IBFLAG
C> - 0, if bit map supplied by user
C> - #, Number of predefined center bit map
C> @param[in] IBMAP Integer array containing user bit map
C> @param[in] IBLEN Length of bit map
C> @param[out] BMS Completed grib bit map section
C> @param[out] LENBMS Length of bit map section
C> @param[out] IER 0 normal exit, 8 = ibmap values are all zero
C> @param[in] IBMAP Integer array containing user bit map.
C> @param[in] IBLEN Length of bit map.
C> @param[out] BMS Completed grib bit map section.
C> @param[out] LENBMS Length of bit map section in bytes.
C> @param[out] IER 0 normal exit, 8 = ibmap values are all zero.
C>
C> @author M. Farley @date 1992-07-01
SUBROUTINE W3FI73 (IBFLAG,IBMAP,IBLEN,BMS,LENBMS,IER)
Expand Down
3 changes: 3 additions & 0 deletions tests/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -16,7 +16,10 @@ endfunction()
# These are the tests.
foreach(kind ${kinds})
create_test(test_summary ${kind})
create_test(test_gbytec ${kind})
create_test(test_w3tagb ${kind})
create_test(test_w3fi71 ${kind})
# create_test(test_w3fi72 ${kind})
create_test(test_w3fi73 ${kind})
create_test(test_w3fi74 ${kind})
endforeach()
186 changes: 186 additions & 0 deletions tests/test_gbytec.F90
Original file line number Diff line number Diff line change
@@ -0,0 +1,186 @@
! This is a test program for the NCEPLIBS-w3emc project.
!
! This program tests the sbyte subroutines. The sbyte code is a
! duplicate of code in the NCEPLIBS-g2 project, and this test was
! developed for the NCEPLIBS-g2 project, and ported to NCEPLIBS-w3emc.
!
! See https://github.com/NOAA-EMC/NCEPLIBS-g2/issues/568 and
! https://github.com/NOAA-EMC/NCEPLIBS-w3emc/issues/202.
!
! Ed Hartnett, 9/27/23
program test_gbytec
implicit none

character*1 :: out(1)
character*1 :: out4(4)
character*1 :: out5(5)
character*1 :: out8(8)
character*1 :: out10(10)
integer, parameter :: n = 1
integer :: in(n)
real :: r_in(n)
integer, parameter :: n2 = 2
integer :: in2(n2)
real :: r_in2(n2)
integer, parameter :: n5 = 5
integer :: in5(n5)
integer :: iskip = 0
integer :: nbits = 8
integer :: nskip = 0
integer :: i
integer :: num

print *, 'Testing sbyte subroutines.'

print *, 'Testing sbytec()...'
in(1) = 3
out(1) = char(0)
call sbytec(out, in, iskip, nbits)
if (ichar(out(1)) .ne. in(1)) stop 10

print *, 'Testing sbytesc()...'
in(1) = 3
out(1) = char(0)
call sbytesc(out, in, iskip, nbits, nskip, n)
if (ichar(out(1)) .ne. in(1)) stop 20

! THIS will pack the numbers 1 and 2 into the first two chars of the
! buffer. The rest of the output buffer will remain zeros.
print *, 'Testing g2_sbytesc() packing 2 values...'
in2(1) = 1
in2(2) = 2
do i = 1, 8
out8(i) = char(0)
end do
nbits = 8
call sbytesc(out8, in2, iskip, nbits, nskip, n2)
do i = 1, 8
if (i .le. 2) then
if (ichar(out8(i)) .ne. in2(i)) stop 30;
else
if (ichar(out8(i)) .ne. 0) stop 31;
endif
end do

! Now pack 5 values into the 5 character array out5.
print *, 'Testing g2_sbytesc() packing 5 values...'
in5(1) = 1
in5(2) = 2
in5(3) = 3
in5(4) = 4
in5(5) = 5
nbits = 8
nskip = 0
do i = 1, 5
out5(i) = char(0)
end do
call sbytesc(out5, in5, iskip, nbits, nskip, n5)
do i = 1, 5
if (ichar(out5(i)) .ne. in5(i)) stop 40;
end do

! Now pack 5 values into the 10 character array out10. Skip every
! other byte in the output.
print *, 'Testing g2_sbytesc() packing 5 values, skipping every other byte...'
nbits = 8
nskip = 0
do i = 1, 10
out10(i) = char(0)
end do
call sbytesc(out10, in5, iskip, nbits, 8, 5)
do i = 1, 10
! print '(z2.2)', out10(i)
if (mod(i, 2) .gt. 0) then
if (ichar(out10(i)) .ne. in5(int(i/2) + 1)) stop 51;
else
if (ichar(out10(i)) .ne. 0) stop 50;
endif
end do

print *, 'Testing sbytec() with iskip of 1...'
in(1) = 1
out(1) = char(0)
call sbytec(out, in, 1, 6)
! print '(z2.2)', out(1)
if (ichar(out(1)) .ne. 2) stop 20

print *, 'Testing g2_sbytesc() with a small array...'
iskip = 0
nbits = 32
nskip = 0
num = 1
in(1) = 1
call sbytesc(out4, in, iskip, nbits, nskip, num)
if (ichar(out4(1)) .ne. 0 .and. ichar(out4(2)) .ne. 0 .and. ichar(out4(3)) .ne. 0 .and. ichar(out4(4)) .ne. 1) stop 50
!print '(z2.2)', out4(1)

! For this test to pass the -fallow-argument-mismatch flag must be
! used, because I am passing in a real array instead of an int array
! for the in parameter. This is how g2_sbytesc() is called in
! addfield.F90.
print *, 'Testing sbytesc() with a real array (size 1) instead of an int array...'
iskip = 0
nbits = 32
nskip = 0
num = 1
r_in(1) = 1
call sbytesc(out4, r_in, iskip, nbits, nskip, num)
! Note that the 32-bit IEEE representation of 1.0 is 3f800000. The
! decimal for 3f is 63, the decimal for 80 is 128.
if (ichar(out4(1)) .ne. 63 .and. ichar(out4(2)) .ne. 128 .and. ichar(out4(3)) .ne. 0 .and. ichar(out4(4)) .ne. 0) stop 50
! print '(z2.2)', out4(1)
! print '(z2.2)', out4(2)
! print '(z2.2)', out4(3)
! print '(z2.2)', out4(4)

! This test is the same as above, but does not require the -fallow-argument-mismatch flag.
print *, 'Testing g2_sbytesc() with a real array (size 1) instead of an int array, using transfer() intrinsic...'
iskip = 0
nbits = 32
nskip = 0
num = 1
r_in(1) = 1
in = transfer(r_in, in)
call sbytesc(out4, in, iskip, nbits, nskip, num)
! Note that the 32-bit IEEE representation of 1.0 is 3f800000. The
! decimal for 3f is 63, the decimal for 80 is 128.
if (ichar(out4(1)) .ne. 63 .and. ichar(out4(2)) .ne. 128 .and. ichar(out4(3)) .ne. 0 .and. ichar(out4(4)) .ne. 0) stop 50
! print '(z2.2)', out4(1)

! For this test to pass the -fallow-argument-mismatch flag must be
! used, because I am passing in a real array instead of an int array
! for the in parameter. This is how g2_sbytesc() is called in
! addfield.F90.
print *, 'Testing sbytesc() with a real array instead of an int array...'
iskip = 0
nbits = 32
nskip = 0
num = 2
r_in2(1) = 1
r_in2(2) = 1
call sbytesc(out8, r_in2, iskip, nbits, nskip, num)
! Note that the 32-bit IEEE representation of 1.0 is 3f800000. The
! decimal for 3f is 63, the decimal for 80 is 128.
if (ichar(out8(1)) .ne. 63 .and. ichar(out8(2)) .ne. 128 .and. ichar(out8(3)) .ne. 0 .and. ichar(out8(4)) .ne. 0) stop 50
if (ichar(out8(5)) .ne. 63 .and. ichar(out8(6)) .ne. 128 .and. ichar(out8(7)) .ne. 0 .and. ichar(out8(8)) .ne. 0) stop 50
! print '(z2.2)', out8(1)

! This test is the same as above, but does not require the -fallow-argument-mismatch flag.
print *, 'Testing sbytesc() with a real array instead of an int array, using transfer() intrinsic...'
iskip = 0
nbits = 32
nskip = 0
num = 2
r_in2(1) = 1
r_in2(2) = 1
in = transfer(r_in2, in2)
call sbytesc(out8, in2, iskip, nbits, nskip, num)
! Note that the 32-bit IEEE representation of 1.0 is 3f800000. The
! decimal for 3f is 63, the decimal for 80 is 128.
if (ichar(out4(1)) .ne. 63 .and. ichar(out4(2)) .ne. 128 .and. ichar(out4(3)) .ne. 0 .and. ichar(out4(4)) .ne. 0) stop 50
if (ichar(out8(5)) .ne. 63 .and. ichar(out8(6)) .ne. 128 .and. ichar(out8(7)) .ne. 0 .and. ichar(out8(8)) .ne. 0) stop 50
! print '(z2.2)', out4(1)

print *, 'SUCCESS!'

end program test_gbytec
55 changes: 55 additions & 0 deletions tests/test_w3fi72.F90
Original file line number Diff line number Diff line change
@@ -0,0 +1,55 @@
! This is a test in the NCEPLIBS-w3emc project.
!
! Test the w3fi72() function.
!
! Ed Hartnett, 2/28/23
program test_w3fi72
implicit none
integer :: i, iret
integer :: kf, nbit
parameter(kf = 489900)
parameter(nbit = 8)
integer maxbit
parameter(maxbit=24)
character pds(40000)
real, dimension(:), allocatable :: fld
character*1, dimension(:), allocatable :: grib
integer igflag, igrid
parameter (igflag = 0)
integer ibm(kf),ipds(200),igds(200),ibds(200)
integer kfo, lgrib, icomp
integer ierr

print *, "Testing w3fi72..."

allocate(fld(kf))
allocate(grib(1000+kf*12))

! Fill up some test data.
do i = 1, 200
ipds(i) = 0
end do
do i = 1, kf
fld(i) = i / 2
end do

! Fill the igds array. This call comes from test_w3fi71.F90.
icomp = 0
igrid = 172
call w3fi71(igrid, igds, ierr)
if (ierr .ne. 0) stop 1

! This call comes from NCEPLIBS-grib_util cnvgrb, putbexn.F90.
call w3fi72(0, fld, 0, nbit, 1, ipds, pds, &
igflag, igrid, igds, icomp, 0, ibm, kf, ibds, &
kfo, grib, lgrib, iret)
print *, kfo, lgrib, iret
if (kfo .ne. 489900) stop 2
! if (lgrib .ne. 489956) stop 3
if (iret .ne. 0) stop 4

deallocate(fld)
deallocate(grib)

print *, "SUCCESS"
end program test_w3fi72
39 changes: 39 additions & 0 deletions tests/test_w3fi73.F90
Original file line number Diff line number Diff line change
@@ -0,0 +1,39 @@
! This is a test in the NCEPLIBS-w3emc project.
!
! Test the w3fi73() subroutine.
!
! Ed Hartnett, 9/27/23
program test_w3fi73
implicit none
integer :: ibflag, iblen, lenbms
integer :: BLEN, BMSLEN
parameter (BLEN = 3)
parameter (BMSLEN = 3)
integer :: ibmap(BLEN), bms(BMSLEN)
integer :: i
integer :: ierr

print*,"Testing w3fi73..."

! Return error code 8 if all ibmap values are 0.
ibflag = 0
do i = 1, BLEN
ibmap(i) = 0
end do
iblen = BLEN
call w3fi73(ibflag, ibmap, iblen, bms, lenbms, ierr)
if (ierr .ne. 8) stop 2

! Return error code 8 if all ibmap values are 0.
ibflag = 0
do i = 1, BLEN
ibmap(i) = 1
end do
iblen = BLEN
call w3fi73(ibflag, ibmap, iblen, bms, lenbms, ierr)
if (ierr .ne. 0) stop 4
if (lenbms .ne. 8) stop 5
! if (bms(1) .ne. 218628096 .or. bms(2) .ne. 14680064) stop 7

print*,"SUCCESS"
end program test_w3fi73
1 change: 0 additions & 1 deletion tests/test_w3fi74.F90
Original file line number Diff line number Diff line change
Expand Up @@ -31,7 +31,6 @@ program test_w3fi74

! Fill the igds array. This call comes from w3if72.f.
icomp = 0
npts = 4
call w3fi74(igds, icomp, gds, lengds, npts, ierr)
if (ierr .ne. 0) stop 1
if (lengds .ne. 32 .or. npts .ne. 489900) stop 2
Expand Down