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