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