FMS  2024.03
Flexible Modeling System
fms_diag_yaml.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 
20 !> @defgroup fms_diag_yaml_mod fms_diag_yaml_mod
21 !> @ingroup diag_manager
22 !! @brief fms_diag_yaml_mod is an integral part of
23 !! diag_manager_mod. Its function is to read the diag_table.yaml to fill in
24 !! the diag_yaml_object
25 
26 !> @file
27 !> @brief File for @ref diag_yaml_mod
28 
29 !> @addtogroup fms_diag_yaml_mod
30 !> @{
31 module fms_diag_yaml_mod
32 #ifdef use_yaml
33 use diag_data_mod, only: diag_null, diag_ocean, diag_all, diag_other, set_base_time, latlon_gridtype, &
34  index_gridtype, null_gridtype, diag_seconds, diag_minutes, diag_hours, diag_days, &
35  diag_months, diag_years, time_average, time_rms, time_max, time_min, time_sum, &
36  time_diurnal, time_power, time_none, r8, i8, r4, i4, diag_not_registered, &
38 use yaml_parser_mod, only: open_and_parse_file, get_value_from_key, get_num_blocks, get_nkeys, &
40 use fms_yaml_output_mod, only: fmsyamloutkeys_type, fmsyamloutvalues_type, write_yaml_from_struct_3, &
41  yaml_out_add_level2key, initialize_key_struct, initialize_val_struct
42 use mpp_mod, only: mpp_error, fatal, note, mpp_pe, mpp_root_pe, stdout
43 use, intrinsic :: iso_c_binding, only : c_ptr, c_null_char
44 use fms_string_utils_mod, only: fms_array_to_pointer, fms_find_my_string, fms_sort_this, fms_find_unique, string, &
46 use platform_mod, only: r4_kind, i4_kind, r8_kind, i8_kind, fms_file_len
47 use fms_mod, only: lowercase
48 
49 implicit none
50 
51 private
52 
53 public :: diag_yaml
58 public :: get_diag_field_ids
59 public :: dump_diag_yaml_obj
60 public :: fms_diag_yaml_out
61 public :: max_subaxes
62 !> @}
63 
64 integer, parameter :: basedate_size = 6
65 integer, parameter :: NUM_SUB_REGION_ARRAY = 8
66 integer, parameter :: MAX_FREQ = 12
67 integer :: MAX_SUBAXES = 0 !< Max number of subaxis, set in diag_yaml_object_init depending on
68  !! what is in the diag yaml
69 
70 
71 !> @brief type to hold an array of sorted diag_fiels
73  character(len=255), allocatable :: var_name(:) !< Array of diag_field
74  type(c_ptr), allocatable :: var_pointer(:) !< Array of pointers
75  integer, allocatable :: diag_field_indices(:) !< Index of the field in the diag_field array
76 end type
77 
78 !> @brief type to hold an array of sorted diag_files
80  character(len=FMS_FILE_LEN), allocatable :: file_name(:) !< Array of diag_field
81  type(c_ptr), allocatable :: file_pointer(:) !< Array of pointers
82  integer, allocatable :: diag_file_indices(:) !< Index of the file in the diag_file array
83 end type
84 
85 !> @brief type to hold the sub region information about a file
87  INTEGER :: grid_type !< Flag indicating the type of region,
88  !! acceptable values are latlon_gridtype, index_gridtype,
89  !! null_gridtype
90  class(*), allocatable :: corners(:,:)!< (x, y) coordinates of the four corner of the region
91  integer :: tile !< Tile number of the sub region
92  !! required if using the "index" grid type
93 
94 end type subregion_type
95 
96 !> @brief type to hold the diag_file information
98  private
99  character (len=:), allocatable :: file_fname !< file name
100  integer :: file_frequnit(max_freq) !< the frequency unit (DIAG_SECONDS,
101  !! DIAG_MINUTES, DIAG_HOURS, DIAG_DAYS,
102  !! DIAG_YEARS)
103  integer :: file_freq(max_freq) !< the frequency of data
104  integer :: file_timeunit !< The unit of time (DIAG_SECONDS,
105  !! DIAG_MINUTES, DIAG_HOURS, DIAG_DAYS,
106  !! DIAG_YEARS)
107  character (len=:), allocatable :: file_unlimdim !< The name of the unlimited dimension
108  type(subregion_type) :: file_sub_region !< type containing info about the subregion
109  integer :: file_new_file_freq(max_freq) !< Frequency for closing the existing file
110  integer :: file_new_file_freq_units(max_freq) !< Time units for creating a new file.
111  !! Required if “new_file_freq” used
112  !! (DIAG_SECONDS, DIAG_MINUTES, &
113  !! DIAG_HOURS, DIAG_DAYS, DIAG_YEARS)
114  character (len=:), allocatable :: file_start_time !< Time to start the file for the
115  !! first time. Requires “new_file_freq”
116  integer :: filename_time !< The time to use when setting the name of
117  !! new files: begin, middle, or end of the
118  !! time_bounds
119  integer :: file_duration(max_freq) !< How long the file should receive data
120  !! after start time in file_duration_units.
121  !! This optional field can only be used if
122  !! the start_time field is present.  If this
123  !! field is absent, then the file duration
124  !! will be equal to the frequency for
125  !! creating new files. NOTE: The
126  !! file_duration_units field must also
127  !! be present if this field is present.
128  integer :: file_duration_units(max_freq) !< The file duration units
129  !! (DIAG_SECONDS, DIAG_MINUTES, &
130  !! DIAG_HOURS, DIAG_DAYS, DIAG_YEARS)
131  integer :: current_new_file_freq_index !< The index of the new_file_freq array
132  !< Need to use `MAX_STR_LEN` because not all filenames/global attributes are the same length
133  character (len=MAX_STR_LEN), allocatable :: file_varlist(:) !< An array of variable names
134  !! within a file
135  character (len=MAX_STR_LEN), allocatable :: file_outlist(:) !< An array of variable output names
136  !! within a file, used to distinguish
137  !! varlist names for yaml output
138  character (len=MAX_STR_LEN), allocatable :: file_global_meta(:,:) !< Array of key(dim=1)
139  !! and values(dim=2) to be
140  !! added as global meta data to
141  !! the file
142  character (len=:), allocatable :: default_var_precision !< The precision for all of the variables in the file
143  !! This may be overriden if the precison was defined
144  !! at the variable level
145  character (len=:), allocatable :: default_var_reduction !< The reduction for all of the variables in the file
146  !! This may be overriden if the reduction was defined at
147  !! the variable level
148  character (len=:), allocatable :: default_var_module !< The module for all of the variables in the file
149  !! This may be overriden if the modules was defined at the
150  !! variable level
151  contains
152 
153  !> All getter functions (functions named get_x(), for member field named x)
154  !! return copies of the member variables unless explicitly noted.
155  procedure, public :: size_file_varlist
156  procedure, public :: get_file_fname
157  procedure, public :: get_file_frequnit
158  procedure, public :: get_file_freq
159  procedure, public :: get_file_timeunit
160  procedure, public :: get_file_unlimdim
161  procedure, public :: get_file_sub_region
162  procedure, public :: get_file_new_file_freq
163  procedure, public :: get_file_new_file_freq_units
164  procedure, public :: get_file_start_time
165  procedure, public :: get_file_duration
166  procedure, public :: get_file_duration_units
167  procedure, public :: get_file_varlist
168  procedure, public :: get_file_global_meta
169  procedure, public :: get_filename_time
170  procedure, public :: is_global_meta
171  !> Has functions to determine if allocatable variables are true. If a variable is not an allocatable
172  !! then is will always return .true.
173  procedure, public :: has_file_fname
174  procedure, public :: has_file_frequnit
175  procedure, public :: has_file_freq
176  procedure, public :: has_file_timeunit
177  procedure, public :: has_file_unlimdim
178  procedure, public :: has_file_sub_region
179  procedure, public :: has_file_new_file_freq
180  procedure, public :: has_file_new_file_freq_units
181  procedure, public :: has_file_start_time
182  procedure, public :: has_file_duration
183  procedure, public :: has_file_duration_units
184  procedure, public :: has_file_varlist
185  procedure, public :: has_file_global_meta
186  procedure, public :: increase_new_file_freq_index
187 end type diagyamlfiles_type
188 
189 !> @brief type to hold the info a diag_field
191  character (len=:), private, allocatable :: var_fname !< The field/diagnostic name
192  character (len=:), private, allocatable :: var_varname !< The name of the variable
193  integer , private, allocatable :: var_reduction !< Reduction to be done on var
194  !! time_average, time_rms, time_max,
195  !! time_min, time_sum, time_diurnal, time_power
196  character (len=:), private, allocatable :: var_module !< The module that th variable is in
197  integer , private, allocatable :: var_kind !< The type/kind of the variable
198  character (len=:), private, allocatable :: var_outname !< Name of the variable as written to the file
199  character (len=:), private, allocatable :: var_longname !< Overwrites the long name of the variable
200  character (len=:), private, allocatable :: var_units !< Overwrites the units
201  real(kind=r4_kind), private :: var_zbounds(2) !< The z axis limits [vert_min, vert_max]
202  integer , private :: n_diurnal !< Number of diurnal samples
203  !! 0 if var_reduction is not "diurnalXX"
204  integer , private :: pow_value !< The power value
205  !! 0 if pow_value is not "powXX"
206  logical , private :: var_file_is_subregional !< true if the file this entry
207  !! belongs to is subregional
208 
209  !< Need to use `MAX_STR_LEN` because not all filenames/global attributes are the same length
210  character (len=MAX_STR_LEN), dimension (:, :), private, allocatable :: var_attributes !< Attributes to overwrite or
211  !! add from diag_yaml
212  character(len=:), allocatable :: var_axes_names !< list of axes names
213  contains
214  !> All getter functions (functions named get_x(), for member field named x)
215  !! return copies of the member variables unless explicitly noted.
216  procedure :: get_var_fname
217  procedure :: get_var_varname
218  procedure :: get_var_reduction
219  procedure :: get_var_module
220  procedure :: get_var_kind
221  procedure :: get_var_outname
222  procedure :: get_var_longname
223  procedure :: get_var_units
224  procedure :: get_var_zbounds
225  procedure :: get_var_attributes
226  procedure :: get_n_diurnal
227  procedure :: get_pow_value
228  procedure :: is_var_attributes
229 
230  procedure :: has_var_fname
231  procedure :: has_var_varname
232  procedure :: has_var_reduction
233  procedure :: has_var_module
234  procedure :: has_var_kind
235  procedure :: has_var_outname
236  procedure :: has_var_longname
237  procedure :: has_var_units
238  procedure :: has_var_zbounds
239  procedure :: has_var_attributes
240  procedure :: has_n_diurnal
241  procedure :: has_pow_value
242  procedure :: add_axis_name
243  procedure :: is_file_subregional
244 
245 end type diagyamlfilesvar_type
246 
247 !> @brief Object that holds the information of the diag_yaml
248 !> @ingroup fms_diag_yaml_mod
250  character(len=:), allocatable, private :: diag_title !< Experiment name
251  integer, private, dimension (basedate_size) :: diag_basedate !< basedate array
252  type(diagyamlfiles_type), allocatable, public, dimension (:) :: diag_files!< History file info
253  type(diagyamlfilesvar_type), allocatable, public, dimension (:) :: diag_fields !< Diag fields info
254  contains
255  procedure :: size_diag_files
256 
257  procedure :: get_title !< Returns the title
258  procedure :: get_basedate !< Returns the basedate array
259  procedure :: get_diag_files !< Returns the diag_files array
260  procedure :: get_diag_fields !< Returns the diag_field array
261  procedure :: get_diag_field_from_id
262 
263  procedure :: has_diag_title
264  procedure :: has_diag_basedate
265  procedure :: has_diag_files
266  procedure :: has_diag_fields
267 
268 end type diagyamlobject_type
269 
270 type (diagyamlobject_type), target :: diag_yaml !< Obj containing the contents of the diag_table.yaml
271 type (varlist_type), save :: variable_list !< List of all the variables in the diag_table.yaml
272 type (filelist_type), save :: file_list !< List of all files in the diag_table.yaml
273 
274 logical, private :: diag_yaml_module_initialized = .false.
275 
276 
277 !> @addtogroup fms_diag_yaml_mod
278 !> @{
279 contains
280 
281 !> @brief gets the diag_yaml module variable
282 !! @return a copy of the diag_yaml module variable
283 function get_diag_yaml_obj() &
284 result(res)
285  type (diagyamlobject_type), pointer :: res
286 
287  res => diag_yaml
288 end function get_diag_yaml_obj
289 
290 !> @brief get the basedate of a diag_yaml type
291 !! @return the basedate as an integer array
292 pure function get_basedate (this) &
293 result(diag_basedate)
294  class(diagyamlobject_type), intent(in) :: this !< The diag_yaml
295  integer, dimension (basedate_size) :: diag_basedate !< Basedate array result to return
296 
297  diag_basedate = this%diag_basedate
298 end function get_basedate
299 
300 !> @brief Find the number of files listed in the diag yaml
301 !! @return the number of files in the diag yaml
302 pure integer function size_diag_files(this)
303  class(diagyamlobject_type), intent(in) :: this !< The diag_yaml
304  if (this%has_diag_files()) then
305  size_diag_files = size(this%diag_files)
306  else
307  size_diag_files = 0
308  endif
309 end function size_diag_files
310 
311 !> @brief get the title of a diag_yaml type
312 !! @return the title of the diag table as an allocated string
313 pure function get_title (this) &
314  result(diag_title)
315  class(diagyamlobject_type), intent(in) :: this !< The diag_yaml
316  character(len=:),allocatable :: diag_title !< Basedate array result to return
317 
318  diag_title = this%diag_title
319 end function get_title
320 
321 !> @brief get the diag_files of a diag_yaml type
322 !! @return the diag_files
323 function get_diag_files(this) &
324 result(diag_files)
325  class(diagyamlobject_type), intent(in) :: this !< The diag_yaml
326  type(diagyamlfiles_type), allocatable, dimension (:) :: diag_files !< History file info
327 
328  diag_files = this%diag_files
329 end function get_diag_files
330 
331 !> @brief Get the diag_field yaml corresponding to a yaml_id
332 !! @return Pointer to the diag_field yaml entry
333 function get_diag_field_from_id(this, yaml_id) &
334  result(diag_field)
335  class(diagyamlobject_type), target, intent(in) :: this !< The diag_yaml
336  integer, intent(in) :: yaml_id !< Yaml id
337 
338  type(diagyamlfilesvar_type), pointer :: diag_field !< Diag fields info
339 
340  if (yaml_id .eq. diag_not_registered) call mpp_error(fatal, &
341  "Diag_manager: The yaml id for this field is not is not set")
342 
343  diag_field => this%diag_fields(variable_list%diag_field_indices(yaml_id))
344 
345 end function get_diag_field_from_id
346 
347 !> @brief get the diag_fields of a diag_yaml type
348 !! @return the diag_fields
349 pure function get_diag_fields(this) &
350 result(diag_fields)
351  class(diagyamlobject_type), intent(in) :: this !< The diag_yaml
352  type(diagyamlfilesvar_type), allocatable, dimension (:) :: diag_fields !< Diag fields info
353 
354  diag_fields = this%diag_fields
355 end function get_diag_fields
356 
357 !> @brief Uses the yaml_parser_mod to read in the diag_table and fill in the
358 !! diag_yaml object
359 subroutine diag_yaml_object_init(diag_subset_output)
360  integer, intent(in) :: diag_subset_output !< DIAG_ALL - Current PE is in the one and only pelist
361  !! DIAG_OTHER - Current PE is not in the ocean pelist
362  !! and there are multiple pelists
363  !! DIAG_OCEAN - Current PE is in the ocean pelist
364  !! and there are multiple pelists
365  integer :: diag_yaml_id !< Id for the diag_table yaml
366  integer :: nfiles !< Number of files in the diag_table yaml
367  integer, allocatable :: diag_file_ids(:) !< Ids of the files in the diag_table yaml
368  integer :: i, j !< For do loops
369  integer :: total_nvars !< The total number of variables in the diag_table yaml
370  integer :: var_count !< The current number of variables added to the diag_yaml obj
371  integer :: file_var_count !< The current number of variables added in the diag_file
372  integer :: nvars !< The number of variables in the current file
373  integer, allocatable :: var_ids(:) !< Ids of the variables in diag_table yaml
374  logical :: is_ocean !< Flag indicating if it is an ocean file
375  logical, allocatable :: ignore(:) !< Flag indicating if the diag_file is going to be ignored
376  integer :: actual_num_files !< The actual number of files that were saved
377  integer :: file_count !! The current number of files added to the diag_yaml obj
378  logical :: write_file !< Flag indicating if the user wants the file to be written
379  logical :: write_var !< Flag indicating if the user wants the variable to be written
380  logical :: allow_averages !< .True. if averages are allowed (the file is not static of you are
381  !! outputing data at every frequency)
382  character(len=:), allocatable :: filename!< Diag file name (for error messages)
383  logical :: is_instantaneous !< .True. if the file is instantaneous (i.e no averaging)
384 
385  if (diag_yaml_module_initialized) return
386 
387  diag_yaml_id = open_and_parse_file("diag_table.yaml")
388 
389  call diag_get_value_from_key(diag_yaml_id, 0, "title", diag_yaml%diag_title)
390  call get_value_from_key(diag_yaml_id, 0, "base_date", diag_yaml%diag_basedate)
391  call set_base_time(diag_yaml%diag_basedate)
392 
393  nfiles = get_num_blocks(diag_yaml_id, "diag_files")
394  allocate(diag_file_ids(nfiles))
395  allocate(ignore(nfiles))
396 
397  call get_block_ids(diag_yaml_id, "diag_files", diag_file_ids)
398 
399  ignore = .false.
400  total_nvars = 0
401  !< If you are on two seperate pelists
402  if(diag_subset_output .ne. diag_all) then
403  do i = 1, nfiles
404  is_ocean = .false.
405  call get_value_from_key(diag_yaml_id, diag_file_ids(i), "is_ocean", is_ocean, is_optional=.true.)
406  !< If you are on the ocean pelist and the file is not an ocean file, skip the file
407  if (diag_subset_output .eq. diag_ocean .and. .not. is_ocean) ignore(i) = .true.
408 
409  !< If you are not on the ocean pelist and the file is ocean, skip the file
410  if(diag_subset_output .eq. diag_other .and. is_ocean) ignore(i) = .true.
411  enddo
412  endif
413 
414  !< Determine how many files are in the diag_yaml, ignoring those with write_file = False
415  actual_num_files = 0
416  do i = 1, nfiles
417  write_file = .true.
418  call get_value_from_key(diag_yaml_id, diag_file_ids(i), "write_file", write_file, is_optional=.true.)
419  if(.not. write_file) ignore(i) = .true.
420 
421  !< If ignoring the file, ignore the fields in that file too!
422  if (.not. ignore(i)) then
423  nvars = get_total_num_vars(diag_yaml_id, diag_file_ids(i))
424  total_nvars = total_nvars + nvars
425  if (nvars .ne. 0) then
426  actual_num_files = actual_num_files + 1
427  else
428  call diag_get_value_from_key(diag_yaml_id, diag_file_ids(i), "file_name", filename)
429  call mpp_error(note, "diag_manager_mod:: the file:"//trim(filename)//" has no variables defined. Ignoring!")
430  if (allocated(filename)) deallocate(filename)
431  ignore(i) = .true.
432  endif
433  endif
434  enddo
435 
436  allocate(diag_yaml%diag_files(actual_num_files))
437  allocate(diag_yaml%diag_fields(total_nvars))
438  allocate(variable_list%var_name(total_nvars))
439  allocate(variable_list%diag_field_indices(total_nvars))
440  allocate(file_list%file_name(actual_num_files))
441  allocate(file_list%diag_file_indices(actual_num_files))
442 
443  var_count = 0
444  file_count = 0
445  !> Loop through the number of nfiles and fill in the diag_yaml obj
446  nfiles_loop: do i = 1, nfiles
447  if(ignore(i)) cycle
448  file_count = file_count + 1
449  call diag_yaml_files_obj_init(diag_yaml%diag_files(file_count))
450  call fill_in_diag_files(diag_yaml_id, diag_file_ids(i), diag_yaml%diag_files(file_count))
451 
452  !> Save the file name in the file_list
453  file_list%file_name(file_count) = trim(diag_yaml%diag_files(file_count)%file_fname)//c_null_char
454  file_list%diag_file_indices(file_count) = file_count
455 
456  nvars = 0
457  nvars = get_num_blocks(diag_yaml_id, "varlist", parent_block_id=diag_file_ids(i))
458  allocate(var_ids(nvars))
459  call get_block_ids(diag_yaml_id, "varlist", var_ids, parent_block_id=diag_file_ids(i))
460  file_var_count = 0
461  allocate(diag_yaml%diag_files(file_count)%file_varlist(get_total_num_vars(diag_yaml_id, diag_file_ids(i))))
462  allocate(diag_yaml%diag_files(file_count)%file_outlist(get_total_num_vars(diag_yaml_id, diag_file_ids(i))))
463  allow_averages = .not. diag_yaml%diag_files(file_count)%file_freq(1) < 1
464  is_instantaneous = .false.
465  nvars_loop: do j = 1, nvars
466  write_var = .true.
467  call get_value_from_key(diag_yaml_id, var_ids(j), "write_var", write_var, is_optional=.true.)
468  if (.not. write_var) cycle
469 
470  var_count = var_count + 1
471  file_var_count = file_var_count + 1
472 
473  !> Save the filename in the diag_field type
474  diag_yaml%diag_fields(var_count)%var_fname = diag_yaml%diag_files(file_count)%file_fname
475 
476  !> initialize axes string
477  diag_yaml%diag_fields(var_count)%var_axes_names = ""
478  diag_yaml%diag_fields(var_count)%var_file_is_subregional = diag_yaml%diag_files(file_count)%has_file_sub_region()
479 
480  call fill_in_diag_fields(diag_yaml_id, diag_yaml%diag_files(file_count), var_ids(j), &
481  diag_yaml%diag_fields(var_count), allow_averages)
482 
483  !> Save the variable name in the diag_file type
484  diag_yaml%diag_files(file_count)%file_varlist(file_var_count) = diag_yaml%diag_fields(var_count)%var_varname
485  if(diag_yaml%diag_fields(var_count)%has_var_outname()) then
486  diag_yaml%diag_files(file_count)%file_outlist(file_var_count) = diag_yaml%diag_fields(var_count)%var_outname
487  else
488  diag_yaml%diag_files(file_count)%file_outlist(file_var_count) = ""
489  endif
490 
491  !> Save the variable name and the module name in the variable_list
492  variable_list%var_name(var_count) = trim(diag_yaml%diag_fields(var_count)%var_varname)//&
493  ":"//trim(diag_yaml%diag_fields(var_count)%var_module)//c_null_char
494  !! The diag_table is not case sensitive (so we are saving it as lowercase)
495  variable_list%var_name(var_count) = lowercase(variable_list%var_name(var_count))
496  variable_list%diag_field_indices(var_count) = var_count
497  enddo nvars_loop
498  deallocate(var_ids)
499  enddo nfiles_loop
500 
501  !> Sort the file list in alphabetical order
502  file_list%file_pointer = fms_array_to_pointer(file_list%file_name)
503  call fms_sort_this(file_list%file_pointer, actual_num_files, file_list%diag_file_indices)
504 
505  variable_list%var_pointer = fms_array_to_pointer(variable_list%var_name)
506  call fms_sort_this(variable_list%var_pointer, total_nvars, variable_list%diag_field_indices)
507 
508  deallocate(diag_file_ids)
509  diag_yaml_module_initialized = .true.
510 end subroutine
511 
512 !> @brief Destroys the diag_yaml object
514  integer :: i !< For do loops
515 
516  do i = 1, size(diag_yaml%diag_files, 1)
517  if(allocated(diag_yaml%diag_files(i)%file_varlist)) deallocate(diag_yaml%diag_files(i)%file_varlist)
518  if(allocated(diag_yaml%diag_files(i)%file_outlist)) deallocate(diag_yaml%diag_files(i)%file_outlist)
519  if(allocated(diag_yaml%diag_files(i)%file_global_meta)) deallocate(diag_yaml%diag_files(i)%file_global_meta)
520  if(allocated(diag_yaml%diag_files(i)%file_sub_region%corners)) &
521  deallocate(diag_yaml%diag_files(i)%file_sub_region%corners)
522  enddo
523  if(allocated(diag_yaml%diag_files)) deallocate(diag_yaml%diag_files)
524 
525  do i = 1, size(diag_yaml%diag_fields, 1)
526  if(allocated(diag_yaml%diag_fields(i)%var_attributes)) deallocate(diag_yaml%diag_fields(i)%var_attributes)
527  enddo
528  if(allocated(diag_yaml%diag_fields)) deallocate(diag_yaml%diag_fields)
529 
530  if(allocated(file_list%file_pointer)) deallocate(file_list%file_pointer)
531  if(allocated(file_list%file_name)) deallocate(file_list%file_name)
532  if(allocated(file_list%diag_file_indices)) deallocate(file_list%diag_file_indices)
533 
534  if(allocated(variable_list%var_pointer)) deallocate(variable_list%var_pointer)
535  if(allocated(variable_list%var_name)) deallocate(variable_list%var_name)
536  if(allocated(variable_list%diag_field_indices)) deallocate(variable_list%diag_field_indices)
537 
538 end subroutine diag_yaml_object_end
539 
540 !> @brief Fills in a diagYamlFiles_type with the contents of a file block in diag_table.yaml
541 subroutine fill_in_diag_files(diag_yaml_id, diag_file_id, yaml_fileobj)
542  integer, intent(in) :: diag_yaml_id !< Id of the diag_table.yaml
543  integer, intent(in) :: diag_file_id !< Id of the file block to read
544  type(diagyamlfiles_type), intent(inout) :: yaml_fileobj !< diagYamlFiles_type obj to read the contents into
545 
546  integer :: nsubregion !< Flag indicating of there any regions (0 or 1)
547  integer :: sub_region_id(1) !< Id of the sub_region block
548  integer :: natt !< Number of global attributes in the current file
549  integer :: global_att_id(1) !< Id of the global attributes block
550  integer :: nkeys !< Number of key/value global attributes pair
551  integer :: j !< For do loops
552 
553  integer, allocatable :: key_ids(:) !< Id of the gloabl atttributes key/value pairs
554  character(len=:), ALLOCATABLE :: grid_type !< grid_type as it is read in from the yaml
555  character(len=:), ALLOCATABLE :: buffer !< buffer to store any *_units as it is read from the yaml
556 
557  yaml_fileobj%file_frequnit = 0
558 
559  call diag_get_value_from_key(diag_yaml_id, diag_file_id, "file_name", yaml_fileobj%file_fname)
560  call diag_get_value_from_key(diag_yaml_id, diag_file_id, "freq", buffer)
561  call parse_key(yaml_fileobj%file_fname, buffer, yaml_fileobj%file_freq, yaml_fileobj%file_frequnit, "freq")
562  deallocate(buffer)
563 
564  call diag_get_value_from_key(diag_yaml_id, diag_file_id, "unlimdim", yaml_fileobj%file_unlimdim)
565  call diag_get_value_from_key(diag_yaml_id, diag_file_id, "time_units", buffer)
566  call set_file_time_units(yaml_fileobj, buffer)
567  deallocate(buffer)
568 
569  call diag_get_value_from_key(diag_yaml_id, diag_file_id, "new_file_freq", buffer, is_optional=.true.)
570  call parse_key(yaml_fileobj%file_fname, buffer, yaml_fileobj%file_new_file_freq, &
571  yaml_fileobj%file_new_file_freq_units, "new_file_freq")
572  deallocate(buffer)
573 
574  call diag_get_value_from_key(diag_yaml_id, diag_file_id, "filename_time", buffer, is_optional=.true.)
575  call set_filename_time(yaml_fileobj, buffer)
576  deallocate(buffer)
577 
578  call diag_get_value_from_key(diag_yaml_id, diag_file_id, "start_time", &
579  yaml_fileobj%file_start_time, is_optional=.true.)
580  call diag_get_value_from_key(diag_yaml_id, diag_file_id, "file_duration", buffer, is_optional=.true.)
581  call parse_key(yaml_fileobj%file_fname, buffer, yaml_fileobj%file_duration, yaml_fileobj%file_duration_units, &
582  "file_duration")
583 
584  nsubregion = 0
585  nsubregion = get_num_blocks(diag_yaml_id, "sub_region", parent_block_id=diag_file_id)
586  if (nsubregion .eq. 1) then
587  max_subaxes = max_subaxes + 1
588  call get_block_ids(diag_yaml_id, "sub_region", sub_region_id, parent_block_id=diag_file_id)
589  call diag_get_value_from_key(diag_yaml_id, sub_region_id(1), "grid_type", grid_type)
590  call get_sub_region(diag_yaml_id, sub_region_id(1), yaml_fileobj%file_sub_region, grid_type, &
591  yaml_fileobj%file_fname)
592  elseif (nsubregion .eq. 0) then
593  yaml_fileobj%file_sub_region%grid_type = null_gridtype
594  else
595  call mpp_error(fatal, "diag_yaml_object_init: file "//trim(yaml_fileobj%file_fname)//" has multiple region blocks")
596  endif
597 
598  natt = 0
599  natt = get_num_blocks(diag_yaml_id, "global_meta", parent_block_id=diag_file_id)
600  if (natt .eq. 1) then
601  call get_block_ids(diag_yaml_id, "global_meta", global_att_id, parent_block_id=diag_file_id)
602  nkeys = get_nkeys(diag_yaml_id, global_att_id(1))
603  allocate(key_ids(nkeys))
604  call get_key_ids(diag_yaml_id, global_att_id(1), key_ids)
605 
606  allocate(yaml_fileobj%file_global_meta(nkeys, 2))
607  do j = 1, nkeys
608  call get_key_name(diag_yaml_id, key_ids(j), yaml_fileobj%file_global_meta(j, 1))
609  call get_key_value(diag_yaml_id, key_ids(j), yaml_fileobj%file_global_meta(j, 2))
610  enddo
611  deallocate(key_ids)
612  elseif (natt .ne. 0) then
613  call mpp_error(fatal, "diag_yaml_object_init: file "//trim(yaml_fileobj%file_fname)//&
614  &" has multiple global_meta blocks")
615  endif
616 
617  call diag_get_value_from_key(diag_yaml_id, diag_file_id, "reduction", yaml_fileobj%default_var_reduction, &
618  is_optional=.true.)
619  call diag_get_value_from_key(diag_yaml_id, diag_file_id, "kind", yaml_fileobj%default_var_precision, &
620  is_optional=.true.)
621  call diag_get_value_from_key(diag_yaml_id, diag_file_id, "module", yaml_fileobj%default_var_module, &
622  is_optional=.true.)
623 end subroutine
624 
625 !> @brief Fills in a diagYamlFilesVar_type with the contents of a variable block in
626 !! diag_table.yaml
627 subroutine fill_in_diag_fields(diag_file_id, yaml_fileobj, var_id, field, allow_averages)
628  integer, intent(in) :: diag_file_id !< Id of the file block in the yaml file
629  type(diagyamlfiles_type), intent(in) :: yaml_fileobj !< The yaml file obj for the variables
630  integer, intent(in) :: var_id !< Id of the variable block in the yaml file
631  type(diagyamlfilesvar_type), intent(inout) :: field !< diagYamlFilesVar_type obj to read the contents into
632  logical, intent(in) :: allow_averages !< .True. if averages are allowed for this file
633 
634  integer :: natt !< Number of attributes in variable
635  integer :: var_att_id(1) !< Id of the variable attribute block
636  integer :: nkeys !< Number of key/value pairs of attributes
637  integer :: j !< For do loops
638 
639  integer, allocatable :: key_ids(:) !< Id of each attribute key/value pair
640  character(len=:), ALLOCATABLE :: buffer !< buffer to store the reduction method as it is read from the yaml
641 
642  call diag_get_value_from_key(diag_file_id, var_id, "var_name", field%var_varname)
643 
644  if (yaml_fileobj%default_var_reduction .eq. "") then
645  !! If there is no default, the reduction method is required
646  call diag_get_value_from_key(diag_file_id, var_id, "reduction", buffer)
647  else
648  call diag_get_value_from_key(diag_file_id, var_id, "reduction", buffer, is_optional=.true.)
649  !! If the reduction was not set for the variable, override it with the default
650  if (trim(buffer) .eq. "") buffer = yaml_fileobj%default_var_reduction
651  endif
652  call set_field_reduction(field, buffer)
653  deallocate(buffer)
654 
655  if (.not. allow_averages) then
656  if (field%var_reduction .ne. time_none) &
657  call mpp_error(fatal, "The file "//field%var_fname//" can only have variables that have none as "//&
658  "the reduction method because the frequency is either -1 or 0. "//&
659  "Check your diag_table.yaml for the field:"//trim(field%var_varname))
660  endif
661 
662  if (yaml_fileobj%default_var_module .eq. "") then
663  call diag_get_value_from_key(diag_file_id, var_id, "module", field%var_module)
664  else
665  call diag_get_value_from_key(diag_file_id, var_id, "module", buffer, is_optional=.true.)
666  !! If the module was set for the variable, override it with the default
667  if (trim(buffer) .eq. "") then
668  field%var_module = yaml_fileobj%default_var_module
669  else
670  field%var_module = trim(buffer)
671  endif
672  deallocate(buffer)
673  endif
674 
675  if (yaml_fileobj%default_var_precision .eq. "") then
676  !! If there is no default, the kind is required
677  call diag_get_value_from_key(diag_file_id, var_id, "kind", buffer)
678  else
679  call diag_get_value_from_key(diag_file_id, var_id, "kind", buffer, is_optional=.true.)
680  !! If the kind was set for the variable, override it with the default
681  if (trim(buffer) .eq. "") buffer = yaml_fileobj%default_var_precision
682  endif
683  call set_field_kind(field, buffer)
684 
685  call diag_get_value_from_key(diag_file_id, var_id, "output_name", field%var_outname, is_optional=.true.)
686  call diag_get_value_from_key(diag_file_id, var_id, "long_name", field%var_longname, is_optional=.true.)
687  !! VAR_UNITS !!
688 
689  natt = 0
690  natt = get_num_blocks(diag_file_id, "attributes", parent_block_id=var_id)
691  if (natt .eq. 1) then
692  call get_block_ids(diag_file_id, "attributes", var_att_id, parent_block_id=var_id)
693  nkeys = get_nkeys(diag_file_id, var_att_id(1))
694  allocate(key_ids(nkeys))
695  call get_key_ids(diag_file_id, var_att_id(1), key_ids)
696 
697  allocate(field%var_attributes(nkeys, 2))
698  do j = 1, nkeys
699  call get_key_name(diag_file_id, key_ids(j), field%var_attributes(j, 1))
700  call get_key_value(diag_file_id, key_ids(j), field%var_attributes(j, 2))
701  enddo
702  deallocate(key_ids)
703  elseif (natt .ne. 0) then
704  call mpp_error(fatal, "diag_yaml_object_init: variable "//trim(field%var_varname)//&
705  " has multiple attribute blocks")
706  endif
707 
708  !> Set the zbounds if they exist
709  field%var_zbounds = diag_null
710  call get_value_from_key(diag_file_id, var_id, "zbounds", field%var_zbounds, is_optional=.true.)
711  if (field%has_var_zbounds()) max_subaxes = max_subaxes + 1
712 end subroutine
713 
714 !> @brief diag_manager wrapper to get_value_from_key to use for allocatable
715 !! string variables
716 subroutine diag_get_value_from_key(diag_file_id, par_id, key_name, value_name, is_optional)
717  integer, intent(in) :: diag_file_id!< Id of the file block in the yaml file
718  integer, intent(in) :: par_id !< Id of the parent block in the yaml file
719  character(len=*), intent(in) :: key_name !< Key to look for in the parent block
720  character(len=:), allocatable :: value_name !< Value of the key
721  logical, intent(in), optional :: is_optional !< Flag indicating if the key is optional
722 
723  character(len=255) :: buffer !< String buffer to read in to
724 
725  buffer = "" !< Needs to be initialized for optional keys that are not present
726  call get_value_from_key(diag_file_id, par_id, trim(key_name), buffer, is_optional= is_optional)
727  allocate(character(len=len_trim(buffer)) :: value_name)
728  value_name = trim(buffer)
729 
730 end subroutine diag_get_value_from_key
731 
732 !> @brief gets the lat/lon of the sub region to use in a diag_table yaml
733 subroutine get_sub_region(diag_yaml_id, sub_region_id, sub_region, grid_type, fname)
734  integer, intent(in) :: diag_yaml_id !< Id of the diag_table yaml file
735  integer, intent(in) :: sub_region_id !< Id of the region block to read from
736  type(subregion_type),intent(inout) :: sub_region !< Type that stores the sub_region
737  character(len=*), intent(in) :: grid_type !< The grid_type as it is read from the file
738  character(len=*), intent(in) :: fname !< filename of the subregion (for error messages)
739 
740  select case (trim(grid_type))
741  case ("latlon")
742  sub_region%grid_type = latlon_gridtype
743  allocate(real(kind=r4_kind) :: sub_region%corners(4,2))
744  case ("index")
745  sub_region%grid_type = index_gridtype
746  allocate(integer(kind=i4_kind) :: sub_region%corners(4,2))
747 
748  call get_value_from_key(diag_yaml_id, sub_region_id, "tile", sub_region%tile, is_optional=.true.)
749  if (sub_region%tile .eq. diag_null) call mpp_error(fatal, &
750  "The tile number is required when defining a "//&
751  "subregion. Check your subregion entry for "//trim(fname))
752  case default
753  call mpp_error(fatal, trim(grid_type)//" is not a valid region type. &
754  &The acceptable values are latlon and index. &
755  &Check your entry for file:"//trim(fname))
756  end select
757 
758  call get_value_from_key(diag_yaml_id, sub_region_id, "corner1", sub_region%corners(1,:))
759  call get_value_from_key(diag_yaml_id, sub_region_id, "corner2", sub_region%corners(2,:))
760  call get_value_from_key(diag_yaml_id, sub_region_id, "corner3", sub_region%corners(3,:))
761  call get_value_from_key(diag_yaml_id, sub_region_id, "corner4", sub_region%corners(4,:))
762 
763 end subroutine get_sub_region
764 
765 !> @brief gets the total number of variables in the diag_table yaml file
766 !! @return total number of variables
767 function get_total_num_vars(diag_yaml_id, diag_file_id) &
768 result(total_nvars)
769 
770  integer, intent(in) :: diag_yaml_id !< Id for the diag_table yaml
771  integer, intent(in) :: diag_file_id !< Id of the file in the diag_table yaml
772  integer :: total_nvars
773 
774  integer :: i !< For do loop
775  integer :: nvars !< Number of variables in a file
776  integer, allocatable :: var_ids(:) !< Id of the variables in the file block of the yaml file
777  logical :: var_write !< Flag indicating if the user wants the variable to be written
778 
779  nvars = get_num_blocks(diag_yaml_id, "varlist", parent_block_id=diag_file_id)
780  allocate(var_ids(nvars))
781  call get_block_ids(diag_yaml_id, "varlist", var_ids, parent_block_id=diag_file_id)
782 
783  !< Loop through all the variables in the diag_file block and only count those that don't have write_var=false
784  total_nvars = 0
785  do i = 1, nvars
786  var_write = .true.
787  call get_value_from_key(diag_yaml_id, var_ids(i), "write_var", var_write, is_optional=.true.)
788  if (var_write) total_nvars = total_nvars + 1
789  end do
790 end function
791 
792 !> @brief This parses the freq, new_file_freq, or file_duration keys which are read in as a comma list
793 subroutine parse_key(filename, buffer, file_freq, file_frequnit, var)
794  character(len=*), intent(in) :: filename !< The name of the file (for error messages)
795  character(len=*), intent(inout) :: buffer !< Buffer that was read in from the yaml
796  integer, intent(out) :: file_freq(:) !< buffer to store the freq, new_file_freq, or
797  !! file_duration after it is parsed
798  integer, intent(out) :: file_frequnit(:) !< buffer to store the freq units, new_file_freq units,
799  !! or file_duration units after it is parsed
800  character(len=*), intent(in) :: var !< Name of the key parsing
801 
802  integer :: j !< location of the ",' in the buffer
803  integer :: k !< location of the " " that seperated the units
804  logical :: finished !< .true. if the parsing is complete
805  integer :: count !< Number of keys that have been parsed
806  character(len=255) :: str !< Member of the comma seperated list
807  character(len=10) :: units !< String to hold the units
808  integer :: err_unit !< Error key
809 
810  if (buffer .eq. "") return
811 
812  finished = .false.
813  j = 0
814  count = 0
815  do while (.not. finished)
816  count = count + 1
817  buffer = buffer(j+1:len_trim(buffer))
818  j = index(buffer, ",")
819  if (j == 0) then
820  !< There is only 1 member in the list
821  j = len_trim(buffer)+1
822  finished = .true.
823  endif
824 
825  str = adjustl(buffer(1:j-1))
826 
827  k = index(str, " ")
828  read(str(1:k-1), *, iostat=err_unit) file_freq(count)
829  units = str(k+1:len_trim(str))
830 
831  if (err_unit .ne. 0) &
832  call mpp_error(fatal, "Error parsing "//trim(var)//". Check your entry for file"//&
833  trim(filename))
834 
835  if (file_freq(count) .lt. -1) &
836  call mpp_error(fatal, trim(var)//" is not valid. &
837  &Check your entry for file:"//trim(filename))
838 
839  if (file_freq(count) .eq. -1 .or. file_freq(count) .eq. 0) then
840  !! The file is static so no need to read the units
841  file_frequnit(count) = diag_days
842  else
843  if (trim(units) .eq. "") &
844  call mpp_error(fatal, trim(var)//" units is required. &
845  &Check your entry for file:"//trim(filename))
846 
847  file_frequnit(count) = set_valid_time_units(units, &
848  trim(var)//" for file:"//trim(filename))
849  endif
850  enddo
851 end subroutine parse_key
852 
853 !> @brief This checks if the time unit in a diag file is valid and sets the integer equivalent
854 subroutine set_file_time_units (yaml_fileobj, file_timeunit)
855  type(diagyamlfiles_type), intent(inout) :: yaml_fileobj !< diagYamlFiles_type obj to checK
856  character(len=*), intent(in) :: file_timeunit !< file_timeunit as it is read from the diag_table
857 
858  yaml_fileobj%file_timeunit = set_valid_time_units(file_timeunit, "timeunit for file:"//trim(yaml_fileobj%file_fname))
859 end subroutine set_file_time_units
860 
861 !> @brief This checks if the filename_time in a diag file is correct and sets the integer equivalent
862 subroutine set_filename_time(yaml_fileobj, filename_time)
863  type(diagyamlfiles_type), intent(inout) :: yaml_fileobj !< diagYamlFiles_type obj to check
864  character(len=*), intent(in) :: filename_time !< filename_time as it is read from the yaml
865 
866  select case (trim(filename_time))
867  case ("")
868  yaml_fileobj%filename_time = middle_time !< This is the default
869  case ("begin")
870  yaml_fileobj%filename_time = begin_time
871  case ("middle")
872  yaml_fileobj%filename_time = middle_time
873  case ("end")
874  yaml_fileobj%filename_time = end_time
875  case default
876  call mpp_error(fatal, trim(filename_time)//" is an invalid filename_time &
877  &The acceptable values are begin, middle, and end. &
878  &Check your entry for file "//trim(yaml_fileobj%file_fname))
879  end select
880 end subroutine set_filename_time
881 
882 !> @brief This checks if the kind of a diag field is valid and sets it
883 subroutine set_field_kind(field, skind)
884  type(diagyamlfilesvar_type), intent(inout) :: field !< diagYamlFilesVar_type obj to read the contents into
885  character(len=*), intent(in) :: skind !< The variable kind as read from diag_yaml
886 
887  select case (trim(skind))
888  case ("r4")
889  field%var_kind = r4
890  case ("r8")
891  field%var_kind = r8
892  case ("i4")
893  field%var_kind = i4
894  case ("i8")
895  field%var_kind = i8
896  case default
897  call mpp_error(fatal, trim(skind)//" is an invalid kind! &
898  &The acceptable values are r4, r8, i4, i8. &
899  &Check your entry for file:"//trim(field%var_varname)//" in file "//trim(field%var_fname))
900  end select
901 
902 end subroutine set_field_kind
903 
904 !> @brief This checks if the reduction of a diag field is valid and sets it
905 !! If the reduction method is diurnalXX or powXX, it gets the number of diurnal sample and the power value
906 subroutine set_field_reduction(field, reduction_method)
907  type(diagyamlfilesvar_type), intent(inout) :: field !< diagYamlFilesVar_type obj to read the contents into
908  character(len=*) , intent(in) :: reduction_method!< reduction method as read from the yaml
909 
910  integer :: n_diurnal !< number of diurnal samples
911  integer :: pow_value !< The power value
912  integer :: ioerror !< io error status after reading in the diurnal samples
913 
914  n_diurnal = 0
915  pow_value = 0
916  ioerror = 0
917  if (index(reduction_method, "diurnal") .ne. 0) then
918  READ (reduction_method(8:len_trim(reduction_method)), fmt=*, iostat=ioerror) n_diurnal
919  if (ioerror .ne. 0) &
920  call mpp_error(fatal, "Error getting the number of diurnal samples from "//trim(reduction_method))
921  if (n_diurnal .le. 0) &
922  call mpp_error(fatal, "Diurnal samples should be greater than 0. &
923  & Check your entry for file:"//trim(field%var_varname)//" in file "//trim(field%var_fname))
924  field%var_reduction = time_diurnal
925  elseif (index(reduction_method, "pow") .ne. 0) then
926  READ (reduction_method(4:len_trim(reduction_method)), fmt=*, iostat=ioerror) pow_value
927  if (ioerror .ne. 0) &
928  call mpp_error(fatal, "Error getting the power value from "//trim(reduction_method))
929  if (pow_value .le. 0) &
930  call mpp_error(fatal, "The power value should be greater than 0. &
931  & Check your entry for file:"//trim(field%var_varname)//" in file "//trim(field%var_fname))
932  field%var_reduction = time_power
933  else
934  select case (reduction_method)
935  case ("none")
936  field%var_reduction = time_none
937  case ("average")
938  field%var_reduction = time_average
939  case ("min")
940  field%var_reduction = time_min
941  case ("max")
942  field%var_reduction = time_max
943  case ("rms")
944  field%var_reduction = time_rms
945  case ("sum")
946  field%var_reduction = time_sum
947  case default
948  call mpp_error(fatal, trim(reduction_method)//" is an invalid reduction method! &
949  &The acceptable values are none, average, pow##, diurnal##, min, max, and rms. &
950  &Check your entry for file:"//trim(field%var_varname)//" in file "//trim(field%var_fname))
951  end select
952  endif
953 
954  field%n_diurnal = n_diurnal
955  field%pow_value = pow_value
956 end subroutine set_field_reduction
957 
958 !> @brief This checks if a time unit is valid and if it is, it assigns the integer equivalent
959 !! @return The integer equivalent to the time units
960 function set_valid_time_units(time_units, error_msg) &
961 result(time_units_int)
962 
963  character(len=*), intent(in) :: time_units !< The time_units as a string
964  character(len=*), intent(in) :: error_msg !< Error message to append
965 
966  integer :: time_units_int !< The integer equivalent of the time_units
967 
968  select case (trim(time_units))
969  case ("seconds")
970  time_units_int = diag_seconds
971  case ("minutes")
972  time_units_int = diag_minutes
973  case ("hours")
974  time_units_int = diag_hours
975  case ("days")
976  time_units_int = diag_days
977  case ("months")
978  time_units_int = diag_months
979  case ("years")
980  time_units_int = diag_years
981  case default
982  time_units_int =diag_null
983  call mpp_error(fatal, trim(error_msg)//" is not valid. Acceptable values are &
984  &seconds, minutes, hours, days, months, years")
985  end select
986 end function set_valid_time_units
987 
988 !!!!!!! YAML FILE INQUIRIES !!!!!!!
989 !> @brief Finds the number of variables in the file_varlist
990 !! @return the size of the diag_files_obj%file_varlist array
991 integer pure function size_file_varlist (this)
992  class(diagyamlfiles_type), intent(in) :: this !< The object being inquiried
993  size_file_varlist = size(this%file_varlist)
994 end function size_file_varlist
995 
996 !> @brief Inquiry for diag_files_obj%file_fname
997 !! @return file_fname of a diag_yaml_file obj
998 pure function get_file_fname (this) &
999 result(res)
1000  class(diagyamlfiles_type), intent(in) :: this !< The object being inquiried
1001  character (len=:), allocatable :: res !< What is returned
1002  res = this%file_fname
1003 end function get_file_fname
1004 !> @brief Inquiry for diag_files_obj%file_frequnit
1005 !! @return file_frequnit of a diag_yaml_file_obj
1006 pure function get_file_frequnit (this) &
1007 result(res)
1008  class(diagyamlfiles_type), intent(in) :: this !< The object being inquiried
1009  integer :: res !< What is returned
1010  res = this%file_frequnit(this%current_new_file_freq_index)
1011 end function get_file_frequnit
1012 !> @brief Inquiry for diag_files_obj%file_freq
1013 !! @return file_freq of a diag_yaml_file_obj
1014 pure function get_file_freq(this) &
1015 result(res)
1016  class(diagyamlfiles_type), intent(in) :: this !< The object being inquiried
1017  integer :: res !< What is returned
1018  res = this%file_freq(this%current_new_file_freq_index)
1019 end function get_file_freq
1020 !> @brief Inquiry for diag_files_obj%file_timeunit
1021 !! @return file_timeunit of a diag_yaml_file_obj
1022 pure function get_file_timeunit (this) &
1023 result(res)
1024  class(diagyamlfiles_type), intent(in) :: this !< The object being inquiried
1025  integer :: res !< What is returned
1026  res = this%file_timeunit
1027 end function get_file_timeunit
1028 !> @brief Inquiry for diag_files_obj%file_unlimdim
1029 !! @return file_unlimdim of a diag_yaml_file_obj
1030 pure function get_file_unlimdim(this) &
1031 result(res)
1032  class(diagyamlfiles_type), intent(in) :: this !< The object being inquiried
1033  character (len=:), allocatable :: res !< What is returned
1034  res = this%file_unlimdim
1035 end function get_file_unlimdim
1036 !> @brief Inquiry for diag_files_obj%file_subregion
1037 !! @return file_sub_region of a diag_yaml_file_obj
1038 function get_file_sub_region (this) &
1039 result(res)
1040  class(diagyamlfiles_type), target, intent(in) :: this !< The object being inquiried
1041  type(subregion_type), pointer :: res !< What is returned
1042  res => this%file_sub_region
1043 end function get_file_sub_region
1044 !> @brief Inquiry for diag_files_obj%file_new_file_freq
1045 !! @return file_new_file_freq of a diag_yaml_file_obj
1046 pure function get_file_new_file_freq(this) &
1047 result(res)
1048  class(diagyamlfiles_type), intent(in) :: this !< The object being inquiried
1049  integer :: res !< What is returned
1050  res = this%file_new_file_freq(this%current_new_file_freq_index)
1051 end function get_file_new_file_freq
1052 !> @brief Inquiry for diag_files_obj%file_new_file_freq_units
1053 !! @return file_new_file_freq_units of a diag_yaml_file_obj
1054 pure function get_file_new_file_freq_units (this) &
1055 result(res)
1056  class(diagyamlfiles_type), intent(in) :: this !< The object being inquiried
1057  integer :: res !< What is returned
1058  res = this%file_new_file_freq_units(this%current_new_file_freq_index)
1059 end function get_file_new_file_freq_units
1060 !> @brief Inquiry for diag_files_obj%file_start_time
1061 !! @return file_start_time of a diag_yaml_file_obj
1062 pure function get_file_start_time (this) &
1063 result(res)
1064  class(diagyamlfiles_type), intent(in) :: this !< The object being inquiried
1065  character (len=:), allocatable :: res !< What is returned
1066  res = this%file_start_time
1067 end function get_file_start_time
1068 !> @brief Inquiry for diag_files_obj%file_duration
1069 !! @return file_duration of a diag_yaml_file_obj
1070 pure function get_file_duration (this) &
1071 result(res)
1072  class(diagyamlfiles_type), intent(in) :: this !< The object being inquiried
1073  integer :: res !< What is returned
1074  res = this%file_duration(this%current_new_file_freq_index)
1075 end function get_file_duration
1076 !> @brief Inquiry for diag_files_obj%file_duration_units
1077 !! @return file_duration_units of a diag_yaml_file_obj
1078 pure function get_file_duration_units (this) &
1079 result(res)
1080  class(diagyamlfiles_type), intent(in) :: this !< The object being inquiried
1081  integer :: res !< What is returned
1082  res = this%file_duration_units(this%current_new_file_freq_index)
1083 end function get_file_duration_units
1084 !> @brief Inquiry for diag_files_obj%file_varlist
1085 !! @return file_varlist of a diag_yaml_file_obj
1086 pure function get_file_varlist (this) &
1087 result(res)
1088  class(diagyamlfiles_type), intent(in) :: this !< The object being inquiried
1089  character (:), allocatable :: res(:) !< What is returned
1090  res = this%file_varlist
1091 end function get_file_varlist
1092 !> @brief Inquiry for diag_files_obj%file_global_meta
1093 !! @return file_global_meta of a diag_yaml_file_obj
1094 pure function get_file_global_meta (this) &
1095 result(res)
1096  class(diagyamlfiles_type), intent(in) :: this !< The object being inquiried
1097  character (len=MAX_STR_LEN), allocatable :: res(:,:) !< What is returned
1098  res = this%file_global_meta
1099 end function get_file_global_meta
1100 !> @brief Get the integer equivalent of the time to use to determine the filename,
1101 !! if using a wildcard file name (i.e ocn%4yr%2mo%2dy%2hr)
1102 !! @return the integer equivalent of the time to use to determine the filename
1103 pure function get_filename_time(this) &
1104 result(res)
1105  class(diagyamlfiles_type), intent(in) :: this !< The object being inquiried
1106  integer :: res !< What is returned
1107  res = this%filename_time
1108 end function
1109 !> @brief Inquiry for whether file_global_meta is allocated
1110 !! @return Flag indicating if file_global_meta is allocated
1111 function is_global_meta(this) &
1112  result(res)
1113  class(diagyamlfiles_type), intent(in) :: this !< The object being inquiried
1114  logical :: res
1115  res = .false.
1116  if (allocated(this%file_global_meta)) &
1117  res = .true.
1118 end function
1119 
1120 !> @brief Increate the current_new_file_freq_index by 1
1122  class(diagyamlfiles_type), intent(inout) :: this !< The file object
1123  this%current_new_file_freq_index = this%current_new_file_freq_index + 1
1124 end subroutine
1125 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1126 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1127 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1128 !!!!!!! VARIABLES ROUTINES AND FUNCTIONS !!!!!!!
1129 
1130 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1131 !!!!!!! YAML VAR INQUIRIES !!!!!!!
1132 !> @brief Inquiry for diag_yaml_files_var_obj%var_fname
1133 !! @return var_fname of a diag_yaml_files_var_obj
1134 pure function get_var_fname (this) &
1135 result(res)
1136  class(diagyamlfilesvar_type), intent(in) :: this !< The object being inquiried
1137  character (len=:), allocatable :: res !< What is returned
1138  res = this%var_fname
1139 end function get_var_fname
1140 !> @brief Inquiry for diag_yaml_files_var_obj%var_varname
1141 !! @return var_varname of a diag_yaml_files_var_obj
1142 pure function get_var_varname (this) &
1143 result(res)
1144  class(diagyamlfilesvar_type), intent(in) :: this !< The object being inquiried
1145  character (len=:), allocatable :: res !< What is returned
1146  res = this%var_varname
1147 end function get_var_varname
1148 !> @brief Inquiry for diag_yaml_files_var_obj%var_reduction
1149 !! @return var_reduction of a diag_yaml_files_var_obj
1150 pure function get_var_reduction (this) &
1151 result(res)
1152  class(diagyamlfilesvar_type), intent(in) :: this !< The object being inquiried
1153  integer, allocatable :: res !< What is returned
1154  res = this%var_reduction
1155 end function get_var_reduction
1156 !> @brief Inquiry for diag_yaml_files_var_obj%var_module
1157 !! @return var_module of a diag_yaml_files_var_obj
1158 pure function get_var_module (this) &
1159 result(res)
1160  class(diagyamlfilesvar_type), intent(in) :: this !< The object being inquiried
1161  character (len=:), allocatable :: res !< What is returned
1162  res = this%var_module
1163 end function get_var_module
1164 !> @brief Inquiry for diag_yaml_files_var_obj%var_kind
1165 !! @return var_kind of a diag_yaml_files_var_obj
1166 pure function get_var_kind (this) &
1167 result(res)
1168  class(diagyamlfilesvar_type), intent(in) :: this !< The object being inquiried
1169  integer, allocatable :: res !< What is returned
1170  res = this%var_kind
1171 end function get_var_kind
1172 !> @brief Inquiry for diag_yaml_files_var_obj%var_outname
1173 !! @return var_outname of a diag_yaml_files_var_obj
1174 pure function get_var_outname (this) &
1175 result(res)
1176  class(diagyamlfilesvar_type), intent(in) :: this !< The object being inquiried
1177  character (len=:), allocatable :: res !< What is returned
1178 
1179  if (this%has_var_outname()) then
1180  res = this%var_outname
1181  else
1182  res = this%var_varname !< If outname is not set, the variable name will be used
1183  endif
1184 end function get_var_outname
1185 !> @brief Inquiry for diag_yaml_files_var_obj%var_longname
1186 !! @return var_longname of a diag_yaml_files_var_obj
1187 pure function get_var_longname (this) &
1188 result(res)
1189  class(diagyamlfilesvar_type), intent(in) :: this !< The object being inquiried
1190  character (len=:), allocatable :: res !< What is returned
1191  res = this%var_longname
1192 end function get_var_longname
1193 !> @brief Inquiry for diag_yaml_files_var_obj%var_units
1194 !! @return var_units of a diag_yaml_files_var_obj
1195 pure function get_var_units (this) &
1196 result(res)
1197  class(diagyamlfilesvar_type), intent(in) :: this !< The object being inquiried
1198  character (len=:), allocatable :: res !< What is returned
1199  res = this%var_units
1200 end function get_var_units
1201 !> @brief Inquiry for diag_yaml_files_var_obj%var_zbounds
1202 !! @return var_zbounds of a diag_yaml_files_var_obj
1203 pure function get_var_zbounds (this) &
1204 result(res)
1205  class(diagyamlfilesvar_type), intent(in) :: this !< The object being inquiried
1206  real(kind=r4_kind) :: res(2) !< What is returned
1207  res = this%var_zbounds
1208 end function get_var_zbounds
1209 !> @brief Inquiry for diag_yaml_files_var_obj%var_attributes
1210 !! @return var_attributes of a diag_yaml_files_var_obj
1211 pure function get_var_attributes(this) &
1212 result(res)
1213  class(diagyamlfilesvar_type), intent(in) :: this !< The object being inquiried
1214  character (len=MAX_STR_LEN), allocatable :: res (:,:) !< What is returned
1215  res = this%var_attributes
1216 end function get_var_attributes
1217 !> @brief Inquiry for diag_yaml_files_var_obj%n_diurnal
1218 !! @return the number of diurnal samples of a diag_yaml_files_var_obj
1219 pure function get_n_diurnal(this) &
1220 result(res)
1221  class(diagyamlfilesvar_type), intent(in) :: this !< The object being inquiried
1222  integer :: res !< What is returned
1223  res = this%n_diurnal
1224 end function get_n_diurnal
1225 !> @brief Inquiry for diag_yaml_files_var_obj%pow_value
1226 !! @return the pow_value of a diag_yaml_files_var_obj
1227 pure function get_pow_value(this) &
1228 result(res)
1229  class(diagyamlfilesvar_type), intent(in) :: this !< The object being inquiried
1230  integer :: res !< What is returned
1231  res = this%pow_value
1232 end function get_pow_value
1233 !> @brief Inquiry for whether var_attributes is allocated
1234 !! @return Flag indicating if var_attributes is allocated
1235 function is_var_attributes(this) &
1236 result(res)
1237  class(diagyamlfilesvar_type), intent(in) :: this !< The object being inquiried
1238  logical :: res
1239  res = .false.
1240  if (allocated(this%var_attributes)) &
1241  res = .true.
1242 end function is_var_attributes
1243 
1244 !> @brief Initializes the non string values of a diagYamlFiles_type to its
1245 !! default values
1247  type(diagyamlfiles_type), intent(out) :: obj !< diagYamlFiles_type object to initialize
1248 
1249  obj%file_freq = diag_null
1250  obj%file_sub_region%tile = diag_null
1251  obj%file_new_file_freq = diag_null
1252  obj%file_duration = diag_null
1253  obj%file_new_file_freq_units = diag_null
1254  obj%file_duration_units = diag_null
1255  obj%current_new_file_freq_index = 1
1256 end subroutine diag_yaml_files_obj_init
1257 
1258 !> @brief Checks if diag_file_obj%file_fname is allocated
1259 !! @return true if diag_file_obj%file_fname is allocated
1260 pure logical function has_file_fname (this)
1261  class(diagyamlfiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1262  has_file_fname = allocated(this%file_fname)
1263 end function has_file_fname
1264 !> @brief Checks if diag_file_obj%file_frequnit is allocated
1265 !! @return true if diag_file_obj%file_frequnit is allocated
1266 pure logical function has_file_frequnit (this)
1267  class(diagyamlfiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1268  has_file_frequnit = this%file_frequnit(this%current_new_file_freq_index) .NE. diag_null
1269 end function has_file_frequnit
1270 !> @brief diag_file_obj%file_freq is on the stack, so the object always has it
1271 !! @return true if diag_file_obj%file_freq is allocated
1272 pure logical function has_file_freq (this)
1273  class(diagyamlfiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1274  has_file_freq = .true.
1275 end function has_file_freq
1276 !> @brief Checks if diag_file_obj%file_timeunit is allocated
1277 !! @return true if diag_file_obj%file_timeunit is allocated
1278 pure logical function has_file_timeunit (this)
1279  class(diagyamlfiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1280  has_file_timeunit = this%file_timeunit .ne. diag_null
1281 end function has_file_timeunit
1282 !> @brief Checks if diag_file_obj%file_unlimdim is allocated
1283 !! @return true if diag_file_obj%file_unlimdim is allocated
1284 pure logical function has_file_unlimdim (this)
1285  class(diagyamlfiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1286  has_file_unlimdim = allocated(this%file_unlimdim)
1287 end function has_file_unlimdim
1288 !> @brief Checks if diag_file_obj%file_write is on the stack, so this will always be true
1289 !! @return true
1290 pure logical function has_file_write (this)
1291  class(diagyamlfiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1292  has_file_write = .true.
1293 end function has_file_write
1294 !> @brief Checks if diag_file_obj%file_sub_region is being used and has the sub region variables allocated
1295 !! @return true if diag_file_obj%file_sub_region sub region variables are allocated
1296 pure logical function has_file_sub_region (this)
1297  class(diagyamlfiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1298  if ( this%file_sub_region%grid_type .eq. latlon_gridtype .or. this%file_sub_region%grid_type .eq. index_gridtype) then
1299  has_file_sub_region = .true.
1300  else
1301  has_file_sub_region = .false.
1302  endif
1303 end function has_file_sub_region
1304 !> @brief diag_file_obj%file_new_file_freq is defined on the stack, so this will return true
1305 !! @return true
1306 pure logical function has_file_new_file_freq (this)
1307  class(diagyamlfiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1308  has_file_new_file_freq = this%file_new_file_freq(this%current_new_file_freq_index) .ne. diag_null
1309 end function has_file_new_file_freq
1310 !> @brief Checks if diag_file_obj%file_new_file_freq_units is allocated
1311 !! @return true if diag_file_obj%file_new_file_freq_units is allocated
1312 pure logical function has_file_new_file_freq_units (this)
1313  class(diagyamlfiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1314  has_file_new_file_freq_units = this%file_new_file_freq_units(this%current_new_file_freq_index) .ne. diag_null
1315 end function has_file_new_file_freq_units
1316 !> @brief Checks if diag_file_obj%file_start_time is allocated
1317 !! @return true if diag_file_obj%file_start_time is allocated
1318 pure logical function has_file_start_time (this)
1319  class(diagyamlfiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1320  has_file_start_time = allocated(this%file_start_time)
1321 end function has_file_start_time
1322 !> @brief diag_file_obj%file_duration is allocated on th stack, so this is always true
1323 !! @return true
1324 pure logical function has_file_duration (this)
1325  class(diagyamlfiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1326  has_file_duration = this%file_duration(this%current_new_file_freq_index) .ne. diag_null
1327 end function has_file_duration
1328 !> @brief diag_file_obj%file_duration_units is on the stack, so this will retrun true
1329 !! @return true
1330 pure logical function has_file_duration_units (this)
1331  class(diagyamlfiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1332  has_file_duration_units = this%file_duration_units(this%current_new_file_freq_index) .ne. diag_null
1333 end function has_file_duration_units
1334 !> @brief Checks if diag_file_obj%file_varlist is allocated
1335 !! @return true if diag_file_obj%file_varlist is allocated
1336 pure logical function has_file_varlist (this)
1337  class(diagyamlfiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1338  has_file_varlist = allocated(this%file_varlist)
1339 end function has_file_varlist
1340 !> @brief Checks if diag_file_obj%file_global_meta is allocated
1341 !! @return true if diag_file_obj%file_global_meta is allocated
1342 pure logical function has_file_global_meta (this)
1343  class(diagyamlfiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1344  has_file_global_meta = allocated(this%file_global_meta)
1345 end function has_file_global_meta
1346 
1347 !> @brief Checks if diag_file_obj%var_fname is allocated
1348 !! @return true if diag_file_obj%var_fname is allocated
1349 pure logical function has_var_fname (this)
1350  class(diagyamlfilesvar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1351  has_var_fname = allocated(this%var_fname)
1352 end function has_var_fname
1353 !> @brief Checks if diag_file_obj%var_varname is allocated
1354 !! @return true if diag_file_obj%var_varname is allocated
1355 pure logical function has_var_varname (this)
1356  class(diagyamlfilesvar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1357  has_var_varname = allocated(this%var_varname)
1358 end function has_var_varname
1359 !> @brief Checks if diag_file_obj%var_reduction is allocated
1360 !! @return true if diag_file_obj%var_reduction is allocated
1361 pure logical function has_var_reduction (this)
1362  class(diagyamlfilesvar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1363  has_var_reduction = allocated(this%var_reduction)
1364 end function has_var_reduction
1365 !> @brief Checks if diag_file_obj%var_module is allocated
1366 !! @return true if diag_file_obj%var_module is allocated
1367 pure logical function has_var_module (this)
1368  class(diagyamlfilesvar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1369  has_var_module = allocated(this%var_module)
1370 end function has_var_module
1371 !> @brief Checks if diag_file_obj%var_kind is allocated
1372 !! @return true if diag_file_obj%var_kind is allocated
1373 pure logical function has_var_kind (this)
1374  class(diagyamlfilesvar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1375  has_var_kind = allocated(this%var_kind)
1376 end function has_var_kind
1377 !> @brief diag_file_obj%var_write is on the stack, so this returns true
1378 !! @return true
1379 pure logical function has_var_write (this)
1380  class(diagyamlfilesvar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1381  has_var_write = .true.
1382 end function has_var_write
1383 !> @brief Checks if diag_file_obj%var_outname is allocated
1384 !! @return true if diag_file_obj%var_outname is allocated
1385 pure logical function has_var_outname (this)
1386  class(diagyamlfilesvar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1387  if (allocated(this%var_outname)) then
1388  if (trim(this%var_outname) .ne. "") then
1389  has_var_outname = .true.
1390  else
1391  has_var_outname = .false.
1392  endif
1393  else
1394  has_var_outname = .true.
1395  endif
1396 end function has_var_outname
1397 !> @brief Checks if diag_file_obj%var_longname is allocated
1398 !! @return true if diag_file_obj%var_longname is allocated
1399 pure logical function has_var_longname (this)
1400  class(diagyamlfilesvar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1401  has_var_longname = allocated(this%var_longname)
1402 end function has_var_longname
1403 !> @brief Checks if diag_file_obj%var_units is allocated
1404 !! @return true if diag_file_obj%var_units is allocated
1405 pure logical function has_var_units (this)
1406  class(diagyamlfilesvar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1407  has_var_units = allocated(this%var_units)
1408 end function has_var_units
1409 !> @brief Checks if diag_file_obj%var_zbounds is allocated
1410 !! @return true if diag_file_obj%var_zbounds is allocated
1411 pure logical function has_var_zbounds (this)
1412  class(diagyamlfilesvar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1413  has_var_zbounds = any(this%var_zbounds .ne. diag_null)
1414 end function has_var_zbounds
1415 !> @brief Checks if diag_file_obj%var_attributes is allocated
1416 !! @return true if diag_file_obj%var_attributes is allocated
1417 pure logical function has_var_attributes (this)
1418  class(diagyamlfilesvar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1419  has_var_attributes = allocated(this%var_attributes)
1420 end function has_var_attributes
1421 !> @brief Checks if diag_file_obj%n_diurnal is set
1422 !! @return true if diag_file_obj%n_diurnal is set
1423 pure logical function has_n_diurnal(this)
1424  class(diagyamlfilesvar_type), intent(in) :: this !< diagYamlvar_type object to inquire
1425  has_n_diurnal = (this%n_diurnal .ne. 0)
1426 end function has_n_diurnal
1427 !> @brief Checks if diag_file_obj%pow_value is set
1428 !! @return true if diag_file_obj%pow_value is set
1429 pure logical function has_pow_value(this)
1430  class(diagyamlfilesvar_type), intent(in) :: this !< diagYamlvar_type object to inquire
1431  has_pow_value = (this%pow_value .ne. 0)
1432 end function has_pow_value
1433 
1434 !> @brief Checks if diag_file_obj%diag_title is allocated
1435 !! @return true if diag_file_obj%diag_title is allocated
1436 pure logical function has_diag_title (this)
1437  class(diagyamlobject_type), intent(in) :: this !< diagYamlObject_type object to inquire
1438  has_diag_title = allocated(this%diag_title)
1439 end function has_diag_title
1440 !> @brief diag_file_obj%diag_basedate is on the stack, so this is always true
1441 !! @return true
1442 pure logical function has_diag_basedate (this)
1443  class(diagyamlobject_type), intent(in) :: this !< diagYamlObject_type object to initialize
1444  has_diag_basedate = .true.
1445 end function has_diag_basedate
1446 !> @brief Checks if diag_file_obj%diag_files is allocated
1447 !! @return true if diag_file_obj%diag_files is allocated
1448 pure logical function has_diag_files (this)
1449  class(diagyamlobject_type), intent(in) :: this !< diagYamlObject_type object to initialize
1450  has_diag_files = allocated(this%diag_files)
1451 end function has_diag_files
1452 !> @brief Checks if diag_file_obj%diag_fields is allocated
1453 !! @return true if diag_file_obj%diag_fields is allocated
1454 pure logical function has_diag_fields (this)
1455  class(diagyamlobject_type), intent(in) :: this !< diagYamlObject_type object to initialize
1456  has_diag_fields = allocated(this%diag_fields)
1457 end function has_diag_fields
1458 
1459 !> @brief Determine the number of unique diag_fields in the diag_yaml_object
1460 !! @return The number of unique diag_fields
1462  result(nfields)
1463  integer :: nfields
1464  nfields = fms_find_unique(variable_list%var_pointer, size(variable_list%var_pointer))
1465 
1466 end function get_num_unique_fields
1467 
1468 !> @brief Determines if a diag_field is in the diag_yaml_object
1469 !! @return Indices of the locations where the field was found
1470 function find_diag_field(diag_field_name, module_name) &
1471 result(indices)
1472 
1473  character(len=*), intent(in) :: diag_field_name !< diag_field name to search for
1474  character(len=*), intent(in) :: module_name !< Name of the module, the variable is in
1475 
1476  integer, allocatable :: indices(:)
1477 
1478  indices = fms_find_my_string(variable_list%var_pointer, size(variable_list%var_pointer), &
1479  & lowercase(trim(diag_field_name))//":"//lowercase(trim(module_name)//c_null_char))
1480 end function find_diag_field
1481 
1482 !> @brief Gets the diag_field entries corresponding to the indices of the sorted variable_list
1483 !! @return Array of diag_fields
1484 function get_diag_fields_entries(indices) &
1485  result(diag_field)
1486 
1487  integer, intent(in) :: indices(:) !< Indices of the field in the sorted variable_list array
1488  type(diagyamlfilesvar_type), dimension (:), allocatable :: diag_field
1489 
1490  integer :: i !< For do loops
1491  integer :: field_id !< Indices of the field in the diag_yaml array
1492 
1493  allocate(diag_field(size(indices)))
1494 
1495  do i = 1, size(indices)
1496  field_id = variable_list%diag_field_indices(indices(i))
1497  diag_field(i) = diag_yaml%diag_fields(field_id)
1498  end do
1499 
1500 end function get_diag_fields_entries
1501 
1502 !> @brief Gets field indices corresponding to the indices (input argument) in the sorted variable_list
1503 !! @return Copy of array of field indices
1504 function get_diag_field_ids(indices) result(field_ids)
1505 
1506  integer, intent(in) :: indices(:) !< Indices of the fields in the sorted variable_list array
1507  integer, allocatable :: field_ids(:)
1508  integer :: i !< For do loop
1509 
1510  allocate(field_ids(size(indices)))
1511 
1512  do i = 1, size(indices)
1513  field_ids(i) = variable_list%diag_field_indices(indices(i))
1514  end do
1515 
1516 end function get_diag_field_ids
1517 
1518 !> @brief Finds the indices of the diag_yaml%diag_files(:) corresponding to fields in variable_list(indices)
1519 !! @return indices of the diag_yaml%diag_files(:)
1520 function get_diag_files_id(indices) &
1521  result(file_id)
1522 
1523  integer, intent(in) :: indices(:) !< Indices of the field in the sorted variable_list
1524  integer, allocatable :: file_id(:)
1525 
1526  integer :: field_id !< Indices of the field in the diag_yaml field array
1527  integer :: i !< For do loops
1528  character(len=FMS_FILE_LEN) :: filename !< Filename of the field
1529  integer, allocatable :: file_indices(:) !< Indices of the file in the sorted variable_list
1530 
1531  allocate(file_id(size(indices)))
1532 
1533  do i = 1, size(indices)
1534  field_id = variable_list%diag_field_indices(indices(i))
1535  !< Get the filename of the field
1536  filename = diag_yaml%diag_fields(field_id)%var_fname
1537 
1538  !< File indice of that file in the array of list of sorted files
1539  file_indices = fms_find_my_string(file_list%file_pointer, size(file_list%file_pointer), &
1540  & trim(filename)//c_null_char)
1541 
1542  if (size(file_indices) .ne. 1) &
1543  & call mpp_error(fatal, "get_diag_files_id: Error getting the correct number of file indices!"//&
1544  " The diag file "//trim(filename)//" was defined "//string(size(file_indices))&
1545  // " times")
1546 
1547  if (file_indices(1) .eq. diag_null) &
1548  & call mpp_error(fatal, "get_diag_files_id: Error finding the file "//trim(filename)//" in the diag_files yaml")
1549 
1550  !< Get the index of the file in the diag_yaml file
1551  file_id(i) = file_list%diag_file_indices(file_indices(1))
1552  end do
1553 
1554 end function get_diag_files_id
1555 
1556 !> Prints out values from diag_yaml object for debugging.
1557 !! Only writes on root.
1558 subroutine dump_diag_yaml_obj( filename )
1559  character(len=*), optional, intent(in) :: filename !< optional name of logfile to write to, otherwise
1560  !! prints to stdout
1561  type(diagyamlfilesvar_type), allocatable :: fields(:)
1562  type(diagyamlfiles_type), pointer :: files(:)
1563  integer :: i, unit_num
1564  if( present(filename)) then
1565  open(newunit=unit_num, file=trim(filename), action='WRITE')
1566  else
1567  unit_num = stdout()
1568  endif
1569  !! TODO write to log
1570  if( mpp_pe() .eq. mpp_root_pe()) then
1571  write(unit_num, *) '**********Dumping diag_yaml object**********'
1572  if( diag_yaml%has_diag_title()) write(unit_num, *) 'Title:', diag_yaml%diag_title
1573  if( diag_yaml%has_diag_basedate()) write(unit_num, *) 'basedate array:', diag_yaml%diag_basedate
1574  write(unit_num, *) 'FILES'
1575  allocate(fields(SIZE(diag_yaml%get_diag_fields())))
1576  files => diag_yaml%diag_files
1577  fields = diag_yaml%get_diag_fields()
1578  do i=1, SIZE(files)
1579  write(unit_num, *) 'File: ', files(i)%get_file_fname()
1580  if(files(i)%has_file_frequnit()) write(unit_num, *) 'file_frequnit:', files(i)%get_file_frequnit()
1581  if(files(i)%has_file_freq()) write(unit_num, *) 'freq:', files(i)%get_file_freq()
1582  if(files(i)%has_file_timeunit()) write(unit_num, *) 'timeunit:', files(i)%get_file_timeunit()
1583  if(files(i)%has_file_unlimdim()) write(unit_num, *) 'unlimdim:', files(i)%get_file_unlimdim()
1584  !if(files(i)%has_file_sub_region()) write(unit_num, *) 'sub_region:', files(i)%get_file_sub_region()
1585  if(files(i)%has_file_new_file_freq()) write(unit_num, *) 'new_file_freq:', files(i)%get_file_new_file_freq()
1586  if(files(i)%has_file_new_file_freq_units()) write(unit_num, *) 'new_file_freq_units:', &
1587  & files(i)%get_file_new_file_freq_units()
1588  if(files(i)%has_file_start_time()) write(unit_num, *) 'start_time:', files(i)%get_file_start_time()
1589  if(files(i)%has_file_duration()) write(unit_num, *) 'duration:', files(i)%get_file_duration()
1590  if(files(i)%has_file_duration_units()) write(unit_num, *) 'duration_units:', files(i)%get_file_duration_units()
1591  if(files(i)%has_file_varlist()) write(unit_num, *) 'varlist:', files(i)%get_file_varlist()
1592  if(files(i)%has_file_global_meta()) write(unit_num, *) 'global_meta:', files(i)%get_file_global_meta()
1593  if(files(i)%is_global_meta()) write(unit_num, *) 'global_meta:', files(i)%is_global_meta()
1594  write(unit_num, *) ''
1595  enddo
1596  write(unit_num, *) 'FIELDS'
1597  do i=1, SIZE(fields)
1598  write(unit_num, *) 'Field: ', fields(i)%get_var_fname()
1599  if(fields(i)%has_var_fname()) write(unit_num, *) 'fname:', fields(i)%get_var_fname()
1600  if(fields(i)%has_var_varname()) write(unit_num, *) 'varname:', fields(i)%get_var_varname()
1601  if(fields(i)%has_var_reduction()) write(unit_num, *) 'reduction:', fields(i)%get_var_reduction()
1602  if(fields(i)%has_var_module()) write(unit_num, *) 'module:', fields(i)%get_var_module()
1603  if(fields(i)%has_var_kind()) write(unit_num, *) 'kind:', fields(i)%get_var_kind()
1604  if(fields(i)%has_var_outname()) write(unit_num, *) 'outname:', fields(i)%get_var_outname()
1605  if(fields(i)%has_var_longname()) write(unit_num, *) 'longname:', fields(i)%get_var_longname()
1606  if(fields(i)%has_var_units()) write(unit_num, *) 'units:', fields(i)%get_var_units()
1607  if(fields(i)%has_var_zbounds()) write(unit_num, *) 'zbounds:', fields(i)%get_var_zbounds()
1608  if(fields(i)%has_var_attributes()) write(unit_num, *) 'attributes:', fields(i)%get_var_attributes()
1609  if(fields(i)%has_n_diurnal()) write(unit_num, *) 'n_diurnal:', fields(i)%get_n_diurnal()
1610  if(fields(i)%has_pow_value()) write(unit_num, *) 'pow_value:', fields(i)%get_pow_value()
1611  if(fields(i)%has_var_attributes()) write(unit_num, *) 'is_var_attributes:', fields(i)%is_var_attributes()
1612  enddo
1613  deallocate(fields)
1614  nullify(files)
1615  if( present(filename)) then
1616  close(unit_num)
1617  endif
1618  endif
1619 end subroutine
1620 
1621 !> Writes an output yaml with all available information on the written files.
1622 !! Will only write with root pe.
1623 !! Global attributes are limited to 16 per file.
1624 subroutine fms_diag_yaml_out()
1625  type(diagyamlfiles_type), pointer :: fileptr !< pointer for individual variables
1626  type(diagyamlfilesvar_type), pointer :: varptr !< pointer for individual variables
1627  type (fmsyamloutkeys_type), allocatable :: keys(:), keys2(:), keys3(:)
1628  type (fmsyamloutvalues_type), allocatable :: vals(:), vals2(:), vals3(:)
1629  integer :: i, j, k
1630  character(len=128) :: tmpstr1, tmpstr2 !< string to store output fields
1631  integer, parameter :: tier1size = 3 !< size of first tier, will always be 3 for basedate, title and diag_files
1632  integer :: tier2size, tier3size !< size of each 'tier'(based one numbers of tabs) in the yaml
1633  integer, allocatable :: tier3each(:) !< tier 3 list sizes corresponding to where they are in the second tier
1634  integer, dimension(basedate_size) :: basedate_loc !< local copy of basedate to loop through
1635  integer :: varnum_i, key3_i, gm
1636  character(len=32), allocatable :: st_vals(:) !< start times for gcc bug
1637 
1638  if( mpp_pe() .ne. mpp_root_pe()) return
1639 
1640  allocate(tier3each(SIZE(diag_yaml%diag_files) * 3))
1641  tier3size = 0; tier3each = 0
1642 
1643  !! allocations for key+val structs
1644  allocate(keys(1))
1645  allocate(vals(1))
1646  allocate(keys2(SIZE(diag_yaml%diag_files)))
1647  allocate(vals2(SIZE(diag_yaml%diag_files)))
1648  allocate(st_vals(SIZE(diag_yaml%diag_files)))
1649  do i=1, SIZE(diag_yaml%diag_files)
1650  call initialize_key_struct(keys2(i))
1651  call initialize_val_struct(vals2(i))
1652  if (allocated(diag_yaml%diag_files(i)%file_varlist) ) then
1653  do j=1, SIZE(diag_yaml%diag_files(i)%file_varlist)
1654  tier3size = tier3size + 1
1655  enddo
1656  endif
1657  tier3size = tier3size + 2
1658  enddo
1659  allocate(keys3(tier3size))
1660  allocate(vals3(tier3size))
1661 
1662  !! tier 1 - title, basedate, diag_files
1663  call initialize_key_struct(keys(1))
1664  call initialize_val_struct(vals(1))
1665  call fms_f2c_string( keys(1)%key1, 'title')
1666  call fms_f2c_string( vals(1)%val1, diag_yaml%diag_title)
1667  call fms_f2c_string( keys(1)%key2, 'base_date')
1668  basedate_loc = diag_yaml%get_basedate()
1669  tmpstr1 = ''; tmpstr2 = ''
1670  tmpstr1 = string(basedate_loc(1))
1671  tmpstr2 = trim(tmpstr1)
1672  do i=2, basedate_size
1673  tmpstr1 = string(basedate_loc(i))
1674  tmpstr2 = trim(tmpstr2) // ' ' // trim(tmpstr1)
1675  enddo
1676  call fms_f2c_string(vals(1)%val2, trim(tmpstr2))
1677  call yaml_out_add_level2key('diag_files', keys(1))
1678  key3_i = 0
1679  !! tier 2 - diag files
1680  do i=1, SIZE(diag_yaml%diag_files)
1681  fileptr => diag_yaml%diag_files(i)
1682 
1683  call fms_f2c_string(keys2(i)%key1, 'file_name')
1684  call fms_f2c_string(keys2(i)%key2, 'freq')
1685  call fms_f2c_string(keys2(i)%key3, 'freq_units')
1686  call fms_f2c_string(keys2(i)%key4, 'time_units')
1687  call fms_f2c_string(keys2(i)%key5, 'unlimdim')
1688  call fms_f2c_string(keys2(i)%key6, 'new_file_freq')
1689  call fms_f2c_string(keys2(i)%key7, 'new_file_freq_units')
1690  call fms_f2c_string(keys2(i)%key8, 'start_time')
1691  call fms_f2c_string(keys2(i)%key9, 'file_duration')
1692  call fms_f2c_string(keys2(i)%key10, 'file_duration_units')
1693 
1694  call fms_f2c_string(vals2(i)%val1, fileptr%file_fname)
1695  call fms_f2c_string(vals2(i)%val5, fileptr%file_unlimdim)
1696  call fms_f2c_string(vals2(i)%val4, get_diag_unit_string((/fileptr%file_timeunit/)))
1697  tmpstr1 = ''
1698  do k=1, SIZE(fileptr%file_freq)
1699  if(fileptr%file_freq(k) .eq. diag_null) exit
1700  tmpstr2 = ''
1701  tmpstr2 = string(fileptr%file_freq(k))
1702  tmpstr1 = trim(tmpstr1)//" "//trim(tmpstr2)
1703  enddo
1704  call fms_f2c_string(vals2(i)%val2, adjustl(tmpstr1))
1705  call fms_f2c_string(vals2(i)%val3, get_diag_unit_string(fileptr%file_frequnit))
1706  tmpstr1 = ''
1707  do k=1, SIZE(fileptr%file_new_file_freq)
1708  if(fileptr%file_new_file_freq(k) .eq. diag_null) exit
1709  tmpstr2 = ''
1710  tmpstr2 = string(fileptr%file_new_file_freq(k))
1711  tmpstr1 = trim(tmpstr1)//" "//trim(tmpstr2)
1712  enddo
1713  call fms_f2c_string(vals2(i)%val6, adjustl(tmpstr1))
1714  call fms_f2c_string(vals2(i)%val7, get_diag_unit_string(fileptr%file_new_file_freq_units))
1715  call fms_f2c_string(vals2(i)%val8, trim(fileptr%get_file_start_time()))
1716  st_vals(i) = fileptr%get_file_start_time()
1717  tmpstr1 = ''
1718  do k=1, SIZE(fileptr%file_duration)
1719  if(fileptr%file_duration(k) .eq. diag_null) exit
1720  tmpstr2 = ''
1721  tmpstr2 = string(fileptr%file_duration(k))
1722  tmpstr1 = trim(tmpstr1)//" "//trim(tmpstr2)
1723  enddo
1724  call fms_f2c_string(vals2(i)%val9, adjustl(tmpstr1))
1725  call fms_f2c_string(vals2(i)%val10, get_diag_unit_string(fileptr%file_duration_units))
1726 
1727  !! tier 3 - varlists, subregion, global metadata
1728  call yaml_out_add_level2key('varlist', keys2(i))
1729  j = 0
1730  if( SIZE(fileptr%file_varlist) .gt. 0) then
1731  do j=1, SIZE(fileptr%file_varlist)
1732  key3_i = key3_i + 1
1733  call initialize_key_struct(keys3(key3_i))
1734  call initialize_val_struct(vals3(key3_i))
1735  !! find the variable object from the list
1736  varptr => null()
1737  do varnum_i=1, SIZE(diag_yaml%diag_fields)
1738  if( trim(diag_yaml%diag_fields(varnum_i)%var_varname ) .eq. trim(fileptr%file_varlist(j)) .and. &
1739  trim(diag_yaml%diag_fields(varnum_i)%var_fname) .eq. trim(fileptr%file_fname)) then
1740  ! if theres a output name, that should match as well
1741  if(diag_yaml%diag_fields(varnum_i)%has_var_outname()) then
1742  if(trim(diag_yaml%diag_fields(varnum_i)%var_outname) .eq. trim(fileptr%file_outlist(j))) then
1743  varptr => diag_yaml%diag_fields(varnum_i)
1744  exit
1745  endif
1746  else
1747  varptr => diag_yaml%diag_fields(varnum_i)
1748  exit
1749  endif
1750  endif
1751  enddo
1752  if( .not. associated(varptr)) call mpp_error(fatal, "diag_yaml_output: could not find variable in list."//&
1753  " var: "// trim(fileptr%file_varlist(j)))
1754  call fms_f2c_string(keys3(key3_i)%key1, 'module')
1755  call fms_f2c_string(keys3(key3_i)%key2, 'var_name')
1756  call fms_f2c_string(keys3(key3_i)%key3, 'reduction')
1757  call fms_f2c_string(keys3(key3_i)%key4, 'kind')
1758  call fms_f2c_string(keys3(key3_i)%key5, 'output_name')
1759  call fms_f2c_string(keys3(key3_i)%key6, 'long_name')
1760  call fms_f2c_string(keys3(key3_i)%key7, 'units')
1761  call fms_f2c_string(keys3(key3_i)%key8, 'zbounds')
1762  call fms_f2c_string(keys3(key3_i)%key9, 'n_diurnal')
1763  call fms_f2c_string(keys3(key3_i)%key10, 'pow_value')
1764  call fms_f2c_string(keys3(key3_i)%key11, 'dimensions')
1765  if (varptr%has_var_module()) call fms_f2c_string(vals3(key3_i)%val1, varptr%var_module)
1766  if (varptr%has_var_varname()) call fms_f2c_string(vals3(key3_i)%val2, varptr%var_varname)
1767  if (varptr%has_var_reduction()) then
1768  call fms_f2c_string(vals3(key3_i)%val3, &
1769  get_diag_reduction_string((/varptr%var_reduction/)))
1770  endif
1771  if (varptr%has_var_outname()) call fms_f2c_string(vals3(key3_i)%val5, varptr%var_outname)
1772  if (varptr%has_var_longname()) call fms_f2c_string(vals3(key3_i)%val6, varptr%var_longname)
1773  if (varptr%has_var_units()) call fms_f2c_string(vals3(key3_i)%val7, varptr%var_units)
1774  if (varptr%has_var_kind()) then
1775  select case(varptr%var_kind)
1776  case(i4)
1777  call fms_f2c_string(vals3(key3_i)%val4, 'i4')
1778  case(i8)
1779  call fms_f2c_string(vals3(key3_i)%val4, 'i8')
1780  case(r4)
1781  call fms_f2c_string(vals3(key3_i)%val4, 'r4')
1782  case(r8)
1783  call fms_f2c_string(vals3(key3_i)%val4, 'r8')
1784  end select
1785  endif
1786 
1787  if( abs(varptr%var_zbounds(1) - real(diag_null, r4_kind)) .gt. 1.0e-5 ) then
1788  tmpstr2 = string(varptr%var_zbounds(1), "F8.2") // ' ' // string(varptr%var_zbounds(2), "F8.2")
1789  call fms_f2c_string(vals3(key3_i)%val8, trim(tmpstr2))
1790  endif
1791 
1792  if( varptr%n_diurnal .gt. 0) then
1793  tmpstr1 = ''; tmpstr1 = string(varptr%n_diurnal)
1794  call fms_f2c_string(vals3(key3_i)%val9, tmpstr1)
1795  endif
1796 
1797  if( varptr%pow_value .gt. 0) then
1798  tmpstr1 = ''; tmpstr1 = string(varptr%pow_value)
1799  call fms_f2c_string(vals3(key3_i)%val10, tmpstr1)
1800  endif
1801 
1802  tmpstr1 = ''; tmpstr1 = varptr%var_axes_names
1803  call fms_f2c_string(vals3(key3_i)%val11, trim(adjustl(tmpstr1)))
1804  enddo
1805  endif
1806 
1807  key3_i = key3_i + 1
1808  tier3each(i*3-2) = j-1 ! j-1 structs to print for varlist keys
1809  tier3each(i*3-1) = 1 ! 1 struct per sub_region key
1810  tier3each(i*3) = 1 ! 1 struct per global metadata key
1811  call initialize_key_struct(keys3(key3_i))
1812  call initialize_val_struct(vals3(key3_i))
1813  !! sub region
1814  call yaml_out_add_level2key('sub_region', keys2(i))
1815  call fms_f2c_string(keys3(key3_i)%key1, 'grid_type')
1816  call fms_f2c_string(keys3(key3_i)%key2, 'tile')
1817  call fms_f2c_string(keys3(key3_i)%key3, 'corner1')
1818  call fms_f2c_string(keys3(key3_i)%key4, 'corner2')
1819  call fms_f2c_string(keys3(key3_i)%key5, 'corner3')
1820  call fms_f2c_string(keys3(key3_i)%key6, 'corner4')
1821 
1822  select case (fileptr%file_sub_region%grid_type)
1823  case(latlon_gridtype)
1824  call fms_f2c_string(vals3(key3_i)%val1, 'latlon')
1825  case(index_gridtype)
1826  call fms_f2c_string(vals3(key3_i)%val1, 'index')
1827  end select
1828  if(fileptr%file_sub_region%tile .ne. diag_null) then
1829  tmpstr1 = ''; tmpstr1 = string(fileptr%file_sub_region%tile)
1830  call fms_f2c_string(vals3(key3_i)%val2, tmpstr1)
1831  endif
1832  if(fileptr%has_file_sub_region()) then
1833  if( allocated(fileptr%file_sub_region%corners)) then
1834  select type (corners => fileptr%file_sub_region%corners)
1835  type is (real(r8_kind))
1836  tmpstr1 = ''; tmpstr1 = string(corners(1,1))
1837  tmpstr2 = ''; tmpstr2 = string(corners(1,2))
1838  call fms_f2c_string(vals3(key3_i)%val3, trim(tmpstr1)//' '//trim(tmpstr2))
1839  tmpstr1 = ''; tmpstr1 = string(corners(2,1))
1840  tmpstr2 = ''; tmpstr2 = string(corners(2,2))
1841  call fms_f2c_string(vals3(key3_i)%val4, trim(tmpstr1)//' '//trim(tmpstr2))
1842  tmpstr1 = ''; tmpstr1 = string(corners(3,1))
1843  tmpstr2 = ''; tmpstr2 = string(corners(3,2))
1844  call fms_f2c_string(vals3(key3_i)%val5, trim(tmpstr1)//' '//trim(tmpstr2))
1845  tmpstr1 = ''; tmpstr1 = string(corners(4,1))
1846  tmpstr2 = ''; tmpstr2 = string(corners(4,2))
1847  call fms_f2c_string(vals3(key3_i)%val6, trim(tmpstr1)//' '//trim(tmpstr2))
1848  type is (real(r4_kind))
1849  tmpstr1 = ''; tmpstr1 = string(corners(1,1))
1850  tmpstr2 = ''; tmpstr2 = string(corners(1,2))
1851  call fms_f2c_string(vals3(key3_i)%val3, trim(tmpstr1)//' '//trim(tmpstr2))
1852  tmpstr1 = ''; tmpstr1 = string(corners(2,1))
1853  tmpstr2 = ''; tmpstr2 = string(corners(2,2))
1854  call fms_f2c_string(vals3(key3_i)%val4, trim(tmpstr1)//' '//trim(tmpstr2))
1855  tmpstr1 = ''; tmpstr1 = string(corners(3,1))
1856  tmpstr2 = ''; tmpstr2 = string(corners(3,2))
1857  call fms_f2c_string(vals3(key3_i)%val5, trim(tmpstr1)//' '//trim(tmpstr2))
1858  tmpstr1 = ''; tmpstr1 = string(corners(4,1))
1859  tmpstr2 = ''; tmpstr2 = string(corners(4,2))
1860  call fms_f2c_string(vals3(key3_i)%val6, trim(tmpstr1)//' '//trim(tmpstr2))
1861  type is (integer(i4_kind))
1862  tmpstr1 = ''; tmpstr1 = string(corners(1,1))
1863  tmpstr2 = ''; tmpstr2 = string(corners(1,2))
1864  call fms_f2c_string(vals3(key3_i)%val3, trim(tmpstr1)//' '//trim(tmpstr2))
1865  tmpstr1 = ''; tmpstr1 = string(corners(2,1))
1866  tmpstr2 = ''; tmpstr2 = string(corners(2,2))
1867  call fms_f2c_string(vals3(key3_i)%val4, trim(tmpstr1)//' '//trim(tmpstr2))
1868  tmpstr1 = ''; tmpstr1 = string(corners(3,1))
1869  tmpstr2 = ''; tmpstr2 = string(corners(3,2))
1870  call fms_f2c_string(vals3(key3_i)%val5, trim(tmpstr1)//' '//trim(tmpstr2))
1871  tmpstr1 = ''; tmpstr1 = string(corners(4,1))
1872  tmpstr2 = ''; tmpstr2 = string(corners(4,2))
1873  call fms_f2c_string(vals3(key3_i)%val6, trim(tmpstr1)//' '//trim(tmpstr2))
1874  type is (integer(i8_kind))
1875  tmpstr1 = ''; tmpstr1 = string(corners(1,1))
1876  tmpstr2 = ''; tmpstr2 = string(corners(1,2))
1877  call fms_f2c_string(vals3(key3_i)%val3, trim(tmpstr1)//' '//trim(tmpstr2))
1878  tmpstr1 = ''; tmpstr1 = string(corners(2,1))
1879  tmpstr2 = ''; tmpstr2 = string(corners(2,2))
1880  call fms_f2c_string(vals3(key3_i)%val4, trim(tmpstr1)//' '//trim(tmpstr2))
1881  tmpstr1 = ''; tmpstr1 = string(corners(3,1))
1882  tmpstr2 = ''; tmpstr2 = string(corners(3,2))
1883  call fms_f2c_string(vals3(key3_i)%val5, trim(tmpstr1)//' '//trim(tmpstr2))
1884  tmpstr1 = ''; tmpstr1 = string(corners(4,1))
1885  tmpstr2 = ''; tmpstr2 = string(corners(4,2))
1886  call fms_f2c_string(vals3(key3_i)%val6, trim(tmpstr1)//' '//trim(tmpstr2))
1887  end select
1888  endif
1889  endif
1890  !! global metadata
1891  key3_i = key3_i + 1
1892  call initialize_key_struct(keys3(key3_i))
1893  call initialize_val_struct(vals3(key3_i))
1894  call yaml_out_add_level2key('global_meta', keys2(i))
1895  if ( fileptr%has_file_global_meta()) then
1896  do gm=1, SIZE(fileptr%file_global_meta, 1)
1897  select case(gm)
1898  case (1)
1899  call fms_f2c_string(keys3(key3_i)%key1, fileptr%file_global_meta(1,1))
1900  call fms_f2c_string(vals3(key3_i)%val1, fileptr%file_global_meta(1,2))
1901  case (2)
1902  call fms_f2c_string(keys3(key3_i)%key2, fileptr%file_global_meta(2,1))
1903  call fms_f2c_string(vals3(key3_i)%val2, fileptr%file_global_meta(2,2))
1904  case (3)
1905  call fms_f2c_string(keys3(key3_i)%key3, fileptr%file_global_meta(3,1))
1906  call fms_f2c_string(vals3(key3_i)%val3, fileptr%file_global_meta(3,2))
1907  case (4)
1908  call fms_f2c_string(keys3(key3_i)%key4, fileptr%file_global_meta(4,1))
1909  call fms_f2c_string(vals3(key3_i)%val4, fileptr%file_global_meta(4,2))
1910  case (5)
1911  call fms_f2c_string(keys3(key3_i)%key5, fileptr%file_global_meta(5,1))
1912  call fms_f2c_string(vals3(key3_i)%val5, fileptr%file_global_meta(5,2))
1913  case (6)
1914  call fms_f2c_string(keys3(key3_i)%key6, fileptr%file_global_meta(6,1))
1915  call fms_f2c_string(vals3(key3_i)%val6, fileptr%file_global_meta(6,2))
1916  case (7)
1917  call fms_f2c_string(keys3(key3_i)%key7, fileptr%file_global_meta(7,1))
1918  call fms_f2c_string(vals3(key3_i)%val7, fileptr%file_global_meta(7,2))
1919  case (8)
1920  call fms_f2c_string(keys3(key3_i)%key8, fileptr%file_global_meta(8,1))
1921  call fms_f2c_string(vals3(key3_i)%val8, fileptr%file_global_meta(8,2))
1922  case (9)
1923  call fms_f2c_string(keys3(key3_i)%key9, fileptr%file_global_meta(9,1))
1924  call fms_f2c_string(vals3(key3_i)%val9, fileptr%file_global_meta(9,2))
1925  case (10)
1926  call fms_f2c_string(keys3(key3_i)%key10, fileptr%file_global_meta(10,1))
1927  call fms_f2c_string(vals3(key3_i)%val10, fileptr%file_global_meta(10,2))
1928  case (11)
1929  call fms_f2c_string(keys3(key3_i)%key11, fileptr%file_global_meta(11,1))
1930  call fms_f2c_string(vals3(key3_i)%val11, fileptr%file_global_meta(11,2))
1931  case (12)
1932  call fms_f2c_string(keys3(key3_i)%key12, fileptr%file_global_meta(12,1))
1933  call fms_f2c_string(vals3(key3_i)%val12, fileptr%file_global_meta(12,2))
1934  case (13)
1935  call fms_f2c_string(keys3(key3_i)%key13, fileptr%file_global_meta(13,1))
1936  call fms_f2c_string(vals3(key3_i)%val13, fileptr%file_global_meta(13,2))
1937  case (14)
1938  call fms_f2c_string(keys3(key3_i)%key14, fileptr%file_global_meta(14,1))
1939  call fms_f2c_string(vals3(key3_i)%val14, fileptr%file_global_meta(14,2))
1940  case (15)
1941  call fms_f2c_string(keys3(key3_i)%key15, fileptr%file_global_meta(15,1))
1942  call fms_f2c_string(vals3(key3_i)%val15, fileptr%file_global_meta(15,2))
1943  case (16)
1944  call fms_f2c_string(keys3(key3_i)%key16, fileptr%file_global_meta(16,1))
1945  call fms_f2c_string(vals3(key3_i)%val16, fileptr%file_global_meta(16,2))
1946  end select
1947  enddo
1948  endif
1949  enddo
1950  tier2size = i
1951 
1952  call write_yaml_from_struct_3( 'diag_out.yaml'//c_null_char, 1, keys, vals, &
1953  SIZE(diag_yaml%diag_files), keys2, vals2, &
1954  tier3size, tier3each, keys3, vals3, &
1955  (/size(diag_yaml%diag_files), 0, 0, 0, 0, 0, 0, 0/))
1956  deallocate( keys, keys2, keys3, vals, vals2, vals3)
1957 
1958 end subroutine
1959 
1960 !> private function for getting unit string from diag_data parameter values
1961 pure function get_diag_unit_string( unit_param )
1962  integer, intent(in) :: unit_param(:) !< diag unit parameter values from diag_data_mod.
1963  !! <br>eg. DIAG_SECONDS, DIAG_MINUTES,DIAG_HOURS, DIAG_DAYS, DIAG_YEARS
1964  character(len=8 * SIZE(unit_param)) :: get_diag_unit_string
1965  character(len=7) :: tmp
1966  integer :: i
1967  get_diag_unit_string = ' '
1968  do i=1, SIZE(unit_param)
1969  select case(unit_param(i))
1970  case (diag_seconds)
1971  tmp = 'seconds'
1972  case (diag_minutes)
1973  tmp = 'minutes'
1974  case (diag_hours)
1975  tmp = 'hours'
1976  case (diag_days)
1977  tmp = 'days'
1978  case (diag_months)
1979  tmp = 'months'
1980  case (diag_years)
1981  tmp = 'years'
1982  case default
1983  exit
1984  end select
1985  get_diag_unit_string = trim(get_diag_unit_string)//" "//trim(tmp)
1986  enddo
1988 end function
1989 
1990 !> private function for getting reduction type string from parameter values
1991 pure function get_diag_reduction_string( reduction_val )
1992  integer, intent(in) :: reduction_val(:) !< reduction types (eg. time_average)
1993  integer :: i
1994  character(len=8 * MAX_FREQ) :: get_diag_reduction_string
1995  character(len=7) :: tmp
1997  do i=1, SIZE(reduction_val)
1998  select case (reduction_val(i))
1999  case (time_none)
2000  tmp = 'none'
2001  case (time_average)
2002  tmp = 'average'
2003  case (time_min)
2004  tmp = 'min'
2005  case (time_max)
2006  tmp = 'max'
2007  case (time_rms)
2008  tmp = 'rms'
2009  case (time_sum)
2010  tmp = 'sum'
2011  case (time_diurnal)
2012  tmp = 'diurnal'
2013  case default
2014  exit
2015  end select
2016  get_diag_reduction_string = trim(get_diag_reduction_string) //" "//trim(tmp)
2017  enddo
2019 end function
2020 
2021 subroutine add_axis_name( this, axis_name )
2022  class(diagyamlfilesvar_type), intent(inout) :: this
2023  character(len=:), allocatable, intent(in) :: axis_name
2024  character(len=:), allocatable :: tmp_str
2025 
2026  this%var_axes_names = trim(axis_name)//" "//trim(this%var_axes_names)
2027 
2028 end subroutine add_axis_name
2029 
2030 pure function is_file_subregional( this ) &
2031  result(res)
2032  class(diagyamlfilesvar_type), intent(in) :: this
2033  logical :: res
2034 
2035  res = this%var_file_is_subregional
2036 end function is_file_subregional
2037 
2038 #endif
2039 end module fms_diag_yaml_mod
2040 !> @}
2041 ! close documentation grouping
integer, parameter max_str_len
Max length for a string.
Definition: diag_data.F90:129
integer, parameter end_time
Use the end of the time average bounds.
Definition: diag_data.F90:128
integer, parameter time_min
The reduction method is min value.
Definition: diag_data.F90:117
subroutine set_base_time(base_time_int)
Set the module variable base_time.
Definition: diag_data.F90:463
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 begin_time
Use the begining of the time average bounds.
Definition: diag_data.F90:126
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 middle_time
Use the middle of the time average bounds.
Definition: diag_data.F90:127
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 pure function size_file_varlist(this)
Finds the number of variables in the file_varlist.
pure logical function has_diag_basedate(this)
diag_file_objdiag_basedate is on the stack, so this is always true
type(diagyamlobject_type) function, pointer, public get_diag_yaml_obj()
gets the diag_yaml module variable
type(diagyamlfilesvar_type) function, dimension(:), allocatable, public get_diag_fields_entries(indices)
Gets the diag_field entries corresponding to the indices of the sorted variable_list.
pure integer function get_filename_time(this)
Get the integer equivalent of the time to use to determine the filename, if using a wildcard file nam...
pure logical function has_file_duration_units(this)
diag_file_objfile_duration_units is on the stack, so this will retrun true
pure logical function has_var_reduction(this)
Checks if diag_file_objvar_reduction is allocated.
pure logical function has_var_attributes(this)
Checks if diag_file_objvar_attributes is allocated.
subroutine set_filename_time(yaml_fileobj, filename_time)
This checks if the filename_time in a diag file is correct and sets the integer equivalent.
subroutine fill_in_diag_fields(diag_file_id, yaml_fileobj, var_id, field, allow_averages)
Fills in a diagYamlFilesVar_type with the contents of a variable block in diag_table....
pure logical function has_n_diurnal(this)
Checks if diag_file_objn_diurnal is set.
subroutine, public diag_yaml_object_end()
Destroys the diag_yaml object.
pure integer function get_pow_value(this)
Inquiry for diag_yaml_files_var_objpow_value.
pure integer function get_file_new_file_freq_units(this)
Inquiry for diag_files_objfile_new_file_freq_units.
pure logical function has_var_varname(this)
Checks if diag_file_objvar_varname is allocated.
subroutine set_file_time_units(yaml_fileobj, file_timeunit)
This checks if the time unit in a diag file is valid and sets the integer equivalent.
pure logical function has_file_timeunit(this)
Checks if diag_file_objfile_timeunit is allocated.
pure integer function size_diag_files(this)
Find the number of files listed in the diag yaml.
pure logical function has_file_start_time(this)
Checks if diag_file_objfile_start_time is allocated.
pure logical function has_file_sub_region(this)
Checks if diag_file_objfile_sub_region is being used and has the sub region variables allocated.
pure integer function, allocatable get_var_reduction(this)
Inquiry for diag_yaml_files_var_objvar_reduction.
pure integer function get_file_duration_units(this)
Inquiry for diag_files_objfile_duration_units.
type(diagyamlfiles_type) function, dimension(:), allocatable get_diag_files(this)
get the diag_files of a diag_yaml type
pure logical function has_file_global_meta(this)
Checks if diag_file_objfile_global_meta is allocated.
pure logical function has_var_longname(this)
Checks if diag_file_objvar_longname is allocated.
pure character(len=:) function, allocatable get_var_module(this)
Inquiry for diag_yaml_files_var_objvar_module.
pure logical function has_file_frequnit(this)
Checks if diag_file_objfile_frequnit is allocated.
logical function is_global_meta(this)
Inquiry for whether file_global_meta is allocated.
pure real(kind=r4_kind) function, dimension(2) get_var_zbounds(this)
Inquiry for diag_yaml_files_var_objvar_zbounds.
pure integer function, dimension(basedate_size) get_basedate(this)
get the basedate of a diag_yaml type
integer function set_valid_time_units(time_units, error_msg)
This checks if a time unit is valid and if it is, it assigns the integer equivalent.
pure logical function has_diag_title(this)
Checks if diag_file_objdiag_title is allocated.
subroutine get_sub_region(diag_yaml_id, sub_region_id, sub_region, grid_type, fname)
gets the lat/lon of the sub region to use in a diag_table yaml
pure character(len=:) function, allocatable get_var_fname(this)
Inquiry for diag_yaml_files_var_objvar_fname.
pure integer function get_file_timeunit(this)
Inquiry for diag_files_objfile_timeunit.
subroutine increase_new_file_freq_index(this)
Increate the current_new_file_freq_index by 1.
subroutine diag_yaml_files_obj_init(obj)
Initializes the non string values of a diagYamlFiles_type to its default values.
integer function, dimension(:), allocatable, public get_diag_field_ids(indices)
Gets field indices corresponding to the indices (input argument) in the sorted variable_list.
pure character(len=8 *size(unit_param)) function get_diag_unit_string(unit_param)
private function for getting unit string from diag_data parameter values
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.
pure logical function has_var_write(this)
diag_file_objvar_write is on the stack, so this returns true
logical function is_var_attributes(this)
Inquiry for whether var_attributes is allocated.
pure character(len=:) function, allocatable get_var_units(this)
Inquiry for diag_yaml_files_var_objvar_units.
subroutine, public dump_diag_yaml_obj(filename)
Prints out values from diag_yaml object for debugging. Only writes on root.
pure integer function, allocatable get_var_kind(this)
Inquiry for diag_yaml_files_var_objvar_kind.
integer function, public get_num_unique_fields()
Determine the number of unique diag_fields in the diag_yaml_object.
pure character(:) function, dimension(:), allocatable get_file_varlist(this)
Inquiry for diag_files_objfile_varlist.
integer function get_total_num_vars(diag_yaml_id, diag_file_id)
gets the total number of variables in the diag_table yaml file
pure integer function get_file_frequnit(this)
Inquiry for diag_files_objfile_frequnit.
type(diagyamlfilesvar_type) function, pointer get_diag_field_from_id(this, yaml_id)
Get the diag_field yaml corresponding to a yaml_id.
pure character(len=:) function, allocatable get_file_fname(this)
Inquiry for diag_files_objfile_fname.
pure logical function has_var_outname(this)
Checks if diag_file_objvar_outname is allocated.
pure logical function has_file_new_file_freq_units(this)
Checks if diag_file_objfile_new_file_freq_units is allocated.
pure character(len=max_str_len) function, dimension(:,:), allocatable get_file_global_meta(this)
Inquiry for diag_files_objfile_global_meta.
pure character(len=:) function, allocatable get_file_start_time(this)
Inquiry for diag_files_objfile_start_time.
pure character(len=:) function, allocatable get_var_varname(this)
Inquiry for diag_yaml_files_var_objvar_varname.
pure integer function get_file_new_file_freq(this)
Inquiry for diag_files_objfile_new_file_freq.
subroutine set_field_reduction(field, reduction_method)
This checks if the reduction of a diag field is valid and sets it If the reduction method is diurnalX...
pure logical function has_diag_fields(this)
Checks if diag_file_objdiag_fields is allocated.
pure integer function get_n_diurnal(this)
Inquiry for diag_yaml_files_var_objn_diurnal.
pure logical function has_file_freq(this)
diag_file_objfile_freq is on the stack, so the object always has it
subroutine set_field_kind(field, skind)
This checks if the kind of a diag field is valid and sets it.
pure logical function has_file_new_file_freq(this)
diag_file_objfile_new_file_freq is defined on the stack, so this will return true
pure type(diagyamlfilesvar_type) function, dimension(:), allocatable get_diag_fields(this)
get the diag_fields of a diag_yaml type
pure logical function has_var_module(this)
Checks if diag_file_objvar_module is allocated.
pure character(len=:) function, allocatable get_file_unlimdim(this)
Inquiry for diag_files_objfile_unlimdim.
pure logical function has_file_varlist(this)
Checks if diag_file_objfile_varlist is allocated.
subroutine, public fms_diag_yaml_out()
Writes an output yaml with all available information on the written files. Will only write with root ...
pure character(len=:) function, allocatable get_var_outname(this)
Inquiry for diag_yaml_files_var_objvar_outname.
pure character(len=:) function, allocatable get_var_longname(this)
Inquiry for diag_yaml_files_var_objvar_longname.
pure character(len=8 *max_freq) function get_diag_reduction_string(reduction_val)
private function for getting reduction type string from parameter values
pure logical function has_diag_files(this)
Checks if diag_file_objdiag_files is allocated.
pure character(len=max_str_len) function, dimension(:,:), allocatable get_var_attributes(this)
Inquiry for diag_yaml_files_var_objvar_attributes.
pure logical function has_var_kind(this)
Checks if diag_file_objvar_kind is allocated.
pure integer function get_file_duration(this)
Inquiry for diag_files_objfile_duration.
subroutine fill_in_diag_files(diag_yaml_id, diag_file_id, yaml_fileobj)
Fills in a diagYamlFiles_type with the contents of a file block in diag_table.yaml.
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 logical function has_var_zbounds(this)
Checks if diag_file_objvar_zbounds is allocated.
pure character(len=:) function, allocatable get_title(this)
get the title of a diag_yaml type
type(subregion_type) function, pointer get_file_sub_region(this)
Inquiry for diag_files_objfile_subregion.
pure integer function get_file_freq(this)
Inquiry for diag_files_objfile_freq.
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 parse_key(filename, buffer, file_freq, file_frequnit, var)
This parses the freq, new_file_freq, or file_duration keys which are read in as a comma list.
pure logical function has_file_write(this)
Checks if diag_file_objfile_write is on the stack, so this will always be true.
pure logical function has_pow_value(this)
Checks if diag_file_objpow_value is set.
pure logical function has_var_units(this)
Checks if diag_file_objvar_units is allocated.
pure logical function has_var_fname(this)
Checks if diag_file_objvar_fname is allocated.
pure logical function has_file_unlimdim(this)
Checks if diag_file_objfile_unlimdim is allocated.
subroutine diag_get_value_from_key(diag_file_id, par_id, key_name, value_name, is_optional)
diag_manager wrapper to get_value_from_key to use for allocatable string variables
pure logical function has_file_fname(this)
Checks if diag_file_objfile_fname is allocated.
pure logical function has_file_duration(this)
diag_file_objfile_duration is allocated on th stack, so this is always true
Object that holds the information of the diag_yaml.
character(:) function, allocatable, public string(v, fmt)
Converts a number or a Boolean value to a string.
type(c_ptr) function, dimension(:), allocatable, public fms_array_to_pointer(my_array)
Converts a character array to an array of c pointers!
integer function, dimension(:), allocatable, public fms_find_my_string(my_pointer, narray, string_to_find)
Searches through a SORTED array of pointers for a string.
subroutine, public fms_f2c_string(dest, str_in)
Copies a Fortran string into a C string and puts c_null_char in any trailing spaces.
subroutine, public initialize_val_struct(yv)
subroutine, public initialize_key_struct(yk)
Keys for the output yaml on a given level corresponding to the struct in yaml_output_functions....
Values for the output yaml on a given level corresponding to the struct in yaml_output_functions....
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_nkeys(file_id, block_id)
Gets the number of key-value pairs in a block.
subroutine, public get_key_name(file_id, key_id, key_name)
Gets the key from a file id.
integer function, public open_and_parse_file(filename)
Opens and parses a yaml file.
subroutine, public get_key_ids(file_id, block_id, key_ids)
Gets the ids of the key-value pairs in a block.
subroutine, public get_block_ids(file_id, block_name, block_ids, parent_block_id)
Gets the the ids of the blocks with block_name in the yaml file If parent_block_id is present,...
integer function, public get_num_blocks(file_id, block_name, parent_block_id)
Determines the number of blocks with block_name in the yaml file If parent_block_id is present,...
subroutine, public get_key_value(file_id, key_id, key_value)
Gets the value from a file id.
Dermine the value of a key from a keyname.
Definition: yaml_parser.F90:57
c function that finds the number of unique strings in an array of c pointers
Sorts an array of pointers (my pointer) of size (p_size) in alphabetical order.
type to hold the diag_file information
type to hold the info a diag_field
type to hold an array of sorted diag_files
type to hold the sub region information about a file
type to hold an array of sorted diag_fiels