FMS  2026.01
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=i4_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 function get_varname (this, filename, to_write) &
810 result(rslt)
811  class(fmsdiagfield_type), intent(in) :: this !< diag object
812  character(len=*), optional, intent(in) :: filename !< Name of the diag_file we are writing to
813  logical, optional, intent(in) :: to_write !< .true. if getting the varname that will be writen to the file
814  character(len=:), allocatable :: rslt
815 
816  integer :: i
817 
818  rslt = this%varname
819 
820  !< If writing the varname can be the outname which is defined in the yaml
821  if (present(to_write)) then
822  if (to_write) then
823  if (.not. present(filename)) then
824  call mpp_error(fatal, "get_varname was called using the to_write optional argument, "//&
825  "but a filename was not provided!")
826  endif
827  ! Loop through all of the file the variable is in and find the outputname of the variable
828  do i = 1, size(this%diag_field)
829  if (trim(filename) .eq. trim(this%diag_field(i)%get_var_fname())) then
830  rslt = this%diag_field(i)%get_var_outname()
831  return
832  endif
833  enddo
834  endif
835  endif
836 
837 end function get_varname
838 
839 !> @brief Gets longname
840 !! @return copy of the variable long name or a single string if there is no long name
841 pure function get_longname (this) &
842 result(rslt)
843  class(fmsdiagfield_type), intent(in) :: this !< diag object
844  character(len=:), allocatable :: rslt
845  if (allocated(this%longname)) then
846  rslt = this%longname
847  else
848  rslt = diag_null_string
849  endif
850 end function get_longname
851 
852 !> @brief Gets standname
853 !! @return copy of the standard name or an empty string if standname is not allocated
854 pure function get_standname (this) &
855 result(rslt)
856  class(fmsdiagfield_type), intent(in) :: this !< diag object
857  character(len=:), allocatable :: rslt
858  if (allocated(this%standname)) then
859  rslt = this%standname
860  else
861  rslt = diag_null_string
862  endif
863 end function get_standname
864 
865 !> @brief Gets units
866 !! @return copy of the units or an empty string if not allocated
867 pure function get_units (this) &
868 result(rslt)
869  class(fmsdiagfield_type), intent(in) :: this !< diag object
870  character(len=:), allocatable :: rslt
871  if (allocated(this%units)) then
872  rslt = this%units
873  else
874  rslt = diag_null_string
875  endif
876 end function get_units
877 
878 !> @brief Gets modname
879 !! @return copy of the module name that the variable is in or an empty string if not allocated
880 pure function get_modname (this) &
881 result(rslt)
882  class(fmsdiagfield_type), intent(in) :: this !< diag object
883  character(len=:), allocatable :: rslt
884  if (allocated(this%modname)) then
885  rslt = this%modname
886  else
887  rslt = diag_null_string
888  endif
889 end function get_modname
890 
891 !> @brief Gets realm
892 !! @return copy of the variables modeling realm or an empty string if not allocated
893 pure function get_realm (this) &
894 result(rslt)
895  class(fmsdiagfield_type), intent(in) :: this !< diag object
896  character(len=:), allocatable :: rslt
897  if (allocated(this%realm)) then
898  rslt = this%realm
899  else
900  rslt = diag_null_string
901  endif
902 end function get_realm
903 
904 !> @brief Gets interp_method
905 !! @return copy of The interpolation method or an empty string if not allocated
906 pure function get_interp_method (this) &
907 result(rslt)
908  class(fmsdiagfield_type), intent(in) :: this !< diag object
909  character(len=:), allocatable :: rslt
910  if (allocated(this%interp_method)) then
911  rslt = this%interp_method
912  else
913  rslt = diag_null_string
914  endif
915 end function get_interp_method
916 
917 !> @brief Gets frequency
918 !! @return copy of the frequency or DIAG_NULL if obj%frequency is not allocated
919 pure function get_frequency (this) &
920 result(rslt)
921  class(fmsdiagfield_type), intent(in) :: this !< diag object
922  integer, allocatable, dimension (:) :: rslt
923  if (allocated(this%frequency)) then
924  allocate (rslt(size(this%frequency)))
925  rslt = this%frequency
926  else
927  allocate (rslt(1))
928  rslt = diag_null
929  endif
930 end function get_frequency
931 
932 !> @brief Gets tile_count
933 !! @return copy of the number of tiles or diag_null if tile_count is not allocated
934 pure function get_tile_count (this) &
935 result(rslt)
936  class(fmsdiagfield_type), intent(in) :: this !< diag object
937  integer :: rslt
938  if (allocated(this%tile_count)) then
939  rslt = this%tile_count
940  else
941  rslt = diag_null
942  endif
943 end function get_tile_count
944 
945 !> @brief Gets area
946 !! @return copy of the area or diag_null if not allocated
947 pure function get_area (this) &
948 result(rslt)
949  class(fmsdiagfield_type), intent(in) :: this !< diag object
950  integer :: rslt
951  if (allocated(this%area)) then
952  rslt = this%area
953  else
954  rslt = diag_null
955  endif
956 end function get_area
957 
958 !> @brief Gets volume
959 !! @return copy of the volume or diag_null if volume is not allocated
960 pure function get_volume (this) &
961 result(rslt)
962  class(fmsdiagfield_type), intent(in) :: this !< diag object
963  integer :: rslt
964  if (allocated(this%volume)) then
965  rslt = this%volume
966  else
967  rslt = diag_null
968  endif
969 end function get_volume
970 
971 !> @brief Gets missing_value
972 !! @return copy of The missing value
973 !! @note Netcdf requires the type of the variable and the type of the missing_value and _Fillvalue to be the same
974 !! var_type is the type of the variable which may not be in the same type as the missing_value in the register call
975 !! For example, if compiling with r8 but the in diag_table.yaml the kind is r4
976 function get_missing_value (this, var_type) &
977 result(rslt)
978  class(fmsdiagfield_type), intent(in) :: this !< diag object
979  integer, intent(in) :: var_type !< The type of the variable as it will writen to the netcdf file
980  !! and the missing value is return as
981 
982  class(*),allocatable :: rslt
983 
984  if (.not. allocated(this%missing_value)) then
985  call mpp_error ("get_missing_value", &
986  "The missing value is not allocated", fatal)
987  endif
988 
989  !< The select types are needed so that the missing_value can be correctly converted and copied as the needed variable
990  !! type
991  select case (var_type)
992  case (r4)
993  allocate (real(kind=r4_kind) :: rslt)
994  select type (miss => this%missing_value)
995  type is (real(kind=r4_kind))
996  select type (rslt)
997  type is (real(kind=r4_kind))
998  rslt = real(miss, kind=r4_kind)
999  end select
1000  type is (real(kind=r8_kind))
1001  select type (rslt)
1002  type is (real(kind=r4_kind))
1003  rslt = real(miss, kind=r4_kind)
1004  end select
1005  end select
1006  case (r8)
1007  allocate (real(kind=r8_kind) :: rslt)
1008  select type (miss => this%missing_value)
1009  type is (real(kind=r4_kind))
1010  select type (rslt)
1011  type is (real(kind=r8_kind))
1012  rslt = real(miss, kind=r8_kind)
1013  end select
1014  type is (real(kind=r8_kind))
1015  select type (rslt)
1016  type is (real(kind=r8_kind))
1017  rslt = real(miss, kind=r8_kind)
1018  end select
1019  end select
1020  end select
1021 
1022 end function get_missing_value
1023 
1024 !> @brief Gets data_range
1025 !! @return copy of the data range
1026 !! @note Netcdf requires the type of the variable and the type of the range to be the same
1027 !! var_type is the type of the variable which may not be in the same type as the range in the register call
1028 !! For example, if compiling with r8 but the in diag_table.yaml the kind is r4
1029 function get_data_range (this, var_type) &
1030 result(rslt)
1031  class(fmsdiagfield_type), intent(in) :: this !< diag object
1032  integer, intent(in) :: var_type !< The type of the variable as it will writen to the netcdf file
1033  !! and the data_range is returned as
1034  class(*),allocatable :: rslt(:)
1035 
1036  if ( .not. allocated(this%data_RANGE)) call mpp_error ("get_data_RANGE", &
1037  "The data_RANGE value is not allocated", fatal)
1038 
1039  !< The select types are needed so that the range can be correctly converted and copied as the needed variable
1040  !! type
1041  select case (var_type)
1042  case (r4)
1043  allocate (real(kind=r4_kind) :: rslt(2))
1044  select type (r => this%data_RANGE)
1045  type is (real(kind=r4_kind))
1046  select type (rslt)
1047  type is (real(kind=r4_kind))
1048  rslt = real(r, kind=r4_kind)
1049  end select
1050  type is (real(kind=r8_kind))
1051  select type (rslt)
1052  type is (real(kind=r4_kind))
1053  rslt = real(r, kind=r4_kind)
1054  end select
1055  end select
1056  case (r8)
1057  allocate (real(kind=r8_kind) :: rslt(2))
1058  select type (r => this%data_RANGE)
1059  type is (real(kind=r4_kind))
1060  select type (rslt)
1061  type is (real(kind=r8_kind))
1062  rslt = real(r, kind=r8_kind)
1063  end select
1064  type is (real(kind=r8_kind))
1065  select type (rslt)
1066  type is (real(kind=r8_kind))
1067  rslt = real(r, kind=r8_kind)
1068  end select
1069  end select
1070  end select
1071 end function get_data_range
1072 
1073 !> @brief Gets axis_ids
1074 !! @return pointer to the axis ids
1075 function get_axis_id (this) &
1076 result(rslt)
1077  class(fmsdiagfield_type), target, intent(in) :: this !< diag object
1078  integer, pointer, dimension(:) :: rslt !< field's axis_ids
1079 
1080  if(allocated(this%axis_ids)) then
1081  rslt => this%axis_ids
1082  else
1083  rslt => null()
1084  endif
1085 end function get_axis_id
1086 
1087 !> @brief Gets field's domain
1088 !! @return pointer to the domain
1089 function get_domain (this) &
1090 result(rslt)
1091  class(fmsdiagfield_type), target, intent(in) :: this !< diag field
1092  class(diagdomain_t), pointer :: rslt !< field's domain
1093 
1094  if (associated(this%domain)) then
1095  rslt => this%domain
1096  else
1097  rslt => null()
1098  endif
1099 
1100 end function get_domain
1101 
1102 !> @brief Gets field's type of domain
1103 !! @return integer defining the type of domain (NO_DOMAIN, TWO_D_DOMAIN, UG_DOMAIN)
1104 pure function get_type_of_domain (this) &
1105 result(rslt)
1106  class(fmsdiagfield_type), target, intent(in) :: this !< diag field
1107  integer :: rslt !< field's domain
1108 
1109  rslt = this%type_of_domain
1110 end function get_type_of_domain
1111 
1112 !> @brief Set the file ids of the files that the field belongs to
1113 subroutine set_file_ids(this, file_ids)
1114  class(fmsdiagfield_type), intent(inout) :: this !< diag field
1115  integer, intent(in) :: file_ids(:) !< File_ids to add
1116 
1117  allocate(this%file_ids(size(file_ids)))
1118  this%file_ids = file_ids
1119 end subroutine set_file_ids
1120 
1121 !> @brief Get the kind of the variable based on the yaml
1122 !! @return A string indicating the kind of the variable (as it is used in fms2_io)
1123 pure function get_var_skind(this, field_yaml) &
1124 result(rslt)
1125  class(fmsdiagfield_type), intent(in) :: this !< diag field
1126  type(diagyamlfilesvar_type), intent(in) :: field_yaml !< The corresponding yaml of the field
1127 
1128  character(len=:), allocatable :: rslt
1129 
1130  integer :: var_kind !< The integer corresponding to the kind of the variable (i4, i8, r4, r8)
1131 
1132  var_kind = field_yaml%get_var_kind()
1133  select case (var_kind)
1134  case (r4)
1135  rslt = "float"
1136  case (r8)
1137  rslt = "double"
1138  case (i4)
1139  rslt = "int"
1140  case (i8)
1141  rslt = "int64"
1142  end select
1143 
1144 end function get_var_skind
1145 
1146 !> @brief Get the multiple_send_data member of the field object
1147 !! @return multiple_send_data of the field
1148 pure function get_multiple_send_data(this) &
1149 result(rslt)
1150  class(fmsdiagfield_type), intent(in) :: this !< diag field
1151  logical :: rslt
1152  rslt = this%multiple_send_data
1153 end function get_multiple_send_data
1154 
1155 !> @brief Determine the long name to write for the field
1156 !! @return Long name to write
1157 pure function get_longname_to_write(this, field_yaml) &
1158 result(rslt)
1159  class(fmsdiagfield_type), intent(in) :: this !< diag field
1160  type(diagyamlfilesvar_type), intent(in) :: field_yaml !< The corresponding yaml of the field
1161 
1162  character(len=:), allocatable :: rslt
1163 
1164  rslt = field_yaml%get_var_longname() !! This is the long name defined in the yaml
1165  if (rslt .eq. "") then !! If the long name is not defined in the yaml, use the long name in the
1166  !! register_diag_field
1167  rslt = this%get_longname()
1168  else
1169  return
1170  endif
1171  if (rslt .eq. "") then !! If the long name is not defined in the yaml and in the register_diag_field
1172  !! use the variable name
1173  rslt = field_yaml%get_var_varname()
1174  endif
1175 end function get_longname_to_write
1176 
1177 !> @brief Determine the dimension names to use when registering the field to fms2_io
1178 subroutine get_dimnames(this, diag_axis, field_yaml, unlim_dimname, dimnames, is_regional, &
1179  file_axis_ids)
1180  class(fmsdiagfield_type), target, intent(inout) :: this !< diag field
1181  class(fmsdiagaxiscontainer_type), target, intent(in) :: diag_axis(:) !< Diag_axis object
1182  type(diagyamlfilesvar_type), intent(in) :: field_yaml !< Field info from diag_table yaml
1183  character(len=*), intent(in) :: unlim_dimname !< The name of unlimited dimension
1184  character(len=120), allocatable, intent(out) :: dimnames(:) !< Array of the dimension names
1185  !! for the field
1186  logical, intent(in) :: is_regional !< Flag indicating if the field is regional
1187  integer, intent(in) :: file_axis_ids(:) !< Ids of the file axis
1188 
1189  integer :: i !< For do loops
1190  integer :: naxis !< Number of axis for the field
1191  class(fmsdiagaxiscontainer_type), pointer :: axis_ptr !diag_axis(this%axis_ids(i), for convenience
1192  character(len=23) :: diurnal_axis_name !< name of the diurnal axis
1193 
1194  if (this%is_static()) then
1195  naxis = size(this%axis_ids)
1196  else
1197  naxis = size(this%axis_ids) + 1 !< Adding 1 more dimension for the unlimited dimension
1198  endif
1199 
1200  if (field_yaml%has_n_diurnal()) then
1201  naxis = naxis + 1 !< Adding 1 more dimension for the diurnal axis
1202  endif
1203 
1204  allocate(dimnames(naxis))
1205 
1206  !< Duplicated do loops for performance
1207  if (field_yaml%has_var_zbounds()) then
1208  do i = 1, size(this%axis_ids)
1209  axis_ptr => diag_axis(this%axis_ids(i))
1210  if (axis_ptr%axis%is_z_axis()) then
1211  call find_z_sub_axis_name(dimnames(i), this%axis_ids(i), file_axis_ids, field_yaml, diag_axis)
1212  else
1213  dimnames(i) = axis_ptr%axis%get_axis_name(is_regional)
1214  endif
1215  enddo
1216  else
1217  do i = 1, size(this%axis_ids)
1218  axis_ptr => diag_axis(this%axis_ids(i))
1219  dimnames(i) = axis_ptr%axis%get_axis_name(is_regional)
1220  enddo
1221  endif
1222 
1223  !< The second to last dimension is always the diurnal axis
1224  if (field_yaml%has_n_diurnal()) then
1225  WRITE (diurnal_axis_name,'(a,i2.2)') 'time_of_day_', field_yaml%get_n_diurnal()
1226  dimnames(naxis - 1) = trim(diurnal_axis_name)
1227  endif
1228 
1229  !< The last dimension is always the unlimited dimensions
1230  if (.not. this%is_static()) dimnames(naxis) = unlim_dimname
1231 
1232 end subroutine get_dimnames
1233 
1234 !> @brief Wrapper for the register_field call. The select types are needed so that the code can go
1235 !! in the correct interface
1236 subroutine register_field_wrap(fms2io_fileobj, varname, vartype, dimensions, chunksizes)
1237  class(fmsnetcdffile_t), INTENT(INOUT) :: fms2io_fileobj!< Fms2_io fileobj to write to
1238  character(len=*), INTENT(IN) :: varname !< Name of the variable
1239  character(len=*), INTENT(IN) :: vartype !< The type of the variable
1240  character(len=*), optional, INTENT(IN) :: dimensions(:) !< The dimension names of the field
1241  integer, optional, INTENT(IN) :: chunksizes(:) !< Chunksize to use, only relevant when using
1242  !! NETCDF-4 and the variable is domain decomposed
1243 
1244  select type(fms2io_fileobj)
1245  type is (fmsnetcdffile_t)
1246  call register_field(fms2io_fileobj, varname, vartype, dimensions)
1247  type is (fmsnetcdfdomainfile_t)
1248  call register_field(fms2io_fileobj, varname, vartype, dimensions, chunksizes=chunksizes)
1249  type is (fmsnetcdfunstructureddomainfile_t)
1250  call register_field(fms2io_fileobj, varname, vartype, dimensions)
1251  end select
1252 end subroutine register_field_wrap
1253 
1254 !> @brief Write the field's metadata to the file
1255 subroutine write_field_metadata(this, fms2io_fileobj, file_id, yaml_id, diag_axis, unlim_dimname, is_regional, &
1256  cell_measures, use_collective_writes, file_axis_ids)
1257  class(fmsdiagfield_type), target, intent(inout) :: this !< diag field
1258  class(fmsnetcdffile_t), INTENT(INOUT) :: fms2io_fileobj!< Fms2_io fileobj to write to
1259  integer, intent(in) :: file_id !< File id of the file to write to
1260  integer, intent(in) :: yaml_id !< Yaml id of the yaml entry of this field
1261  class(fmsdiagaxiscontainer_type), intent(in) :: diag_axis(:) !< Diag_axis object
1262  character(len=*), intent(in) :: unlim_dimname !< The name of the unlimited dimension
1263  logical, intent(in) :: is_regional !< Flag indicating if the field is regional
1264  character(len=*), intent(in) :: cell_measures !< The cell measures attribute to write
1265  logical, intent(in) :: use_collective_writes !< True if using collective writes
1266  !! for this variable
1267  integer, intent(in) :: file_axis_ids(:) !< Ids of all of the axis in thje file
1268 
1269  type(diagyamlfilesvar_type), pointer :: field_yaml !< pointer to the yaml entry
1270  character(len=:), allocatable :: var_name !< Variable name
1271  character(len=:), allocatable :: long_name !< Longname to write
1272  character(len=:), allocatable :: units !< Units of the field to write
1273  character(len=120), allocatable :: dimnames(:) !< Dimension names of the field
1274  character(len=120) :: cell_methods!< Cell methods attribute to write
1275  integer :: i !< For do loops
1276  character (len=MAX_STR_LEN), allocatable :: yaml_field_attributes(:,:) !< Variable attributes defined in the yaml
1277  character(len=:), allocatable :: interp_method_tmp !< temp to hold the name of the interpolation method
1278  integer :: interp_method_len !< length of the above string
1279 
1280  integer, allocatable :: chunksizes(:) !< Chunksizes to use for the variable
1281 
1282  field_yaml => diag_yaml%get_diag_field_from_id(yaml_id)
1283  var_name = field_yaml%get_var_outname()
1284 
1285  if (allocated(this%axis_ids)) then
1286  call this%get_dimnames(diag_axis, field_yaml, unlim_dimname, dimnames, is_regional, file_axis_ids)
1287 
1288  !! Collective writes are only used for 2D+ variables
1289  if ((use_collective_writes .and. size(this%axis_ids) >= 2) .or. field_yaml%has_chunksizes()) then
1290  chunksizes = this%get_chunksizes(diag_axis, field_yaml)
1291  call register_field_wrap(fms2io_fileobj, var_name, this%get_var_skind(field_yaml), dimnames, &
1292  chunksizes = chunksizes)
1293  else
1294  call register_field_wrap(fms2io_fileobj, var_name, this%get_var_skind(field_yaml), dimnames)
1295  endif
1296  else
1297  if (this%is_static()) then
1298  call register_field_wrap(fms2io_fileobj, var_name, this%get_var_skind(field_yaml))
1299  else
1300  !< In this case, the scalar variable is a function of time, so we need to pass in the
1301  !! unlimited dimension as a dimension
1302  call register_field_wrap(fms2io_fileobj, var_name, this%get_var_skind(field_yaml), (/unlim_dimname/))
1303  endif
1304  endif
1305 
1306  long_name = this%get_longname_to_write(field_yaml)
1307  call register_variable_attribute(fms2io_fileobj, var_name, "long_name", long_name, str_len=len_trim(long_name))
1308 
1309  units = this%get_units()
1310  if (units .ne. diag_null_string) &
1311  call register_variable_attribute(fms2io_fileobj, var_name, "units", units, str_len=len_trim(units))
1312 
1313  if (this%has_missing_value()) then
1314  call register_variable_attribute(fms2io_fileobj, var_name, "missing_value", &
1315  this%get_missing_value(field_yaml%get_var_kind()))
1316  call register_variable_attribute(fms2io_fileobj, var_name, "_FillValue", &
1317  this%get_missing_value(field_yaml%get_var_kind()))
1318  else
1319  call register_variable_attribute(fms2io_fileobj, var_name, "missing_value", &
1320  get_default_missing_value(field_yaml%get_var_kind()))
1321  call register_variable_attribute(fms2io_fileobj, var_name, "_FillValue", &
1322  get_default_missing_value(field_yaml%get_var_kind()))
1323  endif
1324 
1325  if (this%has_data_RANGE()) then
1326  call register_variable_attribute(fms2io_fileobj, var_name, "valid_range", &
1327  this%get_data_range(field_yaml%get_var_kind()))
1328  endif
1329 
1330  if (this%has_interp_method()) then
1331  interp_method_tmp = this%interp_method
1332  interp_method_len = len_trim(interp_method_tmp)
1333  call register_variable_attribute(fms2io_fileobj, var_name, "interp_method", interp_method_tmp, &
1334  str_len=interp_method_len)
1335  endif
1336 
1337  cell_methods = ""
1338  !< Check if any of the attributes defined via a "diag_field_add_attribute" call
1339  !! are the cell_methods, if so add to the "cell_methods" variable:
1340  do i = 1, this%num_attributes
1341  call this%attributes(i)%write_metadata(fms2io_fileobj, var_name, &
1342  cell_methods=cell_methods)
1343  enddo
1344 
1345  !< Append the time cell methods based on the variable's reduction
1346  call this%append_time_cell_methods(cell_methods, field_yaml)
1347  if (trim(cell_methods) .ne. "") &
1348  call register_variable_attribute(fms2io_fileobj, var_name, "cell_methods", &
1349  trim(adjustl(cell_methods)), str_len=len_trim(adjustl(cell_methods)))
1350 
1351  !< Write out the cell_measures attribute (i.e Area, Volume)
1352  !! The diag field ids for the Area and Volume are sent in the register call
1353  !! This was defined in file object and passed in here
1354  if (trim(cell_measures) .ne. "") &
1355  call register_variable_attribute(fms2io_fileobj, var_name, "cell_measures", &
1356  trim(adjustl(cell_measures)), str_len=len_trim(adjustl(cell_measures)))
1357 
1358  !< Write out the standard_name (this was defined in the register call)
1359  if (this%has_standname()) &
1360  call register_variable_attribute(fms2io_fileobj, var_name, "standard_name", &
1361  trim(this%get_standname()), str_len=len_trim(this%get_standname()))
1362 
1363  call this%write_coordinate_attribute(fms2io_fileobj, var_name, diag_axis)
1364 
1365  if (field_yaml%has_var_attributes()) then
1366  yaml_field_attributes = field_yaml%get_var_attributes()
1367  do i = 1, size(yaml_field_attributes,1)
1368  call register_variable_attribute(fms2io_fileobj, var_name, trim(yaml_field_attributes(i,1)), &
1369  trim(yaml_field_attributes(i,2)), str_len=len_trim(yaml_field_attributes(i,2)))
1370  enddo
1371  deallocate(yaml_field_attributes)
1372  endif
1373 end subroutine write_field_metadata
1374 
1375 !> @brief Determine the appropriate chunksizes for a diagnostic field based on its axes.
1376 !! For "X" and "Y" axes, the function returns a chunksize equal to the axis size divided by the layout.
1377 !! If the dimension is not evenly divisible by the layout, the function raises an error.
1378 !! For the other axis (i.e z axis) it return a chunksize equal to the axis length
1379 !! For sub-z axes (e.g., layer bounds), a chunksize of 1 is returned for now, as this case is not yet implemented.
1380 !! @return An integer array of chunksizes, one per diagnostic axis, plus one for the unlimited dimension.
1381 function get_chunksizes(this, diag_axis, field_yaml) &
1382  result(chunksizes)
1383 
1384  class(fmsdiagfield_type), target, intent(inout) :: this !< diag field
1385  class(fmsdiagaxiscontainer_type), target, intent(in) :: diag_axis(:) !< Diag_axis object
1386  type(diagyamlfilesvar_type), intent(in) :: field_yaml !< Field info from diag_table yaml
1387 
1388  integer, allocatable :: chunksizes(:)
1389 
1390  integer :: i !< For do loops
1391  integer :: ndim !< Number of spatial dimensions
1392  integer :: dim_size !< Dimensions size for the variable
1393  integer :: layout !< Layout to use for the variable
1394  integer :: specified_chunksizes(max_dimensions) !< Chunksizes specified in the yaml
1395 
1396  ndim = size(this%axis_ids)
1397  allocate(chunksizes(ndim + 1)) !! Adding 1 because of the unlimited dimension
1398  chunksizes = 1
1399  if (field_yaml%has_chunksizes()) then
1400  specified_chunksizes = field_yaml%get_chunksizes()
1401  chunksizes = specified_chunksizes(1:ndim+1)
1402  return
1403  endif
1404 
1405  ! Determine some default chunksizes to use, based on the compute domain
1406  do i = 1, ndim
1407  select type (axis => diag_axis(this%axis_ids(i))%axis)
1408  type is (fmsdiagfullaxis_type)
1409  if (axis%is_x_or_y_axis()) then
1410  call axis%get_dim_size_layout(dim_size, layout)
1411  if (mod(dim_size, layout) == 0) then
1412  chunksizes(i) = dim_size / layout
1413  else
1414  call mpp_error(fatal, "The variable "//field_yaml%get_var_varname()//" has a layout that is not"//&
1415  "evenly divisible by dimension size for axis "//axis%get_axis_name()//"."&
1416  "This may lead to poor performance when using collective writes. "//&
1417  "Consider using a different layout, disabling collective writes, "//&
1418  "or specifying chunksizes manually via the diag table yaml")
1419  endif
1420  else if (axis%is_z_axis() .and. field_yaml%has_var_zbounds()) then
1421  !TODO Handle edge case: chunking for z-axis with layer bounds (not yet implemented)
1422  else
1423  chunksizes(i) = axis%axis_length()
1424  endif
1425  end select
1426  enddo
1427 end function get_chunksizes
1428 
1429 !> @brief Writes the coordinate attribute of a field if any of the field's axis has an
1430 !! auxiliary axis
1431 subroutine write_coordinate_attribute (this, fms2io_fileobj, var_name, diag_axis)
1432  CLASS(fmsdiagfield_type), intent(in) :: this !< The field object
1433  class(fmsnetcdffile_t), INTENT(INOUT) :: fms2io_fileobj!< Fms2_io fileobj to write to
1434  character(len=*), intent(in) :: var_name !< Variable name
1435  class(fmsdiagaxiscontainer_type), intent(in) :: diag_axis(:) !< Diag_axis object
1436 
1437  integer :: i !< For do loops
1438  character(len = 252) :: aux_coord !< Auxuliary axis name
1439 
1440  !> If the variable is a scalar, go away
1441  if (.not. allocated(this%axis_ids)) return
1442 
1443  !> Determine if any of the field's axis has an auxiliary axis and the
1444  !! axis_names as a variable attribute
1445  aux_coord = ""
1446  do i = 1, size(this%axis_ids)
1447  select type (obj => diag_axis(this%axis_ids(i))%axis)
1448  type is (fmsdiagfullaxis_type)
1449  if (obj%has_aux()) then
1450  aux_coord = trim(aux_coord)//" "//obj%get_aux()
1451  endif
1452  end select
1453  enddo
1454 
1455  if (trim(aux_coord) .eq. "") return
1456 
1457  call register_variable_attribute(fms2io_fileobj, var_name, "coordinates", &
1458  trim(adjustl(aux_coord)), str_len=len_trim(adjustl(aux_coord)))
1459 
1460 end subroutine write_coordinate_attribute
1461 
1462 !> @brief Gets a fields data buffer
1463 !! @return a pointer to the data buffer
1464 function get_data_buffer (this) &
1465  result(rslt)
1466  class(fmsdiagfield_type), target, intent(in) :: this !< diag field
1467  class(*),dimension(:,:,:,:), pointer :: rslt !< The field's data buffer
1468 
1469  if (.not. this%data_buffer_is_allocated) &
1470  call mpp_error(fatal, "The input data buffer for the field:"&
1471  //trim(this%varname)//" was never allocated.")
1472 
1473  rslt => this%input_data_buffer%get_buffer()
1474 end function get_data_buffer
1475 
1476 
1477 !> @brief Gets a fields weight buffer
1478 !! @return a pointer to the weight buffer
1479 function get_weight (this) &
1480  result(rslt)
1481  class(fmsdiagfield_type), target, intent(in) :: this !< diag field
1482  type(real(kind=r8_kind)), pointer :: rslt
1483 
1484  if (.not. this%data_buffer_is_allocated) &
1485  call mpp_error(fatal, "The input data buffer for the field:"&
1486  //trim(this%varname)//" was never allocated.")
1487 
1488  rslt => this%input_data_buffer%get_weight()
1489 end function get_weight
1490 
1491 !> Gets the flag telling if the math functions need to be done
1492 !! \return Copy of math_needs_to_be_done flag
1493 pure logical function get_math_needs_to_be_done(this)
1494  class(fmsdiagfield_type), intent(in) :: this !< diag object
1495  get_math_needs_to_be_done = .false.
1496  if (allocated(this%math_needs_to_be_done)) get_math_needs_to_be_done = this%math_needs_to_be_done
1497 end function get_math_needs_to_be_done
1498 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1499 !!!!! Allocation checks
1500 
1501 !!> @brief Checks if obj%diag_field is allocated
1502 !!! @return true if obj%diag_field is allocated
1503 !logical function has_diag_field (obj)
1504 ! class (fmsDiagField_type), intent(in) :: obj !< diag object
1505 ! has_diag_field = allocated(obj%diag_field)
1506 !end function has_diag_field
1507 !> @brief Checks if obj%diag_id is allocated
1508 !! @return true if obj%diag_id is allocated
1509 pure logical function has_diag_id (this)
1510  class(fmsdiagfield_type), intent(in) :: this !< diag object
1511  has_diag_id = allocated(this%diag_id)
1512 end function has_diag_id
1513 
1514 !> @brief Checks if obj%metadata is allocated
1515 !! @return true if obj%metadata is allocated
1516 pure logical function has_attributes (this)
1517  class(fmsdiagfield_type), intent(in) :: this !< diag object
1518  has_attributes = this%num_attributes > 0
1519 end function has_attributes
1520 
1521 !> @brief Checks if obj%static is allocated
1522 !! @return true if obj%static is allocated
1523 pure logical function has_static (this)
1524  class(fmsdiagfield_type), intent(in) :: this !< diag object
1525  has_static = allocated(this%static)
1526 end function has_static
1527 
1528 !> @brief Checks if obj%registered is allocated
1529 !! @return true if obj%registered is allocated
1530 pure logical function has_registered (this)
1531  class(fmsdiagfield_type), intent(in) :: this !< diag object
1532  has_registered = allocated(this%registered)
1533 end function has_registered
1534 
1535 !> @brief Checks if obj%mask_variant is allocated
1536 !! @return true if obj%mask_variant is allocated
1537 pure logical function has_mask_variant (this)
1538  class(fmsdiagfield_type), intent(in) :: this !< diag object
1539  has_mask_variant = allocated(this%mask_variant)
1540 end function has_mask_variant
1541 
1542 !> @brief Checks if obj%local is allocated
1543 !! @return true if obj%local is allocated
1544 pure logical function has_local (this)
1545  class(fmsdiagfield_type), intent(in) :: this !< diag object
1546  has_local = allocated(this%local)
1547 end function has_local
1548 
1549 !> @brief Checks if obj%vartype is allocated
1550 !! @return true if obj%vartype is allocated
1551 pure logical function has_vartype (this)
1552  class(fmsdiagfield_type), intent(in) :: this !< diag object
1553  has_vartype = allocated(this%vartype)
1554 end function has_vartype
1555 
1556 !> @brief Checks if obj%varname is allocated
1557 !! @return true if obj%varname is allocated
1558 pure logical function has_varname (this)
1559  class(fmsdiagfield_type), intent(in) :: this !< diag object
1560  has_varname = allocated(this%varname)
1561 end function has_varname
1562 
1563 !> @brief Checks if obj%longname is allocated
1564 !! @return true if obj%longname is allocated
1565 pure logical function has_longname (this)
1566  class(fmsdiagfield_type), intent(in) :: this !< diag object
1567  has_longname = allocated(this%longname)
1568 end function has_longname
1569 
1570 !> @brief Checks if obj%standname is allocated
1571 !! @return true if obj%standname is allocated
1572 pure logical function has_standname (this)
1573  class(fmsdiagfield_type), intent(in) :: this !< diag object
1574  has_standname = allocated(this%standname)
1575 end function has_standname
1576 
1577 !> @brief Checks if obj%units is allocated
1578 !! @return true if obj%units is allocated
1579 pure logical function has_units (this)
1580  class(fmsdiagfield_type), intent(in) :: this !< diag object
1581  has_units = allocated(this%units)
1582 end function has_units
1583 
1584 !> @brief Checks if obj%modname is allocated
1585 !! @return true if obj%modname is allocated
1586 pure logical function has_modname (this)
1587  class(fmsdiagfield_type), intent(in) :: this !< diag object
1588  has_modname = allocated(this%modname)
1589 end function has_modname
1590 
1591 !> @brief Checks if obj%realm is allocated
1592 !! @return true if obj%realm is allocated
1593 pure logical function has_realm (this)
1594  class(fmsdiagfield_type), intent(in) :: this !< diag object
1595  has_realm = allocated(this%realm)
1596 end function has_realm
1597 
1598 !> @brief Checks if obj%interp_method is allocated
1599 !! @return true if obj%interp_method is allocated
1600 pure logical function has_interp_method (this)
1601  class(fmsdiagfield_type), intent(in) :: this !< diag object
1602  has_interp_method = allocated(this%interp_method)
1603 end function has_interp_method
1604 
1605 !> @brief Checks if obj%frequency is allocated
1606 !! @return true if obj%frequency is allocated
1607 pure logical function has_frequency (this)
1608  class(fmsdiagfield_type), intent(in) :: this !< diag object
1609  has_frequency = allocated(this%frequency)
1610 end function has_frequency
1611 
1612 !> @brief Checks if obj%tile_count is allocated
1613 !! @return true if obj%tile_count is allocated
1614 pure logical function has_tile_count (this)
1615  class(fmsdiagfield_type), intent(in) :: this !< diag object
1616  has_tile_count = allocated(this%tile_count)
1617 end function has_tile_count
1618 
1619 !> @brief Checks if axis_ids of the object is allocated
1620 !! @return true if it is allocated
1621 pure logical function has_axis_ids (this)
1622  class(fmsdiagfield_type), intent(in) :: this !< diag field object
1623  has_axis_ids = allocated(this%axis_ids)
1624 end function has_axis_ids
1625 
1626 !> @brief Checks if obj%area is allocated
1627 !! @return true if obj%area is allocated
1628 pure logical function has_area (this)
1629  class(fmsdiagfield_type), intent(in) :: this !< diag object
1630  has_area = allocated(this%area)
1631 end function has_area
1632 
1633 !> @brief Checks if obj%volume is allocated
1634 !! @return true if obj%volume is allocated
1635 pure logical function has_volume (this)
1636  class(fmsdiagfield_type), intent(in) :: this !< diag object
1637  has_volume = allocated(this%volume)
1638 end function has_volume
1639 
1640 !> @brief Checks if obj%missing_value is allocated
1641 !! @return true if obj%missing_value is allocated
1642 pure logical function has_missing_value (this)
1643  class(fmsdiagfield_type), intent(in) :: this !< diag object
1644  has_missing_value = allocated(this%missing_value)
1645 end function has_missing_value
1646 
1647 !> @brief Checks if obj%data_RANGE is allocated
1648 !! @return true if obj%data_RANGE is allocated
1649 pure logical function has_data_range (this)
1650  class(fmsdiagfield_type), intent(in) :: this !< diag object
1651  has_data_range = allocated(this%data_RANGE)
1652 end function has_data_range
1653 
1654 !> @brief Checks if obj%input_data_buffer is allocated
1655 !! @return true if obj%input_data_buffer is allocated
1656 pure logical function has_input_data_buffer (this)
1657  class(fmsdiagfield_type), intent(in) :: this !< diag object
1658  has_input_data_buffer = allocated(this%input_data_buffer)
1659 end function has_input_data_buffer
1660 
1661 !> @brief Add a attribute to the diag_obj using the diag_field_id
1662 subroutine diag_field_add_attribute(this, att_name, att_value)
1663  class(fmsdiagfield_type), intent (inout) :: this !< The field object
1664  character(len=*), intent(in) :: att_name !< Name of the attribute
1665  class(*), intent(in) :: att_value(:) !< The attribute value to add
1666 
1667  this%num_attributes = this%num_attributes + 1
1668  if (this%num_attributes > max_field_attributes) &
1669  call mpp_error(fatal, "diag_field_add_attribute: Number of attributes exceeds max_field_attributes for field:"&
1670  //trim(this%varname)//". Increase diag_manager_nml:max_field_attributes.")
1671 
1672  call this%attributes(this%num_attributes)%add(att_name, att_value)
1673 end subroutine diag_field_add_attribute
1674 
1675 !> @brief Determine the default missing value to use based on the requested variable type
1676 !! @return The missing value
1677 function get_default_missing_value(var_type) &
1678  result(rslt)
1679 
1680  integer, intent(in) :: var_type !< The type of the variable to return the missing value as
1681  class(*),allocatable :: rslt
1682 
1683  select case(var_type)
1684  case (r4)
1685  allocate(real(kind=r4_kind) :: rslt)
1686  rslt = real(cmor_missing_value, kind=r4_kind)
1687  case (r8)
1688  allocate(real(kind=r8_kind) :: rslt)
1689  rslt = real(cmor_missing_value, kind=r8_kind)
1690  case default
1691  end select
1692 end function
1693 
1694 !> @brief Determines the diag_obj id corresponding to a module name and field_name
1695 !> @return diag_obj id
1696 FUNCTION diag_field_id_from_name(this, module_name, field_name) &
1697  result(diag_field_id)
1698  CLASS(fmsdiagfield_type), INTENT(in) :: this !< The field object
1699  CHARACTER(len=*), INTENT(in) :: module_name !< Module name that registered the variable
1700  CHARACTER(len=*), INTENT(in) :: field_name !< Variable name
1701 
1702  integer :: diag_field_id
1703 
1704  diag_field_id = diag_field_not_found
1705  if (this%get_varname() .eq. trim(field_name) .and. &
1706  this%get_modname() .eq. trim(module_name)) then
1707  diag_field_id = this%get_id()
1708  endif
1709 end function diag_field_id_from_name
1710 
1711 !> @brief Adds the area and volume id to a field object
1712 subroutine add_area_volume(this, area, volume)
1713  CLASS(fmsdiagfield_type), intent(inout) :: this !< The field object
1714  INTEGER, optional, INTENT(in) :: area !< diag ids of area
1715  INTEGER, optional, INTENT(in) :: volume !< diag ids of volume
1716 
1717  if (present(area)) then
1718  if (area > 0) then
1719  this%area = area
1720  else
1721  call mpp_error(fatal, "diag_field_add_cell_measures: the area id is not valid. &
1722  &Verify that the area_id passed in to the field:"//this%varname// &
1723  " is valid and that the field is registered and in the diag_table.yaml")
1724  endif
1725  endif
1726 
1727  if (present(volume)) then
1728  if (volume > 0) then
1729  this%volume = volume
1730  else
1731  call mpp_error(fatal, "diag_field_add_cell_measures: the volume id is not valid. &
1732  &Verify that the volume_id passed in to the field:"//this%varname// &
1733  " is valid and that the field is registered and in the diag_table.yaml")
1734  endif
1735  endif
1736 
1737 end subroutine add_area_volume
1738 
1739 !> @brief Append the time cell meathods based on the variable's reduction
1740 subroutine append_time_cell_methods(this, cell_methods, field_yaml)
1741  class(fmsdiagfield_type), target, intent(inout) :: this !< diag field
1742  character(len=*), intent(inout) :: cell_methods !< The cell methods var to append to
1743  type(diagyamlfilesvar_type), intent(in) :: field_yaml !< The field's yaml
1744 
1745  if (this%static) then
1746  cell_methods = trim(cell_methods)//" time: point "
1747  return
1748  endif
1749 
1750  select case (field_yaml%get_var_reduction())
1751  case (time_none)
1752  cell_methods = trim(cell_methods)//" time: point "
1753  case (time_diurnal)
1754  cell_methods = trim(cell_methods)//" time: mean"
1755  case (time_power)
1756  cell_methods = trim(cell_methods)//" time: mean_pow"//int2str(field_yaml%get_pow_value())
1757  case (time_rms)
1758  cell_methods = trim(cell_methods)//" time: root_mean_square"
1759  case (time_max)
1760  cell_methods = trim(cell_methods)//" time: max"
1761  case (time_min)
1762  cell_methods = trim(cell_methods)//" time: min"
1763  case (time_average)
1764  cell_methods = trim(cell_methods)//" time: mean"
1765  case (time_sum)
1766  cell_methods = trim(cell_methods)//" time: sum"
1767  end select
1768 end subroutine append_time_cell_methods
1769 
1770 !> Dumps any data from a given fmsDiagField_type object
1771 subroutine dump_field_obj (this, unit_num)
1772  class(fmsdiagfield_type), intent(in) :: this
1773  integer, intent(in) :: unit_num !< passed in from dump_diag_obj if log file is being written to
1774  integer :: i
1775 
1776  if( mpp_pe() .eq. mpp_root_pe()) then
1777  if( allocated(this%file_ids)) write(unit_num, *) 'file_ids:' ,this%file_ids
1778  if( allocated(this%diag_id)) write(unit_num, *) 'diag_id:' ,this%diag_id
1779  if( allocated(this%static)) write(unit_num, *) 'static:' ,this%static
1780  if( allocated(this%registered)) write(unit_num, *) 'registered:' ,this%registered
1781  if( allocated(this%mask_variant)) write(unit_num, *) 'mask_variant:' ,this%mask_variant
1782  if( allocated(this%do_not_log)) write(unit_num, *) 'do_not_log:' ,this%do_not_log
1783  if( allocated(this%local)) write(unit_num, *) 'local:' ,this%local
1784  if( allocated(this%vartype)) write(unit_num, *) 'vartype:' ,this%vartype
1785  if( allocated(this%varname)) write(unit_num, *) 'varname:' ,this%varname
1786  if( allocated(this%longname)) write(unit_num, *) 'longname:' ,this%longname
1787  if( allocated(this%standname)) write(unit_num, *) 'standname:' ,this%standname
1788  if( allocated(this%units)) write(unit_num, *) 'units:' ,this%units
1789  if( allocated(this%modname)) write(unit_num, *) 'modname:' ,this%modname
1790  if( allocated(this%realm)) write(unit_num, *) 'realm:' ,this%realm
1791  if( allocated(this%interp_method)) write(unit_num, *) 'interp_method:' ,this%interp_method
1792  if( allocated(this%tile_count)) write(unit_num, *) 'tile_count:' ,this%tile_count
1793  if( allocated(this%axis_ids)) write(unit_num, *) 'axis_ids:' ,this%axis_ids
1794  write(unit_num, *) 'type_of_domain:' ,this%type_of_domain
1795  if( allocated(this%area)) write(unit_num, *) 'area:' ,this%area
1796  if( allocated(this%missing_value)) then
1797  select type(missing_val => this%missing_value)
1798  type is (real(r4_kind))
1799  write(unit_num, *) 'missing_value:', missing_val
1800  type is (real(r8_kind))
1801  write(unit_num, *) 'missing_value:' ,missing_val
1802  type is(integer(i4_kind))
1803  write(unit_num, *) 'missing_value:' ,missing_val
1804  type is(integer(i8_kind))
1805  write(unit_num, *) 'missing_value:' ,missing_val
1806  end select
1807  endif
1808  if( allocated( this%data_RANGE)) then
1809  select type(drange => this%data_RANGE)
1810  type is (real(r4_kind))
1811  write(unit_num, *) 'data_RANGE:' ,drange
1812  type is (real(r8_kind))
1813  write(unit_num, *) 'data_RANGE:' ,drange
1814  type is(integer(i4_kind))
1815  write(unit_num, *) 'data_RANGE:' ,drange
1816  type is(integer(i8_kind))
1817  write(unit_num, *) 'data_RANGE:' ,drange
1818  end select
1819  endif
1820  write(unit_num, *) 'num_attributes:' ,this%num_attributes
1821  if( allocated(this%attributes)) then
1822  do i=1, this%num_attributes
1823  if( allocated(this%attributes(i)%att_value)) then
1824  select type( val => this%attributes(i)%att_value)
1825  type is (real(r8_kind))
1826  write(unit_num, *) 'attribute name', this%attributes(i)%att_name, 'val:', val
1827  type is (real(r4_kind))
1828  write(unit_num, *) 'attribute name', this%attributes(i)%att_name, 'val:', val
1829  type is (integer(i4_kind))
1830  write(unit_num, *) 'attribute name', this%attributes(i)%att_name, 'val:', val
1831  type is (integer(i8_kind))
1832  write(unit_num, *) 'attribute name', this%attributes(i)%att_name, 'val:', val
1833  end select
1834  endif
1835  enddo
1836  endif
1837 
1838  endif
1839 
1840 end subroutine
1841 
1842 !< @brief Get the starting compute domain indices for a set of axis
1843 !! @return compute domain starting indices
1844 function get_starting_compute_domain(axis_ids, diag_axis) &
1845 result(compute_domain)
1846  integer, intent(in) :: axis_ids(:) !< Array of axis ids
1847  class(fmsdiagaxiscontainer_type),intent(in) :: diag_axis(:) !< Array of axis object
1848 
1849  integer :: compute_domain(4)
1850  integer :: a !< For looping through axes
1851  integer :: compute_idx(2) !< Compute domain indices (starting, ending)
1852  logical :: dummy !< Dummy variable for the `get_compute_domain` subroutine
1853 
1854  compute_domain = 1
1855  axis_loop: do a = 1,size(axis_ids)
1856  select type (axis => diag_axis(axis_ids(a))%axis)
1857  type is (fmsdiagfullaxis_type)
1858  call axis%get_compute_domain(compute_idx, dummy)
1859  if ( compute_idx(1) .ne. diag_null) compute_domain(a) = compute_idx(1)
1860  end select
1861  enddo axis_loop
1862 end function get_starting_compute_domain
1863 
1864 !> Get list of field ids
1865 pure function get_file_ids(this)
1866  class(fmsdiagfield_type), intent(in) :: this
1867  integer, allocatable :: get_file_ids(:) !< Ids of the FMS_diag_files the variable
1868  get_file_ids = this%file_ids
1869 end function
1870 
1871 !> @brief Get the mask from the input buffer object
1872 !! @return a pointer to the mask
1873 function get_mask(this)
1874  class(fmsdiagfield_type), target, intent(in) :: this !< input buffer object
1875  logical, pointer :: get_mask(:,:,:,:)
1876  get_mask => this%mask
1877 end function get_mask
1878 
1879 !> @brief If in openmp region, omp_axis should be provided in order to allocate to the given axis lengths.
1880 !! Otherwise mask will be allocated to the size of mask_in
1881 subroutine allocate_mask(this, mask_in, omp_axis)
1882  class(fmsdiagfield_type), target, intent(inout) :: this !< input buffer object
1883  logical, intent(in) :: mask_in(:,:,:,:)
1884  class(fmsdiagaxiscontainer_type), intent(in), optional :: omp_axis(:) !< true if calling from omp region
1885  integer :: axis_num, length(4)
1886  integer, pointer :: id_num
1887  ! if not omp just allocate to whatever is given
1888  if(.not. present(omp_axis)) then
1889  allocate(this%mask(size(mask_in,1), size(mask_in,2), size(mask_in,3), &
1890  size(mask_in,4)))
1891  ! otherwise loop through axis and get sizes
1892  else
1893  length = 1
1894  do axis_num=1, size(this%axis_ids)
1895  id_num => this%axis_ids(axis_num)
1896  select type(axis => omp_axis(id_num)%axis)
1897  type is (fmsdiagfullaxis_type)
1898  length(axis_num) = axis%axis_length()
1899  end select
1900  enddo
1901  allocate(this%mask(length(1), length(2), length(3), length(4)))
1902  endif
1903 end subroutine allocate_mask
1904 
1905 !> Sets previously allocated mask to mask_in at given index ranges
1906 subroutine set_mask(this, mask_in, field_info, is, js, ks, ie, je, ke)
1907  class(fmsdiagfield_type), intent(inout) :: this
1908  logical, intent(in) :: mask_in(:,:,:,:)
1909  character(len=*), intent(in) :: field_info !< Field info to add to error message
1910  integer, optional, intent(in) :: is, js, ks, ie, je, ke
1911  if(present(is)) then
1912  if(is .lt. lbound(this%mask,1) .or. ie .gt. ubound(this%mask,1) .or. &
1913  js .lt. lbound(this%mask,2) .or. je .gt. ubound(this%mask,2) .or. &
1914  ks .lt. lbound(this%mask,3) .or. ke .gt. ubound(this%mask,3)) then
1915  print *, "PE:", int2str(mpp_pe()), "The size of the mask is", &
1916  shape(this%mask), &
1917  "But the indices passed in are is=", int2str(is), " ie=", int2str(ie),&
1918  " js=", int2str(js), " je=", int2str(je), &
1919  " ks=", int2str(ks), " ke=", int2str(ke), &
1920  " ", trim(field_info)
1921  call mpp_error(fatal,"set_mask:: given indices out of bounds for allocated mask")
1922  endif
1923  this%mask(is:ie, js:je, ks:ke, :) = mask_in
1924  else
1925  this%mask = mask_in
1926  endif
1927 end subroutine set_mask
1928 
1929 !> sets halo_present to true
1930 subroutine set_halo_present(this)
1931  class(fmsdiagfield_type), intent(inout) :: this !< field object to modify
1932  this%halo_present = .true.
1933 end subroutine set_halo_present
1934 
1935 !> Getter for halo_present
1936 pure function is_halo_present(this)
1937  class(fmsdiagfield_type), intent(in) :: this !< field object to get from
1938  logical :: is_halo_present
1939  is_halo_present = this%halo_present
1940 end function is_halo_present
1941 
1942 !> Helper routine to find and set the netcdf missing value for a field
1943 !! Always returns r8 due to reduction routine args
1944 !! casts up to r8 from given missing val or default if needed
1945 function find_missing_value(this, missing_val) &
1946  result(res)
1947  class(fmsdiagfield_type), intent(in) :: this !< field object to get missing value for
1948  class(*), allocatable, intent(out) :: missing_val !< outputted netcdf missing value (oriignal type)
1949  real(r8_kind), allocatable :: res !< returned r8 copy of missing_val
1950  integer :: vtype !< temp to hold enumerated variable type
1951 
1952  if(this%has_missing_value()) then
1953  missing_val = this%get_missing_value(this%get_vartype())
1954  else
1955  vtype = this%get_vartype()
1956  if(vtype .eq. r8) then
1957  missing_val = cmor_missing_value
1958  else
1959  missing_val = real(cmor_missing_value, r4_kind)
1960  endif
1961  endif
1962 
1963  select type(missing_val)
1964  type is (real(r8_kind))
1965  res = missing_val
1966  type is (real(r4_kind))
1967  res = real(missing_val, r8_kind)
1968  end select
1969 end function find_missing_value
1970 
1971 !> @returns allocation status of logical mask array
1972 !! this just indicates whether the mask array itself has been alloc'd
1973 !! this is different from @ref has_mask_variant, which is set earlier for whether a mask is being used at all
1974 pure logical function has_mask_allocated(this)
1975  class(fmsdiagfield_type),intent(in) :: this !< field object to check mask allocation for
1976  has_mask_allocated = allocated(this%mask)
1977 end function has_mask_allocated
1978 
1979 !> @brief Determine if the variable is in the file
1980 !! @return .True. if the varibale is in the file
1981 pure function is_variable_in_file(this, file_id) &
1982 result(res)
1983  class(fmsdiagfield_type), intent(in) :: this !< field object to check
1984  integer, intent(in) :: file_id !< File id to check
1985  logical :: res
1986 
1987  integer :: i
1988 
1989  res = .false.
1990  if (any(this%file_ids .eq. file_id)) res = .true.
1991 end function is_variable_in_file
1992 
1993 !> @brief Determine the name of the first file the variable is in
1994 !! @return filename
1995 function get_field_file_name(this) &
1996  result(res)
1997  class(fmsdiagfield_type), intent(in) :: this !< Field object to query
1998  character(len=:), allocatable :: res
1999 
2000  res = this%diag_field(1)%get_var_fname()
2001 end function get_field_file_name
2002 
2003 !> @brief Generate the associated files attribute
2004 subroutine generate_associated_files_att(this, att, start_time)
2005  class(fmsdiagfield_type) , intent(in) :: this !< diag_field_object for the area/volume field
2006  character(len=*), intent(inout) :: att !< associated_files_att
2007  type(time_type), intent(in) :: start_time !< The start_time for the field's file
2008 
2009  character(len=:), allocatable :: field_name !< Name of the area/volume field
2010  character(len=FMS_FILE_LEN) :: file_name !< Name of the file the area/volume field is in!
2011  character(len=128) :: start_date !< Start date to append to the begining of the filename
2012 
2013  integer :: year, month, day, hour, minute, second
2014 
2015  file_name = this%get_field_file_name()
2016  field_name = this%get_varname(to_write = .true., filename=file_name)
2017 
2018  ! Check if the field is already in the associated files attribute (i.e the area can be associated with multiple
2019  ! fields in the file, but it only needs to be added once)
2020  if (index(att, field_name) .ne. 0) return
2021 
2022  if (prepend_date) then
2023  call get_date(start_time, year, month, day, hour, minute, second)
2024  write (start_date, '(1I20.4, 2I2.2)') year, month, day
2025  file_name = trim(adjustl(start_date))//'.'//trim(file_name)
2026  endif
2027 
2028  att = trim(att)//" "//trim(field_name)//": "//trim(file_name)//".nc"
2029 end subroutine generate_associated_files_att
2030 
2031 !> @brief Determines if the compute domain has been divide further into slices (i.e openmp blocks)
2032 !! @return .True. if the compute domain has been divided furter into slices
2033 function check_for_slices(field, diag_axis, var_size) &
2034  result(rslt)
2035  type(fmsdiagfield_type), intent(in) :: field !< Field object
2036  type(fmsdiagaxiscontainer_type), target, intent(in) :: diag_axis(:) !< Array of diag axis
2037  integer, intent(in) :: var_size(:) !< The size of the buffer pass into send_data
2038 
2039  logical :: rslt
2040  integer :: i !< For do loops
2041 
2042  rslt = .false.
2043 
2044  if (.not. field%has_axis_ids()) then
2045  rslt = .false.
2046  return
2047  endif
2048  do i = 1, size(field%axis_ids)
2049  select type (axis_obj => diag_axis(field%axis_ids(i))%axis)
2050  type is (fmsdiagfullaxis_type)
2051  if (axis_obj%axis_length() .ne. var_size(i)) then
2052  rslt = .true.
2053  return
2054  endif
2055  end select
2056  enddo
2057 end function
2058 #endif
2059 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.
subroutine, public find_z_sub_axis_name(dim_name, parent_axis_id, file_axis_id, field_yaml, diag_axis)
Determine the name of the z subaxis by matching the parent axis id and the zbounds in the diag table ...
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