FMS  2025.02.01
Flexible Modeling System
fms_diag_file_object.F90
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 !> @defgroup fms_diag_output_yaml_mod fms_diag_output_yaml_mod
20 !> @ingroup diag_manager
21 !! @brief fms_diag_file_object_mod handles the file objects data, functions, and subroutines.
22 !! @author Tom Robinson
23 !! @description The fmsDiagFile_type contains the information for each history file to be written. It has
24 !! a pointer to the information from the diag yaml, additional metadata that comes from the model, and a
25 !! list of the variables and their variable IDs that are in the file.
26 module fms_diag_file_object_mod
27 #ifdef use_yaml
28 use fms2_io_mod, only: fmsnetcdffile_t, fmsnetcdfunstructureddomainfile_t, fmsnetcdfdomainfile_t, &
29  get_instance_filename, open_file, close_file, get_mosaic_tile_file, unlimited, &
30  register_axis, register_field, register_variable_attribute, write_data, &
31  dimension_exists, register_global_attribute, flush_file
32 use diag_data_mod, only: diag_null, no_domain, max_axes, sub_regional, get_base_time, diag_not_registered, &
33  two_d_domain, ug_domain, prepend_date, diag_days, very_large_file_freq, &
37  middle_time, begin_time, end_time, max_str_len, index_gridtype, latlon_gridtype, &
38  null_gridtype, flush_nc_files, diag_init_time
39 use time_manager_mod, only: time_type, operator(>), operator(/=), operator(==), get_date, get_calendar_type, &
40  valid_calendar_types, operator(>=), date_to_string, &
41  OPERATOR(/), OPERATOR(+), operator(<)
42 use fms_diag_time_utils_mod, only: diag_time_inc, get_time_string, get_date_dif
43 use fms_diag_yaml_mod, only: diag_yaml, diagyamlobject_type, diagyamlfiles_type, subregion_type, diagyamlfilesvar_type
44 use fms_diag_axis_object_mod, only: diagdomain_t, get_domain_and_domain_type, fmsdiagaxis_type, &
49 use fms_diag_field_object_mod, only: fmsdiagfield_type
51 use mpp_mod, only: mpp_get_current_pelist, mpp_npes, mpp_root_pe, mpp_pe, mpp_error, fatal, stdout, &
52  uppercase, lowercase, note, mpp_max
53 use platform_mod, only: fms_file_len
54 use mpp_domains_mod, only: mpp_get_ntile_count, mpp_get_ug_domain_ntiles, mpp_get_io_domain_layout, &
55  mpp_get_io_domain_ug_layout
56 implicit none
57 private
58 
60 public :: fmsdiagfile_type, fms_diag_files_object_init, fms_diag_files_object_initialized
61 
62 logical :: fms_diag_files_object_initialized = .false.
63 
64 integer, parameter :: var_string_len = 25
65 
67  private
68  integer :: id !< The number associated with this file in the larger array of files
69  TYPE(time_type) :: model_time !< The last time data was sent for any of the buffers in this file object
70  TYPE(time_type) :: start_time !< The start time for the file
71  TYPE(time_type) :: last_output !< Time of the last time output was writen
72  TYPE(time_type) :: next_output !< Time of the next write
73  TYPE(time_type) :: next_next_output !< Time of the next next write
74  TYPE(time_type) :: no_more_data !< Time to stop receiving data for this file
75  logical :: done_writing_data!< .True. if finished writing data
76 
77  !< This will be used when using the new_file_freq keys in the diag_table.yaml
78  TYPE(time_type) :: next_close !< Time to close the file
79  logical :: is_file_open !< .True. if the file is opened
80 
81  class(fmsnetcdffile_t), allocatable :: fms2io_fileobj !< fms2_io file object for this history file
82  type(diagyamlfiles_type), pointer :: diag_yaml_file => null() !< Pointer to the diag_yaml_file data
83  integer :: type_of_domain !< The type of domain to use to open the file
84  !! NO_DOMAIN, TWO_D_DOMAIN, UG_DOMAIN, SUB_REGIONAL
85  class(diagdomain_t), pointer :: domain !< The domain to use,
86  !! null if NO_DOMAIN or SUB_REGIONAL
87  character(len=:) , dimension(:), allocatable :: file_metadata_from_model !< File metadata that comes from
88  !! the model.
89  integer, dimension(:), allocatable :: field_ids !< Variable IDs corresponding to file_varlist
90  integer, dimension(:), allocatable :: yaml_ids !< IDs corresponding to the yaml field section
91  logical, dimension(:), private, allocatable :: field_registered !< Array corresponding to `field_ids`, .true.
92  !! if the variable has been registered and
93  !! `field_id` has been set for the variable
94  integer, allocatable :: num_registered_fields !< The number of fields registered
95  !! to the file
96  integer, dimension(:), allocatable :: axis_ids !< Array of axis ids in the file
97  integer :: number_of_axis !< Number of axis in the file
98  integer, dimension(:), allocatable :: buffer_ids !< array of buffer ids associated with the file
99  integer :: number_of_buffers !< Number of buffers that have been added to the file
100  logical, allocatable :: time_ops !< .True. if file contains variables that are time_min, time_max, time_average,
101  !! or time_sum
102  integer :: unlim_dimension_level !< The unlimited dimension level currently being written
103  integer :: num_time_levels !< The number of time levels that were actually written to the file
104  logical :: data_has_been_written !< .True. if data has been written for the current unlimited dimension level
105  logical :: is_static !< .True. if the frequency is -1
106  integer :: nz_subaxis !< The number of Z axis currently added to the file
107 
108  contains
109  procedure, public :: add_field_and_yaml_id
110  procedure, public :: add_buffer_id
111  procedure, public :: is_field_registered
112  procedure, public :: init_diurnal_axis
113  procedure, public :: has_file_metadata_from_model
114  procedure, public :: has_fileobj
115  procedure, public :: has_diag_yaml_file
116  procedure, public :: set_domain_from_axis
117  procedure, public :: set_file_domain
118  procedure, public :: add_axes
119  procedure, public :: add_new_axis
120  procedure, public :: update_write_on_this_pe
121  procedure, public :: get_write_on_this_pe
122  procedure, public :: does_axis_exist
123  procedure, public :: define_new_subaxis
124  procedure, public :: add_start_time
125  procedure, public :: set_file_time_ops
126  procedure, public :: has_field_ids
127  procedure, public :: get_time_ops
128  procedure, public :: get_id
129 ! TODO procedure, public :: get_fileobj ! TODO
130 ! TODO procedure, public :: get_diag_yaml_file ! TODO
131  procedure, public :: get_file_metadata_from_model
132  procedure, public :: get_field_ids
133 ! The following fuctions come will use the yaml inquiry functions
134  procedure, public :: get_file_fname
135  procedure, public :: get_file_frequnit
136  procedure, public :: get_file_freq
137  procedure, public :: get_file_timeunit
138  procedure, public :: get_file_unlimdim
139  procedure, public :: get_file_sub_region
140  procedure, public :: get_file_sub_region_grid_type
141  procedure, public :: get_file_new_file_freq
142  procedure, public :: get_filename_time
143  procedure, public :: get_file_new_file_freq_units
144  procedure, public :: get_file_start_time
145  procedure, public :: get_file_duration
146  procedure, public :: get_file_duration_units
147  procedure, public :: get_file_varlist
148  procedure, public :: get_file_global_meta
149  procedure, public :: is_done_writing_data
150  procedure, public :: has_file_fname
151  procedure, public :: has_file_frequnit
152  procedure, public :: has_file_freq
153  procedure, public :: has_file_timeunit
154  procedure, public :: has_file_unlimdim
155  procedure, public :: has_file_sub_region
156  procedure, public :: has_file_new_file_freq
157  procedure, public :: has_file_new_file_freq_units
158  procedure, public :: has_file_start_time
159  procedure, public :: has_file_duration
160  procedure, public :: has_file_duration_units
161  procedure, public :: has_file_varlist
162  procedure, public :: has_file_global_meta
163  procedure, public :: dump_file_obj
164  procedure, public :: get_buffer_ids
165  procedure, public :: get_number_of_buffers
166  procedure, public :: has_send_data_been_called
167  procedure, public :: check_buffer_times
168 end type fmsdiagfile_type
169 
171  integer, dimension(:), allocatable :: sub_axis_ids !< Array of axis ids in the file
172  logical :: write_on_this_pe !< Flag indicating if the subregion is on the current PE
173  logical :: is_subaxis_defined !< Flag indicating if the subaxes have already been defined
174 end type subregionalfile_type
175 
176 !> \brief A container for fmsDiagFile_type. This is used to create the array of files
178  class(fmsdiagfile_type),allocatable :: fms_diag_file !< The individual file object
179 
180  contains
181  procedure :: is_regional
182  procedure :: is_file_static
183  procedure :: open_diag_file
184  procedure :: write_global_metadata
185  procedure :: write_time_metadata
186  procedure :: write_field_data
187  procedure :: write_axis_metadata
188  procedure :: write_field_metadata
189  procedure :: write_axis_data
190  procedure :: writing_on_this_pe
191  procedure :: check_file_times
192  procedure :: is_time_to_close_file
193  procedure :: write_time_data
194  procedure :: update_next_write
195  procedure :: prepare_for_force_write
196  procedure :: init_unlim_dim
197  procedure :: update_current_new_file_freq_index
198  procedure :: get_unlim_dimension_level
199  procedure :: get_num_time_levels
200  procedure :: get_num_tiles
201  procedure :: get_ndistributedfiles
202  procedure :: flush_diag_file
203  procedure :: get_next_output
204  procedure :: get_next_next_output
205  procedure :: close_diag_file
206  procedure :: set_model_time
207  procedure :: get_model_time
208  procedure :: time_to_start_doing_math
210 
211 !type(fmsDiagFile_type), dimension (:), allocatable, target :: FMS_diag_file !< The array of diag files
212 !class(fmsDiagFileContainer_type),dimension (:), allocatable, target :: FMS_diag_file
213 
214 contains
215 
216 !< @brief Allocates the number of files and sets an ID based for each file
217 !! @return true if there are files allocated in the YAML object
218 logical function fms_diag_files_object_init (files_array)
219  class(fmsdiagfilecontainer_type), allocatable, target, intent(inout) :: files_array (:) !< array of diag files
220  class(fmsdiagfile_type), pointer :: obj => null() !< Pointer for each member of the array
221  integer :: nfiles !< Number of files in the diag yaml
222  integer :: i !< Looping iterator
223  if (diag_yaml%has_diag_files()) then
224  nfiles = diag_yaml%size_diag_files()
225  allocate (files_array(nfiles))
226  set_ids_loop: do i= 1,nfiles
227  !> If the file has a sub_regional, define it as one and allocate the sub_axis_ids array.
228  !! This will be set in a add_axes
229  if (diag_yaml%diag_files(i)%has_file_sub_region()) then
230  allocate(subregionalfile_type :: files_array(i)%FMS_diag_file)
231  obj => files_array(i)%FMS_diag_file
232  select type (obj)
233  type is (subregionalfile_type)
234  allocate(obj%sub_axis_ids(max_axes))
235  obj%sub_axis_ids = diag_null
236  obj%write_on_this_pe = .true.
237  obj%is_subaxis_defined = .false.
238  obj%number_of_axis = 0
239  end select
240  else
241  allocate(fmsdiagfile_type::files_array(i)%FMS_diag_file)
242  obj => files_array(i)%FMS_diag_file
243  endif
244  !!
245  obj%diag_yaml_file => diag_yaml%diag_files(i)
246  obj%id = i
247  allocate(obj%field_ids(diag_yaml%diag_files(i)%size_file_varlist()))
248  allocate(obj%buffer_ids(diag_yaml%diag_files(i)%size_file_varlist()))
249  allocate(obj%yaml_ids(diag_yaml%diag_files(i)%size_file_varlist()))
250  allocate(obj%field_registered(diag_yaml%diag_files(i)%size_file_varlist()))
251  !! Initialize the integer arrays
252  obj%field_ids = diag_not_registered
253  obj%yaml_ids = diag_not_registered
254  obj%buffer_ids = diag_not_registered
255  obj%field_registered = .false.
256  obj%num_registered_fields = 0
257  obj%number_of_buffers = 0
258 
259  !> These will be set in a set_file_domain
260  obj%type_of_domain = no_domain
261  obj%domain => null()
262 
263  !> This will be set in a add_axes
264  allocate(obj%axis_ids(max_axes))
265  obj%number_of_axis = 0
266 
267  !> Set the start_time of the file to the base_time and set up the *_output variables
268  obj%done_writing_data = .false.
269 
270  !! Set this to the time passed in to diag_manager_init
271  !! This will be the base_time if nothing was passed in
272  !! This time is appended to the filename if the prepend_date namelist is .True.
273  if (obj%has_file_start_time()) then
274  obj%start_time = obj%get_file_start_time()
275  else
276  obj%start_time = diag_init_time
277  endif
278  obj%last_output = obj%start_time
279  obj%model_time = diag_init_time
280  obj%next_output = diag_time_inc(obj%start_time, obj%get_file_freq(), obj%get_file_frequnit())
281  obj%next_next_output = diag_time_inc(obj%next_output, obj%get_file_freq(), obj%get_file_frequnit())
282 
283  if (obj%has_file_new_file_freq()) then
284  obj%next_close = diag_time_inc(obj%start_time, obj%get_file_new_file_freq(), &
285  obj%get_file_new_file_freq_units())
286  else
287  if (obj%has_file_duration()) then
288  obj%next_close = diag_time_inc(obj%start_time, obj%get_file_duration(), &
289  obj%get_file_duration_units())
290  else
291  obj%next_close = diag_time_inc(obj%start_time, very_large_file_freq, diag_days)
292  endif
293  endif
294  obj%is_file_open = .false.
295 
296  if(obj%has_file_duration()) then
297  obj%no_more_data = diag_time_inc(obj%start_time, obj%get_file_duration(), &
298  obj%get_file_duration_units())
299  else
300  obj%no_more_data = diag_time_inc(obj%start_time, very_large_file_freq, diag_days)
301  endif
302 
303  obj%unlim_dimension_level = 0
304  obj%num_time_levels = 0
305  obj%is_static = obj%get_file_freq() .eq. -1
306  obj%nz_subaxis = 0
307 
308  nullify(obj)
309  enddo set_ids_loop
310  fms_diag_files_object_init = .true.
311  else
312  fms_diag_files_object_init = .false.
313 ! mpp_error("fms_diag_files_object_init: The diag_table.yaml file has not been correctly parsed.",&
314 ! FATAL)
315  endif
316 end function fms_diag_files_object_init
317 
318 !< @brief Determine if the field corresponding to the field_id was registered to the file
319 !! @return .True. if the field was registed to the file
320 pure logical function is_field_registered(this, field_id)
321  class(fmsdiagfile_type), intent(in) :: this !< The file object
322  integer, intent(in) :: field_id !< Id of the field to check
323 
324  is_field_registered = this%field_registered(field_id)
325 end function is_field_registered
326 
327 !> \brief Adds a field and yaml ID to the file
328 subroutine add_field_and_yaml_id (this, new_field_id, yaml_id)
329  class(fmsdiagfile_type), intent(inout) :: this !< The file object
330  integer, intent(in) :: new_field_id !< The field ID to be added to field_ids
331  integer, intent(in) :: yaml_id !< The yaml_id
332 
333  this%num_registered_fields = this%num_registered_fields + 1
334  if (this%num_registered_fields .le. size(this%field_ids)) then
335  this%field_ids( this%num_registered_fields ) = new_field_id
336  this%yaml_ids( this%num_registered_fields ) = yaml_id
337  this%field_registered( this%num_registered_fields ) = .true.
338  else
339  call mpp_error(fatal, "The file: "//this%get_file_fname()//" has already been assigned its maximum "//&
340  "number of fields.")
341  endif
342 end subroutine add_field_and_yaml_id
343 
344 !> \brief Adds a buffer_id to the file object
345 subroutine add_buffer_id (this, buffer_id)
346  class(fmsdiagfile_type), intent(inout) :: this !< The file object
347  integer, intent(in) :: buffer_id !< Buffer id to add to the file
348 
349  this%number_of_buffers = this%number_of_buffers + 1
350  this%buffer_ids(this%number_of_buffers) = buffer_id
351 
352 end subroutine add_buffer_id
353 
354 !> \brief Initializes a diurnal axis for a fileobj
355 !! \note This is going to be called for every variable in the file, if the variable is not a diurnal variable
356 !! it will do nothing. It only defined a diurnal axis once.
357 subroutine init_diurnal_axis(this, diag_axis, naxis, yaml_id)
358  class(fmsdiagfile_type), intent(inout) :: this !< The file object
359  class(fmsdiagaxiscontainer_type), intent(inout) :: diag_axis(:) !< Array of diag_axis object
360  integer, intent(inout) :: naxis !< Number of diag_axis that heve been defined
361  integer, intent(in) :: yaml_id !< The ID to the variable's yaml
362 
363  integer :: i !< For do loops
364  type(diagyamlfilesvar_type), pointer :: field_yaml !< pointer to the yaml entry
365 
366  field_yaml => diag_yaml%get_diag_field_from_id(yaml_id)
367 
368  !< Go away if the file does not need a diurnal axis
369  if (.not. field_yaml%has_n_diurnal()) return
370 
371  !< Check if the diurnal axis is already defined for this number of diurnal samples
372  do i = 1, this%number_of_axis
373  select type(axis=>diag_axis(this%axis_ids(i))%axis)
374  type is (fmsdiagdiurnalaxis_type)
375  if(field_yaml%get_n_diurnal() .eq. axis%get_diurnal_axis_samples()) return
376  end select
377  end do
378 
379  !< If it is not already defined, define it
380  call define_diurnal_axis(diag_axis, naxis, field_yaml%get_n_diurnal(), .true.)
381  call define_diurnal_axis(diag_axis, naxis, field_yaml%get_n_diurnal(), .false.)
382 
383  !< Add it to the list of axis for the file
384  this%number_of_axis = this%number_of_axis + 1
385  this%axis_ids(this%number_of_axis) = naxis !< This is the diurnal axis edges
386 
387  this%number_of_axis = this%number_of_axis + 1
388  this%axis_ids(this%number_of_axis) = naxis - 1 !< This the diurnal axis
389 
390 end subroutine init_diurnal_axis
391 
392 !> \brief Set the time_ops variable in the diag_file object
393 subroutine set_file_time_ops(this, VarYaml, is_static)
394  class(fmsdiagfile_type), intent(inout) :: this !< The file object
395  type (diagyamlfilesvar_type), intent(in) :: varyaml !< The variable's yaml file
396  logical, intent(in) :: is_static !< Flag indicating if variable is static
397  integer, allocatable :: var_reduct !< temp to hold enumerated reduction type
398 
399  !< Go away if the file is static
400  if (this%is_static) return
401  if (is_static) return
402 
403  ! Set time_ops the first time this subroutine it is called
404  if (.not. allocated(this%time_ops)) then
405  var_reduct = varyaml%get_var_reduction()
406 
407  select case (var_reduct)
409  this%time_ops = .true.
410  case (time_none)
411  this%time_ops = .false.
412  end select
413 
414  return
415  endif
416 
417  if (this%time_ops) then
418  if (varyaml%get_var_reduction() .eq. time_none) &
419  call mpp_error(fatal, "The file: "//this%get_file_fname()//&
420  " has variables that are time averaged and instantaneous")
421  else
422  if (varyaml%get_var_reduction() .ne. time_none) &
423  call mpp_error(fatal, "The file: "//this%get_file_fname()//&
424  " has variables that are time averaged and instantaneous")
425  endif
426 end subroutine set_file_time_ops
427 
428 !> \brief Logical function to determine if the variable file_metadata_from_model has been allocated or associated
429 !! \return .True. if file_metadata_from_model exists .False. if file_metadata_from_model has not been set
430 pure logical function has_file_metadata_from_model (this)
431  class(fmsdiagfile_type), intent(in) :: this !< The file object
432  has_file_metadata_from_model = allocated(this%file_metadata_from_model)
433 end function has_file_metadata_from_model
434 
435 !> \brief Logical function to determine if the variable fileobj has been allocated or associated
436 !! \return .True. if fileobj exists .False. if fileobj has not been set
437 pure logical function has_fileobj (this)
438  class(fmsdiagfile_type), intent(in) :: this !< The file object
439  has_fileobj = allocated(this%fms2io_fileobj)
440 end function has_fileobj
441 
442 !> \brief Logical function to determine if the variable diag_yaml_file has been allocated or associated
443 !! \return .True. if diag_yaml_file exists .False. if diag_yaml has not been set
444 pure logical function has_diag_yaml_file (this)
445  class(fmsdiagfile_type), intent(in) :: this !< The file object
446  has_diag_yaml_file = associated(this%diag_yaml_file)
447 end function has_diag_yaml_file
448 
449 !> \brief Get the time to use to determine the filename, if using a wildcard file name (i.e ocn%4yr%2mo%2dy%2hr)
450 !! \return The time to use when determining the filename
451 function get_filename_time(this) &
452  result(res)
453  class(fmsdiagfile_type), intent(in) :: this !< The file object
454  type(time_type) :: res
455 
456  select case (this%diag_yaml_file%get_filename_time())
457  case (begin_time)
458  res = this%last_output
459  case (middle_time)
460  res = (this%last_output + this%next_close)/2
461  case (end_time)
462  res = this%next_close
463  end select
464 end function get_filename_time
465 
466 !> \brief Logical function to determine if the variable field_ids has been allocated or associated
467 !! \return .True. if field_ids exists .False. if field_ids has not been set
468 pure logical function has_field_ids (this)
469  class(fmsdiagfile_type), intent(in) :: this !< The file object
470  has_field_ids = allocated(this%field_ids)
471 end function has_field_ids
472 
473 !> \brief Returns a copy of the value of id
474 !! \return A copy of id
475 pure function get_id (this) result (res)
476  class(fmsdiagfile_type), intent(in) :: this !< The file object
477  integer :: res
478  res = this%id
479 end function get_id
480 
481 !> \brief Returns a copy of the value of time_ops
482 !! \return A copy of time_ops
483 pure function get_time_ops (this) result (res)
484  class(fmsdiagfile_type), intent(in) :: this !< The file object
485  logical :: res
486 
487  if (.not. allocated(this%time_ops)) then
488  res = .false.
489  else
490  res = this%time_ops
491  endif
492 end function get_time_ops
493 
494 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
495 !! TODO
496 !> \brief Returns a copy of the value of fileobj
497 !! \return A copy of fileobj
498 !pure function get_fileobj (obj) result (res)
499 ! class(fmsDiagFile_type), intent(in) :: obj !< The file object
500 ! class(FmsNetcdfFile_t) :: res
501 ! res = obj%fileobj
502 !end function get_fileobj
503 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
504 !! TODO
505 !!> \brief Returns a copy of the value of diag_yaml_file
506 !!! \return A copy of diag_yaml_file
507 !pure function get_diag_yaml_file (obj) result (res)
508 ! class(fmsDiagFile_type), intent(in) :: obj !< The file object
509 ! type(diagYamlFiles_type) :: res
510 ! res = obj%diag_yaml_file
511 !end function get_diag_yaml_file
512 
513 !> \brief Returns a copy of the value of file_metadata_from_model
514 !! \return A copy of file_metadata_from_model
515 pure function get_file_metadata_from_model (this) result (res)
516  class(fmsdiagfile_type), intent(in) :: this !< The file object
517  character(len=:), dimension(:), allocatable :: res
518  res = this%file_metadata_from_model
519 end function get_file_metadata_from_model
520 
521 !> \brief Returns a copy of the value of field_ids
522 !! \return A copy of field_ids
523 pure function get_field_ids (this) result (res)
524  class(fmsdiagfile_type), intent(in) :: this !< The file object
525  integer, dimension(:), allocatable :: res
526  allocate(res(size(this%field_ids)))
527  res = this%field_ids
528 end function get_field_ids
529 
530 !!!!!!!!! Functions from diag_yaml_file
531 !> \brief Returns a copy of file_fname from the yaml object
532 !! \return Copy of file_fname
533 pure function get_file_fname (this) result(res)
534  class(fmsdiagfile_type), intent(in) :: this !< The file object
535  character (len=:), allocatable :: res
536  res = this%diag_yaml_file%get_file_fname()
537 end function get_file_fname
538 
539 !> \brief Returns a copy of file_frequnit from the yaml object
540 !! \return Copy of file_frequnit
541 pure function get_file_frequnit (this) result(res)
542  class(fmsdiagfile_type), intent(in) :: this !< The file object
543  integer :: res
544  res = this%diag_yaml_file%get_file_frequnit()
545 end function get_file_frequnit
546 
547 !> \brief Returns a copy of file_freq from the yaml object
548 !! \return Copy of file_freq
549 pure function get_file_freq (this) result(res)
550  class(fmsdiagfile_type), intent(in) :: this !< The file object
551  integer :: res
552  res = this%diag_yaml_file%get_file_freq()
553 end function get_file_freq
554 
555 !> \brief Returns a copy of file_timeunit from the yaml object
556 !! \return Copy of file_timeunit
557 pure function get_file_timeunit (this) result(res)
558  class(fmsdiagfile_type), intent(in) :: this !< The file object
559  integer :: res
560  res = this%diag_yaml_file%get_file_timeunit()
561 end function get_file_timeunit
562 
563 !> \brief Returns a copy of file_unlimdim from the yaml object
564 !! \return Copy of file_unlimdim
565 pure function get_file_unlimdim (this) result(res)
566  class(fmsdiagfile_type), intent(in) :: this !< The file object
567  character (len=:), allocatable :: res
568  res = this%diag_yaml_file%get_file_unlimdim()
569 end function get_file_unlimdim
570 
571 !> \brief Returns a copy of file_sub_region from the yaml object
572 !! \return Pointer to file_sub_region
573 function get_file_sub_region (obj) result(res)
574  class(fmsdiagfile_type), target, intent(in) :: obj !< The file object
575  type(subregion_type), pointer :: res
576  res => obj%diag_yaml_file%get_file_sub_region()
577 end function get_file_sub_region
578 
579 !< @brief Query for the subregion grid type (latlon or index)
580 !! @return Pointer to subregion grid type
581 function get_file_sub_region_grid_type(this) &
582  result(res)
583  class(fmsdiagfile_type), intent(in) :: this !< Diag file object
584  integer :: res
585 
586  type(subregion_type), pointer :: subregion !< Subregion type
587 
588  if(this%diag_yaml_file%has_file_sub_region()) then
589  subregion => this%diag_yaml_file%get_file_sub_region()
590  res = subregion%grid_type
591  else
592  res = null_gridtype
593  endif
594 end function get_file_sub_region_grid_type
595 
596 !> \brief Returns a copy of file_new_file_freq from the yaml object
597 !! \return Copy of file_new_file_freq
598 pure function get_file_new_file_freq (this) result(res)
599  class(fmsdiagfile_type), intent(in) :: this !< The file object
600  integer :: res
601  res = this%diag_yaml_file%get_file_new_file_freq()
602 end function get_file_new_file_freq
603 
604 !> \brief Returns a copy of file_new_file_freq_units from the yaml object
605 !! \return Copy of file_new_file_freq_units
606 pure function get_file_new_file_freq_units (this) result(res)
607  class(fmsdiagfile_type), intent(in) :: this !< The file object
608  integer :: res
609  res = this%diag_yaml_file%get_file_new_file_freq_units()
610 end function get_file_new_file_freq_units
611 
612 !> \brief Returns a copy of file_start_time from the yaml object
613 !! \return Copy of file_start_time
614 pure function get_file_start_time (this) result(res)
615  class(fmsdiagfile_type), intent(in) :: this !< The file object
616  type(time_type):: res
617  res = this%diag_yaml_file%get_file_start_time()
618 end function get_file_start_time
619 
620 !> \brief Returns a copy of file_duration from the yaml object
621 !! \return Copy of file_duration
622 pure function get_file_duration (this) result(res)
623  class(fmsdiagfile_type), intent(in) :: this !< The file object
624  integer :: res
625  res = this%diag_yaml_file%get_file_duration()
626 end function get_file_duration
627 
628 !> \brief Returns a copy of file_duration_units from the yaml object
629 !! \return Copy of file_duration_units
630 pure function get_file_duration_units (this) result(res)
631  class(fmsdiagfile_type), intent(in) :: this !< The file object
632  integer :: res
633  res = this%diag_yaml_file%get_file_duration_units()
634 end function get_file_duration_units
635 
636 !> \brief Returns a copy of file_varlist from the yaml object
637 !! \return Copy of file_varlist
638 pure function get_file_varlist (this) result(res)
639  class(fmsdiagfile_type), intent(in) :: this !< The file object
640  character (len=:), allocatable, dimension(:) :: res
641  res = this%diag_yaml_file%get_file_varlist()
642 end function get_file_varlist
643 
644 !> \brief Returns a copy of file_global_meta from the yaml object
645 !! \return Copy of file_global_meta
646 pure function get_file_global_meta (this) result(res)
647  class(fmsdiagfile_type), intent(in) :: this !< The file object
648  character (len=MAX_STR_LEN), allocatable, dimension(:,:) :: res
649  res = this%diag_yaml_file%get_file_global_meta()
650 end function get_file_global_meta
651 
652 !> \brief Determines if done writing data
653 !! \return .True. if done writing data
654 pure function is_done_writing_data (this) result(res)
655  class(fmsdiagfile_type), intent(in) :: this !< The file object
656  logical :: res
657  res = this%done_writing_data
658  if (this%is_file_open) res = .false.
659 end function is_done_writing_data
660 
661 !> \brief Checks if file_fname is allocated in the yaml object
662 !! \return true if file_fname is allocated
663 pure function has_file_fname (this) result(res)
664  class(fmsdiagfile_type), intent(in) :: this !< The file object
665  logical :: res
666  res = this%diag_yaml_file%has_file_fname()
667 end function has_file_fname
668 
669 !> \brief Checks if file_frequnit is allocated in the yaml object
670 !! \return true if file_frequnit is allocated
671 pure function has_file_frequnit (this) result(res)
672  class(fmsdiagfile_type), intent(in) :: this !< The file object
673  logical :: res
674  res = this%diag_yaml_file%has_file_frequnit()
675 end function has_file_frequnit
676 
677 !> \brief Checks if file_freq is allocated in the yaml object
678 !! \return true if file_freq is allocated
679 pure function has_file_freq (this) result(res)
680  class(fmsdiagfile_type), intent(in) :: this !< The file object
681  logical :: res
682  res = this%diag_yaml_file%has_file_freq()
683 end function has_file_freq
684 
685 !> \brief Checks if file_timeunit is allocated in the yaml object
686 !! \return true if file_timeunit is allocated
687 pure function has_file_timeunit (this) result(res)
688  class(fmsdiagfile_type), intent(in) :: this !< The file object
689  logical :: res
690  res = this%diag_yaml_file%has_file_timeunit()
691 end function has_file_timeunit
692 
693 !> \brief Checks if file_unlimdim is allocated in the yaml object
694 !! \return true if file_unlimdim is allocated
695 pure function has_file_unlimdim (this) result(res)
696  class(fmsdiagfile_type), intent(in) :: this !< The file object
697  logical :: res
698  res = this%diag_yaml_file%has_file_unlimdim()
699 end function has_file_unlimdim
700 
701 !> \brief Checks if file_sub_region is allocated in the yaml object
702 !! \return true if file_sub_region is allocated
703 pure function has_file_sub_region (this) result(res)
704  class(fmsdiagfile_type), intent(in) :: this !< The file object
705  logical :: res
706  res = this%diag_yaml_file%has_file_sub_region()
707 end function has_file_sub_region
708 
709 !> \brief Checks if file_new_file_freq is allocated in the yaml object
710 !! \return true if file_new_file_freq is allocated
711 pure function has_file_new_file_freq (this) result(res)
712  class(fmsdiagfile_type), intent(in) :: this !< The file object
713  logical :: res
714  res = this%diag_yaml_file%has_file_new_file_freq()
715 end function has_file_new_file_freq
716 
717 !> \brief Checks if file_new_file_freq_units is allocated in the yaml object
718 !! \return true if file_new_file_freq_units is allocated
719 pure function has_file_new_file_freq_units (this) result(res)
720  class(fmsdiagfile_type), intent(in) :: this !< The file object
721  logical :: res
722  res = this%diag_yaml_file%has_file_new_file_freq_units()
723 end function has_file_new_file_freq_units
724 
725 !> \brief Checks if file_start_time is allocated in the yaml object
726 !! \return true if file_start_time is allocated
727 pure function has_file_start_time (this) result(res)
728  class(fmsdiagfile_type), intent(in) :: this !< The file object
729  logical :: res
730  res = this%diag_yaml_file%has_file_start_time()
731 end function has_file_start_time
732 
733 !> \brief Checks if file_duration is allocated in the yaml object
734 !! \return true if file_duration is allocated
735 pure function has_file_duration (this) result(res)
736  class(fmsdiagfile_type), intent(in) :: this !< The file object
737  logical :: res
738  res = this%diag_yaml_file%has_file_duration()
739 end function has_file_duration
740 
741 !> \brief Checks if file_duration_units is allocated in the yaml object
742 !! \return true if file_duration_units is allocated
743 pure function has_file_duration_units (this) result(res)
744  class(fmsdiagfile_type), intent(in) :: this !< The file object
745  logical :: res
746  res = this%diag_yaml_file%has_file_duration_units()
747 end function has_file_duration_units
748 
749 !> \brief Checks if file_varlist is allocated in the yaml object
750 !! \return true if file_varlist is allocated
751 pure function has_file_varlist (this) result(res)
752  class(fmsdiagfile_type), intent(in) :: this !< The file object
753  logical :: res
754  res = this%diag_yaml_file%has_file_varlist()
755 end function has_file_varlist
756 
757 !> \brief Checks if file_global_meta is allocated in the yaml object
758 !! \return true if file_global_meta is allocated
759 pure function has_file_global_meta (this) result(res)
760  class(fmsdiagfile_type), intent(in) :: this !< The file object
761  logical :: res
762  res = this%diag_yaml_file%has_file_global_meta()
763 end function has_file_global_meta
764 
765 !> @brief Sets the domain and type of domain from the axis IDs
766 subroutine set_domain_from_axis(this, diag_axis, axes)
767  class(fmsdiagfile_type), intent(inout) :: this !< The file object
768  class(fmsdiagaxiscontainer_type), intent(in) :: diag_axis(:) !< Array of diag_axis
769  integer, intent(in) :: axes (:)
770 
771  call get_domain_and_domain_type(diag_axis, axes, this%type_of_domain, this%domain, this%get_file_fname())
772 end subroutine set_domain_from_axis
773 
774 !> @brief Set the domain and the type_of_domain for a file
775 !> @details This subroutine is going to be called once by every variable in the file
776 !! in register_diag_field. It will update the domain and the type_of_domain if needed and verify that
777 !! all the variables are in the same domain
778 subroutine set_file_domain(this, domain, type_of_domain)
779  class(fmsdiagfile_type), intent(inout) :: this !< The file object
780  integer, INTENT(in) :: type_of_domain !< fileobj_type to use
781  CLASS(diagdomain_t), INTENT(in), pointer :: domain !< Domain
782 
783  if (type_of_domain .ne. this%type_of_domain) then
784  !! If the current type_of_domain in the file obj is not the same as the variable calling this subroutine
785 
786  if (type_of_domain .eq. no_domain .or. this%type_of_domain .eq. no_domain) then
787  !! If they are not the same then one of them can be NO_DOMAIN
788  !! (i.e a file can have variables that are not domain decomposed and variables that are)
789 
790  if (type_of_domain .ne. no_domain) then
791  !! Update the file's type_of_domain and domain if needed
792  this%type_of_domain = type_of_domain
793  this%domain => domain
794  endif
795 
796  else
797  !! If they are not the same and of them is not NO_DOMAIN, then crash because the variables don't have the
798  !! same domain (i.e a file has a variable is that in a 2D domain and one that is in a UG domain)
799 
800  call mpp_error(fatal, "The file: "//this%get_file_fname()//" has variables that are not in the same domain")
801  endif
802  endif
803 
804 end subroutine set_file_domain
805 
806 !> @brief Loops through a variable's axis_ids and adds them to the FMSDiagFile object if they don't exist
807 subroutine add_axes(this, axis_ids, diag_axis, naxis, yaml_id, buffer_id, output_buffers)
808  class(fmsdiagfile_type), intent(inout) :: this !< The file object
809  integer, INTENT(in) :: axis_ids(:) !< Array of axes_ids
810  class(fmsdiagaxiscontainer_type), intent(inout) :: diag_axis(:) !< Diag_axis object
811  integer, intent(inout) :: naxis !< Number of axis that have been
812  !! registered
813  integer, intent(in) :: yaml_id !< Yaml id of the field section for
814  !! this var
815  integer, intent(in) :: buffer_id !< ID of the buffer
816  type(fmsdiagoutputbuffer_type), intent(inout) :: output_buffers(:) !< Array of output buffers
817 
818  type(diagyamlfilesvar_type), pointer :: field_yaml !< pointer to the yaml entry
819 
820  integer :: i, j !< For do loops
821  logical :: is_cube_sphere !< Flag indicating if the file's domain is a cubesphere
822  logical :: axis_found !< Flag indicating that the axis was already to the file obj
823  integer, allocatable :: var_axis_ids(:) !< Array of the variable's axis ids
824  integer :: x_y_axis_id(2) !< Ids of the x and y axis
825  integer :: x_or_y !< integer indicating if the axis is x or y
826  logical :: is_x_or_y !< flag indicating if the axis is x or y
827  integer :: subregion_gridtype !< The type of the subregion (latlon or index)
828  logical :: write_on_this_pe !< Flag indicating if the current pe is in the subregion
829 
830  is_cube_sphere = .false.
831  subregion_gridtype = this%get_file_sub_region_grid_type()
832 
833  field_yaml => diag_yaml%get_diag_field_from_id(yaml_id)
834 
835  !< Created a copy here, because if the variable has a z subaxis var_axis_ids will be modified in
836  !! `create_new_z_subaxis` to contain the id of the new z subaxis instead of the parent axis,
837  !! which will be added to the the list of axis in the file object (axis_ids is intent(in),
838  !! which is why the copy was needed)
839  var_axis_ids = axis_ids
840 
841  if (field_yaml%has_var_zbounds()) then
842  call create_new_z_subaxis(field_yaml%get_var_zbounds(), var_axis_ids, diag_axis, naxis, &
843  this%axis_ids, this%number_of_axis, this%nz_subaxis)
844  endif
845 
846  select type(this)
847  type is (subregionalfile_type)
848  if (associated(this%domain)) then
849  if (this%domain%get_ntiles() .eq. 6) is_cube_sphere = .true.
850  endif
851  if (.not. this%get_write_on_this_pe()) return
852  subaxis_defined: if (this%is_subaxis_defined) then
853  do i = 1, size(var_axis_ids)
854  select type (parent_axis => diag_axis(var_axis_ids(i))%axis)
855  type is (fmsdiagfullaxis_type)
856  axis_found = .false.
857  is_x_or_y = parent_axis%is_x_or_y_axis()
858  do j = 1, this%number_of_axis
859  if (is_x_or_y) then
860  if(is_parent_axis(this%axis_ids(j), var_axis_ids(i), diag_axis)) then
861  axis_found = .true.
862  var_axis_ids(i) = this%axis_ids(j) !Set the var_axis_id to the sub axis_id
863  cycle
864  endif
865  elseif (var_axis_ids(i) .eq. this%axis_ids(j)) then
866  axis_found = .true.
867  endif
868  enddo
869 
870  if (.not. axis_found) then
871  if (is_x_or_y) then
872  if (subregion_gridtype .eq. latlon_gridtype .and. is_cube_sphere) &
873  call mpp_error(fatal, "If using the cube sphere and defining the subregion with latlon "//&
874  "the variable need to have the same x and y axis. Please check the variables in the file "//&
875  trim(this%get_file_fname())//" or use indices to define the subregion.")
876 
877  select case (subregion_gridtype)
878  case (index_gridtype)
879  call define_new_subaxis_index(parent_axis, this%get_file_sub_region(), diag_axis, naxis, &
880  i, write_on_this_pe)
881  case (latlon_gridtype)
882  call define_new_subaxis_latlon(diag_axis, var_axis_ids(i:i), naxis, this%get_file_sub_region(), &
883  .false., write_on_this_pe)
884  end select
885  call this%update_write_on_this_pe(write_on_this_pe)
886  if (.not. this%get_write_on_this_pe()) cycle
887  call this%add_new_axis(naxis)
888  var_axis_ids(i) = naxis
889  else
890  call this%add_new_axis(var_axis_ids(i))
891  endif
892  endif
893  type is (fmsdiagsubaxis_type)
894  axis_found = this%does_axis_exist(var_axis_ids(i))
895  if (.not. axis_found) call this%add_new_axis(var_axis_ids(i))
896  end select
897  enddo
898  else
899  x_y_axis_id = diag_null
900  do i = 1, size(var_axis_ids)
901  select type (parent_axis => diag_axis(var_axis_ids(i))%axis)
902  type is (fmsdiagfullaxis_type)
903  if (.not. parent_axis%is_x_or_y_axis(x_or_y)) then
904  axis_found = this%does_axis_exist(var_axis_ids(i))
905  if (.not. axis_found) call this%add_new_axis(var_axis_ids(i))
906  else
907  x_y_axis_id(x_or_y) = var_axis_ids(i)
908  endif
909  type is (fmsdiagsubaxis_type)
910  axis_found = this%does_axis_exist(var_axis_ids(i))
911  if (.not. axis_found) call this%add_new_axis(var_axis_ids(i))
912  end select
913  enddo
914 
915  call this%define_new_subaxis(var_axis_ids, x_y_axis_id, is_cube_sphere, diag_axis, naxis)
916  this%is_subaxis_defined = .true.
917  endif subaxis_defined
918  type is (fmsdiagfile_type)
919  do i = 1, size(var_axis_ids)
920  axis_found = this%does_axis_exist(var_axis_ids(i))
921  if (.not. axis_found) call this%add_new_axis(var_axis_ids(i))
922  enddo
923  end select
924  !> Add the axis to the buffer object
925  call output_buffers(buffer_id)%add_axis_ids(var_axis_ids)
926 end subroutine add_axes
927 
928 !> @brief Adds a new axis the list of axis in the diag file object
929 subroutine add_new_axis(this, var_axis_id)
930  class(fmsdiagfile_type), intent(inout) :: this !< The file object
931  integer, intent(in) :: var_axis_id !< Axis id of the variable
932 
933  this%number_of_axis = this%number_of_axis + 1
934  this%axis_ids(this%number_of_axis) = var_axis_id
935 end subroutine add_new_axis
936 
937 !> @brief This updates write on this pe
938 subroutine update_write_on_this_pe(this, write_on_this_pe)
939  class(fmsdiagfile_type), intent(inout) :: this !< The file object
940  logical, intent(in) :: write_on_this_pe !< .True. if the current PE is in
941  !! subregion
942 
943  select type (this)
944  type is (subregionalfile_type)
945  if (this%write_on_this_pe) this%write_on_this_pe = write_on_this_pe
946  end select
947 end subroutine update_write_on_this_pe
948 
949 !> @brief Query for the write_on_this_pe member of the diag file object
950 !! @return the write_on_this_pe member of the diag file object
951 function get_write_on_this_pe(this) &
952  result(rslt)
953  class(fmsdiagfile_type), intent(inout) :: this !< The file object
954  logical :: rslt
955  rslt = .true.
956  select type (this)
957  type is (subregionalfile_type)
958  rslt= this%write_on_this_pe
959  end select
960 end function get_write_on_this_pe
961 
962 !< @brief Determine if an axis is already in the list of axis for a diag file
963 !! @return .True. if the axis is already in the list of axis for a diag file
964 function does_axis_exist(this, var_axis_id) &
965  result(rslt)
966  class(fmsdiagfile_type), intent(inout) :: this !< The file object
967  integer, intent(in) :: var_axis_id !< Variable axis id to check
968 
969  logical :: rslt
970  integer :: j !< For do loops
971 
972  rslt = .false.
973  do j = 1, this%number_of_axis
974  !> Check if the axis already exists, move on
975  if (var_axis_id .eq. this%axis_ids(j)) then
976  rslt = .true.
977  return
978  endif
979  enddo
980 end function
981 
982 !> @brief Define a new sub axis
983 subroutine define_new_subaxis(this, var_axis_ids, x_y_axis_id, is_cube_sphere, diag_axis, naxis)
984  class(fmsdiagfile_type), intent(inout) :: this !< The file object
985  integer, INTENT(inout) :: var_axis_ids(:) !< Original variable axis ids
986  integer, INTENT(in) :: x_y_axis_id(:) !< The ids of the x and y axis
987  logical, intent(in) :: is_cube_sphere !< .True. if the axis is in the cubesphere
988  integer, intent(inout) :: naxis !< Number of axis current registered
989  class(fmsdiagaxiscontainer_type), intent(inout) :: diag_axis(:) !< Diag_axis object
990 
991  logical :: write_on_this_pe !< .True. if the current PE is in the subregion
992  integer :: i, j !< For do loop
993 
994  select case (this%get_file_sub_region_grid_type())
995  case(latlon_gridtype)
996  call define_new_subaxis_latlon(diag_axis, x_y_axis_id, naxis, this%get_file_sub_region(), is_cube_sphere, &
997  write_on_this_pe)
998  call this%update_write_on_this_pe(write_on_this_pe)
999  if (.not. this%get_write_on_this_pe()) return
1000  call this%add_new_axis(naxis)
1001  call this%add_new_axis(naxis-1)
1002  do j = 1, size(var_axis_ids)
1003  if (x_y_axis_id(1) .eq. var_axis_ids(j)) var_axis_ids(j) = naxis - 1
1004  if (x_y_axis_id(2) .eq. var_axis_ids(j)) var_axis_ids(j) = naxis
1005  enddo
1006  case (index_gridtype)
1007  do i = 1, size(x_y_axis_id)
1008  select type (parent_axis => diag_axis(x_y_axis_id(i))%axis)
1009  type is (fmsdiagfullaxis_type)
1010  call define_new_subaxis_index(parent_axis, this%get_file_sub_region(), diag_axis, naxis, i, &
1011  write_on_this_pe)
1012  call this%update_write_on_this_pe(write_on_this_pe)
1013  if (.not. this%get_write_on_this_pe()) return
1014  call this%add_new_axis(naxis)
1015  do j = 1, size(var_axis_ids)
1016  if (x_y_axis_id(i) .eq. var_axis_ids(j)) var_axis_ids(j) = naxis
1017  enddo
1018  end select
1019  enddo
1020  end select
1021 end subroutine define_new_subaxis
1022 
1023 !> @brief adds the start time to the fileobj
1024 !! @note This should be called from the register field calls. It can be called multiple times (one for each variable)
1025 !! So it needs to make sure that the start_time is the same for each variable. The initial value is the base_time
1026 subroutine add_start_time(this, start_time)
1027  class(fmsdiagfile_type), intent(inout) :: this !< The file object
1028  TYPE(time_type), intent(in) :: start_time !< Start time passed into register_diag_field
1029 
1030  !< If the start_time sent in is equal to the diag_init_time return because
1031  !! this%start_time was already set to the diag_init_time
1032  if (start_time .eq. diag_init_time) return
1033 
1034  !< If the start_time sent is is greater than or equal to the start time already
1035  !! in the diag file obj return because either this%start_time was already updated
1036  !! or the file has start_time defined in the yaml
1037  if (this%start_time >= start_time) return
1038 
1039  if (this%start_time .ne. diag_init_time) then
1040  !> If the this%start_time is not equal to the diag_init_time from the diag_table
1041  !! this%start_time was already updated so make sure it is the same for the current variable
1042  !! or error out
1043  if (this%start_time .ne. start_time)&
1044  call mpp_error(fatal, "The variables associated with the file:"//this%get_file_fname()//" have &
1045  &different start_time")
1046  else
1047  !> If the this%start_time is equal to the diag_init_time,
1048  !! simply update it with the start_time and set up the *_output variables
1049  this%model_time = start_time
1050  this%start_time = start_time
1051  this%last_output = start_time
1052  this%next_output = diag_time_inc(start_time, this%get_file_freq(), this%get_file_frequnit())
1053  this%next_next_output = diag_time_inc(this%next_output, this%get_file_freq(), this%get_file_frequnit())
1054  if (this%has_file_new_file_freq()) then
1055  this%next_close = diag_time_inc(this%start_time, this%get_file_new_file_freq(), &
1056  this%get_file_new_file_freq_units())
1057  else
1058  if (this%is_static) then
1059  ! If the file is static, set the close time to be equal to the start_time, so that it can be closed
1060  ! after the first write!
1061  this%next_close = this%start_time
1062  this%next_next_output = diag_time_inc(this%start_time, very_large_file_freq, diag_days)
1063  else
1064  this%next_close = diag_time_inc(this%start_time, very_large_file_freq, diag_days)
1065  endif
1066  endif
1067 
1068  if(this%has_file_duration()) then
1069  this%no_more_data = diag_time_inc(this%start_time, this%get_file_duration(), &
1070  this%get_file_duration_units())
1071  else
1072  this%no_more_data = diag_time_inc(this%start_time, very_large_file_freq, diag_days)
1073  endif
1074 
1075  endif
1076 
1077 end subroutine
1078 
1079 !> writes out internal values for fmsDiagFile_type object
1080 subroutine dump_file_obj(this, unit_num)
1081  class(fmsdiagfile_type), intent(in) :: this !< the file object
1082  integer, intent(in) :: unit_num !< passed in from dump_diag_obj
1083  !! will either be for new log file or stdout
1084  write( unit_num, *) 'file id:', this%id
1085  write( unit_num, *) 'start time:', date_to_string(this%start_time)
1086  write( unit_num, *) 'last_output', date_to_string(this%last_output)
1087  write( unit_num, *) 'next_output', date_to_string(this%next_output)
1088  write( unit_num, *)'next_next_output', date_to_string(this%next_next_output)
1089  write( unit_num, *)'next_close', date_to_string(this%next_close)
1090 
1091  if( allocated(this%fms2io_fileobj)) write( unit_num, *)'fileobj path', this%fms2io_fileobj%path
1092 
1093  write( unit_num, *)'type_of_domain', this%type_of_domain
1094  if( allocated(this%file_metadata_from_model)) write( unit_num, *) 'file_metadata_from_model', &
1095  this%file_metadata_from_model
1096  if( allocated(this%field_ids)) write( unit_num, *)'field_ids', this%field_ids
1097  if( allocated(this%field_registered)) write( unit_num, *)'field_registered', this%field_registered
1098  if( allocated(this%num_registered_fields)) write( unit_num, *)'num_registered_fields', this%num_registered_fields
1099  if( allocated(this%axis_ids)) write( unit_num, *)'axis_ids', this%axis_ids(1:this%number_of_axis)
1100 
1101 end subroutine
1102 
1103 !> @brief Determine if a file is regional
1104 !! @return Flag indicating if the file is regional or not
1105 logical pure function is_regional(this)
1106  class(fmsdiagfilecontainer_type), intent(in) :: this !< The file object
1107 
1108  select type (wut=>this%FMS_diag_file)
1109  type is (subregionalfile_type)
1110  is_regional = .true.
1111  type is (fmsdiagfile_type)
1112  is_regional = .false.
1113  end select
1114 
1115 end function is_regional
1116 
1117 !> @brief Determine if a file is static
1118 !! @return Flag indicating if the file is static or not
1119 logical pure function is_file_static(this)
1120 class(fmsdiagfilecontainer_type), intent(in) :: this !< The file object
1121 
1122 is_file_static = .false.
1123 
1124 select type (fileptr=>this%FMS_diag_file)
1125 type is (fmsdiagfile_type)
1126  is_file_static = fileptr%is_static
1127 end select
1128 
1129 end function is_file_static
1130 
1131 !< @brief Opens the diag_file if it is time to do so
1132 subroutine open_diag_file(this, time_step, file_is_opened)
1133  class(fmsdiagfilecontainer_type), intent(inout), target :: this !< The file object
1134  TYPE(time_type), intent(in) :: time_step !< Current model step time
1135  logical, intent(out) :: file_is_opened !< .true. if the file was opened in this
1136  !! time
1137 
1138  class(fmsdiagfile_type), pointer :: diag_file !< Diag_file object to open
1139  class(diagdomain_t), pointer :: domain !< The domain used in the file
1140  character(len=:), allocatable :: diag_file_name !< The file name as defined in the yaml
1141  character(len=FMS_FILE_LEN) :: base_name !< The file name as defined in the yaml
1142  !! without the wildcard definition
1143  character(len=FMS_FILE_LEN) :: file_name !< The file name as it will be written to disk
1144  character(len=FMS_FILE_LEN) :: temp_name !< Temp variable to store the file_name
1145  character(len=128) :: start_date !< The start_time as a string that will be added to
1146  !! the begining of the filename (start_date.filename)
1147  character(len=128) :: suffix !< The current time as a string that will be added to
1148  !! the end of filename
1149  integer :: pos !< Index of the filename with the first "%" in the file name
1150  INTEGER :: year !< The year of the start_date
1151  INTEGER :: month !< The month of the start_date
1152  INTEGER :: day !< The day of the start_date
1153  INTEGER :: hour !< The hour of the start_date
1154  INTEGER :: minute !< The minute of the start_date
1155  INTEGER :: second !< The second of the start_date
1156  character(len=4) :: mype_string !< The pe as a string
1157  logical :: is_regional !< Flag indicating if the file is regional
1158  integer, allocatable :: pes(:) !< Array of the pes in the current pelist
1159 
1160  diag_file => this%FMS_diag_file
1161  domain => diag_file%domain
1162 
1163  file_is_opened = .false.
1164  !< Go away if it the file is already open
1165  if (diag_file%is_file_open) return
1166 
1167  is_regional = .false.
1168  !< Figure out what fms2io_fileobj to use!
1169  if (.not. allocated(diag_file%fms2io_fileobj)) then
1170  select type (diag_file)
1171  type is (subregionalfile_type)
1172  !< In this case each PE is going to write its own file
1173  allocate(fmsnetcdffile_t :: diag_file%fms2io_fileobj)
1174  is_regional = .true.
1175  type is (fmsdiagfile_type)
1176  !< Use the type_of_domain to get the correct fms2io_fileobj
1177  select case (diag_file%type_of_domain)
1178  case (no_domain)
1179  allocate(fmsnetcdffile_t :: diag_file%fms2io_fileobj)
1180  case (two_d_domain)
1181  allocate(fmsnetcdfdomainfile_t :: diag_file%fms2io_fileobj)
1182  case (ug_domain)
1183  allocate(fmsnetcdfunstructureddomainfile_t :: diag_file%fms2io_fileobj)
1184  end select
1185  end select
1186  endif
1187 
1188  !< Figure out what to name of the file
1189  diag_file_name = diag_file%get_file_fname()
1190 
1191  !< If using the new_file_freq figure out what the name is based on the current time
1192  if (diag_file%has_file_new_file_freq()) then
1193  !< If using a wildcard file name (i.e ocn%4yr%2mo%2dy%2hr), get the basename (i.e ocn)
1194  pos = index(diag_file_name, '%')
1195  if (pos > 0) base_name = diag_file_name(1:pos-1)
1196  suffix = get_time_string(diag_file_name, diag_file%get_filename_time())
1197  base_name = trim(base_name)//trim(suffix)
1198  else
1199  base_name = trim(diag_file_name)
1200  endif
1201 
1202  !< Add the ens number to the file name (if it exists)
1203  file_name = trim(base_name)
1204  call get_instance_filename(base_name, file_name)
1205 
1206  !< Prepend the file start_time to the file name if prepend_date == .TRUE. in
1207  !! the namelist
1208  IF ( prepend_date ) THEN
1209  call get_date(diag_file%start_time, year, month, day, hour, minute, second)
1210  write (start_date, '(1I20.4, 2I2.2)') year, month, day
1211 
1212  file_name = trim(adjustl(start_date))//'.'//trim(file_name)
1213  END IF
1214 
1215  file_name = trim(file_name)//".nc"
1216 
1217  !< If this is a regional file add the PE and the tile_number to the filename
1218  if (is_regional) then
1219  !< Get the pe number that will be appended to the end of the file
1220  write(mype_string,'(I0.4)') mpp_pe()
1221 
1222  !< Add the tile number if appropriate
1223  select type (domain)
1224  type is (diagdomain2d_t)
1225  temp_name = file_name
1226  call get_mosaic_tile_file(temp_name, file_name, .true., domain%domain2)
1227  end select
1228 
1229  file_name = trim(file_name)//"."//trim(mype_string)
1230  endif
1231 
1232  !< Open the file!
1233  select type (fms2io_fileobj => diag_file%fms2io_fileobj)
1234  type is (fmsnetcdffile_t)
1235  if (is_regional) then
1236  if (.not. open_file(fms2io_fileobj, file_name, "overwrite", pelist=(/mpp_pe()/))) &
1237  &call mpp_error(fatal, "Error opening the file:"//file_name)
1238  call register_global_attribute(fms2io_fileobj, "is_subregional", "True", str_len=4)
1239  else
1240  allocate(pes(mpp_npes()))
1241  call mpp_get_current_pelist(pes)
1242 
1243  if (.not. open_file(fms2io_fileobj, file_name, "overwrite", pelist=pes)) &
1244  &call mpp_error(fatal, "Error opening the file:"//file_name)
1245  endif
1246  type is (fmsnetcdfdomainfile_t)
1247  select type (domain)
1248  type is (diagdomain2d_t)
1249  if (.not. open_file(fms2io_fileobj, file_name, "overwrite", domain%Domain2)) &
1250  &call mpp_error(fatal, "Error opening the file:"//file_name)
1251  end select
1252  type is (fmsnetcdfunstructureddomainfile_t)
1253  select type (domain)
1254  type is (diagdomainug_t)
1255  if (.not. open_file(fms2io_fileobj, file_name, "overwrite", domain%DomainUG)) &
1256  &call mpp_error(fatal, "Error opening the file:"//file_name)
1257  end select
1258  end select
1259 
1260  file_is_opened = .true.
1261  diag_file%is_file_open = file_is_opened
1262  domain => null()
1263  diag_file => null()
1264 end subroutine open_diag_file
1265 
1266 !< @brief Write global attributes in the diag_file
1267 subroutine write_global_metadata(this)
1268  class(fmsdiagfilecontainer_type), intent(inout), target :: this !< The file object
1269 
1270  class(fmsnetcdffile_t), pointer :: fms2io_fileobj !< The fileobj to write to
1271  integer :: i !< For do loops
1272  character (len=MAX_STR_LEN), allocatable :: yaml_file_attributes(:,:) !< Global attributes defined in the yaml
1273 
1274  type(diagyamlfiles_type), pointer :: diag_file_yaml !< The diag_file yaml
1275 
1276  diag_file_yaml => this%FMS_diag_file%diag_yaml_file
1277  fms2io_fileobj => this%FMS_diag_file%fms2io_fileobj
1278 
1279  if (diag_file_yaml%has_file_global_meta()) then
1280  yaml_file_attributes = diag_file_yaml%get_file_global_meta()
1281  do i = 1, size(yaml_file_attributes,1)
1282  call register_global_attribute(fms2io_fileobj, trim(yaml_file_attributes(i,1)), &
1283  trim(yaml_file_attributes(i,2)), str_len=len_trim(yaml_file_attributes(i,2)))
1284  enddo
1285  deallocate(yaml_file_attributes)
1286  endif
1287 end subroutine write_global_metadata
1288 
1289 !< @brief Writes a variable's metadata in the netcdf file
1290 subroutine write_var_metadata(fms2io_fileobj, variable_name, dimensions, long_name, units)
1291  class(fmsnetcdffile_t), intent(inout) :: fms2io_fileobj !< The file object to write into
1292  character(len=*) , intent(in) :: variable_name !< The name of the time variables
1293  character(len=*) , intent(in) :: dimensions(:) !< The dimensions of the variable
1294  character(len=*) , intent(in) :: long_name !< The long_name of the variable
1295  character(len=*) , intent(in) :: units !< The units of the variable
1296 
1297  call register_field(fms2io_fileobj, variable_name, pack_size_str, dimensions)
1298  call register_variable_attribute(fms2io_fileobj, variable_name, "long_name", &
1299  trim(long_name), str_len=len_trim(long_name))
1300  if (trim(units) .ne. no_units) &
1301  call register_variable_attribute(fms2io_fileobj, variable_name, "units", &
1302  trim(units), str_len=len_trim(units))
1303 end subroutine write_var_metadata
1304 
1305 !> \brief Write the time metadata to the diag file
1306 subroutine write_time_metadata(this)
1307  class(fmsdiagfilecontainer_type), intent(inout), target :: this !< The file object
1308 
1309  class(fmsdiagfile_type), pointer :: diag_file !< Diag_file object to open
1310  class(fmsnetcdffile_t), pointer :: fms2io_fileobj !< The fileobj to write to
1311  character(len=50) :: time_units_str !< Time units written as a string
1312  character(len=50) :: calendar !< The calendar name
1313 
1314  character(len=:), allocatable :: time_var_name !< The name of the time variable as it is defined in the yaml
1315  character(len=50) :: dimensions(2) !< Array of dimensions names for the variable
1316 
1317  diag_file => this%FMS_diag_file
1318  fms2io_fileobj => diag_file%fms2io_fileobj
1319 
1320  time_var_name = diag_file%get_file_unlimdim()
1321  call register_axis(fms2io_fileobj, time_var_name, unlimited)
1322 
1323  WRITE(time_units_str, 11) &
1324  trim(time_unit_list(diag_file%get_file_timeunit())), get_base_year(),&
1325  & get_base_month(), get_base_day(), get_base_hour(), get_base_minute(), get_base_second()
1326 11 FORMAT(a, ' since ', i4.4, '-', i2.2, '-', i2.2, ' ', i2.2, ':', i2.2, ':', i2.2)
1327 
1328  dimensions(1) = "nv"
1329  dimensions(2) = trim(time_var_name)
1330 
1331  call write_var_metadata(fms2io_fileobj, time_var_name, dimensions(2:2), &
1332  time_var_name, time_units_str)
1333 
1334  !< Add additional variables to the time variable
1335  call register_variable_attribute(fms2io_fileobj, time_var_name, "axis", "T", str_len=1 )
1336 
1337  !TODO no need to have both attributes, probably?
1338  calendar = valid_calendar_types(get_calendar_type())
1339  call register_variable_attribute(fms2io_fileobj, time_var_name, "calendar_type", &
1340  uppercase(trim(calendar)), str_len=len_trim(calendar))
1341  call register_variable_attribute(fms2io_fileobj, time_var_name, "calendar", &
1342  lowercase(trim(calendar)), str_len=len_trim(calendar))
1343 
1344  if (diag_file%get_time_ops()) then
1345  call register_variable_attribute(fms2io_fileobj, time_var_name, "bounds", &
1346  trim(time_var_name)//"_bnds", str_len=len_trim(time_var_name//"_bnds"))
1347 
1348  !< It is possible that the "nv" "axis" was registered via "diag_axis_init" call
1349  !! so only adding it if it doesn't exist already
1350  if ( .not. dimension_exists(fms2io_fileobj, "nv")) then
1351  call register_axis(fms2io_fileobj, "nv", 2) !< Time bounds need a vertex number
1352  call write_var_metadata(fms2io_fileobj, "nv", dimensions(1:1), &
1353  "vertex number", no_units)
1354  endif
1355  call write_var_metadata(fms2io_fileobj, time_var_name//"_bnds", dimensions, &
1356  trim(time_var_name)//" axis boundaries", time_units_str)
1357  endif
1358 
1359 end subroutine write_time_metadata
1360 
1361 !> \brief Write out the field data to the file
1362 subroutine write_field_data(this, field_obj, buffer_obj, unlim_dim_was_increased)
1363  class(fmsdiagfilecontainer_type), intent(in), target :: this !< The diag file object to write to
1364  type(fmsdiagfield_type), intent(in), target :: field_obj !< The field object to write from
1365  type(fmsdiagoutputbuffer_type), intent(inout), target :: buffer_obj !< The buffer object with the data
1366  logical, intent(inout) :: unlim_dim_was_increased
1367 
1368  class(fmsdiagfile_type), pointer :: diag_file !< Diag_file object to open
1369  class(fmsnetcdffile_t), pointer :: fms2io_fileobj !< Fileobj to write to
1370  logical :: has_diurnal !< indicates if theres a diurnal axis to adjust for
1371 
1372  diag_file => this%FMS_diag_file
1373  fms2io_fileobj => diag_file%fms2io_fileobj
1374 
1375  !< Increase the unlim dimension index for the output buffer and update the output buffer for the file
1376  !! if haven't already
1377  call buffer_obj%increase_unlim_dim()
1378  if (buffer_obj%get_unlim_dim() > diag_file%unlim_dimension_level) then
1379  diag_file%unlim_dimension_level = buffer_obj%get_unlim_dim()
1380  unlim_dim_was_increased = .true.
1381  endif
1382 
1383  !TODO This may be offloaded in the future
1384  if (diag_file%is_static) then
1385  !< Here the file is static so there is no need for the unlimited dimension
1386  !! as a variables are static
1387  call buffer_obj%write_buffer(fms2io_fileobj)
1388  diag_file%data_has_been_written = .true.
1389  else
1390  if (field_obj%is_static()) then
1391  !< If the variable is static, only write it the first time
1392  if (buffer_obj%get_unlim_dim() .eq. 1) then
1393  call buffer_obj%write_buffer(fms2io_fileobj)
1394  diag_file%data_has_been_written = .true.
1395  endif
1396  else
1397  if (unlim_dim_was_increased) diag_file%data_has_been_written = .true.
1398  has_diurnal = buffer_obj%get_diurnal_sample_size() .gt. 1
1399  call buffer_obj%write_buffer(fms2io_fileobj, &
1400  unlim_dim_level=buffer_obj%get_unlim_dim(), is_diurnal=has_diurnal)
1401  endif
1402  endif
1403 
1404 end subroutine write_field_data
1405 
1406 !> \brief Determine if it is time to close the file
1407 !! \return .True. if it is time to close the file
1408 logical function is_time_to_close_file (this, time_step, force_close)
1409  class(fmsdiagfilecontainer_type), intent(in), target :: this !< The file object
1410  TYPE(time_type), intent(in) :: time_step !< Current model step time
1411  logical, intent(in) :: force_close !< if .true. return true
1412 
1413  if (force_close) then
1414  is_time_to_close_file = .true.
1415  elseif (time_step >= this%FMS_diag_file%next_close) then
1416  is_time_to_close_file = .true.
1417  else
1418  if (this%FMS_diag_file%is_static) then
1419  is_time_to_close_file = .true.
1420  else
1421  is_time_to_close_file = .false.
1422  endif
1423  endif
1424 end function
1425 
1426 !> \brief Determine if it is time to start doing mathz
1427 !! \return .True. if it is time to start doing mathz
1428 logical function time_to_start_doing_math (this)
1429  class(fmsdiagfilecontainer_type), intent(in), target :: this !< The file object
1430  time_to_start_doing_math = .false.
1431  if (this%FMS_diag_file%model_time >= this%FMS_diag_file%start_time) then
1432  time_to_start_doing_math = .true.
1433  endif
1434 end function
1435 
1436 !> \brief Determine if it is time to "write" to the file
1437 subroutine check_file_times(this, time_step, output_buffers, diag_fields, do_not_write)
1438  class(fmsdiagfilecontainer_type), intent(inout), target :: this !< The file object
1439  TYPE(time_type), intent(in) :: time_step !< Current model step time
1440  type(fmsdiagoutputbuffer_type), intent(in) :: output_buffers(:) !< Array of output buffer.
1441  !! This is needed for error messages!
1442  type(fmsdiagfield_type), intent(in) :: diag_fields(:) !< Array of diag_fields objects
1443  logical, intent(out) :: do_not_write !< .True. only if this is not a new
1444  !! time step and you are writting
1445  !! at every time step
1446 
1447  do_not_write = .false.
1448  if (time_step > this%FMS_diag_file%next_output) then
1449  if (this%FMS_diag_file%is_static) return
1450  if (time_step > this%FMS_diag_file%next_next_output) then
1451  if (this%FMS_diag_file%get_file_freq() .eq. 0) then
1452  !! If the diag file is being written at every time step
1453  if (time_step .ne. this%FMS_diag_file%next_output) then
1454  !! Only write and update the next_output if it is a new time
1455  call this%FMS_diag_file%check_buffer_times(output_buffers, diag_fields)
1456  this%FMS_diag_file%next_output = time_step
1457  this%FMS_diag_file%next_next_output = time_step
1458  endif
1459  elseif (this%FMS_diag_file%num_registered_fields .eq. 0) then
1460  !! If no variables have been registered, write a dummy time dimension for the first level
1461  !! At least one time level is needed for the combiner to work ...
1462  if (this%FMS_diag_file%unlim_dimension_level .eq. 0) then
1463  call mpp_error(note, this%FMS_diag_file%get_file_fname()//&
1464  ": diag_manager_mod: This file does not have any variables registered. Fill values will be written")
1465  this%FMS_diag_file%data_has_been_written = .true.
1466  this%FMS_diag_file%unlim_dimension_level = 1
1467  endif
1468  else
1469  !! Only fail if send data has actually been called for at least one variable
1470  if (this%FMS_diag_file%has_send_data_been_called(output_buffers, .false.)) &
1471  call mpp_error(fatal, this%FMS_diag_file%get_file_fname()//&
1472  ": diag_manager_mod: You skipped a time_step. Be sure that diag_send_complete is called at every "//&
1473  "time_step needed by the file.")
1474  endif
1475  endif
1476  else
1477  if(this%FMS_diag_file%get_file_freq() .eq. 0) then
1478  do_not_write = .true.
1479  endif
1480  endif
1481 end subroutine check_file_times
1482 
1483 !> \brief Determine if the current PE has data to write
1484 logical function writing_on_this_pe(this)
1485  class(fmsdiagfilecontainer_type), intent(in), target :: this !< The file object
1486 
1487  select type(diag_file => this%FMS_diag_file)
1488  type is (subregionalfile_type)
1489  writing_on_this_pe = diag_file%write_on_this_pe
1490  class default
1491  writing_on_this_pe = .true.
1492  end select
1493 
1494 end function
1495 
1496 !> \brief Write out the time data to the file
1497 subroutine write_time_data(this)
1498  class(fmsdiagfilecontainer_type), intent(in), target :: this !< The file object
1499 
1500  real :: dif !< The time as a real number
1501  class(fmsdiagfile_type), pointer :: diag_file !< Diag_file object to open
1502  class(fmsnetcdffile_t), pointer :: fms2io_fileobj !< The fileobj to write to
1503  TYPE(time_type) :: middle_time !< The middle time of the averaging period
1504 
1505  real :: t1 !< The beginning time of the averaging period
1506  real :: t2 !< The ending time of the averaging period
1507 
1508  diag_file => this%FMS_diag_file
1509  fms2io_fileobj => diag_file%fms2io_fileobj
1510 
1511  !< If data has not been written for the current unlimited dimension leave the subroutine
1512  if (.not. diag_file%data_has_been_written) return
1513 
1514  if (diag_file%get_time_ops()) then
1515  middle_time = (diag_file%last_output+diag_file%next_output)/2
1516  dif = get_date_dif(middle_time, get_base_time(), diag_file%get_file_timeunit())
1517  else
1518  dif = get_date_dif(diag_file%next_output, get_base_time(), diag_file%get_file_timeunit())
1519  endif
1520 
1521  call write_data(fms2io_fileobj, diag_file%get_file_unlimdim(), dif, &
1522  unlim_dim_level=diag_file%unlim_dimension_level)
1523 
1524  if (diag_file%get_time_ops()) then
1525  t1 = get_date_dif(diag_file%last_output, get_base_time(), diag_file%get_file_timeunit())
1526  t2 = get_date_dif(diag_file%next_output, get_base_time(), diag_file%get_file_timeunit())
1527 
1528  call write_data(fms2io_fileobj, trim(diag_file%get_file_unlimdim())//"_bnds", &
1529  (/t1, t2/), unlim_dim_level=diag_file%unlim_dimension_level)
1530 
1531  if (diag_file%unlim_dimension_level .eq. 1) then
1532  call write_data(fms2io_fileobj, "nv", (/1, 2/))
1533  endif
1534  endif
1535 
1536  diag_file%data_has_been_written = .false.
1537 end subroutine write_time_data
1538 
1539 !> \brief Updates the current_new_file_freq_index if using a new_file_freq
1540 subroutine update_current_new_file_freq_index(this, time_step)
1541  class(fmsdiagfilecontainer_type), intent(inout), target :: this !< The file object
1542  TYPE(time_type), intent(in) :: time_step !< Current model step time
1543 
1544  class(fmsdiagfile_type), pointer :: diag_file !< Diag_file object to open
1545 
1546  diag_file => this%FMS_diag_file
1547 
1548  if (time_step >= diag_file%no_more_data) then
1549  call diag_file%diag_yaml_file%increase_new_file_freq_index()
1550 
1551  if (diag_file%has_file_duration()) then
1552  diag_file%no_more_data = diag_time_inc(diag_file%no_more_data, diag_file%get_file_duration(), &
1553  diag_file%get_file_duration_units())
1554  else
1555  !< At this point you are done writing data
1556  diag_file%done_writing_data = .true.
1557  diag_file%no_more_data = diag_time_inc(diag_file%no_more_data, very_large_file_freq, diag_days)
1558  diag_file%next_output = diag_file%no_more_data
1559  diag_file%next_next_output = diag_file%no_more_data
1560  diag_file%last_output = diag_file%no_more_data
1561  endif
1562  endif
1563 
1564  if (diag_file%is_static) diag_file%done_writing_data = .true.
1565 end subroutine update_current_new_file_freq_index
1566 
1567 !> \brief Set up the next_output and next_next_output variable in a file obj
1568 subroutine update_next_write(this, time_step)
1569  class(fmsdiagfilecontainer_type), intent(in), target :: this !< The file object
1570  TYPE(time_type), intent(in) :: time_step !< Current model step time
1571 
1572  class(fmsdiagfile_type), pointer :: diag_file !< Diag_file object to open
1573 
1574  diag_file => this%FMS_diag_file
1575  if (diag_file%is_static) then
1576  diag_file%last_output = diag_file%next_output
1577  diag_file%next_output = diag_time_inc(diag_file%next_output, very_large_file_freq, diag_days)
1578  diag_file%next_next_output = diag_time_inc(diag_file%next_output, very_large_file_freq, diag_days)
1579  else
1580  diag_file%last_output = diag_file%next_output
1581  diag_file%next_output = diag_time_inc(diag_file%next_output, diag_file%get_file_freq(), &
1582  diag_file%get_file_frequnit())
1583  diag_file%next_next_output = diag_time_inc(diag_file%next_output, diag_file%get_file_freq(), &
1584  diag_file%get_file_frequnit())
1585  endif
1586 
1587 end subroutine update_next_write
1588 
1589 !> \brief Prepare the diag file for the force_write
1590 subroutine prepare_for_force_write(this)
1591  class(fmsdiagfilecontainer_type), intent(inout), target :: this !< The file object
1592 
1593  if (this%FMS_diag_file%unlim_dimension_level .eq. 0) then
1594  this%FMS_diag_file%unlim_dimension_level = 1
1595  this%FMS_diag_file%data_has_been_written = .true.
1596  endif
1597 end subroutine prepare_for_force_write
1598 
1599 !> \brief Initialize the unlim dimension in the file and in its buffer objects to 0
1600 subroutine init_unlim_dim(this, output_buffers)
1601  class(fmsdiagfilecontainer_type), intent(inout), target :: this !< The file object
1602  type(fmsdiagoutputbuffer_type), intent(in), target :: output_buffers(:) !< Array of output buffer.
1603 
1604  class(fmsdiagfile_type), pointer :: diag_file !< Diag_file object
1605  type(fmsdiagoutputbuffer_type), pointer :: output_buffer_obj !< Buffer object
1606  integer :: i !< For looping through buffers
1607 
1608  diag_file => this%FMS_diag_file
1609  diag_file%unlim_dimension_level = 0
1610  do i = 1, diag_file%number_of_buffers
1611  output_buffer_obj => output_buffers(diag_file%buffer_ids(i))
1612  call output_buffer_obj%init_buffer_unlim_dim()
1613  enddo
1614 end subroutine init_unlim_dim
1615 
1616 !> \brief Get the number of time levels that were actually written to the file
1617 !! \return Number of time levels that were actually written to the file
1618 function get_num_time_levels(this) &
1619 result(res)
1620  class(fmsdiagfilecontainer_type), intent(in), target :: this !< The file object
1621  integer :: res
1622 
1623  if (this%is_regional()) then
1624  !! If this is a subregional file, num_time_levels will be set to 0 for
1625  !! all PEs that are not in the subregion. If the root pe is not in the subregion
1626  !! then the num_time_levels will not be correct, so this is just getting the max
1627  !! from all PEs
1628  res = this%FMS_diag_file%num_time_levels
1629  call mpp_max(res)
1630  else
1631  res = this%FMS_diag_file%num_time_levels
1632  endif
1633 end function
1634 
1635 !> \brief Get the number of tiles in the file's domain
1636 !! \return Number of tiles in the file's domain
1637 function get_num_tiles(this) &
1638 result(res)
1639  class(fmsdiagfilecontainer_type), intent(in), target :: this !< The file object
1640  integer :: res
1641 
1642  select case(this%FMS_diag_file%type_of_domain)
1643  case (two_d_domain, ug_domain)
1644  select type(domain => this%FMS_diag_file%domain)
1645  type is (diagdomain2d_t)
1646  res = mpp_get_ntile_count(domain%Domain2)
1647  type is (diagdomainug_t)
1648  res = mpp_get_ug_domain_ntiles(domain%DomainUG)
1649  end select
1650  case default
1651  res = 1
1652  end select
1653 end function get_num_tiles
1654 
1655 !> \brief Get number of distributed files that were written
1656 !! \return The number of distributed files that were written
1657 function get_ndistributedfiles(this) &
1658 result(res)
1659  class(fmsdiagfilecontainer_type), intent(in), target :: this !< The file object
1660  integer :: res
1661  integer :: io_layout(2) !< The io_layout
1662 
1663  select case(this%FMS_diag_file%type_of_domain)
1664  case (two_d_domain, ug_domain)
1665  select type(domain => this%FMS_diag_file%domain)
1666  type is (diagdomain2d_t)
1667  io_layout = mpp_get_io_domain_layout(domain%Domain2)
1668  res = io_layout(1) * io_layout(2)
1669  type is (diagdomainug_t)
1670  res = mpp_get_io_domain_ug_layout(domain%DomainUG)
1671  end select
1672  case default
1673  res = 1
1674  end select
1675 end function get_ndistributedfiles
1676 
1677 !> \brief Get the unlimited dimension level that is in the file
1678 !! \return The unlimited dimension
1679 pure function get_unlim_dimension_level(this) &
1680 result(res)
1681  class(fmsdiagfilecontainer_type), intent(in), target :: this !< The file object
1682  integer :: res
1683 
1684  res = this%FMS_diag_file%unlim_dimension_level
1685 end function
1686 
1687 !> \brief Flushes the netcdf file to disk if flush_nc_files is set to .True. in the namelist
1688 subroutine flush_diag_file(this)
1689  class(fmsdiagfilecontainer_type), intent(inout), target :: this !< The file object
1690 
1691  if (flush_nc_files) then
1692  call flush_file(this%FMS_diag_file%fms2io_fileobj)
1693  endif
1694 end subroutine flush_diag_file
1695 
1696 !> \brief Get the next_output for the file object
1697 !! \return The next_output
1698 pure function get_next_output(this) &
1699 result(res)
1700  class(fmsdiagfilecontainer_type), intent(in), target :: this !< The file object
1701  type(time_type) :: res
1702 
1703  res = this%FMS_diag_file%next_output
1704 end function get_next_output
1705 
1706 !> \brief Get the next_output for the file object
1707 !! \return The next_output
1708 pure function get_next_next_output(this) &
1709 result(res)
1710  class(fmsdiagfilecontainer_type), intent(in), target :: this !< The file object
1711  type(time_type) :: res
1712 
1713  res = this%FMS_diag_file%next_next_output
1714  if (this%FMS_diag_file%is_static) then
1715  res = this%FMS_diag_file%no_more_data
1716  endif
1717 end function get_next_next_output
1718 
1719 !< @brief Writes the axis metadata for the file
1720 subroutine write_axis_metadata(this, diag_axis)
1721  class(fmsdiagfilecontainer_type), intent(inout), target :: this !< The file object
1722  class(fmsdiagaxiscontainer_type), intent(in), target :: diag_axis(:) !< Diag_axis object
1723 
1724  class(fmsdiagfile_type), pointer :: diag_file !< Diag_file object to open
1725  class(fmsnetcdffile_t), pointer :: fms2io_fileobj !< The fileobj to write to
1726  integer :: i,k !< For do loops
1727  integer :: parent_axis_id !< Id of the parent_axis
1728  integer :: structured_ids(2) !< Ids of the uncompress axis
1729  integer :: edges_id !< Id of the axis edge
1730 
1731  class(fmsdiagaxiscontainer_type), pointer :: axis_ptr !< pointer to the axis object currently writing
1732  logical :: edges_in_file !< .true. if the edges are already in the file
1733 
1734  diag_file => this%FMS_diag_file
1735  fms2io_fileobj => diag_file%fms2io_fileobj
1736 
1737  do i = 1, diag_file%number_of_axis
1738  edges_in_file = .false.
1739  axis_ptr => diag_axis(diag_file%axis_ids(i))
1740  parent_axis_id = axis_ptr%axis%get_parent_axis_id()
1741 
1742  edges_id = axis_ptr%axis%get_edges_id()
1743  if (edges_id .ne. diag_null) then
1744  !< write the edges if is not in the list of axis in the file, otherwrise ignore
1745  if (any(diag_file%axis_ids(1:diag_file%number_of_axis) .eq. edges_id)) then
1746  edges_in_file = .true.
1747  else
1748  call diag_axis(edges_id)%axis%write_axis_metadata(fms2io_fileobj, .true.)
1749  call diag_file%add_new_axis(edges_id)
1750  endif
1751  endif
1752 
1753  if (parent_axis_id .eq. diag_null) then
1754  call axis_ptr%axis%write_axis_metadata(fms2io_fileobj, edges_in_file)
1755  else
1756  call axis_ptr%axis%write_axis_metadata(fms2io_fileobj, edges_in_file, diag_axis(parent_axis_id)%axis)
1757  endif
1758 
1759  if (axis_ptr%axis%is_unstructured_grid()) then
1760  structured_ids = axis_ptr%axis%get_structured_axis()
1761  do k = 1, size(structured_ids)
1762  call diag_axis(structured_ids(k))%axis%write_axis_metadata(fms2io_fileobj, .false.)
1763  enddo
1764  endif
1765 
1766  enddo
1767 
1768 end subroutine write_axis_metadata
1769 
1770 !< @brief Writes the field metadata for the file
1771 subroutine write_field_metadata(this, diag_field, diag_axis)
1772  class(fmsdiagfilecontainer_type), intent(inout), target :: this !< The file object
1773  class(fmsdiagfield_type) , intent(inout), target :: diag_field(:) !<
1774  class(fmsdiagaxiscontainer_type), intent(in) :: diag_axis(:) !< Diag_axis object
1775 
1776  class(fmsnetcdffile_t), pointer :: fms2io_fileobj !< The fileobj to write to
1777  class(fmsdiagfile_type), pointer :: diag_file !< Diag_file object to open
1778  class(fmsdiagfield_type), pointer :: field_ptr !< diag_field(diag_file%field_ids(i)), for convenience
1779 
1780  integer :: i !< For do loops
1781  logical :: is_regional !< Flag indicating if the field is in a regional file
1782  character(len=255) :: cell_measures !< cell_measures attributes for the field
1783  logical :: need_associated_files !< .True. if the 'associated_files' global attribute is needed
1784  character(len=FMS_FILE_LEN) :: associated_files !< Associated files attribute to add
1785 
1786  is_regional = this%is_regional()
1787 
1788  diag_file => this%FMS_diag_file
1789  fms2io_fileobj => diag_file%fms2io_fileobj
1790 
1791  associated_files = ""
1792  need_associated_files = .false.
1793  do i = 1, size(diag_file%field_ids)
1794  if (.not. diag_file%field_registered(i)) cycle !TODO do something else here
1795  field_ptr => diag_field(diag_file%field_ids(i))
1796 
1797  cell_measures = ""
1798  if (field_ptr%has_area()) then
1799  cell_measures = "area: "//diag_field(field_ptr%get_area())%get_varname(to_write=.true.)
1800 
1801  !! Determine if the area field is already in the file. If it is not create the "associated_files" attribute
1802  !! which contains the file name of the file the area field is in. This is needed for PP/fregrid.
1803  if (.not. diag_field(field_ptr%get_area())%is_variable_in_file(diag_file%id)) then
1804  need_associated_files = .true.
1805  call diag_field(field_ptr%get_area())%generate_associated_files_att(associated_files, diag_file%start_time)
1806  endif
1807  endif
1808 
1809  if (field_ptr%has_volume()) then
1810  cell_measures = trim(cell_measures)//" volume: "//diag_field(field_ptr%get_volume())%get_varname(to_write=.true.)
1811 
1812  !! Determine if the volume field is already in the file. If it is not create the "associated_files" attribute
1813  !! which contains the file name of the file the volume field is in. This is needed for PP/fregrid.
1814  if (.not. diag_field(field_ptr%get_volume())%is_variable_in_file(diag_file%id)) then
1815  need_associated_files = .true.
1816  call diag_field(field_ptr%get_volume())%generate_associated_files_att(associated_files, diag_file%start_time)
1817  endif
1818  endif
1819 
1820  call field_ptr%write_field_metadata(fms2io_fileobj, diag_file%id, diag_file%yaml_ids(i), diag_axis, &
1821  this%FMS_diag_file%get_file_unlimdim(), is_regional, cell_measures)
1822  enddo
1823 
1824  if (need_associated_files) &
1825  call register_global_attribute(fms2io_fileobj, "associated_files", trim(adjustl(associated_files)), &
1826  str_len=len_trim(adjustl(associated_files)))
1827 
1828 end subroutine write_field_metadata
1829 
1830 !< @brief Writes the axis data for the file
1831 subroutine write_axis_data(this, diag_axis)
1832  class(fmsdiagfilecontainer_type), intent(inout), target :: this !< The file object
1833  class(fmsdiagaxiscontainer_type), intent(in) :: diag_axis(:) !< Diag_axis object
1834 
1835  class(fmsdiagfile_type), pointer :: diag_file !< Diag_file object to open
1836  class(fmsnetcdffile_t), pointer :: fms2io_fileobj !< The fileobj to write to
1837  integer :: i, k !< For do loops
1838  integer :: j !< diag_file%axis_ids(i) (for less typing)
1839  integer :: parent_axis_id !< Id of the parent_axis
1840  integer :: structured_ids(2) !< Ids of the uncompress axis
1841 
1842  diag_file => this%FMS_diag_file
1843  fms2io_fileobj => diag_file%fms2io_fileobj
1844 
1845  do i = 1, diag_file%number_of_axis
1846  j = diag_file%axis_ids(i)
1847  parent_axis_id = diag_axis(j)%axis%get_parent_axis_id()
1848  if (parent_axis_id .eq. diag_null) then
1849  call diag_axis(j)%axis%write_axis_data(fms2io_fileobj)
1850  else
1851  call diag_axis(j)%axis%write_axis_data(fms2io_fileobj, diag_axis(parent_axis_id)%axis)
1852  endif
1853 
1854  if (diag_axis(j)%axis%is_unstructured_grid()) then
1855  structured_ids = diag_axis(j)%axis%get_structured_axis()
1856  do k = 1, size(structured_ids)
1857  call diag_axis(structured_ids(k))%axis%write_axis_data(fms2io_fileobj)
1858  enddo
1859  endif
1860  enddo
1861 
1862 end subroutine write_axis_data
1863 
1864 !< @brief Closes the diag_file
1865 subroutine close_diag_file(this, output_buffers, model_end_time, diag_fields)
1866  class(fmsdiagfilecontainer_type), intent(inout), target :: this !< The file object
1867  type(fmsdiagoutputbuffer_type), intent(in) :: output_buffers(:) !< Array of output buffers
1868  !! This is needed for error checking
1869  type(time_type), intent(in) :: model_end_time !< Time that simulation ends
1870  type(fmsdiagfield_type), intent(in), optional :: diag_fields(:) !< Array of diag fields
1871  !! This is needed for error checking
1872 
1873  if (.not. this%FMS_diag_file%is_file_open) return
1874 
1875  !< The select types are needed here because otherwise the code will go to the
1876  !! wrong close_file routine and things will not close propertly
1877  select type( fms2io_fileobj => this%FMS_diag_file%fms2io_fileobj)
1878  type is (fmsnetcdfdomainfile_t)
1879  call close_file(fms2io_fileobj)
1880  type is (fmsnetcdffile_t)
1881  call close_file(fms2io_fileobj)
1882  type is (fmsnetcdfunstructureddomainfile_t)
1883  call close_file(fms2io_fileobj)
1884  end select
1885 
1886  !< Keep track of the number of time levels that were written to the file
1887  !! If the file is using the new_file_freq key, it will be closing the file multiple
1888  !! time, so this ensures that we keep a running count
1889  !! This is going be written to the output yaml, after all the files are closed
1890  this%FMS_diag_file%num_time_levels = this%FMS_diag_file%num_time_levels + &
1891  this%FMS_diag_file%unlim_dimension_level
1892 
1893  !< Reset the unlimited dimension level back to 0, in case the fms2io_fileobj is re-used
1894  this%FMS_diag_file%unlim_dimension_level = 0
1895  this%FMS_diag_file%is_file_open = .false.
1896 
1897  if (this%FMS_diag_file%has_file_new_file_freq()) then
1898  this%FMS_diag_file%next_close = diag_time_inc(this%FMS_diag_file%next_close, &
1899  this%FMS_diag_file%get_file_new_file_freq(), &
1900  this%FMS_diag_file%get_file_new_file_freq_units())
1901  else
1902  this%FMS_diag_file%next_close = model_end_time
1903  endif
1904 
1905  if (this%FMS_diag_file%model_time >= model_end_time) &
1906  this%FMS_diag_file%done_writing_data = .true.
1907  if (this%FMS_diag_file%has_send_data_been_called(output_buffers, .true., diag_fields)) return
1908 end subroutine close_diag_file
1909 
1910 !> \brief Set the model time for the diag file object
1911 subroutine set_model_time(this, model_time)
1912  class(fmsdiagfilecontainer_type), intent(inout) :: this !< The file object
1913  type(time_type), intent(in) :: model_time !< Model time to add
1914 
1915  if (model_time > this%FMS_diag_file%model_time) this%FMS_diag_file%model_time = model_time
1916 end subroutine
1917 
1918 !> \brief Get the model time from the file object
1919 !! \result A pointer to the model time
1920 function get_model_time(this) &
1921  result(rslt)
1922  class(fmsdiagfilecontainer_type), intent(inout), target :: this !< The file object
1923  type(time_type), pointer :: rslt
1924 
1925  rslt => this%FMS_diag_file%model_time
1926 end function get_model_time
1927 
1928 !> \brief Gets the buffer_id list from the file object
1929 pure function get_buffer_ids (this)
1930  class(fmsdiagfile_type), intent(in) :: this !< The file object
1931  integer, allocatable :: get_buffer_ids(:) !< returned buffer ids for this file
1932 
1933  allocate(get_buffer_ids(this%number_of_buffers))
1934  get_buffer_ids = this%buffer_ids(1:this%number_of_buffers)
1935 end function get_buffer_ids
1936 
1937 !> Gets the stored number of buffers from the file object
1938 pure function get_number_of_buffers(this)
1939  class(fmsdiagfile_type), intent(in) :: this !< file object
1940  integer :: get_number_of_buffers !< returned number of buffers
1941  get_number_of_buffers = this%number_of_buffers
1942 end function get_number_of_buffers
1943 
1944 !> Check to ensure that send_data was called at the time step for every output buffer in the file
1945 !! This is only needed when you are output data at every time step
1946 subroutine check_buffer_times(this, output_buffers, diag_fields)
1947  class(fmsdiagfile_type), intent(in) :: this !< file object
1948  type(fmsdiagoutputbuffer_type), intent(in), target :: output_buffers(:) !< Array of output buffers
1949  type(fmsdiagfield_type), intent(in) :: diag_fields(:) !< Array of diag_fields
1950 
1951  integer :: i !< For do loop
1952  type(time_type) :: current_buffer_time !< The buffer time for the current buffer in the do loop
1953  character(len=:), allocatable :: field_name !< The field name (for error messages)
1954  logical :: buffer_time_set !< .True. if current_buffer_time has been set
1955  type(fmsdiagoutputbuffer_type), pointer :: output_buffer_obj !< Pointer to the output buffer
1956 
1957  buffer_time_set = .false.
1958  do i = 1, this%number_of_buffers
1959  output_buffer_obj => output_buffers(this%buffer_ids(i))
1960  if (diag_fields(output_buffer_obj%get_field_id())%is_static()) cycle
1961  if (.not. buffer_time_set) then
1962  current_buffer_time = output_buffer_obj%get_buffer_time()
1963  field_name = output_buffer_obj%get_buffer_name()
1964  buffer_time_set = .true.
1965  else
1966  if (current_buffer_time .ne. output_buffer_obj%get_buffer_time()) &
1967  call mpp_error(fatal, "Send data has not been called at the same time steps for the fields:"//&
1968  field_name//" and "//output_buffer_obj%get_buffer_name()//&
1969  " in file:"//this%get_file_fname())
1970  endif
1971  enddo
1972 end subroutine
1973 
1974 !> @brief Determine if send_data has been called for any fields in the file. Prints out warnings, if indicated
1975 !! @return .True. if send_data has been called for any fields in the file
1976 function has_send_data_been_called(this, output_buffers, print_warnings, diag_fields) &
1977 result(rslt)
1978  class(fmsdiagfile_type), intent(in) :: this !< file object
1979  type(fmsdiagoutputbuffer_type), intent(in), target :: output_buffers(:) !< Array of output buffers
1980  logical, intent(in) :: print_warnings !< .True. if printing warnings
1981  type(fmsdiagfield_type), intent(in), optional :: diag_fields(:) !< Array of diag fields
1982 
1983  logical :: rslt
1984  integer :: i !< For do loops
1985  integer :: field_id !< Field id
1986 
1987  rslt = .false.
1988 
1989  if (print_warnings) then
1990  do i = 1, this%number_of_buffers
1991  if (.not. output_buffers(this%buffer_ids(i))%is_there_data_to_write()) then
1992  field_id = output_buffers(this%buffer_ids(i))%get_field_id()
1993  call mpp_error(note, "Send data was never called for field:"//&
1994  trim(diag_fields(field_id)%get_varname())//" mod: "//trim(diag_fields(field_id)%get_modname())//&
1995  " in file: "//trim(this%get_file_fname())//". Writting FILL VALUES!")
1996  endif
1997  enddo
1998  else
1999  do i = 1, this%number_of_buffers
2000  if (output_buffers(this%buffer_ids(i))%is_there_data_to_write()) then
2001  rslt = .true.
2002  return
2003  endif
2004  enddo
2005  endif
2006 end function has_send_data_been_called
2007 #endif
2008 end module fms_diag_file_object_mod
character(len=8) no_units
String indicating that the variable has no units.
Definition: diag_data.F90:125
integer, parameter sub_regional
This is a file with a sub_region use the FmsNetcdfFile_t fileobj.
Definition: diag_data.F90:104
integer, parameter max_str_len
Max length for a string.
Definition: diag_data.F90:129
integer function get_base_minute()
gets the module variable base_minute
Definition: diag_data.F90:551
integer function get_base_year()
gets the module variable base_year
Definition: diag_data.F90:519
integer function get_base_hour()
gets the module variable base_hour
Definition: diag_data.F90:543
logical flush_nc_files
Control if diag_manager will force a flush of the netCDF file on each write. Note: changing this to ....
Definition: diag_data.F90:365
integer, parameter no_domain
Use the FmsNetcdfFile_t fileobj.
Definition: diag_data.F90:101
character(len=7) avg_name
Name of the average fields.
Definition: diag_data.F90:124
integer, parameter end_time
Use the end of the time average bounds.
Definition: diag_data.F90:128
character(len=6) pack_size_str
Pack size as a string to be used in fms2_io register call set to "double" or "float".
Definition: diag_data.F90:409
type(time_type) function get_base_time()
gets the module variable base_time
Definition: diag_data.F90:511
integer max_axes
Maximum number of independent axes.
Definition: diag_data.F90:361
integer function get_base_day()
gets the module variable base_day
Definition: diag_data.F90:535
type(time_type) diag_init_time
Time diag_manager_init called. If init_time not included in diag_manager_init call,...
Definition: diag_data.F90:417
integer, parameter time_min
The reduction method is min value.
Definition: diag_data.F90:117
integer, parameter ug_domain
Use the FmsNetcdfUnstructuredDomainFile_t fileobj.
Definition: diag_data.F90:103
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
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 begin_time
Use the begining of the time average bounds.
Definition: diag_data.F90:126
integer, parameter time_average
The reduction method is average of values.
Definition: diag_data.F90:120
integer function get_base_month()
gets the module variable base_month
Definition: diag_data.F90:527
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 middle_time
Use the middle of the time average bounds.
Definition: diag_data.F90:127
integer, parameter time_none
There is no reduction method.
Definition: diag_data.F90:116
integer function get_base_second()
gets the module variable base_second
Definition: diag_data.F90:559
integer, parameter time_max
The reduction method is max value.
Definition: diag_data.F90:118
integer, parameter two_d_domain
Use the FmsNetcdfDomainFile_t fileobj.
Definition: diag_data.F90:102
Close a netcdf or domain file opened with open_file or open_virtual_file.
Definition: fms2_io.F90:166
Opens a given netcdf or domain file.
Definition: fms2_io.F90:122
Add a dimension to a given file.
Definition: fms2_io.F90:190
Defines a new field within the given file Example usage:
Definition: fms2_io.F90:212
Write data to a defined field within a file Example usage:
Definition: fms2_io.F90:262
character(len=128) function, public get_time_string(filename, current_time)
This function determines a string based on current time. This string is used as suffix in output file...
real function, public get_date_dif(t2, t1, units)
Return the difference between two times in units.
subroutine, public define_new_subaxis_latlon(diag_axis, axis_ids, naxis, subRegion, is_cube_sphere, write_on_this_pe)
Fill in the subaxis object for a subRegion defined by lat lon.
subroutine, public define_new_subaxis_index(parent_axis, subRegion, diag_axis, naxis, is_x_or_y, write_on_this_pe)
Fill in the subaxis object for a subRegion defined by index.
subroutine, public define_diurnal_axis(diag_axis, naxis, n_diurnal_samples, is_edges)
Defined a new diurnal axis.
subroutine, public create_new_z_subaxis(zbounds, var_axis_ids, diag_axis, naxis, file_axis_id, nfile_axis, nz_subaxis)
Creates a new z subaxis to use.
logical function, public is_parent_axis(axis_id, parent_axis_id, diag_axis)
Determine if the diag_axis(parent_axis_id) is the parent of diag_axis(axis_id)
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.
Object that holds the information of the diag_yaml.
integer function mpp_get_ntile_count(domain)
Returns number of tiles in mosaic.
integer function, dimension(2) mpp_get_io_domain_layout(domain)
Set user stack size.
integer function stdout()
This function returns the current standard fortran unit numbers for output.
Definition: mpp_util.inc:43
integer function mpp_npes()
Returns processor count for current pelist.
Definition: mpp_util.inc:421
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:407
Error handler.
Definition: mpp.F90:382
Reduction operations. Find the max of scalar a from the PEs in pelist result is also automatically br...
Definition: mpp.F90:538
character(len=15) function, public date_to_string(time, err_msg)
Get the a character string that represents the time. The format will be yyyymmdd.hhmmss.
character(len=24) function, public valid_calendar_types(ncal, err_msg)
Returns a character string that describes the calendar type corresponding to the input integer.
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...
integer function, public get_calendar_type()
Returns default calendar type for mapping from time to date.
Type to represent amounts of time. Implemented as seconds and days to allow for larger intervals.
Contains buffer types and routines for the diag manager.
Type to hold the domain info for an axis This type was created to avoid having to send in "Domain",...
Type to hold the unstructured 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.
A container for fmsDiagFile_type. This is used to create the array of files.
type to hold the diag_file information
type to hold the info a diag_field
type to hold the sub region information about a file