Skip to content

Commit

Permalink
+Add the new interface field_checksum
Browse files Browse the repository at this point in the history
  Added the new publicly visible interface field_checksum to return the checksum
that is used to verify the arrays that are written to or read from files, while
handling both the unscaling and rotation of the arrays.  The 5 internal
rotated_field_chksum_real_... functions were renamed to field_checksum_real_...
to reflect their slightly broadened purpose and now have new optional unscale
arguments.  The previous rotated_field_chksum interface has been retained
because it is still being used in SIS2, but this may change in the future.  All
answers are bitwise identical, but there is a new public function interface and
a new optional argument to a preexisting one.
  • Loading branch information
Hallberg-NOAA committed Feb 24, 2025
1 parent 809d56e commit ad0a8b8
Showing 1 changed file with 105 additions and 34 deletions.
139 changes: 105 additions & 34 deletions src/framework/MOM_checksums.F90
Original file line number Diff line number Diff line change
Expand Up @@ -17,7 +17,7 @@ module MOM_checksums

implicit none ; private

public :: chksum0, zchksum, rotated_field_chksum
public :: chksum0, zchksum, rotated_field_chksum, field_checksum
public :: hchksum, Bchksum, uchksum, vchksum, qchksum, is_NaN, chksum
public :: hchksum_pair, uvchksum, Bchksum_pair
public :: MOM_checksums_init
Expand Down Expand Up @@ -84,15 +84,29 @@ module MOM_checksums
module procedure is_NaN_0d, is_NaN_1d, is_NaN_2d, is_NaN_3d
end interface

!> Rotate and compute the checksum of a field
!> Compute the checksum on all elements of a field that may need to be rotated or unscaled.
!! This interface uses the field_chksum function that is used to verify file contents, which
!! may differ from the bitcount function used for other checksums in this module.
interface rotated_field_chksum
module procedure rotated_field_chksum_real_0d
module procedure rotated_field_chksum_real_1d
module procedure rotated_field_chksum_real_2d
module procedure rotated_field_chksum_real_3d
module procedure rotated_field_chksum_real_4d
module procedure field_checksum_real_0d
module procedure field_checksum_real_1d
module procedure field_checksum_real_2d
module procedure field_checksum_real_3d
module procedure field_checksum_real_4d
end interface rotated_field_chksum


!> Compute the checksum on all elements of a field that may need to be rotated or unscaled.
!! This interface uses the field_chksum function that is used to verify file contents, which
!! may differ from the bitcount function used for other checksums in this module.
interface field_checksum
module procedure field_checksum_real_0d
module procedure field_checksum_real_1d
module procedure field_checksum_real_2d
module procedure field_checksum_real_3d
module procedure field_checksum_real_4d
end interface field_checksum

integer, parameter :: bc_modulus = 1000000000 !< Modulus of checksum bitcount
integer, parameter :: default_shift=0 !< The default array shift
logical :: calculateStatistics=.true. !< If true, report min, max and mean.
Expand Down Expand Up @@ -2366,119 +2380,176 @@ function is_NaN_3d(x)

end function is_NaN_3d

! The following set of routines do a checksum across the computational domain of
! a field, with the potential for rotation of this field and masking.
! The following set of routines do a checksum across all elements of a field,
! with the potential for the unscaling and rotation of this field and masking.

!> Compute the field checksum of a scalar.
function rotated_field_chksum_real_0d(field, pelist, mask_val, turns) &
!> Compute the field checksum of a scalar that may need to be unscaled.
!! This uses the field_chksum function that is used to verify file contents, which may differ
!! from the bitcount function used for other checksums in this module.
function field_checksum_real_0d(field, pelist, mask_val, turns, unscale) &
result(chksum)
real, intent(in) :: field !< Input scalar [arbitrary]
real, intent(in) :: field !< Input scalar to be checksummed in arbitrary,
!! possibly rescaled units [A ~> a]
integer, optional, intent(in) :: pelist(:) !< PE list of ranks to checksum
real, optional, intent(in) :: mask_val !< FMS mask value [nondim]
integer, optional, intent(in) :: turns !< Number of quarter turns
real, optional, intent(in) :: unscale !< A factor to convert this array back to
!! unscaled units for checksums [a A-1 ~> 1]
integer(kind=int64) :: chksum !< checksum of scalar

real :: scale_fac ! A local copy of unscale if it is present [a A-1 ~> 1] or 1 otherwise

