FMS  2024.03
Flexible Modeling System
fms_diag_object.F90
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 module fms_diag_object_mod
20 use mpp_mod, only: fatal, note, warning, mpp_error, mpp_pe, mpp_root_pe, stdout
21 use diag_data_mod, only: diag_null, diag_not_found, diag_not_registered, diag_registered_id, &
22  &diag_field_not_found, diag_not_registered, max_axes, two_d_domain, &
23  &get_base_time, null_axis_id, get_var_type, diag_not_registered, &
26 
27  USE time_manager_mod, ONLY: set_time, set_date, get_time, time_type, OPERATOR(>=), OPERATOR(>),&
28  & OPERATOR(<), OPERATOR(==), OPERATOR(/=), OPERATOR(/), OPERATOR(+), ASSIGNMENT(=), get_date, &
30 #ifdef use_yaml
31 use fms_diag_file_object_mod, only: fmsdiagfilecontainer_type, fmsdiagfile_type, fms_diag_files_object_init
32 use fms_diag_field_object_mod, only: fmsdiagfield_type, fms_diag_fields_object_init, get_default_missing_value, &
33  check_for_slices
34 use fms_diag_yaml_mod, only: diag_yaml_object_init, diag_yaml_object_end, find_diag_field, &
36 use fms_diag_axis_object_mod, only: fms_diag_axis_object_init, fmsdiagaxis_type, fmsdiagsubaxis_type, &
41 use fms_mod, only: fms_error_handler
42 use fms_diag_reduction_methods_mod, only: check_indices_order, init_mask, set_weight
43 use constants_mod, only: seconds_per_day
44 #endif
45 USE fms_diag_bbox_mod, ONLY: fmsdiagibounds_type, determine_if_block_is_in_region
46 #if defined(_OPENMP)
47 use omp_lib
48 #endif
49 use mpp_domains_mod, only: domain1d, domain2d, domainug, null_domain2d
50 use fms_string_utils_mod, only: string
51 use platform_mod
52 implicit none
53 private
54 
56 !TODO add container arrays
57 #ifdef use_yaml
58 private
59 !TODO: Remove FMS prefix from variables in this type
60  class(fmsdiagfilecontainer_type), allocatable :: fms_diag_files (:) !< array of diag files
61  type(fmsdiagfield_type), allocatable :: fms_diag_fields(:) !< Array of diag fields
62  type(fmsdiagoutputbuffer_type), allocatable :: fms_diag_output_buffers(:) !< array of output buffer objects
63  !! one for each variable in the diag_table.yaml
64  integer, private :: registered_buffers = 0 !< number of registered buffers, per dimension
65  class(fmsdiagaxiscontainer_type), allocatable :: diag_axis(:) !< Array of diag_axis
66  integer, private :: registered_variables !< Number of registered variables
67  integer, private :: registered_axis !< Number of registered axis
68  logical, private :: initialized=.false. !< True if the fmsDiagObject is initialized
69  logical, private :: files_initialized=.false. !< True if the fmsDiagObject is initialized
70  logical, private :: fields_initialized=.false. !< True if the fmsDiagObject is initialized
71  logical, private :: buffers_initialized=.false. !< True if the fmsDiagObject is initialized
72  logical, private :: axes_initialized=.false. !< True if the fmsDiagObject is initialized
73 #endif
74  contains
75  procedure :: init => fms_diag_object_init
76  procedure :: diag_end => fms_diag_object_end
77  procedure :: fms_register_diag_field_scalar
78  procedure :: fms_register_diag_field_array
79  procedure :: fms_register_static_field
80  procedure :: fms_diag_axis_init
81  procedure :: register => fms_register_diag_field_obj !! Merely initialize fields.
82  procedure :: fms_diag_field_add_attribute
83  procedure :: fms_diag_axis_add_attribute
84  procedure :: fms_get_domain2d
85  procedure :: fms_get_axis_length
86  procedure :: fms_get_diag_field_id_from_name
87  procedure :: fms_get_field_name_from_id
88  procedure :: fms_get_axis_name_from_id
89  procedure :: fms_diag_accept_data
90  procedure :: fms_diag_send_complete
91  procedure :: do_buffer_math
92  procedure :: fms_diag_do_io
93  procedure :: fms_diag_do_reduction
94  procedure :: fms_diag_field_add_cell_measures
95  procedure :: allocate_diag_field_output_buffers
96  procedure :: fms_diag_compare_window
97 #ifdef use_yaml
98  procedure :: get_diag_buffer
99 #endif
100 end type fmsdiagobject_type
101 
102 type (fmsdiagobject_type), target :: fms_diag_object
103 
104 public :: fms_register_diag_field_obj
105 public :: fms_register_diag_field_scalar
106 public :: fms_register_diag_field_array
107 public :: fms_register_static_field
108 public :: fms_diag_field_add_attribute
109 public :: fms_get_diag_field_id_from_name
110 public :: fms_diag_object
111 public :: fmsdiagobject_type
112 integer, private :: registered_variables !< Number of registered variables
113 public :: dump_diag_obj
114 
115 contains
116 
117 !> @brief Initiliazes the fms_diag_object.
118 !! Reads the diag_table.yaml and fills in the yaml object
119 !! Allocates the diag manager object arrays for files, fields, and buffers
120 !! Initializes variables
121 subroutine fms_diag_object_init (this,diag_subset_output, time_init)
122  class(fmsdiagobject_type) :: this !< Diag mediator/controller object
123  integer :: diag_subset_output !< Subset of the diag output?
124  INTEGER, DIMENSION(6), OPTIONAL, INTENT(IN) :: time_init !< Model time diag_manager initialized
125 
126 #ifdef use_yaml
127  if (this%initialized) return
128 
129 ! allocate(diag_objs(get_num_unique_fields()))
130  CALL diag_yaml_object_init(diag_subset_output)
131 
132  !! Doing this here, because the base_time is not set until the yaml is parsed
133  !! if time_init is present, it will be set in diag_manager_init
134  if (.not. present(time_init)) then
136  endif
137 
138  this%axes_initialized = fms_diag_axis_object_init(this%diag_axis)
139  this%files_initialized = fms_diag_files_object_init(this%FMS_diag_files)
140  this%fields_initialized = fms_diag_fields_object_init(this%FMS_diag_fields)
141  this%buffers_initialized =fms_diag_output_buffer_init(this%FMS_diag_output_buffers,SIZE(diag_yaml%get_diag_fields()))
142  this%registered_variables = 0
143  this%registered_axis = 0
144  this%initialized = .true.
145 #else
146  call mpp_error("fms_diag_object_init",&
147  "You must compile with -Duse_yaml to use the option use_modern_diag", fatal)
148 #endif
149 end subroutine fms_diag_object_init
150 
151 !> \description Loops through all files and does one final write.
152 !! Closes all files
153 !! Deallocates all buffers, fields, and files
154 !! Uninitializes the fms_diag_object
155 subroutine fms_diag_object_end (this, time)
156  class(fmsdiagobject_type) :: this
157  TYPE(time_type), INTENT(in) :: time
158 
159  integer :: i
160 #ifdef use_yaml
161  !TODO: loop through files and force write
162  if (.not. this%initialized) return
163 
164  ! write output yaml
165  call fms_diag_yaml_out()
166 
167  call this%do_buffer_math()
168  call this%fms_diag_do_io(end_time=time)
169  !TODO: Deallocate diag object arrays and clean up all memory
170  do i=1, size(this%FMS_diag_output_buffers)
171  call this%FMS_diag_output_buffers(i)%flush_buffer()
172  enddo
173  deallocate(this%FMS_diag_output_buffers)
174  this%axes_initialized = fms_diag_axis_object_end(this%diag_axis)
175  this%initialized = .false.
177 #else
178  call mpp_error(fatal, "You can not call fms_diag_object%end without yaml")
179 #endif
180 end subroutine fms_diag_object_end
181 
182 !> @brief Registers a field.
183 !! @description This to avoid having duplicate code in each of the _scalar, _array and _static register calls
184 !! @return field index to be used in subsequent calls to send_data or DIAG_FIELD_NOT_FOUND if the field is not
185 !! in the diag_table.yaml
186 integer function fms_register_diag_field_obj &
187  (this, modname, varname, axes, init_time, &
188  longname, units, missing_value, varrange, mask_variant, standname, &
189  do_not_log, err_msg, interp_method, tile_count, area, volume, realm, static, &
190  multiple_send_data)
191 
192  class(fmsdiagobject_type),TARGET,INTENT(inout):: this !< Diaj_obj to fill
193  CHARACTER(len=*), INTENT(in) :: modname !< The module name
194  CHARACTER(len=*), INTENT(in) :: varname !< The variable name
195  TYPE(time_type), OPTIONAL, INTENT(in) :: init_time !< Initial time
196  INTEGER, TARGET, OPTIONAL, INTENT(in) :: axes(:) !< The axes indicies
197  CHARACTER(len=*), OPTIONAL, INTENT(in) :: longname !< THe variables long name
198  CHARACTER(len=*), OPTIONAL, INTENT(in) :: units !< The units of the variables
199  CHARACTER(len=*), OPTIONAL, INTENT(in) :: standname !< The variables stanard name
200  class(*), OPTIONAL, INTENT(in) :: missing_value !< Missing value to add as a attribute
201  class(*), OPTIONAL, INTENT(in) :: varrange(2) !< Range to add as a attribute
202  LOGICAL, OPTIONAL, INTENT(in) :: mask_variant !< .True. if mask changes over time
203  LOGICAL, OPTIONAL, INTENT(in) :: do_not_log !< if TRUE, field info is not logged
204  CHARACTER(len=*), OPTIONAL, INTENT(out) :: err_msg !< Error message to be passed back up
205  CHARACTER(len=*), OPTIONAL, INTENT(in) :: interp_method !< The interp method to be used when
206  !! regridding the field in post-processing.
207  !! Valid options are "conserve_order1",
208  !! "conserve_order2", and "none".
209  INTEGER, OPTIONAL, INTENT(in) :: tile_count !< the number of tiles
210  INTEGER, OPTIONAL, INTENT(in) :: area !< diag_field_id of the cell area field
211  INTEGER, OPTIONAL, INTENT(in) :: volume !< diag_field_id of the cell volume field
212  CHARACTER(len=*), OPTIONAL, INTENT(in) :: realm !< String to set as the value to the
213  !! modeling_realm attribute
214  LOGICAL, OPTIONAL, INTENT(in) :: static !< True if the variable is static
215  LOGICAL, OPTIONAL, INTENT(in) :: multiple_send_data !< .True. if send data is called, multiple
216  !! times for the same time
217 
218 #ifdef use_yaml
219 
220  class(fmsdiagfile_type), pointer :: fileptr !< Pointer to the diag_file
221  class(fmsdiagfield_type), pointer :: fieldptr !< Pointer to the diag_field
222  class(fmsdiagoutputbuffer_type), pointer :: bufferptr !< Pointer to the output buffer
223  class(diagyamlfilesvar_type), pointer :: yamlfptr !< Pointer to yaml object to get the reduction method
224  integer, allocatable :: file_ids(:) !< The file IDs for this variable
225  integer :: i !< For do loops
226  integer, allocatable :: diag_field_indices(:) !< indices where the field was found in the yaml
227 #endif
228 #ifndef use_yaml
229 fms_register_diag_field_obj = diag_field_not_found
230 CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
231 #else
232  diag_field_indices = find_diag_field(varname, modname)
233  if (diag_field_indices(1) .eq. diag_null) then
234  !< The field was not found in the table, so return DIAG_FIELD_NOT_FOUND
235  fms_register_diag_field_obj = diag_field_not_found
236  deallocate(diag_field_indices)
237  return
238  endif
239 
240  this%registered_variables = this%registered_variables + 1
241  fms_register_diag_field_obj = this%registered_variables
242 
243  call this%FMS_diag_fields(this%registered_variables)%&
244  &setid(this%registered_variables)
245 
246 !> Use pointers for convenience
247  fieldptr => this%FMS_diag_fields(this%registered_variables)
248 !> Get the file IDs from the field indicies from the yaml
249  file_ids = get_diag_files_id(diag_field_indices)
250  call fieldptr%set_file_ids(file_ids)
251 
252 !> Allocate and initialize member buffer_allocated of this field
253  fieldptr%buffer_allocated = .false.
254  fieldptr%buffer_ids = get_diag_field_ids(diag_field_indices)
255 
256 !> Register the data for the field
257  call fieldptr%register(modname, varname, diag_field_indices, this%diag_axis, &
258  axes=axes, longname=longname, units=units, missing_value=missing_value, varrange= varrange, &
259  mask_variant= mask_variant, standname=standname, do_not_log=do_not_log, err_msg=err_msg, &
260  interp_method=interp_method, tile_count=tile_count, area=area, volume=volume, realm=realm, &
261  static=static, multiple_send_data=multiple_send_data)
262 
263 !> Add the axis information, initial time, and field IDs to the files
264  if (present(axes) .and. present(init_time)) then
265  do i = 1, size(file_ids)
266  fileptr => this%FMS_diag_files(file_ids(i))%FMS_diag_file
267  call fileptr%add_field_and_yaml_id(fieldptr%get_id(), diag_field_indices(i))
268  call fileptr%add_buffer_id(fieldptr%buffer_ids(i))
269  if(fieldptr%get_type_of_domain() .eq. no_domain) then
270  call fileptr%set_file_domain(null(), fieldptr%get_type_of_domain())
271  else
272  call fileptr%set_file_domain(fieldptr%get_domain(), fieldptr%get_type_of_domain())
273  endif
274  call fileptr%init_diurnal_axis(this%diag_axis, this%registered_axis, diag_field_indices(i))
275  call fileptr%add_axes(axes, this%diag_axis, this%registered_axis, diag_field_indices(i), &
276  fieldptr%buffer_ids(i), this%FMS_diag_output_buffers)
277  call fileptr%add_start_time(init_time)
278  call fileptr%set_file_time_ops (fieldptr%diag_field(i), fieldptr%is_static())
279  enddo
280  elseif (present(axes)) then !only axes present
281  do i = 1, size(file_ids)
282  fileptr => this%FMS_diag_files(file_ids(i))%FMS_diag_file
283  call fileptr%add_field_and_yaml_id(fieldptr%get_id(), diag_field_indices(i))
284  call fileptr%add_buffer_id(fieldptr%buffer_ids(i))
285  call fileptr%init_diurnal_axis(this%diag_axis, this%registered_axis, diag_field_indices(i))
286  if(fieldptr%get_type_of_domain() .eq. no_domain) then
287  call fileptr%set_file_domain(null(), fieldptr%get_type_of_domain())
288  else
289  call fileptr%set_file_domain(fieldptr%get_domain(), fieldptr%get_type_of_domain())
290  endif
291  call fileptr%add_axes(axes, this%diag_axis, this%registered_axis, diag_field_indices(i), &
292  fieldptr%buffer_ids(i), this%FMS_diag_output_buffers)
293  call fileptr%set_file_time_ops (fieldptr%diag_field(i), fieldptr%is_static())
294  enddo
295  elseif (present(init_time)) then !only inti time present
296  do i = 1, size(file_ids)
297  fileptr => this%FMS_diag_files(file_ids(i))%FMS_diag_file
298  call fileptr%add_field_and_yaml_id(fieldptr%get_id(), diag_field_indices(i))
299  call fileptr%add_buffer_id(fieldptr%buffer_ids(i))
300  call fileptr%add_start_time(init_time)
301  call fileptr%set_file_time_ops (fieldptr%diag_field(i), fieldptr%is_static())
302  enddo
303  else !no axis or init time present
304  do i = 1, size(file_ids)
305  fileptr => this%FMS_diag_files(file_ids(i))%FMS_diag_file
306  call fileptr%add_field_and_yaml_id(fieldptr%get_id(), diag_field_indices(i))
307  call fileptr%add_buffer_id(fieldptr%buffer_ids(i))
308  call fileptr%set_file_time_ops (fieldptr%diag_field(i), fieldptr%is_static())
309  enddo
310  endif
311 
312  !> Initialize buffer_ids of this field with the diag_field_indices(diag_field_indices)
313 !! of the sorted variable list
314  do i = 1, size(fieldptr%buffer_ids)
315  bufferptr => this%FMS_diag_output_buffers(fieldptr%buffer_ids(i))
316  call bufferptr%set_field_id(this%registered_variables)
317  call bufferptr%set_yaml_id(fieldptr%buffer_ids(i))
318  ! check if diurnal reduction for this buffer and if so set the diurnal sample size
319  yamlfptr => diag_yaml%diag_fields(fieldptr%buffer_ids(i))
320  if( yamlfptr%get_var_reduction() .eq. time_diurnal) then
321  call bufferptr%set_diurnal_sample_size(yamlfptr%get_n_diurnal())
322  endif
323  call bufferptr%init_buffer_time(init_time)
324  call bufferptr%set_next_output(this%FMS_diag_files(file_ids(i))%get_next_output(), &
325  this%FMS_diag_files(file_ids(i))%get_next_next_output(), is_static=fieldptr%is_static())
326  enddo
327 
328  nullify (fileptr)
329  nullify (fieldptr)
330  deallocate(diag_field_indices)
331 #endif
332 end function fms_register_diag_field_obj
333 
334 !> @brief Registers a scalar field
335 !! @return field index to be used in subsequent calls to send_data or DIAG_FIELD_NOT_FOUND if the field is not
336 !! in the diag_table.yaml
337 INTEGER FUNCTION fms_register_diag_field_scalar(this,module_name, field_name, init_time, &
338  & long_name, units, missing_value, var_range, standard_name, do_not_log, err_msg,&
339  & area, volume, realm, multiple_send_data)
340  class(fmsdiagobject_type),TARGET,INTENT(inout):: this !< Diaj_obj to fill
341  CHARACTER(len=*), INTENT(in) :: module_name !< Module where the field comes from
342  CHARACTER(len=*), INTENT(in) :: field_name !< Name of the field
343  TYPE(time_type), OPTIONAL, INTENT(in) :: init_time !< Time to start writing data from
344  CHARACTER(len=*), OPTIONAL, INTENT(in) :: long_name !< Long_name to add as a variable attribute
345  CHARACTER(len=*), OPTIONAL, INTENT(in) :: units !< Units to add as a variable_attribute
346  CHARACTER(len=*), OPTIONAL, INTENT(in) :: standard_name !< Standard_name to name the variable in the file
347  CLASS(*), OPTIONAL, INTENT(in) :: missing_value !< Missing value to add as a variable attribute
348  CLASS(*), OPTIONAL, INTENT(in) :: var_range(:) !< Range to add a variable attribute
349  LOGICAL, OPTIONAL, INTENT(in) :: do_not_log !< If TRUE, field information is not logged
350  CHARACTER(len=*), OPTIONAL, INTENT(out):: err_msg !< Error_msg from call
351  INTEGER, OPTIONAL, INTENT(in) :: area !< Id of the area field
352  INTEGER, OPTIONAL, INTENT(in) :: volume !< Id of the volume field
353  CHARACTER(len=*), OPTIONAL, INTENT(in) :: realm !< String to set as the modeling_realm attribute
354  LOGICAL, OPTIONAL, INTENT(in) :: multiple_send_data !< .True. if send data is called, multiple times
355  !! for the same time
356 
357 #ifndef use_yaml
358 fms_register_diag_field_scalar=diag_field_not_found
359 CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
360 #else
361  fms_register_diag_field_scalar = this%register(&
362  & module_name, field_name, init_time=init_time, &
363  & longname=long_name, units=units, missing_value=missing_value, varrange=var_range, &
364  & standname=standard_name, do_not_log=do_not_log, err_msg=err_msg, &
365  & area=area, volume=volume, realm=realm, multiple_send_data=multiple_send_data)
366 #endif
367 end function fms_register_diag_field_scalar
368 
369 !> @brief Registers an array field
370 !! @return field index to be used in subsequent calls to send_data or DIAG_FIELD_NOT_FOUND if the field is not
371 !! in the diag_table.yaml
372 INTEGER FUNCTION fms_register_diag_field_array(this, module_name, field_name, axes, init_time, &
373  & long_name, units, missing_value, var_range, mask_variant, standard_name, verbose,&
374  & do_not_log, err_msg, interp_method, tile_count, area, volume, realm, &
375  & multiple_send_data)
376  class(fmsdiagobject_type),TARGET,INTENT(inout):: this !< Diaj_obj to fill
377  CHARACTER(len=*), INTENT(in) :: module_name !< Module where the field comes from
378  CHARACTER(len=*), INTENT(in) :: field_name !< Name of the field
379  INTEGER, INTENT(in) :: axes(:) !< Ids corresponding to the variable axis
380  TYPE(time_type), OPTIONAL, INTENT(in) :: init_time !< Time to start writing data from
381  CHARACTER(len=*), OPTIONAL, INTENT(in) :: long_name !< Long_name to add as a variable attribute
382  CHARACTER(len=*), OPTIONAL, INTENT(in) :: units !< Units to add as a variable_attribute
383  CLASS(*), OPTIONAL, INTENT(in) :: missing_value !< Missing value to add as a variable attribute
384  CLASS(*), OPTIONAL, INTENT(in) :: var_range(:) !< Range to add a variable attribute
385  LOGICAL, OPTIONAL, INTENT(in) :: mask_variant !< .True. if mask changes over time
386  CHARACTER(len=*), OPTIONAL, INTENT(in) :: standard_name !< Standard_name to name the variable in the file
387  LOGICAL, OPTIONAL, INTENT(in) :: verbose !< Print more information
388  LOGICAL, OPTIONAL, INTENT(in) :: do_not_log !< If TRUE, field information is not logged
389  CHARACTER(len=*), OPTIONAL, INTENT(out):: err_msg !< Error_msg from call
390  CHARACTER(len=*), OPTIONAL, INTENT(in) :: interp_method !< The interp method to be used when
391  !! regridding the field in post-processing.
392  !! Valid options are "conserve_order1",
393  !! "conserve_order2", and "none".
394  INTEGER, OPTIONAL, INTENT(in) :: tile_count !< The current tile number
395  INTEGER, OPTIONAL, INTENT(in) :: area !< Id of the area field
396  INTEGER, OPTIONAL, INTENT(in) :: volume !< Id of the volume field
397  CHARACTER(len=*), OPTIONAL, INTENT(in) :: realm !< String to set as the modeling_realm attribute
398  LOGICAL, OPTIONAL, INTENT(in) :: multiple_send_data !< .True. if send data is called, multiple times
399  !! for the same time
400 
401 
402 #ifndef use_yaml
403 fms_register_diag_field_array=diag_field_not_found
404 CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
405 #else
406  fms_register_diag_field_array = this%register( &
407  & module_name, field_name, init_time=init_time, &
408  & axes=axes, longname=long_name, units=units, missing_value=missing_value, varrange=var_range, &
409  & mask_variant=mask_variant, standname=standard_name, do_not_log=do_not_log, err_msg=err_msg, &
410  & interp_method=interp_method, tile_count=tile_count, area=area, volume=volume, realm=realm, &
411  & multiple_send_data=multiple_send_data)
412 #endif
413 end function fms_register_diag_field_array
414 
415 !> @brief Return field index for subsequent call to send_data.
416 !! @return field index to be used in subsequent calls to send_data or DIAG_FIELD_NOT_FOUND if the field is not
417 !! in the diag_table.yaml
418 INTEGER FUNCTION fms_register_static_field(this, module_name, field_name, axes, long_name, units,&
419  & missing_value, range, mask_variant, standard_name, DYNAMIC, do_not_log, interp_method,&
420  & tile_count, area, volume, realm)
421  class(fmsdiagobject_type),TARGET,INTENT(inout):: this !< Diaj_obj to fill
422  CHARACTER(len=*), INTENT(in) :: module_name !< Name of the module, the field is on
423  CHARACTER(len=*), INTENT(in) :: field_name !< Name of the field
424  INTEGER, DIMENSION(:), INTENT(in) :: axes !< Axes_id of the field
425  CHARACTER(len=*), OPTIONAL, INTENT(in) :: long_name !< Longname to be added as a attribute
426  CHARACTER(len=*), OPTIONAL, INTENT(in) :: units !< Units to be added as a attribute
427  CHARACTER(len=*), OPTIONAL, INTENT(in) :: standard_name !< Standard name to be added as a attribute
428  CLASS(*), OPTIONAL, INTENT(in) :: missing_value !< Missing value to be added as a attribute
429  CLASS(*), OPTIONAL, INTENT(in) :: range(:) !< Range to be added as a attribute
430  LOGICAL, OPTIONAL, INTENT(in) :: mask_variant !< .True. if mask changes over time
431  LOGICAL, OPTIONAL, INTENT(in) :: dynamic !< Flag indicating if the field is dynamic
432  LOGICAL, OPTIONAL, INTENT(in) :: do_not_log !< if TRUE, field information is not logged
433  CHARACTER(len=*), OPTIONAL, INTENT(in) :: interp_method !< The interp method to be used when
434  !! regridding the field in post-processing
435  !! Valid options are "conserve_order1",
436  !! "conserve_order2", and "none".
437  INTEGER, OPTIONAL, INTENT(in) :: tile_count !! Number of tiles
438  INTEGER, OPTIONAL, INTENT(in) :: area !< Field ID for the area field associated
439  !! with this field
440  INTEGER, OPTIONAL, INTENT(in) :: volume !< Field ID for the volume field associated
441  !! with this field
442  CHARACTER(len=*), OPTIONAL, INTENT(in) :: realm !< String to set as the value to the
443  !! modeling_realm attribute
444 
445 #ifndef use_yaml
446 fms_register_static_field=diag_field_not_found
447 CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
448 #else
449  !TODO The register_static_field interface does not have the capabiliy to register a variable as a "scalar"
450  ! since the axes argument is required, this forced model code to pass in a null_axis_id as an argument
451  if (size(axes) .eq. 1 .and. axes(1) .eq. null_axis_id) then
452  ! If they are passing in the null_axis_ids, ignore the `axes` argument
453  fms_register_static_field = this%register( &
454  & module_name, field_name, &
455  & longname=long_name, units=units, missing_value=missing_value, varrange=range, &
456  & mask_variant=mask_variant, do_not_log=do_not_log, interp_method=interp_method, tile_count=tile_count, &
457  & standname=standard_name, area=area, volume=volume, realm=realm, &
458  & static=.true.)
459  else
460  fms_register_static_field = this%register( &
461  & module_name, field_name, axes=axes, &
462  & longname=long_name, units=units, missing_value=missing_value, varrange=range, &
463  & mask_variant=mask_variant, do_not_log=do_not_log, interp_method=interp_method, tile_count=tile_count, &
464  & standname=standard_name, area=area, volume=volume, realm=realm, &
465  & static=.true.)
466  endif
467 #endif
468 end function fms_register_static_field
469 
470 !> @brief Wrapper for the register_diag_axis subroutine. This is needed to keep the diag_axis_init
471 !! interface the same
472 !> @return Axis id
473 FUNCTION fms_diag_axis_init(this, axis_name, axis_data, units, cart_name, axis_length, long_name, direction,&
474  & set_name, edges, Domain, Domain2, DomainU, aux, req, tile_count, domain_position ) &
475  & result(id)
476 
477  class(fmsdiagobject_type),TARGET,INTENT(inout):: this !< Diaj_obj to fill
478  CHARACTER(len=*), INTENT(in) :: axis_name !< Name of the axis
479  CLASS(*), INTENT(in) :: axis_data(:) !< Array of coordinate values
480  CHARACTER(len=*), INTENT(in) :: units !< Units for the axis
481  CHARACTER(len=1), INTENT(in) :: cart_name !< Cartesian axis ("X", "Y", "Z", "T", "U", "N")
482  integer, intent(in) :: axis_length !< The length of the axis size(axis_data(:))
483  CHARACTER(len=*), INTENT(in), OPTIONAL :: long_name !< Long name for the axis.
484  CHARACTER(len=*), INTENT(in), OPTIONAL :: set_name !< Name of the parent axis, if it is a subaxis
485  INTEGER, INTENT(in), OPTIONAL :: direction !< Indicates the direction of the axis
486  INTEGER, INTENT(in), OPTIONAL :: edges !< Axis ID for the previously defined "edges axis"
487  TYPE(domain1d), INTENT(in), OPTIONAL :: domain !< 1D domain
488  TYPE(domain2d), INTENT(in), OPTIONAL :: domain2 !< 2D domain
489  TYPE(domainug), INTENT(in), OPTIONAL :: domainu !< Unstructured domain
490  CHARACTER(len=*), INTENT(in), OPTIONAL :: aux !< Auxiliary name, can only be <TT>geolon_t</TT>
491  !! or <TT>geolat_t</TT>
492  CHARACTER(len=*), INTENT(in), OPTIONAL :: req !< Required field names.
493  INTEGER, INTENT(in), OPTIONAL :: tile_count !< Number of tiles
494  INTEGER, INTENT(in), OPTIONAL :: domain_position !< Domain position, "NORTH" or "EAST"
495  integer :: id
496 
497 #ifndef use_yaml
498 id = diag_null
499 CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
500 #else
501  CHARACTER(len=:), ALLOCATABLE :: edges_name !< Name of the edges
502 
503  this%registered_axis = this%registered_axis + 1
504 
505  if (this%registered_axis > max_axes) call mpp_error(fatal, &
506  &"diag_axis_init: max_axes exceeded, increase via diag_manager_nml")
507 
508  allocate(fmsdiagfullaxis_type :: this%diag_axis(this%registered_axis)%axis)
509 
510  select type (axis => this%diag_axis(this%registered_axis)%axis )
511  type is (fmsdiagfullaxis_type)
512  if(present(edges)) then
513  if (edges < 0 .or. edges > this%registered_axis) &
514  call mpp_error(fatal, "diag_axit_init: The edge axis has not been defined. &
515  &Call diag_axis_init for the edge axis first")
516  select type (edges_axis => this%diag_axis(edges)%axis)
517  type is (fmsdiagfullaxis_type)
518  edges_name = edges_axis%get_axis_name()
519  call axis%set_edges(edges_name, edges)
520  end select
521  endif
522  call axis%register(axis_name, axis_data, units, cart_name, long_name=long_name, &
523  & direction=direction, set_name=set_name, domain=domain, domain2=domain2, domainu=domainu, aux=aux, &
524  & req=req, tile_count=tile_count, domain_position=domain_position, axis_length=axis_length)
525 
526  id = this%registered_axis
527  call axis%set_axis_id(id)
528  end select
529 #endif
530 end function fms_diag_axis_init
531 
532 !> Accepts data from the send_data functions. If this is in an openmp region with more than
533 !! one thread, the data is buffered in the field object and processed later. If only a single thread
534 !! is being used, then the processing can be done and stored in the buffer object. The hope is that
535 !! the increase in memory footprint related to buffering can be handled by the shared memory of the
536 !! multithreaded case.
537 !! \note If some of the diag manager is offloaded in the future, then it should be treated similarly
538 !! to the multi-threaded option for processing later
539 logical function fms_diag_accept_data (this, diag_field_id, field_data, mask, rmask, &
540  time, is_in, js_in, ks_in, &
541  ie_in, je_in, ke_in, weight, err_msg)
542  class(fmsdiagobject_type),TARGET, INTENT(inout) :: this !< Diaj_obj to fill
543  INTEGER, INTENT(in) :: diag_field_id !< The ID of the diag field
544  CLASS(*), DIMENSION(:,:,:,:), INTENT(in) :: field_data !< The data for the diag_field
545  LOGICAL, allocatable, INTENT(in) :: mask(:,:,:,:) !< Logical mask indicating the grid
546  !! points to mask (null if no mask)
547  CLASS(*), allocatable, INTENT(in) :: rmask(:,:,:,:)!< real mask indicating the grid
548  !! points to mask (null if no mask)
549  CLASS(*), INTENT(in), OPTIONAL :: weight !< The weight used for averaging
550  TYPE (time_type), INTENT(in), OPTIONAL :: time !< The current time
551  INTEGER, INTENT(in), OPTIONAL :: is_in, js_in, ks_in !< Starting indices
552  INTEGER, INTENT(in), OPTIONAL :: ie_in, je_in, ke_in !< Ending indices
553  CHARACTER(len=*), INTENT(out), OPTIONAL :: err_msg !< An error message returned
554 
555  integer :: is, js, ks !< Starting indicies of the field_data
556  integer :: ie, je, ke !< Ending indicies of the field_data
557  integer :: omp_num_threads !< Number of openmp threads
558  integer :: omp_level !< The openmp active level
559  logical :: buffer_the_data !< True if the user selects to buffer the data and run
560  !! the calculationslater. \note This is experimental
561  character(len=128) :: error_string !< Store error text
562  logical :: data_buffer_is_allocated !< .true. if the data buffer is allocated
563  character(len=256) :: field_info !< String holding info about the field to append to the
564  !! error message
565  logical, allocatable, dimension(:,:,:,:) :: oor_mask !< Out of range mask
566  real(kind=r8_kind) :: field_weight !< Weight to use when averaging (it will be converted
567  !! based on the type of field_data when doing the math)
568  type(fmsdiagibounds_type) :: bounds !< Bounds (starting ending indices) for the field
569  logical :: has_halos !< .True. if field_data contains halos
570  logical :: using_blocking !< .True. if field_data is passed in blocks
571 #ifndef use_yaml
572 CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
573 #else
574 
575  !TODO this%FMS_diag_fields(diag_field_id) should be a pointer!
576  field_info = " Check send data call for field:"//trim(this%FMS_diag_fields(diag_field_id)%get_varname())//&
577  " and module:"//trim(this%FMS_diag_fields(diag_field_id)%get_modname())
578 
579  !< Check if time should be present for this field
580  if (.not.this%FMS_diag_fields(diag_field_id)%is_static() .and. .not.present(time)) &
581  call mpp_error(fatal, "Time must be present if the field is not static. "//trim(field_info))
582 
583  !< Set the field_weight. If "weight" is not present it will be set to 1.0_r8_kind
584  field_weight = set_weight(weight)
585 
586  !< Check that the indices are present in the correct combination
587  error_string = check_indices_order(is_in, ie_in, js_in, je_in)
588  if (trim(error_string) .ne. "") call mpp_error(fatal, trim(error_string)//". "//trim(field_info))
589 
590  using_blocking = .false.
591  if ((present(is_in) .and. .not. present(ie_in)) .or. (present(js_in) .and. .not. present(je_in))) &
592  using_blocking = .true.
593 
594  has_halos = .false.
595  if ((present(is_in) .and. present(ie_in)) .or. (present(js_in) .and. present(je_in))) &
596  has_halos = .true.
597 
598  !< If the field has `mask_variant=.true.`, check that mask OR rmask are present
599  if (this%FMS_diag_fields(diag_field_id)%is_mask_variant()) then
600  if (.not. allocated(mask) .and. .not. allocated(rmask)) call mpp_error(fatal, &
601  "The field was registered with mask_variant, but mask or rmask are not present in the send_data call. "//&
602  trim(field_info))
603  endif
604 
605  !< Check that mask and rmask are not both present
606  if (allocated(mask) .and. allocated(rmask)) call mpp_error(fatal, &
607  "mask and rmask are both present in the send_data call. "//&
608  trim(field_info))
609 
610  !< Create the oor_mask based on the "mask" and "rmask" arguments
611  oor_mask = init_mask(rmask, mask, field_data)
612 
613  !> Does the user want to push off calculations until send_diag_complete?
614  buffer_the_data = .false.
615 
616  !> initialize the number of threads and level to be 0
617  omp_num_threads = 0
618  omp_level = 0
619 #if defined(_OPENMP)
620  omp_num_threads = omp_get_num_threads()
621  omp_level = omp_get_level()
622  buffer_the_data = (omp_num_threads > 1 .AND. omp_level > 0)
623 #endif
624 
625  !> Calculate the i,j,k start and end
626  ! If is, js, or ks not present default them to 1
627  is = 1
628  js = 1
629  ks = 1
630  IF ( PRESENT(is_in) ) is = is_in
631  IF ( PRESENT(js_in) ) js = js_in
632  IF ( PRESENT(ks_in) ) ks = ks_in
633  ie = is+SIZE(field_data, 1)-1
634  je = js+SIZE(field_data, 2)-1
635  ke = ks+SIZE(field_data, 3)-1
636  IF ( PRESENT(ie_in) ) ie = ie_in
637  IF ( PRESENT(je_in) ) je = je_in
638  IF ( PRESENT(ke_in) ) ke = ke_in
639 
640  if (.not. buffer_the_data .and. using_blocking) then
641  ! If running with only 1 thread and using blocking, check if the data was sent in blocks
642  ! if it is, then buffer the data
643  buffer_the_data = check_for_slices(this%FMS_diag_fields(diag_field_id), this%diag_axis, &
644  shape(field_data))
645  endif
646 
647  !< If send data is called multiple times, buffer the data
648  !! This is so that the other reduction methods work and just averaging
649  if (this%FMS_diag_fields(diag_field_id)%get_multiple_send_data()) &
650  buffer_the_data = .true.
651 
652  !If this is true, buffer data
653  main_if: if (buffer_the_data) then
654 !> Only 1 thread allocates the output buffer and sets set_math_needs_to_be_done
655 !$omp critical
656 
657  !< These set_* calls need to be done inside an omp_critical to avoid any race conditions
658  !! and allocation issues
659  if(has_halos) call this%FMS_diag_fields(diag_field_id)%set_halo_present()
660 
661  !< Set the variable type based off passed in field data
662  if(.not. this%FMS_diag_fields(diag_field_id)%has_vartype()) &
663  call this%FMS_diag_fields(diag_field_id)%set_type(field_data(1,1,1,1))
664 
665  if (allocated(mask) .or. allocated(rmask)) then
666  call this%FMS_diag_fields(diag_field_id)%set_var_is_masked(.true.)
667  else
668  call this%FMS_diag_fields(diag_field_id)%set_var_is_masked(.false.)
669  endif
670 
671  if (.not. this%FMS_diag_fields(diag_field_id)%is_data_buffer_allocated()) then
672  data_buffer_is_allocated = &
673  this%FMS_diag_fields(diag_field_id)%allocate_data_buffer(field_data, this%diag_axis)
674  if(.not. this%FMS_diag_fields(diag_field_id)%has_mask_allocated()) &
675  call this%FMS_diag_fields(diag_field_id)%allocate_mask(oor_mask, this%diag_axis)
676  endif
677  call this%FMS_diag_fields(diag_field_id)%set_send_data_time(time)
678  call this%FMS_diag_fields(diag_field_id)%set_data_buffer_is_allocated(.true.)
679  call this%FMS_diag_fields(diag_field_id)%set_math_needs_to_be_done(.true.)
680 !$omp end critical
681  call this%FMS_diag_fields(diag_field_id)%set_data_buffer(field_data, oor_mask, field_weight, &
682  is, js, ks, ie, je, ke)
683  fms_diag_accept_data = .true.
684  return
685  else
686 
687  !< At this point if we are no longer in an openmp region or running with 1 thread
688  !! so it is safe to have these set_* calls
689  if(has_halos) call this%FMS_diag_fields(diag_field_id)%set_halo_present()
690 
691  !< Set the variable type based off passed in field data
692  if(.not. this%FMS_diag_fields(diag_field_id)%has_vartype()) &
693  call this%FMS_diag_fields(diag_field_id)%set_type(field_data(1,1,1,1))
694 
695  if (allocated(mask) .or. allocated(rmask)) then
696  call this%FMS_diag_fields(diag_field_id)%set_var_is_masked(.true.)
697  else
698  call this%FMS_diag_fields(diag_field_id)%set_var_is_masked(.false.)
699  endif
700 
701  error_string = bounds%set_bounds(field_data, is, ie, js, je, ks, ke, has_halos)
702  if (trim(error_string) .ne. "") call mpp_error(fatal, trim(error_string)//". "//trim(field_info))
703 
704  call this%allocate_diag_field_output_buffers(field_data, diag_field_id)
705  error_string = this%fms_diag_do_reduction(field_data, diag_field_id, oor_mask, field_weight, &
706  bounds, using_blocking, time=time)
707  if (trim(error_string) .ne. "") call mpp_error(fatal, trim(error_string)//". "//trim(field_info))
708  call this%FMS_diag_fields(diag_field_id)%set_math_needs_to_be_done(.false.)
709  if(.not. this%FMS_diag_fields(diag_field_id)%has_mask_allocated()) &
710  call this%FMS_diag_fields(diag_field_id)%allocate_mask(oor_mask)
711  call this%FMS_diag_fields(diag_field_id)%set_mask(oor_mask, field_info)
712  return
713  end if main_if
714  !> Return false if nothing is done
715  fms_diag_accept_data = .false.
716  return
717 #endif
718 end function fms_diag_accept_data
719 
720 !< @brief Do the math for all the buffers
721 subroutine do_buffer_math(this)
722  class(fmsdiagobject_type), target, intent (inout) :: this !< The diag object
723 
724 #ifdef use_yaml
725  integer :: i !< For do loops
726  integer :: ifile !< For file loops
727  integer :: ifield !< For field loops
728 
729  class(fmsdiagfilecontainer_type), pointer :: diag_file !< Pointer to this%FMS_diag_files(i) (for convenience
730  class(fmsdiagfield_type), pointer :: diag_field !< Pointer to this%FMS_diag_files(i)%diag_field(j)
731  logical :: math !< True if the math functions need to be called using the data buffer,
732  !! False if the math functions were done in accept_data
733  integer, dimension(:), allocatable :: file_field_ids !< Array of field IDs for a file
734  class(*), pointer :: input_data_buffer(:,:,:,:)
735  character(len=128) :: error_string
736  type(fmsdiagibounds_type) :: bounds
737  integer, dimension(:), allocatable :: file_ids !< Array of file IDs for a field
738  logical, parameter :: debug_sc = .false. !< turn on output for debugging
739 
740  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
741  !! In the future, this may be parallelized for offloading
742  ! loop through each field
743  field_loop: do ifield = 1, size(this%FMS_diag_fields)
744  diag_field => this%FMS_diag_fields(ifield)
745  if(.not. diag_field%is_registered()) cycle
746  if(debug_sc) call mpp_error(note, "fms_diag_send_complete:: var: "//diag_field%get_varname())
747  ! get files the field is in
748  allocate (file_ids(size(diag_field%get_file_ids() )))
749  file_ids = diag_field%get_file_ids()
750  math = diag_field%get_math_needs_to_be_done()
751  ! if doing math loop through each file for given field
752  doing_math: if (size(file_ids) .ge. 1 .and. math) then
753  ! Check if buffer alloc'd
754  has_input_buff: if (diag_field%has_input_data_buffer()) then
755  call diag_field%prepare_data_buffer()
756  input_data_buffer => diag_field%get_data_buffer()
757  ! reset bounds, allocate output buffer, and update it with reduction
758  call bounds%reset_bounds_from_array_4D(input_data_buffer)
759  call this%allocate_diag_field_output_buffers(input_data_buffer, ifield)
760  error_string = this%fms_diag_do_reduction(input_data_buffer, ifield, &
761  diag_field%get_mask(), diag_field%get_weight(), &
762  bounds, .false., time=diag_field%get_send_data_time())
763  call diag_field%init_data_buffer()
764  if (trim(error_string) .ne. "") call mpp_error(fatal, "Field:"//trim(diag_field%get_varname()//&
765  " -"//trim(error_string)))
766  else
767  call mpp_error(fatal, "diag_send_complete:: no input buffer allocated for field"//diag_field%get_longname())
768  endif has_input_buff
769  endif doing_math
770  call diag_field%set_math_needs_to_be_done(.false.)
771  !> Clean up, clean up, everybody do your share
772  if (allocated(file_ids)) deallocate(file_ids)
773  if (associated(diag_field)) nullify(diag_field)
774  enddo field_loop
775 #endif
776 end subroutine do_buffer_math
777 
778 !> @brief Loops through all the files, open the file, writes out axis and
779 !! variable metadata and data when necessary.
780 subroutine fms_diag_send_complete(this, time_step)
781  class(fmsdiagobject_type), target, intent (inout) :: this !< The diag object
782  TYPE (time_type), INTENT(in) :: time_step !< The time_step
783 
784 #ifndef use_yaml
785 CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
786 #else
787  call this%do_buffer_math()
788  call this%fms_diag_do_io()
789 #endif
790 
791 end subroutine fms_diag_send_complete
792 
793 !> @brief Loops through all the files, open the file, writes out axis and
794 !! variable metadata and data when necessary.
795 !! TODO: passing in the saved mask from the field obj to diag_reduction_done_wrapper
796 !! for performance
797 subroutine fms_diag_do_io(this, end_time)
798  class(fmsdiagobject_type), target, intent(inout) :: this !< The diag object
799  type(time_type), optional, target, intent(in) :: end_time !< the model end_time
800 #ifdef use_yaml
801  integer :: i !< For do loops
802  class(fmsdiagfilecontainer_type), pointer :: diag_file !< Pointer to this%FMS_diag_files(i) (for convenience)
803  class(fmsdiagoutputbuffer_type), pointer :: diag_buff !< pointer to output buffers iterated in buff_loop
804  class(fmsdiagfield_type), pointer :: diag_field !< pointer to output buffers iterated in buff_loop
805  class(diagyamlfilesvar_type), pointer :: field_yaml !< Pointer to a field from yaml fields
806  TYPE (time_type), pointer :: model_time!< The current model time
807  integer, allocatable :: buff_ids(:) !< ids for output buffers to loop through
808  integer :: ibuff !< buffer index
809  logical :: file_is_opened_this_time_step !< True if the file was opened in this time_step
810  !! If true the metadata will need to be written
811  logical :: force_write !< force the last write if at end of run
812  logical :: finish_writing !< true if finished writing for all the fields
813  logical :: has_mask !< whether we have a mask
814  logical, parameter :: debug_reduct = .false. !< enables debugging output
815  class(*), allocatable :: missing_val !< netcdf missing value for a given field
816  real(r8_kind) :: mval !< r8 copy of missing value
817  character(len=128) :: error_string !< outputted error string from reducti
818  logical :: unlim_dim_was_increased !< .True. if the unlimited dimension index was increased for any of the buffers
819  logical :: do_not_write !< .True. only if this is not a new time step and you are writting at every time step
820 
821  force_write = .false.
822 
823  do i = 1, size(this%FMS_diag_files)
824  diag_file => this%FMS_diag_files(i)
825 
826  !< Go away if the file is a subregional file and the current PE does not have any data for it
827  if (.not. diag_file%writing_on_this_pe()) cycle
828  if (diag_file%FMS_diag_file%is_done_writing_data()) cycle
829 
830  if (present (end_time)) then
831  force_write = .true.
832  model_time => end_time
833  else
834  model_time => diag_file%get_model_time()
835  endif
836 
837  call diag_file%open_diag_file(model_time, file_is_opened_this_time_step)
838  if (file_is_opened_this_time_step) then
839  ! Initialize unlimited dimension in file and the buffer to 0
840  call diag_file%init_unlim_dim(this%FMS_diag_output_buffers)
841 
842  call diag_file%write_global_metadata()
843  call diag_file%write_axis_metadata(this%diag_axis)
844  call diag_file%write_time_metadata()
845  call diag_file%write_field_metadata(this%FMS_diag_fields, this%diag_axis)
846  call diag_file%write_axis_data(this%diag_axis)
847  endif
848 
849  finish_writing = diag_file%is_time_to_write(model_time, this%FMS_diag_output_buffers, &
850  this%FMS_diag_fields, do_not_write)
851  unlim_dim_was_increased = .false.
852 
853  ! finish reduction method if its time to write
854  buff_ids = diag_file%FMS_diag_file%get_buffer_ids()
855  ! loop through the buffers and finish reduction if needed
856  buff_loop: do ibuff=1, SIZE(buff_ids)
857  diag_buff => this%FMS_diag_output_buffers(buff_ids(ibuff))
858  field_yaml => diag_yaml%diag_fields(diag_buff%get_yaml_id())
859  diag_field => this%FMS_diag_fields(diag_buff%get_field_id())
860 
861  ! Go away if there is no data to write
862  if (.not. diag_buff%is_there_data_to_write()) cycle
863 
864  if ( diag_buff%is_time_to_finish_reduction(end_time) .and. .not. do_not_write) then
865  ! sets missing value
866  mval = diag_field%find_missing_value(missing_val)
867  ! time_average and greater values all involve averaging so need to be "finished" before written
868  if( field_yaml%has_var_reduction()) then
869  if( field_yaml%get_var_reduction() .ge. time_average) then
870  if(debug_reduct)call mpp_error(note, "fms_diag_do_io:: finishing reduction for "//diag_field%get_longname())
871  error_string = diag_buff%diag_reduction_done_wrapper( &
872  field_yaml%get_var_reduction(), &
873  mval, diag_field%get_var_is_masked(), diag_field%get_mask_variant())
874  endif
875  endif
876  call diag_file%write_field_data(diag_field, diag_buff, unlim_dim_was_increased)
877  call diag_buff%set_next_output(diag_file%get_next_output(), diag_file%get_next_next_output())
878  endif
879  nullify(diag_buff)
880  nullify(field_yaml)
881  enddo buff_loop
882  deallocate(buff_ids)
883 
884  if (unlim_dim_was_increased) then
885  call diag_file%write_time_data()
886  call diag_file%flush_diag_file()
887  call diag_file%update_next_write(model_time)
888  endif
889 
890  if (finish_writing) then
891  call diag_file%update_current_new_file_freq_index(model_time)
892  if (diag_file%is_time_to_close_file(model_time)) call diag_file%close_diag_file(this%FMS_diag_output_buffers, &
893  diag_fields = this%FMS_diag_fields)
894  else if (force_write) then
895  call diag_file%prepare_for_force_write()
896  call diag_file%write_time_data()
897  call diag_file%close_diag_file(this%FMS_diag_output_buffers, diag_fields = this%FMS_diag_fields)
898  endif
899  enddo
900 #endif
901 end subroutine fms_diag_do_io
902 
903 !> @brief Computes average, min, max, rms error, etc.
904 !! based on the specified reduction method for the field.
905 !> @return Empty string if successful, error message if it fails
906 function fms_diag_do_reduction(this, field_data, diag_field_id, oor_mask, weight, &
907  bounds, using_blocking, time) &
908  result(error_msg)
909  class(fmsdiagobject_type), intent(inout), target:: this !< Diag Object
910  class(*), intent(in) :: field_data(:,:,:,:) !< Field data
911  integer, intent(in) :: diag_field_id !< ID of the input field
912  logical, intent(in), target :: oor_mask(:,:,:,:) !< mask
913  real(kind=r8_kind), intent(in) :: weight !< Must be a updated weight
914  type(fmsdiagibounds_type), intent(in) :: bounds !< Bounds for the field
915  logical, intent(in) :: using_blocking !< .True. if field data is passed
916  !! in blocks
917  type(time_type), intent(in), optional :: time !< Current time
918 
919  character(len=150) :: error_msg !< Error message to check
920  !TODO Mostly everything
921 #ifdef use_yaml
922  type(fmsdiagfield_type), pointer :: field_ptr !< Pointer to the field's object
923  type(fmsdiagoutputbuffer_type), pointer :: buffer_ptr !< Pointer to the field's buffer
924  class(fmsdiagfilecontainer_type), pointer :: file_ptr !< Pointer to the field's file
925  type(diagyamlfilesvar_type), pointer :: field_yaml_ptr !< Pointer to the field's yaml
926 
927  integer :: reduction_method !< Integer representing a reduction method
928  integer :: ids !< For looping through buffer ids
929  integer :: buffer_id !< Id of the buffer
930  integer :: file_id !< File id
931  integer, pointer :: axis_ids(:) !< Axis ids for the buffer
932  logical :: is_subregional !< .True. if the buffer is subregional
933  logical :: reduced_k_range !< .True. is the field is only outputing a section
934  !! of the z dimension
935  type(fmsdiagibounds_type) :: bounds_in !< Starting and ending indices of the input field_data
936  type(fmsdiagibounds_type) :: bounds_out !< Starting and ending indices of the output buffer
937  integer :: i !< For looping through axid ids
938  integer :: sindex !< Starting index of a subregion
939  integer :: eindex !< Ending index of a subregion
940  integer :: compute_idx(2) !< Starting and Ending of the compute domain
941  character(len=1) :: cart_axis !< Cartesian axis of the axis
942  logical :: block_in_subregion !< .True. if the current block is part of the subregion
943  integer :: starting !< Starting index of the subregion relative to the compute domain
944  integer :: ending !< Ending index of the subregion relative to the compute domain
945  real(kind=r8_kind) :: missing_value !< Missing_value for data points that are masked
946  !! This will obtained as r8 and converted to the right type as
947  !! needed. This is to avoid yet another select type ...
948 
949  !TODO mostly everything
950  field_ptr => this%FMS_diag_fields(diag_field_id)
951  if (field_ptr%has_missing_value()) then
952  select type (missing_val => field_ptr%get_missing_value(r8))
953  type is (real(kind=r8_kind))
954  missing_value = missing_val
955  class default
956  call mpp_error(fatal, "The missing value for the field:"//trim(field_ptr%get_varname())//&
957  &" was not allocated to the correct type. This shouldn't have happened")
958  end select
959  else
960  select type (missing_val => get_default_missing_value(r8))
961  type is (real(kind=r8_kind))
962  missing_value = missing_val
963  class default
964  call mpp_error(fatal, "The missing value for the field:"//trim(field_ptr%get_varname())//&
965  &" was not allocated to the correct type. This shouldn't have happened")
966  end select
967  endif
968 
969  buffer_loop: do ids = 1, size(field_ptr%buffer_ids)
970  error_msg = ""
971  buffer_id = this%FMS_diag_fields(diag_field_id)%buffer_ids(ids)
972  file_id = this%FMS_diag_fields(diag_field_id)%file_ids(ids)
973 
974  !< Gather all the objects needed for the buffer
975  field_yaml_ptr => field_ptr%diag_field(ids)
976  buffer_ptr => this%FMS_diag_output_buffers(buffer_id)
977  file_ptr => this%FMS_diag_files(file_id)
978 
979  !< Go away if the file is a subregional file and the current PE does not have any data for it
980  if (.not. file_ptr%writing_on_this_pe()) cycle
981 
982  !< Go away if finished doing math for this buffer
983  if (buffer_ptr%is_done_with_math()) cycle
984 
985  if (present(time)) call file_ptr%set_model_time(time)
986 
987  bounds_out = bounds
988  if (.not. using_blocking) then
989  !< Set output bounds to start at 1:size(buffer_ptr%buffer)
990  call bounds_out%reset_bounds_from_array_4D(buffer_ptr%buffer(:,:,:,:,1))
991  endif
992 
993  bounds_in = bounds
994  if (.not. bounds%has_halos) then
995  !< If field_data does not contain halos, set bounds_in to start at 1:size(field_data)
996  call bounds_in%reset_bounds_from_array_4D(field_data)
997  endif
998 
999  is_subregional = file_ptr%is_regional()
1000  reduced_k_range = field_yaml_ptr%has_var_zbounds()
1001 
1002  !< Reset the bounds based on the reduced k range and subregional
1003  is_subregional_reduced_k_range: if (is_subregional .or. reduced_k_range) then
1004  call buffer_ptr%get_axis_ids(axis_ids)
1005  block_in_subregion = .true.
1006  axis_loops: do i = 1, size(axis_ids)
1007  !< Move on if the block does not have any data for the subregion
1008  if (.not. block_in_subregion) cycle
1009 
1010  select type (diag_axis => this%diag_axis(axis_ids(i))%axis)
1011  type is (fmsdiagsubaxis_type)
1012  sindex = diag_axis%get_starting_index()
1013  eindex = diag_axis%get_ending_index()
1014  compute_idx = diag_axis%get_compute_indices()
1015  starting=sindex-compute_idx(1)+1
1016  ending=eindex-compute_idx(1)+1
1017  if (using_blocking) then
1018  block_in_subregion = determine_if_block_is_in_region(starting, ending, bounds, i)
1019  if (.not. block_in_subregion) cycle
1020 
1021  !< Set bounds_in so that you can the correct section of the data for the block (starting at 1)
1022  call bounds_in%rebase_input(bounds, starting, ending, i)
1023 
1024  !< Set bounds_out to be the correct section relative to the block starting and ending indices
1025  call bounds_out%rebase_output(starting, ending, i)
1026  else
1027  !< Set bounds_in so that only the subregion section of the data will be used (starting at 1)
1028  call bounds_in%update_index(starting, ending, i, .false.)
1029 
1030  !< Set bounds_out to 1:size(subregion) for the PE
1031  call bounds_out%update_index(1, ending-starting+1, i, .true.)
1032  endif
1033  end select
1034  enddo axis_loops
1035  nullify(axis_ids)
1036  !< Move on to the next buffer if the block does not have any data for the subregion
1037  if (.not. block_in_subregion) cycle
1038  endif is_subregional_reduced_k_range
1039 
1040  !< Determine the reduction method for the buffer
1041  reduction_method = field_yaml_ptr%get_var_reduction()
1042  if (present(time)) call buffer_ptr%update_buffer_time(time)
1043  call buffer_ptr%set_send_data_called()
1044  select case(reduction_method)
1045  case (time_none)
1046  error_msg = buffer_ptr%do_time_none_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
1047  bounds_in, bounds_out, missing_value)
1048  if (trim(error_msg) .ne. "") then
1049  return
1050  endif
1051  case (time_min)
1052  error_msg = buffer_ptr%do_time_min_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
1053  bounds_in, bounds_out, missing_value)
1054  if (trim(error_msg) .ne. "") then
1055  return
1056  endif
1057  case (time_max)
1058  error_msg = buffer_ptr%do_time_max_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
1059  bounds_in, bounds_out, missing_value)
1060  if (trim(error_msg) .ne. "") then
1061  return
1062  endif
1063  case (time_sum)
1064  error_msg = buffer_ptr%do_time_sum_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
1065  field_ptr%get_mask_variant(), bounds_in, bounds_out, missing_value, field_ptr%has_missing_value())
1066  if (trim(error_msg) .ne. "") then
1067  return
1068  endif
1069  case (time_average)
1070  error_msg = buffer_ptr%do_time_sum_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
1071  field_ptr%get_mask_variant(), bounds_in, bounds_out, missing_value, field_ptr%has_missing_value(), &
1072  weight=weight)
1073  if (trim(error_msg) .ne. "") then
1074  return
1075  endif
1076  case (time_power)
1077  error_msg = buffer_ptr%do_time_sum_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
1078  field_ptr%get_mask_variant(), bounds_in, bounds_out, missing_value, field_ptr%has_missing_value(), &
1079  weight=weight, pow_value=field_yaml_ptr%get_pow_value())
1080  if (trim(error_msg) .ne. "") then
1081  return
1082  endif
1083  case (time_rms)
1084  error_msg = buffer_ptr%do_time_sum_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
1085  field_ptr%get_mask_variant(), bounds_in, bounds_out, missing_value, field_ptr%has_missing_value(), &
1086  weight=weight, pow_value = 2)
1087  if (trim(error_msg) .ne. "") then
1088  return
1089  endif
1090  case (time_diurnal)
1091  if(.not. present(time)) call mpp_error(fatal, &
1092  "fms_diag_do_reduction:: time must be present when using diurnal reductions")
1093  ! sets the diurnal index for reduction within the buffer object
1094  call buffer_ptr%set_diurnal_section_index(time)
1095  error_msg = buffer_ptr%do_time_sum_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
1096  field_ptr%get_mask_variant(), bounds_in, bounds_out, missing_value, field_ptr%has_missing_value(), &
1097  weight=weight)
1098  if (trim(error_msg) .ne. "") then
1099  return
1100  endif
1101  case default
1102  error_msg = "The reduction method is not supported. "//&
1103  "Only none, min, max, sum, average, power, rms, and diurnal are supported."
1104  end select
1105 
1106  if (field_ptr%is_static() .or. file_ptr%FMS_diag_file%is_done_writing_data()) then
1107  call buffer_ptr%set_done_with_math()
1108  endif
1109  enddo buffer_loop
1110 #else
1111  error_msg = ""
1112  CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
1113 #endif
1114 end function fms_diag_do_reduction
1115 
1116 !> @brief Adds the diag ids of the Area and or Volume of the diag_field_object
1117 subroutine fms_diag_field_add_cell_measures(this, diag_field_id, area, volume)
1118  class(fmsdiagobject_type), intent (inout) :: this !< The diag object
1119  integer, intent(in) :: diag_field_id !< diag_field to add the are and volume to
1120  INTEGER, optional, INTENT(in) :: area !< diag ids of area
1121  INTEGER, optional, INTENT(in) :: volume !< diag ids of volume
1122 
1123 #ifndef use_yaml
1124  CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
1125 #else
1126  call this%FMS_diag_fields(diag_field_id)%add_area_volume(area, volume)
1127 #endif
1128 end subroutine fms_diag_field_add_cell_measures
1129 
1130 !> @brief Add a attribute to the diag_obj using the diag_field_id
1131 subroutine fms_diag_field_add_attribute(this, diag_field_id, att_name, att_value)
1132  class(fmsdiagobject_type), intent (inout) :: this !< The diag object
1133  integer, intent(in) :: diag_field_id !< Id of the axis to add the attribute to
1134  character(len=*), intent(in) :: att_name !< Name of the attribute
1135  class(*), intent(in) :: att_value(:) !< The attribute value to add
1136 #ifndef use_yaml
1137 CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
1138 #else
1139 !TODO: Value for diag not found
1140  if ( diag_field_id .LE. 0 ) THEN
1141  RETURN
1142  else
1143  if (this%FMS_diag_fields(diag_field_id)%is_registered() ) &
1144  call this%FMS_diag_fields(diag_field_id)%add_attribute(att_name, att_value)
1145  endif
1146 #endif
1147 end subroutine fms_diag_field_add_attribute
1148 
1149 !> @brief Add an attribute to an axis
1150 subroutine fms_diag_axis_add_attribute(this, axis_id, att_name, att_value)
1151  class(fmsdiagobject_type), intent (inout) :: this !< The diag object
1152  integer, intent(in) :: axis_id !< Id of the axis to add the attribute to
1153  character(len=*), intent(in) :: att_name !< Name of the attribute
1154  class(*), intent(in) :: att_value(:) !< The attribute value to add
1155 
1156  character(len=20) :: axis_names(2) !< Names of the uncompress axis
1157  character(len=20) :: set_name !< Name of the axis set
1158  integer :: uncmx_ids(2) !< Ids of the uncompress axis
1159  integer :: j !< For do loops
1160 #ifndef use_yaml
1161 CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
1162 #else
1163  if (axis_id < 0 .and. axis_id > this%registered_axis) &
1164  call mpp_error(fatal, "diag_axis_add_attribute: The axis_id is not valid")
1165 
1166  select type (axis => this%diag_axis(axis_id)%axis)
1167  type is (fmsdiagfullaxis_type)
1168  call axis%add_axis_attribute(att_name, att_value)
1169 
1170  !! Axis that are in the "unstructured" domain require a "compress" attribute for the
1171  !! combiner and PP. This attribute is passed in via a diag_axis_add_attribute call in the model code
1172  !! The compress attribute indicates the names of the axis that were compressed
1173  !! For example grid_index:compress = "grid_yt grid_xt"
1174  !! The metadata and the data for these axis also needs to be written to the file
1175  if (trim(att_name) .eq. "compress") then
1176  !< If the attribute is the "compress" attribute, get the axis names,
1177  !! and the ids of the axis and add it to the axis object so it can be written to netcdf files
1178  !! that use this axis
1179  axis_names = parse_compress_att(att_value)
1180  set_name = ""
1181  if (axis%has_set_name()) set_name = axis%get_set_name()
1182  do j = 1, size(axis_names)
1183  uncmx_ids(j) = get_axis_id_from_name(axis_names(j), this%diag_axis, this%registered_axis, set_name)
1184  if (uncmx_ids(j) .eq. diag_null) call mpp_error(fatal, &
1185  &"Error parsing the compress attribute for axis: "//trim(axis%get_axis_name())//&
1186  &". Be sure that the axes in the compress attribute are registered")
1187  enddo
1188  call axis%add_structured_axis_ids(uncmx_ids)
1189  endif
1190  end select
1191 #endif
1192 end subroutine fms_diag_axis_add_attribute
1193 
1194 !> \brief Gets the field_name from the diag_field
1195 !> \returns a copy of the field_name
1196 function fms_get_field_name_from_id (this, field_id) &
1197  result(field_name)
1198 
1199  class(fmsdiagobject_type), intent (in) :: this !< The diag object, the caller
1200  integer, intent (in) :: field_id !< Field id to get the name for
1201  character(len=:), allocatable :: field_name
1202 #ifndef use_yaml
1203  CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
1204 #else
1205  field_name = this%FMS_diag_fields(field_id)%get_varname()
1206 #endif
1207 end function fms_get_field_name_from_id
1208 
1209 !> \brief Gets the diag field ID from the module name and field name.
1210 !> \returns a copy of the ID of the diag field or DIAG_FIELD_NOT_FOUND if the field is not registered
1211 FUNCTION fms_get_diag_field_id_from_name(this, module_name, field_name) &
1212  result(diag_field_id)
1213  class(fmsdiagobject_type), intent (in) :: this !< The diag object, the caller
1214  CHARACTER(len=*), INTENT(in) :: module_name !< Module name that registered the variable
1215  CHARACTER(len=*), INTENT(in) :: field_name !< Variable name
1216  integer :: diag_field_id
1217 
1218 #ifdef use_yaml
1219  integer :: i !< For looping
1220  integer, allocatable :: diag_field_indices(:) !< indices where the field was found in the yaml
1221 
1222  diag_field_id = diag_field_not_found
1223 
1224  !> Loop through fields to find it.
1225  do i=1, this%registered_variables
1226  !< Check if the field was registered, if it was return the diag_field_id
1227  diag_field_id = this%FMS_diag_fields(i)%id_from_name(module_name, field_name)
1228  if(diag_field_id .ne. diag_field_not_found) return
1229  enddo
1230 
1231  !< Check if the field is in the diag_table.yaml. If it is, return DIAG_FIELD_NOT_REGISTERED
1232  !! Otherwsie it will return DIAG_FIELD_NOT_FOUND
1233  diag_field_indices = find_diag_field(field_name, module_name)
1234  if (diag_field_indices(1) .ne. diag_null) then
1235  diag_field_id = diag_not_registered
1236  endif
1237  deallocate(diag_field_indices)
1238 #else
1239  diag_field_id = diag_field_not_found
1240  CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
1241 #endif
1242 END FUNCTION fms_get_diag_field_id_from_name
1243 
1244 #ifdef use_yaml
1245 !> returns the buffer object for the given id
1246 !! actual data comes from %get_buffer_data() on the returned object
1247 function get_diag_buffer(this, bufferid) &
1248 result(rslt)
1249  class(fmsdiagobject_type), intent(in) :: this
1250  integer, intent(in) :: bufferid
1251  type(fmsdiagoutputbuffer_type),allocatable:: rslt
1252  if( (bufferid .gt. ubound(this%FMS_diag_output_buffers, 1)) .or. &
1253  (bufferid .lt. lbound(this%FMS_diag_output_buffers, 1))) &
1254  call mpp_error(fatal, 'get_diag_bufer: invalid bufferid given')
1255  rslt = this%FMS_diag_output_buffers(bufferid)
1256 end function
1257 #endif
1258 
1259 !> @brief Return the 2D domain for the axis IDs given.
1260 !! @return 2D domain for the axis IDs given
1261 type(domain2d) function fms_get_domain2d(this, ids)
1262  class(fmsdiagobject_type), intent (in) :: this !< The diag object
1263  INTEGER, DIMENSION(:), INTENT(in) :: ids !< Axis IDs.
1264 
1265 #ifndef use_yaml
1266 CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
1267 fms_get_domain2d = null_domain2d
1268 #else
1269  INTEGER :: type_of_domain !< The type of domain
1270  CLASS(diagdomain_t), POINTER :: domain !< Diag Domain pointer
1271 
1272  call get_domain_and_domain_type(fms_diag_object%diag_axis, ids, type_of_domain, domain, "get_domain2d")
1273  if (type_of_domain .ne. two_d_domain) &
1274  call mpp_error(fatal, 'diag_axis_mod::get_domain2d- The axis do not correspond to a 2d Domain')
1275  select type(domain)
1276  type is (diagdomain2d_t)
1277  fms_get_domain2d = domain%domain2
1278  end select
1279 #endif
1280 END FUNCTION fms_get_domain2d
1281 
1282  !> @brief Gets the length of the axis based on the axis_id
1283  !> @return Axis_length
1284  integer function fms_get_axis_length(this, axis_id)
1285  class(fmsdiagobject_type), intent (in) :: this !< The diag object
1286  INTEGER, INTENT(in) :: axis_id !< Axis ID of the axis to the length of
1287 
1288 #ifndef use_yaml
1289 CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
1290 fms_get_axis_length = 0
1291 #else
1292 fms_get_axis_length = 0
1293 
1294  if (axis_id < 0 .and. axis_id > this%registered_axis) &
1295  call mpp_error(fatal, "fms_get_axis_length: The axis_id is not valid")
1296 
1297  select type (axis => this%diag_axis(axis_id)%axis)
1298  type is (fmsdiagfullaxis_type)
1299  fms_get_axis_length = axis%axis_length()
1300  type is (fmsdiagsubaxis_type)
1301  fms_get_axis_length = axis%axis_length()
1302  end select
1303 #endif
1304 end function fms_get_axis_length
1305 
1306 !> @brief Gets the name of the axis based on the axis_id
1307  !> @return The axis_name
1308 function fms_get_axis_name_from_id (this, axis_id) &
1309 result(axis_name)
1310  class(fmsdiagobject_type), intent (in) :: this !< The diag object
1311  INTEGER, INTENT(in) :: axis_id !< Axis ID of the axis to the length of
1312 
1313  character (len=:), allocatable :: axis_name
1314 
1315 #ifndef use_yaml
1316 CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
1317 axis_name=" "
1318 #else
1319  if (axis_id < 0 .and. axis_id > this%registered_axis) &
1320  call mpp_error(fatal, "fms_get_axis_length: The axis_id is not valid")
1321 
1322  !! if its a scalar (null axis id) just returns n/a since no axis is defined
1323  if (axis_id .eq. null_axis_id) then
1324  allocate(character(len=3) :: axis_name)
1325  axis_name = "n/a"
1326  return
1327  endif
1328 
1329 
1330  select type (axis => this%diag_axis(axis_id)%axis)
1331  type is (fmsdiagfullaxis_type)
1332  axis_name = axis%get_axis_name()
1333  end select
1334 #endif
1335 end function fms_get_axis_name_from_id
1336 
1337 !> Dumps as much data as it can from the fmsDiagObject_type.
1338 !! Will dump any fields and files as well (see d)
1339 subroutine dump_diag_obj( filename )
1340  character(len=*), intent(in), optional :: filename !< optional filename to print to,
1341  !! otherwise prints to stdout
1342 #ifdef use_yaml
1343  !type(fmsDiagObject_type) :: diag_obj
1344  type(fmsdiagfile_type), pointer :: fileptr !< pointer for traversing file list
1345  type(fmsdiagfield_type), pointer :: fieldptr !< pointer for traversing field list
1346  integer :: i !< do loops
1347  integer :: unit_num !< unit num of opened log file or stdout
1348 
1349  if( present(filename) ) then
1350  open(newunit=unit_num, file=trim(filename), action='WRITE')
1351  else
1352  unit_num = stdout()
1353  endif
1354  if( mpp_pe() .eq. mpp_root_pe()) then
1355  write(unit_num, *) '********** dumping diag object ***********'
1356  write(unit_num, *) 'registered_variables:', fms_diag_object%registered_variables
1357  write(unit_num, *) 'registered_axis:', fms_diag_object%registered_axis
1358  write(unit_num, *) 'initialized:', fms_diag_object%initialized
1359  write(unit_num, *) 'files_initialized:', fms_diag_object%files_initialized
1360  write(unit_num, *) 'fields_initialized:', fms_diag_object%fields_initialized
1361  write(unit_num, *) 'buffers_initialized:', fms_diag_object%buffers_initialized
1362  write(unit_num, *) 'axes_initialized:', fms_diag_object%axes_initialized
1363  write(unit_num, *) 'Files:'
1364  if( fms_diag_object%files_initialized ) then
1365  do i=1, SIZE(fms_diag_object%FMS_diag_files)
1366  write(unit_num, *) 'File num:', i
1367  fileptr => fms_diag_object%FMS_diag_files(i)%FMS_diag_file
1368  call fileptr%dump_file_obj(unit_num)
1369  enddo
1370  else
1371  write(unit_num, *) 'files not initialized'
1372  endif
1373  if( fms_diag_object%fields_initialized) then
1374  do i=1, SIZE(fms_diag_object%FMS_diag_fields)
1375  write(unit_num, *) 'Field num:', i
1376  fieldptr => fms_diag_object%FMS_diag_fields(i)
1377  call fieldptr%dump_field_obj(unit_num)
1378  enddo
1379  else
1380  write(unit_num, *) 'fields not initialized'
1381  endif
1382  if( present(filename) ) close(unit_num)
1383  endif
1384 #else
1385  call mpp_error( fatal, "You can not use the modern diag manager without compiling with -Duse_yaml")
1386 #endif
1387 end subroutine
1388 
1389 !> @brief Allocates the output buffers of the fields corresponding to the registered variable
1390 !! Input arguments are the field and its ID passed to routine fms_diag_accept_data()
1391 subroutine allocate_diag_field_output_buffers(this, field_data, field_id)
1392  class(fmsdiagobject_type), target, intent(inout) :: this !< diag object
1393  class(*), dimension(:,:,:,:), intent(in) :: field_data !< field data
1394  integer, intent(in) :: field_id !< Id of the field data
1395 #ifdef use_yaml
1396  integer :: ndims !< Number of dimensions in the input field data
1397  integer :: buffer_id !< Buffer index of FMS_diag_buffers
1398  integer :: num_diurnal_samples !< Number of diurnal samples from diag_yaml
1399  integer :: axes_length(4) !< Length of each axis
1400  integer :: i, j !< For looping
1401  class(fmsdiagoutputbuffer_type), pointer :: ptr_diag_buffer_obj !< Pointer to the buffer class
1402  class(diagyamlfilesvar_type), pointer :: ptr_diag_field_yaml !< Pointer to a field from yaml fields
1403  integer, pointer :: axis_ids(:) !< Pointer to indices of axes of the field variable
1404  integer :: var_type !< Stores type of the field data (r4, r8, i4, i8, and string) represented as an integer.
1405  character(len=:), allocatable :: var_name !< Field name to initialize output buffers
1406  logical :: is_scalar !< Flag indicating that the variable is a scalar
1407  integer :: yaml_id !< Yaml id for the buffer
1408  integer :: file_id !< File id for the buffer
1409 
1410  if (this%FMS_diag_fields(field_id)%buffer_allocated) return
1411 
1412  ! Determine the type of the field data
1413  var_type = get_var_type(field_data(1, 1, 1, 1))
1414 
1415  ! Get variable/field name
1416  var_name = this%FMS_diag_fields(field_id)%get_varname()
1417 
1418  ! Determine dimensions of the field
1419  is_scalar = this%FMS_diag_fields(field_id)%is_scalar()
1420 
1421  ! Loop over a number of fields/buffers where this variable occurs
1422  do i = 1, size(this%FMS_diag_fields(field_id)%buffer_ids)
1423  buffer_id = this%FMS_diag_fields(field_id)%buffer_ids(i)
1424  file_id = this%FMS_diag_fields(field_id)%file_ids(i)
1425 
1426  !< Go away if the file is a subregional file and the current PE does not have any data for it
1427  if (.not. this%FMS_diag_files(file_id)%writing_on_this_pe()) cycle
1428 
1429  ndims = 0
1430  if (.not. is_scalar) then
1431  call this%FMS_diag_output_buffers(buffer_id)%get_axis_ids(axis_ids)
1432  ndims = size(axis_ids)
1433  endif
1434 
1435  yaml_id = this%FMS_diag_output_buffers(buffer_id)%get_yaml_id()
1436 
1437  ptr_diag_field_yaml => diag_yaml%diag_fields(yaml_id)
1438  num_diurnal_samples = ptr_diag_field_yaml%get_n_diurnal() !< Get number of diurnal samples
1439 
1440  axes_length = 1
1441  do j = 1, ndims
1442  axes_length(j) = this%fms_get_axis_length(axis_ids(j))
1443  enddo
1444 
1445  if (num_diurnal_samples .ne. 0) then
1446  ndims = ndims + 1 !< Add one more dimension for the diurnal axis
1447  endif
1448 
1449  ptr_diag_buffer_obj => this%FMS_diag_output_buffers(buffer_id)
1450  call ptr_diag_buffer_obj%allocate_buffer(field_data(1, 1, 1, 1), ndims, axes_length(1:4), &
1451  this%FMS_diag_fields(field_id)%get_mask_variant(), var_name, num_diurnal_samples)
1452  call ptr_diag_buffer_obj%initialize_buffer(ptr_diag_field_yaml%get_var_reduction(), var_name)
1453 
1454  enddo
1455  nullify(axis_ids)
1456 
1457  this%FMS_diag_fields(field_id)%buffer_allocated = .true.
1458 #else
1459  call mpp_error( fatal, "allocate_diag_field_output_buffers: "//&
1460  "you can not use the modern diag manager without compiling with -Duse_yaml")
1461 #endif
1462 end subroutine allocate_diag_field_output_buffers
1463 
1464 !> @brief Determines if the window defined by the input bounds is a physics window.
1465 !> @return TRUE if the window size is less then the actual field size else FALSE.
1466 function fms_diag_compare_window(this, field, field_id, &
1467  is_in, ie_in, js_in, je_in, ks_in, ke_in) result(is_phys_win)
1468  class(fmsdiagobject_type), intent(in) :: this !< Diag Object
1469  class(*), intent(in) :: field(:,:,:,:) !< Field data
1470  integer, intent(in) :: field_id !< ID of the input field
1471  integer, intent(in) :: is_in, js_in !< Starting field indices for the first 2 dimensions;
1472  !< pass reconditioned indices fis and fjs
1473  !< which are computed elsewhere.
1474  integer, intent(in) :: ie_in, je_in !< Ending field indices for the first 2 dimensions;
1475  !< pass reconditioned indices fie and fje
1476  !< which are computed elsewhere.
1477  integer, intent(in) :: ks_in, ke_in !< Starting and ending indices of the field in 3rd dimension
1478  logical :: is_phys_win !< Return flag
1479 #ifdef use_yaml
1480  integer, pointer :: axis_ids(:)
1481  integer :: total_elements
1482  integer :: i !< For do loop
1483  integer :: field_size
1484  integer, allocatable :: field_shape(:) !< Shape of the field data
1485  integer :: window_size
1486 
1487  !> Determine shape of the field defined by the input bounds
1488  field_shape = shape(field(is_in:ie_in, js_in:je_in, ks_in:ke_in, :))
1489 
1490  window_size = field_shape(1) * field_shape(2) * field_shape(3)
1491 
1492  total_elements = 1
1493  axis_ids => this%FMS_diag_fields(field_id)%get_axis_id()
1494  do i=1, size(axis_ids)
1495  total_elements = total_elements * this%fms_get_axis_length(axis_ids(i))
1496  enddo
1497 
1498  if (total_elements > window_size) then
1499  is_phys_win = .true.
1500  else
1501  is_phys_win = .false.
1502  end if
1503 #else
1504  is_phys_win = .false.
1505  call mpp_error( fatal, "fms_diag_compare_window: "//&
1506  "you can not use the modern diag manager without compiling with -Duse_yaml")
1507 #endif
1508 end function fms_diag_compare_window
1509 
1510 end module fms_diag_object_mod
integer function get_var_type(var)
gets the type of a variable
Definition: diag_data.F90:599
integer, parameter no_domain
Use the FmsNetcdfFile_t fileobj.
Definition: diag_data.F90:101
type(time_type) function get_base_time()
gets the module variable base_time
Definition: diag_data.F90:511
integer max_axes
Maximum number of independent axes.
Definition: diag_data.F90:361
integer, parameter diag_field_not_found
Return value for a diag_field that isn't found in the diag_table.
Definition: diag_data.F90:112
type(time_type) diag_init_time
Time diag_manager_init called. If init_time not included in diag_manager_init call,...
Definition: diag_data.F90:417
integer, parameter time_min
The reduction method is min value.
Definition: diag_data.F90:117
integer, parameter time_diurnal
The reduction method is diurnal.
Definition: diag_data.F90:122
integer, parameter time_power
The reduction method is average with exponents.
Definition: diag_data.F90:123
integer, parameter time_average
The reduction method is average of values.
Definition: diag_data.F90:120
integer, parameter time_sum
The reduction method is sum of values.
Definition: diag_data.F90:119
integer, parameter time_rms
The reudction method is root mean square of values.
Definition: diag_data.F90:121
integer, parameter time_none
There is no reduction method.
Definition: diag_data.F90:116
integer, parameter time_max
The reduction method is max value.
Definition: diag_data.F90:118
integer, parameter r8
Supported type/kind of the variable.
Definition: diag_data.F90:84
integer, parameter two_d_domain
Use the FmsNetcdfDomainFile_t fileobj.
Definition: diag_data.F90:102
logical pure function, public determine_if_block_is_in_region(subregion_start, subregion_end, bounds, dim)
The PEs grid points are divided further into "blocks". This function determines if a block.
pure real(kind=r8_kind) function, public set_weight(weight)
Sets the weight based on the weight passed into send_data (1.0_r8_kind if the weight is not passed in...
logical function, dimension(:,:,:,:), allocatable, public init_mask(rmask, mask, field)
Sets the logical mask based on mask or rmask.
pure character(len=128) function, public check_indices_order(is_in, ie_in, js_in, je_in)
Checks improper combinations of is, ie, js, and je.
subroutine, public diag_yaml_object_end()
Destroys the diag_yaml object.
logical function, public fms_diag_axis_object_end(axis_array)
pure integer function, public get_axis_id_from_name(axis_name, diag_axis, naxis, set_name)
integer function, dimension(:), allocatable, public get_diag_field_ids(indices)
Gets field indices corresponding to the indices (input argument) in the sorted variable_list.
subroutine, public diag_yaml_object_init(diag_subset_output)
Uses the yaml_parser_mod to read in the diag_table and fill in the diag_yaml object.
logical function, public fms_diag_axis_object_init(axis_array)
subroutine, public fms_diag_yaml_out()
Writes an output yaml with all available information on the written files. Will only write with root ...
integer function, dimension(:), allocatable, public find_diag_field(diag_field_name, module_name)
Determines if a diag_field is in the diag_yaml_object.
pure character(len=120) function, dimension(2), public parse_compress_att(compress_att)
integer function, dimension(:), allocatable, public get_diag_files_id(indices)
Finds the indices of the diag_yamldiag_files(:) corresponding to fields in variable_list(indices)
subroutine, public get_domain_and_domain_type(diag_axis, axis_id, domain_type, domain, var_name)
Loop through a variable's axis_id to determine and return the domain type and domain to use.
logical function, public fms_error_handler(routine, message, err_msg)
Facilitates the control of fatal error conditions.
Definition: fms.F90:525
character(:) function, allocatable, public string(v, fmt)
Converts a number or a Boolean value to a string.
One dimensional domain used to manage shared data access between pes.
The domain2D type contains all the necessary information to define the global, compute and data domai...
Domain information for managing data on unstructured grids.
integer function stdout()
This function returns the current standard fortran unit numbers for output.
Definition: mpp_util.inc:43
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:407
Error handler.
Definition: mpp.F90:382
integer function, public get_ticks_per_second()
Returns the number of ticks per second.
subroutine, public get_time(Time, seconds, days, ticks, err_msg)
Returns days and seconds ( < 86400 ) corresponding to a time. err_msg should be checked for any error...
character(len=15) function, public date_to_string(time, err_msg)
Get the a character string that represents the time. The format will be yyyymmdd.hhmmss.
subroutine, public get_date(time, year, month, day, hour, minute, second, tick, err_msg)
Gets the date for different calendar types. Given a time_interval, returns the corresponding date und...
Given an input date in year, month, days, etc., creates a time_type that represents this time interva...
Given some number of seconds and days, returns the corresponding time_type.
Type to represent amounts of time. Implemented as seconds and days to allow for larger intervals.
Contains buffer types and routines for the diag manager.
Type to hold the domain info for an axis This type was created to avoid having to send in "Domain",...
Type to hold the diagnostic axis description.
Type to hold the diag_axis (either subaxis or a full axis)
Type to hold the diagnostic axis description.
Data structure holding a 3D bounding box. It is commonlyused to represent the interval bounds or limi...
Object that holds all variable information.
A container for fmsDiagFile_type. This is used to create the array of files.
type to hold the info a diag_field