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