FMS 2025.01.02-dev
Flexible Modeling System
Loading...
Searching...
No Matches
fms_diag_field_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_field_object_mod
20!> \author Tom Robinson
21!> \email thomas.robinson@noaa.gov
22!! \brief Contains routines for the diag_objects
23!!
24!! \description The diag_manager passes an object back and forth between the diag routines and the users.
25!! The procedures of this object and the types are all in this module. The fms_dag_object is a type
26!! that contains all of the information of the variable. It is extended by a type that holds the
27!! appropriate buffer for the data for manipulation.
28#ifdef use_yaml
29use diag_data_mod, only: prepend_date, diag_null, cmor_missing_value, diag_null_string, max_str_len
30use diag_data_mod, only: r8, r4, i8, i4, string, null_type_int, no_domain
31use diag_data_mod, only: max_field_attributes, fmsdiagattribute_type
32use diag_data_mod, only: diag_null, diag_not_found, diag_not_registered, diag_registered_id, &
35use fms_string_utils_mod, only: int2str=>string
36use mpp_mod, only: fatal, note, warning, mpp_error, mpp_pe, mpp_root_pe
39use fms_diag_axis_object_mod, only: diagdomain_t, get_domain_and_domain_type, fmsdiagaxis_type, &
41use time_manager_mod, ONLY: time_type, get_date
42use fms2_io_mod, only: fmsnetcdffile_t, fmsnetcdfdomainfile_t, fmsnetcdfunstructureddomainfile_t, register_field, &
43 register_variable_attribute
44use fms_diag_input_buffer_mod, only: fmsdiaginputbuffer_t
45!!!set_time, set_date, get_time, time_type, OPERATOR(>=), OPERATOR(>),&
46!!! & OPERATOR(<), OPERATOR(==), OPERATOR(/=), OPERATOR(/), OPERATOR(+), ASSIGNMENT(=), get_date, &
47!!! & get_ticks_per_second
48
49use platform_mod
50use iso_c_binding
51
52implicit none
53
54private
55
56!> \brief Object that holds all variable information
58 type (diagyamlfilesvar_type), allocatable, dimension(:) :: diag_field !< info from diag_table for this variable
59 integer, allocatable, dimension(:) :: file_ids !< Ids of the FMS_diag_files the variable
60 !! belongs to
61 integer, allocatable, private :: diag_id !< unique id for varable
62 integer, allocatable, dimension(:) :: buffer_ids !< index/id for this field's buffers
63 type(fmsdiagattribute_type), allocatable :: attributes(:) !< attributes for the variable
64 integer, private :: num_attributes !< Number of attributes currently added
65 logical, allocatable, private :: static !< true if this is a static var
66 logical, allocatable, private :: scalar !< .True. if the variable is a scalar
67 logical, allocatable, private :: registered !< true when registered
68 logical, allocatable, private :: mask_variant !< true if the mask changes over time
69 logical, allocatable, private :: var_is_masked !< true if the field is masked
70 logical, allocatable, private :: do_not_log !< .true. if no need to log the diag_field
71 logical, allocatable, private :: local !< If the output is local
72 integer, allocatable, private :: vartype !< the type of varaible
73 character(len=:), allocatable, private :: varname !< the name of the variable
74 character(len=:), allocatable, private :: longname !< longname of the variable
75 character(len=:), allocatable, private :: standname !< standard name of the variable
76 character(len=:), allocatable, private :: units !< the units
77 character(len=:), allocatable, private :: modname !< the module
78 character(len=:), allocatable, private :: realm !< String to set as the value
79 !! to the modeling_realm attribute
80 character(len=:), allocatable, private :: interp_method !< The interp method to be used
81 !! when regridding the field in post-processing.
82 !! Valid options are "conserve_order1",
83 !! "conserve_order2", and "none".
84 integer, allocatable, dimension(:), private :: frequency !< specifies the frequency
85 integer, allocatable, private :: tile_count !< The number of tiles
86 integer, allocatable, dimension(:), private :: axis_ids !< variable axis IDs
87 class(diagdomain_t), pointer, private :: domain !< Domain
88 INTEGER , private :: type_of_domain !< The type of domain ("NO_DOMAIN",
89 !! "TWO_D_DOMAIN", or "UG_DOMAIN")
90 integer, allocatable, private :: area, volume !< The Area and Volume
91 class(*), allocatable, private :: missing_value !< The missing fill value
92 class(*), allocatable, private :: data_range(:) !< The range of the variable data
93 type(fmsdiaginputbuffer_t), allocatable :: input_data_buffer !< Input buffer object for when buffering
94 !! data
95 logical, allocatable, private :: multiple_send_data!< .True. if send_data is called multiple
96 !! times for the same model time
97 logical, allocatable, private :: data_buffer_is_allocated !< True if the buffer has
98 !! been allocated
99 logical, allocatable, private :: math_needs_to_be_done !< If true, do math
100 !! functions. False when done.
101 logical, allocatable :: buffer_allocated !< True if a buffer pointed by
102 !! the corresponding index in
103 !! buffer_ids(:) is allocated.
104 logical, allocatable :: mask(:,:,:,:) !< Mask passed in send_data
105 logical :: halo_present = .false. !< set if any halos are used
106 contains
107! procedure :: send_data => fms_send_data !!TODO
108! Get ID functions
109 procedure :: get_id => fms_diag_get_id
110 procedure :: id_from_name => diag_field_id_from_name
111 procedure :: copy => copy_diag_obj
112 procedure :: register => fms_register_diag_field_obj !! Merely initialize fields.
113 procedure :: setid => set_diag_id
114 procedure :: set_type => set_vartype
115 procedure :: set_data_buffer => set_data_buffer
116 procedure :: prepare_data_buffer
117 procedure :: init_data_buffer
118 procedure :: set_data_buffer_is_allocated
121 procedure :: is_data_buffer_allocated
122 procedure :: allocate_data_buffer
123 procedure :: set_math_needs_to_be_done => set_math_needs_to_be_done
124 procedure :: add_attribute => diag_field_add_attribute
125 procedure :: vartype_inq => what_is_vartype
126 procedure :: set_var_is_masked
127 procedure :: get_var_is_masked
128! Check functions
129 procedure :: is_static => diag_obj_is_static
130 procedure :: is_scalar
131 procedure :: is_registered => get_registered
132 procedure :: is_registeredb => diag_obj_is_registered
133 procedure :: is_mask_variant => get_mask_variant
134 procedure :: is_local => get_local
135! Is variable allocated check functions
136!TODO procedure :: has_diag_field
137 procedure :: has_diag_id
138 procedure :: has_attributes
139 procedure :: has_static
140 procedure :: has_registered
141 procedure :: has_mask_variant
142 procedure :: has_local
143 procedure :: has_vartype
144 procedure :: has_varname
145 procedure :: has_longname
146 procedure :: has_standname
147 procedure :: has_units
148 procedure :: has_modname
149 procedure :: has_realm
150 procedure :: has_interp_method
151 procedure :: has_frequency
152 procedure :: has_tile_count
153 procedure :: has_axis_ids
154 procedure :: has_area
155 procedure :: has_volume
156 procedure :: has_missing_value
157 procedure :: has_data_range
158 procedure :: has_input_data_buffer
159! Get functions
160 procedure :: get_attributes
161 procedure :: get_static
162 procedure :: get_registered
163 procedure :: get_mask_variant
164 procedure :: get_local
165 procedure :: get_vartype
166 procedure :: get_varname
167 procedure :: get_longname
168 procedure :: get_standname
169 procedure :: get_units
170 procedure :: get_modname
171 procedure :: get_realm
172 procedure :: get_interp_method
173 procedure :: get_frequency
174 procedure :: get_tile_count
175 procedure :: get_area
176 procedure :: get_volume
177 procedure :: get_missing_value
178 procedure :: get_data_range
179 procedure :: get_axis_id
180 procedure :: get_data_buffer
181 procedure :: get_mask
182 procedure :: get_weight
183 procedure :: dump_field_obj
184 procedure :: get_domain
185 procedure :: get_type_of_domain
186 procedure :: set_file_ids
187 procedure :: get_dimnames
188 procedure :: get_var_skind
189 procedure :: get_longname_to_write
190 procedure :: get_multiple_send_data
191 procedure :: write_field_metadata
192 procedure :: write_coordinate_attribute
193 procedure :: get_math_needs_to_be_done
194 procedure :: add_area_volume
195 procedure :: append_time_cell_methods
196 procedure :: get_file_ids
197 procedure :: set_mask
198 procedure :: allocate_mask
199 procedure :: set_halo_present
200 procedure :: is_halo_present
201 procedure :: find_missing_value
202 procedure :: has_mask_allocated
203 procedure :: is_variable_in_file
204 procedure :: get_field_file_name
205 procedure :: generate_associated_files_att
206end type fmsdiagfield_type
207!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! variables !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
208type(fmsdiagfield_type) :: null_ob
209
210logical,private :: module_is_initialized = .false. !< Flag indicating if the module is initialized
211
212!type(fmsDiagField_type) :: diag_object_placeholder (10)
213!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
214public :: fmsdiagfield_type
215public :: fms_diag_fields_object_init
216public :: null_ob
217public :: fms_diag_field_object_end
218public :: get_default_missing_value
219public :: check_for_slices
220!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
221!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
222 CONTAINS
223!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
224!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
225
226!> @brief Deallocates the array of diag_objs
227subroutine fms_diag_field_object_end (ob)
228 class(fmsdiagfield_type), allocatable, intent(inout) :: ob(:) !< diag field object
229 if (allocated(ob)) deallocate(ob)
230 module_is_initialized = .false.
231end subroutine fms_diag_field_object_end
232!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
233!> \Description Allocates the diad field object array.
234!! Sets the diag_id to the not registered value.
235!! Initializes the number of registered variables to be 0
236logical function fms_diag_fields_object_init(ob)
237 type(fmsdiagfield_type), allocatable, intent(inout) :: ob(:) !< diag field object
238 integer :: i !< For looping
239 allocate(ob(get_num_unique_fields()))
240 do i = 1,size(ob)
241 ob(i)%diag_id = diag_not_registered !null_ob%diag_id
242 ob(i)%registered = .false.
243 enddo
244 module_is_initialized = .true.
245 fms_diag_fields_object_init = .true.
246end function fms_diag_fields_object_init
247!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
248!> \Description Fills in and allocates (when necessary) the values in the diagnostic object
249subroutine fms_register_diag_field_obj &
250 (this, modname, varname, diag_field_indices, diag_axis, axes, &
251 longname, units, missing_value, varrange, mask_variant, standname, &
252 do_not_log, err_msg, interp_method, tile_count, area, volume, realm, static, &
253 multiple_send_data)
254
255 class(fmsdiagfield_type), INTENT(inout) :: this !< Diaj_obj to fill
256 CHARACTER(len=*), INTENT(in) :: modname !< The module name
257 CHARACTER(len=*), INTENT(in) :: varname !< The variable name
258 integer, INTENT(in) :: diag_field_indices(:) !< Array of indices to the field
259 !! in the yaml object
260 class(fmsdiagaxiscontainer_type),intent(in) :: diag_axis(:) !< Array of diag_axis
261 INTEGER, TARGET, OPTIONAL, INTENT(in) :: axes(:) !< The axes indicies
262 CHARACTER(len=*), OPTIONAL, INTENT(in) :: longname !< THe variables long name
263 CHARACTER(len=*), OPTIONAL, INTENT(in) :: units !< The units of the variables
264 CHARACTER(len=*), OPTIONAL, INTENT(in) :: standname !< The variables stanard name
265 class(*), OPTIONAL, INTENT(in) :: missing_value !< Missing value to add as a attribute
266 class(*), OPTIONAL, INTENT(in) :: varrange(2) !< Range to add as a attribute
267 LOGICAL, OPTIONAL, INTENT(in) :: mask_variant !< Mask
268 LOGICAL, OPTIONAL, INTENT(in) :: do_not_log !< if TRUE, field info is not logged
269 CHARACTER(len=*), OPTIONAL, INTENT(out) :: err_msg !< Error message to be passed back up
270 CHARACTER(len=*), OPTIONAL, INTENT(in) :: interp_method !< The interp method to be used when
271 !! regridding the field in post-processing.
272 !! Valid options are "conserve_order1",
273 !! "conserve_order2", and "none".
274 INTEGER, OPTIONAL, INTENT(in) :: tile_count !< the number of tiles
275 INTEGER, OPTIONAL, INTENT(in) :: area !< diag_field_id of the cell area field
276 INTEGER, OPTIONAL, INTENT(in) :: volume !< diag_field_id of the cell volume field
277 CHARACTER(len=*), OPTIONAL, INTENT(in) :: realm !< String to set as the value to the
278 !! modeling_realm attribute
279 LOGICAL, OPTIONAL, INTENT(in) :: static !< Set to true if it is a static field
280 LOGICAL, OPTIONAL, INTENT(in) :: multiple_send_data !< .True. if send data is called, multiple
281 !! times for the same time
282 integer :: i, j !< for looponig over field/axes indices
283 character(len=:), allocatable, target :: a_name_tmp !< axis name tmp
284 type(diagyamlfilesvar_type), pointer :: yaml_var_ptr !< pointer this fields yaml variable entries
285
286!> Fill in information from the register call
287 this%varname = trim(varname)
288 this%modname = trim(modname)
289
290!> Add the yaml info to the diag_object
291 this%diag_field = get_diag_fields_entries(diag_field_indices)
292
293 if (present(static)) then
294 this%static = static
295 else
296 this%static = .false.
297 endif
298
299!> Add axis and domain information
300 if (present(axes)) then
301
302 this%scalar = .false.
303 this%axis_ids = axes
304 call get_domain_and_domain_type(diag_axis, this%axis_ids, this%type_of_domain, this%domain, this%varname)
305
306 ! store dim names for output
307 ! cant use this%diag_field since they are copies
308 do i=1, SIZE(diag_field_indices)
309 yaml_var_ptr => diag_yaml%get_diag_field_from_id(diag_field_indices(i))
310 ! add dim names from axes
311 do j=1, SIZE(axes)
312 a_name_tmp = diag_axis(axes(j))%axis%get_axis_name( yaml_var_ptr%is_file_subregional())
313 if(yaml_var_ptr%has_var_zbounds() .and. a_name_tmp .eq. 'z') &
314 a_name_tmp = trim(a_name_tmp)//"_sub01"
315 call yaml_var_ptr%add_axis_name(a_name_tmp)
316 enddo
317 ! add time_of_day_N dimension if diurnal
318 if(yaml_var_ptr%has_n_diurnal()) then
319 a_name_tmp = "time_of_day_"// int2str(yaml_var_ptr%get_n_diurnal())
320 call yaml_var_ptr%add_axis_name(a_name_tmp)
321 endif
322 ! add time dimension if not static
323 if(.not. this%static) then
324 a_name_tmp = "time"
325 call yaml_var_ptr%add_axis_name(a_name_tmp)
326 endif
327 enddo
328 else
329 !> The variable is a scalar
330 this%scalar = .true.
331 this%type_of_domain = no_domain
332 this%domain => null()
333 ! store dim name for output (just the time if not static and no axes)
334 if(.not. this%static) then
335 do i=1, SIZE(diag_field_indices)
336 a_name_tmp = "time"
337 yaml_var_ptr => diag_yaml%get_diag_field_from_id(diag_field_indices(i))
338 call yaml_var_ptr%add_axis_name(a_name_tmp)
339 enddo
340 endif
341 endif
342 nullify(yaml_var_ptr)
343
344!> get the optional arguments if included and the diagnostic is in the diag table
345 if (present(longname)) this%longname = trim(longname)
346 if (present(standname)) this%standname = trim(standname)
347
348 !> Ignore the units if they are set to "none". This is to reproduce previous diag_manager behavior
349 if (present(units)) then
350 if (trim(units) .ne. "none") this%units = trim(units)
351 endif
352 if (present(realm)) this%realm = trim(realm)
353 if (present(interp_method)) this%interp_method = trim(interp_method)
354
355 if (present(tile_count)) then
356 allocate(this%tile_count)
357 this%tile_count = tile_count
358 endif
359
360 if (present(missing_value)) then
361 select type (missing_value)
362 type is (integer(kind=i4_kind))
363 allocate(integer(kind=i4_kind) :: this%missing_value)
364 this%missing_value = missing_value
365 type is (integer(kind=i8_kind))
366 allocate(integer(kind=i8_kind) :: this%missing_value)
367 this%missing_value = missing_value
368 type is (real(kind=r4_kind))
369 allocate(real(kind=r4_kind) :: this%missing_value)
370 this%missing_value = missing_value
371 type is (real(kind=r8_kind))
372 allocate(real(kind=r8_kind) :: this%missing_value)
373 this%missing_value = missing_value
374 class default
375 call mpp_error("fms_register_diag_field_obj", &
376 "The missing value passed to register a diagnostic is not a r8, r4, i8, or i4",&
377 fatal)
378 end select
379 endif
380
381 if (present(varrange)) then
382 select type (varrange)
383 type is (integer(kind=i4_kind))
384 allocate(integer(kind=i4_kind) :: this%data_range(2))
385 this%data_RANGE = varrange
386 type is (integer(kind=i8_kind))
387 allocate(integer(kind=i8_kind) :: this%data_range(2))
388 this%data_RANGE = varrange
389 type is (real(kind=r4_kind))
390 allocate(integer(kind=r4_kind) :: this%data_range(2))
391 this%data_RANGE = varrange
392 type is (real(kind=r8_kind))
393 allocate(integer(kind=r8_kind) :: this%data_range(2))
394 this%data_RANGE = varrange
395 class default
396 call mpp_error("fms_register_diag_field_obj", &
397 "The varRange passed to register a diagnostic is not a r8, r4, i8, or i4",&
398 fatal)
399 end select
400 endif
401
402 if (present(area)) then
403 if (area < 0) call mpp_error("fms_register_diag_field_obj", &
404 "The area id passed with field_name"//trim(varname)//" has not been registered. &
405 &Check that there is a register_diag_field call for the AREA measure and that is in the &
406 &diag_table.yaml", fatal)
407 allocate(this%area)
408 this%area = area
409 endif
410
411 if (present(volume)) then
412 if (volume < 0) call mpp_error("fms_register_diag_field_obj", &
413 "The volume id passed with field_name"//trim(varname)//" has not been registered. &
414 &Check that there is a register_diag_field call for the VOLUME measure and that is in the &
415 &diag_table.yaml", fatal)
416 allocate(this%volume)
417 this%volume = volume
418 endif
419
420 this%mask_variant = .false.
421 if (present(mask_variant)) then
422 this%mask_variant = mask_variant
423 endif
424
425 if (present(do_not_log)) then
426 allocate(this%do_not_log)
427 this%do_not_log = do_not_log
428 endif
429
430 if (present(multiple_send_data)) then
431 this%multiple_send_data = multiple_send_data
432 else
433 this%multiple_send_data = .false.
434 endif
435
436 !< Allocate space for any additional variable attributes
437 !< These will be fill out when calling `diag_field_add_attribute`
438 allocate(this%attributes(max_field_attributes))
439 this%num_attributes = 0
440 this%registered = .true.
441end subroutine fms_register_diag_field_obj
442!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
443
444!> \brief Sets the diag_id. This can only be done if a variable is unregistered
445subroutine set_diag_id(this , id)
446 class(fmsdiagfield_type) , intent(inout):: this
447 integer :: id
448 if (allocated(this%registered)) then
449 if (this%registered) then
450 call mpp_error("set_diag_id", "The variable"//this%varname//" is already registered", fatal)
451 else
452 this%diag_id = id
453 endif
454 else
455 this%diag_id = id
456 endif
457end subroutine set_diag_id
458
459!> \brief Find the type of the variable and store it in the object
460subroutine set_vartype(objin , var)
461 class(fmsdiagfield_type) , intent(inout):: objin
462 class(*) :: var
463 select type (var)
464 type is (real(kind=8))
465 objin%vartype = r8
466 type is (real(kind=4))
467 objin%vartype = r4
468 type is (integer(kind=8))
469 objin%vartype = i8
470 type is (integer(kind=4))
471 objin%vartype = i4
472 type is (character(*))
473 objin%vartype = string
474 class default
475 objin%vartype = null_type_int
476 call mpp_error("set_vartype", "The variable"//objin%varname//" is not a supported type "// &
477 " r8, r4, i8, i4, or string.", warning)
478 end select
479end subroutine set_vartype
480
481!> @brief Sets the time send data was called last
482subroutine set_send_data_time (this, time)
483 class(fmsdiagfield_type) , intent(inout):: this !< The field object
484 type(time_type), intent(in) :: time !< Current model time
485
486 call this%input_data_buffer%set_send_data_time(time)
487end subroutine set_send_data_time
488
489!> @brief Get the time send data was called last
490!! @result the time send data was called last
491function get_send_data_time(this) &
492 result(rslt)
493 class(fmsdiagfield_type) , intent(in):: this !< The field object
494 type(time_type) :: rslt
495
496 rslt = this%input_data_buffer%get_send_data_time()
497end function get_send_data_time
498
499!> @brief Prepare the input_data_buffer to do the reduction method
500subroutine prepare_data_buffer(this)
501 class(fmsdiagfield_type) , intent(inout):: this !< The field object
502
503 if (.not. this%multiple_send_data) return
504 if (this%mask_variant) return
505 call this%input_data_buffer%prepare_input_buffer_object(this%modname//":"//this%varname)
506end subroutine prepare_data_buffer
507
508!> @brief Initialize the input_data_buffer
509subroutine init_data_buffer(this)
510 class(fmsdiagfield_type) , intent(inout):: this !< The field object
511
512 if (.not. this%multiple_send_data) return
513 if (this%mask_variant) return
514 call this%input_data_buffer%init_input_buffer_object()
515end subroutine init_data_buffer
516
517!> @brief Adds the input data to the buffered data.
518subroutine set_data_buffer (this, input_data, mask, weight, is, js, ks, ie, je, ke)
519 class(fmsdiagfield_type) , intent(inout):: this !< The field object
520 class(*), intent(in) :: input_data(:,:,:,:) !< The input array
521 logical, intent(in) :: mask(:,:,:,:) !< Mask that is passed into
522 !! send_data
523 real(kind=r8_kind), intent(in) :: weight !< The field weight
524 integer, intent(in) :: is, js, ks !< Starting indicies of the field_data relative
525 !! to the compute domain (1 based)
526 integer, intent(in) :: ie, je, ke !< Ending indicies of the field_data relative
527 !! to the compute domain (1 based)
528
529 character(len=128) :: err_msg !< Error msg
530 if (.not.this%data_buffer_is_allocated) &
531 call mpp_error ("set_data_buffer", "The data buffer for the field "//trim(this%varname)//" was unable to be "//&
532 "allocated.", fatal)
533 if (this%multiple_send_data) then
534 err_msg = this%input_data_buffer%update_input_buffer_object(input_data, is, js, ks, ie, je, ke, &
535 mask, this%mask, this%mask_variant, this%var_is_masked)
536 else
537 this%mask(is:ie, js:je, ks:ke, :) = mask
538 err_msg = this%input_data_buffer%set_input_buffer_object(input_data, weight, is, js, ks, ie, je, ke)
539 endif
540 if (trim(err_msg) .ne. "") call mpp_error(fatal, "Field:"//trim(this%varname)//" -"//trim(err_msg))
541
542end subroutine set_data_buffer
543!> Allocates the global data buffer for a given field using a single thread. Returns true when the
544!! buffer is allocated
545logical function allocate_data_buffer(this, input_data, diag_axis)
546 class(fmsdiagfield_type), target, intent(inout):: this !< The field object
547 class(*), dimension(:,:,:,:), intent(in) :: input_data !< The input array
548 class(fmsdiagaxiscontainer_type),intent(in) :: diag_axis(:) !< Array of diag_axis
549
550 character(len=128) :: err_msg !< Error msg
551 err_msg = ""
552
553 allocate(this%input_data_buffer)
554 err_msg = this%input_data_buffer%allocate_input_buffer_object(input_data, this%axis_ids, diag_axis)
555 if (trim(err_msg) .ne. "") then
556 call mpp_error(fatal, "Field:"//trim(this%varname)//" -"//trim(err_msg))
557 return
558 endif
559
560 allocate_data_buffer = .true.
561end function allocate_data_buffer
562!> Sets the flag saying that the math functions need to be done
563subroutine set_math_needs_to_be_done (this, math_needs_to_be_done)
564 class(fmsdiagfield_type) , intent(inout):: this
565 logical, intent (in) :: math_needs_to_be_done !< Flag saying that the math functions need to be done
566 this%math_needs_to_be_done = math_needs_to_be_done
567end subroutine set_math_needs_to_be_done
568
569!> @brief Set the mask_variant to .true.
570subroutine set_var_is_masked(this, is_masked)
571 class(fmsdiagfield_type) , intent(inout):: this !< The diag field object
572 logical, intent (in) :: is_masked !< .True. if the field is masked
573
574 this%var_is_masked = is_masked
575end subroutine set_var_is_masked
576
577!> @brief Queries a field for the var_is_masked variable
578!! @return var_is_masked
579function get_var_is_masked(this) &
580 result(rslt)
581 class(fmsdiagfield_type) , intent(inout):: this !< The diag field object
582 logical :: rslt !< .True. if the field is masked
583
584 rslt = this%var_is_masked
585end function get_var_is_masked
586
587!> @brief Sets the flag saying that the data buffer is allocated
588subroutine set_data_buffer_is_allocated (this, data_buffer_is_allocated)
589 class(fmsdiagfield_type) , intent(inout) :: this !< The field object
590 logical, intent (in) :: data_buffer_is_allocated !< .true. if the
591 !! data buffer is allocated
592 this%data_buffer_is_allocated = data_buffer_is_allocated
593end subroutine set_data_buffer_is_allocated
594
595!> @brief Determine if the data_buffer is allocated
596!! @return logical indicating if the data_buffer is allocated
597pure logical function is_data_buffer_allocated (this)
598 class(fmsdiagfield_type) , intent(in) :: this !< The field object
599
600 is_data_buffer_allocated = .false.
601 if (allocated(this%data_buffer_is_allocated)) is_data_buffer_allocated = this%data_buffer_is_allocated
602
603end function
604!> \brief Prints to the screen what type the diag variable is
605subroutine what_is_vartype(this)
606 class(fmsdiagfield_type) , intent(inout):: this
607 if (.not. allocated(this%vartype)) then
608 call mpp_error("what_is_vartype", "The variable type has not been set prior to this call", warning)
609 return
610 endif
611 select case (this%vartype)
612 case (r8)
613 call mpp_error("what_is_vartype", "The variable type of "//trim(this%varname)//&
614 " is REAL(kind=8)", note)
615 case (r4)
616 call mpp_error("what_is_vartype", "The variable type of "//trim(this%varname)//&
617 " is REAL(kind=4)", note)
618 case (i8)
619 call mpp_error("what_is_vartype", "The variable type of "//trim(this%varname)//&
620 " is INTEGER(kind=8)", note)
621 case (i4)
622 call mpp_error("what_is_vartype", "The variable type of "//trim(this%varname)//&
623 " is INTEGER(kind=4)", note)
624 case (string)
625 call mpp_error("what_is_vartype", "The variable type of "//trim(this%varname)//&
626 " is CHARACTER(*)", note)
627 case (null_type_int)
628 call mpp_error("what_is_vartype", "The variable type of "//trim(this%varname)//&
629 " was not set", warning)
630 case default
631 call mpp_error("what_is_vartype", "The variable type of "//trim(this%varname)//&
632 " is not supported by diag_manager", fatal)
633 end select
634end subroutine what_is_vartype
635!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
636
637!> \brief Copies the calling object into the object that is the argument of the subroutine
638subroutine copy_diag_obj(this , objout)
639 class(fmsdiagfield_type) , intent(in) :: this
640 class(fmsdiagfield_type) , intent(inout) , allocatable :: objout !< The destination of the copy
641select type (objout)
642 class is (fmsdiagfield_type)
643
644 if (allocated(this%registered)) then
645 objout%registered = this%registered
646 else
647 call mpp_error("copy_diag_obj", "You can only copy objects that have been registered",warning)
648 endif
649 objout%diag_id = this%diag_id
650
651 if (allocated(this%attributes)) objout%attributes = this%attributes
652 objout%static = this%static
653 if (allocated(this%frequency)) objout%frequency = this%frequency
654 if (allocated(this%varname)) objout%varname = this%varname
655end select
656end subroutine copy_diag_obj
657!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
658
659!> \brief Returns the ID integer for a variable
660!! \return the diag ID
661pure integer function fms_diag_get_id (this) result(diag_id)
662 class(fmsdiagfield_type) , intent(in) :: this
663!> Check if the diag_object registration has been done
664 if (allocated(this%registered)) then
665 !> Return the diag_id if the variable has been registered
666 diag_id = this%diag_id
667 else
668!> If the variable is not regitered, then return the unregistered value
669 diag_id = diag_not_registered
670 endif
671end function fms_diag_get_id
672
673!> Function to return a character (string) representation of the most basic
674!> object identity info. Intended for debugging and warning. The format produced is:
675!> [this: o.varname(string|?), vartype (string|?), o.registered (T|F|?), diag_id (id|?)].
676!> A questionmark "?" is set in place of the variable that is not yet allocated
677!>TODO: Add diag_id ?
678function fms_diag_obj_as_string_basic(this) result(rslt)
679 class(fmsdiagfield_type), allocatable, intent(in) :: this
680 character(:), allocatable :: rslt
681 character (len=:), allocatable :: registered, vartype, varname, diag_id
682 if ( .not. allocated (this)) then
683 varname = "?"
684 vartype = "?"
685 registered = "?"
686 diag_id = "?"
687 rslt = "[Obj:" // varname // "," // vartype // "," // registered // "," // diag_id // "]"
688 return
689 end if
690
691! if(allocated (this%registered)) then
692! registered = logical_to_cs (this%registered)
693! else
694! registered = "?"
695! end if
696
697! if(allocated (this%diag_id)) then
698! diag_id = int_to_cs (this%diag_id)
699! else
700! diag_id = "?"
701! end if
702
703! if(allocated (this%vartype)) then
704! vartype = int_to_cs (this%vartype)
705! else
706! registered = "?"
707! end if
708
709 if(allocated (this%varname)) then
710 varname = this%varname
711 else
712 registered = "?"
713 end if
714
715 rslt = "[Obj:" // varname // "," // vartype // "," // registered // "," // diag_id // "]"
716
717end function fms_diag_obj_as_string_basic
718
719
720function diag_obj_is_registered (this) result (rslt)
721 class(fmsdiagfield_type), intent(in) :: this
722 logical :: rslt
723 rslt = this%registered
724end function diag_obj_is_registered
725
726function diag_obj_is_static (this) result (rslt)
727 class(fmsdiagfield_type), intent(in) :: this
728 logical :: rslt
729 rslt = .false.
730 if (allocated(this%static)) rslt = this%static
731end function diag_obj_is_static
732
733!> @brief Determine if the field is a scalar
734!! @return .True. if the field is a scalar
735function is_scalar (this) result (rslt)
736 class(fmsdiagfield_type), intent(in) :: this !< diag_field object
737 logical :: rslt
738 rslt = this%scalar
739end function is_scalar
740
741!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
742!! Get functions
743
744!> @brief Gets attributes
745!! @return A pointer to the attributes of the diag_obj, null pointer if there are no attributes
746function get_attributes (this) &
747result(rslt)
748 class(fmsdiagfield_type), target, intent(in) :: this !< diag object
749 type(fmsdiagattribute_type), pointer :: rslt(:)
750
751 rslt => null()
752 if (this%num_attributes > 0 ) rslt => this%attributes
753end function get_attributes
754
755!> @brief Gets static
756!! @return copy of variable static
757pure function get_static (this) &
758result(rslt)
759 class(fmsdiagfield_type), intent(in) :: this !< diag object
760 logical :: rslt
761 rslt = this%static
762end function get_static
763
764!> @brief Gets regisetered
765!! @return copy of registered
766pure function get_registered (this) &
767result(rslt)
768 class(fmsdiagfield_type), intent(in) :: this !< diag object
769 logical :: rslt
770 rslt = this%registered
771end function get_registered
772
773!> @brief Gets mask variant
774!! @return copy of mask variant
775pure function get_mask_variant (this) &
776result(rslt)
777 class(fmsdiagfield_type), intent(in) :: this !< diag object
778 logical :: rslt
779 rslt = .false.
780 if (allocated(this%mask_variant)) rslt = this%mask_variant
781end function get_mask_variant
782
783!> @brief Gets local
784!! @return copy of local
785pure function get_local (this) &
786result(rslt)
787 class(fmsdiagfield_type), intent(in) :: this !< diag object
788 logical :: rslt
789 rslt = this%local
790end function get_local
791
792!> @brief Gets vartype
793!! @return copy of The integer related to the variable type
794pure function get_vartype (this) &
795result(rslt)
796 class(fmsdiagfield_type), intent(in) :: this !< diag object
797 integer :: rslt
798 rslt = this%vartype
799end function get_vartype
800
801!> @brief Gets varname
802!! @return copy of the variable name
803pure function get_varname (this, to_write) &
804result(rslt)
805 class(fmsdiagfield_type), intent(in) :: this !< diag object
806 logical, optional, intent(in) :: to_write !< .true. if getting the varname that will be writen to the file
807 character(len=:), allocatable :: rslt
808 rslt = this%varname
809
810 !< If writing the varname can be the outname which is defined in the yaml
811 if (present(to_write)) then
812 if (to_write) then
813 !TODO this is wrong
814 rslt = this%diag_field(1)%get_var_outname()
815 endif
816 endif
817
818end function get_varname
819
820!> @brief Gets longname
821!! @return copy of the variable long name or a single string if there is no long name
822pure function get_longname (this) &
823result(rslt)
824 class(fmsdiagfield_type), intent(in) :: this !< diag object
825 character(len=:), allocatable :: rslt
826 if (allocated(this%longname)) then
827 rslt = this%longname
828 else
829 rslt = diag_null_string
830 endif
831end function get_longname
832
833!> @brief Gets standname
834!! @return copy of the standard name or an empty string if standname is not allocated
835pure function get_standname (this) &
836result(rslt)
837 class(fmsdiagfield_type), intent(in) :: this !< diag object
838 character(len=:), allocatable :: rslt
839 if (allocated(this%standname)) then
840 rslt = this%standname
841 else
842 rslt = diag_null_string
843 endif
844end function get_standname
845
846!> @brief Gets units
847!! @return copy of the units or an empty string if not allocated
848pure function get_units (this) &
849result(rslt)
850 class(fmsdiagfield_type), intent(in) :: this !< diag object
851 character(len=:), allocatable :: rslt
852 if (allocated(this%units)) then
853 rslt = this%units
854 else
855 rslt = diag_null_string
856 endif
857end function get_units
858
859!> @brief Gets modname
860!! @return copy of the module name that the variable is in or an empty string if not allocated
861pure function get_modname (this) &
862result(rslt)
863 class(fmsdiagfield_type), intent(in) :: this !< diag object
864 character(len=:), allocatable :: rslt
865 if (allocated(this%modname)) then
866 rslt = this%modname
867 else
868 rslt = diag_null_string
869 endif
870end function get_modname
871
872!> @brief Gets realm
873!! @return copy of the variables modeling realm or an empty string if not allocated
874pure function get_realm (this) &
875result(rslt)
876 class(fmsdiagfield_type), intent(in) :: this !< diag object
877 character(len=:), allocatable :: rslt
878 if (allocated(this%realm)) then
879 rslt = this%realm
880 else
881 rslt = diag_null_string
882 endif
883end function get_realm
884
885!> @brief Gets interp_method
886!! @return copy of The interpolation method or an empty string if not allocated
887pure function get_interp_method (this) &
888result(rslt)
889 class(fmsdiagfield_type), intent(in) :: this !< diag object
890 character(len=:), allocatable :: rslt
891 if (allocated(this%interp_method)) then
892 rslt = this%interp_method
893 else
894 rslt = diag_null_string
895 endif
896end function get_interp_method
897
898!> @brief Gets frequency
899!! @return copy of the frequency or DIAG_NULL if obj%frequency is not allocated
900pure function get_frequency (this) &
901result(rslt)
902 class(fmsdiagfield_type), intent(in) :: this !< diag object
903 integer, allocatable, dimension (:) :: rslt
904 if (allocated(this%frequency)) then
905 allocate (rslt(size(this%frequency)))
906 rslt = this%frequency
907 else
908 allocate (rslt(1))
909 rslt = diag_null
910 endif
911end function get_frequency
912
913!> @brief Gets tile_count
914!! @return copy of the number of tiles or diag_null if tile_count is not allocated
915pure function get_tile_count (this) &
916result(rslt)
917 class(fmsdiagfield_type), intent(in) :: this !< diag object
918 integer :: rslt
919 if (allocated(this%tile_count)) then
920 rslt = this%tile_count
921 else
922 rslt = diag_null
923 endif
924end function get_tile_count
925
926!> @brief Gets area
927!! @return copy of the area or diag_null if not allocated
928pure function get_area (this) &
929result(rslt)
930 class(fmsdiagfield_type), intent(in) :: this !< diag object
931 integer :: rslt
932 if (allocated(this%area)) then
933 rslt = this%area
934 else
935 rslt = diag_null
936 endif
937end function get_area
938
939!> @brief Gets volume
940!! @return copy of the volume or diag_null if volume is not allocated
941pure function get_volume (this) &
942result(rslt)
943 class(fmsdiagfield_type), intent(in) :: this !< diag object
944 integer :: rslt
945 if (allocated(this%volume)) then
946 rslt = this%volume
947 else
948 rslt = diag_null
949 endif
950end function get_volume
951
952!> @brief Gets missing_value
953!! @return copy of The missing value
954!! @note Netcdf requires the type of the variable and the type of the missing_value and _Fillvalue to be the same
955!! var_type is the type of the variable which may not be in the same type as the missing_value in the register call
956!! For example, if compiling with r8 but the in diag_table.yaml the kind is r4
957function get_missing_value (this, var_type) &
958result(rslt)
959 class(fmsdiagfield_type), intent(in) :: this !< diag object
960 integer, intent(in) :: var_type !< The type of the variable as it will writen to the netcdf file
961 !! and the missing value is return as
962
963 class(*),allocatable :: rslt
964
965 if (.not. allocated(this%missing_value)) then
966 call mpp_error ("get_missing_value", &
967 "The missing value is not allocated", fatal)
968 endif
969
970 !< The select types are needed so that the missing_value can be correctly converted and copied as the needed variable
971 !! type
972 select case (var_type)
973 case (r4)
974 allocate (real(kind=r4_kind) :: rslt)
975 select type (miss => this%missing_value)
976 type is (real(kind=r4_kind))
977 select type (rslt)
978 type is (real(kind=r4_kind))
979 rslt = real(miss, kind=r4_kind)
980 end select
981 type is (real(kind=r8_kind))
982 select type (rslt)
983 type is (real(kind=r4_kind))
984 rslt = real(miss, kind=r4_kind)
985 end select
986 end select
987 case (r8)
988 allocate (real(kind=r8_kind) :: rslt)
989 select type (miss => this%missing_value)
990 type is (real(kind=r4_kind))
991 select type (rslt)
992 type is (real(kind=r8_kind))
993 rslt = real(miss, kind=r8_kind)
994 end select
995 type is (real(kind=r8_kind))
996 select type (rslt)
997 type is (real(kind=r8_kind))
998 rslt = real(miss, kind=r8_kind)
999 end select
1000 end select
1001 end select
1002
1003end function get_missing_value
1004
1005!> @brief Gets data_range
1006!! @return copy of the data range
1007!! @note Netcdf requires the type of the variable and the type of the range to be the same
1008!! var_type is the type of the variable which may not be in the same type as the range in the register call
1009!! For example, if compiling with r8 but the in diag_table.yaml the kind is r4
1010function get_data_range (this, var_type) &
1011result(rslt)
1012 class(fmsdiagfield_type), intent(in) :: this !< diag object
1013 integer, intent(in) :: var_type !< The type of the variable as it will writen to the netcdf file
1014 !! and the data_range is returned as
1015 class(*),allocatable :: rslt(:)
1016
1017 if ( .not. allocated(this%data_RANGE)) call mpp_error ("get_data_RANGE", &
1018 "The data_RANGE value is not allocated", fatal)
1019
1020 !< The select types are needed so that the range can be correctly converted and copied as the needed variable
1021 !! type
1022 select case (var_type)
1023 case (r4)
1024 allocate (real(kind=r4_kind) :: rslt(2))
1025 select type (r => this%data_RANGE)
1026 type is (real(kind=r4_kind))
1027 select type (rslt)
1028 type is (real(kind=r4_kind))
1029 rslt = real(r, kind=r4_kind)
1030 end select
1031 type is (real(kind=r8_kind))
1032 select type (rslt)
1033 type is (real(kind=r4_kind))
1034 rslt = real(r, kind=r4_kind)
1035 end select
1036 end select
1037 case (r8)
1038 allocate (real(kind=r8_kind) :: rslt(2))
1039 select type (r => this%data_RANGE)
1040 type is (real(kind=r4_kind))
1041 select type (rslt)
1042 type is (real(kind=r8_kind))
1043 rslt = real(r, kind=r8_kind)
1044 end select
1045 type is (real(kind=r8_kind))
1046 select type (rslt)
1047 type is (real(kind=r8_kind))
1048 rslt = real(r, kind=r8_kind)
1049 end select
1050 end select
1051 end select
1052end function get_data_range
1053
1054!> @brief Gets axis_ids
1055!! @return pointer to the axis ids
1056function get_axis_id (this) &
1057result(rslt)
1058 class(fmsdiagfield_type), target, intent(in) :: this !< diag object
1059 integer, pointer, dimension(:) :: rslt !< field's axis_ids
1060
1061 if(allocated(this%axis_ids)) then
1062 rslt => this%axis_ids
1063 else
1064 rslt => null()
1065 endif
1066end function get_axis_id
1067
1068!> @brief Gets field's domain
1069!! @return pointer to the domain
1070function get_domain (this) &
1071result(rslt)
1072 class(fmsdiagfield_type), target, intent(in) :: this !< diag field
1073 class(diagdomain_t), pointer :: rslt !< field's domain
1074
1075 if (associated(this%domain)) then
1076 rslt => this%domain
1077 else
1078 rslt => null()
1079 endif
1080
1081end function get_domain
1082
1083!> @brief Gets field's type of domain
1084!! @return integer defining the type of domain (NO_DOMAIN, TWO_D_DOMAIN, UG_DOMAIN)
1085pure function get_type_of_domain (this) &
1086result(rslt)
1087 class(fmsdiagfield_type), target, intent(in) :: this !< diag field
1088 integer :: rslt !< field's domain
1089
1090 rslt = this%type_of_domain
1091end function get_type_of_domain
1092
1093!> @brief Set the file ids of the files that the field belongs to
1094subroutine set_file_ids(this, file_ids)
1095 class(fmsdiagfield_type), intent(inout) :: this !< diag field
1096 integer, intent(in) :: file_ids(:) !< File_ids to add
1097
1098 allocate(this%file_ids(size(file_ids)))
1099 this%file_ids = file_ids
1100end subroutine set_file_ids
1101
1102!> @brief Get the kind of the variable based on the yaml
1103!! @return A string indicating the kind of the variable (as it is used in fms2_io)
1104pure function get_var_skind(this, field_yaml) &
1105result(rslt)
1106 class(fmsdiagfield_type), intent(in) :: this !< diag field
1107 type(diagyamlfilesvar_type), intent(in) :: field_yaml !< The corresponding yaml of the field
1108
1109 character(len=:), allocatable :: rslt
1110
1111 integer :: var_kind !< The integer corresponding to the kind of the variable (i4, i8, r4, r8)
1112
1113 var_kind = field_yaml%get_var_kind()
1114 select case (var_kind)
1115 case (r4)
1116 rslt = "float"
1117 case (r8)
1118 rslt = "double"
1119 case (i4)
1120 rslt = "int"
1121 case (i8)
1122 rslt = "int64"
1123 end select
1124
1125end function get_var_skind
1126
1127!> @brief Get the multiple_send_data member of the field object
1128!! @return multiple_send_data of the field
1129pure function get_multiple_send_data(this) &
1130result(rslt)
1131 class(fmsdiagfield_type), intent(in) :: this !< diag field
1132 logical :: rslt
1133 rslt = this%multiple_send_data
1134end function get_multiple_send_data
1135
1136!> @brief Determine the long name to write for the field
1137!! @return Long name to write
1138pure function get_longname_to_write(this, field_yaml) &
1139result(rslt)
1140 class(fmsdiagfield_type), intent(in) :: this !< diag field
1141 type(diagyamlfilesvar_type), intent(in) :: field_yaml !< The corresponding yaml of the field
1142
1143 character(len=:), allocatable :: rslt
1144
1145 rslt = field_yaml%get_var_longname() !! This is the long name defined in the yaml
1146 if (rslt .eq. "") then !! If the long name is not defined in the yaml, use the long name in the
1147 !! register_diag_field
1148 rslt = this%get_longname()
1149 else
1150 return
1151 endif
1152 if (rslt .eq. "") then !! If the long name is not defined in the yaml and in the register_diag_field
1153 !! use the variable name
1154 rslt = field_yaml%get_var_varname()
1155 endif
1156end function get_longname_to_write
1157
1158!> @brief Determine the dimension names to use when registering the field to fms2_io
1159subroutine get_dimnames(this, diag_axis, field_yaml, unlim_dimname, dimnames, is_regional)
1160 class(fmsdiagfield_type), target, intent(inout) :: this !< diag field
1161 class(fmsdiagaxiscontainer_type), target, intent(in) :: diag_axis(:) !< Diag_axis object
1162 type(diagyamlfilesvar_type), intent(in) :: field_yaml !< Field info from diag_table yaml
1163 character(len=*), intent(in) :: unlim_dimname !< The name of unlimited dimension
1164 character(len=120), allocatable, intent(out) :: dimnames(:) !< Array of the dimension names
1165 !! for the field
1166 logical, intent(in) :: is_regional !< Flag indicating if the field is regional
1167
1168 integer :: i !< For do loops
1169 integer :: naxis !< Number of axis for the field
1170 class(fmsdiagaxiscontainer_type), pointer :: axis_ptr !diag_axis(this%axis_ids(i), for convenience
1171 character(len=23) :: diurnal_axis_name !< name of the diurnal axis
1172
1173 if (this%is_static()) then
1174 naxis = size(this%axis_ids)
1175 else
1176 naxis = size(this%axis_ids) + 1 !< Adding 1 more dimension for the unlimited dimension
1177 endif
1178
1179 if (field_yaml%has_n_diurnal()) then
1180 naxis = naxis + 1 !< Adding 1 more dimension for the diurnal axis
1181 endif
1182
1183 allocate(dimnames(naxis))
1184
1185 !< Duplicated do loops for performance
1186 if (field_yaml%has_var_zbounds()) then
1187 do i = 1, size(this%axis_ids)
1188 axis_ptr => diag_axis(this%axis_ids(i))
1189 if (axis_ptr%axis%is_z_axis()) then
1190 dimnames(i) = axis_ptr%axis%get_axis_name(is_regional)//"_sub01"
1191 else
1192 dimnames(i) = axis_ptr%axis%get_axis_name(is_regional)
1193 endif
1194 enddo
1195 else
1196 do i = 1, size(this%axis_ids)
1197 axis_ptr => diag_axis(this%axis_ids(i))
1198 dimnames(i) = axis_ptr%axis%get_axis_name(is_regional)
1199 enddo
1200 endif
1201
1202 !< The second to last dimension is always the diurnal axis
1203 if (field_yaml%has_n_diurnal()) then
1204 WRITE (diurnal_axis_name,'(a,i2.2)') 'time_of_day_', field_yaml%get_n_diurnal()
1205 dimnames(naxis - 1) = trim(diurnal_axis_name)
1206 endif
1207
1208 !< The last dimension is always the unlimited dimensions
1209 if (.not. this%is_static()) dimnames(naxis) = unlim_dimname
1210
1211end subroutine get_dimnames
1212
1213!> @brief Wrapper for the register_field call. The select types are needed so that the code can go
1214!! in the correct interface
1215subroutine register_field_wrap(fms2io_fileobj, varname, vartype, dimensions)
1216 class(fmsnetcdffile_t), INTENT(INOUT) :: fms2io_fileobj!< Fms2_io fileobj to write to
1217 character(len=*), INTENT(IN) :: varname !< Name of the variable
1218 character(len=*), INTENT(IN) :: vartype !< The type of the variable
1219 character(len=*), optional, INTENT(IN) :: dimensions(:) !< The dimension names of the field
1220
1221 select type(fms2io_fileobj)
1222 type is (fmsnetcdffile_t)
1223 call register_field(fms2io_fileobj, varname, vartype, dimensions)
1224 type is (fmsnetcdfdomainfile_t)
1225 call register_field(fms2io_fileobj, varname, vartype, dimensions)
1226 type is (fmsnetcdfunstructureddomainfile_t)
1227 call register_field(fms2io_fileobj, varname, vartype, dimensions)
1228 end select
1229end subroutine register_field_wrap
1230
1231!> @brief Write the field's metadata to the file
1232subroutine write_field_metadata(this, fms2io_fileobj, file_id, yaml_id, diag_axis, unlim_dimname, is_regional, &
1233 cell_measures)
1234 class(fmsdiagfield_type), target, intent(inout) :: this !< diag field
1235 class(fmsnetcdffile_t), INTENT(INOUT) :: fms2io_fileobj!< Fms2_io fileobj to write to
1236 integer, intent(in) :: file_id !< File id of the file to write to
1237 integer, intent(in) :: yaml_id !< Yaml id of the yaml entry of this field
1238 class(fmsdiagaxiscontainer_type), intent(in) :: diag_axis(:) !< Diag_axis object
1239 character(len=*), intent(in) :: unlim_dimname !< The name of the unlimited dimension
1240 logical, intent(in) :: is_regional !< Flag indicating if the field is regional
1241 character(len=*), intent(in) :: cell_measures !< The cell measures attribute to write
1242
1243 type(diagyamlfilesvar_type), pointer :: field_yaml !< pointer to the yaml entry
1244 character(len=:), allocatable :: var_name !< Variable name
1245 character(len=:), allocatable :: long_name !< Longname to write
1246 character(len=:), allocatable :: units !< Units of the field to write
1247 character(len=120), allocatable :: dimnames(:) !< Dimension names of the field
1248 character(len=120) :: cell_methods!< Cell methods attribute to write
1249 integer :: i !< For do loops
1250 character (len=MAX_STR_LEN), allocatable :: yaml_field_attributes(:,:) !< Variable attributes defined in the yaml
1251 character(len=:), allocatable :: interp_method_tmp !< temp to hold the name of the interpolation method
1252 integer :: interp_method_len !< length of the above string
1253
1254 field_yaml => diag_yaml%get_diag_field_from_id(yaml_id)
1255 var_name = field_yaml%get_var_outname()
1256
1257 if (allocated(this%axis_ids)) then
1258 call this%get_dimnames(diag_axis, field_yaml, unlim_dimname, dimnames, is_regional)
1259 call register_field_wrap(fms2io_fileobj, var_name, this%get_var_skind(field_yaml), dimnames)
1260 else
1261 if (this%is_static()) then
1262 call register_field_wrap(fms2io_fileobj, var_name, this%get_var_skind(field_yaml))
1263 else
1264 !< In this case, the scalar variable is a function of time, so we need to pass in the
1265 !! unlimited dimension as a dimension
1266 call register_field_wrap(fms2io_fileobj, var_name, this%get_var_skind(field_yaml), (/unlim_dimname/))
1267 endif
1268 endif
1269
1270 long_name = this%get_longname_to_write(field_yaml)
1271 call register_variable_attribute(fms2io_fileobj, var_name, "long_name", long_name, str_len=len_trim(long_name))
1272
1273 units = this%get_units()
1274 if (units .ne. diag_null_string) &
1275 call register_variable_attribute(fms2io_fileobj, var_name, "units", units, str_len=len_trim(units))
1276
1277 if (this%has_missing_value()) then
1278 call register_variable_attribute(fms2io_fileobj, var_name, "missing_value", &
1279 this%get_missing_value(field_yaml%get_var_kind()))
1280 call register_variable_attribute(fms2io_fileobj, var_name, "_FillValue", &
1281 this%get_missing_value(field_yaml%get_var_kind()))
1282 else
1283 call register_variable_attribute(fms2io_fileobj, var_name, "missing_value", &
1284 get_default_missing_value(field_yaml%get_var_kind()))
1285 call register_variable_attribute(fms2io_fileobj, var_name, "_FillValue", &
1286 get_default_missing_value(field_yaml%get_var_kind()))
1287 endif
1288
1289 if (this%has_data_RANGE()) then
1290 call register_variable_attribute(fms2io_fileobj, var_name, "valid_range", &
1291 this%get_data_range(field_yaml%get_var_kind()))
1292 endif
1293
1294 if (this%has_interp_method()) then
1295 interp_method_tmp = this%interp_method
1296 interp_method_len = len_trim(interp_method_tmp)
1297 call register_variable_attribute(fms2io_fileobj, var_name, "interp_method", interp_method_tmp, &
1298 str_len=interp_method_len)
1299 endif
1300
1301 cell_methods = ""
1302 !< Check if any of the attributes defined via a "diag_field_add_attribute" call
1303 !! are the cell_methods, if so add to the "cell_methods" variable:
1304 do i = 1, this%num_attributes
1305 call this%attributes(i)%write_metadata(fms2io_fileobj, var_name, &
1306 cell_methods=cell_methods)
1307 enddo
1308
1309 !< Append the time cell methods based on the variable's reduction
1310 call this%append_time_cell_methods(cell_methods, field_yaml)
1311 if (trim(cell_methods) .ne. "") &
1312 call register_variable_attribute(fms2io_fileobj, var_name, "cell_methods", &
1313 trim(adjustl(cell_methods)), str_len=len_trim(adjustl(cell_methods)))
1314
1315 !< Write out the cell_measures attribute (i.e Area, Volume)
1316 !! The diag field ids for the Area and Volume are sent in the register call
1317 !! This was defined in file object and passed in here
1318 if (trim(cell_measures) .ne. "") &
1319 call register_variable_attribute(fms2io_fileobj, var_name, "cell_measures", &
1320 trim(adjustl(cell_measures)), str_len=len_trim(adjustl(cell_measures)))
1321
1322 !< Write out the standard_name (this was defined in the register call)
1323 if (this%has_standname()) &
1324 call register_variable_attribute(fms2io_fileobj, var_name, "standard_name", &
1325 trim(this%get_standname()), str_len=len_trim(this%get_standname()))
1326
1327 call this%write_coordinate_attribute(fms2io_fileobj, var_name, diag_axis)
1328
1329 if (field_yaml%has_var_attributes()) then
1330 yaml_field_attributes = field_yaml%get_var_attributes()
1331 do i = 1, size(yaml_field_attributes,1)
1332 call register_variable_attribute(fms2io_fileobj, var_name, trim(yaml_field_attributes(i,1)), &
1333 trim(yaml_field_attributes(i,2)), str_len=len_trim(yaml_field_attributes(i,2)))
1334 enddo
1335 deallocate(yaml_field_attributes)
1336 endif
1337end subroutine write_field_metadata
1338
1339!> @brief Writes the coordinate attribute of a field if any of the field's axis has an
1340!! auxiliary axis
1341subroutine write_coordinate_attribute (this, fms2io_fileobj, var_name, diag_axis)
1342 CLASS(fmsdiagfield_type), intent(in) :: this !< The field object
1343 class(fmsnetcdffile_t), INTENT(INOUT) :: fms2io_fileobj!< Fms2_io fileobj to write to
1344 character(len=*), intent(in) :: var_name !< Variable name
1345 class(fmsdiagaxiscontainer_type), intent(in) :: diag_axis(:) !< Diag_axis object
1346
1347 integer :: i !< For do loops
1348 character(len = 252) :: aux_coord !< Auxuliary axis name
1349
1350 !> If the variable is a scalar, go away
1351 if (.not. allocated(this%axis_ids)) return
1352
1353 !> Determine if any of the field's axis has an auxiliary axis and the
1354 !! axis_names as a variable attribute
1355 aux_coord = ""
1356 do i = 1, size(this%axis_ids)
1357 select type (obj => diag_axis(this%axis_ids(i))%axis)
1358 type is (fmsdiagfullaxis_type)
1359 if (obj%has_aux()) then
1360 aux_coord = trim(aux_coord)//" "//obj%get_aux()
1361 endif
1362 end select
1363 enddo
1364
1365 if (trim(aux_coord) .eq. "") return
1366
1367 call register_variable_attribute(fms2io_fileobj, var_name, "coordinates", &
1368 trim(adjustl(aux_coord)), str_len=len_trim(adjustl(aux_coord)))
1369
1370end subroutine write_coordinate_attribute
1371
1372!> @brief Gets a fields data buffer
1373!! @return a pointer to the data buffer
1374function get_data_buffer (this) &
1375 result(rslt)
1376 class(fmsdiagfield_type), target, intent(in) :: this !< diag field
1377 class(*),dimension(:,:,:,:), pointer :: rslt !< The field's data buffer
1378
1379 if (.not. this%data_buffer_is_allocated) &
1380 call mpp_error(fatal, "The input data buffer for the field:"&
1381 //trim(this%varname)//" was never allocated.")
1382
1383 rslt => this%input_data_buffer%get_buffer()
1384end function get_data_buffer
1385
1386
1387!> @brief Gets a fields weight buffer
1388!! @return a pointer to the weight buffer
1389function get_weight (this) &
1390 result(rslt)
1391 class(fmsdiagfield_type), target, intent(in) :: this !< diag field
1392 type(real(kind=r8_kind)), pointer :: rslt
1393
1394 if (.not. this%data_buffer_is_allocated) &
1395 call mpp_error(fatal, "The input data buffer for the field:"&
1396 //trim(this%varname)//" was never allocated.")
1397
1398 rslt => this%input_data_buffer%get_weight()
1399end function get_weight
1400
1401!> Gets the flag telling if the math functions need to be done
1402!! \return Copy of math_needs_to_be_done flag
1403pure logical function get_math_needs_to_be_done(this)
1404 class(fmsdiagfield_type), intent(in) :: this !< diag object
1405 get_math_needs_to_be_done = .false.
1406 if (allocated(this%math_needs_to_be_done)) get_math_needs_to_be_done = this%math_needs_to_be_done
1407end function get_math_needs_to_be_done
1408!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1409!!!!! Allocation checks
1410
1411!!> @brief Checks if obj%diag_field is allocated
1412!!! @return true if obj%diag_field is allocated
1413!logical function has_diag_field (obj)
1414! class (fmsDiagField_type), intent(in) :: obj !< diag object
1415! has_diag_field = allocated(obj%diag_field)
1416!end function has_diag_field
1417!> @brief Checks if obj%diag_id is allocated
1418!! @return true if obj%diag_id is allocated
1419pure logical function has_diag_id (this)
1420 class(fmsdiagfield_type), intent(in) :: this !< diag object
1421 has_diag_id = allocated(this%diag_id)
1422end function has_diag_id
1423
1424!> @brief Checks if obj%metadata is allocated
1425!! @return true if obj%metadata is allocated
1426pure logical function has_attributes (this)
1427 class(fmsdiagfield_type), intent(in) :: this !< diag object
1428 has_attributes = this%num_attributes > 0
1429end function has_attributes
1430
1431!> @brief Checks if obj%static is allocated
1432!! @return true if obj%static is allocated
1433pure logical function has_static (this)
1434 class(fmsdiagfield_type), intent(in) :: this !< diag object
1435 has_static = allocated(this%static)
1436end function has_static
1437
1438!> @brief Checks if obj%registered is allocated
1439!! @return true if obj%registered is allocated
1440pure logical function has_registered (this)
1441 class(fmsdiagfield_type), intent(in) :: this !< diag object
1442 has_registered = allocated(this%registered)
1443end function has_registered
1444
1445!> @brief Checks if obj%mask_variant is allocated
1446!! @return true if obj%mask_variant is allocated
1447pure logical function has_mask_variant (this)
1448 class(fmsdiagfield_type), intent(in) :: this !< diag object
1449 has_mask_variant = allocated(this%mask_variant)
1450end function has_mask_variant
1451
1452!> @brief Checks if obj%local is allocated
1453!! @return true if obj%local is allocated
1454pure logical function has_local (this)
1455 class(fmsdiagfield_type), intent(in) :: this !< diag object
1456 has_local = allocated(this%local)
1457end function has_local
1458
1459!> @brief Checks if obj%vartype is allocated
1460!! @return true if obj%vartype is allocated
1461pure logical function has_vartype (this)
1462 class(fmsdiagfield_type), intent(in) :: this !< diag object
1463 has_vartype = allocated(this%vartype)
1464end function has_vartype
1465
1466!> @brief Checks if obj%varname is allocated
1467!! @return true if obj%varname is allocated
1468pure logical function has_varname (this)
1469 class(fmsdiagfield_type), intent(in) :: this !< diag object
1470 has_varname = allocated(this%varname)
1471end function has_varname
1472
1473!> @brief Checks if obj%longname is allocated
1474!! @return true if obj%longname is allocated
1475pure logical function has_longname (this)
1476 class(fmsdiagfield_type), intent(in) :: this !< diag object
1477 has_longname = allocated(this%longname)
1478end function has_longname
1479
1480!> @brief Checks if obj%standname is allocated
1481!! @return true if obj%standname is allocated
1482pure logical function has_standname (this)
1483 class(fmsdiagfield_type), intent(in) :: this !< diag object
1484 has_standname = allocated(this%standname)
1485end function has_standname
1486
1487!> @brief Checks if obj%units is allocated
1488!! @return true if obj%units is allocated
1489pure logical function has_units (this)
1490 class(fmsdiagfield_type), intent(in) :: this !< diag object
1491 has_units = allocated(this%units)
1492end function has_units
1493
1494!> @brief Checks if obj%modname is allocated
1495!! @return true if obj%modname is allocated
1496pure logical function has_modname (this)
1497 class(fmsdiagfield_type), intent(in) :: this !< diag object
1498 has_modname = allocated(this%modname)
1499end function has_modname
1500
1501!> @brief Checks if obj%realm is allocated
1502!! @return true if obj%realm is allocated
1503pure logical function has_realm (this)
1504 class(fmsdiagfield_type), intent(in) :: this !< diag object
1505 has_realm = allocated(this%realm)
1506end function has_realm
1507
1508!> @brief Checks if obj%interp_method is allocated
1509!! @return true if obj%interp_method is allocated
1510pure logical function has_interp_method (this)
1511 class(fmsdiagfield_type), intent(in) :: this !< diag object
1512 has_interp_method = allocated(this%interp_method)
1513end function has_interp_method
1514
1515!> @brief Checks if obj%frequency is allocated
1516!! @return true if obj%frequency is allocated
1517pure logical function has_frequency (this)
1518 class(fmsdiagfield_type), intent(in) :: this !< diag object
1519 has_frequency = allocated(this%frequency)
1520end function has_frequency
1521
1522!> @brief Checks if obj%tile_count is allocated
1523!! @return true if obj%tile_count is allocated
1524pure logical function has_tile_count (this)
1525 class(fmsdiagfield_type), intent(in) :: this !< diag object
1526 has_tile_count = allocated(this%tile_count)
1527end function has_tile_count
1528
1529!> @brief Checks if axis_ids of the object is allocated
1530!! @return true if it is allocated
1531pure logical function has_axis_ids (this)
1532 class(fmsdiagfield_type), intent(in) :: this !< diag field object
1533 has_axis_ids = allocated(this%axis_ids)
1534end function has_axis_ids
1535
1536!> @brief Checks if obj%area is allocated
1537!! @return true if obj%area is allocated
1538pure logical function has_area (this)
1539 class(fmsdiagfield_type), intent(in) :: this !< diag object
1540 has_area = allocated(this%area)
1541end function has_area
1542
1543!> @brief Checks if obj%volume is allocated
1544!! @return true if obj%volume is allocated
1545pure logical function has_volume (this)
1546 class(fmsdiagfield_type), intent(in) :: this !< diag object
1547 has_volume = allocated(this%volume)
1548end function has_volume
1549
1550!> @brief Checks if obj%missing_value is allocated
1551!! @return true if obj%missing_value is allocated
1552pure logical function has_missing_value (this)
1553 class(fmsdiagfield_type), intent(in) :: this !< diag object
1554 has_missing_value = allocated(this%missing_value)
1555end function has_missing_value
1556
1557!> @brief Checks if obj%data_RANGE is allocated
1558!! @return true if obj%data_RANGE is allocated
1559pure logical function has_data_range (this)
1560 class(fmsdiagfield_type), intent(in) :: this !< diag object
1561 has_data_range = allocated(this%data_RANGE)
1562end function has_data_range
1563
1564!> @brief Checks if obj%input_data_buffer is allocated
1565!! @return true if obj%input_data_buffer is allocated
1566pure logical function has_input_data_buffer (this)
1567 class(fmsdiagfield_type), intent(in) :: this !< diag object
1568 has_input_data_buffer = allocated(this%input_data_buffer)
1569end function has_input_data_buffer
1570
1571!> @brief Add a attribute to the diag_obj using the diag_field_id
1572subroutine diag_field_add_attribute(this, att_name, att_value)
1573 class(fmsdiagfield_type), intent (inout) :: this !< The field object
1574 character(len=*), intent(in) :: att_name !< Name of the attribute
1575 class(*), intent(in) :: att_value(:) !< The attribute value to add
1576
1577 this%num_attributes = this%num_attributes + 1
1578 if (this%num_attributes > max_field_attributes) &
1579 call mpp_error(fatal, "diag_field_add_attribute: Number of attributes exceeds max_field_attributes for field:"&
1580 //trim(this%varname)//". Increase diag_manager_nml:max_field_attributes.")
1581
1582 call this%attributes(this%num_attributes)%add(att_name, att_value)
1583end subroutine diag_field_add_attribute
1584
1585!> @brief Determine the default missing value to use based on the requested variable type
1586!! @return The missing value
1587function get_default_missing_value(var_type) &
1588 result(rslt)
1589
1590 integer, intent(in) :: var_type !< The type of the variable to return the missing value as
1591 class(*),allocatable :: rslt
1592
1593 select case(var_type)
1594 case (r4)
1595 allocate(real(kind=r4_kind) :: rslt)
1596 rslt = real(cmor_missing_value, kind=r4_kind)
1597 case (r8)
1598 allocate(real(kind=r8_kind) :: rslt)
1599 rslt = real(cmor_missing_value, kind=r8_kind)
1600 case default
1601 end select
1602end function
1603
1604!> @brief Determines the diag_obj id corresponding to a module name and field_name
1605!> @return diag_obj id
1606PURE FUNCTION diag_field_id_from_name(this, module_name, field_name) &
1607 result(diag_field_id)
1608 CLASS(fmsdiagfield_type), INTENT(in) :: this !< The field object
1609 CHARACTER(len=*), INTENT(in) :: module_name !< Module name that registered the variable
1610 CHARACTER(len=*), INTENT(in) :: field_name !< Variable name
1611
1612 integer :: diag_field_id
1613
1614 diag_field_id = diag_field_not_found
1615 if (this%get_varname() .eq. trim(field_name) .and. &
1616 this%get_modname() .eq. trim(module_name)) then
1617 diag_field_id = this%get_id()
1618 endif
1619end function diag_field_id_from_name
1620
1621!> @brief Adds the area and volume id to a field object
1622subroutine add_area_volume(this, area, volume)
1623 CLASS(fmsdiagfield_type), intent(inout) :: this !< The field object
1624 INTEGER, optional, INTENT(in) :: area !< diag ids of area
1625 INTEGER, optional, INTENT(in) :: volume !< diag ids of volume
1626
1627 if (present(area)) then
1628 if (area > 0) then
1629 this%area = area
1630 else
1631 call mpp_error(fatal, "diag_field_add_cell_measures: the area id is not valid. &
1632 &Verify that the area_id passed in to the field:"//this%varname// &
1633 " is valid and that the field is registered and in the diag_table.yaml")
1634 endif
1635 endif
1636
1637 if (present(volume)) then
1638 if (volume > 0) then
1639 this%volume = volume
1640 else
1641 call mpp_error(fatal, "diag_field_add_cell_measures: the volume id is not valid. &
1642 &Verify that the volume_id passed in to the field:"//this%varname// &
1643 " is valid and that the field is registered and in the diag_table.yaml")
1644 endif
1645 endif
1646
1647end subroutine add_area_volume
1648
1649!> @brief Append the time cell meathods based on the variable's reduction
1650subroutine append_time_cell_methods(this, cell_methods, field_yaml)
1651 class(fmsdiagfield_type), target, intent(inout) :: this !< diag field
1652 character(len=*), intent(inout) :: cell_methods !< The cell methods var to append to
1653 type(diagyamlfilesvar_type), intent(in) :: field_yaml !< The field's yaml
1654
1655 if (this%static) then
1656 cell_methods = trim(cell_methods)//" time: point "
1657 return
1658 endif
1659
1660 select case (field_yaml%get_var_reduction())
1661 case (time_none)
1662 cell_methods = trim(cell_methods)//" time: point "
1663 case (time_diurnal)
1664 cell_methods = trim(cell_methods)//" time: mean"
1665 case (time_power)
1666 cell_methods = trim(cell_methods)//" time: mean_pow"//int2str(field_yaml%get_pow_value())
1667 case (time_rms)
1668 cell_methods = trim(cell_methods)//" time: root_mean_square"
1669 case (time_max)
1670 cell_methods = trim(cell_methods)//" time: max"
1671 case (time_min)
1672 cell_methods = trim(cell_methods)//" time: min"
1673 case (time_average)
1674 cell_methods = trim(cell_methods)//" time: mean"
1675 case (time_sum)
1676 cell_methods = trim(cell_methods)//" time: sum"
1677 end select
1678end subroutine append_time_cell_methods
1679
1680!> Dumps any data from a given fmsDiagField_type object
1681subroutine dump_field_obj (this, unit_num)
1682 class(fmsdiagfield_type), intent(in) :: this
1683 integer, intent(in) :: unit_num !< passed in from dump_diag_obj if log file is being written to
1684 integer :: i
1685
1686 if( mpp_pe() .eq. mpp_root_pe()) then
1687 if( allocated(this%file_ids)) write(unit_num, *) 'file_ids:' ,this%file_ids
1688 if( allocated(this%diag_id)) write(unit_num, *) 'diag_id:' ,this%diag_id
1689 if( allocated(this%static)) write(unit_num, *) 'static:' ,this%static
1690 if( allocated(this%registered)) write(unit_num, *) 'registered:' ,this%registered
1691 if( allocated(this%mask_variant)) write(unit_num, *) 'mask_variant:' ,this%mask_variant
1692 if( allocated(this%do_not_log)) write(unit_num, *) 'do_not_log:' ,this%do_not_log
1693 if( allocated(this%local)) write(unit_num, *) 'local:' ,this%local
1694 if( allocated(this%vartype)) write(unit_num, *) 'vartype:' ,this%vartype
1695 if( allocated(this%varname)) write(unit_num, *) 'varname:' ,this%varname
1696 if( allocated(this%longname)) write(unit_num, *) 'longname:' ,this%longname
1697 if( allocated(this%standname)) write(unit_num, *) 'standname:' ,this%standname
1698 if( allocated(this%units)) write(unit_num, *) 'units:' ,this%units
1699 if( allocated(this%modname)) write(unit_num, *) 'modname:' ,this%modname
1700 if( allocated(this%realm)) write(unit_num, *) 'realm:' ,this%realm
1701 if( allocated(this%interp_method)) write(unit_num, *) 'interp_method:' ,this%interp_method
1702 if( allocated(this%tile_count)) write(unit_num, *) 'tile_count:' ,this%tile_count
1703 if( allocated(this%axis_ids)) write(unit_num, *) 'axis_ids:' ,this%axis_ids
1704 write(unit_num, *) 'type_of_domain:' ,this%type_of_domain
1705 if( allocated(this%area)) write(unit_num, *) 'area:' ,this%area
1706 if( allocated(this%missing_value)) then
1707 select type(missing_val => this%missing_value)
1708 type is (real(r4_kind))
1709 write(unit_num, *) 'missing_value:', missing_val
1710 type is (real(r8_kind))
1711 write(unit_num, *) 'missing_value:' ,missing_val
1712 type is(integer(i4_kind))
1713 write(unit_num, *) 'missing_value:' ,missing_val
1714 type is(integer(i8_kind))
1715 write(unit_num, *) 'missing_value:' ,missing_val
1716 end select
1717 endif
1718 if( allocated( this%data_RANGE)) then
1719 select type(drange => this%data_RANGE)
1720 type is (real(r4_kind))
1721 write(unit_num, *) 'data_RANGE:' ,drange
1722 type is (real(r8_kind))
1723 write(unit_num, *) 'data_RANGE:' ,drange
1724 type is(integer(i4_kind))
1725 write(unit_num, *) 'data_RANGE:' ,drange
1726 type is(integer(i8_kind))
1727 write(unit_num, *) 'data_RANGE:' ,drange
1728 end select
1729 endif
1730 write(unit_num, *) 'num_attributes:' ,this%num_attributes
1731 if( allocated(this%attributes)) then
1732 do i=1, this%num_attributes
1733 if( allocated(this%attributes(i)%att_value)) then
1734 select type( val => this%attributes(i)%att_value)
1735 type is (real(r8_kind))
1736 write(unit_num, *) 'attribute name', this%attributes(i)%att_name, 'val:', val
1737 type is (real(r4_kind))
1738 write(unit_num, *) 'attribute name', this%attributes(i)%att_name, 'val:', val
1739 type is (integer(i4_kind))
1740 write(unit_num, *) 'attribute name', this%attributes(i)%att_name, 'val:', val
1741 type is (integer(i8_kind))
1742 write(unit_num, *) 'attribute name', this%attributes(i)%att_name, 'val:', val
1743 end select
1744 endif
1745 enddo
1746 endif
1747
1748 endif
1749
1750end subroutine
1751
1752!< @brief Get the starting compute domain indices for a set of axis
1753!! @return compute domain starting indices
1754function get_starting_compute_domain(axis_ids, diag_axis) &
1755result(compute_domain)
1756 integer, intent(in) :: axis_ids(:) !< Array of axis ids
1757 class(fmsdiagaxiscontainer_type),intent(in) :: diag_axis(:) !< Array of axis object
1758
1759 integer :: compute_domain(4)
1760 integer :: a !< For looping through axes
1761 integer :: compute_idx(2) !< Compute domain indices (starting, ending)
1762 logical :: dummy !< Dummy variable for the `get_compute_domain` subroutine
1763
1764 compute_domain = 1
1765 axis_loop: do a = 1,size(axis_ids)
1766 select type (axis => diag_axis(axis_ids(a))%axis)
1767 type is (fmsdiagfullaxis_type)
1768 call axis%get_compute_domain(compute_idx, dummy)
1769 if ( compute_idx(1) .ne. diag_null) compute_domain(a) = compute_idx(1)
1770 end select
1771 enddo axis_loop
1772end function get_starting_compute_domain
1773
1774!> Get list of field ids
1775pure function get_file_ids(this)
1776 class(fmsdiagfield_type), intent(in) :: this
1777 integer, allocatable :: get_file_ids(:) !< Ids of the FMS_diag_files the variable
1778 get_file_ids = this%file_ids
1779end function
1780
1781!> @brief Get the mask from the input buffer object
1782!! @return a pointer to the mask
1783function get_mask(this)
1784 class(fmsdiagfield_type), target, intent(in) :: this !< input buffer object
1785 logical, pointer :: get_mask(:,:,:,:)
1786 get_mask => this%mask
1787end function get_mask
1788
1789!> @brief If in openmp region, omp_axis should be provided in order to allocate to the given axis lengths.
1790!! Otherwise mask will be allocated to the size of mask_in
1791subroutine allocate_mask(this, mask_in, omp_axis)
1792 class(fmsdiagfield_type), target, intent(inout) :: this !< input buffer object
1793 logical, intent(in) :: mask_in(:,:,:,:)
1794 class(fmsdiagaxiscontainer_type), intent(in), optional :: omp_axis(:) !< true if calling from omp region
1795 integer :: axis_num, length(4)
1796 integer, pointer :: id_num
1797 ! if not omp just allocate to whatever is given
1798 if(.not. present(omp_axis)) then
1799 allocate(this%mask(size(mask_in,1), size(mask_in,2), size(mask_in,3), &
1800 size(mask_in,4)))
1801 ! otherwise loop through axis and get sizes
1802 else
1803 length = 1
1804 do axis_num=1, size(this%axis_ids)
1805 id_num => this%axis_ids(axis_num)
1806 select type(axis => omp_axis(id_num)%axis)
1807 type is (fmsdiagfullaxis_type)
1808 length(axis_num) = axis%axis_length()
1809 end select
1810 enddo
1811 allocate(this%mask(length(1), length(2), length(3), length(4)))
1812 endif
1813end subroutine allocate_mask
1814
1815!> Sets previously allocated mask to mask_in at given index ranges
1816subroutine set_mask(this, mask_in, field_info, is, js, ks, ie, je, ke)
1817 class(fmsdiagfield_type), intent(inout) :: this
1818 logical, intent(in) :: mask_in(:,:,:,:)
1819 character(len=*), intent(in) :: field_info !< Field info to add to error message
1820 integer, optional, intent(in) :: is, js, ks, ie, je, ke
1821 if(present(is)) then
1822 if(is .lt. lbound(this%mask,1) .or. ie .gt. ubound(this%mask,1) .or. &
1823 js .lt. lbound(this%mask,2) .or. je .gt. ubound(this%mask,2) .or. &
1824 ks .lt. lbound(this%mask,3) .or. ke .gt. ubound(this%mask,3)) then
1825 print *, "PE:", int2str(mpp_pe()), "The size of the mask is", &
1826 shape(this%mask), &
1827 "But the indices passed in are is=", int2str(is), " ie=", int2str(ie),&
1828 " js=", int2str(js), " je=", int2str(je), &
1829 " ks=", int2str(ks), " ke=", int2str(ke), &
1830 " ", trim(field_info)
1831 call mpp_error(fatal,"set_mask:: given indices out of bounds for allocated mask")
1832 endif
1833 this%mask(is:ie, js:je, ks:ke, :) = mask_in
1834 else
1835 this%mask = mask_in
1836 endif
1837end subroutine set_mask
1838
1839!> sets halo_present to true
1840subroutine set_halo_present(this)
1841 class(fmsdiagfield_type), intent(inout) :: this !< field object to modify
1842 this%halo_present = .true.
1843end subroutine set_halo_present
1844
1845!> Getter for halo_present
1846pure function is_halo_present(this)
1847 class(fmsdiagfield_type), intent(in) :: this !< field object to get from
1848 logical :: is_halo_present
1849 is_halo_present = this%halo_present
1850end function is_halo_present
1851
1852!> Helper routine to find and set the netcdf missing value for a field
1853!! Always returns r8 due to reduction routine args
1854!! casts up to r8 from given missing val or default if needed
1855function find_missing_value(this, missing_val) &
1856 result(res)
1857 class(fmsdiagfield_type), intent(in) :: this !< field object to get missing value for
1858 class(*), allocatable, intent(out) :: missing_val !< outputted netcdf missing value (oriignal type)
1859 real(r8_kind), allocatable :: res !< returned r8 copy of missing_val
1860 integer :: vtype !< temp to hold enumerated variable type
1861
1862 if(this%has_missing_value()) then
1863 missing_val = this%get_missing_value(this%get_vartype())
1864 else
1865 vtype = this%get_vartype()
1866 if(vtype .eq. r8) then
1867 missing_val = cmor_missing_value
1868 else
1869 missing_val = real(cmor_missing_value, r4_kind)
1870 endif
1871 endif
1872
1873 select type(missing_val)
1874 type is (real(r8_kind))
1875 res = missing_val
1876 type is (real(r4_kind))
1877 res = real(missing_val, r8_kind)
1878 end select
1879end function find_missing_value
1880
1881!> @returns allocation status of logical mask array
1882!! this just indicates whether the mask array itself has been alloc'd
1883!! this is different from @ref has_mask_variant, which is set earlier for whether a mask is being used at all
1884pure logical function has_mask_allocated(this)
1885 class(fmsdiagfield_type),intent(in) :: this !< field object to check mask allocation for
1886 has_mask_allocated = allocated(this%mask)
1887end function has_mask_allocated
1888
1889!> @brief Determine if the variable is in the file
1890!! @return .True. if the varibale is in the file
1891pure function is_variable_in_file(this, file_id) &
1892result(res)
1893 class(fmsdiagfield_type), intent(in) :: this !< field object to check
1894 integer, intent(in) :: file_id !< File id to check
1895 logical :: res
1896
1897 integer :: i
1898
1899 res = .false.
1900 if (any(this%file_ids .eq. file_id)) res = .true.
1901end function is_variable_in_file
1902
1903!> @brief Determine the name of the first file the variable is in
1904!! @return filename
1905function get_field_file_name(this) &
1906 result(res)
1907 class(fmsdiagfield_type), intent(in) :: this !< Field object to query
1908 character(len=:), allocatable :: res
1909
1910 res = this%diag_field(1)%get_var_fname()
1911end function get_field_file_name
1912
1913!> @brief Generate the associated files attribute
1914subroutine generate_associated_files_att(this, att, start_time)
1915 class(fmsdiagfield_type) , intent(in) :: this !< diag_field_object for the area/volume field
1916 character(len=*), intent(inout) :: att !< associated_files_att
1917 type(time_type), intent(in) :: start_time !< The start_time for the field's file
1918
1919 character(len=:), allocatable :: field_name !< Name of the area/volume field
1920 character(len=FMS_FILE_LEN) :: file_name !< Name of the file the area/volume field is in!
1921 character(len=128) :: start_date !< Start date to append to the begining of the filename
1922
1923 integer :: year, month, day, hour, minute, second
1924 field_name = this%get_varname(to_write = .true.)
1925
1926 ! Check if the field is already in the associated files attribute (i.e the area can be associated with multiple
1927 ! fields in the file, but it only needs to be added once)
1928 if (index(att, field_name) .ne. 0) return
1929
1930 file_name = this%get_field_file_name()
1931
1932 if (prepend_date) then
1933 call get_date(start_time, year, month, day, hour, minute, second)
1934 write (start_date, '(1I20.4, 2I2.2)') year, month, day
1935 file_name = trim(adjustl(start_date))//'.'//trim(file_name)
1936 endif
1937
1938 att = trim(att)//" "//trim(field_name)//": "//trim(file_name)//".nc"
1939end subroutine generate_associated_files_att
1940
1941!> @brief Determines if the compute domain has been divide further into slices (i.e openmp blocks)
1942!! @return .True. if the compute domain has been divided furter into slices
1943function check_for_slices(field, diag_axis, var_size) &
1944 result(rslt)
1945 type(fmsdiagfield_type), intent(in) :: field !< Field object
1946 type(fmsdiagaxiscontainer_type), target, intent(in) :: diag_axis(:) !< Array of diag axis
1947 integer, intent(in) :: var_size(:) !< The size of the buffer pass into send_data
1948
1949 logical :: rslt
1950 integer :: i !< For do loops
1951
1952 rslt = .false.
1953
1954 if (.not. field%has_axis_ids()) then
1955 rslt = .false.
1956 return
1957 endif
1958 do i = 1, size(field%axis_ids)
1959 select type (axis_obj => diag_axis(field%axis_ids(i))%axis)
1960 type is (fmsdiagfullaxis_type)
1961 if (axis_obj%axis_length() .ne. var_size(i)) then
1962 rslt = .true.
1963 return
1964 endif
1965 end select
1966 enddo
1967end function
1968#endif
1969end module fms_diag_field_object_mod
integer, parameter max_str_len
Max length for a string.
integer, parameter no_domain
Use the FmsNetcdfFile_t fileobj.
character(len=7) avg_name
Name of the average fields.
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, 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.
real(r8_kind), parameter cmor_missing_value
CMOR standard missing value.
logical prepend_date
Should the history file have the start date prepended to the file name. .TRUE. is only supported if t...
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.
logical module_is_initialized
Indicate if diag_manager has been initialized.
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 max_field_attributes
Maximum number of user definable attributes per field. Liptak: Changed from 2 to 4 20170718.
Type to hold the attributes of the field/axis/file.
Defines a new field within the given file Example usage:
Definition fms2_io.F90:212
subroutine set_send_data_time(this, time)
Sets the time send data was called last.
real(kind=r8_kind) function, pointer get_weight(this)
Get the weight from the input buffer object.
type(time_type) function get_send_data_time(this)
Get the time send data was called last.
Type to hold the information needed for the input buffer This is used when set_math_needs_to_be_done ...
type(diagyamlfilesvar_type) function, dimension(:), allocatable, public get_diag_fields_entries(indices)
Gets the diag_field entries corresponding to the indices of the sorted variable_list.
integer function, public get_num_unique_fields()
Determine the number of unique diag_fields in the diag_yaml_object.
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 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)
Error handler.
Definition mpp.F90:382
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...
Type to represent amounts of time. Implemented as seconds and days to allow for larger intervals.
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.
Object that holds all variable information.
type to hold the info a diag_field