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 adds the associated_files attribute #1451

Merged
merged 12 commits into from
Feb 9, 2024
Next Next commit
Fix the associated_files attribute and add a test
  • Loading branch information
uramirez8707 committed Feb 1, 2024
commit c289ea42fe6511556be18bfc98b2a73bc3208f39
33 changes: 31 additions & 2 deletions diag_manager/fms_diag_field_object.F90
Original file line number Diff line number Diff line change
Expand Up @@ -8,7 +8,7 @@ module fms_diag_field_object_mod
!! that contains all of the information of the variable. It is extended by a type that holds the
!! appropriate buffer for the data for manipulation.
#ifdef use_yaml
use diag_data_mod, only: diag_null, CMOR_MISSING_VALUE, diag_null_string, MAX_STR_LEN
use diag_data_mod, only: prepend_date, diag_null, CMOR_MISSING_VALUE, diag_null_string, MAX_STR_LEN
use diag_data_mod, only: r8, r4, i8, i4, string, null_type_int, NO_DOMAIN
use diag_data_mod, only: max_field_attributes, fmsDiagAttribute_type
use diag_data_mod, only: diag_null, diag_not_found, diag_not_registered, diag_registered_id, &
Expand All @@ -20,7 +20,7 @@ module fms_diag_field_object_mod
& find_diag_field, get_num_unique_fields, diag_yaml
use fms_diag_axis_object_mod, only: diagDomain_t, get_domain_and_domain_type, fmsDiagAxis_type, &
& fmsDiagAxisContainer_type, fmsDiagFullAxis_Type
use time_manager_mod, ONLY: time_type
use time_manager_mod, ONLY: time_type, get_date
use fms2_io_mod, only: FmsNetcdfFile_t, FmsNetcdfDomainFile_t, FmsNetcdfUnstructuredDomainFile_t, register_field, &
register_variable_attribute
use fms_diag_input_buffer_mod, only: fmsDiagInputBuffer_t
Expand Down Expand Up @@ -175,6 +175,7 @@ module fms_diag_field_object_mod
procedure :: has_mask_allocated
procedure :: is_variable_in_file
procedure :: get_field_file_name
procedure :: generate_associated_files_att
end type fmsDiagField_type
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! variables !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
type(fmsDiagField_type) :: null_ob
Expand Down Expand Up @@ -1771,5 +1772,33 @@ function get_field_file_name(this) &
res = this%diag_field(1)%get_var_fname()
end function get_field_file_name

!> @brief Generate the associated files attribute
subroutine generate_associated_files_att(this, att, start_time)
class(fmsDiagField_type) , intent(in) :: this !< diag_field_object for the area/volume field
character(len=*), intent(inout) :: att !< associated_files_att
type(time_type), intent(in) :: start_time
Copy link
Contributor Author

Choose a reason for hiding this comment

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

documentation


character(len=:), allocatable :: field_name !< Name of the area/volume field
character(len=MAX_STR_LEN) :: file_name !< Name of the file the area/volume field is in!
character(len=128) :: start_date !< Start date to append to the begining of the filename

integer :: year, month, day, hour, minute, second
field_name = this%get_varname(to_write = .true.)

! Check if the field is already in the associated files attribute (i.e the area can be associated with multiple
! fields in the file, but it only needs to be added once)
if (index(att, field_name) .ne. 0) return

file_name = this%get_field_file_name()

if (prepend_date) then
call get_date(start_time, year, month, day, hour, minute, second)
write (start_date, '(1I20.4, 2I2.2)') year, month, day
file_name = TRIM(adjustl(start_date))//'.'//TRIM(file_name)
endif

att = trim(att)//" "//trim(field_name)//": "//trim(file_name)//".nc"
end subroutine generate_associated_files_att

#endif
end module fms_diag_field_object_mod
17 changes: 8 additions & 9 deletions diag_manager/fms_diag_file_object.F90
Original file line number Diff line number Diff line change
Expand Up @@ -1403,21 +1403,21 @@ subroutine write_field_metadata(this, diag_field, diag_axis)
diag_file => this%FMS_diag_file
fms2io_fileobj => diag_file%fms2io_fileobj

associated_files = ""
need_associated_files = .false.
do i = 1, size(diag_file%field_ids)
if (.not. diag_file%field_registered(i)) cycle !TODO do something else here
field_ptr => diag_field(diag_file%field_ids(i))

