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 and use the new interface field_checksum #841

Merged
merged 3 commits into from
Feb 24, 2025
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
15 changes: 8 additions & 7 deletions src/diagnostics/MOM_sum_output.F90
Original file line number Diff line number Diff line change
Expand Up @@ -4,8 +4,8 @@ module MOM_sum_output
! This file is part of MOM6. See LICENSE.md for the license.

use iso_fortran_env, only : int64
use MOM_checksums, only : is_NaN
use MOM_coms, only : sum_across_PEs, PE_here, root_PE, num_PEs, max_across_PEs, field_chksum
use MOM_checksums, only : is_NaN, field_checksum
use MOM_coms, only : sum_across_PEs, PE_here, root_PE, num_PEs, max_across_PEs
use MOM_coms, only : reproducing_sum, reproducing_sum_EFP, EFP_to_real, real_to_EFP
use MOM_coms, only : EFP_type, operator(+), operator(-), assignment(=), EFP_sum_across_PEs
use MOM_error_handler, only : MOM_error, FATAL, WARNING, NOTE, is_root_pe
Expand Down Expand Up @@ -1439,22 +1439,23 @@ subroutine get_depth_list_checksums(G, US, depth_chksum, area_chksum)
character(len=16), intent(out) :: depth_chksum !< Depth checksum hexstring
character(len=16), intent(out) :: area_chksum !< Area checksum hexstring

! Local variables
real, allocatable :: field(:,:) ! A temporary array with no halos [Z ~> m] or [L2 ~> m2]
integer :: i, j
real, allocatable :: field(:,:) ! A temporary array for output converted to MKS units [m] or [m2]

allocate(field(G%isc:G%iec, G%jsc:G%jec))

! Depth checksum
do j=G%jsc,G%jec ; do i=G%isc,G%iec
field(i,j) = US%Z_to_m*(G%bathyT(i,j) + G%Z_ref)
field(i,j) = G%bathyT(i,j) + G%Z_ref
enddo ; enddo
write(depth_chksum, '(Z16)') field_chksum(field(:,:))
write(depth_chksum, '(Z16)') field_checksum(field(:,:), unscale=US%Z_to_m)

! Area checksum
do j=G%jsc,G%jec ; do i=G%isc,G%iec
field(i,j) = G%mask2dT(i,j) * US%L_to_m**2*G%areaT(i,j)
field(i,j) = G%mask2dT(i,j) * G%areaT(i,j)
enddo ; enddo
write(area_chksum, '(Z16)') field_chksum(field(:,:))
write(area_chksum, '(Z16)') field_checksum(field(:,:), unscale=US%L_to_m**2)

deallocate(field)
end subroutine get_depth_list_checksums
Expand Down
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
22 changes: 11 additions & 11 deletions src/framework/MOM_restart.F90
Original file line number Diff line number Diff line change
Expand Up @@ -5,7 +5,7 @@ module MOM_restart

