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