Skip to content

+Add symmetric_sum functions#1069

Open
Hallberg-NOAA wants to merge 1 commit intoNOAA-GFDL:dev/gfdlfrom
Hallberg-NOAA:add_symmetric_sum
Open

+Add symmetric_sum functions#1069
Hallberg-NOAA wants to merge 1 commit intoNOAA-GFDL:dev/gfdlfrom
Hallberg-NOAA:add_symmetric_sum

Conversation

@Hallberg-NOAA
Copy link
Copy Markdown
Member

Added the new interfaces symmetric_sum() and symmetric_sum_unit_tests() to the MOM_array_transform module, and added a call to symmetric_sum_unit_tests() to unit_tests(). The two functions wrapped by the symmetric_sum() interface are intended to be a general template for rotationally symmetric sums, but they demonstrably work as intended (as demonstrated by the passing unit tests) for any rectangular arrays. Although symmetric_sum() is not used yet in the MOM6 code apart from the tests, it was developed to address problems with failures of rotational symmetry in the diagnostic downscaling routines, and will be used there soon. All answers are bitwise identical, but there are two new publicly visible interface.

@Hallberg-NOAA Hallberg-NOAA added the enhancement New feature or request label Mar 26, 2026
Copy link
Copy Markdown
Member

@adcroft adcroft left a comment

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

A few issues found in the new code — see inline comments. Also, the symmetric_sum_unit_tests function is not reachable via make run.unit (the .testing unit test suite). That path auto-discovers standalone programs in config_src/drivers/unit_tests/. A new file config_src/drivers/unit_tests/test_MOM_array_transform.F90 is needed, following the pattern of e.g. test_MOM_string_functions.F90:

! This file is part of MOM6, the Modular Ocean Model version 6.
! See the LICENSE file for licensing information.
! SPDX-License-Identifier: Apache-2.0

program test_MOM_array_transform

use MOM_array_transform, only : symmetric_sum_unit_tests
use MOM_error_handler, only : set_skip_mpi

character(len=120) :: fail_message

call set_skip_mpi(.true.) ! This unit test is not expecting MPI to be used

if ( symmetric_sum_unit_tests(fail_message) ) then
  write(*,*) "FAILED: "//trim(fail_message)
  stop 1
endif

end program test_MOM_array_transform

No Makefile changes are needed — the target auto-discovers all .F90 files in that directory.

integer :: i, szi, szi_2

szi = size(field, 1)
szi_2 = szi / 2 ! Note that for an add number sz_2 is rounded down.
Copy link
Copy Markdown
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Typo: "add number" should be "odd number".

Copy link
Copy Markdown
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Good catch. It is only when this was pointed out that my brain was not automatically filling in "odd" for "add"!

sum = symmetric_sum_1d(field(:,1))
else
! This is the general case.
! Note that for an add number szi_2 is rounded down.
Copy link
Copy Markdown
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Typo: "add number" should be "odd number".

Copy link
Copy Markdown
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Corrected.

real :: sum !< The rotation dependent sum of the entries in field [A ~> a]

! Local variables
real :: quad_sum(2,2) ! The sums in each of the quadrants [A ~> a]
Copy link
Copy Markdown
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

quad_sum is declared here but never used in naive_sum_2d. Most compilers will warn about this unused variable.

Copy link
Copy Markdown
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Thank you for noting this; it has now been removed.

call rotate_array_real_2d(array, 2, ar_180)
call rotate_array_real_2d(array, 3, ar_270)

symmetric_sum_unit_tests = fail
Copy link
Copy Markdown
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

This assignment is dead code — fail is always .false. at this point, before the test loop runs. The correct final assignment is on line 565. This line should be removed.

Copy link
Copy Markdown
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Good point. This line has been removed.

Copy link
Copy Markdown
Member

@adcroft adcroft left a comment

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

This is a well-crafted PR — the recursive "Union-Jack" summation pattern is elegant, the unit tests are thorough (mixing positive and negative inexact values to amplify floating-point order-of-arithmetic differences is a nice touch), and it is great to see the unit test paradigm being followed with symmetric_sum_unit_tests. The code is clean, well-commented, and the dimensional annotation is consistent throughout.

A few small issues to address before merging:

  1. Missing standalone test driver: symmetric_sum_unit_tests is not reachable via make run.unit — a new config_src/drivers/unit_tests/test_MOM_array_transform.F90 is needed (see the earlier review comment for a ready-to-use template).
  2. Dead code on line 538: the first symmetric_sum_unit_tests = fail assignment (before the test loop) should be removed — see inline comment.
  3. Two typos ("add" → "odd") in comments — see inline comments on lines 384 and 420.
  4. Unused variable quad_sum in naive_sum_2d — see inline comment on line 488.

  Added the new interfaces `symmetric_sum()` and `symmetric_sum_unit_tests()` to
the MOM_array_transform module, and added a call to `symmetric_sum_unit_tests()`
to `unit_tests()`.  The two functions wrapped by the `symmetric_sum()` interface
are intended to be a general template for rotationally symmetric sums, but they
demonstrably work as intended (as demonstrated by the passing unit tests) for
any rectangular arrays.  The new driver in test_MOM_array_transform.F90 tests
these new arrays.  Although `symmetric_sum()` is not used yet in the MOM6 code
apart from the tests, it was developed to address problems with failures of
rotational symmetry in the diagnostic downscaling routines, and will be used
there soon.  All answers are bitwise identical, but there are two new publicly
visible interface.
@Hallberg-NOAA
Copy link
Copy Markdown
Member Author

In addition to correcting the noted errors, I have also added the new test-driver in test_MOM_array_transform.F90. It differs somewhat from what was suggested above in that it uses calls to MOM_error() and MOM_mesg() to handle output in the failing and passing cases, respectively, but is otherwise functionally equivalent to the suggestion.

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment

Labels

enhancement New feature or request

Projects

None yet

Development

Successfully merging this pull request may close these issues.

2 participants