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