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