Skip to content
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
14 changes: 11 additions & 3 deletions diag_manager/fms_diag_axis_object.F90
Original file line number Diff line number Diff line change
Expand Up @@ -458,7 +458,13 @@ subroutine write_axis_data(this, fms2io_fileobj, parent_axis)
if (present(parent_axis)) then
select type(parent_axis)
type is (fmsDiagFullAxis_type)
call write_data(fms2io_fileobj, this%subaxis_name, parent_axis%axis_data(i:j))
! Added this select type so the the data is indexed correctly
select type (vardata => parent_axis%axis_data)
type is (real(kind=r8_kind))
call write_data(fms2io_fileobj, this%subaxis_name, vardata(i:j))
type is (real(kind=r4_kind))
call write_data(fms2io_fileobj, this%subaxis_name, vardata(i:j))
end select
end select
endif
type is (fmsDiagDiurnalAxis_type)
Expand Down Expand Up @@ -881,8 +887,10 @@ subroutine fill_subaxis(this, starting_index, ending_index, axis_id, parent_id,
if (present(nz_subaxis)) nsubaxis = nz_subaxis

this%axis_id = axis_id
this%starting_index = starting_index
this%ending_index = ending_index

! The min and max were added here to support axis that are both increasing and decreasing
this%starting_index = min(starting_index, ending_index)
this%ending_index = max(starting_index, ending_index)
this%parent_axis_id = parent_id
write(nsubaxis_char, '(i2.2)') nsubaxis
this%subaxis_name = trim(parent_axis_name)//"_sub"//nsubaxis_char
Expand Down
49 changes: 44 additions & 5 deletions test_fms/diag_manager/test_multiple_zbounds.F90
Original file line number Diff line number Diff line change
Expand Up @@ -29,9 +29,12 @@ program test_multiple_zbounds
type(time_type) :: Time_step !< Time_step of the simulation
integer :: id_z1 !< Axis id for the z dimension
integer :: id_z2 !< Axis id for the z dimension
integer :: id_z_reverse !< Axis id for the z dimension (decreasing)
integer :: id_var1 !< var id for the first variable
integer :: id_var2 !< var_id for the second variable
integer :: id_var3 !< var_id for the third variable
real, allocatable :: z(:) !< z axis data
real, allocatable :: z_reverse(:) !< z axis data (decreasing)
integer :: nz !< Size of the z dimension
integer :: i
logical :: used !< Dummy argument to send_data
Expand All @@ -42,25 +45,28 @@ program test_multiple_zbounds

nz = 10
allocate(z(nz))
allocate(z_reverse(nz))
do i=1, nz
z(i) = i
z_reverse(i) = nz - i + 1
enddo

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

id_z1 = diag_axis_init('zaxis1', z, 'z', 'z', long_name='Z1')
id_z2 = diag_axis_init('zaxis2', z, 'z', 'z', long_name='Z2')
id_z_reverse = diag_axis_init('zaxis3', z_reverse, 'z_reverse', 'z', long_name='Z3')
id_var1 = register_diag_field ('atmos', 'ua_1', (/id_z1/), Time)
id_var2 = register_diag_field ('atmos', 'ua_2', (/id_z2/), Time)
id_var3 = register_diag_field ('atmos', 'ua_3', (/id_z_reverse/), Time)

call diag_manager_set_time_end(set_date(2,1,2,0,0,0))
do i = 1, 24
Time = Time + Time_step
z = real(i)
used = send_data(id_var1, z, Time)
used = send_data(id_var2, z, Time)

used = send_data(id_var3, z, Time)
call diag_send_complete(Time_step)
end do

Expand All @@ -77,15 +83,25 @@ subroutine check_output()
integer :: SUB1_SIZE = 3
integer :: SUB2_SIZE = 1
character(len=20) :: EXPECTED_DIM_NAMES(2)
real :: EXPECTED_ZSUBAXIS_1(3)
real :: EXPECTED_ZSUBAXIS_2(1)
real :: EXPECTED_ZSUBAXIS_3(3)

EXPECTED_DIM_NAMES(2) = "time"
if (.not. open_file(fileobj, "test_multiple_zbounds.nc", "read")) then
call mpp_error(FATAL, "Unable to open the expected output file: test_var_masks.nc")
endif

call check_dimension(fileobj, "time", EXPECTED_NTIMES)
call check_dimension(fileobj, "zaxis1_sub01", SUB1_SIZE)
call check_dimension(fileobj, "zaxis2_sub02", SUB2_SIZE)

EXPECTED_ZSUBAXIS_1 = (/3., 4., 5./)
call check_dimension(fileobj, "zaxis1_sub01", SUB1_SIZE, EXPECTED_ZSUBAXIS_1)

EXPECTED_ZSUBAXIS_2 = (/1./)
call check_dimension(fileobj, "zaxis2_sub02", SUB2_SIZE, EXPECTED_ZSUBAXIS_2)

EXPECTED_ZSUBAXIS_3 = (/5., 4., 3./)
call check_dimension(fileobj, "zaxis3_sub03", SUB1_SIZE, EXPECTED_ZSUBAXIS_3)

EXPECTED_DIM_NAMES(1) = "zaxis1_sub01"
call check_variable(fileobj, "ua_1", EXPECTED_DIM_NAMES)
Expand Down Expand Up @@ -114,17 +130,40 @@ subroutine check_variable(fileobj, variable_name, expected_dimnames)

end subroutine check_variable

subroutine check_dimension(fileobj, dimension_name, expected_size)
!> @brief Check dimension data
subroutine check_data(err_msg, actual_data, expected_data)
character(len=*), intent(in) :: err_msg !< Error message to append
real, intent(in) :: actual_data(:) !< Dimension data from file
real, intent(in) :: expected_data(:) !< Expected data

integer :: i

do i = 1, size(actual_data)
if (actual_data(i) .ne. expected_data(i)) &
call mpp_error(FATAL, "The data is not expected for "//trim(err_msg))
enddo
end subroutine check_data

subroutine check_dimension(fileobj, dimension_name, expected_size, expected_data)
type(FmsNetcdfFile_t), intent(in) :: fileobj
character(len=*), intent(in) :: dimension_name
integer, intent(in) :: expected_size
real, optional, intent(in) :: expected_data(:)

integer :: dim_size
real, allocatable :: z_data(:)

call get_dimension_size(fileobj, dimension_name, dim_size)
if (dim_size .ne. expected_size) then
call mpp_error(FATAL, trim(dimension_name)//" is not the expected size!")
endif

if (present(expected_data)) then
allocate(z_data(dim_size))
call read_data(fileobj, dimension_name, z_data)
call check_data(dimension_name, z_data, expected_data)
endif


end subroutine
end program test_multiple_zbounds
2 changes: 2 additions & 0 deletions test_fms/diag_manager/test_multiple_zbounds.sh
Original file line number Diff line number Diff line change
Expand Up @@ -47,6 +47,8 @@ diag_files:
zbounds: 3 5
- var_name: ua_2
zbounds: 1 1
- var_name: ua_3
zbounds: 3 5
_EOF

test_expect_success "Test with multiple zbounds limits (modern diag manager)" '
Expand Down
Loading