1 !***********************************************************************
2 !* GNU Lesser General Public License
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 !> @defgroup fms_diag_output_yaml_mod fms_diag_output_yaml_mod
20 !> @ingroup diag_manager
21 !! @brief fms_diag_file_object_mod handles the file objects data, functions, and subroutines.
22 !! @author Tom Robinson
23 !! @description The fmsDiagFile_type contains the information for each history file to be written. It has
24 !! a pointer to the information from the diag yaml, additional metadata that comes from the model, and a
25 !! list of the variables and their variable IDs that are in the file.
26 module fms_diag_file_object_mod
28 use fms2_io_mod, only: FmsNetcdfFile_t, FmsNetcdfUnstructuredDomainFile_t, FmsNetcdfDomainFile_t, &
29 get_instance_filename, open_file, close_file, get_mosaic_tile_file, unlimited, &
30 register_axis, register_field, register_variable_attribute, write_data, &
31 dimension_exists, register_global_attribute
32 use diag_data_mod, only: DIAG_NULL, NO_DOMAIN, max_axes, SUB_REGIONAL, get_base_time, DIAG_NOT_REGISTERED, &
33 TWO_D_DOMAIN, UG_DOMAIN, prepend_date, DIAG_DAYS, VERY_LARGE_FILE_FREQ, &
34 get_base_year, get_base_month, get_base_day, get_base_hour, get_base_minute, &
35 get_base_second, time_unit_list, time_average, time_rms, time_max, time_min, time_sum, &
36 time_diurnal, time_power, time_none, avg_name, no_units, pack_size_str, &
37 middle_time, begin_time, end_time, MAX_STR_LEN, index_gridtype, latlon_gridtype
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, &
44 fmsDiagAxisContainer_type, DIAGDOMAIN2D_T, DIAGDOMAINUG_T, &
45 fmsDiagFullAxis_type, define_diurnal_axis, &
46 fmsDiagDiurnalAxis_type, create_new_z_subaxis, is_parent_axis, &
47 define_new_subaxis_latlon, define_new_subaxis_index, fmsDiagSubAxis_type
48 use fms_diag_field_object_mod, only: fmsDiagField_type
49 use fms_diag_output_buffer_mod, only: fmsDiagOutputBuffer_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
56 public :: fmsDiagFileContainer_type
57 public :: fmsDiagFile_type, fms_diag_files_object_init, fms_diag_files_object_initialized
59 logical :: fms_diag_files_object_initialized = .false.
61 integer, parameter :: var_string_len = 25
63 type :: fmsDiagFile_type
65 integer :: id !< The number associated with this file in the larger array of files
66 TYPE(time_type) :: start_time !< The start time for the file
67 TYPE(time_type) :: last_output !< Time of the last time output was writen
68 TYPE(time_type) :: next_output !< Time of the next write
69 TYPE(time_type) :: next_next_output !< Time of the next next write
70 TYPE(time_type) :: no_more_data !< Time to stop receiving data for this file
71 logical :: done_writing_data!< .True. if finished writing data
73 !< This will be used when using the new_file_freq keys in the diag_table.yaml
74 TYPE(time_type) :: next_close !< Time to close the file
75 logical :: is_file_open !< .True. if the file is opened
77 class(FmsNetcdfFile_t), allocatable :: fms2io_fileobj !< fms2_io file object for this history file
78 type(diagYamlFiles_type), pointer :: diag_yaml_file => null() !< Pointer to the diag_yaml_file data
79 integer :: type_of_domain !< The type of domain to use to open the file
80 !! NO_DOMAIN, TWO_D_DOMAIN, UG_DOMAIN, SUB_REGIONAL
81 class(diagDomain_t), pointer :: domain !< The domain to use,
82 !! null if NO_DOMAIN or SUB_REGIONAL
83 character(len=:) , dimension(:), allocatable :: file_metadata_from_model !< File metadata that comes from
85 integer, dimension(:), allocatable :: field_ids !< Variable IDs corresponding to file_varlist
86 integer, dimension(:), allocatable :: yaml_ids !< IDs corresponding to the yaml field section
87 logical, dimension(:), private, allocatable :: field_registered !< Array corresponding to `field_ids`, .true.
88 !! if the variable has been registered and
89 !! `field_id` has been set for the variable
90 integer, allocatable :: num_registered_fields !< The number of fields registered
92 integer, dimension(:), allocatable :: axis_ids !< Array of axis ids in the file
93 integer :: number_of_axis !< Number of axis in the file
94 integer, dimension(:), allocatable :: buffer_ids !< array of buffer ids associated with the file
95 integer :: number_of_buffers !< Number of buffers that have been added to the file
96 logical :: time_ops !< .True. if file contains variables that are time_min, time_max, time_average or time_sum
97 integer :: unlim_dimension_level !< The unlimited dimension level currently being written
98 logical :: data_has_been_written !< .True. if data has been written for the current unlimited dimension level
99 logical :: is_static !< .True. if the frequency is -1
100 integer :: nz_subaxis !< The number of Z axis currently added to the file
103 procedure, public :: add_field_and_yaml_id
104 procedure, public :: add_buffer_id
105 procedure, public :: is_field_registered
106 procedure, public :: init_diurnal_axis
107 procedure, public :: has_file_metadata_from_model
108 procedure, public :: has_fileobj
109 procedure, public :: has_diag_yaml_file
110 procedure, public :: set_domain_from_axis
111 procedure, public :: set_file_domain
112 procedure, public :: add_axes
113 procedure, public :: add_new_axis
114 procedure, public :: update_write_on_this_pe
115 procedure, public :: get_write_on_this_pe
116 procedure, public :: does_axis_exist
117 procedure, public :: define_new_subaxis
118 procedure, public :: add_start_time
119 procedure, public :: set_file_time_ops
120 procedure, public :: has_field_ids
121 procedure, public :: get_id
122 ! TODO procedure, public :: get_fileobj ! TODO
123 ! TODO procedure, public :: get_diag_yaml_file ! TODO
124 procedure, public :: get_file_metadata_from_model
125 procedure, public :: get_field_ids
126 ! The following fuctions come will use the yaml inquiry functions
127 procedure, public :: get_file_fname
128 procedure, public :: get_file_frequnit
129 procedure, public :: get_file_freq
130 procedure, public :: get_file_timeunit
131 procedure, public :: get_file_unlimdim
132 procedure, public :: get_file_sub_region
133 procedure, public :: get_file_sub_region_grid_type
134 procedure, public :: get_file_new_file_freq
135 procedure, public :: get_filename_time
136 procedure, public :: get_file_new_file_freq_units
137 procedure, public :: get_file_start_time
138 procedure, public :: get_file_duration
139 procedure, public :: get_file_duration_units
140 procedure, public :: get_file_varlist
141 procedure, public :: get_file_global_meta
142 procedure, public :: is_done_writing_data
143 procedure, public :: has_file_fname
144 procedure, public :: has_file_frequnit
145 procedure, public :: has_file_freq
146 procedure, public :: has_file_timeunit
147 procedure, public :: has_file_unlimdim
148 procedure, public :: has_file_sub_region
149 procedure, public :: has_file_new_file_freq
150 procedure, public :: has_file_new_file_freq_units
151 procedure, public :: has_file_start_time
152 procedure, public :: has_file_duration
153 procedure, public :: has_file_duration_units
154 procedure, public :: has_file_varlist
155 procedure, public :: has_file_global_meta
156 procedure, public :: dump_file_obj
157 procedure, public :: get_buffer_ids
158 procedure, public :: get_number_of_buffers
159 procedure, public :: has_send_data_been_called
160 end type fmsDiagFile_type
162 type, extends (fmsDiagFile_type) :: subRegionalFile_type
163 integer, dimension(:), allocatable :: sub_axis_ids !< Array of axis ids in the file
164 logical :: write_on_this_pe !< Flag indicating if the subregion is on the current PE
165 logical :: is_subaxis_defined !< Flag indicating if the subaxes have already been defined
166 end type subRegionalFile_type
168 !> \brief A container for fmsDiagFile_type. This is used to create the array of files
169 type fmsDiagFileContainer_type
170 class (fmsDiagFile_type),allocatable :: FMS_diag_file !< The individual file object
173 procedure :: is_regional
174 procedure :: is_file_static
175 procedure :: open_diag_file
176 procedure :: write_global_metadata
177 procedure :: write_time_metadata
178 procedure :: write_field_data
179 procedure :: write_axis_metadata
180 procedure :: write_field_metadata
181 procedure :: write_axis_data
182 procedure :: writing_on_this_pe
183 procedure :: is_time_to_write
184 procedure :: is_time_to_close_file
185 procedure :: write_time_data
186 procedure :: update_next_write
187 procedure :: update_current_new_file_freq_index
188 procedure :: increase_unlim_dimension_level
189 procedure :: get_unlim_dimension_level
190 procedure :: get_next_output
191 procedure :: get_next_next_output
192 procedure :: close_diag_file
193 end type fmsDiagFileContainer_type
195 !type(fmsDiagFile_type), dimension (:), allocatable, target :: FMS_diag_file !< The array of diag files
196 !class(fmsDiagFileContainer_type),dimension (:), allocatable, target :: FMS_diag_file
200 !< @brief Allocates the number of files and sets an ID based for each file
201 !! @return true if there are files allocated in the YAML object
202 logical function fms_diag_files_object_init (files_array)
203 class(fmsDiagFileContainer_type), allocatable, target, intent(inout) :: files_array (:) !< array of diag files
204 class(fmsDiagFile_type), pointer :: obj => null() !< Pointer for each member of the array
205 integer :: nFiles !< Number of files in the diag yaml
206 integer :: i !< Looping iterator
207 if (diag_yaml%has_diag_files()) then
208 nFiles = diag_yaml%size_diag_files()
209 allocate (files_array(nFiles))
210 set_ids_loop: do i= 1,nFiles
211 !> If the file has a sub_regional, define it as one and allocate the sub_axis_ids array.
212 !! This will be set in a add_axes
213 if (diag_yaml%diag_files(i)%has_file_sub_region()) then
214 allocate(subRegionalFile_type :: files_array(i)%FMS_diag_file)
215 obj => files_array(i)%FMS_diag_file
217 type is (subRegionalFile_type)
218 allocate(obj%sub_axis_ids(max_axes))
219 obj%sub_axis_ids = diag_null
220 obj%write_on_this_pe = .true.
221 obj%is_subaxis_defined = .false.
222 obj%number_of_axis = 0
225 allocate(FmsDiagFile_type::files_array(i)%FMS_diag_file)
226 obj => files_array(i)%FMS_diag_file
229 obj%diag_yaml_file => diag_yaml%diag_files(i)
231 allocate(obj%field_ids(diag_yaml%diag_files(i)%size_file_varlist()))
232 allocate(obj%buffer_ids(diag_yaml%diag_files(i)%size_file_varlist()))
233 allocate(obj%yaml_ids(diag_yaml%diag_files(i)%size_file_varlist()))
234 allocate(obj%field_registered(diag_yaml%diag_files(i)%size_file_varlist()))
235 !! Initialize the integer arrays
236 obj%field_ids = DIAG_NOT_REGISTERED
237 obj%yaml_ids = DIAG_NOT_REGISTERED
238 obj%buffer_ids = DIAG_NOT_REGISTERED
239 obj%field_registered = .FALSE.
240 obj%num_registered_fields = 0
241 obj%number_of_buffers = 0
243 !> These will be set in a set_file_domain
244 obj%type_of_domain = NO_DOMAIN
247 !> This will be set in a add_axes
248 allocate(obj%axis_ids(max_axes))
249 obj%number_of_axis = 0
251 !> Set the start_time of the file to the base_time and set up the *_output variables
252 obj%done_writing_data = .false.
253 obj%start_time = get_base_time()
254 obj%last_output = get_base_time()
255 obj%next_output = diag_time_inc(obj%start_time, obj%get_file_freq(), obj%get_file_frequnit())
256 obj%next_next_output = diag_time_inc(obj%next_output, obj%get_file_freq(), obj%get_file_frequnit())
258 if (obj%has_file_new_file_freq()) then
259 obj%next_close = diag_time_inc(obj%start_time, obj%get_file_new_file_freq(), &
260 obj%get_file_new_file_freq_units())
262 obj%next_close = diag_time_inc(obj%start_time, VERY_LARGE_FILE_FREQ, DIAG_DAYS)
264 obj%is_file_open = .false.
266 if(obj%has_file_duration()) then
267 obj%no_more_data = diag_time_inc(obj%start_time, obj%get_file_duration(), &
268 obj%get_file_duration_units())
270 obj%no_more_data = diag_time_inc(obj%start_time, VERY_LARGE_FILE_FREQ, DIAG_DAYS)
273 obj%time_ops = .false.
274 obj%unlim_dimension_level = 0
275 obj%is_static = obj%get_file_freq() .eq. -1
280 fms_diag_files_object_init = .true.
282 fms_diag_files_object_init = .false.
283 ! mpp_error("fms_diag_files_object_init: The diag_table.yaml file has not been correctly parsed.",&
286 end function fms_diag_files_object_init
288 !< @brief Determine if the field corresponding to the field_id was registered to the file
289 !! @return .True. if the field was registed to the file
290 pure logical function is_field_registered(this, field_id)
291 class(fmsDiagFile_type), intent(in) :: this !< The file object
292 integer, intent(in) :: field_id !< Id of the field to check
294 is_field_registered = this%field_registered(field_id)
295 end function is_field_registered
297 !> \brief Adds a field and yaml ID to the file
298 subroutine add_field_and_yaml_id (this, new_field_id, yaml_id)
299 class(fmsDiagFile_type), intent(inout) :: this !< The file object
300 integer, intent(in) :: new_field_id !< The field ID to be added to field_ids
301 integer, intent(in) :: yaml_id !< The yaml_id
303 this%num_registered_fields = this%num_registered_fields + 1
304 if (this%num_registered_fields .le. size(this%field_ids)) then
305 this%field_ids( this%num_registered_fields ) = new_field_id
306 this%yaml_ids( this%num_registered_fields ) = yaml_id
307 this%field_registered( this%num_registered_fields ) = .true.
309 call mpp_error(FATAL, "The file: "//this%get_file_fname()//" has already been assigned its maximum "//&
312 end subroutine add_field_and_yaml_id
314 !> \brief Adds a buffer_id to the file object
315 subroutine add_buffer_id (this, buffer_id)
316 class(fmsDiagFile_type), intent(inout) :: this !< The file object
317 integer, intent(in) :: buffer_id !< Buffer id to add to the file
319 this%number_of_buffers = this%number_of_buffers + 1
320 this%buffer_ids(this%number_of_buffers) = buffer_id
322 end subroutine add_buffer_id
324 !> \brief Initializes a diurnal axis for a fileobj
325 !! \note This is going to be called for every variable in the file, if the variable is not a diurnal variable
326 !! it will do nothing. It only defined a diurnal axis once.
327 subroutine init_diurnal_axis(this, diag_axis, naxis, yaml_id)
328 class(fmsDiagFile_type), intent(inout) :: this !< The file object
329 class(fmsDiagAxisContainer_type), intent(inout) :: diag_axis(:) !< Array of diag_axis object
330 integer, intent(inout) :: naxis !< Number of diag_axis that heve been defined
331 integer, intent(in) :: yaml_id !< The ID to the variable's yaml
333 integer :: i !< For do loops
334 type(diagYamlFilesVar_type), pointer :: field_yaml !< pointer to the yaml entry
336 field_yaml => diag_yaml%get_diag_field_from_id(yaml_id)
338 !< Go away if the file does not need a diurnal axis
339 if (.not. field_yaml%has_n_diurnal()) return
341 !< Check if the diurnal axis is already defined for this number of diurnal samples
342 do i = 1, this%number_of_axis
343 select type(axis=>diag_axis(this%axis_ids(i))%axis)
344 type is (fmsDiagDiurnalAxis_type)
345 if(field_yaml%get_n_diurnal() .eq. axis%get_diurnal_axis_samples()) return
349 !< If it is not already defined, define it
350 call define_diurnal_axis(diag_axis, naxis, field_yaml%get_n_diurnal(), .true.)
351 call define_diurnal_axis(diag_axis, naxis, field_yaml%get_n_diurnal(), .False.)
353 !< Add it to the list of axis for the file
354 this%number_of_axis = this%number_of_axis + 1
355 this%axis_ids(this%number_of_axis) = naxis !< This is the diurnal axis edges
357 this%number_of_axis = this%number_of_axis + 1
358 this%axis_ids(this%number_of_axis) = naxis - 1 !< This the diurnal axis
360 end subroutine init_diurnal_axis
362 !> \brief Set the time_ops variable in the diag_file object
363 subroutine set_file_time_ops(this, VarYaml, is_static)
364 class(fmsDiagFile_type), intent(inout) :: this !< The file object
365 type (diagYamlFilesVar_type), intent(in) :: VarYaml !< The variable's yaml file
366 logical, intent(in) :: is_static !< Flag indicating if variable is static
368 !< Go away if the file is static
369 if (this%is_static) return
371 if (this%time_ops) then
372 if (is_static) return
373 if (VarYaml%get_var_reduction() .eq. time_none) then
374 call mpp_error(FATAL, "The file: "//this%get_file_fname()//&
375 " has variables that are time averaged and instantaneous")
378 select case (VarYaml%get_var_reduction())
379 case (time_average, time_rms, time_max, time_min, time_sum, time_diurnal, time_power)
380 this%time_ops = .true.
384 end subroutine set_file_time_ops
386 !> \brief Logical function to determine if the variable file_metadata_from_model has been allocated or associated
387 !! \return .True. if file_metadata_from_model exists .False. if file_metadata_from_model has not been set
388 pure logical function has_file_metadata_from_model (this)
389 class(fmsDiagFile_type), intent(in) :: this !< The file object
390 has_file_metadata_from_model = allocated(this%file_metadata_from_model)
391 end function has_file_metadata_from_model
393 !> \brief Logical function to determine if the variable fileobj has been allocated or associated
394 !! \return .True. if fileobj exists .False. if fileobj has not been set
395 pure logical function has_fileobj (this)
396 class(fmsDiagFile_type), intent(in) :: this !< The file object
397 has_fileobj = allocated(this%fms2io_fileobj)
398 end function has_fileobj
400 !> \brief Logical function to determine if the variable diag_yaml_file has been allocated or associated
401 !! \return .True. if diag_yaml_file exists .False. if diag_yaml has not been set
402 pure logical function has_diag_yaml_file (this)
403 class(fmsDiagFile_type), intent(in) :: this !< The file object
404 has_diag_yaml_file = associated(this%diag_yaml_file)
405 end function has_diag_yaml_file
407 !> \brief Get the time to use to determine the filename, if using a wildcard file name (i.e ocn%4yr%2mo%2dy%2hr)
408 !! \return The time to use when determining the filename
409 function get_filename_time(this) &
411 class(fmsDiagFile_type), intent(in) :: this !< The file object
412 type(time_type) :: res
414 select case (this%diag_yaml_file%get_filename_time())
416 res = this%last_output
418 res = (this%last_output + this%next_close)/2
420 res = this%next_close
422 end function get_filename_time
424 !> \brief Logical function to determine if the variable field_ids has been allocated or associated
425 !! \return .True. if field_ids exists .False. if field_ids has not been set
426 pure logical function has_field_ids (this)
427 class(fmsDiagFile_type), intent(in) :: this !< The file object
428 has_field_ids = allocated(this%field_ids)
429 end function has_field_ids
431 !> \brief Returns a copy of the value of id
432 !! \return A copy of id
433 pure function get_id (this) result (res)
434 class(fmsDiagFile_type), intent(in) :: this !< The file object
439 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
441 !> \brief Returns a copy of the value of fileobj
442 !! \return A copy of fileobj
443 !pure function get_fileobj (obj) result (res)
444 ! class(fmsDiagFile_type), intent(in) :: obj !< The file object
445 ! class(FmsNetcdfFile_t) :: res
447 !end function get_fileobj
448 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
450 !!> \brief Returns a copy of the value of diag_yaml_file
451 !!! \return A copy of diag_yaml_file
452 !pure function get_diag_yaml_file (obj) result (res)
453 ! class(fmsDiagFile_type), intent(in) :: obj !< The file object
454 ! type(diagYamlFiles_type) :: res
455 ! res = obj%diag_yaml_file
456 !end function get_diag_yaml_file
458 !> \brief Returns a copy of the value of file_metadata_from_model
459 !! \return A copy of file_metadata_from_model
460 pure function get_file_metadata_from_model (this) result (res)
461 class(fmsDiagFile_type), intent(in) :: this !< The file object
462 character(len=:), dimension(:), allocatable :: res
463 res = this%file_metadata_from_model
464 end function get_file_metadata_from_model
466 !> \brief Returns a copy of the value of field_ids
467 !! \return A copy of field_ids
468 pure function get_field_ids (this) result (res)
469 class(fmsDiagFile_type), intent(in) :: this !< The file object
470 integer, dimension(:), allocatable :: res
471 allocate(res(size(this%field_ids)))
473 end function get_field_ids
475 !!!!!!!!! Functions from diag_yaml_file
476 !> \brief Returns a copy of file_fname from the yaml object
477 !! \return Copy of file_fname
478 pure function get_file_fname (this) result(res)
479 class(fmsDiagFile_type), intent(in) :: this !< The file object
480 character (len=:), allocatable :: res
481 res = this%diag_yaml_file%get_file_fname()
482 end function get_file_fname
484 !> \brief Returns a copy of file_frequnit from the yaml object
485 !! \return Copy of file_frequnit
486 pure function get_file_frequnit (this) result(res)
487 class(fmsDiagFile_type), intent(in) :: this !< The file object
489 res = this%diag_yaml_file%get_file_frequnit()
490 end function get_file_frequnit
492 !> \brief Returns a copy of file_freq from the yaml object
493 !! \return Copy of file_freq
494 pure function get_file_freq (this) result(res)
495 class(fmsDiagFile_type), intent(in) :: this !< The file object
497 res = this%diag_yaml_file%get_file_freq()
498 end function get_file_freq
500 !> \brief Returns a copy of file_timeunit from the yaml object
501 !! \return Copy of file_timeunit
502 pure function get_file_timeunit (this) result(res)
503 class(fmsDiagFile_type), intent(in) :: this !< The file object
505 res = this%diag_yaml_file%get_file_timeunit()
506 end function get_file_timeunit
508 !> \brief Returns a copy of file_unlimdim from the yaml object
509 !! \return Copy of file_unlimdim
510 pure function get_file_unlimdim (this) result(res)
511 class(fmsDiagFile_type), intent(in) :: this !< The file object
512 character (len=:), allocatable :: res
513 res = this%diag_yaml_file%get_file_unlimdim()
514 end function get_file_unlimdim
516 !> \brief Returns a copy of file_sub_region from the yaml object
517 !! \return Copy of file_sub_region
518 function get_file_sub_region (obj) result(res)
519 class(fmsDiagFile_type), intent(in) :: obj !< The file object
520 type(subRegion_type) :: res
521 res = obj%diag_yaml_file%get_file_sub_region()
522 end function get_file_sub_region
524 !< @brief Query for the subregion grid type (latlon or index)
525 !! @return subregion grid type
526 function get_file_sub_region_grid_type(this) &
528 class(fmsDiagFile_type), intent(in) :: this !< Diag file object
531 type(subRegion_type) :: subregion !< Subregion type
533 subregion = this%diag_yaml_file%get_file_sub_region()
534 res = subregion%grid_type
535 end function get_file_sub_region_grid_type
537 !> \brief Returns a copy of file_new_file_freq from the yaml object
538 !! \return Copy of file_new_file_freq
539 pure function get_file_new_file_freq (this) result(res)
540 class(fmsDiagFile_type), intent(in) :: this !< The file object
542 res = this%diag_yaml_file%get_file_new_file_freq()
543 end function get_file_new_file_freq
545 !> \brief Returns a copy of file_new_file_freq_units from the yaml object
546 !! \return Copy of file_new_file_freq_units
547 pure function get_file_new_file_freq_units (this) result(res)
548 class(fmsDiagFile_type), intent(in) :: this !< The file object
550 res = this%diag_yaml_file%get_file_new_file_freq_units()
551 end function get_file_new_file_freq_units
553 !> \brief Returns a copy of file_start_time from the yaml object
554 !! \return Copy of file_start_time
555 pure function get_file_start_time (this) result(res)
556 class(fmsDiagFile_type), intent(in) :: this !< The file object
557 character (len=:), allocatable :: res
558 res = this%diag_yaml_file%get_file_start_time()
559 end function get_file_start_time
561 !> \brief Returns a copy of file_duration from the yaml object
562 !! \return Copy of file_duration
563 pure function get_file_duration (this) result(res)
564 class(fmsDiagFile_type), intent(in) :: this !< The file object
566 res = this%diag_yaml_file%get_file_duration()
567 end function get_file_duration
569 !> \brief Returns a copy of file_duration_units from the yaml object
570 !! \return Copy of file_duration_units
571 pure function get_file_duration_units (this) result(res)
572 class(fmsDiagFile_type), intent(in) :: this !< The file object
574 res = this%diag_yaml_file%get_file_duration_units()
575 end function get_file_duration_units
577 !> \brief Returns a copy of file_varlist from the yaml object
578 !! \return Copy of file_varlist
579 pure function get_file_varlist (this) result(res)
580 class(fmsDiagFile_type), intent(in) :: this !< The file object
581 character (len=:), allocatable, dimension(:) :: res
582 res = this%diag_yaml_file%get_file_varlist()
583 end function get_file_varlist
585 !> \brief Returns a copy of file_global_meta from the yaml object
586 !! \return Copy of file_global_meta
587 pure function get_file_global_meta (this) result(res)
588 class(fmsDiagFile_type), intent(in) :: this !< The file object
589 character (len=MAX_STR_LEN), allocatable, dimension(:,:) :: res
590 res = this%diag_yaml_file%get_file_global_meta()
591 end function get_file_global_meta
593 !> \brief Determines if done writing data
594 !! \return .True. if done writing data
595 pure function is_done_writing_data (this) result(res)
596 class(fmsDiagFile_type), intent(in) :: this !< The file object
598 res = this%done_writing_data
599 end function is_done_writing_data
601 !> \brief Checks if file_fname is allocated in the yaml object
602 !! \return true if file_fname is allocated
603 pure function has_file_fname (this) result(res)
604 class(fmsDiagFile_type), intent(in) :: this !< The file object
606 res = this%diag_yaml_file%has_file_fname()
607 end function has_file_fname
609 !> \brief Checks if file_frequnit is allocated in the yaml object
610 !! \return true if file_frequnit is allocated
611 pure function has_file_frequnit (this) result(res)
612 class(fmsDiagFile_type), intent(in) :: this !< The file object
614 res = this%diag_yaml_file%has_file_frequnit()
615 end function has_file_frequnit
617 !> \brief Checks if file_freq is allocated in the yaml object
618 !! \return true if file_freq is allocated
619 pure function has_file_freq (this) result(res)
620 class(fmsDiagFile_type), intent(in) :: this !< The file object
622 res = this%diag_yaml_file%has_file_freq()
623 end function has_file_freq
625 !> \brief Checks if file_timeunit is allocated in the yaml object
626 !! \return true if file_timeunit is allocated
627 pure function has_file_timeunit (this) result(res)
628 class(fmsDiagFile_type), intent(in) :: this !< The file object
630 res = this%diag_yaml_file%has_file_timeunit()
631 end function has_file_timeunit
633 !> \brief Checks if file_unlimdim is allocated in the yaml object
634 !! \return true if file_unlimdim is allocated
635 pure function has_file_unlimdim (this) result(res)
636 class(fmsDiagFile_type), intent(in) :: this !< The file object
638 res = this%diag_yaml_file%has_file_unlimdim()
639 end function has_file_unlimdim
641 !> \brief Checks if file_sub_region is allocated in the yaml object
642 !! \return true if file_sub_region is allocated
643 pure function has_file_sub_region (this) result(res)
644 class(fmsDiagFile_type), intent(in) :: this !< The file object
646 res = this%diag_yaml_file%has_file_sub_region()
647 end function has_file_sub_region
649 !> \brief Checks if file_new_file_freq is allocated in the yaml object
650 !! \return true if file_new_file_freq is allocated
651 pure function has_file_new_file_freq (this) result(res)
652 class(fmsDiagFile_type), intent(in) :: this !< The file object
654 res = this%diag_yaml_file%has_file_new_file_freq()
655 end function has_file_new_file_freq
657 !> \brief Checks if file_new_file_freq_units is allocated in the yaml object
658 !! \return true if file_new_file_freq_units is allocated
659 pure function has_file_new_file_freq_units (this) result(res)
660 class(fmsDiagFile_type), intent(in) :: this !< The file object
662 res = this%diag_yaml_file%has_file_new_file_freq_units()
663 end function has_file_new_file_freq_units
665 !> \brief Checks if file_start_time is allocated in the yaml object
666 !! \return true if file_start_time is allocated
667 pure function has_file_start_time (this) result(res)
668 class(fmsDiagFile_type), intent(in) :: this !< The file object
670 res = this%diag_yaml_file%has_file_start_time()
671 end function has_file_start_time
673 !> \brief Checks if file_duration is allocated in the yaml object
674 !! \return true if file_duration is allocated
675 pure function has_file_duration (this) result(res)
676 class(fmsDiagFile_type), intent(in) :: this !< The file object
678 res = this%diag_yaml_file%has_file_duration()
679 end function has_file_duration
681 !> \brief Checks if file_duration_units is allocated in the yaml object
682 !! \return true if file_duration_units is allocated
683 pure function has_file_duration_units (this) result(res)
684 class(fmsDiagFile_type), intent(in) :: this !< The file object
686 res = this%diag_yaml_file%has_file_duration_units()
687 end function has_file_duration_units
689 !> \brief Checks if file_varlist is allocated in the yaml object
690 !! \return true if file_varlist is allocated
691 pure function has_file_varlist (this) result(res)
692 class(fmsDiagFile_type), intent(in) :: this !< The file object
694 res = this%diag_yaml_file%has_file_varlist()
695 end function has_file_varlist
697 !> \brief Checks if file_global_meta is allocated in the yaml object
698 !! \return true if file_global_meta is allocated
699 pure function has_file_global_meta (this) result(res)
700 class(fmsDiagFile_type), intent(in) :: this !< The file object
702 res = this%diag_yaml_file%has_file_global_meta()
703 end function has_file_global_meta
705 !> @brief Sets the domain and type of domain from the axis IDs
706 subroutine set_domain_from_axis(this, diag_axis, axes)
707 class(fmsDiagFile_type), intent(inout) :: this !< The file object
708 class(fmsDiagAxisContainer_type), intent(in) :: diag_axis(:) !< Array of diag_axis
709 integer, intent(in) :: axes (:)
711 call get_domain_and_domain_type(diag_axis, axes, this%type_of_domain, this%domain, this%get_file_fname())
712 end subroutine set_domain_from_axis
714 !> @brief Set the domain and the type_of_domain for a file
715 !> @details This subroutine is going to be called once by every variable in the file
716 !! in register_diag_field. It will update the domain and the type_of_domain if needed and verify that
717 !! all the variables are in the same domain
718 subroutine set_file_domain(this, domain, type_of_domain)
719 class(fmsDiagFile_type), intent(inout) :: this !< The file object
720 integer, INTENT(in) :: type_of_domain !< fileobj_type to use
721 CLASS(diagDomain_t), INTENT(in), target :: domain !< Domain
723 if (type_of_domain .ne. this%type_of_domain) then
724 !! If the current type_of_domain in the file obj is not the same as the variable calling this subroutine
726 if (type_of_domain .eq. NO_DOMAIN .or. this%type_of_domain .eq. NO_DOMAIN) then
727 !! If they are not the same then one of them can be NO_DOMAIN
728 !! (i.e a file can have variables that are not domain decomposed and variables that are)
730 if (type_of_domain .ne. NO_DOMAIN) then
731 !! Update the file's type_of_domain and domain if needed
732 this%type_of_domain = type_of_domain
733 this%domain => domain
737 !! If they are not the same and of them is not NO_DOMAIN, then crash because the variables don't have the
738 !! same domain (i.e a file has a variable is that in a 2D domain and one that is in a UG domain)
740 call mpp_error(FATAL, "The file: "//this%get_file_fname()//" has variables that are not in the same domain")
744 end subroutine set_file_domain
746 !> @brief Loops through a variable's axis_ids and adds them to the FMSDiagFile object if they don't exist
747 subroutine add_axes(this, axis_ids, diag_axis, naxis, yaml_id, buffer_id, output_buffers)
748 class(fmsDiagFile_type), intent(inout) :: this !< The file object
749 integer, INTENT(in) :: axis_ids(:) !< Array of axes_ids
750 class(fmsDiagAxisContainer_type), intent(inout) :: diag_axis(:) !< Diag_axis object
751 integer, intent(inout) :: naxis !< Number of axis that have been
753 integer, intent(in) :: yaml_id !< Yaml id of the field section for
755 integer, intent(in) :: buffer_id !< ID of the buffer
756 type(fmsDiagOutputBuffer_type), intent(inout) :: output_buffers(:) !< Array of output buffers
758 type(diagYamlFilesVar_type), pointer :: field_yaml !< pointer to the yaml entry
760 integer :: i, j !< For do loops
761 logical :: is_cube_sphere !< Flag indicating if the file's domain is a cubesphere
762 logical :: axis_found !< Flag indicating that the axis was already to the file obj
763 integer, allocatable :: var_axis_ids(:) !< Array of the variable's axis ids
764 integer :: x_y_axis_id(2) !< Ids of the x and y axis
765 integer :: x_or_y !< integer indicating if the axis is x or y
766 logical :: is_x_or_y !< flag indicating if the axis is x or y
767 integer :: subregion_gridtype !< The type of the subregion (latlon or index)
768 logical :: write_on_this_pe !< Flag indicating if the current pe is in the subregion
770 is_cube_sphere = .false.
771 subregion_gridtype = this%get_file_sub_region_grid_type()
773 field_yaml => diag_yaml%get_diag_field_from_id(yaml_id)
775 !< Created a copy here, because if the variable has a z subaxis var_axis_ids will be modified in
776 !! `create_new_z_subaxis` to contain the id of the new z subaxis instead of the parent axis,
777 !! which will be added to the the list of axis in the file object (axis_ids is intent(in),
778 !! which is why the copy was needed)
779 var_axis_ids = axis_ids
781 if (field_yaml%has_var_zbounds()) then
782 call create_new_z_subaxis(field_yaml%get_var_zbounds(), var_axis_ids, diag_axis, naxis, &
783 this%axis_ids, this%number_of_axis, this%nz_subaxis)
787 type is (subRegionalFile_type)
788 if (associated(this%domain)) then
789 if (this%domain%get_ntiles() .eq. 6) is_cube_sphere = .true.
791 if (.not. this%get_write_on_this_pe()) return
792 subaxis_defined: if (this%is_subaxis_defined) then
793 do i = 1, size(var_axis_ids)
794 select type (parent_axis => diag_axis(var_axis_ids(i))%axis)
795 type is (fmsDiagFullAxis_type)
797 is_x_or_y = parent_axis%is_x_or_y_axis()
798 do j = 1, this%number_of_axis
800 if(is_parent_axis(this%axis_ids(j), var_axis_ids(i), diag_axis)) then
802 var_axis_ids(i) = this%axis_ids(j) !Set the var_axis_id to the sub axis_id
805 elseif (var_axis_ids(i) .eq. this%axis_ids(j)) then
810 if (.not. axis_found) then
812 if (subregion_gridtype .eq. latlon_gridtype .and. is_cube_sphere) &
813 call mpp_error(FATAL, "If using the cube sphere and defining the subregion with latlon "//&
814 "the variable need to have the same x and y axis. Please check the variables in the file "//&
815 trim(this%get_file_fname())//" or use indices to define the subregion.")
817 select case (subregion_gridtype)
818 case (index_gridtype)
819 call define_new_subaxis_index(parent_axis, this%get_file_sub_region(), diag_axis, naxis, &
821 case (latlon_gridtype)
822 call define_new_subaxis_latlon(diag_axis, var_axis_ids(i:i), naxis, this%get_file_sub_region(), &
823 .false., write_on_this_pe)
825 call this%update_write_on_this_pe(write_on_this_pe)
826 if (.not. this%get_write_on_this_pe()) cycle
827 call this%add_new_axis(naxis)
828 var_axis_ids(i) = naxis
830 call this%add_new_axis(var_axis_ids(i))
833 type is (fmsDiagSubAxis_type)
834 axis_found = this%does_axis_exist(var_axis_ids(i))
835 if (.not. axis_found) call this%add_new_axis(var_axis_ids(i))
839 x_y_axis_id = diag_null
840 do i = 1, size(var_axis_ids)
841 select type (parent_axis => diag_axis(var_axis_ids(i))%axis)
842 type is (fmsDiagFullAxis_type)
843 if (.not. parent_axis%is_x_or_y_axis(x_or_y)) then
844 axis_found = this%does_axis_exist(var_axis_ids(i))
845 if (.not. axis_found) call this%add_new_axis(var_axis_ids(i))
847 x_y_axis_id(x_or_y) = var_axis_ids(i)
849 type is (fmsDiagSubAxis_type)
850 axis_found = this%does_axis_exist(var_axis_ids(i))
851 if (.not. axis_found) call this%add_new_axis(var_axis_ids(i))
855 call this%define_new_subaxis(var_axis_ids, x_y_axis_id, is_cube_sphere, diag_axis, naxis)
856 this%is_subaxis_defined = .true.
857 endif subaxis_defined
858 type is (fmsDiagFile_type)
859 do i = 1, size(var_axis_ids)
860 axis_found = this%does_axis_exist(var_axis_ids(i))
861 if (.not. axis_found) call this%add_new_axis(var_axis_ids(i))
864 !> Add the axis to the buffer object
865 call output_buffers(buffer_id)%add_axis_ids(var_axis_ids)
866 end subroutine add_axes
868 !> @brief Adds a new axis the list of axis in the diag file object
869 subroutine add_new_axis(this, var_axis_id)
870 class(fmsDiagFile_type), intent(inout) :: this !< The file object
871 integer, intent(in) :: var_axis_id !< Axis id of the variable
873 this%number_of_axis = this%number_of_axis + 1
874 this%axis_ids(this%number_of_axis) = var_axis_id
875 end subroutine add_new_axis
877 !> @brief This updates write on this pe
878 subroutine update_write_on_this_pe(this, write_on_this_pe)
879 class(fmsDiagFile_type), intent(inout) :: this !< The file object
880 logical, intent(in) :: write_on_this_pe !< .True. if the current PE is in
884 type is (subRegionalFile_type)
885 if (this%write_on_this_pe) this%write_on_this_pe = write_on_this_pe
887 end subroutine update_write_on_this_pe
889 !> @brief Query for the write_on_this_pe member of the diag file object
890 !! @return the write_on_this_pe member of the diag file object
891 function get_write_on_this_pe(this) &
893 class(fmsDiagFile_type), intent(inout) :: this !< The file object
897 type is (subRegionalFile_type)
898 rslt= this%write_on_this_pe
900 end function get_write_on_this_pe
902 !< @brief Determine if an axis is already in the list of axis for a diag file
903 !! @return .True. if the axis is already in the list of axis for a diag file
904 function does_axis_exist(this, var_axis_id) &
906 class(fmsDiagFile_type), intent(inout) :: this !< The file object
907 integer, intent(in) :: var_axis_id !< Variable axis id to check
910 integer :: j !< For do loops
913 do j = 1, this%number_of_axis
914 !> Check if the axis already exists, move on
915 if (var_axis_id .eq. this%axis_ids(j)) then
922 !> @brief Define a new sub axis
923 subroutine define_new_subaxis(this, var_axis_ids, x_y_axis_id, is_cube_sphere, diag_axis, naxis)
924 class(fmsDiagFile_type), intent(inout) :: this !< The file object
925 integer, INTENT(inout) :: var_axis_ids(:) !< Original variable axis ids
926 integer, INTENT(in) :: x_y_axis_id(:) !< The ids of the x and y axis
927 logical, intent(in) :: is_cube_sphere !< .True. if the axis is in the cubesphere
928 integer, intent(inout) :: naxis !< Number of axis current registered
929 class(fmsDiagAxisContainer_type), intent(inout) :: diag_axis(:) !< Diag_axis object
931 logical :: write_on_this_pe !< .True. if the current PE is in the subregion
932 integer :: i, j !< For do loop
934 select case (this%get_file_sub_region_grid_type())
935 case(latlon_gridtype)
936 call define_new_subaxis_latlon(diag_axis, x_y_axis_id, naxis, this%get_file_sub_region(), is_cube_sphere, &
938 call this%update_write_on_this_pe(write_on_this_pe)
939 if (.not. this%get_write_on_this_pe()) return
940 call this%add_new_axis(naxis)
941 call this%add_new_axis(naxis-1)
942 do j = 1, size(var_axis_ids)
943 if (x_y_axis_id(1) .eq. var_axis_ids(j)) var_axis_ids(j) = naxis - 1
944 if (x_y_axis_id(2) .eq. var_axis_ids(j)) var_axis_ids(j) = naxis
946 case (index_gridtype)
947 do i = 1, size(x_y_axis_id)
948 select type (parent_axis => diag_axis(x_y_axis_id(i))%axis)
949 type is (fmsDiagFullAxis_type)
950 call define_new_subaxis_index(parent_axis, this%get_file_sub_region(), diag_axis, naxis, i, &
952 call this%update_write_on_this_pe(write_on_this_pe)
953 if (.not. this%get_write_on_this_pe()) return
954 call this%add_new_axis(naxis)
955 do j = 1, size(var_axis_ids)
956 if (x_y_axis_id(i) .eq. var_axis_ids(j)) var_axis_ids(j) = naxis
961 end subroutine define_new_subaxis
963 !> @brief adds the start time to the fileobj
964 !! @note This should be called from the register field calls. It can be called multiple times (one for each variable)
965 !! So it needs to make sure that the start_time is the same for each variable. The initial value is the base_time
966 subroutine add_start_time(this, start_time, model_time)
967 class(fmsDiagFile_type), intent(inout) :: this !< The file object
968 TYPE(time_type), intent(in) :: start_time !< Start time to add to the fileobj
969 TYPE(time_type), intent(out) :: model_time !< The current model time
970 !! this will be set to the start_time
971 !! at the begining of the run
973 !< If the start_time sent in is equal to the base_time return because
974 !! this%start_time was already set to the base_time
975 if (start_time .eq. get_base_time()) return
977 if (this%start_time .ne. get_base_time()) then
978 !> If the this%start_time is not equal to the base_time from the diag_table
979 !! this%start_time was already updated so make sure it is the same or error out
980 if (this%start_time .ne. start_time)&
981 call mpp_error(FATAL, "The variables associated with the file:"//this%get_file_fname()//" have"&
982 &" different start_time")
984 !> If the this%start_time is equal to the base_time,
985 !! simply update it with the start_time and set up the *_output variables
986 model_time = start_time
987 this%start_time = start_time
988 this%last_output = start_time
989 this%next_output = diag_time_inc(start_time, this%get_file_freq(), this%get_file_frequnit())
990 this%next_next_output = diag_time_inc(this%next_output, this%get_file_freq(), this%get_file_frequnit())
991 if (this%has_file_new_file_freq()) then
992 this%next_close = diag_time_inc(this%start_time, this%get_file_new_file_freq(), &
993 this%get_file_new_file_freq_units())
995 if (this%is_static) then
996 ! If the file is static, set the close time to be equal to the start_time, so that it can be closed
997 ! after the first write!
998 this%next_close = this%start_time
999 this%next_next_output = diag_time_inc(this%start_time, VERY_LARGE_FILE_FREQ, DIAG_DAYS)
1001 this%next_close = diag_time_inc(this%start_time, VERY_LARGE_FILE_FREQ, DIAG_DAYS)
1005 if(this%has_file_duration()) then
1006 this%no_more_data = diag_time_inc(this%start_time, this%get_file_duration(), &
1007 this%get_file_duration_units())
1009 this%no_more_data = diag_time_inc(this%start_time, VERY_LARGE_FILE_FREQ, DIAG_DAYS)
1016 !> writes out internal values for fmsDiagFile_type object
1017 subroutine dump_file_obj(this, unit_num)
1018 class(fmsDiagFile_type), intent(in) :: this !< the file object
1019 integer, intent(in) :: unit_num !< passed in from dump_diag_obj
1020 !! will either be for new log file or stdout
1021 write( unit_num, *) 'file id:', this%id
1022 write( unit_num, *) 'start time:', date_to_string(this%start_time)
1023 write( unit_num, *) 'last_output', date_to_string(this%last_output)
1024 write( unit_num, *) 'next_output', date_to_string(this%next_output)
1025 write( unit_num, *)'next_next_output', date_to_string(this%next_next_output)
1026 write( unit_num, *)'next_close', date_to_string(this%next_close)
1028 if( allocated(this%fms2io_fileobj)) write( unit_num, *)'fileobj path', this%fms2io_fileobj%path
1030 write( unit_num, *)'type_of_domain', this%type_of_domain
1031 if( allocated(this%file_metadata_from_model)) write( unit_num, *) 'file_metadata_from_model', &
1032 this%file_metadata_from_model
1033 if( allocated(this%field_ids)) write( unit_num, *)'field_ids', this%field_ids
1034 if( allocated(this%field_registered)) write( unit_num, *)'field_registered', this%field_registered
1035 if( allocated(this%num_registered_fields)) write( unit_num, *)'num_registered_fields', this%num_registered_fields
1036 if( allocated(this%axis_ids)) write( unit_num, *)'axis_ids', this%axis_ids(1:this%number_of_axis)
1040 !> @brief Determine if a file is regional
1041 !! @return Flag indicating if the file is regional or not
1042 logical pure function is_regional(this)
1043 class(fmsDiagFileContainer_type), intent(in) :: this !< The file object
1045 select type (wut=>this%FMS_diag_file)
1046 type is (subRegionalFile_type)
1047 is_regional = .true.
1048 type is (fmsDiagFile_type)
1049 is_regional = .false.
1052 end function is_regional
1054 !> @brief Determine if a file is static
1055 !! @return Flag indicating if the file is static or not
1056 logical pure function is_file_static(this)
1057 class(fmsDiagFileContainer_type), intent(in) :: this !< The file object
1059 is_file_static = .false.
1061 select type (fileptr=>this%FMS_diag_file)
1062 type is (fmsDiagFile_type)
1063 is_file_static = fileptr%is_static
1066 end function is_file_static
1068 !< @brief Opens the diag_file if it is time to do so
1069 subroutine open_diag_file(this, time_step, file_is_opened)
1070 class(fmsDiagFileContainer_type), intent(inout), target :: this !< The file object
1071 TYPE(time_type), intent(in) :: time_step !< Current model step time
1072 logical, intent(out) :: file_is_opened !< .true. if the file was opened in this
1075 class(fmsDiagFile_type), pointer :: diag_file !< Diag_file object to open
1076 class(diagDomain_t), pointer :: domain !< The domain used in the file
1077 character(len=:), allocatable :: diag_file_name !< The file name as defined in the yaml
1078 character(len=128) :: base_name !< The file name as defined in the yaml
1079 !! without the wildcard definition
1080 character(len=128) :: file_name !< The file name as it will be written to disk
1081 character(len=128) :: temp_name !< Temp variable to store the file_name
1082 character(len=128) :: start_date !< The start_time as a string that will be added to
1083 !! the begining of the filename (start_date.filename)
1084 character(len=128) :: suffix !< The current time as a string that will be added to
1085 !! the end of filename
1086 integer :: pos !< Index of the filename with the first "%" in the file name
1087 INTEGER :: year !< The year of the start_date
1088 INTEGER :: month !< The month of the start_date
1089 INTEGER :: day !< The day of the start_date
1090 INTEGER :: hour !< The hour of the start_date
1091 INTEGER :: minute !< The minute of the start_date
1092 INTEGER :: second !< The second of the start_date
1093 character(len=4) :: mype_string !< The pe as a string
1094 logical :: is_regional !< Flag indicating if the file is regional
1095 integer, allocatable :: pes(:) !< Array of the pes in the current pelist
1097 diag_file => this%FMS_diag_file
1098 domain => diag_file%domain
1100 file_is_opened = .false.
1101 !< Go away if it the file is already open
1102 if (diag_file%is_file_open) return
1104 is_regional = .false.
1105 !< Figure out what fms2io_fileobj to use!
1106 if (.not. allocated(diag_file%fms2io_fileobj)) then
1107 select type (diag_file)
1108 type is (subRegionalFile_type)
1109 !< In this case each PE is going to write its own file
1110 allocate(FmsNetcdfFile_t :: diag_file%fms2io_fileobj)
1111 is_regional = .true.
1112 type is (fmsDiagFile_type)
1113 !< Use the type_of_domain to get the correct fms2io_fileobj
1114 select case (diag_file%type_of_domain)
1116 allocate(FmsNetcdfFile_t :: diag_file%fms2io_fileobj)
1118 allocate(FmsNetcdfDomainFile_t :: diag_file%fms2io_fileobj)
1120 allocate(FmsNetcdfUnstructuredDomainFile_t :: diag_file%fms2io_fileobj)
1125 !< Figure out what to name of the file
1126 diag_file_name = diag_file%get_file_fname()
1128 !< If using the new_file_freq figure out what the name is based on the current time
1129 if (diag_file%has_file_new_file_freq()) then
1130 !< If using a wildcard file name (i.e ocn%4yr%2mo%2dy%2hr), get the basename (i.e ocn)
1131 pos = INDEX(diag_file_name, '%')
1132 if (pos > 0) base_name = diag_file_name(1:pos-1)
1133 suffix = get_time_string(diag_file_name, diag_file%get_filename_time())
1134 base_name = trim(base_name)//trim(suffix)
1136 base_name = trim(diag_file_name)
1139 !< Add the ens number to the file name (if it exists)
1140 file_name = trim(base_name)
1141 call get_instance_filename(base_name, file_name)
1143 !< Prepend the file start_time to the file name if prepend_date == .TRUE. in
1145 IF ( prepend_date ) THEN
1146 call get_date(diag_file%start_time, year, month, day, hour, minute, second)
1147 write (start_date, '(1I20.4, 2I2.2)') year, month, day
1149 file_name = TRIM(adjustl(start_date))//'.'//TRIM(file_name)
1152 file_name = trim(file_name)//".nc"
1154 !< If this is a regional file add the PE and the tile_number to the filename
1155 if (is_regional) then
1156 !< Get the pe number that will be appended to the end of the file
1157 write(mype_string,'(I0.4)') mpp_pe()
1159 !< Add the tile number if appropriate
1160 select type (domain)
1161 type is (DIAGDOMAIN2D_T)
1162 temp_name = file_name
1163 call get_mosaic_tile_file(temp_name, file_name, .true., domain%domain2)
1166 file_name = trim(file_name)//"."//trim(mype_string)
1170 select type (fms2io_fileobj => diag_file%fms2io_fileobj)
1171 type is (FmsNetcdfFile_t)
1172 if (is_regional) then
1173 if (.not. open_file(fms2io_fileobj, file_name, "overwrite", pelist=(/mpp_pe()/))) &
1174 &call mpp_error(FATAL, "Error opening the file:"//file_name)
1175 call register_global_attribute(fms2io_fileobj, "is_subregional", "True", str_len=4)
1177 allocate(pes(mpp_npes()))
1178 call mpp_get_current_pelist(pes)
1180 if (.not. open_file(fms2io_fileobj, file_name, "overwrite", pelist=pes)) &
1181 &call mpp_error(FATAL, "Error opening the file:"//file_name)
1183 type is (FmsNetcdfDomainFile_t)
1184 select type (domain)
1185 type is (diagDomain2d_t)
1186 if (.not. open_file(fms2io_fileobj, file_name, "overwrite", domain%Domain2)) &
1187 &call mpp_error(FATAL, "Error opening the file:"//file_name)
1189 type is (FmsNetcdfUnstructuredDomainFile_t)
1190 select type (domain)
1191 type is (diagDomainUg_t)
1192 if (.not. open_file(fms2io_fileobj, file_name, "overwrite", domain%DomainUG)) &
1193 &call mpp_error(FATAL, "Error opening the file:"//file_name)
1197 file_is_opened = .true.
1198 diag_file%is_file_open = file_is_opened
1201 end subroutine open_diag_file
1203 !< @brief Write global attributes in the diag_file
1204 subroutine write_global_metadata(this)
1205 class(fmsDiagFileContainer_type), intent(inout), target :: this !< The file object
1207 class(FmsNetcdfFile_t), pointer :: fms2io_fileobj !< The fileobj to write to
1208 integer :: i !< For do loops
1209 character (len=MAX_STR_LEN), allocatable :: yaml_file_attributes(:,:) !< Global attributes defined in the yaml
1211 type(diagYamlFiles_type), pointer :: diag_file_yaml !< The diag_file yaml
1213 diag_file_yaml => this%FMS_diag_file%diag_yaml_file
1214 fms2io_fileobj => this%FMS_diag_file%fms2io_fileobj
1216 if (diag_file_yaml%has_file_global_meta()) then
1217 yaml_file_attributes = diag_file_yaml%get_file_global_meta()
1218 do i = 1, size(yaml_file_attributes,1)
1219 call register_global_attribute(fms2io_fileobj, trim(yaml_file_attributes(i,1)), &
1220 trim(yaml_file_attributes(i,2)), str_len=len_trim(yaml_file_attributes(i,2)))
1222 deallocate(yaml_file_attributes)
1224 end subroutine write_global_metadata
1226 !< @brief Writes a variable's metadata in the netcdf file
1227 subroutine write_var_metadata(fms2io_fileobj, variable_name, dimensions, long_name, units)
1228 class(FmsNetcdfFile_t), intent(inout) :: fms2io_fileobj !< The file object to write into
1229 character(len=*) , intent(in) :: variable_name !< The name of the time variables
1230 character(len=*) , intent(in) :: dimensions(:) !< The dimensions of the variable
1231 character(len=*) , intent(in) :: long_name !< The long_name of the variable
1232 character(len=*) , intent(in) :: units !< The units of the variable
1234 call register_field(fms2io_fileobj, variable_name, pack_size_str, dimensions)
1235 call register_variable_attribute(fms2io_fileobj, variable_name, "long_name", &
1236 trim(long_name), str_len=len_trim(long_name))
1237 if (trim(units) .ne. no_units) &
1238 call register_variable_attribute(fms2io_fileobj, variable_name, "units", &
1239 trim(units), str_len=len_trim(units))
1240 end subroutine write_var_metadata
1242 !> \brief Write the time metadata to the diag file
1243 subroutine write_time_metadata(this)
1244 class(fmsDiagFileContainer_type), intent(inout), target :: this !< The file object
1246 class(fmsDiagFile_type), pointer :: diag_file !< Diag_file object to open
1247 class(FmsNetcdfFile_t), pointer :: fms2io_fileobj !< The fileobj to write to
1248 character(len=50) :: time_units_str !< Time units written as a string
1249 character(len=50) :: calendar !< The calendar name
1251 character(len=:), allocatable :: time_var_name !< The name of the time variable as it is defined in the yaml
1252 character(len=50) :: dimensions(2) !< Array of dimensions names for the variable
1254 diag_file => this%FMS_diag_file
1255 fms2io_fileobj => diag_file%fms2io_fileobj
1257 time_var_name = diag_file%get_file_unlimdim()
1258 call register_axis(fms2io_fileobj, time_var_name, unlimited)
1260 WRITE(time_units_str, 11) &
1261 TRIM(time_unit_list(diag_file%get_file_timeunit())), get_base_year(),&
1262 & get_base_month(), get_base_day(), get_base_hour(), get_base_minute(), get_base_second()
1263 11 FORMAT(a, ' since ', i4.4, '-', i2.2, '-', i2.2, ' ', i2.2, ':', i2.2, ':', i2.2)
1265 dimensions(1) = "nv"
1266 dimensions(2) = trim(time_var_name)
1268 call write_var_metadata(fms2io_fileobj, time_var_name, dimensions(2:2), &
1269 time_var_name, time_units_str)
1271 !< Add additional variables to the time variable
1272 call register_variable_attribute(fms2io_fileobj, time_var_name, "axis", "T", str_len=1 )
1274 !TODO no need to have both attributes, probably?
1275 calendar = valid_calendar_types(get_calendar_type())
1276 call register_variable_attribute(fms2io_fileobj, time_var_name, "calendar_type", &
1277 uppercase(trim(calendar)), str_len=len_trim(calendar))
1278 call register_variable_attribute(fms2io_fileobj, time_var_name, "calendar", &
1279 lowercase(trim(calendar)), str_len=len_trim(calendar))
1281 if (diag_file%time_ops) then
1282 call register_variable_attribute(fms2io_fileobj, time_var_name, "bounds", &
1283 trim(time_var_name)//"_bnds", str_len=len_trim(time_var_name//"_bnds"))
1285 !< It is possible that the "nv" "axis" was registered via "diag_axis_init" call
1286 !! so only adding it if it doesn't exist already
1287 if ( .not. dimension_exists(fms2io_fileobj, "nv")) then
1288 call register_axis(fms2io_fileobj, "nv", 2) !< Time bounds need a vertex number
1289 call write_var_metadata(fms2io_fileobj, "nv", dimensions(1:1), &
1290 "vertex number", no_units)
1292 call write_var_metadata(fms2io_fileobj, time_var_name//"_bnds", dimensions, &
1293 trim(time_var_name)//" axis boundaries", time_units_str)
1296 end subroutine write_time_metadata
1298 !> \brief Write out the field data to the file
1299 subroutine write_field_data(this, field_obj, buffer_obj)
1300 class(fmsDiagFileContainer_type), intent(in), target :: this !< The diag file object to write to
1301 type(fmsDiagField_type), intent(in), target :: field_obj !< The field object to write from
1302 type(fmsDiagOutputBuffer_type), intent(inout), target :: buffer_obj !< The buffer object with the data
1304 class(fmsDiagFile_type), pointer :: diag_file !< Diag_file object to open
1305 class(FmsNetcdfFile_t), pointer :: fms2io_fileobj !< Fileobj to write to
1306 logical :: has_diurnal !< indicates if theres a diurnal axis to adjust for
1308 diag_file => this%FMS_diag_file
1309 fms2io_fileobj => diag_file%fms2io_fileobj
1311 !TODO This may be offloaded in the future
1312 if (diag_file%is_static) then
1313 !< Here the file is static so there is no need for the unlimited dimension
1314 !! as a variables are static
1315 call buffer_obj%write_buffer(fms2io_fileobj)
1316 diag_file%data_has_been_written = .true.
1318 if (field_obj%is_static()) then
1319 !< If the variable is static, only write it the first time
1320 if (diag_file%unlim_dimension_level .eq. 1) then
1321 call buffer_obj%write_buffer(fms2io_fileobj)
1322 diag_file%data_has_been_written = .true.
1325 diag_file%data_has_been_written = .true.
1326 has_diurnal = buffer_obj%get_diurnal_sample_size() .gt. 1
1327 call buffer_obj%write_buffer(fms2io_fileobj, &
1328 unlim_dim_level=diag_file%unlim_dimension_level, is_diurnal=has_diurnal)
1332 end subroutine write_field_data
1334 !> \brief Determine if it is time to close the file
1335 !! \return .True. if it is time to close the file
1336 logical function is_time_to_close_file (this, time_step)
1337 class(fmsDiagFileContainer_type), intent(in), target :: this !< The file object
1338 TYPE(time_type), intent(in) :: time_step !< Current model step time
1340 if (time_step >= this%FMS_diag_file%next_close) then
1341 is_time_to_close_file = .true.
1343 is_time_to_close_file = .false.
1347 !> \brief Determine if it is time to "write" to the file
1348 logical function is_time_to_write(this, time_step, output_buffers)
1349 class(fmsDiagFileContainer_type), intent(inout), target :: this !< The file object
1350 TYPE(time_type), intent(in) :: time_step !< Current model step time
1351 type(fmsDiagOutputBuffer_type), intent(in) :: output_buffers(:) !< Array of output buffer.
1352 !! This is needed for error messages!
1354 if (time_step > this%FMS_diag_file%next_output) then
1355 is_time_to_write = .true.
1356 if (this%FMS_diag_file%is_static) return
1357 if (time_step > this%FMS_diag_file%next_next_output) then
1358 if (this%FMS_diag_file%num_registered_fields .eq. 0) then
1359 !! If no variables have been registered, write a dummy time dimension for the first level
1360 !! At least one time level is needed for the combiner to work ...
1361 if (this%FMS_diag_file%unlim_dimension_level .eq. 1) then
1362 call mpp_error(NOTE, this%FMS_diag_file%get_file_fname()//&
1363 ": diag_manager_mod: This file does not have any variables registered. Fill values will be written")
1364 this%FMS_diag_file%data_has_been_written = .true.
1366 is_time_to_write =.false.
1368 !! Only fail if send data has actually been called for at least one variable
1369 if (this%FMS_diag_file%has_send_data_been_called(output_buffers, .false.)) &
1370 call mpp_error(FATAL, this%FMS_diag_file%get_file_fname()//&
1371 ": diag_manager_mod: You skipped a time_step. Be sure that diag_send_complete is called at every "//&
1372 "time_step needed by the file.")
1373 is_time_to_write =.false.
1377 is_time_to_write = .false.
1378 if (this%FMS_diag_file%is_static) then
1379 ! This is to ensure that static files get finished in the begining of the run
1380 if (this%FMS_diag_file%unlim_dimension_level .eq. 1) is_time_to_write = .true.
1383 end function is_time_to_write
1385 !> \brief Determine if the current PE has data to write
1386 logical function writing_on_this_pe(this)
1387 class(fmsDiagFileContainer_type), intent(in), target :: this !< The file object
1389 select type(diag_file => this%FMS_diag_file)
1390 type is (subRegionalFile_type)
1391 writing_on_this_pe = diag_file%write_on_this_pe
1393 writing_on_this_pe = .true.
1398 !> \brief Write out the time data to the file
1399 subroutine write_time_data(this)
1400 class(fmsDiagFileContainer_type), intent(in), target :: this !< The file object
1402 real :: dif !< The time as a real number
1403 class(fmsDiagFile_type), pointer :: diag_file !< Diag_file object to open
1404 class(FmsNetcdfFile_t), pointer :: fms2io_fileobj !< The fileobj to write to
1405 TYPE(time_type) :: middle_time !< The middle time of the averaging period
1407 real :: T1 !< The beginning time of the averaging period
1408 real :: T2 !< The ending time of the averaging period
1410 diag_file => this%FMS_diag_file
1411 fms2io_fileobj => diag_file%fms2io_fileobj
1413 !< If data has not been written for the current unlimited dimension
1414 !! ignore this. The diag_file%unlim_dimension_level .ne. 1 is there to ensure
1415 !! that at least one time level is written (this is needed for the combiner)
1416 if (.not. diag_file%data_has_been_written .and. diag_file%unlim_dimension_level .ne. 1) return
1418 if (diag_file%time_ops) then
1419 middle_time = (diag_file%last_output+diag_file%next_output)/2
1420 dif = get_date_dif(middle_time, get_base_time(), diag_file%get_file_timeunit())
1422 dif = get_date_dif(diag_file%next_output, get_base_time(), diag_file%get_file_timeunit())
1425 call write_data(fms2io_fileobj, diag_file%get_file_unlimdim(), dif, &
1426 unlim_dim_level=diag_file%unlim_dimension_level)
1428 if (diag_file%time_ops) then
1429 T1 = get_date_dif(diag_file%last_output, get_base_time(), diag_file%get_file_timeunit())
1430 T2 = get_date_dif(diag_file%next_output, get_base_time(), diag_file%get_file_timeunit())
1432 call write_data(fms2io_fileobj, trim(diag_file%get_file_unlimdim())//"_bnds", &
1433 (/T1, T2/), unlim_dim_level=diag_file%unlim_dimension_level)
1435 if (diag_file%unlim_dimension_level .eq. 1) then
1436 call write_data(fms2io_fileobj, "nv", (/1, 2/))
1440 end subroutine write_time_data
1442 !> \brief Updates the current_new_file_freq_index if using a new_file_freq
1443 subroutine update_current_new_file_freq_index(this, time_step)
1444 class(fmsDiagFileContainer_type), intent(inout), target :: this !< The file object
1445 TYPE(time_type), intent(in) :: time_step !< Current model step time
1447 class(fmsDiagFile_type), pointer :: diag_file !< Diag_file object to open
1449 diag_file => this%FMS_diag_file
1451 if (time_step >= diag_file%no_more_data) then
1452 call diag_file%diag_yaml_file%increase_new_file_freq_index()
1454 if (diag_file%has_file_duration()) then
1455 diag_file%no_more_data = diag_time_inc(diag_file%no_more_data, diag_file%get_file_duration(), &
1456 diag_file%get_file_duration_units())
1458 !< At this point you are done writing data
1459 diag_file%done_writing_data = .true.
1460 diag_file%no_more_data = diag_time_inc(diag_file%no_more_data, VERY_LARGE_FILE_FREQ, DIAG_DAYS)
1461 diag_file%next_output = diag_file%no_more_data
1462 diag_file%next_next_output = diag_file%no_more_data
1463 diag_file%last_output = diag_file%no_more_data
1464 diag_file%next_close = diag_file%no_more_data
1467 end subroutine update_current_new_file_freq_index
1469 !> \brief Set up the next_output and next_next_output variable in a file obj
1470 subroutine update_next_write(this, time_step)
1471 class(fmsDiagFileContainer_type), intent(in), target :: this !< The file object
1472 TYPE(time_type), intent(in) :: time_step !< Current model step time
1474 class(fmsDiagFile_type), pointer :: diag_file !< Diag_file object to open
1476 diag_file => this%FMS_diag_file
1477 if (diag_file%is_static) then
1478 diag_file%last_output = diag_file%next_output
1479 diag_file%next_output = diag_time_inc(diag_file%next_output, VERY_LARGE_FILE_FREQ, DIAG_DAYS)
1480 diag_file%next_next_output = diag_time_inc(diag_file%next_output, VERY_LARGE_FILE_FREQ, DIAG_DAYS)
1482 diag_file%last_output = diag_file%next_output
1483 diag_file%next_output = diag_time_inc(diag_file%next_output, diag_file%get_file_freq(), &
1484 diag_file%get_file_frequnit())
1485 diag_file%next_next_output = diag_time_inc(diag_file%next_output, diag_file%get_file_freq(), &
1486 diag_file%get_file_frequnit())
1489 end subroutine update_next_write
1491 !> \brief Increase the unlimited dimension level that the file is currently being written to
1492 subroutine increase_unlim_dimension_level(this)
1493 class(fmsDiagFileContainer_type), intent(inout), target :: this !< The file object
1495 this%FMS_diag_file%unlim_dimension_level = this%FMS_diag_file%unlim_dimension_level + 1
1496 this%FMS_diag_file%data_has_been_written = .false.
1497 end subroutine increase_unlim_dimension_level
1499 !> \brief Get the unlimited dimension level that is in the file
1500 !! \return The unlimited dimension
1501 pure function get_unlim_dimension_level(this) &
1503 class(fmsDiagFileContainer_type), intent(in), target :: this !< The file object
1506 res = this%FMS_diag_file%unlim_dimension_level
1509 !> \brief Get the next_output for the file object
1510 !! \return The next_output
1511 pure function get_next_output(this) &
1513 class(fmsDiagFileContainer_type), intent(in), target :: this !< The file object
1514 type(time_type) :: res
1516 res = this%FMS_diag_file%next_output
1517 end function get_next_output
1519 !> \brief Get the next_output for the file object
1520 !! \return The next_output
1521 pure function get_next_next_output(this) &
1523 class(fmsDiagFileContainer_type), intent(in), target :: this !< The file object
1524 type(time_type) :: res
1526 res = this%FMS_diag_file%next_next_output
1527 if (this%FMS_diag_file%is_static) then
1528 res = this%FMS_diag_file%no_more_data
1530 end function get_next_next_output
1532 !< @brief Writes the axis metadata for the file
1533 subroutine write_axis_metadata(this, diag_axis)
1534 class(fmsDiagFileContainer_type), intent(inout), target :: this !< The file object
1535 class(fmsDiagAxisContainer_type), intent(in), target :: diag_axis(:) !< Diag_axis object
1537 class(fmsDiagFile_type), pointer :: diag_file !< Diag_file object to open
1538 class(FmsNetcdfFile_t), pointer :: fms2io_fileobj !< The fileobj to write to
1539 integer :: i,k !< For do loops
1540 integer :: parent_axis_id !< Id of the parent_axis
1541 integer :: structured_ids(2) !< Ids of the uncompress axis
1542 integer :: edges_id !< Id of the axis edge
1544 class(fmsDiagAxisContainer_type), pointer :: axis_ptr !< pointer to the axis object currently writing
1545 logical :: edges_in_file !< .true. if the edges are already in the file
1547 diag_file => this%FMS_diag_file
1548 fms2io_fileobj => diag_file%fms2io_fileobj
1550 do i = 1, diag_file%number_of_axis
1551 edges_in_file = .false.
1552 axis_ptr => diag_axis(diag_file%axis_ids(i))
1553 parent_axis_id = axis_ptr%axis%get_parent_axis_id()
1555 edges_id = axis_ptr%axis%get_edges_id()
1556 if (edges_id .ne. diag_null) then
1557 !< write the edges if is not in the list of axis in the file, otherwrise ignore
1558 if (any(diag_file%axis_ids(1:diag_file%number_of_axis) .eq. edges_id)) then
1559 edges_in_file = .true.
1561 call diag_axis(edges_id)%axis%write_axis_metadata(fms2io_fileobj, .true.)
1562 call diag_file%add_new_axis(edges_id)
1566 if (parent_axis_id .eq. DIAG_NULL) then
1567 call axis_ptr%axis%write_axis_metadata(fms2io_fileobj, edges_in_file)
1569 call axis_ptr%axis%write_axis_metadata(fms2io_fileobj, edges_in_file, diag_axis(parent_axis_id)%axis)
1572 if (axis_ptr%axis%is_unstructured_grid()) then
1573 structured_ids = axis_ptr%axis%get_structured_axis()
1574 do k = 1, size(structured_ids)
1575 call diag_axis(structured_ids(k))%axis%write_axis_metadata(fms2io_fileobj, .false.)
1581 end subroutine write_axis_metadata
1583 !< @brief Writes the field metadata for the file
1584 subroutine write_field_metadata(this, diag_field, diag_axis)
1585 class(fmsDiagFileContainer_type), intent(inout), target :: this !< The file object
1586 class(fmsDiagField_type) , intent(inout), target :: diag_field(:) !<
1587 class(fmsDiagAxisContainer_type), intent(in) :: diag_axis(:) !< Diag_axis object
1589 class(FmsNetcdfFile_t), pointer :: fms2io_fileobj !< The fileobj to write to
1590 class(fmsDiagFile_type), pointer :: diag_file !< Diag_file object to open
1591 class(fmsDiagField_type), pointer :: field_ptr !< diag_field(diag_file%field_ids(i)), for convenience
1593 integer :: i !< For do loops
1594 logical :: is_regional !< Flag indicating if the field is in a regional file
1595 character(len=255) :: cell_measures !< cell_measures attributes for the field
1596 logical :: need_associated_files !< .True. if the 'associated_files' global attribute is needed
1597 character(len=255) :: associated_files !< Associated files attribute to add
1599 is_regional = this%is_regional()
1601 diag_file => this%FMS_diag_file
1602 fms2io_fileobj => diag_file%fms2io_fileobj
1604 associated_files = ""
1605 need_associated_files = .false.
1606 do i = 1, size(diag_file%field_ids)
1607 if (.not. diag_file%field_registered(i)) cycle !TODO do something else here
1608 field_ptr => diag_field(diag_file%field_ids(i))
1611 if (field_ptr%has_area()) then
1612 cell_measures = "area: "//diag_field(field_ptr%get_area())%get_varname(to_write=.true.)
1614 !! Determine if the area field is already in the file. If it is not create the "associated_files" attribute
1615 !! which contains the file name of the file the area field is in. This is needed for PP/fregrid.
1616 if (.not. diag_field(field_ptr%get_area())%is_variable_in_file(diag_file%id)) then
1617 need_associated_files = .true.
1618 call diag_field(field_ptr%get_area())%generate_associated_files_att(associated_files, diag_file%start_time)
1622 if (field_ptr%has_volume()) then
1623 cell_measures = trim(cell_measures)//" volume: "//diag_field(field_ptr%get_volume())%get_varname(to_write=.true.)
1625 !! Determine if the volume field is already in the file. If it is not create the "associated_files" attribute
1626 !! which contains the file name of the file the volume field is in. This is needed for PP/fregrid.
1627 if (.not. diag_field(field_ptr%get_volume())%is_variable_in_file(diag_file%id)) then
1628 need_associated_files = .true.
1629 call diag_field(field_ptr%get_volume())%generate_associated_files_att(associated_files, diag_file%start_time)
1633 call field_ptr%write_field_metadata(fms2io_fileobj, diag_file%id, diag_file%yaml_ids(i), diag_axis, &
1634 this%FMS_diag_file%get_file_unlimdim(), is_regional, cell_measures)
1637 if (need_associated_files) &
1638 call register_global_attribute(fms2io_fileobj, "associated_files", trim(ADJUSTL(associated_files)), &
1639 str_len=len_trim(ADJUSTL(associated_files)))
1641 end subroutine write_field_metadata
1643 !< @brief Writes the axis data for the file
1644 subroutine write_axis_data(this, diag_axis)
1645 class(fmsDiagFileContainer_type), intent(inout), target :: this !< The file object
1646 class(fmsDiagAxisContainer_type), intent(in) :: diag_axis(:) !< Diag_axis object
1648 class(fmsDiagFile_type), pointer :: diag_file !< Diag_file object to open
1649 class(FmsNetcdfFile_t), pointer :: fms2io_fileobj !< The fileobj to write to
1650 integer :: i, k !< For do loops
1651 integer :: j !< diag_file%axis_ids(i) (for less typing)
1652 integer :: parent_axis_id !< Id of the parent_axis
1653 integer :: structured_ids(2) !< Ids of the uncompress axis
1655 diag_file => this%FMS_diag_file
1656 fms2io_fileobj => diag_file%fms2io_fileobj
1658 do i = 1, diag_file%number_of_axis
1659 j = diag_file%axis_ids(i)
1660 parent_axis_id = diag_axis(j)%axis%get_parent_axis_id()
1661 if (parent_axis_id .eq. DIAG_NULL) then
1662 call diag_axis(j)%axis%write_axis_data(fms2io_fileobj)
1664 call diag_axis(j)%axis%write_axis_data(fms2io_fileobj, diag_axis(parent_axis_id)%axis)
1667 if (diag_axis(j)%axis%is_unstructured_grid()) then
1668 structured_ids = diag_axis(j)%axis%get_structured_axis()
1669 do k = 1, size(structured_ids)
1670 call diag_axis(structured_ids(k))%axis%write_axis_data(fms2io_fileobj)
1675 end subroutine write_axis_data
1677 !< @brief Closes the diag_file
1678 subroutine close_diag_file(this, output_buffers, diag_fields)
1679 class(fmsDiagFileContainer_type), intent(inout), target :: this !< The file object
1680 type(fmsDiagOutputBuffer_type), intent(in) :: output_buffers(:) !< Array of output buffers
1681 !! This is needed for error checking
1682 type(fmsDiagField_type), intent(in), optional :: diag_fields(:) !< Array of diag fields
1683 !! This is needed for error checking
1685 if (.not. this%FMS_diag_file%is_file_open) return
1687 !< The select types are needed here because otherwise the code will go to the
1688 !! wrong close_file routine and things will not close propertly
1689 select type( fms2io_fileobj => this%FMS_diag_file%fms2io_fileobj)
1690 type is (FmsNetcdfDomainFile_t)
1691 call close_file(fms2io_fileobj)
1692 type is (FmsNetcdfFile_t)
1693 call close_file(fms2io_fileobj)
1694 type is (FmsNetcdfUnstructuredDomainFile_t)
1695 call close_file(fms2io_fileobj)
1698 !< Reset the unlimited dimension level back to 0, in case the fms2io_fileobj is re-used
1699 this%FMS_diag_file%unlim_dimension_level = 0
1700 this%FMS_diag_file%is_file_open = .false.
1702 if (this%FMS_diag_file%has_file_new_file_freq()) then
1703 this%FMS_diag_file%next_close = diag_time_inc(this%FMS_diag_file%next_close, &
1704 this%FMS_diag_file%get_file_new_file_freq(), &
1705 this%FMS_diag_file%get_file_new_file_freq_units())
1707 this%FMS_diag_file%next_close = diag_time_inc(this%FMS_diag_file%next_close, VERY_LARGE_FILE_FREQ, DIAG_DAYS)
1710 if (this%FMS_diag_file%has_send_data_been_called(output_buffers, .True., diag_fields)) return
1711 end subroutine close_diag_file
1713 !> \brief Gets the buffer_id list from the file object
1714 pure function get_buffer_ids (this)
1715 class(fmsDiagFile_type), intent(in) :: this !< The file object
1716 integer, allocatable :: get_buffer_ids(:) !< returned buffer ids for this file
1718 allocate(get_buffer_ids(this%number_of_buffers))
1719 get_buffer_ids = this%buffer_ids(1:this%number_of_buffers)
1720 end function get_buffer_ids
1722 !> Gets the stored number of buffers from the file object
1723 pure function get_number_of_buffers(this)
1724 class(fmsDiagFile_type), intent(in) :: this !< file object
1725 integer :: get_number_of_buffers !< returned number of buffers
1726 get_number_of_buffers = this%number_of_buffers
1727 end function get_number_of_buffers
1729 !> @brief Determine if send_data has been called for any fields in the file. Prints out warnings, if indicated
1730 !! @return .True. if send_data has been called for any fields in the file
1731 function has_send_data_been_called(this, output_buffers, print_warnings, diag_fields) &
1733 class(fmsDiagFile_type), intent(in) :: this !< file object
1734 type(fmsDiagOutputBuffer_type), intent(in), target :: output_buffers(:) !< Array of output buffers
1735 logical, intent(in) :: print_warnings !< .True. if printing warnings
1736 type(fmsDiagField_type), intent(in), optional :: diag_fields(:) !< Array of diag fields
1739 integer :: i !< For do loops
1740 integer :: field_id !< Field id
1744 if (print_warnings) then
1745 do i = 1, this%number_of_buffers
1746 if (.not. output_buffers(this%buffer_ids(i))%is_there_data_to_write()) then
1747 field_id = output_buffers(this%buffer_ids(i))%get_field_id()
1748 call mpp_error(NOTE, "Send data was never called for field:"//&
1749 trim(diag_fields(field_id)%get_varname())//" mod: "//trim(diag_fields(field_id)%get_modname())//&
1750 " in file: "//trim(this%get_file_fname())//". Writting FILL VALUES!")
1754 do i = 1, this%number_of_buffers
1755 if (output_buffers(this%buffer_ids(i))%is_there_data_to_write()) then
1761 end function has_send_data_been_called
1763 end module fms_diag_file_object_mod