1 !***********************************************************************
2 !* GNU Lesser General Public License
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 !> @author Ryan Mulhall
20 !> @email ryan.mulhall@noaa.gov
21 !! @brief Contains buffer types and routines for the diag manager
23 !! @description Holds buffered data for fmsDiagVars_type objects
24 !! buffer0-5d types extend fmsDiagBuffer_class, and upon allocation
25 !! are added to the module's buffer_lists depending on it's dimension
26 module fms_diag_output_buffer_mod
30 use time_manager_mod, only: time_type, operator(==), operator(>=), get_ticks_per_second, get_time, operator(>)
31 use constants_mod, only: SECONDS_PER_DAY
32 use mpp_mod, only: mpp_error, FATAL, NOTE, mpp_pe, mpp_root_pe
33 use diag_data_mod, only: DIAG_NULL, DIAG_NOT_REGISTERED, i4, i8, r4, r8, get_base_time, MIN_VALUE, MAX_VALUE, EMPTY, &
35 use fms2_io_mod, only: FmsNetcdfFile_t, write_data, FmsNetcdfDomainFile_t, FmsNetcdfUnstructuredDomainFile_t
36 use fms_diag_yaml_mod, only: diag_yaml
37 use fms_diag_bbox_mod, only: fmsDiagIbounds_type
38 use fms_diag_reduction_methods_mod, only: do_time_none, do_time_min, do_time_max, do_time_sum_update, time_update_done
39 use fms_diag_time_utils_mod, only: diag_time_inc
45 !> holds an allocated buffer0-5d object
46 type :: fmsDiagOutputBuffer_type
47 integer :: buffer_id !< index in buffer list
48 integer(i4_kind) :: buffer_type !< set to allocated data type & kind value, one of i4,i8,r4,r8
49 class(*), allocatable :: buffer(:,:,:,:,:) !< 5D numeric data array
50 integer :: ndim !< Number of dimensions for each variable
51 integer, allocatable :: buffer_dims(:) !< holds the size of each dimension in the buffer
52 real(r8_kind), allocatable :: weight_sum(:,:,:,:) !< Weight sum as an array
53 !! (this will be have a size of 1,1,1,1 when not using variable
55 integer, allocatable :: num_elements(:) !< used in time-averaging
56 integer, allocatable :: axis_ids(:) !< Axis ids for the buffer
57 integer :: field_id !< The id of the field the buffer belongs to
58 integer :: yaml_id !< The id of the yaml id the buffer belongs to
59 logical :: done_with_math !< .True. if done doing the math
60 integer :: diurnal_sample_size = -1 !< dirunal sample size as read in from the reduction method
61 !! ie. diurnal24 = sample size of 24
62 integer :: diurnal_section= -1 !< the diurnal section (ie 5th index) calculated from the current model
63 !! time and sample size if using a diurnal reduction
64 logical, allocatable :: send_data_called !< .True. if send_data has been called
65 type(time_type) :: time !< The last time the data was received
66 type(time_type) :: next_output !< The next time to output the data
69 procedure :: add_axis_ids
70 procedure :: get_axis_ids
71 procedure :: set_field_id
72 procedure :: get_field_id
73 procedure :: set_yaml_id
74 procedure :: get_yaml_id
75 procedure :: init_buffer_time
76 procedure :: set_next_output
77 procedure :: update_buffer_time
78 procedure :: is_there_data_to_write
79 procedure :: is_time_to_finish_reduction
80 procedure :: set_send_data_called
81 procedure :: is_done_with_math
82 procedure :: set_done_with_math
83 procedure :: write_buffer
84 !! These are needed because otherwise the write_data calls will go into the wrong interface
85 procedure :: write_buffer_wrapper_netcdf
86 procedure :: write_buffer_wrapper_domain
87 procedure :: write_buffer_wrapper_u
88 procedure :: allocate_buffer
89 procedure :: initialize_buffer
90 procedure :: get_buffer
91 procedure :: flush_buffer
92 procedure :: do_time_none_wrapper
93 procedure :: do_time_min_wrapper
94 procedure :: do_time_max_wrapper
95 procedure :: do_time_sum_wrapper
96 procedure :: diag_reduction_done_wrapper
97 procedure :: get_buffer_dims
98 procedure :: get_diurnal_sample_size
99 procedure :: set_diurnal_sample_size
100 procedure :: set_diurnal_section_index
101 procedure :: get_remapped_diurnal_data
102 end type fmsDiagOutputBuffer_type
105 public :: fmsDiagOutputBuffer_type
108 public :: fms_diag_output_buffer_init
112 !!--------module routines
114 !> Initializes a list of diag buffers
115 !> @returns true if allocation is successfull
116 logical function fms_diag_output_buffer_init(buffobjs, buff_list_size)
117 type(fmsDiagOutputBuffer_type), allocatable, intent(out) :: buffobjs(:) !< an array of buffer container types
119 integer, intent(in) :: buff_list_size !< size of buffer array to allocate
121 if (allocated(buffobjs)) call mpp_error(FATAL,'fms_diag_buffer_init: passed in buffobjs array is already allocated')
122 allocate(buffobjs(buff_list_size))
123 fms_diag_output_buffer_init = allocated(buffobjs)
124 end function fms_diag_output_buffer_init
126 !!--------generic routines for any fmsDiagBuffer_class objects
128 !> Setter for buffer_id for any buffer objects
129 subroutine set_buffer_id(this, id)
130 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< buffer object to set id for
131 integer, intent(in) :: id !< positive integer id to set
134 end subroutine set_buffer_id
136 !> Deallocates data fields from a buffer object.
137 subroutine flush_buffer(this)
138 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< any buffer object
140 this%buffer_id = diag_null
141 this%buffer_type = diag_null
142 this%ndim = diag_null
143 this%field_id = diag_null
144 this%yaml_id = diag_null
145 if (allocated(this%buffer)) deallocate(this%buffer)
146 if (allocated(this%buffer_dims)) deallocate(this%buffer_dims)
147 if (allocated(this%num_elements)) deallocate(this%num_elements)
148 if (allocated(this%axis_ids)) deallocate(this%axis_ids)
149 if (allocated(this%weight_sum)) deallocate(this%weight_sum)
150 end subroutine flush_buffer
152 !> Allocates a 5D buffer to given buff_type.
153 subroutine allocate_buffer(this, buff_type, ndim, buff_sizes, mask_variant, field_name, diurnal_samples)
154 class(fmsDiagOutputBuffer_type), intent(inout), target :: this !< 5D buffer object
155 class(*), intent(in) :: buff_type !< allocates to the type of buff_type
156 integer, intent(in) :: ndim !< Number of dimension
157 integer, intent(in) :: buff_sizes(4) !< dimension buff_sizes
158 logical, intent(in) :: mask_variant !< Mask changes over time
159 character(len=*), intent(in) :: field_name !< field name for error output
160 integer, intent(in) :: diurnal_samples !< number of diurnal samples
162 integer :: n_samples !< number of diurnal samples, defaults to 1
164 n_samples = MAX(1, diurnal_samples)
165 call this%set_diurnal_sample_size(n_samples)
168 if(allocated(this%buffer)) call mpp_error(FATAL, "allocate_buffer: buffer already allocated for field:" // &
170 select type (buff_type)
171 type is (integer(kind=i4_kind))
172 allocate(integer(kind=i4_kind) :: this%buffer(buff_sizes(1),buff_sizes(2),buff_sizes(3),buff_sizes(4), &
174 this%buffer_type = i4
175 type is (integer(kind=i8_kind))
176 allocate(integer(kind=i8_kind) :: this%buffer(buff_sizes(1),buff_sizes(2),buff_sizes(3),buff_sizes(4), &
178 this%buffer_type = i8
179 type is (real(kind=r4_kind))
180 allocate(real(kind=r4_kind) :: this%buffer(buff_sizes(1),buff_sizes(2),buff_sizes(3),buff_sizes(4), &
182 this%buffer_type = r4
183 type is (real(kind=r8_kind))
184 allocate(real(kind=r8_kind) :: this%buffer(buff_sizes(1),buff_sizes(2),buff_sizes(3),buff_sizes(4), &
186 this%buffer_type = r8
188 call mpp_error("allocate_buffer", &
189 "The buff_type value passed to allocate a buffer is not a r8, r4, i8, or i4" // &
190 "for field:" // field_name, FATAL)
192 if (mask_variant) then
193 allocate(this%weight_sum(buff_sizes(1),buff_sizes(2),buff_sizes(3),buff_sizes(4)))
195 allocate(this%weight_sum(1,1,1,1))
197 this%weight_sum = 0.0_r8_kind
199 allocate(this%num_elements(n_samples))
200 this%num_elements = 0
201 this%done_with_math = .false.
202 this%send_data_called = .false.
203 allocate(this%buffer_dims(5))
204 this%buffer_dims(1) = buff_sizes(1)
205 this%buffer_dims(2) = buff_sizes(2)
206 this%buffer_dims(3) = buff_sizes(3)
207 this%buffer_dims(4) = buff_sizes(4)
208 this%buffer_dims(5) = n_samples
209 end subroutine allocate_buffer
211 !> Get routine for 5D buffers.
212 !! Sets the buff_out argument to the integer or real array currently stored in the buffer.
213 subroutine get_buffer (this, buff_out, field_name)
214 class(fmsDiagOutputBuffer_type), intent(in) :: this !< 5d allocated buffer object
215 class(*), allocatable, intent(out) :: buff_out(:,:,:,:,:) !< output of copied buffer data
216 !! must be the same size as the allocated buffer
217 character(len=*), intent(in) :: field_name !< field name for error output
219 integer(i4_kind) :: buff_size(5)!< sizes for allocated buffer
221 if(.not. allocated(this%buffer)) call mpp_error(FATAL, 'get_buffer: buffer not yet allocated for field:' &
223 buff_size(1) = size(this%buffer,1)
224 buff_size(2) = size(this%buffer,2)
225 buff_size(3) = size(this%buffer,3)
226 buff_size(4) = size(this%buffer,4)
227 buff_size(5) = size(this%buffer,5)
229 select type (buff=>this%buffer)
230 type is (real(r4_kind))
231 allocate(real(r4_kind) :: buff_out(buff_size(1), buff_size(2), buff_size(3), buff_size(4), buff_size(5)))
233 type is (real(r8_kind))
234 allocate(real(r8_kind) :: buff_out(buff_size(1), buff_size(2), buff_size(3), buff_size(4), buff_size(5)))
236 type is (integer(i4_kind))
237 allocate(integer(i4_kind) :: buff_out(buff_size(1), buff_size(2), buff_size(3), buff_size(4), buff_size(5)))
239 type is (integer(i8_kind))
240 allocate(integer(i8_kind) :: buff_out(buff_size(1), buff_size(2), buff_size(3), buff_size(4), buff_size(5)))
243 call mpp_error(FATAL, "get_buffer: buffer allocated to invalid type(must be integer or real, kind size 4 or 8)."&
244 //"field name: "// field_name)
248 !> @brief Initializes a buffer based on the reduction method
249 subroutine initialize_buffer (this, reduction_method, field_name)
250 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< allocated 5D buffer object
251 integer, intent(in) :: reduction_method !< The reduction method for the field
252 character(len=*), intent(in) :: field_name !< field name for error output
254 if(.not. allocated(this%buffer)) call mpp_error(FATAL, 'initialize_buffer: field:'// field_name // &
255 'buffer not yet allocated, allocate_buffer() must be called on this object first.')
257 select type(buff => this%buffer)
258 type is(real(r8_kind))
259 select case (reduction_method)
261 buff = real(MIN_VALUE, kind=r8_kind)
263 buff = real(MAX_VALUE, kind=r8_kind)
265 buff = real(EMPTY, kind=r8_kind)
267 type is(real(r4_kind))
268 select case (reduction_method)
270 buff = real(MIN_VALUE, kind=r4_kind)
272 buff = real(MAX_VALUE, kind=r4_kind)
274 buff = real(EMPTY, kind=r4_kind)
276 type is(integer(i8_kind))
277 select case (reduction_method)
279 buff = int(MIN_VALUE, kind=i8_kind)
281 buff = int(MAX_VALUE, kind=i8_kind)
283 buff = int(EMPTY, kind=i8_kind)
285 type is(integer(i4_kind))
286 select case (reduction_method)
288 buff = int(MIN_VALUE, kind=i4_kind)
290 buff = int(MAX_VALUE, kind=i4_kind)
292 buff = int(EMPTY, kind=i4_kind)
295 call mpp_error(FATAL, 'initialize buffer_5d: buffer allocated to invalid data type, this shouldnt happen')
298 end subroutine initialize_buffer
300 !> @brief Adds the axis ids to the buffer object
301 subroutine add_axis_ids(this, axis_ids)
302 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< Buffer object
303 integer, intent(in) :: axis_ids(:) !< Axis ids to add
305 this%axis_ids = axis_ids
308 !> @brief Get the axis_ids for the buffer
309 !! @return Axis_ids, if the buffer doesn't have axis ids it returns diag_null
310 subroutine get_axis_ids(this, res)
311 class(fmsDiagOutputBuffer_type), target, intent(inout) :: this !< Buffer object
312 integer, pointer, intent(out) :: res(:)
314 if (allocated(this%axis_ids)) then
322 !> @brief Get the field id of the buffer
323 !! @return the field id of the buffer
324 function get_field_id(this) &
326 class(fmsDiagOutputBuffer_type), intent(in) :: this !< Buffer object
330 end function get_field_id
332 !> @brief set the field id of the buffer
333 subroutine set_field_id(this, field_id)
334 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< Buffer object
335 integer, intent(in) :: field_id !< field id of the buffer
337 this%field_id = field_id
338 end subroutine set_field_id
340 !> @brief set the field id of the buffer
341 subroutine set_yaml_id(this, yaml_id)
342 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< Buffer object
343 integer, intent(in) :: yaml_id !< yaml id of the buffer
345 this%yaml_id = yaml_id
346 end subroutine set_yaml_id
348 !> @brief inits the buffer time for the buffer
349 subroutine init_buffer_time(this, time)
350 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< Buffer object
351 type(time_type), optional, intent(in) :: time !< time to add to the buffer
353 if (present(time)) then
356 this%time = get_base_time()
358 end subroutine init_buffer_time
360 !> @brief Sets the next output
361 subroutine set_next_output(this, time, is_static)
362 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< Buffer object
363 type(time_type), intent(in) :: time !< time to add to the buffer
364 logical, optional, intent(in) :: is_static !< .True. if the field is static
366 this%next_output = time
367 if (present(is_static)) then
368 !< If the field is static set the next_output to be equal to time
369 !! this should only be used in the init, so next_output will be equal to the the init time
370 if (is_static) this%next_output = this%time
372 end subroutine set_next_output
374 !> @brief Update the buffer time if it is a new time
375 subroutine update_buffer_time(this, time)
376 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< Buffer object
377 type(time_type), intent(in) :: time !< Current model time
379 if (time > this%time) then
382 end subroutine update_buffer_time
384 !> @brief Determine if finished with math
385 !! @return this%done_with_math
386 function is_done_with_math(this) &
388 class(fmsDiagOutputBuffer_type), intent(in) :: this !< Buffer object
391 res = this%done_with_math
392 end function is_done_with_math
394 !> @brief Set done_with_math to .true.
395 subroutine set_done_with_math(this)
396 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< Buffer object
399 this%done_with_math = .true.
400 end subroutine set_done_with_math
402 !> @brief Get the yaml id of the buffer
403 !! @return the yaml id of the buffer
404 function get_yaml_id(this) &
407 class(fmsDiagOutputBuffer_type), intent(in) :: this !< Buffer object
411 end function get_yaml_id
413 !> @brief Write the buffer to the file
414 subroutine write_buffer(this, fms2io_fileobj, unlim_dim_level, is_diurnal)
415 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< buffer object to write
416 class(FmsNetcdfFile_t), intent(in) :: fms2io_fileobj !< fileobj to write to
417 integer, optional, intent(in) :: unlim_dim_level !< unlimited dimension
418 logical, optional, intent(in) :: is_diurnal !< should be set if using diurnal
419 !! reductions so buffer data can be remapped
421 select type(fms2io_fileobj)
422 type is (FmsNetcdfFile_t)
423 call this%write_buffer_wrapper_netcdf(fms2io_fileobj, unlim_dim_level=unlim_dim_level, is_diurnal=is_diurnal)
424 type is (FmsNetcdfDomainFile_t)
425 call this%write_buffer_wrapper_domain(fms2io_fileobj, unlim_dim_level=unlim_dim_level, is_diurnal=is_diurnal)
426 type is (FmsNetcdfUnstructuredDomainFile_t)
427 call this%write_buffer_wrapper_u(fms2io_fileobj, unlim_dim_level=unlim_dim_level, is_diurnal=is_diurnal)
429 call mpp_error(FATAL, "The file "//trim(fms2io_fileobj%path)//" is not one of the accepted types"//&
430 " only FmsNetcdfFile_t, FmsNetcdfDomainFile_t, and FmsNetcdfUnstructuredDomainFile_t are accepted.")
433 call this%initialize_buffer(diag_yaml%diag_fields(this%yaml_id)%get_var_reduction(), &
434 diag_yaml%diag_fields(this%yaml_id)%get_var_outname())
435 !TODO Set the counters back to 0
436 end subroutine write_buffer
438 !> @brief Write the buffer to the FmsNetcdfFile_t fms2io_fileobj
439 subroutine write_buffer_wrapper_netcdf(this, fms2io_fileobj, unlim_dim_level, is_diurnal)
440 class(fmsDiagOutputBuffer_type), intent(in) :: this !< buffer object to write
441 type(FmsNetcdfFile_t), intent(in) :: fms2io_fileobj !< fileobj to write to
442 integer, optional, intent(in) :: unlim_dim_level !< unlimited dimension
443 logical, optional, intent(in) :: is_diurnal !< should be set if using diurnal
444 !! reductions so buffer data can be remapped
445 character(len=:), allocatable :: varname !< name of the variable
446 logical :: using_diurnal !< local copy of is_diurnal if present
447 class(*), allocatable :: buff_ptr(:,:,:,:,:) !< pointer for buffer to write
449 using_diurnal = .false.
450 if( present(is_diurnal) ) using_diurnal = is_diurnal
451 if( using_diurnal ) then
452 call this%get_remapped_diurnal_data(buff_ptr)
454 buff_ptr = this%buffer
457 varname = diag_yaml%diag_fields(this%yaml_id)%get_var_outname()
458 select case(this%ndim)
460 call write_data(fms2io_fileobj, varname, buff_ptr(1,1,1,1,1), unlim_dim_level=unlim_dim_level)
462 call write_data(fms2io_fileobj, varname, buff_ptr(:,1,1,1,1), unlim_dim_level=unlim_dim_level)
464 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,1,1,1), unlim_dim_level=unlim_dim_level)
466 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,1,1), unlim_dim_level=unlim_dim_level)
468 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,1), unlim_dim_level=unlim_dim_level)
470 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,:), unlim_dim_level=unlim_dim_level)
472 end subroutine write_buffer_wrapper_netcdf
474 !> @brief Write the buffer to the FmsNetcdfDomainFile_t fms2io_fileobj
475 subroutine write_buffer_wrapper_domain(this, fms2io_fileobj, unlim_dim_level, is_diurnal)
476 class(fmsDiagOutputBuffer_type), intent(in) :: this !< buffer object to write
477 type(FmsNetcdfDomainFile_t), intent(in) :: fms2io_fileobj !< fileobj to write to
478 integer, optional, intent(in) :: unlim_dim_level !< unlimited dimension
479 logical, optional, intent(in) :: is_diurnal !< should be set if using diurnal
480 !! reductions so buffer data can be remapped
482 character(len=:), allocatable :: varname !< name of the variable
483 logical :: using_diurnal !< local copy of is_diurnal if present
484 class(*), allocatable :: buff_ptr(:,:,:,:,:) !< pointer to buffer to write
486 using_diurnal = .false.
487 if( present(is_diurnal) ) using_diurnal = is_diurnal
488 if( using_diurnal ) then
489 call this%get_remapped_diurnal_data(buff_ptr)
491 buff_ptr = this%buffer
494 varname = diag_yaml%diag_fields(this%yaml_id)%get_var_outname()
495 select case(this%ndim)
497 call write_data(fms2io_fileobj, varname, buff_ptr(1,1,1,1,1), unlim_dim_level=unlim_dim_level)
499 call write_data(fms2io_fileobj, varname, buff_ptr(:,1,1,1,1), unlim_dim_level=unlim_dim_level)
501 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,1,1,1), unlim_dim_level=unlim_dim_level)
503 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,1,1), unlim_dim_level=unlim_dim_level)
505 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,1), unlim_dim_level=unlim_dim_level)
507 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,:), unlim_dim_level=unlim_dim_level)
509 end subroutine write_buffer_wrapper_domain
511 !> @brief Write the buffer to the FmsNetcdfUnstructuredDomainFile_t fms2io_fileobj
512 subroutine write_buffer_wrapper_u(this, fms2io_fileobj, unlim_dim_level, is_diurnal)
513 class(fmsDiagOutputBuffer_type), intent(in) :: this !< buffer object to write
514 type(FmsNetcdfUnstructuredDomainFile_t), intent(in) :: fms2io_fileobj !< fileobj to write to
515 integer, optional, intent(in) :: unlim_dim_level !< unlimited dimension
516 logical, optional, intent(in) :: is_diurnal !< should be set if using diurnal
517 !! reductions so buffer data can be remapped
519 character(len=:), allocatable :: varname !< name of the variable
520 logical :: using_diurnal !< local copy of is_diurnal if present
521 class(*), allocatable :: buff_ptr(:,:,:,:,:) !< pointer for buffer to write
523 using_diurnal = .false.
524 if( present(is_diurnal) ) using_diurnal = is_diurnal
525 if( using_diurnal ) then
526 call this%get_remapped_diurnal_data(buff_ptr)
528 buff_ptr = this%buffer
531 varname = diag_yaml%diag_fields(this%yaml_id)%get_var_outname()
532 select case(this%ndim)
534 call write_data(fms2io_fileobj, varname, buff_ptr(1,1,1,1,1), unlim_dim_level=unlim_dim_level)
536 call write_data(fms2io_fileobj, varname, buff_ptr(:,1,1,1,1), unlim_dim_level=unlim_dim_level)
538 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,1,1,1), unlim_dim_level=unlim_dim_level)
540 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,1,1), unlim_dim_level=unlim_dim_level)
542 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,1), unlim_dim_level=unlim_dim_level)
544 call write_data(fms2io_fileobj, varname, buff_ptr(:,:,:,:,:), unlim_dim_level=unlim_dim_level)
546 end subroutine write_buffer_wrapper_u
548 !> @brief Does the time_none reduction method on the buffer object
549 !! @return Error message if the math was not successful
550 function do_time_none_wrapper(this, field_data, mask, is_masked, bounds_in, bounds_out, missing_value) &
552 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< buffer object to write
553 class(*), intent(in) :: field_data(:,:,:,:) !< Buffer data for current time
554 type(fmsDiagIbounds_type), intent(in) :: bounds_in !< Indicies for the buffer passed in
555 type(fmsDiagIbounds_type), intent(in) :: bounds_out !< Indicies for the output buffer
556 logical, intent(in) :: mask(:,:,:,:) !< Mask for the field
557 logical, intent(in) :: is_masked !< .True. if the field has a mask
558 real(kind=r8_kind), intent(in) :: missing_value !< Missing_value for data points that are masked
559 character(len=50) :: err_msg
561 !TODO This will be expanded for integers
563 select type (output_buffer => this%buffer)
564 type is (real(kind=r8_kind))
565 select type (field_data)
566 type is (real(kind=r8_kind))
567 call do_time_none(output_buffer, field_data, mask, is_masked, bounds_in, bounds_out, missing_value)
569 err_msg="do_time_none_wrapper::the output buffer and the buffer send in are not of the same type (r8_kind)"
571 type is (real(kind=r4_kind))
572 select type (field_data)
573 type is (real(kind=r4_kind))
574 call do_time_none(output_buffer, field_data, mask, is_masked, bounds_in, bounds_out, &
575 real(missing_value, kind=r4_kind))
577 err_msg="do_time_none_wrapper::the output buffer and the buffer send in are not of the same type (r4_kind)"
580 end function do_time_none_wrapper
582 !> @brief Does the time_min reduction method on the buffer object
583 !! @return Error message if the math was not successful
584 function do_time_min_wrapper(this, field_data, mask, is_masked, bounds_in, bounds_out, missing_value) &
586 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< buffer object to write
587 class(*), intent(in) :: field_data(:,:,:,:) !< Buffer data for current time
588 type(fmsDiagIbounds_type), intent(in) :: bounds_in !< Indicies for the buffer passed in
589 type(fmsDiagIbounds_type), intent(in) :: bounds_out !< Indicies for the output buffer
590 logical, intent(in) :: mask(:,:,:,:) !< Mask for the field
591 logical, intent(in) :: is_masked !< .True. if the field has a mask
592 real(kind=r8_kind), intent(in) :: missing_value !< Missing_value for data points that are masked
593 character(len=50) :: err_msg
595 !TODO This will be expanded for integers
597 select type (output_buffer => this%buffer)
598 type is (real(kind=r8_kind))
599 select type (field_data)
600 type is (real(kind=r8_kind))
601 call do_time_min(output_buffer, field_data, mask, is_masked, bounds_in, bounds_out, missing_value)
603 err_msg="do_time_min_wrapper::the output buffer and the buffer send in are not of the same type (r8_kind)"
605 type is (real(kind=r4_kind))
606 select type (field_data)
607 type is (real(kind=r4_kind))
608 call do_time_min(output_buffer, field_data, mask, is_masked, bounds_in, bounds_out, &
609 real(missing_value, kind=r4_kind))
611 err_msg="do_time_min_wrapper::the output buffer and the buffer send in are not of the same type (r4_kind)"
614 end function do_time_min_wrapper
616 !> @brief Does the time_min reduction method on the buffer object
617 !! @return Error message if the math was not successful
618 function do_time_max_wrapper(this, field_data, mask, is_masked, bounds_in, bounds_out, missing_value) &
620 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< buffer object to write
621 class(*), intent(in) :: field_data(:,:,:,:) !< Buffer data for current time
622 type(fmsDiagIbounds_type), intent(in) :: bounds_in !< Indicies for the buffer passed in
623 type(fmsDiagIbounds_type), intent(in) :: bounds_out !< Indicies for the output buffer
624 logical, intent(in) :: mask(:,:,:,:) !< Mask for the field
625 logical, intent(in) :: is_masked !< .True. if the field has a mask
626 real(kind=r8_kind), intent(in) :: missing_value !< Missing_value for data points that are masked
627 character(len=50) :: err_msg
629 !TODO This will be expanded for integers
631 select type (output_buffer => this%buffer)
632 type is (real(kind=r8_kind))
633 select type (field_data)
634 type is (real(kind=r8_kind))
635 call do_time_max(output_buffer, field_data, mask, is_masked, bounds_in, bounds_out, missing_value)
637 err_msg="do_time_max_wrapper::the output buffer and the buffer send in are not of the same type (r8_kind)"
639 type is (real(kind=r4_kind))
640 select type (field_data)
641 type is (real(kind=r4_kind))
642 call do_time_max(output_buffer, field_data, mask, is_masked, bounds_in, bounds_out, &
643 real(missing_value, kind=r4_kind))
645 err_msg="do_time_max_wrapper::the output buffer and the buffer send in are not of the same type (r4_kind)"
648 end function do_time_max_wrapper
650 !> @brief Does the time_sum reduction method on the buffer object
651 !! @return Error message if the math was not successful
652 function do_time_sum_wrapper(this, field_data, mask, is_masked, mask_variant, bounds_in, bounds_out, missing_value, &
653 has_missing_value, pow_value) &
655 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< buffer object to write
656 class(*), intent(in) :: field_data(:,:,:,:) !< Buffer data for current time
657 type(fmsDiagIbounds_type), intent(in) :: bounds_in !< Indicies for the buffer passed in
658 type(fmsDiagIbounds_type), intent(in) :: bounds_out !< Indicies for the output buffer
659 logical, intent(in) :: mask(:,:,:,:) !< Mask for the field
660 logical, intent(in) :: is_masked !< .True. if the field has a mask
661 logical, intent(in) :: mask_variant !< .True. if the mask changes over time
662 real(kind=r8_kind), intent(in) :: missing_value !< Missing_value for data points that are masked
663 logical, intent(in) :: has_missing_value !< .True. if the field was registered with
665 integer, optional, intent(in) :: pow_value !< power value, will calculate field_data^pow
666 !! before adding to buffer should only be
667 !! present if using pow reduction method
668 character(len=150) :: err_msg
670 !TODO This will be expanded for integers
672 select type (output_buffer => this%buffer)
673 type is (real(kind=r8_kind))
674 select type (field_data)
675 type is (real(kind=r8_kind))
676 if (.not. is_masked) then
677 if (any(field_data .eq. missing_value)) &
678 err_msg = "You cannot pass data with missing values without masking them!"
680 call do_time_sum_update(output_buffer, this%weight_sum, field_data, mask, is_masked, mask_variant, &
681 bounds_in, bounds_out, missing_value, this%diurnal_section, &
684 err_msg="do_time_sum_wrapper::the output buffer and the buffer send in are not of the same type (r8_kind)"
686 type is (real(kind=r4_kind))
687 select type (field_data)
688 type is (real(kind=r4_kind))
689 if (.not. is_masked) then
690 if (any(field_data .eq. missing_value)) &
691 err_msg = "You cannot pass data with missing values without masking them!"
693 call do_time_sum_update(output_buffer, this%weight_sum, field_data, mask, is_masked, mask_variant, &
694 bounds_in, bounds_out, real(missing_value, kind=r4_kind), &
695 this%diurnal_section, pow=pow_value)
697 err_msg="do_time_sum_wrapper::the output buffer and the buffer send in are not of the same type (r4_kind)"
700 err_msg="do_time_sum_wrapper::the output buffer is not a valid type, must be real(r8_kind) or real(r4_kind)"
702 end function do_time_sum_wrapper
704 !> Finishes calculations for any reductions that use an average (avg, rms, pow)
705 !! TODO add mask and any other needed args for adjustment, and pass in the adjusted mask
706 !! to time_update_done
707 function diag_reduction_done_wrapper(this, reduction_method, missing_value, has_mask, mask_variant) &
709 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< Updated buffer object
710 integer, intent(in) :: reduction_method !< enumerated reduction type from diag_data
711 real(kind=r8_kind), intent(in) :: missing_value !< missing_value for masked data points
712 logical, intent(in) :: has_mask !< indicates if there was a mask used during buffer updates
713 logical, intent(in) :: mask_variant !< Indicates if the mask changes over time
714 character(len=51) :: err_msg !< error message to return, blank if sucessful
716 if(.not. allocated(this%buffer)) return
719 select type(buff => this%buffer)
720 type is (real(r8_kind))
721 call time_update_done(buff, this%weight_sum, reduction_method, missing_value, has_mask, mask_variant, &
722 this%diurnal_sample_size)
723 type is (real(r4_kind))
724 call time_update_done(buff, this%weight_sum, reduction_method, real(missing_value, r4_kind), has_mask, &
725 mask_variant, this%diurnal_sample_size)
727 this%weight_sum = 0.0_r8_kind
731 !> this leaves out the diurnal index cause its only used for tmp mask allocation
732 pure function get_buffer_dims(this)
733 class(fmsDiagOutputBuffer_type), intent(in) :: this !< buffer object to get from
734 integer :: get_buffer_dims(4)
735 get_buffer_dims = this%buffer_dims(1:4)
738 !> Get diurnal sample size (amount of diurnal sections)
739 pure integer function get_diurnal_sample_size(this)
740 class(fmsDiagOutputBuffer_type), intent(in) :: this !< buffer object to get from
741 get_diurnal_sample_size = this%diurnal_sample_size
742 end function get_diurnal_sample_size
744 !> Set diurnal sample size (amount of diurnal sections)
745 subroutine set_diurnal_sample_size(this, sample_size)
746 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< buffer object to set sample size for
747 integer, intent(in) :: sample_size !< sample size to used to split daily
748 !! data into given amount of sections
749 this%diurnal_sample_size = sample_size
750 end subroutine set_diurnal_sample_size
752 !> Set diurnal section index based off the current time and previously set diurnal_samplesize
753 !! Calculates which diurnal section of daily data the current time is in
754 subroutine set_diurnal_section_index(this, time)
755 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< buffer object to set diurnal index for
756 type(time_type), intent(in) :: time !< current model time
757 integer :: seconds, days, ticks
759 if(this%diurnal_sample_size .lt. 0) call mpp_error(FATAL, "set_diurnal_section_index::"// &
760 " diurnal sample size must be set before trying to set diurnal index for send_data")
762 call get_time(time,seconds,days,ticks) ! get current date
763 ! calculates which diurnal section current time is in for a given amount of diurnal sections(<24)
764 this%diurnal_section = floor( (seconds+real(ticks)/get_ticks_per_second()) &
765 & * this%diurnal_sample_size/SECONDS_PER_DAY) + 1
766 end subroutine set_diurnal_section_index
768 !> Remaps the output buffer array when using the diurnal reduction
769 !! moves the diurnal index to the left-most unused dimension for the io
770 subroutine get_remapped_diurnal_data(this, res)
771 class(fmsDiagOutputBuffer_type), intent(in) :: this !< output buffer object
772 class(*), intent(out), allocatable :: res(:,:,:,:,:) !< resulting remapped data
773 integer :: last_dim !< last dimension thats used
774 integer :: ie, je, ke, ze, de !< ending indices for the new array
775 integer(i4_kind) :: buff_size(5)!< sizes for allocated buffer
777 ! last dim is number of dimensions - 1 for diurnal axis
778 last_dim = this%ndim - 1
779 ! get the bounds of the remapped output array based on # of dims
780 ke = 1; ze = 1; de = 1
781 select case(last_dim)
783 ie = this%buffer_dims(1); je = this%buffer_dims(5)
785 ie = this%buffer_dims(1); je = this%buffer_dims(2)
786 ke = this%buffer_dims(5)
788 ie = this%buffer_dims(1); je = this%buffer_dims(2)
789 ke = this%buffer_dims(3); ze = this%buffer_dims(5)
791 ! no need to remap if 4d
796 select type(buff => this%buffer)
797 type is (real(r8_kind))
798 allocate(real(r8_kind) :: res(1:ie, 1:je, 1:ke, 1:ze, 1:de))
800 type is (real(r8_kind))
801 res(1:ie, 1:je, 1:ke, 1:ze, 1:de) = reshape(buff, SHAPE(res))
803 type is (real(r4_kind))
804 allocate(real(r4_kind) :: res(1:ie, 1:je, 1:ke, 1:ze, 1:de))
806 type is (real(r4_kind))
807 res(1:ie, 1:je, 1:ke, 1:ze, 1:de) = reshape(buff, SHAPE(res))
809 type is (integer(i8_kind))
810 allocate(integer(i8_kind) :: res(1:ie, 1:je, 1:ke, 1:ze, 1:de))
812 type is (integer(i8_kind))
813 res(1:ie, 1:je, 1:ke, 1:ze, 1:de) = reshape(buff, SHAPE(res))
815 type is (integer(i4_kind))
816 allocate(integer(i4_kind) :: res(1:ie, 1:je, 1:ke, 1:ze, 1:de))
818 type is (integer(i4_kind))
819 res(1:ie, 1:je, 1:ke, 1:ze, 1:de) = reshape(buff, SHAPE(res))
823 end subroutine get_remapped_diurnal_data
825 !> @brief Determine if there is any data to write (i.e send_data has been called)
826 !! @return .true. if there is data to write
827 function is_there_data_to_write(this) &
829 class(fmsDiagOutputBuffer_type), intent(in) :: this !< Buffer object
833 if (allocated(this%send_data_called)) then
834 res = this%send_data_called
840 !> @brief Determine if it is time to finish the reduction method
841 !! @return .true. if it is time to finish the reduction method
842 function is_time_to_finish_reduction(this, end_time) &
844 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< Buffer object
845 type(time_type), optional, intent(in) :: end_time !< The time at the end of the run
850 if (this%time >= this%next_output) res = .true.
852 if (present(end_time)) then
853 if (end_time >= this%next_output) res = .true.
855 end function is_time_to_finish_reduction
857 !> @brief Sets send_data_called to .true.
858 subroutine set_send_data_called(this)
859 class(fmsDiagOutputBuffer_type), intent(inout) :: this !< Buffer object
861 this%send_data_called = .true.
862 end subroutine set_send_data_called
864 end module fms_diag_output_buffer_mod