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

Add gbytesc8() and sbytesc8() to handle 64-bit ints, so that we can use a 64-bit int in the index files. #605

Merged
merged 18 commits into from
Feb 7, 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
254 changes: 220 additions & 34 deletions src/g2_gbytesc.F90
Original file line number Diff line number Diff line change
Expand Up @@ -3,17 +3,17 @@
!> string and unpacked array.
!> @author Stephen Gilbert @date 2004-04-27

!> Extract arbitrary size values from a packed bit string, right
!> justifying each value in the unpacked array without skip and
!> interations.
!> Extract one arbitrary size value (up to 32 bits) from a packed bit
!> string, right justifying each value in the unpacked array.
!>
!> This should be used when input array IN has only one element. If in
!> has more elements, use g2_sbytesc().
!>
!> @param[in] in Array input.
!> @param[inout] iout Unpacked array output.
!> @param[in] iskip Initial number of bits to skip.
!> @param[in] nbits Number of bits of each integer in IN to take.
!> @param[in] nbits Number of bits of each integer in IN to take. Must
!> be 32 or less.
!>
!> @author Stephen Gilbert @date 2004-04-27
subroutine g2_gbytec(in, iout, iskip, nbits)
Expand All @@ -25,51 +25,110 @@ subroutine g2_gbytec(in, iout, iskip, nbits)
call g2_gbytesc(in, iout, iskip, nbits, 0, 1)
end subroutine g2_gbytec

!> Put arbitrary size values into a packed bit string, taking the low
!> order bits from each value in the unpacked array without skip and
!> interation.
!> Extract arbitrary size values (up to 32 bits each) from a packed
!> bit string, right justifying each value in the unpacked array.
!>
!> This should be used when input array IN has only one element. If IN
!> has more elements, use G2_SBYTESC().
!>
!> @param[inout] out packed array output
!> @param[in] in unpacked array input
!> @param[in] in array input
!> @param[out] iout unpacked array output
!> @param[in] iskip initial number of bits to skip
!> @param[in] nbits Number of bits of each integer in OUT to fill.
!> @param[in] nbits Number of bits of each integer in IN to take. Must
!> be 32 or less.
!> @param[in] nskip Additional number of bits to skip on each iteration.
!> @param[in] n Number of integers to extract.
!>
!> @author Stephen Gilbert @date 2004-04-27
subroutine g2_sbytec(out, in, iskip, nbits)
subroutine g2_gbytesc(in, iout, iskip, nbits, nskip, n)
implicit none

character*1, intent(inout) :: out(*)
integer, intent(in) :: in(*)
character*1, intent(in) :: in(*)
integer, intent(out) :: iout(*)
integer, intent(in) :: iskip, nbits, nskip, n
integer :: tbit, bitcnt
integer, parameter :: ones(8) = (/ 1, 3, 7, 15, 31, 63, 127, 255 /)

integer :: nbit, i, index, ibit, itmp
integer, external :: mova2i

! nbit is the start position of the field in bits
nbit = iskip
do i = 1, n
bitcnt = nbits
index = nbit / 8 + 1
ibit = mod(nbit, 8)
nbit = nbit + nbits + nskip

! first byte
tbit = min(bitcnt, 8 - ibit)
itmp = iand(mova2i(in(index)), ones(8 - ibit))
if (tbit .ne. 8 - ibit) itmp = ishft(itmp, tbit - 8 + ibit)
index = index + 1
bitcnt = bitcnt - tbit

! now transfer whole bytes
do while (bitcnt .ge. 8)
itmp = ior(ishft(itmp,8), mova2i(in(index)))
bitcnt = bitcnt - 8
index = index + 1
enddo

! get data from last byte
if (bitcnt .gt. 0) then
itmp = ior(ishft(itmp, bitcnt), iand(ishft(mova2i(in(index)), &
- (8 - bitcnt)), ones(bitcnt)))
endif

iout(i) = itmp
enddo

end subroutine g2_gbytesc

!> Extract one arbitrary sized (up to 64-bits) values from a packed bit
!> string, right justifying each value in the unpacked array.
!>
!> This should be used when input array in has only one element. If in
!> has more elements, use g2_sbytesc().
!>
!> @param[in] in Array input.
!> @param[inout] iout Unpacked array output.
!> @param[in] iskip Initial number of bits to skip.
!> @param[in] nbits Number of bits of each integer in IN to take. Must
!> be 64 or less.
!>
!> @author Stephen Gilbert @date 2004-04-27
subroutine g2_gbytec8(in, iout, iskip, nbits)
implicit none

character*1, intent(in) :: in(*)
integer (kind = 8), intent(inout) :: iout(*)
integer, intent(in) :: iskip, nbits
call g2_sbytesc(out, in, iskip, nbits, 0, 1)
end subroutine g2_sbytec
call g2_gbytesc8(in, iout, iskip, nbits, 0, 1)
end subroutine g2_gbytec8

