Skip to content

Commit

Permalink
modern_diag: add the write_data calls (NOAA-GFDL#1320)
Browse files Browse the repository at this point in the history
* add the write_data calls

* fix the test_flexible_time after the update
  • Loading branch information
uramirez8707 authored and rem1776 committed May 1, 2024
1 parent 6a84dd3 commit bd1544a
Show file tree
Hide file tree
Showing 5 changed files with 164 additions and 11 deletions.
2 changes: 1 addition & 1 deletion diag_manager/Makefile.am
Original file line number Diff line number Diff line change
Expand Up @@ -88,7 +88,7 @@ diag_manager_mod.$(FC_MODEXT): diag_axis_mod.$(FC_MODEXT) diag_data_mod.$(FC_MOD
fms_diag_object_container_mod.$(FC_MODEXT) fms_diag_axis_object_mod.$(FC_MODEXT) \
fms_diag_time_reduction_mod.$(FC_MODEXT) fms_diag_outfield_mod.$(FC_MODEXT) \
fms_diag_fieldbuff_update_mod.$(FC_MODEXT)
fms_diag_output_buffer_mod.$(FC_MODEXT): diag_data_mod.$(FC_MODEXT)
fms_diag_output_buffer_mod.$(FC_MODEXT): diag_data_mod.$(FC_MODEXT) fms_diag_yaml_mod.$(FC_MODEXT)

# Mod files are built and then installed as headers.
MODFILES = \
Expand Down
38 changes: 37 additions & 1 deletion diag_manager/fms_diag_file_object.F90
Original file line number Diff line number Diff line change
Expand Up @@ -45,7 +45,7 @@ module fms_diag_file_object_mod
fmsDiagFullAxis_type, define_subaxis, define_diurnal_axis, &
fmsDiagDiurnalAxis_type, create_new_z_subaxis
use fms_diag_field_object_mod, only: fmsDiagField_type
use fms_diag_output_buffer_mod, only: fmsDiagOutputBufferContainer_type
use fms_diag_output_buffer_mod, only: fmsDiagOutputBufferContainer_type, fmsDiagOutputBuffer_class
use mpp_mod, only: mpp_get_current_pelist, mpp_npes, mpp_root_pe, mpp_pe, mpp_error, FATAL, stdout, &
uppercase, lowercase

Expand Down Expand Up @@ -161,6 +161,7 @@ module fms_diag_file_object_mod
procedure :: open_diag_file
procedure :: write_global_metadata
procedure :: write_time_metadata
procedure :: write_field_data
procedure :: write_axis_metadata
procedure :: write_field_metadata
procedure :: write_axis_data
Expand Down Expand Up @@ -1119,6 +1120,41 @@ subroutine write_time_metadata(this)

end subroutine write_time_metadata

!> \brief Write out the field data to the file
subroutine write_field_data(this, field_obj, buffer_obj)
class(fmsDiagFileContainer_type), intent(in), target :: this !< The diag file object to write to
type(fmsDiagField_type), intent(in), target :: field_obj(:) !< The field object to write from
type(fmsDiagOutputBufferContainer_type), intent(in), target :: buffer_obj(:) !< The buffer object with the data

class(fmsDiagFile_type), pointer :: diag_file !< Diag_file object to open
class(FmsNetcdfFile_t), pointer :: fileobj !< Fileobj to write to
integer :: i !< For do loops
integer :: field_id !< The id of the field writing the data from

diag_file => this%FMS_diag_file
fileobj => diag_file%fileobj

!TODO This may be offloaded in the future
if (diag_file%is_static) then
!< Here the file is static so there is no need for the unlimited dimension
!! as a variables are static
do i = 1, diag_file%number_of_buffers
call buffer_obj(diag_file%buffer_ids(i))%write_buffer(fileobj)
enddo
else
do i = 1, diag_file%number_of_buffers
field_id = buffer_obj(diag_file%buffer_ids(i))%get_field_id()
if (field_obj(field_id)%is_static()) then
!< If the variable is static, only write it the first time
if (diag_file%unlim_dimension_level .eq. 1) call buffer_obj(diag_file%buffer_ids(i))%write_buffer(fileobj)
else
call buffer_obj(diag_file%buffer_ids(i))%write_buffer(fileobj, unlim_dim_level=diag_file%unlim_dimension_level)
endif
enddo
endif

end subroutine write_field_data

!> \brief Determine if it is time to close the file
!! \return .True. if it is time to close the file
logical function is_time_to_close_file (this, time_step)
Expand Down
4 changes: 2 additions & 2 deletions diag_manager/fms_diag_object.F90
Original file line number Diff line number Diff line change
Expand Up @@ -227,7 +227,7 @@ integer function fms_register_diag_field_obj &
fieldptr%buffer_ids = get_diag_field_ids(diag_field_indices)
do i = 1, size(fieldptr%buffer_ids)
call this%FMS_diag_output_buffers(fieldptr%buffer_ids(i))%set_field_id(this%registered_variables)
call this%FMS_diag_output_buffers(fieldptr%buffer_ids(i))%set_yaml_id(diag_field_indices(i))
call this%FMS_diag_output_buffers(fieldptr%buffer_ids(i))%set_yaml_id(fieldptr%buffer_ids(i))
enddo

!> Allocate and initialize member buffer_allocated of this field
Expand Down Expand Up @@ -719,7 +719,7 @@ subroutine fms_diag_do_io(this, is_end_of_run)
if (diag_file%is_time_to_write(model_time)) then
call diag_file%increase_unlim_dimension_level()
call diag_file%write_time_data()
!TODO call diag_file%add_variable_data()
call diag_file%write_field_data(this%FMS_diag_fields, this%FMS_diag_output_buffers)
call diag_file%update_next_write(model_time)
call diag_file%update_current_new_file_freq_index(model_time)
if (diag_file%is_time_to_close_file(model_time)) call diag_file%close_diag_file()
Expand Down
113 changes: 111 additions & 2 deletions diag_manager/fms_diag_output_buffer.F90
Original file line number Diff line number Diff line change
Expand Up @@ -24,18 +24,19 @@
!! buffer0-5d types extend fmsDiagBuffer_class, and upon allocation
!! are added to the module's buffer_lists depending on it's dimension
module fms_diag_output_buffer_mod

#ifdef use_yaml
use platform_mod
use iso_c_binding
use time_manager_mod, only: time_type
use mpp_mod, only: mpp_error, FATAL
use diag_data_mod, only: DIAG_NULL, DIAG_NOT_REGISTERED, i4, i8, r4, r8
use fms2_io_mod, only: FmsNetcdfFile_t, write_data, FmsNetcdfDomainFile_t, FmsNetcdfUnstructuredDomainFile_t
use fms_diag_yaml_mod, only: diag_yaml

implicit none

private

#ifdef use_yaml
!> @brief Object that holds buffered data and other diagnostics
!! Abstract to ensure use through its extensions(buffer0-5d types)
type, abstract :: fmsDiagOutputBuffer_class
Expand Down Expand Up @@ -72,6 +73,11 @@ module fms_diag_output_buffer_mod
procedure :: get_field_id
procedure :: set_yaml_id
procedure :: get_yaml_id
procedure :: write_buffer
!! These are needed because otherwise the write_data calls will go into the wrong interface
procedure :: write_buffer_wrapper_netcdf
procedure :: write_buffer_wrapper_domain
procedure :: write_buffer_wrapper_u
end type

!> Scalar buffer type to extend fmsDiagBufferContainer_type
Expand Down Expand Up @@ -1497,5 +1503,108 @@ function get_yaml_id(this) &

res = this%yaml_id
end function get_yaml_id

!> @brief Write the buffer to the file
subroutine write_buffer(this, fileobj, unlim_dim_level)
class(fmsDiagOutputBufferContainer_type), intent(in) :: this !< buffer object to write
class(FmsNetcdfFile_t), intent(in) :: fileobj !< fileobj to write to
integer, optional, intent(in) :: unlim_dim_level !< unlimited dimension

select type(fileobj)
type is (FmsNetcdfFile_t)
call this%write_buffer_wrapper_netcdf(fileobj, unlim_dim_level=unlim_dim_level)
type is (FmsNetcdfDomainFile_t)
call this%write_buffer_wrapper_domain(fileobj, unlim_dim_level=unlim_dim_level)
type is (FmsNetcdfUnstructuredDomainFile_t)
call this%write_buffer_wrapper_u(fileobj, unlim_dim_level=unlim_dim_level)
class default
call mpp_error(FATAL, "The file "//trim(fileobj%path)//" is not one of the accepted types"//&
" only FmsNetcdfFile_t, FmsNetcdfDomainFile_t, and FmsNetcdfUnstructuredDomainFile_t are accepted.")
end select
end subroutine write_buffer

!> @brief Write the buffer to the FmsNetcdfFile_t fileobj
subroutine write_buffer_wrapper_netcdf(this, fileobj, unlim_dim_level)
class(fmsDiagOutputBufferContainer_type), intent(in) :: this !< buffer object to write
type(FmsNetcdfFile_t), intent(in) :: fileobj !< fileobj to write to
integer, optional, intent(in) :: unlim_dim_level !< unlimited dimension

character(len=:), allocatable :: varname !< name of the variable

varname = diag_yaml%diag_fields(this%yaml_id)%get_var_outname()
select type(buffer_obj=>this%diag_buffer_obj)
type is (outputBuffer0d_type)
call write_data(fileobj, varname, buffer_obj%buffer(1), unlim_dim_level=unlim_dim_level)
type is (outputBuffer1d_type)
call write_data(fileobj, varname, buffer_obj%buffer, unlim_dim_level=unlim_dim_level)
type is (outputBuffer2d_type)
call write_data(fileobj, varname, buffer_obj%buffer, unlim_dim_level=unlim_dim_level)
type is (outputBuffer3d_type)
call write_data(fileobj, varname, buffer_obj%buffer, unlim_dim_level=unlim_dim_level)
type is (outputBuffer4d_type)
call write_data(fileobj, varname, buffer_obj%buffer, unlim_dim_level=unlim_dim_level)
type is (outputBuffer5d_type)
call write_data(fileobj, varname, buffer_obj%buffer, unlim_dim_level=unlim_dim_level)
class default
call mpp_error(FATAL, "The field:"//trim(varname)//" does not have a valid buffer object type."//&
" Only 0d, 1d, 2d, 3d, 4d, and 5d buffers are supported.")
end select
end subroutine write_buffer_wrapper_netcdf

!> @brief Write the buffer to the FmsNetcdfDomainFile_t fileobj
subroutine write_buffer_wrapper_domain(this, fileobj, unlim_dim_level)
class(fmsDiagOutputBufferContainer_type), intent(in) :: this !< buffer object to write
type(FmsNetcdfDomainFile_t), intent(in) :: fileobj !< fileobj to write to
integer, optional, intent(in) :: unlim_dim_level !< unlimited dimension

character(len=:), allocatable :: varname !< name of the variable

varname = diag_yaml%diag_fields(this%yaml_id)%get_var_outname()
select type(buffer_obj=>this%diag_buffer_obj)
type is (outputBuffer0d_type)
call write_data(fileobj, varname, buffer_obj%buffer(1), unlim_dim_level=unlim_dim_level)
type is (outputBuffer1d_type)
call write_data(fileobj, varname, buffer_obj%buffer, unlim_dim_level=unlim_dim_level)
type is (outputBuffer2d_type)
call write_data(fileobj, varname, buffer_obj%buffer, unlim_dim_level=unlim_dim_level)
type is (outputBuffer3d_type)
call write_data(fileobj, varname, buffer_obj%buffer, unlim_dim_level=unlim_dim_level)
type is (outputBuffer4d_type)
call write_data(fileobj, varname, buffer_obj%buffer, unlim_dim_level=unlim_dim_level)
type is (outputBuffer5d_type)
call write_data(fileobj, varname, buffer_obj%buffer, unlim_dim_level=unlim_dim_level)
class default
call mpp_error(FATAL, "The field:"//trim(varname)//" does not have a valid buffer object type."//&
" Only 0d, 1d, 2d, 3d, 4d, and 5d buffers are supported.")
end select
end subroutine write_buffer_wrapper_domain

!> @brief Write the buffer to the FmsNetcdfUnstructuredDomainFile_t fileobj
subroutine write_buffer_wrapper_u(this, fileobj, unlim_dim_level)
class(fmsDiagOutputBufferContainer_type), intent(in) :: this !< buffer object to write
type(FmsNetcdfUnstructuredDomainFile_t), intent(in) :: fileobj !< fileobj to write to
integer, optional, intent(in) :: unlim_dim_level !< unlimited dimension

character(len=:), allocatable :: varname !< name of the variable

varname = diag_yaml%diag_fields(this%yaml_id)%get_var_outname()
select type(buffer_obj=>this%diag_buffer_obj)
type is (outputBuffer0d_type)
call write_data(fileobj, varname, buffer_obj%buffer(1), unlim_dim_level=unlim_dim_level)
type is (outputBuffer1d_type)
call write_data(fileobj, varname, buffer_obj%buffer, unlim_dim_level=unlim_dim_level)
type is (outputBuffer2d_type)
call write_data(fileobj, varname, buffer_obj%buffer, unlim_dim_level=unlim_dim_level)
type is (outputBuffer3d_type)
call write_data(fileobj, varname, buffer_obj%buffer, unlim_dim_level=unlim_dim_level)
type is (outputBuffer4d_type)
call write_data(fileobj, varname, buffer_obj%buffer, unlim_dim_level=unlim_dim_level)
type is (outputBuffer5d_type)
call write_data(fileobj, varname, buffer_obj%buffer, unlim_dim_level=unlim_dim_level)
class default
call mpp_error(FATAL, "The field:"//trim(varname)//" does not have a valid buffer object type."//&
" Only 0d, 1d, 2d, 3d, 4d, and 5d buffers are supported.")
end select
end subroutine write_buffer_wrapper_u
#endif
end module fms_diag_output_buffer_mod
18 changes: 13 additions & 5 deletions test_fms/diag_manager/test_flexible_time.F90
Original file line number Diff line number Diff line change
Expand Up @@ -21,15 +21,19 @@
program test_flexible_time
use fms_mod, only: fms_init, fms_end
use time_manager_mod, only: set_date, time_type, increment_date, set_calendar_type, &
JULIAN, set_time
JULIAN, set_time, operator(+)
use diag_manager_mod, only: diag_manager_init, diag_axis_init, register_diag_field, &
diag_manager_set_time_end, diag_send_complete, diag_manager_end
diag_manager_set_time_end, diag_send_complete, diag_manager_end, &
send_data
use mpp_mod, only: FATAL, mpp_error
use platform_mod, only: r8_kind

implicit none

real(kind=r8_kind) :: var_data(2) !< Dummy data
logical :: used !< .True. if send_data was sucessful
type(time_type) :: Time !< Time of the simulation
type(time_type) :: Start_Time !< Start time of the simulation
type(time_type) :: Time_step !< Start time of the simulation
type(time_type) :: End_Time !< End Time of the simulation
integer :: i
integer :: id_z, id_var
Expand All @@ -39,18 +43,22 @@ program test_flexible_time
call diag_manager_init

!< Starting time of the simulation
Start_Time = set_date(2,1,1,3,0,0) !02/01/01 hour 3
Time = set_date(2,1,1,3,0,0) !02/01/01 hour 3

!< Set up a dummy variable
id_z = diag_axis_init('z', (/1. ,2. /), 'point_Z', 'z', long_name='point_Z')
id_var = register_diag_field ('atm_mod', 'var1', (/id_z/), Start_Time, 'Var not domain decomposed', 'mullions')
id_var = register_diag_field ('atm_mod', 'var1', (/id_z/), Time, 'Var not domain decomposed', 'mullions')

!< Set up the end of the simulation (i.e 2 days long)
End_Time = set_date(2,1,3,3,0,0)
call diag_manager_set_time_end(End_Time)

!< Set up the simulation
Time_step = set_time (3600,0) !< 1 hour
do i=1,48
var_data = real(i, kind=r8_kind)
Time = Time + Time_step
used = send_data(id_var, var_data, Time)
call diag_send_complete(set_time(3600,0))
enddo

Expand Down

0 comments on commit bd1544a

Please sign in to comment.