From bfe4cad91008781065863b2e2ead9a1d9a3240be Mon Sep 17 00:00:00 2001 From: uramirez8707 <49168881+uramirez8707@users.noreply.github.com> Date: Tue, 1 Aug 2023 13:18:41 -0400 Subject: [PATCH] Write buffer prep (#1311) * Adds the buffer ids to the diag file object in the register call * add functionality to get the size of the subaxis * Add the field id and the yaml id to the buffer object (with getter and setter functions) --- diag_manager/fms_diag_axis_object.F90 | 11 +++++++ diag_manager/fms_diag_file_object.F90 | 17 +++++++++- diag_manager/fms_diag_object.F90 | 10 ++++++ diag_manager/fms_diag_output_buffer.F90 | 43 +++++++++++++++++++++++++ 4 files changed, 80 insertions(+), 1 deletion(-) diff --git a/diag_manager/fms_diag_axis_object.F90 b/diag_manager/fms_diag_axis_object.F90 index 61555b52e6..8ae2a325b9 100644 --- a/diag_manager/fms_diag_axis_object.F90 +++ b/diag_manager/fms_diag_axis_object.F90 @@ -118,6 +118,7 @@ module fms_diag_axis_object_mod real(kind=r4_kind), allocatable, private :: zbounds(:) !< Bounds of the Z axis contains procedure :: fill_subaxis + procedure :: axis_length END TYPE fmsDiagSubAxis_type !> @brief Type to hold the diurnal axis @@ -755,6 +756,16 @@ subroutine fill_subaxis(this, starting_index, ending_index, axis_id, parent_id, endif end subroutine fill_subaxis + !> @brief Get the axis length of a subaxis + !> @return the axis length + function axis_length(this) & + result(res) + class(fmsDiagSubAxis_type) , INTENT(IN) :: this !< diag_sub_axis obj + integer :: res + + res = this%ending_index - this%starting_index + 1 + end function + !> @brief Get the ntiles in a domain !> @return the number of tiles in a domain function get_ntiles(this) & diff --git a/diag_manager/fms_diag_file_object.F90 b/diag_manager/fms_diag_file_object.F90 index 62f01f7fe0..51a3396eb5 100644 --- a/diag_manager/fms_diag_file_object.F90 +++ b/diag_manager/fms_diag_file_object.F90 @@ -88,14 +88,16 @@ module fms_diag_file_object_mod integer, allocatable :: num_registered_fields !< The number of fields registered !! to the file integer, dimension(:), allocatable :: axis_ids !< Array of axis ids in the file - integer, dimension(:), allocatable :: buffer_ids !< array of buffer ids associated with the file integer :: number_of_axis !< Number of axis in the file + integer, dimension(:), allocatable :: buffer_ids !< array of buffer ids associated with the file + integer :: number_of_buffers !< Number of buffers that have been added to the file logical :: time_ops !< .True. if file contains variables that are time_min, time_max, time_average or time_sum integer :: unlim_dimension_level !< The unlimited dimension level currently being written logical :: is_static !< .True. if the frequency is -1 contains procedure, public :: add_field_and_yaml_id + procedure, public :: add_buffer_id procedure, public :: is_field_registered procedure, public :: init_diurnal_axis procedure, public :: has_file_metadata_from_model @@ -210,13 +212,16 @@ logical function fms_diag_files_object_init (files_array) obj%diag_yaml_file => diag_yaml%diag_files(i) obj%id = i allocate(obj%field_ids(diag_yaml%diag_files(i)%size_file_varlist())) + allocate(obj%buffer_ids(diag_yaml%diag_files(i)%size_file_varlist())) allocate(obj%yaml_ids(diag_yaml%diag_files(i)%size_file_varlist())) allocate(obj%field_registered(diag_yaml%diag_files(i)%size_file_varlist())) !! Initialize the integer arrays obj%field_ids = DIAG_NOT_REGISTERED obj%yaml_ids = DIAG_NOT_REGISTERED + obj%buffer_ids = DIAG_NOT_REGISTERED obj%field_registered = .FALSE. obj%num_registered_fields = 0 + obj%number_of_buffers = 0 !> These will be set in a set_file_domain obj%type_of_domain = NO_DOMAIN @@ -287,6 +292,16 @@ subroutine add_field_and_yaml_id (this, new_field_id, yaml_id) endif end subroutine add_field_and_yaml_id +!> \brief Adds a buffer_id to the file object +subroutine add_buffer_id (this, buffer_id) + class(fmsDiagFile_type), intent(inout) :: this !< The file object + integer, intent(in) :: buffer_id !< Buffer id to add to the file + + this%number_of_buffers = this%number_of_buffers + 1 + this%buffer_ids(this%number_of_buffers) = buffer_id + +end subroutine add_buffer_id + !> \brief Initializes a diurnal axis for a fileobj !! \note This is going to be called for every variable in the file, if the variable is not a diurnal variable !! it will do nothing. It only defined a diurnal axis once. diff --git a/diag_manager/fms_diag_object.F90 b/diag_manager/fms_diag_object.F90 index 76015ae4c4..08d324f615 100644 --- a/diag_manager/fms_diag_object.F90 +++ b/diag_manager/fms_diag_object.F90 @@ -225,6 +225,10 @@ integer function fms_register_diag_field_obj & !> Initialize buffer_ids of this field with the diag_field_indices(diag_field_indices) !! of the sorted variable list 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)) + enddo !> Allocate and initialize member buffer_allocated of this field allocate(fieldptr%buffer_allocated(size(diag_field_indices))) @@ -244,6 +248,7 @@ integer function fms_register_diag_field_obj & do i = 1, size(file_ids) fileptr => this%FMS_diag_files(file_ids(i))%FMS_diag_file call fileptr%add_field_and_yaml_id(fieldptr%get_id(), diag_field_indices(i)) + call fileptr%add_buffer_id(fieldptr%buffer_ids(i)) call fileptr%set_file_domain(fieldptr%get_domain(), fieldptr%get_type_of_domain()) call fileptr%init_diurnal_axis(this%diag_axis, this%registered_axis, diag_field_indices(i)) call fileptr%add_axes(axes, this%diag_axis, this%registered_axis, diag_field_indices(i), & @@ -255,6 +260,7 @@ integer function fms_register_diag_field_obj & do i = 1, size(file_ids) fileptr => this%FMS_diag_files(file_ids(i))%FMS_diag_file call fileptr%add_field_and_yaml_id(fieldptr%get_id(), diag_field_indices(i)) + call fileptr%add_buffer_id(fieldptr%buffer_ids(i)) call fileptr%init_diurnal_axis(this%diag_axis, this%registered_axis, diag_field_indices(i)) call fileptr%set_file_domain(fieldptr%get_domain(), fieldptr%get_type_of_domain()) call fileptr%add_axes(axes, this%diag_axis, this%registered_axis, diag_field_indices(i), & @@ -265,6 +271,7 @@ integer function fms_register_diag_field_obj & do i = 1, size(file_ids) fileptr => this%FMS_diag_files(file_ids(i))%FMS_diag_file call fileptr%add_field_and_yaml_id(fieldptr%get_id(), diag_field_indices(i)) + call fileptr%add_buffer_id(fieldptr%buffer_ids(i)) call fileptr%add_start_time(init_time, this%current_model_time) call fileptr%set_file_time_ops (fieldptr%diag_field(i), fieldptr%is_static()) enddo @@ -272,6 +279,7 @@ integer function fms_register_diag_field_obj & do i = 1, size(file_ids) fileptr => this%FMS_diag_files(file_ids(i))%FMS_diag_file call fileptr%add_field_and_yaml_id(fieldptr%get_id(), diag_field_indices(i)) + call fileptr%add_buffer_id(fieldptr%buffer_ids(i)) call fileptr%set_file_time_ops (fieldptr%diag_field(i), fieldptr%is_static()) enddo endif @@ -891,6 +899,8 @@ integer function fms_get_axis_length(this, axis_id) select type (axis => this%diag_axis(axis_id)%axis) type is (fmsDiagFullAxis_type) fms_get_axis_length = axis%axis_length() + type is (fmsDiagSubAxis_type) + fms_get_axis_length = axis%axis_length() end select #endif end function fms_get_axis_length diff --git a/diag_manager/fms_diag_output_buffer.F90 b/diag_manager/fms_diag_output_buffer.F90 index 5b3d267dfd..c2a29fcedd 100644 --- a/diag_manager/fms_diag_output_buffer.F90 +++ b/diag_manager/fms_diag_output_buffer.F90 @@ -62,10 +62,16 @@ module fms_diag_output_buffer_mod type :: fmsDiagOutputBufferContainer_type class(fmsDiagOutputBuffer_class), allocatable :: diag_buffer_obj !< any 0-5d buffer object integer, allocatable :: axis_ids(:) !< Axis ids for the buffer + integer :: field_id !< The id of the field the buffer belongs to + integer :: yaml_id !< The id of the yaml id the buffer belongs to contains procedure :: add_axis_ids procedure :: get_axis_ids + procedure :: set_field_id + procedure :: get_field_id + procedure :: set_yaml_id + procedure :: get_yaml_id end type !> Scalar buffer type to extend fmsDiagBufferContainer_type @@ -1455,5 +1461,42 @@ function get_axis_ids(this) & endif end function +!> @brief Get the field id of the buffer +!! @return the field id of the buffer +function get_field_id(this) & + result(res) + + class(fmsDiagOutputBufferContainer_type), intent(in) :: this !< Buffer object + integer :: res + + res = this%field_id +end function get_field_id + +!> @brief set the field id of the buffer +subroutine set_field_id(this, field_id) + class(fmsDiagOutputBufferContainer_type), intent(inout) :: this !< Buffer object + integer, intent(in) :: field_id !< field id of the buffer + + this%field_id = field_id +end subroutine set_field_id + +!> @brief set the field id of the buffer +subroutine set_yaml_id(this, yaml_id) + class(fmsDiagOutputBufferContainer_type), intent(inout) :: this !< Buffer object + integer, intent(in) :: yaml_id !< yaml id of the buffer + + this%yaml_id = yaml_id +end subroutine set_yaml_id + +!> @brief Get the yaml id of the buffer +!! @return the yaml id of the buffer +function get_yaml_id(this) & + result(res) + + class(fmsDiagOutputBufferContainer_type), intent(in) :: this !< Buffer object + integer :: res + + res = this%yaml_id +end function get_yaml_id #endif end module fms_diag_output_buffer_mod