!> Extract arbitrary size values from a packed bit string, right
!> justifying each value in the unpacked array with skip and
!> interation options.
!> Extract arbitrary sized (up to 64-bits) values from a packed bit
!> string, right justifying each value in the unpacked array.
!>
!> @param[in] in array input
!> @param[out] iout unpacked array output
!> @param[in] iskip initial number of bits to skip
!> @param[in] nbits Number of bits of each integer in IN to take.
!> @param[in] nbits Number of bits of each integer in IN to take. Must
!> be 64 or less.
!> @param[in] nskip Additional number of bits to skip on each iteration.
!> @param[in] n Number of integers to extract.
!>
!> @author Stephen Gilbert @date 2004-04-27
subroutine g2_gbytesc(in, iout, iskip, nbits, nskip, n)
subroutine g2_gbytesc8(in, iout, iskip, nbits, nskip, n)
implicit none

character*1, intent(in) :: in(*)
integer, intent(out) :: iout(*)
integer (kind = 8), intent(out) :: iout(*)
integer, intent(in) :: iskip, nbits, nskip, n
integer :: tbit, bitcnt
integer, parameter :: ones(8) = (/ 1, 3, 7, 15, 31, 63, 127, 255 /)

integer :: nbit, i, index, ibit, itmp
integer (kind = 8) :: itmp8, itmp8_2, itmp8_3
integer, external :: mova2i
integer (kind = 8), external :: mova2i8

! nbit is the start position of the field in bits
nbit = iskip
Expand All @@ -81,37 +140,65 @@ subroutine g2_gbytesc(in, iout, iskip, nbits, nskip, n)

! first byte
tbit = min(bitcnt, 8 - ibit)
itmp8 = iand(mova2i8(in(index)), int(ones(8 - ibit), kind = 8))
itmp = iand(mova2i(in(index)), ones(8 - ibit))
if (tbit .ne. 8 - ibit) itmp = ishft(itmp, tbit - 8 + ibit)
if (tbit .ne. 8 - ibit) itmp8 = ishft(itmp8, tbit - 8 + ibit)
index = index + 1
bitcnt = bitcnt - tbit

! now transfer whole bytes
do while (bitcnt .ge. 8)
itmp = ior(ishft(itmp,8), mova2i(in(index)))
itmp8 = ior(ishft(itmp8,8), mova2i8(in(index)))
bitcnt = bitcnt - 8
index = index + 1
enddo

! get data from last byte
if (bitcnt .gt. 0) then
itmp = ior(ishft(itmp, bitcnt), iand(ishft(mova2i(in(index)), &
- (8 - bitcnt)), ones(bitcnt)))
itmp = ior(ishft(itmp, bitcnt), iand(ishft(mova2i(in(index)), - (8 - bitcnt)), ones(bitcnt)))
itmp8_2 = ishft(mova2i8(in(index)), int(-(8 - bitcnt), kind(8)))
itmp8_3 = int(ones(bitcnt), kind(8))
itmp8 = ior(ishft(itmp8, bitcnt), iand(itmp8_2, itmp8_3))
endif

iout(i) = itmp
iout(i) = itmp8
enddo

end subroutine g2_gbytesc
end subroutine g2_gbytesc8

!> Put one arbitrary sized (up to 32 bits) values into a packed bit
!> string, taking the low order bits from the value in the unpacked
!> array.
!>
!> This should be used when input array IN has only one element. If IN
!> has more elements, use g2_sbytesc().
!>
!> @param[inout] out packed array output
!> @param[in] in unpacked array input
!> @param[in] iskip initial number of bits to skip
!> @param[in] nbits Number of bits of each integer in OUT to fill.
!>
!> @author Stephen Gilbert @date 2004-04-27
subroutine g2_sbytec(out, in, iskip, nbits)
implicit none

character*1, intent(inout) :: out(*)
integer, intent(in) :: in(*)
integer, intent(in) :: iskip, nbits
call g2_sbytesc(out, in, iskip, nbits, 0, 1)
end subroutine g2_sbytec

!> Put arbitrary size values into a packed bit string, taking the low
!> order bits from each value in the unpacked array with skip and
!> interation options.
!> Put arbitrary size (up to 32 bits each) values into a packed bit
!> string, taking the low order bits from each value in the unpacked
!> array.
!>
!> @param[out] out Packed array output.
!> @param[in] in Unpacked array input.
!> @param[in] iskip Initial number of bits to skip.
!> @param[in] nbits Number of bits of each integer in OUT to fill.
!> @param[in] nbits Number of bits of each integer in OUT to
!> fill. Must be 32 or less.
!> @param[in] nskip Additional number of bits to skip on each iteration.
!> @param[in] n Number of iterations.
!>
Expand All @@ -131,13 +218,15 @@ subroutine g2_sbytesc(out, in, iskip, nbits, nskip, n)
! number bits from zero to ...
! nbit is the last bit of the field to be filled
nbit = iskip + nbits - 1
!print *, 'nbit', nbit, 'nbits ', nbits, 'nskip', nskip, 'n', n
do i = 1, n
itmp = in(i)
bitcnt = nbits
index = nbit / 8 + 1
ibit = mod(nbit, 8)
nbit = nbit + nbits + nskip