use, intrinsic :: iso_fortran_env, only : int64
use MOM_array_transform, only : rotate_array, rotate_vector, rotate_array_pair
use MOM_checksums, only : chksum => rotated_field_chksum
use MOM_checksums, only : chksum => field_checksum
use MOM_domains, only : PE_here, num_PEs, AGRID, BGRID_NE, CGRID_NE
use MOM_error_handler, only : MOM_error, MOM_mesg, FATAL, WARNING, NOTE, is_root_pe, MOM_get_verbosity
use MOM_file_parser, only : get_param, log_param, log_version, param_file_type
Expand Down Expand Up @@ -1738,15 +1738,15 @@ subroutine save_restart(directory, time, G, CS, time_stamped, filename, GV, num_

conv = CS%restart_field(m)%conv
if (associated(CS%var_ptr3d(m)%p)) then
check_val(m-start_var+1,1) = chksum(conv*CS%var_ptr3d(m)%p(isL:ieL,jsL:jeL,:), turns=-turns)
check_val(m-start_var+1,1) = chksum(CS%var_ptr3d(m)%p(isL:ieL,jsL:jeL,:), turns=-turns, unscale=conv)
elseif (associated(CS%var_ptr2d(m)%p)) then
check_val(m-start_var+1,1) = chksum(conv*CS%var_ptr2d(m)%p(isL:ieL,jsL:jeL), turns=-turns)
check_val(m-start_var+1,1) = chksum(CS%var_ptr2d(m)%p(isL:ieL,jsL:jeL), turns=-turns, unscale=conv)
elseif (associated(CS%var_ptr4d(m)%p)) then
check_val(m-start_var+1,1) = chksum(conv*CS%var_ptr4d(m)%p(isL:ieL,jsL:jeL,:,:), turns=-turns)
check_val(m-start_var+1,1) = chksum(CS%var_ptr4d(m)%p(isL:ieL,jsL:jeL,:,:), turns=-turns, unscale=conv)
elseif (associated(CS%var_ptr1d(m)%p)) then
check_val(m-start_var+1,1) = chksum(conv*CS%var_ptr1d(m)%p(:))
check_val(m-start_var+1,1) = chksum(CS%var_ptr1d(m)%p(:), unscale=conv)
elseif (associated(CS%var_ptr0d(m)%p)) then
check_val(m-start_var+1,1) = chksum(conv*CS%var_ptr0d(m)%p, pelist=(/PE_here()/))
check_val(m-start_var+1,1) = chksum(CS%var_ptr0d(m)%p, pelist=(/PE_here()/), unscale=conv)
endif
enddo

Expand Down Expand Up @@ -1924,11 +1924,11 @@ subroutine restore_state(filename, directory, day, G, CS)
! Read a 1d array, which should be invariant to domain decomposition.
call MOM_read_data(unit_path(n), varname, CS%var_ptr1d(m)%p, &
timelevel=1, scale=scale, MOM_Domain=G%Domain)
if (is_there_a_checksum) checksum_data = chksum(conv*CS%var_ptr1d(m)%p(:))
if (is_there_a_checksum) checksum_data = chksum(CS%var_ptr1d(m)%p(:), unscale=conv)
elseif (associated(CS%var_ptr0d(m)%p)) then ! Read a scalar...
call MOM_read_data(unit_path(n), varname, CS%var_ptr0d(m)%p, &
timelevel=1, scale=scale, MOM_Domain=G%Domain)
if (is_there_a_checksum) checksum_data = chksum(conv*CS%var_ptr0d(m)%p, pelist=(/PE_here()/))
if (is_there_a_checksum) checksum_data = chksum(CS%var_ptr0d(m)%p, pelist=(/PE_here()/), unscale=conv)
elseif (associated(CS%var_ptr2d(m)%p)) then ! Read a 2d array.
if (pos /= 0) then
call MOM_read_data(unit_path(n), varname, CS%var_ptr2d(m)%p, &
Expand All @@ -1938,7 +1938,7 @@ subroutine restore_state(filename, directory, day, G, CS)
"MOM_restart does not support 2-d arrays without domain decomposition.")
! call read_data(unit_path(n), varname, CS%var_ptr2d(m)%p,no_domain=.true., timelevel=1)
endif
if (is_there_a_checksum) checksum_data = chksum(conv*CS%var_ptr2d(m)%p(isL:ieL,jsL:jeL))
if (is_there_a_checksum) checksum_data = chksum(CS%var_ptr2d(m)%p(isL:ieL,jsL:jeL), unscale=conv)
elseif (associated(CS%var_ptr3d(m)%p)) then ! Read a 3d array.
if (pos /= 0) then
call MOM_read_data(unit_path(n), varname, CS%var_ptr3d(m)%p, &
Expand All @@ -1948,7 +1948,7 @@ subroutine restore_state(filename, directory, day, G, CS)
"MOM_restart does not support 3-d arrays without domain decomposition.")
! call read_data(unit_path(n), varname, CS%var_ptr3d(m)%p, no_domain=.true., timelevel=1)
endif
if (is_there_a_checksum) checksum_data = chksum(conv*CS%var_ptr3d(m)%p(isL:ieL,jsL:jeL,:))
if (is_there_a_checksum) checksum_data = chksum(CS%var_ptr3d(m)%p(isL:ieL,jsL:jeL,:), unscale=conv)
elseif (associated(CS%var_ptr4d(m)%p)) then ! Read a 4d array.
if (pos /= 0) then
call MOM_read_data(unit_path(n), varname, CS%var_ptr4d(m)%p, &
Expand All @@ -1959,7 +1959,7 @@ subroutine restore_state(filename, directory, day, G, CS)
"MOM_restart does not support 4-d arrays without domain decomposition.")
! call read_data(unit_path(n), varname, CS%var_ptr4d(m)%p, no_domain=.true., timelevel=1)
endif
if (is_there_a_checksum) checksum_data = chksum(conv*CS%var_ptr4d(m)%p(isL:ieL,jsL:jeL,:,:))
if (is_there_a_checksum) checksum_data = chksum(CS%var_ptr4d(m)%p(isL:ieL,jsL:jeL,:,:), unscale=conv)
else
call MOM_error(FATAL, "MOM_restart restore_state: No pointers set for "//trim(varname))
endif
Expand Down