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