if (present(turns)) call MOM_error(FATAL, "Rotation not supported for 0d fields.")

chksum = field_chksum(field, pelist=pelist, mask_val=mask_val)
end function rotated_field_chksum_real_0d
scale_fac = 1.0 ; if (present(unscale)) scale_fac = unscale

chksum = field_chksum(scale_fac*field, pelist=pelist, mask_val=mask_val)
end function field_checksum_real_0d


!> Compute the field checksum of a 1d field.
function rotated_field_chksum_real_1d(field, pelist, mask_val, turns) &
!> Compute the field checksum of an entire 1d array that may need to be unscaled.
!! This uses the field_chksum function that is used to verify file contents, which may differ
!! from the bitcount function used for other checksums in this module.
function field_checksum_real_1d(field, pelist, mask_val, turns, unscale) &
result(chksum)
real, dimension(:), intent(in) :: field !< Input array [arbitrary]
real, dimension(:), intent(in) :: field !< Input array to be checksummed in arbitrary,
!! possibly rescaled units [A ~> a]
integer, optional, intent(in) :: pelist(:) !< PE list of ranks to checksum
real, optional, intent(in) :: mask_val !< FMS mask value [nondim]
integer, optional, intent(in) :: turns !< Number of quarter turns
real, optional, intent(in) :: unscale !< A factor to convert this array back to
!! unscaled units for checksums [a A-1 ~> 1]
integer(kind=int64) :: chksum !< checksum of array

real :: scale_fac ! A local copy of unscale if it is present [a A-1 ~> 1] or 1 otherwise

if (present(turns)) call MOM_error(FATAL, "Rotation not supported for 1d fields.")

chksum = field_chksum(field, pelist=pelist, mask_val=mask_val)
end function rotated_field_chksum_real_1d
scale_fac = 1.0 ; if (present(unscale)) scale_fac = unscale

chksum = field_chksum(scale_fac*field(:), pelist=pelist, mask_val=mask_val)
end function field_checksum_real_1d


!> Compute the field checksum of a rotated 2d field.
function rotated_field_chksum_real_2d(field, pelist, mask_val, turns) &
!> Compute the field checksum of an entire 2d array that may need to be rotated or unscaled.
!! This uses the field_chksum function that is used to verify file contents, which may differ
!! from the bitcount function used for other checksums in this module.
function field_checksum_real_2d(field, pelist, mask_val, turns, unscale) &
result(chksum)
real, dimension(:,:), intent(in) :: field !< Unrotated input field [arbitrary]
real, dimension(:,:), intent(in) :: field !< Unrotated input field to be checksummed in
!! arbitrary, possibly rescaled units [A ~> a]
integer, optional, intent(in) :: pelist(:) !< PE list of ranks to checksum
real, optional, intent(in) :: mask_val !< FMS mask value [nondim]
integer, optional, intent(in) :: turns !< Number of quarter turns
real, optional, intent(in) :: unscale !< A factor to convert this array back to
!! unscaled units for checksums [a A-1 ~> 1]
integer(kind=int64) :: chksum !< checksum of array

! Local variables
real, allocatable :: field_rot(:,:) ! A rotated version of field, with the same units [arbitrary]
integer :: qturns ! The number of quarter turns through which to rotate field
logical :: do_unscale ! If true, unscale the variable before it is checksummed

qturns = 0
if (present(turns)) &
qturns = modulo(turns, 4)

do_unscale = .false. ; if (present(unscale)) do_unscale = (unscale /= 1.0)

if (qturns == 0) then
chksum = field_chksum(field, pelist=pelist, mask_val=mask_val)
if (do_unscale) then
chksum = field_chksum(unscale*field(:,:), pelist=pelist, mask_val=mask_val)
else
chksum = field_chksum(field, pelist=pelist, mask_val=mask_val)
endif
else
call allocate_rotated_array(field, [1,1], qturns, field_rot)
call rotate_array(field, qturns, field_rot)
if (do_unscale) field_rot(:,:) = unscale*field_rot(:,:)
chksum = field_chksum(field_rot, pelist=pelist, mask_val=mask_val)
deallocate(field_rot)
endif
end function rotated_field_chksum_real_2d
end function field_checksum_real_2d

