26module fms_diag_file_object_mod
28use fms2_io_mod,
only: fmsnetcdffile_t, fmsnetcdfunstructureddomainfile_t, fmsnetcdfdomainfile_t, &
31 dimension_exists, register_global_attribute, flush_file
41 OPERATOR(/),
OPERATOR(+),
operator(<)
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
59public ::
fmsdiagfile_type, fms_diag_files_object_init, fms_diag_files_object_initialized
61logical :: fms_diag_files_object_initialized = .false.
63integer,
parameter :: var_string_len = 25
74 logical :: done_writing_data
78 logical :: is_file_open
80 class(fmsnetcdffile_t),
allocatable :: fms2io_fileobj
82 integer :: type_of_domain
86 character(len=:) ,
dimension(:),
allocatable :: file_metadata_from_model
88 integer,
dimension(:),
allocatable :: field_ids
89 integer,
dimension(:),
allocatable :: yaml_ids
90 logical,
dimension(:),
private,
allocatable :: field_registered
93 integer,
allocatable :: num_registered_fields
95 integer,
dimension(:),
allocatable :: axis_ids
96 integer :: number_of_axis
97 integer,
dimension(:),
allocatable :: buffer_ids
98 integer :: number_of_buffers
99 logical,
allocatable :: time_ops
101 integer :: unlim_dimension_level
102 logical :: data_has_been_written
104 integer :: nz_subaxis
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
129 procedure,
public :: get_file_metadata_from_model
130 procedure,
public :: get_field_ids
138 procedure,
public :: get_file_sub_region_grid_type
147 procedure,
public :: is_done_writing_data
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
169 integer,
dimension(:),
allocatable :: sub_axis_ids
170 logical :: write_on_this_pe
171 logical :: is_subaxis_defined
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
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
213logical function fms_diag_files_object_init (files_array)
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
224 if (diag_yaml%diag_files(i)%has_file_sub_region())
then
226 obj => files_array(i)%FMS_diag_file
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
237 obj => files_array(i)%FMS_diag_file
240 obj%diag_yaml_file => diag_yaml%diag_files(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()))
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
260 obj%number_of_axis = 0
263 obj%done_writing_data = .false.
268 if (obj%has_file_start_time())
then
269 obj%start_time = obj%get_file_start_time()
273 obj%last_output = obj%start_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())
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())
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())
286 obj%next_close = diag_time_inc(obj%start_time, very_large_file_freq, diag_days)
289 obj%is_file_open = .false.
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())
295 obj%no_more_data = diag_time_inc(obj%start_time, very_large_file_freq, diag_days)
298 obj%unlim_dimension_level = 0
299 obj%is_static = obj%get_file_freq() .eq. -1
304 fms_diag_files_object_init = .true.
306 fms_diag_files_object_init = .false.
310end function fms_diag_files_object_init
314pure logical function is_field_registered(this, field_id)
316 integer,
intent(in) :: field_id
318 is_field_registered = this%field_registered(field_id)
319end function is_field_registered
322subroutine add_field_and_yaml_id (this, new_field_id, yaml_id)
324 integer,
intent(in) :: new_field_id
325 integer,
intent(in) :: yaml_id
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.
333 call mpp_error(fatal,
"The file: "//this%get_file_fname()//
" has already been assigned its maximum "//&
336end subroutine add_field_and_yaml_id
339subroutine add_buffer_id (this, buffer_id)
341 integer,
intent(in) :: buffer_id
343 this%number_of_buffers = this%number_of_buffers + 1
344 this%buffer_ids(this%number_of_buffers) = buffer_id
346end subroutine add_buffer_id
351subroutine init_diurnal_axis(this, diag_axis, naxis, yaml_id)
354 integer,
intent(inout) :: naxis
355 integer,
intent(in) :: yaml_id
360 field_yaml => diag_yaml%get_diag_field_from_id(yaml_id)
363 if (.not. field_yaml%has_n_diurnal())
return
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
378 this%number_of_axis = this%number_of_axis + 1
379 this%axis_ids(this%number_of_axis) = naxis
381 this%number_of_axis = this%number_of_axis + 1
382 this%axis_ids(this%number_of_axis) = naxis - 1
384end subroutine init_diurnal_axis
387subroutine set_file_time_ops(this, VarYaml, is_static)
390 logical,
intent(in) :: is_static
391 integer,
allocatable :: var_reduct
394 if (this%is_static)
return
395 if (is_static)
return
398 if (.not.
allocated(this%time_ops))
then
399 var_reduct = varyaml%get_var_reduction()
401 select case (var_reduct)
403 this%time_ops = .true.
405 this%time_ops = .false.
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")
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")
420end subroutine set_file_time_ops
424pure logical function has_file_metadata_from_model (this)
426 has_file_metadata_from_model =
allocated(this%file_metadata_from_model)
427end function has_file_metadata_from_model
431pure logical function has_fileobj (this)
433 has_fileobj =
allocated(this%fms2io_fileobj)
434end function has_fileobj
438pure logical function has_diag_yaml_file (this)
440 has_diag_yaml_file =
associated(this%diag_yaml_file)
441end function has_diag_yaml_file
450 select case (this%diag_yaml_file%get_filename_time())
452 res = this%last_output
454 res = (this%last_output + this%next_close)/2
456 res = this%next_close
462pure logical function has_field_ids (this)
464 has_field_ids =
allocated(this%field_ids)
465end function has_field_ids
469pure function get_id (this)
result (res)
477pure function get_time_ops (this)
result (res)
481 if (.not.
allocated(this%time_ops))
then
486end function get_time_ops
509pure function get_file_metadata_from_model (this)
result (res)
511 character(len=:),
dimension(:),
allocatable :: res
512 res = this%file_metadata_from_model
513end function get_file_metadata_from_model
517pure function get_field_ids (this)
result (res)
519 integer,
dimension(:),
allocatable :: res
520 allocate(res(
size(this%field_ids)))
522end function get_field_ids
529 character (len=:),
allocatable :: res
530 res = this%diag_yaml_file%get_file_fname()
538 res = this%diag_yaml_file%get_file_frequnit()
546 res = this%diag_yaml_file%get_file_freq()
554 res = this%diag_yaml_file%get_file_timeunit()
561 character (len=:),
allocatable :: res
562 res = this%diag_yaml_file%get_file_unlimdim()
570 res => obj%diag_yaml_file%get_file_sub_region()
575function get_file_sub_region_grid_type(this) &
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
588end function get_file_sub_region_grid_type
595 res = this%diag_yaml_file%get_file_new_file_freq()
603 res = this%diag_yaml_file%get_file_new_file_freq_units()
611 res = this%diag_yaml_file%get_file_start_time()
619 res = this%diag_yaml_file%get_file_duration()
627 res = this%diag_yaml_file%get_file_duration_units()
634 character (len=:),
allocatable,
dimension(:) :: res
635 res = this%diag_yaml_file%get_file_varlist()
642 character (len=MAX_STR_LEN),
allocatable,
dimension(:,:) :: res
643 res = this%diag_yaml_file%get_file_global_meta()
648pure function is_done_writing_data (this)
result(res)
651 res = this%done_writing_data
652 if (this%is_file_open) res = .false.
653end function is_done_writing_data
660 res = this%diag_yaml_file%has_file_fname()
668 res = this%diag_yaml_file%has_file_frequnit()
676 res = this%diag_yaml_file%has_file_freq()
684 res = this%diag_yaml_file%has_file_timeunit()
692 res = this%diag_yaml_file%has_file_unlimdim()
700 res = this%diag_yaml_file%has_file_sub_region()
708 res = this%diag_yaml_file%has_file_new_file_freq()
716 res = this%diag_yaml_file%has_file_new_file_freq_units()
724 res = this%diag_yaml_file%has_file_start_time()
732 res = this%diag_yaml_file%has_file_duration()
740 res = this%diag_yaml_file%has_file_duration_units()
748 res = this%diag_yaml_file%has_file_varlist()
756 res = this%diag_yaml_file%has_file_global_meta()
760subroutine set_domain_from_axis(this, diag_axis, axes)
763 integer,
intent(in) :: axes (:)
766end subroutine set_domain_from_axis
772subroutine set_file_domain(this, domain, type_of_domain)
774 integer,
INTENT(in) :: type_of_domain
777 if (type_of_domain .ne. this%type_of_domain)
then
786 this%type_of_domain = type_of_domain
787 this%domain => domain
794 call mpp_error(fatal,
"The file: "//this%get_file_fname()//
" has variables that are not in the same domain")
798end subroutine set_file_domain
801subroutine add_axes(this, axis_ids, diag_axis, naxis, yaml_id, buffer_id, output_buffers)
803 integer,
INTENT(in) :: axis_ids(:)
805 integer,
intent(inout) :: naxis
807 integer,
intent(in) :: yaml_id
809 integer,
intent(in) :: buffer_id
815 logical :: is_cube_sphere
816 logical :: axis_found
817 integer,
allocatable :: var_axis_ids(:)
818 integer :: x_y_axis_id(2)
821 integer :: subregion_gridtype
822 logical :: write_on_this_pe
824 is_cube_sphere = .false.
825 subregion_gridtype = this%get_file_sub_region_grid_type()
827 field_yaml => diag_yaml%get_diag_field_from_id(yaml_id)
833 var_axis_ids = axis_ids
835 if (field_yaml%has_var_zbounds())
then
837 this%axis_ids, this%number_of_axis, this%nz_subaxis)
842 if (
associated(this%domain))
then
843 if (this%domain%get_ntiles() .eq. 6) is_cube_sphere = .true.
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)
851 is_x_or_y = parent_axis%is_x_or_y_axis()
852 do j = 1, this%number_of_axis
854 if(
is_parent_axis(this%axis_ids(j), var_axis_ids(i), diag_axis))
then
856 var_axis_ids(i) = this%axis_ids(j)
859 elseif (var_axis_ids(i) .eq. this%axis_ids(j))
then
864 if (.not. axis_found)
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.")
871 select case (subregion_gridtype)
872 case (index_gridtype)
875 case (latlon_gridtype)
877 .false., write_on_this_pe)
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
884 call this%add_new_axis(var_axis_ids(i))
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))
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)
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))
901 x_y_axis_id(x_or_y) = var_axis_ids(i)
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))
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
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))
919 call output_buffers(buffer_id)%add_axis_ids(var_axis_ids)
920end subroutine add_axes
923subroutine add_new_axis(this, var_axis_id)
925 integer,
intent(in) :: var_axis_id
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
932subroutine update_write_on_this_pe(this, write_on_this_pe)
934 logical,
intent(in) :: write_on_this_pe
939 if (this%write_on_this_pe) this%write_on_this_pe = write_on_this_pe
941end subroutine update_write_on_this_pe
945function get_write_on_this_pe(this) &
952 rslt= this%write_on_this_pe
954end function get_write_on_this_pe
958function does_axis_exist(this, var_axis_id) &
961 integer,
intent(in) :: var_axis_id
967 do j = 1, this%number_of_axis
969 if (var_axis_id .eq. this%axis_ids(j))
then
977subroutine define_new_subaxis(this, var_axis_ids, x_y_axis_id, is_cube_sphere, diag_axis, naxis)
979 integer,
INTENT(inout) :: var_axis_ids(:)
980 integer,
INTENT(in) :: x_y_axis_id(:)
981 logical,
intent(in) :: is_cube_sphere
982 integer,
intent(inout) :: naxis
985 logical :: write_on_this_pe
988 select case (this%get_file_sub_region_grid_type())
989 case(latlon_gridtype)
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
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)
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
1015end subroutine define_new_subaxis
1020subroutine add_start_time(this, start_time)
1022 TYPE(
time_type),
intent(in) :: start_time
1031 if (this%start_time >= start_time)
return
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")
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())
1052 if (this%is_static)
then
1055 this%next_close = this%start_time
1056 this%next_next_output = diag_time_inc(this%start_time, very_large_file_freq, diag_days)
1058 this%next_close = diag_time_inc(this%start_time, very_large_file_freq, diag_days)
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())
1066 this%no_more_data = diag_time_inc(this%start_time, very_large_file_freq, diag_days)
1074subroutine dump_file_obj(this, unit_num)
1076 integer,
intent(in) :: unit_num
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)
1085 if(
allocated(this%fms2io_fileobj))
write( unit_num, *)
'fileobj path', this%fms2io_fileobj%path
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)
1099logical pure function is_regional(this)
1102 select type (wut=>this%FMS_diag_file)
1104 is_regional = .true.
1106 is_regional = .false.
1109end function is_regional
1113logical pure function is_file_static(this)
1116is_file_static = .false.
1118select type (fileptr=>this%FMS_diag_file)
1120 is_file_static = fileptr%is_static
1123end function is_file_static
1126subroutine open_diag_file(this, time_step, file_is_opened)
1128 TYPE(time_type),
intent(in) :: time_step
1129 logical,
intent(out) :: file_is_opened
1133 class(diagdomain_t),
pointer :: domain
1134 character(len=:),
allocatable :: diag_file_name
1135 character(len=FMS_FILE_LEN) :: base_name
1137 character(len=FMS_FILE_LEN) :: file_name
1138 character(len=FMS_FILE_LEN) :: temp_name
1139 character(len=128) :: start_date
1141 character(len=128) :: suffix
1150 character(len=4) :: mype_string
1151 logical :: is_regional
1152 integer,
allocatable :: pes(:)
1154 diag_file => this%FMS_diag_file
1155 domain => diag_file%domain
1157 file_is_opened = .false.
1159 if (diag_file%is_file_open)
return
1161 is_regional = .false.
1163 if (.not.
allocated(diag_file%fms2io_fileobj))
then
1164 select type (diag_file)
1167 allocate(fmsnetcdffile_t :: diag_file%fms2io_fileobj)
1168 is_regional = .true.
1171 select case (diag_file%type_of_domain)
1173 allocate(fmsnetcdffile_t :: diag_file%fms2io_fileobj)
1175 allocate(fmsnetcdfdomainfile_t :: diag_file%fms2io_fileobj)
1177 allocate(fmsnetcdfunstructureddomainfile_t :: diag_file%fms2io_fileobj)
1183 diag_file_name = diag_file%get_file_fname()
1186 if (diag_file%has_file_new_file_freq())
then
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)
1193 base_name = trim(diag_file_name)
1197 file_name = trim(base_name)
1198 call get_instance_filename(base_name, file_name)
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
1206 file_name = trim(adjustl(start_date))//
'.'//trim(file_name)
1209 file_name = trim(file_name)//
".nc"
1212 if (is_regional)
then
1214 write(mype_string,
'(I0.4)') mpp_pe()
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)
1223 file_name = trim(file_name)//
"."//trim(mype_string)
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)
1234 allocate(pes(mpp_npes()))
1235 call mpp_get_current_pelist(pes)
1237 if (.not. open_file(fms2io_fileobj, file_name,
"overwrite", pelist=pes)) &
1238 &
call mpp_error(fatal,
"Error opening the file:"//file_name)
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)
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)
1254 file_is_opened = .true.
1255 diag_file%is_file_open = file_is_opened
1258end subroutine open_diag_file
1261subroutine write_global_metadata(this)
1264 class(fmsnetcdffile_t),
pointer :: fms2io_fileobj
1266 character (len=MAX_STR_LEN),
allocatable :: yaml_file_attributes(:,:)
1268 type(diagyamlfiles_type),
pointer :: diag_file_yaml
1270 diag_file_yaml => this%FMS_diag_file%diag_yaml_file
1271 fms2io_fileobj => this%FMS_diag_file%fms2io_fileobj
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)))
1279 deallocate(yaml_file_attributes)
1281end subroutine write_global_metadata
1284subroutine write_var_metadata(fms2io_fileobj, variable_name, dimensions, long_name, units)
1285 class(fmsnetcdffile_t),
intent(inout) :: fms2io_fileobj
1286 character(len=*) ,
intent(in) :: variable_name
1287 character(len=*) ,
intent(in) :: dimensions(:)
1288 character(len=*) ,
intent(in) :: long_name
1289 character(len=*) ,
intent(in) :: units
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
1300subroutine write_time_metadata(this)
1304 class(fmsnetcdffile_t),
pointer :: fms2io_fileobj
1305 character(len=50) :: time_units_str
1306 character(len=50) :: calendar
1308 character(len=:),
allocatable :: time_var_name
1309 character(len=50) :: dimensions(2)
1311 diag_file => this%FMS_diag_file
1312 fms2io_fileobj => diag_file%fms2io_fileobj
1314 time_var_name = diag_file%get_file_unlimdim()
1315 call register_axis(fms2io_fileobj, time_var_name, unlimited)
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)
1322 dimensions(1) =
"nv"
1323 dimensions(2) = trim(time_var_name)
1325 call write_var_metadata(fms2io_fileobj, time_var_name, dimensions(2:2), &
1326 time_var_name, time_units_str)
1329 call register_variable_attribute(fms2io_fileobj, time_var_name,
"axis",
"T", str_len=1 )
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))
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"))
1344 if ( .not. dimension_exists(fms2io_fileobj,
"nv"))
then
1345 call register_axis(fms2io_fileobj,
"nv", 2)
1346 call write_var_metadata(fms2io_fileobj,
"nv", dimensions(1:1), &
1347 "vertex number", no_units)
1349 call write_var_metadata(fms2io_fileobj, time_var_name//
"_bnds", dimensions, &
1350 trim(time_var_name)//
" axis boundaries", time_units_str)
1353end subroutine write_time_metadata
1356subroutine write_field_data(this, field_obj, buffer_obj, unlim_dim_was_increased)
1358 type(fmsdiagfield_type),
intent(in),
target :: field_obj
1359 type(fmsdiagoutputbuffer_type),
intent(inout),
target :: buffer_obj
1360 logical,
intent(inout) :: unlim_dim_was_increased
1363 class(fmsnetcdffile_t),
pointer :: fms2io_fileobj
1364 logical :: has_diurnal
1366 diag_file => this%FMS_diag_file
1367 fms2io_fileobj => diag_file%fms2io_fileobj
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.
1378 if (diag_file%is_static)
then
1381 call buffer_obj%write_buffer(fms2io_fileobj)
1382 diag_file%data_has_been_written = .true.
1384 if (field_obj%is_static())
then
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.
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)
1398end subroutine write_field_data
1402logical function is_time_to_close_file (this, time_step, force_close)
1404 TYPE(time_type),
intent(in) :: time_step
1405 logical,
intent(in) :: force_close
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.
1412 if (this%FMS_diag_file%is_static)
then
1413 is_time_to_close_file = .true.
1415 is_time_to_close_file = .false.
1422logical function time_to_start_doing_math (this)
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.
1431subroutine check_file_times(this, time_step, output_buffers, diag_fields, do_not_write)
1433 TYPE(time_type),
intent(in) :: time_step
1434 type(fmsdiagoutputbuffer_type),
intent(in) :: output_buffers(:)
1436 type(fmsdiagfield_type),
intent(in) :: diag_fields(:)
1437 logical,
intent(out) :: do_not_write
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
1447 if (time_step .ne. this%FMS_diag_file%next_output)
then
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
1453 elseif (this%FMS_diag_file%num_registered_fields .eq. 0)
then
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
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.")
1471 if(this%FMS_diag_file%get_file_freq() .eq. 0)
then
1472 do_not_write = .true.
1475end subroutine check_file_times
1478logical function writing_on_this_pe(this)
1481 select type(diag_file => this%FMS_diag_file)
1483 writing_on_this_pe = diag_file%write_on_this_pe
1485 writing_on_this_pe = .true.
1491subroutine write_time_data(this)
1496 class(fmsnetcdffile_t),
pointer :: fms2io_fileobj
1497 TYPE(time_type) :: middle_time
1502 diag_file => this%FMS_diag_file
1503 fms2io_fileobj => diag_file%fms2io_fileobj
1506 if (.not. diag_file%data_has_been_written)
return
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())
1512 dif = get_date_dif(diag_file%next_output, get_base_time(), diag_file%get_file_timeunit())
1515 call write_data(fms2io_fileobj, diag_file%get_file_unlimdim(), dif, &
1516 unlim_dim_level=diag_file%unlim_dimension_level)
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())
1522 call write_data(fms2io_fileobj, trim(diag_file%get_file_unlimdim())//
"_bnds", &
1523 (/t1, t2/), unlim_dim_level=diag_file%unlim_dimension_level)
1525 if (diag_file%unlim_dimension_level .eq. 1)
then
1526 call write_data(fms2io_fileobj,
"nv", (/1, 2/))
1530 diag_file%data_has_been_written = .false.
1531end subroutine write_time_data
1534subroutine update_current_new_file_freq_index(this, time_step)
1536 TYPE(time_type),
intent(in) :: time_step
1540 diag_file => this%FMS_diag_file
1542 if (time_step >= diag_file%no_more_data)
then
1543 call diag_file%diag_yaml_file%increase_new_file_freq_index()
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())
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
1558 if (diag_file%is_static) diag_file%done_writing_data = .true.
1559end subroutine update_current_new_file_freq_index
1562subroutine update_next_write(this, time_step)
1564 TYPE(time_type),
intent(in) :: time_step
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)
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())
1581end subroutine update_next_write
1584subroutine prepare_for_force_write(this)
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.
1591end subroutine prepare_for_force_write
1594subroutine init_unlim_dim(this, output_buffers)
1596 type(fmsdiagoutputbuffer_type),
intent(in),
target :: output_buffers(:)
1599 type(fmsdiagoutputbuffer_type),
pointer :: output_buffer_obj
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()
1608end subroutine init_unlim_dim
1612pure function get_unlim_dimension_level(this) &
1617 res = this%FMS_diag_file%unlim_dimension_level
1621subroutine flush_diag_file(this)
1624 if (flush_nc_files)
then
1625 call flush_file(this%FMS_diag_file%fms2io_fileobj)
1627end subroutine flush_diag_file
1631pure function get_next_output(this) &
1634 type(time_type) :: res
1636 res = this%FMS_diag_file%next_output
1637end function get_next_output
1641pure function get_next_next_output(this) &
1644 type(time_type) :: res
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
1650end function get_next_next_output
1653subroutine write_axis_metadata(this, diag_axis)
1655 class(fmsdiagaxiscontainer_type),
intent(in),
target :: diag_axis(:)
1658 class(fmsnetcdffile_t),
pointer :: fms2io_fileobj
1660 integer :: parent_axis_id
1661 integer :: structured_ids(2)
1664 class(fmsdiagaxiscontainer_type),
pointer :: axis_ptr
1665 logical :: edges_in_file
1667 diag_file => this%FMS_diag_file
1668 fms2io_fileobj => diag_file%fms2io_fileobj
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()
1675 edges_id = axis_ptr%axis%get_edges_id()
1676 if (edges_id .ne. diag_null)
then
1678 if (any(diag_file%axis_ids(1:diag_file%number_of_axis) .eq. edges_id))
then
1679 edges_in_file = .true.
1681 call diag_axis(edges_id)%axis%write_axis_metadata(fms2io_fileobj, .true.)
1682 call diag_file%add_new_axis(edges_id)
1686 if (parent_axis_id .eq. diag_null)
then
1687 call axis_ptr%axis%write_axis_metadata(fms2io_fileobj, edges_in_file)
1689 call axis_ptr%axis%write_axis_metadata(fms2io_fileobj, edges_in_file, diag_axis(parent_axis_id)%axis)
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.)
1701end subroutine write_axis_metadata
1704subroutine write_field_metadata(this, diag_field, diag_axis)
1706 class(fmsdiagfield_type) ,
intent(inout),
target :: diag_field(:)
1707 class(fmsdiagaxiscontainer_type),
intent(in) :: diag_axis(:)
1709 class(fmsnetcdffile_t),
pointer :: fms2io_fileobj
1711 class(fmsdiagfield_type),
pointer :: field_ptr
1714 logical :: is_regional
1715 character(len=255) :: cell_measures
1716 logical :: need_associated_files
1717 character(len=FMS_FILE_LEN) :: associated_files
1719 is_regional = this%is_regional()
1721 diag_file => this%FMS_diag_file
1722 fms2io_fileobj => diag_file%fms2io_fileobj
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
1728 field_ptr => diag_field(diag_file%field_ids(i))
1731 if (field_ptr%has_area())
then
1732 cell_measures =
"area: "//diag_field(field_ptr%get_area())%get_varname(to_write=.true.)
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)
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.)
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)
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)
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)))
1761end subroutine write_field_metadata
1764subroutine write_axis_data(this, diag_axis)
1766 class(fmsdiagaxiscontainer_type),
intent(in) :: diag_axis(:)
1769 class(fmsnetcdffile_t),
pointer :: fms2io_fileobj
1772 integer :: parent_axis_id
1773 integer :: structured_ids(2)
1775 diag_file => this%FMS_diag_file
1776 fms2io_fileobj => diag_file%fms2io_fileobj
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)
1784 call diag_axis(j)%axis%write_axis_data(fms2io_fileobj, diag_axis(parent_axis_id)%axis)
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)
1795end subroutine write_axis_data
1798subroutine close_diag_file(this, output_buffers, model_end_time, diag_fields)
1800 type(fmsdiagoutputbuffer_type),
intent(in) :: output_buffers(:)
1802 type(time_type),
intent(in) :: model_end_time
1803 type(fmsdiagfield_type),
intent(in),
optional :: diag_fields(:)
1806 if (.not. this%FMS_diag_file%is_file_open)
return
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)
1820 this%FMS_diag_file%unlim_dimension_level = 0
1821 this%FMS_diag_file%is_file_open = .false.
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())
1828 this%FMS_diag_file%next_close = model_end_time
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
1837subroutine set_model_time(this, model_time)
1839 type(time_type),
intent(in) :: model_time
1841 if (model_time > this%FMS_diag_file%model_time) this%FMS_diag_file%model_time = model_time
1846function get_model_time(this) &
1849 type(time_type),
pointer :: rslt
1851 rslt => this%FMS_diag_file%model_time
1852end function get_model_time
1855pure function get_buffer_ids (this)
1857 integer,
allocatable :: get_buffer_ids(:)
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
1864pure function get_number_of_buffers(this)
1866 integer :: get_number_of_buffers
1867 get_number_of_buffers = this%number_of_buffers
1868end function get_number_of_buffers
1872subroutine check_buffer_times(this, output_buffers, diag_fields)
1874 type(fmsdiagoutputbuffer_type),
intent(in),
target :: output_buffers(:)
1875 type(fmsdiagfield_type),
intent(in) :: diag_fields(:)
1878 type(time_type) :: current_buffer_time
1879 character(len=:),
allocatable :: field_name
1880 logical :: buffer_time_set
1881 type(fmsdiagoutputbuffer_type),
pointer :: output_buffer_obj
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.
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())
1902function has_send_data_been_called(this, output_buffers, print_warnings, diag_fields) &
1905 type(fmsdiagoutputbuffer_type),
intent(in),
target :: output_buffers(:)
1906 logical,
intent(in) :: print_warnings
1907 type(fmsdiagfield_type),
intent(in),
optional :: diag_fields(:)
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!")
1925 do i = 1, this%number_of_buffers
1926 if (output_buffers(this%buffer_ids(i))%is_there_data_to_write())
then
1932end function has_send_data_been_called
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.
Opens a given netcdf or domain file.
Add a dimension to a given file.
Defines a new field within the given file Example usage:
Write data to a defined field within a file Example usage:
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.
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 2d domain.
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 diurnal axis.
Type to hold the diagnostic axis description.
Type to hold the subaxis.
Object that holds all variable information.
A container for fmsDiagFile_type. This is used to create the array of files.
holds an allocated buffer0-5d object
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