!print *, 'i', i, 'itmp', itmp, 'bitcnt', bitcnt, 'index', index, 'ibit', ibit, 'nbit', nbit

! make byte aligned
if (ibit .ne. 7) then
tbit = min(bitcnt, ibit + 1)
Expand All @@ -155,17 +244,114 @@ subroutine g2_sbytesc(out, in, iskip, nbits, nskip, n)
! do by bytes
do while (bitcnt .ge. 8)
out(index) = char(iand(itmp, 255))
!print '(z2.2, x, z2.2, x, z2.2)', out(index), itmp, iand(itmp, 255)
itmp = ishft(itmp, -8)
bitcnt = bitcnt - 8
index = index - 1
enddo

! do last byte
! Do left over bits.
if (bitcnt .gt. 0) then
itmp2 = iand(itmp, ones(bitcnt))
!print '(z2.2, x, z2.2)', ones(bitcnt), itmp2
itmp3 = iand(mova2i(out(index)), 255 - ones(bitcnt))
out(index) = char(ior(itmp2, itmp3))
endif
enddo

end subroutine g2_sbytesc

!> Put one arbitrary sized (up to 64 bits) values into a packed bit
!> string, taking the low order bits from each value in the unpacked
!> array.
!>
!> This should be used when input array IN has only one element. If IN
!> has more elements, use g2_sbytesc().
!>
!> @param[inout] out packed array output
!> @param[in] in unpacked array input
!> @param[in] iskip initial number of bits to skip
!> @param[in] nbits Number of bits of each integer in OUT to
!> fill. Must be 64 or less.
!>
!> @author Stephen Gilbert @date 2004-04-27
subroutine g2_sbytec8(out, in, iskip, nbits)
implicit none

character*1, intent(inout) :: out(*)
integer (kind = 8), intent(in) :: in(*)
integer, intent(in) :: iskip, nbits
call g2_sbytesc8(out, in, iskip, nbits, 0, 1)
end subroutine g2_sbytec8

!> Put arbitrary sized (up to 64 bits each) values into a packed bit
!> string, taking the low order bits from each value in the unpacked
!> array.
!>
!> @param[out] out Packed array output.
!> @param[in] in Unpacked array input.
!> @param[in] iskip Initial number of bits to skip.
!> @param[in] nbits Number of bits of each integer in OUT to
!> fill. Must be 64 or less.
!> @param[in] nskip Additional number of bits to skip on each iteration.
!> @param[in] n Number of iterations.
!>
!> @author Stephen Gilbert @date 2004-04-27
subroutine g2_sbytesc8(out, in, iskip, nbits, nskip, n)
implicit none

character*1, intent(out) :: out(*)
integer (kind = 8), intent(in) :: in(n)
integer, intent(in) :: iskip, nbits, nskip, n

integer :: bitcnt, tbit
integer, parameter :: ones(8)=(/ 1, 3, 7, 15, 31, 63, 127, 255/)
integer :: nbit, i, index, ibit, imask, itmp1, itmp2, itmp3
integer (kind = 8) :: itmp8, itmp8_2
integer, external :: mova2i

! number bits from zero to ...
! nbit is the last bit of the field to be filled
nbit = iskip + nbits - 1
print *, 'nbit', nbit, 'nbits ', nbits, 'nskip', nskip, 'n', n
do i = 1, n
itmp8 = in(i)
bitcnt = nbits
index = nbit / 8 + 1
ibit = mod(nbit, 8)
nbit = nbit + nbits + nskip
print *, 'i', i, 'itmp8', itmp8, 'bitcnt', bitcnt, 'index', index, 'ibit', ibit, 'nbit', nbit

! make byte aligned
if (ibit .ne. 7) then
tbit = min(bitcnt, ibit + 1)
imask = ishft(ones(tbit), 7 - ibit)
itmp1 = int(ishft(itmp8, int(7 - ibit, kind(8))), kind(4))
itmp2 = iand(itmp1, imask)
itmp3 = iand(mova2i(out(index)), 255 - imask)
out(index) = char(ior(itmp2, itmp3))
bitcnt = bitcnt - tbit
itmp8 = ishft(itmp8, -tbit)
index = index - 1
endif

! now byte aligned

! Process a byte at a time.
do while (bitcnt .ge. 8)
!print *, bitcnt, iand(itmp8, 255_8)
out(index) = char(iand(itmp8, 255_8))
itmp8 = ishft(itmp8, -8)
bitcnt = bitcnt - 8
index = index - 1
enddo

! Take care of left over bits.
if (bitcnt .gt. 0) then
itmp8_2 = int(ones(bitcnt), kind(8))
itmp2 = int(iand(itmp8, itmp8_2), kind(4))
itmp3 = iand(mova2i(out(index)), 255 - ones(bitcnt))
out(index) = char(ior(itmp2, itmp3))
endif
enddo
end subroutine g2_sbytesc8
Loading
Loading