!> Compute the field checksum of a rotated 3d field.
function rotated_field_chksum_real_3d(field, pelist, mask_val, turns) &
!> Compute the field checksum of an entire 3d array that may need to be rotated or unscaled.
!! This uses the field_chksum function that is used to verify file contents, which may differ
!! from the bitcount function used for other checksums in this module.
function field_checksum_real_3d(field, pelist, mask_val, turns, unscale) &
result(chksum)
real, dimension(:,:,:), intent(in) :: field !< Unrotated input field [arbitrary]
real, dimension(:,:,:), intent(in) :: field !< Unrotated input field to be checksummed in
!! arbitrary, possibly rescaled units [A ~> a]
integer, optional, intent(in) :: pelist(:) !< PE list of ranks to checksum
real, optional, intent(in) :: mask_val !< FMS mask value [nondim]
integer, optional, intent(in) :: turns !< Number of quarter turns
real, optional, intent(in) :: unscale !< A factor to convert this array back to
!! unscaled units for checksums [a A-1 ~> 1]
integer(kind=int64) :: chksum !< checksum of array

! Local variables
real, allocatable :: field_rot(:,:,:) ! A rotated version of field, with the same units [arbitrary]
integer :: qturns ! The number of quarter turns through which to rotate field
logical :: do_unscale ! If true, unscale the variable before it is checksummed

qturns = 0
if (present(turns)) &
qturns = modulo(turns, 4)

do_unscale = .false. ; if (present(unscale)) do_unscale = (unscale /= 1.0)

if (qturns == 0) then
chksum = field_chksum(field, pelist=pelist, mask_val=mask_val)
if (do_unscale) then
chksum = field_chksum(unscale*field(:,:,:), pelist=pelist, mask_val=mask_val)
else
chksum = field_chksum(field, pelist=pelist, mask_val=mask_val)
endif
else
call allocate_rotated_array(field, [1,1,1], qturns, field_rot)
call rotate_array(field, qturns, field_rot)
if (do_unscale) field_rot(:,:,:) = unscale*field_rot(:,:,:)
chksum = field_chksum(field_rot, pelist=pelist, mask_val=mask_val)
deallocate(field_rot)
endif
end function rotated_field_chksum_real_3d
end function field_checksum_real_3d

!> Compute the field checksum of a rotated 4d field.
function rotated_field_chksum_real_4d(field, pelist, mask_val, turns) &
!> Compute the field checksum of an entire 4d array that may need to be rotated or unscaled.
!! This uses the field_chksum function that is used to verify file contents, which may differ
!! from the bitcount function used for other checksums in this module.
function field_checksum_real_4d(field, pelist, mask_val, turns, unscale) &
result(chksum)
real, dimension(:,:,:,:), intent(in) :: field !< Unrotated input field [arbitrary]
real, dimension(:,:,:,:), intent(in) :: field !< Unrotated input field to be checksummed in
!! arbitrary, possibly rescaled units [A ~> a]
integer, optional, intent(in) :: pelist(:) !< PE list of ranks to checksum
real, optional, intent(in) :: mask_val !< FMS mask value [nondim]
integer, optional, intent(in) :: turns !< Number of quarter turns
real, optional, intent(in) :: unscale !< A factor to convert this array back to
!! unscaled units for checksums [a A-1 ~> 1]
integer(kind=int64) :: chksum !< checksum of array

! Local variables
real, allocatable :: field_rot(:,:,:,:) ! A rotated version of field, with the same units [arbitrary]
integer :: qturns ! The number of quarter turns through which to rotate field
logical :: do_unscale ! If true, unscale the variable before it is checksummed

qturns = 0
if (present(turns)) &
qturns = modulo(turns, 4)

do_unscale = .false. ; if (present(unscale)) do_unscale = (unscale /= 1.0)

if (qturns == 0) then
chksum = field_chksum(field, pelist=pelist, mask_val=mask_val)
if (do_unscale) then
chksum = field_chksum(unscale*field(:,:,:,:), pelist=pelist, mask_val=mask_val)
else
chksum = field_chksum(field, pelist=pelist, mask_val=mask_val)
endif
else
call allocate_rotated_array(field, [1,1,1,1], qturns, field_rot)
call rotate_array(field, qturns, field_rot)
if (do_unscale) field_rot(:,:,:,:) = unscale*field_rot(:,:,:,:)
chksum = field_chksum(field_rot, pelist=pelist, mask_val=mask_val)
deallocate(field_rot)
endif
end function rotated_field_chksum_real_4d
end function field_checksum_real_4d


!> Write a message including the checksum of the non-shifted array
Expand Down

0 comments on commit ad0a8b8

Please sign in to comment.