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