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