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