FMS 2025.01.02-dev
Flexible Modeling System
Loading...
Searching...
No Matches
fms_diag_object.F90
1!***********************************************************************
2!* GNU Lesser General Public License
3!*
4!* This file is part of the GFDL Flexible Modeling System (FMS).
5!*
6!* FMS is free software: you can redistribute it and/or modify it under
7!* the terms of the GNU Lesser General Public License as published by
8!* the Free Software Foundation, either version 3 of the License, or (at
9!* your option) any later version.
10!*
11!* FMS is distributed in the hope that it will be useful, but WITHOUT
12!* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13!* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14!* for more details.
15!*
16!* You should have received a copy of the GNU Lesser General Public
17!* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18!***********************************************************************
19module fms_diag_object_mod
20use mpp_mod, only: fatal, note, warning, mpp_error, mpp_pe, mpp_root_pe, stdout
21use diag_data_mod, only: diag_null, diag_not_found, diag_not_registered, diag_registered_id, &
22 &diag_field_not_found, diag_not_registered, max_axes, two_d_domain, &
23 &get_base_time, null_axis_id, get_var_type, diag_not_registered, &
26
27 USE time_manager_mod, ONLY: set_time, set_date, get_time, time_type, OPERATOR(>=), OPERATOR(>),&
28 & OPERATOR(<), OPERATOR(==), OPERATOR(/=), OPERATOR(/), OPERATOR(+), ASSIGNMENT(=), get_date, &
30#ifdef use_yaml
31use fms_diag_file_object_mod, only: fmsdiagfilecontainer_type, fmsdiagfile_type, fms_diag_files_object_init
32use fms_diag_field_object_mod, only: fmsdiagfield_type, fms_diag_fields_object_init, get_default_missing_value, &
33 check_for_slices
34use fms_diag_yaml_mod, only: diag_yaml_object_init, diag_yaml_object_end, find_diag_field, &
36use fms_diag_axis_object_mod, only: fms_diag_axis_object_init, fmsdiagaxis_type, fmsdiagsubaxis_type, &
41use fms_mod, only: fms_error_handler
42use fms_diag_reduction_methods_mod, only: check_indices_order, init_mask, set_weight
43use constants_mod, only: seconds_per_day
44#endif
45USE fms_diag_bbox_mod, ONLY: fmsdiagibounds_type, determine_if_block_is_in_region
46#if defined(_OPENMP)
47use omp_lib
48#endif
49use mpp_domains_mod, only: domain1d, domain2d, domainug, null_domain2d
50use fms_string_utils_mod, only: string
51use platform_mod
52implicit none
53private
54
56!TODO add container arrays
57#ifdef use_yaml
58private
59!TODO: Remove FMS prefix from variables in this type
60 class(fmsdiagfilecontainer_type), allocatable :: fms_diag_files (:) !< array of diag files
61 type(fmsdiagfield_type), allocatable :: fms_diag_fields(:) !< Array of diag fields
62 type(fmsdiagoutputbuffer_type), allocatable :: fms_diag_output_buffers(:) !< array of output buffer objects
63 !! one for each variable in the diag_table.yaml
64 logical, private :: data_was_send !< True if send_data has been successfully called for at least one variable
65 !< diag_send_complete does nothing if it is .false.
66 integer, private :: registered_buffers = 0 !< number of registered buffers, per dimension
67 class(fmsdiagaxiscontainer_type), allocatable :: diag_axis(:) !< Array of diag_axis
68 integer, private :: registered_variables !< Number of registered variables
69 integer, private :: registered_axis !< Number of registered axis
70 logical, private :: initialized=.false. !< True if the fmsDiagObject is initialized
71 logical, private :: files_initialized=.false. !< True if the fmsDiagObject is initialized
72 logical, private :: fields_initialized=.false. !< True if the fmsDiagObject is initialized
73 logical, private :: buffers_initialized=.false. !< True if the fmsDiagObject is initialized
74 logical, private :: axes_initialized=.false. !< True if the fmsDiagObject is initialized
75 type(time_type) :: model_end_time !< The time that the simulation is going to end
76 !! (set by calling diag_manager_set_time_end)
77#endif
78 contains
79 procedure :: init => fms_diag_object_init
80 procedure :: diag_end => fms_diag_object_end
81 procedure :: fms_register_diag_field_scalar
82 procedure :: fms_register_diag_field_array
83 procedure :: fms_register_static_field
84 procedure :: fms_diag_axis_init
85 procedure :: register => fms_register_diag_field_obj !! Merely initialize fields.
86 procedure :: fms_diag_field_add_attribute
87 procedure :: fms_diag_axis_add_attribute
88 procedure :: fms_get_domain2d
89 procedure :: fms_get_axis_length
90 procedure :: fms_get_diag_field_id_from_name
91 procedure :: fms_get_field_name_from_id
92 procedure :: fms_get_axis_name_from_id
93 procedure :: fms_diag_accept_data
94 procedure :: fms_diag_send_complete
95 procedure :: do_buffer_math
96 procedure :: fms_diag_do_io
97 procedure :: fms_diag_do_reduction
98 procedure :: fms_diag_field_add_cell_measures
99 procedure :: allocate_diag_field_output_buffers
100 procedure :: fms_diag_compare_window
101 procedure :: set_time_end
102 procedure :: write_diag_manifest
103#ifdef use_yaml
104 procedure :: get_diag_buffer
105#endif
106end type fmsdiagobject_type
107
108type (fmsdiagobject_type), target :: fms_diag_object
109
110public :: fms_register_diag_field_obj
111public :: fms_register_diag_field_scalar
112public :: fms_register_diag_field_array
113public :: fms_register_static_field
114public :: fms_diag_field_add_attribute
115public :: fms_get_diag_field_id_from_name
116public :: fms_diag_object
117public :: fmsdiagobject_type
118integer, private :: registered_variables !< Number of registered variables
119public :: dump_diag_obj
120
121contains
122
123!> @brief Initiliazes the fms_diag_object.
124!! Reads the diag_table.yaml and fills in the yaml object
125!! Allocates the diag manager object arrays for files, fields, and buffers
126!! Initializes variables
127subroutine fms_diag_object_init (this,diag_subset_output, time_init)
128 class(fmsdiagobject_type) :: this !< Diag mediator/controller object
129 integer :: diag_subset_output !< Subset of the diag output?
130 INTEGER, DIMENSION(6), OPTIONAL, INTENT(IN) :: time_init !< Model time diag_manager initialized
131
132#ifdef use_yaml
133 if (this%initialized) return
134
135! allocate(diag_objs(get_num_unique_fields()))
136 CALL diag_yaml_object_init(diag_subset_output)
137
138 !! Doing this here, because the base_time is not set until the yaml is parsed
139 !! if time_init is present, it will be set in diag_manager_init
140 if (.not. present(time_init)) then
142 endif
143
144 this%axes_initialized = fms_diag_axis_object_init(this%diag_axis)
145 this%files_initialized = fms_diag_files_object_init(this%FMS_diag_files)
146 this%fields_initialized = fms_diag_fields_object_init(this%FMS_diag_fields)
147 this%buffers_initialized =fms_diag_output_buffer_init(this%FMS_diag_output_buffers,SIZE(diag_yaml%get_diag_fields()))
148 this%registered_variables = 0
149 this%registered_axis = 0
150 this%data_was_send = .false.
151 this%initialized = .true.
152#else
153 call mpp_error("fms_diag_object_init",&
154 "You must compile with -Duse_yaml to use the option use_modern_diag", fatal)
155#endif
156end subroutine fms_diag_object_init
157
158!> @brief Writes out the diag manifest file
159subroutine write_diag_manifest(this)
160 class(fmsdiagobject_type) :: this !< Diag object
161
162 integer, allocatable :: ntimes(:) !< Number of times written in each file
163 integer, allocatable :: ntiles(:) !< Number of tiles for each file domain
164 integer, allocatable :: ndistributedfiles(:) !< Number of distributed files
165
166 integer :: i !< For looping through the files
167 integer :: nfiles !< Number of files in the diag object
168
169#ifdef use_yaml
170 nfiles = size(this%FMS_diag_files)
171 allocate(ntimes(nfiles))
172 allocate(ntiles(nfiles))
173 allocate(ndistributedfiles(nfiles))
174
175 do i = 1, size(this%FMS_diag_files)
176 ntimes(i) = this%FMS_diag_files(i)%get_num_time_levels()
177 ntiles(i) = this%FMS_diag_files(i)%get_num_tiles()
178 ndistributedfiles(i) = this%FMS_diag_files(i)%get_ndistributedfiles()
179 enddo
180 call fms_diag_yaml_out(ntimes, ntiles, ndistributedfiles)
181
182#else
183 call mpp_error(fatal, "You must compile with -Duse_yaml to call fms_diag_object%write_diag_manifest!")
184#endif
185
186end subroutine write_diag_manifest
187
188!> \description Loops through all files and does one final write.
189!! Closes all files
190!! Deallocates all buffers, fields, and files
191!! Uninitializes the fms_diag_object
192subroutine fms_diag_object_end (this, time)
193 class(fmsdiagobject_type) :: this
194 TYPE(time_type), INTENT(in) :: time
195
196 integer :: i
197#ifdef use_yaml
198 !TODO: loop through files and force write
199 if (.not. this%initialized) return
200
201 call this%do_buffer_math()
202 call this%fms_diag_do_io(end_time=time)
203
204 call this%write_diag_manifest()
205
206 !TODO: Deallocate diag object arrays and clean up all memory
207 do i=1, size(this%FMS_diag_output_buffers)
208 call this%FMS_diag_output_buffers(i)%flush_buffer()
209 enddo
210 deallocate(this%FMS_diag_output_buffers)
211 this%axes_initialized = fms_diag_axis_object_end(this%diag_axis)
212 this%initialized = .false.
214#else
215 call mpp_error(fatal, "You can not call fms_diag_object%end without yaml")
216#endif
217end subroutine fms_diag_object_end
218
219!> @brief Registers a field.
220!! @description This to avoid having duplicate code in each of the _scalar, _array and _static register calls
221!! @return field index to be used in subsequent calls to send_data or DIAG_FIELD_NOT_FOUND if the field is not
222!! in the diag_table.yaml
223integer function fms_register_diag_field_obj &
224 (this, modname, varname, axes, init_time, &
225 longname, units, missing_value, varrange, mask_variant, standname, &
226 do_not_log, err_msg, interp_method, tile_count, area, volume, realm, static, &
227 multiple_send_data)
228
229 class(fmsdiagobject_type),TARGET,INTENT(inout):: this !< Diaj_obj to fill
230 CHARACTER(len=*), INTENT(in) :: modname !< The module name
231 CHARACTER(len=*), INTENT(in) :: varname !< The variable name
232 TYPE(time_type), OPTIONAL, INTENT(in) :: init_time !< Initial time
233 INTEGER, TARGET, OPTIONAL, INTENT(in) :: axes(:) !< The axes indicies
234 CHARACTER(len=*), OPTIONAL, INTENT(in) :: longname !< THe variables long name
235 CHARACTER(len=*), OPTIONAL, INTENT(in) :: units !< The units of the variables
236 CHARACTER(len=*), OPTIONAL, INTENT(in) :: standname !< The variables stanard name
237 class(*), OPTIONAL, INTENT(in) :: missing_value !< Missing value to add as a attribute
238 class(*), OPTIONAL, INTENT(in) :: varrange(2) !< Range to add as a attribute
239 LOGICAL, OPTIONAL, INTENT(in) :: mask_variant !< .True. if mask changes over time
240 LOGICAL, OPTIONAL, INTENT(in) :: do_not_log !< if TRUE, field info is not logged
241 CHARACTER(len=*), OPTIONAL, INTENT(out) :: err_msg !< Error message to be passed back up
242 CHARACTER(len=*), OPTIONAL, INTENT(in) :: interp_method !< The interp method to be used when
243 !! regridding the field in post-processing.
244 !! Valid options are "conserve_order1",
245 !! "conserve_order2", and "none".
246 INTEGER, OPTIONAL, INTENT(in) :: tile_count !< the number of tiles
247 INTEGER, OPTIONAL, INTENT(in) :: area !< diag_field_id of the cell area field
248 INTEGER, OPTIONAL, INTENT(in) :: volume !< diag_field_id of the cell volume field
249 CHARACTER(len=*), OPTIONAL, INTENT(in) :: realm !< String to set as the value to the
250 !! modeling_realm attribute
251 LOGICAL, OPTIONAL, INTENT(in) :: static !< True if the variable is static
252 LOGICAL, OPTIONAL, INTENT(in) :: multiple_send_data !< .True. if send data is called, multiple
253 !! times for the same time
254
255#ifdef use_yaml
256
257 class(fmsdiagfile_type), pointer :: fileptr !< Pointer to the diag_file
258 class(fmsdiagfield_type), pointer :: fieldptr !< Pointer to the diag_field
259 class(fmsdiagoutputbuffer_type), pointer :: bufferptr !< Pointer to the output buffer
260 class(diagyamlfilesvar_type), pointer :: yamlfptr !< Pointer to yaml object to get the reduction method
261 integer, allocatable :: file_ids(:) !< The file IDs for this variable
262 integer :: i !< For do loops
263 integer, allocatable :: diag_field_indices(:) !< indices where the field was found in the yaml
264 class(diagdomain_t), pointer :: null_diag_domain => null() !< Workaround for a Cray bug which will be fixed in CCE 19
265#endif
266#ifndef use_yaml
267fms_register_diag_field_obj = diag_field_not_found
268CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
269#else
270 diag_field_indices = find_diag_field(varname, modname)
271 if (diag_field_indices(1) .eq. diag_null) then
272 !< The field was not found in the table, so return DIAG_FIELD_NOT_FOUND
273 fms_register_diag_field_obj = diag_field_not_found
274 deallocate(diag_field_indices)
275 return
276 endif
277
278 this%registered_variables = this%registered_variables + 1
279 fms_register_diag_field_obj = this%registered_variables
280
281 call this%FMS_diag_fields(this%registered_variables)%&
282 &setid(this%registered_variables)
283
284!> Use pointers for convenience
285 fieldptr => this%FMS_diag_fields(this%registered_variables)
286!> Get the file IDs from the field indicies from the yaml
287 file_ids = get_diag_files_id(diag_field_indices)
288 call fieldptr%set_file_ids(file_ids)
289
290!> Allocate and initialize member buffer_allocated of this field
291 fieldptr%buffer_allocated = .false.
292 fieldptr%buffer_ids = get_diag_field_ids(diag_field_indices)
293
294!> Register the data for the field
295 call fieldptr%register(modname, varname, diag_field_indices, this%diag_axis, &
296 axes=axes, longname=longname, units=units, missing_value=missing_value, varrange= varrange, &
297 mask_variant= mask_variant, standname=standname, do_not_log=do_not_log, err_msg=err_msg, &
298 interp_method=interp_method, tile_count=tile_count, area=area, volume=volume, realm=realm, &
299 static=static, multiple_send_data=multiple_send_data)
300
301!> Add the axis information, initial time, and field IDs to the files
302 if (present(axes) .and. present(init_time)) then
303 do i = 1, size(file_ids)
304 fileptr => this%FMS_diag_files(file_ids(i))%FMS_diag_file
305 call fileptr%add_field_and_yaml_id(fieldptr%get_id(), diag_field_indices(i))
306 call fileptr%add_buffer_id(fieldptr%buffer_ids(i))
307 if(fieldptr%get_type_of_domain() .eq. no_domain) then
308 call fileptr%set_file_domain(null_diag_domain, fieldptr%get_type_of_domain())
309 else
310 call fileptr%set_file_domain(fieldptr%get_domain(), fieldptr%get_type_of_domain())
311 endif
312 call fileptr%init_diurnal_axis(this%diag_axis, this%registered_axis, diag_field_indices(i))
313 call fileptr%add_axes(axes, this%diag_axis, this%registered_axis, diag_field_indices(i), &
314 fieldptr%buffer_ids(i), this%FMS_diag_output_buffers)
315 call fileptr%add_start_time(init_time)
316 call fileptr%set_file_time_ops (fieldptr%diag_field(i), fieldptr%is_static())
317 enddo
318 elseif (present(axes)) then !only axes present
319 do i = 1, size(file_ids)
320 fileptr => this%FMS_diag_files(file_ids(i))%FMS_diag_file
321 call fileptr%add_field_and_yaml_id(fieldptr%get_id(), diag_field_indices(i))
322 call fileptr%add_buffer_id(fieldptr%buffer_ids(i))
323 call fileptr%init_diurnal_axis(this%diag_axis, this%registered_axis, diag_field_indices(i))
324 if(fieldptr%get_type_of_domain() .eq. no_domain) then
325 call fileptr%set_file_domain(null_diag_domain, fieldptr%get_type_of_domain())
326 else
327 call fileptr%set_file_domain(fieldptr%get_domain(), fieldptr%get_type_of_domain())
328 endif
329 call fileptr%add_axes(axes, this%diag_axis, this%registered_axis, diag_field_indices(i), &
330 fieldptr%buffer_ids(i), this%FMS_diag_output_buffers)
331 call fileptr%set_file_time_ops (fieldptr%diag_field(i), fieldptr%is_static())
332 enddo
333 elseif (present(init_time)) then !only inti time present
334 do i = 1, size(file_ids)
335 fileptr => this%FMS_diag_files(file_ids(i))%FMS_diag_file
336 call fileptr%add_field_and_yaml_id(fieldptr%get_id(), diag_field_indices(i))
337 call fileptr%add_buffer_id(fieldptr%buffer_ids(i))
338 call fileptr%add_start_time(init_time)
339 call fileptr%set_file_time_ops (fieldptr%diag_field(i), fieldptr%is_static())
340 enddo
341 else !no axis or init time present
342 do i = 1, size(file_ids)
343 fileptr => this%FMS_diag_files(file_ids(i))%FMS_diag_file
344 call fileptr%add_field_and_yaml_id(fieldptr%get_id(), diag_field_indices(i))
345 call fileptr%add_buffer_id(fieldptr%buffer_ids(i))
346 call fileptr%set_file_time_ops (fieldptr%diag_field(i), fieldptr%is_static())
347 enddo
348 endif
349
350 !> Initialize buffer_ids of this field with the diag_field_indices(diag_field_indices)
351!! of the sorted variable list
352 do i = 1, size(fieldptr%buffer_ids)
353 bufferptr => this%FMS_diag_output_buffers(fieldptr%buffer_ids(i))
354 call bufferptr%set_field_id(this%registered_variables)
355 call bufferptr%set_yaml_id(fieldptr%buffer_ids(i))
356 ! check if diurnal reduction for this buffer and if so set the diurnal sample size
357 yamlfptr => diag_yaml%diag_fields(fieldptr%buffer_ids(i))
358 if( yamlfptr%get_var_reduction() .eq. time_diurnal) then
359 call bufferptr%set_diurnal_sample_size(yamlfptr%get_n_diurnal())
360 endif
361 call bufferptr%init_buffer_time(init_time)
362 call bufferptr%set_next_output(this%FMS_diag_files(file_ids(i))%get_next_output(), &
363 this%FMS_diag_files(file_ids(i))%get_next_next_output(), is_static=fieldptr%is_static())
364 enddo
365
366 nullify (fileptr)
367 nullify (fieldptr)
368 deallocate(diag_field_indices)
369#endif
370end function fms_register_diag_field_obj
371
372!> @brief Registers a scalar field
373!! @return field index to be used in subsequent calls to send_data or DIAG_FIELD_NOT_FOUND if the field is not
374!! in the diag_table.yaml
375INTEGER FUNCTION fms_register_diag_field_scalar(this,module_name, field_name, init_time, &
376 & long_name, units, missing_value, var_range, standard_name, do_not_log, err_msg,&
377 & area, volume, realm, multiple_send_data)
378 class(fmsdiagobject_type),TARGET,INTENT(inout):: this !< Diaj_obj to fill
379 CHARACTER(len=*), INTENT(in) :: module_name !< Module where the field comes from
380 CHARACTER(len=*), INTENT(in) :: field_name !< Name of the field
381 TYPE(time_type), OPTIONAL, INTENT(in) :: init_time !< Time to start writing data from
382 CHARACTER(len=*), OPTIONAL, INTENT(in) :: long_name !< Long_name to add as a variable attribute
383 CHARACTER(len=*), OPTIONAL, INTENT(in) :: units !< Units to add as a variable_attribute
384 CHARACTER(len=*), OPTIONAL, INTENT(in) :: standard_name !< Standard_name to name the variable in the file
385 CLASS(*), OPTIONAL, INTENT(in) :: missing_value !< Missing value to add as a variable attribute
386 CLASS(*), OPTIONAL, INTENT(in) :: var_range(:) !< Range to add a variable attribute
387 LOGICAL, OPTIONAL, INTENT(in) :: do_not_log !< If TRUE, field information is not logged
388 CHARACTER(len=*), OPTIONAL, INTENT(out):: err_msg !< Error_msg from call
389 INTEGER, OPTIONAL, INTENT(in) :: area !< Id of the area field
390 INTEGER, OPTIONAL, INTENT(in) :: volume !< Id of the volume field
391 CHARACTER(len=*), OPTIONAL, INTENT(in) :: realm !< String to set as the modeling_realm attribute
392 LOGICAL, OPTIONAL, INTENT(in) :: multiple_send_data !< .True. if send data is called, multiple times
393 !! for the same time
394
395#ifndef use_yaml
396fms_register_diag_field_scalar=diag_field_not_found
397CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
398#else
399 fms_register_diag_field_scalar = this%register(&
400 & module_name, field_name, init_time=init_time, &
401 & longname=long_name, units=units, missing_value=missing_value, varrange=var_range, &
402 & standname=standard_name, do_not_log=do_not_log, err_msg=err_msg, &
403 & area=area, volume=volume, realm=realm, multiple_send_data=multiple_send_data)
404#endif
405end function fms_register_diag_field_scalar
406
407!> @brief Registers an array field
408!! @return field index to be used in subsequent calls to send_data or DIAG_FIELD_NOT_FOUND if the field is not
409!! in the diag_table.yaml
410INTEGER FUNCTION fms_register_diag_field_array(this, module_name, field_name, axes, init_time, &
411 & long_name, units, missing_value, var_range, mask_variant, standard_name, verbose,&
412 & do_not_log, err_msg, interp_method, tile_count, area, volume, realm, &
413 & multiple_send_data)
414 class(fmsdiagobject_type),TARGET,INTENT(inout):: this !< Diaj_obj to fill
415 CHARACTER(len=*), INTENT(in) :: module_name !< Module where the field comes from
416 CHARACTER(len=*), INTENT(in) :: field_name !< Name of the field
417 INTEGER, INTENT(in) :: axes(:) !< Ids corresponding to the variable axis
418 TYPE(time_type), OPTIONAL, INTENT(in) :: init_time !< Time to start writing data from
419 CHARACTER(len=*), OPTIONAL, INTENT(in) :: long_name !< Long_name to add as a variable attribute
420 CHARACTER(len=*), OPTIONAL, INTENT(in) :: units !< Units to add as a variable_attribute
421 CLASS(*), OPTIONAL, INTENT(in) :: missing_value !< Missing value to add as a variable attribute
422 CLASS(*), OPTIONAL, INTENT(in) :: var_range(:) !< Range to add a variable attribute
423 LOGICAL, OPTIONAL, INTENT(in) :: mask_variant !< .True. if mask changes over time
424 CHARACTER(len=*), OPTIONAL, INTENT(in) :: standard_name !< Standard_name to name the variable in the file
425 LOGICAL, OPTIONAL, INTENT(in) :: verbose !< Print more information
426 LOGICAL, OPTIONAL, INTENT(in) :: do_not_log !< If TRUE, field information is not logged
427 CHARACTER(len=*), OPTIONAL, INTENT(out):: err_msg !< Error_msg from call
428 CHARACTER(len=*), OPTIONAL, INTENT(in) :: interp_method !< The interp method to be used when
429 !! regridding the field in post-processing.
430 !! Valid options are "conserve_order1",
431 !! "conserve_order2", and "none".
432 INTEGER, OPTIONAL, INTENT(in) :: tile_count !< The current tile number
433 INTEGER, OPTIONAL, INTENT(in) :: area !< Id of the area field
434 INTEGER, OPTIONAL, INTENT(in) :: volume !< Id of the volume field
435 CHARACTER(len=*), OPTIONAL, INTENT(in) :: realm !< String to set as the modeling_realm attribute
436 LOGICAL, OPTIONAL, INTENT(in) :: multiple_send_data !< .True. if send data is called, multiple times
437 !! for the same time
438
439
440#ifndef use_yaml
441fms_register_diag_field_array=diag_field_not_found
442CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
443#else
444 fms_register_diag_field_array = this%register( &
445 & module_name, field_name, init_time=init_time, &
446 & axes=axes, longname=long_name, units=units, missing_value=missing_value, varrange=var_range, &
447 & mask_variant=mask_variant, standname=standard_name, do_not_log=do_not_log, err_msg=err_msg, &
448 & interp_method=interp_method, tile_count=tile_count, area=area, volume=volume, realm=realm, &
449 & multiple_send_data=multiple_send_data)
450#endif
451end function fms_register_diag_field_array
452
453!> @brief Return field index for subsequent call to send_data.
454!! @return field index to be used in subsequent calls to send_data or DIAG_FIELD_NOT_FOUND if the field is not
455!! in the diag_table.yaml
456INTEGER FUNCTION fms_register_static_field(this, module_name, field_name, axes, long_name, units,&
457 & missing_value, range, mask_variant, standard_name, DYNAMIC, do_not_log, interp_method,&
458 & tile_count, area, volume, realm)
459 class(fmsdiagobject_type),TARGET,INTENT(inout):: this !< Diaj_obj to fill
460 CHARACTER(len=*), INTENT(in) :: module_name !< Name of the module, the field is on
461 CHARACTER(len=*), INTENT(in) :: field_name !< Name of the field
462 INTEGER, DIMENSION(:), INTENT(in) :: axes !< Axes_id of the field
463 CHARACTER(len=*), OPTIONAL, INTENT(in) :: long_name !< Longname to be added as a attribute
464 CHARACTER(len=*), OPTIONAL, INTENT(in) :: units !< Units to be added as a attribute
465 CHARACTER(len=*), OPTIONAL, INTENT(in) :: standard_name !< Standard name to be added as a attribute
466 CLASS(*), OPTIONAL, INTENT(in) :: missing_value !< Missing value to be added as a attribute
467 CLASS(*), OPTIONAL, INTENT(in) :: range(:) !< Range to be added as a attribute
468 LOGICAL, OPTIONAL, INTENT(in) :: mask_variant !< .True. if mask changes over time
469 LOGICAL, OPTIONAL, INTENT(in) :: dynamic !< Flag indicating if the field is dynamic
470 LOGICAL, OPTIONAL, INTENT(in) :: do_not_log !< if TRUE, field information is not logged
471 CHARACTER(len=*), OPTIONAL, INTENT(in) :: interp_method !< The interp method to be used when
472 !! regridding the field in post-processing
473 !! Valid options are "conserve_order1",
474 !! "conserve_order2", and "none".
475 INTEGER, OPTIONAL, INTENT(in) :: tile_count !! Number of tiles
476 INTEGER, OPTIONAL, INTENT(in) :: area !< Field ID for the area field associated
477 !! with this field
478 INTEGER, OPTIONAL, INTENT(in) :: volume !< Field ID for the volume field associated
479 !! with this field
480 CHARACTER(len=*), OPTIONAL, INTENT(in) :: realm !< String to set as the value to the
481 !! modeling_realm attribute
482
483#ifndef use_yaml
484fms_register_static_field=diag_field_not_found
485CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
486#else
487 !TODO The register_static_field interface does not have the capabiliy to register a variable as a "scalar"
488 ! since the axes argument is required, this forced model code to pass in a null_axis_id as an argument
489 if (size(axes) .eq. 1 .and. axes(1) .eq. null_axis_id) then
490 ! If they are passing in the null_axis_ids, ignore the `axes` argument
491 fms_register_static_field = this%register( &
492 & module_name, field_name, &
493 & longname=long_name, units=units, missing_value=missing_value, varrange=range, &
494 & mask_variant=mask_variant, do_not_log=do_not_log, interp_method=interp_method, tile_count=tile_count, &
495 & standname=standard_name, area=area, volume=volume, realm=realm, &
496 & static=.true.)
497 else
498 fms_register_static_field = this%register( &
499 & module_name, field_name, axes=axes, &
500 & longname=long_name, units=units, missing_value=missing_value, varrange=range, &
501 & mask_variant=mask_variant, do_not_log=do_not_log, interp_method=interp_method, tile_count=tile_count, &
502 & standname=standard_name, area=area, volume=volume, realm=realm, &
503 & static=.true.)
504 endif
505#endif
506end function fms_register_static_field
507
508!> @brief Wrapper for the register_diag_axis subroutine. This is needed to keep the diag_axis_init
509!! interface the same
510!> @return Axis id
511FUNCTION fms_diag_axis_init(this, axis_name, axis_data, units, cart_name, axis_length, long_name, direction,&
512 & set_name, edges, Domain, Domain2, DomainU, aux, req, tile_count, domain_position ) &
513 & result(id)
514
515 class(fmsdiagobject_type),TARGET,INTENT(inout):: this !< Diaj_obj to fill
516 CHARACTER(len=*), INTENT(in) :: axis_name !< Name of the axis
517 CLASS(*), INTENT(in) :: axis_data(:) !< Array of coordinate values
518 CHARACTER(len=*), INTENT(in) :: units !< Units for the axis
519 CHARACTER(len=1), INTENT(in) :: cart_name !< Cartesian axis ("X", "Y", "Z", "T", "U", "N")
520 integer, intent(in) :: axis_length !< The length of the axis size(axis_data(:))
521 CHARACTER(len=*), INTENT(in), OPTIONAL :: long_name !< Long name for the axis.
522 CHARACTER(len=*), INTENT(in), OPTIONAL :: set_name !< Name of the parent axis, if it is a subaxis
523 INTEGER, INTENT(in), OPTIONAL :: direction !< Indicates the direction of the axis
524 INTEGER, INTENT(in), OPTIONAL :: edges !< Axis ID for the previously defined "edges axis"
525 TYPE(domain1d), INTENT(in), OPTIONAL :: domain !< 1D domain
526 TYPE(domain2d), INTENT(in), OPTIONAL :: domain2 !< 2D domain
527 TYPE(domainug), INTENT(in), OPTIONAL :: domainu !< Unstructured domain
528 CHARACTER(len=*), INTENT(in), OPTIONAL :: aux !< Auxiliary name, can only be <TT>geolon_t</TT>
529 !! or <TT>geolat_t</TT>
530 CHARACTER(len=*), INTENT(in), OPTIONAL :: req !< Required field names.
531 INTEGER, INTENT(in), OPTIONAL :: tile_count !< Number of tiles
532 INTEGER, INTENT(in), OPTIONAL :: domain_position !< Domain position, "NORTH" or "EAST"
533 integer :: id
534
535#ifndef use_yaml
536id = diag_null
537CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
538#else
539 CHARACTER(len=:), ALLOCATABLE :: edges_name !< Name of the edges
540
541 this%registered_axis = this%registered_axis + 1
542
543 if (this%registered_axis > max_axes) call mpp_error(fatal, &
544 &"diag_axis_init: max_axes exceeded, increase via diag_manager_nml")
545
546 allocate(fmsdiagfullaxis_type :: this%diag_axis(this%registered_axis)%axis)
547
548 select type (axis => this%diag_axis(this%registered_axis)%axis )
549 type is (fmsdiagfullaxis_type)
550 if(present(edges)) then
551 if (edges < 0 .or. edges > this%registered_axis) &
552 call mpp_error(fatal, "diag_axit_init: The edge axis has not been defined. &
553 &Call diag_axis_init for the edge axis first")
554 select type (edges_axis => this%diag_axis(edges)%axis)
555 type is (fmsdiagfullaxis_type)
556 edges_name = edges_axis%get_axis_name()
557 call axis%set_edges(edges_name, edges)
558 end select
559 endif
560 call axis%register(axis_name, axis_data, units, cart_name, long_name=long_name, &
561 & direction=direction, set_name=set_name, domain=domain, domain2=domain2, domainu=domainu, aux=aux, &
562 & req=req, tile_count=tile_count, domain_position=domain_position, axis_length=axis_length)
563
564 id = this%registered_axis
565 call axis%set_axis_id(id)
566 end select
567#endif
568end function fms_diag_axis_init
569
570!> Accepts data from the send_data functions. If this is in an openmp region with more than
571!! one thread, the data is buffered in the field object and processed later. If only a single thread
572!! is being used, then the processing can be done and stored in the buffer object. The hope is that
573!! the increase in memory footprint related to buffering can be handled by the shared memory of the
574!! multithreaded case.
575!! \note If some of the diag manager is offloaded in the future, then it should be treated similarly
576!! to the multi-threaded option for processing later
577subroutine fms_diag_accept_data (this, diag_field_id, field_data, mask, rmask, &
578 time, is_in, js_in, ks_in, &
579 ie_in, je_in, ke_in, weight, err_msg)
580 class(fmsdiagobject_type),TARGET, INTENT(inout) :: this !< Diaj_obj to fill
581 INTEGER, INTENT(in) :: diag_field_id !< The ID of the diag field
582 CLASS(*), DIMENSION(:,:,:,:), INTENT(in) :: field_data !< The data for the diag_field
583 LOGICAL, allocatable, INTENT(in) :: mask(:,:,:,:) !< Logical mask indicating the grid
584 !! points to mask (null if no mask)
585 CLASS(*), allocatable, INTENT(in) :: rmask(:,:,:,:)!< real mask indicating the grid
586 !! points to mask (null if no mask)
587 CLASS(*), INTENT(in), OPTIONAL :: weight !< The weight used for averaging
588 TYPE (time_type), INTENT(in), OPTIONAL :: time !< The current time
589 INTEGER, INTENT(in), OPTIONAL :: is_in, js_in, ks_in !< Starting indices
590 INTEGER, INTENT(in), OPTIONAL :: ie_in, je_in, ke_in !< Ending indices
591 CHARACTER(len=*), INTENT(out), OPTIONAL :: err_msg !< An error message returned
592
593 integer :: is, js, ks !< Starting indicies of the field_data
594 integer :: ie, je, ke !< Ending indicies of the field_data
595 integer :: omp_num_threads !< Number of openmp threads
596 integer :: omp_level !< The openmp active level
597 logical :: buffer_the_data !< True if the user selects to buffer the data and run
598 !! the calculationslater. \note This is experimental
599 character(len=128) :: error_string !< Store error text
600 logical :: data_buffer_is_allocated !< .true. if the data buffer is allocated
601 character(len=256) :: field_info !< String holding info about the field to append to the
602 !! error message
603 logical, allocatable, dimension(:,:,:,:) :: oor_mask !< Out of range mask
604 real(kind=r8_kind) :: field_weight !< Weight to use when averaging (it will be converted
605 !! based on the type of field_data when doing the math)
606 type(fmsdiagibounds_type) :: bounds !< Bounds (starting ending indices) for the field
607 logical :: has_halos !< .True. if field_data contains halos
608 logical :: using_blocking !< .True. if field_data is passed in blocks
609#ifndef use_yaml
610CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
611#else
612
613 !TODO this%FMS_diag_fields(diag_field_id) should be a pointer!
614 field_info = " Check send data call for field:"//trim(this%FMS_diag_fields(diag_field_id)%get_varname())//&
615 " and module:"//trim(this%FMS_diag_fields(diag_field_id)%get_modname())
616
617 !< Check if time should be present for this field
618 if (.not.this%FMS_diag_fields(diag_field_id)%is_static() .and. .not.present(time)) &
619 call mpp_error(fatal, "Time must be present if the field is not static. "//trim(field_info))
620
621 !< Set the field_weight. If "weight" is not present it will be set to 1.0_r8_kind
622 field_weight = set_weight(weight)
623
624 !< Check that the indices are present in the correct combination
625 error_string = check_indices_order(is_in, ie_in, js_in, je_in)
626 if (trim(error_string) .ne. "") call mpp_error(fatal, trim(error_string)//". "//trim(field_info))
627
628 using_blocking = .false.
629 if ((present(is_in) .and. .not. present(ie_in)) .or. (present(js_in) .and. .not. present(je_in))) &
630 using_blocking = .true.
631
632 has_halos = .false.
633 if ((present(is_in) .and. present(ie_in)) .or. (present(js_in) .and. present(je_in))) &
634 has_halos = .true.
635
636 !< If the field has `mask_variant=.true.`, check that mask OR rmask are present
637 if (this%FMS_diag_fields(diag_field_id)%is_mask_variant()) then
638 if (.not. allocated(mask) .and. .not. allocated(rmask)) call mpp_error(fatal, &
639 "The field was registered with mask_variant, but mask or rmask are not present in the send_data call. "//&
640 trim(field_info))
641 endif
642
643 !< Check that mask and rmask are not both present
644 if (allocated(mask) .and. allocated(rmask)) call mpp_error(fatal, &
645 "mask and rmask are both present in the send_data call. "//&
646 trim(field_info))
647
648 !< Create the oor_mask based on the "mask" and "rmask" arguments
649 oor_mask = init_mask(rmask, mask, field_data)
650
651 !> Does the user want to push off calculations until send_diag_complete?
652 buffer_the_data = .false.
653
654 !> initialize the number of threads and level to be 0
655 omp_num_threads = 0
656 omp_level = 0
657#if defined(_OPENMP)
658 omp_num_threads = omp_get_num_threads()
659 omp_level = omp_get_level()
660 buffer_the_data = (omp_num_threads > 1 .AND. omp_level > 0)
661#endif
662
663 !> Calculate the i,j,k start and end
664 ! If is, js, or ks not present default them to 1
665 is = 1
666 js = 1
667 ks = 1
668 IF ( PRESENT(is_in) ) is = is_in
669 IF ( PRESENT(js_in) ) js = js_in
670 IF ( PRESENT(ks_in) ) ks = ks_in
671 ie = is+SIZE(field_data, 1)-1
672 je = js+SIZE(field_data, 2)-1
673 ke = ks+SIZE(field_data, 3)-1
674 IF ( PRESENT(ie_in) ) ie = ie_in
675 IF ( PRESENT(je_in) ) je = je_in
676 IF ( PRESENT(ke_in) ) ke = ke_in
677
678 if (.not. buffer_the_data .and. using_blocking) then
679 ! If running with only 1 thread and using blocking, check if the data was sent in blocks
680 ! if it is, then buffer the data
681 buffer_the_data = check_for_slices(this%FMS_diag_fields(diag_field_id), this%diag_axis, &
682 shape(field_data))
683 endif
684
685 !< If send data is called multiple times, buffer the data
686 !! This is so that the other reduction methods work and just averaging
687 if (this%FMS_diag_fields(diag_field_id)%get_multiple_send_data()) &
688 buffer_the_data = .true.
689
690 !If this is true, buffer data
691 main_if: if (buffer_the_data) then
692!> Only 1 thread allocates the output buffer and sets set_math_needs_to_be_done
693!$omp critical
694 !< Let diag_send_complete that there is new data to procress
695 if (.not. this%data_was_send) this%data_was_send = .true.
696
697 !< These set_* calls need to be done inside an omp_critical to avoid any race conditions
698 !! and allocation issues
699 if(has_halos) call this%FMS_diag_fields(diag_field_id)%set_halo_present()
700
701 !< Set the variable type based off passed in field data
702 if(.not. this%FMS_diag_fields(diag_field_id)%has_vartype()) &
703 call this%FMS_diag_fields(diag_field_id)%set_type(field_data(1,1,1,1))
704
705 if (allocated(mask) .or. allocated(rmask)) then
706 call this%FMS_diag_fields(diag_field_id)%set_var_is_masked(.true.)
707 else
708 call this%FMS_diag_fields(diag_field_id)%set_var_is_masked(.false.)
709 endif
710
711 if (.not. this%FMS_diag_fields(diag_field_id)%is_data_buffer_allocated()) then
712 data_buffer_is_allocated = &
713 this%FMS_diag_fields(diag_field_id)%allocate_data_buffer(field_data, this%diag_axis)
714 if(.not. this%FMS_diag_fields(diag_field_id)%has_mask_allocated()) &
715 call this%FMS_diag_fields(diag_field_id)%allocate_mask(oor_mask, this%diag_axis)
716 endif
717 call this%FMS_diag_fields(diag_field_id)%set_send_data_time(time)
718 call this%FMS_diag_fields(diag_field_id)%set_data_buffer_is_allocated(.true.)
719 call this%FMS_diag_fields(diag_field_id)%set_math_needs_to_be_done(.true.)
720!$omp end critical
721 call this%FMS_diag_fields(diag_field_id)%set_data_buffer(field_data, oor_mask, field_weight, &
722 is, js, ks, ie, je, ke)
723 else
724
725 !< Let diag_send_complete that there is new data to procress
726 if (.not. this%data_was_send) this%data_was_send = .true.
727
728 !< At this point if we are no longer in an openmp region or running with 1 thread
729 !! so it is safe to have these set_* calls
730 if(has_halos) call this%FMS_diag_fields(diag_field_id)%set_halo_present()
731
732 !< Set the variable type based off passed in field data
733 if(.not. this%FMS_diag_fields(diag_field_id)%has_vartype()) &
734 call this%FMS_diag_fields(diag_field_id)%set_type(field_data(1,1,1,1))
735
736 if (allocated(mask) .or. allocated(rmask)) then
737 call this%FMS_diag_fields(diag_field_id)%set_var_is_masked(.true.)
738 else
739 call this%FMS_diag_fields(diag_field_id)%set_var_is_masked(.false.)
740 endif
741
742 error_string = bounds%set_bounds(field_data, is, ie, js, je, ks, ke, has_halos)
743 if (trim(error_string) .ne. "") call mpp_error(fatal, trim(error_string)//". "//trim(field_info))
744
745 call this%allocate_diag_field_output_buffers(field_data, diag_field_id)
746 error_string = this%fms_diag_do_reduction(field_data, diag_field_id, oor_mask, field_weight, &
747 bounds, using_blocking, time=time)
748 if (trim(error_string) .ne. "") call mpp_error(fatal, trim(error_string)//". "//trim(field_info))
749 call this%FMS_diag_fields(diag_field_id)%set_math_needs_to_be_done(.false.)
750 if(.not. this%FMS_diag_fields(diag_field_id)%has_mask_allocated()) &
751 call this%FMS_diag_fields(diag_field_id)%allocate_mask(oor_mask)
752 call this%FMS_diag_fields(diag_field_id)%set_mask(oor_mask, field_info)
753 end if main_if
754 !> Return false if nothing is done
755#endif
756end subroutine fms_diag_accept_data
757
758!< @brief Do the math for all the buffers
759subroutine do_buffer_math(this)
760 class(fmsdiagobject_type), target, intent (inout) :: this !< The diag object
761
762#ifdef use_yaml
763 integer :: i !< For do loops
764 integer :: ifile !< For file loops
765 integer :: ifield !< For field loops
766
767 class(fmsdiagfilecontainer_type), pointer :: diag_file !< Pointer to this%FMS_diag_files(i) (for convenience
768 class(fmsdiagfield_type), pointer :: diag_field !< Pointer to this%FMS_diag_files(i)%diag_field(j)
769 logical :: math !< True if the math functions need to be called using the data buffer,
770 !! False if the math functions were done in accept_data
771 integer, dimension(:), allocatable :: file_field_ids !< Array of field IDs for a file
772 class(*), pointer :: input_data_buffer(:,:,:,:)
773 character(len=128) :: error_string
774 type(fmsdiagibounds_type) :: bounds
775 integer, dimension(:), allocatable :: file_ids !< Array of file IDs for a field
776 logical, parameter :: debug_sc = .false. !< turn on output for debugging
777
778 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
779 !! In the future, this may be parallelized for offloading
780 ! loop through each field
781 field_loop: do ifield = 1, size(this%FMS_diag_fields)
782 diag_field => this%FMS_diag_fields(ifield)
783 if(.not. diag_field%is_registered()) cycle
784 if(debug_sc) call mpp_error(note, "fms_diag_send_complete:: var: "//diag_field%get_varname())
785 ! get files the field is in
786 allocate (file_ids(size(diag_field%get_file_ids() )))
787 file_ids = diag_field%get_file_ids()
788 math = diag_field%get_math_needs_to_be_done()
789 ! if doing math loop through each file for given field
790 doing_math: if (size(file_ids) .ge. 1 .and. math) then
791 ! Check if buffer alloc'd
792 has_input_buff: if (diag_field%has_input_data_buffer()) then
793 call diag_field%prepare_data_buffer()
794 input_data_buffer => diag_field%get_data_buffer()
795 ! reset bounds, allocate output buffer, and update it with reduction
796 call bounds%reset_bounds_from_array_4D(input_data_buffer)
797 call this%allocate_diag_field_output_buffers(input_data_buffer, ifield)
798 error_string = this%fms_diag_do_reduction(input_data_buffer, ifield, &
799 diag_field%get_mask(), diag_field%get_weight(), &
800 bounds, .false., time=diag_field%get_send_data_time())
801 call diag_field%init_data_buffer()
802 if (trim(error_string) .ne. "") call mpp_error(fatal, "Field:"//trim(diag_field%get_varname()//&
803 " -"//trim(error_string)))
804 else
805 call mpp_error(fatal, "diag_send_complete:: no input buffer allocated for field"//diag_field%get_longname())
806 endif has_input_buff
807 endif doing_math
808 call diag_field%set_math_needs_to_be_done(.false.)
809 !> Clean up, clean up, everybody do your share
810 if (allocated(file_ids)) deallocate(file_ids)
811 if (associated(diag_field)) nullify(diag_field)
812 enddo field_loop
813#endif
814end subroutine do_buffer_math
815
816!> @brief Loops through all the files, open the file, writes out axis and
817!! variable metadata and data when necessary.
818subroutine fms_diag_send_complete(this, time_step)
819 class(fmsdiagobject_type), target, intent (inout) :: this !< The diag object
820 TYPE (time_type), INTENT(in) :: time_step !< The time_step
821
822#ifndef use_yaml
823CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
824#else
825 !< Go away if there is no new data
826 if (.not. this%data_was_send) return
827
828 call this%do_buffer_math()
829 call this%fms_diag_do_io()
830
831 this%data_was_send = .false.
832#endif
833
834end subroutine fms_diag_send_complete
835
836!> @brief Loops through all the files, open the file, writes out axis and
837!! variable metadata and data when necessary.
838!! TODO: passing in the saved mask from the field obj to diag_reduction_done_wrapper
839!! for performance
840subroutine fms_diag_do_io(this, end_time)
841 class(fmsdiagobject_type), target, intent(inout) :: this !< The diag object
842 type(time_type), optional, target, intent(in) :: end_time !< the model end_time
843#ifdef use_yaml
844 integer :: i !< For do loops
845 class(fmsdiagfilecontainer_type), pointer :: diag_file !< Pointer to this%FMS_diag_files(i) (for convenience)
846 class(fmsdiagoutputbuffer_type), pointer :: diag_buff !< pointer to output buffers iterated in buff_loop
847 class(fmsdiagfield_type), pointer :: diag_field !< pointer to output buffers iterated in buff_loop
848 class(diagyamlfilesvar_type), pointer :: field_yaml !< Pointer to a field from yaml fields
849 TYPE (time_type), pointer :: model_time!< The current model time
850 integer, allocatable :: buff_ids(:) !< ids for output buffers to loop through
851 integer :: ibuff !< buffer index
852 logical :: file_is_opened_this_time_step !< True if the file was opened in this time_step
853 !! If true the metadata will need to be written
854 logical :: force_write !< force the last write if at end of run
855 logical :: finish_writing !< true if finished writing for all the fields
856 logical :: has_mask !< whether we have a mask
857 logical, parameter :: debug_reduct = .false. !< enables debugging output
858 class(*), allocatable :: missing_val !< netcdf missing value for a given field
859 real(r8_kind) :: mval !< r8 copy of missing value
860 character(len=128) :: error_string !< outputted error string from reducti
861 logical :: unlim_dim_was_increased !< .True. if the unlimited dimension index was increased for any of the buffers
862 logical :: do_not_write !< .True. only if this is not a new time step and you are writting at every time step
863
864 force_write = .false.
865
866 do i = 1, size(this%FMS_diag_files)
867 diag_file => this%FMS_diag_files(i)
868
869 !< Go away if the file is a subregional file and the current PE does not have any data for it
870 if (.not. diag_file%writing_on_this_pe()) cycle
871
872 if (present (end_time)) then
873 force_write = .true.
874 model_time => end_time
875 else
876 model_time => diag_file%get_model_time()
877 endif
878 if (diag_file%FMS_diag_file%is_done_writing_data()) cycle
879
880 call diag_file%open_diag_file(model_time, file_is_opened_this_time_step)
881 if (file_is_opened_this_time_step) then
882 ! Initialize unlimited dimension in file and the buffer to 0
883 call diag_file%init_unlim_dim(this%FMS_diag_output_buffers)
884
885 call diag_file%write_global_metadata()
886 call diag_file%write_axis_metadata(this%diag_axis)
887 call diag_file%write_time_metadata()
888 call diag_file%write_field_metadata(this%FMS_diag_fields, this%diag_axis)
889 call diag_file%write_axis_data(this%diag_axis)
890 endif
891
892 call diag_file%check_file_times(model_time, this%FMS_diag_output_buffers, &
893 this%FMS_diag_fields, do_not_write)
894 unlim_dim_was_increased = .false.
895
896 ! finish reduction method if its time to write
897 buff_ids = diag_file%FMS_diag_file%get_buffer_ids()
898 ! loop through the buffers and finish reduction if needed
899 buff_loop: do ibuff=1, SIZE(buff_ids)
900 diag_buff => this%FMS_diag_output_buffers(buff_ids(ibuff))
901 field_yaml => diag_yaml%diag_fields(diag_buff%get_yaml_id())
902 diag_field => this%FMS_diag_fields(diag_buff%get_field_id())
903
904 ! Go away if there is no data to write
905 if (.not. diag_buff%is_there_data_to_write()) cycle
906 if (diag_field%is_static() .and. diag_buff%get_unlim_dim() > 0) cycle
907
908 if ( diag_buff%is_time_to_finish_reduction(end_time) .and. .not. do_not_write) then
909 ! sets missing value
910 mval = diag_field%find_missing_value(missing_val)
911 ! time_average and greater values all involve averaging so need to be "finished" before written
912 if( field_yaml%has_var_reduction()) then
913 if( field_yaml%get_var_reduction() .ge. time_average) then
914 if(debug_reduct)call mpp_error(note, "fms_diag_do_io:: finishing reduction for "//diag_field%get_longname())
915 error_string = diag_buff%diag_reduction_done_wrapper( &
916 field_yaml%get_var_reduction(), &
917 mval, diag_field%get_var_is_masked(), diag_field%get_mask_variant())
918 endif
919 endif
920 call diag_file%write_field_data(diag_field, diag_buff, unlim_dim_was_increased)
921 call diag_buff%set_next_output(diag_file%get_next_output(), diag_file%get_next_next_output())
922 endif
923 nullify(diag_buff)
924 nullify(field_yaml)
925 enddo buff_loop
926 deallocate(buff_ids)
927
928 if (unlim_dim_was_increased) then
929 call diag_file%write_time_data()
930 call diag_file%flush_diag_file()
931 call diag_file%update_next_write(model_time)
932 call diag_file%update_current_new_file_freq_index(model_time)
933 if (diag_file%is_time_to_close_file(model_time, force_write)) &
934 call diag_file%close_diag_file(this%FMS_diag_output_buffers, &
935 this%model_end_time, diag_fields = this%FMS_diag_fields)
936 else if (force_write) then
937 call diag_file%prepare_for_force_write()
938 call diag_file%write_time_data()
939 call diag_file%close_diag_file(this%FMS_diag_output_buffers, &
940 this%model_end_time, diag_fields = this%FMS_diag_fields)
941 endif
942 enddo
943#endif
944end subroutine fms_diag_do_io
945
946!> @brief Computes average, min, max, rms error, etc.
947!! based on the specified reduction method for the field.
948!> @return Empty string if successful, error message if it fails
949function fms_diag_do_reduction(this, field_data, diag_field_id, oor_mask, weight, &
950 bounds, using_blocking, time) &
951 result(error_msg)
952 class(fmsdiagobject_type), intent(inout), target:: this !< Diag Object
953 class(*), intent(in) :: field_data(:,:,:,:) !< Field data
954 integer, intent(in) :: diag_field_id !< ID of the input field
955 logical, intent(in), target :: oor_mask(:,:,:,:) !< mask
956 real(kind=r8_kind), intent(in) :: weight !< Must be a updated weight
957 type(fmsdiagibounds_type), intent(in) :: bounds !< Bounds for the field
958 logical, intent(in) :: using_blocking !< .True. if field data is passed
959 !! in blocks
960 type(time_type), intent(in), optional :: time !< Current time
961
962 character(len=150) :: error_msg !< Error message to check
963 !TODO Mostly everything
964#ifdef use_yaml
965 type(fmsdiagfield_type), pointer :: field_ptr !< Pointer to the field's object
966 type(fmsdiagoutputbuffer_type), pointer :: buffer_ptr !< Pointer to the field's buffer
967 class(fmsdiagfilecontainer_type), pointer :: file_ptr !< Pointer to the field's file
968 type(diagyamlfilesvar_type), pointer :: field_yaml_ptr !< Pointer to the field's yaml
969
970 integer :: reduction_method !< Integer representing a reduction method
971 integer :: ids !< For looping through buffer ids
972 integer :: buffer_id !< Id of the buffer
973 integer :: file_id !< File id
974 integer, pointer :: axis_ids(:) !< Axis ids for the buffer
975 logical :: is_subregional !< .True. if the buffer is subregional
976 logical :: reduced_k_range !< .True. is the field is only outputing a section
977 !! of the z dimension
978 type(fmsdiagibounds_type) :: bounds_in !< Starting and ending indices of the input field_data
979 type(fmsdiagibounds_type) :: bounds_out !< Starting and ending indices of the output buffer
980 integer :: i !< For looping through axid ids
981 integer :: sindex !< Starting index of a subregion
982 integer :: eindex !< Ending index of a subregion
983 integer :: compute_idx(2) !< Starting and Ending of the compute domain
984 character(len=1) :: cart_axis !< Cartesian axis of the axis
985 logical :: block_in_subregion !< .True. if the current block is part of the subregion
986 integer :: starting !< Starting index of the subregion relative to the compute domain
987 integer :: ending !< Ending index of the subregion relative to the compute domain
988 real(kind=r8_kind) :: missing_value !< Missing_value for data points that are masked
989 !! This will obtained as r8 and converted to the right type as
990 !! needed. This is to avoid yet another select type ...
991
992 !TODO mostly everything
993 field_ptr => this%FMS_diag_fields(diag_field_id)
994 if (field_ptr%has_missing_value()) then
995 select type (missing_val => field_ptr%get_missing_value(r8))
996 type is (real(kind=r8_kind))
997 missing_value = missing_val
998 class default
999 call mpp_error(fatal, "The missing value for the field:"//trim(field_ptr%get_varname())//&
1000 &" was not allocated to the correct type. This shouldn't have happened")
1001 end select
1002 else
1003 select type (missing_val => get_default_missing_value(r8))
1004 type is (real(kind=r8_kind))
1005 missing_value = missing_val
1006 class default
1007 call mpp_error(fatal, "The missing value for the field:"//trim(field_ptr%get_varname())//&
1008 &" was not allocated to the correct type. This shouldn't have happened")
1009 end select
1010 endif
1011
1012 buffer_loop: do ids = 1, size(field_ptr%buffer_ids)
1013 error_msg = ""
1014 buffer_id = this%FMS_diag_fields(diag_field_id)%buffer_ids(ids)
1015 file_id = this%FMS_diag_fields(diag_field_id)%file_ids(ids)
1016
1017 !< Gather all the objects needed for the buffer
1018 field_yaml_ptr => field_ptr%diag_field(ids)
1019 buffer_ptr => this%FMS_diag_output_buffers(buffer_id)
1020 file_ptr => this%FMS_diag_files(file_id)
1021
1022 !< Go away if the file is a subregional file and the current PE does not have any data for it
1023 if (.not. file_ptr%writing_on_this_pe()) cycle
1024
1025 !< Go away if finished doing math for this buffer
1026 if (buffer_ptr%is_done_with_math()) cycle
1027
1028 if (present(time)) call file_ptr%set_model_time(time)
1029 if (.not. file_ptr%time_to_start_doing_math()) cycle
1030
1031 bounds_out = bounds
1032 if (.not. using_blocking) then
1033 !< Set output bounds to start at 1:size(buffer_ptr%buffer)
1034 call bounds_out%reset_bounds_from_array_4D(buffer_ptr%buffer(:,:,:,:,1))
1035 endif
1036
1037 bounds_in = bounds
1038 if (.not. bounds%has_halos) then
1039 !< If field_data does not contain halos, set bounds_in to start at 1:size(field_data)
1040 call bounds_in%reset_bounds_from_array_4D(field_data)
1041 endif
1042
1043 is_subregional = file_ptr%is_regional()
1044 reduced_k_range = field_yaml_ptr%has_var_zbounds()
1045
1046 !< Reset the bounds based on the reduced k range and subregional
1047 is_subregional_reduced_k_range: if (is_subregional .or. reduced_k_range) then
1048 call buffer_ptr%get_axis_ids(axis_ids)
1049 block_in_subregion = .true.
1050 axis_loops: do i = 1, size(axis_ids)
1051 !< Move on if the block does not have any data for the subregion
1052 if (.not. block_in_subregion) cycle
1053
1054 select type (diag_axis => this%diag_axis(axis_ids(i))%axis)
1055 type is (fmsdiagsubaxis_type)
1056 sindex = diag_axis%get_starting_index()
1057 eindex = diag_axis%get_ending_index()
1058 compute_idx = diag_axis%get_compute_indices()
1059 starting=sindex-compute_idx(1)+1
1060 ending=eindex-compute_idx(1)+1
1061 if (using_blocking) then
1062 block_in_subregion = determine_if_block_is_in_region(starting, ending, bounds, i)
1063 if (.not. block_in_subregion) cycle
1064
1065 !< Set bounds_in so that you can the correct section of the data for the block (starting at 1)
1066 call bounds_in%rebase_input(bounds, starting, ending, i)
1067
1068 !< Set bounds_out to be the correct section relative to the block starting and ending indices
1069 call bounds_out%rebase_output(starting, ending, i)
1070 else
1071 !< Set bounds_in so that only the subregion section of the data will be used (starting at 1)
1072 call bounds_in%update_index(starting, ending, i, .false.)
1073
1074 !< Set bounds_out to 1:size(subregion) for the PE
1075 call bounds_out%update_index(1, ending-starting+1, i, .true.)
1076 endif
1077 end select
1078 enddo axis_loops
1079 nullify(axis_ids)
1080 !< Move on to the next buffer if the block does not have any data for the subregion
1081 if (.not. block_in_subregion) cycle
1082 endif is_subregional_reduced_k_range
1083
1084 !< Determine the reduction method for the buffer
1085 reduction_method = field_yaml_ptr%get_var_reduction()
1086 if (present(time)) call buffer_ptr%update_buffer_time(time)
1087 call buffer_ptr%set_send_data_called()
1088 select case(reduction_method)
1089 case (time_none)
1090 error_msg = buffer_ptr%do_time_none_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
1091 bounds_in, bounds_out, missing_value)
1092 if (trim(error_msg) .ne. "") then
1093 return
1094 endif
1095 case (time_min)
1096 error_msg = buffer_ptr%do_time_min_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
1097 bounds_in, bounds_out, missing_value)
1098 if (trim(error_msg) .ne. "") then
1099 return
1100 endif
1101 case (time_max)
1102 error_msg = buffer_ptr%do_time_max_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
1103 bounds_in, bounds_out, missing_value)
1104 if (trim(error_msg) .ne. "") then
1105 return
1106 endif
1107 case (time_sum)
1108 error_msg = buffer_ptr%do_time_sum_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
1109 field_ptr%get_mask_variant(), bounds_in, bounds_out, missing_value, field_ptr%has_missing_value())
1110 if (trim(error_msg) .ne. "") then
1111 return
1112 endif
1113 case (time_average)
1114 error_msg = buffer_ptr%do_time_sum_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
1115 field_ptr%get_mask_variant(), bounds_in, bounds_out, missing_value, field_ptr%has_missing_value(), &
1116 weight=weight)
1117 if (trim(error_msg) .ne. "") then
1118 return
1119 endif
1120 case (time_power)
1121 error_msg = buffer_ptr%do_time_sum_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
1122 field_ptr%get_mask_variant(), bounds_in, bounds_out, missing_value, field_ptr%has_missing_value(), &
1123 weight=weight, pow_value=field_yaml_ptr%get_pow_value())
1124 if (trim(error_msg) .ne. "") then
1125 return
1126 endif
1127 case (time_rms)
1128 error_msg = buffer_ptr%do_time_sum_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
1129 field_ptr%get_mask_variant(), bounds_in, bounds_out, missing_value, field_ptr%has_missing_value(), &
1130 weight=weight, pow_value = 2)
1131 if (trim(error_msg) .ne. "") then
1132 return
1133 endif
1134 case (time_diurnal)
1135 if(.not. present(time)) call mpp_error(fatal, &
1136 "fms_diag_do_reduction:: time must be present when using diurnal reductions")
1137 ! sets the diurnal index for reduction within the buffer object
1138 call buffer_ptr%set_diurnal_section_index(time)
1139 error_msg = buffer_ptr%do_time_sum_wrapper(field_data, oor_mask, field_ptr%get_var_is_masked(), &
1140 field_ptr%get_mask_variant(), bounds_in, bounds_out, missing_value, field_ptr%has_missing_value(), &
1141 weight=weight)
1142 if (trim(error_msg) .ne. "") then
1143 return
1144 endif
1145 case default
1146 error_msg = "The reduction method is not supported. "//&
1147 "Only none, min, max, sum, average, power, rms, and diurnal are supported."
1148 end select
1149
1150 if (field_ptr%is_static() .or. file_ptr%FMS_diag_file%is_done_writing_data()) then
1151 call buffer_ptr%set_done_with_math()
1152 endif
1153 enddo buffer_loop
1154#else
1155 error_msg = ""
1156 CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
1157#endif
1158end function fms_diag_do_reduction
1159
1160!> @brief Adds the diag ids of the Area and or Volume of the diag_field_object
1161subroutine fms_diag_field_add_cell_measures(this, diag_field_id, area, volume)
1162 class(fmsdiagobject_type), intent (inout) :: this !< The diag object
1163 integer, intent(in) :: diag_field_id !< diag_field to add the are and volume to
1164 INTEGER, optional, INTENT(in) :: area !< diag ids of area
1165 INTEGER, optional, INTENT(in) :: volume !< diag ids of volume
1166
1167#ifndef use_yaml
1168 CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
1169#else
1170 call this%FMS_diag_fields(diag_field_id)%add_area_volume(area, volume)
1171#endif
1172end subroutine fms_diag_field_add_cell_measures
1173
1174!> @brief Add a attribute to the diag_obj using the diag_field_id
1175subroutine fms_diag_field_add_attribute(this, diag_field_id, att_name, att_value)
1176 class(fmsdiagobject_type), intent (inout) :: this !< The diag object
1177 integer, intent(in) :: diag_field_id !< Id of the axis to add the attribute to
1178 character(len=*), intent(in) :: att_name !< Name of the attribute
1179 class(*), intent(in) :: att_value(:) !< The attribute value to add
1180#ifndef use_yaml
1181CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
1182#else
1183!TODO: Value for diag not found
1184 if ( diag_field_id .LE. 0 ) THEN
1185 RETURN
1186 else
1187 if (this%FMS_diag_fields(diag_field_id)%is_registered() ) &
1188 call this%FMS_diag_fields(diag_field_id)%add_attribute(att_name, att_value)
1189 endif
1190#endif
1191end subroutine fms_diag_field_add_attribute
1192
1193!> @brief Add an attribute to an axis
1194subroutine fms_diag_axis_add_attribute(this, axis_id, att_name, att_value)
1195 class(fmsdiagobject_type), intent (inout) :: this !< The diag object
1196 integer, intent(in) :: axis_id !< Id of the axis to add the attribute to
1197 character(len=*), intent(in) :: att_name !< Name of the attribute
1198 class(*), intent(in) :: att_value(:) !< The attribute value to add
1199
1200 character(len=20) :: axis_names(2) !< Names of the uncompress axis
1201 character(len=20) :: set_name !< Name of the axis set
1202 integer :: uncmx_ids(2) !< Ids of the uncompress axis
1203 integer :: j !< For do loops
1204#ifndef use_yaml
1205CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
1206#else
1207 if (axis_id < 0 .and. axis_id > this%registered_axis) &
1208 call mpp_error(fatal, "diag_axis_add_attribute: The axis_id is not valid")
1209
1210 select type (axis => this%diag_axis(axis_id)%axis)
1211 type is (fmsdiagfullaxis_type)
1212 call axis%add_axis_attribute(att_name, att_value)
1213
1214 !! Axis that are in the "unstructured" domain require a "compress" attribute for the
1215 !! combiner and PP. This attribute is passed in via a diag_axis_add_attribute call in the model code
1216 !! The compress attribute indicates the names of the axis that were compressed
1217 !! For example grid_index:compress = "grid_yt grid_xt"
1218 !! The metadata and the data for these axis also needs to be written to the file
1219 if (trim(att_name) .eq. "compress") then
1220 !< If the attribute is the "compress" attribute, get the axis names,
1221 !! and the ids of the axis and add it to the axis object so it can be written to netcdf files
1222 !! that use this axis
1223 axis_names = parse_compress_att(att_value)
1224 set_name = ""
1225 if (axis%has_set_name()) set_name = axis%get_set_name()
1226 do j = 1, size(axis_names)
1227 uncmx_ids(j) = get_axis_id_from_name(axis_names(j), this%diag_axis, this%registered_axis, set_name)
1228 if (uncmx_ids(j) .eq. diag_null) call mpp_error(fatal, &
1229 &"Error parsing the compress attribute for axis: "//trim(axis%get_axis_name())//&
1230 &". Be sure that the axes in the compress attribute are registered")
1231 enddo
1232 call axis%add_structured_axis_ids(uncmx_ids)
1233 endif
1234 end select
1235#endif
1236end subroutine fms_diag_axis_add_attribute
1237
1238!> \brief Gets the field_name from the diag_field
1239!> \returns a copy of the field_name
1240function fms_get_field_name_from_id (this, field_id) &
1241 result(field_name)
1242
1243 class(fmsdiagobject_type), intent (in) :: this !< The diag object, the caller
1244 integer, intent (in) :: field_id !< Field id to get the name for
1245 character(len=:), allocatable :: field_name
1246#ifndef use_yaml
1247 CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
1248#else
1249 field_name = this%FMS_diag_fields(field_id)%get_varname()
1250#endif
1251end function fms_get_field_name_from_id
1252
1253!> \brief Gets the diag field ID from the module name and field name.
1254!> \returns a copy of the ID of the diag field or DIAG_FIELD_NOT_FOUND if the field is not registered
1255FUNCTION fms_get_diag_field_id_from_name(this, module_name, field_name) &
1256 result(diag_field_id)
1257 class(fmsdiagobject_type), intent (in) :: this !< The diag object, the caller
1258 CHARACTER(len=*), INTENT(in) :: module_name !< Module name that registered the variable
1259 CHARACTER(len=*), INTENT(in) :: field_name !< Variable name
1260 integer :: diag_field_id
1261
1262#ifdef use_yaml
1263 integer :: i !< For looping
1264 integer, allocatable :: diag_field_indices(:) !< indices where the field was found in the yaml
1265
1266 diag_field_id = diag_field_not_found
1267
1268 !> Loop through fields to find it.
1269 do i=1, this%registered_variables
1270 !< Check if the field was registered, if it was return the diag_field_id
1271 diag_field_id = this%FMS_diag_fields(i)%id_from_name(module_name, field_name)
1272 if(diag_field_id .ne. diag_field_not_found) return
1273 enddo
1274
1275 !< Check if the field is in the diag_table.yaml. If it is, return DIAG_FIELD_NOT_REGISTERED
1276 !! Otherwsie it will return DIAG_FIELD_NOT_FOUND
1277 diag_field_indices = find_diag_field(field_name, module_name)
1278 if (diag_field_indices(1) .ne. diag_null) then
1279 diag_field_id = diag_not_registered
1280 endif
1281 deallocate(diag_field_indices)
1282#else
1283 diag_field_id = diag_field_not_found
1284 CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
1285#endif
1286END FUNCTION fms_get_diag_field_id_from_name
1287
1288#ifdef use_yaml
1289!> returns the buffer object for the given id
1290!! actual data comes from %get_buffer_data() on the returned object
1291function get_diag_buffer(this, bufferid) &
1292result(rslt)
1293 class(fmsdiagobject_type), intent(in) :: this
1294 integer, intent(in) :: bufferid
1295 type(fmsdiagoutputbuffer_type),allocatable:: rslt
1296 if( (bufferid .gt. ubound(this%FMS_diag_output_buffers, 1)) .or. &
1297 (bufferid .lt. lbound(this%FMS_diag_output_buffers, 1))) &
1298 call mpp_error(fatal, 'get_diag_bufer: invalid bufferid given')
1299 rslt = this%FMS_diag_output_buffers(bufferid)
1300end function
1301#endif
1302
1303!> @brief Return the 2D domain for the axis IDs given.
1304!! @return 2D domain for the axis IDs given
1305type(domain2d) function fms_get_domain2d(this, ids)
1306 class(fmsdiagobject_type), intent (in) :: this !< The diag object
1307 INTEGER, DIMENSION(:), INTENT(in) :: ids !< Axis IDs.
1308
1309#ifndef use_yaml
1310CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
1311fms_get_domain2d = null_domain2d
1312#else
1313 INTEGER :: type_of_domain !< The type of domain
1314 CLASS(diagdomain_t), POINTER :: domain !< Diag Domain pointer
1315
1316 call get_domain_and_domain_type(fms_diag_object%diag_axis, ids, type_of_domain, domain, "get_domain2d")
1317 if (type_of_domain .ne. two_d_domain) &
1318 call mpp_error(fatal, 'diag_axis_mod::get_domain2d- The axis do not correspond to a 2d Domain')
1319 select type(domain)
1320 type is (diagdomain2d_t)
1321 fms_get_domain2d = domain%domain2
1322 end select
1323#endif
1324END FUNCTION fms_get_domain2d
1325
1326 !> @brief Gets the length of the axis based on the axis_id
1327 !> @return Axis_length
1328 integer function fms_get_axis_length(this, axis_id)
1329 class(fmsdiagobject_type), intent (in) :: this !< The diag object
1330 INTEGER, INTENT(in) :: axis_id !< Axis ID of the axis to the length of
1331
1332#ifndef use_yaml
1333CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
1334fms_get_axis_length = 0
1335#else
1336fms_get_axis_length = 0
1337
1338 if (axis_id < 0 .and. axis_id > this%registered_axis) &
1339 call mpp_error(fatal, "fms_get_axis_length: The axis_id is not valid")
1340
1341 select type (axis => this%diag_axis(axis_id)%axis)
1342 type is (fmsdiagfullaxis_type)
1343 fms_get_axis_length = axis%axis_length()
1344 type is (fmsdiagsubaxis_type)
1345 fms_get_axis_length = axis%axis_length()
1346 end select
1347#endif
1348end function fms_get_axis_length
1349
1350!> @brief Gets the name of the axis based on the axis_id
1351 !> @return The axis_name
1352function fms_get_axis_name_from_id (this, axis_id) &
1353result(axis_name)
1354 class(fmsdiagobject_type), intent (in) :: this !< The diag object
1355 INTEGER, INTENT(in) :: axis_id !< Axis ID of the axis to the length of
1356
1357 character (len=:), allocatable :: axis_name
1358
1359#ifndef use_yaml
1360CALL mpp_error(fatal,"You can not use the modern diag manager without compiling with -Duse_yaml")
1361axis_name=" "
1362#else
1363 if (axis_id < 0 .and. axis_id > this%registered_axis) &
1364 call mpp_error(fatal, "fms_get_axis_length: The axis_id is not valid")
1365
1366 !! if its a scalar (null axis id) just returns n/a since no axis is defined
1367 if (axis_id .eq. null_axis_id) then
1368 allocate(character(len=3) :: axis_name)
1369 axis_name = "n/a"
1370 return
1371 endif
1372
1373
1374 select type (axis => this%diag_axis(axis_id)%axis)
1375 type is (fmsdiagfullaxis_type)
1376 axis_name = axis%get_axis_name()
1377 end select
1378#endif
1379end function fms_get_axis_name_from_id
1380
1381!> Dumps as much data as it can from the fmsDiagObject_type.
1382!! Will dump any fields and files as well (see d)
1383subroutine dump_diag_obj( filename )
1384 character(len=*), intent(in), optional :: filename !< optional filename to print to,
1385 !! otherwise prints to stdout
1386#ifdef use_yaml
1387 !type(fmsDiagObject_type) :: diag_obj
1388 type(fmsdiagfile_type), pointer :: fileptr !< pointer for traversing file list
1389 type(fmsdiagfield_type), pointer :: fieldptr !< pointer for traversing field list
1390 integer :: i !< do loops
1391 integer :: unit_num !< unit num of opened log file or stdout
1392
1393 if( present(filename) ) then
1394 open(newunit=unit_num, file=trim(filename), action='WRITE')
1395 else
1396 unit_num = stdout()
1397 endif
1398 if( mpp_pe() .eq. mpp_root_pe()) then
1399 write(unit_num, *) '********** dumping diag object ***********'
1400 write(unit_num, *) 'registered_variables:', fms_diag_object%registered_variables
1401 write(unit_num, *) 'registered_axis:', fms_diag_object%registered_axis
1402 write(unit_num, *) 'initialized:', fms_diag_object%initialized
1403 write(unit_num, *) 'files_initialized:', fms_diag_object%files_initialized
1404 write(unit_num, *) 'fields_initialized:', fms_diag_object%fields_initialized
1405 write(unit_num, *) 'buffers_initialized:', fms_diag_object%buffers_initialized
1406 write(unit_num, *) 'axes_initialized:', fms_diag_object%axes_initialized
1407 write(unit_num, *) 'Files:'
1408 if( fms_diag_object%files_initialized ) then
1409 do i=1, SIZE(fms_diag_object%FMS_diag_files)
1410 write(unit_num, *) 'File num:', i
1411 fileptr => fms_diag_object%FMS_diag_files(i)%FMS_diag_file
1412 call fileptr%dump_file_obj(unit_num)
1413 enddo
1414 else
1415 write(unit_num, *) 'files not initialized'
1416 endif
1417 if( fms_diag_object%fields_initialized) then
1418 do i=1, SIZE(fms_diag_object%FMS_diag_fields)
1419 write(unit_num, *) 'Field num:', i
1420 fieldptr => fms_diag_object%FMS_diag_fields(i)
1421 call fieldptr%dump_field_obj(unit_num)
1422 enddo
1423 else
1424 write(unit_num, *) 'fields not initialized'
1425 endif
1426 if( present(filename) ) close(unit_num)
1427 endif
1428#else
1429 call mpp_error( fatal, "You can not use the modern diag manager without compiling with -Duse_yaml")
1430#endif
1431end subroutine
1432
1433!> @brief Allocates the output buffers of the fields corresponding to the registered variable
1434!! Input arguments are the field and its ID passed to routine fms_diag_accept_data()
1435subroutine allocate_diag_field_output_buffers(this, field_data, field_id)
1436 class(fmsdiagobject_type), target, intent(inout) :: this !< diag object
1437 class(*), dimension(:,:,:,:), intent(in) :: field_data !< field data
1438 integer, intent(in) :: field_id !< Id of the field data
1439#ifdef use_yaml
1440 integer :: ndims !< Number of dimensions in the input field data
1441 integer :: buffer_id !< Buffer index of FMS_diag_buffers
1442 integer :: num_diurnal_samples !< Number of diurnal samples from diag_yaml
1443 integer :: axes_length(4) !< Length of each axis
1444 integer :: i, j !< For looping
1445 class(fmsdiagoutputbuffer_type), pointer :: ptr_diag_buffer_obj !< Pointer to the buffer class
1446 class(diagyamlfilesvar_type), pointer :: ptr_diag_field_yaml !< Pointer to a field from yaml fields
1447 integer, pointer :: axis_ids(:) !< Pointer to indices of axes of the field variable
1448 integer :: var_type !< Stores type of the field data (r4, r8, i4, i8, and string) represented as an integer.
1449 character(len=:), allocatable :: var_name !< Field name to initialize output buffers
1450 logical :: is_scalar !< Flag indicating that the variable is a scalar
1451 integer :: yaml_id !< Yaml id for the buffer
1452 integer :: file_id !< File id for the buffer
1453
1454 if (this%FMS_diag_fields(field_id)%buffer_allocated) return
1455
1456 ! Determine the type of the field data
1457 var_type = get_var_type(field_data(1, 1, 1, 1))
1458
1459 ! Get variable/field name
1460 var_name = this%FMS_diag_fields(field_id)%get_varname()
1461
1462 ! Determine dimensions of the field
1463 is_scalar = this%FMS_diag_fields(field_id)%is_scalar()
1464
1465 ! Loop over a number of fields/buffers where this variable occurs
1466 do i = 1, size(this%FMS_diag_fields(field_id)%buffer_ids)
1467 buffer_id = this%FMS_diag_fields(field_id)%buffer_ids(i)
1468 file_id = this%FMS_diag_fields(field_id)%file_ids(i)
1469
1470 !< Go away if the file is a subregional file and the current PE does not have any data for it
1471 if (.not. this%FMS_diag_files(file_id)%writing_on_this_pe()) cycle
1472
1473 ndims = 0
1474 if (.not. is_scalar) then
1475 call this%FMS_diag_output_buffers(buffer_id)%get_axis_ids(axis_ids)
1476 ndims = size(axis_ids)
1477 endif
1478
1479 yaml_id = this%FMS_diag_output_buffers(buffer_id)%get_yaml_id()
1480
1481 ptr_diag_field_yaml => diag_yaml%diag_fields(yaml_id)
1482 num_diurnal_samples = ptr_diag_field_yaml%get_n_diurnal() !< Get number of diurnal samples
1483
1484 axes_length = 1
1485 do j = 1, ndims
1486 axes_length(j) = this%fms_get_axis_length(axis_ids(j))
1487 enddo
1488
1489 if (num_diurnal_samples .ne. 0) then
1490 ndims = ndims + 1 !< Add one more dimension for the diurnal axis
1491 endif
1492
1493 ptr_diag_buffer_obj => this%FMS_diag_output_buffers(buffer_id)
1494 call ptr_diag_buffer_obj%allocate_buffer(field_data(1, 1, 1, 1), ndims, axes_length(1:4), &
1495 this%FMS_diag_fields(field_id)%get_mask_variant(), var_name, num_diurnal_samples)
1496 call ptr_diag_buffer_obj%initialize_buffer(ptr_diag_field_yaml%get_var_reduction(), var_name)
1497
1498 enddo
1499 nullify(axis_ids)
1500
1501 this%FMS_diag_fields(field_id)%buffer_allocated = .true.
1502#else
1503 call mpp_error( fatal, "allocate_diag_field_output_buffers: "//&
1504 "you can not use the modern diag manager without compiling with -Duse_yaml")
1505#endif
1506end subroutine allocate_diag_field_output_buffers
1507
1508!> @brief Determines if the window defined by the input bounds is a physics window.
1509!> @return TRUE if the window size is less then the actual field size else FALSE.
1510function fms_diag_compare_window(this, field, field_id, &
1511 is_in, ie_in, js_in, je_in, ks_in, ke_in) result(is_phys_win)
1512 class(fmsdiagobject_type), intent(in) :: this !< Diag Object
1513 class(*), intent(in) :: field(:,:,:,:) !< Field data
1514 integer, intent(in) :: field_id !< ID of the input field
1515 integer, intent(in) :: is_in, js_in !< Starting field indices for the first 2 dimensions;
1516 !< pass reconditioned indices fis and fjs
1517 !< which are computed elsewhere.
1518 integer, intent(in) :: ie_in, je_in !< Ending field indices for the first 2 dimensions;
1519 !< pass reconditioned indices fie and fje
1520 !< which are computed elsewhere.
1521 integer, intent(in) :: ks_in, ke_in !< Starting and ending indices of the field in 3rd dimension
1522 logical :: is_phys_win !< Return flag
1523#ifdef use_yaml
1524 integer, pointer :: axis_ids(:)
1525 integer :: total_elements
1526 integer :: i !< For do loop
1527 integer :: field_size
1528 integer, allocatable :: field_shape(:) !< Shape of the field data
1529 integer :: window_size
1530
1531 !> Determine shape of the field defined by the input bounds
1532 field_shape = shape(field(is_in:ie_in, js_in:je_in, ks_in:ke_in, :))
1533
1534 window_size = field_shape(1) * field_shape(2) * field_shape(3)
1535
1536 total_elements = 1
1537 axis_ids => this%FMS_diag_fields(field_id)%get_axis_id()
1538 do i=1, size(axis_ids)
1539 total_elements = total_elements * this%fms_get_axis_length(axis_ids(i))
1540 enddo
1541
1542 if (total_elements > window_size) then
1543 is_phys_win = .true.
1544 else
1545 is_phys_win = .false.
1546 end if
1547#else
1548 is_phys_win = .false.
1549 call mpp_error( fatal, "fms_diag_compare_window: "//&
1550 "you can not use the modern diag manager without compiling with -Duse_yaml")
1551#endif
1552end function fms_diag_compare_window
1553
1554!> @brief Set the model_end_time in a diag object
1555subroutine set_time_end(this, time_end_in)
1556 class(fmsdiagobject_type), intent(inout) :: this !< Diag Object
1557 type(time_type), intent(in) :: time_end_in !< Time at the end of the simulation
1558#ifdef use_yaml
1559 this%model_end_time = time_end_in
1560#endif
1561end subroutine
1562
1563end module fms_diag_object_mod
integer, parameter no_domain
Use the FmsNetcdfFile_t fileobj.
integer, parameter end_time
Use the end of the time average bounds.
integer max_axes
Maximum number of independent axes.
integer, parameter diag_field_not_found
Return value for a diag_field that isn't found in the diag_table.
integer, parameter string
s is the 19th letter of the alphabet
Definition diag_data.F90:88
integer function get_var_type(var)
gets the type of a variable
type(time_type) diag_init_time
Time diag_manager_init called. If init_time not included in diag_manager_init call,...
integer, parameter time_min
The reduction method is min value.
integer, parameter time_diurnal
The reduction method is diurnal.
integer, parameter time_power
The reduction method is average with exponents.
type(time_type) function get_base_time()
gets the module variable base_time
integer, parameter time_average
The reduction method is average of values.
integer, parameter time_sum
The reduction method is sum of values.
integer, parameter time_rms
The reudction method is root mean square of values.
integer, parameter time_none
There is no reduction method.
integer, parameter time_max
The reduction method is max value.
integer, parameter r8
Supported type/kind of the variable.
Definition diag_data.F90:84
integer, parameter two_d_domain
Use the FmsNetcdfDomainFile_t fileobj.
logical pure function, public determine_if_block_is_in_region(subregion_start, subregion_end, bounds, dim)
The PEs grid points are divided further into "blocks". This function determines if a block.
logical function, dimension(:,:,:,:), allocatable, public init_mask(rmask, mask, field)
Sets the logical mask based on mask or rmask.
pure real(kind=r8_kind) function, public set_weight(weight)
Sets the weight based on the weight passed into send_data (1.0_r8_kind if the weight is not passed in...
pure character(len=128) function, public check_indices_order(is_in, ie_in, js_in, je_in)
Checks improper combinations of is, ie, js, and je.
logical function, public fms_diag_axis_object_end(axis_array)
subroutine, public diag_yaml_object_end()
Destroys the diag_yaml object.
pure character(len=120) function, dimension(2), public parse_compress_att(compress_att)
integer function, dimension(:), allocatable, public get_diag_field_ids(indices)
Gets field indices corresponding to the indices (input argument) in the sorted variable_list.
integer function, dimension(:), allocatable, public find_diag_field(diag_field_name, module_name)
Determines if a diag_field is in the diag_yaml_object.
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.
integer function axis_length(this)
Get the axis length of a subaxis.
subroutine, public get_domain_and_domain_type(diag_axis, axis_id, domain_type, domain, var_name)
Loop through a variable's axis_id to determine and return the domain type and domain to use.
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, 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 ...
logical function, public fms_diag_axis_object_init(axis_array)
pure integer function, public get_axis_id_from_name(axis_name, diag_axis, naxis, set_name)
One dimensional domain used to manage shared data access between pes.
The domain2D type contains all the necessary information to define the global, compute and data domai...
Domain information for managing data on unstructured grids.
Error handler.
Definition mpp.F90:382
integer function, public get_ticks_per_second()
Returns the number of ticks per second.
subroutine, public get_time(time, seconds, days, ticks, err_msg)
Returns days and seconds ( < 86400 ) corresponding to a time. err_msg should be checked for any error...
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.
subroutine, public get_date(time, year, month, day, hour, minute, second, tick, err_msg)
Gets the date for different calendar types. Given a time_interval, returns the corresponding date und...
Given an input date in year, month, days, etc., creates a time_type that represents this time interva...
Given some number of seconds and days, returns the corresponding time_type.
Type to represent amounts of time. Implemented as seconds and days to allow for larger intervals.
Contains buffer types and routines for the diag manager.
Type to hold the domain info for an axis This type was created to avoid having to send in "Domain",...
Type to hold the diagnostic axis description.
Type to hold the diag_axis (either subaxis or a full axis)
Type to hold the diagnostic axis description.
Data structure holding a 3D bounding box. It is commonlyused to represent the interval bounds or limi...
Object that holds all variable information.
A container for fmsDiagFile_type. This is used to create the array of files.
type to hold the info a diag_field