cell_measures = ""
associated_files = ""
need_associated_files = .false.
if (field_ptr%has_area()) then
cell_measures = "area: "//diag_field(field_ptr%get_area())%get_varname(to_write=.true.)

!! Determine if the area field is already in the file. If it is not create the "associated_files" attribute
!! which contains the file name of the file the area field is in. This is needed for PP/fregrid.
if (.not. diag_field(field_ptr%get_area())%is_variable_in_file(diag_file%id)) then
need_associated_files = .true.
associated_files = "area: "//diag_field(field_ptr%get_area())%get_field_file_name()//".nc"
call diag_field(field_ptr%get_area())%generate_associated_files_att(associated_files, diag_file%start_time)
endif
endif

Expand All @@ -1428,19 +1428,18 @@ subroutine write_field_metadata(this, diag_field, diag_axis)
!! which contains the file name of the file the volume field is in. This is needed for PP/fregrid.
if (.not. diag_field(field_ptr%get_volume())%is_variable_in_file(diag_file%id)) then
need_associated_files = .true.
associated_files = trim(associated_files)//&
" volume:"//diag_field(field_ptr%get_volume())%get_field_file_name()//".nc"
call diag_field(field_ptr%get_volume())%generate_associated_files_att(associated_files, diag_file%start_time)
endif
endif

call field_ptr%write_field_metadata(fms2io_fileobj, diag_file%id, diag_file%yaml_ids(i), diag_axis, &
this%FMS_diag_file%get_file_unlimdim(), is_regional, cell_measures)

if (need_associated_files) &
call register_global_attribute(fms2io_fileobj, "associated_files", trim(ADJUSTL(associated_files)), &
str_len=len_trim(ADJUSTL(associated_files)))
enddo

if (need_associated_files) &
call register_global_attribute(fms2io_fileobj, "associated_files", trim(ADJUSTL(associated_files)), &
str_len=len_trim(ADJUSTL(associated_files)))

end subroutine write_field_metadata

!< @brief Writes the axis data for the file
Expand Down
8 changes: 5 additions & 3 deletions test_fms/diag_manager/Makefile.am
Original file line number Diff line number Diff line change
Expand Up @@ -32,7 +32,7 @@ check_PROGRAMS = test_diag_manager test_diag_manager_time \
test_diag_yaml test_diag_ocean test_modern_diag test_diag_buffer \
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_time_pow check_time_rms test_cell_measures

# This is the source code for the test.
test_diag_manager_SOURCES = test_diag_manager.F90
Expand All @@ -53,20 +53,22 @@ check_time_avg_SOURCES = testing_utils.F90 check_time_avg.F90
check_time_diurnal_SOURCES = testing_utils.F90 check_time_diurnal.F90
check_time_pow_SOURCES = testing_utils.F90 check_time_pow.F90
check_time_rms_SOURCES = testing_utils.F90 check_time_rms.F90
test_cell_measures_SOURCES = test_cell_measures.F90

TEST_EXTENSIONS = .sh
SH_LOG_DRIVER = env AM_TAP_AWK='$(AWK)' $(SHELL) \
$(abs_top_srcdir)/test_fms/tap-driver.sh

# 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_time_avg.sh test_time_pow.sh test_time_rms.sh test_time_diurnal.sh test_cell_measures.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_time_sum.sh test_time_avg.sh test_time_pow.sh test_time_rms.sh test_time_diurnal.sh \
test_cell_measures.sh

if USING_YAML
skipflag=""
Expand Down
105 changes: 105 additions & 0 deletions test_fms/diag_manager/test_cell_measures.F90
Original file line number Diff line number Diff line change
@@ -0,0 +1,105 @@
!***********************************************************************
!* 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/>.
!***********************************************************************

!> @brief This program tests the diag_manager with fields with cell measures (area, volume)
program test_cell_measures
use fms_mod, only: fms_init, fms_end
use diag_manager_mod, only: diag_axis_init, register_static_field, diag_send_complete, send_data
use diag_manager_mod, only : register_diag_field
use diag_manager_mod, only: diag_manager_init, diag_manager_end, diag_manager_set_time_end
use time_manager_mod, only: time_type, set_calendar_type, set_date, JULIAN, set_time, OPERATOR(+)
use mpp_mod, only: mpp_error, FATAL
use fms2_io_mod
use platform_mod, only: r4_kind

