amun-code/sources/hash.F90
Grzegorz Kowal 3dbbd0b8b1 HASH: Change kind of local variables.
Local variables remain and offset cannot be larger than the size of
input array, so they can be 4-byte integers.

Signed-off-by: Grzegorz Kowal <grzegorz@amuncode.org>
2020-05-13 08:32:27 -03:00

563 lines
15 KiB
Fortran

!!******************************************************************************
!!
!! This file is part of the AMUN source code, a program to perform
!! Newtonian or relativistic magnetohydrodynamical simulations on uniform or
!! adaptive mesh.
!!
!! Copyright (C) 2012-2020 Yann Collet
!! Copyright (C) 2020 Grzegorz Kowal <grzegorz@amuncode.org>
!!
!! This program is free software: you can redistribute it and/or modify
!! it under the terms of the GNU General Public License as published by
!! the Free Software Foundation, either version 3 of the License, or
!! (at your option) any later version.
!!
!! This program is distributed in the hope that it will be useful,
!! but WITHOUT ANY WARRANTY; without even the implied warranty of
!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
!! GNU General Public License for more details.
!!
!! You should have received a copy of the GNU General Public License
!! along with this program. If not, see <http://www.gnu.org/licenses/>.
!!
!!******************************************************************************
!!
!! module: HASH
!!
!! This module provides 64-bit version of the xxHash64 by Yann Collet.
!! This is a Fortran implementation based on the XXH64 specification
!! published at
!! https://github.com/Cyan4973/xxHash/blob/dev/doc/xxhash_spec.md
!!
!! For additional info, see
!! http://www.xxhash.com or https://github.com/Cyan4973/xxHash
!!
!!******************************************************************************
!
module hash
! module variables are not implicit by default
!
implicit none
! hash parameters
!
integer(kind=8), parameter :: seed = 0_8
integer(kind=8), parameter :: prime1 = -7046029288634856825_8, &
prime2 = -4417276706812531889_8, &
prime3 = 1609587929392839161_8, &
prime4 = -8796714831421723037_8, &
prime5 = 2870177450012600261_8
! by default everything is private
!
private
! declare public subroutines
!
public :: xxh64_integer, xxh64_long, xxh64_double, xxh64_complex
!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
!
contains
!
!===============================================================================
!!
!!*** PRIVATE SUBROUTINES ****************************************************
!!
!===============================================================================
!
!
!===============================================================================
!
! function XXH64_INTEGER:
! ----------------------
!
! Function calculates XXH64 hash for a given integer vector.
!
! Arguments:
!
! n - the size of the input vector;
! input - the input vactor of integer values;
!
!===============================================================================
!
integer(kind=8) function xxh64_integer(n, input) result(hash)
implicit none
! subroutine arguments
!
integer(kind=4) , intent(in) :: n
integer(kind=4), dimension(n), intent(in) :: input
! local variables
!
integer(kind=4) :: remain, offset
! local arrays
!
integer(kind=8), dimension(4) :: lane, chk
!-------------------------------------------------------------------------------
!
hash = 0_8
offset = 1
remain = n
if (remain >= 8) then
lane(1) = seed + prime1 + prime2
lane(2) = seed + prime2
lane(3) = seed + 0_8
lane(4) = seed - prime1
do while (remain >= 8)
chk(1:4) = transfer(input(offset:offset+7), chk)
lane(1) = xxh64_round(lane(1), chk(1))
lane(2) = xxh64_round(lane(2), chk(2))
lane(3) = xxh64_round(lane(3), chk(3))
lane(4) = xxh64_round(lane(4), chk(4))
offset = offset + 8
remain = remain - 8
end do
hash = xxh64_rotl(lane(1), 1) + xxh64_rotl(lane(2), 7) + &
xxh64_rotl(lane(3), 12) + xxh64_rotl(lane(4), 18)
hash = xxh64_merge(hash, lane(1))
hash = xxh64_merge(hash, lane(2))
hash = xxh64_merge(hash, lane(3))
hash = xxh64_merge(hash, lane(4))
else
hash = seed + prime5
end if
hash = hash + 4_8 * n
do while (remain >= 2)
chk(1) = transfer(input(offset:offset+1), chk(1))
hash = ieor(hash, xxh64_round(0_8, chk(1)))
hash = xxh64_rotl(hash, 27)
hash = hash * prime1 + prime4
offset = offset + 2
remain = remain - 2
end do
if (remain == 1) then
chk(1) = transfer((/ input(offset), 0 /), chk(1))
hash = ieor(hash, chk(1) * prime1)
hash = xxh64_rotl(hash, 23)
hash = hash * prime2 + prime3
offset = offset + 1
remain = remain - 1
end if
hash = xxh64_aval(hash)
return
!-------------------------------------------------------------------------------
!
end function xxh64_integer
!
!===============================================================================
!
! function XXH64_LONG:
! -------------------
!
! Function calculates XXH64 hash for a given long integer vector.
!
! Arguments:
!
! n - the size of the input vector;
! input - the input vactor of real values;
!
!===============================================================================
!
integer(kind=8) function xxh64_long(n, input) result(hash)
implicit none
! subroutine arguments
!
integer(kind=4) , intent(in) :: n
integer(kind=8), dimension(n), intent(in) :: input
! local variables
!
integer(kind=4) :: remain, offset
! local arrays
!
integer(kind=8), dimension(4) :: lane, chk
!-------------------------------------------------------------------------------
!
hash = 0_8
offset = 1
remain = n
if (remain >= 4) then
lane(1) = seed + prime1 + prime2
lane(2) = seed + prime2
lane(3) = seed + 0_8
lane(4) = seed - prime1
do while (remain >= 4)
chk(1:4) = transfer(input(offset:offset+3), chk)
lane(1) = xxh64_round(lane(1), chk(1))
lane(2) = xxh64_round(lane(2), chk(2))
lane(3) = xxh64_round(lane(3), chk(3))
lane(4) = xxh64_round(lane(4), chk(4))
offset = offset + 4
remain = remain - 4
end do
hash = xxh64_rotl(lane(1), 1) + xxh64_rotl(lane(2), 7) + &
xxh64_rotl(lane(3), 12) + xxh64_rotl(lane(4), 18)
hash = xxh64_merge(hash, lane(1))
hash = xxh64_merge(hash, lane(2))
hash = xxh64_merge(hash, lane(3))
hash = xxh64_merge(hash, lane(4))
else
hash = seed + prime5
end if
hash = hash + 8_8 * n
do while (remain >= 1)
hash = ieor(hash, xxh64_round(0_8, transfer(input(offset), 0_8)))
hash = xxh64_rotl(hash, 27)
hash = hash * prime1 + prime4
offset = offset + 1
remain = remain - 1
end do
hash = xxh64_aval(hash)
return
!-------------------------------------------------------------------------------
!
end function xxh64_long
!
!===============================================================================
!
! function XXH64_DOUBLE:
! ---------------------
!
! Function calculates XXH64 hash for a given double precision vector.
!
! Arguments:
!
! n - the size of the input vector;
! input - the input vactor of real values;
!
!===============================================================================
!
integer(kind=8) function xxh64_double(n, input) result(hash)
implicit none
! subroutine arguments
!
integer(kind=4) , intent(in) :: n
real(kind=8), dimension(n), intent(in) :: input
! local variables
!
integer(kind=4) :: remain, offset
! local arrays
!
integer(kind=8), dimension(4) :: lane, chk
!-------------------------------------------------------------------------------
!
hash = 0_8
offset = 1
remain = n
if (remain >= 4) then
lane(1) = seed + prime1 + prime2
lane(2) = seed + prime2
lane(3) = seed + 0_8
lane(4) = seed - prime1
do while (remain >= 4)
chk(1:4) = transfer(input(offset:offset+3), chk)
lane(1) = xxh64_round(lane(1), chk(1))
lane(2) = xxh64_round(lane(2), chk(2))
lane(3) = xxh64_round(lane(3), chk(3))
lane(4) = xxh64_round(lane(4), chk(4))
offset = offset + 4
remain = remain - 4
end do
hash = xxh64_rotl(lane(1), 1) + xxh64_rotl(lane(2), 7) + &
xxh64_rotl(lane(3), 12) + xxh64_rotl(lane(4), 18)
hash = xxh64_merge(hash, lane(1))
hash = xxh64_merge(hash, lane(2))
hash = xxh64_merge(hash, lane(3))
hash = xxh64_merge(hash, lane(4))
else
hash = seed + prime5
end if
hash = hash + 8_8 * n
do while (remain >= 1)
hash = ieor(hash, xxh64_round(0_8, transfer(input(offset), 0_8)))
hash = xxh64_rotl(hash, 27)
hash = hash * prime1 + prime4
offset = offset + 1
remain = remain - 1
end do
hash = xxh64_aval(hash)
return
!-------------------------------------------------------------------------------
!
end function xxh64_double
!
!===============================================================================
!
! function XXH64_COMPLEX:
! ----------------------
!
! Function calculates XXH64 hash for a given double precision complex vector.
!
! Arguments:
!
! n - the size of the input vector;
! input - the input vactor of real values;
!
!===============================================================================
!
integer(kind=8) function xxh64_complex(n, input) result(hash)
implicit none
! subroutine arguments
!
integer(kind=4) , intent(in) :: n
complex(kind=8), dimension(n), intent(in) :: input
! local variables
!
integer(kind=4) :: remain, offset
! local arrays
!
integer(kind=8), dimension(4) :: lane, chk
!-------------------------------------------------------------------------------
!
hash = 0_8
offset = 1
remain = n
if (remain >= 2) then
lane(1) = seed + prime1 + prime2
lane(2) = seed + prime2
lane(3) = seed + 0_8
lane(4) = seed - prime1
do while (remain >= 2)
chk(1:4) = transfer(input(offset:offset+1), chk)
lane(1) = xxh64_round(lane(1), chk(1))
lane(2) = xxh64_round(lane(2), chk(2))
lane(3) = xxh64_round(lane(3), chk(3))
lane(4) = xxh64_round(lane(4), chk(4))
offset = offset + 2
remain = remain - 2
end do
hash = xxh64_rotl(lane(1), 1) + xxh64_rotl(lane(2), 7) + &
xxh64_rotl(lane(3), 12) + xxh64_rotl(lane(4), 18)
hash = xxh64_merge(hash, lane(1))
hash = xxh64_merge(hash, lane(2))
hash = xxh64_merge(hash, lane(3))
hash = xxh64_merge(hash, lane(4))
else
hash = seed + prime5
end if
hash = hash + 16_8 * n
if (remain == 1) then
hash = ieor(hash, xxh64_round(0_8, transfer(dreal(input(offset)), 0_8)))
hash = xxh64_rotl(hash, 27)
hash = hash * prime1 + prime4
hash = ieor(hash, xxh64_round(0_8, transfer(dimag(input(offset)), 0_8)))
hash = xxh64_rotl(hash, 27)
hash = hash * prime1 + prime4
offset = offset + 1
remain = remain - 1
end if
hash = xxh64_aval(hash)
return
!-------------------------------------------------------------------------------
!
end function xxh64_complex
!
!===============================================================================
!
! function XXH64_ROUND:
! --------------------
!
! Function processes one stripe of the input data updating
! the correponding lane.
!
! Arguments:
!
! lane - the lane;
! input - the 8-byte data to process;
!
!===============================================================================
!
integer(kind=8) function xxh64_round(lane, input)
implicit none
! subroutine arguments
!
integer(kind=8), intent(in) :: lane, input
!-------------------------------------------------------------------------------
!
xxh64_round = lane + (input * prime2)
xxh64_round = xxh64_rotl(xxh64_round, 31)
xxh64_round = xxh64_round * prime1
return
!-------------------------------------------------------------------------------
!
end function xxh64_round
!
!===============================================================================
!
! function XXH64_MERGE:
! --------------------
!
! Function performs merging of the given lane in to the hash.
!
! Arguments:
!
! hash - the hash to merge to;
! lane - the lane being merged;
!
!===============================================================================
!
integer(kind=8) function xxh64_merge(hash, lane)
implicit none
! subroutine arguments
!
integer(kind=8), intent(in) :: hash, lane
!-------------------------------------------------------------------------------
!
xxh64_merge = ieor(hash, xxh64_round(0_8, lane))
xxh64_merge = xxh64_merge * prime1 + prime4
return
!-------------------------------------------------------------------------------
!
end function xxh64_merge
!
!===============================================================================
!
! function XXH64_AVAL:
! -------------------
!
! Function calculates the final mix of the hash.
!
! Arguments:
!
! hash - the hash to mix;
!
!===============================================================================
!
integer(kind=8) function xxh64_aval(hash)
implicit none
! subroutine arguments
!
integer(kind=8), intent(in) :: hash
!-------------------------------------------------------------------------------
!
xxh64_aval = hash
xxh64_aval = ieor(xxh64_aval, ishft(xxh64_aval, -33)) * prime2
xxh64_aval = ieor(xxh64_aval, ishft(xxh64_aval, -29)) * prime3
xxh64_aval = ieor(xxh64_aval, ishft(xxh64_aval, -32))
return
!-------------------------------------------------------------------------------
!
end function xxh64_aval
!
!===============================================================================
!
! function XXH64_ROTL:
! -------------------
!
! Function calculates the rotation of the input 8-byte word by a given amount.
!
! Arguments:
!
! byte - the byte to be rotates;
! amount - the amount by which rotate the input byte;
!
!===============================================================================
!
integer(kind=8) function xxh64_rotl(byte, amount)
implicit none
! subroutine arguments
!
integer(kind=8), intent(in) :: byte
integer(kind=4), intent(in) :: amount
!-------------------------------------------------------------------------------
!
xxh64_rotl = ior(ishft(byte, amount), ishft(byte, amount - 64))
return
!-------------------------------------------------------------------------------
!
end function xxh64_rotl
!===============================================================================
!
end module hash