1 !***********************************************************************
2 !* GNU Lesser General Public License
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
20 !> @defgroup fms_diag_yaml_mod fms_diag_yaml_mod
21 !> @ingroup diag_manager
22 !! @brief fms_diag_yaml_mod is an integral part of
23 !! diag_manager_mod. Its function is to read the diag_table.yaml to fill in
24 !! the diag_yaml_object
27 !> @brief File for @ref diag_yaml_mod
29 !> @addtogroup fms_diag_yaml_mod
31 module fms_diag_yaml_mod
33 use diag_data_mod, only: DIAG_NULL, DIAG_OCEAN, DIAG_ALL, DIAG_OTHER, set_base_time, latlon_gridtype, &
34 index_gridtype, null_gridtype, DIAG_SECONDS, DIAG_MINUTES, DIAG_HOURS, DIAG_DAYS, &
35 DIAG_MONTHS, DIAG_YEARS, time_average, time_rms, time_max, time_min, time_sum, &
36 time_diurnal, time_power, time_none, r8, i8, r4, i4, DIAG_NOT_REGISTERED, &
37 middle_time, begin_time, end_time, MAX_STR_LEN
38 use yaml_parser_mod, only: open_and_parse_file, get_value_from_key, get_num_blocks, get_nkeys, &
39 get_block_ids, get_key_value, get_key_ids, get_key_name
40 use mpp_mod, only: mpp_error, FATAL, NOTE, mpp_pe, mpp_root_pe, stdout
41 use, intrinsic :: iso_c_binding, only : c_ptr, c_null_char
42 use fms_string_utils_mod, only: fms_array_to_pointer, fms_find_my_string, fms_sort_this, fms_find_unique, string
43 use platform_mod, only: r4_kind, i4_kind
44 use fms_mod, only: lowercase
51 public :: diag_yaml_object_init, diag_yaml_object_end
52 public :: diagYamlObject_type, get_diag_yaml_obj, subRegion_type
53 public :: diagYamlFiles_type, diagYamlFilesVar_type
54 public :: get_num_unique_fields, find_diag_field, get_diag_fields_entries, get_diag_files_id
55 public :: get_diag_field_ids
56 public :: dump_diag_yaml_obj
60 integer, parameter :: basedate_size = 6
61 integer, parameter :: NUM_SUB_REGION_ARRAY = 8
62 integer, parameter :: MAX_FREQ = 12
63 integer :: MAX_SUBAXES = 0 !< Max number of subaxis, set in diag_yaml_object_init depending on
64 !! what is in the diag yaml
67 !> @brief type to hold an array of sorted diag_fiels
69 character(len=255), allocatable :: var_name(:) !< Array of diag_field
70 type(c_ptr), allocatable :: var_pointer(:) !< Array of pointers
71 integer, allocatable :: diag_field_indices(:) !< Index of the field in the diag_field array
74 !> @brief type to hold an array of sorted diag_files
76 character(len=255), allocatable :: file_name(:) !< Array of diag_field
77 type(c_ptr), allocatable :: file_pointer(:) !< Array of pointers
78 integer, allocatable :: diag_file_indices(:) !< Index of the file in the diag_file array
81 !> @brief type to hold the sub region information about a file
83 INTEGER :: grid_type !< Flag indicating the type of region,
84 !! acceptable values are latlon_gridtype, index_gridtype,
86 class(*), allocatable :: corners(:,:)!< (x, y) coordinates of the four corner of the region
87 integer :: tile !< Tile number of the sub region
88 !! required if using the "index" grid type
90 end type subRegion_type
92 !> @brief type to hold the diag_file information
93 type diagYamlFiles_type
95 character (len=:), allocatable :: file_fname !< file name
96 integer :: file_frequnit(MAX_FREQ) !< the frequency unit (DIAG_SECONDS,
97 !! DIAG_MINUTES, DIAG_HOURS, DIAG_DAYS,
99 integer :: file_freq(MAX_FREQ) !< the frequency of data
100 integer :: file_timeunit !< The unit of time (DIAG_SECONDS,
101 !! DIAG_MINUTES, DIAG_HOURS, DIAG_DAYS,
103 character (len=:), allocatable :: file_unlimdim !< The name of the unlimited dimension
104 type(subRegion_type) :: file_sub_region !< type containing info about the subregion
105 integer :: file_new_file_freq(MAX_FREQ) !< Frequency for closing the existing file
106 integer :: file_new_file_freq_units(MAX_FREQ) !< Time units for creating a new file.
107 !! Required if “new_file_freq” used
108 !! (DIAG_SECONDS, DIAG_MINUTES, &
109 !! DIAG_HOURS, DIAG_DAYS, DIAG_YEARS)
110 character (len=:), allocatable :: file_start_time !< Time to start the file for the
111 !! first time. Requires “new_file_freq”
112 integer :: filename_time !< The time to use when setting the name of
113 !! new files: begin, middle, or end of the
115 integer :: file_duration(MAX_FREQ) !< How long the file should receive data
116 !! after start time in file_duration_units.
117 !! This optional field can only be used if
118 !! the start_time field is present. If this
119 !! field is absent, then the file duration
120 !! will be equal to the frequency for
121 !! creating new files. NOTE: The
122 !! file_duration_units field must also
123 !! be present if this field is present.
124 integer :: file_duration_units(MAX_FREQ) !< The file duration units
125 !! (DIAG_SECONDS, DIAG_MINUTES, &
126 !! DIAG_HOURS, DIAG_DAYS, DIAG_YEARS)
127 integer :: current_new_file_freq_index !< The index of the new_file_freq array
128 !< Need to use `MAX_STR_LEN` because not all filenames/global attributes are the same length
129 character (len=MAX_STR_LEN), allocatable :: file_varlist(:) !< An array of variable names
131 character (len=MAX_STR_LEN), allocatable :: file_global_meta(:,:) !< Array of key(dim=1)
132 !! and values(dim=2) to be
133 !! added as global meta data to
137 !> All getter functions (functions named get_x(), for member field named x)
138 !! return copies of the member variables unless explicitly noted.
139 procedure, public :: size_file_varlist
140 procedure, public :: get_file_fname
141 procedure, public :: get_file_frequnit
142 procedure, public :: get_file_freq
143 procedure, public :: get_file_timeunit
144 procedure, public :: get_file_unlimdim
145 procedure, public :: get_file_sub_region
146 procedure, public :: get_file_new_file_freq
147 procedure, public :: get_file_new_file_freq_units
148 procedure, public :: get_file_start_time
149 procedure, public :: get_file_duration
150 procedure, public :: get_file_duration_units
151 procedure, public :: get_file_varlist
152 procedure, public :: get_file_global_meta
153 procedure, public :: get_filename_time
154 procedure, public :: is_global_meta
155 !> Has functions to determine if allocatable variables are true. If a variable is not an allocatable
156 !! then is will always return .true.
157 procedure, public :: has_file_fname
158 procedure, public :: has_file_frequnit
159 procedure, public :: has_file_freq
160 procedure, public :: has_file_timeunit
161 procedure, public :: has_file_unlimdim
162 procedure, public :: has_file_sub_region
163 procedure, public :: has_file_new_file_freq
164 procedure, public :: has_file_new_file_freq_units
165 procedure, public :: has_file_start_time
166 procedure, public :: has_file_duration
167 procedure, public :: has_file_duration_units
168 procedure, public :: has_file_varlist
169 procedure, public :: has_file_global_meta
170 procedure, public :: increase_new_file_freq_index
171 end type diagYamlFiles_type
173 !> @brief type to hold the info a diag_field
174 type diagYamlFilesVar_type
175 character (len=:), private, allocatable :: var_fname !< The field/diagnostic name
176 character (len=:), private, allocatable :: var_varname !< The name of the variable
177 integer , private, allocatable :: var_reduction !< Reduction to be done on var
178 !! time_average, time_rms, time_max,
179 !! time_min, time_sum, time_diurnal, time_power
180 character (len=:), private, allocatable :: var_module !< The module that th variable is in
181 integer , private, allocatable :: var_kind !< The type/kind of the variable
182 character (len=:), private, allocatable :: var_outname !< Name of the variable as written to the file
183 character (len=:), private, allocatable :: var_longname !< Overwrites the long name of the variable
184 character (len=:), private, allocatable :: var_units !< Overwrites the units
185 real(kind=r4_kind), private :: var_zbounds(2) !< The z axis limits [vert_min, vert_max]
186 integer , private :: n_diurnal !< Number of diurnal samples
187 !! 0 if var_reduction is not "diurnalXX"
188 integer , private :: pow_value !< The power value
189 !! 0 if pow_value is not "powXX"
191 !< Need to use `MAX_STR_LEN` because not all filenames/global attributes are the same length
192 character (len=MAX_STR_LEN), dimension (:, :), private, allocatable :: var_attributes !< Attributes to overwrite or
193 !! add from diag_yaml
195 !> All getter functions (functions named get_x(), for member field named x)
196 !! return copies of the member variables unless explicitly noted.
197 procedure :: get_var_fname
198 procedure :: get_var_varname
199 procedure :: get_var_reduction
200 procedure :: get_var_module
201 procedure :: get_var_kind
202 procedure :: get_var_outname
203 procedure :: get_var_longname
204 procedure :: get_var_units
205 procedure :: get_var_zbounds
206 procedure :: get_var_attributes
207 procedure :: get_n_diurnal
208 procedure :: get_pow_value
209 procedure :: is_var_attributes
211 procedure :: has_var_fname
212 procedure :: has_var_varname
213 procedure :: has_var_reduction
214 procedure :: has_var_module
215 procedure :: has_var_kind
216 procedure :: has_var_outname
217 procedure :: has_var_longname
218 procedure :: has_var_units
219 procedure :: has_var_zbounds
220 procedure :: has_var_attributes
221 procedure :: has_n_diurnal
222 procedure :: has_pow_value
224 end type diagYamlFilesVar_type
226 !> @brief Object that holds the information of the diag_yaml
227 !> @ingroup fms_diag_yaml_mod
228 type diagYamlObject_type
229 character(len=:), allocatable, private :: diag_title !< Experiment name
230 integer, private, dimension (basedate_size) :: diag_basedate !< basedate array
231 type(diagYamlFiles_type), allocatable, public, dimension (:) :: diag_files!< History file info
232 type(diagYamlFilesVar_type), allocatable, public, dimension (:) :: diag_fields !< Diag fields info
234 procedure :: size_diag_files
236 procedure :: get_title !< Returns the title
237 procedure :: get_basedate !< Returns the basedate array
238 procedure :: get_diag_files !< Returns the diag_files array
239 procedure :: get_diag_fields !< Returns the diag_field array
240 procedure :: get_diag_field_from_id
242 procedure :: has_diag_title
243 procedure :: has_diag_basedate
244 procedure :: has_diag_files
245 procedure :: has_diag_fields
247 end type diagYamlObject_type
249 type (diagYamlObject_type), target :: diag_yaml !< Obj containing the contents of the diag_table.yaml
250 type (varList_type), save :: variable_list !< List of all the variables in the diag_table.yaml
251 type (fileList_type), save :: file_list !< List of all files in the diag_table.yaml
253 logical, private :: diag_yaml_module_initialized = .false.
256 !> @addtogroup fms_diag_yaml_mod
260 !> @brief gets the diag_yaml module variable
261 !! @return a copy of the diag_yaml module variable
262 function get_diag_yaml_obj() &
264 type (diagYamlObject_type) :: res
267 end function get_diag_yaml_obj
269 !> @brief get the basedate of a diag_yaml type
270 !! @return the basedate as an integer array
271 pure function get_basedate (this) &
272 result (diag_basedate)
273 class (diagYamlObject_type), intent(in) :: this !< The diag_yaml
274 integer, dimension (basedate_size) :: diag_basedate !< Basedate array result to return
276 diag_basedate = this%diag_basedate
277 end function get_basedate
279 !> @brief Find the number of files listed in the diag yaml
280 !! @return the number of files in the diag yaml
281 pure integer function size_diag_files(this)
282 class (diagYamlObject_type), intent(in) :: this !< The diag_yaml
283 if (this%has_diag_files()) then
284 size_diag_files = size(this%diag_files)
288 end function size_diag_files
290 !> @brief get the title of a diag_yaml type
291 !! @return the title of the diag table as an allocated string
292 pure function get_title (this) &
294 class (diagYamlObject_type), intent(in) :: this !< The diag_yaml
295 character(len=:),allocatable :: diag_title !< Basedate array result to return
297 diag_title = this%diag_title
298 end function get_title
300 !> @brief get the diag_files of a diag_yaml type
301 !! @return the diag_files
302 function get_diag_files(this) &
304 class (diagYamlObject_type), intent(in) :: this !< The diag_yaml
305 type(diagYamlFiles_type), allocatable, dimension (:) :: diag_files !< History file info
307 diag_files = this%diag_files
308 end function get_diag_files
310 !> @brief Get the diag_field yaml corresponding to a yaml_id
311 !! @return Pointer to the diag_field yaml entry
312 function get_diag_field_from_id(this, yaml_id) &
314 class (diagYamlObject_type), target, intent(in) :: this !< The diag_yaml
315 integer, intent(in) :: yaml_id !< Yaml id
317 type(diagYamlFilesVar_type), pointer :: diag_field !< Diag fields info
319 if (yaml_id .eq. DIAG_NOT_REGISTERED) call mpp_error(FATAL, &
320 "Diag_manager: The yaml id for this field is not is not set")
322 diag_field => this%diag_fields(variable_list%diag_field_indices(yaml_id))
324 end function get_diag_field_from_id
326 !> @brief get the diag_fields of a diag_yaml type
327 !! @return the diag_fields
328 pure function get_diag_fields(this) &
330 class (diagYamlObject_type), intent(in) :: this !< The diag_yaml
331 type(diagYamlFilesVar_type), allocatable, dimension (:) :: diag_fields !< Diag fields info
333 diag_fields = this%diag_fields
334 end function get_diag_fields
336 !> @brief Uses the yaml_parser_mod to read in the diag_table and fill in the
338 subroutine diag_yaml_object_init(diag_subset_output)
339 integer, intent(in) :: diag_subset_output !< DIAG_ALL - Current PE is in the one and only pelist
340 !! DIAG_OTHER - Current PE is not in the ocean pelist
341 !! and there are multiple pelists
342 !! DIAG_OCEAN - Current PE is in the ocean pelist
343 !! and there are multiple pelists
344 integer :: diag_yaml_id !< Id for the diag_table yaml
345 integer :: nfiles !< Number of files in the diag_table yaml
346 integer, allocatable :: diag_file_ids(:) !< Ids of the files in the diag_table yaml
347 integer :: i, j !< For do loops
348 integer :: total_nvars !< The total number of variables in the diag_table yaml
349 integer :: var_count !< The current number of variables added to the diag_yaml obj
350 integer :: file_var_count !< The current number of variables added in the diag_file
351 integer :: nvars !< The number of variables in the current file
352 integer, allocatable :: var_ids(:) !< Ids of the variables in diag_table yaml
353 logical :: is_ocean !< Flag indicating if it is an ocean file
354 logical, allocatable :: ignore(:) !< Flag indicating if the diag_file is going to be ignored
355 integer :: actual_num_files !< The actual number of files that were saved
356 integer :: file_count !! The current number of files added to the diag_yaml obj
357 logical :: write_file !< Flag indicating if the user wants the file to be written
358 logical :: write_var !< Flag indicating if the user wants the variable to be written
359 character(len=:), allocatable :: filename!< Diag file name (for error messages)
361 if (diag_yaml_module_initialized) return
363 diag_yaml_id = open_and_parse_file("diag_table.yaml")
365 call diag_get_value_from_key(diag_yaml_id, 0, "title", diag_yaml%diag_title)
366 call get_value_from_key(diag_yaml_id, 0, "base_date", diag_yaml%diag_basedate)
367 call set_base_time(diag_yaml%diag_basedate)
369 nfiles = get_num_blocks(diag_yaml_id, "diag_files")
370 allocate(diag_file_ids(nfiles))
371 allocate(ignore(nfiles))
373 call get_block_ids(diag_yaml_id, "diag_files", diag_file_ids)
377 !< If you are on two seperate pelists
378 if(diag_subset_output .ne. DIAG_ALL) then
381 call get_value_from_key(diag_yaml_id, diag_file_ids(i), "is_ocean", is_ocean, is_optional=.true.)
382 !< If you are on the ocean pelist and the file is not an ocean file, skip the file
383 if (diag_subset_output .eq. DIAG_OCEAN .and. .not. is_ocean) ignore(i) = .true.
385 !< If you are not on the ocean pelist and the file is ocean, skip the file
386 if(diag_subset_output .eq. DIAG_OTHER .and. is_ocean) ignore(i) = .true.
390 !< Determine how many files are in the diag_yaml, ignoring those with write_file = False
394 call get_value_from_key(diag_yaml_id, diag_file_ids(i), "write_file", write_file, is_optional=.true.)
395 if(.not. write_file) ignore(i) = .true.
397 !< If ignoring the file, ignore the fields in that file too!
398 if (.not. ignore(i)) then
399 nvars = get_total_num_vars(diag_yaml_id, diag_file_ids(i))
400 total_nvars = total_nvars + nvars
401 if (nvars .ne. 0) then
402 actual_num_files = actual_num_files + 1
404 call diag_get_value_from_key(diag_yaml_id, diag_file_ids(i), "file_name", filename)
405 call mpp_error(NOTE, "diag_manager_mod:: the file:"//trim(filename)//" has no variables defined. Ignoring!")
411 allocate(diag_yaml%diag_files(actual_num_files))
412 allocate(diag_yaml%diag_fields(total_nvars))
413 allocate(variable_list%var_name(total_nvars))
414 allocate(variable_list%diag_field_indices(total_nvars))
415 allocate(file_list%file_name(actual_num_files))
416 allocate(file_list%diag_file_indices(actual_num_files))
420 !> Loop through the number of nfiles and fill in the diag_yaml obj
421 nfiles_loop: do i = 1, nfiles
423 file_count = file_count + 1
424 call diag_yaml_files_obj_init(diag_yaml%diag_files(file_count))
425 call fill_in_diag_files(diag_yaml_id, diag_file_ids(i), diag_yaml%diag_files(file_count))
427 !> Save the file name in the file_list
428 file_list%file_name(file_count) = trim(diag_yaml%diag_files(file_count)%file_fname)//c_null_char
429 file_list%diag_file_indices(file_count) = file_count
432 nvars = get_num_blocks(diag_yaml_id, "varlist", parent_block_id=diag_file_ids(i))
433 allocate(var_ids(nvars))
434 call get_block_ids(diag_yaml_id, "varlist", var_ids, parent_block_id=diag_file_ids(i))
436 allocate(diag_yaml%diag_files(file_count)%file_varlist(get_total_num_vars(diag_yaml_id, diag_file_ids(i))))
437 nvars_loop: do j = 1, nvars
439 call get_value_from_key(diag_yaml_id, var_ids(j), "write_var", write_var, is_optional=.true.)
440 if (.not. write_var) cycle
442 var_count = var_count + 1
443 file_var_count = file_var_count + 1
445 !> Save the filename in the diag_field type
446 diag_yaml%diag_fields(var_count)%var_fname = diag_yaml%diag_files(file_count)%file_fname
448 call fill_in_diag_fields(diag_yaml_id, var_ids(j), diag_yaml%diag_fields(var_count))
450 !> Save the variable name in the diag_file type
451 diag_yaml%diag_files(file_count)%file_varlist(file_var_count) = diag_yaml%diag_fields(var_count)%var_varname
453 !> Save the variable name and the module name in the variable_list
454 variable_list%var_name(var_count) = trim(diag_yaml%diag_fields(var_count)%var_varname)//&
455 ":"//trim(diag_yaml%diag_fields(var_count)%var_module)//c_null_char
456 !! The diag_table is not case sensitive (so we are saving it as lowercase)
457 variable_list%var_name(var_count) = lowercase(variable_list%var_name(var_count))
458 variable_list%diag_field_indices(var_count) = var_count
463 !> Sort the file list in alphabetical order
464 file_list%file_pointer = fms_array_to_pointer(file_list%file_name)
465 call fms_sort_this(file_list%file_pointer, actual_num_files, file_list%diag_file_indices)
467 variable_list%var_pointer = fms_array_to_pointer(variable_list%var_name)
468 call fms_sort_this(variable_list%var_pointer, total_nvars, variable_list%diag_field_indices)
470 deallocate(diag_file_ids)
471 diag_yaml_module_initialized = .true.
474 !> @brief Destroys the diag_yaml object
475 subroutine diag_yaml_object_end()
476 integer :: i !< For do loops
478 do i = 1, size(diag_yaml%diag_files, 1)
479 if(allocated(diag_yaml%diag_files(i)%file_varlist)) deallocate(diag_yaml%diag_files(i)%file_varlist)
480 if(allocated(diag_yaml%diag_files(i)%file_global_meta)) deallocate(diag_yaml%diag_files(i)%file_global_meta)
481 if(allocated(diag_yaml%diag_files(i)%file_sub_region%corners)) &
482 deallocate(diag_yaml%diag_files(i)%file_sub_region%corners)
484 if(allocated(diag_yaml%diag_files)) deallocate(diag_yaml%diag_files)
486 do i = 1, size(diag_yaml%diag_fields, 1)
487 if(allocated(diag_yaml%diag_fields(i)%var_attributes)) deallocate(diag_yaml%diag_fields(i)%var_attributes)
489 if(allocated(diag_yaml%diag_fields)) deallocate(diag_yaml%diag_fields)
491 if(allocated(file_list%file_pointer)) deallocate(file_list%file_pointer)
492 if(allocated(file_list%file_name)) deallocate(file_list%file_name)
493 if(allocated(file_list%diag_file_indices)) deallocate(file_list%diag_file_indices)
495 if(allocated(variable_list%var_pointer)) deallocate(variable_list%var_pointer)
496 if(allocated(variable_list%var_name)) deallocate(variable_list%var_name)
497 if(allocated(variable_list%diag_field_indices)) deallocate(variable_list%diag_field_indices)
499 end subroutine diag_yaml_object_end
501 !> @brief Fills in a diagYamlFiles_type with the contents of a file block in diag_table.yaml
502 subroutine fill_in_diag_files(diag_yaml_id, diag_file_id, yaml_fileobj)
503 integer, intent(in) :: diag_yaml_id !< Id of the diag_table.yaml
504 integer, intent(in) :: diag_file_id !< Id of the file block to read
505 type(diagYamlFiles_type), intent(inout) :: yaml_fileobj !< diagYamlFiles_type obj to read the contents into
507 integer :: nsubregion !< Flag indicating of there any regions (0 or 1)
508 integer :: sub_region_id(1) !< Id of the sub_region block
509 integer :: natt !< Number of global attributes in the current file
510 integer :: global_att_id(1) !< Id of the global attributes block
511 integer :: nkeys !< Number of key/value global attributes pair
512 integer :: j !< For do loops
514 integer, allocatable :: key_ids(:) !< Id of the gloabl atttributes key/value pairs
515 character(len=:), ALLOCATABLE :: grid_type !< grid_type as it is read in from the yaml
516 character(len=:), ALLOCATABLE :: buffer !< buffer to store any *_units as it is read from the yaml
518 call diag_get_value_from_key(diag_yaml_id, diag_file_id, "file_name", yaml_fileobj%file_fname)
519 call diag_get_value_from_key(diag_yaml_id, diag_file_id, "freq", buffer)
520 call parse_key(yaml_fileobj%file_fname, buffer, yaml_fileobj%file_freq, yaml_fileobj%file_frequnit, "freq")
523 call diag_get_value_from_key(diag_yaml_id, diag_file_id, "unlimdim", yaml_fileobj%file_unlimdim)
524 call diag_get_value_from_key(diag_yaml_id, diag_file_id, "time_units", buffer)
525 call set_file_time_units(yaml_fileobj, buffer)
528 call diag_get_value_from_key(diag_yaml_id, diag_file_id, "new_file_freq", buffer, is_optional=.true.)
529 call parse_key(yaml_fileobj%file_fname, buffer, yaml_fileobj%file_new_file_freq, &
530 yaml_fileobj%file_new_file_freq_units, "new_file_freq")
533 call diag_get_value_from_key(diag_yaml_id, diag_file_id, "filename_time", buffer, is_optional=.true.)
534 call set_filename_time(yaml_fileobj, buffer)
537 call diag_get_value_from_key(diag_yaml_id, diag_file_id, "start_time", &
538 yaml_fileobj%file_start_time, is_optional=.true.)
539 call diag_get_value_from_key(diag_yaml_id, diag_file_id, "file_duration", buffer, is_optional=.true.)
540 call parse_key(yaml_fileobj%file_fname, buffer, yaml_fileobj%file_duration, yaml_fileobj%file_duration_units, &
544 nsubregion = get_num_blocks(diag_yaml_id, "sub_region", parent_block_id=diag_file_id)
545 if (nsubregion .eq. 1) then
546 MAX_SUBAXES = MAX_SUBAXES + 1
547 call get_block_ids(diag_yaml_id, "sub_region", sub_region_id, parent_block_id=diag_file_id)
548 call diag_get_value_from_key(diag_yaml_id, sub_region_id(1), "grid_type", grid_type)
549 call get_sub_region(diag_yaml_id, sub_region_id(1), yaml_fileobj%file_sub_region, grid_type, &
550 yaml_fileobj%file_fname)
551 elseif (nsubregion .eq. 0) then
552 yaml_fileobj%file_sub_region%grid_type = null_gridtype
554 call mpp_error(FATAL, "diag_yaml_object_init: file "//trim(yaml_fileobj%file_fname)//" has multiple region blocks")
558 natt = get_num_blocks(diag_yaml_id, "global_meta", parent_block_id=diag_file_id)
559 if (natt .eq. 1) then
560 call get_block_ids(diag_yaml_id, "global_meta", global_att_id, parent_block_id=diag_file_id)
561 nkeys = get_nkeys(diag_yaml_id, global_att_id(1))
562 allocate(key_ids(nkeys))
563 call get_key_ids(diag_yaml_id, global_att_id(1), key_ids)
565 allocate(yaml_fileobj%file_global_meta(nkeys, 2))
567 call get_key_name(diag_yaml_id, key_ids(j), yaml_fileobj%file_global_meta(j, 1))
568 call get_key_value(diag_yaml_id, key_ids(j), yaml_fileobj%file_global_meta(j, 2))
571 elseif (natt .ne. 0) then
572 call mpp_error(FATAL, "diag_yaml_object_init: file "//trim(yaml_fileobj%file_fname)//&
573 &" has multiple global_meta blocks")
578 !> @brief Fills in a diagYamlFilesVar_type with the contents of a variable block in
580 subroutine fill_in_diag_fields(diag_file_id, var_id, field)
581 integer, intent(in) :: diag_file_id !< Id of the file block in the yaml file
582 integer, intent(in) :: var_id !< Id of the variable block in the yaml file
583 type(diagYamlFilesVar_type), intent(inout) :: field !< diagYamlFilesVar_type obj to read the contents into
585 integer :: natt !< Number of attributes in variable
586 integer :: var_att_id(1) !< Id of the variable attribute block
587 integer :: nkeys !< Number of key/value pairs of attributes
588 integer :: j !< For do loops
590 integer, allocatable :: key_ids(:) !< Id of each attribute key/value pair
591 character(len=:), ALLOCATABLE :: buffer !< buffer to store the reduction method as it is read from the yaml
593 call diag_get_value_from_key(diag_file_id, var_id, "var_name", field%var_varname)
594 call diag_get_value_from_key(diag_file_id, var_id, "reduction", buffer)
595 call set_field_reduction(field, buffer)
597 call diag_get_value_from_key(diag_file_id, var_id, "module", field%var_module)
599 call diag_get_value_from_key(diag_file_id, var_id, "kind", buffer)
600 call set_field_kind(field, buffer)
602 call diag_get_value_from_key(diag_file_id, var_id, "output_name", field%var_outname, is_optional=.true.)
603 call diag_get_value_from_key(diag_file_id, var_id, "long_name", field%var_longname, is_optional=.true.)
607 natt = get_num_blocks(diag_file_id, "attributes", parent_block_id=var_id)
608 if (natt .eq. 1) then
609 call get_block_ids(diag_file_id, "attributes", var_att_id, parent_block_id=var_id)
610 nkeys = get_nkeys(diag_file_id, var_att_id(1))
611 allocate(key_ids(nkeys))
612 call get_key_ids(diag_file_id, var_att_id(1), key_ids)
614 allocate(field%var_attributes(nkeys, 2))
616 call get_key_name(diag_file_id, key_ids(j), field%var_attributes(j, 1))
617 call get_key_value(diag_file_id, key_ids(j), field%var_attributes(j, 2))
620 elseif (natt .ne. 0) then
621 call mpp_error(FATAL, "diag_yaml_object_init: variable "//trim(field%var_varname)//" has multiple attribute blocks")
624 !> Set the zbounds if they exist
625 field%var_zbounds = DIAG_NULL
626 call get_value_from_key(diag_file_id, var_id, "zbounds", field%var_zbounds, is_optional=.true.)
627 if (field%has_var_zbounds()) MAX_SUBAXES = MAX_SUBAXES + 1
630 !> @brief diag_manager wrapper to get_value_from_key to use for allocatable
632 subroutine diag_get_value_from_key(diag_file_id, par_id, key_name, value_name, is_optional)
633 integer, intent(in) :: diag_file_id!< Id of the file block in the yaml file
634 integer, intent(in) :: par_id !< Id of the parent block in the yaml file
635 character(len=*), intent(in) :: key_name !< Key to look for in the parent block
636 character(len=:), allocatable :: value_name !< Value of the key
637 logical, intent(in), optional :: is_optional !< Flag indicating if the key is optional
639 character(len=255) :: buffer !< String buffer to read in to
641 buffer = "" !< Needs to be initialized for optional keys that are not present
642 call get_value_from_key(diag_file_id, par_id, trim(key_name), buffer, is_optional= is_optional)
643 allocate(character(len=len_trim(buffer)) :: value_name)
644 value_name = trim(buffer)
646 end subroutine diag_get_value_from_key
648 !> @brief gets the lat/lon of the sub region to use in a diag_table yaml
649 subroutine get_sub_region(diag_yaml_id, sub_region_id, sub_region, grid_type, fname)
650 integer, intent(in) :: diag_yaml_id !< Id of the diag_table yaml file
651 integer, intent(in) :: sub_region_id !< Id of the region block to read from
652 type(subRegion_type),intent(inout) :: sub_region !< Type that stores the sub_region
653 character(len=*), intent(in) :: grid_type !< The grid_type as it is read from the file
654 character(len=*), intent(in) :: fname !< filename of the subregion (for error messages)
656 select case (trim(grid_type))
658 sub_region%grid_type = latlon_gridtype
659 allocate(real(kind=r4_kind) :: sub_region%corners(4,2))
661 sub_region%grid_type = index_gridtype
662 allocate(integer(kind=i4_kind) :: sub_region%corners(4,2))
664 call get_value_from_key(diag_yaml_id, sub_region_id, "tile", sub_region%tile, is_optional=.true.)
665 if (sub_region%tile .eq. DIAG_NULL) call mpp_error(FATAL, &
666 "The tile number is required when defining a "//&
667 "subregion. Check your subregion entry for "//trim(fname))
669 call mpp_error(FATAL, trim(grid_type)//" is not a valid region type. &
670 &The acceptable values are latlon and index. &
671 &Check your entry for file:"//trim(fname))
674 call get_value_from_key(diag_yaml_id, sub_region_id, "corner1", sub_region%corners(1,:))
675 call get_value_from_key(diag_yaml_id, sub_region_id, "corner2", sub_region%corners(2,:))
676 call get_value_from_key(diag_yaml_id, sub_region_id, "corner3", sub_region%corners(3,:))
677 call get_value_from_key(diag_yaml_id, sub_region_id, "corner4", sub_region%corners(4,:))
679 end subroutine get_sub_region
681 !> @brief gets the total number of variables in the diag_table yaml file
682 !! @return total number of variables
683 function get_total_num_vars(diag_yaml_id, diag_file_id) &
686 integer, intent(in) :: diag_yaml_id !< Id for the diag_table yaml
687 integer, intent(in) :: diag_file_id !< Id of the file in the diag_table yaml
688 integer :: total_nvars
690 integer :: i !< For do loop
691 integer :: nvars !< Number of variables in a file
692 integer, allocatable :: var_ids(:) !< Id of the variables in the file block of the yaml file
693 logical :: var_write !< Flag indicating if the user wants the variable to be written
695 nvars = get_num_blocks(diag_yaml_id, "varlist", parent_block_id=diag_file_id)
696 allocate(var_ids(nvars))
697 call get_block_ids(diag_yaml_id, "varlist", var_ids, parent_block_id=diag_file_id)
699 !< Loop through all the variables in the diag_file block and only count those that don't have write_var=false
703 call get_value_from_key(diag_yaml_id, var_ids(i), "write_var", var_write, is_optional=.true.)
704 if (var_write) total_nvars = total_nvars + 1
708 !> @brief This parses the freq, new_file_freq, or file_duration keys which are read in as a comma list
709 subroutine parse_key(filename, buffer, file_freq, file_frequnit, var)
710 character(len=*), intent(in) :: filename !< The name of the file (for error messages)
711 character(len=*), intent(inout) :: buffer !< Buffer that was read in from the yaml
712 integer, intent(out) :: file_freq(:) !< buffer to store the freq, new_file_freq, or
713 !! file_duration after it is parsed
714 integer, intent(out) :: file_frequnit(:) !< buffer to store the freq units, new_file_freq units,
715 !! or file_duration units after it is parsed
716 character(len=*), intent(in) :: var !< Name of the key parsing
718 integer :: j !< location of the ",' in the buffer
719 integer :: k !< location of the " " that seperated the units
720 logical :: finished !< .true. if the parsing is complete
721 integer :: count !< Number of keys that have been parsed
722 character(len=255) :: str !< Member of the comma seperated list
723 character(len=10) :: units !< String to hold the units
724 integer :: err_unit !< Error key
726 if (buffer .eq. "") return
731 do while (.not. finished)
733 buffer = buffer(j+1:len_trim(buffer))
734 j = index(buffer, ",")
736 !< There is only 1 member in the list
737 j = len_trim(buffer)+1
741 str = adjustl(buffer(1:j-1))
744 read(str(1:k-1), *, iostat=err_unit) file_freq(count)
745 units = str(k+1:len_trim(str))
747 if (err_unit .ne. 0) &
748 call mpp_error(FATAL, "Error parsing "//trim(var)//". Check your entry for file"//&
751 if (file_freq(count) .lt. -1) &
752 call mpp_error(FATAL, trim(var)//" is not valid. &
753 &Check your entry for file:"//trim(filename))
755 if (file_freq(count) .eq. -1 .or. file_freq(count) .eq. 0) then
756 !! The file is static so no need to read the units
757 file_frequnit(count) = DIAG_DAYS
759 if (trim(units) .eq. "") &
760 call mpp_error(FATAL, trim(var)//" units is required. &
761 &Check your entry for file:"//trim(filename))
763 file_frequnit(count) = set_valid_time_units(units, &
764 trim(var)//" for file:"//trim(filename))
767 end subroutine parse_key
769 !> @brief This checks if the time unit in a diag file is valid and sets the integer equivalent
770 subroutine set_file_time_units (yaml_fileobj, file_timeunit)
771 type(diagYamlFiles_type), intent(inout) :: yaml_fileobj !< diagYamlFiles_type obj to checK
772 character(len=*), intent(in) :: file_timeunit !< file_timeunit as it is read from the diag_table
774 yaml_fileobj%file_timeunit = set_valid_time_units(file_timeunit, "timeunit for file:"//trim(yaml_fileobj%file_fname))
775 end subroutine set_file_time_units
777 !> @brief This checks if the filename_time in a diag file is correct and sets the integer equivalent
778 subroutine set_filename_time(yaml_fileobj, filename_time)
779 type(diagYamlFiles_type), intent(inout) :: yaml_fileobj !< diagYamlFiles_type obj to check
780 character(len=*), intent(in) :: filename_time !< filename_time as it is read from the yaml
782 select case (trim(filename_time))
784 yaml_fileobj%filename_time = middle_time !< This is the default
786 yaml_fileobj%filename_time = begin_time
788 yaml_fileobj%filename_time = middle_time
790 yaml_fileobj%filename_time = end_time
792 call mpp_error(FATAL, trim(filename_time)//" is an invalid filename_time &
793 &The acceptable values are begin, middle, and end. &
794 &Check your entry for file "//trim(yaml_fileobj%file_fname))
796 end subroutine set_filename_time
798 !> @brief This checks if the kind of a diag field is valid and sets it
799 subroutine set_field_kind(field, skind)
800 type(diagYamlFilesVar_type), intent(inout) :: field !< diagYamlFilesVar_type obj to read the contents into
801 character(len=*), intent(in) :: skind !< The variable kind as read from diag_yaml
803 select case (TRIM(skind))
813 call mpp_error(FATAL, trim(skind)//" is an invalid kind! &
814 &The acceptable values are r4, r8, i4, i8. &
815 &Check your entry for file:"//trim(field%var_varname)//" in file "//trim(field%var_fname))
818 end subroutine set_field_kind
820 !> @brief This checks if the reduction of a diag field is valid and sets it
821 !! If the reduction method is diurnalXX or powXX, it gets the number of diurnal sample and the power value
822 subroutine set_field_reduction(field, reduction_method)
823 type(diagYamlFilesVar_type), intent(inout) :: field !< diagYamlFilesVar_type obj to read the contents into
824 character(len=*) , intent(in) :: reduction_method!< reduction method as read from the yaml
826 integer :: n_diurnal !< number of diurnal samples
827 integer :: pow_value !< The power value
828 integer :: ioerror !< io error status after reading in the diurnal samples
833 if (index(reduction_method, "diurnal") .ne. 0) then
834 READ (reduction_method(8:LEN_TRIM(reduction_method)), FMT=*, IOSTAT=ioerror) n_diurnal
835 if (ioerror .ne. 0) &
836 call mpp_error(FATAL, "Error getting the number of diurnal samples from "//trim(reduction_method))
837 if (n_diurnal .le. 0) &
838 call mpp_error(FATAL, "Diurnal samples should be greater than 0. &
839 & Check your entry for file:"//trim(field%var_varname)//" in file "//trim(field%var_fname))
840 field%var_reduction = time_diurnal
841 elseif (index(reduction_method, "pow") .ne. 0) then
842 READ (reduction_method(4:LEN_TRIM(reduction_method)), FMT=*, IOSTAT=ioerror) pow_value
843 if (ioerror .ne. 0) &
844 call mpp_error(FATAL, "Error getting the power value from "//trim(reduction_method))
845 if (pow_value .le. 0) &
846 call mpp_error(FATAL, "The power value should be greater than 0. &
847 & Check your entry for file:"//trim(field%var_varname)//" in file "//trim(field%var_fname))
848 field%var_reduction = time_power
850 select case (reduction_method)
852 field%var_reduction = time_none
854 field%var_reduction = time_average
856 field%var_reduction = time_min
858 field%var_reduction = time_max
860 field%var_reduction = time_rms
862 field%var_reduction = time_sum
864 call mpp_error(FATAL, trim(reduction_method)//" is an invalid reduction method! &
865 &The acceptable values are none, average, pow##, diurnal##, min, max, and rms. &
866 &Check your entry for file:"//trim(field%var_varname)//" in file "//trim(field%var_fname))
870 field%n_diurnal = n_diurnal
871 field%pow_value = pow_value
872 end subroutine set_field_reduction
874 !> @brief This checks if a time unit is valid and if it is, it assigns the integer equivalent
875 !! @return The integer equivalent to the time units
876 function set_valid_time_units(time_units, error_msg) &
877 result(time_units_int)
879 character(len=*), intent(in) :: time_units !< The time_units as a string
880 character(len=*), intent(in) :: error_msg !< Error message to append
882 integer :: time_units_int !< The integer equivalent of the time_units
884 select case (TRIM(time_units))
886 time_units_int = DIAG_SECONDS
888 time_units_int = DIAG_MINUTES
890 time_units_int = DIAG_HOURS
892 time_units_int = DIAG_DAYS
894 time_units_int = DIAG_MONTHS
896 time_units_int = DIAG_YEARS
898 time_units_int =DIAG_NULL
899 call mpp_error(FATAL, trim(error_msg)//" is not valid. Acceptable values are "&
900 "seconds, minutes, hours, days, months, years")
902 end function set_valid_time_units
904 !!!!!!! YAML FILE INQUIRIES !!!!!!!
905 !> @brief Finds the number of variables in the file_varlist
906 !! @return the size of the diag_files_obj%file_varlist array
907 integer pure function size_file_varlist (this)
908 class (diagYamlFiles_type), intent(in) :: this !< The object being inquiried
909 size_file_varlist = size(this%file_varlist)
910 end function size_file_varlist
912 !> @brief Inquiry for diag_files_obj%file_fname
913 !! @return file_fname of a diag_yaml_file obj
914 pure function get_file_fname (this) &
916 class (diagYamlFiles_type), intent(in) :: this !< The object being inquiried
917 character (len=:), allocatable :: res !< What is returned
918 res = this%file_fname
919 end function get_file_fname
920 !> @brief Inquiry for diag_files_obj%file_frequnit
921 !! @return file_frequnit of a diag_yaml_file_obj
922 pure function get_file_frequnit (this) &
924 class (diagYamlFiles_type), intent(in) :: this !< The object being inquiried
925 integer :: res !< What is returned
926 res = this%file_frequnit(this%current_new_file_freq_index)
927 end function get_file_frequnit
928 !> @brief Inquiry for diag_files_obj%file_freq
929 !! @return file_freq of a diag_yaml_file_obj
930 pure function get_file_freq(this) &
932 class (diagYamlFiles_type), intent(in) :: this !< The object being inquiried
933 integer :: res !< What is returned
934 res = this%file_freq(this%current_new_file_freq_index)
935 end function get_file_freq
936 !> @brief Inquiry for diag_files_obj%file_timeunit
937 !! @return file_timeunit of a diag_yaml_file_obj
938 pure function get_file_timeunit (this) &
940 class (diagYamlFiles_type), intent(in) :: this !< The object being inquiried
941 integer :: res !< What is returned
942 res = this%file_timeunit
943 end function get_file_timeunit
944 !> @brief Inquiry for diag_files_obj%file_unlimdim
945 !! @return file_unlimdim of a diag_yaml_file_obj
946 pure function get_file_unlimdim(this) &
948 class (diagYamlFiles_type), intent(in) :: this !< The object being inquiried
949 character (len=:), allocatable :: res !< What is returned
950 res = this%file_unlimdim
951 end function get_file_unlimdim
952 !> @brief Inquiry for diag_files_obj%file_subregion
953 !! @return file_sub_region of a diag_yaml_file_obj
954 function get_file_sub_region (this) &
956 class (diagYamlFiles_type), intent(in) :: this !< The object being inquiried
957 type(subRegion_type) :: res !< What is returned
958 res = this%file_sub_region
959 end function get_file_sub_region
960 !> @brief Inquiry for diag_files_obj%file_new_file_freq
961 !! @return file_new_file_freq of a diag_yaml_file_obj
962 pure function get_file_new_file_freq(this) &
964 class (diagYamlFiles_type), intent(in) :: this !< The object being inquiried
965 integer :: res !< What is returned
966 res = this%file_new_file_freq(this%current_new_file_freq_index)
967 end function get_file_new_file_freq
968 !> @brief Inquiry for diag_files_obj%file_new_file_freq_units
969 !! @return file_new_file_freq_units of a diag_yaml_file_obj
970 pure function get_file_new_file_freq_units (this) &
972 class (diagYamlFiles_type), intent(in) :: this !< The object being inquiried
973 integer :: res !< What is returned
974 res = this%file_new_file_freq_units(this%current_new_file_freq_index)
975 end function get_file_new_file_freq_units
976 !> @brief Inquiry for diag_files_obj%file_start_time
977 !! @return file_start_time of a diag_yaml_file_obj
978 pure function get_file_start_time (this) &
980 class (diagYamlFiles_type), intent(in) :: this !< The object being inquiried
981 character (len=:), allocatable :: res !< What is returned
982 res = this%file_start_time
983 end function get_file_start_time
984 !> @brief Inquiry for diag_files_obj%file_duration
985 !! @return file_duration of a diag_yaml_file_obj
986 pure function get_file_duration (this) &
988 class (diagYamlFiles_type), intent(in) :: this !< The object being inquiried
989 integer :: res !< What is returned
990 res = this%file_duration(this%current_new_file_freq_index)
991 end function get_file_duration
992 !> @brief Inquiry for diag_files_obj%file_duration_units
993 !! @return file_duration_units of a diag_yaml_file_obj
994 pure function get_file_duration_units (this) &
996 class (diagYamlFiles_type), intent(in) :: this !< The object being inquiried
997 integer :: res !< What is returned
998 res = this%file_duration_units(this%current_new_file_freq_index)
999 end function get_file_duration_units
1000 !> @brief Inquiry for diag_files_obj%file_varlist
1001 !! @return file_varlist of a diag_yaml_file_obj
1002 pure function get_file_varlist (this) &
1004 class (diagYamlFiles_type), intent(in) :: this !< The object being inquiried
1005 character (:), allocatable :: res(:) !< What is returned
1006 res = this%file_varlist
1007 end function get_file_varlist
1008 !> @brief Inquiry for diag_files_obj%file_global_meta
1009 !! @return file_global_meta of a diag_yaml_file_obj
1010 pure function get_file_global_meta (this) &
1012 class (diagYamlFiles_type), intent(in) :: this !< The object being inquiried
1013 character (len=MAX_STR_LEN), allocatable :: res(:,:) !< What is returned
1014 res = this%file_global_meta
1015 end function get_file_global_meta
1016 !> @brief Get the integer equivalent of the time to use to determine the filename,
1017 !! if using a wildcard file name (i.e ocn%4yr%2mo%2dy%2hr)
1018 !! @return the integer equivalent of the time to use to determine the filename
1019 pure function get_filename_time(this) &
1021 class (diagYamlFiles_type), intent(in) :: this !< The object being inquiried
1022 integer :: res !< What is returned
1023 res = this%filename_time
1025 !> @brief Inquiry for whether file_global_meta is allocated
1026 !! @return Flag indicating if file_global_meta is allocated
1027 function is_global_meta(this) &
1029 class (diagYamlFiles_type), intent(in) :: this !< The object being inquiried
1032 if (allocated(this%file_global_meta)) &
1036 !> @brief Increate the current_new_file_freq_index by 1
1037 subroutine increase_new_file_freq_index(this)
1038 class(diagYamlFiles_type), intent(inout) :: this !< The file object
1039 this%current_new_file_freq_index = this%current_new_file_freq_index + 1
1041 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1042 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1043 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1044 !!!!!!! VARIABLES ROUTINES AND FUNCTIONS !!!!!!!
1046 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1047 !!!!!!! YAML VAR INQUIRIES !!!!!!!
1048 !> @brief Inquiry for diag_yaml_files_var_obj%var_fname
1049 !! @return var_fname of a diag_yaml_files_var_obj
1050 pure function get_var_fname (this) &
1052 class (diagYamlFilesVar_type), intent(in) :: this !< The object being inquiried
1053 character (len=:), allocatable :: res !< What is returned
1054 res = this%var_fname
1055 end function get_var_fname
1056 !> @brief Inquiry for diag_yaml_files_var_obj%var_varname
1057 !! @return var_varname of a diag_yaml_files_var_obj
1058 pure function get_var_varname (this) &
1060 class (diagYamlFilesVar_type), intent(in) :: this !< The object being inquiried
1061 character (len=:), allocatable :: res !< What is returned
1062 res = this%var_varname
1063 end function get_var_varname
1064 !> @brief Inquiry for diag_yaml_files_var_obj%var_reduction
1065 !! @return var_reduction of a diag_yaml_files_var_obj
1066 pure function get_var_reduction (this) &
1068 class (diagYamlFilesVar_type), intent(in) :: this !< The object being inquiried
1069 integer, allocatable :: res !< What is returned
1070 res = this%var_reduction
1071 end function get_var_reduction
1072 !> @brief Inquiry for diag_yaml_files_var_obj%var_module
1073 !! @return var_module of a diag_yaml_files_var_obj
1074 pure function get_var_module (this) &
1076 class (diagYamlFilesVar_type), intent(in) :: this !< The object being inquiried
1077 character (len=:), allocatable :: res !< What is returned
1078 res = this%var_module
1079 end function get_var_module
1080 !> @brief Inquiry for diag_yaml_files_var_obj%var_kind
1081 !! @return var_kind of a diag_yaml_files_var_obj
1082 pure function get_var_kind (this) &
1084 class (diagYamlFilesVar_type), intent(in) :: this !< The object being inquiried
1085 integer, allocatable :: res !< What is returned
1087 end function get_var_kind
1088 !> @brief Inquiry for diag_yaml_files_var_obj%var_outname
1089 !! @return var_outname of a diag_yaml_files_var_obj
1090 pure function get_var_outname (this) &
1092 class (diagYamlFilesVar_type), intent(in) :: this !< The object being inquiried
1093 character (len=:), allocatable :: res !< What is returned
1095 if (this%has_var_outname()) then
1096 res = this%var_outname
1098 res = this%var_varname !< If outname is not set, the variable name will be used
1100 end function get_var_outname
1101 !> @brief Inquiry for diag_yaml_files_var_obj%var_longname
1102 !! @return var_longname of a diag_yaml_files_var_obj
1103 pure function get_var_longname (this) &
1105 class (diagYamlFilesVar_type), intent(in) :: this !< The object being inquiried
1106 character (len=:), allocatable :: res !< What is returned
1107 res = this%var_longname
1108 end function get_var_longname
1109 !> @brief Inquiry for diag_yaml_files_var_obj%var_units
1110 !! @return var_units of a diag_yaml_files_var_obj
1111 pure function get_var_units (this) &
1113 class (diagYamlFilesVar_type), intent(in) :: this !< The object being inquiried
1114 character (len=:), allocatable :: res !< What is returned
1115 res = this%var_units
1116 end function get_var_units
1117 !> @brief Inquiry for diag_yaml_files_var_obj%var_zbounds
1118 !! @return var_zbounds of a diag_yaml_files_var_obj
1119 pure function get_var_zbounds (this) &
1121 class (diagYamlFilesVar_type), intent(in) :: this !< The object being inquiried
1122 real(kind=r4_kind) :: res(2) !< What is returned
1123 res = this%var_zbounds
1124 end function get_var_zbounds
1125 !> @brief Inquiry for diag_yaml_files_var_obj%var_attributes
1126 !! @return var_attributes of a diag_yaml_files_var_obj
1127 pure function get_var_attributes(this) &
1129 class (diagYamlFilesVar_type), intent(in) :: this !< The object being inquiried
1130 character (len=MAX_STR_LEN), allocatable :: res (:,:) !< What is returned
1131 res = this%var_attributes
1132 end function get_var_attributes
1133 !> @brief Inquiry for diag_yaml_files_var_obj%n_diurnal
1134 !! @return the number of diurnal samples of a diag_yaml_files_var_obj
1135 pure function get_n_diurnal(this) &
1137 class (diagYamlFilesVar_type), intent(in) :: this !< The object being inquiried
1138 integer :: res !< What is returned
1139 res = this%n_diurnal
1140 end function get_n_diurnal
1141 !> @brief Inquiry for diag_yaml_files_var_obj%pow_value
1142 !! @return the pow_value of a diag_yaml_files_var_obj
1143 pure function get_pow_value(this) &
1145 class (diagYamlFilesVar_type), intent(in) :: this !< The object being inquiried
1146 integer :: res !< What is returned
1147 res = this%pow_value
1148 end function get_pow_value
1149 !> @brief Inquiry for whether var_attributes is allocated
1150 !! @return Flag indicating if var_attributes is allocated
1151 function is_var_attributes(this) &
1153 class (diagYamlFilesVar_type), intent(in) :: this !< The object being inquiried
1156 if (allocated(this%var_attributes)) &
1158 end function is_var_attributes
1160 !> @brief Initializes the non string values of a diagYamlFiles_type to its
1162 subroutine diag_yaml_files_obj_init(obj)
1163 type(diagYamlFiles_type), intent(out) :: obj !< diagYamlFiles_type object to initialize
1165 obj%file_freq = DIAG_NULL
1166 obj%file_sub_region%tile = DIAG_NULL
1167 obj%file_new_file_freq = DIAG_NULL
1168 obj%file_duration = DIAG_NULL
1169 obj%file_new_file_freq_units = DIAG_NULL
1170 obj%file_duration_units = DIAG_NULL
1171 obj%current_new_file_freq_index = 1
1172 end subroutine diag_yaml_files_obj_init
1174 !> @brief Checks if diag_file_obj%file_fname is allocated
1175 !! @return true if diag_file_obj%file_fname is allocated
1176 pure logical function has_file_fname (this)
1177 class(diagYamlFiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1178 has_file_fname = allocated(this%file_fname)
1179 end function has_file_fname
1180 !> @brief Checks if diag_file_obj%file_frequnit is allocated
1181 !! @return true if diag_file_obj%file_frequnit is allocated
1182 pure logical function has_file_frequnit (this)
1183 class(diagYamlFiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1184 has_file_frequnit = this%file_frequnit(this%current_new_file_freq_index) .NE. DIAG_NULL
1185 end function has_file_frequnit
1186 !> @brief diag_file_obj%file_freq is on the stack, so the object always has it
1187 !! @return true if diag_file_obj%file_freq is allocated
1188 pure logical function has_file_freq (this)
1189 class(diagYamlFiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1190 has_file_freq = .true.
1191 end function has_file_freq
1192 !> @brief Checks if diag_file_obj%file_timeunit is allocated
1193 !! @return true if diag_file_obj%file_timeunit is allocated
1194 pure logical function has_file_timeunit (this)
1195 class(diagYamlFiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1196 has_file_timeunit = this%file_timeunit .ne. diag_null
1197 end function has_file_timeunit
1198 !> @brief Checks if diag_file_obj%file_unlimdim is allocated
1199 !! @return true if diag_file_obj%file_unlimdim is allocated
1200 pure logical function has_file_unlimdim (this)
1201 class(diagYamlFiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1202 has_file_unlimdim = allocated(this%file_unlimdim)
1203 end function has_file_unlimdim
1204 !> @brief Checks if diag_file_obj%file_write is on the stack, so this will always be true
1206 pure logical function has_file_write (this)
1207 class(diagYamlFiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1208 has_file_write = .true.
1209 end function has_file_write
1210 !> @brief Checks if diag_file_obj%file_sub_region is being used and has the sub region variables allocated
1211 !! @return true if diag_file_obj%file_sub_region sub region variables are allocated
1212 pure logical function has_file_sub_region (this)
1213 class(diagYamlFiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1214 if ( this%file_sub_region%grid_type .eq. latlon_gridtype .or. this%file_sub_region%grid_type .eq. index_gridtype) then
1215 has_file_sub_region = .true.
1217 has_file_sub_region = .false.
1219 end function has_file_sub_region
1220 !> @brief diag_file_obj%file_new_file_freq is defined on the stack, so this will return true
1222 pure logical function has_file_new_file_freq (this)
1223 class(diagYamlFiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1224 has_file_new_file_freq = this%file_new_file_freq(this%current_new_file_freq_index) .ne. DIAG_NULL
1225 end function has_file_new_file_freq
1226 !> @brief Checks if diag_file_obj%file_new_file_freq_units is allocated
1227 !! @return true if diag_file_obj%file_new_file_freq_units is allocated
1228 pure logical function has_file_new_file_freq_units (this)
1229 class(diagYamlFiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1230 has_file_new_file_freq_units = this%file_new_file_freq_units(this%current_new_file_freq_index) .ne. diag_null
1231 end function has_file_new_file_freq_units
1232 !> @brief Checks if diag_file_obj%file_start_time is allocated
1233 !! @return true if diag_file_obj%file_start_time is allocated
1234 pure logical function has_file_start_time (this)
1235 class(diagYamlFiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1236 has_file_start_time = allocated(this%file_start_time)
1237 end function has_file_start_time
1238 !> @brief diag_file_obj%file_duration is allocated on th stack, so this is always true
1240 pure logical function has_file_duration (this)
1241 class(diagYamlFiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1242 has_file_duration = this%file_duration(this%current_new_file_freq_index) .ne. DIAG_NULL
1243 end function has_file_duration
1244 !> @brief diag_file_obj%file_duration_units is on the stack, so this will retrun true
1246 pure logical function has_file_duration_units (this)
1247 class(diagYamlFiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1248 has_file_duration_units = this%file_duration_units(this%current_new_file_freq_index) .ne. diag_null
1249 end function has_file_duration_units
1250 !> @brief Checks if diag_file_obj%file_varlist is allocated
1251 !! @return true if diag_file_obj%file_varlist is allocated
1252 pure logical function has_file_varlist (this)
1253 class(diagYamlFiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1254 has_file_varlist = allocated(this%file_varlist)
1255 end function has_file_varlist
1256 !> @brief Checks if diag_file_obj%file_global_meta is allocated
1257 !! @return true if diag_file_obj%file_global_meta is allocated
1258 pure logical function has_file_global_meta (this)
1259 class(diagYamlFiles_type), intent(in) :: this !< diagYamlFiles_type object to initialize
1260 has_file_global_meta = allocated(this%file_global_meta)
1261 end function has_file_global_meta
1263 !> @brief Checks if diag_file_obj%var_fname is allocated
1264 !! @return true if diag_file_obj%var_fname is allocated
1265 pure logical function has_var_fname (this)
1266 class(diagYamlFilesVar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1267 has_var_fname = allocated(this%var_fname)
1268 end function has_var_fname
1269 !> @brief Checks if diag_file_obj%var_varname is allocated
1270 !! @return true if diag_file_obj%var_varname is allocated
1271 pure logical function has_var_varname (this)
1272 class(diagYamlFilesVar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1273 has_var_varname = allocated(this%var_varname)
1274 end function has_var_varname
1275 !> @brief Checks if diag_file_obj%var_reduction is allocated
1276 !! @return true if diag_file_obj%var_reduction is allocated
1277 pure logical function has_var_reduction (this)
1278 class(diagYamlFilesVar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1279 has_var_reduction = allocated(this%var_reduction)
1280 end function has_var_reduction
1281 !> @brief Checks if diag_file_obj%var_module is allocated
1282 !! @return true if diag_file_obj%var_module is allocated
1283 pure logical function has_var_module (this)
1284 class(diagYamlFilesVar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1285 has_var_module = allocated(this%var_module)
1286 end function has_var_module
1287 !> @brief Checks if diag_file_obj%var_kind is allocated
1288 !! @return true if diag_file_obj%var_kind is allocated
1289 pure logical function has_var_kind (this)
1290 class(diagYamlFilesVar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1291 has_var_kind = allocated(this%var_kind)
1292 end function has_var_kind
1293 !> @brief diag_file_obj%var_write is on the stack, so this returns true
1295 pure logical function has_var_write (this)
1296 class(diagYamlFilesVar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1297 has_var_write = .true.
1298 end function has_var_write
1299 !> @brief Checks if diag_file_obj%var_outname is allocated
1300 !! @return true if diag_file_obj%var_outname is allocated
1301 pure logical function has_var_outname (this)
1302 class(diagYamlFilesVar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1303 if (allocated(this%var_outname)) then
1304 if (trim(this%var_outname) .ne. "") then
1305 has_var_outname = .true.
1307 has_var_outname = .false.
1310 has_var_outname = .true.
1312 end function has_var_outname
1313 !> @brief Checks if diag_file_obj%var_longname is allocated
1314 !! @return true if diag_file_obj%var_longname is allocated
1315 pure logical function has_var_longname (this)
1316 class(diagYamlFilesVar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1317 has_var_longname = allocated(this%var_longname)
1318 end function has_var_longname
1319 !> @brief Checks if diag_file_obj%var_units is allocated
1320 !! @return true if diag_file_obj%var_units is allocated
1321 pure logical function has_var_units (this)
1322 class(diagYamlFilesVar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1323 has_var_units = allocated(this%var_units)
1324 end function has_var_units
1325 !> @brief Checks if diag_file_obj%var_zbounds is allocated
1326 !! @return true if diag_file_obj%var_zbounds is allocated
1327 pure logical function has_var_zbounds (this)
1328 class(diagYamlFilesVar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1329 has_var_zbounds = any(this%var_zbounds .ne. diag_null)
1330 end function has_var_zbounds
1331 !> @brief Checks if diag_file_obj%var_attributes is allocated
1332 !! @return true if diag_file_obj%var_attributes is allocated
1333 pure logical function has_var_attributes (this)
1334 class(diagYamlFilesVar_type), intent(in) :: this !< diagYamlvar_type object to initialize
1335 has_var_attributes = allocated(this%var_attributes)
1336 end function has_var_attributes
1337 !> @brief Checks if diag_file_obj%n_diurnal is set
1338 !! @return true if diag_file_obj%n_diurnal is set
1339 pure logical function has_n_diurnal(this)
1340 class(diagYamlFilesVar_type), intent(in) :: this !< diagYamlvar_type object to inquire
1341 has_n_diurnal = (this%n_diurnal .ne. 0)
1342 end function has_n_diurnal
1343 !> @brief Checks if diag_file_obj%pow_value is set
1344 !! @return true if diag_file_obj%pow_value is set
1345 pure logical function has_pow_value(this)
1346 class(diagYamlFilesVar_type), intent(in) :: this !< diagYamlvar_type object to inquire
1347 has_pow_value = (this%pow_value .ne. 0)
1348 end function has_pow_value
1350 !> @brief Checks if diag_file_obj%diag_title is allocated
1351 !! @return true if diag_file_obj%diag_title is allocated
1352 pure logical function has_diag_title (this)
1353 class(diagYamlObject_type), intent(in) :: this !< diagYamlObject_type object to inquire
1354 has_diag_title = allocated(this%diag_title)
1355 end function has_diag_title
1356 !> @brief diag_file_obj%diag_basedate is on the stack, so this is always true
1358 pure logical function has_diag_basedate (this)
1359 class(diagYamlObject_type), intent(in) :: this !< diagYamlObject_type object to initialize
1360 has_diag_basedate = .true.
1361 end function has_diag_basedate
1362 !> @brief Checks if diag_file_obj%diag_files is allocated
1363 !! @return true if diag_file_obj%diag_files is allocated
1364 pure logical function has_diag_files (this)
1365 class(diagYamlObject_type), intent(in) :: this !< diagYamlObject_type object to initialize
1366 has_diag_files = allocated(this%diag_files)
1367 end function has_diag_files
1368 !> @brief Checks if diag_file_obj%diag_fields is allocated
1369 !! @return true if diag_file_obj%diag_fields is allocated
1370 pure logical function has_diag_fields (this)
1371 class(diagYamlObject_type), intent(in) :: this !< diagYamlObject_type object to initialize
1372 has_diag_fields = allocated(this%diag_fields)
1373 end function has_diag_fields
1375 !> @brief Determine the number of unique diag_fields in the diag_yaml_object
1376 !! @return The number of unique diag_fields
1377 function get_num_unique_fields() &
1380 nfields = fms_find_unique(variable_list%var_pointer, size(variable_list%var_pointer))
1382 end function get_num_unique_fields
1384 !> @brief Determines if a diag_field is in the diag_yaml_object
1385 !! @return Indices of the locations where the field was found
1386 function find_diag_field(diag_field_name, module_name) &
1389 character(len=*), intent(in) :: diag_field_name !< diag_field name to search for
1390 character(len=*), intent(in) :: module_name !< Name of the module, the variable is in
1392 integer, allocatable :: indices(:)
1394 indices = fms_find_my_string(variable_list%var_pointer, size(variable_list%var_pointer), &
1395 & lowercase(trim(diag_field_name))//":"//lowercase(trim(module_name)//c_null_char))
1396 end function find_diag_field
1398 !> @brief Gets the diag_field entries corresponding to the indices of the sorted variable_list
1399 !! @return Array of diag_fields
1400 function get_diag_fields_entries(indices) &
1403 integer, intent(in) :: indices(:) !< Indices of the field in the sorted variable_list array
1404 type(diagYamlFilesVar_type), dimension (:), allocatable :: diag_field
1406 integer :: i !< For do loops
1407 integer :: field_id !< Indices of the field in the diag_yaml array
1409 allocate(diag_field(size(indices)))
1411 do i = 1, size(indices)
1412 field_id = variable_list%diag_field_indices(indices(i))
1413 diag_field(i) = diag_yaml%diag_fields(field_id)
1416 end function get_diag_fields_entries
1418 !> @brief Gets field indices corresponding to the indices (input argument) in the sorted variable_list
1419 !! @return Copy of array of field indices
1420 function get_diag_field_ids(indices) result(field_ids)
1422 integer, intent(in) :: indices(:) !< Indices of the fields in the sorted variable_list array
1423 integer, allocatable :: field_ids(:)
1424 integer :: i !< For do loop
1426 allocate(field_ids(size(indices)))
1428 do i = 1, size(indices)
1429 field_ids(i) = variable_list%diag_field_indices(indices(i))
1432 end function get_diag_field_ids
1434 !> @brief Finds the indices of the diag_yaml%diag_files(:) corresponding to fields in variable_list(indices)
1435 !! @return indices of the diag_yaml%diag_files(:)
1436 function get_diag_files_id(indices) &
1439 integer, intent(in) :: indices(:) !< Indices of the field in the sorted variable_list
1440 integer, allocatable :: file_id(:)
1442 integer :: field_id !< Indices of the field in the diag_yaml field array
1443 integer :: i !< For do loops
1444 character(len=120) :: filename !< Filename of the field
1445 integer, allocatable :: file_indices(:) !< Indices of the file in the sorted variable_list
1447 allocate(file_id(size(indices)))
1449 do i = 1, size(indices)
1450 field_id = variable_list%diag_field_indices(indices(i))
1451 !< Get the filename of the field
1452 filename = diag_yaml%diag_fields(field_id)%var_fname
1454 !< File indice of that file in the array of list of sorted files
1455 file_indices = fms_find_my_string(file_list%file_pointer, size(file_list%file_pointer), &
1456 & trim(filename)//c_null_char)
1458 if (size(file_indices) .ne. 1) &
1459 & call mpp_error(FATAL, "get_diag_files_id: Error getting the correct number of file indices!"//&
1460 " The diag file "//trim(filename)//" was defined "//string(size(file_indices))&
1463 if (file_indices(1) .eq. diag_null) &
1464 & call mpp_error(FATAL, "get_diag_files_id: Error finding the file "//trim(filename)//" in the diag_files yaml")
1466 !< Get the index of the file in the diag_yaml file
1467 file_id(i) = file_list%diag_file_indices(file_indices(1))
1470 end function get_diag_files_id
1472 !> Prints out values from diag_yaml object for debugging.
1473 !! Only writes on root.
1474 subroutine dump_diag_yaml_obj( filename )
1475 character(len=*), optional, intent(in) :: filename !< optional name of logfile to write to, otherwise
1477 type(diagyamlfilesvar_type), allocatable :: fields(:)
1478 type(diagyamlfiles_type), allocatable :: files(:)
1479 integer :: i, unit_num
1480 if( present(filename)) then
1481 open(newunit=unit_num, file=trim(filename), action='WRITE')
1485 !! TODO write to log
1486 if( mpp_pe() .eq. mpp_root_pe()) then
1487 write(unit_num, *) '**********Dumping diag_yaml object**********'
1488 if( diag_yaml%has_diag_title()) write(unit_num, *) 'Title:', diag_yaml%diag_title
1489 if( diag_yaml%has_diag_basedate()) write(unit_num, *) 'basedate array:', diag_yaml%diag_basedate
1490 write(unit_num, *) 'FILES'
1491 allocate(fields(SIZE(diag_yaml%get_diag_fields())))
1492 allocate(files(SIZE(diag_yaml%get_diag_files())))
1493 files = diag_yaml%get_diag_files()
1494 fields = diag_yaml%get_diag_fields()
1496 write(unit_num, *) 'File: ', files(i)%get_file_fname()
1497 if(files(i)%has_file_frequnit()) write(unit_num, *) 'file_frequnit:', files(i)%get_file_frequnit()
1498 if(files(i)%has_file_freq()) write(unit_num, *) 'freq:', files(i)%get_file_freq()
1499 if(files(i)%has_file_timeunit()) write(unit_num, *) 'timeunit:', files(i)%get_file_timeunit()
1500 if(files(i)%has_file_unlimdim()) write(unit_num, *) 'unlimdim:', files(i)%get_file_unlimdim()
1501 !if(files(i)%has_file_sub_region()) write(unit_num, *) 'sub_region:', files(i)%get_file_sub_region()
1502 if(files(i)%has_file_new_file_freq()) write(unit_num, *) 'new_file_freq:', files(i)%get_file_new_file_freq()
1503 if(files(i)%has_file_new_file_freq_units()) write(unit_num, *) 'new_file_freq_units:', &
1504 & files(i)%get_file_new_file_freq_units()
1505 if(files(i)%has_file_start_time()) write(unit_num, *) 'start_time:', files(i)%get_file_start_time()
1506 if(files(i)%has_file_duration()) write(unit_num, *) 'duration:', files(i)%get_file_duration()
1507 if(files(i)%has_file_duration_units()) write(unit_num, *) 'duration_units:', files(i)%get_file_duration_units()
1508 if(files(i)%has_file_varlist()) write(unit_num, *) 'varlist:', files(i)%get_file_varlist()
1509 if(files(i)%has_file_global_meta()) write(unit_num, *) 'global_meta:', files(i)%get_file_global_meta()
1510 if(files(i)%is_global_meta()) write(unit_num, *) 'global_meta:', files(i)%is_global_meta()
1511 write(unit_num, *) ''
1513 write(unit_num, *) 'FIELDS'
1514 do i=1, SIZE(fields)
1515 write(unit_num, *) 'Field: ', fields(i)%get_var_fname()
1516 if(fields(i)%has_var_fname()) write(unit_num, *) 'fname:', fields(i)%get_var_fname()
1517 if(fields(i)%has_var_varname()) write(unit_num, *) 'varname:', fields(i)%get_var_varname()
1518 if(fields(i)%has_var_reduction()) write(unit_num, *) 'reduction:', fields(i)%get_var_reduction()
1519 if(fields(i)%has_var_module()) write(unit_num, *) 'module:', fields(i)%get_var_module()
1520 if(fields(i)%has_var_kind()) write(unit_num, *) 'kind:', fields(i)%get_var_kind()
1521 if(fields(i)%has_var_outname()) write(unit_num, *) 'outname:', fields(i)%get_var_outname()
1522 if(fields(i)%has_var_longname()) write(unit_num, *) 'longname:', fields(i)%get_var_longname()
1523 if(fields(i)%has_var_units()) write(unit_num, *) 'units:', fields(i)%get_var_units()
1524 if(fields(i)%has_var_zbounds()) write(unit_num, *) 'zbounds:', fields(i)%get_var_zbounds()
1525 if(fields(i)%has_var_attributes()) write(unit_num, *) 'attributes:', fields(i)%get_var_attributes()
1526 if(fields(i)%has_n_diurnal()) write(unit_num, *) 'n_diurnal:', fields(i)%get_n_diurnal()
1527 if(fields(i)%has_pow_value()) write(unit_num, *) 'pow_value:', fields(i)%get_pow_value()
1528 if(fields(i)%has_var_attributes()) write(unit_num, *) 'is_var_attributes:', fields(i)%is_var_attributes()
1530 deallocate(files, fields)
1531 if( present(filename)) then
1538 end module fms_diag_yaml_mod
1540 ! close documentation grouping