implicit none

type(time_type) :: Time !< Time of the simulation
type(time_type) :: Time_step !< Time_step of the simulation
integer :: i !< For looping
integer :: id_axis1 !< Id of axis1
integer :: naxis1 !< Size of axis1
real(kind=r4_kind), allocatable :: axis1_data(:) !< Data for axis1
integer :: id_var1 !< Id of var1
real(kind=r4_kind), allocatable :: var1_data(:) !< Data for "var1"
real(kind=r4_kind), allocatable :: area_data(:) !< Data for the "area"
integer :: id_area !< Id of the "area" field
logical :: used !< Used for send_data call

naxis1 = 10
call fms_init()
call set_calendar_type(JULIAN)
call diag_manager_init()

Time = set_date(2,1,1,0,0,0)
Time_step = set_time (3600,0)
call diag_manager_set_time_end(set_date(2,1,2,0,0,0))

allocate(axis1_data(naxis1))
allocate(var1_data(naxis1))
allocate(area_data(naxis1))
do i = 1, naxis1
axis1_data = real(i, kind=r4_kind)
area_data = real(i/100, kind=r4_kind)
var1_data = real(i*10, kind=r4_kind)
enddo

id_axis1 = diag_axis_init('axis1', axis1_data, 'axis1', 'x')
id_area = register_static_field ('fun_mod', 'area', (/id_axis1/))
id_var1 = register_diag_field ('fun_mod', 'var1', (/id_axis1/), init_time=Time, area=id_area)

used = send_data(id_area, area_data)

do i = 1, 6
Time = Time + Time_step
call diag_send_complete(Time_step)
used = send_data(id_var1, var1_data, Time)
enddo
call diag_manager_end(Time)

call check_output()
call fms_end()

contains
subroutine check_output()
type(FmsNetcdfFile_t) :: fileobj !< FMS2io fileobj
character(len=256) :: buffer !< Buffer to read stuff into

! Check that the static_file.nc was created and it contains the area attribute
if (.not. open_file(fileobj, "static_file.nc", "read")) &
call mpp_error(FATAL, "static_file.nc was not created by the diag manager!")
if (.not. variable_exists(fileobj, "area")) &
call mpp_error(FATAL, "area is not in static_file.nc")
call close_file(fileobj)

! Check that file1.nc exists, that it contains the associated files attribute and it is correct,
! that the var1 exists and it contains the cell_measures attributes
if (.not. open_file(fileobj, "file1.nc", "read")) &
call mpp_error(FATAL, "file1.nc was not created by the diag manager!")

call get_global_attribute(fileobj, "associated_files", buffer)
if (trim(buffer) .ne. "area: static_file.nc") &
call mpp_error(FATAL, "The associated_files global attribute is not the expected result! "//trim(buffer))

call get_variable_attribute(fileobj, "var1", "cell_measures", buffer)
if (trim(buffer) .ne. "area: area") &
call mpp_error(FATAL, "The cell_measures attribute is not the expected result! "//trim(buffer))
call close_file(fileobj)
end subroutine check_output
end program
64 changes: 64 additions & 0 deletions test_fms/diag_manager/test_cell_measures.sh
Original file line number Diff line number Diff line change
@@ -0,0 +1,64 @@
#!/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_diag_manager
base_date: 2 1 1 0 0 0

diag_files:
- file_name: static_file
freq: -1
time_units: hours
unlimdim: time
varlist:
- module: fun_mod
var_name: area
reduction: none
kind: r4
# Here file 1 does not have the "area" variable so the associated files attribute is expected
- file_name: file1
freq: 6 hours
time_units: hours
unlimdim: time
varlist:
- module: fun_mod
var_name: var1
reduction: average
kind: r4
_EOF

# remove any existing files that would result in false passes during checks
rm -f *.nc

my_test_count=1
printf "&diag_manager_nml \n use_modern_diag=.true. \n/" | cat > input.nml
test_expect_success "Running diag_manager with fields with cell measures (area, volume) (test $my_test_count)" '
mpirun -n 1 ../test_cell_measures
'
fi
test_done