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