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

Modern diag manager: Correctly uses the weight argument from send_data #1513

Merged
merged 4 commits into from
May 9, 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
10 changes: 6 additions & 4 deletions diag_manager/fms_diag_object.F90
Original file line number Diff line number Diff line change
Expand Up @@ -1057,21 +1057,22 @@ function fms_diag_do_reduction(this, field_data, diag_field_id, oor_mask, weight
endif
case (time_average)
error_msg = buffer_ptr%do_time_sum_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
field_ptr%get_mask_variant(), bounds_in, bounds_out, missing_value, field_ptr%has_missing_value())
field_ptr%get_mask_variant(), bounds_in, bounds_out, missing_value, field_ptr%has_missing_value(), &
weight=weight)
if (trim(error_msg) .ne. "") then
return
endif
case (time_power)
error_msg = buffer_ptr%do_time_sum_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
field_ptr%get_mask_variant(), bounds_in, bounds_out, missing_value, field_ptr%has_missing_value(), &
pow_value=field_yaml_ptr%get_pow_value())
weight=weight, pow_value=field_yaml_ptr%get_pow_value())
if (trim(error_msg) .ne. "") then
return
endif
case (time_rms)
error_msg = buffer_ptr%do_time_sum_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
field_ptr%get_mask_variant(), bounds_in, bounds_out, missing_value, field_ptr%has_missing_value(), &
pow_value = 2)
weight=weight, pow_value = 2)
if (trim(error_msg) .ne. "") then
return
endif
Expand All @@ -1081,7 +1082,8 @@ function fms_diag_do_reduction(this, field_data, diag_field_id, oor_mask, weight
! sets the diurnal index for reduction within the buffer object
call buffer_ptr%set_diurnal_section_index(time)
error_msg = buffer_ptr%do_time_sum_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
field_ptr%get_mask_variant(), bounds_in, bounds_out, missing_value, field_ptr%has_missing_value())
field_ptr%get_mask_variant(), bounds_in, bounds_out, missing_value, field_ptr%has_missing_value(), &
weight=weight)
if (trim(error_msg) .ne. "") then
return
endif
Expand Down
7 changes: 4 additions & 3 deletions diag_manager/fms_diag_output_buffer.F90
Original file line number Diff line number Diff line change
Expand Up @@ -716,7 +716,7 @@ end function do_time_max_wrapper
!> @brief Does the time_sum reduction method on the buffer object
!! @return Error message if the math was not successful
function do_time_sum_wrapper(this, field_data, mask, is_masked, mask_variant, bounds_in, bounds_out, missing_value, &
has_missing_value, pow_value) &
has_missing_value, pow_value, weight) &
result(err_msg)
class(fmsDiagOutputBuffer_type), intent(inout) :: this !< buffer object to write
class(*), intent(in) :: field_data(:,:,:,:) !< Buffer data for current time
Expand All @@ -731,6 +731,7 @@ function do_time_sum_wrapper(this, field_data, mask, is_masked, mask_variant, bo
integer, optional, intent(in) :: pow_value !< power value, will calculate field_data^pow
!! before adding to buffer should only be
!! present if using pow reduction method
real(kind=r8_kind), optional, intent(in) :: weight !< The weight to use when suming
character(len=150) :: err_msg

!TODO This will be expanded for integers
Expand All @@ -745,7 +746,7 @@ function do_time_sum_wrapper(this, field_data, mask, is_masked, mask_variant, bo
endif
call do_time_sum_update(output_buffer, this%weight_sum, field_data, mask, is_masked, mask_variant, &
bounds_in, bounds_out, missing_value, this%diurnal_section, &
pow=pow_value)
pow=pow_value, weight=weight)
class default
err_msg="do_time_sum_wrapper::the output buffer and the buffer send in are not of the same type (r8_kind)"
end select
Expand All @@ -758,7 +759,7 @@ function do_time_sum_wrapper(this, field_data, mask, is_masked, mask_variant, bo
endif
call do_time_sum_update(output_buffer, this%weight_sum, field_data, mask, is_masked, mask_variant, &
bounds_in, bounds_out, real(missing_value, kind=r4_kind), &
this%diurnal_section, pow=pow_value)
this%diurnal_section, pow=pow_value, weight=weight)
class default
err_msg="do_time_sum_wrapper::the output buffer and the buffer send in are not of the same type (r4_kind)"
end select
Expand Down
2 changes: 1 addition & 1 deletion diag_manager/include/fms_diag_reduction_methods.inc
Original file line number Diff line number Diff line change
Expand Up @@ -247,7 +247,7 @@ subroutine DO_TIME_SUM_UPDATE_(data_out, weight_sum, data_in, mask, is_masked, m
!! will be 1 unless using a diurnal reduction

if(present(weight)) then
weight_scale = weight
weight_scale = real(weight, kind=kindl)
else
weight_scale = 1.0_kindl
endif
Expand Down
9 changes: 6 additions & 3 deletions test_fms/diag_manager/Makefile.am
Original file line number Diff line number Diff line change
Expand Up @@ -33,10 +33,12 @@ check_PROGRAMS = test_diag_manager test_diag_manager_time \
test_flexible_time test_diag_update_buffer test_reduction_methods check_time_none \
check_time_min check_time_max check_time_sum check_time_avg test_diag_diurnal check_time_diurnal \
check_time_pow check_time_rms check_subregional test_cell_measures test_var_masks \
check_var_masks test_multiple_send_data test_diag_out_yaml test_output_every_freq
check_var_masks test_multiple_send_data test_diag_out_yaml test_output_every_freq \
test_dm_weights

# This is the source code for the test.
test_output_every_freq_SOURCES = test_output_every_freq.F90
test_dm_weights_SOURCES = test_dm_weights.F90
test_diag_manager_SOURCES = test_diag_manager.F90
test_diag_manager_time_SOURCES = test_diag_manager_time.F90
test_diag_update_buffer_SOURCES= test_diag_update_buffer.F90
Expand Down Expand Up @@ -69,15 +71,16 @@ SH_LOG_DRIVER = env AM_TAP_AWK='$(AWK)' $(SHELL) \
# Run the test.
TESTS = test_diag_manager2.sh test_time_none.sh test_time_min.sh test_time_max.sh test_time_sum.sh \
test_time_avg.sh test_time_pow.sh test_time_rms.sh test_time_diurnal.sh test_cell_measures.sh \
test_subregional.sh test_var_masks.sh test_multiple_send_data.sh test_output_every_freq.sh
test_subregional.sh test_var_masks.sh test_multiple_send_data.sh test_output_every_freq.sh \
test_dm_weights.sh

testing_utils.mod: testing_utils.$(OBJEXT)

# Copy over other needed files to the srcdir
EXTRA_DIST = test_diag_manager2.sh check_crashes.sh test_time_none.sh test_time_min.sh test_time_max.sh \
test_time_sum.sh test_time_avg.sh test_time_pow.sh test_time_rms.sh test_time_diurnal.sh \
test_cell_measures.sh test_subregional.sh test_var_masks.sh test_multiple_send_data.sh \
test_output_every_freq.sh
test_dm_weights.sh test_output_every_freq.sh

if USING_YAML
skipflag=""
Expand Down
87 changes: 87 additions & 0 deletions test_fms/diag_manager/test_dm_weights.F90
Original file line number Diff line number Diff line change
@@ -0,0 +1,87 @@
program test_dm_weights
use fms_mod, only: fms_init, fms_end
use diag_manager_mod, only: diag_manager_init, diag_axis_init, register_diag_field, send_data, &
diag_send_complete, diag_manager_set_time_end, diag_manager_end
use mpp_mod, only: mpp_error, mpp_pe, mpp_root_pe, FATAL
use time_manager_mod, only: time_type, set_calendar_type, JULIAN, set_time, set_date, operator(+)
use fms2_io_mod

implicit none

type(time_type) :: Time !< Time of the simulation
type(time_type) :: Time_step !< Time_step of the simulation
integer :: nx !< Number of x points
integer :: ny !< Number of y points
integer :: id_x !< Axis id for the x dimension
integer :: id_y !< Axis id for the y dimension
integer :: id_var1 !< Field id for 1st variable
integer :: id_var2 !< Field id for 2nd variable
logical :: used !< Dummy argument to send_data
real, allocatable :: x(:) !< X axis data
real, allocatable :: y(:) !< Y axis_data
real, allocatable :: var1_data(:,:) !< Data for variable 1
integer :: i !< For do loops
integer :: ntimes !< Number of times to run the simulation for

call fms_init()
call set_calendar_type(JULIAN)
call diag_manager_init()

nx = 10
ny = 15
ntimes = 6

allocate(x(nx), y(ny))
allocate(var1_data(nx,ny))
do i=1,nx
x(i) = i
enddo
do i=1,ny
y(i) = -91 + i
enddo

Time = set_date(2,1,1,0,0,0)
Time_step = set_time (3600,0) !< 1 hour

id_x = diag_axis_init('x', x, 'point_E', 'x', long_name='point_E')
id_y = diag_axis_init('y', y, 'point_N', 'y', long_name='point_N')

id_var1 = register_diag_field ('atmos', 'ua', (/id_x, id_y/), Time)

call diag_manager_set_time_end(set_date(2,1,1,ntimes,0,0))
do i = 1, ntimes
Time = Time + Time_step
var1_data = real(i)
used = send_data(id_var1, var1_data, Time, weight=real(i/10.))
call diag_send_complete(Time_step)
enddo

call diag_manager_end(Time)

call check_answers()
call fms_end()

contains

subroutine check_answers()
type(FmsNetcdfFile_t) :: fileobj
integer :: j
real :: ans_var
real :: vardata_out(nx, ny)

if (.not. open_file(fileobj, "test_weights.nc", "read")) &
call mpp_error(FATAL, "unable to open test_var_masks.nc for reading")

ans_var = 0
do j = 1, ntimes
if (mod(j,2) .eq. 0) then
print *, "Checking answers for time = ", j/2
ans_var = (j*(j/10.)+(j-1.)*(j-1.)/10.)/(j/10. + (j-1)/10.)
call read_data(fileobj, "ua", vardata_out, unlim_dim_level=j/2)
if (any(abs(ans_var-vardata_out) > 0.0000001)) &
call mpp_error(FATAL, "The answer is not the expected result!")
endif
enddo
call close_file(fileobj)
end subroutine check_answers
end program test_dm_weights
60 changes: 60 additions & 0 deletions test_fms/diag_manager/test_dm_weights.sh
Original file line number Diff line number Diff line change
@@ -0,0 +1,60 @@
#!/bin/sh

#***********************************************************************
#* GNU Lesser General Public License
#*
#* This file is part of the GFDL Flexible Modeling System (FMS).
#*
#* FMS is free software: you can redistribute it and/or modify it under
#* the terms of the GNU Lesser General Public License as published by
#* the Free Software Foundation, either version 3 of the License, or (at
#* your option) any later version.
#*
#* FMS 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 Lesser General Public
#* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
#***********************************************************************

# Set common test settings.
. ../test-lib.sh

if [ -z "${skipflag}" ]; then
# create and enter directory for in/output files
output_dir

cat <<_EOF > diag_table.yaml
title: test_weights
base_date: 2 1 1 0 0 0
diag_files:
- file_name: test_weights
time_units: hours
unlimdim: time
freq: 2 hours
varlist:
- module: atmos
var_name: ua
reduction: average
kind: r4
_EOF

# remove any existing files that would result in false passes during checks
rm -f *.nc
export OMP_NUM_THREADS=1
my_test_count=1
printf "&diag_manager_nml \n use_modern_diag=.true. \n/" | cat > input.nml
test_expect_success "Running diag_manager with weight passed in no threads (test $my_test_count)" '
mpirun -n 1 ../test_var_masks
'

export OMP_NUM_THREADS=2
my_test_count=`expr $my_test_count + 1`
test_expect_success "Running diag_manager with weight passed in 2 threads (test $my_test_count)" '
mpirun -n 1 ../test_var_masks
'
export OMP_NUM_THREADS=1
fi
test_done
12 changes: 6 additions & 6 deletions test_fms/diag_manager/test_time_diurnal.sh
Original file line number Diff line number Diff line change
Expand Up @@ -83,7 +83,7 @@ test_expect_success "monthly simple diurnal output" '
'

test_expect_success "checking results for diurnal test simple" '
mpirun -n 6 ../check_time_diurnal
mpirun -n 1 ../check_time_diurnal
'

printf "&test_diag_diurnal_nml \n test_case=0 \n mask_case=1 \n / \n" >> input.nml
Expand All @@ -92,7 +92,7 @@ test_expect_success "monthly diurnal output with logical mask" '
mpirun -n 6 ../test_diag_diurnal
'
test_expect_success "checking results for diurnal test with logical mask" '
mpirun -n 6 ../check_time_diurnal
mpirun -n 1 ../check_time_diurnal
'

printf "&test_diag_diurnal_nml \n test_case=0 \n mask_case=2 \n / \n" >> input.nml
Expand All @@ -101,7 +101,7 @@ test_expect_success "monthly diurnal output with real mask" '
mpirun -n 6 ../test_diag_diurnal
'
test_expect_success "checking results for diurnal test with real mask" '
mpirun -n 6 ../check_time_diurnal
mpirun -n 1 ../check_time_diurnal
'

export OMP_NUM_THREADS=2
Expand All @@ -112,7 +112,7 @@ test_expect_success "monthly diurnal output with openmp" '
mpirun -n 6 ../test_diag_diurnal
'
test_expect_success "checking results for diurnal test with openmp" '
mpirun -n 6 ../check_time_diurnal
mpirun -n 1 ../check_time_diurnal
'

printf "&test_diag_diurnal_nml \n test_case=1 \n mask_case=1 \n / \n" >> input.nml
Expand All @@ -121,7 +121,7 @@ test_expect_success "monthly diurnal output with openmp and real mask" '
mpirun -n 6 ../test_diag_diurnal
'
test_expect_success "checking results for diurnal test with openmp and real mask" '
mpirun -n 6 ../check_time_diurnal
mpirun -n 1 ../check_time_diurnal
'

printf "&test_diag_diurnal_nml \n test_case=1 \n mask_case=2 \n / \n" >> input.nml
Expand All @@ -130,7 +130,7 @@ test_expect_success "monthly diurnal output with openmp and logical mask" '
mpirun -n 6 ../test_diag_diurnal
'
test_expect_success "checking results for diurnal test with openmp and logical mask" '
mpirun -n 6 ../check_time_diurnal
mpirun -n 1 ../check_time_diurnal
'

fi
Expand Down
Loading